Panzer  Version of the Day
Panzer_DOF_PointValues_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_DOF_IMPL_HPP
44 #define PANZER_DOF_IMPL_HPP
45 
46 #include <algorithm>
48 #include "Panzer_BasisIRLayout.hpp"
51 #include "Panzer_DOF_Functors.hpp"
52 
53 #include "Intrepid2_FunctionSpaceTools.hpp"
54 
55 namespace panzer {
56 
57 //**********************************************************************
58 //* DOF_PointValues evaluator
59 //**********************************************************************
60 
61 //**********************************************************************
62 // MOST EVALUATION TYPES
63 //**********************************************************************
64 
65 //**********************************************************************
66 template<typename EvalT, typename TRAITS>
68 DOF_PointValues(const Teuchos::ParameterList & p)
69 {
70  const std::string fieldName = p.get<std::string>("Name");
71  basis = p.get< Teuchos::RCP<const PureBasis> >("Basis");
72  Teuchos::RCP<const PointRule> pointRule = p.get< Teuchos::RCP<const PointRule> >("Point Rule");
73  is_vector_basis = basis->isVectorBasis();
74 
75  std::string evalName = fieldName+"_"+pointRule->getName();
76  if(p.isType<bool>("Use DOF Name")) {
77  if(p.get<bool>("Use DOF Name"))
78  evalName = fieldName;
79  }
80 
81  dof_basis = PHX::MDField<const ScalarT,Cell,Point>(fieldName, basis->functional);
82 
83  this->addDependentField(dof_basis);
84 
85  // setup all basis fields that are required
86  Teuchos::RCP<BasisIRLayout> layout = Teuchos::rcp(new BasisIRLayout(basis,*pointRule));
87  basisValues = Teuchos::rcp(new BasisValues2<ScalarT>(basis->name()+"_"+pointRule->getName()+"_"));
88  basisValues->setupArrays(layout,false);
89 
90  // the field manager will allocate all of these field
91  // swap between scalar basis value, or vector basis value
92  if(basis->isScalarBasis()) {
93  dof_ip_scalar = PHX::MDField<ScalarT,Cell,Point>(
94  evalName,
95  pointRule->dl_scalar);
96  this->addEvaluatedField(dof_ip_scalar);
97  constBasisRefScalar_ = basisValues->basis_ref_scalar;
98  constBasisScalar_ = basisValues->basis_scalar;
99  this->addDependentField(constBasisRefScalar_);
100  this->addDependentField(constBasisScalar_);
101  }
102  else if(basis->isVectorBasis()) {
103  dof_ip_vector = PHX::MDField<ScalarT,Cell,Point,Dim>(
104  evalName,
105  pointRule->dl_vector);
106  this->addEvaluatedField(dof_ip_vector);
107  constBasisRefVector_ = basisValues->basis_ref_vector;
108  constBasisVector_ = basisValues->basis_vector;
109  this->addDependentField(constBasisRefVector_);
110  this->addDependentField(constBasisVector_);
111  }
112  else
113  { TEUCHOS_ASSERT(false); }
114 
115  std::string n = "DOF_PointValues: " + dof_basis.fieldTag().name();
116  this->setName(n);
117 }
118 
119 //**********************************************************************
120 template<typename EvalT, typename TRAITS>
122 postRegistrationSetup(typename TRAITS::SetupData /* sd */,
124 {
125  if(!is_vector_basis) {
126  this->utils.setFieldData(basisValues->basis_ref_scalar,fm);
127  this->utils.setFieldData(basisValues->basis_scalar,fm);
128  }
129  else {
130  this->utils.setFieldData(basisValues->basis_ref_vector,fm);
131  this->utils.setFieldData(basisValues->basis_vector,fm);
132  }
133 }
134 
135 //**********************************************************************
136 template<typename EvalT, typename TRAITS>
138 evaluateFields(typename TRAITS::EvalData workset)
139 {
140  if(is_vector_basis) {
141  int spaceDim = basisValues->basis_vector.extent(3);
142  if(spaceDim==3) {
143  dof_functors::EvaluateDOFWithSens_Vector<ScalarT,typename BasisValues2<ScalarT>::Array_CellBasisIPDim,3> functor(dof_basis,dof_ip_vector,basisValues->basis_vector);
144  Kokkos::parallel_for(workset.num_cells,functor);
145  }
146  else {
147  dof_functors::EvaluateDOFWithSens_Vector<ScalarT,typename BasisValues2<ScalarT>::Array_CellBasisIPDim,2> functor(dof_basis,dof_ip_vector,basisValues->basis_vector);
148  Kokkos::parallel_for(workset.num_cells,functor);
149  }
150  }
151  else {
152  dof_functors::EvaluateDOFWithSens_Scalar<ScalarT,typename BasisValues2<ScalarT>::Array_CellBasisIP> functor(dof_basis,dof_ip_scalar,basisValues->basis_scalar);
153  Kokkos::parallel_for(workset.num_cells,functor);
154  }
155 }
156 
157 //**********************************************************************
158 
159 //**********************************************************************
160 // JACOBIAN EVALUATION TYPES
161 //**********************************************************************
162 
163 //**********************************************************************
164 template<typename TRAITS>
166 DOF_PointValues(const Teuchos::ParameterList & p)
167 {
168  const std::string fieldName = p.get<std::string>("Name");
169  basis = p.get< Teuchos::RCP<const PureBasis> >("Basis");
170  Teuchos::RCP<const PointRule> pointRule = p.get< Teuchos::RCP<const PointRule> >("Point Rule");
171  is_vector_basis = basis->isVectorBasis();
172 
173  if(p.isType<Teuchos::RCP<const std::vector<int> > >("Jacobian Offsets Vector")) {
174  const std::vector<int> & offsets = *p.get<Teuchos::RCP<const std::vector<int> > >("Jacobian Offsets Vector");
175 
176  // allocate and copy offsets vector to Kokkos array
177  offsets_array = Kokkos::View<int*,PHX::Device>("offsets",offsets.size());
178  for(std::size_t i=0;i<offsets.size();i++)
179  offsets_array(i) = offsets[i];
180 
181  accelerate_jacobian = true; // short cut for identity matrix
182  }
183  else
184  accelerate_jacobian = false; // don't short cut for identity matrix
185 
186  std::string evalName = fieldName+"_"+pointRule->getName();
187  if(p.isType<bool>("Use DOF Name")) {
188  if(p.get<bool>("Use DOF Name"))
189  evalName = fieldName;
190  }
191 
192  dof_basis = PHX::MDField<const ScalarT,Cell,Point>(fieldName, basis->functional);
193 
194  this->addDependentField(dof_basis);
195 
196  // setup all basis fields that are required
197  Teuchos::RCP<BasisIRLayout> layout = Teuchos::rcp(new BasisIRLayout(basis,*pointRule));
198  basisValues = Teuchos::rcp(new BasisValues2<ScalarT>(basis->name()+"_"+pointRule->getName()+"_"));
199  basisValues->setupArrays(layout,false);
200 
201  // the field manager will allocate all of these field
202  // swap between scalar basis value, or vector basis value
203  if(basis->isScalarBasis()) {
204  dof_ip_scalar = PHX::MDField<ScalarT,Cell,Point>(
205  evalName,
206  pointRule->dl_scalar);
207  this->addEvaluatedField(dof_ip_scalar);
208  constBasisRefScalar_ = basisValues->basis_ref_scalar;
209  constBasisScalar_ = basisValues->basis_scalar;
210  this->addDependentField(constBasisRefScalar_);
211  this->addDependentField(constBasisScalar_);
212  }
213  else if(basis->isVectorBasis()) {
214  dof_ip_vector = PHX::MDField<ScalarT,Cell,Point,Dim>(
215  evalName,
216  pointRule->dl_vector);
217  this->addEvaluatedField(dof_ip_vector);
218  constBasisRefVector_ = basisValues->basis_ref_vector;
219  constBasisVector_ = basisValues->basis_vector;
220  this->addDependentField(constBasisRefVector_);
221  this->addDependentField(constBasisVector_);
222  }
223  else
224  { TEUCHOS_ASSERT(false); }
225 
226  std::string n = "DOF_PointValues: " + dof_basis.fieldTag().name() + " Jacobian";
227  this->setName(n);
228 }
229 
230 //**********************************************************************
231 template<typename TRAITS>
233 postRegistrationSetup(typename TRAITS::SetupData /* sd */,
235 {
236  if(!is_vector_basis) {
237  this->utils.setFieldData(basisValues->basis_ref_scalar,fm);
238  this->utils.setFieldData(basisValues->basis_scalar,fm);
239  }
240  else {
241  this->utils.setFieldData(basisValues->basis_ref_vector,fm);
242  this->utils.setFieldData(basisValues->basis_vector,fm);
243  }
244 }
245 
246 //**********************************************************************
247 template<typename TRAITS>
249 evaluateFields(typename TRAITS::EvalData workset)
250 {
251  if(is_vector_basis) {
252  if(accelerate_jacobian) {
253  int spaceDim = basisValues->basis_vector.extent(3);
254  if(spaceDim==3) {
256  Kokkos::parallel_for(workset.num_cells,functor);
257  }
258  else {
260  Kokkos::parallel_for(workset.num_cells,functor);
261  }
262  }
263  else {
264  int spaceDim = basisValues->basis_vector.extent(3);
265  if(spaceDim==3) {
267  Kokkos::parallel_for(workset.num_cells,functor);
268  }
269  else {
271  Kokkos::parallel_for(workset.num_cells,functor);
272  }
273  }
274  }
275  else {
276  if(accelerate_jacobian) {
278  Kokkos::parallel_for(workset.num_cells,functor);
279  }
280  else {
282  Kokkos::parallel_for(workset.num_cells,functor);
283  }
284  }
285 }
286 
287 }
288 
289 #endif
Teuchos::RCP< BasisValues2< ScalarT > > basisValues
void evaluateFields(typename TRAITS::EvalData d)
PHX::MDField< const ScalarT, BASIS, IP, void, void, void, void, void, void > constBasisRefScalar_
PHX::MDField< ScalarT, Cell, Point, Dim > dof_ip_vector
PHX::MDField< const ScalarT, BASIS, IP, Dim, void, void, void, void, void > constBasisRefVector_
PHX::MDField< const ScalarT, Cell, Point > dof_basis
Teuchos::RCP< const PureBasis > basis
void postRegistrationSetup(typename TRAITS::SetupData d, PHX::FieldManager< TRAITS > &fm)
PHX::MDField< const ScalarT, Cell, BASIS, IP, void, void, void, void, void > constBasisScalar_
PHX::MDField< const ScalarT, Cell, BASIS, IP, Dim, void, void, void, void > constBasisVector_
PHX::MDField< ScalarT, Cell, Point > dof_ip_scalar
DOF_PointValues(const Teuchos::ParameterList &p)
Kokkos::View< const int *, PHX::Device > offsets