11 #ifndef PANZER_NEUMANN_RESIDUAL_IMPL_HPP
12 #define PANZER_NEUMANN_RESIDUAL_IMPL_HPP
19 #include "Intrepid2_FunctionSpaceTools.hpp"
25 template<
typename EvalT,
typename Traits>
30 std::string residual_name = p.
get<std::string>(
"Residual Name");
31 std::string flux_name = p.
get<std::string>(
"Flux Name");
32 std::string normal_name = p.
get<std::string>(
"Normal Name");
33 std::string normal_dot_flux_name = normal_name +
" dot " + flux_name;
47 this->addEvaluatedField(residual);
48 this->addEvaluatedField(normal_dot_flux);
49 this->addDependentField(normal);
50 this->addDependentField(flux);
54 std::string n =
"Neumann Residual Evaluator";
59 template<
typename EvalT,
typename Traits>
66 num_ip = flux.extent(1);
67 num_dim = flux.extent(2);
76 template<
typename EvalT,
typename Traits>
82 residual.deep_copy(
ScalarT(0.0));
83 auto normal_dot_flux_v = normal_dot_flux.get_static_view();
84 auto normal_v = normal.get_static_view();
85 auto flux_v = flux.get_static_view();
87 std::size_t l_num_ip = num_ip, l_num_dim = num_dim;
88 Kokkos::parallel_for(
"NeumannResidual", workset.
num_cells, KOKKOS_LAMBDA (
const index_t cell) {
89 for (std::size_t ip = 0; ip < l_num_ip; ++ip) {
91 normal_dot_flux_v(cell,ip) = 0.0;
92 for (std::size_t dim = 0; dim < l_num_dim; ++dim) {
93 normal_dot_flux_v(cell,ip) += normal_v(cell,ip,dim) * flux_v(cell,ip,dim);
98 auto weighted_basis_scalar = this->wda(workset).bases[basis_index]->weighted_basis_scalar.get_static_view();
99 auto residual_v = residual.get_static_view();
100 Kokkos::parallel_for(
"NeumannResidual", workset.
num_cells, KOKKOS_LAMBDA (
const index_t cell) {
101 for (std::size_t basis = 0; basis < residual_v.extent(1); ++basis) {
102 for (std::size_t qp = 0; qp < l_num_ip; ++qp) {
103 residual_v(cell,basis) += normal_dot_flux_v(cell,qp)*weighted_basis_scalar(cell,basis,qp);
int num_cells
DEPRECATED - use: numCells()
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_