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 // * * * * * * * * * * * * * * * * 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  // Force the construction of the mesh deltaCoeffs which may be needed
55  // for the construction of the derived models and BCs
56  this->mesh_.deltaCoeffs();
57 }
58 
59 
60 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
61 
62 template<class BasicMomentumTransportModel>
65 (
66  const alphaField& alpha,
67  const rhoField& rho,
68  const volVectorField& U,
69  const surfaceScalarField& alphaRhoPhi,
70  const surfaceScalarField& phi,
71  const viscosity& viscosity
72 )
73 {
74  const IOdictionary modelDict
75  (
77  (
78  U.db(),
79  alphaRhoPhi.group()
80  )
81  );
82 
83  if (modelDict.found("laminar"))
84  {
85  const word modelType =
86  modelDict.subDict("laminar").lookupBackwardsCompatible<word>
87  (
88  {"model", "laminarModel"}
89  );
90 
91  Info<< indent
92  << "Selecting laminar stress model " << modelType << endl;
93 
94  typename dictionaryConstructorTable::iterator cstrIter =
95  dictionaryConstructorTablePtr_->find(modelType);
96 
97  if (cstrIter == dictionaryConstructorTablePtr_->end())
98  {
100  << "Unknown laminarModel type "
101  << modelType << nl << nl
102  << "Valid laminarModel types:" << endl
103  << dictionaryConstructorTablePtr_->sortedToc()
104  << exit(FatalError);
105  }
106 
107  Info<< incrIndent;
108 
109  autoPtr<laminarModel> modelPtr
110  (
111  cstrIter()
112  (
113  alpha,
114  rho,
115  U,
116  alphaRhoPhi,
117  phi,
118  viscosity
119  )
120  );
121 
122  Info<< decrIndent;
123 
124  return modelPtr;
125  }
126  else
127  {
128  Info<< indent
129  << "Selecting laminar stress model "
131  << endl;
132 
133  Info<< incrIndent;
134 
135  autoPtr<laminarModel> modelPtr
136  (
138  (
139  alpha,
140  rho,
141  U,
142  alphaRhoPhi,
143  phi,
144  viscosity
145  )
146  );
147 
148  Info<< decrIndent;
149 
150  return modelPtr;
151  }
152 }
153 
154 
155 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
156 
157 template<class BasicMomentumTransportModel>
158 const Foam::dictionary&
160 {
161  return this->subDict("laminar");
162 }
163 
164 
165 template<class BasicMomentumTransportModel>
166 const Foam::dictionary&
168 {
169  return this->laminarDict().optionalSubDict(type() + "Coeffs");
170 }
171 
172 
173 template<class BasicMomentumTransportModel>
175 {
177  {
178  return true;
179  }
180  else
181  {
182  return false;
183  }
184 }
185 
186 
187 template<class BasicMomentumTransportModel>
190 {
191  return volScalarField::New
192  (
193  this->groupName("nut"),
194  this->mesh_,
196  );
197 }
198 
199 
200 template<class BasicMomentumTransportModel>
203 (
204  const label patchi
205 ) const
206 {
207  return tmp<scalarField>
208  (
209  new scalarField(this->mesh_.boundary()[patchi].size(), 0.0)
210  );
211 }
212 
213 
214 template<class BasicMomentumTransportModel>
217 {
218  return volScalarField::New
219  (
220  this->groupName("k"),
221  this->mesh_,
223  );
224 }
225 
226 
227 template<class BasicMomentumTransportModel>
230 {
231  return volScalarField::New
232  (
233  this->groupName("epsilon"),
234  this->mesh_,
236  );
237 }
238 
239 
240 template<class BasicMomentumTransportModel>
243 {
244  return volScalarField::New
245  (
246  this->groupName("omega"),
247  this->mesh_,
249  );
250 }
251 
252 
253 template<class BasicMomentumTransportModel>
256 {
258  (
259  this->groupName("sigma"),
260  this->mesh_,
261  dimensionedSymmTensor(sqr(this->U_.dimensions()), Zero)
262  );
263 }
264 
265 
266 template<class BasicMomentumTransportModel>
268 {
269  BasicMomentumTransportModel::predict();
270 }
271 
272 
273 template<class BasicMomentumTransportModel>
275 {
277 }
278 
279 
280 // ************************************************************************* //
Generic GeometricField class.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, 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:131
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
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:539
BasicMomentumTransportModel::alphaField alphaField
Definition: laminarModel.H:67
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate,.
Definition: laminarModel.C:229
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor [m^2/s^2], i.e. 0 for laminar flow.
Definition: laminarModel.C:255
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:65
virtual void correct()
Predict the laminar viscosity.
Definition: laminarModel.C:274
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy, i.e. 0 for laminar flow.
Definition: laminarModel.C:216
const dictionary & laminarDict() const
Const access to the laminar dictionary.
Definition: laminarModel.C:159
virtual void predict()
Predict the laminar viscosity.
Definition: laminarModel.C:267
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:33
const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Definition: laminarModel.C:167
virtual tmp< volScalarField > nut() const
Return the turbulence viscosity, i.e. 0 for laminar flow.
Definition: laminarModel.C:189
virtual tmp< volScalarField > omega() const
Return the turbulence specific dissipation rate,.
Definition: laminarModel.C:242
virtual bool read()
Read model coefficients if they have changed.
Definition: laminarModel.C:174
BasicMomentumTransportModel::rhoField rhoField
Definition: laminarModel.H:68
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)
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
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:242
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
const dimensionSet dimKinematicViscosity
const dimensionSet dimless
messageStream Info
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:235
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const dimensionSet dimTime
dimensioned< symmTensor > dimensionedSymmTensor
Dimensioned tensor obtained from generic dimensioned type.
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
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