Intrepid2
Intrepid2_OrientationToolsDefCoeffMatrix_HDIV.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
59#ifndef __INTREPID2_ORIENTATIONTOOLS_DEF_COEFF_MATRIX_HDIV_HPP__
60#define __INTREPID2_ORIENTATIONTOOLS_DEF_COEFF_MATRIX_HDIV_HPP__
61
62// disable clang warnings
63#if defined (__clang__) && !defined (__INTEL_COMPILER)
64#pragma clang system_header
65#endif
66
67namespace Intrepid2 {
68
69namespace Impl {
70
71namespace Debug {
72
73#ifdef HAVE_INTREPID2_DEBUG
74template<typename subcellBasisType,
75typename cellBasisType>
76inline
77void
78check_getCoeffMatrix_HDIV(const subcellBasisType& subcellBasis,
79 const cellBasisType& cellBasis,
80 const ordinal_type subcellId,
81 const ordinal_type subcellOrt) {
82 const shards::CellTopology cellTopo = cellBasis.getBaseCellTopology();
83 const shards::CellTopology subcellTopo = subcellBasis.getBaseCellTopology();
84
85 const ordinal_type cellDim = cellTopo.getDimension();
86 const ordinal_type subcellDim = subcellTopo.getDimension();
87
88 INTREPID2_TEST_FOR_EXCEPTION( subcellDim >= cellDim,
89 std::logic_error,
90 ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HDIV): " \
91 "cellDim must be greater than subcellDim.");
92
93 const auto subcellBaseKey = subcellTopo.getBaseKey();
94 const auto cellBaseKey = cellTopo.getBaseKey();
95
96 INTREPID2_TEST_FOR_EXCEPTION( subcellBaseKey != shards::Line<>::key &&
97 subcellBaseKey != shards::Quadrilateral<>::key &&
98 subcellBaseKey != shards::Triangle<>::key,
99 std::logic_error,
100 ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HDIV): " \
101 "subcellBasis must have line, quad, or triangle topology.");
102
103 INTREPID2_TEST_FOR_EXCEPTION( cellBaseKey != shards::Quadrilateral<>::key &&
104 cellBaseKey != shards::Triangle<>::key &&
105 cellBaseKey != shards::Hexahedron<>::key &&
106 cellBaseKey != shards::Wedge<>::key &&
107 cellBaseKey != shards::Tetrahedron<>::key,
108 std::logic_error,
109 ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HDIV): " \
110 "cellBasis must have quad, triangle, hexhedron or tetrahedron topology.");
111
112 //
113 // Function space
114 //
115 {
116 const bool isHDIV = cellBasis.getFunctionSpace() == FUNCTION_SPACE_HDIV;
117 INTREPID2_TEST_FOR_EXCEPTION( !isHDIV,
118 std::logic_error,
119 ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HDIV): "
120 "cellBasis is not HDIV.");
121 {
122 const bool subcellBasisIsHGRAD = subcellBasis.getFunctionSpace() == FUNCTION_SPACE_HGRAD;
123 const bool subcellBasisIsHVOL = subcellBasis.getFunctionSpace() == FUNCTION_SPACE_HVOL;
124 const bool cellIsTri = cellBaseKey == shards::Triangle<>::key;
125 const bool cellIsTet = cellBaseKey == shards::Tetrahedron<>::key;
126 const bool cellIsHex = cellBaseKey == shards::Hexahedron<>::key;
127 const bool cellIsQuad = cellBaseKey == shards::Quadrilateral<>::key;
128
129 switch (subcellDim) {
130 case 1: {
131 //TODO: Hex, QUAD, TET and TRI element should have the same 1d basis
132 if (cellIsHex || cellIsQuad) {
133 INTREPID2_TEST_FOR_EXCEPTION( !(subcellBasisIsHGRAD||subcellBasisIsHVOL),
134 std::logic_error,
135 ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_DIV): "
136 "subcellBasis function space (1d) is not consistent to cellBasis, which should be open line hgrad, order -1.");
137 } else if (cellIsTet || cellIsTri) {
138 INTREPID2_TEST_FOR_EXCEPTION( !subcellBasisIsHVOL,
139 std::logic_error,
140 ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_DIV): "
141 "subcellBasis function space (1d) is not consistent to cellBasis, which should be HVOL line, order -1.");
142 }
143 break;
144 }
145 case 2: {
146 if (subcellBaseKey == shards::Quadrilateral<>::key) {
147 // quad face basis is tensor product of open line basis functions
148 INTREPID2_TEST_FOR_EXCEPTION( !(subcellBasisIsHGRAD||subcellBasisIsHVOL),
149 std::logic_error,
150 ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HDIV): "
151 "subcellBasis function space is not compatible, which should be open line hgrad, order -1.");
152 } else if (subcellBaseKey == shards::Triangle<>::key) {
153 // triangle face basis comes from HVOL basis
154 INTREPID2_TEST_FOR_EXCEPTION( !subcellBasisIsHVOL,
155 std::logic_error,
156 ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HDIV): "
157 "subcellBasis function space is not compatible, which should HVOL, order-1.");
158 }
159 break;
160 }
161 }
162 }
163 }
164}
165#endif
166} //Debug Namespace
167
168template<typename OutputViewType,
169typename subcellBasisHostType,
170typename cellBasisHostType>
171inline
172void
174getCoeffMatrix_HDIV(OutputViewType &output,
177 const ordinal_type subcellId,
178 const ordinal_type subcellOrt) {
179
180#ifdef HAVE_INTREPID2_DEBUG
181 Debug::check_getCoeffMatrix_HDIV(subcellBasis,cellBasis,subcellId,subcellOrt);
182#endif
183
184 using value_type = typename OutputViewType::non_const_value_type;
185 using host_device_type = Kokkos::HostSpace::device_type;
186
187 //
188 // Collocation points
189 //
190 const shards::CellTopology cellTopo = cellBasis.getBaseCellTopology();
191 const shards::CellTopology subcellTopo = subcellBasis.getBaseCellTopology();
192 const ordinal_type cellDim = cellTopo.getDimension();
193 const ordinal_type subcellDim = subcellTopo.getDimension();
194 const auto subcellBaseKey = subcellTopo.getBaseKey();
195 const ordinal_type numCellBasis = cellBasis.getCardinality();
196 const ordinal_type numSubcellBasis = subcellBasis.getCardinality();
197 const ordinal_type ndofSubcell = cellBasis.getDofCount(subcellDim,subcellId);
198
199 //
200 // Reference points
201 //
202
203 //use lattice to compute reference subcell points \xi_j
204 auto latticeDegree = (subcellBaseKey == shards::Triangle<>::key) ?
205 cellBasis.getDegree()+2 : cellBasis.getDegree()+1;
206
207 // Reference points \xi_j on the subcell
208 Kokkos::DynRankView<value_type,host_device_type> refPtsSubcell("refPtsSubcell", ndofSubcell, subcellDim);
210 INTREPID2_TEST_FOR_EXCEPTION( latticeSize != ndofSubcell,
211 std::logic_error,
212 ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HDIV): " \
213 "Lattice size should be equal to the onber of subcell internal DoFs");
214 PointTools::getLattice(refPtsSubcell, subcellTopo, latticeDegree, 1);//, POINTTYPE_WARPBLEND);
215
216 // evaluate values on the modified cell
218
219 // refPtsCell = F_s (\eta_o (refPtsSubcell))
220 Kokkos::DynRankView<value_type,host_device_type> refPtsCell("refPtsCell", ndofSubcell, cellDim);
221 // map points from the subcell manifold into the cell one
223
224 //computing normal to the subcell accounting for orientation
225 Kokkos::DynRankView<value_type,host_device_type> tangentsAndNormal("trJacobianF", cellDim, cellDim );
227 auto sideNormal = Kokkos::subview(tangentsAndNormal, cellDim-1, Kokkos::ALL());
228
229
230 //
231 // Basis evaluation on the collocation points
232 //
233
234 // cellBasisValues = \psi_k(F_s (\eta_o (\xi_j)))
235 Kokkos::DynRankView<value_type,host_device_type> cellBasisValues("cellBasisValues", numCellBasis, ndofSubcell, cellDim);
236 cellBasis.getValues(cellBasisValues, refPtsCell, OPERATOR_VALUE);
237
238 // subcellBasisValues = \phi_i (\xi_j)
239 Kokkos::DynRankView<value_type,host_device_type> subCellValues("subCellValues", numSubcellBasis, ndofSubcell);
240 subcellBasis.getValues(subCellValues, refPtsSubcell, OPERATOR_VALUE);
241
242 //
243 // Compute Psi_jk = \psi_k(F_s (\eta_o (\xi_j))) \cdot (n_s det(J_\eta))
244 // and Phi_ji = \phi_i (\xi_j), and solve
245 // Psi A^T = Phi
246 //
247
248 // construct Psi and Phi matrices. LAPACK wants left layout
249 Kokkos::DynRankView<value_type,Kokkos::LayoutLeft,host_device_type>
250 PsiMat("PsiMat", ndofSubcell, ndofSubcell),
251 PhiMat("PhiMat", ndofSubcell, ndofSubcell);
252
253 auto cellTagToOrdinal = cellBasis.getAllDofOrdinal();
254 auto subcellTagToOrdinal = subcellBasis.getAllDofOrdinal();
255
256 for (ordinal_type i=0;i<ndofSubcell;++i) {
257 const ordinal_type ic = cellTagToOrdinal(subcellDim, subcellId, i);
258 const ordinal_type isc = subcellTagToOrdinal(subcellDim, 0, i);
259 for (ordinal_type j=0;j<ndofSubcell;++j) {
261 for (ordinal_type k=0;k<cellDim;++k)
263 }
264 }
265
266 auto RefMat = PsiMat;
267 auto OrtMat = PhiMat;
268
269 // solve the system
270 {
271 Teuchos::LAPACK<ordinal_type,value_type> lapack;
272 ordinal_type info = 0;
273 Kokkos::DynRankView<ordinal_type,host_device_type> pivVec("pivVec", ndofSubcell);
274
276 PsiMat.data(),
277 PsiMat.stride_1(),
278 pivVec.data(),
279 PhiMat.data(),
280 PhiMat.stride_1(),
281 &info);
282
283 if (info) {
284 std::stringstream ss;
285 ss << ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HDIV): "
286 << "LAPACK return with error code: "
287 << info;
288 INTREPID2_TEST_FOR_EXCEPTION( true, std::runtime_error, ss.str().c_str() );
289 }
290
291 //After solving the system w/ LAPACK, Phi contains A^T
292
293 // transpose and clean up numerical noise (for permutation matrices)
294 const double eps = tolerence();
295 for (ordinal_type i=0;i<ndofSubcell;++i) {
296 auto intmatii = std::round(PhiMat(i,i));
297 PhiMat(i,i) = (std::abs(PhiMat(i,i) - intmatii) < eps) ? intmatii : PhiMat(i,i);
298 for (ordinal_type j=i+1;j<ndofSubcell;++j) {
299 auto matij = PhiMat(i,j);
300
301 auto intmatji = std::round(PhiMat(j,i));
302 PhiMat(i,j) = (std::abs(PhiMat(j,i) - intmatji) < eps) ? intmatji : PhiMat(j,i);
303
304 auto intmatij = std::round(matij);
305 PhiMat(j,i) = (std::abs(matij - intmatij) < eps) ? intmatij : matij;
306 }
307 }
308
309 // Print Matrix A
310 /*
311 {
312 std::cout << "|";
313 for (ordinal_type i=0;i<ndofSubcell;++i) {
314 for (ordinal_type j=0;j<ndofSubcell;++j) {
315 std::cout << PhiMat(i,j) << " ";
316 }
317 std::cout << "| ";
318 }
319 std::cout <<std::endl;
320 }
321 */
322
323 }
324
325 {
326 // move the data to original device memory
327 const Kokkos::pair<ordinal_type,ordinal_type> range(0, ndofSubcell);
328 auto suboutput = Kokkos::subview(output, range, range);
329 auto tmp = Kokkos::create_mirror_view_and_copy(typename OutputViewType::device_type::memory_space(), PhiMat);
330 Kokkos::deep_copy(suboutput, tmp);
331 }
332}
333}
334
335}
336#endif
KOKKOS_FORCEINLINE_FUNCTION constexpr std::enable_if<!std::is_pod< T >::value, typenameScalarTraits< T >::scalar_type >::type get_scalar_value(const T &obj)
functions returning the scalar value. for pod types, they return the input object itself....
static KOKKOS_INLINE_FUNCTION void getRefSideTangentsAndNormal(TanNormViewType tangentsAndNormal, const ParamViewType subCellParametrization, const unsigned subcellTopoKey, const ordinal_type subCellOrd, const ordinal_type ort)
Computes the (oriented) tangents and normal of the side subCell The normals are only defined for side...
static KOKKOS_INLINE_FUNCTION void mapSubcellCoordsToRefCell(coordsViewType cellCoords, const subcellCoordsViewType subCellCoords, const ParamViewType subcellParametrization, const unsigned subcellTopoKey, const ordinal_type subCellOrd, const ordinal_type ort)
Maps points defined on the subCell manifold into the parent Cell.
static void getCoeffMatrix_HDIV(OutputViewType &output, const subcellBasisHostType &subcellBasis, const cellBasisHostType &cellBasis, const ordinal_type subcellId, const ordinal_type subcellOrt)
Compute orientation matrix for HDIV basis for a given subcell and its reference basis.
static ordinal_type getLatticeSize(const shards::CellTopology cellType, const ordinal_type order, const ordinal_type offset=0)
Computes the number of points in a lattice of a given order on a simplex (currently disabled for othe...
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...
static ConstViewType get(const ordinal_type subcellDim, const unsigned parentCellKey)
Returns a Kokkos view with the coefficients of the parametrization maps for the edges or faces of a r...