multiphaseSystem.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "multiphaseSystem.H"
27 #include "alphaContactAngleFvPatchScalarField.H"
29 #include "Time.H"
30 #include "subCycle.H"
31 #include "MULES.H"
32 #include "surfaceInterpolate.H"
33 #include "fvcGrad.H"
34 #include "fvcSnGrad.H"
35 #include "fvcDiv.H"
36 #include "fvcFlux.H"
37 #include "fvcAverage.H"
38 
39 // * * * * * * * * * * * * * * * Static Member Data * * * * * * * * * * * * //
40 
41 const Foam::scalar Foam::multiphaseSystem::convertToRad =
43 
44 
45 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46 
47 void Foam::multiphaseSystem::calcAlphas()
48 {
49  scalar level = 0.0;
50  alphas_ == 0.0;
51 
52  forAllIter(PtrDictionary<phaseModel>, phases_, iter)
53  {
54  alphas_ += level*iter();
55  level += 1.0;
56  }
57 }
58 
59 
60 void Foam::multiphaseSystem::solveAlphas()
61 {
62  PtrList<surfaceScalarField> alphaPhiCorrs(phases_.size());
63  int phasei = 0;
64 
65  forAllIter(PtrDictionary<phaseModel>, phases_, iter)
66  {
67  phaseModel& phase = iter();
68  volScalarField& alpha1 = phase;
69 
70  alphaPhiCorrs.set
71  (
72  phasei,
74  (
75  "phi" + alpha1.name() + "Corr",
76  fvc::flux
77  (
78  phi_,
79  phase,
80  "div(phi," + alpha1.name() + ')'
81  )
82  )
83  );
84 
85  surfaceScalarField& alphaPhiCorr = alphaPhiCorrs[phasei];
86 
87  forAllIter(PtrDictionary<phaseModel>, phases_, iter2)
88  {
89  phaseModel& phase2 = iter2();
90  volScalarField& alpha2 = phase2;
91 
92  if (&phase2 == &phase) continue;
93 
94  surfaceScalarField phir(phase.phi() - phase2.phi());
95 
97  (
98  cAlphas_.find(interfacePair(phase, phase2))
99  );
100 
101  if (cAlpha != cAlphas_.end())
102  {
104  (
105  (mag(phi_) + mag(phir))/mesh_.magSf()
106  );
107 
108  phir += min(cAlpha()*phic, max(phic))*nHatf(phase, phase2);
109  }
110 
111  word phirScheme
112  (
113  "div(phir," + alpha2.name() + ',' + alpha1.name() + ')'
114  );
115 
116  alphaPhiCorr += fvc::flux
117  (
118  -fvc::flux(-phir, phase2, phirScheme),
119  phase,
120  phirScheme
121  );
122  }
123 
124  phase.correctInflowOutflow(alphaPhiCorr);
125 
127  (
128  1.0/mesh_.time().deltaT().value(),
129  geometricOneField(),
130  phase,
131  phi_,
132  alphaPhiCorr,
133  zeroField(),
134  zeroField(),
135  oneField(),
136  zeroField(),
137  true
138  );
139 
140  phasei++;
141  }
142 
143  MULES::limitSum(alphaPhiCorrs);
144 
145  volScalarField sumAlpha
146  (
147  IOobject
148  (
149  "sumAlpha",
150  mesh_.time().timeName(),
151  mesh_
152  ),
153  mesh_,
154  dimensionedScalar("sumAlpha", dimless, 0)
155  );
156 
157  phasei = 0;
158 
159  forAllIter(PtrDictionary<phaseModel>, phases_, iter)
160  {
161  phaseModel& phase = iter();
162 
163  surfaceScalarField& alphaPhi = alphaPhiCorrs[phasei];
164  alphaPhi += upwind<scalar>(mesh_, phi_).flux(phase);
165  phase.correctInflowOutflow(alphaPhi);
166 
168  (
169  geometricOneField(),
170  phase,
171  alphaPhi
172  );
173 
174  phase.alphaPhi() = alphaPhi;
175 
176  Info<< phase.name() << " volume fraction, min, max = "
177  << phase.weightedAverage(mesh_.V()).value()
178  << ' ' << min(phase).value()
179  << ' ' << max(phase).value()
180  << endl;
181 
182  sumAlpha += phase;
183 
184  phasei++;
185  }
186 
187  Info<< "Phase-sum volume fraction, min, max = "
188  << sumAlpha.weightedAverage(mesh_.V()).value()
189  << ' ' << min(sumAlpha).value()
190  << ' ' << max(sumAlpha).value()
191  << endl;
192 
193  // Correct the sum of the phase-fractions to avoid 'drift'
194  volScalarField sumCorr(1.0 - sumAlpha);
195  forAllIter(PtrDictionary<phaseModel>, phases_, iter)
196  {
197  phaseModel& phase = iter();
198  volScalarField& alpha = phase;
199  alpha += alpha*sumCorr;
200  }
201 
202  calcAlphas();
203 }
204 
205 
206 Foam::tmp<Foam::surfaceVectorField> Foam::multiphaseSystem::nHatfv
207 (
208  const volScalarField& alpha1,
209  const volScalarField& alpha2
210 ) const
211 {
212  /*
213  // Cell gradient of alpha
214  volVectorField gradAlpha =
215  alpha2*fvc::grad(alpha1) - alpha1*fvc::grad(alpha2);
216 
217  // Interpolated face-gradient of alpha
218  surfaceVectorField gradAlphaf = fvc::interpolate(gradAlpha);
219  */
220 
221  surfaceVectorField gradAlphaf
222  (
224  - fvc::interpolate(alpha1)*fvc::interpolate(fvc::grad(alpha2))
225  );
226 
227  // Face unit interface normal
228  return gradAlphaf/(mag(gradAlphaf) + deltaN_);
229 }
230 
231 
232 Foam::tmp<Foam::surfaceScalarField> Foam::multiphaseSystem::nHatf
233 (
234  const volScalarField& alpha1,
235  const volScalarField& alpha2
236 ) const
237 {
238  // Face unit interface normal flux
239  return nHatfv(alpha1, alpha2) & mesh_.Sf();
240 }
241 
242 
243 // Correction for the boundary condition on the unit normal nHat on
244 // walls to produce the correct contact angle.
245 
246 // The dynamic contact angle is calculated from the component of the
247 // velocity on the direction of the interface, parallel to the wall.
248 
249 void Foam::multiphaseSystem::correctContactAngle
250 (
251  const phaseModel& phase1,
252  const phaseModel& phase2,
253  surfaceVectorField::Boundary& nHatb
254 ) const
255 {
256  const volScalarField::Boundary& gbf
257  = phase1.boundaryField();
258 
259  const fvBoundaryMesh& boundary = mesh_.boundary();
260 
261  forAll(boundary, patchi)
262  {
263  if (isA<alphaContactAngleFvPatchScalarField>(gbf[patchi]))
264  {
265  const alphaContactAngleFvPatchScalarField& acap =
266  refCast<const alphaContactAngleFvPatchScalarField>(gbf[patchi]);
267 
268  vectorField& nHatPatch = nHatb[patchi];
269 
270  vectorField AfHatPatch
271  (
272  mesh_.Sf().boundaryField()[patchi]
273  /mesh_.magSf().boundaryField()[patchi]
274  );
275 
276  alphaContactAngleFvPatchScalarField::thetaPropsTable::
277  const_iterator tp =
278  acap.thetaProps().find(interfacePair(phase1, phase2));
279 
280  if (tp == acap.thetaProps().end())
281  {
283  << "Cannot find interface " << interfacePair(phase1, phase2)
284  << "\n in table of theta properties for patch "
285  << acap.patch().name()
286  << exit(FatalError);
287  }
288 
289  bool matched = (tp.key().first() == phase1.name());
290 
291  scalar theta0 = convertToRad*tp().theta0(matched);
292  scalarField theta(boundary[patchi].size(), theta0);
293 
294  scalar uTheta = tp().uTheta();
295 
296  // Calculate the dynamic contact angle if required
297  if (uTheta > small)
298  {
299  scalar thetaA = convertToRad*tp().thetaA(matched);
300  scalar thetaR = convertToRad*tp().thetaR(matched);
301 
302  // Calculated the component of the velocity parallel to the wall
303  vectorField Uwall
304  (
305  phase1.U().boundaryField()[patchi].patchInternalField()
306  - phase1.U().boundaryField()[patchi]
307  );
308  Uwall -= (AfHatPatch & Uwall)*AfHatPatch;
309 
310  // Find the direction of the interface parallel to the wall
311  vectorField nWall
312  (
313  nHatPatch - (AfHatPatch & nHatPatch)*AfHatPatch
314  );
315 
316  // Normalise nWall
317  nWall /= (mag(nWall) + small);
318 
319  // Calculate Uwall resolved normal to the interface parallel to
320  // the interface
321  scalarField uwall(nWall & Uwall);
322 
323  theta += (thetaA - thetaR)*tanh(uwall/uTheta);
324  }
325 
326 
327  // Reset nHatPatch to correspond to the contact angle
328 
329  scalarField a12(nHatPatch & AfHatPatch);
330 
331  scalarField b1(cos(theta));
332 
333  scalarField b2(nHatPatch.size());
334 
335  forAll(b2, facei)
336  {
337  b2[facei] = cos(acos(a12[facei]) - theta[facei]);
338  }
339 
340  scalarField det(1.0 - a12*a12);
341 
342  scalarField a((b1 - a12*b2)/det);
343  scalarField b((b2 - a12*b1)/det);
344 
345  nHatPatch = a*AfHatPatch + b*nHatPatch;
346 
347  nHatPatch /= (mag(nHatPatch) + deltaN_.value());
348  }
349  }
350 }
351 
352 
353 Foam::tmp<Foam::volScalarField> Foam::multiphaseSystem::K
354 (
355  const phaseModel& phase1,
356  const phaseModel& phase2
357 ) const
358 {
359  tmp<surfaceVectorField> tnHatfv = nHatfv(phase1, phase2);
360 
361  correctContactAngle(phase1, phase2, tnHatfv.ref().boundaryFieldRef());
362 
363  // Simple expression for curvature
364  return -fvc::div(tnHatfv & mesh_.Sf());
365 }
366 
367 
368 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
369 
371 (
372  const volVectorField& U,
373  const surfaceScalarField& phi
374 )
375 :
377  (
378  IOobject
379  (
380  "transportProperties",
381  U.time().constant(),
382  U.db(),
385  )
386  ),
387 
388  phases_(lookup("phases"), phaseModel::iNew(U.mesh())),
389 
390  mesh_(U.mesh()),
391  phi_(phi),
392 
393  alphas_
394  (
395  IOobject
396  (
397  "alphas",
398  mesh_.time().timeName(),
399  mesh_,
400  IOobject::NO_READ,
402  ),
403  mesh_,
404  dimensionedScalar("alphas", dimless, 0.0)
405  ),
406 
407  sigmas_(lookup("sigmas")),
408  dimSigma_(1, 0, -2, 0, 0),
409  cAlphas_(lookup("interfaceCompression")),
410  Cvms_(lookup("virtualMass")),
411  deltaN_
412  (
413  "deltaN",
414  1e-8/pow(average(mesh_.V()), 1.0/3.0)
415  )
416 {
417  calcAlphas();
418  alphas_.write();
419 
420  interfaceDictTable dragModelsDict(lookup("drag"));
421 
422  forAllConstIter(interfaceDictTable, dragModelsDict, iter)
423  {
424  dragModels_.insert
425  (
426  iter.key(),
428  (
429  iter(),
430  *phases_.lookup(iter.key().first()),
431  *phases_.lookup(iter.key().second())
432  ).ptr()
433  );
434  }
435 
436  forAllConstIter(PtrDictionary<phaseModel>, phases_, iter1)
437  {
438  const phaseModel& phase1 = iter1();
439 
440  forAllConstIter(PtrDictionary<phaseModel>, phases_, iter2)
441  {
442  const phaseModel& phase2 = iter2();
443 
444  if (&phase2 != &phase1)
445  {
446  scalarCoeffSymmTable::const_iterator sigma
447  (
448  sigmas_.find(interfacePair(phase1, phase2))
449  );
450 
451  if (sigma != sigmas_.end())
452  {
453  scalarCoeffSymmTable::const_iterator cAlpha
454  (
455  cAlphas_.find(interfacePair(phase1, phase2))
456  );
457 
458  if (cAlpha == cAlphas_.end())
459  {
461  << "Compression coefficient not specified for "
462  "phase pair ("
463  << phase1.name() << ' ' << phase2.name()
464  << ") for which a surface tension "
465  "coefficient is specified"
466  << endl;
467  }
468  }
469  }
470  }
471  }
472 }
473 
474 
475 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
476 
478 {
479  PtrDictionary<phaseModel>::const_iterator iter = phases_.begin();
480 
481  tmp<volScalarField> trho = iter()*iter().rho();
482  volScalarField& rho = trho.ref();
483 
484  for (++iter; iter != phases_.end(); ++iter)
485  {
486  rho += iter()*iter().rho();
487  }
488 
489  return trho;
490 }
491 
492 
494 Foam::multiphaseSystem::rho(const label patchi) const
495 {
496  PtrDictionary<phaseModel>::const_iterator iter = phases_.begin();
497 
498  tmp<scalarField> trho = iter().boundaryField()[patchi]*iter().rho().value();
499  scalarField& rho = trho.ref();
500 
501  for (++iter; iter != phases_.end(); ++iter)
502  {
503  rho += iter().boundaryField()[patchi]*iter().rho().value();
504  }
505 
506  return trho;
507 }
508 
509 
511 {
512  PtrDictionary<phaseModel>::const_iterator iter = phases_.begin();
513 
514  tmp<volScalarField> tmu = iter()*(iter().rho()*iter().nu());
515  volScalarField& mu = tmu.ref();
516 
517  for (++iter; iter != phases_.end(); ++iter)
518  {
519  mu += iter()*(iter().rho()*iter().nu());
520  }
521 
522  return tmu/rho();
523 }
524 
525 
527 Foam::multiphaseSystem::nu(const label patchi) const
528 {
529  PtrDictionary<phaseModel>::const_iterator iter = phases_.begin();
530 
531  tmp<scalarField> tmu =
532  iter().boundaryField()[patchi]
533  *(iter().rho().value()*iter().nu().value());
534  scalarField& mu = tmu.ref();
535 
536  for (++iter; iter != phases_.end(); ++iter)
537  {
538  mu +=
539  iter().boundaryField()[patchi]
540  *(iter().rho().value()*iter().nu().value());
541  }
542 
543  return tmu/rho(patchi);
544 }
545 
546 
548 (
549  const phaseModel& phase
550 ) const
551 {
552  tmp<volScalarField> tCvm
553  (
554  new volScalarField
555  (
556  IOobject
557  (
558  "Cvm",
559  mesh_.time().timeName(),
560  mesh_
561  ),
562  mesh_,
564  (
565  "Cvm",
566  dimensionSet(1, -3, 0, 0, 0),
567  0
568  )
569  )
570  );
571 
572  forAllConstIter(PtrDictionary<phaseModel>, phases_, iter)
573  {
574  const phaseModel& phase2 = iter();
575 
576  if (&phase2 != &phase)
577  {
578  scalarCoeffTable::const_iterator Cvm
579  (
580  Cvms_.find(interfacePair(phase, phase2))
581  );
582 
583  if (Cvm != Cvms_.end())
584  {
585  tCvm.ref() += Cvm()*phase2.rho()*phase2;
586  }
587  else
588  {
589  Cvm = Cvms_.find(interfacePair(phase2, phase));
590 
591  if (Cvm != Cvms_.end())
592  {
593  tCvm.ref() += Cvm()*phase.rho()*phase2;
594  }
595  }
596  }
597  }
598 
599  return tCvm;
600 }
601 
602 
604 (
605  const phaseModel& phase
606 ) const
607 {
608  tmp<volVectorField> tSvm
609  (
610  new volVectorField
611  (
612  IOobject
613  (
614  "Svm",
615  mesh_.time().timeName(),
616  mesh_
617  ),
618  mesh_,
620  (
621  "Svm",
622  dimensionSet(1, -2, -2, 0, 0),
623  Zero
624  )
625  )
626  );
627 
628  forAllConstIter(PtrDictionary<phaseModel>, phases_, iter)
629  {
630  const phaseModel& phase2 = iter();
631 
632  if (&phase2 != &phase)
633  {
634  scalarCoeffTable::const_iterator Cvm
635  (
636  Cvms_.find(interfacePair(phase, phase2))
637  );
638 
639  if (Cvm != Cvms_.end())
640  {
641  tSvm.ref() += Cvm()*phase2.rho()*phase2*phase2.DDtU();
642  }
643  else
644  {
645  Cvm = Cvms_.find(interfacePair(phase2, phase));
646 
647  if (Cvm != Cvms_.end())
648  {
649  tSvm.ref() += Cvm()*phase.rho()*phase2*phase2.DDtU();
650  }
651  }
652  }
653  }
654 
655  volVectorField::Boundary& SvmBf =
656  tSvm.ref().boundaryFieldRef();
657 
658  // Remove virtual mass at fixed-flux boundaries
659  forAll(phase.phi().boundaryField(), patchi)
660  {
661  if
662  (
663  isA<fixedValueFvsPatchScalarField>
664  (
665  phase.phi().boundaryField()[patchi]
666  )
667  )
668  {
669  SvmBf[patchi] = Zero;
670  }
671  }
672 
673  return tSvm;
674 }
675 
676 
679 {
680  autoPtr<dragCoeffFields> dragCoeffsPtr(new dragCoeffFields);
681 
682  forAllConstIter(dragModelTable, dragModels_, iter)
683  {
684  const dragModel& dm = *iter();
685 
686  volScalarField* Kptr =
687  (
688  max
689  (
690  // fvc::average(dm.phase1()*dm.phase2()),
691  // fvc::average(dm.phase1())*fvc::average(dm.phase2()),
692  dm.phase1()*dm.phase2(),
693  dm.residualPhaseFraction()
694  )
695  *dm.K
696  (
697  max
698  (
699  mag(dm.phase1().U() - dm.phase2().U()),
700  dm.residualSlip()
701  )
702  )
703  ).ptr();
704 
705  volScalarField::Boundary& Kbf = Kptr->boundaryFieldRef();
706 
707  // Remove drag at fixed-flux boundaries
708  forAll(dm.phase1().phi().boundaryField(), patchi)
709  {
710  if
711  (
712  isA<fixedValueFvsPatchScalarField>
713  (
714  dm.phase1().phi().boundaryField()[patchi]
715  )
716  )
717  {
718  Kbf[patchi] = 0.0;
719  }
720  }
721 
722  dragCoeffsPtr().insert(iter.key(), Kptr);
723  }
724 
725  return dragCoeffsPtr;
726 }
727 
728 
730 (
731  const phaseModel& phase,
732  const dragCoeffFields& dragCoeffs
733 ) const
734 {
735  tmp<volScalarField> tdragCoeff
736  (
737  new volScalarField
738  (
739  IOobject
740  (
741  "dragCoeff",
742  mesh_.time().timeName(),
743  mesh_
744  ),
745  mesh_,
747  (
748  "dragCoeff",
749  dimensionSet(1, -3, -1, 0, 0),
750  0
751  )
752  )
753  );
754 
755  dragModelTable::const_iterator dmIter = dragModels_.begin();
756  dragCoeffFields::const_iterator dcIter = dragCoeffs.begin();
757  for
758  (
759  ;
760  dmIter != dragModels_.end() && dcIter != dragCoeffs.end();
761  ++dmIter, ++dcIter
762  )
763  {
764  if
765  (
766  &phase == &dmIter()->phase1()
767  || &phase == &dmIter()->phase2()
768  )
769  {
770  tdragCoeff.ref() += *dcIter();
771  }
772  }
773 
774  return tdragCoeff;
775 }
776 
777 
779 (
780  const phaseModel& phase1
781 ) const
782 {
783  tmp<surfaceScalarField> tSurfaceTension
784  (
786  (
787  IOobject
788  (
789  "surfaceTension",
790  mesh_.time().timeName(),
791  mesh_
792  ),
793  mesh_,
795  (
796  "surfaceTension",
797  dimensionSet(1, -2, -2, 0, 0),
798  0
799  )
800  )
801  );
802 
803  forAllConstIter(PtrDictionary<phaseModel>, phases_, iter)
804  {
805  const phaseModel& phase2 = iter();
806 
807  if (&phase2 != &phase1)
808  {
809  scalarCoeffSymmTable::const_iterator sigma
810  (
811  sigmas_.find(interfacePair(phase1, phase2))
812  );
813 
814  if (sigma != sigmas_.end())
815  {
816  tSurfaceTension.ref() +=
817  dimensionedScalar("sigma", dimSigma_, sigma())
818  *fvc::interpolate(K(phase1, phase2))*
819  (
820  fvc::interpolate(phase2)*fvc::snGrad(phase1)
821  - fvc::interpolate(phase1)*fvc::snGrad(phase2)
822  );
823  }
824  }
825  }
826 
827  return tSurfaceTension;
828 }
829 
830 
833 {
834  tmp<volScalarField> tnearInt
835  (
836  new volScalarField
837  (
838  IOobject
839  (
840  "nearInterface",
841  mesh_.time().timeName(),
842  mesh_
843  ),
844  mesh_,
845  dimensionedScalar("nearInterface", dimless, 0.0)
846  )
847  );
848 
849  forAllConstIter(PtrDictionary<phaseModel>, phases_, iter)
850  {
851  tnearInt.ref() =
852  max(tnearInt(), pos0(iter() - 0.01)*pos0(0.99 - iter()));
853  }
854 
855  return tnearInt;
856 }
857 
858 
860 {
861  forAllIter(PtrDictionary<phaseModel>, phases_, iter)
862  {
863  iter().correct();
864  }
865 
866  const Time& runTime = mesh_.time();
867 
868  const dictionary& alphaControls = mesh_.solverDict("alpha");
869  label nAlphaSubCycles(readLabel(alphaControls.lookup("nAlphaSubCycles")));
870 
871  if (nAlphaSubCycles > 1)
872  {
873  dimensionedScalar totalDeltaT = runTime.deltaT();
874 
875  PtrList<volScalarField> alpha0s(phases_.size());
876  PtrList<surfaceScalarField> alphaPhiSums(phases_.size());
877 
878  int phasei = 0;
879  forAllIter(PtrDictionary<phaseModel>, phases_, iter)
880  {
881  phaseModel& phase = iter();
883 
884  alpha0s.set
885  (
886  phasei,
887  new volScalarField(alpha.oldTime())
888  );
889 
890  alphaPhiSums.set
891  (
892  phasei,
894  (
895  IOobject
896  (
897  "phiSum" + alpha.name(),
898  runTime.timeName(),
899  mesh_
900  ),
901  mesh_,
902  dimensionedScalar("0", dimensionSet(0, 3, -1, 0, 0), 0)
903  )
904  );
905 
906  phasei++;
907  }
908 
909  for
910  (
911  subCycleTime alphaSubCycle
912  (
913  const_cast<Time&>(runTime),
915  );
916  !(++alphaSubCycle).end();
917  )
918  {
919  solveAlphas();
920 
921  int phasei = 0;
922  forAllIter(PtrDictionary<phaseModel>, phases_, iter)
923  {
924  alphaPhiSums[phasei] += iter().alphaPhi()/nAlphaSubCycles;
925  phasei++;
926  }
927  }
928 
929  phasei = 0;
930  forAllIter(PtrDictionary<phaseModel>, phases_, iter)
931  {
932  phaseModel& phase = iter();
933  volScalarField& alpha = phase;
934 
935  phase.alphaPhi() = alphaPhiSums[phasei];
936 
937  // Correct the time index of the field
938  // to correspond to the global time
939  alpha.timeIndex() = runTime.timeIndex();
940 
941  // Reset the old-time field value
942  alpha.oldTime() = alpha0s[phasei];
943  alpha.oldTime().timeIndex() = runTime.timeIndex();
944 
945  phasei++;
946  }
947  }
948  else
949  {
950  solveAlphas();
951  }
952 }
953 
954 
956 {
957  if (regIOobject::read())
958  {
959  bool readOK = true;
960 
961  PtrList<entry> phaseData(lookup("phases"));
962  label phasei = 0;
963 
964  forAllIter(PtrDictionary<phaseModel>, phases_, iter)
965  {
966  readOK &= iter().read(phaseData[phasei++].dict());
967  }
968 
969  lookup("sigmas") >> sigmas_;
970  lookup("interfaceCompression") >> cAlphas_;
971  lookup("virtualMass") >> Cvms_;
972 
973  return readOK;
974  }
975  else
976  {
977  return false;
978  }
979 }
980 
981 
982 // ************************************************************************* //
void limit(const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phi, surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su, const PsiMaxType &psiMax, const PsiMinType &psiMin, const bool returnCorr)
dimensionedScalar tanh(const dimensionedScalar &ds)
dictionary dict
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
dimensionedScalar acos(const dimensionedScalar &ds)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
label phasei
Definition: pEqn.H:27
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:110
virtual bool read()
Read object.
error FatalError
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:453
IOobject(const word &name, const fileName &instance, const objectRegistry &registry, readOption r=NO_READ, writeOption w=NO_WRITE, bool registerObject=true)
Construct from name, instance, registry, io options.
Definition: IOobject.C:209
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
tmp< volScalarField > trho
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
autoPtr< multiphaseSystem::dragCoeffFields > dragCoeffs(fluid.dragCoeffs())
Calculate the snGrad of the given volField.
dimensionedScalar det(const dimensionedSphericalTensor &dt)
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
CGAL::Exact_predicates_exact_constructions_kernel K
tmp< volScalarField > dragCoeff(const phaseModel &phase, const dragCoeffFields &dragCoeffs) const
Return the sum of the drag coefficients for the given phase.
Area-weighted average a surfaceField creating a volField.
const surfaceScalarField & alphaPhi() const
Definition: phaseModel.H:200
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:142
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
tmp< surfaceScalarField > surfaceTension(const phaseModel &phase) const
stressControl lookup("compactNormalStress") >> compactNormalStress
tmp< volScalarField > nearInterface() const
Indicator of the proximity of the interface.
static autoPtr< dragModel > New(const dictionary &interfaceDict, const phaseModel &phase1, const phaseModel &phase2)
const phaseModel & phase() const
Return the phase.
dimensionedScalar cos(const dimensionedScalar &ds)
const dictionary & alphaControls
Definition: alphaControls.H:1
void limitSum(UPtrList< scalarField > &phiPsiCorrs)
Definition: MULES.C:30
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Calculate the gradient of the given field.
surfaceScalarField phir(fvc::flux(UdmModel.Udm()))
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Calculate the face-flux of the given field.
friend class const_iterator
Declare friendship with the const_iterator.
Definition: HashTable.H:197
tmp< volVectorField > U() const
Return the mixture velocity.
const Type & value() const
Return const reference to value.
tmp< volVectorField > Svm(const phaseModel &phase) const
Return the virtual-mass source for the given phase.
word timeName
Definition: getTimeIndex.H:3
static const zero Zero
Definition: zero.H:97
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
faceListList boundary(nPatches)
void solve()
Solve for the mixture phase-fractions.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
label readLabel(Istream &is)
Definition: label.H:64
phaseModel & phase1
tmp< volScalarField > rho() const
Return the mixture density.
Calculate the divergence of the given field.
dimensionedScalar pos0(const dimensionedScalar &ds)
tmp< volScalarField > nu() const
Return the mixture laminar viscosity.
bool read()
Read base transportProperties dictionary.
autoPtr< dragCoeffFields > dragCoeffs() const
Return the drag coefficients for all of the interfaces.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
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.
label timeIndex() const
Return the time index of the field.
const dimensionedScalar mu
Atomic mass unit.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
IOdictionary(const IOobject &)
Construct given an IOobject.
Definition: IOdictionary.C:30
label patchi
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
surfaceScalarField phic(mixture.cAlpha() *mag(phi/mesh.magSf()))
#define WarningInFunction
Report a warning using Foam::Warning.
const surfaceScalarField & phi() const
Return the mixture flux.
Definition: phaseSystemI.H:111
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const Time & time() const
Return time.
Definition: IOobject.C:367
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
Field< vector > vectorField
Specialisation of Field<T> for vector.
MULES: Multidimensional universal limiter for explicit solution.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:98
phaseModel & phase2
A class for managing temporary objects.
Definition: PtrList.H:53
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Definition: fvcFlux.C:32
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:361
surfaceScalarField alphaPhi(phi.name()+alpha1.name(), fvc::flux(phi, alpha1, alphaScheme))
multiphaseSystem(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
const fvMesh & mesh() const
Return the mesh.
Definition: phaseSystemI.H:28
tmp< volScalarField > Cvm(const phaseModel &phase) const
Return the virtual-mass coefficient for the given phase.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
friend class const_iterator
Definition: LPtrList.H:87
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576
label nAlphaSubCycles(readLabel(alphaControls.lookup("nAlphaSubCycles")))