WARPXM v1.10.0
Loading...
Searching...
No Matches
artificial_viscosity_limiter_historical.h
Go to the documentation of this file.
1#pragma once
2
3// Wm includes
5
6namespace wxm
7{
8namespace apps
9{
10namespace five_moment
11{
12namespace limiters
13{
14
27{
28public:
33
38
43 void setup(const WxCryptSet& wxc) override;
44
50 const std::vector<int>& getInputVariableIndexes(int flag) const override
51 {
52 return _input_variables;
53 }
54
60 const std::vector<int>& getAuxiliaryVariableIndexes(int flag) const override
61 {
62 return _aux_variables;
63 }
64
70 const std::vector<int>& getOutputVariableIndexes(int flag) const override
71 {
72 return _output_variables;
73 }
74
75
86 const real* q_r,
87 const real* aux_l,
88 const real* aux_r,
89 const solverVariables_t* pFV,
90 real* numericalFlux) const override;
91
100 const real* aux,
101 const solverVariables_t* pSV,
102 std::vector<std::vector<real>>& internalFlux) const override;
103
104 // specific calls
105
116 const real* q_r,
117 const real* aux_l,
118 const real* aux_r,
119 const solverVariables_t* pFV,
120 real* numericalFlux) const;
121
132 const real* q_r,
133 const real* aux_l,
134 const real* aux_r,
135 const solverVariables_t* pFV,
136 real* numericalFlux) const;
137
138protected:
139 // viscous multiplier of the limiter
141 // gradient calculation method ("ip" or "ldg")
142 std::string _gradient_method;
143 // diffusion factor multiplier ("fabs" or "squared")
144 std::string _diffusion_factor;
145
146 std::vector<int> _input_variables;
147 std::vector<int> _aux_variables;
148 std::vector<int> _output_variables;
149};
150
151} // namespace limiters
152} // namespace five_moment
153} // namespace apps
154} // namespace wxm
Base Class for physics applications.
Definition: wmapplication.h:93
WxCryptSet extends WxCrypt by providing, in addition to name-value pairs, an set of named WxCryptSets...
Definition: wxcryptset.h:35
Class used implementing artificial viscosity limiter based on velocity divergence Adds a diffusion te...
Definition: artificial_viscosity_limiter_historical.h:27
void setup(const WxCryptSet &wxc) override
Setup function.
real _coefficient
Definition: artificial_viscosity_limiter_historical.h:140
real numerical_flux_ip(const real *q_l, const real *q_r, const real *aux_l, const real *aux_r, const solverVariables_t *pFV, real *numericalFlux) const
Numerical flux call as determined by interior penalty method.
const std::vector< int > & getAuxiliaryVariableIndexes(int flag) const override
Grab auxiliary variable indices.
Definition: artificial_viscosity_limiter_historical.h:60
std::string _diffusion_factor
Definition: artificial_viscosity_limiter_historical.h:144
real numerical_flux(const real *q_l, const real *q_r, const real *aux_l, const real *aux_r, const solverVariables_t *pFV, real *numericalFlux) const override
Numerical flux redefinition from parent.
std::string _gradient_method
Definition: artificial_viscosity_limiter_historical.h:142
std::vector< int > _output_variables
Definition: artificial_viscosity_limiter_historical.h:148
const std::vector< int > & getOutputVariableIndexes(int flag) const override
Grab output variable indices.
Definition: artificial_viscosity_limiter_historical.h:70
std::vector< int > _input_variables
Definition: artificial_viscosity_limiter_historical.h:146
std::vector< int > _aux_variables
Definition: artificial_viscosity_limiter_historical.h:147
const std::vector< int > & getInputVariableIndexes(int flag) const override
Grab input variable indices.
Definition: artificial_viscosity_limiter_historical.h:50
real numerical_flux_ldg(const real *q_l, const real *q_r, const real *aux_l, const real *aux_r, const solverVariables_t *pFV, real *numericalFlux) const
Numerical flux call as determined by local discontinuous galerkin method.
real internal_flux(const real *q, const real *aux, const solverVariables_t *pSV, std::vector< std::vector< real > > &internalFlux) const override
Internal flux redefinition from parent.
list apps
Definition: shock_tube.py:33
Base namespace for everything not included in the global namespace.
Definition: field_source.h:8
Definition: wmapplication.h:38
#define real
Definition: wmoclunstructuredreconstruction.h:11