59 #ifndef __INTREPID2_ORIENTATIONTOOLS_DEF_COEFF_MATRIX_HCURL_HPP__
60 #define __INTREPID2_ORIENTATIONTOOLS_DEF_COEFF_MATRIX_HCURL_HPP__
63 #if defined (__clang__) && !defined (__INTEL_COMPILER)
64 #pragma clang system_header
71 template<
typename outputViewType,
72 typename subcellBasisType,
73 typename cellBasisType>
78 const subcellBasisType subcellBasis,
79 const cellBasisType cellBasis,
80 const ordinal_type subcellId,
81 const ordinal_type subcellOrt) {
82 typedef typename outputViewType::execution_space space_type;
83 typedef typename outputViewType::value_type value_type;
87 Kokkos::Impl::is_space<space_type>::host_mirror_space::execution_space host_space_type;
89 typedef Kokkos::DynRankView<value_type,host_space_type> DynRankViewHostType;
95 const shards::CellTopology cellTopo = cellBasis.getBaseCellTopology();
96 const shards::CellTopology subcellTopo = subcellBasis.getBaseCellTopology();
98 const ordinal_type cellDim = cellTopo.getDimension();
99 const ordinal_type subcellDim = subcellTopo.getDimension();
101 INTREPID2_TEST_FOR_EXCEPTION( subcellDim >= cellDim,
103 ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HCURL): " \
104 "cellDim must be greater than subcellDim.");
106 const auto subcellBaseKey = subcellTopo.getBaseKey();
107 const auto cellBaseKey = cellTopo.getBaseKey();
109 INTREPID2_TEST_FOR_EXCEPTION( subcellBaseKey != shards::Line<>::key &&
110 subcellBaseKey != shards::Quadrilateral<>::key &&
111 subcellBaseKey != shards::Triangle<>::key,
113 ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HCURL): " \
114 "subcellBasis must have line, quad, or triangle topology.");
116 INTREPID2_TEST_FOR_EXCEPTION( cellBaseKey != shards::Quadrilateral<>::key &&
117 cellBaseKey != shards::Triangle<>::key &&
118 cellBaseKey != shards::Hexahedron<>::key &&
119 cellBaseKey != shards::Tetrahedron<>::key,
121 ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HCURL): " \
122 "cellBasis must have quad, triangle, hexhedron or tetrahedron topology.");
128 const std::string cellBasisName(cellBasis.getName());
129 if (cellBasisName.find(
"HCURL") != std::string::npos) {
130 const std::string subcellBasisName(subcellBasis.getName());
132 switch (subcellDim) {
135 if ((cellBasisName.find(
"HEX") != std::string::npos) || (cellBasisName.find(
"QUAD") != std::string::npos)) {
136 INTREPID2_TEST_FOR_EXCEPTION( subcellBasisName.find(
"HGRAD") == std::string::npos,
138 ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HCURL): "
139 "subcellBasis function space (1d) is not consistent to cellBasis.");
140 }
else if ((cellBasisName.find(
"TET") != std::string::npos) || (cellBasisName.find(
"TRI") != std::string::npos)) {
141 INTREPID2_TEST_FOR_EXCEPTION( subcellBasisName.find(
"HVOL") == std::string::npos,
143 ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HCURL): "
144 "subcellBasis function space (1d) is not consistent to cellBasis.");
149 INTREPID2_TEST_FOR_EXCEPTION( subcellBasisName.find(
"HCURL") == std::string::npos,
151 ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HCURL): "
152 "subcellBasis function space (2d) is not consistent to cellBasis.");
159 const ordinal_type numCellBasis = cellBasis.getCardinality();
160 const ordinal_type numSubcellBasis = subcellBasis.getCardinality();
162 const ordinal_type ordSubcell = cellBasis.getDofOrdinal(subcellDim, subcellId, 0);
163 INTREPID2_TEST_FOR_EXCEPTION( ordSubcell == -1,
165 ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HCURL): " \
166 "Invalid subcellId returns -1 ordSubcell.");
168 const ordinal_type ndofSubcell = cellBasis.getDofTag(ordSubcell)(3);
171 DynRankViewHostType refPtsSubcell(
"refPtsSubcell", ndofSubcell, subcellDim);
172 DynRankViewHostType subcellDofCoords(
"subcellDofCoords", numSubcellBasis, subcellDim);
173 subcellBasis.getDofCoords(subcellDofCoords);
174 for(ordinal_type i=0; i<ndofSubcell; ++i)
175 for(ordinal_type d=0; d <subcellDim; ++d)
176 refPtsSubcell(i,d) = subcellDofCoords(subcellBasis.getDofOrdinal(subcellDim, 0, i),d);
179 DynRankViewHostType ortPtsSubcell(
"ortPtsSubcell", ndofSubcell, subcellDim);
186 DynRankViewHostType refPtsCell(
"refPtsCell", ndofSubcell, cellDim);
198 DynRankViewHostType refValues(
"refValues", numCellBasis, ndofSubcell, cellDim);
199 cellBasis.getValues(refValues, refPtsCell, OPERATOR_VALUE);
202 DynRankViewHostType subcellTangents(
"subcellTangents", numSubcellBasis, subcellDim);
203 DynRankViewHostType ortJacobian(
"ortJacobian", subcellDim, subcellDim);
210 DynRankViewHostType jacobianF(
"jacobianF", cellDim, subcellDim );
211 switch (subcellBaseKey) {
212 case shards::Line<>::key: {
213 auto lineDofCoeffs = Kokkos::subview(subcellTangents, Kokkos::ALL(),0);
214 subcellBasis.getDofCoeffs(lineDofCoeffs);
215 auto edgeTan = Kokkos::subview(jacobianF, Kokkos::ALL(),0);
217 if((cellBaseKey == shards::Triangle<>::key) || (cellBaseKey == shards::Tetrahedron<>::key))
218 for(ordinal_type i=0; i<cellDim; ++i)
222 case shards::Quadrilateral<>::key:
223 case shards::Triangle<>::key: {
224 subcellBasis.getDofCoeffs(subcellTangents);
226 auto faceTanU = Kokkos::subview(jacobianF, Kokkos::ALL(), 0);
227 auto faceTanV = Kokkos::subview(jacobianF, Kokkos::ALL(), 1);
232 INTREPID2_TEST_FOR_EXCEPTION(
true, std::runtime_error,
"Should not come here" );
236 Kokkos::View<value_type**,Kokkos::LayoutLeft,host_space_type>
237 refMat(
"refMat", ndofSubcell, ndofSubcell),
238 ortMat(
"ortMat", ndofSubcell, ndofSubcell);
239 for (ordinal_type i=0;i<ndofSubcell;++i) {
240 const ordinal_type iout = cellBasis.getDofOrdinal(subcellDim, subcellId, i);
241 for (ordinal_type j=0;j<ndofSubcell;++j) {
243 const ordinal_type jsc = subcellBasis.getDofOrdinal(subcellDim, 0, j);
244 for (ordinal_type k=0;k<subcellDim;++k)
245 for (ordinal_type l=0;l<subcellDim;++l)
246 for (ordinal_type d=0; d<cellDim; ++d)
247 tmp += refValues(iout,j,d)*jacobianF(d,l)*ortJacobian(l,k)*subcellTangents(jsc,k);
249 ortMat(j,i) = (i==j);
259 Teuchos::LAPACK<ordinal_type,value_type> lapack;
260 ordinal_type info = 0;
261 Kokkos::View<value_type**,Kokkos::LayoutLeft,host_space_type> pivVec(
"pivVec", ndofSubcell, 1);
263 lapack.GESV(ndofSubcell, ndofSubcell,
266 (ordinal_type*)pivVec.data(),
272 std::stringstream ss;
273 ss <<
">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HCURL): "
274 <<
"LAPACK return with error code: "
276 INTREPID2_TEST_FOR_EXCEPTION(
true, std::runtime_error, ss.str().c_str() );
282 const double eps = threshold();
283 for (ordinal_type i=0;i<ndofSubcell;++i)
284 for (ordinal_type j=i;j<ndofSubcell;++j) {
285 auto intOrtMat = std::round(ortMat(i,j));
286 ortMat(i,j) = (std::abs(ortMat(i,j) - std::round(ortMat(i,j))) < eps) ? intOrtMat : ortMat(i,j);
293 const Kokkos::pair<ordinal_type,ordinal_type> range(0, ndofSubcell);
294 Kokkos::deep_copy(Kokkos::subview(output, range, range),
295 Kokkos::subview(ortMat, range, range));