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 // Panzer: A partial differential equation assembly
4 // engine for strongly coupled complex multiphysics systems
5 //
6 // Copyright 2011 NTESS and the Panzer contributors.
7 // SPDX-License-Identifier: BSD-3-Clause
8 // *****************************************************************************
9 // @HEADER
10 
11 #ifndef PANZER_EVALUATOR_SCALAR_IMPL_HPP
12 #define PANZER_EVALUATOR_SCALAR_IMPL_HPP
13 
14 #include "Intrepid2_FunctionSpaceTools.hpp"
17 #include "Kokkos_ViewFactory.hpp"
18 #include "Phalanx_DataLayout_MDALayout.hpp"
19 
20 namespace panzer {
21 
22 //**********************************************************************
23 template<typename EvalT, typename Traits>
26  const Teuchos::ParameterList& p) : quad_index(0)
27 {
29  p.validateParameters(*valid_params);
30 
33 
35  integral = PHX::MDField<ScalarT>( p.get<std::string>("Integral Name"), dl_cell);
36  scalar = PHX::MDField<const ScalarT,Cell,IP>( p.get<std::string>("Integrand Name"), ir->dl_scalar);
37 
38  this->addEvaluatedField(integral);
39  this->addDependentField(scalar);
40 
41  multiplier = 1.0;
42  if(p.isType<double>("Multiplier"))
43  multiplier = p.get<double>("Multiplier");
44 
45  if (p.isType<Teuchos::RCP<const std::vector<std::string> > >("Field Multipliers")) {
46  const std::vector<std::string>& field_multiplier_names =
47  *(p.get<Teuchos::RCP<const std::vector<std::string> > >("Field Multipliers"));
48 
49  field_multipliers_h = typename PHX::View<PHX::UnmanagedView<const ScalarT**>* >::HostMirror("FieldMultipliersHost", field_multiplier_names.size());
50  field_multipliers = PHX::View<PHX::UnmanagedView<const ScalarT**>* >("FieldMultipliersDevice", field_multiplier_names.size());
51 
52  int cnt=0;
53  for (std::vector<std::string>::const_iterator name =
54  field_multiplier_names.begin();
55  name != field_multiplier_names.end(); ++name, ++cnt) {
57  this->addDependentField(tmp_field.fieldTag(),field_multipliers_h(cnt));
58  }
59  }
60 
61  std::string n = "Integrator_Scalar: " + integral.fieldTag().name();
62  this->setName(n);
63 }
64 
65 //**********************************************************************
66 template<typename EvalT, typename Traits>
67 void
70  typename Traits::SetupData sd,
71  PHX::FieldManager<Traits>& /* fm */)
72 {
73  num_qp = scalar.extent(1);
74  tmp = Kokkos::createDynRankView(scalar.get_static_view(),"tmp", scalar.extent(0), num_qp);
75  quad_index = panzer::getIntegrationRuleIndex(quad_order,(*sd.worksets_)[0], this->wda);
76 
77  Kokkos::deep_copy(field_multipliers, field_multipliers_h);
78 }
79 
80 
81 //**********************************************************************
82 template<typename EvalT, typename Traits>
83 void
86  typename Traits::EvalData workset)
87 {
88 
89  const auto wm = this->wda(workset).int_rules[quad_index]->weighted_measure;
90 
91  auto l_tmp = tmp;
92  auto l_multiplier = multiplier;
93  auto l_scalar = scalar;
94  auto l_field_multipliers = field_multipliers;
95  auto l_num_qp = num_qp;
96  auto l_integral = integral.get_static_view();
97 
98  Kokkos::parallel_for("IntegratorScalar", workset.num_cells, KOKKOS_LAMBDA (int cell) {
99  for (std::size_t qp = 0; qp < l_num_qp; ++qp) {
100  l_tmp(cell,qp) = l_multiplier * l_scalar(cell,qp);
101  for (size_t f=0;f<l_field_multipliers.extent(0); ++f)
102  l_tmp(cell,qp) *= l_field_multipliers(f)(cell,qp);
103  }
104  l_integral(cell) = 0.0;
105  for (std::size_t qp = 0; qp < l_num_qp; ++qp) {
106  l_integral(cell) += l_tmp(cell, qp)*wm(cell, qp);
107  }
108  } );
109  Kokkos::fence();
110 }
111 
112 //**********************************************************************
113 template<typename EvalT, typename TRAITS>
116 {
118  p->set<std::string>("Integral Name", "?");
119  p->set<std::string>("Integrand Name", "?");
120 
122  p->set("IR", ir);
123  p->set<double>("Multiplier", 1.0);
124 
126  p->set("Field Multipliers", fms);
127  return p;
128 }
129 
130 //**********************************************************************
131 
132 }
133 
134 #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.
int num_cells
DEPRECATED - use: numCells()
T & get(const std::string &name, T def_value)
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
std::vector< int >::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)
PHX::View< PHX::UnmanagedView< const ScalarT ** > * >::HostMirror field_multipliers_h
Teuchos::RCP< PHX::DataLayout > dl_scalar
Data layout for scalar fields.
PHX::View< PHX::UnmanagedView< const ScalarT ** > * > field_multipliers
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
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_