DAkEpsilon.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 kEpsilon model
8 
9  This file is modified from OpenFOAM's source code
10  src/TurbulenceModels/turbulenceModels/RAS/kEpsilon/kEpsilon.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 DAkEpsilon_H
34 #define DAkEpsilon_H
35 
36 #include "DATurbulenceModel.H"
37 #include "addToRunTimeSelectionTable.H"
38 
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43 
44 /*---------------------------------------------------------------------------*\
45  Class DAkEpsilon Declaration
46 \*---------------------------------------------------------------------------*/
47 
49  : public DATurbulenceModel
50 {
51 
52 protected:
54 
55  dimensionedScalar Cmu_;
56  dimensionedScalar C1_;
57  dimensionedScalar C2_;
58  dimensionedScalar C3_;
59  dimensionedScalar sigmak_;
60  dimensionedScalar sigmaEps_;
62 
64 
65  tmp<fvScalarMatrix> kSource() const;
66  tmp<fvScalarMatrix> epsilonSource() const;
67  //- Return the effective diffusivity for k
68  tmp<volScalarField> DkEff() const;
69  //- Return the effective diffusivity for epsilon
70  tmp<volScalarField> DepsilonEff() const;
71  //- Return the turbulence kinetic energy
72  virtual tmp<volScalarField> k() const
73  {
74  return k_;
75  }
76 
77  //- Return the turbulence kinetic energy dissipation rate
78  virtual tmp<volScalarField> epsilon() const
79  {
80  return epsilon_;
81  }
83 
85 
86  volScalarField& epsilon_;
87  volScalarField epsilonRes_;
88  volScalarField& k_;
89  volScalarField kRes_;
91 
96  scalarList epsilonNearWall_;
97 
99  label solveTurbState_ = 0;
100 
103 
104 public:
105  TypeName("kEpsilon");
106  // Constructors
107 
108  //- Construct from components
109  DAkEpsilon(
110  const word modelType,
111  const fvMesh& mesh,
112  const DAOption& daOption);
113 
114  //- Destructor
115  virtual ~DAkEpsilon()
116  {
117  }
118 
119  // Member functions
120 
122  virtual void correctModelStates(wordList& modelStates) const;
123 
125  virtual void correctNut();
126 
128  virtual void correctBoundaryConditions();
129 
131  virtual void updateIntermediateVariables();
132 
134  virtual void correctStateResidualModelCon(List<List<word>>& stateCon) const;
135 
137  virtual void addModelResidualCon(HashTable<List<List<word>>>& allCon) const;
138 
140  virtual void calcResiduals(const dictionary& options);
141 
143  virtual void correct();
144 
146  void saveEpsilonNearWall();
147 
149  void setEpsilonNearWall();
150 
153 };
154 
155 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
156 
157 } // End namespace Foam
158 
159 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
160 
161 #endif
162 
163 // ************************************************************************* //
Foam::DAkEpsilon::DkEff
tmp< volScalarField > DkEff() const
Definition: DAkEpsilon.C:170
Foam::DAkEpsilon::epsilonNearWall_
scalarList epsilonNearWall_
Definition: DAkEpsilon.H:96
Foam::DAkEpsilon::Cmu_
dimensionedScalar Cmu_
Definition: DAkEpsilon.H:55
Foam::DAkEpsilon::TypeName
TypeName("kEpsilon")
Foam::DAkEpsilon::addModelResidualCon
virtual void addModelResidualCon(HashTable< List< List< word >>> &allCon) const
add the model residual connectivity to stateCon
Definition: DAkEpsilon.C:392
Foam::DAkEpsilon::epsilon_
volScalarField & epsilon_
Definition: DAkEpsilon.H:86
Foam::DAOption
Definition: DAOption.H:29
Foam::DAkEpsilon::correctBoundaryConditions
virtual void correctBoundaryConditions()
update turbulence variable boundary values
Definition: DAkEpsilon.C:241
Foam::DAkEpsilon::~DAkEpsilon
virtual ~DAkEpsilon()
Definition: DAkEpsilon.H:115
daOption
DAOption daOption(mesh, pyOptions_)
Foam::DAkEpsilon::correct
virtual void correct()
solve the residual equations and update the state
Definition: DAkEpsilon.C:484
Foam::DAkEpsilon::calcResiduals
virtual void calcResiduals(const dictionary &options)
compute the turbulence residuals
Definition: DAkEpsilon.C:504
Foam::DAkEpsilon::printInterval_
label printInterval_
time step interval to print residual
Definition: DAkEpsilon.H:102
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
Foam::DAkEpsilon::k
virtual tmp< volScalarField > k() const
Definition: DAkEpsilon.H:72
Foam::DAkEpsilon::correctModelStates
virtual void correctModelStates(wordList &modelStates) const
update the turbulence state for DAStateInfo::regStates_
Definition: DAkEpsilon.C:187
Foam::DAkEpsilon::sigmaEps_
dimensionedScalar sigmaEps_
Definition: DAkEpsilon.H:60
Foam::DAkEpsilon::epsilon
virtual tmp< volScalarField > epsilon() const
Definition: DAkEpsilon.H:78
Foam::DAkEpsilon::epsilonRes_
volScalarField epsilonRes_
Definition: DAkEpsilon.H:87
Foam
Definition: multiFreqScalarFvPatchField.C:144
Foam::DAkEpsilon::C1_
dimensionedScalar C1_
Definition: DAkEpsilon.H:56
Foam::DAkEpsilon
Definition: DAkEpsilon.H:48
Foam::DAkEpsilon::kSource
tmp< fvScalarMatrix > kSource() const
Definition: DAkEpsilon.C:152
Foam::DAkEpsilon::correctEpsilonBoundaryConditions
void correctEpsilonBoundaryConditions()
specially treatment to correct epsilon BC
Definition: DAkEpsilon.C:253
Foam::DAkEpsilon::setEpsilonNearWall
void setEpsilonNearWall()
set epsilonNearWall_ to near wall epsilon values
Definition: DAkEpsilon.C:304
Foam::DAkEpsilon::C3_
dimensionedScalar C3_
Definition: DAkEpsilon.H:58
Foam::DAkEpsilon::correctNut
virtual void correctNut()
update nut based on other turbulence variables and update the BCs
Definition: DAkEpsilon.C:223
Foam::DAkEpsilon::k_
volScalarField & k_
Definition: DAkEpsilon.H:88
Foam::DAkEpsilon::epsilonSource
tmp< fvScalarMatrix > epsilonSource() const
Definition: DAkEpsilon.C:161
Foam::DAkEpsilon::kRes_
volScalarField kRes_
Definition: DAkEpsilon.H:89
Foam::DATurbulenceModel
Definition: DATurbulenceModel.H:61
Foam::DAkEpsilon::correctStateResidualModelCon
virtual void correctStateResidualModelCon(List< List< word >> &stateCon) const
update the original variable connectivity for the adjoint state residuals in stateCon
Definition: DAkEpsilon.C:341
Foam::DAkEpsilon::DAkEpsilon
DAkEpsilon(const word modelType, const fvMesh &mesh, const DAOption &daOption)
Definition: DAkEpsilon.C:41
Foam::DAkEpsilon::saveEpsilonNearWall
void saveEpsilonNearWall()
save near wall epsilon values to epsilonNearWall_
Definition: DAkEpsilon.C:280
Foam::DAkEpsilon::updateIntermediateVariables
virtual void updateIntermediateVariables()
update any intermediate variables that are dependent on state variables and are used in calcResiduals
Definition: DAkEpsilon.C:330
Foam::DAkEpsilon::solveTurbState_
label solveTurbState_
whether to solve for turb states
Definition: DAkEpsilon.H:99
Foam::DAkEpsilon::sigmak_
dimensionedScalar sigmak_
Definition: DAkEpsilon.H:59
Foam::DAkEpsilon::DepsilonEff
tmp< volScalarField > DepsilonEff() const
Definition: DAkEpsilon.C:178
Foam::DAkEpsilon::C2_
dimensionedScalar C2_
Definition: DAkEpsilon.H:57
DATurbulenceModel.H