Intrepid2
Intrepid2_OrientationToolsDefCoeffMatrix_HVOL.hpp
Go to the documentation of this file.
1// @HEADER
2// ************************************************************************
3//
4// Intrepid2 Package
5// Copyright (2007) Sandia Corporation
6//
7// Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8// license for use of this work by or on behalf of the U.S. Government.
9//
10// Redistribution and use in source and binary forms, with or without
11// modification, are permitted provided that the following conditions are
12// met:
13//
14// 1. Redistributions of source code must retain the above copyright
15// notice, this list of conditions and the following disclaimer.
16//
17// 2. Redistributions in binary form must reproduce the above copyright
18// notice, this list of conditions and the following disclaimer in the
19// documentation and/or other materials provided with the distribution.
20//
21// 3. Neither the name of the Corporation nor the names of the
22// contributors may be used to endorse or promote products derived from
23// this software without specific prior written permission.
24//
25// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36//
37// Questions? Contact Kyungjoo Kim (kyukim@sandia.gov), or
38// Mauro Perego (mperego@sandia.gov)
39//
40// ************************************************************************
41// @HEADER
42
43
56#ifndef __INTREPID2_ORIENTATIONTOOLS_DEF_COEFF_MATRIX_HVOL_HPP__
57#define __INTREPID2_ORIENTATIONTOOLS_DEF_COEFF_MATRIX_HVOL_HPP__
58
59// disable clang warnings
60#if defined (__clang__) && !defined (__INTEL_COMPILER)
61#pragma clang system_header
62#endif
63
64namespace Intrepid2 {
65
66namespace Impl {
67namespace Debug {
68
69#ifdef HAVE_INTREPID2_DEBUG
70template<typename cellBasisType>
71inline
72void
73check_getCoeffMatrix_HVOL(const cellBasisType& cellBasis,
74 const ordinal_type cellOrt) {
75
76 // populate points on a subcell and map to subcell
77 const shards::CellTopology cellTopo = cellBasis.getBaseCellTopology();
78 const ordinal_type cellDim = cellTopo.getDimension();
79
80 INTREPID2_TEST_FOR_EXCEPTION( cellDim > 2,
81 std::logic_error,
82 ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HVOL): " \
83 "HVOL orientation supported only for (side) cells with dimension less than 3.");
84
85 const auto cellBaseKey = cellTopo.getBaseKey();
86
87 INTREPID2_TEST_FOR_EXCEPTION( cellBaseKey != shards::Line<>::key &&
88 cellBaseKey != shards::Quadrilateral<>::key &&
89 cellBaseKey != shards::Triangle<>::key,
90 std::logic_error,
91 ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HVOL): " \
92 "cellBasis must have line, quad, or triangle topology.");
93
94
95 //
96 // Function space
97 //
98
99 {
100 INTREPID2_TEST_FOR_EXCEPTION( cellBasis.getFunctionSpace() != FUNCTION_SPACE_HVOL,
101 std::logic_error,
102 ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HVOL): " \
103 "cellBasis function space is not HVOL.");
104 }
105 }
106#endif
107} // Debug namespace
108
109template<typename OutputViewType,
110typename cellBasisHostType>
111inline
112void
114getCoeffMatrix_HVOL(OutputViewType &output,
116 const ordinal_type cellOrt) {
117
118#ifdef HAVE_INTREPID2_DEBUG
119 Debug::check_getCoeffMatrix_HVOL(cellBasis,cellOrt);
120#endif
121
122 using host_device_type = typename Kokkos::HostSpace::device_type;
123 using value_type = typename OutputViewType::non_const_value_type;
124
125 //
126 // Topology
127 //
128
129 const shards::CellTopology cellTopo = cellBasis.getBaseCellTopology();
130 const ordinal_type cellDim = cellTopo.getDimension();
131 const auto cellBaseKey = cellTopo.getBaseKey();
132 const ordinal_type cardinality = cellBasis.getCardinality();
133
134 //
135 // Reference points
136 //
137
138 // Reference points \xi_j on the subcell
139 Kokkos::DynRankView<value_type,host_device_type> refPtsCell("refPtsCell", cardinality, cellDim),refPtsCellNotOriented("refPtsCellNotOriented", cardinality, cellDim);
140
141 ordinal_type latticeOffset(1);
142
143 // this work for line and quadrilateral topologies
144 ordinal_type latticeOrder = (cellTopo.getBaseKey() == shards::Triangle<>::key) ?
145 cellBasis.getDegree() + 3 * latticeOffset : // triangle
146 cellBasis.getDegree() + 2 * latticeOffset; // line and quad
147
149
150 // map the points into the parent, cell accounting for orientation
152
153 //
154 // Bases evaluation on the reference points
155 //
156
157 // cellBasisValues = \psi_k(\eta_o (\xi_j))
158 Kokkos::DynRankView<value_type,host_device_type> cellBasisValues("cellBasisValues", cardinality, cardinality);
159
160 // basisValues = \phi_i (\xi_j)
161 Kokkos::DynRankView<value_type,host_device_type> nonOrientedBasisValues("subcellBasisValues", cardinality, cardinality);
162
163 cellBasis.getValues(cellBasisValues, refPtsCell, OPERATOR_VALUE);
164 cellBasis.getValues(nonOrientedBasisValues, refPtsCellNotOriented, OPERATOR_VALUE);
165
166 //
167 // Compute Psi_jk = \psi_k(\eta_o (\xi_j)) det(J_\eta) and Phi_ji = \phi_i (\xi_j),
168 // and solve
169 // Psi A^T = Phi
170 //
171
172 // construct Psi and Phi matrices. LAPACK wants left layout
173 Kokkos::DynRankView<value_type,Kokkos::LayoutLeft,host_device_type>
174 PsiMat("PsiMat", cardinality, cardinality),
175 PhiMat("PhiMat", cardinality, cardinality);
176
177 auto cellTagToOrdinal = cellBasis.getAllDofOrdinal();
178
179 Kokkos::DynRankView<value_type,host_device_type> jac("jacobian",cellDim,cellDim);
181 value_type jacDet(0);
182 if(cellDim == 2) {
183 jacDet = jac(0,0)*jac(1,1)-jac(0,1)*jac(1,0);
184 } else { //celldim == 1
185 jacDet = jac(0,0);
186 }
187
188 for (ordinal_type i=0;i<cardinality;++i) {
189 const ordinal_type ic = cellTagToOrdinal(cellDim, 0, i);
190 for (ordinal_type j=0;j<cardinality;++j) {
193 }
194 }
195
196 // Solve the system
197 {
198 Teuchos::LAPACK<ordinal_type,value_type> lapack;
199 ordinal_type info = 0;
200
201 Kokkos::DynRankView<ordinal_type,Kokkos::LayoutLeft,host_device_type> pivVec("pivVec", cardinality);
203 PsiMat.data(),
204 PsiMat.stride_1(),
205 pivVec.data(),
206 PhiMat.data(),
207 PhiMat.stride_1(),
208 &info);
209
210 if (info) {
211 std::stringstream ss;
212 ss << ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HVOL): "
213 << "LAPACK return with error code: "
214 << info;
215 INTREPID2_TEST_FOR_EXCEPTION( true, std::runtime_error, ss.str().c_str() );
216 }
217
218 //After solving the system w/ LAPACK, Phi contains A^T
219
220 // transpose B and clean up numerical noise (for permutation matrices)
221 const double eps = tolerence();
222 for (ordinal_type i=0;i<cardinality;++i) {
223 auto intmatii = std::round(PhiMat(i,i));
224 PhiMat(i,i) = (std::abs(PhiMat(i,i) - intmatii) < eps) ? intmatii : PhiMat(i,i);
225 for (ordinal_type j=i+1;j<cardinality;++j) {
226 auto matij = PhiMat(i,j);
227
228 auto intmatji = std::round(PhiMat(j,i));
229 PhiMat(i,j) = (std::abs(PhiMat(j,i) - intmatji) < eps) ? intmatji : PhiMat(j,i);
230
231 auto intmatij = std::round(matij);
232 PhiMat(j,i) = (std::abs(matij - intmatij) < eps) ? intmatij : matij;
233 }
234 }
235
236 }
237
238 // Print A Matrix
239 /*
240 {
241 std::cout << "Ort: " << cellOrt << ": |";
242 for (ordinal_type i=0;i<cardinality;++i) {
243 for (ordinal_type j=0;j<cardinality;++j) {
244 std::cout << PhiMat(i,j) << " ";
245 }
246 std::cout << "| ";
247 }
248 std::cout <<std::endl;
249 }
250 */
251
252 {
253 // move the data to original device memory
254 const Kokkos::pair<ordinal_type,ordinal_type> range(0, cardinality);
255 auto suboutput = Kokkos::subview(output, range, range);
256 auto tmp = Kokkos::create_mirror_view_and_copy(typename OutputViewType::device_type::memory_space(), PhiMat);
257 Kokkos::deep_copy(suboutput, tmp);
258 }
259}
260}
261
262}
263#endif
static void getCoeffMatrix_HVOL(OutputViewType &output, const cellBasisHostType &cellBasis, const ordinal_type cellOrt)
Compute orientation matrix for HVOL basis for a given (2D or 1D) cell and its reference basis....
static void getJacobianOfOrientationMap(JacobianViewType jacobian, const shards::CellTopology cellTopo, const ordinal_type cellOrt)
Computes jacobian of the parameterization maps of 1- and 2-subcells with orientation.
static void mapToModifiedReference(outPointViewType outPoints, const refPointViewType refPoints, const shards::CellTopology cellTopo, const ordinal_type cellOrt=0)
Computes modified parameterization maps of 1- and 2-subcells with orientation.
static void getLattice(Kokkos::DynRankView< pointValueType, pointProperties... > points, const shards::CellTopology cellType, const ordinal_type order, const ordinal_type offset=0, const EPointType pointType=POINTTYPE_EQUISPACED)
Computes a lattice of points of a given order on a reference simplex, quadrilateral or hexahedron (cu...