Panzer  Version of the Day
Panzer_Integrator_DivBasisTimesScalar_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_Integrator_DivBasisTimesScalar_impl_hpp__
44 #define __Panzer_Integrator_DivBasisTimesScalar_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"
63 
64 namespace panzer
65 {
67  //
68  // Main Constructor
69  //
71  template<typename EvalT, typename Traits>
75  const std::string& resName,
76  const std::string& valName,
77  const panzer::BasisIRLayout& basis,
78  const panzer::IntegrationRule& ir,
79  const double& multiplier, /* = 1 */
80  const std::vector<std::string>& fmNames /* =
81  std::vector<std::string>() */)
82  :
83  evalStyle_(evalStyle),
84  multiplier_(multiplier),
85  basisName_(basis.name()),
86  use_shared_memory(false)
87  {
88  using PHX::View;
89  using panzer::BASIS;
90  using panzer::Cell;
92  using panzer::IP;
93  using panzer::PureBasis;
94  using PHX::MDField;
95  using std::invalid_argument;
96  using std::logic_error;
97  using std::string;
98  using Teuchos::RCP;
99 
100  // Ensure the input makes sense.
101  TEUCHOS_TEST_FOR_EXCEPTION(resName == "", invalid_argument, "Error: " \
102  "Integrator_DivBasisTimesScalar called with an empty residual name.")
103  TEUCHOS_TEST_FOR_EXCEPTION(valName == "", invalid_argument, "Error: " \
104  "Integrator_DivBasisTimesScalar called with an empty value name.")
105  RCP<const PureBasis> tmpBasis = basis.getBasis();
106  TEUCHOS_TEST_FOR_EXCEPTION(not tmpBasis->supportsDiv(), logic_error,
107  "Error: Integrator_DivBasisTimesScalar: Basis of type \""
108  << tmpBasis->name() << "\" does not support DIV.")
109  TEUCHOS_TEST_FOR_EXCEPTION(not tmpBasis->requiresOrientations(), logic_error,
110  "Error: Integration_DivBasisTimesScalar: Basis of type \""
111  << tmpBasis->name() << "\" should require orientations.")
112 
113  // Create the field for the scalar quantity we're integrating.
114  scalar_ = MDField<const ScalarT, Cell, IP>(valName, ir.dl_scalar);
115  this->addDependentField(scalar_);
116 
117  // Create the field that we're either contributing to or evaluating
118  // (storing).
119  field_ = MDField<ScalarT, Cell, BASIS>(resName, basis.functional);
121  this->addContributedField(field_);
122  else // if (evalStyle_ == EvaluatorStyle::EVALUATES)
123  this->addEvaluatedField(field_);
124 
125  // Add the dependent field multipliers, if there are any.
126  int i(0);
127  fieldMults_.resize(fmNames.size());
129  View<View<const ScalarT**>*>("BasisTimesScalar::KokkosFieldMultipliers",
130  fmNames.size());
131  for (const auto& name : fmNames)
132  {
133  fieldMults_[i++] = MDField<const ScalarT, Cell, IP>(name, ir.dl_scalar);
134  this->addDependentField(fieldMults_[i - 1]);
135  } // end loop over the field multipliers
136 
137  // Set the name of this object.
138  string n("Integrator_DivBasisTimesScalar (");
140  n += "CONTRIBUTES";
141  else // if (evalStyle_ == EvaluatorStyle::EVALUATES)
142  n += "EVALUATES";
143  n += "): " + field_.fieldTag().name();
144  this->setName(n);
145  } // end of Main Constructor
146 
148  //
149  // ParameterList Constructor
150  //
152  template<typename EvalT, typename Traits>
155  const Teuchos::ParameterList& p)
156  :
159  p.get<std::string>("Residual Name"),
160  p.get<std::string>("Value Name"),
161  (*p.get<Teuchos::RCP<panzer::BasisIRLayout>>("Basis")),
162  (*p.get<Teuchos::RCP<panzer::IntegrationRule>>("IR")),
163  p.get<double>("Multiplier"),
164  p.isType<Teuchos::RCP<const std::vector<std::string>>>
165  ("Field Multipliers") ?
166  (*p.get<Teuchos::RCP<const std::vector<std::string>>>
167  ("Field Multipliers")) : std::vector<std::string>())
168  {
169  using Teuchos::ParameterList;
170  using Teuchos::RCP;
171 
172  // Ensure that the input ParameterList didn't contain any bogus entries.
173  RCP<ParameterList> validParams = this->getValidParameters();
174  p.validateParameters(*validParams);
175  } // end of ParameterList Constructor
176 
178  //
179  // postRegistrationSetup()
180  //
182  template<typename EvalT, typename Traits>
183  void
186  typename Traits::SetupData sd,
187  PHX::FieldManager<Traits>& /* fm */)
188  {
190  using panzer::getBasisIndex;
191 
192  // Get the PHX::Views of the field multipliers.
193  for (size_t i(0); i < fieldMults_.size(); ++i)
194  kokkosFieldMults_(i) = fieldMults_[i].get_static_view();
195 
196  // Determine the index in the Workset bases for our particular basis name.
197  basisIndex_ = getBasisIndex(basisName_, (*sd.worksets_)[0], this->wda);
198  } // end of postRegistrationSetup()
199 
201  //
202  // No shared memory operator()()
203  //
205  template<typename EvalT, typename Traits>
206  template<int NUM_FIELD_MULT>
207  KOKKOS_INLINE_FUNCTION
208  void
211  const FieldMultTag<NUM_FIELD_MULT>& /* tag */,
212  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
213  {
215  const int cell = team.league_rank();
216 
217  // Initialize the evaluated field.
218  const int numQP(scalar_.extent(1)), numBases(basis_.extent(1));
219  if (evalStyle_ == EvaluatorStyle::EVALUATES) {
220  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (const int& basis) {
221  field_(cell, basis) = 0.0;
222  });
223  }
224 
225  // The following if-block is for the sake of optimization depending on the
226  // number of field multipliers.
227  ScalarT tmp;
228  if (NUM_FIELD_MULT == 0)
229  {
230  // Loop over the quadrature points, scale the integrand by the
231  // multiplier, and then perform the actual integration, looping over the
232  // bases.
233  for (int qp(0); qp < numQP; ++qp)
234  {
235  tmp = multiplier_ * scalar_(cell, qp);
236  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (const int& basis) {
237  field_(cell, basis) += basis_(cell, basis, qp) * tmp;
238  });
239  } // end loop over the quadrature points
240  }
241  else if (NUM_FIELD_MULT == 1)
242  {
243  // Loop over the quadrature points, scale the integrand by the multiplier
244  // and the single field multiplier, and then perform the actual
245  // integration, looping over the bases.
246  for (int qp(0); qp < numQP; ++qp)
247  {
248  tmp = multiplier_ * scalar_(cell, qp) * kokkosFieldMults_(0)(cell, qp);
249  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (const int& basis) {
250  field_(cell, basis) += basis_(cell, basis, qp) * tmp;
251  });
252  } // end loop over the quadrature points
253  }
254  else
255  {
256  // Loop over the quadrature points and pre-multiply all the field
257  // multipliers together, scale the integrand by the multiplier and the
258  // combination of field multipliers, and then perform the actual
259  // integration, looping over the bases.
260  const int numFieldMults(kokkosFieldMults_.extent(0));
261  for (int qp(0); qp < numQP; ++qp)
262  {
263  ScalarT fieldMultsTotal(1);
264  for (int fm(0); fm < numFieldMults; ++fm)
265  fieldMultsTotal *= kokkosFieldMults_(fm)(cell, qp);
266  tmp = multiplier_ * scalar_(cell, qp) * fieldMultsTotal;
267  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (const int& basis) {
268  field_(cell, basis) += basis_(cell, basis, qp) * tmp;
269  });
270  } // end loop over the quadrature points
271  } // end if (NUM_FIELD_MULT == something)
272  } // end of operator()
273 
275  //
276  // Shared memory operator()()
277  //
279  template<typename EvalT, typename Traits>
280  template<int NUM_FIELD_MULT>
281  KOKKOS_INLINE_FUNCTION
282  void
285  const SharedFieldMultTag<NUM_FIELD_MULT>& /* tag */,
286  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
287  {
289  const int cell = team.league_rank();
290  const int numQP = scalar_.extent(1);
291  const int numBases = basis_.extent(1);
292  const int fadSize = Kokkos::dimension_scalar(field_.get_view());
293 
294  scratch_view tmp;
295  scratch_view tmp_field;
296  if (Sacado::IsADType<ScalarT>::value) {
297  tmp = scratch_view(team.team_shmem(),1,fadSize);
298  tmp_field = scratch_view(team.team_shmem(),numBases,fadSize);
299  }
300  else {
301  tmp = scratch_view(team.team_shmem(),1);
302  tmp_field = scratch_view(team.team_shmem(),numBases);
303  }
304 
305  // Initialize the evaluated field.
306  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (const int& basis) {
307  tmp_field(basis) = 0.0;
308  });
309 
310  // The following if-block is for the sake of optimization depending on the
311  // number of field multipliers.
312  if (NUM_FIELD_MULT == 0)
313  {
314  // Loop over the quadrature points, scale the integrand by the
315  // multiplier, and then perform the actual integration, looping over the
316  // bases.
317  for (int qp(0); qp < numQP; ++qp)
318  {
319  team.team_barrier();
320  tmp(0) = multiplier_ * scalar_(cell, qp);
321  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (const int& basis) {
322  tmp_field(basis) += basis_(cell, basis, qp) * tmp(0);
323  });
324  } // end loop over the quadrature points
325  }
326  else if (NUM_FIELD_MULT == 1)
327  {
328  // Loop over the quadrature points, scale the integrand by the multiplier
329  // and the single field multiplier, and then perform the actual
330  // integration, looping over the bases.
331  for (int qp(0); qp < numQP; ++qp)
332  {
333  team.team_barrier();
334  tmp(0) = multiplier_ * scalar_(cell, qp) * kokkosFieldMults_(0)(cell, qp);
335  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (const int& basis) {
336  tmp_field(basis) += basis_(cell, basis, qp) * tmp(0);
337  });
338  } // end loop over the quadrature points
339  }
340  else
341  {
342  // Loop over the quadrature points and pre-multiply all the field
343  // multipliers together, scale the integrand by the multiplier and the
344  // combination of field multipliers, and then perform the actual
345  // integration, looping over the bases.
346  const int numFieldMults = kokkosFieldMults_.extent(0);
347  for (int qp(0); qp < numQP; ++qp)
348  {
349  team.team_barrier();
350  ScalarT fieldMultsTotal(1); // need shared mem here
351  for (int fm(0); fm < numFieldMults; ++fm)
352  fieldMultsTotal *= kokkosFieldMults_(fm)(cell, qp);
353  tmp(0) = multiplier_ * scalar_(cell, qp) * fieldMultsTotal;
354  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (const int& basis) {
355  tmp_field(basis) += basis_(cell, basis, qp) * tmp(0);
356  });
357  } // end loop over the quadrature points
358  } // end if (NUM_FIELD_MULT == something)
359 
360  // Put final values into target field
361  if (evalStyle_ == EvaluatorStyle::EVALUATES) {
362  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),[&] (const int& basis) {
363  field_(cell,basis) = tmp_field(basis);
364  });
365  }
366  else { // Contributed
367  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),[&] (const int& basis) {
368  field_(cell,basis) += tmp_field(basis);
369  });
370  }
371 
372  } // end of operator()
373 
375  //
376  // evaluateFields()
377  //
379  template<typename EvalT, typename Traits>
380  void
383  typename Traits::EvalData workset)
384  {
385  using Kokkos::parallel_for;
386  using Kokkos::TeamPolicy;
387 
388  // Grab the basis information.
389  basis_ = this->wda(workset).bases[basisIndex_]->weighted_div_basis;
390 
391  use_shared_memory = panzer::HP::inst().useSharedMemory<ScalarT>();
392 
393  if (use_shared_memory) {
394  int bytes;
395  if (Sacado::IsADType<ScalarT>::value) {
396  const int fadSize = Kokkos::dimension_scalar(field_.get_view());
397  bytes = scratch_view::shmem_size(1,fadSize) + scratch_view::shmem_size(basis_.extent(1),fadSize);
398  }
399  else
400  bytes = scratch_view::shmem_size(1) + scratch_view::shmem_size(basis_.extent(1));
401 
402  // The following if-block is for the sake of optimization depending on the
403  // number of field multipliers. The parallel_fors will loop over the cells
404  // in the Workset and execute operator()() above.
405  if (fieldMults_.size() == 0) {
406  auto policy = panzer::HP::inst().teamPolicy<ScalarT,SharedFieldMultTag<0>,PHX::Device>(workset.num_cells).set_scratch_size(0,Kokkos::PerTeam(bytes));
407  parallel_for(policy, *this, this->getName());
408  } else if (fieldMults_.size() == 1) {
409  auto policy = panzer::HP::inst().teamPolicy<ScalarT,SharedFieldMultTag<1>,PHX::Device>(workset.num_cells).set_scratch_size(0,Kokkos::PerTeam(bytes));
410  parallel_for(policy, *this, this->getName());
411  } else {
412  auto policy = panzer::HP::inst().teamPolicy<ScalarT,SharedFieldMultTag<-1>,PHX::Device>(workset.num_cells).set_scratch_size(0,Kokkos::PerTeam(bytes));
413  parallel_for(policy, *this, this->getName());
414  }
415  }
416  else {
417  // The following if-block is for the sake of optimization depending on the
418  // number of field multipliers. The parallel_fors will loop over the cells
419  // in the Workset and execute operator()() above.
420  if (fieldMults_.size() == 0) {
421  auto policy = panzer::HP::inst().teamPolicy<ScalarT,FieldMultTag<0>,PHX::Device>(workset.num_cells);
422  parallel_for(policy, *this, this->getName());
423  } else if (fieldMults_.size() == 1) {
424  auto policy = panzer::HP::inst().teamPolicy<ScalarT,FieldMultTag<1>,PHX::Device>(workset.num_cells);
425  parallel_for(policy, *this, this->getName());
426  } else {
427  auto policy = panzer::HP::inst().teamPolicy<ScalarT,FieldMultTag<-1>,PHX::Device>(workset.num_cells);
428  parallel_for(policy, *this, this->getName());
429  }
430  }
431  } // end of evaluateFields()
432 
434  //
435  // getValidParameters()
436  //
438  template<typename EvalT, typename TRAITS>
439  Teuchos::RCP<Teuchos::ParameterList>
441  {
442  using panzer::BasisIRLayout;
444  using std::string;
445  using std::vector;
446  using Teuchos::ParameterList;
447  using Teuchos::RCP;
448  using Teuchos::rcp;
449 
450  // Create a ParameterList with all the valid keys we support.
451  RCP<ParameterList> p = rcp(new ParameterList);
452  p->set<string>("Residual Name", "?");
453  p->set<string>("Value Name", "?");
454  p->set<string>("Test Field Name", "?");
455  RCP<BasisIRLayout> basis;
456  p->set("Basis", basis);
457  RCP<IntegrationRule> ir;
458  p->set("IR", ir);
459  p->set<double>("Multiplier", 1.0);
460  RCP<const vector<string>> fms;
461  p->set("Field Multipliers", fms);
462  return p;
463  } // end of getValidParameters()
464 
465 } // end of namespace panzer
466 
467 #endif // __Panzer_Integrator_DivBasisTimesScalar_impl_hpp__
int num_cells
DEPRECATED - use: numCells()
void evaluateFields(typename Traits::EvalData d)
Evaluate Fields.
Kokkos::TeamPolicy< TeamPolicyProperties... > teamPolicy(const int &league_size)
Returns a TeamPolicy for hierarchic parallelism.
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< 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 ( ).
PHX::MDField< const ScalarT, Cell, IP > scalar_
A field representing the scalar function we&#39;re integrating ( ).
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.
Teuchos::RCP< PHX::DataLayout > dl_scalar
Data layout for scalar fields.
panzer::EvaluatorStyle evalStyle
The EvaluatorStyle of the parent Integrator_CurlBasisDotVector object.
This empty struct allows us to optimize operator()() depending on the number of field multipliers...
KOKKOS_INLINE_FUNCTION void operator()(const FieldMultTag< NUM_FIELD_MULT > &tag, const Kokkos::TeamPolicy< PHX::exec_space >::member_type &team) const
Perform the integration. Main memory version.
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< ScalarT, Cell, BASIS > field_
A field to which we&#39;ll contribute, or in which we&#39;ll store, the result of computing this integral...
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
Get Valid Parameters.
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
Post-Registration Setup.
Integrator_DivBasisTimesScalar(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.
static HP & inst()
Private ctor.
Kokkos::View< ScalarT *,typename PHX::DevLayout< ScalarT >::type, typename PHX::exec_space::scratch_memory_space, Kokkos::MemoryUnmanaged > scratch_view
Type for shared memory.
Description and data layouts associated with a particular basis.
This empty struct allows us to optimize operator()() depending on the number of field multipliers...
Teuchos::RCP< PHX::DataLayout > functional
<Cell,Basis>
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_
std::vector< PHX::MDField< const ScalarT, Cell, IP > > fieldMults_
The (possibly empty) list of fields that are multipliers out in front of the integral ( ...