Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_CellAverage.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_hpp__
12#define __Panzer_CellAverage_hpp__
13
14#include <string>
15#include "Panzer_Dimension.hpp"
16#include "Phalanx_Evaluator_Macros.hpp"
17#include "Phalanx_MDField.hpp"
18#include "Kokkos_DynRankView.hpp"
19
21
22namespace panzer {
23
37template<typename EvalT, typename Traits>
39 :
40 public panzer::EvaluatorWithBaseImpl<Traits>,
41 public PHX::EvaluatorDerived<EvalT, Traits>
42{
43 public:
44
46 const Teuchos::ParameterList& p);
47
48 void
50 typename Traits::SetupData d,
52
53 void
55 typename Traits::EvalData d);
56
57 private:
58
59 using ScalarT = typename EvalT::ScalarT;
60
61 PHX::MDField<ScalarT,Cell> average; // result
62
63 PHX::MDField<const ScalarT,Cell,IP> scalar; // function to be integrated
64
65 std::vector<PHX::MDField<const ScalarT,Cell,IP> > field_multipliers;
66 double multiplier;
67
68 std::size_t num_qp;
69 std::size_t quad_index;
71
72public:
73 // for testing purposes
74 const PHX::FieldTag & getFieldTag() const
75 { return average.fieldTag(); }
76
77private:
78 Teuchos::RCP<Teuchos::ParameterList> getValidParameters() const;
79
80}; // end of class CellAverage
81
82
94template <typename EvalT,typename Traits>
95Teuchos::RCP<PHX::Evaluator<Traits> > cellAverageEvaluator(const std::string & averageName,
96 const std::string & fieldName,
97 const Teuchos::RCP<const panzer::IntegrationRule> & ir)
98{
99 using Teuchos::RCP;
100 using Teuchos::rcp;
101 using Teuchos::rcp_const_cast;
102
103 Teuchos::ParameterList input;
104 input.set("Average Name",averageName);
105 input.set("Field Name",fieldName);
106 input.set("IR",rcp_const_cast<panzer::IntegrationRule>(ir));
107
108 return rcp(new CellAverage<EvalT,Traits>(input));
109}
110
111}
112
113#endif
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
const PHX::FieldTag & getFieldTag() const
PHX::MDField< ScalarT, Cell > average
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
Wrapper to PHX::EvaluatorWithBaseImpl that implements Panzer-specific helpers.
Teuchos::RCP< PHX::Evaluator< Traits > > cellAverageEvaluator(const std::string &averageName, const std::string &fieldName, const Teuchos::RCP< const panzer::IntegrationRule > &ir)