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 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_.lookupOrDefault<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(dict.get<word>("laminarModel"));
116 
117  Info<< "Selecting laminar stress model " << modelType << endl;
118 
119  auto cstrIter = dictionaryConstructorTablePtr_->cfind(modelType);
120 
121  if (!cstrIter.found())
122  {
124  (
125  dict,
126  "laminarModel",
127  modelType,
128  *dictionaryConstructorTablePtr_
129  ) << exit(FatalIOError);
130  }
131 
132  return autoPtr<laminarModel>
133  (
134  cstrIter()
135  (
136  alpha,
137  rho,
138  U,
139  alphaRhoPhi,
140  phi,
141  transport, propertiesName)
142  );
143  }
144  else
145  {
146  Info<< "Selecting laminar stress model "
148 
149  return autoPtr<laminarModel>
150  (
152  (
153  alpha,
154  rho,
155  U,
156  alphaRhoPhi,
157  phi,
158  transport,
159  propertiesName
160  )
161  );
162  }
163 }
164 
165 
166 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
167 
168 template<class BasicTurbulenceModel>
170 {
172  {
173  laminarDict_ <<= this->subDict("laminar");
174 
175  coeffDict_ <<= laminarDict_.optionalSubDict(type() + "Coeffs");
176 
177  return true;
178  }
179 
180  return false;
181 }
182 
183 
184 template<class BasicTurbulenceModel>
187 {
189  (
190  IOobject
191  (
192  IOobject::groupName("nut", this->alphaRhoPhi_.group()),
193  this->runTime_.timeName(),
194  this->mesh_,
195  IOobject::NO_READ,
196  IOobject::NO_WRITE,
197  false
198  ),
199  this->mesh_,
201  );
202 }
203 
204 
205 template<class BasicTurbulenceModel>
208 (
209  const label patchi
210 ) const
211 {
212  return tmp<scalarField>
213  (
214  new scalarField(this->mesh_.boundary()[patchi].size(), Zero)
215  );
216 }
217 
218 
219 template<class BasicTurbulenceModel>
222 {
223  return tmp<volScalarField>
224  (
225  new volScalarField
226  (
227  IOobject::groupName("nuEff", this->alphaRhoPhi_.group()), this->nu()
228  )
229  );
230 }
231 
232 
233 template<class BasicTurbulenceModel>
236 (
237  const label patchi
238 ) const
239 {
240  return this->nu(patchi);
241 }
242 
243 
244 template<class BasicTurbulenceModel>
247 {
249  (
250  IOobject
251  (
252  IOobject::groupName("k", this->alphaRhoPhi_.group()),
253  this->runTime_.timeName(),
254  this->mesh_,
255  IOobject::NO_READ,
256  IOobject::NO_WRITE,
257  false
258  ),
259  this->mesh_,
260  dimensionedScalar(sqr(this->U_.dimensions()), Zero)
261  );
262 }
263 
264 
265 template<class BasicTurbulenceModel>
268 {
270  (
271  IOobject
272  (
273  IOobject::groupName("epsilon", this->alphaRhoPhi_.group()),
274  this->runTime_.timeName(),
275  this->mesh_,
276  IOobject::NO_READ,
277  IOobject::NO_WRITE,
278  false
279  ),
280  this->mesh_,
281  dimensionedScalar(sqr(this->U_.dimensions())/dimTime, Zero)
282  );
283 }
284 
285 
286 template<class BasicTurbulenceModel>
289 {
291  (
292  IOobject
293  (
294  IOobject::groupName("R", this->alphaRhoPhi_.group()),
295  this->runTime_.timeName(),
296  this->mesh_,
297  IOobject::NO_READ,
298  IOobject::NO_WRITE,
299  false
300  ),
301  this->mesh_,
302  dimensionedSymmTensor(sqr(this->U_.dimensions()), Zero)
303  );
304 }
305 
306 
307 template<class BasicTurbulenceModel>
309 {
311 }
312 
313 
314 // ************************************************************************* //
Foam::dictionary::findDict
dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX)
Find and return a sub-dictionary pointer if present.
Definition: dictionary.C:503
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:104
Foam::Switch
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:70
Foam::laminarModel::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: laminarModel.H:84
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
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:59
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
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:337
Foam::laminarModels::Stokes
Turbulence model for Stokes flow.
Definition: Stokes.H:58
rho
rho
Definition: readInitialConditions.H:96
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:54
FatalIOErrorInLookup
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:380
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
correct
fvOptions correct(rho)
Foam::laminarModel::read
virtual bool read()
Read model coefficients if they have changed.
Definition: laminarModel.C:169
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::laminarModel::printCoeffs
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: laminarModel.C:35
Foam::blockMeshTools::read
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
Definition: blockMeshTools.C:33
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:121
Foam::laminarModel::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: laminarModel.H:86
Foam::laminarModel::R
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor, i.e. 0 for laminar flow.
Definition: laminarModel.C:288
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:308
Foam::laminarModel::k
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy, i.e. 0 for laminar flow.
Definition: laminarModel.C:246
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:221
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:267
Foam::GeometricField< vector, fvPatchField, volMesh >
Foam::laminarModel::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: laminarModel.H:85
Foam::laminarModel::nut
virtual tmp< volScalarField > nut() const
Return the turbulence viscosity, i.e. 0 for laminar flow.
Definition: laminarModel.C:186