Panzer  Version of the Day
Panzer_Integrator_GradBasisTimesScalar_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_EVALUATOR_GRADBASISTIMESSCALAR_IMPL_HPP
44 #define PANZER_EVALUATOR_GRADBASISTIMESSCALAR_IMPL_HPP
45 
47 //
48 // Include Files
49 //
51 
52 // Intrepid2
53 #include "Intrepid2_FunctionSpaceTools.hpp"
54 
55 // Kokkos
56 #include "Kokkos_ViewFactory.hpp"
57 
58 // Panzer
59 #include "Panzer_BasisIRLayout.hpp"
62 
63 namespace panzer
64 {
66  //
67  // Constructor
68  //
70  template<typename EvalT, typename Traits>
74  const std::vector<std::string>& resNames,
75  const std::string& scalar,
76  const panzer::BasisIRLayout& basis,
77  const panzer::IntegrationRule& ir,
78  const double& multiplier, /* = 1 */
79  const std::vector<std::string>& fmNames) /* =
80  std::vector<std::string>() */
81  :
82  evalStyle_(evalStyle),
83  multiplier_(multiplier),
84  numDims_(ir.dl_vector->extent(2)),
85  basisName_(basis.name())
86  {
87  using Kokkos::View;
88  using panzer::BASIS;
89  using panzer::Cell;
91  using panzer::IP;
92  using panzer::PureBasis;
93  using PHX::Device;
94  using PHX::DevLayout;
95  using PHX::MDField;
96  using std::invalid_argument;
97  using std::logic_error;
98  using std::string;
99  using Teuchos::RCP;
100 
101  // Ensure the input makes sense.
102  for (const auto& name : resNames)
103  TEUCHOS_TEST_FOR_EXCEPTION(name == "", invalid_argument, "Error: " \
104  "Integrator_GradBasisTimesScalar called with an empty residual name.")
105  TEUCHOS_TEST_FOR_EXCEPTION(scalar == "", invalid_argument, "Error: " \
106  "Integrator_GradBasisTimesScalar called with an empty scalar name.")
107  TEUCHOS_TEST_FOR_EXCEPTION(numDims_ != static_cast<int>(resNames.size()),
108  logic_error, "Error: Integrator_GradBasisTimesScalar: The number of " \
109  "residuals must equal the number of gradient components (mesh " \
110  "dimensions).")
111  RCP<const PureBasis> tmpBasis = basis.getBasis();
112  TEUCHOS_TEST_FOR_EXCEPTION(not tmpBasis->supportsGrad(), logic_error,
113  "Error: Integrator_GradBasisTimesScalar: Basis of type \""
114  << tmpBasis->name() << "\" does not support the gradient operation.")
115 
116  // Create the field or the scalar function we're integrating.
117  scalar_ = MDField<const ScalarT, Cell, IP>(scalar, ir.dl_scalar);
118  this->addDependentField(scalar_);
119 
120  // Create the fields that we're either contributing to or evaluating
121  // (storing).
122  fields_ = Kokkos::View<PHX::MDField<ScalarT, Cell, BASIS>*>("Integrator_GradBasisTimesScalar",resNames.size());
123  {
124  int i=0;
125  for (const auto& name : resNames)
126  fields_(i++) = MDField<ScalarT, Cell, BASIS>(name, basis.functional);
127  }
128  for (std::size_t i=0; i<fields_.extent(0); ++i) {
129  const auto& field = fields_(i);
131  this->addContributedField(field);
132  else // if (evalStyle_ == EvaluatorStyle::EVALUATES)
133  this->addEvaluatedField(field);
134  }
135 
136  // Add the dependent field multipliers, if there are any.
137  int i = 0;
138  fieldMults_.resize(fmNames.size());
140  typename DevLayout<ScalarT>::type, Device>*>(
141  "GradBasisTimesScalar::KokkosFieldMultipliers", fmNames.size());
142  for (const auto& name : fmNames)
143  {
144  fieldMults_[i++] = MDField<const ScalarT, Cell, IP>(name, ir.dl_scalar);
145  this->addDependentField(fieldMults_[i - 1]);
146  } // end loop over the field multipliers
147 
148  // Set the name of this object.
149  string n("Integrator_GradBasisTimesScalar (");
151  n += "CONTRIBUTES";
152  else // if (evalStyle_ == EvaluatorStyle::EVALUATES)
153  n += "EVALUATES";
154  n += "): {";
155  for (size_t j=0; j < fields_.extent(0) - 1; ++j)
156  n += fields_[j].fieldTag().name() + ", ";
157  n += fields_(fields_.extent(0)-1).fieldTag().name() + "}";
158  this->setName(n);
159  } // end of Constructor
160 
162  //
163  // ParameterList Constructor
164  //
166  template<typename EvalT, typename Traits>
169  const Teuchos::ParameterList& p)
170  :
173  p.get<const std::vector<std::string>>("Residual Names"),
174  p.get<std::string>("Scalar Name"),
175  (*p.get<Teuchos::RCP<panzer::BasisIRLayout>>("Basis")),
176  (*p.get<Teuchos::RCP<panzer::IntegrationRule>>("IR")),
177  p.get<double>("Multiplier"),
178  p.isType<Teuchos::RCP<const std::vector<std::string>>>
179  ("Field Multipliers") ?
180  (*p.get<Teuchos::RCP<const std::vector<std::string>>>
181  ("Field Multipliers")) : std::vector<std::string>())
182  {
183  using Teuchos::ParameterList;
184  using Teuchos::RCP;
185 
186  // Ensure that the input ParameterList didn't contain any bogus entries.
187  RCP<ParameterList> validParams = this->getValidParameters();
188  p.validateParameters(*validParams);
189  } // end of ParameterList Constructor
190 
192  //
193  // postRegistrationSetup()
194  //
196  template<typename EvalT, typename Traits>
197  void
200  typename Traits::SetupData sd,
201  PHX::FieldManager<Traits>& /* fm */)
202  {
204  using panzer::getBasisIndex;
205  using PHX::Device;
206 
207  // Get the Kokkos::Views of the field multipliers.
208  for (size_t i(0); i < fieldMults_.size(); ++i)
209  kokkosFieldMults_(i) = fieldMults_[i].get_static_view();
210  Device().fence();
211 
212  // Determine the index in the Workset bases for our particular basis name.
213  basisIndex_ = getBasisIndex(basisName_, (*sd.worksets_)[0], this->wda);
214  } // end of postRegistrationSetup()
215 
216  template<typename EvalT, typename Traits>
217  template<int NUM_FIELD_MULT>
218  KOKKOS_INLINE_FUNCTION
219  void
222  const FieldMultTag<NUM_FIELD_MULT>& /* tag */,
223  const size_t& cell) const
224  {
226 
227  // Initialize the evaluated fields.
228  const int numQP(scalar_.extent(1)), numBases(fields_[0].extent(1));
229  if (evalStyle_ == EvaluatorStyle::EVALUATES)
230  for (int dim(0); dim < numDims_; ++dim)
231  for (int basis(0); basis < numBases; ++basis)
232  fields_[dim](cell, basis) = 0.0;
233 
234  // The following if-block is for the sake of optimization depending on the
235  // number of field multipliers.
236  ScalarT tmp;
237  if (NUM_FIELD_MULT == 0)
238  {
239  // Loop over the quadrature points, scale the integrand by the
240  // multiplier, and then perform the actual integration, looping over the
241  // bases and dimensions.
242  for (int qp(0); qp < numQP; ++qp)
243  {
244  tmp = multiplier_ * scalar_(cell, qp);
245  for (int basis(0); basis < numBases; ++basis)
246  for (int dim(0); dim < numDims_; ++dim)
247  fields_[dim](cell, basis) += basis_(cell, basis, qp, dim) * tmp;
248  } // end loop over the quadrature points
249  }
250  else if (NUM_FIELD_MULT == 1)
251  {
252  // Loop over the quadrature points, scale the integrand by the multiplier
253  // and the single field multiplier, and then perform the actual
254  // integration, looping over the bases and dimensions.
255  for (int qp(0); qp < numQP; ++qp)
256  {
257  tmp = multiplier_ * scalar_(cell, qp) * kokkosFieldMults_(0)(cell, qp);
258  for (int basis(0); basis < numBases; ++basis)
259  for (int dim(0); dim < numDims_; ++dim)
260  fields_[dim](cell, basis) += basis_(cell, basis, qp, dim) * tmp;
261  } // end loop over the quadrature points
262  }
263  else
264  {
265  // Loop over the quadrature points and pre-multiply all the field
266  // multipliers together, scale the integrand by the multiplier and
267  // the combination of the field multipliers, and then perform the actual
268  // integration, looping over the bases and dimensions.
269  const int numFieldMults(kokkosFieldMults_.extent(0));
270  for (int qp(0); qp < numQP; ++qp)
271  {
272  ScalarT fieldMultsTotal(1);
273  for (int fm(0); fm < numFieldMults; ++fm)
274  fieldMultsTotal *= kokkosFieldMults_(fm)(cell, qp);
275  tmp = multiplier_ * scalar_(cell, qp) * fieldMultsTotal;
276  for (int basis(0); basis < numBases; ++basis)
277  for (int dim(0); dim < numDims_; ++dim)
278  fields_[dim](cell, basis) += basis_(cell, basis, qp, dim) * tmp;
279  } // end loop over the quadrature points
280  } // end if (NUM_FIELD_MULT == something)
281  } // end of operator()()
282 
284  //
285  // evaluateFields()
286  //
288  template<typename EvalT, typename Traits>
289  void
292  typename Traits::EvalData workset)
293  {
294  using Kokkos::parallel_for;
295  using Kokkos::RangePolicy;
296 
297  // Grab the basis information.
298  basis_ = this->wda(workset).bases[basisIndex_]->weighted_grad_basis;
299 
300  // The following if-block is for the sake of optimization depending on the
301  // number of field multipliers. The parallel_fors will loop over the cells
302  // in the Workset and execute operator()() above.
303  if (fieldMults_.size() == 0)
304  parallel_for(RangePolicy<FieldMultTag<0>>(0, workset.num_cells), *this);
305  else if (fieldMults_.size() == 1)
306  parallel_for(RangePolicy<FieldMultTag<1>>(0, workset.num_cells), *this);
307  else
308  parallel_for(RangePolicy<FieldMultTag<-1>>(0, workset.num_cells), *this);
309  } // end of evaluateFields()
310 
312  //
313  // getValidParameters()
314  //
316  template<typename EvalT, typename TRAITS>
317  Teuchos::RCP<Teuchos::ParameterList>
319  getValidParameters() const
320  {
321  using panzer::BasisIRLayout;
323  using std::string;
324  using std::vector;
325  using Teuchos::ParameterList;
326  using Teuchos::RCP;
327  using Teuchos::rcp;
328 
329  // Create a ParameterList with all the valid keys we support.
330  RCP<ParameterList> p = rcp(new ParameterList);
331 
332  RCP<const vector<string>> resNames;
333  p->set("Residual Names", resNames);
334  p->set<string>("Scalar Name", "?");
335  RCP<BasisIRLayout> basis;
336  p->set("Basis", basis);
337  RCP<IntegrationRule> ir;
338  p->set("IR", ir);
339  p->set<double>("Multiplier", 1.0);
340  RCP<const vector<string>> fms;
341  p->set("Field Multipliers", fms);
342 
343  return p;
344  } // end of getValidParameters()
345 
346 } // end of namespace panzer
347 
348 #endif // PANZER_EVALUATOR_GRADBASISTIMESSCALAR_IMPL_HPP
View
panzer::EvaluatorStyle evalStyle
The EvaluatorStyle of the parent Integrator_CurlBasisDotVector object.
double multiplier
The scalar multiplier out in front of the integral ( ).
PHX::MDField< ScalarT, panzer::Cell, panzer::BASIS > field
A field to which we'll contribute, or in which we'll store, the result of computing this integral.
Teuchos::RCP< PHX::DataLayout > functional
<Cell,Basis>
Teuchos::RCP< const PureBasis > getBasis() const
int numDims_
The number of dimensions associated with the gradient.
PHX::MDField< const ScalarT, Cell, IP > scalar_
A field representing the scalar function we're integrating ( ).
Integrator_GradBasisTimesScalar(const panzer::EvaluatorStyle &evalStyle, const std::vector< std::string > &resNames, const std::string &scalar, const panzer::BasisIRLayout &basis, const panzer::IntegrationRule &ir, const double &multiplier=1, const std::vector< std::string > &fmNames=std::vector< std::string >())
Main Constructor.
KOKKOS_INLINE_FUNCTION void operator()(const FieldMultTag< NUM_FIELD_MULT > &tag, const std::size_t &cell) const
Perform the integration.
Kokkos::View< PHX::MDField< ScalarT, Cell, BASIS > * > fields_
The fields to which we'll contribute, or in which we'll store, the result of computing this integral.
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
Post-Registration Setup.
std::vector< PHX::MDField< const ScalarT, Cell, IP > > fieldMults_
The scalar multiplier out in front of the integral ( ).
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
Get Valid Parameters.
void evaluateFields(typename Traits::EvalData d)
Evaluate 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...
const panzer::EvaluatorStyle evalStyle_
An enum determining the behavior of this Evaluator.
Teuchos::RCP< PHX::DataLayout > dl_scalar
Data layout for scalar fields.
Description and data layouts associated with a particular basis.
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.
Kokkos::DynRankView< typename InputArray::value_type, PHX::Device > createDynRankView(const InputArray &a, const std::string &name, const DimensionPack... dims)
Wrapper to simplify Panzer use of Sacado ViewFactory.
EvaluatorStyle
An indication of how an Evaluator will behave.
This empty struct allows us to optimize operator()() depending on the number of field multipliers.
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_