11#ifndef __Panzer_Integerator_BasisTimesTensorTimesVector_impl_hpp__
12#define __Panzer_Integerator_BasisTimesTensorTimesVector_impl_hpp__
32 template<
typename EvalT,
typename Traits>
36 const std::string& resName,
37 const std::string& valName,
40 const std::string& tensorName)
43 useDescriptors_(false),
44 basisName_(basis.name())
54 using std::invalid_argument;
55 using std::logic_error;
60 TEUCHOS_TEST_FOR_EXCEPTION(resName ==
"", invalid_argument,
"Error: " \
61 "Integrator_BasisTimesTensorTimesVector called with an empty residual name.")
62 TEUCHOS_TEST_FOR_EXCEPTION(valName ==
"", invalid_argument,
"Error: " \
63 "Integrator_BasisTimesTensorTimesVector called with an empty value name.")
64 RCP<const PureBasis> tmpBasis = basis.
getBasis();
65 TEUCHOS_TEST_FOR_EXCEPTION(not tmpBasis->isVectorBasis(), logic_error,
66 "Error: Integrator_BasisTimesTensorTimesVector: Basis of type \""
67 << tmpBasis->name() <<
"\" is not a vector basis.");
68 TEUCHOS_TEST_FOR_EXCEPTION(not tmpBasis->requiresOrientations(),
69 logic_error,
"Integrator_BasisTimesTensorTimesVector: Basis of type \""
70 << tmpBasis->name() <<
"\" does not require orientations. This seems " \
71 "very strange, so I'm failing.");
75 this->addDependentField(
vector_);
81 this->addContributedField(
field_);
83 this->addEvaluatedField(
field_);
86 tensor_ = MDField<const ScalarT, Cell, IP, Dim, Dim>(tensorName, ir.
dl_tensor);
87 this->addDependentField(
tensor_);
90 string n(
"Integrator_BasisTimesTensorTimesVector (");
95 n +=
", " + print<EvalT>() +
"): " +
field_.fieldTag().name();
104 template<
typename EvalT,
typename Traits>
108 const PHX::FieldTag& resTag,
109 const PHX::FieldTag& valTag,
112 const PHX::FieldTag& tensorTag)
115 useDescriptors_(true),
127 using std::invalid_argument;
128 using std::logic_error;
133 TEUCHOS_TEST_FOR_EXCEPTION(not ((
bd_.
getType() ==
"HCurl") or
134 (
bd_.
getType() ==
"HDiv")), logic_error,
"Error: " \
135 "Integrator_BasisTimesTensorTimesVector: Basis of type \"" <<
bd_.
getType()
136 <<
"\" is not a vector basis.")
140 this->addDependentField(
vector_);
146 this->addContributedField(
field_);
148 this->addEvaluatedField(
field_);
152 this->addDependentField(
tensor_);
155 string n(
"Integrator_BasisTimesTensorTimesVector (");
160 n +=
", " + print<EvalT>() +
"): " +
field_.fieldTag().name();
169 template<
typename EvalT,
typename Traits>
172 const Teuchos::ParameterList& p)
176 p.get<
std::string>(
"Residual Name"),
177 p.get<
std::string>(
"Value Name"),
180 p.get<
std::string>(
"Tensor Name"))
182 using Teuchos::ParameterList;
187 p.validateParameters(*validParams);
195 template<
typename EvalT,
typename Traits>
206 kokkosTensor_ = tensor_.get_static_view();
210 numQP_ = vector_.extent(1);
211 numDim_ = vector_.extent(2);
214 if (not useDescriptors_)
217 if (Sacado::IsADType<ScalarT>::value) {
218 const auto fadSize = Kokkos::dimension_scalar(field_.get_view());
219 tmp_ = PHX::View<ScalarT*>(
"GradBasisDotVector::tmp_",field_.extent(0),fadSize);
221 tmp_ = PHX::View<ScalarT*>(
"GradBasisDotVector::tmp_",field_.extent(0));
231 template<
typename EvalT,
typename Traits>
232 KOKKOS_INLINE_FUNCTION
236 const size_t& cell)
const
241 const int numBases(basis_.extent(1));
243 for (
int basis(0); basis < numBases; ++basis)
244 field_(cell, basis) = 0.0;
249 for (
int qp(0); qp < numQP_; ++qp)
251 for (
int dim(0); dim < numDim_; ++dim)
254 for (
int dim2(0); dim2 < numDim_; ++dim2)
255 tmp_(cell) += kokkosTensor_(cell, qp, dim, dim2) * vector_(cell, qp, dim2);
256 for (
int basis(0); basis < numBases; ++basis)
257 field_(cell, basis) += basis_(cell, basis, qp, dim) * tmp_(cell);
267 template<
typename EvalT,
typename Traits>
273 using Kokkos::parallel_for;
274 using Kokkos::RangePolicy;
279 *this->wda(workset).bases[basisIndex_];
283 parallel_for(RangePolicy<>(0, workset.
num_cells), *
this);
291 template<
typename EvalT,
typename Traits>
292 Teuchos::RCP<Teuchos::ParameterList>
300 using Teuchos::ParameterList;
305 RCP<ParameterList> p = rcp(
new ParameterList);
306 p->set<
string>(
"Residual Name",
"?");
307 p->set<
string>(
"Value Name",
"?");
308 RCP<BasisIRLayout> basis;
309 p->set(
"Basis", basis);
310 RCP<IntegrationRule> ir;
312 p->set<std::string>(
"Tensor Name",
"?");
panzer::EvaluatorStyle evalStyle
The EvaluatorStyle of the parent Integrator_CurlBasisDotVector object.
const std::string & getType() const
Get type of basis.
Teuchos::RCP< PHX::DataLayout > functional
<Cell,Basis>
Teuchos::RCP< const PureBasis > getBasis() const
ConstArray_CellBasisIP getBasisValues(const bool weighted, const bool cache=true, const bool force=false) const
Get the basis values evaluated at mesh points.
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.
void postRegistrationSetup(typename Traits::SetupData sd, PHX::FieldManager< Traits > &fm)
Post-Registration Setup.
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
Get Valid Parameters.
KOKKOS_INLINE_FUNCTION void operator()(const std::size_t &cell) const
Perform the integration.
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_BasisTimesTensorTimesVector(const panzer::EvaluatorStyle &evalStyle, const std::string &resName, const std::string &valName, const panzer::BasisIRLayout &basis, const panzer::IntegrationRule &ir, const std::string &tensorName)
Main Constructor.
void evaluateFields(typename Traits::EvalData workset)
Evaluate Fields.
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( ).
Teuchos::RCP< PHX::DataLayout > dl_vector
Data layout for vector fields.
Teuchos::RCP< PHX::DataLayout > dl_tensor
Data layout for rank-2 tensor 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.
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_