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 PHX::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_ = PHX::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());
139  kokkosFieldMults_ = View<View<const ScalarT**>*>(
140  "GradBasisTimesScalar::KokkosFieldMultipliers", fmNames.size());
141  for (const auto& name : fmNames)
142  {
143  fieldMults_[i++] = MDField<const ScalarT, Cell, IP>(name, ir.dl_scalar);
144  this->addDependentField(fieldMults_[i - 1]);
145  } // end loop over the field multipliers
146 
147  // Set the name of this object.
148  string n("Integrator_GradBasisTimesScalar (");
150  n += "CONTRIBUTES";
151  else // if (evalStyle_ == EvaluatorStyle::EVALUATES)
152  n += "EVALUATES";
153  n += "): {";
154  for (size_t j=0; j < fields_.extent(0) - 1; ++j)
155  n += fields_[j].fieldTag().name() + ", ";
156  n += fields_(fields_.extent(0)-1).fieldTag().name() + "}";
157  this->setName(n);
158  } // end of Constructor
159 
161  //
162  // ParameterList Constructor
163  //
165  template<typename EvalT, typename Traits>
168  const Teuchos::ParameterList& p)
169  :
172  p.get<const std::vector<std::string>>("Residual Names"),
173  p.get<std::string>("Scalar Name"),
174  (*p.get<Teuchos::RCP<panzer::BasisIRLayout>>("Basis")),
175  (*p.get<Teuchos::RCP<panzer::IntegrationRule>>("IR")),
176  p.get<double>("Multiplier"),
177  p.isType<Teuchos::RCP<const std::vector<std::string>>>
178  ("Field Multipliers") ?
179  (*p.get<Teuchos::RCP<const std::vector<std::string>>>
180  ("Field Multipliers")) : std::vector<std::string>())
181  {
182  using Teuchos::ParameterList;
183  using Teuchos::RCP;
184 
185  // Ensure that the input ParameterList didn't contain any bogus entries.
186  RCP<ParameterList> validParams = this->getValidParameters();
187  p.validateParameters(*validParams);
188  } // end of ParameterList Constructor
189 
191  //
192  // postRegistrationSetup()
193  //
195  template<typename EvalT, typename Traits>
196  void
199  typename Traits::SetupData sd,
200  PHX::FieldManager<Traits>& /* fm */)
201  {
203  using panzer::getBasisIndex;
204 
205  // Get the PHX::Views of the field multipliers.
206  for (size_t i(0); i < fieldMults_.size(); ++i)
207  kokkosFieldMults_(i) = fieldMults_[i].get_static_view();
208 
209  // Determine the index in the Workset bases for our particular basis name.
210  basisIndex_ = getBasisIndex(basisName_, (*sd.worksets_)[0], this->wda);
211  } // end of postRegistrationSetup()
212 
213  template<typename EvalT, typename Traits>
214  template<int NUM_FIELD_MULT>
215  KOKKOS_INLINE_FUNCTION
216  void
219  const FieldMultTag<NUM_FIELD_MULT>& /* tag */,
220  const size_t& cell) const
221  {
223 
224  // Initialize the evaluated fields.
225  const int numQP(scalar_.extent(1)), numBases(fields_[0].extent(1));
226  if (evalStyle_ == EvaluatorStyle::EVALUATES)
227  for (int dim(0); dim < numDims_; ++dim)
228  for (int basis(0); basis < numBases; ++basis)
229  fields_[dim](cell, basis) = 0.0;
230 
231  // The following if-block is for the sake of optimization depending on the
232  // number of field multipliers.
233  ScalarT tmp;
234  if (NUM_FIELD_MULT == 0)
235  {
236  // Loop over the quadrature points, scale the integrand by the
237  // multiplier, and then perform the actual integration, looping over the
238  // bases and dimensions.
239  for (int qp(0); qp < numQP; ++qp)
240  {
241  tmp = multiplier_ * scalar_(cell, qp);
242  for (int basis(0); basis < numBases; ++basis)
243  for (int dim(0); dim < numDims_; ++dim)
244  fields_[dim](cell, basis) += basis_(cell, basis, qp, dim) * tmp;
245  } // end loop over the quadrature points
246  }
247  else if (NUM_FIELD_MULT == 1)
248  {
249  // Loop over the quadrature points, scale the integrand by the multiplier
250  // and the single field multiplier, and then perform the actual
251  // integration, looping over the bases and dimensions.
252  for (int qp(0); qp < numQP; ++qp)
253  {
254  tmp = multiplier_ * scalar_(cell, qp) * kokkosFieldMults_(0)(cell, qp);
255  for (int basis(0); basis < numBases; ++basis)
256  for (int dim(0); dim < numDims_; ++dim)
257  fields_[dim](cell, basis) += basis_(cell, basis, qp, dim) * tmp;
258  } // end loop over the quadrature points
259  }
260  else
261  {
262  // Loop over the quadrature points and pre-multiply all the field
263  // multipliers together, scale the integrand by the multiplier and
264  // the combination of the field multipliers, and then perform the actual
265  // integration, looping over the bases and dimensions.
266  const int numFieldMults(kokkosFieldMults_.extent(0));
267  for (int qp(0); qp < numQP; ++qp)
268  {
269  ScalarT fieldMultsTotal(1);
270  for (int fm(0); fm < numFieldMults; ++fm)
271  fieldMultsTotal *= kokkosFieldMults_(fm)(cell, qp);
272  tmp = multiplier_ * scalar_(cell, qp) * fieldMultsTotal;
273  for (int basis(0); basis < numBases; ++basis)
274  for (int dim(0); dim < numDims_; ++dim)
275  fields_[dim](cell, basis) += basis_(cell, basis, qp, dim) * tmp;
276  } // end loop over the quadrature points
277  } // end if (NUM_FIELD_MULT == something)
278  } // end of operator()()
279 
281  //
282  // evaluateFields()
283  //
285  template<typename EvalT, typename Traits>
286  void
289  typename Traits::EvalData workset)
290  {
291  using Kokkos::parallel_for;
292  using Kokkos::RangePolicy;
293 
294  // Grab the basis information.
295  basis_ = this->wda(workset).bases[basisIndex_]->weighted_grad_basis;
296 
297  // The following if-block is for the sake of optimization depending on the
298  // number of field multipliers. The parallel_fors will loop over the cells
299  // in the Workset and execute operator()() above.
300  if (fieldMults_.size() == 0)
301  parallel_for(RangePolicy<FieldMultTag<0>>(0, workset.num_cells), *this);
302  else if (fieldMults_.size() == 1)
303  parallel_for(RangePolicy<FieldMultTag<1>>(0, workset.num_cells), *this);
304  else
305  parallel_for(RangePolicy<FieldMultTag<-1>>(0, workset.num_cells), *this);
306  } // end of evaluateFields()
307 
309  //
310  // getValidParameters()
311  //
313  template<typename EvalT, typename TRAITS>
314  Teuchos::RCP<Teuchos::ParameterList>
317  {
318  using panzer::BasisIRLayout;
320  using std::string;
321  using std::vector;
322  using Teuchos::ParameterList;
323  using Teuchos::RCP;
324  using Teuchos::rcp;
325 
326  // Create a ParameterList with all the valid keys we support.
327  RCP<ParameterList> p = rcp(new ParameterList);
328 
329  RCP<const vector<string>> resNames;
330  p->set("Residual Names", resNames);
331  p->set<string>("Scalar Name", "?");
332  RCP<BasisIRLayout> basis;
333  p->set("Basis", basis);
334  RCP<IntegrationRule> ir;
335  p->set("IR", ir);
336  p->set<double>("Multiplier", 1.0);
337  RCP<const vector<string>> fms;
338  p->set("Field Multipliers", fms);
339 
340  return p;
341  } // end of getValidParameters()
342 
343 } // end of namespace panzer
344 
345 #endif // PANZER_EVALUATOR_GRADBASISTIMESSCALAR_IMPL_HPP
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
Post-Registration Setup.
PHX::MDField< const ScalarT, Cell, IP > scalar_
A field representing the scalar function we&#39;re integrating ( ).
PHX::View< PHX::MDField< ScalarT, Cell, BASIS > * > fields_
The fields to which we&#39;ll contribute, or in which we&#39;ll store, the result of computing this integral...
int num_cells
DEPRECATED - use: numCells()
PHX::View< PHX::View< const ScalarT ** > * > kokkosFieldMults_
The PHX::View representation of the (possibly empty) list of fields that are multipliers out in front...
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
Get Valid Parameters.
Teuchos::RCP< const PureBasis > getBasis() const
EvaluatorStyle
An indication of how an Evaluator will behave.
const panzer::EvaluatorStyle evalStyle_
An enum determining the behavior of this Evaluator.
double multiplier
The scalar multiplier out in front of the integral ( ).
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.
void evaluateFields(typename Traits::EvalData d)
Evaluate Fields.
Teuchos::RCP< PHX::DataLayout > dl_scalar
Data layout for scalar fields.
std::vector< PHX::MDField< const ScalarT, Cell, IP > > fieldMults_
The scalar multiplier out in front of the integral ( ).
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.
int numDims_
The number of dimensions associated with the gradient.
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...
KOKKOS_INLINE_FUNCTION void operator()(const FieldMultTag< NUM_FIELD_MULT > &tag, const std::size_t &cell) const
Perform the integration.
Description and data layouts associated with a particular basis.
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.
Teuchos::RCP< PHX::DataLayout > functional
<Cell,Basis>
This empty struct allows us to optimize operator()() depending on the number of field multipliers...
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_