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-2023 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 {
56  auto tphi = surfaceScalarField::New
57  (
58  "phi",
60  fvc::interpolate(phaseModels[0])*phaseModels[0].phi()
61  );
62 
63  for (label phasei=1; phasei<phaseModels.size(); ++phasei)
64  {
65  tphi.ref() +=
66  fvc::interpolate(phaseModels[phasei])*phaseModels[phasei].phi();
67  }
68 
69  return tphi;
70 }
71 
72 
74 (
75  const dictTable& modelDicts
76 )
77 {
78  forAllConstIters(modelDicts, iter)
79  {
80  const phasePairKey& key = iter.key();
81 
82  // pair already exists
83  if (phasePairs_.found(key))
84  {}
85 
86  // new ordered pair
87  else if (key.ordered())
88  {
89  phasePairs_.insert
90  (
91  key,
92  autoPtr<phasePair>
93  (
94  new orderedPhasePair
95  (
96  phaseModels_[key.first()],
97  phaseModels_[key.second()]
98  )
99  )
100  );
101  }
102 
103  // new unordered pair
104  else
105  {
106  phasePairs_.insert
107  (
108  key,
109  autoPtr<phasePair>
110  (
111  new phasePair
112  (
113  phaseModels_[key.first()],
114  phaseModels_[key.second()]
115  )
116  )
117  );
118  }
119  }
120 }
121 
122 
123 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
124 
126 (
127  const fvMesh& mesh
128 )
129 :
130  IOdictionary
131  (
132  IOobject
133  (
134  "phaseProperties",
135  mesh.time().constant(),
136  mesh,
137  IOobject::READ_MODIFIED,
138  IOobject::NO_WRITE,
139  IOobject::REGISTER
140  )
141  ),
142 
143  mesh_(mesh),
144 
145  phaseModels_(lookup("phases"), phaseModel::iNew(*this)),
146 
147  phi_(calcPhi(phaseModels_)),
148 
149  dpdt_
150  (
151  IOobject
152  (
153  "dpdt",
154  mesh.time().timeName(),
155  mesh
156  ),
157  mesh,
159  ),
160 
161  MRF_(mesh_)
162 {
163  // Groupings
164  label movingPhasei = 0;
165  label stationaryPhasei = 0;
166  label anisothermalPhasei = 0;
167  label multiComponentPhasei = 0;
169  {
170  phaseModel& phase = phaseModels_[phasei];
171  movingPhasei += !phase.stationary();
172  stationaryPhasei += phase.stationary();
173  anisothermalPhasei += !phase.isothermal();
174  multiComponentPhasei += !phase.pure();
175  }
176  movingPhaseModels_.resize(movingPhasei);
177  stationaryPhaseModels_.resize(stationaryPhasei);
178  anisothermalPhaseModels_.resize(anisothermalPhasei);
179  multiComponentPhaseModels_.resize(multiComponentPhasei);
180 
181  movingPhasei = 0;
182  stationaryPhasei = 0;
183  anisothermalPhasei = 0;
184  multiComponentPhasei = 0;
186  {
187  phaseModel& phase = phaseModels_[phasei];
188  if (!phase.stationary())
189  {
190  movingPhaseModels_.set(movingPhasei ++, &phase);
191  }
192  if (phase.stationary())
193  {
194  stationaryPhaseModels_.set(stationaryPhasei ++, &phase);
195  }
196  if (!phase.isothermal())
197  {
198  anisothermalPhaseModels_.set(anisothermalPhasei ++, &phase);
199  }
200  if (!phase.pure())
201  {
202  multiComponentPhaseModels_.set(multiComponentPhasei ++, &phase);
203  }
204  }
205 
206  // Write phi
208 
209  // Blending methods
210  forAllConstIter(dictionary, subDict("blending"), iter)
211  {
213  (
214  iter().keyword(),
216  (
217  iter().keyword(),
218  iter().dict(),
219  phaseModels_.toc()
220  )
221  );
222  }
223 
224  // Sub-models
227 
228  // Update motion fields
230 }
231 
232 
233 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
236 {}
237 
238 
239 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
240 
242 {
243  auto phasei = movingPhaseModels_.cbegin();
244 
246 
247  for (++phasei; phasei != movingPhaseModels_.cend(); ++phasei)
248  {
249  trho.ref() += phasei()*phasei().rho();
250  }
251 
252  if (stationaryPhaseModels_.empty())
253  {
254  return trho;
255  }
256 
257  phasei = movingPhaseModels_.cbegin();
258 
260  for (++phasei; phasei != movingPhaseModels_.cend(); ++phasei)
261  {
262  alpha += phasei();
263  }
264 
265  return trho/alpha;
266 }
267 
268 
270 {
271  auto phasei = movingPhaseModels_.cbegin();
272 
273  tmp<volVectorField> tU(phasei()*phasei().U());
274 
275  for (++phasei; phasei != movingPhaseModels_.cend(); ++phasei)
276  {
277  tU.ref() += phasei()*phasei().U();
278  }
279 
280  if (stationaryPhaseModels_.empty())
281  {
282  return tU;
283  }
284 
285  phasei = movingPhaseModels_.cbegin();
286 
288  for (++phasei; phasei != movingPhaseModels_.cend(); ++phasei)
289  {
290  alpha += phasei();
291  }
293 
294  return tU/alpha;
295 }
296 
297 
300 {
301  if (aspectRatioModels_.found(key))
302  {
303  return aspectRatioModels_[key]->E();
304  }
305 
306  return volScalarField::New
307  (
308  IOobject::scopedName(aspectRatioModel::typeName, "E"),
310  this->mesh_,
311  dimensionedScalar("one", dimless, 1)
312  );
313 }
314 
315 
318 {
319  if (surfaceTensionModels_.found(key))
320  {
321  return surfaceTensionModels_[key]->sigma();
322  }
323 
324  return volScalarField::New
325  (
327  (
328  reactingMultiphaseEuler::surfaceTensionModel::typeName, "sigma"
329  ),
331  this->mesh_,
333  (
336  )
337  );
338 }
339 
340 
342 (
343  const phasePairKey& key
344 ) const
345 {
346  return volScalarField::New
347  (
348  IOobject::groupName("dmdt", phasePairs_[key]->name()),
350  this->mesh_,
352  );
353 }
354 
355 
357 {
358  PtrList<volScalarField> dmdts(this->phaseModels_.size());
359 
360  return dmdts;
361 }
362 
363 
365 {}
366 
367 
369 {
370  for (phaseModel& phase : phaseModels_)
371  {
372  phase.correct();
373  }
374 }
375 
376 
378 {
379  bool updateDpdt = false;
380 
381  for (phaseModel& phase : phaseModels_)
382  {
383  phase.correctKinematics();
384 
385  updateDpdt = updateDpdt || phase.thermo().dpdt();
386  }
387 
388  // Update the pressure time-derivative if required
389  if (updateDpdt)
390  {
391  dpdt_ = fvc::ddt(phaseModels_.cbegin()().thermo().p());
392  }
393 }
394 
395 
397 {
398  for (phaseModel& phase : phaseModels_)
399  {
400  phase.correctThermo();
401  }
402 }
403 
404 
406 {
407  for (phaseModel& phase : phaseModels_)
408  {
409  phase.correctTurbulence();
410  }
411 }
412 
413 
415 {
416  for (phaseModel& phase : phaseModels_)
417  {
418  phase.correctEnergyTransport();
419  }
420 }
421 
422 
424 {
425  if (regIOobject::read())
426  {
427  bool readOK = true;
428 
429  for (phaseModel& phase : phaseModels_)
430  {
431  readOK &= phase.read();
432  }
433 
434  // models ...
435 
436  return readOK;
437  }
438 
439  return false;
440 }
441 
442 
444 {
445  if (fv::localEulerDdt::enabled(vf.mesh()))
446  {
447  return fv::localEulerDdt::localRDeltaT(vf.mesh())*vf;
448  }
449  else
450  {
451  return vf/vf.mesh().time().deltaT();
452  }
453 }
454 
455 
457 {
458  if (fv::localEulerDdt::enabled(sf.mesh()))
459  {
460  return fv::localEulerDdt::localRDeltaTf(sf.mesh())*sf;
461  }
462  else
463  {
464  return sf/sf.mesh().time().deltaT();
465  }
466 }
467 
468 
469 // ************************************************************************* //
virtual void correctKinematics()
Correct the kinematics.
Definition: phaseSystem.C:370
virtual void correctEnergyTransport()
Correct the energy transport e.g. alphat.
Definition: phaseSystem.C:407
dictionary dict
void generatePairs(const dictTable &modelDicts)
Generate pairs.
Definition: phaseSystem.C:67
surfaceTensionModelTable surfaceTensionModels_
Surface tension models.
Definition: phaseSystem.H:185
static const surfaceScalarField & localRDeltaTf(const fvMesh &mesh)
Return the reciprocal of the local face time-step.
Definition: localEulerDdt.C:53
label phasei
Definition: pEqn.H:27
virtual bool read()
Read object.
virtual tmp< volScalarField > dmdt(const phasePairKey &key) const
Return the mass transfer rate for a pair.
Definition: phaseSystem.C:335
tmp< volScalarField > trho
phaseModelPartialList stationaryPhaseModels_
Stationary phase models.
Definition: phaseSystem.H:142
const dimensionSet dimless
Dimensionless.
void correct()
Correct the phase properties.
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:478
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:441
phaseModelList phaseModels_
Phase models.
Definition: phaseSystem.H:132
virtual void correctTurbulence()
Correct the turbulence.
Definition: phaseSystem.C:398
phaseModelPartialList multiComponentPhaseModels_
Multi-component phase models.
Definition: phaseSystem.H:152
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
virtual volScalarField & p()
Pressure [Pa].
Definition: basicThermo.C:607
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:40
dictionary()
Default construct, a top-level empty dictionary.
Definition: dictionary.C:68
tmp< volScalarField > rho() const
Return the mixture density.
Definition: phaseSystem.C:234
virtual void correctThermo()
Correct the thermodynamics.
Definition: phaseSystem.C:389
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition: HashTableI.H:152
Calculate the first temporal derivative.
psiReactionThermo & thermo
Definition: createFields.H:28
tmp< volScalarField > sigma(const phasePairKey &key) const
Return the surface tension coefficient for a pair.
Definition: phaseSystem.C:310
word timeName
Definition: getTimeIndex.H:3
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition: IOobjectI.H:40
writeOption writeOpt() const noexcept
Get the write option.
dynamicFvMesh & mesh
static bool enabled(const fvMesh &mesh)
Return true if LTS is enabled.
Definition: localEulerDdt.C:32
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
phaseSystem(const fvMesh &mesh)
Construct from fvMesh.
Definition: phaseSystem.C:119
static autoPtr< blendingMethod > New(const word &modelName, const dictionary &dict, const wordList &phaseNames)
virtual bool read()
Read base phaseProperties dictionary.
Definition: phaseSystem.C:416
A class for handling words, derived from Foam::string.
Definition: word.H:63
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
tmp< volScalarField > E(const phasePairKey &key) const
Return the aspect-ratio for a pair.
Definition: phaseSystem.C:292
tmp< volScalarField > byDt(const volScalarField &vf)
Definition: phaseSystem.C:436
static const volScalarField & localRDeltaT(const fvMesh &mesh)
Return the reciprocal of the local time-step.
Definition: localEulerDdt.C:41
void resize(const label newLen)
Change the size of the list. Any new entries are nullptr.
Definition: UPtrListI.H:251
const dimensionSet dimPressure
phaseModelPartialList anisothermalPhaseModels_
Anisothermal phase models.
Definition: phaseSystem.H:147
virtual void solve()
Solve for the phase fractions.
Definition: phaseSystem.C:357
static tmp< GeometricField< scalar, fvsPatchField, surfaceMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvsPatchField< scalar >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
An ordered or unorder pair of phase names. Typically specified as follows.
Definition: phasePairKey.H:61
phaseModelPartialList movingPhaseModels_
Moving phase models.
Definition: phaseSystem.H:137
const Mesh & mesh() const noexcept
Return mesh.
aspectRatioModelTable aspectRatioModels_
Aspect ratio models.
Definition: phaseSystem.H:190
defineTypeNameAndDebug(combustionModel, 0)
static const dimensionSet dimSigma
Coefficient dimensions.
static const word propertiesName
Default name of the phase properties dictionary.
Definition: phaseSystem.H:340
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
void generatePairsAndSubModels(const word &modelName, HashTable< autoPtr< modelType >, phasePairKey, phasePairKey::hash > &models)
Generate pairs and sub-model tables.
U
Definition: pEqn.H:72
const surfaceScalarField & phi() const
Return the mixture flux.
Definition: phaseSystemI.H:106
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
auto key(const Type &t) -> typename std::enable_if< std::is_enum< Type >::value, typename std::underlying_type< Type >::type >::type
Definition: foamGltfBase.H:103
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
virtual void correct()
Correct the fluid properties other than those listed below.
Definition: phaseSystem.C:361
blendingMethodTable blendingMethods_
Blending methods.
Definition: phaseSystem.H:177
Automatically write from objectRegistry::writeObject()
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
tmp< volVectorField > U() const
Return the mixture velocity.
Definition: phaseSystem.C:262
tmp< surfaceScalarField > calcPhi(const phaseModelList &phaseModels) const
Calculate and return the mixture flux.
Definition: phaseSystem.C:45
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:50
virtual ~phaseSystem()
Destructor.
Definition: phaseSystem.C:228
virtual PtrList< volScalarField > dmdts() const
Return the mass transfer rates for each phase.
Definition: phaseSystem.C:349
surfaceScalarField phi_
Total volumetric flux.
Definition: phaseSystem.H:162
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
wordList toc() const
The table of contents (as a sorted list)
Do not request registration (bool: false)
Namespace for OpenFOAM.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie...
Definition: UPtrList.H:366
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127