Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_GatherNormals_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_GATHER_NORMALS_IMPL_HPP
44 #define PANZER_GATHER_NORMALS_IMPL_HPP
45 
46 #include "Teuchos_Assert.hpp"
47 #include "Phalanx_DataLayout.hpp"
48 
49 #include "Panzer_PureBasis.hpp"
50 #include "Kokkos_ViewFactory.hpp"
51 
52 #include "Intrepid2_Kernels.hpp"
53 #include "Intrepid2_CellTools.hpp"
54 #include "Intrepid2_OrientationTools.hpp"
55 
56 #include "Teuchos_FancyOStream.hpp"
57 
58 template<typename EvalT,typename Traits>
61  const Teuchos::ParameterList& p)
62 {
63  dof_name = (p.get< std::string >("DOF Name"));
64 
65  if(p.isType< Teuchos::RCP<PureBasis> >("Basis"))
66  basis = p.get< Teuchos::RCP<PureBasis> >("Basis");
67  else
68  basis = p.get< Teuchos::RCP<const PureBasis> >("Basis");
69 
70  pointRule = p.get<Teuchos::RCP<const PointRule> >("Point Rule");
71 
72  Teuchos::RCP<PHX::DataLayout> basis_layout = basis->functional;
73  Teuchos::RCP<PHX::DataLayout> vector_layout_vector = basis->functional_grad;
74 
75  // some sanity checks
76  TEUCHOS_ASSERT(basis->isVectorBasis());
77 
78  // setup the orientation field
79  std::string orientationFieldName = basis->name() + " Orientation";
80  // setup all fields to be evaluated and constructed
81  pointValues = panzer::PointValues2<double> (pointRule->getName()+"_",false);
82  pointValues.setupArrays(pointRule);
83 
84  // the field manager will allocate all of these field
85  constJac_ = pointValues.jac;
86  this->addDependentField(constJac_);
87 
88  gatherFieldNormals = PHX::MDField<ScalarT,Cell,NODE,Dim>(dof_name+"_Normals",vector_layout_vector);
89  this->addEvaluatedField(gatherFieldNormals);
90 
91  this->setName("Gather Normals");
92 }
93 
94 // **********************************************************************
95 template<typename EvalT,typename Traits>
99 {
100  orientations = d.orientations_;
101  this->utils.setFieldData(pointValues.jac,fm);
102  faceNormal = Kokkos::createDynRankView(gatherFieldNormals.get_static_view(),
103  "faceNormal",
104  gatherFieldNormals.extent(0),
105  gatherFieldNormals.extent(1),
106  gatherFieldNormals.extent(2));
107 }
108 
109 // **********************************************************************
110 template<typename EvalT,typename Traits>
113 {
114 
115  if(workset.num_cells<=0)
116  return;
117 
118  const shards::CellTopology & parentCell = *basis->getCellTopology();
119  int cellDim = parentCell.getDimension();
120  int numFaces = gatherFieldNormals.extent(1);
121 
122  // allocate space that is sized correctly for AD
123  auto refEdges = Kokkos::createDynRankView(gatherFieldNormals.get_static_view(),"ref_edges", 2, cellDim);
124  auto phyEdges = Kokkos::createDynRankView(gatherFieldNormals.get_static_view(),"phy_edges", 2, cellDim);
125 
126  const WorksetDetails & details = workset;
127  const auto worksetJacobians = pointValues.jac.get_view();
128 
129  // Loop over workset faces and edge points
130  for(index_t c=0;c<workset.num_cells;c++) {
131  int faceOrts[6] = {};
132  orientations->at(details.cell_local_ids[c]).getFaceOrientation(faceOrts, numFaces);
133 
134  for(int pt = 0; pt < numFaces; pt++) {
135  auto ortEdgeTan_U = Kokkos::subview(refEdges, 0, Kokkos::ALL());
136  auto ortEdgeTan_V = Kokkos::subview(refEdges, 1, Kokkos::ALL());
137 
138  // Apply parent cell Jacobian to ref. edge tangent
139  Intrepid2::Orientation::getReferenceFaceTangents(ortEdgeTan_U,
140  ortEdgeTan_V,
141  pt,
142  parentCell,
143  faceOrts[pt]);
144 
145  auto phyEdgeTan_U = Kokkos::subview(phyEdges, 0, Kokkos::ALL());
146  auto phyEdgeTan_V = Kokkos::subview(phyEdges, 1, Kokkos::ALL());
147  auto J = Kokkos::subview(worksetJacobians, c, pt, Kokkos::ALL(), Kokkos::ALL());
148 
149  Intrepid2::Kernels::Serial::matvec_product(phyEdgeTan_U, J, ortEdgeTan_U);
150  Intrepid2::Kernels::Serial::matvec_product(phyEdgeTan_V, J, ortEdgeTan_V);
151 
152  // take the cross product of the two vectors
153  gatherFieldNormals(c,pt,0) = (phyEdgeTan_U(1)*phyEdgeTan_V(2) - phyEdgeTan_U(2)*phyEdgeTan_V(1));
154  gatherFieldNormals(c,pt,1) = (phyEdgeTan_U(2)*phyEdgeTan_V(0) - phyEdgeTan_U(0)*phyEdgeTan_V(2));
155  gatherFieldNormals(c,pt,2) = (phyEdgeTan_U(0)*phyEdgeTan_V(1) - phyEdgeTan_U(1)*phyEdgeTan_V(0));
156  }
157  }
158 
159 }
160 
161 #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 setupArrays(const Teuchos::RCP< const panzer::PointRule > &pr)
Sizes/allocates memory for arrays.
void evaluateFields(typename Traits::EvalData d)
bool isType(const std::string &name) const
#define TEUCHOS_ASSERT(assertion_test)
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &vm)