Panzer Version of the Day
Loading...
Searching...
No Matches
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
62
63namespace 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,
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_host_.resize(resNames.size());
141 fields_ = OuterView("Integrator_GradBasisCrossVector::fields_", resNames.size());
142 {
143 int i=0;
144 for (const auto& name : resNames)
145 fields_host_[i++] = MDField<ScalarT, Cell, BASIS>(name, basis.functional);
146 } // end loop over resNames
147
148 for (std::size_t i=0; i< fields_.extent(0); ++i) {
149 const auto& field = fields_host_[i];
151 this->addContributedField(field);
152 else // if (evalStyle_ == EvaluatorStyle::EVALUATES)
153 this->addEvaluatedField(field);
154 }
155
156 // Add the dependent field multipliers, if there are any.
157 int i = 0;
158 fieldMults_.resize(fmNames.size());
159 kokkosFieldMults_ = PHX::View<PHX::UnmanagedView<const ScalarT**>*>("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_host_.size() - 1; ++j)
174 n += resNames[j] + ", ";
175 n += resNames[resNames.size()-1] + "}";
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 :
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 {
204 using Teuchos::ParameterList;
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,
223 {
224 using Kokkos::createDynRankView;
226
227 // Get the PHX::Views of the fields.
228 auto fields_host_mirror_ = Kokkos::create_mirror_view(fields_);
229 for (size_t i=0; i < fields_host_.size(); ++i) {
230 fields_host_mirror_(i) = fields_host_[i].get_static_view();
231 }
232 Kokkos::deep_copy(fields_,fields_host_mirror_);
233
234 // Get the PHX::Views of the field multipliers.
235 auto field_mults_host_mirror_ = Kokkos::create_mirror_view(kokkosFieldMults_);
236 for (size_t i=0; i < fieldMults_.size(); ++i)
237 field_mults_host_mirror_(i) = fieldMults_[i].get_static_view();
238 Kokkos::deep_copy(kokkosFieldMults_,field_mults_host_mirror_);
239
240 // Determine the index in the Workset bases for our particular basis name.
241 basisIndex_ = getBasisIndex(basisName_, (*sd.worksets_)[0], this->wda);
242 } // end of postRegistrationSetup()
243
245 //
246 // operator()()
247 //
249 template<typename EvalT, typename Traits>
250 template<int NUM_FIELD_MULT>
251 KOKKOS_INLINE_FUNCTION
252 void
255 const FieldMultTag<NUM_FIELD_MULT>& /* tag */,
256 const size_t& cell) const
257 {
259
260 // Initialize the evaluated fields.
261 const int numBases(fields_[0].extent(1)), numQP(vector_.extent(1));
262 if (evalStyle_ == EvaluatorStyle::EVALUATES)
263 for (int dim(0); dim < numDims_; ++dim)
264 for (int basis(0); basis < numBases; ++basis)
265 fields_[dim](cell, basis) = 0.0;
266
267 // The following if-block is for the sake of optimization depending on the
268 // number of field multipliers.
269 ScalarT tmp[3];
270 const int X(0), Y(1), Z(2);
271 if (NUM_FIELD_MULT == 0)
272 {
273 if (numGradDims_ == 1)
274 {
275 for (int qp(0); qp < numQP; ++qp)
276 {
277 tmp[Y] = multiplier_ * vector_(cell, qp, Y);
278 tmp[Z] = multiplier_ * vector_(cell, qp, Z);
279 for (int basis(0); basis < numBases; ++basis)
280 {
281 fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
282 fields_[Z](cell, basis) += -tmp[Y] * basis_(cell, basis, qp, X);
283 } // end loop over the bases
284 } // end loop over the quadrature points
285 }
286 else if (numGradDims_ == 2)
287 {
288 for (int qp(0); qp < numQP; ++qp)
289 {
290 tmp[X] = multiplier_ * vector_(cell, qp, X);
291 tmp[Y] = multiplier_ * vector_(cell, qp, Y);
292 tmp[Z] = multiplier_ * vector_(cell, qp, Z);
293 for (int basis(0); basis < numBases; ++basis)
294 {
295 fields_[X](cell, basis) += -tmp[Z] * basis_(cell, basis, qp, Y);
296 fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
297 fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
298 tmp[Y] * basis_(cell, basis, qp, X);
299 } // end loop over the bases
300 } // end loop over the quadrature points
301 }
302 else if (numGradDims_ == 3)
303 {
304 for (int qp(0); qp < numQP; ++qp)
305 {
306 tmp[X] = multiplier_ * vector_(cell, qp, X);
307 tmp[Y] = multiplier_ * vector_(cell, qp, Y);
308 tmp[Z] = multiplier_ * vector_(cell, qp, Z);
309 for (int basis(0); basis < numBases; ++basis)
310 {
311 fields_[X](cell, basis) += tmp[Y] * basis_(cell, basis, qp, Z) -
312 tmp[Z] * basis_(cell, basis, qp, Y);
313 fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X) -
314 tmp[X] * basis_(cell, basis, qp, Z);
315 fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
316 tmp[Y] * basis_(cell, basis, qp, X);
317 } // end loop over the bases
318 } // end loop over the quadrature points
319 } // end if (numGradDims_ == something)
320 }
321 else if (NUM_FIELD_MULT == 1)
322 {
323 if (numGradDims_ == 1)
324 {
325 for (int qp(0); qp < numQP; ++qp)
326 {
327 tmp[Y] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
328 vector_(cell, qp, Y);
329 tmp[Z] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
330 vector_(cell, qp, Z);
331 for (int basis(0); basis < numBases; ++basis)
332 {
333 fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
334 fields_[Z](cell, basis) += -tmp[Y] * basis_(cell, basis, qp, X);
335 } // end loop over the bases
336 } // end loop over the quadrature points
337 }
338 else if (numGradDims_ == 2)
339 {
340 for (int qp(0); qp < numQP; ++qp)
341 {
342 tmp[X] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
343 vector_(cell, qp, X);
344 tmp[Y] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
345 vector_(cell, qp, Y);
346 tmp[Z] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
347 vector_(cell, qp, Z);
348 for (int basis(0); basis < numBases; ++basis)
349 {
350 fields_[X](cell, basis) += -tmp[Z] * basis_(cell, basis, qp, Y);
351 fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
352 fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
353 tmp[Y] * basis_(cell, basis, qp, X);
354 } // end loop over the bases
355 } // end loop over the quadrature points
356 }
357 else if (numGradDims_ == 3)
358 {
359 for (int qp(0); qp < numQP; ++qp)
360 {
361 tmp[X] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
362 vector_(cell, qp, X);
363 tmp[Y] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
364 vector_(cell, qp, Y);
365 tmp[Z] = multiplier_ * kokkosFieldMults_(0)(cell, qp) *
366 vector_(cell, qp, Z);
367 for (int basis(0); basis < numBases; ++basis)
368 {
369 fields_[X](cell, basis) += tmp[Y] * basis_(cell, basis, qp, Z) -
370 tmp[Z] * basis_(cell, basis, qp, Y);
371 fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X) -
372 tmp[X] * basis_(cell, basis, qp, Z);
373 fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
374 tmp[Y] * basis_(cell, basis, qp, X);
375 } // end loop over the bases
376 } // end loop over the quadrature points
377 } // end if (numGradDims_ == something)
378 }
379 else
380 {
381 const int numFieldMults(kokkosFieldMults_.extent(0));
382 if (numGradDims_ == 1)
383 {
384 for (int qp(0); qp < numQP; ++qp)
385 {
386 tmp[Y] = multiplier_ * vector_(cell, qp, Y);
387 tmp[Z] = multiplier_ * vector_(cell, qp, Z);
388 for (int fm(0); fm < numFieldMults; ++fm)
389 {
390 tmp[Y] *= kokkosFieldMults_(fm)(cell, qp);
391 tmp[Z] *= kokkosFieldMults_(fm)(cell, qp);
392 } // end loop over the field multipliers
393 for (int basis(0); basis < numBases; ++basis)
394 {
395 fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
396 fields_[Z](cell, basis) += -tmp[Y] * basis_(cell, basis, qp, X);
397 } // end loop over the bases
398 } // end loop over the quadrature points
399 }
400 else if (numGradDims_ == 2)
401 {
402 for (int qp(0); qp < numQP; ++qp)
403 {
404 tmp[X] = multiplier_ * vector_(cell, qp, X);
405 tmp[Y] = multiplier_ * vector_(cell, qp, Y);
406 tmp[Z] = multiplier_ * vector_(cell, qp, Z);
407 for (int fm(0); fm < numFieldMults; ++fm)
408 {
409 tmp[X] *= kokkosFieldMults_(fm)(cell, qp);
410 tmp[Y] *= kokkosFieldMults_(fm)(cell, qp);
411 tmp[Z] *= kokkosFieldMults_(fm)(cell, qp);
412 } // end loop over the field multipliers
413 for (int basis(0); basis < numBases; ++basis)
414 {
415 fields_[X](cell, basis) += -tmp[Z] * basis_(cell, basis, qp, Y);
416 fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X);
417 fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
418 tmp[Y] * basis_(cell, basis, qp, X);
419 } // end loop over the bases
420 } // end loop over the quadrature points
421 }
422 else if (numGradDims_ == 3)
423 {
424 for (int qp(0); qp < numQP; ++qp)
425 {
426 tmp[X] = multiplier_ * vector_(cell, qp, X);
427 tmp[Y] = multiplier_ * vector_(cell, qp, Y);
428 tmp[Z] = multiplier_ * vector_(cell, qp, Z);
429 for (int fm(0); fm < numFieldMults; ++fm)
430 {
431 tmp[X] *= kokkosFieldMults_(fm)(cell, qp);
432 tmp[Y] *= kokkosFieldMults_(fm)(cell, qp);
433 tmp[Z] *= kokkosFieldMults_(fm)(cell, qp);
434 } // end loop over the field multipliers
435 for (int basis(0); basis < numBases; ++basis)
436 {
437 fields_[X](cell, basis) += tmp[Y] * basis_(cell, basis, qp, Z) -
438 tmp[Z] * basis_(cell, basis, qp, Y);
439 fields_[Y](cell, basis) += tmp[Z] * basis_(cell, basis, qp, X) -
440 tmp[X] * basis_(cell, basis, qp, Z);
441 fields_[Z](cell, basis) += tmp[X] * basis_(cell, basis, qp, Y) -
442 tmp[Y] * basis_(cell, basis, qp, X);
443 } // end loop over the bases
444 } // end loop over the quadrature points
445 } // end if (numGradDims_ == something)
446 } // end if (NUM_FIELD_MULT == something)
447 } // end of operator()()
448
450 //
451 // evaluateFields()
452 //
454 template<typename EvalT, typename Traits>
455 void
458 typename Traits::EvalData workset)
459 {
460 using Kokkos::parallel_for;
461 using Kokkos::RangePolicy;
462
463 // Grab the basis information.
464 basis_ = this->wda(workset).bases[basisIndex_]->weighted_grad_basis;
465
466 // The following if-block is for the sake of optimization depending on the
467 // number of field multipliers. The parallel_fors will loop over the cells
468 // in the Workset and execute operator()() above.
469 if (fieldMults_.size() == 0)
470 parallel_for(RangePolicy<FieldMultTag<0>>(0, workset.num_cells), *this);
471 else if (fieldMults_.size() == 1)
472 parallel_for(RangePolicy<FieldMultTag<1>>(0, workset.num_cells), *this);
473 else
474 parallel_for(RangePolicy<FieldMultTag<-1>>(0, workset.num_cells), *this);
475 } // end of evaluateFields()
476
478 //
479 // getValidParameters()
480 //
482 template<typename EvalT, typename TRAITS>
483 Teuchos::RCP<Teuchos::ParameterList>
486 {
489 using PHX::DataLayout;
490 using std::string;
491 using std::vector;
492 using Teuchos::ParameterList;
493 using Teuchos::RCP;
494 using Teuchos::rcp;
495
496 // Create a ParameterList with all the valid keys we support.
497 RCP<ParameterList> p = rcp(new ParameterList);
498
499 RCP<const vector<string>> resNames;
500 p->set("Residual Names", resNames);
501 p->set<string>("Vector Name", "?");
502 RCP<BasisIRLayout> basis;
503 p->set("Basis", basis);
504 RCP<IntegrationRule> ir;
505 p->set("IR", ir);
506 p->set<double>("Multiplier", 1.0);
507 RCP<const vector<string>> fms;
508 p->set("Field Multipliers", fms);
509 RCP<DataLayout> vecDL;
510 p->set("Data Layout Vector", vecDL);
511
512 return p;
513 } // end of getValidParameters()
514
515} // end of namespace panzer
516
517#endif // PANZER_EVALUATOR_GRADBASISCROSSVECTOR_IMPL_HPP
panzer::EvaluatorStyle evalStyle
The EvaluatorStyle of the parent Integrator_CurlBasisDotVector object.
double multiplier
The scalar multiplier out in front of the integral ( ).
PHX::MDField< ScalarT, panzer::Cell, panzer::BASIS > field
A field to which we'll contribute, or in which we'll store, the result of computing this integral.
Teuchos::RCP< PHX::DataLayout > functional
<Cell,Basis>
Teuchos::RCP< const PureBasis > getBasis() const
PHX::MDField< const ScalarT, Cell, IP, Dim > vector_
A field representing the vector-valued function we're integrating ( ).
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
Get Valid Parameters.
std::vector< PHX::MDField< const ScalarT, Cell, IP > > fieldMults_
The (possibly empty) list of fields that are multipliers out in front of the integral ( ,...
void evaluateFields(typename Traits::EvalData d)
Evaluate Fields.
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
Post-Registration Setup.
PHX::View< PHX::UnmanagedView< const ScalarT ** > * > kokkosFieldMults_
The PHX::View representation of the (possibly empty) list of fields that are multipliers out in front...
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.
const panzer::EvaluatorStyle evalStyle_
An enum determining the behavior of this Evaluator.
int numDims_
The number of dimensions associated with the vector.
KOKKOS_INLINE_FUNCTION void operator()(const FieldMultTag< NUM_FIELD_MULT > &tag, const std::size_t &cell) const
Perform the integration.
std::vector< PHX::MDField< ScalarT, Cell, BASIS > > fields_host_
The fields to which we'll contribute, or in which we'll store, the result of computing this integral.
int numGradDims_
The number of dimensions associated with the gradient.
Teuchos::RCP< PHX::DataLayout > dl_vector
Data layout for vector fields.
Teuchos::RCP< PHX::DataLayout > dl_scalar
Data layout for scalar fields.
int num_cells
DEPRECATED - use: numCells()
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.
EvaluatorStyle
An indication of how an Evaluator will behave.
This empty struct allows us to optimize operator()() depending on the number of field multipliers.
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_