43 #ifndef PANZER_EVALUATOR_GRADBASISDOTVECTOR_IMPL_HPP 44 #define PANZER_EVALUATOR_GRADBASISDOTVECTOR_IMPL_HPP 46 #include "Intrepid2_FunctionSpaceTools.hpp" 50 #include "Kokkos_ViewFactory.hpp" 52 #define PANZER_USE_FAST_QUAD 1 75 flux = PHX::MDField<const ScalarT,Cell,IP,Dim>( p.get<std::string>(
"Flux Name"),
vector);
82 "Integrator_GradBasisDotVector: Basis of type \"" <<
basis->
name() <<
"\" does not support GRAD");
86 "Integrator_GradBasisDotVector: Dimension of space exceeds dimension of vector.");
89 this->addDependentField(
flux);
93 if (p.isType<
Teuchos::RCP<
const std::vector<std::string> > >(
"Field Multipliers"))
95 const std::vector<std::string>& field_multiplier_names =
98 for (std::vector<std::string>::const_iterator name = field_multiplier_names.begin();
99 name != field_multiplier_names.end(); ++name)
107 this->addDependentField(
field);
111 "Integrator_GradBasisDotVector: " +
residual.fieldTag().name();
119 this->utils.setFieldData(
residual,fm);
120 this->utils.setFieldData(
flux,fm);
123 this->utils.setFieldData(
field,fm);
139 Kokkos::deep_copy(
residual.get_static_view(), ScalarT(0.0));
141 #if PANZER_USE_FAST_QUAD 143 for (
int i=0; i <
flux.extent_int(0); ++i)
144 for (
int j=0; j <
flux.extent_int(1); ++j)
145 for (
int k=0; k <
num_dim; ++k)
152 PHX::MDField<const ScalarT,Cell,IP> field_data =
field;
154 for (index_t cell = 0; cell < workset.num_cells; ++cell) {
155 for (std::size_t qp = 0; qp <
num_qp; ++qp) {
156 ScalarT tmpVar = field_data(cell,qp);
158 for (std::size_t dim = 0; dim <
num_dim; ++dim)
159 tmp(cell,qp,dim) *= tmpVar;
168 for (index_t cell = 0; cell < workset.num_cells; ++cell) {
170 for (std::size_t qp = 0; qp <
num_qp; ++qp) {
171 for (std::size_t dim = 0; dim <
num_dim; ++dim) {
184 for (index_t cell = 0; cell < workset.num_cells; ++cell)
186 for (std::size_t qp = 0; qp <
num_qp; ++qp)
188 ScalarT tmpVar = 1.0;
191 tmpVar = tmpVar * (*
field)(cell,qp);
193 for (std::size_t dim = 0; dim <
num_dim; ++dim)
198 if(workset.num_cells>0)
199 Intrepid2::FunctionSpaceTools::
201 (this->wda(workset).bases[
basis_index])->weighted_grad_basis,
202 Intrepid2::COMP_CPP);
PHX::MDField< ScalarT > residual
Evaluates a Dirichlet BC residual corresponding to a field value.
std::vector< std::string >::size_type getBasisIndex(std::string basis_name, panzer::Workset &workset, WorksetDetailsAccessor &wda)
Returns the index in the workset bases for a particular BasisIRLayout name.
std::string name() const
A unique key that is the combination of the basis type and basis order.
PHX::MDField< ScalarT, Cell > tmp
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
PHX::MDField< ScalarT > flux
bool supportsGrad() const
PHX::MDField< ScalarT > vector
Kokkos::DynRankView< typename InputArray::value_type, PHX::Device > createDynRankView(const InputArray &a, const std::string &name, const DimensionPack... dims)
Wrapper to simplify Panzer use of Sacado ViewFactory.
Teuchos::RCP< PHX::DataLayout > dl_scalar
Data layout for scalar fields.
PHX_EVALUATOR_CTOR(BasisValues_Evaluator, p)
PHX::MDField< const ScalarT, Cell, IP > field
Array_CellBasisIPDim weighted_grad_basis
PHX_EVALUATE_FIELDS(BasisValues_Evaluator, workset)
Teuchos::RCP< const panzer::PureBasis > basis
Interpolates basis DOF values to IP DOF values.
std::vector< PHX::MDField< const ScalarT, Cell, IP > > field_multipliers
PHX_POST_REGISTRATION_SETUP(BasisValues_Evaluator, sd, fm)