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 Kokkos::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**,typename PHX::DevLayout<ScalarT>::type,PHX::Device>*>("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  using PHX::Device;
192 
193  // Get the Kokkos::Views of the field multipliers.
194  for (size_t i(0); i < fieldMults_.size(); ++i)
195  kokkosFieldMults_(i) = fieldMults_[i].get_static_view();
196  Device().fence();
197 
198  // Determine the index in the Workset bases for our particular basis name.
199  basisIndex_ = getBasisIndex(basisName_, (*sd.worksets_)[0], this->wda);
200  } // end of postRegistrationSetup()
201 
203  //
204  // No shared memory operator()()
205  //
207  template<typename EvalT, typename Traits>
208  template<int NUM_FIELD_MULT>
209  KOKKOS_INLINE_FUNCTION
210  void
213  const FieldMultTag<NUM_FIELD_MULT>& /* tag */,
214  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
215  {
217  const int cell = team.league_rank();
218 
219  // Initialize the evaluated field.
220  const int numQP(scalar_.extent(1)), numBases(basis_.extent(1));
221  if (evalStyle_ == EvaluatorStyle::EVALUATES) {
222  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (const int& basis) {
223  field_(cell, basis) = 0.0;
224  });
225  }
226 
227  // The following if-block is for the sake of optimization depending on the
228  // number of field multipliers.
229  ScalarT tmp;
230  if (NUM_FIELD_MULT == 0)
231  {
232  // Loop over the quadrature points, scale the integrand by the
233  // multiplier, and then perform the actual integration, looping over the
234  // bases.
235  for (int qp(0); qp < numQP; ++qp)
236  {
237  tmp = multiplier_ * scalar_(cell, qp);
238  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (const int& basis) {
239  field_(cell, basis) += basis_(cell, basis, qp) * tmp;
240  });
241  } // end loop over the quadrature points
242  }
243  else if (NUM_FIELD_MULT == 1)
244  {
245  // Loop over the quadrature points, scale the integrand by the multiplier
246  // and the single field multiplier, and then perform the actual
247  // integration, looping over the bases.
248  for (int qp(0); qp < numQP; ++qp)
249  {
250  tmp = multiplier_ * scalar_(cell, qp) * kokkosFieldMults_(0)(cell, qp);
251  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (const int& basis) {
252  field_(cell, basis) += basis_(cell, basis, qp) * tmp;
253  });
254  } // end loop over the quadrature points
255  }
256  else
257  {
258  // Loop over the quadrature points and pre-multiply all the field
259  // multipliers together, scale the integrand by the multiplier and the
260  // combination of field multipliers, and then perform the actual
261  // integration, looping over the bases.
262  const int numFieldMults(kokkosFieldMults_.extent(0));
263  for (int qp(0); qp < numQP; ++qp)
264  {
265  ScalarT fieldMultsTotal(1);
266  for (int fm(0); fm < numFieldMults; ++fm)
267  fieldMultsTotal *= kokkosFieldMults_(fm)(cell, qp);
268  tmp = multiplier_ * scalar_(cell, qp) * fieldMultsTotal;
269  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (const int& basis) {
270  field_(cell, basis) += basis_(cell, basis, qp) * tmp;
271  });
272  } // end loop over the quadrature points
273  } // end if (NUM_FIELD_MULT == something)
274  } // end of operator()
275 
277  //
278  // Shared memory operator()()
279  //
281  template<typename EvalT, typename Traits>
282  template<int NUM_FIELD_MULT>
283  KOKKOS_INLINE_FUNCTION
284  void
287  const SharedFieldMultTag<NUM_FIELD_MULT>& /* tag */,
288  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
289  {
291  const int cell = team.league_rank();
292  const int numQP = scalar_.extent(1);
293  const int numBases = basis_.extent(1);
294  const int fadSize = Kokkos::dimension_scalar(field_.get_view());
295 
296  scratch_view tmp;
297  scratch_view tmp_field;
298  if (Sacado::IsADType<ScalarT>::value) {
299  tmp = scratch_view(team.team_shmem(),1,fadSize);
300  tmp_field = scratch_view(team.team_shmem(),numBases,fadSize);
301  }
302  else {
303  tmp = scratch_view(team.team_shmem(),1);
304  tmp_field = scratch_view(team.team_shmem(),numBases);
305  }
306 
307  // Initialize the evaluated field.
308  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (const int& basis) {
309  tmp_field(basis) = 0.0;
310  });
311 
312  // The following if-block is for the sake of optimization depending on the
313  // number of field multipliers.
314  if (NUM_FIELD_MULT == 0)
315  {
316  // Loop over the quadrature points, scale the integrand by the
317  // multiplier, and then perform the actual integration, looping over the
318  // bases.
319  for (int qp(0); qp < numQP; ++qp)
320  {
321  team.team_barrier();
322  tmp(0) = multiplier_ * scalar_(cell, qp);
323  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (const int& basis) {
324  tmp_field(basis) += basis_(cell, basis, qp) * tmp(0);
325  });
326  } // end loop over the quadrature points
327  }
328  else if (NUM_FIELD_MULT == 1)
329  {
330  // Loop over the quadrature points, scale the integrand by the multiplier
331  // and the single field multiplier, and then perform the actual
332  // integration, looping over the bases.
333  for (int qp(0); qp < numQP; ++qp)
334  {
335  team.team_barrier();
336  tmp(0) = multiplier_ * scalar_(cell, qp) * kokkosFieldMults_(0)(cell, qp);
337  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (const int& basis) {
338  tmp_field(basis) += basis_(cell, basis, qp) * tmp(0);
339  });
340  } // end loop over the quadrature points
341  }
342  else
343  {
344  // Loop over the quadrature points and pre-multiply all the field
345  // multipliers together, scale the integrand by the multiplier and the
346  // combination of field multipliers, and then perform the actual
347  // integration, looping over the bases.
348  const int numFieldMults = kokkosFieldMults_.extent(0);
349  for (int qp(0); qp < numQP; ++qp)
350  {
351  team.team_barrier();
352  ScalarT fieldMultsTotal(1); // need shared mem here
353  for (int fm(0); fm < numFieldMults; ++fm)
354  fieldMultsTotal *= kokkosFieldMults_(fm)(cell, qp);
355  tmp(0) = multiplier_ * scalar_(cell, qp) * fieldMultsTotal;
356  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),KOKKOS_LAMBDA (const int& basis) {
357  tmp_field(basis) += basis_(cell, basis, qp) * tmp(0);
358  });
359  } // end loop over the quadrature points
360  } // end if (NUM_FIELD_MULT == something)
361 
362  // Put final values into target field
363  if (evalStyle_ == EvaluatorStyle::EVALUATES) {
364  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),[&] (const int& basis) {
365  field_(cell,basis) = tmp_field(basis);
366  });
367  }
368  else { // Contributed
369  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),[&] (const int& basis) {
370  field_(cell,basis) += tmp_field(basis);
371  });
372  }
373 
374  } // end of operator()
375 
377  //
378  // evaluateFields()
379  //
381  template<typename EvalT, typename Traits>
382  void
385  typename Traits::EvalData workset)
386  {
387  using Kokkos::parallel_for;
388  using Kokkos::TeamPolicy;
389 
390  // Grab the basis information.
391  basis_ = this->wda(workset).bases[basisIndex_]->weighted_div_basis;
392 
393  use_shared_memory = panzer::HP::inst().useSharedMemory<ScalarT>();
394 
395  if (use_shared_memory) {
396  int bytes;
397  if (Sacado::IsADType<ScalarT>::value) {
398  const int fadSize = Kokkos::dimension_scalar(field_.get_view());
399  bytes = scratch_view::shmem_size(1,fadSize) + scratch_view::shmem_size(basis_.extent(1),fadSize);
400  }
401  else
402  bytes = scratch_view::shmem_size(1) + scratch_view::shmem_size(basis_.extent(1));
403 
404  // The following if-block is for the sake of optimization depending on the
405  // number of field multipliers. The parallel_fors will loop over the cells
406  // in the Workset and execute operator()() above.
407  if (fieldMults_.size() == 0) {
408  auto policy = panzer::HP::inst().teamPolicy<ScalarT,SharedFieldMultTag<0>,PHX::Device>(workset.num_cells).set_scratch_size(0,Kokkos::PerTeam(bytes));
409  parallel_for(policy, *this, this->getName());
410  } else if (fieldMults_.size() == 1) {
411  auto policy = panzer::HP::inst().teamPolicy<ScalarT,SharedFieldMultTag<1>,PHX::Device>(workset.num_cells).set_scratch_size(0,Kokkos::PerTeam(bytes));
412  parallel_for(policy, *this, this->getName());
413  } else {
414  auto policy = panzer::HP::inst().teamPolicy<ScalarT,SharedFieldMultTag<-1>,PHX::Device>(workset.num_cells).set_scratch_size(0,Kokkos::PerTeam(bytes));
415  parallel_for(policy, *this, this->getName());
416  }
417  }
418  else {
419  // The following if-block is for the sake of optimization depending on the
420  // number of field multipliers. The parallel_fors will loop over the cells
421  // in the Workset and execute operator()() above.
422  if (fieldMults_.size() == 0) {
423  auto policy = panzer::HP::inst().teamPolicy<ScalarT,FieldMultTag<0>,PHX::Device>(workset.num_cells);
424  parallel_for(policy, *this, this->getName());
425  } else if (fieldMults_.size() == 1) {
426  auto policy = panzer::HP::inst().teamPolicy<ScalarT,FieldMultTag<1>,PHX::Device>(workset.num_cells);
427  parallel_for(policy, *this, this->getName());
428  } else {
429  auto policy = panzer::HP::inst().teamPolicy<ScalarT,FieldMultTag<-1>,PHX::Device>(workset.num_cells);
430  parallel_for(policy, *this, this->getName());
431  }
432  }
433  } // end of evaluateFields()
434 
436  //
437  // getValidParameters()
438  //
440  template<typename EvalT, typename TRAITS>
441  Teuchos::RCP<Teuchos::ParameterList>
443  {
444  using panzer::BasisIRLayout;
446  using std::string;
447  using std::vector;
448  using Teuchos::ParameterList;
449  using Teuchos::RCP;
450  using Teuchos::rcp;
451 
452  // Create a ParameterList with all the valid keys we support.
453  RCP<ParameterList> p = rcp(new ParameterList);
454  p->set<string>("Residual Name", "?");
455  p->set<string>("Value Name", "?");
456  p->set<string>("Test Field Name", "?");
457  RCP<BasisIRLayout> basis;
458  p->set("Basis", basis);
459  RCP<IntegrationRule> ir;
460  p->set("IR", ir);
461  p->set<double>("Multiplier", 1.0);
462  RCP<const vector<string>> fms;
463  p->set("Field Multipliers", fms);
464  return p;
465  } // end of getValidParameters()
466 
467 } // end of namespace panzer
468 
469 #endif // __Panzer_Integrator_DivBasisTimesScalar_impl_hpp__
panzer::EvaluatorStyle evalStyle
The EvaluatorStyle of the parent Integrator_CurlBasisDotVector object.
double multiplier
The scalar multiplier out in front of the integral ( ).
Teuchos::RCP< PHX::DataLayout > functional
<Cell,Basis>
Teuchos::RCP< const PureBasis > getBasis() const
static HP & inst()
Private ctor.
Kokkos::TeamPolicy< TeamPolicyProperties... > teamPolicy(const int &league_size)
Returns a TeamPolicy for hierarchic parallelism.
PHX::MDField< const ScalarT, Cell, IP > scalar_
A field representing the scalar function we're integrating ( ).
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
Get Valid Parameters.
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.
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.
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
Post-Registration Setup.
Kokkos::View< ScalarT *,typename PHX::DevLayout< ScalarT >::type, typename PHX::exec_space::scratch_memory_space, Kokkos::MemoryUnmanaged > scratch_view
Type for shared memory.
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 evaluateFields(typename Traits::EvalData d)
Evaluate Fields.
const panzer::EvaluatorStyle evalStyle_
An enum determining the behavior of this Evaluator.
PHX::MDField< ScalarT, Cell, BASIS > field_
A field to which we'll contribute, or in which we'll store, the result of computing this integral.
std::vector< PHX::MDField< const ScalarT, Cell, IP > > fieldMults_
The (possibly empty) list of fields that are multipliers out in front of the integral ( ,...
Teuchos::RCP< PHX::DataLayout > dl_scalar
Data layout for scalar fields.
Description and data layouts associated with a particular basis.
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.
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.
EvaluatorStyle
An indication of how an Evaluator will behave.
This empty struct allows us to optimize operator()() depending on the number of field multipliers.
This empty struct allows us to optimize operator()() depending on the number of field multipliers.
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_