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"
56 
57 namespace panzer
58 {
60  //
61  // Main Constructor
62  //
64  template<typename EvalT, typename Traits>
68  const std::string& resName,
69  const std::string& valName,
70  const panzer::BasisIRLayout& basis,
71  const panzer::IntegrationRule& ir,
72  const double& multiplier, /* = 1 */
73  const std::vector<std::string>& fmNames /* =
74  std::vector<std::string>() */)
75  :
76  evalStyle_(evalStyle),
77  useDescriptors_(false),
78  multiplier_(multiplier),
79  basisName_(basis.name())
80  {
81  using Kokkos::View;
82  using panzer::BASIS;
83  using panzer::Cell;
85  using panzer::IP;
86  using panzer::PureBasis;
87  using PHX::MDField;
88  using PHX::typeAsString;
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_BasisTimesVector called with an empty residual name.")
97  TEUCHOS_TEST_FOR_EXCEPTION(valName == "", invalid_argument, "Error: " \
98  "Integrator_BasisTimesVector called with an empty value name.")
99  RCP<const PureBasis> tmpBasis = basis.getBasis();
100  TEUCHOS_TEST_FOR_EXCEPTION(not tmpBasis->isVectorBasis(), logic_error,
101  "Error: Integrator_BasisTimesVector: Basis of type \""
102  << tmpBasis->name() << "\" is not a vector basis.");
103  TEUCHOS_TEST_FOR_EXCEPTION(not tmpBasis->requiresOrientations(),
104  logic_error, "Integrator_BasisTimesVector: Basis of type \""
105  << tmpBasis->name() << "\" does not require orientations. This seems " \
106  "very strange, so I'm failing.");
107 
108  // Create the field for the vector-valued quantity we're integrating.
109  vector_ = MDField<const ScalarT, Cell, IP, Dim>(valName, ir.dl_vector);
110  this->addDependentField(vector_);
111 
112  // Create the field that we're either contributing to or evaluating
113  // (storing).
114  field_ = MDField<ScalarT, Cell, BASIS>(resName, basis.functional);
115  if (evalStyle == EvaluatorStyle::CONTRIBUTES)
116  this->addContributedField(field_);
117  else // if (evalStyle == EvaluatorStyle::EVALUATES)
118  this->addEvaluatedField(field_);
119 
120  // Add the dependent field multipliers, if there are any.
121  int i(0);
122  fieldMults_.resize(fmNames.size());
124  View<View<const ScalarT**,typename PHX::DevLayout<ScalarT>::type,PHX::Device>*>("BasisTimesVector::KokkosFieldMultipliers",
125  fmNames.size());
126  for (const auto& name : fmNames)
127  {
128  fieldMults_[i++] = MDField<const ScalarT, Cell, IP>(name, ir.dl_scalar);
129  this->addDependentField(fieldMults_[i - 1]);
130  } // end loop over the field multipliers
131 
132  // Set the name of this object.
133  string n("Integrator_BasisTimesVector (");
134  if (evalStyle == EvaluatorStyle::CONTRIBUTES)
135  n += "Cont";
136  else // if (evalStyle == EvaluatorStyle::EVALUATES)
137  n += "Eval";
138  n += ", " + typeAsString<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 Kokkos::View;
166  using panzer::BASIS;
167  using panzer::Cell;
169  using panzer::IP;
170  using panzer::PureBasis;
171  using PHX::MDField;
172  using PHX::typeAsString;
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  // Create the field for the vector-valued quantity we're integrating.
185  vector_ = valTag;
186  this->addDependentField(vector_);
187 
188  // Create the field that we're either contributing to or evaluating
189  // (storing).
190  field_ = resTag;
191  if (evalStyle == EvaluatorStyle::CONTRIBUTES)
192  this->addContributedField(field_);
193  else // if (evalStyle == EvaluatorStyle::EVALUATES)
194  this->addEvaluatedField(field_);
195 
196  // Add the dependent field multipliers, if there are any.
197  int i(0);
198  fieldMults_.resize(multipliers.size());
200  View<View<const ScalarT**,typename PHX::DevLayout<ScalarT>::type,PHX::Device>*>("BasisTimesVector::KokkosFieldMultipliers",
201  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  if (evalStyle == EvaluatorStyle::CONTRIBUTES)
211  n += "Cont";
212  else // if (evalStyle == EvaluatorStyle::EVALUATES)
213  n += "Eval";
214  n += ", " + typeAsString<EvalT>() + "): " + field_.fieldTag().name();
215  this->setName(n);
216  } // end of Main Constructor
217 
219  //
220  // ParameterList Constructor
221  //
223  template<typename EvalT, typename Traits>
226  const Teuchos::ParameterList& p)
227  :
229  panzer::EvaluatorStyle::EVALUATES,
230  p.get<std::string>("Residual Name"),
231  p.get<std::string>("Value Name"),
232  (*p.get<Teuchos::RCP<panzer::BasisIRLayout>>("Basis")),
233  (*p.get<Teuchos::RCP<panzer::IntegrationRule>>("IR")),
234  p.get<double>("Multiplier"),
235  p.isType<Teuchos::RCP<const std::vector<std::string>>>
236  ("Field Multipliers") ?
237  (*p.get<Teuchos::RCP<const std::vector<std::string>>>
238  ("Field Multipliers")) : std::vector<std::string>())
239  {
241  using Teuchos::RCP;
242 
243  // Ensure that the input ParameterList didn't contain any bogus entries.
244  RCP<ParameterList> validParams = this->getValidParameters();
245  p.validateParameters(*validParams);
246  } // end of ParameterList Constructor
247 
249  //
250  // postRegistrationSetup()
251  //
253  template<typename EvalT, typename Traits>
254  void
257  typename Traits::SetupData sd,
258  PHX::FieldManager<Traits>& /* fm */)
259  {
260  using panzer::getBasisIndex;
261  using std::size_t;
262 
263  // Get the Kokkos::Views of the field multipliers.
264  for (size_t i(0); i < fieldMults_.size(); ++i)
265  kokkosFieldMults_(i) = fieldMults_[i].get_static_view();
266 
267  // Determine the number of quadrature points and the dimensionality of the
268  // vector that we're integrating.
269  numQP_ = vector_.extent(1);
270  numDim_ = vector_.extent(2);
271 
272  // Determine the index in the Workset bases for our particular basis name.
273  if (not useDescriptors_)
274  basisIndex_ = getBasisIndex(basisName_, (*sd.worksets_)[0], this->wda);
275  } // end of postRegistrationSetup()
276 
278  //
279  // operator()()
280  //
282  template<typename EvalT, typename Traits>
283  template<int NUM_FIELD_MULT>
284  KOKKOS_INLINE_FUNCTION
285  void
288  const FieldMultTag<NUM_FIELD_MULT>& /* tag */,
289  const size_t& cell) const
290  {
292 
293  // Initialize the evaluated field.
294  const int numBases(basis_.extent(1));
295  if (evalStyle_ == EvaluatorStyle::EVALUATES)
296  for (int basis(0); basis < numBases; ++basis)
297  field_(cell, basis) = 0.0;
298 
299  // The following if-block is for the sake of optimization depending on the
300  // number of field multipliers.
301  ScalarT tmp;
302  if (NUM_FIELD_MULT == 0)
303  {
304  // Loop over the quadrature points and dimensions of our vector fields,
305  // scale the integrand by the multiplier, and then perform the
306  // integration, looping over the bases.
307  for (int qp(0); qp < numQP_; ++qp)
308  {
309  for (int dim(0); dim < numDim_; ++dim)
310  {
311  tmp = multiplier_ * vector_(cell, qp, dim);
312  for (int basis(0); basis < numBases; ++basis)
313  field_(cell, basis) += basis_(cell, basis, qp, dim) * tmp;
314  } // end loop over the dimensions of the vector field
315  } // end loop over the quadrature points
316  }
317  else if (NUM_FIELD_MULT == 1)
318  {
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  {
324  for (int dim(0); dim < numDim_; ++dim)
325  {
326  tmp = multiplier_ * vector_(cell, qp, dim) *
327  kokkosFieldMults_(0)(cell, qp);
328  for (int basis(0); basis < numBases; ++basis)
329  field_(cell, basis) += basis_(cell, basis, qp, dim) * tmp;
330  } // end loop over the dimensions of the vector field
331  } // end loop over the quadrature points
332  }
333  else
334  {
335  // Loop over the quadrature points and pre-multiply all the field
336  // multipliers together. Then loop over the dimensions of our vector
337  // fields, scale the integrand by the multiplier and the combination of
338  // the field multipliers, and then perform the actual integration,
339  // looping over the bases.
340  const int numFieldMults(kokkosFieldMults_.extent(0));
341  for (int qp(0); qp < numQP_; ++qp)
342  {
343  ScalarT fieldMultsTotal(1);
344  for (int fm(0); fm < numFieldMults; ++fm)
345  fieldMultsTotal *= kokkosFieldMults_(fm)(cell, qp);
346  for (int dim(0); dim < numDim_; ++dim)
347  {
348  tmp = multiplier_ * vector_(cell, qp, dim) * fieldMultsTotal;
349  for (int basis(0); basis < numBases; ++basis)
350  field_(cell, basis) += basis_(cell, basis, qp, dim) * tmp;
351  } // end loop over the dimensions of the vector field
352  } // end loop over the quadrature points
353  } // end if (NUM_FIELD_MULT == something)
354  } // end of operator()()
355 
357  //
358  // evaluateFields()
359  //
361  template<typename EvalT, typename Traits>
362  void
365  typename Traits::EvalData workset)
366  {
367  using Kokkos::parallel_for;
368  using Kokkos::RangePolicy;
369 
370  // Grab the basis information.
371  const panzer::BasisValues2<double>& bv = useDescriptors_ ?
372  this->wda(workset).getBasisValues(bd_,id_) :
373  *this->wda(workset).bases[basisIndex_];
374  basis_ = bv.weighted_basis_vector;
375 
376  // The following if-block is for the sake of optimization depending on the
377  // number of field multipliers. The parallel_fors will loop over the cells
378  // in the Workset and execute operator()() above.
379  if (fieldMults_.size() == 0)
380  parallel_for(RangePolicy<FieldMultTag<0>>(0, workset.num_cells), *this);
381  else if (fieldMults_.size() == 1)
382  parallel_for(RangePolicy<FieldMultTag<1>>(0, workset.num_cells), *this);
383  else
384  parallel_for(RangePolicy<FieldMultTag<-1>>(0, workset.num_cells), *this);
385  } // end of evaluateFields()
386 
388  //
389  // getValidParameters()
390  //
392  template<typename EvalT, typename Traits>
396  {
397  using panzer::BasisIRLayout;
399  using std::string;
400  using std::vector;
402  using Teuchos::RCP;
403  using Teuchos::rcp;
404 
405  // Create a ParameterList with all the valid keys we support.
406  RCP<ParameterList> p = rcp(new ParameterList);
407  p->set<string>("Residual Name", "?");
408  p->set<string>("Value Name", "?");
409  RCP<BasisIRLayout> basis;
410  p->set("Basis", basis);
412  p->set("IR", ir);
413  p->set<double>("Multiplier", 1.0);
415  p->set("Field Multipliers", fms);
416  return p;
417  } // end of getValidParameters()
418 
419 } // end of namespace panzer
420 
421 #endif // __Panzer_Integerator_BasisTimesVector_impl_hpp__
void operator()(const FieldMultTag< NUM_FIELD_MULT > &tag, const std::size_t &cell) const
Perform the integration.
#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 ( ).
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.
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
const std::vector< std::pair< int, LocalOrdinal > > &pid_and_lid const
Teuchos::RCP< PHX::DataLayout > dl_vector
Data layout for vector fields.
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...
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.
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_