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 -------------------------------------------------------------------------------
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 "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  (
63  IOobject
64  (
65  "phaseProperties",
66  mesh.time().constant(),
67  mesh,
68  IOobject::MUST_READ_IF_MODIFIED,
69  IOobject::NO_WRITE
70  )
71  ),
72 
73  mesh_(mesh),
74 
75  phase1_
76  (
77  *this,
78  *this,
79  this->get<wordList>("phases")[0]
80  ),
81 
82  phase2_
83  (
84  *this,
85  *this,
86  this->get<wordList>("phases")[1]
87  ),
88 
89  phi_
90  (
91  IOobject
92  (
93  "phi",
94  mesh.time().timeName(),
95  mesh,
96  IOobject::NO_READ,
97  IOobject::AUTO_WRITE
98  ),
99  this->calcPhi()
100  ),
101 
102  dgdt_
103  (
104  IOobject
105  (
106  "dgdt",
107  mesh.time().timeName(),
108  mesh,
109  IOobject::READ_IF_PRESENT,
110  IOobject::AUTO_WRITE
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  (
152  new orderedPhasePair
153  (
154  phase1_,
155  phase2_,
156  g,
157  sigmaTable,
158  aspectRatioTable
159  )
160  );
161 
162  pair2In1_.reset
163  (
164  new orderedPhasePair
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 * * * * * * * * * * * * * * * //
279 {}
280 
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 
296 Foam::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 
305 {
306  return drag_->K();
307 }
308 
309 
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  const Time& runTime = mesh_.time();
355 
356  volScalarField& alpha1 = phase1_;
357  volScalarField& alpha2 = phase2_;
358 
359  const surfaceScalarField& phi1 = phase1_.phi();
360  const surfaceScalarField& phi2 = phase2_.phi();
361 
362  const dictionary& alphaControls = mesh_.solverDict
363  (
364  alpha1.name()
365  );
366 
367  label nAlphaSubCycles(alphaControls.get<label>("nAlphaSubCycles"));
368  label nAlphaCorr(alphaControls.get<label>("nAlphaCorr"));
369 
370  word alphaScheme("div(phi," + alpha1.name() + ')');
371  word alpharScheme("div(phir," + alpha1.name() + ')');
372 
374 
375 
376  surfaceScalarField phic("phic", phi_);
377  surfaceScalarField phir("phir", phi1 - phi2);
378 
379  tmp<surfaceScalarField> alpha1alpha2f;
380 
381  if (pPrimeByA_.valid())
382  {
383  alpha1alpha2f =
384  fvc::interpolate(max(alpha1, scalar(0)))
385  *fvc::interpolate(max(alpha2, scalar(0)));
386 
387  surfaceScalarField phiP
388  (
389  pPrimeByA_()*fvc::snGrad(alpha1, "bounded")*mesh_.magSf()
390  );
391 
392  phir += phiP;
393  }
394 
395  for (int acorr=0; acorr<nAlphaCorr; acorr++)
396  {
398  (
399  IOobject
400  (
401  "Sp",
402  runTime.timeName(),
403  mesh_
404  ),
405  mesh_,
406  dimensionedScalar("Sp", dgdt_.dimensions(), 0.0)
407  );
408 
410  (
411  IOobject
412  (
413  "Su",
414  runTime.timeName(),
415  mesh_
416  ),
417  // Divergence term is handled explicitly to be
418  // consistent with the explicit transport solution
419  fvc::div(phi_)*min(alpha1, scalar(1))
420  );
421 
422  forAll(dgdt_, celli)
423  {
424  if (dgdt_[celli] > 0.0)
425  {
426  Sp[celli] -= dgdt_[celli]/max(1.0 - alpha1[celli], 1e-4);
427  Su[celli] += dgdt_[celli]/max(1.0 - alpha1[celli], 1e-4);
428  }
429  else if (dgdt_[celli] < 0.0)
430  {
431  Sp[celli] += dgdt_[celli]/max(alpha1[celli], 1e-4);
432  }
433  }
434 
435  surfaceScalarField alphaPhic1
436  (
437  fvc::flux
438  (
439  phic,
440  alpha1,
441  alphaScheme
442  )
443  + fvc::flux
444  (
445  -fvc::flux(-phir, scalar(1) - alpha1, alpharScheme),
446  alpha1,
448  )
449  );
450 
451  phase1_.correctInflowOutflow(alphaPhic1);
452 
453  if (nAlphaSubCycles > 1)
454  {
455  for
456  (
457  subCycle<volScalarField> alphaSubCycle(alpha1, nAlphaSubCycles);
458  !(++alphaSubCycle).end();
459  )
460  {
461  surfaceScalarField alphaPhic10(alphaPhic1);
462 
464  (
465  geometricOneField(),
466  alpha1,
467  phi_,
468  alphaPhic10,
469  (alphaSubCycle.index()*Sp)(),
470  (Su - (alphaSubCycle.index() - 1)*Sp*alpha1)(),
471  UniformField<scalar>(phase1_.alphaMax()),
472  zeroField()
473  );
474 
475  if (alphaSubCycle.index() == 1)
476  {
477  phase1_.alphaPhi() = alphaPhic10;
478  }
479  else
480  {
481  phase1_.alphaPhi() += alphaPhic10;
482  }
483  }
484 
485  phase1_.alphaPhi() /= nAlphaSubCycles;
486  }
487  else
488  {
490  (
491  geometricOneField(),
492  alpha1,
493  phi_,
494  alphaPhic1,
495  Sp,
496  Su,
497  UniformField<scalar>(phase1_.alphaMax()),
498  zeroField()
499  );
500 
501  phase1_.alphaPhi() = alphaPhic1;
502  }
503 
504  if (pPrimeByA_.valid())
505  {
506  fvScalarMatrix alpha1Eqn
507  (
509  - fvm::laplacian(alpha1alpha2f()*pPrimeByA_(), alpha1, "bounded")
510  );
511 
512  alpha1Eqn.relax();
513  alpha1Eqn.solve();
514 
515  phase1_.alphaPhi() += alpha1Eqn.flux();
516  }
517 
518  phase1_.alphaRhoPhi() =
519  fvc::interpolate(phase1_.rho())*phase1_.alphaPhi();
520 
521  phase2_.alphaPhi() = phi_ - phase1_.alphaPhi();
522  phase2_.correctInflowOutflow(phase2_.alphaPhi());
523  phase2_.alphaRhoPhi() =
524  fvc::interpolate(phase2_.rho())*phase2_.alphaPhi();
525 
526  Info<< alpha1.name() << " volume fraction = "
527  << alpha1.weightedAverage(mesh_.V()).value()
528  << " Min(" << alpha1.name() << ") = " << min(alpha1).value()
529  << " Max(" << alpha1.name() << ") = " << max(alpha1).value()
530  << endl;
531 
532  // Ensure the phase-fractions are bounded
534 
535  alpha2 = scalar(1) - alpha1;
536  }
537 }
538 
539 
541 {
542  phase1_.correct();
543  phase2_.correct();
544 }
545 
546 
548 {
549  phase1_.turbulence().correct();
550  phase2_.turbulence().correct();
551 }
552 
553 
555 {
556  if (regIOobject::read())
557  {
558  bool readOK = true;
559 
560  readOK &= phase1_.read(*this);
561  readOK &= phase2_.read(*this);
562 
563  // models ...
564 
565  return readOK;
566  }
567 
568  return false;
569 }
570 
571 
573 {
574  return pair_->sigma();
575 }
576 
577 
578 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
virtual ~twoPhaseSystem()
Destructor.
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionaryI.H:53
phaseModel & phase1_
Phase model 1.
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
void clamp_range(const dimensioned< MinMax< Type >> &range)
Clamp field values (in-place) to the specified range.
dictionary dict
tmp< volScalarField > Kd() const
Return the drag coefficient.
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:37
tmp< volVectorField > F() const
Return the combined force (lift + wall-lubrication)
tmp< volVectorField > U() const
Return the mixture velocity.
zeroField Su
Definition: alphaSuSp.H:1
void correct()
Correct two-phase properties other than turbulence.
virtual bool read()
Read object.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
tmp< volScalarField > rho() const
Return the mixture density.
tmp< surfaceScalarField > Ff() const
Return the combined face-force (lift + wall-lubrication)
bool found(const Key &key) const
Same as contains()
Definition: HashTable.H:1354
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
tmp< volScalarField > sigma() const
Return the surface tension coefficient.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:42
word alpharScheme("div(phirb,alpha)")
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
Calculate the matrix for the laplacian of the field.
vtkSmartPointer< vtkFloatArray > zeroField(const word &name, const label size)
Create named field initialized to zero.
Calculate the snGrad of the given volField.
const volScalarField & alpha2
const dimensionSet dimless
Dimensionless.
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.
surfaceScalarField & phi1
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:478
DimensionedField< scalar, volMesh > Internal
The internal field type from which this GeometricField is derived.
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:441
tmp< volScalarField > D() const
Return the turbulent diffusivity.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:40
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
virtual void solve()
Solve for the phase fractions.
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.
word timeName
Definition: getTimeIndex.H:3
tmp< surfaceScalarField > Vmf() const
Return the face virtual mass coefficient.
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
phaseModel & phase2_
Phase model 2.
dynamicFvMesh & mesh
const dictionary & alphaControls
Definition: alphaControls.H:1
static autoPtr< blendingMethod > New(const word &modelName, const dictionary &dict, const wordList &phaseNames)
surfaceScalarField phir(fvc::flux(UdmModel.Udm()))
Calculate the curl of the given volField by constructing the Hodge-dual of the symmetric part of the ...
const Time & time() const noexcept
Return time registry.
Calculate the face-flux of the given field.
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:41
Calculate the matrix for the first temporal derivative.
tmp< volScalarField > Vm() const
Return the virtual mass coefficient.
twoPhaseSystem(const fvMesh &)
Construct from fvMesh.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
surfaceScalarField & phi2
bool read()
Read base phaseProperties dictionary.
Calculate the divergence of the given field.
const uniformDimensionedVectorField & g
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
tmp< volScalarField > Kh() const
Return the heat transfer coefficient.
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:201
tmp< surfaceScalarField > Kdf() const
Return the face drag coefficient.
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.
void correctTurbulence()
Correct two-phase turbulence.
Represents 0/1 range or concept. Used for tagged dispatch or clamping.
Definition: pTraits.H:65
surfaceScalarField phic(mixture.cAlpha() *mag(alphaPhic/mesh.magSf()))
label nAlphaSubCycles(alphaControls.get< label >("nAlphaSubCycles"))
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:41
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
void correctBoundaryConditions()
Correct boundary field.
messageStream Info
Information stream (stdout output on master, null elsewhere)
MULES: Multidimensional universal limiter for explicit solution.
dimensioned< Type > weightedAverage(const DimensionedField< scalar, GeoMesh > &weightField) const
Calculate and return weighted average.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
label nAlphaCorr(alphaControls.get< label >("nAlphaCorr"))
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition: phasePair.H:49
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:40
zeroField Sp
Definition: alphaSuSp.H:2
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
const volScalarField & alpha1