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 // 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_INTEGRATOR_BASISTIMESSCALAR_IMPL_HPP
12 #define PANZER_INTEGRATOR_BASISTIMESSCALAR_IMPL_HPP
13 
15 //
16 // Include Files
17 //
19 
20 // Panzer
21 #include "Panzer_BasisIRLayout.hpp"
24 
25 namespace panzer
26 {
28  //
29  // Main Constructor
30  //
32  template<typename EvalT, typename Traits>
36  const std::string& resName,
37  const std::string& valName,
38  const panzer::BasisIRLayout& basis,
39  const panzer::IntegrationRule& ir,
40  const double& multiplier, /* = 1 */
41  const std::vector<std::string>& fmNames /* =
42  std::vector<std::string>() */)
43  :
44  evalStyle_(evalStyle),
45  multiplier_(multiplier),
46  basisName_(basis.name())
47  {
48  using PHX::View;
49  using panzer::BASIS;
50  using panzer::Cell;
52  using panzer::IP;
53  using panzer::PureBasis;
54  using PHX::MDField;
55  using std::invalid_argument;
56  using std::logic_error;
57  using std::string;
58  using Teuchos::RCP;
59 
60  // Ensure the input makes sense.
61  TEUCHOS_TEST_FOR_EXCEPTION(resName == "", invalid_argument, "Error: " \
62  "Integrator_BasisTimesScalar called with an empty residual name.")
63  TEUCHOS_TEST_FOR_EXCEPTION(valName == "", invalid_argument, "Error: " \
64  "Integrator_BasisTimesScalar called with an empty value name.")
65  RCP<const PureBasis> tmpBasis = basis.getBasis();
66  TEUCHOS_TEST_FOR_EXCEPTION(not tmpBasis->isScalarBasis(), logic_error,
67  "Error: Integrator_BasisTimesScalar: Basis of type \""
68  << tmpBasis->name() << "\" is not a scalar basis.")
69 
70  // Create the field for the scalar quantity we're integrating.
71  scalar_ = MDField<const ScalarT, Cell, IP>(valName, ir.dl_scalar);
72  this->addDependentField(scalar_);
73 
74  // Create the field that we're either contributing to or evaluating
75  // (storing).
76  field_ = MDField<ScalarT, Cell, BASIS>(resName, basis.functional);
77  if (evalStyle == EvaluatorStyle::CONTRIBUTES)
78  this->addContributedField(field_);
79  else // if (evalStyle == EvaluatorStyle::EVALUATES)
80  this->addEvaluatedField(field_);
81 
82  // Add the dependent field multipliers, if there are any.
83  int i(0);
84  fieldMults_.resize(fmNames.size());
86  View<Kokkos::View<const ScalarT**, typename PHX::DevLayout<ScalarT>::type, Kokkos::MemoryUnmanaged>*>("BasisTimesScalar::KokkosFieldMultipliers",
87  fmNames.size());
88  for (const auto& name : fmNames)
89  {
90  fieldMults_[i++] = MDField<const ScalarT, Cell, IP>(name, ir.dl_scalar);
91  this->addDependentField(fieldMults_[i - 1]);
92  } // end loop over the field multipliers
93 
94  // Set the name of this object.
95  string n("Integrator_BasisTimesScalar (");
97  n += "CONTRIBUTES";
98  else // if (evalStyle_ == EvaluatorStyle::EVALUATES)
99  n += "EVALUATES";
100  n += "): " + field_.fieldTag().name();
101  this->setName(n);
102  } // end of Main Constructor
103 
105  //
106  // ParameterList Constructor
107  //
109  template<typename EvalT, typename Traits>
112  const panzer::EvaluatorStyle es)
113  :
115  p.get<std::string>("Residual Name"),
116  p.get<std::string>("Value Name"),
117  (*p.get<Teuchos::RCP<panzer::BasisIRLayout>>("Basis")),
118  (*p.get<Teuchos::RCP<panzer::IntegrationRule>>("IR")),
119  p.get<double>("Multiplier"),
120  p.isType<Teuchos::RCP<const std::vector<std::string>>>
121  ("Field Multipliers") ?
122  (*p.get<Teuchos::RCP<const std::vector<std::string>>>
123  ("Field Multipliers")) : std::vector<std::string>())
124  {
126  using Teuchos::RCP;
127 
128  // Ensure that the input ParameterList didn't contain any bogus entries.
129  RCP<ParameterList> validParams = this->getValidParameters();
130  p.validateParameters(*validParams);
131  } // end of ParameterList Constructor
132 
134  //
135  // postRegistrationSetup()
136  //
138  template<typename EvalT, typename Traits>
139  void
142  typename Traits::SetupData sd,
143  PHX::FieldManager<Traits>& /* fm */)
144  {
145  using panzer::getBasisIndex;
146  using std::size_t;
147 
148  auto kokkosFieldMults_h = Kokkos::create_mirror_view(kokkosFieldMults_);
149 
150  // Get the PHX::Views of the field multipliers.
151  for (size_t i(0); i < fieldMults_.size(); ++i)
152  kokkosFieldMults_h(i) = fieldMults_[i].get_static_view();
153 
154  Kokkos::deep_copy(kokkosFieldMults_, kokkosFieldMults_h);
155 
156  // Determine the index in the Workset bases for our particular basis name.
157  basisIndex_ = getBasisIndex(basisName_, (*sd.worksets_)[0], this->wda);
158 
159  // Allocate temporary
160  if (Sacado::IsADType<ScalarT>::value) {
161  const auto fadSize = Kokkos::dimension_scalar(field_.get_view());
162  tmp_ = PHX::View<ScalarT*>("IntegratorBasisTimesScalar::tmp_",field_.extent(0),fadSize);
163  if (fieldMults_.size() > 1)
164  tmp2_ = PHX::View<ScalarT*>("IntegratorBasisTimesScalar::tmp_",field_.extent(0),fadSize);
165  } else {
166  tmp_ = PHX::View<ScalarT*>("IntegratorBasisTimesScalar::tmp_",field_.extent(0));
167  if (fieldMults_.size() > 1)
168  tmp2_ = PHX::View<ScalarT*>("IntegratorBasisTimesScalar::tmp_",field_.extent(0));
169  }
170  } // end of postRegistrationSetup()
171 
173  //
174  // operator()()
175  //
177  template<typename EvalT, typename Traits>
178  template<int NUM_FIELD_MULT>
179  KOKKOS_INLINE_FUNCTION
180  void
183  const FieldMultTag<NUM_FIELD_MULT>& /* tag */,
184  const size_t& cell) const
185  {
187 
188  // Initialize the evaluated field.
189  const int numQP(scalar_.extent(1)), numBases(basis_.extent(1));
190  if (evalStyle_ == EvaluatorStyle::EVALUATES)
191  for (int basis(0); basis < numBases; ++basis)
192  field_(cell, basis) = 0.0;
193 
194  // The following if-block is for the sake of optimization depending on the
195  // number of field multipliers.
196  if (NUM_FIELD_MULT == 0)
197  {
198  // Loop over the quadrature points, scale the integrand by the
199  // multiplier, and then perform the actual integration, looping over the
200  // bases.
201  for (int qp(0); qp < numQP; ++qp)
202  {
203  tmp_(cell) = multiplier_ * scalar_(cell, qp);
204  for (int basis(0); basis < numBases; ++basis)
205  field_(cell, basis) += basis_(cell, basis, qp) * tmp_(cell);
206  } // end loop over the quadrature points
207  }
208  else if (NUM_FIELD_MULT == 1)
209  {
210  // Loop over the quadrature points, scale the integrand by the multiplier
211  // and the single field multiplier, and then perform the actual
212  // integration, looping over the bases.
213  for (int qp(0); qp < numQP; ++qp)
214  {
215  tmp_(cell) = multiplier_ * scalar_(cell, qp) * kokkosFieldMults_(0)(cell, qp);
216  for (int basis(0); basis < numBases; ++basis)
217  field_(cell, basis) += basis_(cell, basis, qp) * tmp_(cell);
218  } // end loop over the quadrature points
219  }
220  else
221  {
222  // Loop over the quadrature points and pre-multiply all the field
223  // multipliers together, scale the integrand by the multiplier and
224  // the combination of the field multipliers, and then perform the actual
225  // integration, looping over the bases.
226  const int numFieldMults(kokkosFieldMults_.extent(0));
227  for (int qp(0); qp < numQP; ++qp)
228  {
229  tmp2_(cell) = 1.0;
230  for (int fm(0); fm < numFieldMults; ++fm)
231  tmp2_(cell) *= kokkosFieldMults_(fm)(cell, qp);
232  tmp_(cell) = multiplier_ * scalar_(cell, qp) * tmp2_(cell);
233  for (int basis(0); basis < numBases; ++basis)
234  field_(cell, basis) += basis_(cell, basis, qp) * tmp_(cell);
235  } // end loop over the quadrature points
236  } // end if (NUM_FIELD_MULT == something)
237  } // end of operator()()
238 
240  //
241  // evaluateFields()
242  //
244  template<typename EvalT, typename Traits>
245  void
248  typename Traits::EvalData workset)
249  {
250  using Kokkos::parallel_for;
251  using Kokkos::RangePolicy;
252 
253  // Grab the basis information.
254  basis_ = this->wda(workset).bases[basisIndex_]->weighted_basis_scalar;
255 
256  // The following if-block is for the sake of optimization depending on the
257  // number of field multipliers. The parallel_fors will loop over the cells
258  // in the Workset and execute operator()() above.
259  if (fieldMults_.size() == 0)
260  parallel_for(RangePolicy<FieldMultTag<0>>(0, workset.num_cells), *this);
261  else if (fieldMults_.size() == 1)
262  parallel_for(RangePolicy<FieldMultTag<1>>(0, workset.num_cells), *this);
263  else
264  parallel_for(RangePolicy<FieldMultTag<-1>>(0, workset.num_cells), *this);
265  } // end of evaluateFields()
266 
268  //
269  // getValidParameters()
270  //
272  template<typename EvalT, typename TRAITS>
276  {
277  using panzer::BasisIRLayout;
279  using std::string;
280  using std::vector;
282  using Teuchos::RCP;
283  using Teuchos::rcp;
284 
285  // Create a ParameterList with all the valid keys we support.
286  RCP<ParameterList> p = rcp(new ParameterList);
287  p->set<string>("Residual Name", "?");
288  p->set<string>("Value Name", "?");
289  RCP<BasisIRLayout> basis;
290  p->set("Basis", basis);
292  p->set("IR", ir);
293  p->set<double>("Multiplier", 1.0);
295  p->set("Field Multipliers", fms);
296  return p;
297  } // end of getValidParameters()
298 
299 } // end of namespace panzer
300 
301 #endif // PANZER_INTEGRATOR_BASISTIMESSCALAR_IMPL_HPP
int num_cells
DEPRECATED - use: numCells()
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.
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 ( ).
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
PHX::View< Kokkos::View< const ScalarT **, typename PHX::DevLayout< ScalarT >::type, Kokkos::MemoryUnmanaged > * > kokkosFieldMults_
The PHX::View representation of the (possibly empty) list of fields that are multipliers out in front...
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_