Panzer  Version of the Day
Panzer_Integrator_BasisTimesScalar_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_EVALUATOR_BASISTIMESSCALAR_IMPL_HPP
44 #define PANZER_EVALUATOR_BASISTIMESSCALAR_IMPL_HPP
45 
46 #include "Intrepid2_FunctionSpaceTools.hpp"
48 #include "Panzer_BasisIRLayout.hpp"
50 #include "Kokkos_ViewFactory.hpp"
51 
52 #define PANZER_USE_FAST_QUAD 1
53 // #define PANZER_USE_FAST_QUAD 0
54 
55 namespace panzer {
56 
57 //**********************************************************************
58 PHX_EVALUATOR_CTOR(Integrator_BasisTimesScalar,p) :
59  residual( p.get<std::string>("Residual Name"),
60  p.get< Teuchos::RCP<panzer::BasisIRLayout> >("Basis")->functional),
61  scalar( p.get<std::string>("Value Name"),
62  p.get< Teuchos::RCP<panzer::IntegrationRule> >("IR")->dl_scalar),
63  basis_name(p.get< Teuchos::RCP<panzer::BasisIRLayout> >("Basis")->name())
64 {
65  using Teuchos::RCP;
66 
68  p.validateParameters(*valid_params);
69 
71  = p.get< Teuchos::RCP<BasisIRLayout> >("Basis")->getBasis();
72 
73  // Verify that this basis supports the gradient operation
74  TEUCHOS_TEST_FOR_EXCEPTION(!basis->isScalarBasis(),std::logic_error,
75  "Integrator_BasisTimesScalar: Basis of type \"" << basis->name() << "\" is not "
76  "a scalar basis");
77 
78  this->addEvaluatedField(residual);
79  this->addDependentField(scalar);
80 
81  multiplier = p.get<double>("Multiplier");
82 
83 
84  // build field multpliers if vector is nonnull (defaults to null)
85  if(p.isType<Teuchos::RCP<const std::vector<std::string> > >("Field Multipliers")) {
86  RCP<const std::vector<std::string> > field_multiplier_names =
87  p.get<Teuchos::RCP<const std::vector<std::string> > >("Field Multipliers");
88  if(field_multiplier_names!=Teuchos::null) {
89  for (std::vector<std::string>::const_iterator name =
90  field_multiplier_names->begin();
91  name != field_multiplier_names->end(); ++name) {
92  PHX::MDField<ScalarT,Cell,IP> tmp_field(*name, p.get< Teuchos::RCP<panzer::IntegrationRule> >("IR")->dl_scalar);
93  field_multipliers.push_back(tmp_field);
94  }
95  }
96  }
97 
98  // add dependent field multiplers
99  for (typename std::vector<PHX::MDField<ScalarT,Cell,IP> >::iterator field = field_multipliers.begin();
100  field != field_multipliers.end(); ++field)
101  this->addDependentField(*field);
102 
103  std::string n = "Integrator_BasisTimesScalar: " + residual.fieldTag().name();
104  this->setName(n);
105 }
106 
107 //**********************************************************************
108 PHX_POST_REGISTRATION_SETUP(Integrator_BasisTimesScalar,sd,fm)
109 {
110  this->utils.setFieldData(residual,fm);
111  this->utils.setFieldData(scalar,fm);
112 
113  for (typename std::vector<PHX::MDField<ScalarT,Cell,IP> >::iterator field = field_multipliers.begin();
114  field != field_multipliers.end(); ++field)
115  this->utils.setFieldData(*field,fm);
116 
117  num_nodes = residual.dimension(1);
118  num_qp = scalar.dimension(1);
119 
120  basis_index = panzer::getBasisIndex(basis_name, (*sd.worksets_)[0], this->wda);
121 
122  tmp = Kokkos::createDynRankView(scalar.get_static_view(),"tmp",scalar.dimension(0), num_qp);
123 }
124 
125 //**********************************************************************
126 PHX_EVALUATE_FIELDS(Integrator_BasisTimesScalar,workset)
127 {
128  // for (int i=0; i < residual.size(); ++i)
129  // residual[i] = 0.0;
130  // Irina modified
131  residual.deep_copy(ScalarT(0.0));
132 
133 #if PANZER_USE_FAST_QUAD
134  // do a scaled copy
135  // Irina modified
136  // for (int i=0; i < scalar.size(); ++i)
137  // tmp[i] = multiplier * scalar[i];
138  for (int i=0; i < scalar.extent_int(0); ++i)
139  for (int j=0; j < scalar.extent_int(1); ++j)
140  tmp(i,j) = multiplier * scalar(i,j);
141 
142  for (typename std::vector<PHX::MDField<ScalarT,Cell,IP> >::iterator field = field_multipliers.begin();
143  field != field_multipliers.end(); ++field) {
144  const PHX::MDField<ScalarT,Cell,IP> & field_data = *field;
145 
146  //Irina modified
147  //for (int i=0; i < field_data.size(); ++i)
148  // tmp[i] *= field_data[i];
149  for (int i=0; i < scalar.extent_int(0); ++i)
150  for (int j=0; j < scalar.extent_int(1); ++j)
151  tmp(i,j) = tmp(i,j) * field_data(i,j);
152  }
153 
154  // const Kokkos::DynRankView<double,PHX::Device> & weighted_basis = this->wda(workset).bases[basis_index]->weighted_basis;
155  const BasisValues2<double> & bv = *this->wda(workset).bases[basis_index];
156 
157  for (index_t cell = 0; cell < workset.num_cells; ++cell) {
158  for (std::size_t basis = 0; basis < num_nodes; ++basis) {
159  for (std::size_t qp = 0; qp < num_qp; ++qp) {
160  residual(cell,basis) += tmp(cell,qp)*bv.weighted_basis_scalar(cell,basis,qp);
161  }
162  }
163  }
164 
165 #else
166  for (index_t cell = 0; cell < workset.num_cells; ++cell) {
167  for (std::size_t qp = 0; qp < num_qp; ++qp) {
168  tmp(cell,qp) = multiplier * scalar(cell,qp);
169  for (typename std::vector<PHX::MDField<ScalarT,Cell,IP> >::iterator field = field_multipliers.begin();
170  field != field_multipliers.end(); ++field)
171  tmp(cell,qp) = tmp(cell,qp) * (*field)(cell,qp);
172  }
173  }
174 
175  if(workset.num_cells>0)
176  Intrepid2::FunctionSpaceTools::
177  integrate<ScalarT>(residual, tmp,
178  (this->wda(workset).bases[basis_index])->weighted_basis,
179  Intrepid2::COMP_CPP);
180 #endif
181 }
182 
183 //**********************************************************************
184 template<typename EvalT, typename TRAITS>
187 {
189  p->set<std::string>("Residual Name", "?");
190  p->set<std::string>("Value Name", "?");
192  p->set("Basis", basis);
194  p->set("IR", ir);
195  p->set<double>("Multiplier", 1.0);
197  p->set("Field Multipliers", fms);
198  return p;
199 }
200 
201 //**********************************************************************
202 
203 }
204 
205 #endif
206 
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 name() const
A unique key that is the combination of the basis type and basis order.
std::string basis_name
PHX::MDField< ScalarT, Cell > tmp
std::size_t num_qp
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Array_CellBasisIP weighted_basis_scalar
T * get() const
PHX::MDField< ScalarT > vector
Kokkos::DynRankView< typename InputArray::value_type, PHX::Device > createDynRankView(const InputArray &a, const std::string &name, const DimensionPack... dims)
Wrapper to simplify Panzer use of Sacado ViewFactory.
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Teuchos::RCP< PHX::DataLayout > dl_scalar
Data layout for scalar fields.
PHX_EVALUATOR_CTOR(BasisValues_Evaluator, p)
PHX::MDField< const ScalarT, Cell, IP > field
bool isScalarBasis() const
double multiplier
PHX::MDField< const ScalarT, Cell, IP > scalar
PHX_EVALUATE_FIELDS(BasisValues_Evaluator, workset)
Teuchos::RCP< const panzer::PureBasis > basis
Interpolates basis DOF values to IP DOF values.
std::vector< PHX::MDField< const ScalarT, Cell, IP > > field_multipliers
PHX_POST_REGISTRATION_SETUP(BasisValues_Evaluator, sd, fm)