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-2018 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 BasicTurbulenceModel>
33 {
34  if (printCoeffs_)
35  {
36  Info<< coeffDict_.dictName() << coeffDict_ << endl;
37  }
38 }
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
43 template<class BasicTurbulenceModel>
45 (
46  const word& type,
47  const alphaField& alpha,
48  const rhoField& rho,
49  const volVectorField& U,
50  const surfaceScalarField& alphaRhoPhi,
51  const surfaceScalarField& phi,
52  const transportModel& transport,
53  const word& propertiesName
54 )
55 :
56  BasicTurbulenceModel
57  (
58  type,
59  alpha,
60  rho,
61  U,
62  alphaRhoPhi,
63  phi,
64  transport,
65  propertiesName
66  ),
67 
68  laminarDict_(this->subOrEmptyDict("laminar")),
69  printCoeffs_(laminarDict_.lookupOrDefault<Switch>("printCoeffs", false)),
70  coeffDict_(laminarDict_.optionalSubDict(type + "Coeffs"))
71 {
72  // Force the construction of the mesh deltaCoeffs which may be needed
73  // for the construction of the derived models and BCs
74  this->mesh_.deltaCoeffs();
75 }
76 
77 
78 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
79 
80 template<class BasicTurbulenceModel>
83 (
84  const alphaField& alpha,
85  const rhoField& rho,
86  const volVectorField& U,
87  const surfaceScalarField& alphaRhoPhi,
88  const surfaceScalarField& phi,
89  const transportModel& transport,
90  const word& propertiesName
91 )
92 {
93  IOdictionary modelDict
94  (
95  IOobject
96  (
97  IOobject::groupName(propertiesName, alphaRhoPhi.group()),
98  U.time().constant(),
99  U.db(),
100  IOobject::MUST_READ_IF_MODIFIED,
101  IOobject::NO_WRITE,
102  false
103  )
104  );
105 
106  if (modelDict.found("laminar"))
107  {
108  // get model name, but do not register the dictionary
109  // otherwise it is registered in the database twice
110  const word modelType
111  (
112  modelDict.subDict("laminar").lookup("laminarModel")
113  );
114 
115  Info<< "Selecting laminar stress model " << modelType << endl;
116 
117  typename dictionaryConstructorTable::iterator cstrIter =
118  dictionaryConstructorTablePtr_->find(modelType);
119 
120  if (cstrIter == dictionaryConstructorTablePtr_->end())
121  {
123  << "Unknown laminarModel type "
124  << modelType << nl << nl
125  << "Valid laminarModel types:" << endl
126  << dictionaryConstructorTablePtr_->sortedToc()
127  << exit(FatalError);
128  }
129 
130  return autoPtr<laminarModel>
131  (
132  cstrIter()
133  (
134  alpha,
135  rho,
136  U,
137  alphaRhoPhi,
138  phi,
139  transport, propertiesName)
140  );
141  }
142  else
143  {
144  Info<< "Selecting laminar stress model "
146 
147  return autoPtr<laminarModel>
148  (
150  (
151  alpha,
152  rho,
153  U,
154  alphaRhoPhi,
155  phi,
156  transport,
157  propertiesName
158  )
159  );
160  }
161 }
162 
163 
164 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
165 
166 template<class BasicTurbulenceModel>
168 {
170  {
171  laminarDict_ <<= this->subDict("laminar");
172 
173  coeffDict_ <<= laminarDict_.optionalSubDict(type() + "Coeffs");
174 
175  return true;
176  }
177  else
178  {
179  return false;
180  }
181 }
182 
183 
184 template<class BasicTurbulenceModel>
187 {
188  return tmp<volScalarField>
189  (
190  new volScalarField
191  (
192  IOobject
193  (
194  IOobject::groupName("nut", this->alphaRhoPhi_.group()),
195  this->runTime_.timeName(),
196  this->mesh_,
197  IOobject::NO_READ,
198  IOobject::NO_WRITE,
199  false
200  ),
201  this->mesh_,
202  dimensionedScalar("nut", dimViscosity, 0.0)
203  )
204  );
205 }
206 
207 
208 template<class BasicTurbulenceModel>
211 (
212  const label patchi
213 ) const
214 {
215  return tmp<scalarField>
216  (
217  new scalarField(this->mesh_.boundary()[patchi].size(), 0.0)
218  );
219 }
220 
221 
222 template<class BasicTurbulenceModel>
225 {
226  return tmp<volScalarField>
227  (
228  new volScalarField
229  (
230  IOobject::groupName("nuEff", this->alphaRhoPhi_.group()), this->nu()
231  )
232  );
233 }
234 
235 
236 template<class BasicTurbulenceModel>
239 (
240  const label patchi
241 ) const
242 {
243  return this->nu(patchi);
244 }
245 
246 
247 template<class BasicTurbulenceModel>
250 {
251  return tmp<volScalarField>
252  (
253  new volScalarField
254  (
255  IOobject
256  (
257  IOobject::groupName("k", this->alphaRhoPhi_.group()),
258  this->runTime_.timeName(),
259  this->mesh_,
260  IOobject::NO_READ,
261  IOobject::NO_WRITE,
262  false
263  ),
264  this->mesh_,
265  dimensionedScalar("k", sqr(this->U_.dimensions()), 0.0)
266  )
267  );
268 }
269 
270 
271 template<class BasicTurbulenceModel>
274 {
275  return tmp<volScalarField>
276  (
277  new volScalarField
278  (
279  IOobject
280  (
281  IOobject::groupName("epsilon", this->alphaRhoPhi_.group()),
282  this->runTime_.timeName(),
283  this->mesh_,
284  IOobject::NO_READ,
285  IOobject::NO_WRITE,
286  false
287  ),
288  this->mesh_,
290  (
291  "epsilon", sqr(this->U_.dimensions())/dimTime, 0.0
292  )
293  )
294  );
295 }
296 
297 
298 template<class BasicTurbulenceModel>
301 {
303  (
305  (
306  IOobject
307  (
308  IOobject::groupName("R", this->alphaRhoPhi_.group()),
309  this->runTime_.timeName(),
310  this->mesh_,
311  IOobject::NO_READ,
312  IOobject::NO_WRITE,
313  false
314  ),
315  this->mesh_,
317  (
318  "R", sqr(this->U_.dimensions()), Zero
319  )
320  )
321  );
322 }
323 
324 
325 template<class BasicTurbulenceModel>
327 {
329 }
330 
331 
332 // ************************************************************************* //
static word group(const word &name)
Return group (extension part of name)
Definition: IOobject.C:176
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:58
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:431
Templated abstract base class for laminar transport models.
Definition: laminarModel.H:49
type
Types of root.
Definition: Roots.H:52
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
surfaceScalarField & phi
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity, i.e. the laminar viscosity.
Definition: laminarModel.C:224
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimViscosity
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none.
Definition: Switch.H:60
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
virtual bool read()
Read model coefficients if they have changed.
Definition: laminarModel.C:167
BasicTurbulenceModel::transportModel transportModel
Definition: laminarModel.H:89
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:52
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:692
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:83
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
BasicTurbulenceModel::rhoField rhoField
Definition: laminarModel.H:88
BasicTurbulenceModel::alphaField alphaField
Definition: laminarModel.H:87
A class for handling words, derived from string.
Definition: word.H:59
const word & constant() const
Return constant name.
Definition: TimePaths.H:124
static const zero Zero
Definition: zero.H:97
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate,.
Definition: laminarModel.C:273
volScalarField scalarField(fieldObject, mesh)
virtual tmp< volScalarField > nut() const
Return the turbulence viscosity, i.e. 0 for laminar flow.
Definition: laminarModel.C:186
static const char nl
Definition: Ostream.H:265
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(psi *p+alphal *rhol0+((alphav *psiv+alphal *psil) - psi) *pSat, rhoMin);# 1 "/home/ubuntu/OpenFOAM-6/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:72
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy, i.e. 0 for laminar flow.
Definition: laminarModel.C:249
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: laminarModel.C:32
label patchi
U
Definition: pEqn.H:72
virtual void correct()
Correct the laminar transport.
Definition: laminarModel.C:326
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const Time & time() const
Return time.
Definition: IOobject.C:367
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor, i.e. 0 for laminar flow.
Definition: laminarModel.C:300
messageStream Info
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
Turbulence model for Stokes flow.
Definition: Stokes.H:52
A class for managing temporary objects.
Definition: PtrList.H:53
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:361
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
volScalarField & nu
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576