DAkOmegaFieldInversionOmega.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 kOmega model
8 
9  This file is modified from OpenFOAM's source code
10  src/TurbulenceModels/turbulenceModels/RAS/kOmega/kOmega.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 DAkOmegaFieldInversionOmega_H
34 #define DAkOmegaFieldInversionOmega_H
35 
36 #include "DATurbulenceModel.H"
37 #include "addToRunTimeSelectionTable.H"
38 
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43 
44 /*---------------------------------------------------------------------------*\
45  Class DAkOmegaFieldInversionOmega Declaration
46 \*---------------------------------------------------------------------------*/
47 
49  : public DATurbulenceModel
50 {
51 
52 protected:
54 
55  dimensionedScalar Cmu_;
56  dimensionedScalar beta_;
57  dimensionedScalar gamma_;
58  dimensionedScalar alphaK_;
59  dimensionedScalar alphaOmega_;
61 
63 
64  //- Return the effective diffusivity for k
65  tmp<volScalarField> DkEff() const
66  {
67  return tmp<volScalarField>
68  (
69  new volScalarField
70  (
71  "DkEff",
72  alphaK_*this->nut_ + this->nu()
73  )
74  );
75  }
76 
77  //- Return the effective diffusivity for omega
78  tmp<volScalarField> DomegaEff() const
79  {
80  return tmp<volScalarField>
81  (
82  new volScalarField
83  (
84  "DomegaEff",
85  alphaOmega_*this->nut_ + this->nu()
86  )
87  );
88  }
90 
92 
93  volScalarField& omega_;
94  volScalarField omegaRes_;
95  volScalarField& k_;
96  volScalarField kRes_;
98 
103  scalarList omegaNearWall_;
104 
106  label solveTurbState_ = 0;
107 
110 
112  volScalarField& betaFieldInversion_;
113 
115  volScalarField surfaceFriction_;
116 
118  volScalarField surfaceFrictionData_;
119 
121  volScalarField pData_;
122 
124  volVectorField UData_;
125 
127  volScalarField USingleComponentData_;
128 
129 public:
130  TypeName("kOmegaFieldInversionOmega");
131  // Constructors
132 
133  //- Construct from components
135  const word modelType,
136  const fvMesh& mesh,
137  const DAOption& daOption);
138 
139  //- Destructor
141  {
142  }
143 
144  // Member functions
145 
147  virtual void correctModelStates(wordList& modelStates) const;
148 
150  virtual void correctNut();
151 
153  virtual void correctBoundaryConditions();
154 
156  virtual void updateIntermediateVariables();
157 
159  virtual void correctStateResidualModelCon(List<List<word>>& stateCon) const;
160 
162  virtual void addModelResidualCon(HashTable<List<List<word>>>& allCon) const;
163 
165  virtual void calcResiduals(const dictionary& options);
166 
168  virtual void correct();
169 
171  void saveOmegaNearWall();
172 
174  void setOmegaNearWall();
175 
178 };
179 
180 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
181 
182 } // End namespace Foam
183 
184 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
185 
186 #endif
187 
188 // ************************************************************************* //
Foam::DAkOmegaFieldInversionOmega::UData_
volVectorField UData_
reference field (e.g. velocity for DNS)
Definition: DAkOmegaFieldInversionOmega.H:124
Foam::DAkOmegaFieldInversionOmega::betaFieldInversion_
volScalarField & betaFieldInversion_
A beta field multiplying to the production term.
Definition: DAkOmegaFieldInversionOmega.H:112
Foam::DAkOmegaFieldInversionOmega::updateIntermediateVariables
virtual void updateIntermediateVariables()
update any intermediate variables that are dependent on state variables and are used in calcResiduals
Definition: DAkOmegaFieldInversionOmega.C:305
Foam::DAkOmegaFieldInversionOmega::DomegaEff
tmp< volScalarField > DomegaEff() const
Definition: DAkOmegaFieldInversionOmega.H:78
Foam::DAOption
Definition: DAOption.H:29
Foam::DAkOmegaFieldInversionOmega::pData_
volScalarField pData_
the reference field for surfacePressure
Definition: DAkOmegaFieldInversionOmega.H:121
Foam::DAkOmegaFieldInversionOmega::gamma_
dimensionedScalar gamma_
Definition: DAkOmegaFieldInversionOmega.H:57
Foam::DAkOmegaFieldInversionOmega::saveOmegaNearWall
void saveOmegaNearWall()
save near wall omega values to omegaNearWall_
Definition: DAkOmegaFieldInversionOmega.C:255
daOption
DAOption daOption(mesh, pyOptions_)
Foam::DAkOmegaFieldInversionOmega::printInterval_
label printInterval_
time step interval to print residual
Definition: DAkOmegaFieldInversionOmega.H:109
Foam::DAkOmegaFieldInversionOmega::correctModelStates
virtual void correctModelStates(wordList &modelStates) const
update the turbulence state for DAStateInfo::regStates_
Definition: DAkOmegaFieldInversionOmega.C:162
Foam::DAkOmegaFieldInversionOmega::omegaNearWall_
scalarList omegaNearWall_
Definition: DAkOmegaFieldInversionOmega.H:103
Foam::DAkOmegaFieldInversionOmega::omega_
volScalarField & omega_
Definition: DAkOmegaFieldInversionOmega.H:93
Foam::DAkOmegaFieldInversionOmega::surfaceFrictionData_
volScalarField surfaceFrictionData_
the reference field for surfaceFriction
Definition: DAkOmegaFieldInversionOmega.H:118
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
Foam::DAkOmegaFieldInversionOmega::Cmu_
dimensionedScalar Cmu_
Definition: DAkOmegaFieldInversionOmega.H:55
Foam::DAkOmegaFieldInversionOmega::surfaceFriction_
volScalarField surfaceFriction_
a surface friction 'field' when using skin friction data for field inversion
Definition: DAkOmegaFieldInversionOmega.H:115
Foam::DAkOmegaFieldInversionOmega::kRes_
volScalarField kRes_
Definition: DAkOmegaFieldInversionOmega.H:96
Foam::DAkOmegaFieldInversionOmega::DkEff
tmp< volScalarField > DkEff() const
Definition: DAkOmegaFieldInversionOmega.H:65
Foam::DATurbulenceModel::nut_
volScalarField & nut_
turbulence viscosity
Definition: DATurbulenceModel.H:83
Foam::DAkOmegaFieldInversionOmega::alphaK_
dimensionedScalar alphaK_
Definition: DAkOmegaFieldInversionOmega.H:58
Foam
Definition: multiFreqScalarFvPatchField.C:144
Foam::DAkOmegaFieldInversionOmega::DAkOmegaFieldInversionOmega
DAkOmegaFieldInversionOmega(const word modelType, const fvMesh &mesh, const DAOption &daOption)
Definition: DAkOmegaFieldInversionOmega.C:41
Foam::DAkOmegaFieldInversionOmega::calcResiduals
virtual void calcResiduals(const dictionary &options)
compute the turbulence residuals
Definition: DAkOmegaFieldInversionOmega.C:527
Foam::DAkOmegaFieldInversionOmega::k_
volScalarField & k_
Definition: DAkOmegaFieldInversionOmega.H:95
Foam::DAkOmegaFieldInversionOmega::omegaRes_
volScalarField omegaRes_
Definition: DAkOmegaFieldInversionOmega.H:94
Foam::DAkOmegaFieldInversionOmega::correct
virtual void correct()
solve the residual equations and update the state
Definition: DAkOmegaFieldInversionOmega.C:507
Foam::DATurbulenceModel::nu
tmp< volScalarField > nu() const
get the nu field
Definition: DATurbulenceModel.C:262
Foam::DAkOmegaFieldInversionOmega::correctStateResidualModelCon
virtual void correctStateResidualModelCon(List< List< word >> &stateCon) const
update the original variable connectivity for the adjoint state residuals in stateCon
Definition: DAkOmegaFieldInversionOmega.C:316
Foam::DAkOmegaFieldInversionOmega::USingleComponentData_
volScalarField USingleComponentData_
the reference pressure field data
Definition: DAkOmegaFieldInversionOmega.H:127
Foam::DAkOmegaFieldInversionOmega::correctBoundaryConditions
virtual void correctBoundaryConditions()
update turbulence variable boundary values
Definition: DAkOmegaFieldInversionOmega.C:216
Foam::DAkOmegaFieldInversionOmega::alphaOmega_
dimensionedScalar alphaOmega_
Definition: DAkOmegaFieldInversionOmega.H:59
Foam::DATurbulenceModel
Definition: DATurbulenceModel.H:61
Foam::DAkOmegaFieldInversionOmega::beta_
dimensionedScalar beta_
Definition: DAkOmegaFieldInversionOmega.H:56
Foam::DAkOmegaFieldInversionOmega
Definition: DAkOmegaFieldInversionOmega.H:48
Foam::DAkOmegaFieldInversionOmega::setOmegaNearWall
void setOmegaNearWall()
set omegaNearWall_ to near wall omega values
Definition: DAkOmegaFieldInversionOmega.C:279
Foam::DAkOmegaFieldInversionOmega::solveTurbState_
label solveTurbState_
whether to solve for turb states
Definition: DAkOmegaFieldInversionOmega.H:106
Foam::DAkOmegaFieldInversionOmega::addModelResidualCon
virtual void addModelResidualCon(HashTable< List< List< word >>> &allCon) const
add the model residual connectivity to stateCon
Definition: DAkOmegaFieldInversionOmega.C:415
Foam::DAkOmegaFieldInversionOmega::correctNut
virtual void correctNut()
update nut based on other turbulence variables and update the BCs
Definition: DAkOmegaFieldInversionOmega.C:198
Foam::DAkOmegaFieldInversionOmega::TypeName
TypeName("kOmegaFieldInversionOmega")
DATurbulenceModel.H
Foam::DAkOmegaFieldInversionOmega::correctOmegaBoundaryConditions
void correctOmegaBoundaryConditions()
specially treatment to correct epsilon BC
Definition: DAkOmegaFieldInversionOmega.C:228
Foam::DAkOmegaFieldInversionOmega::~DAkOmegaFieldInversionOmega
virtual ~DAkOmegaFieldInversionOmega()
Definition: DAkOmegaFieldInversionOmega.H:140