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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2015-2018 OpenFOAM Foundation
9  Copyright (C) 2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "phaseSystem.H"
30 #include "surfaceTensionModel.H"
31 #include "aspectRatioModel.H"
32 #include "surfaceInterpolate.H"
33 #include "fvcDdt.H"
34 #include "localEulerDdtScheme.H"
35 
36 #include "dragModel.H"
37 #include "BlendedInterfacialModel.H"
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43  defineTypeNameAndDebug(phaseSystem, 0);
44 }
45 
46 const Foam::word Foam::phaseSystem::propertiesName("phaseProperties");
47 
48 
49 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
50 
52 (
53  const phaseModelList& phaseModels
54 ) const
55 {
57  (
59  (
60  "phi",
61  fvc::interpolate(phaseModels[0])*phaseModels[0].phi()
62  )
63  );
64 
65  for (label phasei=1; phasei<phaseModels.size(); ++phasei)
66  {
67  tphi.ref() +=
68  fvc::interpolate(phaseModels[phasei])*phaseModels[phasei].phi();
69  }
70 
71  return tphi;
72 }
73 
74 
76 (
77  const dictTable& modelDicts
78 )
79 {
80  forAllConstIters(modelDicts, iter)
81  {
82  const phasePairKey& key = iter.key();
83 
84  // pair already exists
85  if (phasePairs_.found(key))
86  {}
87 
88  // new ordered pair
89  else if (key.ordered())
90  {
91  phasePairs_.insert
92  (
93  key,
95  (
97  (
98  phaseModels_[key.first()],
99  phaseModels_[key.second()]
100  )
101  )
102  );
103  }
104 
105  // new unordered pair
106  else
107  {
108  phasePairs_.insert
109  (
110  key,
112  (
113  new phasePair
114  (
115  phaseModels_[key.first()],
116  phaseModels_[key.second()]
117  )
118  )
119  );
120  }
121  }
122 }
123 
124 
125 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
126 
128 (
129  const fvMesh& mesh
130 )
131 :
133  (
134  IOobject
135  (
136  "phaseProperties",
137  mesh.time().constant(),
138  mesh,
139  IOobject::MUST_READ_IF_MODIFIED,
140  IOobject::NO_WRITE
141  )
142  ),
143 
144  mesh_(mesh),
145 
146  phaseModels_(lookup("phases"), phaseModel::iNew(*this)),
147 
148  phi_(calcPhi(phaseModels_)),
149 
150  dpdt_
151  (
152  IOobject
153  (
154  "dpdt",
155  mesh.time().timeName(),
156  mesh
157  ),
158  mesh,
160  ),
161 
162  MRF_(mesh_)
163 {
164  // Groupings
165  label movingPhasei = 0;
166  label stationaryPhasei = 0;
167  label anisothermalPhasei = 0;
168  label multiComponentPhasei = 0;
169  forAll(phaseModels_, phasei)
170  {
171  phaseModel& phase = phaseModels_[phasei];
172  movingPhasei += !phase.stationary();
173  stationaryPhasei += phase.stationary();
174  anisothermalPhasei += !phase.isothermal();
175  multiComponentPhasei += !phase.pure();
176  }
177  movingPhaseModels_.resize(movingPhasei);
178  stationaryPhaseModels_.resize(stationaryPhasei);
179  anisothermalPhaseModels_.resize(anisothermalPhasei);
180  multiComponentPhaseModels_.resize(multiComponentPhasei);
181 
182  movingPhasei = 0;
183  stationaryPhasei = 0;
184  anisothermalPhasei = 0;
185  multiComponentPhasei = 0;
186  forAll(phaseModels_, phasei)
187  {
188  phaseModel& phase = phaseModels_[phasei];
189  if (!phase.stationary())
190  {
191  movingPhaseModels_.set(movingPhasei ++, &phase);
192  }
193  if (phase.stationary())
194  {
195  stationaryPhaseModels_.set(stationaryPhasei ++, &phase);
196  }
197  if (!phase.isothermal())
198  {
199  anisothermalPhaseModels_.set(anisothermalPhasei ++, &phase);
200  }
201  if (!phase.pure())
202  {
203  multiComponentPhaseModels_.set(multiComponentPhasei ++, &phase);
204  }
205  }
206 
207  // Write phi
208  phi_.writeOpt() = IOobject::AUTO_WRITE;
209 
210  // Blending methods
211  forAllConstIter(dictionary, subDict("blending"), iter)
212  {
213  blendingMethods_.insert
214  (
215  iter().keyword(),
217  (
218  iter().keyword(),
219  iter().dict(),
220  phaseModels_.toc()
221  )
222  );
223  }
224 
225  // Sub-models
226  generatePairsAndSubModels("surfaceTension", surfaceTensionModels_);
227  generatePairsAndSubModels("aspectRatio", aspectRatioModels_);
228 
229  // Update motion fields
230  correctKinematics();
231 }
232 
233 
234 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
235 
237 {}
238 
239 
240 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
241 
243 {
244  auto phasei = movingPhaseModels_.cbegin();
245 
247 
248  for (++phasei; phasei != movingPhaseModels_.cend(); ++phasei)
249  {
250  trho.ref() += phasei()*phasei().rho();
251  }
252 
253  if (stationaryPhaseModels_.empty())
254  {
255  return trho;
256  }
257 
258  phasei = movingPhaseModels_.cbegin();
259 
261  for (++phasei; phasei != movingPhaseModels_.cend(); ++phasei)
262  {
263  alpha += phasei();
264  }
265 
266  return trho/alpha;
267 }
268 
269 
271 {
272  auto phasei = movingPhaseModels_.cbegin();
273 
275 
276  for (++phasei; phasei != movingPhaseModels_.cend(); ++phasei)
277  {
278  tU.ref() += phasei()*phasei().U();
279  }
280 
281  if (stationaryPhaseModels_.empty())
282  {
283  return tU;
284  }
285 
286  phasei = movingPhaseModels_.cbegin();
287 
289  for (++phasei; phasei != movingPhaseModels_.cend(); ++phasei)
290  {
291  alpha += phasei();
292  }
293 
294 
295  return tU/alpha;
296 }
297 
298 
301 {
302  if (aspectRatioModels_.found(key))
303  {
304  return aspectRatioModels_[key]->E();
305  }
306  else
307  {
308  return volScalarField::New
309  (
310  aspectRatioModel::typeName + ":E",
311  this->mesh_,
312  dimensionedScalar("one", dimless, 1)
313  );
314  }
315 }
316 
317 
320 {
321  if (surfaceTensionModels_.found(key))
322  {
323  return surfaceTensionModels_[key]->sigma();
324  }
325  else
326  {
327  return volScalarField::New
328  (
329  surfaceTensionModel::typeName + ":sigma",
330  this->mesh_,
332  );
333  }
334 }
335 
336 
338 (
339  const phasePairKey& key
340 ) const
341 {
342  return volScalarField::New
343  (
344  IOobject::groupName("dmdt", phasePairs_[key]->name()),
345  this->mesh_,
347  );
348 }
349 
350 
352 {
353  PtrList<volScalarField> dmdts(this->phaseModels_.size());
354 
355  return dmdts;
356 }
357 
358 
360 {}
361 
362 
364 {
365  for (phaseModel& phase : phaseModels_)
366  {
367  phase.correct();
368  }
369 }
370 
371 
373 {
374  bool updateDpdt = false;
375 
376  for (phaseModel& phase : phaseModels_)
377  {
378  phase.correctKinematics();
379 
380  updateDpdt = updateDpdt || phase.thermo().dpdt();
381  }
382 
383  // Update the pressure time-derivative if required
384  if (updateDpdt)
385  {
386  dpdt_ = fvc::ddt(phaseModels_.cbegin()().thermo().p());
387  }
388 }
389 
390 
392 {
393  for (phaseModel& phase : phaseModels_)
394  {
395  phase.correctThermo();
396  }
397 }
398 
399 
401 {
402  for (phaseModel& phase : phaseModels_)
403  {
404  phase.correctTurbulence();
405  }
406 }
407 
408 
410 {
411  for (phaseModel& phase : phaseModels_)
412  {
413  phase.correctEnergyTransport();
414  }
415 }
416 
417 
419 {
420  if (regIOobject::read())
421  {
422  bool readOK = true;
423 
424  for (phaseModel& phase : phaseModels_)
425  {
426  readOK &= phase.read();
427  }
428 
429  // models ...
430 
431  return readOK;
432  }
433 
434  return false;
435 }
436 
437 
439 {
440  if (fv::localEulerDdt::enabled(vf.mesh()))
441  {
442  return fv::localEulerDdt::localRDeltaT(vf.mesh())*vf;
443  }
444  else
445  {
446  return vf/vf.mesh().time().deltaT();
447  }
448 }
449 
450 
452 {
453  if (fv::localEulerDdt::enabled(sf.mesh()))
454  {
455  return fv::localEulerDdt::localRDeltaTf(sf.mesh())*sf;
456  }
457  else
458  {
459  return sf/sf.mesh().time().deltaT();
460  }
461 }
462 
463 
464 // ************************************************************************* //
Foam::Pair::second
const T & second() const noexcept
Return second element, which is also the last element.
Definition: PairI.H:122
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
Foam::phaseModel
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:57
Foam::dimPressure
const dimensionSet dimPressure
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::phasePair
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition: phasePair.H:51
Foam::phase::read
bool read(const dictionary &phaseDict)
Read base transportProperties dictionary.
Foam::fv::localEulerDdt::localRDeltaT
static const volScalarField & localRDeltaT(const fvMesh &mesh)
Return the reciprocal of the local time-step.
Definition: localEulerDdt.C:48
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::phaseSystem::calcPhi
tmp< surfaceScalarField > calcPhi(const phaseModelList &phaseModels) const
Calculate and return the mixture flux.
Definition: phaseSystem.C:52
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::basicThermo::p
virtual volScalarField & p()
Pressure [Pa].
Definition: basicThermo.C:512
Foam::phase
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:54
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::phaseSystem::solve
virtual void solve()
Solve for the phase fractions.
Definition: phaseSystem.C:359
Foam::phaseSystem::sigma
tmp< volScalarField > sigma(const phasePairKey &key) const
Return the surface tension coefficient for a pair.
Definition: phaseSystem.C:319
Foam::dimDensity
const dimensionSet dimDensity
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
thermo
psiReactionThermo & thermo
Definition: createFields.H:28
Foam::surfaceTensionModel::dimSigma
static const dimensionSet dimSigma
Coefficient dimensions.
Definition: surfaceTensionModel.H:90
Foam::phaseSystem::phaseSystem
phaseSystem(const fvMesh &mesh)
Construct from fvMesh.
Definition: phaseSystem.C:128
Foam::regIOobject::read
virtual bool read()
Read object.
Definition: regIOobjectRead.C:202
Foam::fv::localEulerDdt::localRDeltaTf
static const surfaceScalarField & localRDeltaTf(const fvMesh &mesh)
Return the reciprocal of the local face time-step.
Definition: localEulerDdt.C:60
phasei
label phasei
Definition: pEqn.H:27
Foam::phaseSystem::U
tmp< volVectorField > U() const
Return the mixture velocity.
Definition: phaseSystem.C:270
Foam::phaseSystem::~phaseSystem
virtual ~phaseSystem()
Destructor.
Definition: phaseSystem.C:236
Foam::phaseSystem::dmdts
virtual PtrList< volScalarField > dmdts() const
Return the mass transfer rates for each phase.
Definition: phaseSystem.C:351
rho
rho
Definition: readInitialConditions.H:96
forAllConstIter
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:338
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::GeometricField< scalar, fvPatchField, volMesh >::New
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &ds, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return tmp field from name, mesh, dimensions and patch type.
Definition: GeometricFieldNew.C:34
Foam::PtrListDictionary< phaseModel >
Foam::fv::localEulerDdt::enabled
static bool enabled(const fvMesh &mesh)
Return true if LTS is enabled.
Definition: localEulerDdt.C:39
Foam::phase::correct
void correct()
Correct the phase properties.
trho
tmp< volScalarField > trho
Definition: setRegionSolidFields.H:4
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:258
localEulerDdtScheme.H
Foam::label
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:62
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::cellModeller::lookup
const cellModel * lookup(const word &modelName)
Deprecated(2017-11) equivalent to cellModel::ptr static method.
Definition: cellModeller.H:46
Foam::phaseSystem::correctKinematics
virtual void correctKinematics()
Correct the kinematics.
Definition: phaseSystem.C:372
Foam::orderedPhasePair
Definition: orderedPhasePair.H:49
Foam::phaseSystem::correct
virtual void correct()
Correct the fluid properties other than those listed below.
Definition: phaseSystem.C:363
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:65
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::phasePairKey
Definition: phasePairKey.H:59
Foam::phaseModel::iNew
Return a pointer to a new phase created on freestore.
Definition: phaseModel.H:129
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::phaseSystem::correctEnergyTransport
virtual void correctEnergyTransport()
Correct the energy transport e.g. alphat.
Definition: phaseSystem.C:409
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::phaseSystem::propertiesName
static const word propertiesName
Default name of the phase properties dictionary.
Definition: phaseSystem.H:287
Foam::phaseSystem::dmdt
virtual tmp< volScalarField > dmdt(const phasePairKey &key) const
Return the mass transfer rate for a pair.
Definition: phaseSystem.C:338
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::phasePairKey::ordered
bool ordered() const
Return the ordered flag.
Definition: phasePairKey.C:60
phaseSystem.H
Foam::HashTable
A HashTable similar to std::unordered_map.
Definition: HashTable.H:105
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::autoPtr< phasePair >
U
U
Definition: pEqn.H:72
Foam::phaseSystem::correctTurbulence
virtual void correctTurbulence()
Correct the turbulence.
Definition: phaseSystem.C:400
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
Foam::IOobject::groupName
static word groupName(StringType name, const word &group)
Create dot-delimited name.group.
Foam::phaseSystem::correctThermo
virtual void correctThermo()
Correct the thermodynamics.
Definition: phaseSystem.C:391
fvcDdt.H
Calculate the first temporal derivative.
Foam::fvc::ddt
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:47
Foam::byDt
tmp< volScalarField > byDt(const volScalarField &vf)
Definition: phaseSystem.C:438
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::phaseSystem::rho
tmp< volScalarField > rho() const
Return the mixture density.
Definition: phaseSystem.C:242
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::phaseSystem::read
virtual bool read()
Read base phaseProperties dictionary.
Definition: phaseSystem.C:418
Foam::fac::interpolate
static tmp< GeometricField< Type, faePatchField, edgeMesh > > interpolate(const GeometricField< Type, faPatchField, areaMesh > &tvf, const edgeScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Foam::phaseSystem::E
tmp< volScalarField > E(const phasePairKey &key) const
Return the aspect-ratio for a pair.
Definition: phaseSystem.C:300
Foam::phaseSystem::generatePairs
void generatePairs(const dictTable &modelDicts)
Generate pairs.
Definition: phaseSystem.C:76