Loading...
Searching...
No Matches
twoPhaseSystem.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) 2013-2018 OpenFOAM Foundation
9 Copyright (C) 2020 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 "twoPhaseSystem.H"
31#include "BlendedInterfacialModel.H"
32#include "virtualMassModel.H"
33#include "heatTransferModel.H"
34#include "liftModel.H"
35#include "wallLubricationModel.H"
36#include "turbulentDispersionModel.H"
37#include "fvMatrix.H"
38#include "surfaceInterpolate.H"
39#include "MULES.H"
40#include "subCycle.H"
41#include "fvcDdt.H"
42#include "fvcDiv.H"
43#include "fvcSnGrad.H"
44#include "fvcFlux.H"
45#include "fvcCurl.H"
46#include "fvmDdt.H"
47#include "fvmLaplacian.H"
49#include "blendingMethod.H"
50#include "HashPtrTable.H"
51#include "UniformField.H"
52
53// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54
56(
57 const fvMesh& mesh,
58 const dimensionedVector& g
59)
60:
62 (
64 (
65 "phaseProperties",
66 mesh.time().constant(),
67 mesh,
70 )
71 ),
72
73 mesh_(mesh),
74
76 (
77 *this,
78 *this,
79 this->get<wordList>("phases")[0]
80 ),
81
83 (
84 *this,
85 *this,
86 this->get<wordList>("phases")[1]
87 ),
88
89 phi_
90 (
92 (
93 "phi",
94 mesh.time().timeName(),
95 mesh,
98 ),
99 this->calcPhi()
100 ),
101
102 dgdt_
103 (
105 (
106 "dgdt",
107 mesh.time().timeName(),
108 mesh,
111 ),
112 mesh,
114 )
115{
116 phase2_.volScalarField::operator=(scalar(1) - phase1_);
117
118
119 // Blending
120 forAllConstIter(dictionary, subDict("blending"), iter)
121 {
122 blendingMethods_.insert
123 (
124 iter().dict().dictName(),
126 (
127 iter().dict(),
128 this->get<wordList>("phases")
129 )
130 );
131 }
132
133
134 // Pairs
135
136 phasePair::scalarTable sigmaTable(lookup("sigma"));
137 phasePair::dictTable aspectRatioTable(lookup("aspectRatio"));
138
139 pair_.reset
140 (
141 new phasePair
142 (
143 phase1_,
144 phase2_,
145 g,
146 sigmaTable
147 )
148 );
149
150 pair1In2_.reset
151 (
153 (
154 phase1_,
155 phase2_,
156 g,
157 sigmaTable,
158 aspectRatioTable
159 )
160 );
161
162 pair2In1_.reset
163 (
165 (
166 phase2_,
167 phase1_,
168 g,
169 sigmaTable,
170 aspectRatioTable
171 )
172 );
173
174
175 // Models
176
177 drag_.reset
178 (
180 (
181 lookup("drag"),
182 (
183 blendingMethods_.found("drag")
184 ? blendingMethods_["drag"]
185 : blendingMethods_["default"]
186 ),
187 pair_,
188 pair1In2_,
189 pair2In1_,
190 false // Do not zero drag coefficient at fixed-flux BCs
191 )
192 );
193
194 virtualMass_.reset
195 (
197 (
198 lookup("virtualMass"),
199 (
200 blendingMethods_.found("virtualMass")
201 ? blendingMethods_["virtualMass"]
202 : blendingMethods_["default"]
203 ),
204 pair_,
205 pair1In2_,
206 pair2In1_
207 )
208 );
209
210 heatTransfer_.reset
211 (
213 (
214 lookup("heatTransfer"),
215 (
216 blendingMethods_.found("heatTransfer")
217 ? blendingMethods_["heatTransfer"]
218 : blendingMethods_["default"]
219 ),
220 pair_,
221 pair1In2_,
222 pair2In1_
223 )
224 );
225
226 lift_.reset
227 (
229 (
230 lookup("lift"),
231 (
232 blendingMethods_.found("lift")
233 ? blendingMethods_["lift"]
234 : blendingMethods_["default"]
235 ),
236 pair_,
237 pair1In2_,
238 pair2In1_
239 )
240 );
241
242 wallLubrication_.reset
243 (
245 (
246 lookup("wallLubrication"),
247 (
248 blendingMethods_.found("wallLubrication")
249 ? blendingMethods_["wallLubrication"]
250 : blendingMethods_["default"]
251 ),
252 pair_,
253 pair1In2_,
254 pair2In1_
255 )
256 );
257
258 turbulentDispersion_.reset
259 (
261 (
262 lookup("turbulentDispersion"),
263 (
264 blendingMethods_.found("turbulentDispersion")
265 ? blendingMethods_["turbulentDispersion"]
266 : blendingMethods_["default"]
267 ),
268 pair_,
269 pair1In2_,
270 pair2In1_
271 )
272 );
273}
274
275
276// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
281
282// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
285{
286 return phase1_*phase1_.thermo().rho() + phase2_*phase2_.thermo().rho();
287}
288
289
291{
292 return phase1_*phase1_.U() + phase2_*phase2_.U();
293}
294
295
296Foam::tmp<Foam::surfaceScalarField> Foam::twoPhaseSystem::calcPhi() const
297{
298 return
299 fvc::interpolate(phase1_)*phase1_.phi()
300 + fvc::interpolate(phase2_)*phase2_.phi();
301}
302
303
304Foam::tmp<Foam::volScalarField> Foam::twoPhaseSystem::Kd() const
305{
306 return drag_->K();
307}
308
309
310Foam::tmp<Foam::surfaceScalarField> Foam::twoPhaseSystem::Kdf() const
311{
312 return drag_->Kf();
313}
314
317{
318 return virtualMass_->K();
319}
320
323{
324 return virtualMass_->Kf();
325}
326
329{
330 return heatTransfer_->K();
331}
332
335{
336 return lift_->F<vector>() + wallLubrication_->F<vector>();
337}
338
341{
342 return lift_->Ff() + wallLubrication_->Ff();
343}
344
345
347{
348 return turbulentDispersion_->D();
349}
350
351
353{
354 volScalarField& alpha1 = phase1_;
355 volScalarField& alpha2 = phase2_;
356
357 const surfaceScalarField& phi1 = phase1_.phi();
358 const surfaceScalarField& phi2 = phase2_.phi();
359
360 const dictionary& alphaControls = mesh_.solverDict
361 (
362 alpha1.name()
363 );
364
365 label nAlphaSubCycles(alphaControls.get<label>("nAlphaSubCycles"));
366 label nAlphaCorr(alphaControls.get<label>("nAlphaCorr"));
367
368 word alphaScheme("div(phi," + alpha1.name() + ')');
369 word alpharScheme("div(phir," + alpha1.name() + ')');
370
371 alpha1.correctBoundaryConditions();
372
373
374 surfaceScalarField phic("phic", phi_);
375 surfaceScalarField phir("phir", phi1 - phi2);
376
377 tmp<surfaceScalarField> alpha1alpha2f;
378
379 if (pPrimeByA_.valid())
380 {
381 alpha1alpha2f =
382 fvc::interpolate(max(alpha1, scalar(0)))
383 *fvc::interpolate(max(alpha2, scalar(0)));
384
386 (
387 pPrimeByA_()*fvc::snGrad(alpha1, "bounded")*mesh_.magSf()
388 );
389
390 phir += phiP;
391 }
392
393 for (int acorr=0; acorr<nAlphaCorr; acorr++)
394 {
396 (
397 mesh_.newIOobject("Sp"),
398 mesh_,
399 dimensionedScalar(dgdt_.dimensions(), Zero)
400 );
401
403 (
404 mesh_.newIOobject("Su"),
405 // Divergence term is handled explicitly to be
406 // consistent with the explicit transport solution
407 fvc::div(phi_)*min(alpha1, scalar(1))
408 );
409
410 forAll(dgdt_, celli)
411 {
412 if (dgdt_[celli] > 0.0)
413 {
414 Sp[celli] -= dgdt_[celli]/max(1.0 - alpha1[celli], 1e-4);
415 Su[celli] += dgdt_[celli]/max(1.0 - alpha1[celli], 1e-4);
416 }
417 else if (dgdt_[celli] < 0.0)
418 {
419 Sp[celli] += dgdt_[celli]/max(alpha1[celli], 1e-4);
420 }
421 }
422
423 surfaceScalarField alphaPhic1
424 (
426 (
427 phic,
428 alpha1,
429 alphaScheme
430 )
431 + fvc::flux
432 (
433 -fvc::flux(-phir, scalar(1) - alpha1, alpharScheme),
434 alpha1,
436 )
437 );
438
439 phase1_.correctInflowOutflow(alphaPhic1);
440
441 if (nAlphaSubCycles > 1)
442 {
443 for
444 (
446 !(++alphaSubCycle).end();
447 )
448 {
449 surfaceScalarField alphaPhic10(alphaPhic1);
450
452 (
454 alpha1,
455 phi_,
456 alphaPhic10,
457 (alphaSubCycle.index()*Sp)(),
458 (Su - (alphaSubCycle.index() - 1)*Sp*alpha1)(),
459 UniformField<scalar>(phase1_.alphaMax()),
460 zeroField()
461 );
462
463 if (alphaSubCycle.index() == 1)
464 {
465 phase1_.alphaPhi() = alphaPhic10;
466 }
467 else
468 {
469 phase1_.alphaPhi() += alphaPhic10;
470 }
471 }
472
473 phase1_.alphaPhi() /= nAlphaSubCycles;
474 }
475 else
476 {
478 (
480 alpha1,
481 phi_,
482 alphaPhic1,
483 Sp,
484 Su,
485 UniformField<scalar>(phase1_.alphaMax()),
486 zeroField()
487 );
488
489 phase1_.alphaPhi() = alphaPhic1;
490 }
491
492 if (pPrimeByA_.valid())
493 {
494 fvScalarMatrix alpha1Eqn
495 (
497 - fvm::laplacian(alpha1alpha2f()*pPrimeByA_(), alpha1, "bounded")
498 );
499
500 alpha1Eqn.relax();
501 alpha1Eqn.solve();
502
503 phase1_.alphaPhi() += alpha1Eqn.flux();
504 }
505
506 phase1_.alphaRhoPhi() =
507 fvc::interpolate(phase1_.rho())*phase1_.alphaPhi();
508
509 phase2_.alphaPhi() = phi_ - phase1_.alphaPhi();
510 phase2_.correctInflowOutflow(phase2_.alphaPhi());
511 phase2_.alphaRhoPhi() =
512 fvc::interpolate(phase2_.rho())*phase2_.alphaPhi();
513
514 Info<< alpha1.name() << " volume fraction = "
515 << alpha1.weightedAverage(mesh_.V()).value()
516 << " Min(" << alpha1.name() << ") = " << min(alpha1).value()
517 << " Max(" << alpha1.name() << ") = " << max(alpha1).value()
518 << endl;
519
520 // Ensure the phase-fractions are bounded
521 alpha1.clamp_range(zero_one{});
522
523 alpha2 = scalar(1) - alpha1;
524 }
525}
526
527
529{
530 phase1_.correct();
531 phase2_.correct();
532}
533
534
536{
537 phase1_.turbulence().correct();
538 phase2_.turbulence().correct();
539}
540
541
543{
544 if (regIOobject::read())
545 {
546 bool readOK = true;
547
548 readOK &= phase1_.read(*this);
549 readOK &= phase2_.read(*this);
550
551 // models ...
552
553 return readOK;
554 }
555
556 return false;
557}
558
559
561{
562 return pair_->sigma();
563}
564
565
566// ************************************************************************* //
surfaceScalarField phic(mixture.cAlpha() *mag(alphaPhic/mesh.magSf()))
MULES: Multidimensional universal limiter for explicit solution.
const uniformDimensionedVectorField & g
const volScalarField & alpha1
surfaceScalarField & phi2
const volScalarField & alpha2
surfaceScalarField & phi1
dimensioned< Type > weightedAverage(const DimensionedField< scalar, GeoMesh > &weights, const label comm=UPstream::worldComm) const
Return the global weighted average.
DimensionedField< scalar, volMesh > Internal
void clamp_range(const dimensioned< MinMax< Type > > &range)
Clamp field values (in-place) to the specified range.
IOdictionary(const IOobject &io, const dictionary *fallback=nullptr)
Construct given an IOobject and optional fallback dictionary content.
@ NO_READ
Nothing to be read.
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
A class representing the concept of a uniform field which stores only the single value and providing ...
static autoPtr< blendingMethod > New(const word &modelName, const dictionary &dict, const wordList &phaseNames)
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition dictionary.C:441
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return an entry data stream. FatalIOError if not found, or not a stream.
Definition dictionary.C:367
word dictName() const
The local dictionary name (final part of scoped name).
Definition dictionaryI.H:53
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition phasePair.H:52
HashTable< scalar, phasePairKey, phasePairKey::hash > scalarTable
Scalar hash table.
Definition phasePair.H:65
HashTable< dictionary, phasePairKey, phasePairKey::hash > dictTable
Dictionary hash table.
Definition phasePair.H:59
virtual bool read()
Read object.
bool end() const
Return true if the number of sub-cycles has been reached.
Perform a subCycleTime on a field.
Definition subCycle.H:137
A class for managing temporary objects.
Definition tmp.H:75
tmp< surfaceScalarField > Ff() const
Return the combined face-force (lift + wall-lubrication).
tmp< volVectorField > F() const
Return the combined force (lift + wall-lubrication).
tmp< volScalarField > D() const
Return the turbulent diffusivity.
tmp< volScalarField > Kd() const
Return the drag coefficient.
void correct()
Correct two-phase properties other than turbulence.
virtual ~twoPhaseSystem()
Destructor.
phaseModel & phase1_
Phase model 1.
tmp< surfaceScalarField > Kdf() const
Return the face drag coefficient.
tmp< surfaceScalarField > Vmf() const
Return the face virtual mass coefficient.
tmp< volVectorField > U() const
Return the mixture velocity.
phaseModel & phase2_
Phase model 2.
void correctTurbulence()
Correct two-phase turbulence.
const fvMesh & mesh() const
Return the mesh.
tmp< volScalarField > sigma() const
Return the surface tension coefficient.
tmp< volScalarField > rho() const
Return the mixture density.
twoPhaseSystem(const fvMesh &)
Construct from fvMesh.
tmp< volScalarField > Kh() const
Return the heat transfer coefficient.
virtual void solve()
Solve for the phase fractions.
bool read()
Read base phaseProperties dictionary.
tmp< volScalarField > Vm() const
Return the virtual mass coefficient.
A class for handling words, derived from Foam::string.
Definition word.H:66
A class representing the concept of a field of 0 used to avoid unnecessary manipulations for objects ...
Definition zeroField.H:50
Represents 0/1 range or concept. Used for tagged dispatch or clamping.
Definition pTraits.H:51
surfaceScalarField phir(fvc::flux(UdmModel.Udm()))
word alpharScheme("div(phirb,alpha)")
Calculate the curl of the given volField by constructing the Hodge-dual of the symmetric part of the ...
Calculate the first temporal derivative.
Calculate the divergence of the given field.
Calculate the face-flux of the given field.
Calculate the snGrad of the given volField.
Calculate the matrix for the first temporal derivative.
Calculate the matrix for the laplacian of the field.
zeroField Su
Definition alphaSuSp.H:1
zeroField Sp
Definition alphaSuSp.H:2
word timeName
Definition getTimeIndex.H:3
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
Different types of constants.
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, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition fvcSnGrad.C:40
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcDiv.C:42
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition fvcDdt.C:40
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Definition fvcFlux.C:27
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition fvmDdt.C:41
List< word > wordList
List of word.
Definition fileName.H:60
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
const dimensionSet dimless
Dimensionless.
fvMatrix< scalar > fvScalarMatrix
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
void Sp(fvMatrix< typename Expr::value_type > &m, const Expr2 &mult, const Expr &expression)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
void Su(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
dictionary dict
volScalarField & e
const dictionary & alphaControls
label nAlphaSubCycles(alphaControls.get< label >("nAlphaSubCycles"))
label nAlphaCorr(alphaControls.get< label >("nAlphaCorr"))
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object.
Definition stdFoam.H:356