Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 template<typename EvalT, typename Traits>
58  const Teuchos::ParameterList& p) : quad_index(-1)
59 {
61  p.validateParameters(*valid_params);
62 
65 
66  Teuchos::RCP<PHX::DataLayout> dl_cell = Teuchos::rcp(new PHX::MDALayout<Cell>(ir->dl_scalar->extent(0)));
67  integral = PHX::MDField<ScalarT>( p.get<std::string>("Integral Name"), dl_cell);
68  scalar = PHX::MDField<const ScalarT,Cell,IP>( p.get<std::string>("Integrand Name"), ir->dl_scalar);
69 
70  this->addEvaluatedField(integral);
71  this->addDependentField(scalar);
72 
73  multiplier = 1.0;
74  if(p.isType<double>("Multiplier"))
75  multiplier = p.get<double>("Multiplier");
76 
77  if (p.isType<Teuchos::RCP<const std::vector<std::string> > >("Field Multipliers")) {
78  const std::vector<std::string>& field_multiplier_names =
79  *(p.get<Teuchos::RCP<const std::vector<std::string> > >("Field Multipliers"));
80 
81  for (std::vector<std::string>::const_iterator name =
82  field_multiplier_names.begin();
83  name != field_multiplier_names.end(); ++name) {
84  PHX::MDField<const ScalarT,Cell,IP> tmp_field(*name, p.get< Teuchos::RCP<panzer::IntegrationRule> >("IR")->dl_scalar);
85  field_multipliers.push_back(tmp_field);
86  }
87  }
88 
89  for (typename std::vector<PHX::MDField<const ScalarT,Cell,IP> >::iterator field = field_multipliers.begin();
90  field != field_multipliers.end(); ++field)
91  this->addDependentField(*field);
92 
93  std::string n = "Integrator_Scalar: " + integral.fieldTag().name();
94  this->setName(n);
95 }
96 
97 //**********************************************************************
98 template<typename EvalT, typename Traits>
99 void
102  typename Traits::SetupData sd,
103  PHX::FieldManager<Traits>& /* fm */)
104 {
105  num_qp = scalar.extent(1);
106  tmp = Kokkos::createDynRankView(scalar.get_static_view(),"tmp", scalar.extent(0), num_qp);
107  quad_index = panzer::getIntegrationRuleIndex(quad_order,(*sd.worksets_)[0], this->wda);
108 }
109 
110 //**********************************************************************
111 template<typename EvalT, typename Traits>
112 void
115  typename Traits::EvalData workset)
116 {
117 /*
118  for (index_t cell = 0; cell < workset.num_cells; ++cell)
119  integral(cell) = 0.0;
120 */
121 
122  for (index_t cell = 0; cell < workset.num_cells; ++cell) {
123  for (std::size_t qp = 0; qp < num_qp; ++qp) {
124  tmp(cell,qp) = multiplier * scalar(cell,qp);
125  for (typename std::vector<PHX::MDField<const ScalarT,Cell,IP> >::iterator field = field_multipliers.begin();
126  field != field_multipliers.end(); ++field)
127  tmp(cell,qp) = tmp(cell,qp) * (*field)(cell,qp);
128  }
129  }
130 
131  // Switch to Kokkos means we can no longer use intrepid. There is a
132  // hard coded call to blas that grab a reference to the first
133  // element. You can't grab refrerences to memory anymore with
134  // kokkos. When Intrepid2 is fixed, we can switch this back.
135  /*
136  if(workset.num_cells>0)
137  Intrepid2::FunctionSpaceTools::
138  integrate<ScalarT>(integral, tmp,
139  (this->wda(workset).int_rules[quad_index])->weighted_measure,
140  Intrepid2::COMP_CPP);
141  */
142 
143  // NOTE: this is not portable to GPUs. Need to remove all uses of
144  // intrepid field container for MDFields. This is rather involved
145  // since we need to change the Worksets.
146 
147  // const Kokkos::DynRankView<double,PHX::Device>& rightFields = (this->wda(workset).int_rules[quad_index])->weighted_measure;
148  const IntegrationValues2<double> & iv = *this->wda(workset).int_rules[quad_index];
149 
150  int numPoints = tmp.extent(1);
151 
152  for(index_t cl = 0; cl < workset.num_cells; cl++) {
153  integral(cl) = tmp(cl, 0)*iv.weighted_measure(cl, 0);
154  for(int qp = 1; qp < numPoints; qp++)
155  // integral(cl) += tmp(cl, qp)*rightFields(cl, qp);
156  integral(cl) += tmp(cl, qp)*iv.weighted_measure(cl, qp);
157  } // cl-loop
158 }
159 
160 //**********************************************************************
161 template<typename EvalT, typename TRAITS>
164 {
166  p->set<std::string>("Integral Name", "?");
167  p->set<std::string>("Integrand Name", "?");
168 
170  p->set("IR", ir);
171  p->set<double>("Multiplier", 1.0);
172 
174  p->set("Field Multipliers", fms);
175  return p;
176 }
177 
178 //**********************************************************************
179 
180 }
181 
182 #endif
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.
T & get(const std::string &name, T def_value)
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
int numPoints
std::vector< std::string >::size_type getIntegrationRuleIndex(int ir_degree, const panzer::Workset &workset, WorksetDetailsAccessor &wda)
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
double multiplier
The scalar multiplier out in front of the integral ( ).
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Teuchos::RCP< PHX::DataLayout > dl_scalar
Data layout for scalar fields.
Integrator_Scalar(const Teuchos::ParameterList &p)
PHX::MDField< const ScalarT, Cell, IP > scalar
void validateParameters(ParameterList const &validParamList, int const depth=1000, EValidateUsed const validateUsed=VALIDATE_USED_ENABLED, EValidateDefaults const validateDefaults=VALIDATE_DEFAULTS_ENABLED) const
PHX::MDField< ScalarT, panzer::Cell, panzer::BASIS > field
A field to which we&#39;ll contribute, or in which we&#39;ll store, the result of computing this integral...
std::vector< PHX::MDField< const ScalarT, Cell, IP > > field_multipliers
bool isType(const std::string &name) const
void evaluateFields(typename Traits::EvalData d)
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_