laminarModel.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) 2016-2024 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 "laminarModel.H"
27 #include "Stokes.H"
28 
29 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
30 
31 template<class BasicMomentumTransportModel>
33 (
34  const word& type
35 )
36 {
37  if (printCoeffs_)
38  {
39  Info<< coeffDict_.dictName() << coeffDict_ << endl;
40  }
41 }
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
46 template<class BasicMomentumTransportModel>
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 viscosity& viscosity
56 )
57 :
58  BasicMomentumTransportModel
59  (
60  type,
61  alpha,
62  rho,
63  U,
64  alphaRhoPhi,
65  phi,
66  viscosity
67  ),
68 
69  laminarDict_(this->subOrEmptyDict("laminar")),
70  printCoeffs_(laminarDict_.lookupOrDefault<Switch>("printCoeffs", false)),
71  coeffDict_(laminarDict_.optionalSubDict(type + "Coeffs"))
72 {
73  // Force the construction of the mesh deltaCoeffs which may be needed
74  // for the construction of the derived models and BCs
75  this->mesh_.deltaCoeffs();
76 }
77 
78 
79 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
80 
81 template<class BasicMomentumTransportModel>
84 (
85  const alphaField& alpha,
86  const rhoField& rho,
87  const volVectorField& U,
88  const surfaceScalarField& alphaRhoPhi,
89  const surfaceScalarField& phi,
90  const viscosity& viscosity
91 )
92 {
93  const IOdictionary modelDict
94  (
96  (
97  U.db(),
98  alphaRhoPhi.group()
99  )
100  );
101 
102  if (modelDict.found("laminar"))
103  {
104  const word modelType =
105  modelDict.subDict("laminar").lookupBackwardsCompatible<word>
106  (
107  {"model", "laminarModel"}
108  );
109 
110  Info<< "Selecting laminar stress model " << modelType << endl;
111 
112  typename dictionaryConstructorTable::iterator cstrIter =
113  dictionaryConstructorTablePtr_->find(modelType);
114 
115  if (cstrIter == dictionaryConstructorTablePtr_->end())
116  {
118  << "Unknown laminarModel type "
119  << modelType << nl << nl
120  << "Valid laminarModel types:" << endl
121  << dictionaryConstructorTablePtr_->sortedToc()
122  << exit(FatalError);
123  }
124 
125  return autoPtr<laminarModel>
126  (
127  cstrIter()
128  (
129  alpha,
130  rho,
131  U,
132  alphaRhoPhi,
133  phi,
134  viscosity
135  )
136  );
137  }
138  else
139  {
140  Info<< "Selecting laminar stress model "
142  << endl;
143 
144  return autoPtr<laminarModel>
145  (
147  (
148  alpha,
149  rho,
150  U,
151  alphaRhoPhi,
152  phi,
153  viscosity
154  )
155  );
156  }
157 }
158 
159 
160 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
161 
162 template<class BasicMomentumTransportModel>
164 {
166  {
167  laminarDict_ <<= this->subDict("laminar");
168 
169  coeffDict_ <<= laminarDict_.optionalSubDict(type() + "Coeffs");
170 
171  return true;
172  }
173  else
174  {
175  return false;
176  }
177 }
178 
179 
180 template<class BasicMomentumTransportModel>
183 {
184  return volScalarField::New
185  (
186  this->groupName("nut"),
187  this->mesh_,
189  );
190 }
191 
192 
193 template<class BasicMomentumTransportModel>
196 (
197  const label patchi
198 ) const
199 {
200  return tmp<scalarField>
201  (
202  new scalarField(this->mesh_.boundary()[patchi].size(), 0.0)
203  );
204 }
205 
206 
207 template<class BasicMomentumTransportModel>
210 {
211  return volScalarField::New
212  (
213  this->groupName("k"),
214  this->mesh_,
216  );
217 }
218 
219 
220 template<class BasicMomentumTransportModel>
223 {
224  return volScalarField::New
225  (
226  this->groupName("epsilon"),
227  this->mesh_,
229  );
230 }
231 
232 
233 template<class BasicMomentumTransportModel>
236 {
237  return volScalarField::New
238  (
239  this->groupName("omega"),
240  this->mesh_,
242  );
243 }
244 
245 
246 template<class BasicMomentumTransportModel>
249 {
251  (
252  this->groupName("sigma"),
253  this->mesh_,
254  dimensionedSymmTensor(sqr(this->U_.dimensions()), Zero)
255  );
256 }
257 
258 
259 template<class BasicMomentumTransportModel>
261 {
262  BasicMomentumTransportModel::predict();
263 }
264 
265 
266 template<class BasicMomentumTransportModel>
268 {
270 }
271 
272 
273 // ************************************************************************* //
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
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:134
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
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:843
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:721
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:509
BasicMomentumTransportModel::alphaField alphaField
Definition: laminarModel.H:76
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate,.
Definition: laminarModel.C:222
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: laminarModel.C:33
static autoPtr< laminarModel > 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 laminar model.
Definition: laminarModel.C:84
virtual void correct()
Predict the laminar viscosity.
Definition: laminarModel.C:267
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy, i.e. 0 for laminar flow.
Definition: laminarModel.C:209
virtual void predict()
Predict the laminar viscosity.
Definition: laminarModel.C:260
laminarModel(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: laminarModel.C:48
virtual tmp< volSymmTensorField > sigma() const
Return the stress tensor [m^2/s^2], i.e. 0 for laminar flow.
Definition: laminarModel.C:248
virtual tmp< volScalarField > nut() const
Return the turbulence viscosity, i.e. 0 for laminar flow.
Definition: laminarModel.C:182
virtual tmp< volScalarField > omega() const
Return the turbulence specific dissipation rate,.
Definition: laminarModel.C:235
virtual bool read()
Read model coefficients if they have changed.
Definition: laminarModel.C:163
BasicMomentumTransportModel::rhoField rhoField
Definition: laminarModel.H:77
Momentum transport model for Stokes flow.
Definition: Stokes.H:55
static typeIOobject< IOdictionary > readModelDict(const objectRegistry &obr, const word &group, bool registerObject=false)
A class for managing temporary objects.
Definition: tmp.H:55
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
label patchi
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, const SuType &Su)
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
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:59
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
const dimensionSet dimKinematicViscosity
const dimensionSet dimless
messageStream Info
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const dimensionSet dimTime
dimensioned< symmTensor > dimensionedSymmTensor
Dimensioned tensor obtained from generic dimensioned type.
const dimensionSet dimVelocity
error FatalError
static const char nl
Definition: Ostream.H:266
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488