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 // Panzer: A partial differential equation assembly
4 // engine for strongly coupled complex multiphysics systems
5 //
6 // Copyright 2011 NTESS and the Panzer contributors.
7 // SPDX-License-Identifier: BSD-3-Clause
8 // *****************************************************************************
9 // @HEADER
10 
11 #ifndef __Panzer_Integrator_GradBasisDotVector_impl_hpp__
12 #define __Panzer_Integrator_GradBasisDotVector_impl_hpp__
13 
15 //
16 // Include Files
17 //
19 
20 // Panzer
21 #include "Panzer_BasisIRLayout.hpp"
25 
26 namespace panzer
27 {
29  //
30  // Main Constructor
31  //
33  template<typename EvalT, typename Traits>
37  const std::string& resName,
38  const std::string& fluxName,
39  const panzer::BasisIRLayout& basis,
40  const panzer::IntegrationRule& ir,
41  const double& multiplier, /* = 1 */
42  const std::vector<std::string>& fmNames, /* =
43  std::vector<std::string>() */
44  const Teuchos::RCP<PHX::DataLayout>& vecDL /* = Teuchos::null */)
45  :
46  evalStyle_(evalStyle),
47  multiplier_(multiplier),
48  basisName_(basis.name())
49  {
50  using PHX::View;
51  using panzer::BASIS;
52  using panzer::Cell;
54  using panzer::IP;
55  using PHX::DataLayout;
56  using PHX::MDField;
57  using std::invalid_argument;
58  using std::logic_error;
59  using std::string;
60  using Teuchos::RCP;
61 
62  // Ensure the input makes sense.
63  TEUCHOS_TEST_FOR_EXCEPTION(resName == "", invalid_argument, "Error: " \
64  "Integrator_GradBasisDotVector called with an empty residual name.")
65  TEUCHOS_TEST_FOR_EXCEPTION(fluxName == "", invalid_argument, "Error: " \
66  "Integrator_GradBasisDotVector called with an empty flux name.")
67  RCP<const PureBasis> tmpBasis = basis.getBasis();
68  TEUCHOS_TEST_FOR_EXCEPTION(not tmpBasis->supportsGrad(), logic_error,
69  "Error: Integrator_GradBasisDotVector: Basis of type \""
70  << tmpBasis->name() << "\" does not support the gradient operator.")
71  RCP<DataLayout> tmpVecDL = ir.dl_vector;
72  if (not vecDL.is_null())
73  {
74  tmpVecDL = vecDL;
76  tmpVecDL->extent(2) < ir.dl_vector->extent(2), logic_error,
77  "Integrator_GradBasisDotVector: Dimension of space exceeds " \
78  "dimension of Vector Data Layout.");
79  } // end if (not vecDL.is_null())
80 
81  // Create the field for the vector-valued function we're integrating.
82  vector_ = MDField<const ScalarT, Cell, IP, Dim>(fluxName, tmpVecDL);
83  this->addDependentField(vector_);
84 
85  // Create the field that we're either contributing to or evaluating
86  // (storing).
87  field_ = MDField<ScalarT, Cell, BASIS>(resName, basis.functional);
89  this->addContributedField(field_);
90  else // if (evalStyle_ == EvaluatorStyle::EVALUATES)
91  this->addEvaluatedField(field_);
92 
93  // Add the dependent field multipliers, if there are any.
94  int i(0);
95  fieldMults_.resize(fmNames.size());
96  kokkosFieldMults_ = PHX::View<PHX::UnmanagedView<const ScalarT**>*>("GradBasisDotVector::KokkosFieldMultipliers",
97  fmNames.size());
98  for (const auto& name : fmNames)
99  {
100  fieldMults_[i++] = MDField<const ScalarT, Cell, IP>(name, ir.dl_scalar);
101  this->addDependentField(fieldMults_[i - 1]);
102  } // end loop over the field multipliers
103 
104  // Set the name of this object.
105  string n("Integrator_GradBasisDotVector (");
107  n += "CONTRIBUTES";
108  else // if (evalStyle_ == EvaluatorStyle::EVALUATES)
109  n += "EVALUATES";
110  n += "): " + field_.fieldTag().name();
111  this->setName(n);
112  } // end of Main Constructor
113 
115  //
116  // ParameterList Constructor
117  //
119  template<typename EvalT, typename Traits>
122  const panzer::EvaluatorStyle es)
123  :
125  p.get<std::string>("Residual Name"),
126  p.get<std::string>("Flux Name"),
127  (*p.get<Teuchos::RCP<panzer::BasisIRLayout>>("Basis")),
128  (*p.get<Teuchos::RCP<panzer::IntegrationRule>>("IR")),
129  p.get<double>("Multiplier"),
130  p.isType<Teuchos::RCP<const std::vector<std::string>>>
131  ("Field Multipliers") ?
132  (*p.get<Teuchos::RCP<const std::vector<std::string>>>
133  ("Field Multipliers")) : std::vector<std::string>(),
134  p.isType<Teuchos::RCP<PHX::DataLayout>>("Vector Data Layout") ?
135  p.get<Teuchos::RCP<PHX::DataLayout>>("Vector Data Layout") :
136  Teuchos::null)
137  {
139  using Teuchos::RCP;
140 
141  // Ensure that the input ParameterList didn't contain any bogus entries.
142  RCP<ParameterList> validParams = this->getValidParameters();
143  p.validateParameters(*validParams);
144  } // end of ParameterList Constructor
145 
147  //
148  // postRegistrationSetup()
149  //
151  template<typename EvalT, typename Traits>
152  void
155  typename Traits::SetupData sd,
156  PHX::FieldManager<Traits>& /* fm */)
157  {
158  using panzer::getBasisIndex;
159  using std::size_t;
160 
161  // Get the PHX::Views of the field multipliers.
162  auto field_mults_host_mirror = Kokkos::create_mirror_view(kokkosFieldMults_);
163  for (size_t i(0); i < fieldMults_.size(); ++i)
164  field_mults_host_mirror(i) = fieldMults_[i].get_static_view();
165  Kokkos::deep_copy(kokkosFieldMults_,field_mults_host_mirror);
166 
167  // Determine the index in the Workset bases for our particular basis name.
168  basisIndex_ = getBasisIndex(basisName_, (*sd.worksets_)[0], this->wda);
169 
170  // Allocate temporary if not using shared memory
171  bool use_shared_memory = panzer::HP::inst().useSharedMemory<ScalarT>();
172  if (!use_shared_memory) {
173  if (Sacado::IsADType<ScalarT>::value) {
174  const auto fadSize = Kokkos::dimension_scalar(field_.get_view());
175  tmp_ = PHX::View<ScalarT*>("GradBasisDotVector::tmp_",field_.extent(0),fadSize);
176  } else {
177  tmp_ = PHX::View<ScalarT*>("GradBasisDotVector::tmp_",field_.extent(0));
178  }
179  }
180  } // end of postRegistrationSetup()
181 
183  //
184  // operator()() NO shared memory
185  //
187  template<typename EvalT, typename Traits>
188  template<int NUM_FIELD_MULT>
189  KOKKOS_INLINE_FUNCTION
190  void
193  const FieldMultTag<NUM_FIELD_MULT>& /* tag */,
194  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
195  {
197  const int cell = team.league_rank();
198 
199  // Initialize the evaluated field.
200  const int numQP(vector_.extent(1)), numDim(vector_.extent(2)),
201  numBases(basis_.extent(1));
202  if (evalStyle_ == EvaluatorStyle::EVALUATES)
203  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (const int basis) {
204  field_(cell, basis) = 0.0;
205  });
206 
207  for (int qp(0); qp < numQP; ++qp) {
208  for (int dim(0); dim < numDim; ++dim) {
209  if constexpr (NUM_FIELD_MULT == 0) {
210  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (const int basis) {
211  field_(cell, basis) += basis_(cell, basis, qp, dim) * multiplier_ * vector_(cell, qp, dim);
212  });
213  }
214  else if constexpr (NUM_FIELD_MULT == 1) {
215  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (const int basis) {
216  field_(cell, basis) += basis_(cell, basis, qp, dim) * multiplier_ * vector_(cell, qp, dim)
217  * kokkosFieldMults_(0)(cell, qp);
218  });
219  }
220  else if constexpr (NUM_FIELD_MULT == 2) {
221  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (const int basis) {
222  field_(cell, basis) += basis_(cell, basis, qp, dim) * multiplier_ * vector_(cell, qp, dim)
223  * kokkosFieldMults_(0)(cell, qp) * kokkosFieldMults_(1)(cell, qp);
224  });
225  }
226  else if constexpr (NUM_FIELD_MULT == 3) {
227  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (const int basis) {
228  field_(cell, basis) += basis_(cell, basis, qp, dim) * multiplier_ * vector_(cell, qp, dim)
229  * kokkosFieldMults_(0)(cell, qp) * kokkosFieldMults_(1)(cell, qp) * kokkosFieldMults_(2)(cell, qp);
230  });
231  }
232  else {
233  Kokkos::abort("Panzer_Integrator_GradBasisDotVector: NUM_FIELD_MULT out of range!");
234  }
235  } // end dim loop
236  } // end qp loop
237  } // end of operator()()
238 
240  //
241  // operator()() With shared memory
242  //
244  template<typename EvalT, typename Traits>
245  template<int NUM_FIELD_MULT>
246  KOKKOS_INLINE_FUNCTION
247  void
250  const SharedFieldMultTag<NUM_FIELD_MULT>& /* tag */,
251  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
252  {
254  const int cell = team.league_rank();
255  const int numQP = vector_.extent(1);
256  const int numDim = vector_.extent(2);
257  const int numBases = basis_.extent(1);
258  const int fadSize = Kokkos::dimension_scalar(field_.get_view());
259 
260  scratch_view tmp_field;
261  if (Sacado::IsADType<ScalarT>::value) {
262  tmp_field = scratch_view(team.team_shmem(),numBases,fadSize);
263  }
264  else {
265  tmp_field = scratch_view(team.team_shmem(),numBases);
266  }
267 
268  // Initialize the evaluated field.
269  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (const int basis) {
270  tmp_field(basis) = 0.0;
271  });
272 
273  // The following if-block is for the sake of optimization depending on the
274  // number of field multipliers.
275  for (int qp(0); qp < numQP; ++qp) {
276  for (int dim(0); dim < numDim; ++dim) {
277  if constexpr (NUM_FIELD_MULT == 0) {
278  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (const int basis) {
279  tmp_field(basis) += basis_(cell, basis, qp, dim) * multiplier_ * vector_(cell, qp, dim);
280  });
281  }
282  else if constexpr (NUM_FIELD_MULT == 1) {
283  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (const int basis) {
284  tmp_field(basis) += basis_(cell, basis, qp, dim) * multiplier_ * vector_(cell, qp, dim)
285  * kokkosFieldMults_(0)(cell, qp);
286  });
287  }
288  else if constexpr (NUM_FIELD_MULT == 2) {
289  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (const int basis) {
290  tmp_field(basis) += basis_(cell, basis, qp, dim) * multiplier_ * vector_(cell, qp, dim)
291  * kokkosFieldMults_(0)(cell, qp) * kokkosFieldMults_(1)(cell, qp);
292  });
293  }
294  else if constexpr (NUM_FIELD_MULT == 3) {
295  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (const int basis) {
296  tmp_field(basis) += basis_(cell, basis, qp, dim) * multiplier_ * vector_(cell, qp, dim)
297  * kokkosFieldMults_(0)(cell, qp) * kokkosFieldMults_(1)(cell, qp) * kokkosFieldMults_(2)(cell, qp);
298  });
299  }
300  else {
301  Kokkos::abort("Panzer_Integrator_GradBasisDotVector: NUM_FIELD_MULT out of range!");
302  }
303  } // end loop over the dimensions of the vector field
304  } // end loop over the quadrature points
305 
306  // Put final values into target field
307  if (evalStyle_ == EvaluatorStyle::EVALUATES) {
308  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),[&] (const int basis) {
309  field_(cell,basis) = tmp_field(basis);
310  });
311  }
312  else { // Contributed
313  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),[&] (const int basis) {
314  field_(cell,basis) += tmp_field(basis);
315  });
316  }
317 
318  } // end of operator()()
319 
321  //
322  // evaluateFields()
323  //
325  template<typename EvalT, typename Traits>
326  void
329  typename Traits::EvalData workset)
330  {
331  using Kokkos::parallel_for;
332  using Kokkos::TeamPolicy;
333 
334  // Grab the basis information.
335  basis_ = this->wda(workset).bases[basisIndex_]->weighted_grad_basis;
336 
337  bool use_shared_memory = panzer::HP::inst().useSharedMemory<ScalarT>();
338  if (use_shared_memory) {
339  int bytes;
340  if (Sacado::IsADType<ScalarT>::value) {
341  const int fadSize = Kokkos::dimension_scalar(field_.get_view());
342  bytes = scratch_view::shmem_size(basis_.extent(1),fadSize);
343  }
344  else
345  bytes = scratch_view::shmem_size(basis_.extent(1));
346 
347  // The following if-block is for the sake of optimization depending on the
348  // number of field multipliers. The parallel_fors will loop over the cells
349  // in the Workset and execute operator()() above.
350  if (fieldMults_.size() == 0) {
351  auto policy = panzer::HP::inst().teamPolicy<ScalarT,SharedFieldMultTag<0>,PHX::Device>(workset.num_cells).set_scratch_size(0,Kokkos::PerTeam(bytes));
352  parallel_for(this->getName(), policy, *this);
353  } else if (fieldMults_.size() == 1) {
354  auto policy = panzer::HP::inst().teamPolicy<ScalarT,SharedFieldMultTag<1>,PHX::Device>(workset.num_cells).set_scratch_size(0,Kokkos::PerTeam(bytes));
355  parallel_for(this->getName(), policy, *this);
356  } else if (fieldMults_.size() == 2) {
357  auto policy = panzer::HP::inst().teamPolicy<ScalarT,SharedFieldMultTag<2>,PHX::Device>(workset.num_cells).set_scratch_size(0,Kokkos::PerTeam(bytes));
358  parallel_for(this->getName(), policy, *this);
359  } else if (fieldMults_.size() == 3) {
360  auto policy = panzer::HP::inst().teamPolicy<ScalarT,SharedFieldMultTag<3>,PHX::Device>(workset.num_cells).set_scratch_size(0,Kokkos::PerTeam(bytes));
361  parallel_for(this->getName(), policy, *this);
362  } else {
363  TEUCHOS_TEST_FOR_EXCEPTION(fieldMults_.size() > 3,std::runtime_error,
364  "ERROR: Panzer_Integrator_GradBasisDotVector supports up to three field multipliers! User requested "
365  << fieldMults_.size() << "!");
366  }
367  }
368  else {
369  // The following if-block is for the sake of optimization depending on the
370  // number of field multipliers. The parallel_fors will loop over the cells
371  // in the Workset and execute operator()() above.
372  if (fieldMults_.size() == 0) {
373  auto policy = panzer::HP::inst().teamPolicy<ScalarT,FieldMultTag<0>,PHX::Device>(workset.num_cells);
374  parallel_for(this->getName(), policy, *this);
375  } else if (fieldMults_.size() == 1) {
376  auto policy = panzer::HP::inst().teamPolicy<ScalarT,FieldMultTag<1>,PHX::Device>(workset.num_cells);
377  parallel_for(this->getName(), policy, *this);
378  } else if (fieldMults_.size() == 2) {
379  auto policy = panzer::HP::inst().teamPolicy<ScalarT,FieldMultTag<2>,PHX::Device>(workset.num_cells);
380  parallel_for(this->getName(), policy, *this);
381  } else if (fieldMults_.size() == 3) {
382  auto policy = panzer::HP::inst().teamPolicy<ScalarT,FieldMultTag<3>,PHX::Device>(workset.num_cells);
383  parallel_for(this->getName(), policy, *this);
384  } else {
385  TEUCHOS_TEST_FOR_EXCEPTION(fieldMults_.size() > 3,std::runtime_error,
386  "ERROR: Panzer_Integrator_GradBasisDotVector supports up to three field multipliers! User requested "
387  << fieldMults_.size() << "!");
388  }
389  }
390  } // end of evaluateFields()
391 
393  //
394  // getValidParameters()
395  //
397  template<typename EvalT, typename TRAITS>
401  {
402  using panzer::BasisIRLayout;
404  using PHX::DataLayout;
405  using std::string;
406  using std::vector;
408  using Teuchos::RCP;
409  using Teuchos::rcp;
410 
411  // Create a ParameterList with all the valid keys we support.
412  RCP<ParameterList> p = rcp(new ParameterList);
413  p->set<string>("Residual Name", "?");
414  p->set<string>("Flux Name", "?");
415  RCP<BasisIRLayout> basis;
416  p->set("Basis", basis);
418  p->set("IR", ir);
419  p->set<double>("Multiplier", 1.0);
421  p->set("Field Multipliers", fms);
422  RCP<DataLayout> vecDL;
423  p->set("Vector Data Layout", vecDL);
424  return p;
425  } // end of getValidParameters()
426 
427 } // end of namespace panzer
428 
429 #endif // __Panzer_Integrator_GradBasisDotVector_impl_hpp__
int num_cells
DEPRECATED - use: numCells()
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.
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.
PHX::View< PHX::UnmanagedView< const ScalarT ** > * > kokkosFieldMults_
The PHX::View representation of the (possibly empty) list of fields that are multipliers out in front...
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;
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_
bool is_null() const