LESModel.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 | Copyright (C) 2013-2015 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 "LESModel.H"
27 
28 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
29 
30 template<class BasicTurbulenceModel>
32 {
33  if (printCoeffs_)
34  {
35  Info<< type << "Coeffs" << coeffDict_ << endl;
36  }
37 }
38 
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
42 template<class BasicTurbulenceModel>
44 (
45  const word& type,
46  const alphaField& alpha,
47  const rhoField& rho,
48  const volVectorField& U,
49  const surfaceScalarField& alphaRhoPhi,
50  const surfaceScalarField& phi,
51  const transportModel& transport,
52  const word& propertiesName
53 )
54 :
55  BasicTurbulenceModel
56  (
57  type,
58  alpha,
59  rho,
60  U,
61  alphaRhoPhi,
62  phi,
63  transport,
64  propertiesName
65  ),
66 
67  LESDict_(this->subOrEmptyDict("LES")),
68  turbulence_(LESDict_.lookup("turbulence")),
69  printCoeffs_(LESDict_.lookupOrDefault<Switch>("printCoeffs", false)),
70  coeffDict_(LESDict_.subOrEmptyDict(type + "Coeffs")),
71 
72  kMin_
73  (
75  (
76  "kMin",
77  LESDict_,
79  SMALL
80  )
81  ),
82 
83  delta_
84  (
86  (
87  IOobject::groupName("delta", U.group()),
88  *this,
89  LESDict_
90  )
91  )
92 {
93  // Force the construction of the mesh deltaCoeffs which may be needed
94  // for the construction of the derived models and BCs
95  this->mesh_.deltaCoeffs();
96 }
97 
98 
99 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
100 
101 template<class BasicTurbulenceModel>
104 (
105  const alphaField& alpha,
106  const rhoField& rho,
107  const volVectorField& U,
108  const surfaceScalarField& alphaRhoPhi,
109  const surfaceScalarField& phi,
110  const transportModel& transport,
111  const word& propertiesName
112 )
113 {
114  // get model name, but do not register the dictionary
115  // otherwise it is registered in the database twice
116  const word modelType
117  (
119  (
120  IOobject
121  (
122  IOobject::groupName(propertiesName, U.group()),
123  U.time().constant(),
124  U.db(),
125  IOobject::MUST_READ_IF_MODIFIED,
126  IOobject::NO_WRITE,
127  false
128  )
129  ).subDict("LES").lookup("LESModel")
130  );
131 
132  Info<< "Selecting LES turbulence model " << modelType << endl;
133 
134  typename dictionaryConstructorTable::iterator cstrIter =
135  dictionaryConstructorTablePtr_->find(modelType);
136 
137  if (cstrIter == dictionaryConstructorTablePtr_->end())
138  {
140  (
141  "LESModel::New"
142  "("
143  "const volScalarField&, "
144  "const volVectorField&, "
145  "const surfaceScalarField&, "
146  "transportModel&, "
147  "const word&"
148  ")"
149  ) << "Unknown LESModel type "
150  << modelType << nl << nl
151  << "Valid LESModel types:" << endl
152  << dictionaryConstructorTablePtr_->sortedToc()
153  << exit(FatalError);
154  }
155 
156  return autoPtr<LESModel>
157  (
158  cstrIter()(alpha, rho, U, alphaRhoPhi, phi, transport, propertiesName)
159  );
160 }
161 
162 
163 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
164 
165 template<class BasicTurbulenceModel>
167 {
169  {
170  LESDict_ <<= this->subDict("LES");
171  LESDict_.lookup("turbulence") >> turbulence_;
172 
173  if (const dictionary* dictPtr = LESDict_.subDictPtr(type() + "Coeffs"))
174  {
175  coeffDict_ <<= *dictPtr;
176  }
177 
178  delta_().read(LESDict_);
179 
180  kMin_.readIfPresent(LESDict_);
181 
182  return true;
183  }
184  else
185  {
186  return false;
187  }
188 }
189 
190 
191 template<class BasicTurbulenceModel>
193 {
194  delta_().correct();
196 }
197 
198 
199 // ************************************************************************* //
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/dm2/henry/OpenFOAM/OpenFOAM-3.0.x/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:74
static autoPtr< LESModel > 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 LES model.
Definition: LESModel.C:104
A class for handling words, derived from string.
Definition: word.H:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Templated abstract base class for LES SGS models.
Definition: LESModel.H:56
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
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: LESModel.C:31
messageStream Info
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:638
virtual bool read()
Read model coefficients if they have changed.
Definition: LESModel.C:166
const dictionary * subDictPtr(const word &) const
Find and return a sub-dictionary pointer if present.
Definition: dictionary.C:623
BasicTurbulenceModel::rhoField rhoField
Definition: LESModel.H:104
const word & constant() const
Return constant name.
Definition: TimePaths.H:124
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Generic dimensioned Type class.
static const char nl
Definition: Ostream.H:260
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
BasicTurbulenceModel::transportModel transportModel
Definition: LESModel.H:105
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: LESModel.C:192
BasicTurbulenceModel::alphaField alphaField
Definition: LESModel.H:103
word group() const
Return group (extension part of name)
Definition: IOobject.C:263
surfaceScalarField & phi
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:452
error FatalError
const Time & time() const
Return time.
Definition: IOobject.C:251
bool read(const char *, int32_t &)
Definition: int32IO.C:87
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:245
const dimensionSet dimVelocity
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
U
Definition: pEqn.H:82