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