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 // Panzer: A partial differential equation assembly
4 // engine for strongly coupled complex multiphysics systems
5 //
6 // Copyright 2011 NTESS and the Panzer contributors.
7 // SPDX-License-Identifier: BSD-3-Clause
8 // *****************************************************************************
9 // @HEADER
10 
11 #ifndef PANZER_EVALUATOR_GRADBASISTIMESSCALAR_IMPL_HPP
12 #define PANZER_EVALUATOR_GRADBASISTIMESSCALAR_IMPL_HPP
13 
15 //
16 // Include Files
17 //
19 
20 // Intrepid2
21 #include "Intrepid2_FunctionSpaceTools.hpp"
22 
23 // Kokkos
24 #include "Kokkos_ViewFactory.hpp"
25 
26 // Panzer
27 #include "Panzer_BasisIRLayout.hpp"
30 
31 namespace panzer
32 {
34  //
35  // Constructor
36  //
38  template<typename EvalT, typename Traits>
42  const std::vector<std::string>& resNames,
43  const std::string& scalar,
44  const panzer::BasisIRLayout& basis,
45  const panzer::IntegrationRule& ir,
46  const double& multiplier, /* = 1 */
47  const std::vector<std::string>& fmNames) /* =
48  std::vector<std::string>() */
49  :
50  evalStyle_(evalStyle),
51  multiplier_(multiplier),
52  numDims_(ir.dl_vector->extent(2)),
53  basisName_(basis.name())
54  {
55  using PHX::View;
56  using panzer::BASIS;
57  using panzer::Cell;
59  using panzer::IP;
60  using panzer::PureBasis;
61  using PHX::Device;
62  using PHX::DevLayout;
63  using PHX::MDField;
64  using std::invalid_argument;
65  using std::logic_error;
66  using std::string;
67  using Teuchos::RCP;
68 
69  // Ensure the input makes sense.
70  for (const auto& name : resNames)
71  TEUCHOS_TEST_FOR_EXCEPTION(name == "", invalid_argument, "Error: " \
72  "Integrator_GradBasisTimesScalar called with an empty residual name.")
73  TEUCHOS_TEST_FOR_EXCEPTION(scalar == "", invalid_argument, "Error: " \
74  "Integrator_GradBasisTimesScalar called with an empty scalar name.")
75  TEUCHOS_TEST_FOR_EXCEPTION(numDims_ != static_cast<int>(resNames.size()),
76  logic_error, "Error: Integrator_GradBasisTimesScalar: The number of " \
77  "residuals must equal the number of gradient components (mesh " \
78  "dimensions).")
79  RCP<const PureBasis> tmpBasis = basis.getBasis();
80  TEUCHOS_TEST_FOR_EXCEPTION(not tmpBasis->supportsGrad(), logic_error,
81  "Error: Integrator_GradBasisTimesScalar: Basis of type \""
82  << tmpBasis->name() << "\" does not support the gradient operation.")
83 
84  // Create the field or the scalar function we're integrating.
85  scalar_ = MDField<const ScalarT, Cell, IP>(scalar, ir.dl_scalar);
86  this->addDependentField(scalar_);
87 
88  // Create the fields that we're either contributing to or evaluating
89  // (storing).
90  fields_host_.resize(resNames.size());
91  fields_ = OuterView("Integrator_GradBasisCrossVector::fields_", resNames.size());
92  {
93  int i=0;
94  for (const auto& name : resNames)
96  }
97  for (std::size_t i=0; i<fields_.extent(0); ++i) {
98  const auto& field = fields_host_[i];
100  this->addContributedField(field);
101  else // if (evalStyle_ == EvaluatorStyle::EVALUATES)
102  this->addEvaluatedField(field);
103  }
104 
105  // Add the dependent field multipliers, if there are any.
106  int i = 0;
107  fieldMults_.resize(fmNames.size());
108  kokkosFieldMults_ = PHX::View<PHX::UnmanagedView<const ScalarT**>*>(
109  "GradBasisTimesScalar::KokkosFieldMultipliers", fmNames.size());
110  for (const auto& name : fmNames)
111  {
112  fieldMults_[i++] = MDField<const ScalarT, Cell, IP>(name, ir.dl_scalar);
113  this->addDependentField(fieldMults_[i - 1]);
114  } // end loop over the field multipliers
115 
116  // Set the name of this object.
117  string n("Integrator_GradBasisTimesScalar (");
119  n += "CONTRIBUTES";
120  else // if (evalStyle_ == EvaluatorStyle::EVALUATES)
121  n += "EVALUATES";
122  n += "): {";
123  for (size_t j=0; j < fields_host_.size() - 1; ++j)
124  n += resNames[j] + ", ";
125  n += resNames[resNames.size()-1] + "}";
126  this->setName(n);
127  } // end of Constructor
128 
130  //
131  // ParameterList Constructor
132  //
134  template<typename EvalT, typename Traits>
137  const Teuchos::ParameterList& p)
138  :
140  panzer::EvaluatorStyle::EVALUATES,
141  p.get<const std::vector<std::string>>("Residual Names"),
142  p.get<std::string>("Scalar Name"),
143  (*p.get<Teuchos::RCP<panzer::BasisIRLayout>>("Basis")),
144  (*p.get<Teuchos::RCP<panzer::IntegrationRule>>("IR")),
145  p.get<double>("Multiplier"),
146  p.isType<Teuchos::RCP<const std::vector<std::string>>>
147  ("Field Multipliers") ?
148  (*p.get<Teuchos::RCP<const std::vector<std::string>>>
149  ("Field Multipliers")) : std::vector<std::string>())
150  {
152  using Teuchos::RCP;
153 
154  // Ensure that the input ParameterList didn't contain any bogus entries.
155  RCP<ParameterList> validParams = this->getValidParameters();
156  p.validateParameters(*validParams);
157  } // end of ParameterList Constructor
158 
160  //
161  // postRegistrationSetup()
162  //
164  template<typename EvalT, typename Traits>
165  void
168  typename Traits::SetupData sd,
169  PHX::FieldManager<Traits>& /* fm */)
170  {
172  using panzer::getBasisIndex;
173 
174  // Get the PHX::Views of the fields.
175  auto fields_host_mirror = Kokkos::create_mirror_view(fields_);
176  for (size_t i(0); i < fields_host_.size(); ++i)
177  fields_host_mirror(i) = fields_host_[i].get_static_view();
178  Kokkos::deep_copy(fields_,fields_host_mirror);
179 
180  // Get the PHX::Views of the field multipliers.
181  auto field_mults_host_mirror = Kokkos::create_mirror_view(kokkosFieldMults_);
182  for (size_t i(0); i < fieldMults_.size(); ++i)
183  field_mults_host_mirror(i) = fieldMults_[i].get_static_view();
184  Kokkos::deep_copy(kokkosFieldMults_,field_mults_host_mirror);
185 
186  // Determine the index in the Workset bases for our particular basis name.
187  basisIndex_ = getBasisIndex(basisName_, (*sd.worksets_)[0], this->wda);
188  } // end of postRegistrationSetup()
189 
190  template<typename EvalT, typename Traits>
191  template<int NUM_FIELD_MULT>
192  KOKKOS_INLINE_FUNCTION
193  void
196  const FieldMultTag<NUM_FIELD_MULT>& /* tag */,
197  const size_t& cell) const
198  {
200 
201  // Initialize the evaluated fields.
202  const int numQP(scalar_.extent(1)), numBases(fields_[0].extent(1));
203  if (evalStyle_ == EvaluatorStyle::EVALUATES)
204  for (int dim(0); dim < numDims_; ++dim)
205  for (int basis(0); basis < numBases; ++basis)
206  fields_[dim](cell, basis) = 0.0;
207 
208  // The following if-block is for the sake of optimization depending on the
209  // number of field multipliers.
210  ScalarT tmp;
211  if (NUM_FIELD_MULT == 0)
212  {
213  // Loop over the quadrature points, scale the integrand by the
214  // multiplier, and then perform the actual integration, looping over the
215  // bases and dimensions.
216  for (int qp(0); qp < numQP; ++qp)
217  {
218  tmp = multiplier_ * scalar_(cell, qp);
219  for (int basis(0); basis < numBases; ++basis)
220  for (int dim(0); dim < numDims_; ++dim)
221  fields_[dim](cell, basis) += basis_(cell, basis, qp, dim) * tmp;
222  } // end loop over the quadrature points
223  }
224  else if (NUM_FIELD_MULT == 1)
225  {
226  // Loop over the quadrature points, scale the integrand by the multiplier
227  // and the single field multiplier, and then perform the actual
228  // integration, looping over the bases and dimensions.
229  for (int qp(0); qp < numQP; ++qp)
230  {
231  tmp = multiplier_ * scalar_(cell, qp) * kokkosFieldMults_(0)(cell, qp);
232  for (int basis(0); basis < numBases; ++basis)
233  for (int dim(0); dim < numDims_; ++dim)
234  fields_[dim](cell, basis) += basis_(cell, basis, qp, dim) * tmp;
235  } // end loop over the quadrature points
236  }
237  else
238  {
239  // Loop over the quadrature points and pre-multiply all the field
240  // multipliers together, scale the integrand by the multiplier and
241  // the combination of the field multipliers, and then perform the actual
242  // integration, looping over the bases and dimensions.
243  const int numFieldMults(kokkosFieldMults_.extent(0));
244  for (int qp(0); qp < numQP; ++qp)
245  {
246  ScalarT fieldMultsTotal(1);
247  for (int fm(0); fm < numFieldMults; ++fm)
248  fieldMultsTotal *= kokkosFieldMults_(fm)(cell, qp);
249  tmp = multiplier_ * scalar_(cell, qp) * fieldMultsTotal;
250  for (int basis(0); basis < numBases; ++basis)
251  for (int dim(0); dim < numDims_; ++dim)
252  fields_[dim](cell, basis) += basis_(cell, basis, qp, dim) * 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_grad_basis;
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 
306  RCP<const vector<string>> resNames;
307  p->set("Residual Names", resNames);
308  p->set<string>("Scalar Name", "?");
309  RCP<BasisIRLayout> basis;
310  p->set("Basis", basis);
312  p->set("IR", ir);
313  p->set<double>("Multiplier", 1.0);
315  p->set("Field Multipliers", fms);
316 
317  return p;
318  } // end of getValidParameters()
319 
320 } // end of namespace panzer
321 
322 #endif // PANZER_EVALUATOR_GRADBASISTIMESSCALAR_IMPL_HPP
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 ( ).
int num_cells
DEPRECATED - use: numCells()
#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.
PHX::View< PHX::UnmanagedView< const ScalarT ** > * > kokkosFieldMults_
The PHX::View representation of the (possibly empty) list of fields that are multipliers out in front...
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 ( ).
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.
std::vector< PHX::MDField< ScalarT, Cell, BASIS > > fields_host_
The fields to which we&#39;ll contribute, or in which we&#39;ll store, the result of computing this integral...
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
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...
const std::vector< std::pair< int, LocalOrdinal > > &pid_and_lid const
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_