43#ifndef __Panzer_Integerator_BasisTimesTensorTimesVector_decl_hpp__
44#define __Panzer_Integerator_BasisTimesTensorTimesVector_decl_hpp__
56#include "Kokkos_DynRankView.hpp"
63#include "Phalanx_Evaluator_Derived.hpp"
64#include "Phalanx_MDField.hpp"
80 template<
typename EvalT,
typename Traits>
84 public PHX::EvaluatorDerived<EvalT, Traits>
119 const std::string& resName,
120 const std::string& valName,
123 const std::string& tensorName );
163 const Teuchos::ParameterList& p);
197 const PHX::FieldTag& resTag,
198 const PHX::FieldTag& valTag,
201 const PHX::FieldTag& tensorTag);
244 KOKKOS_INLINE_FUNCTION
245 void operator()(
const std::size_t& cell)
const;
259 Teuchos::RCP<Teuchos::ParameterList>
297 PHX::MDField<ScalarT, panzer::Cell, panzer::BASIS>
field_;
303 PHX::MDField<const ScalarT, panzer::Cell, panzer::IP, panzer::Dim>
309 PHX::MDField<const ScalarT, panzer::Cell, panzer::IP, panzer::Dim, panzer::Dim>
panzer::EvaluatorStyle evalStyle
The EvaluatorStyle of the parent Integrator_CurlBasisDotVector object.
Wrapper to PHX::EvaluatorWithBaseImpl that implements Panzer-specific helpers.
IntegrationDescriptor id_
The IntegrationDescriptor for the quadrature to use.
void postRegistrationSetup(typename Traits::SetupData sd, PHX::FieldManager< Traits > &fm)
Post-Registration Setup.
int numDim_
The dimensionality of our vector-valued fields.
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
Get Valid Parameters.
PHX::View< ScalarT * > tmp_
Scratch space for caching temporary values in the kokkos kernel.
std::size_t basisIndex_
The index in the Workset bases for our particular BasisIRLayout name.
typename EvalT::ScalarT ScalarT
The scalar type.
KOKKOS_INLINE_FUNCTION void operator()(const std::size_t &cell) const
Perform the integration.
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.
PHX::View< const ScalarT **** > kokkosTensor_
The PHX::View representation of the tensor fields that are multiplied out in front of the integral ( ...
void evaluateFields(typename Traits::EvalData workset)
Evaluate Fields.
std::string basisName_
The name of the basis we're using.
bool useDescriptors_
A flag indicating whether or not to use the descriptor interface.
int numQP_
The number of quadrature points for each cell.
PHX::MDField< const double, panzer::Cell, panzer::BASIS, panzer::IP, panzer::Dim > basis_
The vector basis information necessary for integration.
PHX::MDField< const ScalarT, panzer::Cell, panzer::IP, panzer::Dim > vector_
A field representing the vector-valued function we're integrating ( ).
BasisDescriptor bd_
The BasisDescriptor for the basis to use.
PHX::MDField< const ScalarT, panzer::Cell, panzer::IP, panzer::Dim, panzer::Dim > tensor_
The tensor field( ).
EvaluatorStyle
An indication of how an Evaluator will behave.