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<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 
98  this->addDependentField(basisValues->basis_ref_scalar);
99  this->addDependentField(basisValues->basis_scalar);
100  }
101  else if(basis->isVectorBasis()) {
102  dof_ip_vector = PHX::MDField<ScalarT,Cell,Point,Dim>(
103  evalName,
104  pointRule->dl_vector);
105  this->addEvaluatedField(dof_ip_vector);
106 
107  this->addDependentField(basisValues->basis_ref_vector);
108  this->addDependentField(basisValues->basis_vector);
109  }
110  else
111  { TEUCHOS_ASSERT(false); }
112 
113  std::string n = "DOF_PointValues: " + dof_basis.fieldTag().name();
114  this->setName(n);
115 }
116 
117 //**********************************************************************
118 template<typename EvalT, typename TRAITS>
120 postRegistrationSetup(typename TRAITS::SetupData sd,
122 {
123  this->utils.setFieldData(dof_basis,fm);
124 
125  if(!is_vector_basis) {
126  this->utils.setFieldData(dof_ip_scalar,fm);
127 
128  // setup the pointers for the basis values data structure
129  this->utils.setFieldData(basisValues->basis_ref_scalar,fm);
130  this->utils.setFieldData(basisValues->basis_scalar,fm);
131  }
132  else {
133  this->utils.setFieldData(dof_ip_vector,fm);
134 
135  // setup the pointers for the basis values data structure
136  this->utils.setFieldData(basisValues->basis_ref_vector,fm);
137  this->utils.setFieldData(basisValues->basis_vector,fm);
138  }
139 }
140 
141 //**********************************************************************
142 template<typename EvalT, typename TRAITS>
144 evaluateFields(typename TRAITS::EvalData workset)
145 {
146  // evaluateDOF_withSens(dof_basis,dof_ip,dof_orientation,is_vector_basis,workset.num_cells,basisValues.basis);
147 
148  if(is_vector_basis) {
149  int spaceDim = basisValues->basis_vector.dimension(3);
150  if(spaceDim==3) {
152  Kokkos::parallel_for(workset.num_cells,functor);
153  }
154  else {
156  Kokkos::parallel_for(workset.num_cells,functor);
157  }
158  }
159  else {
161  Kokkos::parallel_for(workset.num_cells,functor);
162  }
163 }
164 
165 //**********************************************************************
166 
167 //**********************************************************************
168 // JACOBIAN EVALUATION TYPES
169 //**********************************************************************
170 
171 //**********************************************************************
172 template<typename TRAITS>
174 DOF_PointValues(const Teuchos::ParameterList & p)
175 {
176  const std::string fieldName = p.get<std::string>("Name");
177  basis = p.get< Teuchos::RCP<const PureBasis> >("Basis");
178  Teuchos::RCP<const PointRule> pointRule = p.get< Teuchos::RCP<const PointRule> >("Point Rule");
179  is_vector_basis = basis->isVectorBasis();
180 
181  if(p.isType<Teuchos::RCP<const std::vector<int> > >("Jacobian Offsets Vector")) {
182  const std::vector<int> & offsets = *p.get<Teuchos::RCP<const std::vector<int> > >("Jacobian Offsets Vector");
183 
184  // allocate and copy offsets vector to Kokkos array
185  offsets_array = Kokkos::View<int*,PHX::Device>("offsets",offsets.size());
186  for(std::size_t i=0;i<offsets.size();i++)
187  offsets_array(i) = offsets[i];
188 
189  accelerate_jacobian = true; // short cut for identity matrix
190  }
191  else
192  accelerate_jacobian = false; // don't short cut for identity matrix
193 
194  std::string evalName = fieldName+"_"+pointRule->getName();
195  if(p.isType<bool>("Use DOF Name")) {
196  if(p.get<bool>("Use DOF Name"))
197  evalName = fieldName;
198  }
199 
200  dof_basis = PHX::MDField<ScalarT,Cell,Point>(fieldName, basis->functional);
201 
202  this->addDependentField(dof_basis);
203 
204  // setup all basis fields that are required
205  Teuchos::RCP<BasisIRLayout> layout = Teuchos::rcp(new BasisIRLayout(basis,*pointRule));
206  basisValues = Teuchos::rcp(new BasisValues2<ScalarT>(basis->name()+"_"+pointRule->getName()+"_"));
207  basisValues->setupArrays(layout,false);
208 
209  // the field manager will allocate all of these field
210  // swap between scalar basis value, or vector basis value
211  if(basis->isScalarBasis()) {
212  dof_ip_scalar = PHX::MDField<ScalarT,Cell,Point>(
213  evalName,
214  pointRule->dl_scalar);
215  this->addEvaluatedField(dof_ip_scalar);
216 
217  this->addDependentField(basisValues->basis_ref_scalar);
218  this->addDependentField(basisValues->basis_scalar);
219  }
220  else if(basis->isVectorBasis()) {
221  dof_ip_vector = PHX::MDField<ScalarT,Cell,Point,Dim>(
222  evalName,
223  pointRule->dl_vector);
224  this->addEvaluatedField(dof_ip_vector);
225 
226  this->addDependentField(basisValues->basis_ref_vector);
227  this->addDependentField(basisValues->basis_vector);
228  }
229  else
230  { TEUCHOS_ASSERT(false); }
231 
232  std::string n = "DOF_PointValues: " + dof_basis.fieldTag().name() + " Jacobian";
233  this->setName(n);
234 }
235 
236 //**********************************************************************
237 template<typename TRAITS>
239 postRegistrationSetup(typename TRAITS::SetupData sd,
241 {
242  this->utils.setFieldData(dof_basis,fm);
243 
244  if(!is_vector_basis) {
245  this->utils.setFieldData(dof_ip_scalar,fm);
246 
247  // setup the pointers for the basis values data structure
248  this->utils.setFieldData(basisValues->basis_ref_scalar,fm);
249  this->utils.setFieldData(basisValues->basis_scalar,fm);
250  }
251  else {
252  this->utils.setFieldData(dof_ip_vector,fm);
253 
254  // setup the pointers for the basis values data structure
255  this->utils.setFieldData(basisValues->basis_ref_vector,fm);
256  this->utils.setFieldData(basisValues->basis_vector,fm);
257  }
258 }
259 
260 //**********************************************************************
261 template<typename TRAITS>
263 evaluateFields(typename TRAITS::EvalData workset)
264 {
265  if(is_vector_basis) {
266  if(accelerate_jacobian) {
267  int spaceDim = basisValues->basis_vector.dimension(3);
268  if(spaceDim==3) {
270  Kokkos::parallel_for(workset.num_cells,functor);
271  }
272  else {
274  Kokkos::parallel_for(workset.num_cells,functor);
275  }
276  }
277  else {
278  int spaceDim = basisValues->basis_vector.dimension(3);
279  if(spaceDim==3) {
281  Kokkos::parallel_for(workset.num_cells,functor);
282  }
283  else {
285  Kokkos::parallel_for(workset.num_cells,functor);
286  }
287  }
288  }
289  else {
290  if(accelerate_jacobian) {
292  Kokkos::parallel_for(workset.num_cells,functor);
293  }
294  else {
296  Kokkos::parallel_for(workset.num_cells,functor);
297  }
298  }
299 }
300 
301 }
302 
303 #endif
Teuchos::RCP< BasisValues2< ScalarT > > basisValues
void evaluateFields(typename TRAITS::EvalData d)
PHX::MDField< ScalarT, Cell, Point, Dim > dof_ip_vector
PHX::MDField< ScalarT, Cell, Point > dof_basis
Teuchos::RCP< const panzer::PointRule > pointRule
Teuchos::RCP< const PureBasis > basis
void postRegistrationSetup(typename TRAITS::SetupData d, PHX::FieldManager< TRAITS > &fm)
Teuchos::RCP< BasisValues2< ScalarT > > basisValues
Teuchos::RCP< const panzer::PureBasis > basis
Interpolates basis DOF values to IP DOF values.
PHX::MDField< ScalarT, Cell, Point > dof_ip_scalar
DOF_PointValues(const Teuchos::ParameterList &p)
Kokkos::View< const int *, PHX::Device > offsets