Panzer  Version of the Day
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 template<typename EvalT,typename Traits>
66  const Teuchos::ParameterList& p)
67 {
68  dof_name = (p.get< std::string >("DOF Name"));
69 
70  if(p.isType< Teuchos::RCP<PureBasis> >("Basis"))
71  basis = p.get< Teuchos::RCP<PureBasis> >("Basis");
72  else
73  basis = p.get< Teuchos::RCP<const PureBasis> >("Basis");
74 
75  quad_degree = 0;
76  if(p.isType<int>("Quadrature Order"))
77  quad_degree = p.get<int>("Quadrature Order");
78 
79  Teuchos::RCP<PHX::DataLayout> basis_layout = basis->functional;
80  Teuchos::RCP<PHX::DataLayout> vector_layout = basis->functional_grad;
81 
82  // some sanity checks
83  TEUCHOS_ASSERT(basis->isVectorBasis());
84 
85  result = PHX::MDField<ScalarT,Cell,BASIS>(dof_name,basis_layout);
86  this->addEvaluatedField(result);
87 
88  normals = PHX::MDField<const ScalarT,Cell,BASIS,Dim>(dof_name+"_Normals",vector_layout);
89  this->addDependentField(normals);
90 
91  if(quad_degree > 0){
92  const shards::CellTopology & parentCell = *basis->getCellTopology();
93  Intrepid2::DefaultCubatureFactory quadFactory;
94  Teuchos::RCP< Intrepid2::Cubature<PHX::exec_space,double,double> > quadRule
95  = quadFactory.create<PHX::exec_space,double,double>(parentCell.getCellTopologyData(2,0), quad_degree);
96  int numQPoints = quadRule->getNumPoints();
97 
98  vector_values.resize(numQPoints);
99  for (int qp(0); qp < numQPoints; ++qp)
100  {
101  vector_values[qp] = PHX::MDField<const ScalarT,Cell,BASIS,Dim>(dof_name+"_Vector"+"_qp_"+std::to_string(qp),vector_layout);
102  this->addDependentField(vector_values[qp]);
103  }
104 
105  gatherFieldNormals = PHX::MDField<ScalarT,Cell,NODE,Dim>(dof_name+"_Normals",basis->functional_grad);
106  this->addEvaluatedField(gatherFieldNormals);
107 
108  } else {
109  vector_values.resize(1);
110  vector_values[0] = PHX::MDField<const ScalarT,Cell,BASIS,Dim>(dof_name+"_Vector",vector_layout);
111  this->addDependentField(vector_values[0]);
112  }
113 
114  this->setName("Project To Faces");
115 }
116 
117 // **********************************************************************
118 template<typename EvalT,typename Traits>
121  PHX::FieldManager<Traits>& /* fm */)
122 {
123  orientations = d.orientations_;
124 
125  num_pts = result.extent(1);
126  num_dim = vector_values[0].extent(2);
127 
128  TEUCHOS_ASSERT(result.extent(1) == normals.extent(1));
129  TEUCHOS_ASSERT(vector_values[0].extent(2) == normals.extent(2));
130 }
131 
132 // **********************************************************************
133 template<typename EvalT,typename Traits>
136 {
137 
138  // The coefficients being calculated here in the projection to the face basis
139  // are defined as the integral over the face of the field dotted with the face
140  // normal vector. For a first-order face basis, single point integration is
141  // adequate, so the cubature here just provides the proper weighting.
142  // For higher order, a distinction between "cell" and Gauss points will need
143  // to be made so the field is appropriately projected.
144  const shards::CellTopology & parentCell = *basis->getCellTopology();
145  Intrepid2::DefaultCubatureFactory quadFactory;
146  Teuchos::RCP< Intrepid2::Cubature<PHX::exec_space,double,double> > faceQuad;
147 
148  // One point quadrature if higher order quadrature not requested
149  if (quad_degree == 0){
150 
151  // Collect the reference face information. For now, do nothing with the quadPts.
152  const unsigned num_faces = parentCell.getFaceCount();
153  std::vector<double> refFaceWt(num_faces, 0.0);
154  for (unsigned f=0; f<num_faces; f++) {
155  faceQuad = quadFactory.create<PHX::exec_space,double,double>(parentCell.getCellTopologyData(2,f), 1);
156  const int numQPoints = faceQuad->getNumPoints();
157  Kokkos::DynRankView<double,PHX::Device> quadWts("quadWts",numQPoints);
158  Kokkos::DynRankView<double,PHX::Device> quadPts("quadPts",numQPoints,num_dim);
159  faceQuad->getCubature(quadPts,quadWts);
160  for (int q=0; q<numQPoints; q++)
161  refFaceWt[f] += quadWts(q);
162  }
163 
164 
165  // Loop over the faces of the workset cells.
166  for (index_t cell = 0; cell < workset.num_cells; ++cell) {
167  for (int p = 0; p < num_pts; ++p) {
168  result(cell,p) = ScalarT(0.0);
169  for (int dim = 0; dim < num_dim; ++dim)
170  result(cell,p) += vector_values[0](cell,p,dim) * normals(cell,p,dim);
171  result(cell,p) *= refFaceWt[p];
172  }
173  }
174 
175  } else {
176  PHX::MDField<double,Cell,panzer::NODE,Dim> vertex_coords = workset.cell_vertex_coordinates;
177  int subcell_dim = 2;
178 
179  int cellDim = parentCell.getDimension();
180  int numFaces = Teuchos::as<int>(parentCell.getFaceCount());
181 
182  // allocate space that is sized correctly for AD
183  auto refEdges = Kokkos::createDynRankView(result.get_static_view(),"ref_edges", 2, cellDim);
184  auto phyEdges = Kokkos::createDynRankView(result.get_static_view(),"phy_edges", 2, cellDim);
185 
186  const WorksetDetails & details = workset;
187 
188  // Loop over the faces of the workset cells
189  for (index_t cell = 0; cell < workset.num_cells; ++cell) {
190 
191  // get nodal coordinates for this cell
192  Kokkos::DynRankView<double,PHX::Device> physicalNodes("physicalNodes",1,vertex_coords.extent(1),num_dim);
193  for (int point(0); point < vertex_coords.extent_int(1); ++point)
194  {
195  for (int ict(0); ict < num_dim; ict++)
196  physicalNodes(0,point,ict) = vertex_coords(cell,point,ict);
197  }
198 
199  int faceOrts[6] = {};
200  orientations->at(details.cell_local_ids[cell]).getFaceOrientation(faceOrts, numFaces);
201 
202  // loop over faces
203  for (int p = 0; p < num_pts; ++p){
204  result(cell,p) = ScalarT(0.0);
205 
206  auto ortEdgeTan_U = Kokkos::subview(refEdges, 0, Kokkos::ALL());
207  auto ortEdgeTan_V = Kokkos::subview(refEdges, 1, Kokkos::ALL());
208 
209  // Apply parent cell Jacobian to ref. edge tangent
210  Intrepid2::Orientation::getReferenceFaceTangents(ortEdgeTan_U,
211  ortEdgeTan_V,
212  p,
213  parentCell,
214  faceOrts[p]);
215 
216  // get quad weights/pts on reference 2d cell
217  const shards::CellTopology & subcell = parentCell.getCellTopologyData(subcell_dim,p);
218  faceQuad = quadFactory.create<PHX::exec_space,double,double>(subcell, quad_degree);
219  TEUCHOS_ASSERT(
220  faceQuad->getNumPoints() == static_cast<int>(vector_values.size()));
221  Kokkos::DynRankView<double,PHX::Device> quadWts("quadWts",faceQuad->getNumPoints());
222  Kokkos::DynRankView<double,PHX::Device> quadPts("quadPts",faceQuad->getNumPoints(),subcell_dim);
223  faceQuad->getCubature(quadPts,quadWts);
224 
225  // map 2d quad pts to reference cell (3d)
226  Kokkos::DynRankView<double,PHX::Device> refQuadPts("refQuadPts",faceQuad->getNumPoints(),num_dim);
227  Intrepid2::CellTools<PHX::exec_space>::mapToReferenceSubcell(refQuadPts, quadPts, subcell_dim, p, parentCell);
228 
229  // Calculate side jacobian
230  Kokkos::DynRankView<double,PHX::Device> jacobianSide("jacobianSide", 1, faceQuad->getNumPoints(), num_dim, num_dim);
231  Intrepid2::CellTools<PHX::exec_space>::setJacobian(jacobianSide, refQuadPts, physicalNodes, parentCell);
232 
233  // Calculate weighted measure at quadrature points
234  Kokkos::DynRankView<double,PHX::Device> weighted_measure("weighted_measure",1,faceQuad->getNumPoints());
235  Kokkos::DynRankView<double,PHX::Device> scratch_space("scratch_space",jacobianSide.span());
236  Intrepid2::FunctionSpaceTools<PHX::exec_space>::computeFaceMeasure(weighted_measure, jacobianSide, quadWts, p, parentCell, scratch_space);
237 
238  // loop over quadrature points
239  for (int qp = 0; qp < faceQuad->getNumPoints(); ++qp) {
240 
241  auto phyEdgeTan_U = Kokkos::subview(phyEdges, 0, Kokkos::ALL());
242  auto phyEdgeTan_V = Kokkos::subview(phyEdges, 1, Kokkos::ALL());
243  auto J = Kokkos::subview(jacobianSide, 0, qp, Kokkos::ALL(), Kokkos::ALL());
244 
245  Intrepid2::Kernels::Serial::matvec_product(phyEdgeTan_U, J, ortEdgeTan_U);
246  Intrepid2::Kernels::Serial::matvec_product(phyEdgeTan_V, J, ortEdgeTan_V);
247 
248  // normal = TanU x TanV
249  std::vector<ScalarT> normal(3,0.0);
250  normal[0] = (phyEdgeTan_U(1)*phyEdgeTan_V(2) - phyEdgeTan_U(2)*phyEdgeTan_V(1));
251  normal[1] = (phyEdgeTan_U(2)*phyEdgeTan_V(0) - phyEdgeTan_U(0)*phyEdgeTan_V(2));
252  normal[2] = (phyEdgeTan_U(0)*phyEdgeTan_V(1) - phyEdgeTan_U(1)*phyEdgeTan_V(0));
253 
254  // compute the magnitude of the normal vector
255  ScalarT nnorm(0.0);
256  for(int dim = 0; dim < num_dim; ++dim){
257  nnorm += normal[dim]*normal[dim];
258  }
259  nnorm = std::sqrt(nnorm);
260 
261  // integrate vector dot n
262  // normalize n since jacobian information is factored into both weighted measure and normal
263  for (int dim = 0; dim < num_dim; ++dim)
264  result(cell,p) += weighted_measure(0,qp) * vector_values[qp](cell,p,dim) * normal[dim] / nnorm;
265  }
266  }
267 
268  }
269 
270  } // end else (high order quad)
271 
272 }
273 
274 #endif
Teuchos::RCP< const std::vector< Intrepid2::Orientation > > orientations_
void evaluateFields(typename Traits::EvalData d)
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &vm)
CellCoordArray cell_vertex_coordinates
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.
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.