43 #ifndef PANZER_EVALUATOR_GRADBASISTIMESSCALAR_IMPL_HPP 44 #define PANZER_EVALUATOR_GRADBASISTIMESSCALAR_IMPL_HPP 46 #include "Intrepid2_FunctionSpaceTools.hpp" 50 #include "Kokkos_ViewFactory.hpp" 55 template<
typename EvalT,
typename Traits>
60 _num_quadrature_points(0),
73 TEUCHOS_TEST_FOR_EXCEPTION(!basis->supportsGrad(),std::logic_error,
"Integrator_GradBasisTimesScalar: Basis of type \"" << basis->name() <<
"\" does not support GRAD");
79 for(
const auto & name : names){
80 PHX::MDField<ScalarT,Cell,BASIS> res(name, basis_layout->functional);
89 _scalar = PHX::MDField<const ScalarT,Cell,IP>(p.
get<std::string>(
"Scalar Name"), ir->
dl_scalar);
97 for (
const std::string & name : field_multiplier_names){
104 this->addEvaluatedField(residual);
108 this->addDependentField(
_scalar);
110 this->addDependentField(
field);
114 this->setName(
"Integrator_GradBasisTimesScalar: " +
_scalar.fieldTag().name());
118 template<
typename EvalT,
typename Traits>
125 _num_basis_nodes = _residuals[0].extent(1);
126 _num_quadrature_points = _scalar.extent(1);
135 template<
typename EvalT,
typename Traits>
143 for (
int i(0); i < static_cast<int>(
_num_dims); ++i)
144 Kokkos::deep_copy(_residuals[i].get_static_view(),
ScalarT(0.0));
147 for (
int i=0; i < _scalar.extent_int(0); ++i)
148 for (
int j=0; j < _scalar.extent_int(1); ++j)
149 _tmp(i,j) = _multiplier * _scalar(i,j);
152 for(
auto & field_data : _field_multipliers){
153 for (index_t cell = 0; cell < workset.
num_cells; ++cell) {
154 for (std::size_t qp = 0; qp < _num_quadrature_points; ++qp) {
155 _tmp(cell,qp) *= field_data(cell,qp);
162 for (index_t cell = 0; cell < workset.
num_cells; ++cell){
163 for (std::size_t basis = 0; basis < _num_basis_nodes; ++basis){
164 for (std::size_t qp = 0; qp < _num_quadrature_points; ++qp){
165 for(std::size_t dim = 0; dim <
_num_dims; ++dim){
typename EvalT::ScalarT ScalarT
PHX::MDField< const ScalarT, Cell, IP > _scalar
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
T & get(const std::string &name, T def_value)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
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.
void evaluateFields(typename Traits::EvalData d)
Teuchos::RCP< PHX::DataLayout > dl_scalar
Data layout for scalar fields.
bool isType(const std::string &name) 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.
std::vector< PHX::MDField< ScalarT, Cell, BASIS > > _residuals
Array_CellBasisIPDim weighted_grad_basis
Integrator_GradBasisTimesScalar(const Teuchos::ParameterList &p)
Teuchos::RCP< PHX::DataLayout > dl_vector
Data layout for vector fields.
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::vector< PHX::MDField< const ScalarT, Cell, IP > > _field_multipliers
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_