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 bool cellBasisIsHCURL = cellBasis.getFunctionSpace() == FUNCTION_SPACE_HCURL;
 
  129         if (cellBasisIsHCURL) {
 
  130           const bool subcellBasisIsHGRAD = subcellBasis.getFunctionSpace() == FUNCTION_SPACE_HGRAD;
 
  131           const bool subcellBasisIsHVOL  = subcellBasis.getFunctionSpace() == FUNCTION_SPACE_HVOL;
 
  132           const bool subcellBasisIsHCURL = subcellBasis.getFunctionSpace() == FUNCTION_SPACE_HCURL;
 
  133           const bool cellIsTri  = cellBaseKey == shards::Triangle<>::key;
 
  134           const bool cellIsTet  = cellBaseKey == shards::Tetrahedron<>::key;
 
  135           const bool cellIsHex  = cellBaseKey == shards::Hexahedron<>::key;
 
  136           const bool cellIsQuad = cellBaseKey == shards::Quadrilateral<>::key;
 
  138           switch (subcellDim) {
 
  141             if (cellIsHex || cellIsQuad) {
 
  142               INTREPID2_TEST_FOR_EXCEPTION( !subcellBasisIsHGRAD,
 
  144                                           ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HCURL): "  
  145                                           "subcellBasis function space (1d) is not consistent to cellBasis.");
 
  146             } 
else if (cellIsTet || cellIsTri) {
 
  147               INTREPID2_TEST_FOR_EXCEPTION( !subcellBasisIsHVOL,
 
  149                                           ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HCURL): "  
  150                                           "subcellBasis function space (1d) is not consistent to cellBasis.");
 
  155             INTREPID2_TEST_FOR_EXCEPTION( !subcellBasisIsHCURL,
 
  157                                           ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HCURL): "  
  158                                           "subcellBasis function space (2d) is not consistent to cellBasis.");
 
  165       const ordinal_type numCellBasis = cellBasis.getCardinality();
 
  166       const ordinal_type numSubcellBasis = subcellBasis.getCardinality();
 
  168       const ordinal_type ordSubcell = cellBasis.getDofOrdinal(subcellDim, subcellId, 0);
 
  169       INTREPID2_TEST_FOR_EXCEPTION( ordSubcell == -1,
 
  171                                     ">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HCURL): " \
 
  172                                     "Invalid subcellId returns -1 ordSubcell.");
 
  174       const ordinal_type ndofSubcell = cellBasis.getDofTag(ordSubcell)(3);
 
  177       DynRankViewHostType refPtsSubcell(
"refPtsSubcell", ndofSubcell, subcellDim);
 
  178       DynRankViewHostType subcellDofCoords(
"subcellDofCoords", numSubcellBasis, subcellDim);
 
  179       subcellBasis.getDofCoords(subcellDofCoords);
 
  180       for(ordinal_type i=0; i<ndofSubcell; ++i)
 
  181         for(ordinal_type d=0; d <subcellDim; ++d)
 
  182           refPtsSubcell(i,d) = subcellDofCoords(subcellBasis.getDofOrdinal(subcellDim, 0, i),d);
 
  185       DynRankViewHostType ortPtsSubcell(
"ortPtsSubcell", ndofSubcell, subcellDim);
 
  192       DynRankViewHostType refPtsCell(
"refPtsCell", ndofSubcell, cellDim);
 
  204       DynRankViewHostType refValues(
"refValues", numCellBasis, ndofSubcell, cellDim);
 
  205       cellBasis.getValues(refValues, refPtsCell, OPERATOR_VALUE);
 
  208       DynRankViewHostType subcellTangents(
"subcellTangents", numSubcellBasis, subcellDim);
 
  209       DynRankViewHostType ortJacobian(
"ortJacobian", subcellDim, subcellDim);
 
  216       DynRankViewHostType jacobianF(
"jacobianF", cellDim, subcellDim );
 
  217       switch (subcellBaseKey) {
 
  218       case shards::Line<>::key: {
 
  219         auto lineDofCoeffs = Kokkos::subview(subcellTangents, Kokkos::ALL(),0);
 
  220         subcellBasis.getDofCoeffs(lineDofCoeffs);
 
  221         auto edgeTan = Kokkos::subview(jacobianF, Kokkos::ALL(),0);
 
  223         if((cellBaseKey == shards::Triangle<>::key) || (cellBaseKey == shards::Tetrahedron<>::key))
 
  224           for(ordinal_type i=0; i<cellDim; ++i)
 
  228       case shards::Quadrilateral<>::key:
 
  229       case shards::Triangle<>::key: {
 
  230         subcellBasis.getDofCoeffs(subcellTangents);
 
  232         auto faceTanU = Kokkos::subview(jacobianF, Kokkos::ALL(), 0);
 
  233         auto faceTanV = Kokkos::subview(jacobianF, Kokkos::ALL(), 1);
 
  238               INTREPID2_TEST_FOR_EXCEPTION( 
true, std::runtime_error, 
"Should not come here" );
 
  242       Kokkos::View<value_type**,Kokkos::LayoutLeft,host_space_type> 
 
  243         refMat(
"refMat", ndofSubcell, ndofSubcell),
 
  244         ortMat(
"ortMat", ndofSubcell, ndofSubcell);
 
  245       for (ordinal_type i=0;i<ndofSubcell;++i) {
 
  246         const ordinal_type iout = cellBasis.getDofOrdinal(subcellDim, subcellId, i);
 
  247         for (ordinal_type j=0;j<ndofSubcell;++j) {
 
  249           const ordinal_type jsc = subcellBasis.getDofOrdinal(subcellDim, 0, j);
 
  250           for (ordinal_type k=0;k<subcellDim;++k)
 
  251             for (ordinal_type l=0;l<subcellDim;++l)
 
  252             for (ordinal_type d=0; d<cellDim; ++d)
 
  253                 tmp +=  refValues(iout,j,d)*jacobianF(d,l)*ortJacobian(l,k)*subcellTangents(jsc,k);
 
  255           ortMat(j,i) = (i==j);  
 
  265         Teuchos::LAPACK<ordinal_type,value_type> lapack;
 
  266         ordinal_type info = 0;
 
  267         Kokkos::View<value_type**,Kokkos::LayoutLeft,host_space_type> pivVec(
"pivVec", ndofSubcell, 1);
 
  269         lapack.GESV(ndofSubcell, ndofSubcell,
 
  272                             (ordinal_type*)pivVec.data(),
 
  278           std::stringstream ss;
 
  279           ss << 
">>> ERROR (Intrepid::OrientationTools::getCoeffMatrix_HCURL): " 
  280              << 
"LAPACK return with error code: " 
  282           INTREPID2_TEST_FOR_EXCEPTION( 
true, std::runtime_error, ss.str().c_str() );
 
  288           const double eps = threshold();
 
  289           for (ordinal_type i=0;i<ndofSubcell;++i)
 
  290             for (ordinal_type j=i;j<ndofSubcell;++j) {
 
  291               auto intOrtMat = std::round(ortMat(i,j));
 
  292               ortMat(i,j) = (std::abs(ortMat(i,j) - std::round(ortMat(i,j))) < eps) ? intOrtMat : ortMat(i,j);
 
  299         const Kokkos::pair<ordinal_type,ordinal_type> range(0, ndofSubcell);
 
  300         Kokkos::deep_copy(Kokkos::subview(output, range, range),
 
  301                           Kokkos::subview(ortMat, range, range));