Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_Dirichlet_Residual_EdgeBasis_impl.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Panzer: A partial differential equation assembly
4// engine for strongly coupled complex multiphysics systems
5//
6// Copyright 2011 NTESS and the Panzer contributors.
7// SPDX-License-Identifier: BSD-3-Clause
8// *****************************************************************************
9// @HEADER
10
11#ifndef PANZER_DIRICHLET_RESIDUAL_EDGEBASIS_IMPL_HPP
12#define PANZER_DIRICHLET_RESIDUAL_EDGEBASIS_IMPL_HPP
13
14#include <cstddef>
15#include <string>
16#include <vector>
17
18#include "Intrepid2_Kernels.hpp"
19#include "Intrepid2_CellTools.hpp"
20#include "Intrepid2_OrientationTools.hpp"
21
22#include "Phalanx_Print.hpp"
23
25#include "Kokkos_ViewFactory.hpp"
26
27namespace panzer {
28
29//**********************************************************************
30template<typename EvalT, typename Traits>
33 const Teuchos::ParameterList& p)
34{
35 std::string residual_name = p.get<std::string>("Residual Name");
36
37 basis = p.get<Teuchos::RCP<const panzer::PureBasis> >("Basis");
38 pointRule = p.get<Teuchos::RCP<const panzer::PointRule> >("Point Rule");
39
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");
43
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;
47
48 // some sanity checks
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));
56
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);
60
61 // setup all basis fields that are required
62
63 // setup all fields to be evaluated and constructed
64 pointValues = PointValues2<double>(pointRule->getName()+"_",false);
65 pointValues.setupArrays(pointRule);
66
67 // the field manager will allocate all of these field
68 this->addNonConstDependentField(pointValues.jac);
69
70 this->addEvaluatedField(residual);
71 this->addDependentField(dof);
72 this->addDependentField(value);
73
74 std::string n = "Dirichlet Residual Edge Basis Evaluator";
75 this->setName(n);
76}
77
78//**********************************************************************
79template<typename EvalT, typename Traits>
80void
83 typename Traits::SetupData sd,
85{
86 orientations = sd.orientations_;
87 this->utils.setFieldData(pointValues.jac,fm);
88
89 const auto cellTopo = *basis->getCellTopology();
90 const int edgeDim = 1;
91 const int faceDim = 2;
92 if(cellTopo.getDimension() > edgeDim)
93 edgeParam = Intrepid2::RefSubcellParametrization<Kokkos::HostSpace>::get(edgeDim, cellTopo.getKey());
94
95 if(cellTopo.getDimension() > faceDim)
96 faceParam = Intrepid2::RefSubcellParametrization<Kokkos::HostSpace>::get(faceDim, cellTopo.getKey());
97}
98
99//**********************************************************************
100template<typename EvalT, typename Traits>
101void
104 typename Traits::EvalData workset)
105{
106 const int numCells = workset.num_cells;
107 if(numCells <= 0)
108 return;
109 else {
110 residual.deep_copy(ScalarT(0.0));
111
112 // dofs are already oriented but tangent directions are not oriented
113
114 const int subcellDim = workset.subcell_dim;
115 const int subcellOrd = this->wda(workset).subcell_index;
116
117 const auto cellTopo = *basis->getCellTopology();
118 const auto worksetJacobians = pointValues.jac.get_view();
119
120 const int cellDim = cellTopo.getDimension();
121 const int edgeDim = 1;
122 const int faceDim = 2;
123
124 auto intrepid_basis = basis->getIntrepid2Basis();
125 const WorksetDetails & details = workset;
126
127 //const bool is_normalize = true;
128 auto work = Kokkos::create_mirror_view(Kokkos::createDynRankView(residual.get_static_view(),"work", 4, cellDim));
129
130 // compute residual
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) {
138 case 1: { // 2D element Tri and Quad
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());
142
143 const int ndofsEdge = intrepid_basis->getDofCount(1, subcellOrd);
144 const int numEdges = cellTopo.getEdgeCount();
145 /* */ int edgeOrts[4] = {};
146
147 for(index_t c=0;c<workset.num_cells;c++) {
148 orientations->at(details.cell_local_ids[c]).getEdgeOrientation(edgeOrts, numEdges);
149
150 Intrepid2::Impl::OrientationTools::getRefSubcellTangents(work,
151 edgeParam,
152 cellTopo.getKey(edgeDim,subcellOrd),
153 subcellOrd,
154 edgeOrts[subcellOrd]);
155
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);
160
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);
163 }
164 }
165 }
166 }
167 break;
168 }
169 case 2: { // 3D element Tet and Hex
170 const int numEdges = cellTopo.getEdgeCount();
171 const int numFaces = cellTopo.getFaceCount();
172
173 {
174
175 auto ortEdgeTan = Kokkos::subview(work, 0, Kokkos::ALL());
176 auto phyEdgeTan = Kokkos::subview(work, 1, Kokkos::ALL());
177
178 const int numEdgesOfFace= cellTopo.getEdgeCount(2, subcellOrd);
179
180 int edgeOrts[12] = {};
181 for(index_t c=0;c<workset.num_cells;c++) {
182 for (int i=0;i<numEdgesOfFace;++i) {
183
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);
187
188 Intrepid2::Impl::OrientationTools::getRefSubcellTangents(work,
189 edgeParam,
190 cellTopo.getKey(edgeDim,edgeOrd),
191 edgeOrd,
192 edgeOrts[edgeOrd]);
193
194 {
195 auto J = Kokkos::subview(worksetJacobians_h, c, b, Kokkos::ALL(), Kokkos::ALL());
196 Intrepid2::Kernels::Serial::matvec_product(phyEdgeTan, J, ortEdgeTan);
197
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);
200 }
201 }
202 }
203 }
204 }
205
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());
211
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);
215
216 Intrepid2::Impl::OrientationTools::getRefSubcellTangents(work,
217 faceParam,
218 cellTopo.getKey(faceDim,subcellOrd),
219 subcellOrd,
220 faceOrts[subcellOrd]);
221
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);
226
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);
230 }
231 }
232 }
233 }
234
235 break;
236 }
237 }
238 Kokkos::deep_copy(residual.get_static_view(), residual_h);
239 }
240
241}
242
243//**********************************************************************
244
245}
246
247#endif
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
int subcell_dim
DEPRECATED - use: getSubcellDimension()
int num_cells
DEPRECATED - use: numCells()
Teuchos::RCP< const std::vector< Intrepid2::Orientation > > orientations_