Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_WeakDirichlet_Residual_impl.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Panzer: A partial differential equation assembly
4// engine for strongly coupled complex multiphysics systems
5//
6// Copyright 2011 NTESS and the Panzer contributors.
7// SPDX-License-Identifier: BSD-3-Clause
8// *****************************************************************************
9// @HEADER
10
11#ifndef PANZER_WEAKDIRICHLET_RESIDUAL_IMPL_HPP
12#define PANZER_WEAKDIRICHLET_RESIDUAL_IMPL_HPP
13
14#include <cstddef>
15#include <string>
16#include <vector>
19#include "Intrepid2_FunctionSpaceTools.hpp"
20#include "Teuchos_RCP.hpp"
21
22namespace panzer {
23
24//**********************************************************************
25template<typename EvalT, typename Traits>
28 const Teuchos::ParameterList& p)
29{
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;
34 std::string dof_name = p.get<std::string>("DOF Name");
35 std::string value_name = p.get<std::string>("Value Name");
36 std::string sigma_name = p.get<std::string>("Sigma Name");
37
38 const Teuchos::RCP<const panzer::PureBasis> basis =
39 p.get< Teuchos::RCP<const panzer::PureBasis> >("Basis");
40
41 const Teuchos::RCP<const panzer::IntegrationRule> ir =
42 p.get< Teuchos::RCP<const panzer::IntegrationRule> >("IR");
43
44
45 residual = PHX::MDField<ScalarT>(residual_name, basis->functional);
46 normal_dot_flux_plus_pen = PHX::MDField<ScalarT>(normal_dot_flux_name, ir->dl_scalar);
47 flux = PHX::MDField<const ScalarT>(flux_name, ir->dl_vector);
48 normal = PHX::MDField<const ScalarT>(normal_name, ir->dl_vector);
49 dof = PHX::MDField<const ScalarT>(dof_name, ir->dl_scalar);
50 value = PHX::MDField<const ScalarT>(value_name, ir->dl_scalar);
51 sigma = PHX::MDField<const ScalarT>(sigma_name, ir->dl_scalar);
52
53 this->addEvaluatedField(residual);
54 this->addEvaluatedField(normal_dot_flux_plus_pen);
55 this->addDependentField(normal);
56 this->addDependentField(flux);
57 this->addDependentField(dof);
58 this->addDependentField(value);
59 this->addDependentField(sigma);
60
61 basis_name = panzer::basisIRLayout(basis,*ir)->name();
62
63 std::string n = "Weak Dirichlet Residual Evaluator";
64 this->setName(n);
65}
66
67//**********************************************************************
68template<typename EvalT, typename Traits>
69void
72 typename Traits::SetupData sd,
74{
75 num_ip = flux.extent(1);
76 num_dim = flux.extent(2);
77
78 TEUCHOS_ASSERT(flux.extent(1) == normal.extent(1));
79 TEUCHOS_ASSERT(flux.extent(2) == normal.extent(2));
80
81 basis_index = panzer::getBasisIndex(basis_name, (*sd.worksets_)[0], this->wda);
82}
83
84//**********************************************************************
85template<typename EvalT, typename Traits>
86void
89 typename Traits::EvalData workset)
90{
91 for (index_t cell = 0; cell < workset.num_cells; ++cell) {
92 for (std::size_t ip = 0; ip < num_ip; ++ip) {
93 normal_dot_flux_plus_pen(cell,ip) = ScalarT(0.0);
94 for (std::size_t dim = 0; dim < num_dim; ++dim) {
95 normal_dot_flux_plus_pen(cell,ip) += normal(cell,ip,dim) * flux(cell,ip,dim);
96 }
97 normal_dot_flux_plus_pen(cell,ip) += sigma(cell,ip) * (dof(cell,ip) - value(cell,ip));
98 }
99 }
100
101 if(workset.num_cells>0)
102 Intrepid2::FunctionSpaceTools<PHX::exec_space>::
103 integrate(residual.get_view(),
104 normal_dot_flux_plus_pen.get_view(),
105 (this->wda(workset).bases[basis_index])->weighted_basis_scalar.get_view());
106
107}
108
109//**********************************************************************
110
111}
112
113#endif
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
WeakDirichletResidual(const Teuchos::ParameterList &p)
void evaluateFields(typename Traits::EvalData d)
int num_cells
DEPRECATED - use: numCells()
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.
Teuchos::RCP< panzer::BasisIRLayout > basisIRLayout(std::string basis_type, const int basis_order, const PointRule &pt_rule)
Nonmember constructor.
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_