55 #ifndef __INTREPID2_ORIENTATIONTOOLS_DEF_COEFF_MATRIX_HGRAD_HPP__
56 #define __INTREPID2_ORIENTATIONTOOLS_DEF_COEFF_MATRIX_HGRAD_HPP__
59 #if defined (__clang__) && !defined (__INTEL_COMPILER)
60 #pragma clang system_header
67 template<
typename OutputViewType,
68 typename subcellBasisType,
69 typename cellBasisType>
74 const subcellBasisType subcellBasis,
75 const cellBasisType cellBasis,
76 const ordinal_type subcellId,
77 const ordinal_type subcellOrt) {
78 typedef typename OutputViewType::execution_space space_type;
79 typedef typename OutputViewType::value_type value_type;
83 Kokkos::Impl::is_space<space_type>::host_mirror_space::execution_space host_space_type;
85 typedef Kokkos::DynRankView<value_type,host_space_type> DynRankViewHostType;
93 const shards::CellTopology cellTopo = cellBasis.getBaseCellTopology();
94 const shards::CellTopology subcellTopo = subcellBasis.getBaseCellTopology();
96 const ordinal_type cellDim = cellTopo.getDimension();
97 const ordinal_type subcellDim = subcellTopo.getDimension();
99 INTREPID2_TEST_FOR_EXCEPTION( subcellDim >= cellDim,
101 ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HGRAD): " \
102 "cellDim must be greater than subcellDim.");
104 const auto subcellBaseKey = subcellTopo.getBaseKey();
106 INTREPID2_TEST_FOR_EXCEPTION( subcellBaseKey != shards::Line<>::key &&
107 subcellBaseKey != shards::Quadrilateral<>::key &&
108 subcellBaseKey != shards::Triangle<>::key,
110 ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HGRAD): " \
111 "subcellBasis must have line, quad, or triangle topology.");
119 const bool cellBasisIsHGRAD = cellBasis.getFunctionSpace() == FUNCTION_SPACE_HGRAD;
120 const bool subcellBasisIsHGRAD = subcellBasis.getFunctionSpace() == FUNCTION_SPACE_HGRAD;
121 if (cellBasisIsHGRAD) {
122 INTREPID2_TEST_FOR_EXCEPTION( !subcellBasisIsHGRAD,
124 ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HGRAD): " \
125 "subcellBasis function space is not consistent to cellBasis.");
128 INTREPID2_TEST_FOR_EXCEPTION( subcellBasis.getDegree() != cellBasis.getDegree(),
130 ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HGRAD): " \
131 "subcellBasis has a different polynomial degree from cellBasis' degree.");
139 const ordinal_type numCellBasis = cellBasis.getCardinality();
140 const ordinal_type numSubcellBasis = subcellBasis.getCardinality();
142 const ordinal_type ordSubcell = cellBasis.getDofOrdinal(subcellDim, subcellId, 0);
143 INTREPID2_TEST_FOR_EXCEPTION( ordSubcell == -1,
145 ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HGRAD): " \
146 "Invalid subcellId returns -1 ordSubcell.");
148 const ordinal_type ndofSubcell = cellBasis.getDofTag(ordSubcell)(3);
151 DynRankViewHostType refPtsSubcell(
"refPtsSubcell", ndofSubcell, subcellDim);
152 DynRankViewHostType subcellDofCoords(
"subcellDofCoords", numSubcellBasis, subcellDim);
153 subcellBasis.getDofCoords(subcellDofCoords);
154 for(ordinal_type i=0; i<ndofSubcell; ++i)
155 for(ordinal_type d=0; d <subcellDim; ++d)
156 refPtsSubcell(i,d) = subcellDofCoords(subcellBasis.getDofOrdinal(subcellDim, 0, i),d);
159 DynRankViewHostType ortPtsSubcell(
"ortPtsSubcell", ndofSubcell, subcellDim);
166 DynRankViewHostType refPtsCell(
"refPtsCell", ndofSubcell, cellDim);
177 DynRankViewHostType refValues(
"refValues", numCellBasis, ndofSubcell);
178 cellBasis.getValues(refValues, refPtsCell, OPERATOR_VALUE);
186 Kokkos::View<value_type**,Kokkos::LayoutLeft,host_space_type>
187 refMat(
"refMat", ndofSubcell, ndofSubcell),
188 ortMat(
"ortMat", ndofSubcell, ndofSubcell),
189 pivVec(
"pivVec", ndofSubcell, 1);
191 for (ordinal_type i=0;i<ndofSubcell;++i) {
192 const ordinal_type iref = cellBasis.getDofOrdinal(subcellDim, subcellId, i);
193 for (ordinal_type j=0;j<ndofSubcell;++j) {
194 refMat(i,j) = refValues(iref,j);
195 ortMat(i,j) = (i==j);
201 Teuchos::LAPACK<ordinal_type,value_type> lapack;
202 ordinal_type info = 0;
204 lapack.GESV(ndofSubcell, ndofSubcell,
207 (ordinal_type*)pivVec.data(),
213 std::stringstream ss;
214 ss <<
">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HGRAD): "
215 <<
"LAPACK return with error code: "
217 INTREPID2_TEST_FOR_EXCEPTION(
true, std::runtime_error, ss.str().c_str() );
222 const double eps = threshold();
223 for (ordinal_type i=0;i<ndofSubcell;++i)
224 for (ordinal_type j=i;j<ndofSubcell;++j) {
225 auto intOrtMat = std::round(ortMat(i,j));
226 ortMat(i,j) = (std::abs(ortMat(i,j) - std::round(ortMat(i,j))) < eps) ? intOrtMat : ortMat(i,j);
233 const Kokkos::pair<ordinal_type,ordinal_type> range(0, ndofSubcell);
234 Kokkos::deep_copy(Kokkos::subview(output, range, range), ortMat);