Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_CellExtreme_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_CellExtreme_impl_hpp__
12 #define __Panzer_CellExtreme_impl_hpp__
13 
14 #include <limits>
15 
16 #include "Intrepid2_FunctionSpaceTools.hpp"
19 #include "Phalanx_DataLayout_MDALayout.hpp"
20 
21 namespace panzer {
22 
23 //**********************************************************************
24 template<typename EvalT, typename Traits>
27  const Teuchos::ParameterList& p) : quad_index(0)
28 {
30  p.validateParameters(*valid_params);
31 
32  // setup default value
33  use_max = true;
34  if(p.isType<bool>("Use Max"))
35  use_max = p.get<bool>("Use Max");
36 
38  quad_order = ir->cubature_degree;
39 
40  Teuchos::RCP<PHX::DataLayout> dl_cell = Teuchos::rcp(new PHX::MDALayout<Cell>(ir->dl_scalar->extent(0)));
41  extreme = PHX::MDField<ScalarT>( p.get<std::string>("Extreme Name"), dl_cell);
42  scalar = PHX::MDField<const ScalarT,Cell,IP>( p.get<std::string>("Field Name"), ir->dl_scalar);
43 
44  this->addEvaluatedField(extreme);
45  this->addDependentField(scalar);
46 
47  multiplier = 1.0;
48  if(p.isType<double>("Multiplier"))
49  multiplier = p.get<double>("Multiplier");
50 
51  if (p.isType<Teuchos::RCP<const std::vector<std::string> > >("Field Multipliers")) {
52  const std::vector<std::string>& field_multiplier_names =
53  *(p.get<Teuchos::RCP<const std::vector<std::string> > >("Field Multipliers"));
54 
55  for (std::vector<std::string>::const_iterator name =
56  field_multiplier_names.begin();
57  name != field_multiplier_names.end(); ++name) {
59  field_multipliers.push_back(tmp_field);
60  }
61  }
62 
63  for (typename std::vector<PHX::MDField<const ScalarT,Cell,IP> >::iterator field = field_multipliers.begin();
64  field != field_multipliers.end(); ++field)
65  this->addDependentField(*field);
66 
67  std::string n = "CellExtreme: " + extreme.fieldTag().name();
68  this->setName(n);
69 }
70 
71 //**********************************************************************
72 template<typename EvalT, typename Traits>
73 void
76  typename Traits::SetupData sd,
77  PHX::FieldManager<Traits>& /* fm */)
78 {
79  num_qp = scalar.extent(1);
80  quad_index = panzer::getIntegrationRuleIndex(quad_order,(*sd.worksets_)[0], this->wda);
81 }
82 
83 //**********************************************************************
84 template<typename EvalT, typename Traits>
85 void
88  typename Traits::EvalData workset)
89 {
90  double mult = this->multiplier;
91  std::size_t num_pt = this->num_qp;
92  bool extreme_max = this->use_max;
93  auto scalar_view = scalar.get_view();
94  auto extreme_view = extreme.get_view();
95 
96  // compute field weighted with multipliers
97  PHX::View<ScalarT**> mult_field("Multiplied Field", workset.num_cells, scalar.extent(1));
98 
99  // initialize to scalar value
100  Kokkos::parallel_for("Initialize Muliplier Field", workset.num_cells, KOKKOS_LAMBDA( const int cell) {
101  for (std::size_t qp = 0; qp < num_pt; ++qp) {
102  mult_field(cell, qp) = mult * scalar_view(cell, qp);
103  }
104  });
105 
106  // multiply by field values
107  for (std::size_t field_num = 0; field_num < field_multipliers.size(); ++field_num)
108  {
109  auto field = field_multipliers[field_num];
110  Kokkos::parallel_for("CellExtreme: Multiply Fields", workset.num_cells, KOKKOS_LAMBDA( const int cell) {
111  for (std::size_t qp = 0; qp < num_pt; ++qp) {
112  mult_field(cell, qp) *= field(cell, qp);
113  }
114  });
115  }
116 
117  // take extreme over points in each cell
118  if (extreme_max) {
119  Kokkos::parallel_for ("CellExtreme (max)", workset.num_cells, KOKKOS_LAMBDA( const int cell) {
120  for (std::size_t qp = 0; qp < num_pt; ++qp) {
121  auto current = mult_field(cell, qp);
122 
123  // take first point
124  if (qp == 0)
125  extreme_view(cell) = current;
126  else
127  extreme_view(cell) = extreme_view(cell)<current ? current : extreme_view(cell);
128  }
129  });
130  } else {
131  Kokkos::parallel_for ("CellExtreme (min)", workset.num_cells, KOKKOS_LAMBDA( const int cell) {
132  for (std::size_t qp = 0; qp < num_pt; ++qp) {
133  auto current = mult_field(cell, qp);
134 
135  // take first point
136  if (qp == 0)
137  extreme_view(cell) = current;
138  else // use_min
139  extreme_view(cell) = extreme_view(cell)>current ? current : extreme_view(cell);
140  }
141  });
142  }
143 }
144 
145 //**********************************************************************
146 template<typename EvalT, typename TRAITS>
149 {
151  p->set<std::string>("Extreme Name", "?");
152  p->set<std::string>("Field Name", "?");
153 
155  p->set("IR", ir);
156  p->set<bool>("Use Max", true);
157  p->set<double>("Multiplier", 1.0);
158 
160  p->set("Field Multipliers", fms);
161  return p;
162 }
163 
164 //**********************************************************************
165 
166 }
167 
168 #endif
int num_cells
DEPRECATED - use: numCells()
void evaluateFields(typename Traits::EvalData d)
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)
CellExtreme(const Teuchos::ParameterList &p)
std::vector< PHX::MDField< const ScalarT, Cell, IP > > field_multipliers
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)
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
PHX::MDField< ScalarT > extreme
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...
bool isType(const std::string &name) const
PHX::MDField< const ScalarT, Cell, IP > scalar
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_