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_.clear();
141  for (const auto& name : resNames)
142  {
143  MDField<ScalarT, Cell, BASIS> res(name, basis.functional);
144  fields_.push_back(res);
145  } // end loop over resNames
146  for (const auto& field : fields_)
148  this->addContributedField(field);
149  else // if (evalStyle_ == EvaluatorStyle::EVALUATES)
150  this->addEvaluatedField(field);
151 
152  // Add the dependent field multipliers, if there are any.
153  int i = 0;
154  fieldMults_.resize(fmNames.size());
156  typename DevLayout<ScalarT>::type, Device>*>(
157  "GradBasisCrossVector::KokkosFieldMultipliers", fmNames.size());
158  for (const auto& name : fmNames)
159  {
160  fieldMults_[i++] = MDField<const ScalarT, Cell, IP>(name, ir.dl_scalar);
161  this->addDependentField(fieldMults_[i - 1]);
162  } // end loop over the field multipliers
163 
164  // Set the name of this object.
165  string n("Integrator_GradBasisCrossVector (");
167  n += "CONTRIBUTES";
168  else // if (evalStyle_ == EvaluatorStyle::EVALUATES)
169  n += "EVALUATES";
170  n += "): {";
171  for (size_t j=0; j < fields_.size() - 1; ++j)
172  n += fields_[j].fieldTag().name() + ", ";
173  n += fields_.back().fieldTag().name() + "}";
174  this->setName(n);
175  } // end of Constructor
176 
178  //
179  // ParameterList Constructor
180  //
182  template<typename EvalT, typename Traits>
185  const Teuchos::ParameterList& p)
186  :
188  panzer::EvaluatorStyle::EVALUATES,
189  p.get<const std::vector<std::string>>("Residual Names"),
190  p.get<std::string>("Vector Name"),
191  (*p.get<Teuchos::RCP<panzer::BasisIRLayout>>("Basis")),
192  (*p.get<Teuchos::RCP<panzer::IntegrationRule>>("IR")),
193  p.get<double>("Multiplier"),
194  p.isType<Teuchos::RCP<const std::vector<std::string>>>
195  ("Field Multipliers") ?
196  (*p.get<Teuchos::RCP<const std::vector<std::string>>>
197  ("Field Multipliers")) : std::vector<std::string>(),
198  p.isType<Teuchos::RCP<PHX::DataLayout>>("Data Layout Vector") ?
199  p.get<Teuchos::RCP<PHX::DataLayout>>("Data Layout Vector") :
200  Teuchos::null)
201  {
203  using Teuchos::RCP;
204 
205  // Ensure that the input ParameterList didn't contain any bogus entries.
206  RCP<ParameterList> validParams = this->getValidParameters();
207  p.validateParameters(*validParams);
208  } // end of ParameterList Constructor
209 
211  //
212  // postRegistrationSetup()
213  //
215  template<typename EvalT, typename Traits>
216  void
219  typename Traits::SetupData sd,
220  PHX::FieldManager<Traits>& /* fm */)
221  {
223  using panzer::getBasisIndex;
224  using PHX::Device;
225 
226  // Get the Kokkos::Views of the field multipliers.
227  for (size_t i(0); i < fieldMults_.size(); ++i)
228  kokkosFieldMults_(i) = fieldMults_[i].get_static_view();
229  Device::fence();
230 
231  // Determine the index in the Workset bases for our particular basis name.
232  basisIndex_ = getBasisIndex(basisName_, (*sd.worksets_)[0], this->wda);
233  } // end of postRegistrationSetup()
234 
236  //
237  // operator()()
238  //
240  template<typename EvalT, typename Traits>
241  template<int NUM_FIELD_MULT>
242  KOKKOS_INLINE_FUNCTION
243  void
246  const FieldMultTag<NUM_FIELD_MULT>& /* tag */,
247  const size_t& cell) const
248  {
250 
251  // Initialize the evaluated fields.
252  const int numBases(fields_[0].extent(1)), numQP(vector_.extent(1));
253  if (evalStyle_ == EvaluatorStyle::EVALUATES)
254  for (int dim(0); dim < numDims_; ++dim)
255  for (int basis(0); basis < numBases; ++basis)
256  fields_[dim](cell, basis) = 0.0;
257 
258  // The following if-block is for the sake of optimization depending on the
259  // number of field multipliers.
260  ScalarT tmp[3];
261  const int X(0), Y(1), Z(2);
262  if (NUM_FIELD_MULT == 0)
263  {
264  if (numGradDims_ == 1)
265  {
266  for (int qp(0); qp < numQP; ++qp)
267  {
268  tmp[Y] = multiplier_ * vector_(cell, qp, Y);
269  tmp[Z] = multiplier_ * vector_(cell, qp, Z);
270  for (int basis(0); basis < numBases; ++basis)
271  {
272  fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
273  fields_[Z](cell, basis) += -tmp[Y] * basis_(cell, basis, qp, X);
274  } // end loop over the bases
275  } // end loop over the quadrature points
276  }
277  else if (numGradDims_ == 2)
278  {
279  for (int qp(0); qp < numQP; ++qp)
280  {
281  tmp[X] = multiplier_ * vector_(cell, qp, X);
282  tmp[Y] = multiplier_ * vector_(cell, qp, Y);
283  tmp[Z] = multiplier_ * vector_(cell, qp, Z);
284  for (int basis(0); basis < numBases; ++basis)
285  {
286  fields_[X](cell, basis) += -tmp[Z] * basis_(cell, basis, qp, Y);
287  fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
288  fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
289  tmp[Y] * basis_(cell, basis, qp, X);
290  } // end loop over the bases
291  } // end loop over the quadrature points
292  }
293  else if (numGradDims_ == 3)
294  {
295  for (int qp(0); qp < numQP; ++qp)
296  {
297  tmp[X] = multiplier_ * vector_(cell, qp, X);
298  tmp[Y] = multiplier_ * vector_(cell, qp, Y);
299  tmp[Z] = multiplier_ * vector_(cell, qp, Z);
300  for (int basis(0); basis < numBases; ++basis)
301  {
302  fields_[X](cell, basis) += tmp[Y] * basis_(cell, basis, qp, Z) -
303  tmp[Z] * basis_(cell, basis, qp, Y);
304  fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X) -
305  tmp[X] * basis_(cell, basis, qp, Z);
306  fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
307  tmp[Y] * basis_(cell, basis, qp, X);
308  } // end loop over the bases
309  } // end loop over the quadrature points
310  } // end if (numGradDims_ == something)
311  }
312  else if (NUM_FIELD_MULT == 1)
313  {
314  if (numGradDims_ == 1)
315  {
316  for (int qp(0); qp < numQP; ++qp)
317  {
318  tmp[Y] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
319  vector_(cell, qp, Y);
320  tmp[Z] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
321  vector_(cell, qp, Z);
322  for (int basis(0); basis < numBases; ++basis)
323  {
324  fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
325  fields_[Z](cell, basis) += -tmp[Y] * basis_(cell, basis, qp, X);
326  } // end loop over the bases
327  } // end loop over the quadrature points
328  }
329  else if (numGradDims_ == 2)
330  {
331  for (int qp(0); qp < numQP; ++qp)
332  {
333  tmp[X] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
334  vector_(cell, qp, X);
335  tmp[Y] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
336  vector_(cell, qp, Y);
337  tmp[Z] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
338  vector_(cell, qp, Z);
339  for (int basis(0); basis < numBases; ++basis)
340  {
341  fields_[X](cell, basis) += -tmp[Z] * basis_(cell, basis, qp, Y);
342  fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
343  fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
344  tmp[Y] * basis_(cell, basis, qp, X);
345  } // end loop over the bases
346  } // end loop over the quadrature points
347  }
348  else if (numGradDims_ == 3)
349  {
350  for (int qp(0); qp < numQP; ++qp)
351  {
352  tmp[X] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
353  vector_(cell, qp, X);
354  tmp[Y] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
355  vector_(cell, qp, Y);
356  tmp[Z] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
357  vector_(cell, qp, Z);
358  for (int basis(0); basis < numBases; ++basis)
359  {
360  fields_[X](cell, basis) += tmp[Y] * basis_(cell, basis, qp, Z) -
361  tmp[Z] * basis_(cell, basis, qp, Y);
362  fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X) -
363  tmp[X] * basis_(cell, basis, qp, Z);
364  fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
365  tmp[Y] * basis_(cell, basis, qp, X);
366  } // end loop over the bases
367  } // end loop over the quadrature points
368  } // end if (numGradDims_ == something)
369  }
370  else
371  {
372  const int numFieldMults(kokkosFieldMults_.extent(0));
373  if (numGradDims_ == 1)
374  {
375  for (int qp(0); qp < numQP; ++qp)
376  {
377  tmp[Y] = multiplier_ * vector_(cell, qp, Y);
378  tmp[Z] = multiplier_ * vector_(cell, qp, Z);
379  for (int fm(0); fm < numFieldMults; ++fm)
380  {
381  tmp[Y] *= kokkosFieldMults_(fm)(cell, qp);
382  tmp[Z] *= kokkosFieldMults_(fm)(cell, qp);
383  } // end loop over the field multipliers
384  for (int basis(0); basis < numBases; ++basis)
385  {
386  fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
387  fields_[Z](cell, basis) += -tmp[Y] * basis_(cell, basis, qp, X);
388  } // end loop over the bases
389  } // end loop over the quadrature points
390  }
391  else if (numGradDims_ == 2)
392  {
393  for (int qp(0); qp < numQP; ++qp)
394  {
395  tmp[X] = multiplier_ * vector_(cell, qp, X);
396  tmp[Y] = multiplier_ * vector_(cell, qp, Y);
397  tmp[Z] = multiplier_ * vector_(cell, qp, Z);
398  for (int fm(0); fm < numFieldMults; ++fm)
399  {
400  tmp[X] *= kokkosFieldMults_(fm)(cell, qp);
401  tmp[Y] *= kokkosFieldMults_(fm)(cell, qp);
402  tmp[Z] *= kokkosFieldMults_(fm)(cell, qp);
403  } // end loop over the field multipliers
404  for (int basis(0); basis < numBases; ++basis)
405  {
406  fields_[X](cell, basis) += -tmp[Z] * basis_(cell, basis, qp, Y);
407  fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
408  fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
409  tmp[Y] * basis_(cell, basis, qp, X);
410  } // end loop over the bases
411  } // end loop over the quadrature points
412  }
413  else if (numGradDims_ == 3)
414  {
415  for (int qp(0); qp < numQP; ++qp)
416  {
417  tmp[X] = multiplier_ * vector_(cell, qp, X);
418  tmp[Y] = multiplier_ * vector_(cell, qp, Y);
419  tmp[Z] = multiplier_ * vector_(cell, qp, Z);
420  for (int fm(0); fm < numFieldMults; ++fm)
421  {
422  tmp[X] *= kokkosFieldMults_(fm)(cell, qp);
423  tmp[Y] *= kokkosFieldMults_(fm)(cell, qp);
424  tmp[Z] *= kokkosFieldMults_(fm)(cell, qp);
425  } // end loop over the field multipliers
426  for (int basis(0); basis < numBases; ++basis)
427  {
428  fields_[X](cell, basis) += tmp[Y] * basis_(cell, basis, qp, Z) -
429  tmp[Z] * basis_(cell, basis, qp, Y);
430  fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X) -
431  tmp[X] * basis_(cell, basis, qp, Z);
432  fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
433  tmp[Y] * basis_(cell, basis, qp, X);
434  } // end loop over the bases
435  } // end loop over the quadrature points
436  } // end if (numGradDims_ == something)
437  } // end if (NUM_FIELD_MULT == something)
438  } // end of operator()()
439 
441  //
442  // evaluateFields()
443  //
445  template<typename EvalT, typename Traits>
446  void
449  typename Traits::EvalData workset)
450  {
451  using Kokkos::parallel_for;
452  using Kokkos::RangePolicy;
453 
454  // Grab the basis information.
455  basis_ = this->wda(workset).bases[basisIndex_]->weighted_grad_basis;
456 
457  // The following if-block is for the sake of optimization depending on the
458  // number of field multipliers. The parallel_fors will loop over the cells
459  // in the Workset and execute operator()() above.
460  if (fieldMults_.size() == 0)
461  parallel_for(RangePolicy<FieldMultTag<0>>(0, workset.num_cells), *this);
462  else if (fieldMults_.size() == 1)
463  parallel_for(RangePolicy<FieldMultTag<1>>(0, workset.num_cells), *this);
464  else
465  parallel_for(RangePolicy<FieldMultTag<-1>>(0, workset.num_cells), *this);
466  } // end of evaluateFields()
467 
469  //
470  // getValidParameters()
471  //
473  template<typename EvalT, typename TRAITS>
477  {
478  using panzer::BasisIRLayout;
480  using PHX::DataLayout;
481  using std::string;
482  using std::vector;
484  using Teuchos::RCP;
485  using Teuchos::rcp;
486 
487  // Create a ParameterList with all the valid keys we support.
488  RCP<ParameterList> p = rcp(new ParameterList);
489 
490  RCP<const vector<string>> resNames;
491  p->set("Residual Names", resNames);
492  p->set<string>("Vector Name", "?");
493  RCP<BasisIRLayout> basis;
494  p->set("Basis", basis);
496  p->set("IR", ir);
497  p->set<double>("Multiplier", 1.0);
499  p->set("Field Multipliers", fms);
500  RCP<DataLayout> vecDL;
501  p->set("Data Layout Vector", vecDL);
502 
503  return p;
504  } // end of getValidParameters()
505 
506 } // end of namespace panzer
507 
508 #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.
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
const std::vector< std::pair< int, LocalOrdinal > > &pid_and_lid 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.
std::vector< 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...
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