DATurbulenceModel.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v3
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 "DARegDbNearWallDist.H"
34 #include "DAMacroFunctions.H"
35 #include "nutkWallFunctionFvPatchScalarField.H" // for yPlus
36 #include "wallFvPatch.H" // for yPlus
37 
38 #ifdef IncompressibleFlow
39 #include "singlePhaseTransportModel.H"
40 #include "turbulentTransportModel.H"
43 #endif
44 
45 #ifdef CompressibleFlow
46 #include "fluidThermo.H"
47 #include "turbulentFluidThermoModel.H"
49 #include "DARegDbFluidThermo.H"
50 #endif
51 
52 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 
54 namespace Foam
55 {
56 
57 /*---------------------------------------------------------------------------*\
58  Class DATurbulenceModel Declaration
59 \*---------------------------------------------------------------------------*/
60 
62  : public regIOobject
63 {
64 
65 private:
66  //- Disallow default bitwise copy construct
68 
69  //- Disallow default bitwise assignment
70  void operator=(const DATurbulenceModel&);
71 
72 protected:
74  const fvMesh& mesh_;
75 
78 
80  const dictionary& allOptions_;
81 
83  volScalarField& nut_;
84 
86  volVectorField& U_;
87 
89  surfaceScalarField& phi_;
90 
92  volScalarField phase_;
93 
95  surfaceScalarField& phaseRhoPhi_;
96 
97 #ifdef IncompressibleFlow
98  const DARegDbSinglePhaseTransportModel& daRegDbTransport_;
101  const singlePhaseTransportModel& laminarTransport_;
103  const DARegDbTurbulenceModelIncompressible& daRegDbTurbIncomp_;
105  const incompressible::turbulenceModel& turbulence_;
107  volScalarField rho_;
108 #endif
109 
110 #ifdef CompressibleFlow
111  const DARegDbFluidThermo& daRegDbThermo_;
114  const fluidThermo& thermo_;
116  const DARegDbTurbulenceModelCompressible& daRegDbTurbComp_;
118  const compressible::turbulenceModel& turbulence_;
120  volScalarField& rho_;
121 #endif
122 
124  IOdictionary turbDict_;
125 
127  dictionary coeffDict_;
128 
130  dimensionedScalar kMin_;
131 
133  dimensionedScalar epsilonMin_;
134 
136  dimensionedScalar omegaMin_;
137 
139  dimensionedScalar nuTildaMin_;
140 
142  scalar Pr_;
143 
145  scalar Prt_ = -9999.0;
146 
147 public:
148  //- Runtime type information
149  TypeName("DATurbulenceModel");
150 
151  // Declare run-time constructor selection table
153  autoPtr,
155  dictionary,
156  (const word modelType,
157  const fvMesh& mesh,
158  const DAOption& daOption),
159  (modelType, mesh, daOption));
160 
161  // Constructors
162 
163  //- Construct from components
165  const word modelType,
166  const fvMesh& mesh,
167  const DAOption& daOption);
168 
169  // Selectors
170 
171  //- Return a reference to the selected model
172  static autoPtr<DATurbulenceModel> New(
173  const word modelType,
174  const fvMesh& mesh,
175  const DAOption& daOption);
176 
177  //- Destructor
179  {
180  }
181 
182  // Members
183 
185  void correctWallDist();
186 
188  void correctAlphat();
189 
191  virtual void correctNut() = 0;
192 
194  virtual void correctModelStates(wordList& modelStates) const = 0;
195 
197  virtual void correctBoundaryConditions() = 0;
198 
200  virtual void updateIntermediateVariables() = 0;
201 
203  virtual void correctStateResidualModelCon(List<List<word>>& stateCon) const = 0;
204 
206  virtual void addModelResidualCon(HashTable<List<List<word>>>& allCon) const = 0;
207 
209  virtual void calcResiduals(const dictionary& options) = 0;
210 
212  virtual void correct() = 0;
213 
215  virtual void getTurbProdTerm(scalarList& prodTerm) const;
216 
218  tmp<volSymmTensorField> devRhoReff() const;
219 
221  tmp<fvVectorMatrix> divDevRhoReff(volVectorField& U);
222 
224  tmp<fvVectorMatrix> divDevReff(volVectorField& U);
225 
227  tmp<volScalarField> nuEff() const;
228 
230  tmp<volScalarField> getNut()
231  {
232  return nut_;
233  }
234 
236  tmp<volScalarField> alphaEff();
237 
239  tmp<volScalarField> nu() const;
240 
242  tmp<volScalarField> getAlpha() const;
243 
245  tmp<volScalarField> getRho()
246  {
247  return rho_;
248  }
249 
251  tmp<volScalarField> getPhase()
252  {
253  return phase_;
254  }
255 
257  scalar getPrt()
258  {
259  if (Prt_ > 0)
260  {
261  return Prt_;
262  }
263  else
264  {
265  // alphat not found in Db so Prt is not set
266  FatalErrorIn("getPrt") << exit(FatalError);
267  return -9999.0;
268  }
269  }
270 
271 #ifdef CompressibleFlow
272  const fluidThermo& getThermo() const;
273 #endif
274 
276  tmp<Foam::volScalarField> getMu() const;
277 
279  bool writeData(Ostream& os) const;
280 
282  void printYPlus() const;
283 
284  label isPrintTime(
285  const Time& runTime,
286  const label printInterval) const;
287 
289  virtual void invTranProdNuTildaEqn(
290  const volScalarField& mySource,
291  volScalarField& pseudoNuTilda);
292 
293  virtual void constructPseudoNuTildaEqn();
294 
295  virtual void rhsSolvePseudoNuTildaEqn(const volScalarField& nuTildaSource);
296 
298  virtual void calcLduResidualTurb(volScalarField& nuTildaRes);
299 };
300 
301 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
302 
303 } // End namespace Foam
304 
305 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
306 
307 #endif
308 
309 // ************************************************************************* //
Foam::DARegDbTurbulenceModelCompressible
Definition: DARegDbTurbulenceModelCompressible.H:43
DAOption.H
Foam::DATurbulenceModel::nuEff
tmp< volScalarField > nuEff() const
return effective viscosity
Definition: DATurbulenceModel.C:223
U
U
Definition: pEqnPimpleDyM.H:60
DARegDbTurbulenceModelIncompressible.H
Foam::DATurbulenceModel::constructPseudoNuTildaEqn
virtual void constructPseudoNuTildaEqn()
Definition: DATurbulenceModel.C:521
Foam::DATurbulenceModel::phi_
surfaceScalarField & phi_
face flux
Definition: DATurbulenceModel.H:89
Foam::DATurbulenceModel::divDevRhoReff
tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U)
divDev terms
Definition: DATurbulenceModel.C:324
Foam::DATurbulenceModel::divDevReff
tmp< fvVectorMatrix > divDevReff(volVectorField &U)
divDev terms
Definition: DATurbulenceModel.C:347
Foam::DATurbulenceModel::kMin_
dimensionedScalar kMin_
Lower limit of k.
Definition: DATurbulenceModel.H:130
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
DARegDbTurbulenceModelCompressible.H
daOption
DAOption daOption(mesh, pyOptions_)
pseudoNuTilda
volScalarField & pseudoNuTilda
Definition: pseudoEqns.H:34
Foam::DATurbulenceModel::nuTildaMin_
dimensionedScalar nuTildaMin_
Lower limit for nuTilda.
Definition: DATurbulenceModel.H:139
Foam::DATurbulenceModel::getPrt
scalar getPrt()
get the turbulent Prandtl number
Definition: DATurbulenceModel.H:257
Foam::DATurbulenceModel::correctNut
virtual void correctNut()=0
update nut based on other turbulence variables and update the BCs
Foam::DATurbulenceModel::phaseRhoPhi_
surfaceScalarField & phaseRhoPhi_
phase*phi*density field
Definition: DATurbulenceModel.H:95
Foam::DATurbulenceModel::allOptions_
const dictionary & allOptions_
all DAFoam options
Definition: DATurbulenceModel.H:80
Foam::DATurbulenceModel::New
static autoPtr< DATurbulenceModel > New(const word modelType, const fvMesh &mesh, const DAOption &daOption)
Definition: DATurbulenceModel.C:164
DAMacroFunctions.H
Foam::DATurbulenceModel::Prt_
scalar Prt_
turbulent Prandtl number
Definition: DATurbulenceModel.H:145
Foam::DATurbulenceModel::getPhase
tmp< volScalarField > getPhase()
get the phase field
Definition: DATurbulenceModel.H:251
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:509
Foam::DATurbulenceModel::declareRunTimeSelectionTable
declareRunTimeSelectionTable(autoPtr, DATurbulenceModel, dictionary,(const word modelType, const fvMesh &mesh, const DAOption &daOption),(modelType, mesh, daOption))
Foam::DATurbulenceModel::getMu
tmp< Foam::volScalarField > getMu() const
get mu
Definition: DATurbulenceModel.C:289
Foam::DATurbulenceModel::correct
virtual void correct()=0
solve the residual equations and update the state
Foam::DATurbulenceModel::calcResiduals
virtual void calcResiduals(const dictionary &options)=0
compute the turbulence residuals
Foam::DATurbulenceModel::getAlpha
tmp< volScalarField > getAlpha() const
get alpha field
Definition: DATurbulenceModel.C:280
Foam::DATurbulenceModel::getNut
tmp< volScalarField > getNut()
get the nut field
Definition: DATurbulenceModel.H:230
Foam::DATurbulenceModel::correctWallDist
void correctWallDist()
update wall distance for d_. Note: y_ will be automatically updated in mesh_ object
Definition: DATurbulenceModel.C:358
Foam::DATurbulenceModel::U_
volVectorField & U_
velocity
Definition: DATurbulenceModel.H:86
Foam::DATurbulenceModel::turbDict_
IOdictionary turbDict_
turbulence model property dict
Definition: DATurbulenceModel.H:124
Foam::DATurbulenceModel::nut_
volScalarField & nut_
turbulence viscosity
Definition: DATurbulenceModel.H:83
Foam
Definition: multiFreqScalarFvPatchField.C:144
Foam::DATurbulenceModel::invTranProdNuTildaEqn
virtual void invTranProdNuTildaEqn(const volScalarField &mySource, volScalarField &pseudoNuTilda)
Inverse transpose product, M_nuTilda^(-T)
Definition: DATurbulenceModel.C:495
Foam::DATurbulenceModel::devRhoReff
tmp< volSymmTensorField > devRhoReff() const
dev terms
Definition: DATurbulenceModel.C:306
Foam::DARegDbTurbulenceModelIncompressible
Definition: DARegDbTurbulenceModelIncompressible.H:43
Foam::DATurbulenceModel::~DATurbulenceModel
virtual ~DATurbulenceModel()
Definition: DATurbulenceModel.H:178
Foam::DATurbulenceModel::getRho
tmp< volScalarField > getRho()
get the density field
Definition: DATurbulenceModel.H:245
Foam::DATurbulenceModel::rhsSolvePseudoNuTildaEqn
virtual void rhsSolvePseudoNuTildaEqn(const volScalarField &nuTildaSource)
Definition: DATurbulenceModel.C:533
Foam::DATurbulenceModel::nu
tmp< volScalarField > nu() const
get the nu field
Definition: DATurbulenceModel.C:262
Foam::DARegDbFluidThermo
Definition: DARegDbFluidThermo.H:40
Foam::DATurbulenceModel::addModelResidualCon
virtual void addModelResidualCon(HashTable< List< List< word >>> &allCon) const =0
add the model residual connectivity to stateCon
Foam::DATurbulenceModel::correctAlphat
void correctAlphat()
update alphat
Definition: DATurbulenceModel.C:206
DARegDbFluidThermo.H
Foam::DATurbulenceModel
Definition: DATurbulenceModel.H:61
Foam::DATurbulenceModel::daOption_
const DAOption & daOption_
DAOption object.
Definition: DATurbulenceModel.H:77
Foam::DATurbulenceModel::getTurbProdTerm
virtual void getTurbProdTerm(scalarList &prodTerm) const
return the value of the production term from the Spalart Allmaras model
Definition: DATurbulenceModel.C:483
Foam::DATurbulenceModel::mesh_
const fvMesh & mesh_
fvMesh
Definition: DATurbulenceModel.H:74
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::isPrintTime
label isPrintTime(const Time &runTime, const label printInterval) const
Definition: DATurbulenceModel.C:464
Foam::DATurbulenceModel::alphaEff
tmp< volScalarField > alphaEff()
return effective thermal diffusivity
Definition: DATurbulenceModel.C:236
Foam::DATurbulenceModel::omegaMin_
dimensionedScalar omegaMin_
Lower limit for omega.
Definition: DATurbulenceModel.H:136
Foam::DATurbulenceModel::printYPlus
void printYPlus() const
print the min max and mean yPlus to screen
Definition: DATurbulenceModel.C:385
Foam::DATurbulenceModel::coeffDict_
dictionary coeffDict_
turbulence model parameters dict
Definition: DATurbulenceModel.H:127
runTime
Time & runTime
Definition: createRefsHeatTransfer.H:1
DARegDbNearWallDist.H
Foam::DATurbulenceModel::Pr_
scalar Pr_
Prandtl number.
Definition: DATurbulenceModel.H:142
Foam::DARegDbSinglePhaseTransportModel
Definition: DARegDbSinglePhaseTransportModel.H:40
Foam::DATurbulenceModel::phase_
volScalarField phase_
phase field
Definition: DATurbulenceModel.H:92
DARegDbSinglePhaseTransportModel.H
Foam::DATurbulenceModel::epsilonMin_
dimensionedScalar epsilonMin_
Lower limit of epsilon.
Definition: DATurbulenceModel.H:133
Foam::DATurbulenceModel::writeData
bool writeData(Ostream &os) const
this is a virtual function for regIOobject
Definition: DATurbulenceModel.C:200