Intrepid2
Intrepid2_CellTools.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
16#ifndef __INTREPID2_CELLTOOLS_HPP__
17#define __INTREPID2_CELLTOOLS_HPP__
18
19#include "Intrepid2_ConfigDefs.hpp"
20
21#include "Shards_CellTopology.hpp"
22#include "Shards_BasicTopologies.hpp"
23
24#include "Teuchos_RCP.hpp"
25
26#include "Intrepid2_Types.hpp"
27#include "Intrepid2_Utils.hpp"
28
30
31#include "Intrepid2_Basis.hpp"
32
35
38
41
45
48
51
54
55#include "Intrepid2_Data.hpp"
57
58namespace Intrepid2 {
59
60 //============================================================================================//
61 // //
62 // CellTools //
63 // //
64 //============================================================================================//
65
75 namespace FunctorCellTools {
76
79 template<typename OutputViewType,
80 typename InputViewType,
81 typename ImplBasis,
82 EOperator Operator>
83 struct F_getValues {
84 OutputViewType output_;
85 const InputViewType input_;
86
87 KOKKOS_INLINE_FUNCTION
88 F_getValues(OutputViewType output,
89 const InputViewType input)
90 : output_(output),
91 input_(input) {}
92
93 KOKKOS_INLINE_FUNCTION
94 void
95 operator()(const ordinal_type cl, const ordinal_type pt) const {
96 auto in = Kokkos::subview(input_, cl, pt, Kokkos::ALL()); // D
97 if constexpr(Operator == OPERATOR_VALUE) {
98 auto out = Kokkos::subview(output_, cl, Kokkos::ALL(), pt); // N
99 ImplBasis::template Serial<OPERATOR_VALUE>::getValues(out, in);
100 }
101 else if constexpr (Operator == OPERATOR_GRAD) {
102 auto out = Kokkos::subview(output_, cl, Kokkos::ALL(), pt, Kokkos::ALL()); // N, D
103 ImplBasis::template Serial<OPERATOR_GRAD>::getValues(out, in);
104 } else {
105 static_assert((Operator != OPERATOR_VALUE) && (Operator != OPERATOR_GRAD), "Invalid template parameter. Only OPERATOR_VALUE and OPERATOR_GRAD are supported.");
106 }
107 }
108 };
109 }
110
111 template<typename DeviceType>
112 class CellTools {
113 using ExecSpaceType = typename DeviceType::execution_space;
114 using MemSpaceType = typename DeviceType::memory_space;
115 public:
116
121 inline
122 static bool
123 hasReferenceCell( const shards::CellTopology cellTopo ) {
125 }
126
130 template<typename outputValueType,
131 typename pointValueType>
132 static Teuchos::RCP<Basis<DeviceType,outputValueType,pointValueType> >
133 createHGradBasis( const shards::CellTopology cellTopo ) {
134 Teuchos::RCP<Basis<DeviceType,outputValueType,pointValueType> > r_val;
135
136 switch (cellTopo.getKey()) {
137 case shards::Line<2>::key: r_val = Teuchos::rcp(new Basis_HGRAD_LINE_C1_FEM <DeviceType,outputValueType,pointValueType>()); break;
138 case shards::Triangle<3>::key: r_val = Teuchos::rcp(new Basis_HGRAD_TRI_C1_FEM <DeviceType,outputValueType,pointValueType>()); break;
139 case shards::Quadrilateral<4>::key: r_val = Teuchos::rcp(new Basis_HGRAD_QUAD_C1_FEM <DeviceType,outputValueType,pointValueType>()); break;
140 case shards::Tetrahedron<4>::key: r_val = Teuchos::rcp(new Basis_HGRAD_TET_C1_FEM <DeviceType,outputValueType,pointValueType>()); break;
141 case shards::Hexahedron<8>::key: r_val = Teuchos::rcp(new Basis_HGRAD_HEX_C1_FEM <DeviceType,outputValueType,pointValueType>()); break;
142 case shards::Wedge<6>::key: r_val = Teuchos::rcp(new Basis_HGRAD_WEDGE_C1_FEM <DeviceType,outputValueType,pointValueType>()); break;
143 case shards::Pyramid<5>::key: r_val = Teuchos::rcp(new Basis_HGRAD_PYR_C1_FEM <DeviceType,outputValueType,pointValueType>()); break;
144
145 case shards::Line<3>::key: r_val = Teuchos::rcp(new Basis_HGRAD_LINE_C2_FEM <DeviceType,outputValueType,pointValueType>()); break;
146 case shards::Triangle<6>::key: r_val = Teuchos::rcp(new Basis_HGRAD_TRI_C2_FEM <DeviceType,outputValueType,pointValueType>()); break;
147 case shards::Quadrilateral<8>::key: r_val = Teuchos::rcp(new Basis_HGRAD_QUAD_I2_FEM <DeviceType,outputValueType,pointValueType>()); break;
148 case shards::Quadrilateral<9>::key: r_val = Teuchos::rcp(new Basis_HGRAD_QUAD_C2_FEM <DeviceType,outputValueType,pointValueType>()); break;
149 case shards::Tetrahedron<10>::key: r_val = Teuchos::rcp(new Basis_HGRAD_TET_C2_FEM <DeviceType,outputValueType,pointValueType>()); break;
150 case shards::Tetrahedron<11>::key: r_val = Teuchos::rcp(new Basis_HGRAD_TET_COMP12_FEM<DeviceType,outputValueType,pointValueType>()); break;
151 case shards::Hexahedron<20>::key: r_val = Teuchos::rcp(new Basis_HGRAD_HEX_I2_FEM <DeviceType,outputValueType,pointValueType>()); break;
152 case shards::Hexahedron<27>::key: r_val = Teuchos::rcp(new Basis_HGRAD_HEX_C2_FEM <DeviceType,outputValueType,pointValueType>()); break;
153 case shards::Wedge<15>::key: r_val = Teuchos::rcp(new Basis_HGRAD_WEDGE_I2_FEM <DeviceType,outputValueType,pointValueType>()); break;
154 case shards::Wedge<18>::key: r_val = Teuchos::rcp(new Basis_HGRAD_WEDGE_C2_FEM <DeviceType,outputValueType,pointValueType>()); break;
155 case shards::Pyramid<13>::key: r_val = Teuchos::rcp(new Basis_HGRAD_PYR_I2_FEM <DeviceType,outputValueType,pointValueType>()); break;
156
157 case shards::Beam<2>::key: r_val = Teuchos::rcp(new Basis_HGRAD_LINE_C1_FEM <DeviceType,outputValueType,pointValueType>()); break;
158 case shards::Beam<3>::key: r_val = Teuchos::rcp(new Basis_HGRAD_LINE_C2_FEM <DeviceType,outputValueType,pointValueType>()); break;
159 case shards::ShellLine<2>::key: r_val = Teuchos::rcp(new Basis_HGRAD_LINE_C1_FEM <DeviceType,outputValueType,pointValueType>()); break;
160 case shards::ShellLine<3>::key: r_val = Teuchos::rcp(new Basis_HGRAD_LINE_C2_FEM <DeviceType,outputValueType,pointValueType>()); break;
161 case shards::ShellTriangle<3>::key: r_val = Teuchos::rcp(new Basis_HGRAD_TRI_C1_FEM <DeviceType,outputValueType,pointValueType>()); break;
162 case shards::ShellTriangle<6>::key: r_val = Teuchos::rcp(new Basis_HGRAD_TRI_C2_FEM <DeviceType,outputValueType,pointValueType>()); break;
163 case shards::ShellQuadrilateral<4>::key: r_val = Teuchos::rcp(new Basis_HGRAD_QUAD_C1_FEM <DeviceType,outputValueType,pointValueType>()); break;
164 case shards::ShellQuadrilateral<8>::key: r_val = Teuchos::rcp(new Basis_HGRAD_QUAD_I2_FEM <DeviceType,outputValueType,pointValueType>()); break;
165 case shards::ShellQuadrilateral<9>::key: r_val = Teuchos::rcp(new Basis_HGRAD_QUAD_C2_FEM <DeviceType,outputValueType,pointValueType>()); break;
166 default: {
167 INTREPID2_TEST_FOR_EXCEPTION( true, std::invalid_argument,
168 ">>> ERROR (Intrepid2::CellTools::createHGradBasis): Cell topology not supported.");
169 }
170 }
171 return r_val;
172 }
173
174
175
183 template<EOperator Operator,
184 typename OutputViewType,
185 typename InputViewType
186 >
187 static void
188 getHGradValues( OutputViewType output,
189 const InputViewType input,
191
192 using outputValueType = typename OutputViewType::value_type;
193 using pointValueType = typename InputViewType::value_type;
194 using range_policy_type = Kokkos::MDRangePolicy<typename DeviceType::execution_space, Kokkos::Rank<2>, Kokkos::IndexType<ordinal_type> >;
195 range_policy_type policy( { 0, 0 }, { input.extent(0), input.extent(1) } );
196
197 #define KOKKOS_BASIS_PARALLEL_FOR(BasisName) \
198 do { \
199 using FunctorType = FunctorCellTools::F_getValues<OutputViewType, InputViewType, Intrepid2::Impl::BasisName, Operator>; \
200 Kokkos::parallel_for(policy, FunctorType(output, input)); \
201 } while (0)
202
204 KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_LINE_C1_FEM);
206 KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_TRI_C1_FEM);
208 KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_QUAD_C1_FEM);
210 KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_TET_C1_FEM);
212 KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_HEX_C1_FEM);
214 KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_WEDGE_C1_FEM);
216 KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_PYR_C1_FEM);
217
219 KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_LINE_C2_FEM);
221 KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_TRI_C2_FEM);
223 KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_QUAD_DEG2_FEM<true>);
225 KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_QUAD_DEG2_FEM<false>);
227 KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_TET_C2_FEM);
229 KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_TET_COMP12_FEM);
231 KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_HEX_DEG2_FEM<true>);
233 KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_HEX_DEG2_FEM<false>);
235 KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_WEDGE_DEG2_FEM<true>);
237 KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_WEDGE_DEG2_FEM<false>);
239 KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_PYR_I2_FEM);
240 else {
241 INTREPID2_TEST_FOR_EXCEPTION( true, std::invalid_argument,
242 ">>> ERROR (Intrepid2::CellTools::getHGradValues): Basis type not supported.");
243 }
244 }
245
246
254 template<EOperator Operator,
255 typename OutputViewType,
256 typename InputViewType
257 >
258 static void
259 getHGradValues( OutputViewType output,
260 const InputViewType input,
261 const shards::CellTopology cellTopo) {
262
263 using range_policy_type = Kokkos::MDRangePolicy<typename DeviceType::execution_space, Kokkos::Rank<2>, Kokkos::IndexType<ordinal_type> >;
264 range_policy_type policy( { 0, 0 }, { input.extent(0), input.extent(1) } );
265
266 #define KOKKOS_BASIS_PARALLEL_FOR(BasisName) \
267 do { \
268 using FunctorType = FunctorCellTools::F_getValues<OutputViewType, InputViewType, Intrepid2::Impl::BasisName, Operator>; \
269 Kokkos::parallel_for(policy, FunctorType(output, input)); \
270 } while (0)
271
272 switch (cellTopo.getKey()) {
273 case shards::Line<2>::key: KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_LINE_C1_FEM); break;
274 case shards::Triangle<3>::key: KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_TRI_C1_FEM); break;
275 case shards::Quadrilateral<4>::key: KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_QUAD_C1_FEM); break;
276 case shards::Tetrahedron<4>::key: KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_TET_C1_FEM); break;
277 case shards::Hexahedron<8>::key: KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_HEX_C1_FEM); break;
278 case shards::Wedge<6>::key: KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_WEDGE_C1_FEM); break;
279 case shards::Pyramid<5>::key: KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_PYR_C1_FEM); break;
280
281 case shards::Line<3>::key: KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_LINE_C2_FEM); break;
282 case shards::Triangle<6>::key: KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_TRI_C2_FEM); break;
283 case shards::Quadrilateral<8>::key: KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_QUAD_DEG2_FEM<true>); break;
284 case shards::Quadrilateral<9>::key: KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_QUAD_DEG2_FEM<false>); break;
285 case shards::Tetrahedron<10>::key: KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_TET_C2_FEM); break;
286 case shards::Tetrahedron<11>::key: KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_TET_COMP12_FEM); break;
287 case shards::Hexahedron<20>::key: KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_HEX_DEG2_FEM<true>); break;
288 case shards::Hexahedron<27>::key: KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_HEX_DEG2_FEM<false>); break;
289 case shards::Wedge<15>::key: KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_WEDGE_DEG2_FEM<true>); break;
290 case shards::Wedge<18>::key: KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_WEDGE_DEG2_FEM<false>); break;
291 case shards::Pyramid<13>::key: KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_PYR_I2_FEM); break;
292
293 case shards::Beam<2>::key: KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_LINE_C1_FEM); break;
294 case shards::Beam<3>::key: KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_LINE_C2_FEM); break;
295 case shards::ShellLine<2>::key: KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_LINE_C1_FEM); break;
296 case shards::ShellLine<3>::key: KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_LINE_C2_FEM); break;
297 case shards::ShellTriangle<3>::key: KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_TRI_C1_FEM); break;
298 case shards::ShellTriangle<6>::key: KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_TRI_C2_FEM); break;
299 case shards::ShellQuadrilateral<4>::key: KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_QUAD_C1_FEM); break;
300 case shards::ShellQuadrilateral<8>::key: KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_QUAD_DEG2_FEM<true>); break;
301 case shards::ShellQuadrilateral<9>::key: KOKKOS_BASIS_PARALLEL_FOR(Basis_HGRAD_QUAD_DEG2_FEM<false>); break;
302
303 default: {
304 INTREPID2_TEST_FOR_EXCEPTION( true, std::invalid_argument,
305 ">>> ERROR (Intrepid2::CellTools::getHGradValues): Cell topology not supported.");
306 }
307 }
308 }
309
310public:
311
314 CellTools() = default;
315
318 ~CellTools() = default;
319
320 //============================================================================================//
321 // //
322 // Jacobian, inverse Jacobian and Jacobian determinant //
323 // //
324 //============================================================================================//
325
361 template<typename JacobianViewType,
362 typename PointViewType,
363 typename WorksetType,
364 typename HGradBasisType>
365 static void
366 setJacobian( JacobianViewType jacobian,
367 const PointViewType points,
368 const WorksetType worksetCell,
369 const Teuchos::RCP<HGradBasisType> basis,
370 const int startCell=0, const int endCell=-1);
371
406 template<typename JacobianViewType,
407 typename BasisGradientsType,
408 typename WorksetType>
409 static void
410 setJacobian( JacobianViewType jacobian,
411 const WorksetType worksetCell,
412 const BasisGradientsType gradients,
413 const int startCell=0, const int endCell=-1);
414
449 template<typename JacobianViewType,
450 typename PointViewType,
451 typename WorksetCellViewType>
452 static void
453 setJacobian( JacobianViewType jacobian,
454 const PointViewType points,
455 const WorksetCellViewType worksetCell,
456 const shards::CellTopology cellTopo ) {
457 using nonConstPointValueType = typename PointViewType::non_const_value_type;
458 auto basis = createHGradBasis<nonConstPointValueType,nonConstPointValueType>(cellTopo);
459 setJacobian(jacobian,
460 points,
461 worksetCell,
462 basis);
463 }
464
465
466
484 template<typename OutViewType,
485 typename TanViewType,
486 typename InViewType>
487 static void
489 OutViewType result,
490 const TanViewType tangents,
491 const InViewType residual,
492 const typename InViewType::value_type scaling);
493
494
495
506 template<typename JacobianInvViewType,
507 typename JacobianViewType>
508 static void
509 setJacobianInv( JacobianInvViewType jacobianInv,
510 const JacobianViewType jacobian );
511
522 template<typename JacobianDetViewType,
523 typename JacobianViewType>
524 static void
525 setJacobianDet( JacobianDetViewType jacobianDet,
526 const JacobianViewType jacobian );
527
533 template<class PointScalar>
535
541 template<class PointScalar>
543
549 template<class PointScalar>
550 static void setJacobianDet( Data<PointScalar,DeviceType> & jacobianDet,
551 const Data<PointScalar,DeviceType> & jacobian);
552
558 template<class PointScalar>
559 static void setJacobianDetInv( Data<PointScalar,DeviceType> & jacobianDetInv,
560 const Data<PointScalar,DeviceType> & jacobian);
561
567 template<class PointScalar>
568 static void setJacobianInv( Data<PointScalar,DeviceType> & jacobianInv,
569 const Data<PointScalar,DeviceType> & jacobian);
570
577 template<class PointScalar>
578 static void setJacobianDividedByDet( Data<PointScalar,DeviceType> & jacobianDividedByDet,
579 const Data<PointScalar,DeviceType> & jacobian,
580 const Data<PointScalar,DeviceType> & jacobianDetInv);
581
582 //============================================================================================//
583 // //
584 // Node information //
585 // //
586 //============================================================================================//
587
588 // the node information can be used inside of kokkos functor and needs kokkos inline and
589 // exception should be an abort. for now, let's not decorate
590
598 template<typename cellCenterValueType, class ...cellCenterProperties>
599 static void
600 getReferenceCellCenter( Kokkos::DynRankView<cellCenterValueType,cellCenterProperties...> cellCenter,
601 const shards::CellTopology cell );
602
611 template<typename cellVertexValueType, class ...cellVertexProperties>
612 static void
613 getReferenceVertex( Kokkos::DynRankView<cellVertexValueType,cellVertexProperties...> cellVertex,
614 const shards::CellTopology cell,
615 const ordinal_type vertexOrd );
616
617
632 template<typename subcellVertexValueType, class ...subcellVertexProperties>
633 static void
634 getReferenceSubcellVertices( Kokkos::DynRankView<subcellVertexValueType,subcellVertexProperties...> subcellVertices,
635 const ordinal_type subcellDim,
636 const ordinal_type subcellOrd,
637 const shards::CellTopology parentCell );
638
639
640
656 template<typename cellNodeValueType, class ...cellNodeProperties>
657 static void
658 getReferenceNode( Kokkos::DynRankView<cellNodeValueType,cellNodeProperties...> cellNode,
659 const shards::CellTopology cell,
660 const ordinal_type nodeOrd );
661
662
663
677 template<typename SubcellNodeViewType>
678 static void
679 getReferenceSubcellNodes( SubcellNodeViewType subcellNodes,
680 const ordinal_type subcellDim,
681 const ordinal_type subcellOrd,
682 const shards::CellTopology parentCell );
683
709 template<typename RefEdgeTangentViewType>
710 static void
711 getReferenceEdgeTangent( RefEdgeTangentViewType refEdgeTangent,
712 const ordinal_type edgeOrd,
713 const shards::CellTopology parentCell );
714
751 template<typename RefFaceTanViewType>
752 static void
753 getReferenceFaceTangents( RefFaceTanViewType refFaceTanU,
754 RefFaceTanViewType refFaceTanV,
755 const ordinal_type faceOrd,
756 const shards::CellTopology parentCell );
757
819 template<typename RefSideNormalViewType>
820 static void
821 getReferenceSideNormal( RefSideNormalViewType refSideNormal,
822 const ordinal_type sideOrd,
823 const shards::CellTopology parentCell );
824
863 template<typename RefFaceNormalViewType>
864 static void
865 getReferenceFaceNormal( RefFaceNormalViewType refFaceNormal,
866 const ordinal_type faceOrd,
867 const shards::CellTopology parentCell );
868
898 template<typename edgeTangentValueType, class ...edgeTangentProperties,
899 typename worksetJacobianValueType, class ...worksetJacobianProperties>
900 static void
901 getPhysicalEdgeTangents( Kokkos::DynRankView<edgeTangentValueType,edgeTangentProperties...> edgeTangents,
902 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
903 const ordinal_type worksetEdgeOrd,
904 const shards::CellTopology parentCell );
905
906
918 template<typename edgeTangentValueType, class ...edgeTangentProperties,
919 typename worksetJacobianValueType, class ...worksetJacobianProperties,
920 typename edgeOrdValueType, class ...edgeOrdProperties>
921 static void
922 getPhysicalEdgeTangents( Kokkos::DynRankView<edgeTangentValueType,edgeTangentProperties...> edgeTangents,
923 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
924 const Kokkos::DynRankView<edgeOrdValueType,edgeOrdProperties...> worksetEdgeOrds,
925 const shards::CellTopology parentCell );
926
966 template<typename faceTanValueType, class ...faceTanProperties,
967 typename worksetJacobianValueType, class ...worksetJacobianProperties>
968 static void
969 getPhysicalFaceTangents( Kokkos::DynRankView<faceTanValueType,faceTanProperties...> faceTanU,
970 Kokkos::DynRankView<faceTanValueType,faceTanProperties...> faceTanV,
971 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
972 const ordinal_type worksetFaceOrd,
973 const shards::CellTopology parentCell );
974
975
988 template<typename faceTanValueType, class ...faceTanProperties,
989 typename worksetJacobianValueType, class ...worksetJacobianProperties,
990 typename faceOrdValueType, class ...faceOrdProperties>
991 static void
992 getPhysicalFaceTangents( Kokkos::DynRankView<faceTanValueType,faceTanProperties...> faceTanU,
993 Kokkos::DynRankView<faceTanValueType,faceTanProperties...> faceTanV,
994 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
995 const Kokkos::DynRankView<faceOrdValueType,faceOrdProperties...> worksetFaceOrds,
996 const shards::CellTopology parentCell );
997
998
999
1060 template<typename sideNormalValueType, class ...sideNormalProperties,
1061 typename worksetJacobianValueType, class ...worksetJacobianProperties>
1062 static void
1063 getPhysicalSideNormals( Kokkos::DynRankView<sideNormalValueType,sideNormalProperties...> sideNormals,
1064 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
1065 const ordinal_type worksetSideOrd,
1066 const shards::CellTopology parentCell );
1067
1068
1079 template<typename sideNormalValueType, class ...sideNormalProperties,
1080 typename worksetJacobianValueType, class ...worksetJacobianProperties,
1081 typename edgeOrdValueType, class ...edgeOrdProperties>
1082 static void
1083 getPhysicalSideNormals( Kokkos::DynRankView<sideNormalValueType,sideNormalProperties...> sideNormals,
1084 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
1085 const Kokkos::DynRankView<edgeOrdValueType,edgeOrdProperties...> worksetSideOrds,
1086 const shards::CellTopology parentCell );
1087
1126 template<typename faceNormalValueType, class ...faceNormalProperties,
1127 typename worksetJacobianValueType, class ...worksetJacobianProperties>
1128 static void
1129 getPhysicalFaceNormals( Kokkos::DynRankView<faceNormalValueType,faceNormalProperties...> faceNormals,
1130 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
1131 const ordinal_type worksetFaceOrd,
1132 const shards::CellTopology parentCell );
1133
1134
1146 template<typename faceNormalValueType, class ...faceNormalProperties,
1147 typename worksetJacobianValueType, class ...worksetJacobianProperties,
1148 typename faceOrdValueType, class ...faceOrdProperties>
1149 static void
1150 getPhysicalFaceNormals( Kokkos::DynRankView<faceNormalValueType,faceNormalProperties...> faceNormals,
1151 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
1152 const Kokkos::DynRankView<faceOrdValueType,faceOrdProperties...> worksetFaceOrds,
1153 const shards::CellTopology parentCell );
1154
1155 //============================================================================================//
1156 // //
1157 // Reference-to-physical frame mapping and its inverse //
1158 // //
1159 //============================================================================================//
1160
1199 template<typename PhysPointValueType,
1200 typename RefPointValueType,
1201 typename WorksetType,
1202 typename HGradBasisPtrType>
1203 static void
1204 mapToPhysicalFrame( PhysPointValueType physPoints,
1205 const RefPointValueType refPoints,
1206 const WorksetType worksetCell,
1207 const HGradBasisPtrType basis );
1208
1254 template<typename PhysPointViewType,
1255 typename RefPointViewType,
1256 typename WorksetCellViewType>
1257 static void
1258 mapToPhysicalFrame( PhysPointViewType physPoints,
1259 const RefPointViewType refPoints,
1260 const WorksetCellViewType worksetCell,
1261 const shards::CellTopology cellTopo );
1262
1263
1314 template<typename refSubcellViewType,
1315 typename paramPointViewType>
1316 static void
1317 mapToReferenceSubcell( refSubcellViewType refSubcellPoints,
1318 const paramPointViewType paramPoints,
1319 const ordinal_type subcellDim,
1320 const ordinal_type subcellOrd,
1321 const shards::CellTopology parentCell );
1322
1323
1330 template<typename refSubcellViewType, typename paramPointViewType>
1331 static void
1332 mapToReferenceSubcell( refSubcellViewType refSubcellPoints,
1333 const paramPointViewType paramPoints,
1334 const typename RefSubcellParametrization<DeviceType>::ConstViewType subcellParametrization,
1335 const ordinal_type subcellOrd);
1336
1342 template<typename refSubcellViewType, typename paramPointViewType, typename ordViewType>
1343 static void
1344 mapToReferenceSubcell( refSubcellViewType refSubcellPoints,
1345 const paramPointViewType paramPoints,
1346 const typename RefSubcellParametrization<DeviceType>::ConstViewType subcellParametrization,
1347 const ordViewType subcellOrd);
1348
1349
1350 //============================================================================================//
1351 // //
1352 // Physical-to-reference frame mapping and its inverse //
1353 // //
1354 //============================================================================================//
1355
1356
1402 template<typename refPointValueType, class ...refPointProperties,
1403 typename physPointValueType, class ...physPointProperties,
1404 typename worksetCellValueType, class ...worksetCellProperties>
1405 static bool
1406 mapToReferenceFrame( Kokkos::DynRankView<refPointValueType,refPointProperties...> refPoints,
1407 const Kokkos::DynRankView<physPointValueType,physPointProperties...> physPoints,
1408 const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
1409 const shards::CellTopology cellTopo,
1410 const physPointValueType shellThickness = -1.0 );
1411
1412
1413
1441 template<typename refPointValueType, class ...refPointProperties,
1442 typename physPointValueType, class ...physPointProperties,
1443 typename worksetCellValueType, class ...worksetCellProperties>
1444
1445 static bool
1446 mapToReferenceFrame( Kokkos::DynRankView<refPointValueType,refPointProperties...> refPoints,
1447 const Kokkos::DynRankView<physPointValueType,physPointProperties...> physPoints,
1448 const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
1450 const physPointValueType shellThickness = -1.0 );
1451
1478 template<typename refPointValueType, class ...refPointProperties,
1479 typename initGuessValueType, class ...initGuessProperties,
1480 typename physPointValueType, class ...physPointProperties,
1481 typename worksetCellValueType, class ...worksetCellProperties,
1482 typename HGradBasisPtrType>
1483 [[deprecated("Deprecated, use mapToReferenceFrame instead with refPoints initialized with initGuess.")]]
1484 static bool
1485 mapToReferenceFrameInitGuess( Kokkos::DynRankView<refPointValueType,refPointProperties...> refPoints,
1486 const Kokkos::DynRankView<initGuessValueType,initGuessProperties...> initGuess,
1487 const Kokkos::DynRankView<physPointValueType,physPointProperties...> physPoints,
1488 const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
1489 const HGradBasisPtrType basis,
1490 const physPointValueType shellThickness = -1.0 ){
1491 Kokkos::deep_copy(refPoints,initGuess);
1492 return mapToReferenceFrame(refPoints,physPoints,worksetCell,basis,shellThickness);
1493 }
1494
1495
1526 template<typename refPointValueType, class ...refPointProperties,
1527 typename initGuessValueType, class ...initGuessProperties,
1528 typename physPointValueType, class ...physPointProperties,
1529 typename worksetCellValueType, class ...worksetCellProperties>
1530 static bool
1531 mapToReferenceFrameInitGuess( Kokkos::DynRankView<refPointValueType,refPointProperties...> refPoints,
1532 const Kokkos::DynRankView<initGuessValueType,initGuessProperties...> initGuess,
1533 const Kokkos::DynRankView<physPointValueType,physPointProperties...> physPoints,
1534 const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
1535 const shards::CellTopology cellTopo,
1536 const physPointValueType shellThickness = -1.0 );
1537
1538
1539 //============================================================================================//
1540 // //
1541 // Control Volume Coordinates //
1542 // //
1543 //============================================================================================//
1544
1618 template<typename subcvCoordValueType, class ...subcvCoordProperties,
1619 typename cellCoordValueType, class ...cellCoordProperties>
1620 static void
1621 getSubcvCoords( Kokkos::DynRankView<subcvCoordValueType,subcvCoordProperties...> subcvCoords,
1622 const Kokkos::DynRankView<cellCoordValueType,cellCoordProperties...> cellCoords,
1623 const shards::CellTopology primaryCell );
1624
1625 //============================================================================================//
1626 // //
1627 // Inclusion tests //
1628 // //
1629 //============================================================================================//
1630
1631
1645 template<typename PointViewType>
1646 [[deprecated("Deprecated, use the templated ParametricDistance struct directly.")]]
1647 static typename ScalarTraits<typename PointViewType::value_type>::scalar_type
1648 parametricDistance( const PointViewType refPoint,
1649 const shards::CellTopology cellTopo);
1650
1660 template<typename PointViewType>
1661 [[deprecated("Deprecated, use checkPointwiseInclusion, or the templated ParametricDistance struct directly.")]]
1662 static bool
1663 checkPointInclusion( const PointViewType refPoint,
1664 const shards::CellTopology cellTopo,
1665 const typename ScalarTraits<typename PointViewType::value_type>::scalar_type thres =
1666 threshold<typename ScalarTraits<typename PointViewType::value_type>::scalar_type>() );
1667
1668
1669
1679 template<unsigned cellTopologyKey,
1680 typename OutputViewType,
1681 typename InputViewType>
1682 static void checkPointwiseInclusion( OutputViewType inCell,
1683 const InputViewType refPoints,
1684 const typename ScalarTraits<typename InputViewType::value_type>::scalar_type thresh =
1685 threshold<typename ScalarTraits<typename InputViewType::value_type>::scalar_type>());
1686
1687
1688
1697 template<typename InCellViewType,
1698 typename PointViewType>
1699 static void checkPointwiseInclusion( InCellViewType inCell,
1700 const PointViewType refPoints,
1701 const shards::CellTopology cellTopo,
1702 const typename ScalarTraits<typename PointViewType::value_type>::scalar_type thres =
1703 threshold<typename ScalarTraits<typename PointViewType::value_type>::scalar_type>() );
1704
1718 template<typename inCellValueType, class ...inCellProperties,
1719 typename pointValueType, class ...pointProperties,
1720 typename cellWorksetValueType, class ...cellWorksetProperties>
1721 static bool checkPointwiseInclusion( Kokkos::DynRankView<inCellValueType,inCellProperties...> inCell,
1722 const Kokkos::DynRankView<pointValueType,pointProperties...> points,
1723 const Kokkos::DynRankView<cellWorksetValueType,cellWorksetProperties...> cellWorkset,
1724 const shards::CellTopology cellTopo,
1725 const typename ScalarTraits<pointValueType>::scalar_type thres =
1726 threshold<typename ScalarTraits<pointValueType>::scalar_type>(),
1727 const typename ScalarTraits<pointValueType>::scalar_type shellThickness = -1.0 );
1728
1729
1730 // //============================================================================================//
1731 // // //
1732 // // Debug //
1733 // // //
1734 // //============================================================================================//
1735
1736
1737 // /** \brief Prints the reference space coordinates of the vertices of the specified subcell
1738 // \param subcellDim [in] - dimension of the subcell where points are mapped to
1739 // \param subcellOrd [in] - subcell ordinal
1740 // \param parentCell [in] - cell topology of the parent cell.
1741 // */
1742 // static void printSubcellVertices(const int subcellDim,
1743 // const int subcellOrd,
1744 // const shards::CellTopology & parentCell);
1745
1746
1747
1748 // /** \brief Prints the nodes of a subcell from a cell workset
1749
1750 // */
1751 // template<class ArrayCell>
1752 // static void printWorksetSubcell(const ArrayCell & cellWorkset,
1753 // const shards::CellTopology & parentCell,
1754 // const int& pCellOrd,
1755 // const int& subcellDim,
1756 // const int& subcellOrd,
1757 // const int& fieldWidth = 3);
1758 };
1759
1760 //============================================================================================//
1761 // //
1762 // Validation of input/output arguments for CellTools methods //
1763 // //
1764 //============================================================================================//
1765
1772 template<typename jacobianViewType,
1773 typename PointViewType,
1774 typename worksetCellViewType>
1775 static void
1776 CellTools_setJacobianArgs( const jacobianViewType jacobian,
1777 const PointViewType points,
1778 const worksetCellViewType worksetCell,
1779 const shards::CellTopology cellTopo );
1780
1785 template<typename jacobianInvViewType,
1786 typename jacobianViewType>
1787 static void
1788 CellTools_setJacobianInvArgs( const jacobianInvViewType jacobianInv,
1789 const jacobianViewType jacobian );
1790
1791
1796 template<typename jacobianDetViewType,
1797 typename jacobianViewType>
1798 static void
1799 CellTools_setJacobianDetArgs( const jacobianDetViewType jacobianDet,
1800 const jacobianViewType jacobian );
1801
1802
1809 template<typename physPointViewType,
1810 typename refPointViewType,
1811 typename worksetCellViewType>
1812 static void
1813 CellTools_mapToPhysicalFrameArgs( const physPointViewType physPoints,
1814 const refPointViewType refPoints,
1815 const worksetCellViewType worksetCell,
1816 const shards::CellTopology cellTopo );
1817
1818
1825 template<typename refPointViewType,
1826 typename physPointViewType,
1827 typename worksetCellViewType>
1828 static void
1829 CellTools_mapToReferenceFrameArgs( const refPointViewType refPoints,
1830 const physPointViewType physPoints,
1831 const worksetCellViewType worksetCell,
1832 const shards::CellTopology cellTopo );
1833
1834
1835
1843 template<typename refPointViewType,
1844 typename initGuessViewType,
1845 typename physPointViewType,
1846 typename worksetCellViewType>
1847 static void
1848 CellTools_mapToReferenceFrameInitGuess( const refPointViewType refPoints,
1849 const initGuessViewType initGuess,
1850 const physPointViewType physPoints,
1851 const worksetCellViewType worksetCell,
1852 const shards::CellTopology cellTopo );
1853
1854}
1855
1857
1860
1864
1866
1868
1869
1870#endif
Header file for the abstract base class Intrepid2::Basis.
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
Header file for the classes: Intrepid2::RefSubcellParametrization, Intrepid2::RefCellNodes,...
static void CellTools_mapToReferenceFrameInitGuess(const refPointViewType refPoints, const initGuessViewType initGuess, const physPointViewType physPoints, const worksetCellViewType worksetCell, const shards::CellTopology cellTopo)
Validates arguments to Intrepid2::CellTools::mapToReferenceFrame with user-defined initial guess.
static void CellTools_setJacobianArgs(const jacobianViewType jacobian, const PointViewType points, const worksetCellViewType worksetCell, const shards::CellTopology cellTopo)
Validates arguments to Intrepid2::CellTools::setJacobian.
Definition file for the control volume functions of the Intrepid2::CellTools class.
Definition file for point inclusion functions of the Intrepid2::CellTools class.
Definition file for the Jacobian functions in the Intrepid2::CellTools class.
Definition file for node data and subcell functions of the Intrepid2::CellTools class.
Definition file for the physical to reference mappings in the Intrepid2::CellTools class.
Definition file for the reference to physical mappings in the Intrepid2::CellTools class.
Definition file for the validate arguments functions of the Intrepid2::CellTools class.
void CellTools_setJacobianInvArgs(const jacobianInvViewType jacobianInv, const jacobianViewType jacobian)
Validates arguments to Intrepid2::CellTools::setJacobianInv.
void CellTools_mapToReferenceFrameArgs(const refPointViewType refPoints, const physPointViewType physPoints, const worksetCellViewType worksetCell, const shards::CellTopology cellTopo)
Validates arguments to Intrepid2::CellTools::mapToReferenceFrame with default initial guess.
void CellTools_mapToPhysicalFrameArgs(const physPointViewType physPoints, const refPointViewType refPoints, const worksetCellViewType worksetCell, const shards::CellTopology cellTopo)
Validates arguments to Intrepid2::CellTools::mapToPhysicalFrame.
void CellTools_setJacobianDetArgs(const jacobianDetViewType jacobianDet, const jacobianViewType jacobian)
Validates arguments to Intrepid2::CellTools::setJacobianDet.
Header file with additional documentation for the Intrepid2::CellTools class.
Defines the Data class, a wrapper around a Kokkos::View that allows data that is constant or repeatin...
Header file for the Intrepid2::Basis_HGRAD_HEX_C1_FEM class.
Header file for the Intrepid2::Basis_HGRAD_HEX_C2_FEM class.
Header file for the Intrepid2::Basis_HGRAD_LINE_C1_FEM class.
Header file for the Intrepid2::Basis_HGRAD_LINE_C2_FEM class.
Header file for the Intrepid2::Basis_HGRAD_PYR_C1_FEM class.
Header file for the Intrepid2::Basis_HGRAD_PYR_I2_FEM class.
Header file for the Intrepid2::Basis_HGRAD_QUAD_C1_FEM class.
Header file for the Intrepid2::Basis_HGRAD_QUAD_C2_FEM class.
Header file for the Intrepid2::Basis_HGRAD_TET_C1_FEM class.
Header file for the Intrepid2::Basis_HGRAD_TET_C2_FEM class.
Header file for the Intrepid2::Basis_HGRAD_TET_COMP12_FEM class.
Header file for the Intrepid2::Basis_HGRAD_TRI_C1_FEM class.
Header file for the Intrepid2::Basis_HGRAD_TRI_C2_FEM class.
Header file for the Intrepid2::Basis_HGRAD_WEDGE_C1_FEM class.
Header file for the Intrepid2::Basis_HGRAD_WEDGE_C2_FEM class.
Header file for Intrepid2::RealSpaceTools class providing basic linear algebra functionality in 1D,...
Contains definitions of custom data types in Intrepid2.
Header function for Intrepid2::Util class and other utility functions.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Hexahedron cell.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Hexahedron cell.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Line cell.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Line cell.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Pyramid cell.
Implementation of an H(grad)-compatible FEM basis of degree 2 on a Pyramid cell.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Quadrilateral cell.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Quadrilateral cell.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Tetrahedron cell.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Tetrahedron cell.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Tetrahedron cell.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Triangle cell.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Triangle cell.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Wedge cell.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Wedge cell.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
A stateless class for operations on cell data. Provides methods for:
static ScalarTraits< typenamePointViewType::value_type >::scalar_type parametricDistance(const PointViewType refPoint, const shards::CellTopology cellTopo)
Compute the parametric distance of the point in the specified reference cell. The parametric distance...
static Teuchos::RCP< Basis< DeviceType, outputValueType, pointValueType > > createHGradBasis(const shards::CellTopology cellTopo)
Generates default HGrad basis based on cell topology.
static void mapToReferenceSubcell(refSubcellViewType refSubcellPoints, const paramPointViewType paramPoints, const typename RefSubcellParametrization< DeviceType >::ConstViewType subcellParametrization, const ordinal_type subcellOrd)
Computes parameterization maps of 1- and 2-subcells of reference cells.
static void setJacobian(JacobianViewType jacobian, const PointViewType points, const WorksetCellViewType worksetCell, const shards::CellTopology cellTopo)
Computes the Jacobian matrix DF of the reference-to-physical frame map F.
static void getReferenceNode(Kokkos::DynRankView< cellNodeValueType, cellNodeProperties... > cellNode, const shards::CellTopology cell, const ordinal_type nodeOrd)
Retrieves the Cartesian coordinates of a reference cell node.
static void setJacobianDividedByDet(Data< PointScalar, DeviceType > &jacobianDividedByDet, const Data< PointScalar, DeviceType > &jacobian, const Data< PointScalar, DeviceType > &jacobianDetInv)
Multiplies the Jacobian with shape (C,P,D,D) by the reciprocals of the determinants,...
static Data< PointScalar, DeviceType > allocateJacobianDet(const Data< PointScalar, DeviceType > &jacobian)
Allocates and returns a Data container suitable for storing determinants corresponding to the Jacobia...
static void getHGradValues(OutputViewType output, const InputViewType input, const Basis< DeviceType, typename OutputViewType::value_type, typename InputViewType::value_type > *basis)
Computes the value or gradient of HGRAD basis functions at input points.
static void getPhysicalFaceTangents(Kokkos::DynRankView< faceTanValueType, faceTanProperties... > faceTanU, Kokkos::DynRankView< faceTanValueType, faceTanProperties... > faceTanV, const Kokkos::DynRankView< worksetJacobianValueType, worksetJacobianProperties... > worksetJacobians, const ordinal_type worksetFaceOrd, const shards::CellTopology parentCell)
Computes non-normalized tangent vector pairs to physical faces in a face workset ; (see Subcell works...
static void setJacobianDet(JacobianDetViewType jacobianDet, const JacobianViewType jacobian)
Computes the determinant of the Jacobian matrix DF of the reference-to-physical frame map F.
static void mapToPhysicalFrame(PhysPointValueType physPoints, const RefPointValueType refPoints, const WorksetType worksetCell, const HGradBasisPtrType basis)
Computes F, the reference-to-physical frame map.
static void setJacobianInv(JacobianInvViewType jacobianInv, const JacobianViewType jacobian)
Computes the inverse of the Jacobian matrix DF of the reference-to-physical frame map F.
static void getSubcvCoords(Kokkos::DynRankView< subcvCoordValueType, subcvCoordProperties... > subcvCoords, const Kokkos::DynRankView< cellCoordValueType, cellCoordProperties... > cellCoords, const shards::CellTopology primaryCell)
Computes coordinates of sub-control volumes in each primary cell.
static bool checkPointInclusion(const PointViewType refPoint, const shards::CellTopology cellTopo, const typename ScalarTraits< typename PointViewType::value_type >::scalar_type thres=threshold< typename ScalarTraits< typename PointViewType::value_type >::scalar_type >())
Checks if a point belongs to a reference cell.
static bool mapToReferenceFrame(Kokkos::DynRankView< refPointValueType, refPointProperties... > refPoints, const Kokkos::DynRankView< physPointValueType, physPointProperties... > physPoints, const Kokkos::DynRankView< worksetCellValueType, worksetCellProperties... > worksetCell, const shards::CellTopology cellTopo, const physPointValueType shellThickness=-1.0)
Computes , the inverse of the reference-to-physical frame map using a default initial guess.
static void getPhysicalFaceNormals(Kokkos::DynRankView< faceNormalValueType, faceNormalProperties... > faceNormals, const Kokkos::DynRankView< worksetJacobianValueType, worksetJacobianProperties... > worksetJacobians, const ordinal_type worksetFaceOrd, const shards::CellTopology parentCell)
Computes non-normalized normal vectors to physical faces in a face workset ; (see Subcell worksets fo...
static bool hasReferenceCell(const shards::CellTopology cellTopo)
Checks if a cell topology has reference cell.
static void mapToReferenceSubcell(refSubcellViewType refSubcellPoints, const paramPointViewType paramPoints, const ordinal_type subcellDim, const ordinal_type subcellOrd, const shards::CellTopology parentCell)
Computes parameterization maps of 1- and 2-subcells of reference cells.
static void getReferenceFaceNormal(RefFaceNormalViewType refFaceNormal, const ordinal_type faceOrd, const shards::CellTopology parentCell)
Computes constant normal vectors to faces of 3D reference cell.
static Data< PointScalar, DeviceType > allocateJacobianInv(const Data< PointScalar, DeviceType > &jacobian)
Allocates and returns a Data container suitable for storing inverses corresponding to the Jacobians i...
static void setJacobian(JacobianViewType jacobian, const PointViewType points, const WorksetType worksetCell, const Teuchos::RCP< HGradBasisType > basis, const int startCell=0, const int endCell=-1)
Computes the Jacobian matrix DF of the reference-to-physical frame map F.
static void scaledResidualNormalComponent(OutViewType result, const TanViewType tangents, const InViewType residual, const typename InViewType::value_type scaling)
Computes the scaled normal component of a d-dimensional resudial vecor with respect to a (d-1) dimens...
CellTools()=default
Default constructor.
static void checkPointwiseInclusion(InCellViewType inCell, const PointViewType refPoints, const shards::CellTopology cellTopo, const typename ScalarTraits< typename PointViewType::value_type >::scalar_type thres=threshold< typename ScalarTraits< typename PointViewType::value_type >::scalar_type >())
Checks every point in multiple sets indexed by a cell ordinal for inclusion in the reference cell of ...
~CellTools()=default
Destructor.
static void getPhysicalEdgeTangents(Kokkos::DynRankView< edgeTangentValueType, edgeTangentProperties... > edgeTangents, const Kokkos::DynRankView< worksetJacobianValueType, worksetJacobianProperties... > worksetJacobians, const ordinal_type worksetEdgeOrd, const shards::CellTopology parentCell)
Computes non-normalized tangent vectors to physical edges in an edge workset ; (see Subcell worksets ...
static void mapToReferenceSubcell(refSubcellViewType refSubcellPoints, const paramPointViewType paramPoints, const typename RefSubcellParametrization< DeviceType >::ConstViewType subcellParametrization, const ordViewType subcellOrd)
Computes parameterization maps of 1- and 2-subcells of reference cells.
static void getReferenceCellCenter(Kokkos::DynRankView< cellCenterValueType, cellCenterProperties... > cellCenter, const shards::CellTopology cell)
Computes the Cartesian coordinates of reference cell barycenter.
static void getReferenceSubcellNodes(SubcellNodeViewType subcellNodes, const ordinal_type subcellDim, const ordinal_type subcellOrd, const shards::CellTopology parentCell)
Retrieves the Cartesian coordinates of all nodes of a reference subcell.
static void getHGradValues(OutputViewType output, const InputViewType input, const shards::CellTopology cellTopo)
Computes the value or gradient of HGRAD basis functions at input points.
static void getReferenceSubcellVertices(Kokkos::DynRankView< subcellVertexValueType, subcellVertexProperties... > subcellVertices, const ordinal_type subcellDim, const ordinal_type subcellOrd, const shards::CellTopology parentCell)
Retrieves the Cartesian coordinates of all vertices of a reference subcell.
static void checkPointwiseInclusion(OutputViewType inCell, const InputViewType refPoints, const typename ScalarTraits< typename InputViewType::value_type >::scalar_type thresh=threshold< typename ScalarTraits< typename InputViewType::value_type >::scalar_type >())
Checks every point for inclusion in the reference cell of a given topology. The points can belong to ...
static void getReferenceSideNormal(RefSideNormalViewType refSideNormal, const ordinal_type sideOrd, const shards::CellTopology parentCell)
Computes constant normal vectors to sides of 2D or 3D reference cells.
static bool mapToReferenceFrameInitGuess(Kokkos::DynRankView< refPointValueType, refPointProperties... > refPoints, const Kokkos::DynRankView< initGuessValueType, initGuessProperties... > initGuess, const Kokkos::DynRankView< physPointValueType, physPointProperties... > physPoints, const Kokkos::DynRankView< worksetCellValueType, worksetCellProperties... > worksetCell, const HGradBasisPtrType basis, const physPointValueType shellThickness=-1.0)
[Deprecated] Computation of , the inverse of the reference-to-physical frame map using user-supplied ...
static void setJacobianDetInv(Data< PointScalar, DeviceType > &jacobianDetInv, const Data< PointScalar, DeviceType > &jacobian)
Computes reciprocals of determinants corresponding to the Jacobians in the Data container provided.
static void getReferenceVertex(Kokkos::DynRankView< cellVertexValueType, cellVertexProperties... > cellVertex, const shards::CellTopology cell, const ordinal_type vertexOrd)
Retrieves the Cartesian coordinates of a reference cell vertex.
static void getReferenceEdgeTangent(RefEdgeTangentViewType refEdgeTangent, const ordinal_type edgeOrd, const shards::CellTopology parentCell)
Computes constant tangent vectors to edges of 2D or 3D reference cells.
static void getReferenceFaceTangents(RefFaceTanViewType refFaceTanU, RefFaceTanViewType refFaceTanV, const ordinal_type faceOrd, const shards::CellTopology parentCell)
Computes pairs of constant tangent vectors to faces of a 3D reference cells.
static void getPhysicalSideNormals(Kokkos::DynRankView< sideNormalValueType, sideNormalProperties... > sideNormals, const Kokkos::DynRankView< worksetJacobianValueType, worksetJacobianProperties... > worksetJacobians, const ordinal_type worksetSideOrd, const shards::CellTopology parentCell)
Computes non-normalized normal vectors to physical sides in a side workset .
Wrapper around a Kokkos::View that allows data that is constant or repeating in various logical dimen...
static bool isSupported(const unsigned cellTopoKey)
Checks if a cell topology has a reference parametrization.
Computes the values or gradients of basis functions at input points.