Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_Interface_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_INTERFACE_RESIDUAL_IMPL_HPP
12#define PANZER_INTERFACE_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
35 const Teuchos::RCP<const panzer::PureBasis> basis =
36 p.get< Teuchos::RCP<const panzer::PureBasis> >("Basis");
37
38 const Teuchos::RCP<const panzer::IntegrationRule> ir =
39 p.get< Teuchos::RCP<const panzer::IntegrationRule> >("IR");
40
41
42 residual = PHX::MDField<ScalarT>(residual_name, basis->functional);
43 normal_dot_flux = PHX::MDField<ScalarT>(normal_dot_flux_name, ir->dl_scalar);
44 flux = PHX::MDField<const ScalarT>(flux_name, ir->dl_vector);
45 normal = PHX::MDField<const ScalarT>(normal_name, ir->dl_vector);
46
47 this->addEvaluatedField(residual);
48 this->addEvaluatedField(normal_dot_flux);
49 this->addDependentField(normal);
50 this->addDependentField(flux);
51
52 basis_name = panzer::basisIRLayout(basis,*ir)->name();
53
54 std::string n = "Interface Residual Evaluator";
55 this->setName(n);
56}
57
58//**********************************************************************
59template<typename EvalT, typename Traits>
60void
63 typename Traits::SetupData sd,
65{
66 num_ip = flux.extent(1);
67 num_dim = flux.extent(2);
68
69 TEUCHOS_ASSERT(flux.extent(1) == normal.extent(1));
70 TEUCHOS_ASSERT(flux.extent(2) == normal.extent(2));
71
72 basis_index = panzer::getBasisIndex(basis_name, (*sd.worksets_)[0], this->wda);
73}
74
75//**********************************************************************
76template<typename EvalT, typename Traits>
77void
80 typename Traits::EvalData workset)
81{
82 residual.deep_copy(ScalarT(0.0));
83
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);
89 }
90 }
91 }
92
93 // const Kokkos::DynRankView<double,PHX::Device> & weighted_basis = this->wda(workset).bases[basis_index]->weighted_basis;
94 const Teuchos::RCP<const BasisValues2<double> > bv = this->wda(workset).bases[basis_index];
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);
99 }
100 }
101 }
102
103 if(workset.num_cells>0)
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());
108}
109
110//**********************************************************************
111
112}
113
114#endif
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
InterfaceResidual(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_