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-2021 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  rhoThermo::composite(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 
129 
131 {
132  PtrDictionary<phaseModel>::iterator phasei = phases_.begin();
133 
134  rho_ = phasei()*phasei().thermo().rho();
135  psi_ = phasei()*phasei().thermo().psi();
136  mu_ = phasei()*phasei().thermo().mu();
137  alpha_ = phasei()*phasei().thermo().alpha();
138 
139  for (++phasei; phasei != phases_.end(); ++phasei)
140  {
141  rho_ += phasei()*phasei().thermo().rho();
142  psi_ += phasei()*phasei().thermo().psi();
143  mu_ += phasei()*phasei().thermo().mu();
144  alpha_ += phasei()*phasei().thermo().alpha();
145  }
146 
147  forAllIter(PtrDictionary<phaseModel>, phases_, phasei)
148  {
149  phasei().Alpha() = phasei()*phasei().thermo().rho()/rho_;
150  }
151 }
152 
153 
155 {
156  forAllIter(PtrDictionary<phaseModel>, phases_, phasei)
157  {
158  phasei().thermo().rho() += phasei().thermo().psi()*dp;
159  }
160 }
161 
162 
164 {
165  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
166 
167  word name = phasei().thermo().thermoName();
168 
169  for (++ phasei; phasei != phases_.end(); ++ phasei)
170  {
171  name += ',' + phasei().thermo().thermoName();
172  }
173 
174  return name;
175 }
176 
177 
179 {
180  bool ico = true;
181 
182  forAllConstIter(PtrDictionary<phaseModel>, phases_, phase)
183  {
184  ico &= phase().thermo().incompressible();
185  }
186 
187  return ico;
188 }
189 
190 
192 {
193  bool iso = true;
194 
195  forAllConstIter(PtrDictionary<phaseModel>, phases_, phase)
196  {
197  iso &= phase().thermo().incompressible();
198  }
199 
200  return iso;
201 }
202 
203 
205 (
206  const volScalarField& p,
207  const volScalarField& T
208 ) const
209 {
210  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
211 
212  tmp<volScalarField> the(phasei().Alpha()*phasei().thermo().he(p, T));
213 
214  for (++phasei; phasei != phases_.end(); ++phasei)
215  {
216  the.ref() += phasei().Alpha()*phasei().thermo().he(p, T);
217  }
218 
219  return the;
220 }
221 
222 
224 (
225  const scalarField& T,
226  const labelList& cells
227 ) const
228 {
229  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
230 
231  tmp<scalarField> the
232  (
233  scalarField(phasei(), cells)*phasei().thermo().he(T, cells)
234  );
235 
236  for (++phasei; phasei != phases_.end(); ++phasei)
237  {
238  the.ref() +=
239  scalarField(phasei(), cells)*phasei().thermo().he(T, cells);
240  }
241 
242  return the;
243 }
244 
245 
247 (
248  const scalarField& T,
249  const label patchi
250 ) const
251 {
252  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
253 
254  tmp<scalarField> the
255  (
256  phasei().Alpha().boundaryField()[patchi]*phasei().thermo().he(T, patchi)
257  );
258 
259  for (++phasei; phasei != phases_.end(); ++phasei)
260  {
261  the.ref() +=
262  phasei().Alpha().boundaryField()[patchi]
263  *phasei().thermo().he(T, patchi);
264  }
265 
266  return the;
267 }
268 
269 
271 {
272  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
273 
274  tmp<volScalarField> ths(phasei().Alpha()*phasei().thermo().hs());
275 
276  for (++phasei; phasei != phases_.end(); ++phasei)
277  {
278  ths.ref() += phasei().Alpha()*phasei().thermo().hs();
279  }
280 
281  return ths;
282 }
283 
284 
286 (
287  const volScalarField& p,
288  const volScalarField& T
289 ) const
290 {
291  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
292 
293  tmp<volScalarField> ths(phasei().Alpha()*phasei().thermo().hs(p, T));
294 
295  for (++phasei; phasei != phases_.end(); ++phasei)
296  {
297  ths.ref() += phasei().Alpha()*phasei().thermo().hs(p, T);
298  }
299 
300  return ths;
301 }
302 
303 
305 (
306  const scalarField& T,
307  const labelList& cells
308 ) const
309 {
310  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
311 
312  tmp<scalarField> ths
313  (
314  scalarField(phasei(), cells)*phasei().thermo().hs(T, cells)
315  );
316 
317  for (++phasei; phasei != phases_.end(); ++phasei)
318  {
319  ths.ref() +=
320  scalarField(phasei(), cells)*phasei().thermo().hs(T, cells);
321  }
322 
323  return ths;
324 }
325 
326 
328 (
329  const scalarField& T,
330  const label patchi
331 ) const
332 {
333  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
334 
335  tmp<scalarField> ths
336  (
337  phasei().Alpha().boundaryField()[patchi]*phasei().thermo().hs(T, patchi)
338  );
339 
340  for (++phasei; phasei != phases_.end(); ++phasei)
341  {
342  ths.ref() +=
343  phasei().Alpha().boundaryField()[patchi]
344  *phasei().thermo().hs(T, patchi);
345  }
346 
347  return ths;
348 }
349 
350 
352 {
353  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
354 
355  tmp<volScalarField> tha(phasei().Alpha()*phasei().thermo().ha());
356 
357  for (++phasei; phasei != phases_.end(); ++phasei)
358  {
359  tha.ref() += phasei().Alpha()*phasei().thermo().ha();
360  }
361 
362  return tha;
363 }
364 
365 
367 (
368  const volScalarField& p,
369  const volScalarField& T
370 ) const
371 {
372  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
373 
374  tmp<volScalarField> tha(phasei().Alpha()*phasei().thermo().ha(p, T));
375 
376  for (++phasei; phasei != phases_.end(); ++phasei)
377  {
378  tha.ref() += phasei().Alpha()*phasei().thermo().ha(p, T);
379  }
380 
381  return tha;
382 }
383 
384 
386 (
387  const scalarField& T,
388  const labelList& cells
389 ) const
390 {
391  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
392 
393  tmp<scalarField> tha
394  (
395  scalarField(phasei(), cells)*phasei().thermo().ha(T, cells)
396  );
397 
398  for (++phasei; phasei != phases_.end(); ++phasei)
399  {
400  tha.ref() +=
401  scalarField(phasei(), cells)*phasei().thermo().ha(T, cells);
402  }
403 
404  return tha;
405 }
406 
407 
409 (
410  const scalarField& T,
411  const label patchi
412 ) const
413 {
414  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
415 
416  tmp<scalarField> tha
417  (
418  phasei().Alpha().boundaryField()[patchi]*phasei().thermo().ha(T, patchi)
419  );
420 
421  for (++phasei; phasei != phases_.end(); ++phasei)
422  {
423  tha.ref() +=
424  phasei().Alpha().boundaryField()[patchi]
425  *phasei().thermo().ha(T, patchi);
426  }
427 
428  return tha;
429 }
430 
431 
433 {
434  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
435 
436  tmp<volScalarField> thc(phasei().Alpha()*phasei().thermo().hc());
437 
438  for (++phasei; phasei != phases_.end(); ++phasei)
439  {
440  thc.ref() += phasei().Alpha()*phasei().thermo().hc();
441  }
442 
443  return thc;
444 }
445 
446 
448 (
449  const volScalarField& h,
450  const volScalarField& p,
451  const volScalarField& T0
452 ) const
453 {
455  return T0;
456 }
457 
458 
460 (
461  const scalarField& h,
462  const scalarField& T0,
463  const labelList& cells
464 ) const
465 {
467  return T0;
468 }
469 
470 
472 (
473  const scalarField& h,
474  const scalarField& T0,
475  const label patchi
476 ) const
477 {
479  return T0;
480 }
481 
482 
484 {
485  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
486 
487  tmp<volScalarField> tCp(phasei().Alpha()*phasei().thermo().Cp());
488 
489  for (++phasei; phasei != phases_.end(); ++phasei)
490  {
491  tCp.ref() += phasei().Alpha()*phasei().thermo().Cp();
492  }
493 
494  return tCp;
495 }
496 
497 
499 (
500  const scalarField& T,
501  const label patchi
502 ) const
503 {
504  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
505 
506  tmp<scalarField> tCp
507  (
508  phasei().Alpha().boundaryField()[patchi]*phasei().thermo().Cp(T, patchi)
509  );
510 
511  for (++phasei; phasei != phases_.end(); ++phasei)
512  {
513  tCp.ref() +=
514  phasei().Alpha().boundaryField()[patchi]
515  *phasei().thermo().Cp(T, patchi);
516  }
517 
518  return tCp;
519 }
520 
521 
523 {
524  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
525 
526  tmp<volScalarField> tCv(phasei().Alpha()*phasei().thermo().Cv());
527 
528  for (++phasei; phasei != phases_.end(); ++phasei)
529  {
530  tCv.ref() += phasei().Alpha()*phasei().thermo().Cv();
531  }
532 
533  return tCv;
534 }
535 
536 
538 (
539  const scalarField& T,
540  const label patchi
541 ) const
542 {
543  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
544 
545  tmp<scalarField> tCv
546  (
547  phasei().Alpha().boundaryField()[patchi]*phasei().thermo().Cv(T, patchi)
548  );
549 
550  for (++phasei; phasei != phases_.end(); ++phasei)
551  {
552  tCv.ref() +=
553  phasei().Alpha().boundaryField()[patchi]
554  *phasei().thermo().Cv(T, patchi);
555  }
556 
557  return tCv;
558 }
559 
560 
562 {
563  return Cp()/Cv();
564 }
565 
566 
568 (
569  const scalarField& T,
570  const label patchi
571 ) const
572 {
573  return Cp(T, patchi)/Cv(T, patchi);
574 }
575 
576 
578 {
579  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
580 
581  tmp<volScalarField> tCpv(phasei().Alpha()*phasei().thermo().Cpv());
582 
583  for (++phasei; phasei != phases_.end(); ++phasei)
584  {
585  tCpv.ref() += phasei().Alpha()*phasei().thermo().Cpv();
586  }
587 
588  return tCpv;
589 }
590 
591 
593 (
594  const scalarField& T,
595  const label patchi
596 ) const
597 {
598  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
599 
600  tmp<scalarField> tCpv
601  (
602  phasei().Alpha().boundaryField()[patchi]
603  *phasei().thermo().Cpv(T, patchi)
604  );
605 
606  for (++phasei; phasei != phases_.end(); ++phasei)
607  {
608  tCpv.ref() +=
609  phasei().boundaryField()[patchi]
610  *phasei().thermo().Cpv(T, patchi);
611  }
612 
613  return tCpv;
614 }
615 
616 
618 {
619  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
620 
621  volScalarField rW(phasei().Alpha()/phasei().thermo().W());
622 
623  for (++phasei; phasei != phases_.end(); ++phasei)
624  {
625  rW += phasei().Alpha()/phasei().thermo().W();
626  }
627 
628  return 1/rW;
629 }
630 
631 
633 (
634  const label patchi
635 ) const
636 {
637  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
638 
639  scalarField rW
640  (
641  phasei().Alpha().boundaryField()[patchi]/phasei().thermo().W(patchi)
642  );
643 
644  for (++phasei; phasei != phases_.end(); ++phasei)
645  {
646  rW +=
647  phasei().Alpha().boundaryField()[patchi]
648  /phasei().thermo().W(patchi);
649  }
650 
651  return 1/rW;
652 }
653 
654 
656 {
657  return mu()/rho();
658 }
659 
660 
662 (
663  const label patchi
664 ) const
665 {
666  return mu(patchi)/rho(patchi);
667 }
668 
669 
671 {
672  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
673 
674  tmp<volScalarField> tkappa(phasei()*phasei().thermo().kappa());
675 
676  for (++phasei; phasei != phases_.end(); ++phasei)
677  {
678  tkappa.ref() += phasei()*phasei().thermo().kappa();
679  }
680 
681  return tkappa;
682 }
683 
684 
686 (
687  const label patchi
688 ) const
689 {
690  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
691 
692  tmp<scalarField> tkappa
693  (
694  phasei().boundaryField()[patchi]*phasei().thermo().kappa(patchi)
695  );
696 
697  for (++phasei; phasei != phases_.end(); ++phasei)
698  {
699  tkappa.ref() +=
700  phasei().boundaryField()[patchi]*phasei().thermo().kappa(patchi);
701  }
702 
703  return tkappa;
704 }
705 
706 
708 {
709  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
710 
711  tmp<volScalarField> talphaEff(phasei()*phasei().thermo().alphahe());
712 
713  for (++phasei; phasei != phases_.end(); ++phasei)
714  {
715  talphaEff.ref() += phasei()*phasei().thermo().alphahe();
716  }
717 
718  return talphaEff;
719 }
720 
721 
723 (
724  const label patchi
725 ) const
726 {
727  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
728 
729  tmp<scalarField> talphaEff
730  (
731  phasei().boundaryField()[patchi]
732  *phasei().thermo().alphahe(patchi)
733  );
734 
735  for (++phasei; phasei != phases_.end(); ++phasei)
736  {
737  talphaEff.ref() +=
738  phasei().boundaryField()[patchi]
739  *phasei().thermo().alphahe(patchi);
740  }
741 
742  return talphaEff;
743 }
744 
745 
747 (
748  const volScalarField& alphat
749 ) const
750 {
751  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
752 
753  tmp<volScalarField> tkappaEff(phasei()*phasei().thermo().kappaEff(alphat));
754 
755  for (++phasei; phasei != phases_.end(); ++phasei)
756  {
757  tkappaEff.ref() += phasei()*phasei().thermo().kappaEff(alphat);
758  }
759 
760  return tkappaEff;
761 }
762 
763 
765 (
766  const scalarField& alphat,
767  const label patchi
768 ) const
769 {
770  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
771 
772  tmp<scalarField> tkappaEff
773  (
774  phasei().boundaryField()[patchi]
775  *phasei().thermo().kappaEff(alphat, patchi)
776  );
777 
778  for (++phasei; phasei != phases_.end(); ++phasei)
779  {
780  tkappaEff.ref() +=
781  phasei().boundaryField()[patchi]
782  *phasei().thermo().kappaEff(alphat, patchi);
783  }
784 
785  return tkappaEff;
786 }
787 
788 
790 (
791  const volScalarField& alphat
792 ) const
793 {
794  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
795 
796  tmp<volScalarField> talphaEff(phasei()*phasei().thermo().alphaEff(alphat));
797 
798  for (++phasei; phasei != phases_.end(); ++phasei)
799  {
800  talphaEff.ref() += phasei()*phasei().thermo().alphaEff(alphat);
801  }
802 
803  return talphaEff;
804 }
805 
806 
808 (
809  const scalarField& alphat,
810  const label patchi
811 ) const
812 {
813  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
814 
815  tmp<scalarField> talphaEff
816  (
817  phasei().boundaryField()[patchi]
818  *phasei().thermo().alphaEff(alphat, patchi)
819  );
820 
821  for (++phasei; phasei != phases_.end(); ++phasei)
822  {
823  talphaEff.ref() +=
824  phasei().boundaryField()[patchi]
825  *phasei().thermo().alphaEff(alphat, patchi);
826  }
827 
828  return talphaEff;
829 }
830 
831 
833 {
834  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
835 
836  tmp<volScalarField> trCv(phasei()/phasei().thermo().Cv());
837 
838  for (++phasei; phasei != phases_.end(); ++phasei)
839  {
840  trCv.ref() += phasei()/phasei().thermo().Cv();
841  }
842 
843  return trCv;
844 }
845 
846 
849 {
850  tmp<surfaceScalarField> tstf
851  (
853  (
854  "surfaceTensionForce",
855  mesh_,
856  dimensionedScalar(dimensionSet(1, -2, -2, 0, 0), 0)
857  )
858  );
859 
860  surfaceScalarField& stf = tstf.ref();
861 
862  forAllConstIter(PtrDictionary<phaseModel>, phases_, phase1)
863  {
864  const phaseModel& alpha1 = phase1();
865 
867  ++phase2;
868 
869  for (; phase2 != phases_.end(); ++phase2)
870  {
871  const phaseModel& alpha2 = phase2();
872 
873  sigmaTable::const_iterator sigma =
874  sigmas_.find(interfacePair(alpha1, alpha2));
875 
876  if (sigma == sigmas_.end())
877  {
879  << "Cannot find interface " << interfacePair(alpha1, alpha2)
880  << " in list of sigma values"
881  << exit(FatalError);
882  }
883 
884  stf += dimensionedScalar(dimSigma_, sigma())
885  *fvc::interpolate(K(alpha1, alpha2))*
886  (
887  fvc::interpolate(alpha2)*fvc::snGrad(alpha1)
888  - fvc::interpolate(alpha1)*fvc::snGrad(alpha2)
889  );
890  }
891  }
892 
893  return tstf;
894 }
895 
896 
898 {
899  const Time& runTime = mesh_.time();
900 
901  const dictionary& alphaControls = mesh_.solverDict("alpha");
902  label nAlphaSubCycles(alphaControls.lookup<label>("nAlphaSubCycles"));
903  scalar cAlpha(alphaControls.lookup<scalar>("cAlpha"));
904 
905  volScalarField& alpha = phases_.first();
906 
907  if (nAlphaSubCycles > 1)
908  {
909  surfaceScalarField rhoPhiSum(0.0*rhoPhi_);
910  dimensionedScalar totalDeltaT = runTime.deltaT();
911 
912  for
913  (
914  subCycle<volScalarField> alphaSubCycle(alpha, nAlphaSubCycles);
915  !(++alphaSubCycle).end();
916  )
917  {
918  solveAlphas(cAlpha);
919  rhoPhiSum += (runTime.deltaT()/totalDeltaT)*rhoPhi_;
920  }
921 
922  rhoPhi_ = rhoPhiSum;
923  }
924  else
925  {
926  solveAlphas(cAlpha);
927  }
928 
929  correct();
930 }
931 
932 
933 Foam::tmp<Foam::surfaceVectorField> Foam::multiphaseMixtureThermo::nHatfv
934 (
935  const volScalarField& alpha1,
936  const volScalarField& alpha2
937 ) const
938 {
939  /*
940  // Cell gradient of alpha
941  volVectorField gradAlpha =
942  alpha2*fvc::grad(alpha1) - alpha1*fvc::grad(alpha2);
943 
944  // Interpolated face-gradient of alpha
945  surfaceVectorField gradAlphaf = fvc::interpolate(gradAlpha);
946  */
947 
948  surfaceVectorField gradAlphaf
949  (
951  - fvc::interpolate(alpha1)*fvc::interpolate(fvc::grad(alpha2))
952  );
953 
954  // Face unit interface normal
955  return gradAlphaf/(mag(gradAlphaf) + deltaN_);
956 }
957 
958 
959 Foam::tmp<Foam::surfaceScalarField> Foam::multiphaseMixtureThermo::nHatf
960 (
961  const volScalarField& alpha1,
962  const volScalarField& alpha2
963 ) const
964 {
965  // Face unit interface normal flux
966  return nHatfv(alpha1, alpha2) & mesh_.Sf();
967 }
968 
969 
970 // Correction for the boundary condition on the unit normal nHat on
971 // walls to produce the correct contact angle.
972 
973 // The dynamic contact angle is calculated from the component of the
974 // velocity on the direction of the interface, parallel to the wall.
975 
976 void Foam::multiphaseMixtureThermo::correctContactAngle
977 (
978  const phaseModel& alpha1,
979  const phaseModel& alpha2,
980  surfaceVectorField::Boundary& nHatb
981 ) const
982 {
983  const volScalarField::Boundary& a1bf = alpha1.boundaryField();
984  const volScalarField::Boundary& a2bf = alpha2.boundaryField();
985 
986  const fvBoundaryMesh& boundary = mesh_.boundary();
987 
988  forAll(boundary, patchi)
989  {
990  if
991  (
992  isA<alphaContactAngleFvPatchScalarField>(a1bf[patchi])
993  || isA<alphaContactAngleFvPatchScalarField>(a2bf[patchi])
994  )
995  {
996  if
997  (
998  isA<alphaContactAngleFvPatchScalarField>(a1bf[patchi])
999  && isA<alphaContactAngleFvPatchScalarField>(a2bf[patchi])
1000  )
1001  {
1003  << "alphaContactAngle boundary condition "
1004  "specified on patch " << boundary[patchi].name()
1005  << " for both " << alpha1.name() << " and " << alpha2.name()
1006  << nl << "which may be inconsistent."
1007  << exit(FatalError);
1008  }
1009 
1010  const alphaContactAngleFvPatchScalarField& acap =
1011  isA<alphaContactAngleFvPatchScalarField>(a1bf[patchi])
1012  ? refCast<const alphaContactAngleFvPatchScalarField>(a1bf[patchi])
1013  : refCast<const alphaContactAngleFvPatchScalarField>(a2bf[patchi])
1014  ;
1015 
1016  vectorField& nHatPatch = nHatb[patchi];
1017 
1018  vectorField AfHatPatch
1019  (
1020  mesh_.Sf().boundaryField()[patchi]
1021  /mesh_.magSf().boundaryField()[patchi]
1022  );
1023 
1024  alphaContactAngleFvPatchScalarField::thetaPropsTable::
1025  const_iterator tp =
1026  acap.thetaProps().find(interfacePair(alpha1, alpha2));
1027 
1028  if (tp == acap.thetaProps().end())
1029  {
1031  << "Cannot find interface " << interfacePair(alpha1, alpha2)
1032  << "\n in table of theta properties for patch "
1033  << acap.patch().name()
1034  << exit(FatalError);
1035  }
1036 
1037  const bool matched = (tp.key().first() == alpha1.name());
1038 
1039  const scalar theta0 = degToRad(tp().theta0(matched));
1040 
1041  scalarField theta(boundary[patchi].size(), theta0);
1042 
1043  const scalar uTheta = tp().uTheta();
1044 
1045  // Calculate the dynamic contact angle if required
1046  if (uTheta > small)
1047  {
1048  const scalar thetaA = degToRad(tp().thetaA(matched));
1049  const scalar thetaR = degToRad(tp().thetaR(matched));
1050 
1051  // Calculated the component of the velocity parallel to the wall
1052  vectorField Uwall
1053  (
1054  U_.boundaryField()[patchi].patchInternalField()
1055  - U_.boundaryField()[patchi]
1056  );
1057  Uwall -= (AfHatPatch & Uwall)*AfHatPatch;
1058 
1059  // Find the direction of the interface parallel to the wall
1060  vectorField nWall
1061  (
1062  nHatPatch - (AfHatPatch & nHatPatch)*AfHatPatch
1063  );
1064 
1065  // Normalise nWall
1066  nWall /= (mag(nWall) + small);
1067 
1068  // Calculate Uwall resolved normal to the interface parallel to
1069  // the interface
1070  const scalarField uwall(nWall & Uwall);
1071 
1072  theta += (thetaA - thetaR)*tanh(uwall/uTheta);
1073  }
1074 
1075 
1076  // Reset nHatPatch to correspond to the contact angle
1077 
1078  const scalarField a12(nHatPatch & AfHatPatch);
1079 
1080  const scalarField b1(cos(theta));
1081 
1082  scalarField b2(nHatPatch.size());
1083 
1084  forAll(b2, facei)
1085  {
1086  b2[facei] = cos(acos(a12[facei]) - theta[facei]);
1087  }
1088 
1089  const scalarField det(1.0 - a12*a12);
1090 
1091  const scalarField a((b1 - a12*b2)/det);
1092  const scalarField b((b2 - a12*b1)/det);
1093 
1094  nHatPatch = a*AfHatPatch + b*nHatPatch;
1095 
1096  nHatPatch /= (mag(nHatPatch) + deltaN_.value());
1097  }
1098  }
1099 }
1100 
1101 
1102 Foam::tmp<Foam::volScalarField> Foam::multiphaseMixtureThermo::K
1103 (
1104  const phaseModel& alpha1,
1105  const phaseModel& alpha2
1106 ) const
1107 {
1108  tmp<surfaceVectorField> tnHatfv = nHatfv(alpha1, alpha2);
1109 
1110  correctContactAngle(alpha1, alpha2, tnHatfv.ref().boundaryFieldRef());
1111 
1112  // Simple expression for curvature
1113  return -fvc::div(tnHatfv & mesh_.Sf());
1114 }
1115 
1116 
1119 {
1120  tmp<volScalarField> tnearInt
1121  (
1123  (
1124  "nearInterface",
1125  mesh_,
1127  )
1128  );
1129 
1130  forAllConstIter(PtrDictionary<phaseModel>, phases_, phase)
1131  {
1132  tnearInt.ref() =
1133  max(tnearInt(), pos0(phase() - 0.01)*pos0(0.99 - phase()));
1134  }
1135 
1136  return tnearInt;
1137 }
1138 
1139 
1140 void Foam::multiphaseMixtureThermo::solveAlphas
1141 (
1142  const scalar cAlpha
1143 )
1144 {
1145  static label nSolves=-1;
1146  nSolves++;
1147 
1148  word alphaScheme("div(phi,alpha)");
1149  word alpharScheme("div(phirb,alpha)");
1150 
1151  surfaceScalarField phic(mag(phi_/mesh_.magSf()));
1152  phic = min(cAlpha*phic, max(phic));
1153 
1154  UPtrList<const volScalarField> alphas(phases_.size());
1155  PtrList<surfaceScalarField> alphaPhis(phases_.size());
1156  int phasei = 0;
1157 
1158  forAllIter(PtrDictionary<phaseModel>, phases_, phase)
1159  {
1160  const phaseModel& alpha = phase();
1161 
1162  alphas.set(phasei, &alpha);
1163 
1164  alphaPhis.set
1165  (
1166  phasei,
1167  new surfaceScalarField
1168  (
1169  phi_.name() + alpha.name(),
1170  fvc::flux
1171  (
1172  phi_,
1173  alpha,
1174  alphaScheme
1175  )
1176  )
1177  );
1178 
1179  surfaceScalarField& alphaPhi = alphaPhis[phasei];
1180 
1181  forAllIter(PtrDictionary<phaseModel>, phases_, phase2)
1182  {
1183  phaseModel& alpha2 = phase2();
1184 
1185  if (&alpha2 == &alpha) continue;
1186 
1187  surfaceScalarField phir(phic*nHatf(alpha, alpha2));
1188 
1189  alphaPhi += fvc::flux
1190  (
1191  -fvc::flux(-phir, alpha2, alpharScheme),
1192  alpha,
1193  alpharScheme
1194  );
1195  }
1196 
1197  // Limit alphaPhi for each phase
1198  MULES::limit
1199  (
1200  1.0/mesh_.time().deltaT().value(),
1201  geometricOneField(),
1202  alpha,
1203  phi_,
1204  alphaPhi,
1205  zeroField(),
1206  zeroField(),
1207  oneField(),
1208  zeroField(),
1209  false
1210  );
1211 
1212  phasei++;
1213  }
1214 
1215  MULES::limitSum(alphas, alphaPhis, phi_);
1216 
1217  rhoPhi_ = dimensionedScalar(dimensionSet(1, 0, -1, 0, 0), 0);
1218 
1219  volScalarField sumAlpha
1220  (
1221  IOobject
1222  (
1223  "sumAlpha",
1224  mesh_.time().timeName(),
1225  mesh_
1226  ),
1227  mesh_,
1229  );
1230 
1232 
1233  phasei = 0;
1234 
1235  forAllIter(PtrDictionary<phaseModel>, phases_, phase)
1236  {
1237  phaseModel& alpha = phase();
1238 
1239  surfaceScalarField& alphaPhi = alphaPhis[phasei];
1240 
1242  (
1243  IOobject
1244  (
1245  "Sp",
1246  mesh_.time().timeName(),
1247  mesh_
1248  ),
1249  mesh_,
1250  dimensionedScalar(alpha.dgdt().dimensions(), 0)
1251  );
1252 
1254  (
1255  IOobject
1256  (
1257  "Su",
1258  mesh_.time().timeName(),
1259  mesh_
1260  ),
1261  // Divergence term is handled explicitly to be
1262  // consistent with the explicit transport solution
1263  divU.v()*min(alpha.v(), scalar(1))
1264  );
1265 
1266  {
1267  const scalarField& dgdt = alpha.dgdt();
1268 
1269  forAll(dgdt, celli)
1270  {
1271  if (dgdt[celli] < 0.0 && alpha[celli] > 0.0)
1272  {
1273  Sp[celli] += dgdt[celli]*alpha[celli];
1274  Su[celli] -= dgdt[celli]*alpha[celli];
1275  }
1276  else if (dgdt[celli] > 0.0 && alpha[celli] < 1.0)
1277  {
1278  Sp[celli] -= dgdt[celli]*(1.0 - alpha[celli]);
1279  }
1280  }
1281  }
1282 
1283  forAllConstIter(PtrDictionary<phaseModel>, phases_, phase2)
1284  {
1285  const phaseModel& alpha2 = phase2();
1286 
1287  if (&alpha2 == &alpha) continue;
1288 
1289  const scalarField& dgdt2 = alpha2.dgdt();
1290 
1291  forAll(dgdt2, celli)
1292  {
1293  if (dgdt2[celli] > 0.0 && alpha2[celli] < 1.0)
1294  {
1295  Sp[celli] -= dgdt2[celli]*(1.0 - alpha2[celli]);
1296  Su[celli] += dgdt2[celli]*alpha[celli];
1297  }
1298  else if (dgdt2[celli] < 0.0 && alpha2[celli] > 0.0)
1299  {
1300  Sp[celli] += dgdt2[celli]*alpha2[celli];
1301  }
1302  }
1303  }
1304 
1306  (
1307  geometricOneField(),
1308  alpha,
1309  alphaPhi,
1310  Sp,
1311  Su
1312  );
1313 
1314  rhoPhi_ += fvc::interpolate(alpha.thermo().rho())*alphaPhi;
1315 
1316  Info<< alpha.name() << " volume fraction, min, max = "
1317  << alpha.weightedAverage(mesh_.V()).value()
1318  << ' ' << min(alpha).value()
1319  << ' ' << max(alpha).value()
1320  << endl;
1321 
1322  sumAlpha += alpha;
1323 
1324  phasei++;
1325  }
1326 
1327  Info<< "Phase-sum volume fraction, min, max = "
1328  << sumAlpha.weightedAverage(mesh_.V()).value()
1329  << ' ' << min(sumAlpha).value()
1330  << ' ' << max(sumAlpha).value()
1331  << endl;
1332 
1333  calcAlphas();
1334 }
1335 
1336 
1337 // ************************************************************************* //
scalar Cv(const scalar p, const scalar T) const
Definition: HtoEthermo.H:2
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
fluidReactionThermo & thermo
Definition: createFields.H:28
volScalarField dgdt(alpha1 *fvc::div(phi))
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
label phasei
Definition: pEqn.H:27
const word & name() const
Return const access to the source name.
Definition: fvModelI.H:28
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:323
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Unit conversion functions.
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:106
virtual tmp< volScalarField > hc() const
Enthalpy of formation [J/kg].
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
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
const dimensionedScalar b
Wien displacement law constant: default SI units: [m K].
Definition: createFields.H:27
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Calculate the snGrad of the given volField.
surfaceScalarField phir(fvc::flux(UdmModel.Udm()))
const dimensionSet dimless
dimensionedScalar det(const dimensionedSphericalTensor &dt)
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.
const dimensionedScalar h
Planck constant.
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 void correct()
Solve the film and update the sources.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
void correctRho(const volScalarField &dp)
Update densities for given pressure change.
const dimensionSet dimTime
void solve()
Solve for the mixture phase-fractions.
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
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].
Calculate the gradient of the given field.
phic
Definition: correctPhic.H:2
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.
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)
static const char nl
Definition: Ostream.H:260
const dimensionedScalar mu
Atomic mass unit.
thermo he()
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 > &)
tmp< volScalarField > divU
Definition: alphaSuSp.H:9
const dimensionSet dimMass
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.
const tmp< volScalarField::Internal > & Sp
Definition: alphaSuSp.H:7
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
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
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 > &)
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
const tmp< volScalarField::Internal > & Su
Definition: alphaSuSp.H:6
virtual tmp< volScalarField > Cv() const
Heat capacity at constant volume [J/kg/K].
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:370
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 void correctThermo()
Correct the thermodynamics of each phase.
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.
virtual tmp< volScalarField > THE(const volScalarField &h, const volScalarField &p, const volScalarField &T0) const
Temperature from enthalpy/internal energy.
scalar T0
Definition: createFields.H:22
virtual bool incompressible() const
Return true if the equation of state is incompressible.
virtual void correct()
Update properties.