multiphaseMixtureThermo.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) 2013-2020 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 "unitConversion.H"
29 #include "Time.H"
30 #include "subCycle.H"
31 #include "MULES.H"
32 #include "fvcDiv.H"
33 #include "fvcGrad.H"
34 #include "fvcSnGrad.H"
35 #include "fvcFlux.H"
36 #include "fvcMeshPhi.H"
37 #include "surfaceInterpolate.H"
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43  defineTypeNameAndDebug(multiphaseMixtureThermo, 0);
44 }
45 
46 
47 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48 
49 void Foam::multiphaseMixtureThermo::calcAlphas()
50 {
51  scalar level = 0.0;
52  alphas_ == 0.0;
53 
54  forAllIter(PtrDictionary<phaseModel>, phases_, phase)
55  {
56  alphas_ += level*phase();
57  level += 1.0;
58  }
59 }
60 
61 
62 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
63 
65 (
66  const volVectorField& U,
67  const surfaceScalarField& phi
68 )
69 :
70  psiThermo(U.mesh(), word::null),
71  phases_(lookup("phases"), phaseModel::iNew(p_, T_)),
72 
73  mesh_(U.mesh()),
74  U_(U),
75  phi_(phi),
76 
77  rhoPhi_
78  (
79  IOobject
80  (
81  "rhoPhi",
82  mesh_.time().timeName(),
83  mesh_,
84  IOobject::NO_READ,
85  IOobject::NO_WRITE
86  ),
87  mesh_,
89  ),
90 
91  alphas_
92  (
93  IOobject
94  (
95  "alphas",
96  mesh_.time().timeName(),
97  mesh_,
98  IOobject::NO_READ,
99  IOobject::AUTO_WRITE
100  ),
101  mesh_,
103  ),
104 
105  sigmas_(lookup("sigmas")),
106  dimSigma_(1, 0, -2, 0, 0),
107  deltaN_
108  (
109  "deltaN",
110  1e-8/pow(average(mesh_.V()), 1.0/3.0)
111  )
112 {
113  calcAlphas();
114  alphas_.write();
115  correct();
116 }
117 
118 
119 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
120 
122 {
123  forAllIter(PtrDictionary<phaseModel>, phases_, phasei)
124  {
125  phasei().correct();
126  }
127 
128  PtrDictionary<phaseModel>::iterator phasei = phases_.begin();
129 
130  psi_ = phasei()*phasei().thermo().psi();
131  mu_ = phasei()*phasei().thermo().mu();
132  alpha_ = phasei()*phasei().thermo().alpha();
133 
134  for (++phasei; phasei != phases_.end(); ++phasei)
135  {
136  psi_ += phasei()*phasei().thermo().psi();
137  mu_ += phasei()*phasei().thermo().mu();
138  alpha_ += phasei()*phasei().thermo().alpha();
139  }
140 }
141 
142 
144 {
145  forAllIter(PtrDictionary<phaseModel>, phases_, phasei)
146  {
147  phasei().thermo().rho() += phasei().thermo().psi()*dp;
148  }
149 }
150 
151 
153 {
154  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
155 
156  word name = phasei().thermo().thermoName();
157 
158  for (++ phasei; phasei != phases_.end(); ++ phasei)
159  {
160  name += ',' + phasei().thermo().thermoName();
161  }
162 
163  return name;
164 }
165 
166 
168 {
169  bool ico = true;
170 
171  forAllConstIter(PtrDictionary<phaseModel>, phases_, phase)
172  {
173  ico &= phase().thermo().incompressible();
174  }
175 
176  return ico;
177 }
178 
179 
181 {
182  bool iso = true;
183 
184  forAllConstIter(PtrDictionary<phaseModel>, phases_, phase)
185  {
186  iso &= phase().thermo().incompressible();
187  }
188 
189  return iso;
190 }
191 
192 
194 (
195  const volScalarField& p,
196  const volScalarField& T
197 ) const
198 {
199  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
200 
201  tmp<volScalarField> the(phasei()*phasei().thermo().he(p, T));
202 
203  for (++phasei; phasei != phases_.end(); ++phasei)
204  {
205  the.ref() += phasei()*phasei().thermo().he(p, T);
206  }
207 
208  return the;
209 }
210 
211 
213 (
214  const scalarField& T,
215  const labelList& cells
216 ) const
217 {
218  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
219 
220  tmp<scalarField> the
221  (
222  scalarField(phasei(), cells)*phasei().thermo().he(T, cells)
223  );
224 
225  for (++phasei; phasei != phases_.end(); ++phasei)
226  {
227  the.ref() +=
228  scalarField(phasei(), cells)*phasei().thermo().he(T, cells);
229  }
230 
231  return the;
232 }
233 
234 
236 (
237  const scalarField& T,
238  const label patchi
239 ) const
240 {
241  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
242 
243  tmp<scalarField> the
244  (
245  phasei().boundaryField()[patchi]*phasei().thermo().he(T, patchi)
246  );
247 
248  for (++phasei; phasei != phases_.end(); ++phasei)
249  {
250  the.ref() +=
251  phasei().boundaryField()[patchi]*phasei().thermo().he(T, patchi);
252  }
253 
254  return the;
255 }
256 
257 
259 {
260  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
261 
262  tmp<volScalarField> ths(phasei()*phasei().thermo().hs());
263 
264  for (++phasei; phasei != phases_.end(); ++phasei)
265  {
266  ths.ref() += phasei()*phasei().thermo().hs();
267  }
268 
269  return ths;
270 }
271 
272 
274 (
275  const volScalarField& p,
276  const volScalarField& T
277 ) const
278 {
279  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
280 
281  tmp<volScalarField> ths(phasei()*phasei().thermo().hs(p, T));
282 
283  for (++phasei; phasei != phases_.end(); ++phasei)
284  {
285  ths.ref() += phasei()*phasei().thermo().hs(p, T);
286  }
287 
288  return ths;
289 }
290 
291 
293 (
294  const scalarField& T,
295  const labelList& cells
296 ) const
297 {
298  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
299 
300  tmp<scalarField> ths
301  (
302  scalarField(phasei(), cells)*phasei().thermo().hs(T, cells)
303  );
304 
305  for (++phasei; phasei != phases_.end(); ++phasei)
306  {
307  ths.ref() +=
308  scalarField(phasei(), cells)*phasei().thermo().hs(T, cells);
309  }
310 
311  return ths;
312 }
313 
314 
316 (
317  const scalarField& T,
318  const label patchi
319 ) const
320 {
321  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
322 
323  tmp<scalarField> ths
324  (
325  phasei().boundaryField()[patchi]*phasei().thermo().hs(T, patchi)
326  );
327 
328  for (++phasei; phasei != phases_.end(); ++phasei)
329  {
330  ths.ref() +=
331  phasei().boundaryField()[patchi]*phasei().thermo().hs(T, patchi);
332  }
333 
334  return ths;
335 }
336 
337 
339 {
340  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
341 
342  tmp<volScalarField> tha(phasei()*phasei().thermo().ha());
343 
344  for (++phasei; phasei != phases_.end(); ++phasei)
345  {
346  tha.ref() += phasei()*phasei().thermo().ha();
347  }
348 
349  return tha;
350 }
351 
352 
354 (
355  const volScalarField& p,
356  const volScalarField& T
357 ) const
358 {
359  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
360 
361  tmp<volScalarField> tha(phasei()*phasei().thermo().ha(p, T));
362 
363  for (++phasei; phasei != phases_.end(); ++phasei)
364  {
365  tha.ref() += phasei()*phasei().thermo().ha(p, T);
366  }
367 
368  return tha;
369 }
370 
371 
373 (
374  const scalarField& T,
375  const labelList& cells
376 ) const
377 {
378  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
379 
380  tmp<scalarField> tha
381  (
382  scalarField(phasei(), cells)*phasei().thermo().ha(T, cells)
383  );
384 
385  for (++phasei; phasei != phases_.end(); ++phasei)
386  {
387  tha.ref() +=
388  scalarField(phasei(), cells)*phasei().thermo().ha(T, cells);
389  }
390 
391  return tha;
392 }
393 
394 
396 (
397  const scalarField& T,
398  const label patchi
399 ) const
400 {
401  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
402 
403  tmp<scalarField> tha
404  (
405  phasei().boundaryField()[patchi]*phasei().thermo().ha(T, patchi)
406  );
407 
408  for (++phasei; phasei != phases_.end(); ++phasei)
409  {
410  tha.ref() +=
411  phasei().boundaryField()[patchi]*phasei().thermo().ha(T, patchi);
412  }
413 
414  return tha;
415 }
416 
417 
419 {
420  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
421 
422  tmp<volScalarField> thc(phasei()*phasei().thermo().hc());
423 
424  for (++phasei; phasei != phases_.end(); ++phasei)
425  {
426  thc.ref() += phasei()*phasei().thermo().hc();
427  }
428 
429  return thc;
430 }
431 
432 
434 (
435  const scalarField& h,
436  const scalarField& T0,
437  const labelList& cells
438 ) const
439 {
441  return T0;
442 }
443 
444 
446 (
447  const scalarField& h,
448  const scalarField& T0,
449  const label patchi
450 ) const
451 {
453  return T0;
454 }
455 
456 
458 {
459  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
460 
461  tmp<volScalarField> trho(phasei()*phasei().thermo().rho());
462 
463  for (++phasei; phasei != phases_.end(); ++phasei)
464  {
465  trho.ref() += phasei()*phasei().thermo().rho();
466  }
467 
468  return trho;
469 }
470 
471 
473 (
474  const label patchi
475 ) const
476 {
477  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
478 
479  tmp<scalarField> trho
480  (
481  phasei().boundaryField()[patchi]*phasei().thermo().rho(patchi)
482  );
483 
484  for (++phasei; phasei != phases_.end(); ++phasei)
485  {
486  trho.ref() +=
487  phasei().boundaryField()[patchi]*phasei().thermo().rho(patchi);
488  }
489 
490  return trho;
491 }
492 
493 
495 {
496  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
497 
498  tmp<volScalarField> tCp(phasei()*phasei().thermo().Cp());
499 
500  for (++phasei; phasei != phases_.end(); ++phasei)
501  {
502  tCp.ref() += phasei()*phasei().thermo().Cp();
503  }
504 
505  return tCp;
506 }
507 
508 
510 (
511  const scalarField& T,
512  const label patchi
513 ) const
514 {
515  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
516 
517  tmp<scalarField> tCp
518  (
519  phasei().boundaryField()[patchi]*phasei().thermo().Cp(T, patchi)
520  );
521 
522  for (++phasei; phasei != phases_.end(); ++phasei)
523  {
524  tCp.ref() +=
525  phasei().boundaryField()[patchi]*phasei().thermo().Cp(T, patchi);
526  }
527 
528  return tCp;
529 }
530 
531 
533 {
534  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
535 
536  tmp<volScalarField> tCv(phasei()*phasei().thermo().Cv());
537 
538  for (++phasei; phasei != phases_.end(); ++phasei)
539  {
540  tCv.ref() += phasei()*phasei().thermo().Cv();
541  }
542 
543  return tCv;
544 }
545 
546 
548 (
549  const scalarField& T,
550  const label patchi
551 ) const
552 {
553  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
554 
555  tmp<scalarField> tCv
556  (
557  phasei().boundaryField()[patchi]*phasei().thermo().Cv(T, patchi)
558  );
559 
560  for (++phasei; phasei != phases_.end(); ++phasei)
561  {
562  tCv.ref() +=
563  phasei().boundaryField()[patchi]*phasei().thermo().Cv(T, patchi);
564  }
565 
566  return tCv;
567 }
568 
569 
571 {
572  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
573 
574  tmp<volScalarField> tgamma(phasei()*phasei().thermo().gamma());
575 
576  for (++phasei; phasei != phases_.end(); ++phasei)
577  {
578  tgamma.ref() += phasei()*phasei().thermo().gamma();
579  }
580 
581  return tgamma;
582 }
583 
584 
586 (
587  const scalarField& T,
588  const label patchi
589 ) const
590 {
591  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
592 
593  tmp<scalarField> tgamma
594  (
595  phasei().boundaryField()[patchi]*phasei().thermo().gamma(T, patchi)
596  );
597 
598  for (++phasei; phasei != phases_.end(); ++phasei)
599  {
600  tgamma.ref() +=
601  phasei().boundaryField()[patchi]
602  *phasei().thermo().gamma(T, patchi);
603  }
604 
605  return tgamma;
606 }
607 
608 
610 {
611  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
612 
613  tmp<volScalarField> tCpv(phasei()*phasei().thermo().Cpv());
614 
615  for (++phasei; phasei != phases_.end(); ++phasei)
616  {
617  tCpv.ref() += phasei()*phasei().thermo().Cpv();
618  }
619 
620  return tCpv;
621 }
622 
623 
625 (
626  const scalarField& T,
627  const label patchi
628 ) const
629 {
630  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
631 
632  tmp<scalarField> tCpv
633  (
634  phasei().boundaryField()[patchi]*phasei().thermo().Cpv(T, patchi)
635  );
636 
637  for (++phasei; phasei != phases_.end(); ++phasei)
638  {
639  tCpv.ref() +=
640  phasei().boundaryField()[patchi]
641  *phasei().thermo().Cpv(T, patchi);
642  }
643 
644  return tCpv;
645 }
646 
647 
649 {
650  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
651 
652  tmp<volScalarField> tCpByCpv(phasei()*phasei().thermo().CpByCpv());
653 
654  for (++phasei; phasei != phases_.end(); ++phasei)
655  {
656  tCpByCpv.ref() += phasei()*phasei().thermo().CpByCpv();
657  }
658 
659  return tCpByCpv;
660 }
661 
662 
664 (
665  const scalarField& T,
666  const label patchi
667 ) const
668 {
669  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
670 
671  tmp<scalarField> tCpByCpv
672  (
673  phasei().boundaryField()[patchi]*phasei().thermo().CpByCpv(T, patchi)
674  );
675 
676  for (++phasei; phasei != phases_.end(); ++phasei)
677  {
678  tCpByCpv.ref() +=
679  phasei().boundaryField()[patchi]
680  *phasei().thermo().CpByCpv(T, patchi);
681  }
682 
683  return tCpByCpv;
684 }
685 
686 
688 {
689  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
690 
691  tmp<volScalarField> tW(phasei()*phasei().thermo().W());
692 
693  for (++phasei; phasei != phases_.end(); ++phasei)
694  {
695  tW.ref() += phasei()*phasei().thermo().W();
696  }
697 
698  return tW;
699 }
700 
701 
703 (
704  const label patchi
705 ) const
706 {
707  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
708 
709  tmp<scalarField> tW
710  (
711  phasei().boundaryField()[patchi]*phasei().thermo().W(patchi)
712  );
713 
714  for (++phasei; phasei != phases_.end(); ++phasei)
715  {
716  tW.ref() +=
717  phasei().boundaryField()[patchi]*phasei().thermo().W(patchi);
718  }
719 
720  return tW;
721 }
722 
723 
725 {
726  return mu()/rho();
727 }
728 
729 
731 (
732  const label patchi
733 ) const
734 {
735  return mu(patchi)/rho(patchi);
736 }
737 
738 
740 {
741  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
742 
743  tmp<volScalarField> tkappa(phasei()*phasei().thermo().kappa());
744 
745  for (++phasei; phasei != phases_.end(); ++phasei)
746  {
747  tkappa.ref() += phasei()*phasei().thermo().kappa();
748  }
749 
750  return tkappa;
751 }
752 
753 
755 (
756  const label patchi
757 ) const
758 {
759  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
760 
761  tmp<scalarField> tkappa
762  (
763  phasei().boundaryField()[patchi]*phasei().thermo().kappa(patchi)
764  );
765 
766  for (++phasei; phasei != phases_.end(); ++phasei)
767  {
768  tkappa.ref() +=
769  phasei().boundaryField()[patchi]*phasei().thermo().kappa(patchi);
770  }
771 
772  return tkappa;
773 }
774 
775 
777 {
778  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
779 
780  tmp<volScalarField> talphaEff(phasei()*phasei().thermo().alphahe());
781 
782  for (++phasei; phasei != phases_.end(); ++phasei)
783  {
784  talphaEff.ref() += phasei()*phasei().thermo().alphahe();
785  }
786 
787  return talphaEff;
788 }
789 
790 
792 (
793  const label patchi
794 ) const
795 {
796  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
797 
798  tmp<scalarField> talphaEff
799  (
800  phasei().boundaryField()[patchi]
801  *phasei().thermo().alphahe(patchi)
802  );
803 
804  for (++phasei; phasei != phases_.end(); ++phasei)
805  {
806  talphaEff.ref() +=
807  phasei().boundaryField()[patchi]
808  *phasei().thermo().alphahe(patchi);
809  }
810 
811  return talphaEff;
812 }
813 
814 
816 (
817  const volScalarField& alphat
818 ) const
819 {
820  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
821 
822  tmp<volScalarField> tkappaEff(phasei()*phasei().thermo().kappaEff(alphat));
823 
824  for (++phasei; phasei != phases_.end(); ++phasei)
825  {
826  tkappaEff.ref() += phasei()*phasei().thermo().kappaEff(alphat);
827  }
828 
829  return tkappaEff;
830 }
831 
832 
834 (
835  const scalarField& alphat,
836  const label patchi
837 ) const
838 {
839  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
840 
841  tmp<scalarField> tkappaEff
842  (
843  phasei().boundaryField()[patchi]
844  *phasei().thermo().kappaEff(alphat, patchi)
845  );
846 
847  for (++phasei; phasei != phases_.end(); ++phasei)
848  {
849  tkappaEff.ref() +=
850  phasei().boundaryField()[patchi]
851  *phasei().thermo().kappaEff(alphat, patchi);
852  }
853 
854  return tkappaEff;
855 }
856 
857 
859 (
860  const volScalarField& alphat
861 ) const
862 {
863  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
864 
865  tmp<volScalarField> talphaEff(phasei()*phasei().thermo().alphaEff(alphat));
866 
867  for (++phasei; phasei != phases_.end(); ++phasei)
868  {
869  talphaEff.ref() += phasei()*phasei().thermo().alphaEff(alphat);
870  }
871 
872  return talphaEff;
873 }
874 
875 
877 (
878  const scalarField& alphat,
879  const label patchi
880 ) const
881 {
882  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
883 
884  tmp<scalarField> talphaEff
885  (
886  phasei().boundaryField()[patchi]
887  *phasei().thermo().alphaEff(alphat, patchi)
888  );
889 
890  for (++phasei; phasei != phases_.end(); ++phasei)
891  {
892  talphaEff.ref() +=
893  phasei().boundaryField()[patchi]
894  *phasei().thermo().alphaEff(alphat, patchi);
895  }
896 
897  return talphaEff;
898 }
899 
900 
902 {
903  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
904 
905  tmp<volScalarField> trCv(phasei()/phasei().thermo().Cv());
906 
907  for (++phasei; phasei != phases_.end(); ++phasei)
908  {
909  trCv.ref() += phasei()/phasei().thermo().Cv();
910  }
911 
912  return trCv;
913 }
914 
915 
918 {
919  tmp<surfaceScalarField> tstf
920  (
922  (
923  "surfaceTensionForce",
924  mesh_,
925  dimensionedScalar(dimensionSet(1, -2, -2, 0, 0), 0)
926  )
927  );
928 
929  surfaceScalarField& stf = tstf.ref();
930 
931  forAllConstIter(PtrDictionary<phaseModel>, phases_, phase1)
932  {
933  const phaseModel& alpha1 = phase1();
934 
936  ++phase2;
937 
938  for (; phase2 != phases_.end(); ++phase2)
939  {
940  const phaseModel& alpha2 = phase2();
941 
942  sigmaTable::const_iterator sigma =
943  sigmas_.find(interfacePair(alpha1, alpha2));
944 
945  if (sigma == sigmas_.end())
946  {
948  << "Cannot find interface " << interfacePair(alpha1, alpha2)
949  << " in list of sigma values"
950  << exit(FatalError);
951  }
952 
953  stf += dimensionedScalar(dimSigma_, sigma())
954  *fvc::interpolate(K(alpha1, alpha2))*
955  (
956  fvc::interpolate(alpha2)*fvc::snGrad(alpha1)
957  - fvc::interpolate(alpha1)*fvc::snGrad(alpha2)
958  );
959  }
960  }
961 
962  return tstf;
963 }
964 
965 
967 {
968  const Time& runTime = mesh_.time();
969 
970  const dictionary& alphaControls = mesh_.solverDict("alpha");
971  label nAlphaSubCycles(alphaControls.lookup<label>("nAlphaSubCycles"));
972  scalar cAlpha(alphaControls.lookup<scalar>("cAlpha"));
973 
974  volScalarField& alpha = phases_.first();
975 
976  if (nAlphaSubCycles > 1)
977  {
978  surfaceScalarField rhoPhiSum(0.0*rhoPhi_);
979  dimensionedScalar totalDeltaT = runTime.deltaT();
980 
981  for
982  (
983  subCycle<volScalarField> alphaSubCycle(alpha, nAlphaSubCycles);
984  !(++alphaSubCycle).end();
985  )
986  {
987  solveAlphas(cAlpha);
988  rhoPhiSum += (runTime.deltaT()/totalDeltaT)*rhoPhi_;
989  }
990 
991  rhoPhi_ = rhoPhiSum;
992  }
993  else
994  {
995  solveAlphas(cAlpha);
996  }
997 }
998 
999 
1000 Foam::tmp<Foam::surfaceVectorField> Foam::multiphaseMixtureThermo::nHatfv
1001 (
1002  const volScalarField& alpha1,
1003  const volScalarField& alpha2
1004 ) const
1005 {
1006  /*
1007  // Cell gradient of alpha
1008  volVectorField gradAlpha =
1009  alpha2*fvc::grad(alpha1) - alpha1*fvc::grad(alpha2);
1010 
1011  // Interpolated face-gradient of alpha
1012  surfaceVectorField gradAlphaf = fvc::interpolate(gradAlpha);
1013  */
1014 
1015  surfaceVectorField gradAlphaf
1016  (
1018  - fvc::interpolate(alpha1)*fvc::interpolate(fvc::grad(alpha2))
1019  );
1020 
1021  // Face unit interface normal
1022  return gradAlphaf/(mag(gradAlphaf) + deltaN_);
1023 }
1024 
1025 
1026 Foam::tmp<Foam::surfaceScalarField> Foam::multiphaseMixtureThermo::nHatf
1027 (
1028  const volScalarField& alpha1,
1029  const volScalarField& alpha2
1030 ) const
1031 {
1032  // Face unit interface normal flux
1033  return nHatfv(alpha1, alpha2) & mesh_.Sf();
1034 }
1035 
1036 
1037 // Correction for the boundary condition on the unit normal nHat on
1038 // walls to produce the correct contact angle.
1039 
1040 // The dynamic contact angle is calculated from the component of the
1041 // velocity on the direction of the interface, parallel to the wall.
1042 
1043 void Foam::multiphaseMixtureThermo::correctContactAngle
1044 (
1045  const phaseModel& alpha1,
1046  const phaseModel& alpha2,
1047  surfaceVectorField::Boundary& nHatb
1048 ) const
1049 {
1050  const volScalarField::Boundary& gbf
1051  = alpha1.boundaryField();
1052 
1053  const fvBoundaryMesh& boundary = mesh_.boundary();
1054 
1055  forAll(boundary, patchi)
1056  {
1057  if (isA<alphaContactAngleFvPatchScalarField>(gbf[patchi]))
1058  {
1059  const alphaContactAngleFvPatchScalarField& acap =
1060  refCast<const alphaContactAngleFvPatchScalarField>(gbf[patchi]);
1061 
1062  vectorField& nHatPatch = nHatb[patchi];
1063 
1064  vectorField AfHatPatch
1065  (
1066  mesh_.Sf().boundaryField()[patchi]
1067  /mesh_.magSf().boundaryField()[patchi]
1068  );
1069 
1070  alphaContactAngleFvPatchScalarField::thetaPropsTable::
1071  const_iterator tp =
1072  acap.thetaProps().find(interfacePair(alpha1, alpha2));
1073 
1074  if (tp == acap.thetaProps().end())
1075  {
1077  << "Cannot find interface " << interfacePair(alpha1, alpha2)
1078  << "\n in table of theta properties for patch "
1079  << acap.patch().name()
1080  << exit(FatalError);
1081  }
1082 
1083  bool matched = (tp.key().first() == alpha1.name());
1084 
1085  scalar theta0 = degToRad(tp().theta0(matched));
1086  scalarField theta(boundary[patchi].size(), theta0);
1087 
1088  scalar uTheta = tp().uTheta();
1089 
1090  // Calculate the dynamic contact angle if required
1091  if (uTheta > small)
1092  {
1093  scalar thetaA = degToRad(tp().thetaA(matched));
1094  scalar thetaR = degToRad(tp().thetaR(matched));
1095 
1096  // Calculated the component of the velocity parallel to the wall
1097  vectorField Uwall
1098  (
1099  U_.boundaryField()[patchi].patchInternalField()
1100  - U_.boundaryField()[patchi]
1101  );
1102  Uwall -= (AfHatPatch & Uwall)*AfHatPatch;
1103 
1104  // Find the direction of the interface parallel to the wall
1105  vectorField nWall
1106  (
1107  nHatPatch - (AfHatPatch & nHatPatch)*AfHatPatch
1108  );
1109 
1110  // Normalise nWall
1111  nWall /= (mag(nWall) + small);
1112 
1113  // Calculate Uwall resolved normal to the interface parallel to
1114  // the interface
1115  scalarField uwall(nWall & Uwall);
1116 
1117  theta += (thetaA - thetaR)*tanh(uwall/uTheta);
1118  }
1119 
1120 
1121  // Reset nHatPatch to correspond to the contact angle
1122 
1123  scalarField a12(nHatPatch & AfHatPatch);
1124 
1125  scalarField b1(cos(theta));
1126 
1127  scalarField b2(nHatPatch.size());
1128 
1129  forAll(b2, facei)
1130  {
1131  b2[facei] = cos(acos(a12[facei]) - theta[facei]);
1132  }
1133 
1134  scalarField det(1.0 - a12*a12);
1135 
1136  scalarField a((b1 - a12*b2)/det);
1137  scalarField b((b2 - a12*b1)/det);
1138 
1139  nHatPatch = a*AfHatPatch + b*nHatPatch;
1140 
1141  nHatPatch /= (mag(nHatPatch) + deltaN_.value());
1142  }
1143  }
1144 }
1145 
1146 
1147 Foam::tmp<Foam::volScalarField> Foam::multiphaseMixtureThermo::K
1148 (
1149  const phaseModel& alpha1,
1150  const phaseModel& alpha2
1151 ) const
1152 {
1153  tmp<surfaceVectorField> tnHatfv = nHatfv(alpha1, alpha2);
1154 
1155  correctContactAngle(alpha1, alpha2, tnHatfv.ref().boundaryFieldRef());
1156 
1157  // Simple expression for curvature
1158  return -fvc::div(tnHatfv & mesh_.Sf());
1159 }
1160 
1161 
1164 {
1165  tmp<volScalarField> tnearInt
1166  (
1168  (
1169  "nearInterface",
1170  mesh_,
1172  )
1173  );
1174 
1175  forAllConstIter(PtrDictionary<phaseModel>, phases_, phase)
1176  {
1177  tnearInt.ref() =
1178  max(tnearInt(), pos0(phase() - 0.01)*pos0(0.99 - phase()));
1179  }
1180 
1181  return tnearInt;
1182 }
1183 
1184 
1185 void Foam::multiphaseMixtureThermo::solveAlphas
1186 (
1187  const scalar cAlpha
1188 )
1189 {
1190  static label nSolves=-1;
1191  nSolves++;
1192 
1193  word alphaScheme("div(phi,alpha)");
1194  word alpharScheme("div(phirb,alpha)");
1195 
1196  surfaceScalarField phic(mag(phi_/mesh_.magSf()));
1197  phic = min(cAlpha*phic, max(phic));
1198 
1199  PtrList<surfaceScalarField> alphaPhiCorrs(phases_.size());
1200  int phasei = 0;
1201 
1202  forAllIter(PtrDictionary<phaseModel>, phases_, phase)
1203  {
1204  phaseModel& alpha = phase();
1205 
1206  alphaPhiCorrs.set
1207  (
1208  phasei,
1209  new surfaceScalarField
1210  (
1211  phi_.name() + alpha.name(),
1212  fvc::flux
1213  (
1214  phi_,
1215  alpha,
1216  alphaScheme
1217  )
1218  )
1219  );
1220 
1221  surfaceScalarField& alphaPhiCorr = alphaPhiCorrs[phasei];
1222 
1223  forAllIter(PtrDictionary<phaseModel>, phases_, phase2)
1224  {
1225  phaseModel& alpha2 = phase2();
1226 
1227  if (&alpha2 == &alpha) continue;
1228 
1229  surfaceScalarField phir(phic*nHatf(alpha, alpha2));
1230 
1231  alphaPhiCorr += fvc::flux
1232  (
1233  -fvc::flux(-phir, alpha2, alpharScheme),
1234  alpha,
1235  alpharScheme
1236  );
1237  }
1238 
1239  MULES::limit
1240  (
1241  1.0/mesh_.time().deltaT().value(),
1242  geometricOneField(),
1243  alpha,
1244  phi_,
1245  alphaPhiCorr,
1246  zeroField(),
1247  zeroField(),
1248  oneField(),
1249  zeroField(),
1250  true
1251  );
1252 
1253  phasei++;
1254  }
1255 
1256  MULES::limitSum(alphaPhiCorrs);
1257 
1258  rhoPhi_ = dimensionedScalar(dimensionSet(1, 0, -1, 0, 0), 0);
1259 
1260  volScalarField sumAlpha
1261  (
1262  IOobject
1263  (
1264  "sumAlpha",
1265  mesh_.time().timeName(),
1266  mesh_
1267  ),
1268  mesh_,
1270  );
1271 
1272 
1274 
1275 
1276  phasei = 0;
1277 
1278  forAllIter(PtrDictionary<phaseModel>, phases_, phase)
1279  {
1280  phaseModel& alpha = phase();
1281 
1282  surfaceScalarField& alphaPhi = alphaPhiCorrs[phasei];
1283  alphaPhi += upwind<scalar>(mesh_, phi_).flux(alpha);
1284 
1286  (
1287  IOobject
1288  (
1289  "Sp",
1290  mesh_.time().timeName(),
1291  mesh_
1292  ),
1293  mesh_,
1294  dimensionedScalar(alpha.dgdt().dimensions(), 0)
1295  );
1296 
1298  (
1299  IOobject
1300  (
1301  "Su",
1302  mesh_.time().timeName(),
1303  mesh_
1304  ),
1305  // Divergence term is handled explicitly to be
1306  // consistent with the explicit transport solution
1307  divU*min(alpha, scalar(1))
1308  );
1309 
1310  {
1311  const scalarField& dgdt = alpha.dgdt();
1312 
1313  forAll(dgdt, celli)
1314  {
1315  if (dgdt[celli] < 0.0 && alpha[celli] > 0.0)
1316  {
1317  Sp[celli] += dgdt[celli]*alpha[celli];
1318  Su[celli] -= dgdt[celli]*alpha[celli];
1319  }
1320  else if (dgdt[celli] > 0.0 && alpha[celli] < 1.0)
1321  {
1322  Sp[celli] -= dgdt[celli]*(1.0 - alpha[celli]);
1323  }
1324  }
1325  }
1326 
1327  forAllConstIter(PtrDictionary<phaseModel>, phases_, phase2)
1328  {
1329  const phaseModel& alpha2 = phase2();
1330 
1331  if (&alpha2 == &alpha) continue;
1332 
1333  const scalarField& dgdt2 = alpha2.dgdt();
1334 
1335  forAll(dgdt2, celli)
1336  {
1337  if (dgdt2[celli] > 0.0 && alpha2[celli] < 1.0)
1338  {
1339  Sp[celli] -= dgdt2[celli]*(1.0 - alpha2[celli]);
1340  Su[celli] += dgdt2[celli]*alpha[celli];
1341  }
1342  else if (dgdt2[celli] < 0.0 && alpha2[celli] > 0.0)
1343  {
1344  Sp[celli] += dgdt2[celli]*alpha2[celli];
1345  }
1346  }
1347  }
1348 
1350  (
1351  geometricOneField(),
1352  alpha,
1353  alphaPhi,
1354  Sp,
1355  Su
1356  );
1357 
1358  rhoPhi_ += fvc::interpolate(alpha.thermo().rho())*alphaPhi;
1359 
1360  Info<< alpha.name() << " volume fraction, min, max = "
1361  << alpha.weightedAverage(mesh_.V()).value()
1362  << ' ' << min(alpha).value()
1363  << ' ' << max(alpha).value()
1364  << endl;
1365 
1366  sumAlpha += alpha;
1367 
1368  phasei++;
1369  }
1370 
1371  Info<< "Phase-sum volume fraction, min, max = "
1372  << sumAlpha.weightedAverage(mesh_.V()).value()
1373  << ' ' << min(sumAlpha).value()
1374  << ' ' << max(sumAlpha).value()
1375  << endl;
1376 
1377  calcAlphas();
1378 }
1379 
1380 
1381 // ************************************************************************* //
volScalarField & he
Definition: YEEqn.H:50
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)
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:434
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
const dimensionedScalar & kappa
Coulomb constant: default SI units: [N.m2/C2].
zeroField Su
Definition: alphaSuSp.H:1
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
label phasei
Definition: pEqn.H:27
error FatalError
virtual tmp< volScalarField > gamma() const
Gamma = Cp/Cv [].
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
Unit conversion functions.
const dimensionedScalar & sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
tmp< volScalarField > tCv
const volScalarField & Cv
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
virtual tmp< volScalarField > hc() const
Enthalpy of formation [J/kg].
static tmp< GeometricField< scalar, fvsPatchField, surfaceMesh > > New(const word &name, const Internal &, const PtrList< fvsPatchField< scalar >> &)
Return a temporary field constructed from name,.
virtual tmp< volScalarField > ha() const
Absolute enthalpy [J/kg].
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
tmp< volScalarField > trho
Calculate the snGrad of the given volField.
surfaceScalarField phir(fvc::flux(UdmModel.Udm()))
dimensionedScalar det(const dimensionedSphericalTensor &dt)
rhoReactionThermo & thermo
Definition: createFields.H:28
friend class iterator
Definition: LPtrList.H:84
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:58
T & first()
Return the first element of the list.
Definition: UListI.H:114
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
CGAL::Exact_predicates_exact_constructions_kernel K
scalar degToRad(const scalar deg)
Conversion from degrees to radians.
virtual tmp< volScalarField > Cp() const
Heat capacity at constant pressure [J/kg/K].
virtual tmp< volScalarField > CpByCpv() const
Heat capacity ratio [].
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
void correctRho(const volScalarField &dp)
Update densities for given pressure change.
void solve()
Solve for the mixture phase-fractions.
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
dimensionedScalar cos(const dimensionedScalar &ds)
void limitSum(UPtrList< scalarField > &phiPsiCorrs)
Definition: MULES.C:30
virtual tmp< volScalarField > alphahe() const
Thermal diffusivity for energy of mixture [kg/m/s].
virtual tmp< volScalarField > rho() const
Density [kg/m^3].
const dimensionedScalar & b
Wien displacement law constant: default SI units: [m K].
Definition: createFields.H:27
Calculate the gradient of the given field.
tmp< volScalarField > rCv() const
Return the phase-averaged reciprocal Cv.
A class for handling words, derived from string.
Definition: word.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Calculate the face-flux of the given field.
virtual tmp< scalarField > THE(const scalarField &h, const scalarField &T0, const labelList &cells) const
Temperature from enthalpy/internal energy for cell-set.
const Type & value() const
Return const reference to value.
surfaceScalarField alphaPhi(phi.name()+alpha1.name(), fvc::flux(phi, alpha1, alphaScheme))
word timeName
Definition: getTimeIndex.H:3
List< label > labelList
A List of labels.
Definition: labelList.H:56
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
faceListList boundary(nPatches)
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
Calculate the divergence of the given field.
dimensionedScalar pos0(const dimensionedScalar &ds)
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-8/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:68
const label nAlphaSubCycles(alphaControls.lookup< label >("nAlphaSubCycles"))
defineTypeNameAndDebug(combustionModel, 0)
virtual tmp< volScalarField > nu() const
Kinematic viscosity of mixture [m^2/s].
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
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.
Calculate the mesh motion flux and convert fluxes from absolute to relative and back.
Internal & ref()
Return a reference to the dimensioned internal field.
const word alphaScheme(mesh.divScheme(divAlphaName)[1].wordToken())
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
virtual tmp< volScalarField > kappaEff(const volScalarField &alphat) const
Effective thermal diffusivity of mixture [W/m/K].
multiphaseMixtureThermo(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
label patchi
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:188
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
const dimensionedScalar & mu
Atomic mass unit.
virtual tmp< volScalarField > W() const
Molecular weight [kg/kmol].
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const scalarList W(::W(thermo))
scalar Cp(const scalar p, const scalar T) const
Definition: EtoHthermo.H:2
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
word alpharScheme("div(phirb,alpha)")
virtual tmp< volScalarField > hs() const
Sensible enthalpy [J/kg].
tmp< volScalarField > nearInterface() const
Indicator of the proximity of the interface.
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
const dimensionedScalar & h
Planck constant.
Field< vector > vectorField
Specialisation of Field<T> for vector.
MULES: Multidimensional universal limiter for explicit solution.
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
virtual tmp< volScalarField > Cv() const
Heat capacity at constant volume [J/kg/K].
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
zeroField divU
Definition: alphaSuSp.H:3
A class for managing temporary objects.
Definition: PtrList.H:53
tmp< surfaceScalarField > surfaceTensionForce() const
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
virtual tmp< volScalarField > alphaEff(const volScalarField &alphat) const
Effective thermal diffusivity of mixture [W/m/K].
virtual tmp< volScalarField > Cpv() const
Heat capacity at constant pressure/volume [J/kg/K].
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:366
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
virtual tmp< volScalarField > kappa() const
Thermal diffusivity for temperature of mixture [W/m/K].
virtual bool isochoric() const
Return true if the equation of state is isochoric.
virtual word thermoName() const
Return the name of the thermo physics.
virtual volScalarField & he()
Enthalpy/Internal energy [J/kg].
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
Namespace for OpenFOAM.
scalar T0
Definition: createFields.H:22
virtual bool incompressible() const
Return true if the equation of state is incompressible.
virtual void correct()
Update properties.
zeroField Sp
Definition: alphaSuSp.H:2