59#ifndef __INTREPID2_ORIENTATIONTOOLS_DEF_COEFF_MATRIX_HDIV_HPP__
60#define __INTREPID2_ORIENTATIONTOOLS_DEF_COEFF_MATRIX_HDIV_HPP__
63#if defined (__clang__) && !defined (__INTEL_COMPILER)
64#pragma clang system_header
73#ifdef HAVE_INTREPID2_DEBUG
74template<
typename subcellBasisType,
75typename cellBasisType>
78check_getCoeffMatrix_HDIV(
const subcellBasisType& subcellBasis,
79 const cellBasisType& cellBasis,
80 const ordinal_type subcellId,
81 const ordinal_type subcellOrt) {
82 const shards::CellTopology cellTopo = cellBasis.getBaseCellTopology();
83 const shards::CellTopology subcellTopo = subcellBasis.getBaseCellTopology();
85 const ordinal_type cellDim = cellTopo.getDimension();
86 const ordinal_type subcellDim = subcellTopo.getDimension();
88 INTREPID2_TEST_FOR_EXCEPTION( subcellDim >= cellDim,
90 ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HDIV): " \
91 "cellDim must be greater than subcellDim.");
93 const auto subcellBaseKey = subcellTopo.getBaseKey();
94 const auto cellBaseKey = cellTopo.getBaseKey();
96 INTREPID2_TEST_FOR_EXCEPTION( subcellBaseKey != shards::Line<>::key &&
97 subcellBaseKey != shards::Quadrilateral<>::key &&
98 subcellBaseKey != shards::Triangle<>::key,
100 ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HDIV): " \
101 "subcellBasis must have line, quad, or triangle topology.");
103 INTREPID2_TEST_FOR_EXCEPTION( cellBaseKey != shards::Quadrilateral<>::key &&
104 cellBaseKey != shards::Triangle<>::key &&
105 cellBaseKey != shards::Hexahedron<>::key &&
106 cellBaseKey != shards::Wedge<>::key &&
107 cellBaseKey != shards::Tetrahedron<>::key,
109 ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HDIV): " \
110 "cellBasis must have quad, triangle, hexhedron or tetrahedron topology.");
116 const bool isHDIV = cellBasis.getFunctionSpace() == FUNCTION_SPACE_HDIV;
117 INTREPID2_TEST_FOR_EXCEPTION( !isHDIV,
119 ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HDIV): "
120 "cellBasis is not HDIV.");
122 const bool subcellBasisIsHGRAD = subcellBasis.getFunctionSpace() == FUNCTION_SPACE_HGRAD;
123 const bool subcellBasisIsHVOL = subcellBasis.getFunctionSpace() == FUNCTION_SPACE_HVOL;
124 const bool cellIsTri = cellBaseKey == shards::Triangle<>::key;
125 const bool cellIsTet = cellBaseKey == shards::Tetrahedron<>::key;
126 const bool cellIsHex = cellBaseKey == shards::Hexahedron<>::key;
127 const bool cellIsQuad = cellBaseKey == shards::Quadrilateral<>::key;
129 switch (subcellDim) {
132 if (cellIsHex || cellIsQuad) {
133 INTREPID2_TEST_FOR_EXCEPTION( !(subcellBasisIsHGRAD||subcellBasisIsHVOL),
135 ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_DIV): "
136 "subcellBasis function space (1d) is not consistent to cellBasis, which should be open line hgrad, order -1.");
137 }
else if (cellIsTet || cellIsTri) {
138 INTREPID2_TEST_FOR_EXCEPTION( !subcellBasisIsHVOL,
140 ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_DIV): "
141 "subcellBasis function space (1d) is not consistent to cellBasis, which should be HVOL line, order -1.");
146 if (subcellBaseKey == shards::Quadrilateral<>::key) {
148 INTREPID2_TEST_FOR_EXCEPTION( !(subcellBasisIsHGRAD||subcellBasisIsHVOL),
150 ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HDIV): "
151 "subcellBasis function space is not compatible, which should be open line hgrad, order -1.");
152 }
else if (subcellBaseKey == shards::Triangle<>::key) {
154 INTREPID2_TEST_FOR_EXCEPTION( !subcellBasisIsHVOL,
156 ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HDIV): "
157 "subcellBasis function space is not compatible, which should HVOL, order-1.");
168template<
typename OutputViewType,
169typename subcellBasisHostType,
170typename cellBasisHostType>
180#ifdef HAVE_INTREPID2_DEBUG
184 using value_type =
typename OutputViewType::non_const_value_type;
192 const ordinal_type cellDim =
cellTopo.getDimension();
212 ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HDIV): " \
213 "Lattice size should be equal to the onber of subcell internal DoFs");
225 Kokkos::DynRankView<value_type,host_device_type>
tangentsAndNormal(
"trJacobianF", cellDim, cellDim );
249 Kokkos::DynRankView<value_type,Kokkos::LayoutLeft,host_device_type>
261 for (ordinal_type
k=0;
k<cellDim;++
k)
271 Teuchos::LAPACK<ordinal_type,value_type>
lapack;
272 ordinal_type
info = 0;
273 Kokkos::DynRankView<ordinal_type,host_device_type>
pivVec(
"pivVec",
ndofSubcell);
284 std::stringstream
ss;
285 ss <<
">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HDIV): "
286 <<
"LAPACK return with error code: "
288 INTREPID2_TEST_FOR_EXCEPTION(
true, std::runtime_error,
ss.str().c_str() );
294 const double eps = tolerence();
329 auto tmp = Kokkos::create_mirror_view_and_copy(
typename OutputViewType::device_type::memory_space(),
PhiMat);
KOKKOS_FORCEINLINE_FUNCTION constexpr std::enable_if<!std::is_pod< T >::value, typenameScalarTraits< T >::scalar_type >::type get_scalar_value(const T &obj)
functions returning the scalar value. for pod types, they return the input object itself....
static ConstViewType get(const ordinal_type subcellDim, const unsigned parentCellKey)
Returns a Kokkos view with the coefficients of the parametrization maps for the edges or faces of a r...