phaseSystem.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) 2015-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 "phaseSystem.H"
27 #include "surfaceTensionModel.H"
28 #include "aspectRatioModel.H"
29 #include "surfaceInterpolate.H"
30 #include "fvcDdt.H"
31 #include "localEulerDdtScheme.H"
32 
33 #include "dragModel.H"
34 #include "BlendedInterfacialModel.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40  defineTypeNameAndDebug(phaseSystem, 0);
41 }
42 
43 const Foam::word Foam::phaseSystem::propertiesName("phaseProperties");
44 
45 
46 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
47 
49 (
50  const phaseModelList& phaseModels
51 ) const
52 {
53  tmp<surfaceScalarField> tmpPhi
54  (
56  (
57  "phi",
58  fvc::interpolate(phaseModels[0])*phaseModels[0].phi()
59  )
60  );
61 
62  for (label phasei=1; phasei<phaseModels.size(); phasei++)
63  {
64  tmpPhi.ref() +=
65  fvc::interpolate(phaseModels[phasei])*phaseModels[phasei].phi();
66  }
67 
68  return tmpPhi;
69 }
70 
71 
73 (
74  const dictTable& modelDicts
75 )
76 {
77  forAllConstIter(dictTable, modelDicts, iter)
78  {
79  const phasePairKey& key = iter.key();
80 
81  // pair already exists
82  if (phasePairs_.found(key))
83  {}
84 
85  // new ordered pair
86  else if (key.ordered())
87  {
88  phasePairs_.insert
89  (
90  key,
91  autoPtr<phasePair>
92  (
93  new orderedPhasePair
94  (
95  phaseModels_[key.first()],
96  phaseModels_[key.second()]
97  )
98  )
99  );
100  }
101 
102  // new unordered pair
103  else
104  {
105  phasePairs_.insert
106  (
107  key,
108  autoPtr<phasePair>
109  (
110  new phasePair
111  (
112  phaseModels_[key.first()],
113  phaseModels_[key.second()]
114  )
115  )
116  );
117  }
118  }
119 }
120 
121 
122 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
123 
125 (
126  const fvMesh& mesh
127 )
128 :
129  IOdictionary
130  (
131  IOobject
132  (
133  "phaseProperties",
134  mesh.time().constant(),
135  mesh,
136  IOobject::MUST_READ_IF_MODIFIED,
137  IOobject::NO_WRITE
138  )
139  ),
140 
141  mesh_(mesh),
142 
143  phaseModels_(lookup("phases"), phaseModel::iNew(*this)),
144 
145  phi_(calcPhi(phaseModels_)),
146 
147  dpdt_
148  (
149  IOobject
150  (
151  "dpdt",
152  mesh.time().timeName(),
153  mesh
154  ),
155  mesh,
157  ),
158 
159  MRF_(mesh_)
160 {
161  // Groupings
162  label movingPhasei = 0;
163  label stationaryPhasei = 0;
164  label anisothermalPhasei = 0;
165  label multiComponentPhasei = 0;
166  forAll(phaseModels_, phasei)
167  {
168  phaseModel& phase = phaseModels_[phasei];
169  movingPhasei += !phase.stationary();
170  stationaryPhasei += phase.stationary();
171  anisothermalPhasei += !phase.isothermal();
172  multiComponentPhasei += !phase.pure();
173  }
174  movingPhaseModels_.resize(movingPhasei);
175  stationaryPhaseModels_.resize(stationaryPhasei);
176  anisothermalPhaseModels_.resize(anisothermalPhasei);
177  multiComponentPhaseModels_.resize(multiComponentPhasei);
178 
179  movingPhasei = 0;
180  stationaryPhasei = 0;
181  anisothermalPhasei = 0;
182  multiComponentPhasei = 0;
183  forAll(phaseModels_, phasei)
184  {
185  phaseModel& phase = phaseModels_[phasei];
186  if (!phase.stationary())
187  {
188  movingPhaseModels_.set(movingPhasei ++, &phase);
189  }
190  if (phase.stationary())
191  {
192  stationaryPhaseModels_.set(stationaryPhasei ++, &phase);
193  }
194  if (!phase.isothermal())
195  {
196  anisothermalPhaseModels_.set(anisothermalPhasei ++, &phase);
197  }
198  if (!phase.pure())
199  {
200  multiComponentPhaseModels_.set(multiComponentPhasei ++, &phase);
201  }
202  }
203 
204  // Write phi
205  phi_.writeOpt() = IOobject::AUTO_WRITE;
206 
207  // Blending methods
208  forAllConstIter(dictionary, subDict("blending"), iter)
209  {
210  blendingMethods_.insert
211  (
212  iter().keyword(),
214  (
215  iter().keyword(),
216  iter().dict(),
217  phaseModels_.toc()
218  )
219  );
220  }
221 
222  // Sub-models
223  generatePairsAndSubModels("surfaceTension", surfaceTensionModels_);
224  generatePairsAndSubModels("aspectRatio", aspectRatioModels_);
225 
226  // Update motion fields
227  correctKinematics();
228 }
229 
230 
231 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
232 
234 {}
235 
236 
237 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
238 
240 {
241  const label nMovingPhases = movingPhaseModels_.size();
242 
243  tmp<volScalarField> rho(movingPhaseModels_[0]*movingPhaseModels_[0].rho());
244  for (label movingPhasei = 1; movingPhasei < nMovingPhases; ++ movingPhasei)
245  {
246  rho.ref() +=
247  movingPhaseModels_[movingPhasei]
248  *movingPhaseModels_[movingPhasei].rho();
249  }
250 
251  if (stationaryPhaseModels_.empty())
252  {
253  return rho;
254  }
255 
256  volScalarField alpha(movingPhaseModels_[0]);
257  for (label movingPhasei = 1; movingPhasei < nMovingPhases; ++ movingPhasei)
258  {
259  alpha += movingPhaseModels_[movingPhasei];
260  }
261 
262  return rho/alpha;
263 }
264 
265 
267 {
268  const label nMovingPhases = movingPhaseModels_.size();
269 
270  tmp<volVectorField> U(movingPhaseModels_[0]*movingPhaseModels_[0].U());
271  for (label movingPhasei = 1; movingPhasei < nMovingPhases; ++ movingPhasei)
272  {
273  U.ref() +=
274  movingPhaseModels_[movingPhasei]
275  *movingPhaseModels_[movingPhasei].U();
276  }
277 
278  if (stationaryPhaseModels_.empty())
279  {
280  return U;
281  }
282 
283  volScalarField alpha(movingPhaseModels_[0]);
284  for (label movingPhasei = 1; movingPhasei < nMovingPhases; ++ movingPhasei)
285  {
286  alpha += movingPhaseModels_[movingPhasei];
287  }
288 
289  return U/alpha;
290 }
291 
292 
294 Foam::phaseSystem::E(const phasePairKey& key) const
295 {
296  if (aspectRatioModels_.found(key))
297  {
298  return aspectRatioModels_[key]->E();
299  }
300  else
301  {
302  return tmp<volScalarField>
303  (
304  new volScalarField
305  (
306  IOobject
307  (
308  aspectRatioModel::typeName + ":E",
309  this->mesh_.time().timeName(),
310  this->mesh_,
313  false
314  ),
315  this->mesh_,
316  dimensionedScalar("zero", dimless, 1)
317  )
318  );
319  }
320 }
321 
322 
324 Foam::phaseSystem::sigma(const phasePairKey& key) const
325 {
326  if (surfaceTensionModels_.found(key))
327  {
328  return surfaceTensionModels_[key]->sigma();
329  }
330  else
331  {
332  return tmp<volScalarField>
333  (
334  new volScalarField
335  (
336  IOobject
337  (
338  surfaceTensionModel::typeName + ":sigma",
339  this->mesh_.time().timeName(),
340  this->mesh_,
343  false
344  ),
345  this->mesh_,
347  )
348  );
349  }
350 }
351 
352 
354 (
355  const phasePairKey& key
356 ) const
357 {
358  return tmp<volScalarField>
359  (
360  new volScalarField
361  (
362  IOobject
363  (
364  IOobject::groupName("dmdt", phasePairs_[key]->name()),
365  this->mesh_.time().timeName(),
366  this->mesh_
367  ),
368  this->mesh_,
370  )
371  );
372 }
373 
374 
376 {
377  PtrList<volScalarField> dmdts(this->phaseModels_.size());
378 
379  return dmdts.xfer();
380 }
381 
382 
384 {}
385 
386 
388 {
389  forAll(phaseModels_, phasei)
390  {
391  phaseModels_[phasei].correct();
392  }
393 }
394 
395 
397 {
398  bool updateDpdt = false;
399 
400  forAll(phaseModels_, phasei)
401  {
402  phaseModels_[phasei].correctKinematics();
403 
404  updateDpdt = updateDpdt || phaseModels_[phasei].thermo().dpdt();
405  }
406 
407  // Update the pressure time-derivative if required
408  if (updateDpdt)
409  {
410  dpdt_ = fvc::ddt(phaseModels_.begin()().thermo().p());
411  }
412 }
413 
414 
416 {
417  forAll(phaseModels_, phasei)
418  {
419  phaseModels_[phasei].correctThermo();
420  }
421 }
422 
423 
425 {
426  forAll(phaseModels_, phasei)
427  {
428  phaseModels_[phasei].correctTurbulence();
429  }
430 }
431 
432 
434 {
435  forAll(phaseModels_, phasei)
436  {
437  phaseModels_[phasei].correctEnergyTransport();
438  }
439 }
440 
441 
443 {
444  if (regIOobject::read())
445  {
446  bool readOK = true;
447 
448  forAll(phaseModels_, phasei)
449  {
450  readOK &= phaseModels_[phasei].read();
451  }
452 
453  // models ...
454 
455  return readOK;
456  }
457  else
458  {
459  return false;
460  }
461 }
462 
463 
465 {
466  if (fv::localEulerDdt::enabled(vf.mesh()))
467  {
468  return fv::localEulerDdt::localRDeltaT(vf.mesh())*vf;
469  }
470  else
471  {
472  return vf/vf.mesh().time().deltaT();
473  }
474 }
475 
476 
478 {
479  if (fv::localEulerDdt::enabled(sf.mesh()))
480  {
481  return fv::localEulerDdt::localRDeltaTf(sf.mesh())*sf;
482  }
483  else
484  {
485  return sf/sf.mesh().time().deltaT();
486  }
487 }
488 
489 
490 // ************************************************************************* //
A simple container for copying or transferring objects of type <T>.
Definition: Xfer.H:85
static const word propertiesName
Default name of the phase properties dictionary.
Definition: phaseSystem.H:285
dictionary dict
void generatePairs(const dictTable &modelDicts)
Generate pairs.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
tmp< volScalarField > E(const phasePairKey &key) const
Return the aspect-ratio for a pair.
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
static const surfaceScalarField & localRDeltaTf(const fvMesh &mesh)
Return the reciprocal of the local face time-step.
Definition: localEulerDdt.C:58
surfaceScalarField & phi
label phasei
Definition: pEqn.H:27
virtual bool read()
Read object.
virtual void solve()
Solve for the phase fractions.
tmp< volScalarField > sigma(const phasePairKey &key) const
Return the surface tension coefficient for a pair.
virtual tmp< volScalarField > dmdt(const phasePairKey &key) const
Return the mass transfer rate for a pair.
rhoReactionThermo & thermo
Definition: createFields.H:28
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.
virtual void correctKinematics()
Correct the kinematics.
stressControl lookup("compactNormalStress") >> compactNormalStress
static bool enabled(const fvMesh &mesh)
Return true if LTS is enabled.
Definition: localEulerDdt.C:37
phaseSystem(const fvMesh &mesh)
Construct from fvMesh.
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
tmp< volVectorField > U() const
Return the mixture velocity.
static autoPtr< blendingMethod > New(const word &modelName, const dictionary &dict, const wordList &phaseNames)
virtual void correctThermo()
Correct the thermodynamics.
word timeName
Definition: getTimeIndex.H:3
static const volScalarField & localRDeltaT(const fvMesh &mesh)
Return the reciprocal of the local time-step.
Definition: localEulerDdt.C:46
const dimensionSet dimPressure
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
virtual void correct()
Correct the fluid properties other than those listed below.
tmp< volScalarField > rho() const
Return the mixture density.
const Mesh & mesh() const
Return mesh.
defineTypeNameAndDebug(combustionModel, 0)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
tmp< volScalarField > byDt(const volScalarField &vf)
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.
const dimensionSet dimDensity
virtual void correctEnergyTransport()
Correct the energy transport e.g. alphat.
tmp< surfaceScalarField > calcPhi(const phaseModelList &phaseModels) const
Calculate and return the mixture flux.
U
Definition: pEqn.H:72
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
static const dimensionSet dimSigma
Surface tension coefficient dimensions.
virtual ~phaseSystem()
Destructor.
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:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Namespace for OpenFOAM.
virtual Xfer< PtrList< volScalarField > > dmdts() const
Return the mass transfer rates for each phase.