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  // -----------------------------------------------------------------------------
260  // Coefficient Matrix
261  //
262  //
263 
273  template<typename OutputViewType,
274  typename subcellBasisType,
275  typename cellBasisType>
276  inline
277  static void
278  getCoeffMatrix_HGRAD(OutputViewType &output,
279  const subcellBasisType subcellBasis,
280  const cellBasisType cellBasis,
281  const ordinal_type subcellId,
282  const ordinal_type subcellOrt);
283 
293  template<typename OutputViewType,
294  typename subcellBasisType,
295  typename cellBasisType>
296  inline
297  static void
298  getCoeffMatrix_HCURL(OutputViewType &output,
299  const subcellBasisType subcellBasis,
300  const cellBasisType cellBasis,
301  const ordinal_type subcellId,
302  const ordinal_type subcellOrt);
303 
313  template<typename OutputViewType,
314  typename subcellBasisType,
315  typename cellBasisType>
316  inline
317  static void
318  getCoeffMatrix_HDIV(OutputViewType &output,
319  const subcellBasisType subcellBasis,
320  const cellBasisType cellBasis,
321  const ordinal_type subcellId,
322  const ordinal_type subcellOrt);
323  };
324  }
325 
329  template<typename ExecSpaceType>
331  public:
334  typedef Kokkos::View<double****,ExecSpaceType> CoeffMatrixDataViewType;
335 
336  //
339  static std::map<std::pair<std::string,ordinal_type>,CoeffMatrixDataViewType> ortCoeffData;
340 
341  private:
342 
343  template<typename BasisType>
344  inline
345  static CoeffMatrixDataViewType createCoeffMatrixInternal(const BasisType* basis);
346 
347  //
348  // High order elements transformation matrices for HGRAD
349  //
350 
351  template<typename BasisType>
352  inline
353  static void init_HGRAD_QUAD(CoeffMatrixDataViewType matData,
354  BasisType const *cellBasis);
355 
356  template<typename BasisType>
357  inline
358  static void init_HGRAD_HEX(CoeffMatrixDataViewType matData,
359  BasisType const *cellBasis);
360 
361  template<typename BasisType>
362  inline
363  static void init_HGRAD_TRI(CoeffMatrixDataViewType matData,
364  BasisType const *cellBasis);
365 
366  template<typename BasisType>
367  inline
368  static void init_HGRAD_TET(CoeffMatrixDataViewType matData,
369  BasisType const *cellBasis);
370 
371  //
372  // High order elements transformation matrices for HCURL
373  //
374 
375  template<typename BasisType>
376  inline
377  static void init_HCURL_QUAD(CoeffMatrixDataViewType matData,
378  BasisType const *cellBasis);
379 
380  template<typename BasisType>
381  inline
382  static void init_HCURL_HEX(CoeffMatrixDataViewType matData,
383  BasisType const *cellBasis);
384 
385  template<typename BasisType>
386  inline
387  static void init_HCURL_TRI(CoeffMatrixDataViewType matData,
388  BasisType const *cellBasis);
389 
390  template<typename BasisType>
391  inline
392  static void init_HCURL_TET(CoeffMatrixDataViewType matData,
393  BasisType const *cellBasis);
394 
395  //
396  // High order elements transformation matrices for HDIV
397  //
398 
399  template<typename BasisType>
400  inline
401  static void init_HDIV_QUAD(CoeffMatrixDataViewType matData,
402  BasisType const *cellBasis);
403 
404 
405  template<typename BasisType>
406  inline
407  static void init_HDIV_HEX(CoeffMatrixDataViewType matData,
408  BasisType const *cellBasis);
409 
410  template<typename BasisType>
411  inline
412  static void init_HDIV_TRI(CoeffMatrixDataViewType matData,
413  BasisType const *cellBasis);
414 
415  template<typename BasisType>
416  inline
417  static void init_HDIV_TET(CoeffMatrixDataViewType matData,
418  BasisType const *cellBasis);
419 
420  //
421  // I1 element specialization
422  //
423 
424  inline
425  static void init_EDGE_ELEMENT_I1_FEM(CoeffMatrixDataViewType matData,
426  const ordinal_type edgeId);
427 
428  inline
429  static void init_TRI_FACE_ELEMENT_I1_FEM(CoeffMatrixDataViewType matData,
430  const ordinal_type faceId);
431 
432  inline
433  static void init_QUAD_FACE_ELEMENT_I1_FEM(CoeffMatrixDataViewType matData,
434  const ordinal_type faceId);
435 
436  public:
437 
441  template<typename BasisType>
442  inline
443  static CoeffMatrixDataViewType createCoeffMatrix(const BasisType* basis);
444 
447  inline
448  static void clearCoeffMatrix();
449 
453  template<typename ptViewType>
454  KOKKOS_INLINE_FUNCTION
455  static bool
456  isLeftHandedCell(const ptViewType pts);
457 
463  template<typename elemOrtValueType, class ...elemOrtProperties,
464  typename elemNodeValueType, class ...elemNodeProperties>
465  inline
466  static void
467  getOrientation( Kokkos::DynRankView<elemOrtValueType,elemOrtProperties...> elemOrts,
468  const Kokkos::DynRankView<elemNodeValueType,elemNodeProperties...> elemNodes,
469  const shards::CellTopology cellTopo);
470 
477  template<typename outputValueType, class ...outputProperties,
478  typename inputValueType, class ...inputProperties,
479  typename ortValueType, class ...ortProperties,
480  typename BasisType>
481  inline
482  static void
483  modifyBasisByOrientation( Kokkos::DynRankView<outputValueType,outputProperties...> output,
484  const Kokkos::DynRankView<inputValueType, inputProperties...> input,
485  const Kokkos::DynRankView<ortValueType, ortProperties...> orts,
486  const BasisType* basis);
487 
488  template<typename ExecutionSpaceType, typename dofCoordsValueType,
489  typename dofCoeffsValueType, typename BasisType>
490  inline
491  static void
492  getSubCellBases(
495  const BasisType& basis,
496  EPointType pointType = POINTTYPE_EQUISPACED);
497 
498  template<typename dofCoordsValueType, class ...dofCoordsProperties,
499  typename dofCoeffsValueType, class ...dofCoeffsProperties,
500  typename ortValueType, class ...ortProperties,
501  typename BasisType>
502  inline
503  static void
504  getModifiedDofs( Kokkos::DynRankView<dofCoordsValueType,dofCoordsProperties...> orientedDofCoords,
505  Kokkos::DynRankView<dofCoeffsValueType,dofCoeffsProperties...> orientedDofCoeffs,
506  const Kokkos::DynRankView<ortValueType, ortProperties...> orts,
507  const BasisType& basis,
508  EPointType pointType = POINTTYPE_EQUISPACED);
509 
510  };
511 
512  template<typename T>
513  std::map<std::pair<std::string,ordinal_type>, typename OrientationTools<T>::CoeffMatrixDataViewType>
515 }
516 
517 // include templated function definitions
524 
525 #endif
526 
527 
528 
529 
530 
531  // class CoeffMatrix {
532  // private:
533  // ordinal_type _m, _n;
534 
535  // Kokkos::View<size_type*,ExecSpaceType> _ap; //!< pointers to column index and values
536  // Kokkos::View<ordinal_type*,ExecSpaceType> _aj; //!< column index compressed format
537  // Kokkos::View<double*,ExecSpaceType> _ax; //!< values
538 
539  // inline
540  // void createInternalArrays(const ordinal_type m,
541  // const ordinal_type n,
542  // const size_type nnz) {
543  // _m = m;
544  // _n = n;
545 
546  // _ap = Kokkos::View<size_type*, SpT>("OrientationTools::CoeffMatrix::RowPtrArray", m+1);
547  // _aj = Kokkos::View<ordinal_type*,SpT>("OrientationTools::CoeffMatrix::ColsArray", nnz);
548  // _ax = Kokkos::View<double*, SpT>("OrientationTools::CoeffMAtrix::ValuesArray", nnz);
549  // }
550 
551  // public:
552  // KOKKOS_INLINE_FUNCTION
553  // CoeffMatrix()
554  // : _m(0), _n(0), _ap(), _aj(), _ax() { }
555 
556  // KOKKOS_INLINE_FUNCTION
557  // CoeffMatrix(const CoeffMatrix &b) = default;
558 
559  // KOKKOS_INLINE_FUNCTION
560  // ordinal_type NumRows() const {
561  // return _m;
562  // }
563 
564  // KOKKOS_INLINE_FUNCTION
565  // ordinal_type NumCols() const {
566  // return _n;
567  // }
568 
569  // KOKKOS_INLINE_FUNCTION
570  // size_type RowPtr(const ordinal_type i) const {
571  // return _ap(i);
572  // }
573 
574  // KOKKOS_INLINE_FUNCTION
575  // Kokkos::View<ordinal_type*,ExecSpaceType> ColsInRow(const ordinal_type i) const {
576  // return Kokkos::subview(_aj, Kokkos::pair<ordinal_type,ordinal_type>(_ap(i), _ap(i+1)));
577  // }
578 
579  // KOKKOS_INLINE_FUNCTION
580  // Kokkos::View<double*,ExecSpaceType> ValuesInRow(const ordinal_type i) const {
581  // return Kokkos::subview(_ax, Kokkos::pair<ordinal_type,ordinal_type>(_ap(i), _ap(i+1)));
582  // }
583 
584  // KOKKOS_INLINE_FUNCTION
585  // ordinal_type NumNonZerosInRow(const ordinal_type i) const {
586  // return (_ap(i+1) - _ap(i));
587  // }
588  // };
589 
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_HDIV(OutputViewType &output, const subcellBasisType subcellBasis, const cellBasisType cellBasis, const ordinal_type subcellId, const ordinal_type subcellOrt)
Compute coefficient matrix for HDIV by collocating point values.
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 void getCoeffMatrix_HCURL(OutputViewType &output, const subcellBasisType subcellBasis, const cellBasisType cellBasis, const ordinal_type subcellId, const ordinal_type subcellOrt)
Compute coefficient matrix for HCURL by collocating point values.
Definition file for the Intrepid2::OrientationTools class.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
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 KOKKOS_INLINE_FUNCTION bool isLeftHandedCell(const ptViewType pts)
Check if left-handed. If an element is alinged left, it is an error.
static void modifyBasisByOrientation(Kokkos::DynRankView< outputValueType, outputProperties...> output, const Kokkos::DynRankView< inputValueType, inputProperties...> input, const Kokkos::DynRankView< ortValueType, ortProperties...> orts, const BasisType *basis)
Modify basis due to orientation.
Kokkos::View< double ****, ExecSpaceType > CoeffMatrixDataViewType
subcell ordinal, orientation, matrix m x n
Creation of orientation matrix A of a face or edge for HDIV elements.
Header function for Intrepid2::Util class and other utility functions.
static void clearCoeffMatrix()
Clear coefficient matrix.
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 getOrientation(Kokkos::DynRankView< elemOrtValueType, elemOrtProperties...> elemOrts, const Kokkos::DynRankView< elemNodeValueType, elemNodeProperties...> elemNodes, const shards::CellTopology cellTopo)
Compute orientations of cells in a workset.
Header file for the Intrepid2::Basis_HDIV_TRI_I1_FEM class.
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.
Creation of orientation matrix A of a face or edge for HCURL elements.
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.
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.
Definition file for functions that modify points due to orientation in the Intrepid2::Impl::Orientati...
static void getCoeffMatrix_HGRAD(OutputViewType &output, const subcellBasisType subcellBasis, const cellBasisType cellBasis, const ordinal_type subcellId, const ordinal_type subcellOrt)
Compute coefficient matrix for HGRAD by collocating point values.
Header file for the Intrepid2::Basis_HCURL_QUAD_I1_FEM class.
Header file for the Intrepid2::Basis_HCURL_WEDGE_I1_FEM class.
Tools to compute orientations for degrees-of-freedom.
Header file for the Intrepid2::Basis_HCURL_TET_I1_FEM class.
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.
Header file for the Intrepid2::Basis_HCURL_TRI_I1_FEM class.
Header file for the Intrepid2::Basis_HVOL_TRI_Cn_FEM class.
static CoeffMatrixDataViewType createCoeffMatrix(const BasisType *basis)
Create 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.