Panzer  Version of the Day
Panzer_GatherTangents_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_TANGENTS_IMPL_HPP
44 #define PANZER_GATHER_TANGENTS_IMPL_HPP
45 
46 #include "Teuchos_Assert.hpp"
47 #include "Phalanx_DataLayout.hpp"
48 
49 #include "Panzer_PureBasis.hpp"
51 #include "Kokkos_ViewFactory.hpp"
52 
53 #include "Teuchos_FancyOStream.hpp"
54 
55 template<typename EvalT,typename Traits>
58  const Teuchos::ParameterList& p)
59 {
60  dof_name = (p.get< std::string >("DOF Name"));
61 
62  if(p.isType< Teuchos::RCP<PureBasis> >("Basis"))
63  basis = p.get< Teuchos::RCP<PureBasis> >("Basis");
64  else
65  basis = p.get< Teuchos::RCP<const PureBasis> >("Basis");
66 
67  pointRule = p.get<Teuchos::RCP<const PointRule> >("Point Rule");
68 
70  Teuchos::RCP<PHX::DataLayout> vector_layout_vector = basis->functional_grad;
71 
72  // some sanity checks
74 
75  // setup the orientation field
76  std::string orientationFieldName = basis->name() + " Orientation";
77  dof_orientation = PHX::MDField<const ScalarT,Cell,NODE>(orientationFieldName, basis_layout);
78 
79  // setup all basis fields that are required
81 
82  // setup all fields to be evaluated and constructed
84 
85  // the field manager will allocate all of these field
86  this->addDependentField(dof_orientation);
87  this->addDependentField(pointValues.jac);
88 
89  gatherFieldTangents = PHX::MDField<ScalarT,Cell,NODE,Dim>(dof_name+"_Tangents",vector_layout_vector);
90  this->addEvaluatedField(gatherFieldTangents);
91 
92  this->setName("Gather Tangents");
93 }
94 
95 // **********************************************************************
96 template<typename EvalT,typename Traits>
100 {
101  // setup the field data object
102  this->utils.setFieldData(gatherFieldTangents,fm);
103  this->utils.setFieldData(dof_orientation,fm);
104  this->utils.setFieldData(pointValues.jac,fm);
105 
106  edgeTan = Kokkos::createDynRankView(gatherFieldTangents.get_static_view(),"edgeTan",gatherFieldTangents.dimension(0),gatherFieldTangents.dimension(1),gatherFieldTangents.dimension(2));
107 }
108 
109 // **********************************************************************
110 template<typename EvalT,typename Traits>
113 {
114 
115  if(workset.num_cells<=0)
116  return;
117  else {
118  const shards::CellTopology & parentCell = *basis->getCellTopology();
119  int cellDim = parentCell.getDimension();
120  int numEdges = gatherFieldTangents.dimension(1);
121 
122  refEdgeTan = Kokkos::createDynRankView(gatherFieldTangents.get_static_view(),"refEdgeTan",numEdges,cellDim);
123 
124  for(int i=0;i<numEdges;i++) {
125  Kokkos::DynRankView<double,PHX::Device> refEdgeTan_local("refEdgeTan_local",cellDim);
126  Intrepid2::CellTools<double>::getReferenceEdgeTangent(refEdgeTan_local, i, parentCell);
127 
128  for(int d=0;d<cellDim;d++)
129  refEdgeTan(i,d) = refEdgeTan_local(d);
130  }
131 
132  // Loop over workset faces and edge points
133  for(index_t c=0;c<workset.num_cells;c++) {
134  for(int pt = 0; pt < numEdges; pt++) {
135 
136  // Apply parent cell Jacobian to ref. edge tangent
137  for(int i = 0; i < cellDim; i++) {
138  edgeTan(c, pt, i) = 0.0;
139  for(int j = 0; j < cellDim; j++){
140  edgeTan(c, pt, i) += pointValues.jac(c, pt, i, j)*refEdgeTan(pt,j);
141  }// for j
142  }// for i
143  }// for pt
144  }// for pCell
145 
146  // Multiply tangent by orientation
147  for(index_t c=0;c<workset.num_cells;c++) {
148  for(int b=0;b<gatherFieldTangents.extent_int(1);b++) {
149  for(int d=0;d<gatherFieldTangents.extent_int(2);d++) {
150  gatherFieldTangents(c,b,d) = edgeTan(c,b,d)*dof_orientation(c,b);
151  }
152  }
153  }
154  }
155 
156 }
157 
158 #endif
std::string name() const
A unique key that is the combination of the basis type and basis order.
T & get(const std::string &name, T def_value)
PointValues2< ScalarT, PHX::MDField > pointValues
Interpolates basis DOF values to IP DOF values.
bool isVectorBasis() const
Kokkos::DynRankView< ScalarT, PHX::Device > edgeTan
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.
PHX::MDField< const ScalarT, Cell, BASIS > dof_orientation
Teuchos::RCP< const panzer::PointRule > pointRule
Teuchos::RCP< PHX::DataLayout > functional_grad
<Cell,Basis,Dim>
bool isType(const std::string &name) const
Kokkos::DynRankView< ScalarT, PHX::Device > refEdgeTan
Teuchos::RCP< const shards::CellTopology > getCellTopology() const
Teuchos::RCP< const panzer::PureBasis > basis
Interpolates basis DOF values to IP DOF values.
Array< Scalar, Cell, IP, Dim, Dim, void, void, void, void > jac
void setupArrays(const Teuchos::RCP< const panzer::PointRule > &pr, const ArrayFactory &af)
Sizes/allocates memory for arrays.
const std::string & getName() const
#define TEUCHOS_ASSERT(assertion_test)
Teuchos::RCP< PHX::DataLayout > functional
<Cell,Basis> or <Cell,Basis>
void evaluateFields(typename Traits::EvalData d)
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &vm)