Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_Integrator_GradBasisCrossVector_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_GRADBASISCROSSVECTOR_IMPL_HPP
44 #define PANZER_EVALUATOR_GRADBASISCROSSVECTOR_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& vecName,
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  const Teuchos::RCP<PHX::DataLayout>& vecDL /* = Teuchos::null */)
82  :
83  evalStyle_(evalStyle),
84  multiplier_(multiplier),
85  numDims_(resNames.size()),
86  numGradDims_(ir.dl_vector->extent(2)),
87  basisName_(basis.name())
88  {
89  using Kokkos::View;
90  using panzer::BASIS;
91  using panzer::Cell;
93  using panzer::IP;
94  using PHX::DataLayout;
95  using PHX::Device;
96  using PHX::DevLayout;
97  using PHX::MDField;
98  using std::invalid_argument;
99  using std::logic_error;
100  using std::size_t;
101  using std::string;
102  using Teuchos::RCP;
103 
104  // Ensure the input makes sense.
105  TEUCHOS_TEST_FOR_EXCEPTION(numDims_ != 3, invalid_argument, "Error: " \
106  "Integrator_GradBasisCrossVector called with the number of residual " \
107  "names not equal to three.")
108  for (const auto& name : resNames)
109  TEUCHOS_TEST_FOR_EXCEPTION(name == "", invalid_argument, "Error: " \
110  "Integrator_GradBasisCrossVector called with an empty residual name.")
111  TEUCHOS_TEST_FOR_EXCEPTION(vecName == "", invalid_argument, "Error: " \
112  "Integrator_GradBasisCrossVector called with an empty vector name.")
113  RCP<const PureBasis> tmpBasis = basis.getBasis();
114  TEUCHOS_TEST_FOR_EXCEPTION(not tmpBasis->supportsGrad(), logic_error,
115  "Error: Integrator_GradBasisCrossVector: Basis of type \""
116  << tmpBasis->name() << "\" does not support the gradient operator.")
117  RCP<DataLayout> tmpVecDL = ir.dl_vector;
118  if (not vecDL.is_null())
119  {
120  tmpVecDL = vecDL;
122  tmpVecDL->extent(2) < ir.dl_vector->extent(2), logic_error,
123  "Error: Integrator_GradBasisCrossVector: Dimension of space " \
124  "exceeds dimension of Vector Data Layout.")
126  static_cast<int>(vecDL->extent(2)), logic_error, "Error: " \
127  "Integrator_GradBasisCrossVector: The vector must be the same " \
128  "length as the number of residuals.")
129  } // end if (not vecDL.is_null())
131  "Error: Integrator_GradBasisCrossVector: The vector must have at " \
132  "least as many components as there are dimensions in the mesh.")
133 
134  // Create the field for the vector-valued function we're integrating.
135  vector_ = MDField<const ScalarT, Cell, IP, Dim>(vecName, tmpVecDL);
136  this->addDependentField(vector_);
137 
138  // Create the fields that we're either contributing to or evaluating
139  // (storing).
140  fields_ = Kokkos::View<PHX::MDField<ScalarT, Cell, BASIS>*>("Integrator_GradBasisCrossVector::fields_", resNames.size());
141 
142  {
143  int i=0;
144  for (const auto& name : resNames)
145  fields_(i++) = MDField<ScalarT, Cell, BASIS>(name, basis.functional);
146  } // end loop over resNames
147  for (std::size_t i=0; i< fields_.extent(0); ++i) {
148  const auto& field = fields_(i);
150  this->addContributedField(field);
151  else // if (evalStyle_ == EvaluatorStyle::EVALUATES)
152  this->addEvaluatedField(field);
153  }
154  // Add the dependent field multipliers, if there are any.
155  int i = 0;
156  fieldMults_.resize(fmNames.size());
158  typename DevLayout<ScalarT>::type, Device>*>(
159  "GradBasisCrossVector::KokkosFieldMultipliers", fmNames.size());
160  for (const auto& name : fmNames)
161  {
162  fieldMults_[i++] = MDField<const ScalarT, Cell, IP>(name, ir.dl_scalar);
163  this->addDependentField(fieldMults_[i - 1]);
164  } // end loop over the field multipliers
165 
166  // Set the name of this object.
167  string n("Integrator_GradBasisCrossVector (");
169  n += "CONTRIBUTES";
170  else // if (evalStyle_ == EvaluatorStyle::EVALUATES)
171  n += "EVALUATES";
172  n += "): {";
173  for (size_t j=0; j < fields_.extent(0) - 1; ++j)
174  n += fields_[j].fieldTag().name() + ", ";
175  n += fields_(fields_.extent(0)-1).fieldTag().name() + "}";
176  this->setName(n);
177  } // end of Constructor
178 
180  //
181  // ParameterList Constructor
182  //
184  template<typename EvalT, typename Traits>
187  const Teuchos::ParameterList& p)
188  :
190  panzer::EvaluatorStyle::EVALUATES,
191  p.get<const std::vector<std::string>>("Residual Names"),
192  p.get<std::string>("Vector Name"),
193  (*p.get<Teuchos::RCP<panzer::BasisIRLayout>>("Basis")),
194  (*p.get<Teuchos::RCP<panzer::IntegrationRule>>("IR")),
195  p.get<double>("Multiplier"),
196  p.isType<Teuchos::RCP<const std::vector<std::string>>>
197  ("Field Multipliers") ?
198  (*p.get<Teuchos::RCP<const std::vector<std::string>>>
199  ("Field Multipliers")) : std::vector<std::string>(),
200  p.isType<Teuchos::RCP<PHX::DataLayout>>("Data Layout Vector") ?
201  p.get<Teuchos::RCP<PHX::DataLayout>>("Data Layout Vector") :
202  Teuchos::null)
203  {
205  using Teuchos::RCP;
206 
207  // Ensure that the input ParameterList didn't contain any bogus entries.
208  RCP<ParameterList> validParams = this->getValidParameters();
209  p.validateParameters(*validParams);
210  } // end of ParameterList Constructor
211 
213  //
214  // postRegistrationSetup()
215  //
217  template<typename EvalT, typename Traits>
218  void
221  typename Traits::SetupData sd,
222  PHX::FieldManager<Traits>& /* fm */)
223  {
225  using panzer::getBasisIndex;
226  using PHX::Device;
227 
228  // Get the Kokkos::Views of the field multipliers.
229  for (size_t i(0); i < fieldMults_.size(); ++i)
230  kokkosFieldMults_(i) = fieldMults_[i].get_static_view();
231  Device().fence();
232 
233  // Determine the index in the Workset bases for our particular basis name.
234  basisIndex_ = getBasisIndex(basisName_, (*sd.worksets_)[0], this->wda);
235  } // end of postRegistrationSetup()
236 
238  //
239  // operator()()
240  //
242  template<typename EvalT, typename Traits>
243  template<int NUM_FIELD_MULT>
244  KOKKOS_INLINE_FUNCTION
245  void
248  const FieldMultTag<NUM_FIELD_MULT>& /* tag */,
249  const size_t& cell) const
250  {
252 
253  // Initialize the evaluated fields.
254  const int numBases(fields_[0].extent(1)), numQP(vector_.extent(1));
255  if (evalStyle_ == EvaluatorStyle::EVALUATES)
256  for (int dim(0); dim < numDims_; ++dim)
257  for (int basis(0); basis < numBases; ++basis)
258  fields_[dim](cell, basis) = 0.0;
259 
260  // The following if-block is for the sake of optimization depending on the
261  // number of field multipliers.
262  ScalarT tmp[3];
263  const int X(0), Y(1), Z(2);
264  if (NUM_FIELD_MULT == 0)
265  {
266  if (numGradDims_ == 1)
267  {
268  for (int qp(0); qp < numQP; ++qp)
269  {
270  tmp[Y] = multiplier_ * vector_(cell, qp, Y);
271  tmp[Z] = multiplier_ * vector_(cell, qp, Z);
272  for (int basis(0); basis < numBases; ++basis)
273  {
274  fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
275  fields_[Z](cell, basis) += -tmp[Y] * basis_(cell, basis, qp, X);
276  } // end loop over the bases
277  } // end loop over the quadrature points
278  }
279  else if (numGradDims_ == 2)
280  {
281  for (int qp(0); qp < numQP; ++qp)
282  {
283  tmp[X] = multiplier_ * vector_(cell, qp, X);
284  tmp[Y] = multiplier_ * vector_(cell, qp, Y);
285  tmp[Z] = multiplier_ * vector_(cell, qp, Z);
286  for (int basis(0); basis < numBases; ++basis)
287  {
288  fields_[X](cell, basis) += -tmp[Z] * basis_(cell, basis, qp, Y);
289  fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
290  fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
291  tmp[Y] * basis_(cell, basis, qp, X);
292  } // end loop over the bases
293  } // end loop over the quadrature points
294  }
295  else if (numGradDims_ == 3)
296  {
297  for (int qp(0); qp < numQP; ++qp)
298  {
299  tmp[X] = multiplier_ * vector_(cell, qp, X);
300  tmp[Y] = multiplier_ * vector_(cell, qp, Y);
301  tmp[Z] = multiplier_ * vector_(cell, qp, Z);
302  for (int basis(0); basis < numBases; ++basis)
303  {
304  fields_[X](cell, basis) += tmp[Y] * basis_(cell, basis, qp, Z) -
305  tmp[Z] * basis_(cell, basis, qp, Y);
306  fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X) -
307  tmp[X] * basis_(cell, basis, qp, Z);
308  fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
309  tmp[Y] * basis_(cell, basis, qp, X);
310  } // end loop over the bases
311  } // end loop over the quadrature points
312  } // end if (numGradDims_ == something)
313  }
314  else if (NUM_FIELD_MULT == 1)
315  {
316  if (numGradDims_ == 1)
317  {
318  for (int qp(0); qp < numQP; ++qp)
319  {
320  tmp[Y] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
321  vector_(cell, qp, Y);
322  tmp[Z] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
323  vector_(cell, qp, Z);
324  for (int basis(0); basis < numBases; ++basis)
325  {
326  fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
327  fields_[Z](cell, basis) += -tmp[Y] * basis_(cell, basis, qp, X);
328  } // end loop over the bases
329  } // end loop over the quadrature points
330  }
331  else if (numGradDims_ == 2)
332  {
333  for (int qp(0); qp < numQP; ++qp)
334  {
335  tmp[X] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
336  vector_(cell, qp, X);
337  tmp[Y] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
338  vector_(cell, qp, Y);
339  tmp[Z] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
340  vector_(cell, qp, Z);
341  for (int basis(0); basis < numBases; ++basis)
342  {
343  fields_[X](cell, basis) += -tmp[Z] * basis_(cell, basis, qp, Y);
344  fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
345  fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
346  tmp[Y] * basis_(cell, basis, qp, X);
347  } // end loop over the bases
348  } // end loop over the quadrature points
349  }
350  else if (numGradDims_ == 3)
351  {
352  for (int qp(0); qp < numQP; ++qp)
353  {
354  tmp[X] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
355  vector_(cell, qp, X);
356  tmp[Y] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
357  vector_(cell, qp, Y);
358  tmp[Z] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
359  vector_(cell, qp, Z);
360  for (int basis(0); basis < numBases; ++basis)
361  {
362  fields_[X](cell, basis) += tmp[Y] * basis_(cell, basis, qp, Z) -
363  tmp[Z] * basis_(cell, basis, qp, Y);
364  fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X) -
365  tmp[X] * basis_(cell, basis, qp, Z);
366  fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
367  tmp[Y] * basis_(cell, basis, qp, X);
368  } // end loop over the bases
369  } // end loop over the quadrature points
370  } // end if (numGradDims_ == something)
371  }
372  else
373  {
374  const int numFieldMults(kokkosFieldMults_.extent(0));
375  if (numGradDims_ == 1)
376  {
377  for (int qp(0); qp < numQP; ++qp)
378  {
379  tmp[Y] = multiplier_ * vector_(cell, qp, Y);
380  tmp[Z] = multiplier_ * vector_(cell, qp, Z);
381  for (int fm(0); fm < numFieldMults; ++fm)
382  {
383  tmp[Y] *= kokkosFieldMults_(fm)(cell, qp);
384  tmp[Z] *= kokkosFieldMults_(fm)(cell, qp);
385  } // end loop over the field multipliers
386  for (int basis(0); basis < numBases; ++basis)
387  {
388  fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
389  fields_[Z](cell, basis) += -tmp[Y] * basis_(cell, basis, qp, X);
390  } // end loop over the bases
391  } // end loop over the quadrature points
392  }
393  else if (numGradDims_ == 2)
394  {
395  for (int qp(0); qp < numQP; ++qp)
396  {
397  tmp[X] = multiplier_ * vector_(cell, qp, X);
398  tmp[Y] = multiplier_ * vector_(cell, qp, Y);
399  tmp[Z] = multiplier_ * vector_(cell, qp, Z);
400  for (int fm(0); fm < numFieldMults; ++fm)
401  {
402  tmp[X] *= kokkosFieldMults_(fm)(cell, qp);
403  tmp[Y] *= kokkosFieldMults_(fm)(cell, qp);
404  tmp[Z] *= kokkosFieldMults_(fm)(cell, qp);
405  } // end loop over the field multipliers
406  for (int basis(0); basis < numBases; ++basis)
407  {
408  fields_[X](cell, basis) += -tmp[Z] * basis_(cell, basis, qp, Y);
409  fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
410  fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
411  tmp[Y] * basis_(cell, basis, qp, X);
412  } // end loop over the bases
413  } // end loop over the quadrature points
414  }
415  else if (numGradDims_ == 3)
416  {
417  for (int qp(0); qp < numQP; ++qp)
418  {
419  tmp[X] = multiplier_ * vector_(cell, qp, X);
420  tmp[Y] = multiplier_ * vector_(cell, qp, Y);
421  tmp[Z] = multiplier_ * vector_(cell, qp, Z);
422  for (int fm(0); fm < numFieldMults; ++fm)
423  {
424  tmp[X] *= kokkosFieldMults_(fm)(cell, qp);
425  tmp[Y] *= kokkosFieldMults_(fm)(cell, qp);
426  tmp[Z] *= kokkosFieldMults_(fm)(cell, qp);
427  } // end loop over the field multipliers
428  for (int basis(0); basis < numBases; ++basis)
429  {
430  fields_[X](cell, basis) += tmp[Y] * basis_(cell, basis, qp, Z) -
431  tmp[Z] * basis_(cell, basis, qp, Y);
432  fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X) -
433  tmp[X] * basis_(cell, basis, qp, Z);
434  fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
435  tmp[Y] * basis_(cell, basis, qp, X);
436  } // end loop over the bases
437  } // end loop over the quadrature points
438  } // end if (numGradDims_ == something)
439  } // end if (NUM_FIELD_MULT == something)
440  } // end of operator()()
441 
443  //
444  // evaluateFields()
445  //
447  template<typename EvalT, typename Traits>
448  void
451  typename Traits::EvalData workset)
452  {
453  using Kokkos::parallel_for;
454  using Kokkos::RangePolicy;
455 
456  // Grab the basis information.
457  basis_ = this->wda(workset).bases[basisIndex_]->weighted_grad_basis;
458 
459  // The following if-block is for the sake of optimization depending on the
460  // number of field multipliers. The parallel_fors will loop over the cells
461  // in the Workset and execute operator()() above.
462  if (fieldMults_.size() == 0)
463  parallel_for(RangePolicy<FieldMultTag<0>>(0, workset.num_cells), *this);
464  else if (fieldMults_.size() == 1)
465  parallel_for(RangePolicy<FieldMultTag<1>>(0, workset.num_cells), *this);
466  else
467  parallel_for(RangePolicy<FieldMultTag<-1>>(0, workset.num_cells), *this);
468  } // end of evaluateFields()
469 
471  //
472  // getValidParameters()
473  //
475  template<typename EvalT, typename TRAITS>
479  {
480  using panzer::BasisIRLayout;
482  using PHX::DataLayout;
483  using std::string;
484  using std::vector;
486  using Teuchos::RCP;
487  using Teuchos::rcp;
488 
489  // Create a ParameterList with all the valid keys we support.
490  RCP<ParameterList> p = rcp(new ParameterList);
491 
492  RCP<const vector<string>> resNames;
493  p->set("Residual Names", resNames);
494  p->set<string>("Vector Name", "?");
495  RCP<BasisIRLayout> basis;
496  p->set("Basis", basis);
498  p->set("IR", ir);
499  p->set<double>("Multiplier", 1.0);
501  p->set("Field Multipliers", fms);
502  RCP<DataLayout> vecDL;
503  p->set("Data Layout Vector", vecDL);
504 
505  return p;
506  } // end of getValidParameters()
507 
508 } // end of namespace panzer
509 
510 #endif // PANZER_EVALUATOR_GRADBASISCROSSVECTOR_IMPL_HPP
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.
std::vector< PHX::MDField< const ScalarT, Cell, IP > > fieldMults_
The (possibly empty) list of fields that are multipliers out in front of the integral ( ...
Integrator_GradBasisCrossVector(const panzer::EvaluatorStyle &evalStyle, const std::vector< std::string > &resNames, const std::string &vecName, const panzer::BasisIRLayout &basis, const panzer::IntegrationRule &ir, const double &multiplier=1, const std::vector< std::string > &fmNames=std::vector< std::string >(), const Teuchos::RCP< PHX::DataLayout > &vecDL=Teuchos::null)
Main Constructor.
#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...
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...
EvaluatorStyle
An indication of how an Evaluator will behave.
const panzer::EvaluatorStyle evalStyle_
An enum determining the behavior of this Evaluator.
Kokkos::View< 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::RCP< Teuchos::ParameterList > getValidParameters() const
Get Valid Parameters.
Teuchos::RCP< const PureBasis > getBasis() const
PHX::MDField< const ScalarT, Cell, IP, Dim > vector_
A field representing the vector-valued function we&#39;re integrating ( ).
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.
void validateParameters(ParameterList const &validParamList, int const depth=1000, EValidateUsed const validateUsed=VALIDATE_USED_ENABLED, EValidateDefaults const validateDefaults=VALIDATE_DEFAULTS_ENABLED) const
Teuchos::RCP< PHX::DataLayout > dl_vector
Data layout for vector fields.
KOKKOS_INLINE_FUNCTION void operator()(const FieldMultTag< NUM_FIELD_MULT > &tag, const std::size_t &cell) const
Perform the integration.
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...
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
Post-Registration Setup.
const std::vector< std::pair< int, LocalOrdinal > > &pid_and_lid const
int numDims_
The number of dimensions associated with the vector.
int numGradDims_
The number of dimensions associated with the gradient.
Teuchos::RCP< PHX::DataLayout > functional
&lt;Cell,Basis&gt;
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_
bool is_null() const