Intrepid2
Intrepid2_OrientationTools.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Intrepid2 Package
4//
5// Copyright 2007 NTESS and the Intrepid2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
15#ifndef __INTREPID2_ORIENTATIONTOOLS_HPP__
16#define __INTREPID2_ORIENTATIONTOOLS_HPP__
17
18#include "Intrepid2_ConfigDefs.hpp"
19#include "Intrepid2_Types.hpp"
20#include "Intrepid2_Utils.hpp"
21
22#include "Shards_CellTopology.hpp"
23#include "Shards_BasicTopologies.hpp"
24
26
27#include "Intrepid2_Basis.hpp"
28
29// -- HGRAD family
33
36
37// -- HCURL family
40
44
45// -- HDIV family
48
52
53// -- Lower order family
56
59
63
67
70
71#include "Teuchos_LAPACK.hpp"
72
73
74namespace Intrepid2 {
75
76 namespace Impl {
77
82 public:
83
84 // -----------------------------------------------------------------------------
85 // Point modification
86 //
87 //
88
95 template<typename ValueType>
96 KOKKOS_INLINE_FUNCTION
97 static void
98 getModifiedLinePoint(ValueType &ot,
99 const ValueType pt,
100 const ordinal_type ort);
101
110 template<typename ValueType>
111 KOKKOS_INLINE_FUNCTION
112 static void
114 ValueType &ot1,
115 const ValueType pt0,
116 const ValueType pt1,
117 const ordinal_type ort);
118
127 template<typename ValueType>
128 KOKKOS_INLINE_FUNCTION
129 static void
131 ValueType &ot1,
132 const ValueType pt0,
133 const ValueType pt1,
134 const ordinal_type ort);
135
143 template<typename outPointViewType,
144 typename refPointViewType>
145 inline
146 static void
147 mapToModifiedReference(outPointViewType outPoints,
148 const refPointViewType refPoints,
149 const shards::CellTopology cellTopo,
150 const ordinal_type cellOrt = 0);
151
159 template<typename outPointViewType,
160 typename refPointViewType>
161 KOKKOS_INLINE_FUNCTION
162 static void
163 mapToModifiedReference(outPointViewType outPoints,
164 const refPointViewType refPoints,
165 const unsigned cellTopoKey,
166 const ordinal_type cellOrt = 0);
167
168
174 template<typename JacobianViewType>
175 KOKKOS_INLINE_FUNCTION
176 static void
177 getLineJacobian(JacobianViewType jacobian, const ordinal_type ort);
178
184 template<typename JacobianViewType>
185 KOKKOS_INLINE_FUNCTION
186 static void
187 getTriangleJacobian(JacobianViewType jacobian, const ordinal_type ort);
188
194 template<typename JacobianViewType>
195 KOKKOS_INLINE_FUNCTION
196 static void
197 getQuadrilateralJacobian(JacobianViewType jacobian, const ordinal_type ort);
198
199
206 template<typename JacobianViewType>
207 inline
208 static void
209 getJacobianOfOrientationMap(JacobianViewType jacobian,
210 const shards::CellTopology cellTopo,
211 const ordinal_type cellOrt);
212
219 template<typename JacobianViewType>
220 KOKKOS_INLINE_FUNCTION
221 static void
222 getJacobianOfOrientationMap(JacobianViewType jacobian,
223 const unsigned cellTopoKey,
224 const ordinal_type cellOrt);
225
226
234 template<typename TanViewType, typename ParamViewType>
235 KOKKOS_INLINE_FUNCTION
236 static void getRefSubcellTangents(TanViewType tangents,
237 const ParamViewType subCellParametrization,
238 const unsigned subcellTopoKey,
239 const ordinal_type subCellOrd,
240 const ordinal_type ort);
241
242
253 template<typename TanNormViewType, typename ParamViewType>
254 KOKKOS_INLINE_FUNCTION
255 static void getRefSideTangentsAndNormal(TanNormViewType tangentsAndNormal,
256 const ParamViewType subCellParametrization,
257 const unsigned subcellTopoKey,
258 const ordinal_type subCellOrd,
259 const ordinal_type ort);
260
261
271 template<typename coordsViewType, typename subcellCoordsViewType, typename ParamViewType>
272 KOKKOS_INLINE_FUNCTION
273 static void mapSubcellCoordsToRefCell(coordsViewType cellCoords,
274 const subcellCoordsViewType subCellCoords,
275 const ParamViewType subcellParametrization,
276 const unsigned subcellTopoKey,
277 const ordinal_type subCellOrd,
278 const ordinal_type ort);
279
280 // -----------------------------------------------------------------------------
281 // Coefficient Matrix
282 //
283 //
284
296 template<typename OutputViewType,
297 typename subcellBasisHostType,
298 typename cellBasisHostType>
299 inline
300 static void
301 getCoeffMatrix_HGRAD(OutputViewType &output,
302 const subcellBasisHostType& subcellBasis,
303 const cellBasisHostType& cellBasis,
304 const ordinal_type subcellId,
305 const ordinal_type subcellOrt,
306 const bool inverse = false);
307
318 template<typename OutputViewType,
319 typename subcellBasisHostType,
320 typename cellBasisHostType>
321 inline
322 static void
323 getCoeffMatrix_HCURL(OutputViewType &output,
324 const subcellBasisHostType& subcellBasis,
325 const cellBasisHostType& cellBasis,
326 const ordinal_type subcellId,
327 const ordinal_type subcellOrt,
328 const bool inverse = false);
329
330
341 template<typename OutputViewType,
342 typename subcellBasisHostType,
343 typename cellBasisHostType>
344 inline
345 static void
346 getCoeffMatrix_HDIV(OutputViewType &output,
347 const subcellBasisHostType& subcellBasis,
348 const cellBasisHostType& cellBasis,
349 const ordinal_type subcellId,
350 const ordinal_type subcellOrt,
351 const bool inverse = false);
352
353
362 template<typename OutputViewType,
363 typename cellBasisHostType>
364 inline
365 static void
366 getCoeffMatrix_HVOL(OutputViewType &output,
367 const cellBasisHostType& cellBasis,
368 const ordinal_type cellOrt,
369 const bool inverse = false);
370
371 };
372 }
373
377 template<typename DeviceType>
379 public:
380
383 using CoeffMatrixDataViewType = Kokkos::View<double****,DeviceType>;
384
385 //
388 using KeyType = std::pair<const std::string,ordinal_type>;
389 using OrtCoeffDataType = std::map<KeyType,Kokkos::View<double****,DeviceType> >;
390
391 static OrtCoeffDataType ortCoeffData;
392 static OrtCoeffDataType ortInvCoeffData;
393
394 private:
395
396 template<typename BasisHostType>
397 inline
398 static CoeffMatrixDataViewType createCoeffMatrixInternal(const BasisHostType* basis, const bool invTrans = false);
399
400
403 template<typename BasisHostType>
404 inline
406 BasisHostType const *cellBasis,
407 const bool inverse = false);
408
411 template<typename BasisHostType>
412 inline
414 BasisHostType const *cellBasis,
415 const bool inverse = false);
416
419 template<typename BasisHostType>
420 inline
422 BasisHostType const *cellBasis,
423 const bool inverse = false);
424
427 template<typename BasisHostType>
428 inline
430 BasisHostType const *cellBasis,
431 const bool inverse = false);
432
433 public:
434
438 template<typename BasisType>
439 inline
440 static CoeffMatrixDataViewType createCoeffMatrix(const BasisType* basis);
441
445 template<typename BasisType>
446 inline
447 static CoeffMatrixDataViewType createInvCoeffMatrix(const BasisType* basis);
448
451 inline
452 static void clearCoeffMatrix();
453
460 template<typename elemOrtValueType, class ...elemOrtProperties,
461 typename elemNodeValueType, class ...elemNodeProperties>
462 inline
463 static void
464 getOrientation(Kokkos::DynRankView<elemOrtValueType,elemOrtProperties...> elemOrts,
465 const Kokkos::DynRankView<elemNodeValueType,elemNodeProperties...> elemNodes,
466 const shards::CellTopology cellTopo,
467 bool isSide = false);
468
469
477 template<typename outputValueType, class ...outputProperties,
478 typename inputValueType, class ...inputProperties,
479 typename OrientationViewType,
480 typename BasisType>
481 inline
482 static void
483 modifyBasisByOrientation(Kokkos::DynRankView<outputValueType,outputProperties...> output,
484 const Kokkos::DynRankView<inputValueType, inputProperties...> input,
485 const OrientationViewType orts,
486 const BasisType * basis,
487 const bool transpose = false);
488
495 template<typename outputValueType, class ...outputProperties,
496 typename inputValueType, class ...inputProperties,
497 typename OrientationViewType,
498 typename BasisType>
499 inline
500 static void
501 modifyBasisByOrientationTranspose(Kokkos::DynRankView<outputValueType,outputProperties...> output,
502 const Kokkos::DynRankView<inputValueType, inputProperties...> input,
503 const OrientationViewType orts,
504 const BasisType * basis);
505
506
514 template<typename outputValueType, class ...outputProperties,
515 typename inputValueType, class ...inputProperties,
516 typename OrientationViewType,
517 typename BasisType>
518 inline
519 static void
520 modifyBasisByOrientationInverse(Kokkos::DynRankView<outputValueType,outputProperties...> output,
521 const Kokkos::DynRankView<inputValueType, inputProperties...> input,
522 const OrientationViewType orts,
523 const BasisType * basis,
524 const bool transpose = false);
525
526
534 template<typename outputValueType, class ...outputProperties,
535 typename inputValueType, class ...inputProperties,
536 typename OrientationViewType,
537 typename BasisTypeLeft,
538 typename BasisTypeRight>
539 inline
540 static void
541 modifyMatrixByOrientation(Kokkos::DynRankView<outputValueType,outputProperties...> output,
542 const Kokkos::DynRankView<inputValueType, inputProperties...> input,
543 const OrientationViewType orts,
544 const BasisTypeLeft* basisLeft,
545 const BasisTypeRight* basisRight);
546 };
547
548 //Definition of static members
549 template<typename T>
550 typename OrientationTools<T>::OrtCoeffDataType OrientationTools<T>::ortCoeffData;
551
552 template<typename T>
553 typename OrientationTools<T>::OrtCoeffDataType OrientationTools<T>::ortInvCoeffData;
554}
555
556// include templated function definitions
559#include "Intrepid2_OrientationToolsDefCoeffMatrix_HCURL.hpp"
564
565#endif
Header file for the abstract base class Intrepid2::Basis.
Header file for the Intrepid2::Basis_HCURL_HEX_I1_FEM class.
Header file for the Intrepid2::Basis_HCURL_HEX_In_FEM class.
Header file for the Intrepid2::Basis_HCURL_QUAD_I1_FEM class.
Header file for the Intrepid2::Basis_HCURL_QUAD_In_FEM class.
Header file for the Intrepid2::Basis_HCURL_TET_I1_FEM class.
Header file for the Intrepid2::Basis_HCURL_TET_In_FEM class.
Header file for the Intrepid2::Basis_HCURL_TRI_I1_FEM class.
Header file for the Intrepid2::Basis_HCURL_TRI_In_FEM class.
Header file for the Intrepid2::Basis_HCURL_WEDGE_I1_FEM class.
Header file for the Intrepid2::Basis_HDIV_HEX_I1_FEM class.
Header file for the Intrepid2::Basis_HDIV_HEX_In_FEM class.
Header file for the Intrepid2::Basis_HDIV_QUAD_I1_FEM class.
Header file for the Intrepid2::Basis_HDIV_QUAD_In_FEM class.
Header file for the Intrepid2::Basis_HDIV_TET_I1_FEM class.
Header file for the Intrepid2::Basis_HDIV_TET_In_FEM class.
Header file for the Intrepid2::Basis_HDIV_TRI_I1_FEM class.
Header file for the Intrepid2::Basis_HDIV_TRI_In_FEM class.
Header file for the Intrepid2::Basis_HDIV_WEDGE_I1_FEM class.
Header file for the Intrepid2::Basis_HGRAD_HEX_Cn_FEM class.
Header file for the Intrepid2::Basis_HGRAD_LINE_Cn_FEM class.
Header file for the Intrepid2::Basis_HGRAD_QUAD_Cn_FEM class.
Header file for the Intrepid2::Basis_HGRAD_TET_Cn_FEM class.
Header file for the Intrepid2::Basis_HGRAD_TRI_Cn_FEM class.
Header file for the Intrepid2::Basis_HVOL_LINE_Cn_FEM class.
Header file for the Intrepid2::Basis_HVOL_TRI_Cn_FEM class.
Stateless classes that act as factories for two families of hierarchical bases. HierarchicalBasisFami...
Stateless class that acts as a factory for a family of nodal bases (hypercube topologies only at this...
Creation of orientation matrix A of a face or edge for HDIV elements.
Creation of orientation matrix A of a face or edge for HGRAD elements.
Creation of orientation matrix A of a face or edge for HGRAD elements.
Definition file for matrix data in the Intrepid2::OrientationTools class.
Definition file for the Intrepid2::OrientationTools class.
Definition file for functions that modify points due to orientation in the Intrepid2::Impl::Orientati...
Header file for Intrepid2::PointTools class to provide utilities for barycentric coordinates,...
Contains definitions of custom data types in Intrepid2.
Header function for Intrepid2::Util class and other utility functions.
Tools to compute orientations for degrees-of-freedom.
static void getCoeffMatrix_HGRAD(OutputViewType &output, const subcellBasisHostType &subcellBasis, const cellBasisHostType &cellBasis, const ordinal_type subcellId, const ordinal_type subcellOrt, const bool inverse=false)
Compute orientation matrix for HGRAD basis for a given subcell and its reference basis.
static KOKKOS_INLINE_FUNCTION void getJacobianOfOrientationMap(JacobianViewType jacobian, const unsigned cellTopoKey, const ordinal_type cellOrt)
Computes jacobian of the parameterization maps of 1- and 2-subcells with orientation.
static KOKKOS_INLINE_FUNCTION void getModifiedTrianglePoint(ValueType &ot0, ValueType &ot1, const ValueType pt0, const ValueType pt1, const ordinal_type ort)
Computes modified point for triangle.
static void getCoeffMatrix_HDIV(OutputViewType &output, const subcellBasisHostType &subcellBasis, const cellBasisHostType &cellBasis, const ordinal_type subcellId, const ordinal_type subcellOrt, const bool inverse=false)
Compute orientation matrix for HDIV basis for a given subcell and its reference basis.
static KOKKOS_INLINE_FUNCTION void getModifiedQuadrilateralPoint(ValueType &ot0, ValueType &ot1, const ValueType pt0, const ValueType pt1, const ordinal_type ort)
Computes modified point for quadrilateral.
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 KOKKOS_INLINE_FUNCTION void getRefSubcellTangents(TanViewType tangents, const ParamViewType subCellParametrization, const unsigned subcellTopoKey, const ordinal_type subCellOrd, const ordinal_type ort)
Computes the (oriented) subCell tangents.
static KOKKOS_INLINE_FUNCTION void getQuadrilateralJacobian(JacobianViewType jacobian, const ordinal_type ort)
Computes Jacobian of orientation map for quadrilateral.
static KOKKOS_INLINE_FUNCTION void getLineJacobian(JacobianViewType jacobian, const ordinal_type ort)
Computes Jacobian of orientation map for line segment.
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 getCoeffMatrix_HVOL(OutputViewType &output, const cellBasisHostType &cellBasis, const ordinal_type cellOrt, const bool inverse=false)
Compute orientation matrix for HVOL basis for a given (2D or 1D) cell and its reference basis....
static void getCoeffMatrix_HCURL(OutputViewType &output, const subcellBasisHostType &subcellBasis, const cellBasisHostType &cellBasis, const ordinal_type subcellId, const ordinal_type subcellOrt, const bool inverse=false)
Compute orientation matrix for HCURL basis for a given subcell and its reference basis.
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 accounting for orientation.
static KOKKOS_INLINE_FUNCTION void getModifiedLinePoint(ValueType &ot, const ValueType pt, const ordinal_type ort)
Computes modified point for line segment.
static KOKKOS_INLINE_FUNCTION void getTriangleJacobian(JacobianViewType jacobian, const ordinal_type ort)
Computes Jacobian of orientation map for triangle.
Tools to compute orientations for degrees-of-freedom.
static void init_HVOL(CoeffMatrixDataViewType matData, BasisHostType const *cellBasis, const bool inverse=false)
Compute orientation matrix for HVOL basis.
static CoeffMatrixDataViewType createCoeffMatrix(const BasisType *basis)
Create coefficient matrix.
std::pair< const std::string, ordinal_type > KeyType
key :: basis name, order, value :: matrix data view type
static void modifyBasisByOrientationInverse(Kokkos::DynRankView< outputValueType, outputProperties... > output, const Kokkos::DynRankView< inputValueType, inputProperties... > input, const OrientationViewType orts, const BasisType *basis, const bool transpose=false)
Modify basis due to orientation, applying the inverse of the operator applied in modifyBasisByOrienta...
static void init_HCURL(CoeffMatrixDataViewType matData, BasisHostType const *cellBasis, const bool inverse=false)
Compute orientation matrix for HCURL basis.
static CoeffMatrixDataViewType createInvCoeffMatrix(const BasisType *basis)
Create inverse of coefficient matrix.
static void clearCoeffMatrix()
Clear coefficient matrix.
static void modifyBasisByOrientation(Kokkos::DynRankView< outputValueType, outputProperties... > output, const Kokkos::DynRankView< inputValueType, inputProperties... > input, const OrientationViewType orts, const BasisType *basis, const bool transpose=false)
Modify basis due to orientation.
Kokkos::View< double ****, DeviceType > CoeffMatrixDataViewType
subcell ordinal, orientation, matrix m x n
static void init_HDIV(CoeffMatrixDataViewType matData, BasisHostType const *cellBasis, const bool inverse=false)
Compute orientation matrix for HDIV basis.
static void modifyMatrixByOrientation(Kokkos::DynRankView< outputValueType, outputProperties... > output, const Kokkos::DynRankView< inputValueType, inputProperties... > input, const OrientationViewType orts, const BasisTypeLeft *basisLeft, const BasisTypeRight *basisRight)
Modify an assembled (C,F1,F2) matrix according to orientation of the cells.
static void getOrientation(Kokkos::DynRankView< elemOrtValueType, elemOrtProperties... > elemOrts, const Kokkos::DynRankView< elemNodeValueType, elemNodeProperties... > elemNodes, const shards::CellTopology cellTopo, bool isSide=false)
Compute orientations of cells in a workset.
static void init_HGRAD(CoeffMatrixDataViewType matData, BasisHostType const *cellBasis, const bool inverse=false)
Compute orientation matrix for HGRAD basis.
static void modifyBasisByOrientationTranspose(Kokkos::DynRankView< outputValueType, outputProperties... > output, const Kokkos::DynRankView< inputValueType, inputProperties... > input, const OrientationViewType orts, const BasisType *basis)
Modify basis due to orientation, applying the transpose of the operator applied in modifyBasisByOrien...