multiphaseMixtureThermo.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) 2013-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 
27 #include "alphaContactAngleFvPatchScalarField.H"
28 #include "Time.H"
29 #include "subCycle.H"
30 #include "MULES.H"
31 #include "fvcDiv.H"
32 #include "fvcGrad.H"
33 #include "fvcSnGrad.H"
34 #include "fvcFlux.H"
35 #include "fvcMeshPhi.H"
36 #include "surfaceInterpolate.H"
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42  defineTypeNameAndDebug(multiphaseMixtureThermo, 0);
43 }
44 
45 
48 
49 
50 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 
53 {
54  scalar level = 0.0;
55  alphas_ == 0.0;
56 
57  forAllIter(PtrDictionary<phaseModel>, phases_, phase)
58  {
59  alphas_ += level*phase();
60  level += 1.0;
61  }
62 
64 }
65 
66 
67 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
68 
70 (
71  const volVectorField& U,
72  const surfaceScalarField& phi
73 )
74 :
75  psiThermo(U.mesh(), word::null),
76  phases_(lookup("phases"), phaseModel::iNew(p_, T_)),
77 
78  mesh_(U.mesh()),
79  U_(U),
80  phi_(phi),
81 
82  rhoPhi_
83  (
84  IOobject
85  (
86  "rhoPhi",
87  mesh_.time().timeName(),
88  mesh_,
89  IOobject::NO_READ,
90  IOobject::NO_WRITE
91  ),
92  mesh_,
93  dimensionedScalar("rhoPhi", dimMass/dimTime, 0.0)
94  ),
95 
96  alphas_
97  (
98  IOobject
99  (
100  "alphas",
101  mesh_.time().timeName(),
102  mesh_,
103  IOobject::NO_READ,
104  IOobject::AUTO_WRITE
105  ),
106  mesh_,
107  dimensionedScalar("alphas", dimless, 0.0),
108  zeroGradientFvPatchScalarField::typeName
109  ),
110 
111  sigmas_(lookup("sigmas")),
112  dimSigma_(1, 0, -2, 0, 0),
113  deltaN_
114  (
115  "deltaN",
116  1e-8/pow(average(mesh_.V()), 1.0/3.0)
117  )
118 {
119  calcAlphas();
120  alphas_.write();
121  correct();
122 }
123 
124 
125 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
126 
128 {
129  forAllIter(PtrDictionary<phaseModel>, phases_, phasei)
130  {
131  phasei().correct();
132  }
133 
134  PtrDictionary<phaseModel>::iterator phasei = phases_.begin();
135 
136  psi_ = phasei()*phasei().thermo().psi();
137  mu_ = phasei()*phasei().thermo().mu();
138  alpha_ = phasei()*phasei().thermo().alpha();
139 
140  for (++phasei; phasei != phases_.end(); ++phasei)
141  {
142  psi_ += phasei()*phasei().thermo().psi();
143  mu_ += phasei()*phasei().thermo().mu();
144  alpha_ += phasei()*phasei().thermo().alpha();
145  }
146 }
147 
148 
150 {
151  forAllIter(PtrDictionary<phaseModel>, phases_, phasei)
152  {
153  phasei().thermo().rho() += phasei().thermo().psi()*dp;
154  }
155 }
156 
157 
159 {
160  bool ico = true;
161 
162  forAllConstIter(PtrDictionary<phaseModel>, phases_, phase)
163  {
164  ico &= phase().thermo().incompressible();
165  }
166 
167  return ico;
168 }
169 
170 
172 {
173  bool iso = true;
174 
175  forAllConstIter(PtrDictionary<phaseModel>, phases_, phase)
176  {
177  iso &= phase().thermo().incompressible();
178  }
179 
180  return iso;
181 }
182 
183 
185 (
186  const volScalarField& p,
187  const volScalarField& T
188 ) const
189 {
190  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
191 
192  tmp<volScalarField> the(phasei()*phasei().thermo().he(p, T));
193 
194  for (++phasei; phasei != phases_.end(); ++phasei)
195  {
196  the() += phasei()*phasei().thermo().he(p, T);
197  }
198 
199  return the;
200 }
201 
202 
204 (
205  const scalarField& p,
206  const scalarField& T,
207  const labelList& cells
208 ) const
209 {
210  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
211 
212  tmp<scalarField> the
213  (
215  );
216 
217  for (++phasei; phasei != phases_.end(); ++phasei)
218  {
219  the() += scalarField(phasei(), cells)*phasei().thermo().he(p, T, cells);
220  }
221 
222  return the;
223 }
224 
225 
227 (
228  const scalarField& p,
229  const scalarField& T,
230  const label patchi
231 ) const
232 {
233  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
234 
235  tmp<scalarField> the
236  (
238  );
239 
240  for (++phasei; phasei != phases_.end(); ++phasei)
241  {
242  the() +=
243  phasei().boundaryField()[patchi]*phasei().thermo().he(p, T, patchi);
244  }
245 
246  return the;
247 }
248 
249 
251 {
252  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
253 
254  tmp<volScalarField> thc(phasei()*phasei().thermo().hc());
255 
256  for (++phasei; phasei != phases_.end(); ++phasei)
257  {
258  thc() += phasei()*phasei().thermo().hc();
259  }
260 
261  return thc;
262 }
263 
264 
266 (
267  const scalarField& h,
268  const scalarField& p,
269  const scalarField& T0,
270  const labelList& cells
271 ) const
272 {
274  return T0;
275 }
276 
277 
279 (
280  const scalarField& h,
281  const scalarField& p,
282  const scalarField& T0,
283  const label patchi
284 ) const
285 {
287  return T0;
288 }
289 
290 
292 {
293  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
294 
295  tmp<volScalarField> trho(phasei()*phasei().thermo().rho());
296 
297  for (++phasei; phasei != phases_.end(); ++phasei)
298  {
299  trho() += phasei()*phasei().thermo().rho();
300  }
301 
302  return trho;
303 }
304 
305 
307 (
308  const label patchi
309 ) const
310 {
311  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
312 
313  tmp<scalarField> trho
314  (
316  );
317 
318  for (++phasei; phasei != phases_.end(); ++phasei)
319  {
320  trho() +=
321  phasei().boundaryField()[patchi]*phasei().thermo().rho(patchi);
322  }
323 
324  return trho;
325 }
326 
327 
329 {
330  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
331 
332  tmp<volScalarField> tCp(phasei()*phasei().thermo().Cp());
333 
334  for (++phasei; phasei != phases_.end(); ++phasei)
335  {
336  tCp() += phasei()*phasei().thermo().Cp();
337  }
338 
339  return tCp;
340 }
341 
342 
344 (
345  const scalarField& p,
346  const scalarField& T,
347  const label patchi
348 ) const
349 {
350  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
351 
352  tmp<scalarField> tCp
353  (
355  );
356 
357  for (++phasei; phasei != phases_.end(); ++phasei)
358  {
359  tCp() +=
360  phasei().boundaryField()[patchi]*phasei().thermo().Cp(p, T, patchi);
361  }
362 
363  return tCp;
364 }
365 
366 
368 {
369  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
370 
371  tmp<volScalarField> tCv(phasei()*phasei().thermo().Cv());
372 
373  for (++phasei; phasei != phases_.end(); ++phasei)
374  {
375  tCv() += phasei()*phasei().thermo().Cv();
376  }
377 
378  return tCv;
379 }
380 
381 
383 (
384  const scalarField& p,
385  const scalarField& T,
386  const label patchi
387 ) const
388 {
389  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
390 
391  tmp<scalarField> tCv
392  (
394  );
395 
396  for (++phasei; phasei != phases_.end(); ++phasei)
397  {
398  tCv() +=
399  phasei().boundaryField()[patchi]*phasei().thermo().Cv(p, T, patchi);
400  }
401 
402  return tCv;
403 }
404 
405 
407 {
408  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
409 
410  tmp<volScalarField> tgamma(phasei()*phasei().thermo().gamma());
411 
412  for (++phasei; phasei != phases_.end(); ++phasei)
413  {
414  tgamma() += phasei()*phasei().thermo().gamma();
415  }
416 
417  return tgamma;
418 }
419 
420 
422 (
423  const scalarField& p,
424  const scalarField& T,
425  const label patchi
426 ) const
427 {
428  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
429 
430  tmp<scalarField> tgamma
431  (
432  phasei().boundaryField()[patchi]*phasei().thermo().gamma(p, T, patchi)
433  );
434 
435  for (++phasei; phasei != phases_.end(); ++phasei)
436  {
437  tgamma() +=
438  phasei().boundaryField()[patchi]
439  *phasei().thermo().gamma(p, T, patchi);
440  }
441 
442  return tgamma;
443 }
444 
445 
447 {
448  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
449 
450  tmp<volScalarField> tCpv(phasei()*phasei().thermo().Cpv());
451 
452  for (++phasei; phasei != phases_.end(); ++phasei)
453  {
454  tCpv() += phasei()*phasei().thermo().Cpv();
455  }
456 
457  return tCpv;
458 }
459 
460 
462 (
463  const scalarField& p,
464  const scalarField& T,
465  const label patchi
466 ) const
467 {
468  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
469 
470  tmp<scalarField> tCpv
471  (
472  phasei().boundaryField()[patchi]*phasei().thermo().Cpv(p, T, patchi)
473  );
474 
475  for (++phasei; phasei != phases_.end(); ++phasei)
476  {
477  tCpv() +=
478  phasei().boundaryField()[patchi]
479  *phasei().thermo().Cpv(p, T, patchi);
480  }
481 
482  return tCpv;
483 }
484 
485 
487 {
488  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
489 
490  tmp<volScalarField> tCpByCpv(phasei()*phasei().thermo().CpByCpv());
491 
492  for (++phasei; phasei != phases_.end(); ++phasei)
493  {
494  tCpByCpv() += phasei()*phasei().thermo().CpByCpv();
495  }
496 
497  return tCpByCpv;
498 }
499 
500 
502 (
503  const scalarField& p,
504  const scalarField& T,
505  const label patchi
506 ) const
507 {
508  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
509 
510  tmp<scalarField> tCpByCpv
511  (
512  phasei().boundaryField()[patchi]*phasei().thermo().CpByCpv(p, T, patchi)
513  );
514 
515  for (++phasei; phasei != phases_.end(); ++phasei)
516  {
517  tCpByCpv() +=
518  phasei().boundaryField()[patchi]
519  *phasei().thermo().CpByCpv(p, T, patchi);
520  }
521 
522  return tCpByCpv;
523 }
524 
525 
527 {
528  return mu()/rho();
529 }
530 
531 
533 (
534  const label patchi
535 ) const
536 {
537  return mu(patchi)/rho(patchi);
538 }
539 
540 
542 {
543  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
544 
545  tmp<volScalarField> tkappa(phasei()*phasei().thermo().kappa());
546 
547  for (++phasei; phasei != phases_.end(); ++phasei)
548  {
549  tkappa() += phasei()*phasei().thermo().kappa();
550  }
551 
552  return tkappa;
553 }
554 
555 
557 (
558  const label patchi
559 ) const
560 {
561  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
562 
563  tmp<scalarField> tkappa
564  (
566  );
567 
568  for (++phasei; phasei != phases_.end(); ++phasei)
569  {
570  tkappa() +=
571  phasei().boundaryField()[patchi]*phasei().thermo().kappa(patchi);
572  }
573 
574  return tkappa;
575 }
576 
577 
579 (
580  const volScalarField& alphat
581 ) const
582 {
583  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
584 
585  tmp<volScalarField> tkappaEff(phasei()*phasei().thermo().kappaEff(alphat));
586 
587  for (++phasei; phasei != phases_.end(); ++phasei)
588  {
589  tkappaEff() += phasei()*phasei().thermo().kappaEff(alphat);
590  }
591 
592  return tkappaEff;
593 }
594 
595 
597 (
598  const scalarField& alphat,
599  const label patchi
600 ) const
601 {
602  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
603 
604  tmp<scalarField> tkappaEff
605  (
607  *phasei().thermo().kappaEff(alphat, patchi)
608  );
609 
610  for (++phasei; phasei != phases_.end(); ++phasei)
611  {
612  tkappaEff() +=
613  phasei().boundaryField()[patchi]
614  *phasei().thermo().kappaEff(alphat, patchi);
615  }
616 
617  return tkappaEff;
618 }
619 
620 
622 (
623  const volScalarField& alphat
624 ) const
625 {
626  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
627 
628  tmp<volScalarField> talphaEff(phasei()*phasei().thermo().alphaEff(alphat));
629 
630  for (++phasei; phasei != phases_.end(); ++phasei)
631  {
632  talphaEff() += phasei()*phasei().thermo().alphaEff(alphat);
633  }
634 
635  return talphaEff;
636 }
637 
638 
640 (
641  const scalarField& alphat,
642  const label patchi
643 ) const
644 {
645  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
646 
647  tmp<scalarField> talphaEff
648  (
651  );
652 
653  for (++phasei; phasei != phases_.end(); ++phasei)
654  {
655  talphaEff() +=
656  phasei().boundaryField()[patchi]
657  *phasei().thermo().alphaEff(alphat, patchi);
658  }
659 
660  return talphaEff;
661 }
662 
663 
665 {
666  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
667 
668  tmp<volScalarField> trCv(phasei()/phasei().thermo().Cv());
669 
670  for (++phasei; phasei != phases_.end(); ++phasei)
671  {
672  trCv() += phasei()/phasei().thermo().Cv();
673  }
674 
675  return trCv;
676 }
677 
678 
681 {
682  tmp<surfaceScalarField> tstf
683  (
685  (
686  IOobject
687  (
688  "surfaceTensionForce",
689  mesh_.time().timeName(),
690  mesh_
691  ),
692  mesh_,
694  (
695  "surfaceTensionForce",
696  dimensionSet(1, -2, -2, 0, 0),
697  0.0
698  )
699  )
700  );
701 
702  surfaceScalarField& stf = tstf();
703 
704  forAllConstIter(PtrDictionary<phaseModel>, phases_, phase1)
705  {
706  const phaseModel& alpha1 = phase1();
707 
708  PtrDictionary<phaseModel>::const_iterator phase2 = phase1;
709  ++phase2;
710 
711  for (; phase2 != phases_.end(); ++phase2)
712  {
713  const phaseModel& alpha2 = phase2();
714 
715  sigmaTable::const_iterator sigma =
716  sigmas_.find(interfacePair(alpha1, alpha2));
717 
718  if (sigma == sigmas_.end())
719  {
721  << "Cannot find interface " << interfacePair(alpha1, alpha2)
722  << " in list of sigma values"
723  << exit(FatalError);
724  }
725 
726  stf += dimensionedScalar("sigma", dimSigma_, sigma())
728  (
731  );
732  }
733  }
734 
735  return tstf;
736 }
737 
738 
740 {
741  const Time& runTime = mesh_.time();
742 
743  const dictionary& alphaControls = mesh_.solverDict("alpha");
744  label nAlphaSubCycles(readLabel(alphaControls.lookup("nAlphaSubCycles")));
745  scalar cAlpha(readScalar(alphaControls.lookup("cAlpha")));
746 
747  volScalarField& alpha = phases_.first();
748 
749  if (nAlphaSubCycles > 1)
750  {
751  surfaceScalarField rhoPhiSum(0.0*rhoPhi_);
752  dimensionedScalar totalDeltaT = runTime.deltaT();
753 
754  for
755  (
756  subCycle<volScalarField> alphaSubCycle(alpha, nAlphaSubCycles);
757  !(++alphaSubCycle).end();
758  )
759  {
760  solveAlphas(cAlpha);
761  rhoPhiSum += (runTime.deltaT()/totalDeltaT)*rhoPhi_;
762  }
763 
764  rhoPhi_ = rhoPhiSum;
765  }
766  else
767  {
768  solveAlphas(cAlpha);
769  }
770 }
771 
772 
774 (
775  const volScalarField& alpha1,
776  const volScalarField& alpha2
777 ) const
778 {
779  /*
780  // Cell gradient of alpha
781  volVectorField gradAlpha =
782  alpha2*fvc::grad(alpha1) - alpha1*fvc::grad(alpha2);
783 
784  // Interpolated face-gradient of alpha
785  surfaceVectorField gradAlphaf = fvc::interpolate(gradAlpha);
786  */
787 
788  surfaceVectorField gradAlphaf
789  (
792  );
793 
794  // Face unit interface normal
795  return gradAlphaf/(mag(gradAlphaf) + deltaN_);
796 }
797 
798 
800 (
801  const volScalarField& alpha1,
802  const volScalarField& alpha2
803 ) const
804 {
805  // Face unit interface normal flux
806  return nHatfv(alpha1, alpha2) & mesh_.Sf();
807 }
808 
809 
810 // Correction for the boundary condition on the unit normal nHat on
811 // walls to produce the correct contact angle.
812 
813 // The dynamic contact angle is calculated from the component of the
814 // velocity on the direction of the interface, parallel to the wall.
815 
817 (
818  const phaseModel& alpha1,
819  const phaseModel& alpha2,
820  surfaceVectorField::GeometricBoundaryField& nHatb
821 ) const
822 {
823  const volScalarField::GeometricBoundaryField& gbf
824  = alpha1.boundaryField();
825 
826  const fvBoundaryMesh& boundary = mesh_.boundary();
827 
829  {
830  if (isA<alphaContactAngleFvPatchScalarField>(gbf[patchi]))
831  {
832  const alphaContactAngleFvPatchScalarField& acap =
833  refCast<const alphaContactAngleFvPatchScalarField>(gbf[patchi]);
834 
835  vectorField& nHatPatch = nHatb[patchi];
836 
837  vectorField AfHatPatch
838  (
839  mesh_.Sf().boundaryField()[patchi]
840  /mesh_.magSf().boundaryField()[patchi]
841  );
842 
843  alphaContactAngleFvPatchScalarField::thetaPropsTable::
844  const_iterator tp =
845  acap.thetaProps().find(interfacePair(alpha1, alpha2));
846 
847  if (tp == acap.thetaProps().end())
848  {
850  << "Cannot find interface " << interfacePair(alpha1, alpha2)
851  << "\n in table of theta properties for patch "
852  << acap.patch().name()
853  << exit(FatalError);
854  }
855 
856  bool matched = (tp.key().first() == alpha1.name());
857 
858  scalar theta0 = convertToRad*tp().theta0(matched);
859  scalarField theta(boundary[patchi].size(), theta0);
860 
861  scalar uTheta = tp().uTheta();
862 
863  // Calculate the dynamic contact angle if required
864  if (uTheta > SMALL)
865  {
866  scalar thetaA = convertToRad*tp().thetaA(matched);
867  scalar thetaR = convertToRad*tp().thetaR(matched);
868 
869  // Calculated the component of the velocity parallel to the wall
870  vectorField Uwall
871  (
872  U_.boundaryField()[patchi].patchInternalField()
873  - U_.boundaryField()[patchi]
874  );
875  Uwall -= (AfHatPatch & Uwall)*AfHatPatch;
876 
877  // Find the direction of the interface parallel to the wall
878  vectorField nWall
879  (
880  nHatPatch - (AfHatPatch & nHatPatch)*AfHatPatch
881  );
882 
883  // Normalise nWall
884  nWall /= (mag(nWall) + SMALL);
885 
886  // Calculate Uwall resolved normal to the interface parallel to
887  // the interface
888  scalarField uwall(nWall & Uwall);
889 
890  theta += (thetaA - thetaR)*tanh(uwall/uTheta);
891  }
892 
893 
894  // Reset nHatPatch to correspond to the contact angle
895 
896  scalarField a12(nHatPatch & AfHatPatch);
897 
898  scalarField b1(cos(theta));
899 
900  scalarField b2(nHatPatch.size());
901 
902  forAll(b2, facei)
903  {
904  b2[facei] = cos(acos(a12[facei]) - theta[facei]);
905  }
906 
907  scalarField det(1.0 - a12*a12);
908 
909  scalarField a((b1 - a12*b2)/det);
910  scalarField b((b2 - a12*b1)/det);
911 
912  nHatPatch = a*AfHatPatch + b*nHatPatch;
913 
914  nHatPatch /= (mag(nHatPatch) + deltaN_.value());
915  }
916  }
917 }
918 
919 
921 (
922  const phaseModel& alpha1,
923  const phaseModel& alpha2
924 ) const
925 {
926  tmp<surfaceVectorField> tnHatfv = nHatfv(alpha1, alpha2);
927 
928  correctContactAngle(alpha1, alpha2, tnHatfv().boundaryField());
929 
930  // Simple expression for curvature
931  return -fvc::div(tnHatfv & mesh_.Sf());
932 }
933 
934 
937 {
938  tmp<volScalarField> tnearInt
939  (
940  new volScalarField
941  (
942  IOobject
943  (
944  "nearInterface",
945  mesh_.time().timeName(),
946  mesh_
947  ),
948  mesh_,
949  dimensionedScalar("nearInterface", dimless, 0.0)
950  )
951  );
952 
953  forAllConstIter(PtrDictionary<phaseModel>, phases_, phase)
954  {
955  tnearInt() = max(tnearInt(), pos(phase() - 0.01)*pos(0.99 - phase()));
956  }
957 
958  return tnearInt;
959 }
960 
961 
963 (
964  const scalar cAlpha
965 )
966 {
967  static label nSolves=-1;
968  nSolves++;
969 
970  word alphaScheme("div(phi,alpha)");
971  word alpharScheme("div(phirb,alpha)");
972 
973  surfaceScalarField phic(mag(phi_/mesh_.magSf()));
974  phic = min(cAlpha*phic, max(phic));
975 
976  PtrList<surfaceScalarField> alphaPhiCorrs(phases_.size());
977  int phasei = 0;
978 
979  forAllIter(PtrDictionary<phaseModel>, phases_, phase)
980  {
981  phaseModel& alpha = phase();
982 
983  alphaPhiCorrs.set
984  (
985  phasei,
987  (
988  phi_.name() + alpha.name(),
989  fvc::flux
990  (
991  phi_,
992  alpha,
993  alphaScheme
994  )
995  )
996  );
997 
998  surfaceScalarField& alphaPhiCorr = alphaPhiCorrs[phasei];
999 
1000  forAllIter(PtrDictionary<phaseModel>, phases_, phase2)
1001  {
1002  phaseModel& alpha2 = phase2();
1003 
1004  if (&alpha2 == &alpha) continue;
1005 
1007 
1008  alphaPhiCorr += fvc::flux
1009  (
1011  alpha,
1012  alpharScheme
1013  );
1014  }
1015 
1016  MULES::limit
1017  (
1018  1.0/mesh_.time().deltaT().value(),
1019  geometricOneField(),
1020  alpha,
1021  phi_,
1022  alphaPhiCorr,
1023  zeroField(),
1024  zeroField(),
1025  1,
1026  0,
1027  true
1028  );
1029 
1030  phasei++;
1031  }
1032 
1033  MULES::limitSum(alphaPhiCorrs);
1034 
1035  rhoPhi_ = dimensionedScalar("0", dimensionSet(1, 0, -1, 0, 0), 0);
1036 
1037  volScalarField sumAlpha
1038  (
1039  IOobject
1040  (
1041  "sumAlpha",
1042  mesh_.time().timeName(),
1043  mesh_
1044  ),
1045  mesh_,
1046  dimensionedScalar("sumAlpha", dimless, 0)
1047  );
1048 
1049 
1051 
1052 
1053  phasei = 0;
1054 
1055  forAllIter(PtrDictionary<phaseModel>, phases_, phase)
1056  {
1057  phaseModel& alpha = phase();
1058 
1059  surfaceScalarField& alphaPhi = alphaPhiCorrs[phasei];
1060  alphaPhi += upwind<scalar>(mesh_, phi_).flux(alpha);
1061 
1063  (
1064  IOobject
1065  (
1066  "Sp",
1067  mesh_.time().timeName(),
1068  mesh_
1069  ),
1070  mesh_,
1071  dimensionedScalar("Sp", alpha.dgdt().dimensions(), 0.0)
1072  );
1073 
1075  (
1076  IOobject
1077  (
1078  "Su",
1079  mesh_.time().timeName(),
1080  mesh_
1081  ),
1082  // Divergence term is handled explicitly to be
1083  // consistent with the explicit transport solution
1084  divU*min(alpha, scalar(1))
1085  );
1086 
1087  {
1088  const scalarField& dgdt = alpha.dgdt();
1089 
1090  forAll(dgdt, celli)
1091  {
1092  if (dgdt[celli] < 0.0 && alpha[celli] > 0.0)
1093  {
1094  Sp[celli] += dgdt[celli]*alpha[celli];
1095  Su[celli] -= dgdt[celli]*alpha[celli];
1096  }
1097  else if (dgdt[celli] > 0.0 && alpha[celli] < 1.0)
1098  {
1099  Sp[celli] -= dgdt[celli]*(1.0 - alpha[celli]);
1100  }
1101  }
1102  }
1103 
1104  forAllConstIter(PtrDictionary<phaseModel>, phases_, phase2)
1105  {
1106  const phaseModel& alpha2 = phase2();
1107 
1108  if (&alpha2 == &alpha) continue;
1109 
1110  const scalarField& dgdt2 = alpha2.dgdt();
1111 
1112  forAll(dgdt2, celli)
1113  {
1114  if (dgdt2[celli] > 0.0 && alpha2[celli] < 1.0)
1115  {
1116  Sp[celli] -= dgdt2[celli]*(1.0 - alpha2[celli]);
1117  Su[celli] += dgdt2[celli]*alpha[celli];
1118  }
1119  else if (dgdt2[celli] < 0.0 && alpha2[celli] > 0.0)
1120  {
1121  Sp[celli] += dgdt2[celli]*alpha2[celli];
1122  }
1123  }
1124  }
1125 
1127  (
1128  geometricOneField(),
1129  alpha,
1130  alphaPhi,
1131  Sp,
1132  Su
1133  );
1134 
1135  rhoPhi_ += fvc::interpolate(alpha.thermo().rho())*alphaPhi;
1136 
1137  Info<< alpha.name() << " volume fraction, min, max = "
1138  << alpha.weightedAverage(mesh_.V()).value()
1139  << ' ' << min(alpha).value()
1140  << ' ' << max(alpha).value()
1141  << endl;
1142 
1143  sumAlpha += alpha;
1144 
1145  phasei++;
1146  }
1147 
1148  Info<< "Phase-sum volume fraction, min, max = "
1149  << sumAlpha.weightedAverage(mesh_.V()).value()
1150  << ' ' << min(sumAlpha).value()
1151  << ' ' << max(sumAlpha).value()
1152  << endl;
1153 
1154  calcAlphas();
1155 }
1156 
1157 
1158 // ************************************************************************* //
Foam::multiphaseMixtureThermo::nu
virtual tmp< volScalarField > nu() const
Kinematic viscosity of mixture [m^2/s].
Foam::multiphaseMixtureThermo::correctContactAngle
void correctContactAngle(const phaseModel &alpha1, const phaseModel &alpha2, surfaceVectorField::GeometricBoundaryField &nHatb) const
Foam::fvc::snGrad
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
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
p
p
Definition: pEqn.H:62
forAllIter
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:431
Foam::multiphaseMixtureThermo::CpByCpv
virtual tmp< volScalarField > CpByCpv() const
Heat capacity ratio [].
Foam::labelList
List< label > labelList
A List of labels.
Definition: labelList.H:56
Foam::constant::physicoChemical::mu
const dimensionedScalar mu
Atomic mass unit.
Definition: createFields.H:13
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
fvcMeshPhi.H
Calculate the mesh motion flux and convert fluxes from absolute to relative and back.
talphaEff
Info<< "Creating turbulence model\n"<< endl;tmp< volScalarField > talphaEff
Definition: setAlphaEff.H:2
subCycle.H
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:216
alphat
alphat
Definition: TEqn.H:2
Foam::multiphaseMixtureThermo::convertToRad
static const scalar convertToRad
Conversion factor for degrees into radians.
Definition: multiphaseMixtureThermo.H:147
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::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
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::multiphaseMixtureThermo::alphas_
volScalarField alphas_
Definition: multiphaseMixtureThermo.H:135
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::multiphaseMixtureThermo::Cpv
virtual tmp< volScalarField > Cpv() const
Heat capacity at constant pressure/volume [J/kg/K].
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:261
Foam::multiphaseMixtureThermo::isochoric
virtual bool isochoric() const
Return true if the equation of state is isochoric.
phir
surfaceScalarField phir(phic *interface.nHatf())
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::constant::electromagnetic::kappa
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
Foam::multiphaseMixtureThermo::nearInterface
tmp< volScalarField > nearInterface() const
Indicator of the proximity of the interface.
Foam::multiphaseMixtureThermo::surfaceTensionForce
tmp< surfaceScalarField > surfaceTensionForce() const
thermo
rhoThermo & thermo
Definition: setRegionFluidFields.H:3
U
U
Definition: pEqn.H:46
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
K
CGAL::Exact_predicates_exact_constructions_kernel K
Definition: CGALTriangulation3DKernel.H:56
Foam::multiphaseMixtureThermo::kappaEff
virtual tmp< volScalarField > kappaEff(const volScalarField &alphat) const
Effective thermal diffusivity of mixture [J/m/s/K].
NotImplemented
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:365
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
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
Foam::fvc::Su
tmp< GeometricField< Type, fvPatchField, volMesh > > Su(const GeometricField< Type, fvPatchField, volMesh > &su, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:44
nAlphaSubCycles
label nAlphaSubCycles(readLabel(alphaControls.lookup("nAlphaSubCycles")))
Foam::tanh
dimensionedScalar tanh(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:263
Foam::constant::universal::h
const dimensionedScalar h
Planck constant.
Definition: createFields.H:6
Foam::Info
messageStream Info
Foam::MULES::limitSum
void limitSum(UPtrList< scalarField > &phiPsiCorrs)
Definition: MULES.C:51
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:55
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::multiphaseMixtureThermo::Cv
virtual tmp< volScalarField > Cv() const
Heat capacity at constant volume [J/kg/K].
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Foam::multiphaseMixtureThermo::rho
virtual tmp< volScalarField > rho() const
Density [kg/m^3].
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:73
alpha2
alpha2
Definition: alphaEqn.H:112
Foam::multiphaseMixtureThermo::phases_
PtrDictionary< phaseModel > phases_
Dictionary of phases.
Definition: multiphaseMixtureThermo.H:127
Foam::FatalError
error FatalError
Foam::multiphaseMixtureThermo::he
virtual volScalarField & he()
Enthalpy/Internal energy [J/kg].
Definition: multiphaseMixtureThermo.H:253
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
boundaryField
cellIbMask *cellIbMaskExt *faceIbMask *cellIbMask boundaryField().evaluateCoupled()
surfaceInterpolate.H
Surface Interpolation.
Foam::twoPhaseMixtureThermo::correct
virtual void correct()
Update properties.
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
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
Foam::multiphaseMixtureThermo::calcAlphas
void calcAlphas()
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
phic
phic
Definition: alphaEqnsSubCycle.H:5
alphaEff
const volScalarField & alphaEff
Definition: setAlphaEff.H:93
alpha1
volScalarField & alpha1
Definition: createFields.H:15
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
multiphaseMixtureThermo.H
Foam::multiphaseMixtureThermo::hc
virtual tmp< volScalarField > hc() const
Chemical enthalpy [J/kg].
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 > &)
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
trho
tmp< volScalarField > trho
Definition: setRegionSolidFields.H:4
Foam::multiphaseMixtureThermo::rCv
tmp< volScalarField > rCv() const
Return the phase-averaged reciprocal Cv.
alpharScheme
word alpharScheme("div(phirb,alpha)")
he
volScalarField & he
Definition: YEEqn.H:56
phasei
label phasei
Definition: pEqn.H:37
Foam::multiphaseMixtureThermo::THE
virtual tmp< scalarField > THE(const scalarField &h, const scalarField &p, const scalarField &T0, const labelList &cells) const
Temperature from enthalpy/internal energy for cell-set.
Foam::readScalar
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:52
Foam::multiphaseMixtureThermo::K
tmp< volScalarField > K(const phaseModel &alpha1, const phaseModel &alpha2) const
divU
volScalarField divU(fvc::div(fvc::absolute(phi, U)))
Foam::acos
dimensionedScalar acos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:259
fvcGrad.H
Calculate the gradient of the given field.
Foam::multiphaseMixtureThermo::Cp
virtual tmp< volScalarField > Cp() const
Heat capacity at constant pressure [J/kg/K].
fvcFlux.H
Calculate the face-flux of the given field.
Foam::multiphaseMixtureThermo::alphaEff
virtual tmp< volScalarField > alphaEff(const volScalarField &alphat) const
Effective thermal diffusivity of mixture [J/m/s/K].
patchi
label patchi
Definition: getPatchFieldScalar.H:1
Foam::multiphaseMixtureThermo::multiphaseMixtureThermo
multiphaseMixtureThermo(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::readLabel
label readLabel(Istream &is)
Definition: label.H:64
Foam::multiphaseMixtureThermo::kappa
virtual tmp< volScalarField > kappa() const
Thermal diffusivity for temperature of mixture [J/m/s/K].
Foam::constant::physicoChemical::sigma
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
phase2
phaseModel & phase2
Definition: createFields.H:13
Foam::multiphaseMixtureThermo::gamma
virtual tmp< volScalarField > gamma() const
Gamma = Cp/Cv [].
Foam::constant::mathematical::pi
const scalar pi(M_PI)
Foam::GeometricField::DimensionedInternalField
DimensionedField< Type, GeoMesh > DimensionedInternalField
Definition: GeometricField.H:100
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::det
dimensionedScalar det(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:60
Foam::multiphaseMixtureThermo::solve
void solve()
Solve for the mixture phase-fractions.
Foam::multiphaseMixtureThermo::incompressible
virtual bool incompressible() const
Return true if the equation of state is incompressible.
Foam::multiphaseMixtureThermo::correct
virtual void correct()
Update properties.
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::fvc::Sp
tmp< GeometricField< Type, fvPatchField, volMesh > > Sp(const volScalarField &sp, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:67
Foam::surfaceVectorField
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
Definition: surfaceFieldsFwd.H:55
Foam::dimensioned::dimensions
const dimensionSet & dimensions() const
Return const reference to dimensions.
Definition: dimensionedType.C:248
MULES.H
MULES: Multidimensional universal limiter for explicit solution.
Foam::multiphaseMixtureThermo::nHatfv
tmp< surfaceVectorField > nHatfv(const volScalarField &alpha1, const volScalarField &alpha2) const
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::multiphaseMixtureThermo::solveAlphas
void solveAlphas(const scalar cAlpha)
Foam::fvc::absolute
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:187
Foam::multiphaseMixtureThermo::correctRho
void correctRho(const volScalarField &dp)
Update densities for given pressure change.
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::multiphaseMixtureThermo::nHatf
tmp< surfaceScalarField > nHatf(const volScalarField &alpha1, const volScalarField &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