Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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_.clear();
123  for (const auto& name : resNames)
124  {
125  MDField<ScalarT, Cell, BASIS> res(name, basis.functional);
126  fields_.push_back(res);
127  } // end loop over resNames
128  for (const auto& field : fields_)
130  this->addContributedField(field);
131  else // if (evalStyle_ == EvaluatorStyle::EVALUATES)
132  this->addEvaluatedField(field);
133 
134  // Add the dependent field multipliers, if there are any.
135  int i = 0;
136  fieldMults_.resize(fmNames.size());
138  typename DevLayout<ScalarT>::type, Device>*>(
139  "GradBasisTimesScalar::KokkosFieldMultipliers", fmNames.size());
140  for (const auto& name : fmNames)
141  {
142  fieldMults_[i++] = MDField<const ScalarT, Cell, IP>(name, ir.dl_scalar);
143  this->addDependentField(fieldMults_[i - 1]);
144  } // end loop over the field multipliers
145 
146  // Set the name of this object.
147  string n("Integrator_GradBasisTimesScalar (");
149  n += "CONTRIBUTES";
150  else // if (evalStyle_ == EvaluatorStyle::EVALUATES)
151  n += "EVALUATES";
152  n += "): {";
153  for (size_t j=0; j < fields_.size() - 1; ++j)
154  n += fields_[j].fieldTag().name() + ", ";
155  n += fields_.back().fieldTag().name() + "}";
156  this->setName(n);
157  } // end of Constructor
158 
160  //
161  // ParameterList Constructor
162  //
164  template<typename EvalT, typename Traits>
167  const Teuchos::ParameterList& p)
168  :
170  panzer::EvaluatorStyle::EVALUATES,
171  p.get<const std::vector<std::string>>("Residual Names"),
172  p.get<std::string>("Scalar Name"),
173  (*p.get<Teuchos::RCP<panzer::BasisIRLayout>>("Basis")),
174  (*p.get<Teuchos::RCP<panzer::IntegrationRule>>("IR")),
175  p.get<double>("Multiplier"),
176  p.isType<Teuchos::RCP<const std::vector<std::string>>>
177  ("Field Multipliers") ?
178  (*p.get<Teuchos::RCP<const std::vector<std::string>>>
179  ("Field Multipliers")) : std::vector<std::string>())
180  {
182  using Teuchos::RCP;
183 
184  // Ensure that the input ParameterList didn't contain any bogus entries.
185  RCP<ParameterList> validParams = this->getValidParameters();
186  p.validateParameters(*validParams);
187  } // end of ParameterList Constructor
188 
190  //
191  // postRegistrationSetup()
192  //
194  template<typename EvalT, typename Traits>
195  void
198  typename Traits::SetupData sd,
199  PHX::FieldManager<Traits>& /* fm */)
200  {
202  using panzer::getBasisIndex;
203  using PHX::Device;
204 
205  // Get the Kokkos::Views of the field multipliers.
206  for (size_t i(0); i < fieldMults_.size(); ++i)
207  kokkosFieldMults_(i) = fieldMults_[i].get_static_view();
208  Device::fence();
209 
210  // Determine the index in the Workset bases for our particular basis name.
211  basisIndex_ = getBasisIndex(basisName_, (*sd.worksets_)[0], this->wda);
212  } // end of postRegistrationSetup()
213 
214  template<typename EvalT, typename Traits>
215  template<int NUM_FIELD_MULT>
216  KOKKOS_INLINE_FUNCTION
217  void
220  const FieldMultTag<NUM_FIELD_MULT>& /* tag */,
221  const size_t& cell) const
222  {
224 
225  // Initialize the evaluated fields.
226  const int numQP(scalar_.extent(1)), numBases(fields_[0].extent(1));
227  if (evalStyle_ == EvaluatorStyle::EVALUATES)
228  for (int dim(0); dim < numDims_; ++dim)
229  for (int basis(0); basis < numBases; ++basis)
230  fields_[dim](cell, basis) = 0.0;
231 
232  // The following if-block is for the sake of optimization depending on the
233  // number of field multipliers.
234  ScalarT tmp;
235  if (NUM_FIELD_MULT == 0)
236  {
237  // Loop over the quadrature points, scale the integrand by the
238  // multiplier, and then perform the actual integration, looping over the
239  // bases and dimensions.
240  for (int qp(0); qp < numQP; ++qp)
241  {
242  tmp = multiplier_ * scalar_(cell, qp);
243  for (int basis(0); basis < numBases; ++basis)
244  for (int dim(0); dim < numDims_; ++dim)
245  fields_[dim](cell, basis) += basis_(cell, basis, qp, dim) * tmp;
246  } // end loop over the quadrature points
247  }
248  else if (NUM_FIELD_MULT == 1)
249  {
250  // Loop over the quadrature points, scale the integrand by the multiplier
251  // and the single field multiplier, and then perform the actual
252  // integration, looping over the bases and dimensions.
253  for (int qp(0); qp < numQP; ++qp)
254  {
255  tmp = multiplier_ * scalar_(cell, qp) * kokkosFieldMults_(0)(cell, qp);
256  for (int basis(0); basis < numBases; ++basis)
257  for (int dim(0); dim < numDims_; ++dim)
258  fields_[dim](cell, basis) += basis_(cell, basis, qp, dim) * tmp;
259  } // end loop over the quadrature points
260  }
261  else
262  {
263  // Loop over the quadrature points and pre-multiply all the field
264  // multipliers together, scale the integrand by the multiplier and
265  // the combination of the field multipliers, and then perform the actual
266  // integration, looping over the bases and dimensions.
267  const int numFieldMults(kokkosFieldMults_.extent(0));
268  for (int qp(0); qp < numQP; ++qp)
269  {
270  ScalarT fieldMultsTotal(1);
271  for (int fm(0); fm < numFieldMults; ++fm)
272  fieldMultsTotal *= kokkosFieldMults_(fm)(cell, qp);
273  tmp = multiplier_ * scalar_(cell, qp) * fieldMultsTotal;
274  for (int basis(0); basis < numBases; ++basis)
275  for (int dim(0); dim < numDims_; ++dim)
276  fields_[dim](cell, basis) += basis_(cell, basis, qp, dim) * tmp;
277  } // end loop over the quadrature points
278  } // end if (NUM_FIELD_MULT == something)
279  } // end of operator()()
280 
282  //
283  // evaluateFields()
284  //
286  template<typename EvalT, typename Traits>
287  void
290  typename Traits::EvalData workset)
291  {
292  using Kokkos::parallel_for;
293  using Kokkos::RangePolicy;
294 
295  // Grab the basis information.
296  basis_ = this->wda(workset).bases[basisIndex_]->weighted_grad_basis;
297 
298  // The following if-block is for the sake of optimization depending on the
299  // number of field multipliers. The parallel_fors will loop over the cells
300  // in the Workset and execute operator()() above.
301  if (fieldMults_.size() == 0)
302  parallel_for(RangePolicy<FieldMultTag<0>>(0, workset.num_cells), *this);
303  else if (fieldMults_.size() == 1)
304  parallel_for(RangePolicy<FieldMultTag<1>>(0, workset.num_cells), *this);
305  else
306  parallel_for(RangePolicy<FieldMultTag<-1>>(0, workset.num_cells), *this);
307  } // end of evaluateFields()
308 
310  //
311  // getValidParameters()
312  //
314  template<typename EvalT, typename TRAITS>
318  {
319  using panzer::BasisIRLayout;
321  using std::string;
322  using std::vector;
324  using Teuchos::RCP;
325  using Teuchos::rcp;
326 
327  // Create a ParameterList with all the valid keys we support.
328  RCP<ParameterList> p = rcp(new ParameterList);
329 
330  RCP<const vector<string>> resNames;
331  p->set("Residual Names", resNames);
332  p->set<string>("Scalar Name", "?");
333  RCP<BasisIRLayout> basis;
334  p->set("Basis", basis);
336  p->set("IR", ir);
337  p->set<double>("Multiplier", 1.0);
339  p->set("Field Multipliers", fms);
340 
341  return p;
342  } // end of getValidParameters()
343 
344 } // end of namespace panzer
345 
346 #endif // PANZER_EVALUATOR_GRADBASISTIMESSCALAR_IMPL_HPP
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...
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
Post-Registration Setup.
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.
PHX::MDField< const ScalarT, Cell, IP > scalar_
A field representing the scalar function we&#39;re integrating ( ).
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
EvaluatorStyle
An indication of how an Evaluator will behave.
Teuchos::RCP< const PureBasis > getBasis() const
const panzer::EvaluatorStyle evalStyle_
An enum determining the behavior of this Evaluator.
KOKKOS_INLINE_FUNCTION void operator()(const FieldMultTag< NUM_FIELD_MULT > &tag, const std::size_t &cell) const
Perform the integration.
double multiplier
The scalar multiplier out in front of the integral ( ).
std::vector< 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...
void evaluateFields(typename Traits::EvalData d)
Evaluate Fields.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
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.
void validateParameters(ParameterList const &validParamList, int const depth=1000, EValidateUsed const validateUsed=VALIDATE_USED_ENABLED, EValidateDefaults const validateDefaults=VALIDATE_DEFAULTS_ENABLED) const
const std::vector< std::pair< int, LocalOrdinal > > &pid_and_lid const
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...
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
&lt;Cell,Basis&gt;
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
Get Valid Parameters.
This empty struct allows us to optimize operator()() depending on the number of field multipliers...
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_