Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_CellAverage_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_CellAverage_impl_hpp__
12#define __Panzer_CellAverage_impl_hpp__
13
14#include <limits>
15
16#include "Intrepid2_FunctionSpaceTools.hpp"
19#include "Phalanx_DataLayout_MDALayout.hpp"
20
21namespace panzer {
22
23//**********************************************************************
24template<typename EvalT, typename Traits>
27 const Teuchos::ParameterList& p) : quad_index(0)
28{
29 Teuchos::RCP<Teuchos::ParameterList> valid_params = this->getValidParameters();
30 p.validateParameters(*valid_params);
31
32 Teuchos::RCP<panzer::IntegrationRule> ir = p.get< Teuchos::RCP<panzer::IntegrationRule> >("IR");
33 quad_order = ir->cubature_degree;
34
35 Teuchos::RCP<PHX::DataLayout> dl_cell = Teuchos::rcp(new PHX::MDALayout<Cell>(ir->dl_scalar->extent(0)));
36 average = PHX::MDField<ScalarT,Cell>( p.get<std::string>("Average Name"), dl_cell);
37 scalar = PHX::MDField<const ScalarT,Cell,IP>( p.get<std::string>("Field Name"), ir->dl_scalar);
38
39 this->addEvaluatedField(average);
40 this->addDependentField(scalar);
41
42 multiplier = 1.0;
43 if(p.isType<double>("Multiplier"))
44 multiplier = p.get<double>("Multiplier");
45
46 if (p.isType<Teuchos::RCP<const std::vector<std::string> > >("Field Multipliers")) {
47 const std::vector<std::string>& field_multiplier_names =
48 *(p.get<Teuchos::RCP<const std::vector<std::string> > >("Field Multipliers"));
49
50 for (std::vector<std::string>::const_iterator name =
51 field_multiplier_names.begin();
52 name != field_multiplier_names.end(); ++name) {
53 PHX::MDField<const ScalarT,Cell,IP> tmp_field(*name, p.get< Teuchos::RCP<panzer::IntegrationRule> >("IR")->dl_scalar);
54 field_multipliers.push_back(tmp_field);
55 }
56 }
57
58 for (typename std::vector<PHX::MDField<const ScalarT,Cell,IP> >::iterator field = field_multipliers.begin();
59 field != field_multipliers.end(); ++field)
60 this->addDependentField(*field);
61
62 std::string n = "CellAverage: " + average.fieldTag().name();
63 this->setName(n);
64}
65
66//**********************************************************************
67template<typename EvalT, typename Traits>
68void
71 typename Traits::SetupData sd,
73{
74 num_qp = scalar.extent(1);
75 quad_index = panzer::getIntegrationRuleIndex(quad_order,(*sd.worksets_)[0], this->wda);
76}
77
78//**********************************************************************
79template<typename EvalT, typename Traits>
80void
83 typename Traits::EvalData workset)
84{
85 for (index_t cell = 0; cell < workset.num_cells; ++cell) {
86
87 // start with no average
88 average(cell) = 0.0;
89
90 // loop over quadrture points, compute simple average
91 for (std::size_t qp = 0; qp < num_qp; ++qp) {
92 ScalarT current= multiplier * scalar(cell,qp);
93 for (typename std::vector<PHX::MDField<const ScalarT,Cell,IP> >::iterator field = field_multipliers.begin();
94 field != field_multipliers.end(); ++field)
95 current *= (*field)(cell,qp);
96
97 // take first quad point value
98 average(cell) += current/num_qp;
99 }
100 }
101}
102
103//**********************************************************************
104template<typename EvalT, typename TRAITS>
105Teuchos::RCP<Teuchos::ParameterList>
107{
108 Teuchos::RCP<Teuchos::ParameterList> p = Teuchos::rcp(new Teuchos::ParameterList);
109 p->set<std::string>("Average Name", "?");
110 p->set<std::string>("Field Name", "?");
111
112 Teuchos::RCP<panzer::IntegrationRule> ir;
113 p->set("IR", ir);
114 p->set<double>("Multiplier", 1.0);
115
116 Teuchos::RCP<const std::vector<std::string> > fms;
117 p->set("Field Multipliers", fms);
118 return p;
119}
120
121//**********************************************************************
122
123}
124
125#endif
double multiplier
The scalar multiplier out in front of the integral ( ).
PHX::MDField< ScalarT, panzer::Cell, panzer::BASIS > field
A field to which we'll contribute, or in which we'll store, the result of computing this integral.
typename EvalT::ScalarT ScalarT
void evaluateFields(typename Traits::EvalData d)
PHX::MDField< const ScalarT, Cell, IP > scalar
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
std::vector< PHX::MDField< const ScalarT, Cell, IP > > field_multipliers
CellAverage(const Teuchos::ParameterList &p)
PHX::MDField< ScalarT, Cell > average
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
int num_cells
DEPRECATED - use: numCells()
std::vector< int >::size_type getIntegrationRuleIndex(int ir_degree, const panzer::Workset &workset, WorksetDetailsAccessor &wda)
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_