laminarModel.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2016-2017 OpenFOAM Foundation
9  Copyright (C) 2019-2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "laminarModel.H"
30 #include "Stokes.H"
31 
32 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
33 
34 template<class BasicTurbulenceModel>
36 {
37  if (printCoeffs_)
38  {
39  Info<< coeffDict_.dictName() << coeffDict_ << endl;
40  }
41 }
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
46 template<class BasicTurbulenceModel>
48 (
49  const word& type,
50  const alphaField& alpha,
51  const rhoField& rho,
52  const volVectorField& U,
53  const surfaceScalarField& alphaRhoPhi,
54  const surfaceScalarField& phi,
55  const transportModel& transport,
56  const word& propertiesName
57 )
58 :
59  BasicTurbulenceModel
60  (
61  type,
62  alpha,
63  rho,
64  U,
65  alphaRhoPhi,
66  phi,
67  transport,
68  propertiesName
69  ),
70 
71  laminarDict_(this->subOrEmptyDict("laminar")),
72  printCoeffs_(laminarDict_.getOrDefault<Switch>("printCoeffs", false)),
73  coeffDict_(laminarDict_.optionalSubDict(type + "Coeffs"))
74 {
75  // Force the construction of the mesh deltaCoeffs which may be needed
76  // for the construction of the derived models and BCs
77  this->mesh_.deltaCoeffs();
78 }
79 
80 
81 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
82 
83 template<class BasicTurbulenceModel>
86 (
87  const alphaField& alpha,
88  const rhoField& rho,
89  const volVectorField& U,
90  const surfaceScalarField& alphaRhoPhi,
91  const surfaceScalarField& phi,
92  const transportModel& transport,
93  const word& propertiesName
94 )
95 {
96  const IOdictionary modelDict
97  (
98  IOobject
99  (
100  IOobject::groupName(propertiesName, alphaRhoPhi.group()),
101  U.time().constant(),
102  U.db(),
103  IOobject::MUST_READ_IF_MODIFIED,
104  IOobject::NO_WRITE,
105  false // Do not register
106  )
107  );
108 
109  const dictionary* dictptr = modelDict.findDict("laminar");
110 
111  if (dictptr)
112  {
113  const dictionary& dict = *dictptr;
114 
115  const word modelType
116  (
117  // laminarModel -> model (after v2006)
118  dict.getCompat<word>("model", {{"laminarModel", -2006}})
119  );
120 
121  Info<< "Selecting laminar stress model " << modelType << endl;
122 
123  auto* ctorPtr = dictionaryConstructorTable(modelType);
124 
125  if (!ctorPtr)
126  {
128  (
129  dict,
130  "laminar model",
131  modelType,
132  *dictionaryConstructorTablePtr_
133  ) << exit(FatalIOError);
134  }
135 
136  return autoPtr<laminarModel>
137  (
138  ctorPtr
139  (
140  alpha,
141  rho,
142  U,
143  alphaRhoPhi,
144  phi,
145  transport, propertiesName)
146  );
147  }
148  else
149  {
150  Info<< "Selecting laminar stress model "
152 
153  return autoPtr<laminarModel>
154  (
156  (
157  alpha,
158  rho,
159  U,
160  alphaRhoPhi,
161  phi,
162  transport,
163  propertiesName
164  )
165  );
166  }
167 }
168 
169 
170 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
171 
172 template<class BasicTurbulenceModel>
174 {
176  {
177  laminarDict_ <<= this->subDict("laminar");
178 
179  coeffDict_ <<= laminarDict_.optionalSubDict(type() + "Coeffs");
180 
181  return true;
182  }
183 
184  return false;
185 }
186 
187 
188 template<class BasicTurbulenceModel>
191 {
193  (
194  IOobject
195  (
196  IOobject::groupName("nut", this->alphaRhoPhi_.group()),
197  this->runTime_.timeName(),
198  this->mesh_,
199  IOobject::NO_READ,
200  IOobject::NO_WRITE,
201  false
202  ),
203  this->mesh_,
205  );
206 }
207 
208 
209 template<class BasicTurbulenceModel>
212 (
213  const label patchi
214 ) const
215 {
216  return tmp<scalarField>
217  (
218  new scalarField(this->mesh_.boundary()[patchi].size(), Zero)
219  );
220 }
221 
222 
223 template<class BasicTurbulenceModel>
226 {
227  return tmp<volScalarField>
228  (
229  new volScalarField
230  (
231  IOobject::groupName("nuEff", this->alphaRhoPhi_.group()), this->nu()
232  )
233  );
234 }
235 
236 
237 template<class BasicTurbulenceModel>
240 (
241  const label patchi
242 ) const
243 {
244  return this->nu(patchi);
245 }
246 
247 
248 template<class BasicTurbulenceModel>
251 {
253  (
254  IOobject
255  (
256  IOobject::groupName("k", this->alphaRhoPhi_.group()),
257  this->runTime_.timeName(),
258  this->mesh_,
259  IOobject::NO_READ,
260  IOobject::NO_WRITE,
261  false
262  ),
263  this->mesh_,
264  dimensionedScalar(sqr(this->U_.dimensions()), Zero)
265  );
266 }
267 
268 
269 template<class BasicTurbulenceModel>
272 {
274  (
275  IOobject
276  (
277  IOobject::groupName("epsilon", this->alphaRhoPhi_.group()),
278  this->runTime_.timeName(),
279  this->mesh_,
280  IOobject::NO_READ,
281  IOobject::NO_WRITE,
282  false
283  ),
284  this->mesh_,
285  dimensionedScalar(sqr(this->U_.dimensions())/dimTime, Zero)
286  );
287 }
288 
289 
290 template<class BasicTurbulenceModel>
293 {
295  (
296  IOobject
297  (
298  IOobject::groupName("omega", this->alphaRhoPhi_.group()),
299  this->runTime_.timeName(),
300  this->mesh_,
301  IOobject::NO_READ,
302  IOobject::NO_WRITE,
303  false
304  ),
305  this->mesh_,
307  );
308 }
309 
310 
311 template<class BasicTurbulenceModel>
314 {
316  (
317  IOobject
318  (
319  IOobject::groupName("R", this->alphaRhoPhi_.group()),
320  this->runTime_.timeName(),
321  this->mesh_,
322  IOobject::NO_READ,
323  IOobject::NO_WRITE,
324  false
325  ),
326  this->mesh_,
327  dimensionedSymmTensor(sqr(this->U_.dimensions()), Zero)
328  );
329 }
330 
331 
332 template<class BasicTurbulenceModel>
334 {
336 }
337 
338 
339 // ************************************************************************* //
Foam::dictionary::findDict
dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX)
Find and return a sub-dictionary pointer if present.
Definition: dictionaryI.H:127
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::Switch
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:77
Foam::laminarModel< BasicMomentumTransportModel >::alphaField
BasicMomentumTransportModel ::alphaField alphaField
Definition: laminarModel.H:85
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::dimensionedSymmTensor
dimensioned< symmTensor > dimensionedSymmTensor
Dimensioned tensor obtained from generic dimensioned type.
Definition: dimensionedSymmTensor.H:50
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::laminarModel::laminarModel
laminarModel(const laminarModel &)=delete
No copy construct.
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::laminarModels::Stokes
Turbulence model for Stokes flow.
Definition: Stokes.H:58
rho
rho
Definition: readInitialConditions.H:88
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
Foam::blockMeshTools::read
void read(Istream &, label &val, const dictionary &)
In-place read with dictionary lookup.
Definition: blockMeshTools.C:57
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
FatalIOErrorInLookup
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:478
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
correct
fvOptions correct(rho)
Foam::laminarModel::read
virtual bool read()
Read model coefficients if they have changed.
Definition: laminarModel.C:173
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::laminarModel::printCoeffs
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: laminarModel.C:35
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dimViscosity
const dimensionSet dimViscosity
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::laminarModel< BasicMomentumTransportModel >::transportModel
BasicMomentumTransportModel ::transportModel transportModel
Definition: laminarModel.H:87
Foam::laminarModel::R
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor, i.e. 0 for laminar flow.
Definition: laminarModel.C:313
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Stokes.H
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
U
U
Definition: pEqn.H:72
Foam::laminarModel::correct
virtual void correct()
Correct the laminar transport.
Definition: laminarModel.C:333
Foam::laminarModel::k
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy, i.e. 0 for laminar flow.
Definition: laminarModel.C:250
laminarModel.H
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::laminarModel::New
static autoPtr< laminarModel > New(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName)
Return a reference to the selected laminar model.
Definition: laminarModel.C:86
Foam::laminarModel::nuEff
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity, i.e. the laminar viscosity.
Definition: laminarModel.C:225
Foam::roots::type
type
Types of root.
Definition: Roots.H:54
Foam::laminarModel::epsilon
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate,.
Definition: laminarModel.C:271
Foam::GeometricField< vector, fvPatchField, volMesh >
Foam::laminarModel< BasicMomentumTransportModel >::rhoField
BasicMomentumTransportModel ::rhoField rhoField
Definition: laminarModel.H:86
Foam::laminarModel::nut
virtual tmp< volScalarField > nut() const
Return the turbulence viscosity, i.e. 0 for laminar flow.
Definition: laminarModel.C:190
Foam::laminarModel::omega
virtual tmp< volScalarField > omega() const
Return the specific dissipation rate, i.e. 0 for laminar flow.
Definition: laminarModel.C:292
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189