33 const Teuchos::ParameterList& p)
35 std::string residual_name = p.get<std::string>(
"Residual Name");
37 basis = p.get<Teuchos::RCP<const panzer::PureBasis> >(
"Basis");
38 pointRule = p.get<Teuchos::RCP<const panzer::PointRule> >(
"Point Rule");
40 std::string field_name = p.get<std::string>(
"DOF Name");
41 std::string dof_name = field_name+
"_"+pointRule->getName();
42 std::string value_name = p.get<std::string>(
"Value Name");
44 Teuchos::RCP<PHX::DataLayout> basis_layout = basis->functional;
45 Teuchos::RCP<PHX::DataLayout> vector_layout_dof = pointRule->dl_vector;
46 Teuchos::RCP<PHX::DataLayout> vector_layout_vector = basis->functional_grad;
49 TEUCHOS_ASSERT(basis->isVectorBasis());
50 TEUCHOS_ASSERT(basis_layout->extent(0)==vector_layout_dof->extent(0));
51 TEUCHOS_ASSERT(basis_layout->extent(1)==vector_layout_dof->extent(1));
52 TEUCHOS_ASSERT(Teuchos::as<unsigned>(basis->dimension())==vector_layout_dof->extent(2));
53 TEUCHOS_ASSERT(vector_layout_vector->extent(0)==vector_layout_dof->extent(0));
54 TEUCHOS_ASSERT(vector_layout_vector->extent(1)==vector_layout_dof->extent(1));
55 TEUCHOS_ASSERT(vector_layout_vector->extent(2)==vector_layout_dof->extent(2));
57 residual = PHX::MDField<ScalarT,Cell,BASIS>(residual_name, basis_layout);
58 dof = PHX::MDField<const ScalarT,Cell,Point,Dim>(dof_name, vector_layout_dof);
59 value = PHX::MDField<const ScalarT,Cell,Point,Dim>(value_name, vector_layout_vector);
65 pointValues.setupArrays(pointRule);
68 this->addNonConstDependentField(pointValues.jac);
70 this->addEvaluatedField(residual);
71 this->addDependentField(dof);
72 this->addDependentField(value);
74 std::string n =
"Dirichlet Residual Edge Basis Evaluator";
110 residual.deep_copy(
ScalarT(0.0));
115 const int subcellOrd = this->wda(workset).subcell_index;
117 const auto cellTopo = *basis->getCellTopology();
118 const auto worksetJacobians = pointValues.jac.get_view();
120 const int cellDim = cellTopo.getDimension();
121 const int edgeDim = 1;
122 const int faceDim = 2;
124 auto intrepid_basis = basis->getIntrepid2Basis();
128 auto work = Kokkos::create_mirror_view(Kokkos::createDynRankView(residual.get_static_view(),
"work", 4, cellDim));
131 auto residual_h = Kokkos::create_mirror_view(residual.get_static_view());
132 auto dof_h = Kokkos::create_mirror_view(dof.get_static_view());
133 auto value_h = Kokkos::create_mirror_view(value.get_static_view());
134 Kokkos::deep_copy(dof_h, dof.get_static_view());
135 Kokkos::deep_copy(value_h, value.get_static_view());
136 auto worksetJacobians_h = Kokkos::create_mirror_view_and_copy(Kokkos::HostSpace(),worksetJacobians);
137 switch (subcellDim) {
139 if (intrepid_basis->getDofCount(1, subcellOrd)) {
140 auto ortEdgeTan = Kokkos::subview(work, 0, Kokkos::ALL());
141 auto phyEdgeTan = Kokkos::subview(work, 1, Kokkos::ALL());
143 const int ndofsEdge = intrepid_basis->getDofCount(1, subcellOrd);
144 const int numEdges = cellTopo.getEdgeCount();
145 int edgeOrts[4] = {};
147 for(index_t c=0;c<workset.
num_cells;c++) {
148 orientations->at(
details.cell_local_ids[c]).getEdgeOrientation(edgeOrts, numEdges);
150 Intrepid2::Impl::OrientationTools::getRefSubcellTangents(work,
152 cellTopo.getKey(edgeDim,subcellOrd),
154 edgeOrts[subcellOrd]);
156 for (
int i=0;i<ndofsEdge;++i) {
157 const int b = intrepid_basis->getDofOrdinal(1, subcellOrd, i);
158 auto J = Kokkos::subview(worksetJacobians_h, c, b, Kokkos::ALL(), Kokkos::ALL());
159 Intrepid2::Kernels::Serial::matvec_product(phyEdgeTan, J, ortEdgeTan);
161 for(
int d=0;d<cellDim;d++) {
162 residual_h(c,b) += (dof_h(c,b,d)-value_h(c,b,d))*phyEdgeTan(d);
170 const int numEdges = cellTopo.getEdgeCount();
171 const int numFaces = cellTopo.getFaceCount();
175 auto ortEdgeTan = Kokkos::subview(work, 0, Kokkos::ALL());
176 auto phyEdgeTan = Kokkos::subview(work, 1, Kokkos::ALL());
178 const int numEdgesOfFace= cellTopo.getEdgeCount(2, subcellOrd);
180 int edgeOrts[12] = {};
181 for(index_t c=0;c<workset.
num_cells;c++) {
182 for (
int i=0;i<numEdgesOfFace;++i) {
184 const int edgeOrd = Intrepid2::Orientation::getEdgeOrdinalOfFace(i, subcellOrd, cellTopo);
185 const int b = edgeOrd;
186 orientations->at(
details.cell_local_ids[c]).getEdgeOrientation(edgeOrts, numEdges);
188 Intrepid2::Impl::OrientationTools::getRefSubcellTangents(work,
190 cellTopo.getKey(edgeDim,edgeOrd),
195 auto J = Kokkos::subview(worksetJacobians_h, c, b, Kokkos::ALL(), Kokkos::ALL());
196 Intrepid2::Kernels::Serial::matvec_product(phyEdgeTan, J, ortEdgeTan);
198 for(
int d=0;d<dof.extent_int(2);d++) {
199 residual_h(c,b) += (dof_h(c,b,d)-value_h(c,b,d))*phyEdgeTan(d);
206 if (intrepid_basis->getDofCount(2, subcellOrd)) {
207 auto ortFaceTanU = Kokkos::subview(work, 0, Kokkos::ALL());
208 auto ortFaceTanV = Kokkos::subview(work, 1, Kokkos::ALL());
209 auto phyFaceTanU = Kokkos::subview(work, 2, Kokkos::ALL());
210 auto phyFaceTanV = Kokkos::subview(work, 3, Kokkos::ALL());
212 int faceOrts[6] = {};
213 for(index_t c=0;c<workset.
num_cells;c++) {
214 orientations->at(
details.cell_local_ids[c]).getFaceOrientation(faceOrts, numFaces);
216 Intrepid2::Impl::OrientationTools::getRefSubcellTangents(work,
218 cellTopo.getKey(faceDim,subcellOrd),
220 faceOrts[subcellOrd]);
222 for(
int b=0;b<dof.extent_int(1);b++) {
223 auto J = Kokkos::subview(worksetJacobians, c, b, Kokkos::ALL(), Kokkos::ALL());
224 Intrepid2::Kernels::Serial::matvec_product(phyFaceTanU, J, ortFaceTanU);
225 Intrepid2::Kernels::Serial::matvec_product(phyFaceTanV, J, ortFaceTanV);
227 for(
int d=0;d<dof.extent_int(2);d++) {
228 residual_h(c,b) += (dof_h(c,b,d)-value_h(c,b,d))*phyFaceTanU(d);
229 residual_h(c,b) += (dof_h(c,b,d)-value_h(c,b,d))*phyFaceTanV(d);
238 Kokkos::deep_copy(residual.get_static_view(), residual_h);