Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 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  std::string dof_name = p.get<std::string>("DOF Name");
67  std::string value_name = p.get<std::string>("Value Name");
68  std::string sigma_name = p.get<std::string>("Sigma Name");
69 
72 
75 
76 
77  residual = PHX::MDField<ScalarT>(residual_name, basis->functional);
78  normal_dot_flux_plus_pen = PHX::MDField<ScalarT>(normal_dot_flux_name, ir->dl_scalar);
79  flux = PHX::MDField<const ScalarT>(flux_name, ir->dl_vector);
80  normal = PHX::MDField<const ScalarT>(normal_name, ir->dl_vector);
81  dof = PHX::MDField<const ScalarT>(dof_name, ir->dl_scalar);
82  value = PHX::MDField<const ScalarT>(value_name, ir->dl_scalar);
83  sigma = PHX::MDField<const ScalarT>(sigma_name, ir->dl_scalar);
84 
85  this->addEvaluatedField(residual);
86  this->addEvaluatedField(normal_dot_flux_plus_pen);
87  this->addDependentField(normal);
88  this->addDependentField(flux);
89  this->addDependentField(dof);
90  this->addDependentField(value);
91  this->addDependentField(sigma);
92 
93  basis_name = panzer::basisIRLayout(basis,*ir)->name();
94 
95  std::string n = "Weak Dirichlet Residual Evaluator";
96  this->setName(n);
97 }
98 
99 //**********************************************************************
100 template<typename EvalT, typename Traits>
101 void
104  typename Traits::SetupData sd,
105  PHX::FieldManager<Traits>& /* fm */)
106 {
107  num_ip = flux.extent(1);
108  num_dim = flux.extent(2);
109 
110  TEUCHOS_ASSERT(flux.extent(1) == normal.extent(1));
111  TEUCHOS_ASSERT(flux.extent(2) == normal.extent(2));
112 
113  basis_index = panzer::getBasisIndex(basis_name, (*sd.worksets_)[0], this->wda);
114 }
115 
116 //**********************************************************************
117 template<typename EvalT, typename Traits>
118 void
121  typename Traits::EvalData workset)
122 {
123  for (index_t cell = 0; cell < workset.num_cells; ++cell) {
124  for (std::size_t ip = 0; ip < num_ip; ++ip) {
125  normal_dot_flux_plus_pen(cell,ip) = ScalarT(0.0);
126  for (std::size_t dim = 0; dim < num_dim; ++dim) {
127  normal_dot_flux_plus_pen(cell,ip) += normal(cell,ip,dim) * flux(cell,ip,dim);
128  }
129  normal_dot_flux_plus_pen(cell,ip) += sigma(cell,ip) * (dof(cell,ip) - value(cell,ip));
130  }
131  }
132 
133  if(workset.num_cells>0)
134  Intrepid2::FunctionSpaceTools<PHX::exec_space>::
135  integrate(residual.get_view(),
136  normal_dot_flux_plus_pen.get_view(),
137  (this->wda(workset).bases[basis_index])->weighted_basis_scalar.get_view());
138 
139 }
140 
141 //**********************************************************************
142 
143 }
144 
145 #endif
WeakDirichletResidual(const Teuchos::ParameterList &p)
T & get(const std::string &name, T def_value)
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
Teuchos::RCP< panzer::BasisIRLayout > basisIRLayout(std::string basis_type, const int basis_order, const PointRule &pt_rule)
Nonmember constructor.
void evaluateFields(typename Traits::EvalData d)
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< const std::vector< panzer::Workset > > worksets_