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 //
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_Integerator_BasisTimesVector_impl_hpp__
44 #define __Panzer_Integerator_BasisTimesVector_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& valName,
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  :
77  evalStyle_(evalStyle),
78  useDescriptors_(false),
79  multiplier_(multiplier),
80  basisName_(basis.name())
81  {
82  using PHX::View;
83  using panzer::BASIS;
84  using panzer::Cell;
86  using panzer::IP;
87  using panzer::PureBasis;
88  using PHX::MDField;
89  using PHX::print;
90  using std::invalid_argument;
91  using std::logic_error;
92  using std::string;
93  using Teuchos::RCP;
94 
95  // Ensure the input makes sense.
96  TEUCHOS_TEST_FOR_EXCEPTION(resName == "", invalid_argument, "Error: " \
97  "Integrator_BasisTimesVector called with an empty residual name.")
98  TEUCHOS_TEST_FOR_EXCEPTION(valName == "", invalid_argument, "Error: " \
99  "Integrator_BasisTimesVector called with an empty value name.")
100  RCP<const PureBasis> tmpBasis = basis.getBasis();
101  TEUCHOS_TEST_FOR_EXCEPTION(not tmpBasis->isVectorBasis(), logic_error,
102  "Error: Integrator_BasisTimesVector: Basis of type \""
103  << tmpBasis->name() << "\" is not a vector basis.");
104  TEUCHOS_TEST_FOR_EXCEPTION(not tmpBasis->requiresOrientations(),
105  logic_error, "Integrator_BasisTimesVector: Basis of type \""
106  << tmpBasis->name() << "\" does not require orientations. This seems " \
107  "very strange, so I'm failing.");
108  TEUCHOS_TEST_FOR_EXCEPTION(fmNames.size() > 3,std::runtime_error,
109  "ERROR: Integrator_BasisTimesVector only supports up to three multipliers!");
110 
111  // Create the field for the vector-valued quantity we're integrating.
112  vector_ = MDField<const ScalarT, Cell, IP, Dim>(valName, ir.dl_vector);
113  this->addDependentField(vector_);
114 
115  // Create the field that we're either contributing to or evaluating
116  // (storing).
117  field_ = MDField<ScalarT, Cell, BASIS>(resName, basis.functional);
118  if (evalStyle == EvaluatorStyle::CONTRIBUTES)
119  this->addContributedField(field_);
120  else // if (evalStyle == EvaluatorStyle::EVALUATES)
121  this->addEvaluatedField(field_);
122 
123  // Add the dependent field multipliers, if there are any.
124  int i(0);
125  fieldMults_.resize(fmNames.size());
126  for (const auto& name : fmNames) {
127  fieldMults_[i++] = MDField<const ScalarT, Cell, IP>(name, ir.dl_scalar);
128  this->addDependentField(fieldMults_[i - 1]);
129  } // end loop over the field multipliers
130 
131  // Set the name of this object.
132  string n("Integrator_BasisTimesVector<");
133  n += std::to_string(fmNames.size()) + ">(";
134  if (evalStyle == EvaluatorStyle::CONTRIBUTES)
135  n += "Cont";
136  else // if (evalStyle == EvaluatorStyle::EVALUATES)
137  n += "Eval";
138  n += ", " + print<EvalT>() + "): " + field_.fieldTag().name();
139  this->setName(n);
140  } // end of Main Constructor
141 
143  //
144  // Descriptor Constructor
145  //
147  template<typename EvalT, typename Traits>
151  const PHX::FieldTag& resTag,
152  const PHX::FieldTag& valTag,
153  const BasisDescriptor& bd,
154  const IntegrationDescriptor& id,
155  const double& multiplier, /* = 1 */
156  const std::vector<PHX::FieldTag>& multipliers /* =
157  std::vector<PHX::FieldTag>() */)
158  :
159  evalStyle_(evalStyle),
160  useDescriptors_(true),
161  bd_(bd),
162  id_(id),
163  multiplier_(multiplier)
164  {
165  using PHX::View;
166  using panzer::BASIS;
167  using panzer::Cell;
169  using panzer::IP;
170  using panzer::PureBasis;
171  using PHX::MDField;
172  using PHX::print;
173  using std::invalid_argument;
174  using std::logic_error;
175  using std::string;
176  using Teuchos::RCP;
177 
178  // Ensure the input makes sense.
179  TEUCHOS_TEST_FOR_EXCEPTION(not ((bd_.getType() == "HCurl") or
180  (bd_.getType() == "HDiv")), logic_error, "Error: " \
181  "Integrator_BasisTimesVector: Basis of type \"" << bd_.getType()
182  << "\" is not a vector basis.")
183 
184  TEUCHOS_TEST_FOR_EXCEPTION(multipliers.size() > 3,std::runtime_error,
185  "ERROR: Integrator_BasisTimesVector only supports up to three multipliers!");
186 
187  // Create the field for the vector-valued quantity we're integrating.
188  vector_ = valTag;
189  this->addDependentField(vector_);
190 
191  // Create the field that we're either contributing to or evaluating
192  // (storing).
193  field_ = resTag;
194  if (evalStyle == EvaluatorStyle::CONTRIBUTES)
195  this->addContributedField(field_);
196  else // if (evalStyle == EvaluatorStyle::EVALUATES)
197  this->addEvaluatedField(field_);
198 
199  // Add the dependent field multipliers, if there are any.
200  int i(0);
201  fieldMults_.resize(multipliers.size());
202  for (const auto& fm : multipliers)
203  {
204  fieldMults_[i++] = fm;
205  this->addDependentField(fieldMults_[i - 1]);
206  } // end loop over the field multipliers
207 
208  // Set the name of this object.
209  string n("Integrator_BasisTimesVector<");
210  n += std::to_string(multipliers.size()) + ">(";
211  if (evalStyle == EvaluatorStyle::CONTRIBUTES)
212  n += "Cont";
213  else // if (evalStyle == EvaluatorStyle::EVALUATES)
214  n += "Eval";
215  n += ", " + print<EvalT>() + "): " + field_.fieldTag().name();
216  this->setName(n);
217  } // end of Main Constructor
218 
220  //
221  // ParameterList Constructor
222  //
224  template<typename EvalT, typename Traits>
227  const Teuchos::ParameterList& p)
228  :
230  panzer::EvaluatorStyle::EVALUATES,
231  p.get<std::string>("Residual Name"),
232  p.get<std::string>("Value Name"),
233  (*p.get<Teuchos::RCP<panzer::BasisIRLayout>>("Basis")),
234  (*p.get<Teuchos::RCP<panzer::IntegrationRule>>("IR")),
235  p.get<double>("Multiplier"),
236  p.isType<Teuchos::RCP<const std::vector<std::string>>>
237  ("Field Multipliers") ?
238  (*p.get<Teuchos::RCP<const std::vector<std::string>>>
239  ("Field Multipliers")) : std::vector<std::string>())
240  {
242  using Teuchos::RCP;
243 
244  // Ensure that the input ParameterList didn't contain any bogus entries.
245  RCP<ParameterList> validParams = this->getValidParameters();
246  p.validateParameters(*validParams);
247  } // end of ParameterList Constructor
248 
250  //
251  // postRegistrationSetup()
252  //
254  template<typename EvalT, typename Traits>
255  void
258  typename Traits::SetupData sd,
259  PHX::FieldManager<Traits>& /* fm */)
260  {
261  using panzer::getBasisIndex;
262  using std::size_t;
263 
264  // Get the PHX::Views of the field multipliers.
265  for (size_t i=0; i < fieldMults_.size(); ++i)
266  kokkosFieldMults_[i] = fieldMults_[i].get_static_view();
267 
268  // Determine the number of quadrature points and the dimensionality of the
269  // vector that we're integrating.
270  numQP_ = vector_.extent(1);
271  numDim_ = vector_.extent(2);
272 
273  // Determine the index in the Workset bases for our particular basis name.
274  if (not useDescriptors_)
275  basisIndex_ = getBasisIndex(basisName_, (*sd.worksets_)[0], this->wda);
276 
277  } // end of postRegistrationSetup()
278 
280  //
281  // operator()()
282  //
284  template<typename EvalT, typename Traits>
285  template<int NUM_FIELD_MULT>
286  KOKKOS_INLINE_FUNCTION
287  void
290  const FieldMultTag<NUM_FIELD_MULT>& /* tag */,
291  const typename Kokkos::TeamPolicy<FieldMultTag<NUM_FIELD_MULT>,PHX::exec_space>::member_type& team) const
292  {
294  const int cell = team.league_rank();
295  const int numBases = basis_.extent(1);
296 
297  // Initialize the evaluated field.
298  if (evalStyle_ == EvaluatorStyle::EVALUATES) {
299  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (const int basis) {
300  field_(cell, basis) = 0.0;
301  });
302  }
303 
304  // The following if-block is for the sake of optimization depending on the
305  // number of field multipliers.
306  if constexpr (NUM_FIELD_MULT == 0) {
307  // Loop over the quadrature points and dimensions of our vector fields,
308  // scale the integrand by the multiplier, and then perform the
309  // integration, looping over the bases.
310  for (int qp(0); qp < numQP_; ++qp) {
311  for (int dim(0); dim < numDim_; ++dim) {
312  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (const int basis) {
313  field_(cell, basis) += basis_(cell, basis, qp, dim) * multiplier_ * vector_(cell, qp, dim);
314  });
315  } // end loop over the dimensions of the vector field
316  } // end loop over the quadrature points
317  }
318  else if constexpr (NUM_FIELD_MULT == 1) {
319  // Loop over the quadrature points and dimensions of our vector fields,
320  // scale the integrand by the multiplier and the single field multiplier,
321  // and then perform the actual integration, looping over the bases.
322  for (int qp(0); qp < numQP_; ++qp) {
323  for (int dim(0); dim < numDim_; ++dim) {
324  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (const int basis) {
325  field_(cell, basis) += basis_(cell, basis, qp, dim) * multiplier_ * vector_(cell, qp, dim)
326  * kokkosFieldMults_[0](cell, qp);
327  });
328  } // end loop over the dimensions of the vector field
329  } // end loop over the quadrature points
330  }
331  else if constexpr (NUM_FIELD_MULT == 2) {
332  // Loop over the quadrature points and dimensions of our vector fields,
333  // scale the integrand by the multiplier and the single field multiplier,
334  // and then perform the actual integration, looping over the bases.
335  for (int qp(0); qp < numQP_; ++qp) {
336  for (int dim(0); dim < numDim_; ++dim) {
337  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (const int basis) {
338  field_(cell, basis) += basis_(cell, basis, qp, dim) * multiplier_ * vector_(cell, qp, dim)
339  * kokkosFieldMults_[0](cell, qp) * kokkosFieldMults_[1](cell, qp);
340  });
341  } // end loop over the dimensions of the vector field
342  } // end loop over the quadrature points
343  }
344  else if constexpr (NUM_FIELD_MULT == 3) {
345  // Loop over the quadrature points and dimensions of our vector fields,
346  // scale the integrand by the multiplier and the single field multiplier,
347  // and then perform the actual integration, looping over the bases.
348  for (int qp(0); qp < numQP_; ++qp) {
349  for (int dim(0); dim < numDim_; ++dim) {
350  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (const int basis) {
351  field_(cell, basis) += basis_(cell, basis, qp, dim) * multiplier_ * vector_(cell, qp, dim)
352  * kokkosFieldMults_[0](cell, qp) * kokkosFieldMults_[1](cell, qp) * kokkosFieldMults_[2](cell, qp);
353  });
354  } // end loop over the dimensions of the vector field
355  } // end loop over the quadrature points
356  }
357  else {
358  Kokkos::abort("Panzer_Integrator_BasisTimesVector: NUM_FIELD_MULT out of range!");
359  }
360 
361  } // end of operator()()
362 
364  //
365  // evaluateFields()
366  //
368  template<typename EvalT, typename Traits>
369  void
372  typename Traits::EvalData workset)
373  {
374  using Kokkos::parallel_for;
375  using Kokkos::RangePolicy;
376 
377  // Grab the basis information.
378  const panzer::BasisValues2<double>& bv = useDescriptors_ ?
379  this->wda(workset).getBasisValues(bd_,id_) :
380  *this->wda(workset).bases[basisIndex_];
382  basis_ = useDescriptors_ ? bv.getVectorBasisValues(true) : Array(bv.weighted_basis_vector);
383 
384  // The following if-block is for the sake of optimization depending on the
385  // number of field multipliers. The parallel_fors will loop over the cells
386  // in the Workset and execute operator()() above.
387  if (fieldMults_.size() == 0) {
388  auto policy = panzer::HP::inst().teamPolicy<ScalarT,FieldMultTag<0>,PHX::exec_space>(workset.num_cells);
389  parallel_for("Panzer_Integrator_BasisTimesVector<0>", policy, *this);
390  }
391  else if (fieldMults_.size() == 1) {
392  auto policy = panzer::HP::inst().teamPolicy<ScalarT,FieldMultTag<1>,PHX::exec_space>(workset.num_cells);
393  parallel_for("Panzer_Integrator_BasisTimesVector<1>", policy, *this);
394  }
395  else if (fieldMults_.size() == 2) {
396  auto policy = panzer::HP::inst().teamPolicy<ScalarT,FieldMultTag<2>,PHX::exec_space>(workset.num_cells);
397  parallel_for("Panzer_Integrator_BasisTimesVector<2>", policy, *this);
398  }
399  else if (fieldMults_.size() == 3) {
400  auto policy = panzer::HP::inst().teamPolicy<ScalarT,FieldMultTag<3>,PHX::exec_space>(workset.num_cells);
401  parallel_for("Panzer_Integrator_BasisTimesVector<3>", policy, *this);
402  }
403  } // end of evaluateFields()
404 
406  //
407  // getValidParameters()
408  //
410  template<typename EvalT, typename Traits>
414  {
415  using panzer::BasisIRLayout;
417  using std::string;
418  using std::vector;
420  using Teuchos::RCP;
421  using Teuchos::rcp;
422 
423  // Create a ParameterList with all the valid keys we support.
424  RCP<ParameterList> p = rcp(new ParameterList);
425  p->set<string>("Residual Name", "?");
426  p->set<string>("Value Name", "?");
427  RCP<BasisIRLayout> basis;
428  p->set("Basis", basis);
430  p->set("IR", ir);
431  p->set<double>("Multiplier", 1.0);
433  p->set("Field Multipliers", fms);
434  return p;
435  } // end of getValidParameters()
436 
437 } // end of namespace panzer
438 
439 #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_