Intrepid2
Intrepid2_OrientationTools.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 
48 #ifndef __INTREPID2_ORIENTATIONTOOLS_HPP__
49 #define __INTREPID2_ORIENTATIONTOOLS_HPP__
50 
51 #include "Intrepid2_ConfigDefs.hpp"
52 #include "Intrepid2_Types.hpp"
53 #include "Intrepid2_Utils.hpp"
54 
55 #include "Shards_CellTopology.hpp"
56 #include "Shards_BasicTopologies.hpp"
57 
58 #include "Intrepid2_PointTools.hpp"
59 
60 #include "Intrepid2_Basis.hpp"
61 
62 // -- HGRAD family
66 
69 
70 // -- HCURL family
73 
77 
78 // -- HDIV family
81 
85 
86 // -- Lower order family
89 
92 
96 
100 
103 
104 #include "Teuchos_LAPACK.hpp"
105 
106 
107 namespace Intrepid2 {
108 
109  namespace Impl {
110 
115  public:
116 
117  // -----------------------------------------------------------------------------
118  // Point modification
119  //
120  //
121 
128  template<typename ValueType>
129  KOKKOS_INLINE_FUNCTION
130  static void
131  getModifiedLinePoint(ValueType &ot,
132  const ValueType pt,
133  const ordinal_type ort);
134 
143  template<typename ValueType>
144  KOKKOS_INLINE_FUNCTION
145  static void
146  getModifiedTrianglePoint(ValueType &ot0,
147  ValueType &ot1,
148  const ValueType pt0,
149  const ValueType pt1,
150  const ordinal_type ort);
151 
160  template<typename ValueType>
161  KOKKOS_INLINE_FUNCTION
162  static void
163  getModifiedQuadrilateralPoint(ValueType &ot0,
164  ValueType &ot1,
165  const ValueType pt0,
166  const ValueType pt1,
167  const ordinal_type ort);
168 
176  template<typename outPointViewType,
177  typename refPointViewType>
178  inline
179  static void
180  mapToModifiedReference(outPointViewType outPoints,
181  const refPointViewType refPoints,
182  const shards::CellTopology cellTopo,
183  const ordinal_type cellOrt = 0);
184 
192  template<typename outPointViewType,
193  typename refPointViewType>
194  KOKKOS_INLINE_FUNCTION
195  static void
196  mapToModifiedReference(outPointViewType outPoints,
197  const refPointViewType refPoints,
198  const unsigned cellTopoKey,
199  const ordinal_type cellOrt = 0);
200 
201 
207  template<typename JacobianViewType>
208  KOKKOS_INLINE_FUNCTION
209  static void
210  getLineJacobian(JacobianViewType jacobian, const ordinal_type ort);
211 
217  template<typename JacobianViewType>
218  KOKKOS_INLINE_FUNCTION
219  static void
220  getTriangleJacobian(JacobianViewType jacobian, const ordinal_type ort);
221 
227  template<typename JacobianViewType>
228  KOKKOS_INLINE_FUNCTION
229  static void
230  getQuadrilateralJacobian(JacobianViewType jacobian, const ordinal_type ort);
231 
232 
239  template<typename JacobianViewType>
240  inline
241  static void
242  getJacobianOfOrientationMap(JacobianViewType jacobian,
243  const shards::CellTopology cellTopo,
244  const ordinal_type cellOrt);
245 
252  template<typename JacobianViewType>
253  KOKKOS_INLINE_FUNCTION
254  static void
255  getJacobianOfOrientationMap(JacobianViewType jacobian,
256  const unsigned cellTopoKey,
257  const ordinal_type cellOrt);
258 
259 
267  template<typename TanViewType, typename ParamViewType>
268  KOKKOS_INLINE_FUNCTION
269  static void getRefSubcellTangents(TanViewType tangents,
270  const ParamViewType subCellParametrization,
271  const unsigned subcellTopoKey,
272  const ordinal_type subCellOrd,
273  const ordinal_type ort);
274 
275 
286  template<typename TanNormViewType, typename ParamViewType>
287  KOKKOS_INLINE_FUNCTION
288  static void getRefSideTangentsAndNormal(TanNormViewType tangentsAndNormal,
289  const ParamViewType subCellParametrization,
290  const unsigned subcellTopoKey,
291  const ordinal_type subCellOrd,
292  const ordinal_type ort);
293 
294 
304  template<typename coordsViewType, typename subcellCoordsViewType, typename ParamViewType>
305  KOKKOS_INLINE_FUNCTION
306  static void mapSubcellCoordsToRefCell(coordsViewType cellCoords,
307  const subcellCoordsViewType subCellCoords,
308  const ParamViewType subcellParametrization,
309  const unsigned subcellTopoKey,
310  const ordinal_type subCellOrd,
311  const ordinal_type ort);
312 
313  // -----------------------------------------------------------------------------
314  // Coefficient Matrix
315  //
316  //
317 
329  template<typename OutputViewType,
330  typename subcellBasisHostType,
331  typename cellBasisHostType>
332  inline
333  static void
334  getCoeffMatrix_HGRAD(OutputViewType &output,
335  const subcellBasisHostType& subcellBasis,
336  const cellBasisHostType& cellBasis,
337  const ordinal_type subcellId,
338  const ordinal_type subcellOrt,
339  const bool inverse = false);
340 
351  template<typename OutputViewType,
352  typename subcellBasisHostType,
353  typename cellBasisHostType>
354  inline
355  static void
356  getCoeffMatrix_HCURL(OutputViewType &output,
357  const subcellBasisHostType& subcellBasis,
358  const cellBasisHostType& cellBasis,
359  const ordinal_type subcellId,
360  const ordinal_type subcellOrt,
361  const bool inverse = false);
362 
363 
374  template<typename OutputViewType,
375  typename subcellBasisHostType,
376  typename cellBasisHostType>
377  inline
378  static void
379  getCoeffMatrix_HDIV(OutputViewType &output,
380  const subcellBasisHostType& subcellBasis,
381  const cellBasisHostType& cellBasis,
382  const ordinal_type subcellId,
383  const ordinal_type subcellOrt,
384  const bool inverse = false);
385 
386 
395  template<typename OutputViewType,
396  typename cellBasisHostType>
397  inline
398  static void
399  getCoeffMatrix_HVOL(OutputViewType &output,
400  const cellBasisHostType& cellBasis,
401  const ordinal_type cellOrt,
402  const bool inverse = false);
403 
404  };
405  }
406 
410  template<typename DeviceType>
412  public:
413 
416  typedef Kokkos::View<double****,DeviceType> CoeffMatrixDataViewType;
417 
418  //
421  static std::map<std::pair<std::string,ordinal_type>,CoeffMatrixDataViewType> ortCoeffData;
422  static std::map<std::pair<std::string,ordinal_type>,CoeffMatrixDataViewType> ortInvCoeffData;
423 
424  private:
425 
426  template<typename BasisHostType>
427  inline
428  static CoeffMatrixDataViewType createCoeffMatrixInternal(const BasisHostType* basis, const bool invTrans = false);
429 
430 
433  template<typename BasisHostType>
434  inline
435  static void init_HGRAD(CoeffMatrixDataViewType matData,
436  BasisHostType const *cellBasis,
437  const bool inverse = false);
438 
441  template<typename BasisHostType>
442  inline
443  static void init_HCURL(CoeffMatrixDataViewType matData,
444  BasisHostType const *cellBasis,
445  const bool inverse = false);
446 
449  template<typename BasisHostType>
450  inline
451  static void init_HDIV(CoeffMatrixDataViewType matData,
452  BasisHostType const *cellBasis,
453  const bool inverse = false);
454 
457  template<typename BasisHostType>
458  inline
459  static void init_HVOL(CoeffMatrixDataViewType matData,
460  BasisHostType const *cellBasis,
461  const bool inverse = false);
462 
463  public:
464 
468  template<typename BasisType>
469  inline
470  static CoeffMatrixDataViewType createCoeffMatrix(const BasisType* basis);
471 
475  template<typename BasisType>
476  inline
477  static CoeffMatrixDataViewType createInvCoeffMatrix(const BasisType* basis);
478 
481  inline
482  static void clearCoeffMatrix();
483 
490  template<typename elemOrtValueType, class ...elemOrtProperties,
491  typename elemNodeValueType, class ...elemNodeProperties>
492  inline
493  static void
494  getOrientation(Kokkos::DynRankView<elemOrtValueType,elemOrtProperties...> elemOrts,
495  const Kokkos::DynRankView<elemNodeValueType,elemNodeProperties...> elemNodes,
496  const shards::CellTopology cellTopo,
497  bool isSide = false);
498 
499 
507  template<typename outputValueType, class ...outputProperties,
508  typename inputValueType, class ...inputProperties,
509  typename OrientationViewType,
510  typename BasisType>
511  inline
512  static void
513  modifyBasisByOrientation(Kokkos::DynRankView<outputValueType,outputProperties...> output,
514  const Kokkos::DynRankView<inputValueType, inputProperties...> input,
515  const OrientationViewType orts,
516  const BasisType * basis,
517  const bool transpose = false);
518 
525  template<typename outputValueType, class ...outputProperties,
526  typename inputValueType, class ...inputProperties,
527  typename OrientationViewType,
528  typename BasisType>
529  inline
530  static void
531  modifyBasisByOrientationTranspose(Kokkos::DynRankView<outputValueType,outputProperties...> output,
532  const Kokkos::DynRankView<inputValueType, inputProperties...> input,
533  const OrientationViewType orts,
534  const BasisType * basis);
535 
536 
544  template<typename outputValueType, class ...outputProperties,
545  typename inputValueType, class ...inputProperties,
546  typename OrientationViewType,
547  typename BasisType>
548  inline
549  static void
550  modifyBasisByOrientationInverse(Kokkos::DynRankView<outputValueType,outputProperties...> output,
551  const Kokkos::DynRankView<inputValueType, inputProperties...> input,
552  const OrientationViewType orts,
553  const BasisType * basis,
554  const bool transpose = false);
555 
556 
564  template<typename outputValueType, class ...outputProperties,
565  typename inputValueType, class ...inputProperties,
566  typename OrientationViewType,
567  typename BasisTypeLeft,
568  typename BasisTypeRight>
569  inline
570  static void
571  modifyMatrixByOrientation(Kokkos::DynRankView<outputValueType,outputProperties...> output,
572  const Kokkos::DynRankView<inputValueType, inputProperties...> input,
573  const OrientationViewType orts,
574  const BasisTypeLeft* basisLeft,
575  const BasisTypeRight* basisRight);
576  };
577 
578  template<typename T>
579  std::map<std::pair<std::string,ordinal_type>, typename OrientationTools<T>::CoeffMatrixDataViewType>
581 
582  template<typename T>
583  std::map<std::pair<std::string,ordinal_type>, typename OrientationTools<T>::CoeffMatrixDataViewType>
585 }
586 
587 // include templated function definitions
590 #include "Intrepid2_OrientationToolsDefCoeffMatrix_HCURL.hpp"
595 
596 #endif
static KOKKOS_INLINE_FUNCTION void getModifiedLinePoint(ValueType &ot, const ValueType pt, const ordinal_type ort)
Computes modified point for line segment.
Stateless classes that act as factories for two families of hierarchical bases. HierarchicalBasisFami...
Header file for the Intrepid2::Basis_HGRAD_LINE_Cn_FEM class.
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.
Header file for the Intrepid2::Basis_HDIV_TET_I1_FEM class.
Header file for the Intrepid2::Basis_HCURL_HEX_In_FEM class.
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.
Header file for the Intrepid2::Basis_HDIV_HEX_I1_FEM class.
Header file for the Intrepid2::Basis_HDIV_WEDGE_I1_FEM class.
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.
Kokkos::View< double ****, DeviceType > CoeffMatrixDataViewType
subcell ordinal, orientation, matrix m x n
Definition file for the Intrepid2::OrientationTools class.
static std::map< std::pair< std::string, ordinal_type >, CoeffMatrixDataViewType > ortCoeffData
key :: basis name, order, value :: matrix data view type
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.
Header file for the Intrepid2::Basis_HDIV_HEX_In_FEM class.
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.
Creation of orientation matrix A of a face or edge for HDIV elements.
static void clearCoeffMatrix()
Clear coefficient matrix.
Header function for Intrepid2::Util class and other utility functions.
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 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.
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_HDIV_QUAD_I1_FEM class.
Header file for the Intrepid2::Basis_HCURL_TET_In_FEM class.
Header file for the Intrepid2::Basis_HGRAD_TET_Cn_FEM class.
static void init_HDIV(CoeffMatrixDataViewType matData, BasisHostType const *cellBasis, const bool inverse=false)
Compute orientation matrix for HDIV basis.
Header file for the Intrepid2::Basis_HDIV_TRI_I1_FEM class.
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.
Header file for the Intrepid2::Basis_HDIV_QUAD_In_FEM class.
Header file for the Intrepid2::Basis_HCURL_QUAD_In_FEM class.
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 KOKKOS_INLINE_FUNCTION void getTriangleJacobian(JacobianViewType jacobian, const ordinal_type ort)
Computes Jacobian of orientation map for triangle.
static KOKKOS_INLINE_FUNCTION void getLineJacobian(JacobianViewType jacobian, const ordinal_type ort)
Computes Jacobian of orientation map for line segment.
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...
Header file for the Intrepid2::Basis_HDIV_TRI_In_FEM class.
Contains definitions of custom data types in Intrepid2.
Header file for the Intrepid2::Basis_HGRAD_QUAD_Cn_FEM class.
Creation of orientation matrix A of a face or edge for HGRAD elements.
Header file for the Intrepid2::Basis_HDIV_TET_In_FEM class.
Creation of orientation matrix A of a face or edge for HGRAD elements.
static void init_HGRAD(CoeffMatrixDataViewType matData, BasisHostType const *cellBasis, const bool inverse=false)
Compute orientation matrix for HGRAD basis.
Definition file for functions that modify points due to orientation in the Intrepid2::Impl::Orientati...
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.
Header file for the Intrepid2::Basis_HCURL_QUAD_I1_FEM class.
Header file for the Intrepid2::Basis_HCURL_WEDGE_I1_FEM class.
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 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...
Tools to compute orientations for degrees-of-freedom.
Header file for the Intrepid2::Basis_HCURL_TET_I1_FEM class.
static CoeffMatrixDataViewType createCoeffMatrix(const BasisType *basis)
Create coefficient matrix.
static KOKKOS_INLINE_FUNCTION void getQuadrilateralJacobian(JacobianViewType jacobian, const ordinal_type ort)
Computes Jacobian of orientation map for quadrilateral.
Header file for the Intrepid2::Basis_HCURL_HEX_I1_FEM class.
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.
Definition file for matrix data in the Intrepid2::OrientationTools class.
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...
Header file for the Intrepid2::Basis_HCURL_TRI_I1_FEM class.
static void init_HCURL(CoeffMatrixDataViewType matData, BasisHostType const *cellBasis, const bool inverse=false)
Compute orientation matrix for HCURL basis.
Header file for the Intrepid2::Basis_HVOL_TRI_Cn_FEM class.
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_HVOL(CoeffMatrixDataViewType matData, BasisHostType const *cellBasis, const bool inverse=false)
Compute orientation matrix for HVOL basis.
static CoeffMatrixDataViewType createInvCoeffMatrix(const BasisType *basis)
Create inverse of coefficient matrix.
Stateless class that acts as a factory for a family of nodal bases (hypercube topologies only at this...
Header file for the abstract base class Intrepid2::Basis.
Header file for the Intrepid2::Basis_HCURL_TRI_In_FEM class.
Header file for Intrepid2::PointTools class to provide utilities for barycentric coordinates, equispaced lattices, and warp-blend point distrubtions.
Tools to compute orientations for degrees-of-freedom.
Header file for the Intrepid2::Basis_HGRAD_HEX_Cn_FEM class.