43 #ifndef PANZER_NEUMANN_RESIDUAL_IMPL_HPP
44 #define PANZER_NEUMANN_RESIDUAL_IMPL_HPP
51 #include "Intrepid2_FunctionSpaceTools.hpp"
57 template<
typename EvalT,
typename Traits>
62 std::string residual_name = p.
get<std::string>(
"Residual Name");
63 std::string flux_name = p.
get<std::string>(
"Flux Name");
64 std::string normal_name = p.
get<std::string>(
"Normal Name");
65 std::string normal_dot_flux_name = normal_name +
" dot " + flux_name;
74 residual = PHX::MDField<ScalarT>(residual_name, basis->
functional);
75 normal_dot_flux = PHX::MDField<ScalarT>(normal_dot_flux_name, ir->dl_scalar);
76 flux = PHX::MDField<const ScalarT>(flux_name, ir->dl_vector);
77 normal = PHX::MDField<const ScalarT>(normal_name, ir->dl_vector);
79 this->addEvaluatedField(residual);
80 this->addEvaluatedField(normal_dot_flux);
81 this->addDependentField(normal);
82 this->addDependentField(flux);
86 std::string n =
"Neumann Residual Evaluator";
91 template<
typename EvalT,
typename Traits>
98 num_ip = flux.extent(1);
99 num_dim = flux.extent(2);
108 template<
typename EvalT,
typename Traits>
114 residual.deep_copy(
ScalarT(0.0));
116 for (index_t cell = 0; cell < workset.
num_cells; ++cell) {
117 for (std::size_t ip = 0; ip < num_ip; ++ip) {
118 normal_dot_flux(cell,ip) =
ScalarT(0.0);
119 for (std::size_t dim = 0; dim < num_dim; ++dim) {
120 normal_dot_flux(cell,ip) += normal(cell,ip,dim) * flux(cell,ip,dim);
127 for (index_t cell = 0; cell < workset.
num_cells; ++cell) {
128 for (std::size_t basis = 0; basis < residual.extent(1); ++basis) {
129 for (std::size_t qp = 0; qp < num_ip; ++qp) {
130 residual(cell,basis) += normal_dot_flux(cell,qp)*bv->weighted_basis_scalar(cell,basis,qp);
136 Intrepid2::FunctionSpaceTools<PHX::exec_space>::
137 integrate<ScalarT>(residual.get_view(),
138 normal_dot_flux.get_view(),
139 (this->wda(workset).bases[basis_index])->weighted_basis_scalar.get_view());
T & get(const std::string &name, T def_value)
void evaluateFields(typename Traits::EvalData d)
Teuchos::RCP< panzer::BasisIRLayout > basisIRLayout(std::string basis_type, const int basis_order, const PointRule &pt_rule)
Nonmember constructor.
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
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
#define TEUCHOS_ASSERT(assertion_test)
Teuchos::RCP< PHX::DataLayout > functional
<Cell,Basis> or <Cell,Basis>
NeumannResidual(const Teuchos::ParameterList &p)
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_