43 #ifndef PANZER_EVALUATOR_TRANSIENT_BASISTIMESSCALAR_IMPL_HPP
44 #define PANZER_EVALUATOR_TRANSIENT_BASISTIMESSCALAR_IMPL_HPP
46 #include "Intrepid2_FunctionSpaceTools.hpp"
50 #include "Kokkos_ViewFactory.hpp"
55 template<
typename EvalT,
typename Traits>
59 residual( p.get<std::string>(
"Residual Name"),
61 scalar( p.get<std::string>(
"Value Name"),
73 "Integrator_TransientBasisTimesScalar: Basis of type \"" << basis->name() <<
"\" is not a "
77 this->addDependentField(
scalar);
83 const std::vector<std::string>& field_multiplier_names =
86 for (std::vector<std::string>::const_iterator name =
87 field_multiplier_names.begin();
88 name != field_multiplier_names.end(); ++name) {
94 for (
typename std::vector<PHX::MDField<const ScalarT,Cell,IP> >::iterator
field =
field_multipliers.begin();
96 this->addDependentField(*
field);
98 std::string n =
"Integrator_TransientBasisTimesScalar: " +
residual.fieldTag().name();
103 template<
typename EvalT,
typename Traits>
110 this->utils.setFieldData(residual,fm);
111 this->utils.setFieldData(scalar,fm);
113 num_nodes = residual.extent(1);
114 num_qp = scalar.extent(1);
122 template<
typename EvalT,
typename Traits>
133 Kokkos::deep_copy (residual.get_static_view(),
ScalarT(0.0));
135 for (index_t cell = 0; cell < workset.
num_cells; ++cell) {
136 for (std::size_t qp = 0; qp < num_qp; ++qp) {
138 for (
typename std::vector<PHX::MDField<const ScalarT,Cell,IP> >::iterator
field = field_multipliers.begin();
140 tmp(cell,qp) = tmp(cell,qp) * (*field)(cell,qp);
145 Intrepid2::FunctionSpaceTools<PHX::exec_space>::
146 integrate<ScalarT>(residual.get_view(),
148 (this->wda(workset).bases[basis_index])->weighted_basis_scalar.get_view());
153 template<
typename EvalT,
typename TRAITS>
158 p->
set<std::string>(
"Residual Name",
"?");
159 p->
set<std::string>(
"Value Name",
"?");
161 p->
set(
"Basis", basis);
164 p->
set<
double>(
"Multiplier", 1.0);
166 p->
set(
"Field Multipliers", fms);
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.
T & get(const std::string &name, T def_value)
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
std::vector< PHX::MDField< const ScalarT, Cell, IP > > field_multipliers
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
bool evaluate_transient_terms
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
double multiplier
The scalar multiplier out in front of the integral ( ).
PHX::MDField< ScalarT, Cell, BASIS > residual
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
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.
typename EvalT::ScalarT ScalarT
Integrator_TransientBasisTimesScalar(const Teuchos::ParameterList &p)
void validateParameters(ParameterList const &validParamList, int const depth=1000, EValidateUsed const validateUsed=VALIDATE_USED_ENABLED, EValidateDefaults const validateDefaults=VALIDATE_DEFAULTS_ENABLED) const
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...
bool isType(const std::string &name) const
PHX::MDField< const ScalarT, Cell, IP > scalar
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_
void evaluateFields(typename Traits::EvalData d)