Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
11License
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\*---------------------------------------------------------------------------*/
29#include "laminarModel.H"
30#include "Stokes.H"
31
32// * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
33
34template<class BasicTurbulenceModel>
36{
37 if (printCoeffs_)
38 {
39 Info<< coeffDict_.dictName() << coeffDict_ << endl;
40 }
41}
42
43
44// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45
46template<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();
79
80
81// * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
82
83template<class BasicTurbulenceModel>
86(
87 const alphaField& alpha,
88 const rhoField& rho,
89 const volVectorField& U,
90 const surfaceScalarField& alphaRhoPhi,
92 const transportModel& transport,
93 const word& propertiesName
94)
95{
96 const IOdictionary modelDict
97 (
99 (
100 IOobject::groupName(propertiesName, alphaRhoPhi.group()),
101 U.time().constant(),
102 U.db(),
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
137 (
138 ctorPtr
139 (
140 alpha,
141 rho,
142 U,
143 alphaRhoPhi,
144 phi,
145 transport, propertiesName)
146 );
148 else
149 {
150 Info<< "Selecting laminar stress model "
152
154 (
156 (
157 alpha,
158 rho,
159 U,
160 alphaRhoPhi,
161 phi,
162 transport,
163 propertiesName
164 )
165 );
167}
168
169
170// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
171
172template<class BasicTurbulenceModel>
174{
175 if (BasicTurbulenceModel::read())
176 {
177 laminarDict_ <<= this->subDict("laminar");
178
179 coeffDict_ <<= laminarDict_.optionalSubDict(type() + "Coeffs");
180
181 return true;
182 }
184 return false;
185}
187
188template<class BasicTurbulenceModel>
191{
193 (
194 IOobject::groupName("nut", this->alphaRhoPhi_.group()),
196 this->mesh_,
198 );
199}
200
201
202template<class BasicTurbulenceModel>
205(
206 const label patchi
207) const
209 return tmp<scalarField>::New(this->mesh_.boundary()[patchi].size(), Zero);
210}
211
212
213template<class BasicTurbulenceModel>
216{
219 IOobject::groupName("nuEff", this->alphaRhoPhi_.group()),
221 this->nu()
222 );
224
225
226template<class BasicTurbulenceModel>
229(
230 const label patchi
231) const
233 return this->nu(patchi);
234}
235
236
237template<class BasicTurbulenceModel>
240{
242 (
243 IOobject::groupName("k", this->alphaRhoPhi_.group()),
245 this->mesh_,
246 dimensionedScalar(sqr(this->U_.dimensions()), Zero)
247 );
248}
249
250
251template<class BasicTurbulenceModel>
254{
256 (
257 IOobject::groupName("epsilon", this->alphaRhoPhi_.group()),
259 this->mesh_,
260 dimensionedScalar(sqr(this->U_.dimensions())/dimTime, Zero)
261 );
262}
263
264
265template<class BasicTurbulenceModel>
268{
270 (
271 IOobject::groupName("omega", this->alphaRhoPhi_.group()),
273 this->mesh_,
275 );
276}
277
278
279template<class BasicTurbulenceModel>
282{
284 (
285 IOobject::groupName("R", this->alphaRhoPhi_.group()),
287 this->mesh_,
288 dimensionedSymmTensor(sqr(this->U_.dimensions()), Zero)
289 );
290}
291
292
293template<class BasicTurbulenceModel>
295{
296 BasicTurbulenceModel::correct();
297}
298
299
300// ************************************************************************* //
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
@ NO_REGISTER
Do not request registration (bool: false).
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static word group(const word &name)
Return group (extension part of name).
Definition IOobject.C:300
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition Switch.H:81
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary pointer if present (and it is a dictionary) otherwise return nullptr...
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate,.
BasicTurbulenceModel::alphaField alphaField
dictionary coeffDict_
Model coefficients dictionary.
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor, i.e. 0 for laminar flow.
virtual void printCoeffs(const word &type)
Print model coefficients.
BasicTurbulenceModel::rhoField rhoField
virtual void correct()
Correct the laminar transport.
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy, i.e. 0 for laminar flow.
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.
dictionary laminarDict_
laminar coefficients dictionary
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity, i.e. the laminar viscosity.
laminarModel(const laminarModel &)=delete
No copy construct.
BasicTurbulenceModel::transportModel transportModel
virtual tmp< volScalarField > nut() const
Return the turbulence viscosity, i.e. 0 for laminar flow.
virtual tmp< volScalarField > omega() const
Return the specific dissipation rate, i.e. 0 for laminar flow.
Switch printCoeffs_
Flag to print the model coeffs at run-time.
virtual bool read()
Read model coefficients if they have changed.
Turbulence model for Stokes flow.
Definition Stokes.H:57
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
Base-class for all transport models used by the incompressible turbulence models.
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition error.H:637
const dimensionSet dimViscosity
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
dimensioned< symmTensor > dimensionedSymmTensor
Dimensioned tensor obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
volScalarField & nu
volScalarField & alpha
dictionary dict