Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_Integrator_BasisTimesVector_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_Integerator_BasisTimesVector_impl_hpp__
12 #define __Panzer_Integerator_BasisTimesVector_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& valName,
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  :
45  evalStyle_(evalStyle),
46  useDescriptors_(false),
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 panzer::PureBasis;
56  using PHX::MDField;
57  using PHX::print;
58  using std::invalid_argument;
59  using std::logic_error;
60  using std::string;
61  using Teuchos::RCP;
62 
63  // Ensure the input makes sense.
64  TEUCHOS_TEST_FOR_EXCEPTION(resName == "", invalid_argument, "Error: " \
65  "Integrator_BasisTimesVector called with an empty residual name.")
66  TEUCHOS_TEST_FOR_EXCEPTION(valName == "", invalid_argument, "Error: " \
67  "Integrator_BasisTimesVector called with an empty value name.")
68  RCP<const PureBasis> tmpBasis = basis.getBasis();
69  TEUCHOS_TEST_FOR_EXCEPTION(not tmpBasis->isVectorBasis(), logic_error,
70  "Error: Integrator_BasisTimesVector: Basis of type \""
71  << tmpBasis->name() << "\" is not a vector basis.");
72  TEUCHOS_TEST_FOR_EXCEPTION(not tmpBasis->requiresOrientations(),
73  logic_error, "Integrator_BasisTimesVector: Basis of type \""
74  << tmpBasis->name() << "\" does not require orientations. This seems " \
75  "very strange, so I'm failing.");
76  TEUCHOS_TEST_FOR_EXCEPTION(fmNames.size() > 3,std::runtime_error,
77  "ERROR: Integrator_BasisTimesVector only supports up to three multipliers!");
78 
79  // Create the field for the vector-valued quantity we're integrating.
80  vector_ = MDField<const ScalarT, Cell, IP, Dim>(valName, ir.dl_vector);
81  this->addDependentField(vector_);
82 
83  // Create the field that we're either contributing to or evaluating
84  // (storing).
85  field_ = MDField<ScalarT, Cell, BASIS>(resName, basis.functional);
86  if (evalStyle == EvaluatorStyle::CONTRIBUTES)
87  this->addContributedField(field_);
88  else // if (evalStyle == EvaluatorStyle::EVALUATES)
89  this->addEvaluatedField(field_);
90 
91  // Add the dependent field multipliers, if there are any.
92  int i(0);
93  fieldMults_.resize(fmNames.size());
94  for (const auto& name : fmNames) {
95  fieldMults_[i++] = MDField<const ScalarT, Cell, IP>(name, ir.dl_scalar);
96  this->addDependentField(fieldMults_[i - 1]);
97  } // end loop over the field multipliers
98 
99  // Set the name of this object.
100  string n("Integrator_BasisTimesVector<");
101  n += std::to_string(fmNames.size()) + ">(";
102  if (evalStyle == EvaluatorStyle::CONTRIBUTES)
103  n += "Cont";
104  else // if (evalStyle == EvaluatorStyle::EVALUATES)
105  n += "Eval";
106  n += ", " + print<EvalT>() + "): " + field_.fieldTag().name();
107  this->setName(n);
108  } // end of Main Constructor
109 
111  //
112  // Descriptor Constructor
113  //
115  template<typename EvalT, typename Traits>
119  const PHX::FieldTag& resTag,
120  const PHX::FieldTag& valTag,
121  const BasisDescriptor& bd,
122  const IntegrationDescriptor& id,
123  const double& multiplier, /* = 1 */
124  const std::vector<PHX::FieldTag>& multipliers /* =
125  std::vector<PHX::FieldTag>() */)
126  :
127  evalStyle_(evalStyle),
128  useDescriptors_(true),
129  bd_(bd),
130  id_(id),
131  multiplier_(multiplier)
132  {
133  using PHX::View;
134  using panzer::BASIS;
135  using panzer::Cell;
137  using panzer::IP;
138  using panzer::PureBasis;
139  using PHX::MDField;
140  using PHX::print;
141  using std::invalid_argument;
142  using std::logic_error;
143  using std::string;
144  using Teuchos::RCP;
145 
146  // Ensure the input makes sense.
147  TEUCHOS_TEST_FOR_EXCEPTION(not ((bd_.getType() == "HCurl") or
148  (bd_.getType() == "HDiv")), logic_error, "Error: " \
149  "Integrator_BasisTimesVector: Basis of type \"" << bd_.getType()
150  << "\" is not a vector basis.")
151 
152  TEUCHOS_TEST_FOR_EXCEPTION(multipliers.size() > 3,std::runtime_error,
153  "ERROR: Integrator_BasisTimesVector only supports up to three multipliers!");
154 
155  // Create the field for the vector-valued quantity we're integrating.
156  vector_ = valTag;
157  this->addDependentField(vector_);
158 
159  // Create the field that we're either contributing to or evaluating
160  // (storing).
161  field_ = resTag;
162  if (evalStyle == EvaluatorStyle::CONTRIBUTES)
163  this->addContributedField(field_);
164  else // if (evalStyle == EvaluatorStyle::EVALUATES)
165  this->addEvaluatedField(field_);
166 
167  // Add the dependent field multipliers, if there are any.
168  int i(0);
169  fieldMults_.resize(multipliers.size());
170  for (const auto& fm : multipliers)
171  {
172  fieldMults_[i++] = fm;
173  this->addDependentField(fieldMults_[i - 1]);
174  } // end loop over the field multipliers
175 
176  // Set the name of this object.
177  string n("Integrator_BasisTimesVector<");
178  n += std::to_string(multipliers.size()) + ">(";
179  if (evalStyle == EvaluatorStyle::CONTRIBUTES)
180  n += "Cont";
181  else // if (evalStyle == EvaluatorStyle::EVALUATES)
182  n += "Eval";
183  n += ", " + print<EvalT>() + "): " + field_.fieldTag().name();
184  this->setName(n);
185  } // end of Main Constructor
186 
188  //
189  // ParameterList Constructor
190  //
192  template<typename EvalT, typename Traits>
195  const Teuchos::ParameterList& p)
196  :
198  panzer::EvaluatorStyle::EVALUATES,
199  p.get<std::string>("Residual Name"),
200  p.get<std::string>("Value Name"),
201  (*p.get<Teuchos::RCP<panzer::BasisIRLayout>>("Basis")),
202  (*p.get<Teuchos::RCP<panzer::IntegrationRule>>("IR")),
203  p.get<double>("Multiplier"),
204  p.isType<Teuchos::RCP<const std::vector<std::string>>>
205  ("Field Multipliers") ?
206  (*p.get<Teuchos::RCP<const std::vector<std::string>>>
207  ("Field Multipliers")) : std::vector<std::string>())
208  {
210  using Teuchos::RCP;
211 
212  // Ensure that the input ParameterList didn't contain any bogus entries.
213  RCP<ParameterList> validParams = this->getValidParameters();
214  p.validateParameters(*validParams);
215  } // end of ParameterList Constructor
216 
218  //
219  // postRegistrationSetup()
220  //
222  template<typename EvalT, typename Traits>
223  void
226  typename Traits::SetupData sd,
227  PHX::FieldManager<Traits>& /* fm */)
228  {
229  using panzer::getBasisIndex;
230  using std::size_t;
231 
232  // Get the PHX::Views of the field multipliers.
233  for (size_t i=0; i < fieldMults_.size(); ++i)
234  kokkosFieldMults_[i] = fieldMults_[i].get_static_view();
235 
236  // Determine the number of quadrature points and the dimensionality of the
237  // vector that we're integrating.
238  numQP_ = vector_.extent(1);
239  numDim_ = vector_.extent(2);
240 
241  // Determine the index in the Workset bases for our particular basis name.
242  if (not useDescriptors_)
243  basisIndex_ = getBasisIndex(basisName_, (*sd.worksets_)[0], this->wda);
244 
245  } // end of postRegistrationSetup()
246 
248  //
249  // operator()()
250  //
252  template<typename EvalT, typename Traits>
253  template<int NUM_FIELD_MULT>
254  KOKKOS_INLINE_FUNCTION
255  void
258  const FieldMultTag<NUM_FIELD_MULT>& /* tag */,
259  const typename Kokkos::TeamPolicy<FieldMultTag<NUM_FIELD_MULT>,PHX::exec_space>::member_type& team) const
260  {
262  const int cell = team.league_rank();
263  const int numBases = basis_.extent(1);
264 
265  // Initialize the evaluated field.
266  if (evalStyle_ == EvaluatorStyle::EVALUATES) {
267  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (const int basis) {
268  field_(cell, basis) = 0.0;
269  });
270  }
271 
272  for (int qp(0); qp < numQP_; ++qp) {
273  for (int dim(0); dim < numDim_; ++dim) {
274  if constexpr (NUM_FIELD_MULT == 0) {
275  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (const int basis) {
276  field_(cell, basis) += basis_(cell, basis, qp, dim) * multiplier_ * vector_(cell, qp, dim);
277  });
278  }
279  else if constexpr (NUM_FIELD_MULT == 1) {
280  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (const int basis) {
281  field_(cell, basis) += basis_(cell, basis, qp, dim) * multiplier_ * vector_(cell, qp, dim)
282  * kokkosFieldMults_[0](cell, qp);
283  });
284  }
285  else if constexpr (NUM_FIELD_MULT == 2) {
286  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (const int basis) {
287  field_(cell, basis) += basis_(cell, basis, qp, dim) * multiplier_ * vector_(cell, qp, dim)
288  * kokkosFieldMults_[0](cell, qp) * kokkosFieldMults_[1](cell, qp);
289  });
290  }
291  else if constexpr (NUM_FIELD_MULT == 3) {
292  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (const int basis) {
293  field_(cell, basis) += basis_(cell, basis, qp, dim) * multiplier_ * vector_(cell, qp, dim)
294  * kokkosFieldMults_[0](cell, qp) * kokkosFieldMults_[1](cell, qp) * kokkosFieldMults_[2](cell, qp);
295  });
296  }
297  else {
298  Kokkos::abort("Panzer_Integrator_BasisTimesVector: NUM_FIELD_MULT out of range!");
299  }
300  } // end loop over the dimensions of the vector field
301  } // end loop over the quadrature points
302  } // end of operator()()
303 
305  //
306  // evaluateFields()
307  //
309  template<typename EvalT, typename Traits>
310  void
313  typename Traits::EvalData workset)
314  {
315  using Kokkos::parallel_for;
316  using Kokkos::RangePolicy;
317 
318  // Grab the basis information.
319  const panzer::BasisValues2<double>& bv = useDescriptors_ ?
320  this->wda(workset).getBasisValues(bd_,id_) :
321  *this->wda(workset).bases[basisIndex_];
323  basis_ = useDescriptors_ ? bv.getVectorBasisValues(true) : Array(bv.weighted_basis_vector);
324 
325  // The following if-block is for the sake of optimization depending on the
326  // number of field multipliers. The parallel_fors will loop over the cells
327  // in the Workset and execute operator()() above.
328  if (fieldMults_.size() == 0) {
329  auto policy = panzer::HP::inst().teamPolicy<ScalarT,FieldMultTag<0>,PHX::exec_space>(workset.num_cells);
330  parallel_for("Panzer_Integrator_BasisTimesVector<0>", policy, *this);
331  }
332  else if (fieldMults_.size() == 1) {
333  auto policy = panzer::HP::inst().teamPolicy<ScalarT,FieldMultTag<1>,PHX::exec_space>(workset.num_cells);
334  parallel_for("Panzer_Integrator_BasisTimesVector<1>", policy, *this);
335  }
336  else if (fieldMults_.size() == 2) {
337  auto policy = panzer::HP::inst().teamPolicy<ScalarT,FieldMultTag<2>,PHX::exec_space>(workset.num_cells);
338  parallel_for("Panzer_Integrator_BasisTimesVector<2>", policy, *this);
339  }
340  else if (fieldMults_.size() == 3) {
341  auto policy = panzer::HP::inst().teamPolicy<ScalarT,FieldMultTag<3>,PHX::exec_space>(workset.num_cells);
342  parallel_for("Panzer_Integrator_BasisTimesVector<3>", policy, *this);
343  }
344  } // end of evaluateFields()
345 
347  //
348  // getValidParameters()
349  //
351  template<typename EvalT, typename Traits>
355  {
356  using panzer::BasisIRLayout;
358  using std::string;
359  using std::vector;
361  using Teuchos::RCP;
362  using Teuchos::rcp;
363 
364  // Create a ParameterList with all the valid keys we support.
365  RCP<ParameterList> p = rcp(new ParameterList);
366  p->set<string>("Residual Name", "?");
367  p->set<string>("Value Name", "?");
368  RCP<BasisIRLayout> basis;
369  p->set("Basis", basis);
371  p->set("IR", ir);
372  p->set<double>("Multiplier", 1.0);
374  p->set("Field Multipliers", fms);
375  return p;
376  } // end of getValidParameters()
377 
378 } // end of namespace panzer
379 
380 #endif // __Panzer_Integerator_BasisTimesVector_impl_hpp__
int num_cells
DEPRECATED - use: numCells()
Kokkos::TeamPolicy< TeamPolicyProperties...> teamPolicy(const int &league_size)
Returns a TeamPolicy for hierarchic parallelism.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
BasisDescriptor bd_
The BasisDescriptor for the basis to use.
EvaluatorStyle
An indication of how an Evaluator will behave.
Teuchos::RCP< const PureBasis > getBasis() const
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
Get Valid Parameters.
double multiplier
The scalar multiplier out in front of the integral ( ).
KOKKOS_INLINE_FUNCTION void operator()(const FieldMultTag< NUM_FIELD_MULT > &tag, const typename Kokkos::TeamPolicy< FieldMultTag< NUM_FIELD_MULT >, PHX::exec_space >::member_type &team) const
Perform the integration.
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...
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
const std::string & getType() const
Get type of basis.
Teuchos::RCP< PHX::DataLayout > dl_scalar
Data layout for scalar fields.
void evaluateFields(typename Traits::EvalData workset)
Evaluate Fields.
ConstArray_CellBasisIP getBasisValues(const bool weighted, const bool cache=true, const bool force=false) const
Get the basis values evaluated at mesh points.
panzer::EvaluatorStyle evalStyle
The EvaluatorStyle of the parent Integrator_CurlBasisDotVector object.
typename EvalT::ScalarT ScalarT
The scalar type.
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.
This empty struct allows us to optimize operator()() depending on the number of field multipliers...
void postRegistrationSetup(typename Traits::SetupData sd, PHX::FieldManager< Traits > &fm)
Post-Registration Setup.
void validateParameters(ParameterList const &validParamList, int const depth=1000, EValidateUsed const validateUsed=VALIDATE_USED_ENABLED, EValidateDefaults const validateDefaults=VALIDATE_DEFAULTS_ENABLED) const
Array_CellBasisIPDim weighted_basis_vector
ConstArray_CellBasisIPDim getVectorBasisValues(const bool weighted, const bool cache=true, const bool force=false) const
Get the vector basis values evaluated at mesh points.
Teuchos::RCP< PHX::DataLayout > dl_vector
Data layout for vector fields.
static HP & inst()
Private ctor.
Integrator_BasisTimesVector(const panzer::EvaluatorStyle &evalStyle, const std::string &resName, const std::string &valName, const panzer::BasisIRLayout &basis, const panzer::IntegrationRule &ir, const double &multiplier=1, const std::vector< std::string > &fmNames=std::vector< std::string >())
Main Constructor.
const std::vector< std::pair< int, LocalOrdinal > > &pid_and_lid const
PHX::MDField< const ScalarT, panzer::Cell, panzer::IP, panzer::Dim > vector_
A field representing the vector-valued function we&#39;re integrating ( ).
Description and data layouts associated with a particular basis.
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 ( ...
Teuchos::RCP< PHX::DataLayout > functional
&lt;Cell,Basis&gt;
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_