Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_ProjectToFaces_impl.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Panzer: A partial differential equation assembly
5 // engine for strongly coupled complex multiphysics systems
6 // Copyright (2011) Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Roger P. Pawlowski (rppawlo@sandia.gov) and
39 // Eric C. Cyr (eccyr@sandia.gov)
40 // ***********************************************************************
41 // @HEADER
42 
43 #ifndef PANZER_PROJECT_TO_FACES_IMPL_HPP
44 #define PANZER_PROJECT_TO_FACES_IMPL_HPP
45 
46 #include "Teuchos_Assert.hpp"
47 #include "Phalanx_DataLayout.hpp"
48 
49 #include "Intrepid2_Cubature.hpp"
50 #include "Intrepid2_DefaultCubatureFactory.hpp"
51 #include "Intrepid2_FunctionSpaceTools.hpp"
52 #include "Intrepid2_Kernels.hpp"
53 #include "Intrepid2_CellTools.hpp"
54 #include "Intrepid2_OrientationTools.hpp"
55 
57 #include "Panzer_PureBasis.hpp"
59 #include "Kokkos_ViewFactory.hpp"
60 
61 #include "Teuchos_FancyOStream.hpp"
62 
63 #include <cstring>
64 
65 template<typename EvalT,typename Traits>
68  : quad_degree(1),
69  use_fast_method_on_rectangular_hex_mesh(false)
70 {
71  dof_name = (p.get< std::string >("DOF Name"));
72 
73  if(p.isType< Teuchos::RCP<PureBasis> >("Basis"))
74  basis = p.get< Teuchos::RCP<PureBasis> >("Basis");
75  else
76  basis = p.get< Teuchos::RCP<const PureBasis> >("Basis");
77 
78  if(p.isType<int>("Quadrature Order"))
79  quad_degree = p.get<int>("Quadrature Order");
80 
81  if(p.isType<bool>("Use Fast Method for Rectangular Hex Mesh"))
82  use_fast_method_on_rectangular_hex_mesh = p.get<bool>("Use Fast Method for Rectangular Hex Mesh");
83 
84  Teuchos::RCP<PHX::DataLayout> basis_layout = basis->functional;
85  Teuchos::RCP<PHX::DataLayout> vector_layout = basis->functional_grad;
86 
87  // some sanity checks
88  TEUCHOS_ASSERT(basis->isVectorBasis());
89 
90  result = PHX::MDField<ScalarT,Cell,BASIS>(dof_name,basis_layout);
91  this->addEvaluatedField(result);
92 
93  normals = PHX::MDField<const ScalarT,Cell,BASIS,Dim>(dof_name+"_Normals",vector_layout);
94  this->addDependentField(normals);
95 
97  const shards::CellTopology & parentCell = *basis->getCellTopology();
98  Intrepid2::DefaultCubatureFactory quadFactory;
100  = quadFactory.create<PHX::exec_space,double,double>(parentCell.getCellTopologyData(2,0), quad_degree);
101  int numQPoints = quadRule->getNumPoints();
102 
103  vector_values.resize(numQPoints);
104  for (int qp(0); qp < numQPoints; ++qp)
105  {
106  vector_values[qp] = PHX::MDField<const ScalarT,Cell,BASIS,Dim>(dof_name+"_Vector"+"_qp_"+std::to_string(qp),vector_layout);
107  this->addDependentField(vector_values[qp]);
108  }
109 
110  gatherFieldNormals = PHX::MDField<ScalarT,Cell,NODE,Dim>(dof_name+"_Normals",basis->functional_grad);
111  this->addEvaluatedField(gatherFieldNormals);
112 
113  } else {
114  TEUCHOS_ASSERT(quad_degree == 1); // One pt quadrature for fast method
115  TEUCHOS_ASSERT(std::strstr(basis->getCellTopology()->getBaseName(),"Hexahedron") != nullptr);
116 
117  vector_values.resize(1);
118  vector_values[0] = PHX::MDField<const ScalarT,Cell,BASIS,Dim>(dof_name+"_Vector",vector_layout);
119  this->addDependentField(vector_values[0]);
120  }
121 
122  if (use_fast_method_on_rectangular_hex_mesh)
123  this->setName("Project To Faces (Fast Rectangular Hex)");
124  else
125  this->setName("Project To Faces");
126 }
127 
128 // **********************************************************************
129 template<typename EvalT,typename Traits>
132  PHX::FieldManager<Traits>& /* fm */)
133 {
134  orientations = d.orientations_;
135 
136  num_pts = result.extent(1);
137  num_dim = vector_values[0].extent(2);
138 
139  TEUCHOS_ASSERT(result.extent(1) == normals.extent(1));
140  TEUCHOS_ASSERT(vector_values[0].extent(2) == normals.extent(2));
141 
142 
143  // Cache off vectors and keep intrepid2 function calls on host (by
144  // switching memory space to ProjectionSpace). This is done to
145  // prevent many small UVM allocations of temporaries in the
146  // intrepid2 functions.
147  if (use_fast_method_on_rectangular_hex_mesh) {
148  // Reminder to fix the fast method in the future
149  // quadWts = Kokkos::DynRankView<double,PHX::Device>("quadWts",numQPoints);
150  // quadPts = Kokkos::DynRankView<double,PHX::Device>("quadPts",numQPoints,num_dim);
151  } else {
152  const shards::CellTopology & parentCell = *basis->getCellTopology();
153  const int cellDim = parentCell.getDimension();
154  using HostDRVType = Kokkos::DynRankView<ScalarT,ProjectionSpace>;
155  refEdges = Kokkos::createDynRankViewWithType<HostDRVType>(result.get_static_view(),"ref_edges", 2, cellDim);
156  phyEdges = Kokkos::createDynRankViewWithType<HostDRVType>(result.get_static_view(),"phy_edges", 2, cellDim);
157 
158  PHX::MDField<double,Cell,panzer::NODE,Dim> vertex_coords = (*d.worksets_)[0].cell_vertex_coordinates;
159  physicalNodes = Kokkos::DynRankView<double,ProjectionSpace>("physicalNodes",1,vertex_coords.extent(1),num_dim);
160 
161  const int subcell_dim = 2;
162  Intrepid2::DefaultCubatureFactory quadFactory;
163  int maxNumFaceQP = 0;
164  faceQuads.resize(num_pts);
165  for (int p = 0; p < num_pts; ++p){
166  const shards::CellTopology & subcell = parentCell.getCellTopologyData(subcell_dim,p);
167  faceQuads[p] = quadFactory.create<ProjectionSpace::execution_space,double,double>(subcell, quad_degree);
168  TEUCHOS_ASSERT(faceQuads[p]->getNumPoints() == static_cast<int>(vector_values.size()));
169  maxNumFaceQP = std::max(maxNumFaceQP,faceQuads[p]->getNumPoints());
170  }
171  quadWts = Kokkos::DynRankView<double,ProjectionSpace>("quadWts",maxNumFaceQP);
172  quadPts = Kokkos::DynRankView<double,ProjectionSpace>("quadPts",maxNumFaceQP,subcell_dim);
173  refQuadPts = Kokkos::DynRankView<double,ProjectionSpace>("refQuadPts",maxNumFaceQP,num_dim);
174  jacobianSide = Kokkos::DynRankView<double,ProjectionSpace>("jacobianSide",1,maxNumFaceQP,num_dim,num_dim);
175  weighted_measure = Kokkos::DynRankView<double,ProjectionSpace>("weighted_measure",1,maxNumFaceQP);
176  scratch_space = Kokkos::DynRankView<double,ProjectionSpace>("scratch_space",jacobianSide.span());
177  }
178 }
179 
180 // **********************************************************************
181 template<typename EvalT,typename Traits>
184 {
185 
186  // The coefficients being calculated here in the projection to the face basis
187  // are defined as the integral over the face of the field dotted with the face
188  // normal vector. For a first-order face basis, single point integration is
189  // adequate, so the cubature here just provides the proper weighting.
190  // For higher order, a distinction between "cell" and Gauss points will need
191  // to be made so the field is appropriately projected.
192  const shards::CellTopology & parentCell = *basis->getCellTopology();
193  Intrepid2::DefaultCubatureFactory quadFactory;
194 
195  // Fast Method: One point quadrature on hex mesh. Assumes that the
196  // face area can be computed from two edge normals. This is only
197  // true if the mesh is square or rectangular. Will not work for
198  // paved meshes.
199  if (use_fast_method_on_rectangular_hex_mesh){
201 
202  // Collect the reference face information. For now, do nothing with the quadPts.
203  const unsigned num_faces = parentCell.getFaceCount();
204  std::vector<double> refFaceWt(num_faces, 0.0);
205  for (unsigned f=0; f<num_faces; f++) {
206  faceQuad = quadFactory.create<PHX::exec_space,double,double>(parentCell.getCellTopologyData(2,f), 1);
207  const int numQPoints = faceQuad->getNumPoints();
208  Kokkos::DynRankView<double,PHX::Device> quadWtsFast("quadWts",numQPoints);
209  Kokkos::DynRankView<double,PHX::Device> quadPtsFast("quadPts",numQPoints,num_dim);
210  faceQuad->getCubature(quadPtsFast,quadWtsFast);
211  for (int q=0; q<numQPoints; q++)
212  refFaceWt[f] += quadWtsFast(q);
213  }
214 
215  // Loop over the faces of the workset cells.
216  for (index_t cell = 0; cell < workset.num_cells; ++cell) {
217  for (int p = 0; p < num_pts; ++p) {
218  result(cell,p) = ScalarT(0.0);
219  for (int dim = 0; dim < num_dim; ++dim)
220  result(cell,p) += vector_values[0](cell,p,dim) * normals(cell,p,dim);
221  result(cell,p) *= refFaceWt[p];
222  }
223  }
224 
225  } else {
226  PHX::MDField<double,Cell,panzer::NODE,Dim> vertex_coords = workset.cell_vertex_coordinates;
227  const int subcell_dim = 2;
228  const int numFaces = Teuchos::as<int>(parentCell.getFaceCount());
229  const WorksetDetails & details = workset;
230 
231  // Loop over the faces of the workset cells
232  for (index_t cell = 0; cell < workset.num_cells; ++cell) {
233 
234  // get nodal coordinates for this cell
235  for (int point(0); point < vertex_coords.extent_int(1); ++point)
236  {
237  for (int ict(0); ict < num_dim; ict++)
238  physicalNodes(0,point,ict) = vertex_coords(cell,point,ict);
239  }
240 
241  int faceOrts[6] = {};
242  orientations->at(details.cell_local_ids[cell]).getFaceOrientation(faceOrts, numFaces);
243 
244  // loop over faces
245  for (int p = 0; p < num_pts; ++p){
246  result(cell,p) = ScalarT(0.0);
247 
248  auto ortEdgeTan_U = Kokkos::subview(refEdges, 0, Kokkos::ALL());
249  auto ortEdgeTan_V = Kokkos::subview(refEdges, 1, Kokkos::ALL());
250 
251  // Apply parent cell Jacobian to ref. edge tangent
252  Intrepid2::Orientation::getReferenceFaceTangents(ortEdgeTan_U,
253  ortEdgeTan_V,
254  p,
255  parentCell,
256  faceOrts[p]);
257 
258  // get quad weights/pts on reference 2d cell
259  const shards::CellTopology & subcell = parentCell.getCellTopologyData(subcell_dim,p);
260  const auto& faceQuad = faceQuads[p];
261  faceQuad->getCubature(quadPts,quadWts);
262 
263  // map 2d quad pts to reference cell (3d)
264  Intrepid2::CellTools<ProjectionSpace::execution_space>::mapToReferenceSubcell(refQuadPts, quadPts, subcell_dim, p, parentCell);
265 
266  // Calculate side jacobian
267  Intrepid2::CellTools<ProjectionSpace::execution_space>::setJacobian(jacobianSide, refQuadPts, physicalNodes, parentCell);
268 
269  // Calculate weighted measure at quadrature points
270  Intrepid2::FunctionSpaceTools<ProjectionSpace::execution_space>::computeFaceMeasure(weighted_measure, jacobianSide, quadWts, p, parentCell, scratch_space);
271 
272  // loop over quadrature points
273  for (int qp = 0; qp < faceQuad->getNumPoints(); ++qp) {
274 
275  auto phyEdgeTan_U = Kokkos::subview(phyEdges, 0, Kokkos::ALL());
276  auto phyEdgeTan_V = Kokkos::subview(phyEdges, 1, Kokkos::ALL());
277  auto J = Kokkos::subview(jacobianSide, 0, qp, Kokkos::ALL(), Kokkos::ALL());
278 
279  Intrepid2::Kernels::Serial::matvec_product(phyEdgeTan_U, J, ortEdgeTan_U);
280  Intrepid2::Kernels::Serial::matvec_product(phyEdgeTan_V, J, ortEdgeTan_V);
281 
282  // normal = TanU x TanV
283  ScalarT normal[3];
284  normal[0] = (phyEdgeTan_U(1)*phyEdgeTan_V(2) - phyEdgeTan_U(2)*phyEdgeTan_V(1));
285  normal[1] = (phyEdgeTan_U(2)*phyEdgeTan_V(0) - phyEdgeTan_U(0)*phyEdgeTan_V(2));
286  normal[2] = (phyEdgeTan_U(0)*phyEdgeTan_V(1) - phyEdgeTan_U(1)*phyEdgeTan_V(0));
287 
288  // compute the magnitude of the normal vector
289  ScalarT nnorm(0.0);
290  for(int dim = 0; dim < num_dim; ++dim){
291  nnorm += normal[dim]*normal[dim];
292  }
293  nnorm = std::sqrt(nnorm);
294 
295  // integrate vector dot n
296  // normalize n since jacobian information is factored into both weighted measure and normal
297  for (int dim = 0; dim < num_dim; ++dim)
298  result(cell,p) += weighted_measure(0,qp) * vector_values[qp](cell,p,dim) * normal[dim] / nnorm;
299  }
300  }
301 
302  }
303 
304  } // end else (high order quad)
305 
306 }
307 
308 #endif
Teuchos::RCP< const std::vector< Intrepid2::Orientation > > orientations_
PHX::MDField< ScalarT, Cell, NODE, Dim > gatherFieldNormals
bool use_fast_method_on_rectangular_hex_mesh
If true, a fast algorithm can be used, but this requires a rectangular/square structured hex mesh...
PHX::MDField< const ScalarT, Cell, BASIS, Dim > normals
void evaluateFields(typename Traits::EvalData d)
T & get(const std::string &name, T def_value)
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &vm)
CellCoordArray cell_vertex_coordinates
PHX::MDField< ScalarT, Cell, BASIS > result
PHX::MDField< ScalarT, panzer::Cell, panzer::IP > result
A field that will be used to build up the result of the integral we&#39;re performing.
Teuchos::RCP< const PureBasis > basis
bool isType(const std::string &name) const
#define TEUCHOS_ASSERT(assertion_test)
std::vector< PHX::MDField< const ScalarT, Cell, BASIS, Dim > > vector_values
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_