DAkOmegaSST.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v3
5 
6  Description:
7  Child class for the kOmegaSST model
8 
9  This file is modified from OpenFOAM's source code
10  src/TurbulenceModels/turbulenceModels/RAS/kOmegaSST/kOmegaSST.H
11 
12  OpenFOAM: The Open Source CFD Toolbox
13 
14  Copyright (C): 2011-2016 OpenFOAM Foundation
15 
16  OpenFOAM License:
17 
18  OpenFOAM is free software: you can redistribute it and/or modify it
19  under the terms of the GNU General Public License as published by
20  the Free Software Foundation, either version 3 of the License, or
21  (at your option) any later version.
22 
23  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
24  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
25  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
26  for more details.
27 
28  You should have received a copy of the GNU General Public License
29  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
30 
31 \*---------------------------------------------------------------------------*/
32 
33 #ifndef DAkOmegaSST_H
34 #define DAkOmegaSST_H
35 
36 #include "DATurbulenceModel.H"
37 #include "addToRunTimeSelectionTable.H"
38 
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43 
44 /*---------------------------------------------------------------------------*\
45  Class DAkOmegaSST Declaration
46 \*---------------------------------------------------------------------------*/
47 
49  : public DATurbulenceModel
50 {
51 
52 protected:
54 
55  dimensionedScalar alphaK1_;
56  dimensionedScalar alphaK2_;
57 
58  dimensionedScalar alphaOmega1_;
59  dimensionedScalar alphaOmega2_;
60 
61  dimensionedScalar gamma1_;
62  dimensionedScalar gamma2_;
63 
64  dimensionedScalar beta1_;
65  dimensionedScalar beta2_;
66 
67  dimensionedScalar betaStar_;
68 
69  dimensionedScalar a1_;
70  dimensionedScalar b1_;
71  dimensionedScalar c1_;
72 
73  Switch F3_;
75 
77 
78  tmp<volScalarField> F1(const volScalarField& CDkOmega) const;
79  tmp<volScalarField> F2() const;
80  tmp<volScalarField> F3() const;
81  tmp<volScalarField> F23() const;
82 
83  tmp<volScalarField> blend(
84  const volScalarField& F1,
85  const dimensionedScalar& psi1,
86  const dimensionedScalar& psi2) const
87  {
88  return F1 * (psi1 - psi2) + psi2;
89  }
90 
91  tmp<volScalarField::Internal> blend(
92  const volScalarField::Internal& F1,
93  const dimensionedScalar& psi1,
94  const dimensionedScalar& psi2) const
95  {
96  return F1 * (psi1 - psi2) + psi2;
97  }
98 
99  tmp<volScalarField> alphaK(const volScalarField& F1) const
100  {
101  return blend(F1, alphaK1_, alphaK2_);
102  }
103 
104  tmp<volScalarField> alphaOmega(const volScalarField& F1) const
105  {
106  return blend(F1, alphaOmega1_, alphaOmega2_);
107  }
108 
109  tmp<volScalarField::Internal> beta(
110  const volScalarField::Internal& F1) const
111  {
112  return blend(F1, beta1_, beta2_);
113  }
114 
115  tmp<volScalarField::Internal> gamma(
116  const volScalarField::Internal& F1) const
117  {
118  return blend(F1, gamma1_, gamma2_);
119  }
120 
121  //- Return the effective diffusivity for k
122  tmp<volScalarField> DkEff(const volScalarField& F1) const
123  {
124  return tmp<volScalarField>(
125  new volScalarField("DkEff", alphaK(F1) * nut_ + this->nu()));
126  }
127 
128  //- Return the effective diffusivity for omega
129  tmp<volScalarField> DomegaEff(const volScalarField& F1) const
130  {
131  return tmp<volScalarField>(
132  new volScalarField(
133  "DomegaEff",
134  alphaOmega(F1) * nut_ + this->nu()));
135  }
136 
137  //- Return k production rate
138  tmp<volScalarField::Internal> Pk(
139  const volScalarField::Internal& G) const;
140 
141  //- Return epsilon/k which for standard RAS is betaStar*omega
142  tmp<volScalarField::Internal> epsilonByk(
143  const volScalarField& F1,
144  const volTensorField& gradU) const;
145 
146  //- Return G/nu
147  tmp<volScalarField::Internal> GbyNu(
148  const volScalarField::Internal& GbyNu0,
149  const volScalarField::Internal& F2,
150  const volScalarField::Internal& S2) const;
151 
152  tmp<fvScalarMatrix> kSource() const;
153 
154  tmp<fvScalarMatrix> omegaSource() const;
155 
156  tmp<fvScalarMatrix> Qsas(
157  const volScalarField::Internal& S2,
158  const volScalarField::Internal& gamma,
159  const volScalarField::Internal& beta) const;
161 
163 
164  volScalarField& omega_;
165  volScalarField omegaRes_;
166  volScalarField& k_;
167  volScalarField kRes_;
169 
171  const volScalarField& y_;
172 
177  scalarList omegaNearWall_;
178 
180  label solveTurbState_ = 0;
181 
184 
185 public:
186  TypeName("kOmegaSST");
187  // Constructors
188 
189  //- Construct from components
190  DAkOmegaSST(
191  const word modelType,
192  const fvMesh& mesh,
193  const DAOption& daOption);
194 
195  //- Destructor
196  virtual ~DAkOmegaSST()
197  {
198  }
199 
200  // Member functions
201 
203  virtual void correctModelStates(wordList& modelStates) const;
204 
206  virtual void correctNut();
207 
209  virtual void correctBoundaryConditions();
210 
212  virtual void updateIntermediateVariables();
213 
215  virtual void correctStateResidualModelCon(List<List<word>>& stateCon) const;
216 
218  virtual void addModelResidualCon(HashTable<List<List<word>>>& allCon) const;
219 
221  virtual void calcResiduals(const dictionary& options);
222 
224  virtual void correct();
225 
227  void saveOmegaNearWall();
228 
230  void setOmegaNearWall();
231 
234 };
235 
236 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
237 
238 } // End namespace Foam
239 
240 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
241 
242 #endif
243 
244 // ************************************************************************* //
Foam::DAkOmegaSST::c1_
dimensionedScalar c1_
Definition: DAkOmegaSST.H:71
Foam::DAkOmegaSST::printInterval_
label printInterval_
time step interval to print residual
Definition: DAkOmegaSST.H:183
Foam::DAkOmegaSST::epsilonByk
tmp< volScalarField::Internal > epsilonByk(const volScalarField &F1, const volTensorField &gradU) const
Definition: DAkOmegaSST.C:251
Foam::DAkOmegaSST::omegaNearWall_
scalarList omegaNearWall_
Definition: DAkOmegaSST.H:177
Foam::DAkOmegaSST::correctModelStates
virtual void correctModelStates(wordList &modelStates) const
update the turbulence state for DAStateInfo::regStates_
Definition: DAkOmegaSST.C:289
Foam::DAkOmegaSST::beta2_
dimensionedScalar beta2_
Definition: DAkOmegaSST.H:65
Foam::DAkOmegaSST::TypeName
TypeName("kOmegaSST")
Foam::DAkOmegaSST::blend
tmp< volScalarField > blend(const volScalarField &F1, const dimensionedScalar &psi1, const dimensionedScalar &psi2) const
Definition: DAkOmegaSST.H:83
Foam::DAkOmegaSST::betaStar_
dimensionedScalar betaStar_
Definition: DAkOmegaSST.H:67
Foam::DAkOmegaSST::kSource
tmp< fvScalarMatrix > kSource() const
Definition: DAkOmegaSST.C:258
Foam::DAOption
Definition: DAOption.H:29
Foam::DAkOmegaSST::F1
tmp< volScalarField > F1(const volScalarField &CDkOmega) const
Definition: DAkOmegaSST.C:181
Foam::DAkOmegaSST::Qsas
tmp< fvScalarMatrix > Qsas(const volScalarField::Internal &S2, const volScalarField::Internal &gamma, const volScalarField::Internal &beta) const
Definition: DAkOmegaSST.C:276
Foam::DAkOmegaSST::k_
volScalarField & k_
Definition: DAkOmegaSST.H:166
daOption
DAOption daOption(mesh, pyOptions_)
Foam::DAkOmegaSST::Pk
tmp< volScalarField::Internal > Pk(const volScalarField::Internal &G) const
Definition: DAkOmegaSST.C:245
Foam::DAkOmegaSST::DAkOmegaSST
DAkOmegaSST(const word modelType, const fvMesh &mesh, const DAOption &daOption)
Definition: DAkOmegaSST.C:41
Foam::DAkOmegaSST::alphaK
tmp< volScalarField > alphaK(const volScalarField &F1) const
Definition: DAkOmegaSST.H:99
Foam::DAkOmegaSST::F3_
Switch F3_
Definition: DAkOmegaSST.H:73
Foam::DAkOmegaSST
Definition: DAkOmegaSST.H:48
Foam::DAkOmegaSST::omegaSource
tmp< fvScalarMatrix > omegaSource() const
Definition: DAkOmegaSST.C:267
Foam::DAkOmegaSST::solveTurbState_
label solveTurbState_
whether to solve for turb states
Definition: DAkOmegaSST.H:180
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
Foam::DAkOmegaSST::omega_
volScalarField & omega_
Definition: DAkOmegaSST.H:164
Foam::DAkOmegaSST::saveOmegaNearWall
void saveOmegaNearWall()
save near wall omega values to omegaNearWall_
Definition: DAkOmegaSST.C:386
Foam::DAkOmegaSST::correctOmegaBoundaryConditions
void correctOmegaBoundaryConditions()
specially treatment to correct epsilon BC
Definition: DAkOmegaSST.C:359
Foam::DAkOmegaSST::gamma2_
dimensionedScalar gamma2_
Definition: DAkOmegaSST.H:62
Foam::DAkOmegaSST::F23
tmp< volScalarField > F23() const
Definition: DAkOmegaSST.C:222
Foam::DAkOmegaSST::DomegaEff
tmp< volScalarField > DomegaEff(const volScalarField &F1) const
Definition: DAkOmegaSST.H:129
Foam::DAkOmegaSST::GbyNu
tmp< volScalarField::Internal > GbyNu(const volScalarField::Internal &GbyNu0, const volScalarField::Internal &F2, const volScalarField::Internal &S2) const
Definition: DAkOmegaSST.C:234
Foam::DAkOmegaSST::~DAkOmegaSST
virtual ~DAkOmegaSST()
Definition: DAkOmegaSST.H:196
Foam::DATurbulenceModel::nut_
volScalarField & nut_
turbulence viscosity
Definition: DATurbulenceModel.H:83
Foam::DAkOmegaSST::alphaK2_
dimensionedScalar alphaK2_
Definition: DAkOmegaSST.H:56
Foam::DAkOmegaSST::a1_
dimensionedScalar a1_
Definition: DAkOmegaSST.H:69
Foam::DAkOmegaSST::blend
tmp< volScalarField::Internal > blend(const volScalarField::Internal &F1, const dimensionedScalar &psi1, const dimensionedScalar &psi2) const
Definition: DAkOmegaSST.H:91
Foam::DAkOmegaSST::y_
const volScalarField & y_
3D wall distance
Definition: DAkOmegaSST.H:171
Foam
Definition: multiFreqScalarFvPatchField.C:144
Foam::DAkOmegaSST::alphaOmega1_
dimensionedScalar alphaOmega1_
Definition: DAkOmegaSST.H:58
Foam::DAkOmegaSST::gamma
tmp< volScalarField::Internal > gamma(const volScalarField::Internal &F1) const
Definition: DAkOmegaSST.H:115
Foam::DAkOmegaSST::gamma1_
dimensionedScalar gamma1_
Definition: DAkOmegaSST.H:61
Foam::DATurbulenceModel::nu
tmp< volScalarField > nu() const
get the nu field
Definition: DATurbulenceModel.C:262
Foam::DAkOmegaSST::F2
tmp< volScalarField > F2() const
Definition: DAkOmegaSST.C:200
Foam::DAkOmegaSST::b1_
dimensionedScalar b1_
Definition: DAkOmegaSST.H:70
Foam::DAkOmegaSST::alphaOmega
tmp< volScalarField > alphaOmega(const volScalarField &F1) const
Definition: DAkOmegaSST.H:104
Foam::DATurbulenceModel
Definition: DATurbulenceModel.H:61
Foam::DAkOmegaSST::updateIntermediateVariables
virtual void updateIntermediateVariables()
update any intermediate variables that are dependent on state variables and are used in calcResiduals
Definition: DAkOmegaSST.C:436
Foam::DAkOmegaSST::DkEff
tmp< volScalarField > DkEff(const volScalarField &F1) const
Definition: DAkOmegaSST.H:122
Foam::DAkOmegaSST::beta1_
dimensionedScalar beta1_
Definition: DAkOmegaSST.H:64
Foam::DAkOmegaSST::correctStateResidualModelCon
virtual void correctStateResidualModelCon(List< List< word >> &stateCon) const
update the original variable connectivity for the adjoint state residuals in stateCon
Definition: DAkOmegaSST.C:447
Foam::DAkOmegaSST::alphaK1_
dimensionedScalar alphaK1_
Definition: DAkOmegaSST.H:55
Foam::DAkOmegaSST::calcResiduals
virtual void calcResiduals(const dictionary &options)
compute the turbulence residuals
Definition: DAkOmegaSST.C:658
Foam::DAkOmegaSST::setOmegaNearWall
void setOmegaNearWall()
set omegaNearWall_ to near wall omega values
Definition: DAkOmegaSST.C:410
Foam::DAkOmegaSST::kRes_
volScalarField kRes_
Definition: DAkOmegaSST.H:167
Foam::DAkOmegaSST::alphaOmega2_
dimensionedScalar alphaOmega2_
Definition: DAkOmegaSST.H:59
Foam::DAkOmegaSST::correctBoundaryConditions
virtual void correctBoundaryConditions()
update turbulence variable boundary values
Definition: DAkOmegaSST.C:347
Foam::DAkOmegaSST::F3
tmp< volScalarField > F3() const
Definition: DAkOmegaSST.C:212
DATurbulenceModel.H
Foam::DAkOmegaSST::correct
virtual void correct()
solve the residual equations and update the state
Definition: DAkOmegaSST.C:638
Foam::DAkOmegaSST::addModelResidualCon
virtual void addModelResidualCon(HashTable< List< List< word >>> &allCon) const
add the model residual connectivity to stateCon
Definition: DAkOmegaSST.C:546
Foam::DAkOmegaSST::beta
tmp< volScalarField::Internal > beta(const volScalarField::Internal &F1) const
Definition: DAkOmegaSST.H:109
Foam::DAkOmegaSST::correctNut
virtual void correctNut()
update nut based on other turbulence variables and update the BCs
Definition: DAkOmegaSST.C:325
Foam::DAkOmegaSST::omegaRes_
volScalarField omegaRes_
Definition: DAkOmegaSST.H:165