Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_ProjectToEdges_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_EDGES_IMPL_HPP
44 #define PANZER_PROJECT_TO_EDGES_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_OrientationTools.hpp"
53 
54 #include "Panzer_PureBasis.hpp"
56 #include "Kokkos_ViewFactory.hpp"
57 
58 #include "Teuchos_FancyOStream.hpp"
59 
60 template<typename EvalT,typename Traits>
63  const Teuchos::ParameterList& p)
64 {
65  dof_name = (p.get< std::string >("DOF Name"));
66 
67  if(p.isType< Teuchos::RCP<PureBasis> >("Basis"))
68  basis = p.get< Teuchos::RCP<PureBasis> >("Basis");
69  else
70  basis = p.get< Teuchos::RCP<const PureBasis> >("Basis");
71 
72  quad_degree = 0;
73  if(p.isType<int>("Quadrature Order"))
74  quad_degree = p.get<int>("Quadrature Order");
75 
76  Teuchos::RCP<PHX::DataLayout> basis_layout = basis->functional;
77  Teuchos::RCP<PHX::DataLayout> vector_layout = basis->functional_grad;
78 
79  // some sanity checks
80  TEUCHOS_ASSERT(basis->isVectorBasis());
81 
82  result = PHX::MDField<ScalarT,Cell,BASIS>(dof_name,basis_layout);
83  this->addEvaluatedField(result);
84 
85  tangents = PHX::MDField<const ScalarT,Cell,BASIS,Dim>(dof_name+"_Tangents",vector_layout);
86  this->addDependentField(tangents);
87 
88  if(quad_degree > 0){
89  const shards::CellTopology & parentCell = *basis->getCellTopology();
90  Intrepid2::DefaultCubatureFactory quadFactory;
92  = quadFactory.create<PHX::exec_space,double,double>(parentCell.getCellTopologyData(1,0), quad_degree);
93  int numQPoints = quadRule->getNumPoints();
94 
95  vector_values.resize(numQPoints);
96  for (int qp(0); qp < numQPoints; ++qp)
97  {
98  vector_values[qp] = PHX::MDField<const ScalarT,Cell,BASIS,Dim>(dof_name+"_Vector"+"_qp_"+std::to_string(qp),vector_layout);
99  this->addDependentField(vector_values[qp]);
100  }
101 
102  // setup the orientation field
103  std::string orientationFieldName = basis->name() + " Orientation";
104  dof_orientation = PHX::MDField<const ScalarT,Cell,NODE>(orientationFieldName,basis_layout);
105  this->addDependentField(dof_orientation);
106 
107  gatherFieldTangents = PHX::MDField<ScalarT,Cell,NODE,Dim>(dof_name+"_Tangents",basis->functional_grad);
108  this->addEvaluatedField(gatherFieldTangents);
109 
110  } else {
111  vector_values.resize(1);
112  vector_values[0] = PHX::MDField<const ScalarT,Cell,BASIS,Dim>(dof_name+"_Vector",vector_layout);
113  this->addDependentField(vector_values[0]);
114  }
115 
116  this->setName("Project To Edges");
117 }
118 
119 // **********************************************************************
120 template<typename EvalT,typename Traits>
123  PHX::FieldManager<Traits>& /* fm */)
124 {
125  orientations = d.orientations_;
126 
127  num_pts = vector_values[0].extent(1);
128  num_dim = vector_values[0].extent(2);
129 
130  TEUCHOS_ASSERT(vector_values[0].extent(1) == tangents.extent(1));
131  TEUCHOS_ASSERT(vector_values[0].extent(2) == tangents.extent(2));
132 }
133 
134 // **********************************************************************
135 template<typename EvalT,typename Traits>
138 {
139  const shards::CellTopology & parentCell = *basis->getCellTopology();
140  const int intDegree = basis->order();
141  TEUCHOS_ASSERT(intDegree == 1);
142  Intrepid2::DefaultCubatureFactory quadFactory;
144 
145  // One point quadrature if higher order quadrature not requested
146  if (quad_degree == 0){
147 
148  // this should be edge cubature and collecting weights are always 2.
149  // // Collect the reference edge information. For now, do nothing with the quadPts.
150  // const unsigned num_edges = parentCell.getEdgeCount();
151  // std::vector<double> refEdgeWt(num_edges, 0.0);
152  // for (unsigned e=0; e<num_edges; e++) {
153  // edgeQuad = quadFactory.create<PHX::exec_space,double,double>(parentCell.getCellTopologyData(1,e), intDegree);
154  // const int numQPoints = edgeQuad->getNumPoints();
155  // Kokkos::DynRankView<double,PHX::Device> quadWts("quadWts",numQPoints);
156  // Kokkos::DynRankView<double,PHX::Device> quadPts("quadPts",numQPoints,num_dim);
157  // edgeQuad->getCubature(quadPts,quadWts);
158  // for (int q=0; q<numQPoints; q++)
159  // refEdgeWt[e] += quadWts(q);
160  // }
161 
162  Kokkos::DynRankView<double,PHX::Device> v0("v0", num_dim), v1("v1", num_dim);
163  const unsigned num_edges = parentCell.getEdgeCount();
164  std::vector<double> refEdgeWt(num_edges, 0.0);
165  for (unsigned e=0; e<num_edges; e++) {
166  const auto v0_id = parentCell.getNodeMap(1, e, 0);
167  const auto v1_id = parentCell.getNodeMap(1, e, 1);
168  Intrepid2::CellTools<PHX::exec_space>::getReferenceVertex(v0, parentCell, v0_id);
169  Intrepid2::CellTools<PHX::exec_space>::getReferenceVertex(v1, parentCell, v1_id);
170 
171  double norm = 0.0;
172  for (int d=0;d<num_dim;++d)
173  norm += (v0(d) - v1(d))*(v0(d) - v1(d));
174 
175  refEdgeWt[e] = sqrt(norm);
176  }
177 
178  // Loop over the edges of the workset cells.
179  for (index_t cell = 0; cell < workset.num_cells; ++cell) {
180  for (int p = 0; p < num_pts; ++p) {
181  result(cell,p) = ScalarT(0.0);
182  for (int dim = 0; dim < num_dim; ++dim)
183  result(cell,p) += vector_values[0](cell,p,dim) * tangents(cell,p,dim);
184  result(cell,p) *= refEdgeWt[p];
185  }
186  }
187 
188  } else {
189 
190  TEUCHOS_ASSERT(false); // this doesn't work since we modified the way orientations are handled
191 
192  PHX::MDField<double,Cell,panzer::NODE,Dim> vertex_coords = workset.cell_vertex_coordinates;
193  int subcell_dim = 1;
194 
195  // to compute tangents at qps (copied from GatherTangents)
196  int numEdges = gatherFieldTangents.extent(1);
197  Kokkos::DynRankView<ScalarT,typename PHX::DevLayout<ScalarT>::type,PHX::Device> refEdgeTan = Kokkos::createDynRankView(gatherFieldTangents.get_static_view(),"refEdgeTan",numEdges,num_dim);
198  for(int i=0;i<numEdges;i++) {
199  Kokkos::DynRankView<double,typename PHX::DevLayout<ScalarT>::type,PHX::Device> refEdgeTan_local("refEdgeTan_local",num_dim);
200  Intrepid2::CellTools<PHX::exec_space>::getReferenceEdgeTangent(refEdgeTan_local, i, parentCell);
201 
202  for(int d=0;d<num_dim;d++)
203  refEdgeTan(i,d) = refEdgeTan_local(d);
204  }
205 
206  // Loop over the faces of the workset cells
207  for (index_t cell = 0; cell < workset.num_cells; ++cell) {
208 
209  // get nodal coordinates for this cell
210  Kokkos::DynRankView<double,PHX::Device> physicalNodes("physicalNodes",1,vertex_coords.extent(1),num_dim);
211  for (int point(0); point < vertex_coords.extent_int(1); ++point)
212  {
213  for (int ict(0); ict < num_dim; ict++)
214  physicalNodes(0, point, ict) = vertex_coords(cell, point, ict);
215  }
216 
217  // loop over edges
218  for (int p = 0; p < num_pts; ++p){
219  result(cell,p) = ScalarT(0.0);
220 
221  // get quad weights/pts on reference 2d cell
222  const shards::CellTopology & subcell = parentCell.getCellTopologyData(subcell_dim,p);
223  edgeQuad = quadFactory.create<PHX::exec_space,double,double>(subcell, quad_degree);
225  edgeQuad->getNumPoints() == static_cast<int>(vector_values.size()));
226  Kokkos::DynRankView<double,PHX::Device> quadWts("quadWts",edgeQuad->getNumPoints());
227  Kokkos::DynRankView<double,PHX::Device> quadPts("quadPts",edgeQuad->getNumPoints(),subcell_dim);
228  edgeQuad->getCubature(quadPts,quadWts);
229 
230  // map 1d quad pts to reference cell (3d)
231  Kokkos::DynRankView<double,PHX::Device> refQuadPts("refQuadPts",edgeQuad->getNumPoints(),num_dim);
232  Intrepid2::CellTools<PHX::exec_space>::mapToReferenceSubcell(refQuadPts, quadPts, subcell_dim, p, parentCell);
233 
234 
235  // Calculate side jacobian
236  Kokkos::DynRankView<double,PHX::Device> jacobianSide("jacobianSide", 1, edgeQuad->getNumPoints(), num_dim, num_dim);
237  Intrepid2::CellTools<PHX::exec_space>::setJacobian(jacobianSide, refQuadPts, physicalNodes, parentCell);
238 
239  // Calculate weighted measure at quadrature points
240  Kokkos::DynRankView<double,PHX::Device> weighted_measure("weighted_measure",1,edgeQuad->getNumPoints());
241  Kokkos::DynRankView<double,PHX::Device> scratch_space("scratch_space",jacobianSide.span());
242  Intrepid2::FunctionSpaceTools<PHX::exec_space>::computeEdgeMeasure(weighted_measure, jacobianSide, quadWts, p, parentCell,scratch_space);
243 
244  // loop over quadrature points
245  for (int qp = 0; qp < edgeQuad->getNumPoints(); ++qp) {
246 
247  // get normal vector at quad points
248  std::vector<ScalarT> edgeTan(num_dim);
249  for(int i = 0; i < 3; i++) {
250  edgeTan[i] = (Sacado::ScalarValue<ScalarT>::eval(jacobianSide(0,qp,i,0)*refEdgeTan(p,0))
251  + Sacado::ScalarValue<ScalarT>::eval(jacobianSide(0,qp,i,1)*refEdgeTan(p,1))
252  + Sacado::ScalarValue<ScalarT>::eval(jacobianSide(0,qp,i,2)*refEdgeTan(p,2)))
253  * dof_orientation(cell,p);
254  }
255 
256  // compute the magnitude of the tangent vector
257  ScalarT tnorm(0.0);
258  for(int dim = 0; dim < num_dim; ++dim){
259  tnorm += edgeTan[dim]*edgeTan[dim];
260  }
261  tnorm = std::sqrt(tnorm);
262 
263  // integrate vector dot t
264  // normalize t since jacobian information is factored into both weighted measure and tangent
265  for (int dim = 0; dim < num_dim; ++dim)
266  result(cell,p) += weighted_measure(0,qp) * vector_values[qp](cell,p,dim) * edgeTan[dim] / tnorm;
267 
268  }
269  }
270 
271  }
272 
273  }
274 
275 
276 }
277 
278 #endif
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)
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.
void evaluateFields(typename Traits::EvalData d)
bool isType(const std::string &name) const
#define TEUCHOS_ASSERT(assertion_test)