Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_Neumann_Residual_impl.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Panzer: A partial differential equation assembly
5 // engine for strongly coupled complex multiphysics systems
6 // Copyright (2011) Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Roger P. Pawlowski (rppawlo@sandia.gov) and
39 // Eric C. Cyr (eccyr@sandia.gov)
40 // ***********************************************************************
41 // @HEADER
42 
43 #ifndef PANZER_NEUMANN_RESIDUAL_IMPL_HPP
44 #define PANZER_NEUMANN_RESIDUAL_IMPL_HPP
45 
46 #include <cstddef>
47 #include <string>
48 #include <vector>
49 #include "Panzer_BasisIRLayout.hpp"
51 #include "Intrepid2_FunctionSpaceTools.hpp"
52 #include "Teuchos_RCP.hpp"
53 
54 namespace panzer {
55 
56 //**********************************************************************
57 template<typename EvalT, typename Traits>
60  const Teuchos::ParameterList& p)
61 {
62  std::string residual_name = p.get<std::string>("Residual Name");
63  std::string flux_name = p.get<std::string>("Flux Name");
64  std::string normal_name = p.get<std::string>("Normal Name");
65  std::string normal_dot_flux_name = normal_name + " dot " + flux_name;
66 
69 
72 
73 
74  residual = PHX::MDField<ScalarT>(residual_name, basis->functional);
75  normal_dot_flux = PHX::MDField<ScalarT>(normal_dot_flux_name, ir->dl_scalar);
76  flux = PHX::MDField<const ScalarT>(flux_name, ir->dl_vector);
77  normal = PHX::MDField<const ScalarT>(normal_name, ir->dl_vector);
78 
79  this->addEvaluatedField(residual);
80  this->addEvaluatedField(normal_dot_flux);
81  this->addDependentField(normal);
82  this->addDependentField(flux);
83 
84  basis_name = panzer::basisIRLayout(basis,*ir)->name();
85 
86  std::string n = "Neumann Residual Evaluator";
87  this->setName(n);
88 }
89 
90 //**********************************************************************
91 template<typename EvalT, typename Traits>
92 void
95  typename Traits::SetupData sd,
96  PHX::FieldManager<Traits>& /* fm */)
97 {
98  num_ip = flux.extent(1);
99  num_dim = flux.extent(2);
100 
101  TEUCHOS_ASSERT(flux.extent(1) == normal.extent(1));
102  TEUCHOS_ASSERT(flux.extent(2) == normal.extent(2));
103 
104  basis_index = panzer::getBasisIndex(basis_name, (*sd.worksets_)[0], this->wda);
105 }
106 
107 //**********************************************************************
108 template<typename EvalT, typename Traits>
109 void
112  typename Traits::EvalData workset)
113 {
114  residual.deep_copy(ScalarT(0.0));
115  auto normal_dot_flux_v = normal_dot_flux.get_static_view();
116  auto normal_v = normal.get_static_view();
117  auto flux_v = flux.get_static_view();
118 
119  std::size_t l_num_ip = num_ip, l_num_dim = num_dim;
120  Kokkos::parallel_for("NeumannResidual", workset.num_cells, KOKKOS_LAMBDA (const index_t cell) {
121  for (std::size_t ip = 0; ip < l_num_ip; ++ip) {
122  // normal_dot_flux_v(cell,ip) = ScalarT(0.0);
123  normal_dot_flux_v(cell,ip) = 0.0;
124  for (std::size_t dim = 0; dim < l_num_dim; ++dim) {
125  normal_dot_flux_v(cell,ip) += normal_v(cell,ip,dim) * flux_v(cell,ip,dim);
126  }
127  }
128  });
129 
130  auto weighted_basis_scalar = this->wda(workset).bases[basis_index]->weighted_basis_scalar.get_static_view();
131  auto residual_v = residual.get_static_view();
132  Kokkos::parallel_for("NeumannResidual", workset.num_cells, KOKKOS_LAMBDA (const index_t cell) {
133  for (std::size_t basis = 0; basis < residual_v.extent(1); ++basis) {
134  for (std::size_t qp = 0; qp < l_num_ip; ++qp) {
135  residual_v(cell,basis) += normal_dot_flux_v(cell,qp)*weighted_basis_scalar(cell,basis,qp);
136  }
137  }
138  });
139 
140  // Intrepid2 integrate calls are broken for sacado scalar types. The
141  // temporaries outside of views use new/delete on host, which is not
142  // supported by HIP. Need to roll our own integrate call for now.
143  // if(workset.num_cells>0)
144  // Intrepid2::FunctionSpaceTools<PHX::exec_space>::
145  // integrate<ScalarT>(residual.get_view(),
146  // normal_dot_flux.get_view(),
147  // (this->wda(workset).bases[basis_index])->weighted_basis_scalar.get_view());
148 
149  Kokkos::fence();
150 
151 }
152 
153 //**********************************************************************
154 
155 }
156 
157 #endif
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.
#define TEUCHOS_ASSERT(assertion_test)
Teuchos::RCP< PHX::DataLayout > functional
&lt;Cell,Basis&gt; or &lt;Cell,Basis&gt;
NeumannResidual(const Teuchos::ParameterList &p)
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_