11 #ifndef PANZER_INTERFACE_RESIDUAL_IMPL_HPP
12 #define PANZER_INTERFACE_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 =
"Interface 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));
84 for (index_t cell = 0; cell < workset.
num_cells; ++cell) {
85 for (std::size_t ip = 0; ip < num_ip; ++ip) {
86 normal_dot_flux(cell,ip) =
ScalarT(0.0);
87 for (std::size_t dim = 0; dim < num_dim; ++dim) {
88 normal_dot_flux(cell,ip) += normal(cell,ip,dim) * flux(cell,ip,dim);
95 for (index_t cell = 0; cell < workset.
num_cells; ++cell) {
96 for (std::size_t basis = 0; basis < residual.extent(1); ++basis) {
97 for (std::size_t qp = 0; qp < num_ip; ++qp) {
98 residual(cell,basis) += normal_dot_flux(cell,qp)*bv->weighted_basis_scalar(cell,basis,qp);
104 Intrepid2::FunctionSpaceTools<PHX::exec_space>::
105 integrate(residual.get_view(),
106 normal_dot_flux.get_view(),
107 (this->wda(workset).bases[basis_index])->weighted_basis_scalar.get_view());
int num_cells
DEPRECATED - use: numCells()
T & get(const std::string &name, T def_value)
typename EvalT::ScalarT ScalarT
Teuchos::RCP< panzer::BasisIRLayout > basisIRLayout(std::string basis_type, const int basis_order, const PointRule &pt_rule)
Nonmember constructor.
InterfaceResidual(const Teuchos::ParameterList &p)
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.
void evaluateFields(typename Traits::EvalData d)
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
#define TEUCHOS_ASSERT(assertion_test)
Teuchos::RCP< PHX::DataLayout > functional
<Cell,Basis> or <Cell,Basis>
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_