Panzer  Version of the Day
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 PHX::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;
121  TEUCHOS_TEST_FOR_EXCEPTION(
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.")
125  TEUCHOS_TEST_FOR_EXCEPTION(numDims_ !=
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())
130  TEUCHOS_TEST_FOR_EXCEPTION(numGradDims_ > numDims_, logic_error,
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_ = PHX::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());
157  kokkosFieldMults_ = View<View<const ScalarT**>*>("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_.extent(0) - 1; ++j)
172  n += fields_[j].fieldTag().name() + ", ";
173  n += fields_(fields_.extent(0)-1).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  :
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  {
202  using Teuchos::ParameterList;
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 
225  // Get the PHX::Views of the field multipliers.
226  for (size_t i(0); i < fieldMults_.size(); ++i)
227  kokkosFieldMults_(i) = fieldMults_[i].get_static_view();
228 
229  // Determine the index in the Workset bases for our particular basis name.
230  basisIndex_ = getBasisIndex(basisName_, (*sd.worksets_)[0], this->wda);
231  } // end of postRegistrationSetup()
232 
234  //
235  // operator()()
236  //
238  template<typename EvalT, typename Traits>
239  template<int NUM_FIELD_MULT>
240  KOKKOS_INLINE_FUNCTION
241  void
244  const FieldMultTag<NUM_FIELD_MULT>& /* tag */,
245  const size_t& cell) const
246  {
248 
249  // Initialize the evaluated fields.
250  const int numBases(fields_[0].extent(1)), numQP(vector_.extent(1));
251  if (evalStyle_ == EvaluatorStyle::EVALUATES)
252  for (int dim(0); dim < numDims_; ++dim)
253  for (int basis(0); basis < numBases; ++basis)
254  fields_[dim](cell, basis) = 0.0;
255 
256  // The following if-block is for the sake of optimization depending on the
257  // number of field multipliers.
258  ScalarT tmp[3];
259  const int X(0), Y(1), Z(2);
260  if (NUM_FIELD_MULT == 0)
261  {
262  if (numGradDims_ == 1)
263  {
264  for (int qp(0); qp < numQP; ++qp)
265  {
266  tmp[Y] = multiplier_ * vector_(cell, qp, Y);
267  tmp[Z] = multiplier_ * vector_(cell, qp, Z);
268  for (int basis(0); basis < numBases; ++basis)
269  {
270  fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
271  fields_[Z](cell, basis) += -tmp[Y] * basis_(cell, basis, qp, X);
272  } // end loop over the bases
273  } // end loop over the quadrature points
274  }
275  else if (numGradDims_ == 2)
276  {
277  for (int qp(0); qp < numQP; ++qp)
278  {
279  tmp[X] = multiplier_ * vector_(cell, qp, X);
280  tmp[Y] = multiplier_ * vector_(cell, qp, Y);
281  tmp[Z] = multiplier_ * vector_(cell, qp, Z);
282  for (int basis(0); basis < numBases; ++basis)
283  {
284  fields_[X](cell, basis) += -tmp[Z] * basis_(cell, basis, qp, Y);
285  fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
286  fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
287  tmp[Y] * basis_(cell, basis, qp, X);
288  } // end loop over the bases
289  } // end loop over the quadrature points
290  }
291  else if (numGradDims_ == 3)
292  {
293  for (int qp(0); qp < numQP; ++qp)
294  {
295  tmp[X] = multiplier_ * vector_(cell, qp, X);
296  tmp[Y] = multiplier_ * vector_(cell, qp, Y);
297  tmp[Z] = multiplier_ * vector_(cell, qp, Z);
298  for (int basis(0); basis < numBases; ++basis)
299  {
300  fields_[X](cell, basis) += tmp[Y] * basis_(cell, basis, qp, Z) -
301  tmp[Z] * basis_(cell, basis, qp, Y);
302  fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X) -
303  tmp[X] * basis_(cell, basis, qp, Z);
304  fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
305  tmp[Y] * basis_(cell, basis, qp, X);
306  } // end loop over the bases
307  } // end loop over the quadrature points
308  } // end if (numGradDims_ == something)
309  }
310  else if (NUM_FIELD_MULT == 1)
311  {
312  if (numGradDims_ == 1)
313  {
314  for (int qp(0); qp < numQP; ++qp)
315  {
316  tmp[Y] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
317  vector_(cell, qp, Y);
318  tmp[Z] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
319  vector_(cell, qp, Z);
320  for (int basis(0); basis < numBases; ++basis)
321  {
322  fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
323  fields_[Z](cell, basis) += -tmp[Y] * basis_(cell, basis, qp, X);
324  } // end loop over the bases
325  } // end loop over the quadrature points
326  }
327  else if (numGradDims_ == 2)
328  {
329  for (int qp(0); qp < numQP; ++qp)
330  {
331  tmp[X] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
332  vector_(cell, qp, X);
333  tmp[Y] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
334  vector_(cell, qp, Y);
335  tmp[Z] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
336  vector_(cell, qp, Z);
337  for (int basis(0); basis < numBases; ++basis)
338  {
339  fields_[X](cell, basis) += -tmp[Z] * basis_(cell, basis, qp, Y);
340  fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
341  fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
342  tmp[Y] * basis_(cell, basis, qp, X);
343  } // end loop over the bases
344  } // end loop over the quadrature points
345  }
346  else if (numGradDims_ == 3)
347  {
348  for (int qp(0); qp < numQP; ++qp)
349  {
350  tmp[X] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
351  vector_(cell, qp, X);
352  tmp[Y] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
353  vector_(cell, qp, Y);
354  tmp[Z] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
355  vector_(cell, qp, Z);
356  for (int basis(0); basis < numBases; ++basis)
357  {
358  fields_[X](cell, basis) += tmp[Y] * basis_(cell, basis, qp, Z) -
359  tmp[Z] * basis_(cell, basis, qp, Y);
360  fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X) -
361  tmp[X] * basis_(cell, basis, qp, Z);
362  fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
363  tmp[Y] * basis_(cell, basis, qp, X);
364  } // end loop over the bases
365  } // end loop over the quadrature points
366  } // end if (numGradDims_ == something)
367  }
368  else
369  {
370  const int numFieldMults(kokkosFieldMults_.extent(0));
371  if (numGradDims_ == 1)
372  {
373  for (int qp(0); qp < numQP; ++qp)
374  {
375  tmp[Y] = multiplier_ * vector_(cell, qp, Y);
376  tmp[Z] = multiplier_ * vector_(cell, qp, Z);
377  for (int fm(0); fm < numFieldMults; ++fm)
378  {
379  tmp[Y] *= kokkosFieldMults_(fm)(cell, qp);
380  tmp[Z] *= kokkosFieldMults_(fm)(cell, qp);
381  } // end loop over the field multipliers
382  for (int basis(0); basis < numBases; ++basis)
383  {
384  fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
385  fields_[Z](cell, basis) += -tmp[Y] * basis_(cell, basis, qp, X);
386  } // end loop over the bases
387  } // end loop over the quadrature points
388  }
389  else if (numGradDims_ == 2)
390  {
391  for (int qp(0); qp < numQP; ++qp)
392  {
393  tmp[X] = multiplier_ * vector_(cell, qp, X);
394  tmp[Y] = multiplier_ * vector_(cell, qp, Y);
395  tmp[Z] = multiplier_ * vector_(cell, qp, Z);
396  for (int fm(0); fm < numFieldMults; ++fm)
397  {
398  tmp[X] *= kokkosFieldMults_(fm)(cell, qp);
399  tmp[Y] *= kokkosFieldMults_(fm)(cell, qp);
400  tmp[Z] *= kokkosFieldMults_(fm)(cell, qp);
401  } // end loop over the field multipliers
402  for (int basis(0); basis < numBases; ++basis)
403  {
404  fields_[X](cell, basis) += -tmp[Z] * basis_(cell, basis, qp, Y);
405  fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
406  fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
407  tmp[Y] * basis_(cell, basis, qp, X);
408  } // end loop over the bases
409  } // end loop over the quadrature points
410  }
411  else if (numGradDims_ == 3)
412  {
413  for (int qp(0); qp < numQP; ++qp)
414  {
415  tmp[X] = multiplier_ * vector_(cell, qp, X);
416  tmp[Y] = multiplier_ * vector_(cell, qp, Y);
417  tmp[Z] = multiplier_ * vector_(cell, qp, Z);
418  for (int fm(0); fm < numFieldMults; ++fm)
419  {
420  tmp[X] *= kokkosFieldMults_(fm)(cell, qp);
421  tmp[Y] *= kokkosFieldMults_(fm)(cell, qp);
422  tmp[Z] *= kokkosFieldMults_(fm)(cell, qp);
423  } // end loop over the field multipliers
424  for (int basis(0); basis < numBases; ++basis)
425  {
426  fields_[X](cell, basis) += tmp[Y] * basis_(cell, basis, qp, Z) -
427  tmp[Z] * basis_(cell, basis, qp, Y);
428  fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X) -
429  tmp[X] * basis_(cell, basis, qp, Z);
430  fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
431  tmp[Y] * basis_(cell, basis, qp, X);
432  } // end loop over the bases
433  } // end loop over the quadrature points
434  } // end if (numGradDims_ == something)
435  } // end if (NUM_FIELD_MULT == something)
436  } // end of operator()()
437 
439  //
440  // evaluateFields()
441  //
443  template<typename EvalT, typename Traits>
444  void
447  typename Traits::EvalData workset)
448  {
449  using Kokkos::parallel_for;
450  using Kokkos::RangePolicy;
451 
452  // Grab the basis information.
453  basis_ = this->wda(workset).bases[basisIndex_]->weighted_grad_basis;
454 
455  // The following if-block is for the sake of optimization depending on the
456  // number of field multipliers. The parallel_fors will loop over the cells
457  // in the Workset and execute operator()() above.
458  if (fieldMults_.size() == 0)
459  parallel_for(RangePolicy<FieldMultTag<0>>(0, workset.num_cells), *this);
460  else if (fieldMults_.size() == 1)
461  parallel_for(RangePolicy<FieldMultTag<1>>(0, workset.num_cells), *this);
462  else
463  parallel_for(RangePolicy<FieldMultTag<-1>>(0, workset.num_cells), *this);
464  } // end of evaluateFields()
465 
467  //
468  // getValidParameters()
469  //
471  template<typename EvalT, typename TRAITS>
472  Teuchos::RCP<Teuchos::ParameterList>
475  {
476  using panzer::BasisIRLayout;
478  using PHX::DataLayout;
479  using std::string;
480  using std::vector;
481  using Teuchos::ParameterList;
482  using Teuchos::RCP;
483  using Teuchos::rcp;
484 
485  // Create a ParameterList with all the valid keys we support.
486  RCP<ParameterList> p = rcp(new ParameterList);
487 
488  RCP<const vector<string>> resNames;
489  p->set("Residual Names", resNames);
490  p->set<string>("Vector Name", "?");
491  RCP<BasisIRLayout> basis;
492  p->set("Basis", basis);
493  RCP<IntegrationRule> ir;
494  p->set("IR", ir);
495  p->set<double>("Multiplier", 1.0);
496  RCP<const vector<string>> fms;
497  p->set("Field Multipliers", fms);
498  RCP<DataLayout> vecDL;
499  p->set("Data Layout Vector", vecDL);
500 
501  return p;
502  } // end of getValidParameters()
503 
504 } // end of namespace panzer
505 
506 #endif // PANZER_EVALUATOR_GRADBASISCROSSVECTOR_IMPL_HPP
int num_cells
DEPRECATED - use: numCells()
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.
This empty struct allows us to optimize operator()() depending on the number of field multipliers...
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
Get Valid Parameters.
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.
void evaluateFields(typename Traits::EvalData d)
Evaluate Fields.
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 ( ).
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.
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< 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::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...
KOKKOS_INLINE_FUNCTION void operator()(const FieldMultTag< NUM_FIELD_MULT > &tag, const std::size_t &cell) const
Perform the integration.
Teuchos::RCP< PHX::DataLayout > dl_vector
Data layout for vector fields.
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.
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
<Cell,Basis>
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_