Panzer  Version of the Day
Panzer_WeakDirichlet_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_WEAKDIRICHLET_RESIDUAL_IMPL_HPP
44 #define PANZER_WEAKDIRICHLET_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(WeakDirichletResidual,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  std::string dof_name = p.get<std::string>("DOF Name");
64  std::string value_name = p.get<std::string>("Value Name");
65  std::string sigma_name = p.get<std::string>("Sigma Name");
66 
67  const Teuchos::RCP<const panzer::PureBasis> basis =
68  p.get< Teuchos::RCP<const panzer::PureBasis> >("Basis");
69 
70  const Teuchos::RCP<const panzer::IntegrationRule> ir =
71  p.get< Teuchos::RCP<const panzer::IntegrationRule> >("IR");
72 
73 
74  residual = PHX::MDField<ScalarT>(residual_name, basis->functional);
75  normal_dot_flux_plus_pen = PHX::MDField<ScalarT>(normal_dot_flux_name, ir->dl_scalar);
76  flux = PHX::MDField<ScalarT>(flux_name, ir->dl_vector);
77  normal = PHX::MDField<ScalarT>(normal_name, ir->dl_vector);
78  dof = PHX::MDField<ScalarT>(dof_name, ir->dl_scalar);
79  value = PHX::MDField<ScalarT>(value_name, ir->dl_scalar);
80  sigma = PHX::MDField<ScalarT>(sigma_name, ir->dl_scalar);
81 
82  this->addEvaluatedField(residual);
83  this->addEvaluatedField(normal_dot_flux_plus_pen);
84  this->addDependentField(normal);
85  this->addDependentField(flux);
86  this->addDependentField(dof);
87  this->addDependentField(value);
88  this->addDependentField(sigma);
89 
90  basis_name = panzer::basisIRLayout(basis,*ir)->name();
91 
92  std::string n = "Weak Dirichlet Residual Evaluator";
93  this->setName(n);
94 }
95 
96 //**********************************************************************
97 PHX_POST_REGISTRATION_SETUP(WeakDirichletResidual,sd,fm)
98 {
99  this->utils.setFieldData(residual,fm);
100  this->utils.setFieldData(normal_dot_flux_plus_pen,fm);
101  this->utils.setFieldData(flux,fm);
102  this->utils.setFieldData(normal,fm);
103  this->utils.setFieldData(dof,fm);
104  this->utils.setFieldData(value,fm);
105  this->utils.setFieldData(sigma,fm);
106 
107  num_ip = flux.dimension(1);
108  num_dim = flux.dimension(2);
109 
110  TEUCHOS_ASSERT(flux.dimension(1) == normal.dimension(1));
111  TEUCHOS_ASSERT(flux.dimension(2) == normal.dimension(2));
112 
113  basis_index = panzer::getBasisIndex(basis_name, (*sd.worksets_)[0], this->wda);
114 }
115 
116 //**********************************************************************
117 PHX_EVALUATE_FIELDS(WeakDirichletResidual,workset)
118 {
119  for (index_t cell = 0; cell < workset.num_cells; ++cell) {
120  for (std::size_t ip = 0; ip < num_ip; ++ip) {
121  normal_dot_flux_plus_pen(cell,ip) = ScalarT(0.0);
122  for (std::size_t dim = 0; dim < num_dim; ++dim) {
123  normal_dot_flux_plus_pen(cell,ip) += normal(cell,ip,dim) * flux(cell,ip,dim);
124  }
125  normal_dot_flux_plus_pen(cell,ip) += sigma(cell,ip) * (dof(cell,ip) - value(cell,ip));
126  }
127  }
128 
129  if(workset.num_cells>0)
130  Intrepid2::FunctionSpaceTools::
131  integrate<ScalarT>(residual, normal_dot_flux_plus_pen,
132  (this->wda(workset).bases[basis_index])->weighted_basis_scalar,
133  Intrepid2::COMP_CPP);
134 
135 }
136 
137 //**********************************************************************
138 
139 }
140 
141 #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.
PHX::MDField< ScalarT > sigma
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::MDField< const ScalarT > dof
PHX_EVALUATOR_CTOR(BasisValues_Evaluator, p)
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)
PHX::MDField< ScalarT > normal_dot_flux_plus_pen