Panzer  Version of the Day
Panzer_Integrator_BasisTimesScalar_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_BASISTIMESSCALAR_IMPL_HPP
44 #define PANZER_INTEGRATOR_BASISTIMESSCALAR_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  multiplier_(multiplier),
78  basisName_(basis.name())
79  {
80  using Kokkos::View;
81  using panzer::BASIS;
82  using panzer::Cell;
84  using panzer::IP;
85  using panzer::PureBasis;
86  using PHX::MDField;
87  using PHX::typeAsString;
88  using std::invalid_argument;
89  using std::logic_error;
90  using std::string;
91  using Teuchos::RCP;
92 
93  // Ensure the input makes sense.
94  TEUCHOS_TEST_FOR_EXCEPTION(resName == "", invalid_argument, "Error: " \
95  "Integrator_BasisTimesScalar called with an empty residual name.")
96  TEUCHOS_TEST_FOR_EXCEPTION(valName == "", invalid_argument, "Error: " \
97  "Integrator_BasisTimesScalar called with an empty value name.")
98  RCP<const PureBasis> tmpBasis = basis.getBasis();
99  TEUCHOS_TEST_FOR_EXCEPTION(not tmpBasis->isScalarBasis(), logic_error,
100  "Error: Integrator_BasisTimesScalar: Basis of type \""
101  << tmpBasis->name() << "\" is not a scalar basis.")
102 
103  // Create the field for the scalar quantity we're integrating.
104  scalar_ = MDField<const ScalarT, Cell, IP>(valName, ir.dl_scalar);
105  this->addDependentField(scalar_);
106 
107  // Create the field that we're either contributing to or evaluating
108  // (storing).
109  field_ = MDField<ScalarT, Cell, BASIS>(resName, basis.functional);
111  this->addContributedField(field_);
112  else // if (evalStyle == EvaluatorStyle::EVALUATES)
113  this->addEvaluatedField(field_);
114 
115  // Add the dependent field multipliers, if there are any.
116  int i(0);
117  fieldMults_.resize(fmNames.size());
119  View<View<const ScalarT**>*>("BasisTimesScalar::KokkosFieldMultipliers",
120  fmNames.size());
121  for (const auto& name : fmNames)
122  {
123  fieldMults_[i++] = MDField<const ScalarT, Cell, IP>(name, ir.dl_scalar);
124  this->addDependentField(fieldMults_[i - 1]);
125  } // end loop over the field multipliers
126 
127  // Set the name of this object.
128  string n("Integrator_BasisTimesScalar (");
130  n += "CONTRIBUTES";
131  else // if (evalStyle_ == EvaluatorStyle::EVALUATES)
132  n += "EVALUATES";
133  n += "): " + field_.fieldTag().name();
134  this->setName(n);
135  } // end of Main Constructor
136 
138  //
139  // ParameterList Constructor
140  //
142  template<typename EvalT, typename Traits>
145  const Teuchos::ParameterList& p)
146  :
149  p.get<std::string>("Residual Name"),
150  p.get<std::string>("Value Name"),
151  (*p.get<Teuchos::RCP<panzer::BasisIRLayout>>("Basis")),
152  (*p.get<Teuchos::RCP<panzer::IntegrationRule>>("IR")),
153  p.get<double>("Multiplier"),
154  p.isType<Teuchos::RCP<const std::vector<std::string>>>
155  ("Field Multipliers") ?
156  (*p.get<Teuchos::RCP<const std::vector<std::string>>>
157  ("Field Multipliers")) : std::vector<std::string>())
158  {
160  using Teuchos::RCP;
161 
162  // Ensure that the input ParameterList didn't contain any bogus entries.
163  RCP<ParameterList> validParams = this->getValidParameters();
164  p.validateParameters(*validParams);
165  } // end of ParameterList Constructor
166 
168  //
169  // postRegistrationSetup()
170  //
172  template<typename EvalT, typename Traits>
173  void
176  typename Traits::SetupData sd,
177  PHX::FieldManager<Traits>& /* fm */)
178  {
179  using panzer::getBasisIndex;
180  using std::size_t;
181 
182  // Get the Kokkos::Views of the field multipliers.
183  for (size_t i(0); i < fieldMults_.size(); ++i)
184  kokkosFieldMults_(i) = fieldMults_[i].get_static_view();
185 
186  // Determine the number of nodes and quadrature points.
187  numQP_ = scalar_.extent(1);
188 
189  // Determine the index in the Workset bases for our particular basis name.
190  basisIndex_ = getBasisIndex(basisName_, (*sd.worksets_)[0], this->wda);
191  } // end of postRegistrationSetup()
192 
194  //
195  // operator()()
196  //
198  template<typename EvalT, typename Traits>
199  template<int NUM_FIELD_MULT>
200  KOKKOS_INLINE_FUNCTION
201  void
204  const FieldMultTag<NUM_FIELD_MULT>& /* tag */,
205  const size_t& cell) const
206  {
208 
209  // Initialize the evaluated field.
210  const int numBases(basis_.extent(1));
211  if (evalStyle_ == EvaluatorStyle::EVALUATES)
212  for (int basis(0); basis < numBases; ++basis)
213  field_(cell, basis) = 0.0;
214 
215  // The following if-block is for the sake of optimization depending on the
216  // number of field multipliers.
217  ScalarT tmp;
218  if (NUM_FIELD_MULT == 0)
219  {
220  // Loop over the quadrature points, scale the integrand by the
221  // multiplier, and then perform the actual integration, looping over the
222  // bases.
223  for (int qp(0); qp < numQP_; ++qp)
224  {
225  tmp = multiplier_ * scalar_(cell, qp);
226  for (int basis(0); basis < numBases; ++basis)
227  field_(cell, basis) += basis_(cell, basis, qp) * tmp;
228  } // end loop over the quadrature points
229  }
230  else if (NUM_FIELD_MULT == 1)
231  {
232  // Loop over the quadrature points, scale the integrand by the multiplier
233  // and the single field multiplier, and then perform the actual
234  // integration, looping over the bases.
235  for (int qp(0); qp < numQP_; ++qp)
236  {
237  tmp = multiplier_ * scalar_(cell, qp) * kokkosFieldMults_(0)(cell, qp);
238  for (int basis(0); basis < numBases; ++basis)
239  field_(cell, basis) += basis_(cell, basis, qp) * tmp;
240  } // end loop over the quadrature points
241  }
242  else
243  {
244  // Loop over the quadrature points and pre-multiply all the field
245  // multipliers together, scale the integrand by the multiplier and
246  // the combination of the field multipliers, and then perform the actual
247  // integration, looping over the bases.
248  const int numFieldMults(kokkosFieldMults_.extent(0));
249  for (int qp(0); qp < numQP_; ++qp)
250  {
251  ScalarT fieldMultsTotal(1);
252  for (int fm(0); fm < numFieldMults; ++fm)
253  fieldMultsTotal *= kokkosFieldMults_(fm)(cell, qp);
254  tmp = multiplier_ * scalar_(cell, qp) * fieldMultsTotal;
255  for (int basis(0); basis < numBases; ++basis)
256  field_(cell, basis) += basis_(cell, basis, qp) * tmp;
257  } // end loop over the quadrature points
258  } // end if (NUM_FIELD_MULT == something)
259  } // end of operator()()
260 
262  //
263  // evaluateFields()
264  //
266  template<typename EvalT, typename Traits>
267  void
270  typename Traits::EvalData workset)
271  {
272  using Kokkos::parallel_for;
273  using Kokkos::RangePolicy;
274 
275  // Grab the basis information.
276  basis_ = this->wda(workset).bases[basisIndex_]->weighted_basis_scalar;
277 
278  // The following if-block is for the sake of optimization depending on the
279  // number of field multipliers. The parallel_fors will loop over the cells
280  // in the Workset and execute operator()() above.
281  if (fieldMults_.size() == 0)
282  parallel_for(RangePolicy<FieldMultTag<0>>(0, workset.num_cells), *this);
283  else if (fieldMults_.size() == 1)
284  parallel_for(RangePolicy<FieldMultTag<1>>(0, workset.num_cells), *this);
285  else
286  parallel_for(RangePolicy<FieldMultTag<-1>>(0, workset.num_cells), *this);
287  } // end of evaluateFields()
288 
290  //
291  // getValidParameters()
292  //
294  template<typename EvalT, typename TRAITS>
298  {
299  using panzer::BasisIRLayout;
301  using std::string;
302  using std::vector;
304  using Teuchos::RCP;
305  using Teuchos::rcp;
306 
307  // Create a ParameterList with all the valid keys we support.
308  RCP<ParameterList> p = rcp(new ParameterList);
309  p->set<string>("Residual Name", "?");
310  p->set<string>("Value Name", "?");
311  RCP<BasisIRLayout> basis;
312  p->set("Basis", basis);
314  p->set("IR", ir);
315  p->set<double>("Multiplier", 1.0);
317  p->set("Field Multipliers", fms);
318  return p;
319  } // end of getValidParameters()
320 
321 } // end of namespace panzer
322 
323 #endif // PANZER_INTEGRATOR_BASISTIMESSCALAR_IMPL_HPP
const panzer::EvaluatorStyle evalStyle_
An enum determining the behavior of this Evaluator.
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)
This empty struct allows us to optimize operator()() depending on the number of field multipliers...
Teuchos::RCP< const PureBasis > getBasis() const
Integrator_BasisTimesScalar(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.
EvaluatorStyle
An indication of how an Evaluator will behave.
void validateParameters(ParameterList const &validParamList, int const depth=1000, EValidateUsed const validateUsed=VALIDATE_USED_ENABLED, EValidateDefaults const validateDefaults=VALIDATE_DEFAULTS_ENABLED) const
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.
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
Get Valid Parameters.
Kokkos::View< Kokkos::View< const ScalarT ** > * > kokkosFieldMults_
The Kokkos::View representation of the (possibly empty) list of fields that are multipliers out in fr...
panzer::EvaluatorStyle evalStyle
The EvaluatorStyle of the parent Integrator_CurlBasisDotVector object.
void operator()(const FieldMultTag< NUM_FIELD_MULT > &tag, const std::size_t &cell) const
Perform the integration.
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.
PHX::MDField< const ScalarT, panzer::Cell, panzer::IP > scalar_
A field representing the scalar function we&#39;re integrating ( ).
typename EvalT::ScalarT ScalarT
The scalar type.
void evaluateFields(typename Traits::EvalData workset)
Evaluate Fields.
void postRegistrationSetup(typename Traits::SetupData sd, PHX::FieldManager< Traits > &fm)
Post-Registration Setup.
Description and data layouts associated with a particular basis.
Teuchos::RCP< PHX::DataLayout > functional
<Cell,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< const std::vector< panzer::Workset > > worksets_