11 #ifndef PANZER_INTEGRATOR_BASISTIMESSCALAR_HPP
12 #define PANZER_INTEGRATOR_BASISTIMESSCALAR_HPP
24 #include "Kokkos_DynRankView.hpp"
31 #include "Phalanx_Evaluator_Derived.hpp"
32 #include "Phalanx_MDField.hpp"
47 template<
typename EvalT,
typename Traits>
91 const std::string& resName,
92 const std::string& valName,
96 const std::vector<std::string>& fmNames =
97 std::vector<std::string>());
177 template<
int NUM_FIELD_MULT>
197 template<
int NUM_FIELD_MULT>
198 KOKKOS_INLINE_FUNCTION
202 const std::size_t& cell)
const;
256 std::vector<PHX::MDField<const ScalarT, panzer::Cell, panzer::IP>>
fieldMults_;
263 PHX::View<Kokkos::View<const ScalarT**, typename PHX::DevLayout<ScalarT>::type, Kokkos::MemoryUnmanaged>*>
kokkosFieldMults_;
295 #endif // PANZER_INTEGRATOR_BASISTIMESSCALAR_HPP
PHX::View< ScalarT * > tmp2_
For storing temporaries, one value per thread.
int numQP_
The number of quadrature points for each cell.
const panzer::EvaluatorStyle evalStyle_
An enum determining the behavior of this Evaluator.
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...
std::string basisName_
The name of the basis we're using.
This empty struct allows us to optimize operator()() depending on the number of field multipliers...
PHX::View< ScalarT * > tmp_
For storing temporaries, one value per thread.
Integrator_BasisTimesScalar(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.
EvaluatorStyle
An indication of how an Evaluator will behave.
KOKKOS_INLINE_FUNCTION void operator()(const FieldMultTag< NUM_FIELD_MULT > &tag, const std::size_t &cell) const
Perform the integration.
double multiplier
The scalar multiplier out in front of the integral ( ).
Wrapper to PHX::EvaluatorWithBaseImpl that implements Panzer-specific helpers.
panzer::EvaluatorStyle evalStyle
The EvaluatorStyle of the parent Integrator_CurlBasisDotVector object.
PHX::MDField< const ScalarT, panzer::Cell, panzer::IP > scalar_
A field representing the scalar function we're integrating ( ).
typename EvalT::ScalarT ScalarT
The scalar type.
void evaluateFields(typename Traits::EvalData workset)
Evaluate Fields.
std::size_t basisIndex_
The index in the Workset bases for our particular BasisIRLayout name.
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...
double multiplier_
The scalar multiplier out in front of the integral ( ).
void postRegistrationSetup(typename Traits::SetupData sd, PHX::FieldManager< Traits > &fm)
Post-Registration Setup.
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 ( ...
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
Get Valid Parameters.
PHX::MDField< double, panzer::Cell, panzer::BASIS, panzer::IP > basis_
The scalar basis information necessary for integration.