Panzer  Version of the Day
Panzer_Interface_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_INTERFACE_RESIDUAL_IMPL_HPP
44 #define PANZER_INTERFACE_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 PHX_EVALUATOR_CTOR(InterfaceResidual,p)
58 {
59  std::string residual_name = p.get<std::string>("Residual Name");
60  std::string flux_name = p.get<std::string>("Flux Name");
61  std::string normal_name = p.get<std::string>("Normal Name");
62  std::string normal_dot_flux_name = normal_name + " dot " + flux_name;
63 
64  const Teuchos::RCP<const panzer::PureBasis> basis =
65  p.get< Teuchos::RCP<const panzer::PureBasis> >("Basis");
66 
67  const Teuchos::RCP<const panzer::IntegrationRule> ir =
68  p.get< Teuchos::RCP<const panzer::IntegrationRule> >("IR");
69 
70 
71  residual = PHX::MDField<ScalarT>(residual_name, basis->functional);
72  normal_dot_flux = PHX::MDField<ScalarT>(normal_dot_flux_name, ir->dl_scalar);
73  flux = PHX::MDField<ScalarT>(flux_name, ir->dl_vector);
74  normal = PHX::MDField<ScalarT>(normal_name, ir->dl_vector);
75 
76  this->addEvaluatedField(residual);
77  this->addEvaluatedField(normal_dot_flux);
78  this->addDependentField(normal);
79  this->addDependentField(flux);
80 
81  basis_name = panzer::basisIRLayout(basis,*ir)->name();
82 
83  std::string n = "Interface Residual Evaluator";
84  this->setName(n);
85 }
86 
87 //**********************************************************************
88 PHX_POST_REGISTRATION_SETUP(InterfaceResidual,sd,fm)
89 {
90  this->utils.setFieldData(residual,fm);
91  this->utils.setFieldData(normal_dot_flux,fm);
92  this->utils.setFieldData(flux,fm);
93  this->utils.setFieldData(normal,fm);
94 
95  num_ip = flux.dimension(1);
96  num_dim = flux.dimension(2);
97 
98  TEUCHOS_ASSERT(flux.dimension(1) == normal.dimension(1));
99  TEUCHOS_ASSERT(flux.dimension(2) == normal.dimension(2));
100 
101  basis_index = panzer::getBasisIndex(basis_name, (*sd.worksets_)[0], this->wda);
102 }
103 
104 //**********************************************************************
105 PHX_EVALUATE_FIELDS(InterfaceResidual,workset)
106 {
107  residual.deep_copy(ScalarT(0.0));
108 
109  for (index_t cell = 0; cell < workset.num_cells; ++cell) {
110  for (std::size_t ip = 0; ip < num_ip; ++ip) {
111  normal_dot_flux(cell,ip) = ScalarT(0.0);
112  for (std::size_t dim = 0; dim < num_dim; ++dim) {
113  normal_dot_flux(cell,ip) += normal(cell,ip,dim) * flux(cell,ip,dim);
114  }
115  }
116  }
117 
118  // const Kokkos::DynRankView<double,PHX::Device> & weighted_basis = this->wda(workset).bases[basis_index]->weighted_basis;
119  const Teuchos::RCP<const BasisValues2<double> > bv = this->wda(workset).bases[basis_index];
120  for (index_t cell = 0; cell < workset.num_cells; ++cell) {
121  for (std::size_t basis = 0; basis < residual.dimension(1); ++basis) {
122  for (std::size_t qp = 0; qp < num_ip; ++qp) {
123  residual(cell,basis) += normal_dot_flux(cell,qp)*bv->weighted_basis_scalar(cell,basis,qp);
124  }
125  }
126  }
127 
128  if(workset.num_cells>0)
129  Intrepid2::FunctionSpaceTools::
130  integrate<ScalarT>(residual, normal_dot_flux,
131  (this->wda(workset).bases[basis_index])->weighted_basis_scalar,
132  Intrepid2::COMP_CPP);
133 }
134 
135 //**********************************************************************
136 
137 }
138 
139 #endif
PHX::MDField< ScalarT > residual
Evaluates a Dirichlet BC residual corresponding to a field value.
std::size_t basis_index
std::vector< std::string >::size_type getBasisIndex(std::string basis_name, panzer::Workset &workset, WorksetDetailsAccessor &wda)
Returns the index in the workset bases for a particular BasisIRLayout name.
std::string basis_name
PHX::MDField< ScalarT > flux
PHX::MDField< ScalarT > normal
Teuchos::RCP< panzer::BasisIRLayout > basisIRLayout(std::string basis_type, const int basis_order, const PointRule &pt_rule)
Nonmember constructor.
PHX_EVALUATOR_CTOR(BasisValues_Evaluator, p)
PHX::MDField< ScalarT > normal_dot_flux
PHX_EVALUATE_FIELDS(BasisValues_Evaluator, workset)
Teuchos::RCP< const panzer::PureBasis > basis
Interpolates basis DOF values to IP DOF values.
PHX_POST_REGISTRATION_SETUP(BasisValues_Evaluator, sd, fm)