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-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
41namespace Foam
42{
44}
45
46const 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 (
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 (
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;
170 {
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;
187 {
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
209
210 // Blending methods
211 forAllConstIter(dictionary, subDict("blending"), iter)
212 {
214 (
215 iter().keyword(),
217 (
218 iter().keyword(),
219 iter().dict(),
221 )
222 );
223 }
224
225 // Sub-models
228
229 // Update motion fields
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 {
309 (
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 {
328 (
330 this->mesh_,
332 (
334 0
335 )
336 );
337 }
338}
339
340
342(
343 const phasePairKey& key
344) const
345{
347 (
348 IOobject::groupName("dmdt", phasePairs_[key]->name()),
349 this->mesh_,
351 );
352}
353
354
356{
357 PtrList<volScalarField> dmdts(this->phaseModels_.size());
358
359 return dmdts;
360}
361
362
364{}
365
366
368{
369 for (phaseModel& phase : phaseModels_)
370 {
371 phase.correct();
372 }
373}
374
375
377{
378 bool updateDpdt = false;
379
380 for (phaseModel& phase : phaseModels_)
381 {
382 phase.correctKinematics();
383
384 updateDpdt = updateDpdt || phase.thermo().dpdt();
385 }
386
387 // Update the pressure time-derivative if required
388 if (updateDpdt)
389 {
390 dpdt_ = fvc::ddt(phaseModels_.cbegin()().thermo().p());
391 }
392}
393
394
396{
397 for (phaseModel& phase : phaseModels_)
398 {
399 phase.correctThermo();
400 }
401}
402
403
405{
406 for (phaseModel& phase : phaseModels_)
407 {
408 phase.correctTurbulence();
409 }
410}
411
412
414{
415 for (phaseModel& phase : phaseModels_)
416 {
417 phase.correctEnergyTransport();
418 }
419}
420
421
423{
424 if (regIOobject::read())
425 {
426 bool readOK = true;
427
428 for (phaseModel& phase : phaseModels_)
429 {
430 readOK &= phase.read();
431 }
432
433 // models ...
434
435 return readOK;
436 }
437
438 return false;
439}
440
441
443{
445 {
446 return fv::localEulerDdt::localRDeltaT(vf.mesh())*vf;
447 }
448 else
449 {
450 return vf/vf.mesh().time().deltaT();
451 }
452}
453
454
456{
458 {
459 return fv::localEulerDdt::localRDeltaTf(sf.mesh())*sf;
460 }
461 else
462 {
463 return sf/sf.mesh().time().deltaT();
464 }
465}
466
467
468// ************************************************************************* //
tmp< volScalarField > dmdt() const
Return the blended mass transfer rate.
wordList toc() const
Return the table of contents (as a sorted list)
const Mesh & mesh() const
Return mesh.
A HashTable similar to std::unordered_map.
Definition: HashTable.H:123
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition: HashTableI.H:180
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
writeOption writeOpt() const noexcept
The write option.
Definition: IOobjectI.H:179
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
const T * set(const label i) const
Definition: UPtrList.H:248
void resize(const label newLen)
Change the size of the list.
Definition: UPtrListI.H:190
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
virtual bool enabled() const
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:460
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
constant condensation/saturation model.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:61
An ordered or unorder pair of phase names. Typically specified as follows.
Definition: phasePairKey.H:68
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition: phasePair.H:56
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:76
phaseModelPartialList stationaryPhaseModels_
Stationary phase models.
Definition: phaseSystem.H:138
virtual void correctKinematics()
Correct the kinematics.
Definition: phaseSystem.C:376
virtual void correctEnergyTransport()
Correct the energy transport e.g. alphat.
Definition: phaseSystem.C:413
static const word propertiesName
Default name of the phase properties dictionary.
Definition: phaseSystem.H:290
void generatePairs(const dictTable &modelDicts)
Generate pairs.
Definition: phaseSystem.C:76
phaseModelList phaseModels_
Phase models.
Definition: phaseSystem.H:132
phaseModelPartialList movingPhaseModels_
Moving phase models.
Definition: phaseSystem.H:135
const surfaceScalarField & phi() const
Return the mixture flux.
Definition: phaseSystemI.H:113
virtual void correctThermo()
Correct the thermodynamics.
Definition: phaseSystem.C:395
virtual void correct()
Correct the fluid properties other than those listed below.
Definition: phaseSystem.C:367
tmp< volScalarField > sigma(const phasePairKey &key) const
Return the surface tension coefficient for a pair.
Definition: phaseSystem.C:319
phaseModelPartialList multiComponentPhaseModels_
Multi-component phase models.
Definition: phaseSystem.H:144
phaseModelPartialList anisothermalPhaseModels_
Anisothermal phase models.
Definition: phaseSystem.H:141
tmp< volVectorField > U() const
Return the mixture velocity.
Definition: phaseSystem.C:270
tmp< surfaceScalarField > calcPhi(const phaseModelList &phaseModels) const
Calculate and return the mixture flux.
Definition: phaseSystem.C:52
aspectRatioModelTable aspectRatioModels_
Aspect ratio models.
Definition: phaseSystem.H:168
surfaceTensionModelTable surfaceTensionModels_
Surface tension models.
Definition: phaseSystem.H:165
virtual void correctTurbulence()
Correct the turbulence.
Definition: phaseSystem.C:404
void generatePairsAndSubModels(const word &modelName, HashTable< autoPtr< modelType >, phasePairKey, phasePairKey::hash > &models)
Generate pairs and sub-model tables.
virtual PtrList< volScalarField > dmdts() const
Return the mass transfer rates for each phase.
Definition: phaseSystem.C:355
surfaceScalarField phi_
Total volumetric flux.
Definition: phaseSystem.H:150
virtual ~phaseSystem()
Destructor.
Definition: phaseSystem.C:236
tmp< volScalarField > rho() const
Return the mixture density.
Definition: phaseSystem.C:242
blendingMethodTable blendingMethods_
Blending methods.
Definition: phaseSystem.H:159
tmp< volScalarField > E(const phasePairKey &key) const
Return the aspect-ratio for a pair.
Definition: phaseSystem.C:300
virtual void solve()
Solve for the phase fractions.
Definition: phaseSystem.C:363
virtual bool read()
Read base phaseProperties dictionary.
Definition: phaseSystem.C:422
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:57
void correct()
Correct the phase properties.
bool read(const dictionary &phaseDict)
Read base transportProperties dictionary.
Lookup type of boundary radiation properties.
Definition: lookup.H:66
static const dimensionSet dimSigma
Coefficient dimensions.
virtual bool read()
Read object.
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
U
Definition: pEqn.H:72
dynamicFvMesh & mesh
Calculate the first temporal derivative.
word timeName
Definition: getTimeIndex.H:3
label phasei
Definition: pEqn.H:27
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.
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:47
Namespace for OpenFOAM.
const dimensionSet dimPressure
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimless
Dimensionless.
tmp< volScalarField > byDt(const volScalarField &vf)
Definition: phaseSystem.C:442
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
const dimensionSet dimDensity
volScalarField & alpha
dictionary dict
tmp< volScalarField > trho
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition: stdFoam.H:278
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:381
static const char *const typeName
The type name used in ensight case files.