43 #ifndef PANZER_DIRICHLET_RESIDUAL_FACEBASIS_IMPL_HPP
44 #define PANZER_DIRICHLET_RESIDUAL_FACEBASIS_IMPL_HPP
50 #include "Intrepid2_CellTools.hpp"
52 #include "Kokkos_ViewFactory.hpp"
57 template<
typename EvalT,
typename Traits>
62 std::string residual_name = p.
get<std::string>(
"Residual Name");
67 std::string field_name = p.
get<std::string>(
"DOF Name");
68 std::string dof_name = field_name+
"_"+pointRule->getName();
69 std::string value_name = p.
get<std::string>(
"Value Name");
77 TEUCHOS_ASSERT(basis_layout->extent(0)==vector_layout_dof->extent(0));
78 TEUCHOS_ASSERT(basis_layout->extent(1)==vector_layout_dof->extent(1));
79 TEUCHOS_ASSERT(basis->dimension()==vector_layout_dof->extent_int(2));
80 TEUCHOS_ASSERT(vector_layout_vector->extent(0)==vector_layout_dof->extent(0));
81 TEUCHOS_ASSERT(vector_layout_vector->extent(1)==vector_layout_dof->extent(1));
82 TEUCHOS_ASSERT(vector_layout_vector->extent(2)==vector_layout_dof->extent(2));
84 residual = PHX::MDField<ScalarT,Cell,BASIS>(residual_name, basis_layout);
85 dof = PHX::MDField<const ScalarT,Cell,Point,Dim>(dof_name, vector_layout_dof);
86 value = PHX::MDField<const ScalarT,Cell,Point,Dim>(value_name, vector_layout_vector);
93 constJac_ = pointValues.jac;
94 this->addDependentField(constJac_);
97 this->addEvaluatedField(residual);
98 this->addDependentField(dof);
99 this->addDependentField(value);
101 std::string n =
"Dirichlet Residual Face Basis Evaluator";
106 template<
typename EvalT,
typename Traits>
114 this->utils.setFieldData(pointValues.jac,fm);
115 faceNormal =
Kokkos::createDynRankView(residual.get_static_view(),
"faceNormal",dof.extent(0),dof.extent(1),dof.extent(2));
119 template<
typename EvalT,
typename Traits>
126 const shards::CellTopology & parentCell = *basis->getCellTopology();
127 const int cellDim = parentCell.getDimension();
130 const int subcellDim = 2;
131 const int subcellOrd = this->wda(workset).subcell_index;
133 const int numFaces = parentCell.getSubcellCount(subcellDim);
134 const int numFaceDofs = dof.extent_int(1);
141 Intrepid2::CellTools<PHX::exec_space>::getPhysicalFaceNormals(faceNormal,
142 pointValues.jac.get_view(),
146 const auto subcellTopo = shards::CellTopology(parentCell.getBaseCellTopologyData(subcellDim, subcellOrd));
147 TEUCHOS_ASSERT(subcellTopo.getBaseKey() == shards::Triangle<>::key ||
148 subcellTopo.getBaseKey() == shards::Quadrilateral<>::key);
152 int faceOrts[6] = {};
153 for(index_t c=0;c<workset.
num_cells;c++) {
154 const auto ort = orientations->at(details.cell_local_ids[c]);
155 ort.getFaceOrientation(faceOrts, numFaces);
158 const double ortVal = faceOrts[subcellOrd] <
static_cast<int>(subcellTopo.getVertexCount()) ? 1.0 : -1.0;
159 for(
int b=0;b<numFaceDofs;b++) {
161 for(
int d=0;d<cellDim;d++)
162 residual(c,b) += (dof(c,b,d)-value(c,b,d))*faceNormal(c,b,d);
163 residual(c,b) *= ortVal;
Kokkos::DynRankView< typename InputArray::value_type, PHX::Device > createDynRankView(const InputArray &a, const std::string &name, const DimensionPack...dims)
Wrapper to simplify Panzer use of Sacado ViewFactory.
Teuchos::RCP< const std::vector< Intrepid2::Orientation > > orientations_
T & get(const std::string &name, T def_value)
typename EvalT::ScalarT ScalarT
DirichletResidual_FaceBasis(const Teuchos::ParameterList &p)
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
void setupArrays(const Teuchos::RCP< const panzer::PointRule > &pr)
Sizes/allocates memory for arrays.
void evaluateFields(typename Traits::EvalData d)
#define TEUCHOS_ASSERT(assertion_test)