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
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
107namespace 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
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
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
328 template<typename OutputViewType,
329 typename subcellBasisHostType,
330 typename cellBasisHostType>
331 inline
332 static void
333 getCoeffMatrix_HGRAD(OutputViewType &output,
334 const subcellBasisHostType& subcellBasis,
335 const cellBasisHostType& cellBasis,
336 const ordinal_type subcellId,
337 const ordinal_type subcellOrt);
338
349 template<typename OutputViewType,
350 typename subcellBasisHostType,
351 typename cellBasisHostType>
352 inline
353 static void
354 getCoeffMatrix_HCURL(OutputViewType &output,
355 const subcellBasisHostType& subcellBasis,
356 const cellBasisHostType& cellBasis,
357 const ordinal_type subcellId,
358 const ordinal_type subcellOrt);
359
360
371 template<typename OutputViewType,
372 typename subcellBasisHostType,
373 typename cellBasisHostType>
374 inline
375 static void
376 getCoeffMatrix_HDIV(OutputViewType &output,
377 const subcellBasisHostType& subcellBasis,
378 const cellBasisHostType& cellBasis,
379 const ordinal_type subcellId,
380 const ordinal_type subcellOrt);
381 };
382 }
383
387 template<typename DeviceType>
389 public:
390
393 typedef Kokkos::View<double****,DeviceType> CoeffMatrixDataViewType;
394
395 //
398 static std::map<std::pair<std::string,ordinal_type>,CoeffMatrixDataViewType> ortCoeffData;
399
400 private:
401
402 template<typename BasisHostType>
403 inline
404 static CoeffMatrixDataViewType createCoeffMatrixInternal(const BasisHostType* basis);
405
406
409 template<typename BasisHostType>
410 inline
412 BasisHostType const *cellBasis);
413
416 template<typename BasisHostType>
417 inline
419 BasisHostType const *cellBasis);
420
423 template<typename BasisHostType>
424 inline
426 BasisHostType const *cellBasis);
427
428 public:
429
433 template<typename BasisType>
434 inline
435 static CoeffMatrixDataViewType createCoeffMatrix(const BasisType* basis);
436
439 inline
440 static void clearCoeffMatrix();
441
447 template<typename elemOrtValueType, class ...elemOrtProperties,
448 typename elemNodeValueType, class ...elemNodeProperties>
449 inline
450 static void
451 getOrientation(Kokkos::DynRankView<elemOrtValueType,elemOrtProperties...> elemOrts,
452 const Kokkos::DynRankView<elemNodeValueType,elemNodeProperties...> elemNodes,
453 const shards::CellTopology cellTopo);
454
461 template<typename outputValueType, class ...outputProperties,
462 typename inputValueType, class ...inputProperties,
463 typename OrientationViewType,
464 typename BasisType>
465 inline
466 static void
467 modifyBasisByOrientation(Kokkos::DynRankView<outputValueType,outputProperties...> output,
468 const Kokkos::DynRankView<inputValueType, inputProperties...> input,
469 const OrientationViewType orts,
470 const BasisType * basis);
471 };
472
473 template<typename T>
474 std::map<std::pair<std::string,ordinal_type>, typename OrientationTools<T>::CoeffMatrixDataViewType>
476}
477
478// include templated function definitions
481#include "Intrepid2_OrientationToolsDefCoeffMatrix_HCURL.hpp"
485
486#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.
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 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 void getCoeffMatrix_HGRAD(OutputViewType &output, const subcellBasisHostType &subcellBasis, const cellBasisHostType &cellBasis, const ordinal_type subcellId, const ordinal_type subcellOrt)
Compute orientation matrix for HGRAD basis for a given subcell and its reference basis.
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 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_HCURL(OutputViewType &output, const subcellBasisHostType &subcellBasis, const cellBasisHostType &cellBasis, const ordinal_type subcellId, const ordinal_type subcellOrt)
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.
static KOKKOS_INLINE_FUNCTION void getModifiedLinePoint(ValueType &ot, const ValueType pt, const ordinal_type ort)
Computes modified point for line segment.
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 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 CoeffMatrixDataViewType createCoeffMatrix(const BasisType *basis)
Create coefficient matrix.
static std::map< std::pair< std::string, ordinal_type >, CoeffMatrixDataViewType > ortCoeffData
key :: basis name, order, value :: matrix data view type
static void init_HDIV(CoeffMatrixDataViewType matData, BasisHostType const *cellBasis)
Compute orientation matrix for HDIV basis.
Kokkos::View< double ****, DeviceType > CoeffMatrixDataViewType
subcell ordinal, orientation, matrix m x n
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)
Modify basis due to orientation.
static void init_HCURL(CoeffMatrixDataViewType matData, BasisHostType const *cellBasis)
Compute orientation matrix for HCURL basis.
static void init_HGRAD(CoeffMatrixDataViewType matData, BasisHostType const *cellBasis)
Compute orientation matrix for HGRAD basis.
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.