Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_Integrator_GradBasisDotVector_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_Integrator_GradBasisDotVector_impl_hpp__
44 #define __Panzer_Integrator_GradBasisDotVector_impl_hpp__
45 
47 //
48 // Include Files
49 //
51 
52 // Panzer
53 #include "Panzer_BasisIRLayout.hpp"
57 
58 namespace panzer
59 {
61  //
62  // Main Constructor
63  //
65  template<typename EvalT, typename Traits>
69  const std::string& resName,
70  const std::string& fluxName,
71  const panzer::BasisIRLayout& basis,
72  const panzer::IntegrationRule& ir,
73  const double& multiplier, /* = 1 */
74  const std::vector<std::string>& fmNames, /* =
75  std::vector<std::string>() */
76  const Teuchos::RCP<PHX::DataLayout>& vecDL /* = Teuchos::null */)
77  :
78  evalStyle_(evalStyle),
79  multiplier_(multiplier),
80  basisName_(basis.name())
81  {
82  using Kokkos::View;
83  using panzer::BASIS;
84  using panzer::Cell;
86  using panzer::IP;
87  using PHX::DataLayout;
88  using PHX::MDField;
89  using std::invalid_argument;
90  using std::logic_error;
91  using std::string;
92  using Teuchos::RCP;
93 
94  // Ensure the input makes sense.
95  TEUCHOS_TEST_FOR_EXCEPTION(resName == "", invalid_argument, "Error: " \
96  "Integrator_GradBasisDotVector called with an empty residual name.")
97  TEUCHOS_TEST_FOR_EXCEPTION(fluxName == "", invalid_argument, "Error: " \
98  "Integrator_GradBasisDotVector called with an empty flux name.")
99  RCP<const PureBasis> tmpBasis = basis.getBasis();
100  TEUCHOS_TEST_FOR_EXCEPTION(not tmpBasis->supportsGrad(), logic_error,
101  "Error: Integrator_GradBasisDotVector: Basis of type \""
102  << tmpBasis->name() << "\" does not support the gradient operator.")
103  RCP<DataLayout> tmpVecDL = ir.dl_vector;
104  if (not vecDL.is_null())
105  {
106  tmpVecDL = vecDL;
108  tmpVecDL->extent(2) < ir.dl_vector->extent(2), logic_error,
109  "Integrator_GradBasisDotVector: Dimension of space exceeds " \
110  "dimension of Vector Data Layout.");
111  } // end if (not vecDL.is_null())
112 
113  // Create the field for the vector-valued function we're integrating.
114  vector_ = MDField<const ScalarT, Cell, IP, Dim>(fluxName, tmpVecDL);
115  this->addDependentField(vector_);
116 
117  // Create the field that we're either contributing to or evaluating
118  // (storing).
119  field_ = MDField<ScalarT, Cell, BASIS>(resName, basis.functional);
121  this->addContributedField(field_);
122  else // if (evalStyle_ == EvaluatorStyle::EVALUATES)
123  this->addEvaluatedField(field_);
124 
125  // Add the dependent field multipliers, if there are any.
126  int i(0);
127  fieldMults_.resize(fmNames.size());
129  View<View<const ScalarT**,typename PHX::DevLayout<ScalarT>::type,PHX::Device>*>("GradBasisDotVector::KokkosFieldMultipliers",
130  fmNames.size());
131  for (const auto& name : fmNames)
132  {
133  fieldMults_[i++] = MDField<const ScalarT, Cell, IP>(name, ir.dl_scalar);
134  this->addDependentField(fieldMults_[i - 1]);
135  } // end loop over the field multipliers
136 
137  // Set the name of this object.
138  string n("Integrator_GradBasisDotVector (");
140  n += "CONTRIBUTES";
141  else // if (evalStyle_ == EvaluatorStyle::EVALUATES)
142  n += "EVALUATES";
143  n += "): " + field_.fieldTag().name();
144  this->setName(n);
145  } // end of Main Constructor
146 
148  //
149  // ParameterList Constructor
150  //
152  template<typename EvalT, typename Traits>
155  const Teuchos::ParameterList& p)
156  :
158  panzer::EvaluatorStyle::EVALUATES,
159  p.get<std::string>("Residual Name"),
160  p.get<std::string>("Flux Name"),
161  (*p.get<Teuchos::RCP<panzer::BasisIRLayout>>("Basis")),
162  (*p.get<Teuchos::RCP<panzer::IntegrationRule>>("IR")),
163  p.get<double>("Multiplier"),
164  p.isType<Teuchos::RCP<const std::vector<std::string>>>
165  ("Field Multipliers") ?
166  (*p.get<Teuchos::RCP<const std::vector<std::string>>>
167  ("Field Multipliers")) : std::vector<std::string>(),
168  p.isType<Teuchos::RCP<PHX::DataLayout>>("Vector Data Layout") ?
169  p.get<Teuchos::RCP<PHX::DataLayout>>("Vector Data Layout") :
170  Teuchos::null)
171  {
173  using Teuchos::RCP;
174 
175  // Ensure that the input ParameterList didn't contain any bogus entries.
176  RCP<ParameterList> validParams = this->getValidParameters();
177  p.validateParameters(*validParams);
178  } // end of ParameterList Constructor
179 
181  //
182  // postRegistrationSetup()
183  //
185  template<typename EvalT, typename Traits>
186  void
189  typename Traits::SetupData sd,
190  PHX::FieldManager<Traits>& /* fm */)
191  {
192  using panzer::getBasisIndex;
193  using std::size_t;
194  using PHX::Device;
195 
196  // Get the Kokkos::Views of the field multipliers.
197  for (size_t i(0); i < fieldMults_.size(); ++i)
198  kokkosFieldMults_(i) = fieldMults_[i].get_static_view();
199  Device().fence();
200 
201  // Determine the index in the Workset bases for our particular basis name.
202  basisIndex_ = getBasisIndex(basisName_, (*sd.worksets_)[0], this->wda);
203  } // end of postRegistrationSetup()
204 
206  //
207  // operator()() NO shared memory
208  //
210  template<typename EvalT, typename Traits>
211  template<int NUM_FIELD_MULT>
212  KOKKOS_INLINE_FUNCTION
213  void
216  const FieldMultTag<NUM_FIELD_MULT>& /* tag */,
217  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
218  {
220  const int cell = team.league_rank();
221 
222  // Initialize the evaluated field.
223  const int numQP(vector_.extent(1)), numDim(vector_.extent(2)),
224  numBases(basis_.extent(1));
225  if (evalStyle_ == EvaluatorStyle::EVALUATES)
226  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (const int& basis) {
227  field_(cell, basis) = 0.0;
228  });
229 
230  // The following if-block is for the sake of optimization depending on the
231  // number of field multipliers.
232  ScalarT tmp;
233  if (NUM_FIELD_MULT == 0)
234  {
235  // Loop over the quadrature points and dimensions of our vector fields,
236  // scale the integrand by the multiplier, and then perform the
237  // integration, looping over the bases.
238  for (int qp(0); qp < numQP; ++qp)
239  {
240  for (int dim(0); dim < numDim; ++dim)
241  {
242  tmp = multiplier_ * vector_(cell, qp, dim);
243  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (const int& basis) {
244  field_(cell, basis) += basis_(cell, basis, qp, dim) * tmp;
245  });
246  } // end loop over the dimensions of the vector field
247  } // end loop over the quadrature points
248  }
249  else if (NUM_FIELD_MULT == 1)
250  {
251  // Loop over the quadrature points and dimensions of our vector fields,
252  // scale the integrand by the multiplier and the single field multiplier,
253  // and then perform the actual integration, looping over the bases.
254  for (int qp(0); qp < numQP; ++qp)
255  {
256  for (int dim(0); dim < numDim; ++dim)
257  {
258  tmp = multiplier_ * vector_(cell, qp, dim) *
259  kokkosFieldMults_(0)(cell, qp);
260  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (const int& basis) {
261  field_(cell, basis) += basis_(cell, basis, qp, dim) * tmp;
262  });
263  } // end loop over the dimensions of the vector field
264  } // end loop over the quadrature points
265  }
266  else
267  {
268  // Loop over the quadrature points and pre-multiply all the field
269  // multipliers together. Then loop over the dimensions of our vector
270  // fields, scale the integrand by the multiplier and the combination of
271  // the field multipliers, and then perform the actual integration,
272  // looping over the bases.
273  const int numFieldMults(kokkosFieldMults_.extent(0));
274  for (int qp(0); qp < numQP; ++qp)
275  {
276  ScalarT fieldMultsTotal(1);
277  for (int fm(0); fm < numFieldMults; ++fm)
278  fieldMultsTotal *= kokkosFieldMults_(fm)(cell, qp);
279  for (int dim(0); dim < numDim; ++dim)
280  {
281  tmp = multiplier_ * vector_(cell, qp, dim) * fieldMultsTotal;
282  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (const int& basis) {
283  field_(cell, basis) += basis_(cell, basis, qp, dim) * tmp;
284  });
285  } // end loop over the dimensions of the vector field
286  } // end loop over the quadrature points
287  } // end if (NUM_FIELD_MULT == something)
288  } // end of operator()()
289 
291  //
292  // operator()() With shared memory
293  //
295  template<typename EvalT, typename Traits>
296  template<int NUM_FIELD_MULT>
297  KOKKOS_INLINE_FUNCTION
298  void
301  const SharedFieldMultTag<NUM_FIELD_MULT>& /* tag */,
302  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
303  {
305  const int cell = team.league_rank();
306  const int numQP = vector_.extent(1);
307  const int numDim = vector_.extent(2);
308  const int numBases = basis_.extent(1);
309  const int fadSize = Kokkos::dimension_scalar(field_.get_view());
310 
311  scratch_view tmp;
312  scratch_view tmp_field;
313  if (Sacado::IsADType<ScalarT>::value) {
314  tmp = scratch_view(team.team_shmem(),1,fadSize);
315  tmp_field = scratch_view(team.team_shmem(),numBases,fadSize);
316  }
317  else {
318  tmp = scratch_view(team.team_shmem(),1);
319  tmp_field = scratch_view(team.team_shmem(),numBases);
320  }
321 
322  // Initialize the evaluated field.
323  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (const int& basis) {
324  tmp_field(basis) = 0.0;
325  });
326 
327  // The following if-block is for the sake of optimization depending on the
328  // number of field multipliers.
329  if (NUM_FIELD_MULT == 0)
330  {
331  // Loop over the quadrature points and dimensions of our vector fields,
332  // scale the integrand by the multiplier, and then perform the
333  // integration, looping over the bases.
334  for (int qp(0); qp < numQP; ++qp)
335  {
336  for (int dim(0); dim < numDim; ++dim)
337  {
338  team.team_barrier();
339  tmp(0) = multiplier_ * vector_(cell, qp, dim);
340  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (const int& basis) {
341  tmp_field(basis) += basis_(cell, basis, qp, dim) * tmp(0);
342  });
343  } // end loop over the dimensions of the vector field
344  } // end loop over the quadrature points
345  }
346  else if (NUM_FIELD_MULT == 1)
347  {
348  // Loop over the quadrature points and dimensions of our vector fields,
349  // scale the integrand by the multiplier and the single field multiplier,
350  // and then perform the actual integration, looping over the bases.
351  for (int qp(0); qp < numQP; ++qp)
352  {
353  for (int dim(0); dim < numDim; ++dim)
354  {
355  team.team_barrier();
356  tmp(0) = multiplier_ * vector_(cell, qp, dim) *
357  kokkosFieldMults_(0)(cell, qp);
358  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (const int& basis) {
359  tmp_field(basis) += basis_(cell, basis, qp, dim) * tmp(0);
360  });
361  } // end loop over the dimensions of the vector field
362  } // end loop over the quadrature points
363  }
364  else
365  {
366  // Loop over the quadrature points and pre-multiply all the field
367  // multipliers together. Then loop over the dimensions of our vector
368  // fields, scale the integrand by the multiplier and the combination of
369  // the field multipliers, and then perform the actual integration,
370  // looping over the bases.
371  const int numFieldMults(kokkosFieldMults_.extent(0));
372  for (int qp(0); qp < numQP; ++qp)
373  {
374  ScalarT fieldMultsTotal(1); // need shared mem here
375  for (int fm(0); fm < numFieldMults; ++fm)
376  fieldMultsTotal *= kokkosFieldMults_(fm)(cell, qp);
377  for (int dim(0); dim < numDim; ++dim)
378  {
379  team.team_barrier();
380  tmp(0) = multiplier_ * vector_(cell, qp, dim) * fieldMultsTotal;
381  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (const int& basis) {
382  tmp_field(basis) += basis_(cell, basis, qp, dim) * tmp(0);
383  });
384  } // end loop over the dimensions of the vector field
385  } // end loop over the quadrature points
386  } // end if (NUM_FIELD_MULT == something)
387 
388  // Put final values into target field
389  if (evalStyle_ == EvaluatorStyle::EVALUATES) {
390  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),[&] (const int& basis) {
391  field_(cell,basis) = tmp_field(basis);
392  });
393  }
394  else { // Contributed
395  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),[&] (const int& basis) {
396  field_(cell,basis) += tmp_field(basis);
397  });
398  }
399 
400  } // end of operator()()
401 
403  //
404  // evaluateFields()
405  //
407  template<typename EvalT, typename Traits>
408  void
411  typename Traits::EvalData workset)
412  {
413  using Kokkos::parallel_for;
414  using Kokkos::TeamPolicy;
415 
416  // Grab the basis information.
417  basis_ = this->wda(workset).bases[basisIndex_]->weighted_grad_basis;
418 
419  bool use_shared_memory = panzer::HP::inst().useSharedMemory<ScalarT>();
420  if (use_shared_memory) {
421  int bytes;
422  if (Sacado::IsADType<ScalarT>::value) {
423  const int fadSize = Kokkos::dimension_scalar(field_.get_view());
424  bytes = scratch_view::shmem_size(1,fadSize) + scratch_view::shmem_size(basis_.extent(1),fadSize);
425  }
426  else
427  bytes = scratch_view::shmem_size(1) + scratch_view::shmem_size(basis_.extent(1));
428 
429  // The following if-block is for the sake of optimization depending on the
430  // number of field multipliers. The parallel_fors will loop over the cells
431  // in the Workset and execute operator()() above.
432  if (fieldMults_.size() == 0) {
433  auto policy = panzer::HP::inst().teamPolicy<ScalarT,SharedFieldMultTag<0>,PHX::Device>(workset.num_cells).set_scratch_size(0,Kokkos::PerTeam(bytes));
434  parallel_for(policy, *this, this->getName());
435  } else if (fieldMults_.size() == 1) {
436  auto policy = panzer::HP::inst().teamPolicy<ScalarT,SharedFieldMultTag<1>,PHX::Device>(workset.num_cells).set_scratch_size(0,Kokkos::PerTeam(bytes));
437  parallel_for(policy, *this, this->getName());
438  } else {
439  auto policy = panzer::HP::inst().teamPolicy<ScalarT,SharedFieldMultTag<-1>,PHX::Device>(workset.num_cells).set_scratch_size(0,Kokkos::PerTeam(bytes));
440  parallel_for(policy, *this, this->getName());
441  }
442  }
443  else {
444  // The following if-block is for the sake of optimization depending on the
445  // number of field multipliers. The parallel_fors will loop over the cells
446  // in the Workset and execute operator()() above.
447  if (fieldMults_.size() == 0) {
448  auto policy = panzer::HP::inst().teamPolicy<ScalarT,FieldMultTag<0>,PHX::Device>(workset.num_cells);
449  parallel_for(policy, *this, this->getName());
450  } else if (fieldMults_.size() == 1) {
451  auto policy = panzer::HP::inst().teamPolicy<ScalarT,FieldMultTag<1>,PHX::Device>(workset.num_cells);
452  parallel_for(policy, *this, this->getName());
453  } else {
454  auto policy = panzer::HP::inst().teamPolicy<ScalarT,FieldMultTag<-1>,PHX::Device>(workset.num_cells);
455  parallel_for(policy, *this, this->getName());
456  }
457  }
458  } // end of evaluateFields()
459 
461  //
462  // getValidParameters()
463  //
465  template<typename EvalT, typename TRAITS>
469  {
470  using panzer::BasisIRLayout;
472  using PHX::DataLayout;
473  using std::string;
474  using std::vector;
476  using Teuchos::RCP;
477  using Teuchos::rcp;
478 
479  // Create a ParameterList with all the valid keys we support.
480  RCP<ParameterList> p = rcp(new ParameterList);
481  p->set<string>("Residual Name", "?");
482  p->set<string>("Flux Name", "?");
483  RCP<BasisIRLayout> basis;
484  p->set("Basis", basis);
486  p->set("IR", ir);
487  p->set<double>("Multiplier", 1.0);
489  p->set("Field Multipliers", fms);
490  RCP<DataLayout> vecDL;
491  p->set("Vector Data Layout", vecDL);
492  return p;
493  } // end of getValidParameters()
494 
495 } // end of namespace panzer
496 
497 #endif // __Panzer_Integrator_GradBasisDotVector_impl_hpp__
PHX::MDField< const ScalarT, panzer::Cell, panzer::IP, panzer::Dim > vector_
A field representing the vector-valued function we&#39;re integrating ( ).
Integrator_GradBasisDotVector(const panzer::EvaluatorStyle &evalStyle, const std::string &resName, const std::string &fluxName, const panzer::BasisIRLayout &basis, const panzer::IntegrationRule &ir, const double &multiplier=1, const std::vector< std::string > &fmNames=std::vector< std::string >(), const Teuchos::RCP< PHX::DataLayout > &vecDL=Teuchos::null)
Main Constructor.
Kokkos::TeamPolicy< TeamPolicyProperties...> teamPolicy(const int &league_size)
Returns a TeamPolicy for hierarchic parallelism.
PHX::MDField< ScalarT, panzer::Cell, panzer::BASIS > field_
A field to which we&#39;ll contribute, or in which we&#39;ll store, the result of computing this integral...
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
Get Valid Parameters.
KOKKOS_INLINE_FUNCTION void operator()(const FieldMultTag< NUM_FIELD_MULT > &tag, const Kokkos::TeamPolicy< PHX::exec_space >::member_type &team) const
Perform the integration.
This empty struct allows us to optimize operator()() depending on the number of field multipliers...
EvaluatorStyle
An indication of how an Evaluator will behave.
const panzer::EvaluatorStyle evalStyle_
An enum determining the behavior of this Evaluator.
Kokkos::View< ScalarT *,typename PHX::DevLayout< ScalarT >::type, typename PHX::exec_space::scratch_memory_space, Kokkos::MemoryUnmanaged > scratch_view
Type for shared memory.
Teuchos::RCP< const PureBasis > getBasis() const
std::vector< PHX::MDField< const ScalarT, panzer::Cell, panzer::IP > > fieldMults_
The (possibly empty) list of fields that are multipliers out in front of the integral ( ...
double multiplier
The scalar multiplier out in front of the integral ( ).
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Teuchos::RCP< PHX::DataLayout > dl_scalar
Data layout for scalar fields.
panzer::EvaluatorStyle evalStyle
The EvaluatorStyle of the parent Integrator_CurlBasisDotVector object.
std::vector< std::string >::size_type getBasisIndex(std::string basis_name, const panzer::Workset &workset, WorksetDetailsAccessor &wda)
Returns the index in the workset bases for a particular BasisIRLayout name.
void evaluateFields(typename Traits::EvalData d)
Evaluate Fields.
void validateParameters(ParameterList const &validParamList, int const depth=1000, EValidateUsed const validateUsed=VALIDATE_USED_ENABLED, EValidateDefaults const validateDefaults=VALIDATE_DEFAULTS_ENABLED) const
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
Post-Registration Setup.
Teuchos::RCP< PHX::DataLayout > dl_vector
Data layout for vector fields.
static HP & inst()
Private ctor.
const std::vector< std::pair< int, LocalOrdinal > > &pid_and_lid const
This empty struct allows us to optimize operator()() depending on the number of field multipliers...
Teuchos::RCP< PHX::DataLayout > functional
&lt;Cell,Basis&gt;
Kokkos::View< Kokkos::View< const ScalarT **, typename PHX::DevLayout< ScalarT >::type, PHX::Device > * > kokkosFieldMults_
The Kokkos::View representation of the (possibly empty) list of fields that are multipliers out in fr...
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_
bool is_null() const