Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 std::invalid_argument;
88  using std::logic_error;
89  using std::string;
90  using Teuchos::RCP;
91 
92  // Ensure the input makes sense.
93  TEUCHOS_TEST_FOR_EXCEPTION(resName == "", invalid_argument, "Error: " \
94  "Integrator_BasisTimesScalar called with an empty residual name.")
95  TEUCHOS_TEST_FOR_EXCEPTION(valName == "", invalid_argument, "Error: " \
96  "Integrator_BasisTimesScalar called with an empty value name.")
97  RCP<const PureBasis> tmpBasis = basis.getBasis();
98  TEUCHOS_TEST_FOR_EXCEPTION(not tmpBasis->isScalarBasis(), logic_error,
99  "Error: Integrator_BasisTimesScalar: Basis of type \""
100  << tmpBasis->name() << "\" is not a scalar basis.")
101 
102  // Create the field for the scalar quantity we're integrating.
103  scalar_ = MDField<const ScalarT, Cell, IP>(valName, ir.dl_scalar);
104  this->addDependentField(scalar_);
105 
106  // Create the field that we're either contributing to or evaluating
107  // (storing).
108  field_ = MDField<ScalarT, Cell, BASIS>(resName, basis.functional);
109  if (evalStyle == EvaluatorStyle::CONTRIBUTES)
110  this->addContributedField(field_);
111  else // if (evalStyle == EvaluatorStyle::EVALUATES)
112  this->addEvaluatedField(field_);
113 
114  // Add the dependent field multipliers, if there are any.
115  int i(0);
116  fieldMults_.resize(fmNames.size());
118  View<View<const ScalarT**,typename PHX::DevLayout<ScalarT>::type,PHX::Device>*>("BasisTimesScalar::KokkosFieldMultipliers",
119  fmNames.size());
120  for (const auto& name : fmNames)
121  {
122  fieldMults_[i++] = MDField<const ScalarT, Cell, IP>(name, ir.dl_scalar);
123  this->addDependentField(fieldMults_[i - 1]);
124  } // end loop over the field multipliers
125 
126  // Set the name of this object.
127  string n("Integrator_BasisTimesScalar (");
129  n += "CONTRIBUTES";
130  else // if (evalStyle_ == EvaluatorStyle::EVALUATES)
131  n += "EVALUATES";
132  n += "): " + field_.fieldTag().name();
133  this->setName(n);
134  } // end of Main Constructor
135 
137  //
138  // ParameterList Constructor
139  //
141  template<typename EvalT, typename Traits>
144  const Teuchos::ParameterList& p)
145  :
147  panzer::EvaluatorStyle::EVALUATES,
148  p.get<std::string>("Residual Name"),
149  p.get<std::string>("Value Name"),
150  (*p.get<Teuchos::RCP<panzer::BasisIRLayout>>("Basis")),
151  (*p.get<Teuchos::RCP<panzer::IntegrationRule>>("IR")),
152  p.get<double>("Multiplier"),
153  p.isType<Teuchos::RCP<const std::vector<std::string>>>
154  ("Field Multipliers") ?
155  (*p.get<Teuchos::RCP<const std::vector<std::string>>>
156  ("Field Multipliers")) : std::vector<std::string>())
157  {
159  using Teuchos::RCP;
160 
161  // Ensure that the input ParameterList didn't contain any bogus entries.
162  RCP<ParameterList> validParams = this->getValidParameters();
163  p.validateParameters(*validParams);
164  } // end of ParameterList Constructor
165 
167  //
168  // postRegistrationSetup()
169  //
171  template<typename EvalT, typename Traits>
172  void
175  typename Traits::SetupData sd,
176  PHX::FieldManager<Traits>& /* fm */)
177  {
178  using panzer::getBasisIndex;
179  using std::size_t;
180 
181  // Get the Kokkos::Views of the field multipliers.
182  for (size_t i(0); i < fieldMults_.size(); ++i)
183  kokkosFieldMults_(i) = fieldMults_[i].get_static_view();
184 
185  // Determine the index in the Workset bases for our particular basis name.
186  basisIndex_ = getBasisIndex(basisName_, (*sd.worksets_)[0], this->wda);
187  } // end of postRegistrationSetup()
188 
190  //
191  // operator()()
192  //
194  template<typename EvalT, typename Traits>
195  template<int NUM_FIELD_MULT>
196  KOKKOS_INLINE_FUNCTION
197  void
200  const FieldMultTag<NUM_FIELD_MULT>& /* tag */,
201  const size_t& cell) const
202  {
204 
205  // Initialize the evaluated field.
206  const int numQP(scalar_.extent(1)), numBases(basis_.extent(1));
207  if (evalStyle_ == EvaluatorStyle::EVALUATES)
208  for (int basis(0); basis < numBases; ++basis)
209  field_(cell, basis) = 0.0;
210 
211  // The following if-block is for the sake of optimization depending on the
212  // number of field multipliers.
213  ScalarT tmp;
214  if (NUM_FIELD_MULT == 0)
215  {
216  // Loop over the quadrature points, scale the integrand by the
217  // multiplier, and then perform the actual integration, looping over the
218  // bases.
219  for (int qp(0); qp < numQP; ++qp)
220  {
221  tmp = multiplier_ * scalar_(cell, qp);
222  for (int basis(0); basis < numBases; ++basis)
223  field_(cell, basis) += basis_(cell, basis, qp) * tmp;
224  } // end loop over the quadrature points
225  }
226  else if (NUM_FIELD_MULT == 1)
227  {
228  // Loop over the quadrature points, scale the integrand by the multiplier
229  // and the single field multiplier, and then perform the actual
230  // integration, looping over the bases.
231  for (int qp(0); qp < numQP; ++qp)
232  {
233  tmp = multiplier_ * scalar_(cell, qp) * kokkosFieldMults_(0)(cell, qp);
234  for (int basis(0); basis < numBases; ++basis)
235  field_(cell, basis) += basis_(cell, basis, qp) * tmp;
236  } // end loop over the quadrature points
237  }
238  else
239  {
240  // Loop over the quadrature points and pre-multiply all the field
241  // multipliers together, scale the integrand by the multiplier and
242  // the combination of the field multipliers, and then perform the actual
243  // integration, looping over the bases.
244  const int numFieldMults(kokkosFieldMults_.extent(0));
245  for (int qp(0); qp < numQP; ++qp)
246  {
247  ScalarT fieldMultsTotal(1);
248  for (int fm(0); fm < numFieldMults; ++fm)
249  fieldMultsTotal *= kokkosFieldMults_(fm)(cell, qp);
250  tmp = multiplier_ * scalar_(cell, qp) * fieldMultsTotal;
251  for (int basis(0); basis < numBases; ++basis)
252  field_(cell, basis) += basis_(cell, basis, qp) * tmp;
253  } // end loop over the quadrature points
254  } // end if (NUM_FIELD_MULT == something)
255  } // end of operator()()
256 
258  //
259  // evaluateFields()
260  //
262  template<typename EvalT, typename Traits>
263  void
266  typename Traits::EvalData workset)
267  {
268  using Kokkos::parallel_for;
269  using Kokkos::RangePolicy;
270 
271  // Grab the basis information.
272  basis_ = this->wda(workset).bases[basisIndex_]->weighted_basis_scalar;
273 
274  // The following if-block is for the sake of optimization depending on the
275  // number of field multipliers. The parallel_fors will loop over the cells
276  // in the Workset and execute operator()() above.
277  if (fieldMults_.size() == 0)
278  parallel_for(RangePolicy<FieldMultTag<0>>(0, workset.num_cells), *this);
279  else if (fieldMults_.size() == 1)
280  parallel_for(RangePolicy<FieldMultTag<1>>(0, workset.num_cells), *this);
281  else
282  parallel_for(RangePolicy<FieldMultTag<-1>>(0, workset.num_cells), *this);
283  } // end of evaluateFields()
284 
286  //
287  // getValidParameters()
288  //
290  template<typename EvalT, typename TRAITS>
294  {
295  using panzer::BasisIRLayout;
297  using std::string;
298  using std::vector;
300  using Teuchos::RCP;
301  using Teuchos::rcp;
302 
303  // Create a ParameterList with all the valid keys we support.
304  RCP<ParameterList> p = rcp(new ParameterList);
305  p->set<string>("Residual Name", "?");
306  p->set<string>("Value Name", "?");
307  RCP<BasisIRLayout> basis;
308  p->set("Basis", basis);
310  p->set("IR", ir);
311  p->set<double>("Multiplier", 1.0);
313  p->set("Field Multipliers", fms);
314  return p;
315  } // end of getValidParameters()
316 
317 } // end of namespace panzer
318 
319 #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...
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.
KOKKOS_INLINE_FUNCTION void operator()(const FieldMultTag< NUM_FIELD_MULT > &tag, const std::size_t &cell) const
Perform the integration.
Teuchos::RCP< const PureBasis > getBasis() 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.
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...
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.
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 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 sd, PHX::FieldManager< Traits > &fm)
Post-Registration Setup.
const std::vector< std::pair< int, LocalOrdinal > > &pid_and_lid const
Description and data layouts associated with a particular basis.
Teuchos::RCP< PHX::DataLayout > functional
&lt;Cell,Basis&gt;
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< Teuchos::ParameterList > getValidParameters() const
Get Valid Parameters.
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_