Panzer  Version of the Day
Panzer_Integrator_Scalar_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_SCALAR_IMPL_HPP
44 #define PANZER_EVALUATOR_SCALAR_IMPL_HPP
45 
46 #include "Intrepid2_FunctionSpaceTools.hpp"
49 #include "Kokkos_ViewFactory.hpp"
50 #include "Phalanx_DataLayout_MDALayout.hpp"
51 
52 namespace panzer {
53 
54 //**********************************************************************
55 PHX_EVALUATOR_CTOR(Integrator_Scalar,p) : quad_index(-1)
56 {
58  p.validateParameters(*valid_params);
59 
62 
63  Teuchos::RCP<PHX::DataLayout> dl_cell = Teuchos::rcp(new PHX::MDALayout<Cell>(ir->dl_scalar->dimension(0)));
64  integral = PHX::MDField<ScalarT>( p.get<std::string>("Integral Name"), dl_cell);
65  scalar = PHX::MDField<ScalarT,Cell,IP>( p.get<std::string>("Integrand Name"), ir->dl_scalar);
66 
67  this->addEvaluatedField(integral);
68  this->addDependentField(scalar);
69 
70  multiplier = 1.0;
71  if(p.isType<double>("Multiplier"))
72  multiplier = p.get<double>("Multiplier");
73 
74  if (p.isType<Teuchos::RCP<const std::vector<std::string> > >("Field Multipliers")) {
75  const std::vector<std::string>& field_multiplier_names =
76  *(p.get<Teuchos::RCP<const std::vector<std::string> > >("Field Multipliers"));
77 
78  for (std::vector<std::string>::const_iterator name =
79  field_multiplier_names.begin();
80  name != field_multiplier_names.end(); ++name) {
81  PHX::MDField<ScalarT,Cell,IP> tmp_field(*name, p.get< Teuchos::RCP<panzer::IntegrationRule> >("IR")->dl_scalar);
82  field_multipliers.push_back(tmp_field);
83  }
84  }
85 
86  for (typename std::vector<PHX::MDField<ScalarT,Cell,IP> >::iterator field = field_multipliers.begin();
87  field != field_multipliers.end(); ++field)
88  this->addDependentField(*field);
89 
90  std::string n = "Integrator_Scalar: " + integral.fieldTag().name();
91  this->setName(n);
92 }
93 
94 //**********************************************************************
95 PHX_POST_REGISTRATION_SETUP(Integrator_Scalar,sd,fm)
96 {
97  this->utils.setFieldData(integral,fm);
98  this->utils.setFieldData(scalar,fm);
99 
100  for (typename std::vector<PHX::MDField<ScalarT,Cell,IP> >::iterator field = field_multipliers.begin();
101  field != field_multipliers.end(); ++field)
102  this->utils.setFieldData(*field,fm);
103 
104  num_qp = scalar.dimension(1);
105 
106  tmp = Kokkos::createDynRankView(scalar.get_static_view(),"tmp", scalar.dimension(0), num_qp);
107 
108  quad_index = panzer::getIntegrationRuleIndex(quad_order,(*sd.worksets_)[0], this->wda);
109 }
110 
111 //**********************************************************************
112 PHX_EVALUATE_FIELDS(Integrator_Scalar,workset)
113 {
114 /*
115  for (index_t cell = 0; cell < workset.num_cells; ++cell)
116  integral(cell) = 0.0;
117 */
118 
119  for (index_t cell = 0; cell < workset.num_cells; ++cell) {
120  for (std::size_t qp = 0; qp < num_qp; ++qp) {
121  tmp(cell,qp) = multiplier * scalar(cell,qp);
122  for (typename std::vector<PHX::MDField<ScalarT,Cell,IP> >::iterator field = field_multipliers.begin();
123  field != field_multipliers.end(); ++field)
124  tmp(cell,qp) = tmp(cell,qp) * (*field)(cell,qp);
125  }
126  }
127 
128  // Switch to Kokkos means we can no longer use intrepid. There is a
129  // hard coded call to blas that grab a reference to the first
130  // element. You can't grab refrerences to memory anymore with
131  // kokkos. When Intrepid2 is fixed, we can switch this back.
132  /*
133  if(workset.num_cells>0)
134  Intrepid2::FunctionSpaceTools::
135  integrate<ScalarT>(integral, tmp,
136  (this->wda(workset).int_rules[quad_index])->weighted_measure,
137  Intrepid2::COMP_CPP);
138  */
139 
140  // NOTE: this is not portable to GPUs. Need to remove all uses of
141  // intrepid field container for MDFields. This is rather involved
142  // since we need to change the Worksets.
143 
144  // const Kokkos::DynRankView<double,PHX::Device>& rightFields = (this->wda(workset).int_rules[quad_index])->weighted_measure;
145  const IntegrationValues2<double> & iv = *this->wda(workset).int_rules[quad_index];
146 
147  int numPoints = tmp.dimension(1);
148 
149  for(index_t cl = 0; cl < workset.num_cells; cl++) {
150  integral(cl) = tmp(cl, 0)*iv.weighted_measure(cl, 0);
151  for(int qp = 1; qp < numPoints; qp++)
152  // integral(cl) += tmp(cl, qp)*rightFields(cl, qp);
153  integral(cl) += tmp(cl, qp)*iv.weighted_measure(cl, qp);
154  } // cl-loop
155 }
156 
157 //**********************************************************************
158 template<typename EvalT, typename TRAITS>
161 {
163  p->set<std::string>("Integral Name", "?");
164  p->set<std::string>("Integrand Name", "?");
165 
167  p->set("IR", ir);
168  p->set<double>("Multiplier", 1.0);
169 
171  p->set("Field Multipliers", fms);
172  return p;
173 }
174 
175 //**********************************************************************
176 
177 }
178 
179 #endif
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)
int numPoints
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::MDField< ScalarT > integral
std::size_t quad_index
PHX_EVALUATOR_CTOR(BasisValues_Evaluator, p)
PHX::MDField< const ScalarT, Cell, IP > field
double multiplier
PHX::MDField< const ScalarT, Cell, IP > scalar
PHX_EVALUATE_FIELDS(BasisValues_Evaluator, workset)
std::vector< std::string >::size_type getIntegrationRuleIndex(int ir_degree, panzer::Workset &workset, WorksetDetailsAccessor &wda)
std::vector< PHX::MDField< const ScalarT, Cell, IP > > field_multipliers
PHX_POST_REGISTRATION_SETUP(BasisValues_Evaluator, sd, fm)