All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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-2020 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 transportModel& transport
56 )
57 :
58  BasicMomentumTransportModel
59  (
60  type,
61  alpha,
62  rho,
63  U,
64  alphaRhoPhi,
65  phi,
66  transport
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 transportModel& transport
91 )
92 {
93  const IOdictionary modelDict
94  (
95  momentumTransportModel::readModelDict
96  (
97  U.db(),
98  alphaRhoPhi.group()
99  )
100  );
101 
102  if (modelDict.found("laminar"))
103  {
104  const word modelType
105  (
106  modelDict.subDict("laminar").found("model")
107  ? modelDict.subDict("laminar").lookup("model")
108  : modelDict.subDict("laminar").lookup("laminarModel")
109  );
110 
111  Info<< "Selecting laminar stress model " << modelType << endl;
112 
113  typename dictionaryConstructorTable::iterator cstrIter =
114  dictionaryConstructorTablePtr_->find(modelType);
115 
116  if (cstrIter == dictionaryConstructorTablePtr_->end())
117  {
119  << "Unknown laminarModel type "
120  << modelType << nl << nl
121  << "Valid laminarModel types:" << endl
122  << dictionaryConstructorTablePtr_->sortedToc()
123  << exit(FatalError);
124  }
125 
126  return autoPtr<laminarModel>
127  (
128  cstrIter()
129  (
130  alpha,
131  rho,
132  U,
133  alphaRhoPhi,
134  phi,
135  transport
136  )
137  );
138  }
139  else
140  {
141  Info<< "Selecting laminar stress model "
143  << endl;
144 
145  return autoPtr<laminarModel>
146  (
148  (
149  alpha,
150  rho,
151  U,
152  alphaRhoPhi,
153  phi,
154  transport
155  )
156  );
157  }
158 }
159 
160 
161 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
162 
163 template<class BasicMomentumTransportModel>
165 {
167  {
168  laminarDict_ <<= this->subDict("laminar");
169 
170  coeffDict_ <<= laminarDict_.optionalSubDict(type() + "Coeffs");
171 
172  return true;
173  }
174  else
175  {
176  return false;
177  }
178 }
179 
180 
181 template<class BasicMomentumTransportModel>
184 {
185  return volScalarField::New
186  (
187  IOobject::groupName("nut", this->alphaRhoPhi_.group()),
188  this->mesh_,
190  );
191 }
192 
193 
194 template<class BasicMomentumTransportModel>
197 (
198  const label patchi
199 ) const
200 {
201  return tmp<scalarField>
202  (
203  new scalarField(this->mesh_.boundary()[patchi].size(), 0.0)
204  );
205 }
206 
207 
208 template<class BasicMomentumTransportModel>
211 {
212  return volScalarField::New
213  (
214  IOobject::groupName("nuEff", this->alphaRhoPhi_.group()),
215  this->nu()
216  );
217 }
218 
219 
220 template<class BasicMomentumTransportModel>
223 (
224  const label patchi
225 ) const
226 {
227  return this->nu(patchi);
228 }
229 
230 
231 template<class BasicMomentumTransportModel>
234 {
235  return volScalarField::New
236  (
237  IOobject::groupName("k", this->alphaRhoPhi_.group()),
238  this->mesh_,
239  dimensionedScalar(sqr(this->U_.dimensions()), 0)
240  );
241 }
242 
243 
244 template<class BasicMomentumTransportModel>
247 {
248  return volScalarField::New
249  (
250  IOobject::groupName("epsilon", this->alphaRhoPhi_.group()),
251  this->mesh_,
252  dimensionedScalar(sqr(this->U_.dimensions())/dimTime, 0)
253  );
254 }
255 
256 
257 template<class BasicMomentumTransportModel>
260 {
262  (
263  IOobject::groupName("R", this->alphaRhoPhi_.group()),
264  this->mesh_,
265  dimensionedSymmTensor(sqr(this->U_.dimensions()), Zero)
266  );
267 }
268 
269 
270 template<class BasicMomentumTransportModel>
272 {
274 }
275 
276 
277 // ************************************************************************* //
static word group(const word &name)
Return group (extension part of name)
Definition: IOobject.C:144
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:667
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
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:210
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimViscosity
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
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/any.
Definition: Switch.H:60
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
BasicMomentumTransportModel::rhoField rhoField
Definition: laminarModel.H:77
BasicMomentumTransportModel::alphaField alphaField
Definition: laminarModel.H:76
BasicMomentumTransportModel::transportModel transportModel
Definition: laminarModel.H:78
virtual bool read()
Read model coefficients if they have changed.
Definition: laminarModel.C:164
phi
Definition: pEqn.H:104
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:934
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
A class for handling words, derived from string.
Definition: word.H:59
static const zero Zero
Definition: zero.H:97
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate,.
Definition: laminarModel.C:246
volScalarField scalarField(fieldObject, mesh)
virtual tmp< volScalarField > nut() const
Return the turbulence viscosity, i.e. 0 for laminar flow.
Definition: laminarModel.C:183
static const char nl
Definition: Ostream.H:260
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-8/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:68
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy, i.e. 0 for laminar flow.
Definition: laminarModel.C:233
dimensioned< symmTensor > dimensionedSymmTensor
Dimensioned tensor obtained from generic dimensioned type.
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: laminarModel.C:33
virtual tmp< volSymmTensorField > sigma() const
Return the stress tensor [m^2/s^2], i.e. 0 for laminar flow.
Definition: laminarModel.C:259
laminarModel(const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)
Construct from components.
Definition: laminarModel.C:48
label patchi
U
Definition: pEqn.H:72
virtual void correct()
Correct the laminar transport.
Definition: laminarModel.C:271
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
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
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
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
A class for managing temporary objects.
Definition: PtrList.H:53
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:322
static autoPtr< laminarModel > New(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)
Return a reference to the selected laminar model.
Definition: laminarModel.C:84
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:812