43#ifndef __Panzer_Integerator_BasisTimesVector_impl_hpp__
44#define __Panzer_Integerator_BasisTimesVector_impl_hpp__
64 template<
typename EvalT,
typename Traits>
68 const std::string& resName,
69 const std::string& valName,
72 const double& multiplier,
73 const std::vector<std::string>& fmNames
89 using std::invalid_argument;
90 using std::logic_error;
95 TEUCHOS_TEST_FOR_EXCEPTION(resName ==
"", invalid_argument,
"Error: " \
96 "Integrator_BasisTimesVector called with an empty residual name.")
97 TEUCHOS_TEST_FOR_EXCEPTION(valName ==
"", invalid_argument,
"Error: " \
98 "Integrator_BasisTimesVector called with an empty value name.")
99 RCP<const PureBasis> tmpBasis = basis.
getBasis();
100 TEUCHOS_TEST_FOR_EXCEPTION(not tmpBasis->isVectorBasis(), logic_error,
101 "Error: Integrator_BasisTimesVector: Basis of type \""
102 << tmpBasis->name() <<
"\" is not a vector basis.");
103 TEUCHOS_TEST_FOR_EXCEPTION(not tmpBasis->requiresOrientations(),
104 logic_error,
"Integrator_BasisTimesVector: Basis of type \""
105 << tmpBasis->name() <<
"\" does not require orientations. This seems " \
106 "very strange, so I'm failing.");
110 this->addDependentField(
vector_);
116 this->addContributedField(
field_);
118 this->addEvaluatedField(
field_);
126 for (
const auto& name : fmNames)
133 string n(
"Integrator_BasisTimesVector (");
138 n +=
", " + print<EvalT>() +
"): " +
field_.fieldTag().name();
147 template<
typename EvalT,
typename Traits>
151 const PHX::FieldTag& resTag,
152 const PHX::FieldTag& valTag,
155 const double& multiplier,
156 const std::vector<PHX::FieldTag>& multipliers
173 using std::invalid_argument;
174 using std::logic_error;
179 TEUCHOS_TEST_FOR_EXCEPTION(not ((
bd_.getType() ==
"HCurl") or
180 (
bd_.getType() ==
"HDiv")), logic_error,
"Error: " \
181 "Integrator_BasisTimesVector: Basis of type \"" <<
bd_.getType()
182 <<
"\" is not a vector basis.")
186 this->addDependentField(
vector_);
192 this->addContributedField(
field_);
194 this->addEvaluatedField(
field_);
202 for (
const auto& fm : multipliers)
209 string n(
"Integrator_BasisTimesVector (");
214 n +=
", " + print<EvalT>() +
"): " +
field_.fieldTag().name();
223 template<
typename EvalT,
typename Traits>
226 const Teuchos::ParameterList& p)
230 p.get<
std::string>(
"Residual Name"),
231 p.get<
std::string>(
"Value Name"),
234 p.get<double>(
"Multiplier"),
236 (
"Field Multipliers") ?
238 (
"Field Multipliers")) :
std::vector<
std::string>())
240 using Teuchos::ParameterList;
245 p.validateParameters(*validParams);
253 template<
typename EvalT,
typename Traits>
266 kokkosFieldMults_h(i) =
fieldMults_[i].get_static_view();
279 if (Sacado::IsADType<ScalarT>::value) {
280 const auto fadSize = Kokkos::dimension_scalar(
field_.get_static_view());
281 tmp_ = PHX::View<ScalarT*>(
"panzer::Integrator::BasisTimesVector::tmp_",
field_.extent(0),fadSize);
283 tmp_ = PHX::View<ScalarT*>(
"panzer::Integrator::BasisTimesVector::tmp_",
field_.extent(0));
293 template<
typename EvalT,
typename Traits>
294 template<
int NUM_FIELD_MULT>
295 KOKKOS_INLINE_FUNCTION
300 const size_t& cell)
const
305 const int numBases(
basis_.extent(1));
307 for (
int basis(0); basis < numBases; ++basis)
308 field_(cell, basis) = 0.0;
312 if (NUM_FIELD_MULT == 0)
317 for (
int qp(0); qp <
numQP_; ++qp)
319 for (
int dim(0); dim <
numDim_; ++dim)
321 for (
int basis(0); basis < numBases; ++basis)
326 else if (NUM_FIELD_MULT == 1)
331 for (
int qp(0); qp <
numQP_; ++qp)
333 for (
int dim(0); dim <
numDim_; ++dim)
335 for (
int basis(0); basis < numBases; ++basis)
348 for (
int qp(0); qp <
numQP_; ++qp)
351 for (
int fm(0); fm < numFieldMults; ++fm)
353 for (
int dim(0); dim <
numDim_; ++dim)
355 for (
int basis(0); basis < numBases; ++basis)
367 template<
typename EvalT,
typename Traits>
373 using Kokkos::parallel_for;
374 using Kokkos::RangePolicy;
378 this->
wda(workset).getBasisValues(
bd_,
id_) :
399 template<
typename EvalT,
typename Traits>
400 Teuchos::RCP<Teuchos::ParameterList>
408 using Teuchos::ParameterList;
413 RCP<ParameterList> p = rcp(
new ParameterList);
414 p->set<
string>(
"Residual Name",
"?");
415 p->set<
string>(
"Value Name",
"?");
416 RCP<BasisIRLayout> basis;
417 p->set(
"Basis", basis);
418 RCP<IntegrationRule> ir;
420 p->set<
double>(
"Multiplier", 1.0);
421 RCP<const vector<string> > fms;
422 p->set(
"Field Multipliers", fms);
Teuchos::RCP< PHX::DataLayout > functional
<Cell,Basis>
Teuchos::RCP< const PureBasis > getBasis() const
PHX::MDField< const Scalar, Cell, BASIS, IP, Dim > ConstArray_CellBasisIPDim
Array_CellBasisIPDim weighted_basis_vector
ConstArray_CellBasisIPDim getVectorBasisValues(const bool weighted, const bool cache=true, const bool force=false) const
Get the vector basis values evaluated at mesh points.
WorksetDetailsAccessor wda
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.
Integrator_BasisTimesVector(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.
PHX::MDField< const ScalarT, panzer::Cell, panzer::IP, panzer::Dim > vector_
A field representing the vector-valued function we're integrating ( ).
bool useDescriptors_
A flag indicating whether or not to use the descriptor interface.
BasisDescriptor bd_
The BasisDescriptor for the basis to use.
PHX::MDField< const double, panzer::Cell, panzer::BASIS, panzer::IP, panzer::Dim > basis_
The vector basis information necessary for integration.
PHX::View< ScalarT * > tmp_
Temporary for caching field multipliers.
void postRegistrationSetup(typename Traits::SetupData sd, PHX::FieldManager< Traits > &fm)
Post-Registration Setup.
std::size_t basisIndex_
The index in the Workset bases for our particular BasisIRLayout name.
KOKKOS_INLINE_FUNCTION void operator()(const FieldMultTag< NUM_FIELD_MULT > &tag, const std::size_t &cell) const
Perform the integration.
const panzer::EvaluatorStyle evalStyle_
An enum determining the behavior of this Evaluator.
std::vector< PHX::MDField< const ScalarT, panzer::Cell, panzer::IP > > fieldMults_
The (possibly empty) list of fields that are multipliers out in front of the integral ( ,...
PHX::View< Kokkos::View< const ScalarT **, typename PHX::DevLayout< ScalarT >::type, Kokkos::MemoryUnmanaged > * > kokkosFieldMults_
The PHX::View representation of the (possibly empty) list of fields that are multipliers out in front...
int numQP_
The number of quadrature points for each cell.
std::string basisName_
The name of the basis we're using.
IntegrationDescriptor id_
The IntegrationDescriptor for the quadrature to use.
int numDim_
The dimensionality of our vector-valued fields.
double multiplier_
The scalar multiplier out in front of the integral ( ).
void evaluateFields(typename Traits::EvalData workset)
Evaluate Fields.
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
Get Valid Parameters.
Teuchos::RCP< PHX::DataLayout > dl_vector
Data layout for vector fields.
Teuchos::RCP< PHX::DataLayout > dl_scalar
Data layout for scalar fields.
Description and data layouts associated with a particular basis.
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_
const panzer::Workset & EvalData