DASpalartAllmarasFv3.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 SpalartAllmarasFv3 model used in OpenFOAM 2.4 and before
8 
9  This file is modified from OpenFOAM's source code
10  src/TurbulenceModels/turbulenceModels/RAS/SpalartAllmaras/SpalartAllmaras.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 DASpalartAllmarasFv3_H
34 #define DASpalartAllmarasFv3_H
35 
36 #include "DATurbulenceModel.H"
37 #include "addToRunTimeSelectionTable.H"
38 
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43 
44 /*---------------------------------------------------------------------------*\
45  Class DASpalartAllmarasFv3 Declaration
46 \*---------------------------------------------------------------------------*/
47 
49  : public DATurbulenceModel
50 {
51 
52 protected:
54 
55  dimensionedScalar sigmaNut_;
56  dimensionedScalar kappa_;
57  dimensionedScalar Cb1_;
58  dimensionedScalar Cb2_;
59  dimensionedScalar Cw1_;
60  dimensionedScalar Cw2_;
61  dimensionedScalar Cw3_;
62  dimensionedScalar Cv1_;
63  dimensionedScalar Cv2_;
65 
67 
68  tmp<volScalarField> chi() const;
69 
70  tmp<volScalarField> fv1(const volScalarField& chi) const;
71 
72  tmp<volScalarField> fv2(
73  const volScalarField& chi,
74  const volScalarField& fv1) const;
75 
76  tmp<volScalarField> fv3(
77  const volScalarField& chi,
78  const volScalarField& fv1) const;
79 
80  tmp<volScalarField> fw(const volScalarField& Stilda) const;
81 
83 
85 
86  volScalarField& nuTilda_;
87  volScalarField nuTildaRes_;
89 
91  volScalarField pseudoNuTilda_;
92  fvScalarMatrix pseudoNuTildaEqn_;
93 
95  const volScalarField& y_;
96 
98  label solveTurbState_ = 0;
99 
102 
103  // fvSolutions parameters
104  scalar relaxNuTildaEqn_ = 1.0;
105  dictionary solverDictNuTilda_;
106 
107 public:
108  TypeName("SpalartAllmarasFv3");
109  // Constructors
110 
111  //- Construct from components
113  const word modelType,
114  const fvMesh& mesh,
115  const DAOption& daOption);
116 
117  //- Destructor
119  {
120  }
121 
122  // Member functions
123 
125  tmp<volScalarField> DnuTildaEff() const;
126 
128  virtual void correctModelStates(wordList& modelStates) const;
129 
131  virtual void correctNut();
132 
134  virtual void correctBoundaryConditions();
135 
137  virtual void updateIntermediateVariables();
138 
140  virtual void correctStateResidualModelCon(List<List<word>>& stateCon) const;
141 
143  virtual void addModelResidualCon(HashTable<List<List<word>>>& allCon) const;
144 
146  virtual void calcResiduals(const dictionary& options);
147 
149  virtual void correct();
150 
152  virtual void invTranProdNuTildaEqn(
153  const volScalarField& mySource,
154  volScalarField& pseudoNuTilda);
155 
157  virtual void constructPseudoNuTildaEqn();
158 
160  virtual void rhsSolvePseudoNuTildaEqn(const volScalarField& nuTildaSource);
161 
163  virtual void calcLduResidualTurb(volScalarField& nuTildaRes);
164 };
165 
166 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
167 
168 } // End namespace Foam
169 
170 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
171 
172 #endif
173 
174 // ************************************************************************* //
Foam::DASpalartAllmarasFv3::solverDictNuTilda_
dictionary solverDictNuTilda_
Definition: DASpalartAllmarasFv3.H:105
Foam::DASpalartAllmarasFv3::y_
const volScalarField & y_
3D wall distance
Definition: DASpalartAllmarasFv3.H:95
Foam::DASpalartAllmarasFv3::fv1
tmp< volScalarField > fv1(const volScalarField &chi) const
Definition: DASpalartAllmarasFv3.C:157
Foam::DASpalartAllmarasFv3::calcResiduals
virtual void calcResiduals(const dictionary &options)
compute the turbulence residuals
Definition: DASpalartAllmarasFv3.C:433
Foam::DAOption
Definition: DAOption.H:29
Foam::DASpalartAllmarasFv3::nuTilda_
volScalarField & nuTilda_
Definition: DASpalartAllmarasFv3.H:86
Foam::DASpalartAllmarasFv3::~DASpalartAllmarasFv3
virtual ~DASpalartAllmarasFv3()
Definition: DASpalartAllmarasFv3.H:118
Foam::DASpalartAllmarasFv3::Cw3_
dimensionedScalar Cw3_
Definition: DASpalartAllmarasFv3.H:61
daOption
DAOption daOption(mesh, pyOptions_)
pseudoNuTilda
volScalarField & pseudoNuTilda
Definition: pseudoEqns.H:34
Foam::DASpalartAllmarasFv3::solveTurbState_
label solveTurbState_
whether to solve for turb states
Definition: DASpalartAllmarasFv3.H:98
Foam::DASpalartAllmarasFv3::Cw1_
dimensionedScalar Cw1_
Definition: DASpalartAllmarasFv3.H:59
Foam::DASpalartAllmarasFv3::correctModelStates
virtual void correctModelStates(wordList &modelStates) const
update the turbulence state for DAStateInfo::regStates_
Definition: DASpalartAllmarasFv3.C:209
Foam::DASpalartAllmarasFv3::relaxNuTildaEqn_
scalar relaxNuTildaEqn_
Definition: DASpalartAllmarasFv3.H:104
Foam::DASpalartAllmarasFv3::nuTildaRes_
volScalarField nuTildaRes_
Definition: DASpalartAllmarasFv3.H:87
Foam::DASpalartAllmarasFv3::Cb1_
dimensionedScalar Cb1_
Definition: DASpalartAllmarasFv3.H:57
Foam::DASpalartAllmarasFv3::Cv1_
dimensionedScalar Cv1_
Definition: DASpalartAllmarasFv3.H:62
Foam::DASpalartAllmarasFv3::correctNut
virtual void correctNut()
update nut based on other turbulence variables and update the BCs
Definition: DASpalartAllmarasFv3.C:243
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
Foam::DASpalartAllmarasFv3::fv3
tmp< volScalarField > fv3(const volScalarField &chi, const volScalarField &fv1) const
Definition: DASpalartAllmarasFv3.C:171
Foam::DASpalartAllmarasFv3::pseudoNuTildaEqn_
fvScalarMatrix pseudoNuTildaEqn_
Definition: DASpalartAllmarasFv3.H:92
Foam::DASpalartAllmarasFv3::addModelResidualCon
virtual void addModelResidualCon(HashTable< List< List< word >>> &allCon) const
add the model residual connectivity to stateCon
Definition: DASpalartAllmarasFv3.C:335
Foam::DASpalartAllmarasFv3::TypeName
TypeName("SpalartAllmarasFv3")
Foam::DASpalartAllmarasFv3::calcLduResidualTurb
virtual void calcLduResidualTurb(volScalarField &nuTildaRes)
calculate the turbulence residual using LDU matrix
Definition: DASpalartAllmarasFv3.C:665
Foam::DASpalartAllmarasFv3::kappa_
dimensionedScalar kappa_
Definition: DASpalartAllmarasFv3.H:56
Foam::DASpalartAllmarasFv3::updateIntermediateVariables
virtual void updateIntermediateVariables()
update any intermediate variables that are dependent on state variables and are used in calcResiduals
Definition: DASpalartAllmarasFv3.C:274
Foam::DASpalartAllmarasFv3::fw
tmp< volScalarField > fw(const volScalarField &Stilda) const
Definition: DASpalartAllmarasFv3.C:184
Foam::DASpalartAllmarasFv3::constructPseudoNuTildaEqn
virtual void constructPseudoNuTildaEqn()
construct pseudoNuTildaEqn_, which is nuTildaEqn with swapped upper and lower arrays
Definition: DASpalartAllmarasFv3.C:595
Foam::DASpalartAllmarasFv3::sigmaNut_
dimensionedScalar sigmaNut_
Definition: DASpalartAllmarasFv3.H:55
Foam::DASpalartAllmarasFv3::invTranProdNuTildaEqn
virtual void invTranProdNuTildaEqn(const volScalarField &mySource, volScalarField &pseudoNuTilda)
Inverse transpose product, A_nuTilda^(-T)
Definition: DASpalartAllmarasFv3.C:520
Foam::DASpalartAllmarasFv3::fv2
tmp< volScalarField > fv2(const volScalarField &chi, const volScalarField &fv1) const
Definition: DASpalartAllmarasFv3.C:164
Foam
Definition: multiFreqScalarFvPatchField.C:144
Foam::DASpalartAllmarasFv3::Cb2_
dimensionedScalar Cb2_
Definition: DASpalartAllmarasFv3.H:58
Foam::DATurbulenceModel
Definition: DATurbulenceModel.H:61
Foam::DASpalartAllmarasFv3::chi
tmp< volScalarField > chi() const
Definition: DASpalartAllmarasFv3.C:152
Foam::DASpalartAllmarasFv3::pseudoNuTilda_
volScalarField pseudoNuTilda_
pseudoNuTilda_ and pseudoNuTildaEqn_ for solving adjoint equation
Definition: DASpalartAllmarasFv3.H:91
Foam::DASpalartAllmarasFv3
Definition: DASpalartAllmarasFv3.H:48
Foam::DASpalartAllmarasFv3::Cv2_
dimensionedScalar Cv2_
Definition: DASpalartAllmarasFv3.H:63
Foam::DASpalartAllmarasFv3::DASpalartAllmarasFv3
DASpalartAllmarasFv3(const word modelType, const fvMesh &mesh, const DAOption &daOption)
Definition: DASpalartAllmarasFv3.C:41
Foam::DASpalartAllmarasFv3::Cw2_
dimensionedScalar Cw2_
Definition: DASpalartAllmarasFv3.H:60
Foam::DASpalartAllmarasFv3::printInterval_
label printInterval_
time step interval to print residual
Definition: DASpalartAllmarasFv3.H:101
Foam::DASpalartAllmarasFv3::correctBoundaryConditions
virtual void correctBoundaryConditions()
update turbulence variable boundary values
Definition: DASpalartAllmarasFv3.C:263
Foam::DASpalartAllmarasFv3::DnuTildaEff
tmp< volScalarField > DnuTildaEff() const
Return the effective diffusivity for nuTilda.
Definition: DASpalartAllmarasFv3.C:202
Foam::DASpalartAllmarasFv3::rhsSolvePseudoNuTildaEqn
virtual void rhsSolvePseudoNuTildaEqn(const volScalarField &nuTildaSource)
solve pseudoNuTildaEqn_ with overwritten rhs
Definition: DASpalartAllmarasFv3.C:633
DATurbulenceModel.H
Foam::DASpalartAllmarasFv3::correct
virtual void correct()
solve the residual equations and update the state
Definition: DASpalartAllmarasFv3.C:413
Foam::DASpalartAllmarasFv3::correctStateResidualModelCon
virtual void correctStateResidualModelCon(List< List< word >> &stateCon) const
update the original variable connectivity for the adjoint state residuals in stateCon
Definition: DASpalartAllmarasFv3.C:285