43 #ifndef __Panzer_Integrator_GradBasisDotTensorTimesVector_impl_hpp__
44 #define __Panzer_Integrator_GradBasisDotTensorTimesVector_impl_hpp__
65 template<
typename EvalT,
typename Traits>
69 const std::string& resName,
70 const std::string& fluxName,
73 const std::string& tensorName,
76 evalStyle_(evalStyle),
77 basisName_(basis.name())
86 using std::invalid_argument;
87 using std::logic_error;
93 "Integrator_GradBasisDotTensorTimesVector called with an empty residual name.")
95 "Integrator_GradBasisDotTensorTimesVector called with an empty flux name.")
98 "Error: Integrator_GradBasisDotTensorTimesVector: Basis of type \""
99 << tmpBasis->name() <<
"\" does not support the gradient operator.")
105 tmpVecDL->extent(2) < ir.
dl_vector->extent(2), logic_error,
106 "Integrator_GradBasisDotTensorTimesVector: Dimension of space exceeds " \
107 "dimension of Vector Data Layout.");
111 vector_ = MDField<const ScalarT, Cell, IP, Dim>(fluxName, tmpVecDL);
112 this->addDependentField(
vector_);
118 this->addContributedField(
field_);
120 this->addEvaluatedField(
field_);
123 tensor_ = MDField<const ScalarT, Cell, IP, Dim, Dim>(tensorName, ir.
dl_tensor);
124 this->addDependentField(
tensor_);
127 string n(
"Integrator_GradBasisDotTensorTimesVector (");
132 n +=
"): " +
field_.fieldTag().name();
141 template<
typename EvalT,
typename Traits>
148 p.get<std::string>(
"Residual Name"),
149 p.get<std::string>(
"Flux Name"),
152 p.get<std::string>(
"Tensor Name"),
153 p.isType<Teuchos::
RCP<PHX::DataLayout>>(
"Vector Data Layout") ?
154 p.get<Teuchos::
RCP<PHX::DataLayout>>(
"Vector Data Layout") :
170 template<
typename EvalT,
typename Traits>
181 kokkosTensor_ = tensor_.get_static_view();
188 if (!use_shared_memory) {
189 if (Sacado::IsADType<ScalarT>::value) {
190 const auto fadSize = Kokkos::dimension_scalar(field_.get_view());
191 tmp_ = PHX::View<ScalarT*>(
"GradBasisDotTensorTimesVector::tmp_",field_.extent(0),fadSize);
193 tmp_ = PHX::View<ScalarT*>(
"GradBasisDotTensorTimesVector::tmp_",field_.extent(0));
203 template<
typename EvalT,
typename Traits>
204 KOKKOS_INLINE_FUNCTION
209 const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const
212 const int cell = team.league_rank();
215 const int numQP(vector_.extent(1)), numDim(vector_.extent(2)),
216 numBases(basis_.extent(1));
218 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (
const int basis) {
219 field_(cell, basis) = 0.0;
226 for (
int qp(0); qp < numQP; ++qp)
228 for (
int dim(0); dim < numDim; ++dim)
231 for (
int dim2(0); dim2 < numDim; ++dim2)
232 tmp_(cell) += kokkosTensor_(cell, qp, dim, dim2) * vector_(cell, qp, dim2);
233 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (
const int basis) {
234 field_(cell, basis) += basis_(cell, basis, qp, dim) * tmp_(cell);
246 template<
typename EvalT,
typename Traits>
247 KOKKOS_INLINE_FUNCTION
252 const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const
255 const int cell = team.league_rank();
256 const int numQP = vector_.extent(1);
257 const int numDim = vector_.extent(2);
258 const int numBases = basis_.extent(1);
259 const int fadSize = Kokkos::dimension_scalar(field_.get_view());
263 if (Sacado::IsADType<ScalarT>::value) {
265 tmp_field =
scratch_view(team.team_shmem(),numBases,fadSize);
273 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (
const int basis) {
274 tmp_field(basis) = 0.0;
281 for (
int qp(0); qp < numQP; ++qp)
283 for (
int dim(0); dim < numDim; ++dim)
286 for (
int dim2(0); dim2 < numDim; ++dim2)
287 tmp_(cell) += kokkosTensor_(cell, qp, dim, dim2) * vector_(cell, qp, dim2);
288 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases), [&] (
const int basis) {
289 tmp_field(basis) += basis_(cell, basis, qp, dim) * tmp_(cell);
297 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),[&] (
const int basis) {
298 field_(cell,basis) = tmp_field(basis);
302 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,numBases),[&] (
const int basis) {
303 field_(cell,basis) += tmp_field(basis);
314 template<
typename EvalT,
typename Traits>
320 using Kokkos::parallel_for;
321 using Kokkos::TeamPolicy;
324 basis_ = this->wda(workset).bases[basisIndex_]->weighted_grad_basis;
327 if (use_shared_memory) {
329 if (Sacado::IsADType<ScalarT>::value) {
330 const int fadSize = Kokkos::dimension_scalar(field_.get_view());
331 bytes = scratch_view::shmem_size(1,fadSize) + scratch_view::shmem_size(basis_.extent(1),fadSize);
334 bytes = scratch_view::shmem_size(1) + scratch_view::shmem_size(basis_.extent(1));
337 parallel_for(this->getName(), policy, *
this);
341 parallel_for(this->getName(), policy, *
this);
350 template<
typename EvalT,
typename TRAITS>
366 p->set<
string>(
"Residual Name",
"?");
367 p->set<
string>(
"Flux Name",
"?");
369 p->set(
"Basis", basis);
372 p->set<std::string>(
"Tensor Name",
"?");
374 p->set(
"Vector Data Layout", vecDL);
380 #endif // __Panzer_Integrator_GradBasisDotTensorTimesVector_impl_hpp__
int num_cells
DEPRECATED - use: numCells()
This empty struct allows us to optimize operator()() depending on the number of field multipliers...
Kokkos::TeamPolicy< TeamPolicyProperties...> teamPolicy(const int &league_size)
Returns a TeamPolicy for hierarchic parallelism.
Teuchos::RCP< PHX::DataLayout > dl_tensor
Data layout for rank-2 tensor fields.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
KOKKOS_INLINE_FUNCTION void operator()(const FieldMultTag &tag, const Kokkos::TeamPolicy< PHX::exec_space >::member_type &team) const
Perform the integration.
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.
typename EvalT::ScalarT ScalarT
The scalar type.
Teuchos::RCP< const PureBasis > getBasis() const
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
Post-Registration Setup.
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
Get Valid Parameters.
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_GradBasisDotTensorTimesVector(const panzer::EvaluatorStyle &evalStyle, const std::string &resName, const std::string &fluxName, const panzer::BasisIRLayout &basis, const panzer::IntegrationRule &ir, const std::string &tensorName, const Teuchos::RCP< PHX::DataLayout > &vecDL=Teuchos::null)
Main Constructor.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
panzer::EvaluatorStyle evalStyle
The EvaluatorStyle of the parent Integrator_CurlBasisDotVector object.
bool useSharedMemory() const
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::MDField< const ScalarT, panzer::Cell, panzer::IP, panzer::Dim > vector_
A field representing the vector-valued function we're integrating ( ).
void validateParameters(ParameterList const &validParamList, int const depth=1000, EValidateUsed const validateUsed=VALIDATE_USED_ENABLED, EValidateDefaults const validateDefaults=VALIDATE_DEFAULTS_ENABLED) const
const panzer::EvaluatorStyle evalStyle_
An enum determining the behavior of this Evaluator.
This empty struct allows us to optimize operator()() depending on the number of field multipliers...
void evaluateFields(typename Traits::EvalData d)
Evaluate Fields.
Teuchos::RCP< PHX::DataLayout > dl_vector
Data layout for vector fields.
static HP & inst()
Private ctor.
Teuchos::RCP< PHX::DataLayout > functional
<Cell,Basis>
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_