Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_Neumann_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_NEUMANN_RESIDUAL_IMPL_HPP
12#define PANZER_NEUMANN_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 = "Neumann 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 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();
86
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) {
90 // normal_dot_flux_v(cell,ip) = ScalarT(0.0);
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);
94 }
95 }
96 });
97
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);
104 }
105 }
106 });
107
108 // Intrepid2 integrate calls are broken for sacado scalar types. The
109 // temporaries outside of views use new/delete on host, which is not
110 // supported by HIP. Need to roll our own integrate call for now.
111 // if(workset.num_cells>0)
112 // Intrepid2::FunctionSpaceTools<PHX::exec_space>::
113 // integrate<ScalarT>(residual.get_view(),
114 // normal_dot_flux.get_view(),
115 // (this->wda(workset).bases[basis_index])->weighted_basis_scalar.get_view());
116
117 Kokkos::fence();
118
119}
120
121//**********************************************************************
122
123}
124
125#endif
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
NeumannResidual(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_