DAkOmegaSSTFieldInversion.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 with a betaFieldInversion field
8  multiplying to the omega transport equation. This betaFieldInversion term
9  can then be trained to improve the kOmegaSST model.
10 
11  This file is modified from OpenFOAM's source code
12  src/TurbulenceModels/turbulenceModels/RAS/kOmega/kOmega.H
13 
14  OpenFOAM: The Open Source CFD Toolbox
15 
16  Copyright (C): 2011-2016 OpenFOAM Foundation
17 
18  OpenFOAM License:
19 
20  OpenFOAM is free software: you can redistribute it and/or modify it
21  under the terms of the GNU General Public License as published by
22  the Free Software Foundation, either version 3 of the License, or
23  (at your option) any later version.
24 
25  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
26  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
27  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
28  for more details.
29 
30  You should have received a copy of the GNU General Public License
31  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef DAkOmegaSSTFieldInversion_H
36 #define DAkOmegaSSTFieldInversion_H
37 
38 #include "DATurbulenceModel.H"
39 #include "addToRunTimeSelectionTable.H"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 
46 /*---------------------------------------------------------------------------*\
47  Class DAkOmegaSSTFieldInversion Declaration
48 \*---------------------------------------------------------------------------*/
49 
51  : public DATurbulenceModel
52 {
53 
54 protected:
56 
57  dimensionedScalar alphaK1_;
58  dimensionedScalar alphaK2_;
59 
60  dimensionedScalar alphaOmega1_;
61  dimensionedScalar alphaOmega2_;
62 
63  dimensionedScalar gamma1_;
64  dimensionedScalar gamma2_;
65 
66  dimensionedScalar beta1_;
67  dimensionedScalar beta2_;
68 
69  dimensionedScalar betaStar_;
70 
71  dimensionedScalar a1_;
72  dimensionedScalar b1_;
73  dimensionedScalar c1_;
74 
75  Switch F3_;
77 
79 
80  tmp<volScalarField> F1(const volScalarField& CDkOmega) const;
81  tmp<volScalarField> F2() const;
82  tmp<volScalarField> F3() const;
83  tmp<volScalarField> F23() const;
84 
85  tmp<volScalarField> blend(
86  const volScalarField& F1,
87  const dimensionedScalar& psi1,
88  const dimensionedScalar& psi2) const
89  {
90  return F1 * (psi1 - psi2) + psi2;
91  }
92 
93  tmp<volScalarField::Internal> blend(
94  const volScalarField::Internal& F1,
95  const dimensionedScalar& psi1,
96  const dimensionedScalar& psi2) const
97  {
98  return F1 * (psi1 - psi2) + psi2;
99  }
100 
101  tmp<volScalarField> alphaK(const volScalarField& F1) const
102  {
103  return blend(F1, alphaK1_, alphaK2_);
104  }
105 
106  tmp<volScalarField> alphaOmega(const volScalarField& F1) const
107  {
108  return blend(F1, alphaOmega1_, alphaOmega2_);
109  }
110 
111  tmp<volScalarField::Internal> beta(
112  const volScalarField::Internal& F1) const
113  {
114  return blend(F1, beta1_, beta2_);
115  }
116 
117  tmp<volScalarField::Internal> gamma(
118  const volScalarField::Internal& F1) const
119  {
120  return blend(F1, gamma1_, gamma2_);
121  }
122 
123  //- Return the effective diffusivity for k
124  tmp<volScalarField> DkEff(const volScalarField& F1) const
125  {
126  return tmp<volScalarField>(
127  new volScalarField("DkEff", alphaK(F1) * nut_ + this->nu()));
128  }
129 
130  //- Return the effective diffusivity for omega
131  tmp<volScalarField> DomegaEff(const volScalarField& F1) const
132  {
133  return tmp<volScalarField>(
134  new volScalarField(
135  "DomegaEff",
136  alphaOmega(F1) * nut_ + this->nu()));
137  }
138 
139  //- Return k production rate
140  tmp<volScalarField::Internal> Pk(
141  const volScalarField::Internal& G) const;
142 
143  //- Return epsilon/k which for standard RAS is betaStar*omega
144  tmp<volScalarField::Internal> epsilonByk(
145  const volScalarField& F1,
146  const volTensorField& gradU) const;
147 
148  //- Return G/nu
149  tmp<volScalarField::Internal> GbyNu(
150  const volScalarField::Internal& GbyNu0,
151  const volScalarField::Internal& F2,
152  const volScalarField::Internal& S2) const;
153 
154  tmp<fvScalarMatrix> kSource() const;
155 
156  tmp<fvScalarMatrix> omegaSource() const;
157 
158  tmp<fvScalarMatrix> Qsas(
159  const volScalarField::Internal& S2,
160  const volScalarField::Internal& gamma,
161  const volScalarField::Internal& beta) const;
163 
165 
166  volScalarField& omega_;
167  volScalarField omegaRes_;
168  volScalarField omegaResRef_;
169  volScalarField omegaResPartDeriv_;
170  volScalarField omegaRef_;
171  volScalarField& k_;
172  volScalarField kRes_;
173  volScalarField kResRef_;
174  volScalarField kResPartDeriv_;
175  volScalarField kRef_;
177 
179  volScalarField& betaFieldInversion_;
180 
182  volScalarField surfaceFriction_;
183 
185  volScalarField surfaceFrictionData_;
186 
188  volScalarField pData_;
189 
191  volVectorField UData_;
192 
194  volScalarField USingleComponentData_;
195 
200  //word eqnToModify_;
201 
203  //scalar betaFieldInversionWeight_;
204 
206  const volScalarField& y_;
207 
212  scalarList omegaNearWall_;
213 
215  label solveTurbState_ = 0;
216 
219 
220 public:
221  TypeName("kOmegaSSTFieldInversion");
222  // Constructors
223 
224  //- Construct from components
226  const word modelType,
227  const fvMesh& mesh,
228  const DAOption& daOption);
229 
230  //- Destructor
232  {
233  }
234 
235  // Member functions
236 
238  virtual void correctModelStates(wordList& modelStates) const;
239 
241  virtual void correctNut();
242 
244  virtual void correctBoundaryConditions();
245 
247  virtual void updateIntermediateVariables();
248 
250  virtual void correctStateResidualModelCon(List<List<word>>& stateCon) const;
251 
253  virtual void addModelResidualCon(HashTable<List<List<word>>>& allCon) const;
254 
256  virtual void calcResiduals(const dictionary& options);
257 
259  virtual void correct();
260 
262  void saveOmegaNearWall();
263 
265  void setOmegaNearWall();
266 
269 };
270 
271 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
272 
273 } // End namespace Foam
274 
275 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
276 
277 #endif
278 
279 // ************************************************************************* //
Foam::DAkOmegaSSTFieldInversion::UData_
volVectorField UData_
reference field (e.g. velocity for DNS)
Definition: DAkOmegaSSTFieldInversion.H:191
Foam::DAkOmegaSSTFieldInversion::USingleComponentData_
volScalarField USingleComponentData_
the reference pressure field data
Definition: DAkOmegaSSTFieldInversion.H:194
Foam::DAkOmegaSSTFieldInversion::DkEff
tmp< volScalarField > DkEff(const volScalarField &F1) const
Definition: DAkOmegaSSTFieldInversion.H:124
Foam::DAkOmegaSSTFieldInversion::printInterval_
label printInterval_
time step interval to print residual
Definition: DAkOmegaSSTFieldInversion.H:218
Foam::DAkOmegaSSTFieldInversion::correctStateResidualModelCon
virtual void correctStateResidualModelCon(List< List< word >> &stateCon) const
update the original variable connectivity for the adjoint state residuals in stateCon
Definition: DAkOmegaSSTFieldInversion.C:508
Foam::DAkOmegaSSTFieldInversion::solveTurbState_
label solveTurbState_
whether to solve for turb states
Definition: DAkOmegaSSTFieldInversion.H:215
Foam::DAkOmegaSSTFieldInversion::blend
tmp< volScalarField::Internal > blend(const volScalarField::Internal &F1, const dimensionedScalar &psi1, const dimensionedScalar &psi2) const
Definition: DAkOmegaSSTFieldInversion.H:93
Foam::DAkOmegaSSTFieldInversion::~DAkOmegaSSTFieldInversion
virtual ~DAkOmegaSSTFieldInversion()
Definition: DAkOmegaSSTFieldInversion.H:231
Foam::DAkOmegaSSTFieldInversion::F3
tmp< volScalarField > F3() const
Definition: DAkOmegaSSTFieldInversion.C:273
Foam::DAkOmegaSSTFieldInversion::alphaK2_
dimensionedScalar alphaK2_
Definition: DAkOmegaSSTFieldInversion.H:58
Foam::DAkOmegaSSTFieldInversion::TypeName
TypeName("kOmegaSSTFieldInversion")
Foam::DAkOmegaSSTFieldInversion::alphaOmega
tmp< volScalarField > alphaOmega(const volScalarField &F1) const
Definition: DAkOmegaSSTFieldInversion.H:106
Foam::DAkOmegaSSTFieldInversion::gamma
tmp< volScalarField::Internal > gamma(const volScalarField::Internal &F1) const
Definition: DAkOmegaSSTFieldInversion.H:117
Foam::DAkOmegaSSTFieldInversion::correctModelStates
virtual void correctModelStates(wordList &modelStates) const
update the turbulence state for DAStateInfo::regStates_
Definition: DAkOmegaSSTFieldInversion.C:350
Foam::DAkOmegaSSTFieldInversion::Qsas
tmp< fvScalarMatrix > Qsas(const volScalarField::Internal &S2, const volScalarField::Internal &gamma, const volScalarField::Internal &beta) const
Definition: DAkOmegaSSTFieldInversion.C:337
Foam::DAkOmegaSSTFieldInversion::omegaResPartDeriv_
volScalarField omegaResPartDeriv_
Definition: DAkOmegaSSTFieldInversion.H:169
Foam::DAkOmegaSSTFieldInversion::correctOmegaBoundaryConditions
void correctOmegaBoundaryConditions()
specially treatment to correct epsilon BC
Definition: DAkOmegaSSTFieldInversion.C:420
Foam::DAkOmegaSSTFieldInversion::setOmegaNearWall
void setOmegaNearWall()
set omegaNearWall_ to near wall omega values
Definition: DAkOmegaSSTFieldInversion.C:471
Foam::DAkOmegaSSTFieldInversion::k_
volScalarField & k_
Definition: DAkOmegaSSTFieldInversion.H:171
Foam::DAkOmegaSSTFieldInversion::epsilonByk
tmp< volScalarField::Internal > epsilonByk(const volScalarField &F1, const volTensorField &gradU) const
Definition: DAkOmegaSSTFieldInversion.C:312
Foam::DAkOmegaSSTFieldInversion::kRes_
volScalarField kRes_
Definition: DAkOmegaSSTFieldInversion.H:172
Foam::DAkOmegaSSTFieldInversion::betaFieldInversion_
volScalarField & betaFieldInversion_
A beta field multiplying to the production term.
Definition: DAkOmegaSSTFieldInversion.H:179
Foam::DAOption
Definition: DAOption.H:29
Foam::DAkOmegaSSTFieldInversion::Pk
tmp< volScalarField::Internal > Pk(const volScalarField::Internal &G) const
Definition: DAkOmegaSSTFieldInversion.C:306
Foam::DAkOmegaSSTFieldInversion::kRef_
volScalarField kRef_
Definition: DAkOmegaSSTFieldInversion.H:175
daOption
DAOption daOption(mesh, pyOptions_)
Foam::DAkOmegaSSTFieldInversion::gamma2_
dimensionedScalar gamma2_
Definition: DAkOmegaSSTFieldInversion.H:64
Foam::DAkOmegaSSTFieldInversion::updateIntermediateVariables
virtual void updateIntermediateVariables()
update any intermediate variables that are dependent on state variables and are used in calcResiduals
Definition: DAkOmegaSSTFieldInversion.C:497
Foam::DAkOmegaSSTFieldInversion::correctBoundaryConditions
virtual void correctBoundaryConditions()
update turbulence variable boundary values
Definition: DAkOmegaSSTFieldInversion.C:408
Foam::DAkOmegaSSTFieldInversion::GbyNu
tmp< volScalarField::Internal > GbyNu(const volScalarField::Internal &GbyNu0, const volScalarField::Internal &F2, const volScalarField::Internal &S2) const
Definition: DAkOmegaSSTFieldInversion.C:295
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
Foam::DAkOmegaSSTFieldInversion::alphaOmega2_
dimensionedScalar alphaOmega2_
Definition: DAkOmegaSSTFieldInversion.H:61
Foam::DAkOmegaSSTFieldInversion::omega_
volScalarField & omega_
Definition: DAkOmegaSSTFieldInversion.H:166
Foam::DAkOmegaSSTFieldInversion::c1_
dimensionedScalar c1_
Definition: DAkOmegaSSTFieldInversion.H:73
Foam::DAkOmegaSSTFieldInversion::omegaRes_
volScalarField omegaRes_
Definition: DAkOmegaSSTFieldInversion.H:167
Foam::DAkOmegaSSTFieldInversion::beta1_
dimensionedScalar beta1_
Definition: DAkOmegaSSTFieldInversion.H:66
Foam::DAkOmegaSSTFieldInversion::beta2_
dimensionedScalar beta2_
Definition: DAkOmegaSSTFieldInversion.H:67
Foam::DAkOmegaSSTFieldInversion::addModelResidualCon
virtual void addModelResidualCon(HashTable< List< List< word >>> &allCon) const
add the model residual connectivity to stateCon
Definition: DAkOmegaSSTFieldInversion.C:607
Foam::DAkOmegaSSTFieldInversion::a1_
dimensionedScalar a1_
Definition: DAkOmegaSSTFieldInversion.H:71
Foam::DAkOmegaSSTFieldInversion::pData_
volScalarField pData_
the reference field for surfacePressure
Definition: DAkOmegaSSTFieldInversion.H:188
Foam::DAkOmegaSSTFieldInversion::alphaOmega1_
dimensionedScalar alphaOmega1_
Definition: DAkOmegaSSTFieldInversion.H:60
Foam::DAkOmegaSSTFieldInversion::kResRef_
volScalarField kResRef_
Definition: DAkOmegaSSTFieldInversion.H:173
Foam::DAkOmegaSSTFieldInversion::kSource
tmp< fvScalarMatrix > kSource() const
Definition: DAkOmegaSSTFieldInversion.C:319
Foam::DATurbulenceModel::nut_
volScalarField & nut_
turbulence viscosity
Definition: DATurbulenceModel.H:83
Foam::DAkOmegaSSTFieldInversion::y_
const volScalarField & y_
Weight of betaFieldInversion for the turbulent transport equations if both are modified.
Definition: DAkOmegaSSTFieldInversion.H:206
Foam::DAkOmegaSSTFieldInversion::blend
tmp< volScalarField > blend(const volScalarField &F1, const dimensionedScalar &psi1, const dimensionedScalar &psi2) const
Definition: DAkOmegaSSTFieldInversion.H:85
Foam::DAkOmegaSSTFieldInversion::calcResiduals
virtual void calcResiduals(const dictionary &options)
compute the turbulence residuals
Definition: DAkOmegaSSTFieldInversion.C:719
Foam::DAkOmegaSSTFieldInversion::surfaceFriction_
volScalarField surfaceFriction_
a surface friction 'field' when using skin friction data for field inversion
Definition: DAkOmegaSSTFieldInversion.H:182
Foam::DAkOmegaSSTFieldInversion::omegaResRef_
volScalarField omegaResRef_
Definition: DAkOmegaSSTFieldInversion.H:168
Foam
Definition: multiFreqScalarFvPatchField.C:144
Foam::DAkOmegaSSTFieldInversion::alphaK
tmp< volScalarField > alphaK(const volScalarField &F1) const
Definition: DAkOmegaSSTFieldInversion.H:101
Foam::DAkOmegaSSTFieldInversion::betaStar_
dimensionedScalar betaStar_
Definition: DAkOmegaSSTFieldInversion.H:69
Foam::DAkOmegaSSTFieldInversion::beta
tmp< volScalarField::Internal > beta(const volScalarField::Internal &F1) const
Definition: DAkOmegaSSTFieldInversion.H:111
Foam::DAkOmegaSSTFieldInversion::F3_
Switch F3_
Definition: DAkOmegaSSTFieldInversion.H:75
Foam::DATurbulenceModel::nu
tmp< volScalarField > nu() const
get the nu field
Definition: DATurbulenceModel.C:262
Foam::DAkOmegaSSTFieldInversion::F2
tmp< volScalarField > F2() const
Definition: DAkOmegaSSTFieldInversion.C:261
Foam::DAkOmegaSSTFieldInversion::alphaK1_
dimensionedScalar alphaK1_
Definition: DAkOmegaSSTFieldInversion.H:57
Foam::DATurbulenceModel
Definition: DATurbulenceModel.H:61
Foam::DAkOmegaSSTFieldInversion::DomegaEff
tmp< volScalarField > DomegaEff(const volScalarField &F1) const
Definition: DAkOmegaSSTFieldInversion.H:131
Foam::DAkOmegaSSTFieldInversion::omegaSource
tmp< fvScalarMatrix > omegaSource() const
Definition: DAkOmegaSSTFieldInversion.C:328
Foam::DAkOmegaSSTFieldInversion::gamma1_
dimensionedScalar gamma1_
Definition: DAkOmegaSSTFieldInversion.H:63
Foam::DAkOmegaSSTFieldInversion::b1_
dimensionedScalar b1_
Definition: DAkOmegaSSTFieldInversion.H:72
Foam::DAkOmegaSSTFieldInversion::kResPartDeriv_
volScalarField kResPartDeriv_
Definition: DAkOmegaSSTFieldInversion.H:174
Foam::DAkOmegaSSTFieldInversion::DAkOmegaSSTFieldInversion
DAkOmegaSSTFieldInversion(const word modelType, const fvMesh &mesh, const DAOption &daOption)
Definition: DAkOmegaSSTFieldInversion.C:41
Foam::DAkOmegaSSTFieldInversion::F23
tmp< volScalarField > F23() const
Definition: DAkOmegaSSTFieldInversion.C:283
Foam::DAkOmegaSSTFieldInversion::F1
tmp< volScalarField > F1(const volScalarField &CDkOmega) const
Definition: DAkOmegaSSTFieldInversion.C:242
Foam::DAkOmegaSSTFieldInversion::omegaNearWall_
scalarList omegaNearWall_
Definition: DAkOmegaSSTFieldInversion.H:212
Foam::DAkOmegaSSTFieldInversion
Definition: DAkOmegaSSTFieldInversion.H:50
Foam::DAkOmegaSSTFieldInversion::omegaRef_
volScalarField omegaRef_
Definition: DAkOmegaSSTFieldInversion.H:170
Foam::DAkOmegaSSTFieldInversion::surfaceFrictionData_
volScalarField surfaceFrictionData_
the reference field for surfaceFriction
Definition: DAkOmegaSSTFieldInversion.H:185
DATurbulenceModel.H
Foam::DAkOmegaSSTFieldInversion::correct
virtual void correct()
solve the residual equations and update the state
Definition: DAkOmegaSSTFieldInversion.C:699
Foam::DAkOmegaSSTFieldInversion::saveOmegaNearWall
void saveOmegaNearWall()
save near wall omega values to omegaNearWall_
Definition: DAkOmegaSSTFieldInversion.C:447
Foam::DAkOmegaSSTFieldInversion::correctNut
virtual void correctNut()
update nut based on other turbulence variables and update the BCs
Definition: DAkOmegaSSTFieldInversion.C:386