DAkOmegaSSTLM.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 kOmegaSSTLM model
8 
9  This file is modified from OpenFOAM's source code
10  src/TurbulenceModels/turbulenceModels/RAS/kOmegaSSTLM/kOmegaSSTLM.C
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 DAkOmegaSSTLM_H
34 #define DAkOmegaSSTLM_H
35 
36 #include "DATurbulenceModel.H"
37 #include "addToRunTimeSelectionTable.H"
38 
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43 
44 /*---------------------------------------------------------------------------*\
45  Class DAkOmegaSSTLM 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  dimensionedScalar ca1_;
79  dimensionedScalar ca2_;
80  dimensionedScalar ce1_;
81  dimensionedScalar ce2_;
82  dimensionedScalar cThetat_;
83  dimensionedScalar sigmaThetat_;
84  scalar lambdaErr_; // Convergence criterion for the lambda/thetat loop
85  label maxLambdaIter_; // Maximum number of iterations to converge the lambda/thetat loop
86  const dimensionedScalar deltaU_; // Stabilization for division by the magnitude of the velocity
88 
90 
91  tmp<volScalarField> F1SST(const volScalarField& CDkOmega) const;
92  tmp<volScalarField> F2() const;
93  tmp<volScalarField> F3() const;
94  tmp<volScalarField> F23() const;
95 
96  tmp<volScalarField> blend(
97  const volScalarField& F1,
98  const dimensionedScalar& psi1,
99  const dimensionedScalar& psi2) const
100  {
101  return F1 * (psi1 - psi2) + psi2;
102  }
103 
104  tmp<volScalarField::Internal> blend(
105  const volScalarField::Internal& F1,
106  const dimensionedScalar& psi1,
107  const dimensionedScalar& psi2) const
108  {
109  return F1 * (psi1 - psi2) + psi2;
110  }
111 
112  tmp<volScalarField> alphaK(const volScalarField& F1) const
113  {
114  return blend(F1, alphaK1_, alphaK2_);
115  }
116 
117  tmp<volScalarField> alphaOmega(const volScalarField& F1) const
118  {
119  return blend(F1, alphaOmega1_, alphaOmega2_);
120  }
121 
122  tmp<volScalarField::Internal> beta(
123  const volScalarField::Internal& F1) const
124  {
125  return blend(F1, beta1_, beta2_);
126  }
127 
128  tmp<volScalarField::Internal> gamma(
129  const volScalarField::Internal& F1) const
130  {
131  return blend(F1, gamma1_, gamma2_);
132  }
133 
134  //- Return k production rate
135  tmp<volScalarField::Internal> PkSST(
136  const volScalarField::Internal& G) const;
137 
138  //- Return epsilon/k which for standard RAS is betaStar*omega
139  tmp<volScalarField::Internal> epsilonBykSST(
140  const volScalarField& F1,
141  const volTensorField& gradU) const;
142 
143  //- Return G/nu
144  tmp<volScalarField::Internal> GbyNu(
145  const volScalarField::Internal& GbyNu0,
146  const volScalarField::Internal& F2,
147  const volScalarField::Internal& S2) const;
148 
149  tmp<fvScalarMatrix> kSource() const;
150 
151  tmp<fvScalarMatrix> omegaSource() const;
152 
153  tmp<fvScalarMatrix> Qsas(
154  const volScalarField::Internal& S2,
155  const volScalarField::Internal& gamma,
156  const volScalarField::Internal& beta) const;
158 
159  // Modified form of the k-omega SST F1 function
160  tmp<volScalarField> F1(const volScalarField& CDkOmega) const;
161 
162  // Modified form of the k-omega SST k production rate
163  tmp<volScalarField::Internal> Pk(
164  const volScalarField::Internal& G) const;
165 
166  // Modified form of the k-omega SST epsilon/k
167  tmp<volScalarField::Internal> epsilonByk(
168  const volScalarField& F1,
169  const volTensorField& gradU) const;
170 
171  // Freestream blending-function
172  tmp<volScalarField::Internal> Fthetat(
173  const volScalarField::Internal& Us,
174  const volScalarField::Internal& Omega,
175  const volScalarField::Internal& nu) const;
176 
177  // Empirical correlation for critical Reynolds number where the
178  // intermittency first starts to increase in the boundary layer
179  tmp<volScalarField::Internal> ReThetac() const;
180 
181  // Empirical correlation that controls the length of the
182  // transition region
183  tmp<volScalarField::Internal> Flength(
184  const volScalarField::Internal& nu) const;
185 
186  // Transition onset location control function
187  tmp<volScalarField::Internal> Fonset(
188  const volScalarField::Internal& Rev,
189  const volScalarField::Internal& ReThetac,
190  const volScalarField::Internal& RT) const;
191 
192  //- Return the transition onset momentum-thickness Reynolds number
193  // (based on freestream conditions)
194  tmp<volScalarField::Internal> ReThetat0(
195  const volScalarField::Internal& Us,
196  const volScalarField::Internal& dUsds,
197  const volScalarField::Internal& nu) const;
198 
200 
201  volScalarField& omega_;
202  volScalarField omegaRes_;
203  volScalarField& k_;
204  volScalarField kRes_;
205  volScalarField& ReThetat_; // Transition onset momentum-thickness Reynolds number
206  volScalarField ReThetatRes_;
207  volScalarField& gammaInt_; // Intermittency
208  volScalarField gammaIntRes_;
210 
212  volScalarField::Internal& gammaIntEff_;
213 
215  const volScalarField& y_;
216 
221  scalarList omegaNearWall_;
222 
224  label solveTurbState_ = 0;
225 
228 
229 public:
230  TypeName("kOmegaSSTLM");
231  // Constructors
232 
233  //- Construct from components
235  const word modelType,
236  const fvMesh& mesh,
237  const DAOption& daOption);
238 
239  //- Destructor
240  virtual ~DAkOmegaSSTLM()
241  {
242  }
243 
244  // Member functions
245 
246  //- Return the effective diffusivity for k
247  tmp<volScalarField> DkEff(const volScalarField& F1) const
248  {
249  return tmp<volScalarField>(
250  new volScalarField("DkEff", alphaK(F1) * nut_ + this->nu()));
251  }
252 
253  //- Return the effective diffusivity for omega
254  tmp<volScalarField> DomegaEff(const volScalarField& F1) const
255  {
256  return tmp<volScalarField>(
257  new volScalarField(
258  "DomegaEff",
259  alphaOmega(F1) * nut_ + this->nu()));
260  }
261 
262  //- Return the effective diffusivity for transition onset
263  // momentum-thickness Reynolds number
264  tmp<volScalarField> DReThetatEff() const
265  {
266  return tmp<volScalarField>(
267  new volScalarField(
268  "DReThetatEff",
269  sigmaThetat_ * (nut_ + this->nu())));
270  }
271  //- Return the effective diffusivity for intermittency
272  tmp<volScalarField> DgammaIntEff() const
273  {
274  return tmp<volScalarField>(
275  new volScalarField(
276  "DgammaIntEff",
277  nut_ + this->nu()));
278  }
279 
281  virtual void correctModelStates(wordList& modelStates) const;
282 
284  virtual void correctNut();
285 
287  virtual void correctBoundaryConditions();
288 
290  virtual void updateIntermediateVariables();
291 
293  virtual void correctStateResidualModelCon(List<List<word>>& stateCon) const;
294 
296  virtual void addModelResidualCon(HashTable<List<List<word>>>& allCon) const;
297 
299  virtual void calcResiduals(const dictionary& options);
300 
302  virtual void correct();
303 
305  void saveOmegaNearWall();
306 
308  void setOmegaNearWall();
309 
312 };
313 
314 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
315 
316 } // End namespace Foam
317 
318 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
319 
320 #endif
321 
322 // ************************************************************************* //
Foam::DAkOmegaSSTLM::omega_
volScalarField & omega_
Definition: DAkOmegaSSTLM.H:201
Foam::DAkOmegaSSTLM::ca2_
dimensionedScalar ca2_
Definition: DAkOmegaSSTLM.H:79
Foam::DAkOmegaSSTLM::alphaOmega1_
dimensionedScalar alphaOmega1_
Definition: DAkOmegaSSTLM.H:58
Foam::DAkOmegaSSTLM::omegaNearWall_
scalarList omegaNearWall_
Definition: DAkOmegaSSTLM.H:221
Foam::DAkOmegaSSTLM::ce2_
dimensionedScalar ce2_
Definition: DAkOmegaSSTLM.H:81
Foam::DAkOmegaSSTLM::epsilonByk
tmp< volScalarField::Internal > epsilonByk(const volScalarField &F1, const volTensorField &gradU) const
Definition: DAkOmegaSSTLM.C:366
Foam::DAkOmegaSSTLM::gammaIntEff_
volScalarField::Internal & gammaIntEff_
Effective intermittency.
Definition: DAkOmegaSSTLM.H:212
Foam::DAkOmegaSSTLM::setOmegaNearWall
void setOmegaNearWall()
set omegaNearWall_ to near wall omega values
Definition: DAkOmegaSSTLM.C:723
Foam::DAkOmegaSSTLM::y_
const volScalarField & y_
3D wall distance
Definition: DAkOmegaSSTLM.H:215
Foam::DAkOmegaSSTLM::addModelResidualCon
virtual void addModelResidualCon(HashTable< List< List< word >>> &allCon) const
add the model residual connectivity to stateCon
Definition: DAkOmegaSSTLM.C:864
Foam::DAkOmegaSSTLM::ca1_
dimensionedScalar ca1_
Definition: DAkOmegaSSTLM.H:78
Foam::DAkOmegaSSTLM::beta
tmp< volScalarField::Internal > beta(const volScalarField::Internal &F1) const
Definition: DAkOmegaSSTLM.H:122
Foam::DAkOmegaSSTLM::ce1_
dimensionedScalar ce1_
Definition: DAkOmegaSSTLM.H:80
Foam::DAOption
Definition: DAOption.H:29
Foam::DAkOmegaSSTLM::F2
tmp< volScalarField > F2() const
Definition: DAkOmegaSSTLM.C:262
daOption
DAOption daOption(mesh, pyOptions_)
Foam::DAkOmegaSSTLM::DomegaEff
tmp< volScalarField > DomegaEff(const volScalarField &F1) const
Definition: DAkOmegaSSTLM.H:254
Foam::DAkOmegaSSTLM::beta2_
dimensionedScalar beta2_
Definition: DAkOmegaSSTLM.H:65
Foam::DAkOmegaSSTLM::F23
tmp< volScalarField > F23() const
Definition: DAkOmegaSSTLM.C:284
Foam::DAkOmegaSSTLM::calcResiduals
virtual void calcResiduals(const dictionary &options)
compute the turbulence residuals
Definition: DAkOmegaSSTLM.C:1008
Foam::DAkOmegaSSTLM::a1_
dimensionedScalar a1_
Definition: DAkOmegaSSTLM.H:69
Foam::DAkOmegaSSTLM::Pk
tmp< volScalarField::Internal > Pk(const volScalarField::Internal &G) const
Definition: DAkOmegaSSTLM.C:360
Foam::DAkOmegaSSTLM::omegaRes_
volScalarField omegaRes_
Definition: DAkOmegaSSTLM.H:202
Foam::DAkOmegaSSTLM::F3_
Switch F3_
Definition: DAkOmegaSSTLM.H:73
Foam::DAkOmegaSSTLM::DReThetatEff
tmp< volScalarField > DReThetatEff() const
Definition: DAkOmegaSSTLM.H:264
Foam::DAkOmegaSSTLM::sigmaThetat_
dimensionedScalar sigmaThetat_
Definition: DAkOmegaSSTLM.H:83
Foam::DAkOmegaSSTLM::~DAkOmegaSSTLM
virtual ~DAkOmegaSSTLM()
Definition: DAkOmegaSSTLM.H:240
Foam::DAkOmegaSSTLM::alphaOmega
tmp< volScalarField > alphaOmega(const volScalarField &F1) const
Definition: DAkOmegaSSTLM.H:117
Foam::DAkOmegaSSTLM::updateIntermediateVariables
virtual void updateIntermediateVariables()
update any intermediate variables that are dependent on state variables and are used in calcResiduals
Definition: DAkOmegaSSTLM.C:749
Foam::DAkOmegaSSTLM
Definition: DAkOmegaSSTLM.H:48
Foam::DAkOmegaSSTLM::correctNut
virtual void correctNut()
update nut based on other turbulence variables and update the BCs
Definition: DAkOmegaSSTLM.C:635
Foam::DAkOmegaSSTLM::correctBoundaryConditions
virtual void correctBoundaryConditions()
update turbulence variable boundary values
Definition: DAkOmegaSSTLM.C:657
Foam::DAkOmegaSSTLM::ReThetat_
volScalarField & ReThetat_
Definition: DAkOmegaSSTLM.H:205
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
Foam::DAkOmegaSSTLM::ReThetac
tmp< volScalarField::Internal > ReThetac() const
Definition: DAkOmegaSSTLM.C:396
Foam::DAkOmegaSSTLM::c1_
dimensionedScalar c1_
Definition: DAkOmegaSSTLM.H:71
Foam::DATurbulenceModel::nut_
volScalarField & nut_
turbulence viscosity
Definition: DATurbulenceModel.H:83
Foam::DAkOmegaSSTLM::gamma
tmp< volScalarField::Internal > gamma(const volScalarField::Internal &F1) const
Definition: DAkOmegaSSTLM.H:128
Foam::DAkOmegaSSTLM::beta1_
dimensionedScalar beta1_
Definition: DAkOmegaSSTLM.H:64
Foam::DAkOmegaSSTLM::alphaK2_
dimensionedScalar alphaK2_
Definition: DAkOmegaSSTLM.H:56
Foam::DAkOmegaSSTLM::maxLambdaIter_
label maxLambdaIter_
Definition: DAkOmegaSSTLM.H:85
Foam::DAkOmegaSSTLM::solveTurbState_
label solveTurbState_
whether to solve for turb states
Definition: DAkOmegaSSTLM.H:224
Foam::DAkOmegaSSTLM::gammaInt_
volScalarField & gammaInt_
Definition: DAkOmegaSSTLM.H:207
Foam::DAkOmegaSSTLM::printInterval_
label printInterval_
time step interval to print residual
Definition: DAkOmegaSSTLM.H:227
Foam
Definition: multiFreqScalarFvPatchField.C:144
Foam::DAkOmegaSSTLM::correct
virtual void correct()
solve the residual equations and update the state
Definition: DAkOmegaSSTLM.C:988
Foam::DAkOmegaSSTLM::Fthetat
tmp< volScalarField::Internal > Fthetat(const volScalarField::Internal &Us, const volScalarField::Internal &Omega, const volScalarField::Internal &nu) const
Definition: DAkOmegaSSTLM.C:374
Foam::DAkOmegaSSTLM::blend
tmp< volScalarField > blend(const volScalarField &F1, const dimensionedScalar &psi1, const dimensionedScalar &psi2) const
Definition: DAkOmegaSSTLM.H:96
Foam::DAkOmegaSSTLM::Qsas
tmp< fvScalarMatrix > Qsas(const volScalarField::Internal &S2, const volScalarField::Internal &gamma, const volScalarField::Internal &beta) const
Definition: DAkOmegaSSTLM.C:338
Foam::DAkOmegaSSTLM::Fonset
tmp< volScalarField::Internal > Fonset(const volScalarField::Internal &Rev, const volScalarField::Internal &ReThetac, const volScalarField::Internal &RT) const
Definition: DAkOmegaSSTLM.C:481
Foam::DAkOmegaSSTLM::alphaK
tmp< volScalarField > alphaK(const volScalarField &F1) const
Definition: DAkOmegaSSTLM.H:112
Foam::DATurbulenceModel::nu
tmp< volScalarField > nu() const
get the nu field
Definition: DATurbulenceModel.C:262
Foam::DAkOmegaSSTLM::deltaU_
const dimensionedScalar deltaU_
Definition: DAkOmegaSSTLM.H:86
Foam::DAkOmegaSSTLM::alphaOmega2_
dimensionedScalar alphaOmega2_
Definition: DAkOmegaSSTLM.H:59
Foam::DAkOmegaSSTLM::gamma2_
dimensionedScalar gamma2_
Definition: DAkOmegaSSTLM.H:62
Foam::DAkOmegaSSTLM::saveOmegaNearWall
void saveOmegaNearWall()
save near wall omega values to omegaNearWall_
Definition: DAkOmegaSSTLM.C:699
Foam::DAkOmegaSSTLM::omegaSource
tmp< fvScalarMatrix > omegaSource() const
Definition: DAkOmegaSSTLM.C:329
Foam::DATurbulenceModel
Definition: DATurbulenceModel.H:61
Foam::DAkOmegaSSTLM::F1
tmp< volScalarField > F1(const volScalarField &CDkOmega) const
Definition: DAkOmegaSSTLM.C:351
Foam::DAkOmegaSSTLM::kRes_
volScalarField kRes_
Definition: DAkOmegaSSTLM.H:204
Foam::DAkOmegaSSTLM::DkEff
tmp< volScalarField > DkEff(const volScalarField &F1) const
Definition: DAkOmegaSSTLM.H:247
Foam::DAkOmegaSSTLM::ReThetatRes_
volScalarField ReThetatRes_
Definition: DAkOmegaSSTLM.H:206
Foam::DAkOmegaSSTLM::DAkOmegaSSTLM
DAkOmegaSSTLM(const word modelType, const fvMesh &mesh, const DAOption &daOption)
Definition: DAkOmegaSSTLM.C:41
Foam::DAkOmegaSSTLM::DgammaIntEff
tmp< volScalarField > DgammaIntEff() const
Definition: DAkOmegaSSTLM.H:272
Foam::DAkOmegaSSTLM::ReThetat0
tmp< volScalarField::Internal > ReThetat0(const volScalarField::Internal &Us, const volScalarField::Internal &dUsds, const volScalarField::Internal &nu) const
Definition: DAkOmegaSSTLM.C:499
Foam::DAkOmegaSSTLM::gammaIntRes_
volScalarField gammaIntRes_
Definition: DAkOmegaSSTLM.H:208
Foam::DAkOmegaSSTLM::epsilonBykSST
tmp< volScalarField::Internal > epsilonBykSST(const volScalarField &F1, const volTensorField &gradU) const
Definition: DAkOmegaSSTLM.C:313
Foam::DAkOmegaSSTLM::Flength
tmp< volScalarField::Internal > Flength(const volScalarField::Internal &nu) const
Definition: DAkOmegaSSTLM.C:427
Foam::DAkOmegaSSTLM::lambdaErr_
scalar lambdaErr_
Definition: DAkOmegaSSTLM.H:84
Foam::DAkOmegaSSTLM::b1_
dimensionedScalar b1_
Definition: DAkOmegaSSTLM.H:70
Foam::DAkOmegaSSTLM::betaStar_
dimensionedScalar betaStar_
Definition: DAkOmegaSSTLM.H:67
Foam::DAkOmegaSSTLM::kSource
tmp< fvScalarMatrix > kSource() const
Definition: DAkOmegaSSTLM.C:320
Foam::DAkOmegaSSTLM::PkSST
tmp< volScalarField::Internal > PkSST(const volScalarField::Internal &G) const
Definition: DAkOmegaSSTLM.C:307
Foam::DAkOmegaSSTLM::correctStateResidualModelCon
virtual void correctStateResidualModelCon(List< List< word >> &stateCon) const
update the original variable connectivity for the adjoint state residuals in stateCon
Definition: DAkOmegaSSTLM.C:763
Foam::DAkOmegaSSTLM::blend
tmp< volScalarField::Internal > blend(const volScalarField::Internal &F1, const dimensionedScalar &psi1, const dimensionedScalar &psi2) const
Definition: DAkOmegaSSTLM.H:104
Foam::DAkOmegaSSTLM::correctOmegaBoundaryConditions
void correctOmegaBoundaryConditions()
specially treatment to correct epsilon BC
Definition: DAkOmegaSSTLM.C:672
Foam::DAkOmegaSSTLM::k_
volScalarField & k_
Definition: DAkOmegaSSTLM.H:203
Foam::DAkOmegaSSTLM::TypeName
TypeName("kOmegaSSTLM")
Foam::DAkOmegaSSTLM::alphaK1_
dimensionedScalar alphaK1_
Definition: DAkOmegaSSTLM.H:55
Foam::DAkOmegaSSTLM::F3
tmp< volScalarField > F3() const
Definition: DAkOmegaSSTLM.C:274
DATurbulenceModel.H
Foam::DAkOmegaSSTLM::cThetat_
dimensionedScalar cThetat_
Definition: DAkOmegaSSTLM.H:82
Foam::DAkOmegaSSTLM::GbyNu
tmp< volScalarField::Internal > GbyNu(const volScalarField::Internal &GbyNu0, const volScalarField::Internal &F2, const volScalarField::Internal &S2) const
Definition: DAkOmegaSSTLM.C:296
Foam::DAkOmegaSSTLM::gamma1_
dimensionedScalar gamma1_
Definition: DAkOmegaSSTLM.H:61
Foam::DAkOmegaSSTLM::correctModelStates
virtual void correctModelStates(wordList &modelStates) const
update the turbulence state for DAStateInfo::regStates_
Definition: DAkOmegaSSTLM.C:597
Foam::DAkOmegaSSTLM::F1SST
tmp< volScalarField > F1SST(const volScalarField &CDkOmega) const
Definition: DAkOmegaSSTLM.C:243