phaseSystem.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) 2015-2016 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 "phaseSystem.H"
27 #include "surfaceTensionModel.H"
28 #include "aspectRatioModel.H"
29 #include "surfaceInterpolate.H"
30 #include "fvcDdt.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(phaseSystem, 0);
37 }
38 
39 const Foam::word Foam::phaseSystem::propertiesName("phaseProperties");
40 
41 
42 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
43 
45 (
46  const phaseModelList& phaseModels
47 ) const
48 {
49  tmp<surfaceScalarField> tmpPhi
50  (
52  (
53  "phi",
54  fvc::interpolate(phaseModels[0])*phaseModels[0].phi()
55  )
56  );
57 
58  for (label phasei=1; phasei<phaseModels.size(); phasei++)
59  {
60  tmpPhi.ref() +=
61  fvc::interpolate(phaseModels[phasei])*phaseModels[phasei].phi();
62  }
63 
64  return tmpPhi;
65 }
66 
67 
69 (
70  const dictTable& modelDicts
71 )
72 {
73  forAllConstIter(dictTable, modelDicts, iter)
74  {
75  const phasePairKey& key = iter.key();
76 
77  // pair already exists
78  if (phasePairs_.found(key))
79  {
80  // do nothing ...
81  }
82 
83  // new ordered pair
84  else if (key.ordered())
85  {
86  phasePairs_.insert
87  (
88  key,
89  autoPtr<phasePair>
90  (
91  new orderedPhasePair
92  (
93  phaseModels_[key.first()],
94  phaseModels_[key.second()]
95  )
96  )
97  );
98  }
99 
100  // new unordered pair
101  else
102  {
103  phasePairs_.insert
104  (
105  key,
106  autoPtr<phasePair>
107  (
108  new phasePair
109  (
110  phaseModels_[key.first()],
111  phaseModels_[key.second()]
112  )
113  )
114  );
115  }
116  }
117 }
118 
119 
120 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
121 
123 (
124  const fvMesh& mesh
125 )
126 :
127  IOdictionary
128  (
129  IOobject
130  (
131  "phaseProperties",
132  mesh.time().constant(),
133  mesh,
134  IOobject::MUST_READ_IF_MODIFIED,
135  IOobject::NO_WRITE
136  )
137  ),
138 
139  mesh_(mesh),
140 
141  phaseModels_(lookup("phases"), phaseModel::iNew(*this)),
142 
143  phi_(calcPhi(phaseModels_)),
144 
145  dpdt_
146  (
147  IOobject
148  (
149  "dpdt",
150  mesh.time().timeName(),
151  mesh
152  ),
153  mesh,
155  ),
156 
157  MRF_(mesh_)
158 {
159  phi_.writeOpt() = IOobject::AUTO_WRITE;
160 
161  // Blending methods
162  forAllConstIter(dictionary, subDict("blending"), iter)
163  {
164  blendingMethods_.insert
165  (
166  iter().dict().dictName(),
168  (
169  iter().dict(),
170  phaseModels_.toc()
171  )
172  );
173  }
174 
175  // Sub-models
176  generatePairsAndSubModels("surfaceTension", surfaceTensionModels_);
177  generatePairsAndSubModels("aspectRatio", aspectRatioModels_);
178 
179  correctKinematics();
180 }
181 
182 
183 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
184 
186 {}
187 
188 
189 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
190 
192 {
193  tmp<volScalarField> tmpRho
194  (
195  phaseModels_[0]*phaseModels_[0].rho()
196  );
197 
198  for (label phasei=1; phasei<phaseModels_.size(); phasei++)
199  {
200  tmpRho.ref() += phaseModels_[phasei]*phaseModels_[phasei].rho();
201  }
202 
203  return tmpRho;
204 }
205 
206 
208 {
209  tmp<volVectorField> tmpU
210  (
211  phaseModels_[0]*phaseModels_[0].U()
212  );
213 
214  for (label phasei=1; phasei<phaseModels_.size(); phasei++)
215  {
216  tmpU.ref() += phaseModels_[phasei]*phaseModels_[phasei].U();
217  }
218 
219  return tmpU;
220 }
221 
222 
224 Foam::phaseSystem::E(const phasePairKey& key) const
225 {
226  if (aspectRatioModels_.found(key))
227  {
228  return aspectRatioModels_[key]->E();
229  }
230  else
231  {
232  return tmp<volScalarField>
233  (
234  new volScalarField
235  (
236  IOobject
237  (
238  aspectRatioModel::typeName + ":E",
239  this->mesh_.time().timeName(),
240  this->mesh_,
243  false
244  ),
245  this->mesh_,
246  dimensionedScalar("zero", dimless, 1)
247  )
248  );
249  }
250 }
251 
252 
254 Foam::phaseSystem::sigma(const phasePairKey& key) const
255 {
256  if (surfaceTensionModels_.found(key))
257  {
258  return surfaceTensionModels_[key]->sigma();
259  }
260  else
261  {
262  return tmp<volScalarField>
263  (
264  new volScalarField
265  (
266  IOobject
267  (
268  surfaceTensionModel::typeName + ":sigma",
269  this->mesh_.time().timeName(),
270  this->mesh_,
273  false
274  ),
275  this->mesh_,
277  )
278  );
279  }
280 }
281 
282 
284 {}
285 
286 
288 {
289  forAll(phaseModels_, phasei)
290  {
291  phaseModels_[phasei].correct();
292  }
293 }
294 
295 
297 {
298  bool updateDpdt = false;
299 
300  forAll(phaseModels_, phasei)
301  {
302  phaseModels_[phasei].correctKinematics();
303 
304  updateDpdt = updateDpdt || phaseModels_[phasei].thermo().dpdt();
305  }
306 
307  // Update the pressure time-derivative if required
308  if (updateDpdt)
309  {
310  dpdt_ = fvc::ddt(phaseModels_.begin()().thermo().p());
311  }
312 }
313 
314 
316 {
317  forAll(phaseModels_, phasei)
318  {
319  phaseModels_[phasei].correctThermo();
320  }
321 }
322 
323 
325 {
326  forAll(phaseModels_, phasei)
327  {
328  phaseModels_[phasei].correctTurbulence();
329  }
330 }
331 
332 
334 {
335  forAll(phaseModels_, phasei)
336  {
337  phaseModels_[phasei].correctEnergyTransport();
338  }
339 }
340 
341 
343 {
344  if (regIOobject::read())
345  {
346  bool readOK = true;
347 
348  forAll(phaseModels_, phasei)
349  {
350  readOK &= phaseModels_[phasei].read();
351  }
352 
353  // models ...
354 
355  return readOK;
356  }
357  else
358  {
359  return false;
360  }
361 }
362 
363 
364 // ************************************************************************* //
surfaceScalarField & phi
static const word propertiesName
Default name of the phase properties dictionary.
Definition: phaseSystem.H:264
dictionary dict
void generatePairs(const dictTable &modelDicts)
Generate pairs.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
U
Definition: pEqn.H:83
virtual void correctTurbulence()
Correct the turbulence.
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
label phasei
Definition: pEqn.H:27
virtual bool read()
Read object.
tmp< volVectorField > U() const
Return the mixture velocity.
virtual void solve()
Solve for the phase fractions.
tmp< volScalarField > E(const phasePairKey &key) const
Return the aspect-ratio.
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Calculate the first temporal derivative.
psiReactionThermo & thermo
Definition: createFields.H:31
virtual void correctKinematics()
Correct the kinematics.
stressControl lookup("compactNormalStress") >> compactNormalStress
phaseSystem(const fvMesh &mesh)
Construct from fvMesh.
static const dimensionSet dimSigma
Coefficient dimensions.
A class for handling words, derived from string.
Definition: word.H:59
virtual void correctThermo()
Correct the thermodynamics.
word timeName
Definition: getTimeIndex.H:3
const dimensionSet dimPressure
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
static autoPtr< blendingMethod > New(const dictionary &dict, const wordList &phaseNames)
virtual void correct()
Correct the fluid properties other than the thermo and turbulence.
tmp< volScalarField > rho() const
Return the mixture density.
defineTypeNameAndDebug(combustionModel, 0)
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
virtual void correctEnergyTransport()
Correct the energy transport e.g. alphat.
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
virtual ~phaseSystem()
Destructor.
word dictName("noiseDict")
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual bool read()
Read base phaseProperties dictionary.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
phasePairKey()
Construct null.
A class for managing temporary objects.
Definition: PtrList.H:54
tmp< surfaceScalarField > calcPhi(const phaseModelList &phaseModels) const
Calculate and return the mixture flux.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
tmp< volScalarField > sigma(const phasePairKey &key) const
Return the surface tension coefficient.
Namespace for OpenFOAM.