DATurbulenceModel.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v4
5 
6  Description:
7  Augmented turbulence model for the adjoint method, including residual
8  calculation functions, etc
9 
10  NOTE 1:
11  Instead of inheriting from the OpenFOAM turbulence implementation, in
12  RASModel's children, we re-write all the correspndong functions for each
13  turbulence model. This is to avoid using template classes and template
14  functions for all the other classes in DAFoam. The downside is that we
15  need to update all the RASModel's children when upgrading to a new version
16  of OpenFOAM. Hopefully, the turbulence model part does not change too much
17  from version to version so the modification will be minimal.
18 
19  NOTE 2:
20  This function will lookup turbulence model object in fvMesh, so make sure
21  you initialize a turbulence model BEFORE initialzing DATurbulenceModel
22 
23 \*---------------------------------------------------------------------------*/
24 
25 #ifndef DATurbulenceModel_H
26 #define DATurbulenceModel_H
27 
28 #include "runTimeSelectionTables.H"
29 #include "wallDist.H"
30 #include "nearWallDist.H"
31 #include "DAUtility.H"
32 #include "DAOption.H"
33 #include "DAMacroFunctions.H"
34 #include "DAGlobalVar.H"
35 #include "nutkWallFunctionFvPatchScalarField.H" // for yPlus
36 #include "wallFvPatch.H" // for yPlus
37 
38 #include "singlePhaseTransportModel.H"
39 #include "turbulentTransportModel.H"
40 #include "fluidThermo.H"
41 #include "turbulentFluidThermoModel.H"
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 namespace Foam
46 {
47 
48 /*---------------------------------------------------------------------------*\
49  Class DATurbulenceModel Declaration
50 \*---------------------------------------------------------------------------*/
51 
53  : public regIOobject
54 {
55 
56 private:
57  //- Disallow default bitwise copy construct
59 
60  //- Disallow default bitwise assignment
61  void operator=(const DATurbulenceModel&);
62 
63 protected:
65  const fvMesh& mesh_;
66 
69 
71  const dictionary& allOptions_;
72 
75 
77  volScalarField& nut_;
78 
80  volVectorField& U_;
81 
83  surfaceScalarField& phi_;
84 
86  volScalarField phase_;
87 
89  surfaceScalarField& phaseRhoPhi_;
90 
92  volScalarField rhoOne_;
93 
95  word turbModelType_ = "None";
96 
98  IOdictionary turbDict_;
99 
101  dictionary coeffDict_;
102 
104  dimensionedScalar kMin_;
105 
107  dimensionedScalar epsilonMin_;
108 
110  dimensionedScalar omegaMin_;
111 
113  dimensionedScalar nuTildaMin_;
114 
116  scalar Pr_;
117 
119  scalar Prt_ = -9999.0;
120 
121 public:
122  //- Runtime type information
123  TypeName("DATurbulenceModel");
124 
125  // Declare run-time constructor selection table
127  autoPtr,
129  dictionary,
130  (const word modelType,
131  const fvMesh& mesh,
132  const DAOption& daOption),
133  (modelType, mesh, daOption));
134 
135  // Constructors
136 
137  //- Construct from components
139  const word modelType,
140  const fvMesh& mesh,
141  const DAOption& daOption);
142 
143  // Selectors
144 
145  //- Return a reference to the selected model
146  static autoPtr<DATurbulenceModel> New(
147  const word modelType,
148  const fvMesh& mesh,
149  const DAOption& daOption);
150 
151  //- Destructor
153  {
154  }
155 
156  // Members
157 
159  void correctWallDist();
160 
162  void correctAlphat();
163 
165  virtual void correctNut() = 0;
166 
168  virtual void correctModelStates(wordList& modelStates) const = 0;
169 
171  virtual void correctBoundaryConditions() = 0;
172 
174  virtual void updateIntermediateVariables() = 0;
175 
177  virtual void correctStateResidualModelCon(List<List<word>>& stateCon) const = 0;
178 
180  virtual void addModelResidualCon(HashTable<List<List<word>>>& allCon) const = 0;
181 
183  virtual void calcResiduals(const dictionary& options) = 0;
184 
186  virtual void correct(label printToScreen) = 0;
187 
189  virtual void getTurbProdTerm(volScalarField& prodTerm) const;
190 
192  virtual void getTurbProdOverDestruct(volScalarField& PoD) const;
193 
195  virtual void getTurbConvOverProd(volScalarField& CoP) const;
196 
198  tmp<volSymmTensorField> devRhoReff() const;
199 
201  tmp<fvVectorMatrix> divDevRhoReff(volVectorField& U);
202 
204  tmp<fvVectorMatrix> divDevReff(volVectorField& U);
205 
207  tmp<volScalarField> nuEff() const;
208 
210  tmp<volScalarField> nut()
211  {
212  return nut_;
213  }
214 
216  tmp<volScalarField> alphaEff();
217 
219  tmp<volScalarField> nu() const;
220 
222  tmp<volScalarField> alpha() const;
223 
225  tmp<volScalarField> rho() const;
226 
228  dimensionSet rhoDimensions() const;
229 
231  tmp<volScalarField> phase()
232  {
233  return phase_;
234  }
235 
237  scalar prt()
238  {
239  if (Prt_ > 0)
240  {
241  return Prt_;
242  }
243  else
244  {
245  // alphat not found in Db so Prt is not set
246  FatalErrorIn("getPrt") << exit(FatalError);
247  return -9999.0;
248  }
249  }
250 
252  tmp<Foam::volScalarField> mu() const;
253 
255  bool writeData(Ostream& os) const;
256 
258  void printYPlus(const label printToScreen) const;
259 
260  word getTurbModelType() const
261  {
262  return turbModelType_;
263  }
264 
265  label isPrintTime(
266  const Time& runTime,
267  const label printInterval) const;
268 
270  virtual void invTranProdNuTildaEqn(
271  const volScalarField& mySource,
272  volScalarField& pseudoNuTilda);
273 
274  virtual void constructPseudoNuTildaEqn();
275 
276  virtual void rhsSolvePseudoNuTildaEqn(const volScalarField& nuTildaSource);
277 
279  virtual void calcLduResidualTurb(volScalarField& nuTildaRes);
280 
282  virtual void getFvMatrixFields(
283  const word varName,
284  scalarField& diag,
285  scalarField& upper,
286  scalarField& lower);
287 
289  virtual void solveAdjointFP(
290  const word varName,
291  const scalarList& rhs,
292  scalarList& dPsi);
293 };
294 
295 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
296 
297 } // End namespace Foam
298 
299 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
300 
301 #endif
302 
303 // ************************************************************************* //
Foam::DATurbulenceModel::getTurbProdTerm
virtual void getTurbProdTerm(volScalarField &prodTerm) const
return the value of the production term from the turbulence model
Definition: DATurbulenceModel.C:540
U
volVectorField & U
Definition: createRefsPimpleDyM.H:7
DAOption.H
Foam::DATurbulenceModel::nuEff
tmp< volScalarField > nuEff() const
return effective viscosity
Definition: DATurbulenceModel.C:220
Foam::DATurbulenceModel::constructPseudoNuTildaEqn
virtual void constructPseudoNuTildaEqn()
Definition: DATurbulenceModel.C:602
Foam::DATurbulenceModel::phi_
surfaceScalarField & phi_
face flux
Definition: DATurbulenceModel.H:83
Foam::DATurbulenceModel::divDevRhoReff
tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U)
divDev terms
Definition: DATurbulenceModel.C:389
Foam::DATurbulenceModel::divDevReff
tmp< fvVectorMatrix > divDevReff(volVectorField &U)
divDev terms
Definition: DATurbulenceModel.C:411
Foam::DATurbulenceModel::kMin_
dimensionedScalar kMin_
Lower limit of k.
Definition: DATurbulenceModel.H:104
Foam::DATurbulenceModel::updateIntermediateVariables
virtual void updateIntermediateVariables()=0
update any intermediate variables that are dependent on state variables and are used in calcResiduals
Foam::DAOption
Definition: DAOption.H:29
Foam::DATurbulenceModel::correctModelStates
virtual void correctModelStates(wordList &modelStates) const =0
update the turbulence state for DAStateInfo::regStates_
Foam::DATurbulenceModel::TypeName
TypeName("DATurbulenceModel")
DAUtility.H
pseudoNuTilda
volScalarField & pseudoNuTilda
Definition: pseudoEqns.H:34
Foam::DATurbulenceModel::nuTildaMin_
dimensionedScalar nuTildaMin_
Lower limit for nuTilda.
Definition: DATurbulenceModel.H:113
Foam::DATurbulenceModel::correctNut
virtual void correctNut()=0
update nut based on other turbulence variables and update the BCs
Foam::DATurbulenceModel::phase
tmp< volScalarField > phase()
get the phase field
Definition: DATurbulenceModel.H:231
DAGlobalVar.H
Foam::DATurbulenceModel::phaseRhoPhi_
surfaceScalarField & phaseRhoPhi_
phase*phi*density field
Definition: DATurbulenceModel.H:89
Foam::DATurbulenceModel::allOptions_
const dictionary & allOptions_
all DAFoam options
Definition: DATurbulenceModel.H:71
Foam::DATurbulenceModel::New
static autoPtr< DATurbulenceModel > New(const word modelType, const fvMesh &mesh, const DAOption &daOption)
Definition: DATurbulenceModel.C:161
Foam::DATurbulenceModel::rho
tmp< volScalarField > rho() const
get the density field
Definition: DATurbulenceModel.C:294
DAMacroFunctions.H
Foam::DATurbulenceModel::prt
scalar prt()
get the turbulent Prandtl number
Definition: DATurbulenceModel.H:237
Foam::DATurbulenceModel::mu
tmp< Foam::volScalarField > mu() const
get mu
Definition: DATurbulenceModel.C:349
Foam::DATurbulenceModel::Prt_
scalar Prt_
turbulent Prandtl number
Definition: DATurbulenceModel.H:119
Foam::DATurbulenceModel::correctBoundaryConditions
virtual void correctBoundaryConditions()=0
update turbulence variable boundary values
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
Foam::DATurbulenceModel::calcLduResidualTurb
virtual void calcLduResidualTurb(volScalarField &nuTildaRes)
calculate the turbulence residual using LDU matrix
Definition: DATurbulenceModel.C:590
Foam::DATurbulenceModel::declareRunTimeSelectionTable
declareRunTimeSelectionTable(autoPtr, DATurbulenceModel, dictionary,(const word modelType, const fvMesh &mesh, const DAOption &daOption),(modelType, mesh, daOption))
Foam::DATurbulenceModel::nut
tmp< volScalarField > nut()
get the nut field
Definition: DATurbulenceModel.H:210
Foam::DATurbulenceModel::calcResiduals
virtual void calcResiduals(const dictionary &options)=0
compute the turbulence residuals
Foam::DATurbulenceModel::correctWallDist
void correctWallDist()
update wall distance for d_. Note: y_ will be automatically updated in mesh_ object
Definition: DATurbulenceModel.C:422
Foam::DATurbulenceModel::daGlobalVar_
DAGlobalVar & daGlobalVar_
global variable
Definition: DATurbulenceModel.H:74
Foam::DATurbulenceModel::U_
volVectorField & U_
velocity
Definition: DATurbulenceModel.H:80
Foam::DATurbulenceModel::turbDict_
IOdictionary turbDict_
turbulence model property dict
Definition: DATurbulenceModel.H:98
Foam::DATurbulenceModel::nut_
volScalarField & nut_
turbulence viscosity
Definition: DATurbulenceModel.H:77
Foam
Definition: checkGeometry.C:32
Foam::DATurbulenceModel::getTurbConvOverProd
virtual void getTurbConvOverProd(volScalarField &CoP) const
return the value of the convective over production term from the turbulence model
Definition: DATurbulenceModel.C:564
Foam::DATurbulenceModel::invTranProdNuTildaEqn
virtual void invTranProdNuTildaEqn(const volScalarField &mySource, volScalarField &pseudoNuTilda)
Inverse transpose product, M_nuTilda^(-T)
Definition: DATurbulenceModel.C:576
Foam::DATurbulenceModel::devRhoReff
tmp< volSymmTensorField > devRhoReff() const
dev terms
Definition: DATurbulenceModel.C:371
Foam::DATurbulenceModel::rhoDimensions
dimensionSet rhoDimensions() const
return the dimension of rho
Definition: DATurbulenceModel.C:321
Foam::DATurbulenceModel::rhoOne_
volScalarField rhoOne_
an uniform rho field with ones
Definition: DATurbulenceModel.H:92
Foam::DATurbulenceModel::~DATurbulenceModel
virtual ~DATurbulenceModel()
Definition: DATurbulenceModel.H:152
Foam::DATurbulenceModel::rhsSolvePseudoNuTildaEqn
virtual void rhsSolvePseudoNuTildaEqn(const volScalarField &nuTildaSource)
Definition: DATurbulenceModel.C:614
Foam::DATurbulenceModel::nu
tmp< volScalarField > nu() const
get the nu field
Definition: DATurbulenceModel.C:267
Foam::DATurbulenceModel::getTurbModelType
word getTurbModelType() const
Definition: DATurbulenceModel.H:260
Foam::DATurbulenceModel::addModelResidualCon
virtual void addModelResidualCon(HashTable< List< List< word >>> &allCon) const =0
add the model residual connectivity to stateCon
Foam::DATurbulenceModel::printYPlus
void printYPlus(const label printToScreen) const
print the min max and mean yPlus to screen
Definition: DATurbulenceModel.C:438
Foam::DATurbulenceModel::correctAlphat
void correctAlphat()
update alphat
Definition: DATurbulenceModel.C:203
Foam::DATurbulenceModel
Definition: DATurbulenceModel.H:52
Foam::DATurbulenceModel::getFvMatrixFields
virtual void getFvMatrixFields(const word varName, scalarField &diag, scalarField &upper, scalarField &lower)
return the diag(), upper(), and lower() scalarFields from the turbulence model's fvMatrix
Definition: DATurbulenceModel.C:627
Foam::DATurbulenceModel::daOption_
const DAOption & daOption_
DAOption object.
Definition: DATurbulenceModel.H:68
Foam::DATurbulenceModel::correct
virtual void correct(label printToScreen)=0
solve the residual equations and update the state
Foam::DATurbulenceModel::mesh_
const fvMesh & mesh_
fvMesh
Definition: DATurbulenceModel.H:65
Foam::DATurbulenceModel::correctStateResidualModelCon
virtual void correctStateResidualModelCon(List< List< word >> &stateCon) const =0
update the original variable connectivity for the adjoint state residuals in stateCon
Foam::DATurbulenceModel::solveAdjointFP
virtual void solveAdjointFP(const word varName, const scalarList &rhs, scalarList &dPsi)
solve the fvMatrixT field with given rhs and solution
Definition: DATurbulenceModel.C:639
Foam::DATurbulenceModel::isPrintTime
label isPrintTime(const Time &runTime, const label printInterval) const
Definition: DATurbulenceModel.C:521
Foam::DATurbulenceModel::alphaEff
tmp< volScalarField > alphaEff()
return effective thermal diffusivity
Definition: DATurbulenceModel.C:233
Foam::DATurbulenceModel::omegaMin_
dimensionedScalar omegaMin_
Lower limit for omega.
Definition: DATurbulenceModel.H:110
Foam::DATurbulenceModel::coeffDict_
dictionary coeffDict_
turbulence model parameters dict
Definition: DATurbulenceModel.H:101
Foam::DATurbulenceModel::alpha
tmp< volScalarField > alpha() const
get alpha field
Definition: DATurbulenceModel.C:340
runTime
Time & runTime
Definition: createRefsHeatTransfer.H:1
Foam::DAGlobalVar
Definition: DAGlobalVar.H:26
Foam::DATurbulenceModel::Pr_
scalar Pr_
Prandtl number.
Definition: DATurbulenceModel.H:116
Foam::DATurbulenceModel::getTurbProdOverDestruct
virtual void getTurbProdOverDestruct(volScalarField &PoD) const
return the ratio of the production over destruction term from the turbulence model
Definition: DATurbulenceModel.C:552
Foam::DATurbulenceModel::turbModelType_
word turbModelType_
whether the turbulence model is incompressible or compressible
Definition: DATurbulenceModel.H:95
Foam::DATurbulenceModel::phase_
volScalarField phase_
phase field
Definition: DATurbulenceModel.H:86
Foam::DATurbulenceModel::epsilonMin_
dimensionedScalar epsilonMin_
Lower limit of epsilon.
Definition: DATurbulenceModel.H:107
Foam::DATurbulenceModel::writeData
bool writeData(Ostream &os) const
this is a virtual function for regIOobject
Definition: DATurbulenceModel.C:197