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>
58 const Teuchos::ParameterList& p):
60 _num_quadrature_points(0),
65 Teuchos::RCP<const BasisIRLayout> basis_layout = p.get< Teuchos::RCP<BasisIRLayout> >(
"Basis").getConst();
66 Teuchos::RCP<const panzer::IntegrationRule> ir = p.get< Teuchos::RCP<panzer::IntegrationRule> >(
"IR").getConst();
69 Teuchos::RCP<const PureBasis> basis = basis_layout->getBasis();
73 TEUCHOS_TEST_FOR_EXCEPTION(!basis->supportsGrad(),std::logic_error,
"Integrator_GradBasisTimesScalar: Basis of type \"" << basis->name() <<
"\" does not support GRAD");
77 if (p.isType<Teuchos::RCP<
const std::vector<std::string> > >(
"Residual Names")){
78 const std::vector<std::string> & names = *(p.get<Teuchos::RCP<const std::vector<std::string> > >(
"Residual Names"));
79 for(
const auto & name : names){
80 PHX::MDField<ScalarT,Cell,BASIS> res(name, basis_layout->functional);
86 TEUCHOS_TEST_FOR_EXCEPTION(
_num_dims !=
_residuals.size(),std::logic_error,
"Number of residuals must equal number of gradient components (mesh dimensions).");
89 _scalar = PHX::MDField<const ScalarT,Cell,IP>(p.get<std::string>(
"Scalar Name"), ir->dl_scalar);
95 if (p.isType<Teuchos::RCP<
const std::vector<std::string> > >(
"Field Multipliers")){
96 const std::vector<std::string> & field_multiplier_names = *(p.get<Teuchos::RCP<const std::vector<std::string> > >(
"Field Multipliers"));
97 for (
const std::string & name : field_multiplier_names){
98 _field_multipliers.push_back(PHX::MDField<const ScalarT,Cell,IP>(name, ir->dl_scalar));
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){