LESModel.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2013-2025 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "LESModel.H"
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
31 template<class BasicMomentumTransportModel>
33 (
34  const word& type,
35  const alphaField& alpha,
36  const rhoField& rho,
37  const volVectorField& U,
38  const surfaceScalarField& alphaRhoPhi,
39  const surfaceScalarField& phi,
40  const viscosity& viscosity
41 )
42 :
43  BasicMomentumTransportModel
44  (
45  type,
46  alpha,
47  rho,
48  U,
49  alphaRhoPhi,
50  phi,
51  viscosity
52  ),
53 
54  turbulence_(LESDict().template lookupOrDefault<Switch>("turbulence", true)),
55  kMin_("kMin", sqr(dimVelocity), LESDict(), small),
56  nutMaxCoeff_("nutMaxCoeff", dimless, LESDict(), 1e5),
57 
58  viscosityModel_
59  (
60  coeffDict().found("viscosityModel")
61  ? laminarModels::generalisedNewtonianViscosityModel::New
62  (
63  coeffDict(),
64  viscosity,
65  U
66  )
67  : autoPtr<laminarModels::generalisedNewtonianViscosityModel>
68  (
69  new laminarModels::generalisedNewtonianViscosityModels::Newtonian
70  (
71  coeffDict(),
72  viscosity,
73  U
74  )
75  )
76  ),
77 
78  delta_
79  (
80  LESdelta::New
81  (
82  this->groupName("delta"),
83  *this,
84  LESDict()
85  )
86  )
87 {
88  // Force the construction of the mesh deltaCoeffs which may be needed
89  // for the construction of the derived models and BCs
90  this->mesh_.deltaCoeffs();
91 }
92 
93 
94 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
95 
96 template<class BasicMomentumTransportModel>
99 (
100  const alphaField& alpha,
101  const rhoField& rho,
102  const volVectorField& U,
103  const surfaceScalarField& alphaRhoPhi,
104  const surfaceScalarField& phi,
105  const viscosity& viscosity
106 )
107 {
108  const IOdictionary modelDict
109  (
111  (
112  U.db(),
113  alphaRhoPhi.group()
114  )
115  );
116 
117  const word modelType =
118  modelDict.subDict("LES").lookupBackwardsCompatible<word>
119  (
120  {"model", "LESModel"}
121  );
122 
123  Info<< indent
124  << "Selecting LES turbulence model " << modelType << endl;
125 
126  typename dictionaryConstructorTable::iterator cstrIter =
127  dictionaryConstructorTablePtr_->find(modelType);
128 
129  if (cstrIter == dictionaryConstructorTablePtr_->end())
130  {
132  << "Unknown LESModel type "
133  << modelType << nl << nl
134  << "Valid LESModel types:" << endl
135  << dictionaryConstructorTablePtr_->sortedToc()
136  << exit(FatalError);
137  }
138 
139  Info<< incrIndent;
140 
141  autoPtr<LESModel> modelPtr
142  (
143  cstrIter()(alpha, rho, U, alphaRhoPhi, phi, viscosity)
144  );
145 
146  Info<< decrIndent;
147 
148  return modelPtr;
149 }
150 
151 
152 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
153 
154 template<class BasicMomentumTransportModel>
155 const Foam::dictionary&
157 {
158  return this->subDict("LES");
159 }
160 
161 
162 template<class BasicMomentumTransportModel>
163 const Foam::dictionary&
165 {
166  return this->LESDict().optionalSubDict(type() + "Coeffs");
167 }
168 
169 
170 template<class BasicMomentumTransportModel>
172 {
174  {
175  turbulence_.readIfPresent("turbulence", LESDict());
176  delta_().read(LESDict());
177 
178  kMin_.readIfPresent(LESDict());
179  nutMaxCoeff_.readIfPresent(LESDict());
180 
181  return true;
182  }
183  else
184  {
185  return false;
186  }
187 }
188 
189 
190 template<class BasicMomentumTransportModel>
192 {
193  viscosityModel_->correct();
194  delta_().correct();
196 }
197 
198 
199 // ************************************************************************* //
bool found
Generic GeometricField class.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
static word group(const word &name)
Return group (extension part of name)
Definition: IOobject.C:131
BasicMomentumTransportModel::alphaField alphaField
Definition: LESModel.H:97
const dictionary & LESDict() const
Const access to the LES dictionary.
Definition: LESModel.C:156
LESModel(const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity)
Construct from components.
Definition: LESModel.C:33
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: LESModel.C:191
const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Definition: LESModel.C:164
static autoPtr< LESModel > New(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity)
Return a reference to the selected LES model.
Definition: LESModel.C:99
virtual bool read()
Read model coefficients if they have changed.
Definition: LESModel.C:171
BasicMomentumTransportModel::rhoField rhoField
Definition: LESModel.H:98
Abstract base class for LES deltas.
Definition: LESdelta.H:51
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
Definition: dictionary.C:927
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:849
ITstream & lookupBackwardsCompatible(const wordList &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream, trying a list of keywords.
Definition: dictionary.C:751
static typeIOobject< IOdictionary > readModelDict(const objectRegistry &obr, const word &group, bool registerObject=false)
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiCorr, const SpType &Sp)
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:242
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
const dimensionSet dimless
messageStream Info
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:235
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
const dimensionSet dimVelocity
error FatalError
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:228
static const char nl
Definition: Ostream.H:267
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488