Loading...
Searching...
No Matches
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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2013-2017 OpenFOAM Foundation
9 Copyright (C) 2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
30#include "alphaContactAngleFvPatchScalarField.H"
31#include "Time.H"
32#include "subCycle.H"
33#include "MULES.H"
34#include "fvcDiv.H"
35#include "fvcGrad.H"
36#include "fvcSnGrad.H"
37#include "fvcFlux.H"
38#include "fvcMeshPhi.H"
39#include "surfaceInterpolate.H"
40#include "unitConversion.H"
41
42// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43
44namespace Foam
45{
47}
48
49
50// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51
52void Foam::multiphaseMixtureThermo::calcAlphas()
53{
54 scalar level = 0.0;
55 alphas_ == 0.0;
56
57 for (const phaseModel& phase : phases_)
58 {
59 alphas_ += level * phase;
60 level += 1.0;
61 }
62}
63
64
65// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
66
68(
69 const volVectorField& U,
71)
72:
73 psiThermo(U.mesh(), word::null),
74 phases_(lookup("phases"), phaseModel::iNew(p_, T_)),
75
76 mesh_(U.mesh()),
77 U_(U),
78 phi_(phi),
79
80 rhoPhi_
81 (
83 (
84 "rhoPhi",
85 mesh_.time().timeName(),
86 mesh_,
87 IOobject::NO_READ,
88 IOobject::NO_WRITE
89 ),
90 mesh_,
92 ),
93
94 alphas_
95 (
97 (
98 "alphas",
99 mesh_.time().timeName(),
100 mesh_,
101 IOobject::NO_READ,
102 IOobject::AUTO_WRITE,
103 IOobject::REGISTER
104 ),
105 mesh_,
107 ),
108
109 sigmas_(lookup("sigmas")),
110 dimSigma_(1, 0, -2, 0, 0),
111 deltaN_
112 (
113 "deltaN",
114 1e-8/cbrt(average(mesh_.V()))
115 )
116{
117 rhoPhi_.setOriented();
118 calcAlphas();
119 alphas_.write();
120 correct();
121}
122
123
124// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
125
127{
128 for (phaseModel& phase : phases_)
129 {
130 phase.correct();
131 }
132
133 auto phasei = phases_.cbegin();
134
135 psi_ = phasei()*phasei().thermo().psi();
136 mu_ = phasei()*phasei().thermo().mu();
137 alpha_ = phasei()*phasei().thermo().alpha();
138
139 for (++phasei; phasei != phases_.cend(); ++phasei)
140 {
141 psi_ += phasei()*phasei().thermo().psi();
142 mu_ += phasei()*phasei().thermo().mu();
143 alpha_ += phasei()*phasei().thermo().alpha();
144 }
145}
146
147
149{
150 for (phaseModel& phase : phases_)
151 {
152 phase.thermo().rho() += phase.thermo().psi()*dp;
153 }
154}
155
156
158{
159 auto phasei = phases_.cbegin();
160
161 word name = phasei().thermo().thermoName();
162
163 for (++phasei; phasei != phases_.cend(); ++phasei)
164 {
165 name += ',' + phasei().thermo().thermoName();
166 }
167
168 return name;
169}
170
171
173{
174 for (const phaseModel& phase : phases_)
175 {
176 if (!phase.thermo().incompressible())
177 {
178 return false;
179 }
180 }
181
182 return true;
183}
184
185
187{
188 for (const phaseModel& phase : phases_)
189 {
190 if (!phase.thermo().isochoric())
191 {
192 return false;
193 }
194 }
195
196 return true;
197}
198
199
200Foam::tmp<Foam::volScalarField> Foam::multiphaseMixtureThermo::he
201(
202 const volScalarField& p,
203 const volScalarField& T
204) const
205{
206 auto phasei = phases_.cbegin();
207
209
210 for (++phasei; phasei != phases_.cend(); ++phasei)
211 {
212 the.ref() += phasei()*phasei().thermo().he(p, T);
213 }
214
215 return the;
216}
217
218
219Foam::tmp<Foam::scalarField> Foam::multiphaseMixtureThermo::he
220(
221 const scalarField& p,
222 const scalarField& T,
223 const labelList& cells
224) const
225{
226 auto phasei = phases_.cbegin();
227
229 (
231 );
232
233 for (++phasei; phasei != phases_.cend(); ++phasei)
234 {
235 the.ref() +=
236 scalarField(phasei(), cells)*phasei().thermo().he(p, T, cells);
237 }
238
239 return the;
240}
241
242
243Foam::tmp<Foam::scalarField> Foam::multiphaseMixtureThermo::he
244(
245 const scalarField& p,
246 const scalarField& T,
247 const label patchi
248) const
249{
250 auto phasei = phases_.cbegin();
251
253 (
254 phasei().boundaryField()[patchi]*phasei().thermo().he(p, T, patchi)
255 );
256
257 for (++phasei; phasei != phases_.cend(); ++phasei)
258 {
259 the.ref() +=
260 phasei().boundaryField()[patchi]*phasei().thermo().he(p, T, patchi);
261 }
262
263 return the;
264}
265
266
267Foam::tmp<Foam::volScalarField> Foam::multiphaseMixtureThermo::hc() const
268{
269 auto phasei = phases_.cbegin();
270
271 tmp<volScalarField> thc(phasei()*phasei().thermo().hc());
272
273 for (++phasei; phasei != phases_.cend(); ++phasei)
274 {
275 thc.ref() += phasei()*phasei().thermo().hc();
276 }
277
278 return thc;
279}
280
281
282Foam::tmp<Foam::scalarField> Foam::multiphaseMixtureThermo::THE
283(
284 const scalarField& h,
285 const scalarField& p,
286 const scalarField& T0,
287 const labelList& cells
288) const
289{
291 return T0;
292}
293
294
295Foam::tmp<Foam::scalarField> Foam::multiphaseMixtureThermo::THE
296(
297 const scalarField& h,
298 const scalarField& p,
299 const scalarField& T0,
300 const label patchi
301) const
302{
304 return T0;
305}
306
307
308Foam::tmp<Foam::volScalarField> Foam::multiphaseMixtureThermo::rho() const
309{
310 auto phasei = phases_.cbegin();
311
313
314 for (++phasei; phasei != phases_.cend(); ++phasei)
315 {
316 trho.ref() += phasei()*phasei().thermo().rho();
317 }
318
319 return trho;
320}
321
322
323Foam::tmp<Foam::scalarField> Foam::multiphaseMixtureThermo::rho
324(
325 const label patchi
326) const
327{
328 auto phasei = phases_.cbegin();
329
331 (
332 phasei().boundaryField()[patchi]*phasei().thermo().rho(patchi)
333 );
334
335 for (++phasei; phasei != phases_.cend(); ++phasei)
336 {
337 trho.ref() +=
338 phasei().boundaryField()[patchi]*phasei().thermo().rho(patchi);
339 }
340
341 return trho;
342}
343
344
345Foam::tmp<Foam::volScalarField> Foam::multiphaseMixtureThermo::Cp() const
346{
347 auto phasei = phases_.cbegin();
348
350
351 for (++phasei; phasei != phases_.cend(); ++phasei)
352 {
353 tCp.ref() += phasei()*phasei().thermo().Cp();
354 }
355
356 return tCp;
357}
358
359
360Foam::tmp<Foam::scalarField> Foam::multiphaseMixtureThermo::Cp
361(
362 const scalarField& p,
363 const scalarField& T,
364 const label patchi
365) const
366{
367 auto phasei = phases_.cbegin();
368
370 (
371 phasei().boundaryField()[patchi]*phasei().thermo().Cp(p, T, patchi)
372 );
373
374 for (++phasei; phasei != phases_.cend(); ++phasei)
375 {
376 tCp.ref() +=
377 phasei().boundaryField()[patchi]*phasei().thermo().Cp(p, T, patchi);
378 }
379
380 return tCp;
381}
382
383
384Foam::tmp<Foam::volScalarField> Foam::multiphaseMixtureThermo::Cv() const
385{
386 auto phasei = phases_.cbegin();
387
389
390 for (++phasei; phasei != phases_.cend(); ++phasei)
391 {
392 tCv.ref() += phasei()*phasei().thermo().Cv();
393 }
394
395 return tCv;
396}
397
398
399Foam::tmp<Foam::scalarField> Foam::multiphaseMixtureThermo::Cv
400(
401 const scalarField& p,
402 const scalarField& T,
403 const label patchi
404) const
405{
406 auto phasei = phases_.cbegin();
407
409 (
410 phasei().boundaryField()[patchi]*phasei().thermo().Cv(p, T, patchi)
411 );
412
413 for (++phasei; phasei != phases_.cend(); ++phasei)
414 {
415 tCv.ref() +=
416 phasei().boundaryField()[patchi]*phasei().thermo().Cv(p, T, patchi);
417 }
418
419 return tCv;
420}
421
422
423Foam::tmp<Foam::volScalarField> Foam::multiphaseMixtureThermo::gamma() const
424{
425 auto phasei = phases_.cbegin();
426
428
429 for (++phasei; phasei != phases_.cend(); ++phasei)
430 {
431 tgamma.ref() += phasei()*phasei().thermo().gamma();
432 }
433
434 return tgamma;
435}
436
437
438Foam::tmp<Foam::scalarField> Foam::multiphaseMixtureThermo::gamma
439(
440 const scalarField& p,
441 const scalarField& T,
442 const label patchi
443) const
444{
445 auto phasei = phases_.cbegin();
446
447 tmp<scalarField> tgamma
448 (
449 phasei().boundaryField()[patchi]*phasei().thermo().gamma(p, T, patchi)
450 );
451
452 for (++phasei; phasei != phases_.cend(); ++phasei)
453 {
454 tgamma.ref() +=
455 phasei().boundaryField()[patchi]
456 *phasei().thermo().gamma(p, T, patchi);
457 }
458
459 return tgamma;
460}
461
462
463Foam::tmp<Foam::volScalarField> Foam::multiphaseMixtureThermo::Cpv() const
464{
465 auto phasei = phases_.cbegin();
466
467 tmp<volScalarField> tCpv(phasei()*phasei().thermo().Cpv());
468
469 for (++phasei; phasei != phases_.cend(); ++phasei)
470 {
471 tCpv.ref() += phasei()*phasei().thermo().Cpv();
472 }
473
474 return tCpv;
475}
476
477
478Foam::tmp<Foam::scalarField> Foam::multiphaseMixtureThermo::Cpv
479(
480 const scalarField& p,
481 const scalarField& T,
482 const label patchi
483) const
484{
485 auto phasei = phases_.cbegin();
486
488 (
489 phasei().boundaryField()[patchi]*phasei().thermo().Cpv(p, T, patchi)
490 );
491
492 for (++phasei; phasei != phases_.cend(); ++phasei)
493 {
494 tCpv.ref() +=
495 phasei().boundaryField()[patchi]
496 *phasei().thermo().Cpv(p, T, patchi);
497 }
498
499 return tCpv;
500}
501
502
503Foam::tmp<Foam::volScalarField> Foam::multiphaseMixtureThermo::CpByCpv() const
504{
505 auto phasei = phases_.cbegin();
506
507 tmp<volScalarField> tCpByCpv(phasei()*phasei().thermo().CpByCpv());
508
509 for (++phasei; phasei != phases_.cend(); ++phasei)
510 {
511 tCpByCpv.ref() += phasei()*phasei().thermo().CpByCpv();
512 }
513
514 return tCpByCpv;
515}
516
517
518Foam::tmp<Foam::scalarField> Foam::multiphaseMixtureThermo::CpByCpv
519(
520 const scalarField& p,
521 const scalarField& T,
522 const label patchi
523) const
524{
525 auto phasei = phases_.cbegin();
526
527 tmp<scalarField> tCpByCpv
528 (
529 phasei().boundaryField()[patchi]*phasei().thermo().CpByCpv(p, T, patchi)
530 );
531
532 for (++phasei; phasei != phases_.cend(); ++phasei)
533 {
534 tCpByCpv.ref() +=
535 phasei().boundaryField()[patchi]
536 *phasei().thermo().CpByCpv(p, T, patchi);
537 }
538
539 return tCpByCpv;
540}
541
542
543Foam::tmp<Foam::volScalarField> Foam::multiphaseMixtureThermo::W() const
544{
545 auto phasei = phases_.cbegin();
546
548
549 for (++phasei; phasei != phases_.cend(); ++phasei)
550 {
551 tW.ref() += phasei()*phasei().thermo().W();
552 }
553
554 return tW;
555}
556
557
558Foam::tmp<Foam::volScalarField> Foam::multiphaseMixtureThermo::nu() const
559{
560 return mu()/rho();
561}
562
563
564Foam::tmp<Foam::scalarField> Foam::multiphaseMixtureThermo::nu
565(
566 const label patchi
567) const
568{
569 return mu(patchi)/rho(patchi);
570}
571
572
573Foam::tmp<Foam::volScalarField> Foam::multiphaseMixtureThermo::kappa() const
574{
575 auto phasei = phases_.cbegin();
576
578
579 for (++phasei; phasei != phases_.cend(); ++phasei)
580 {
581 tkappa.ref() += phasei()*phasei().thermo().kappa();
582 }
583
584 return tkappa;
585}
586
587
588Foam::tmp<Foam::scalarField> Foam::multiphaseMixtureThermo::kappa
589(
590 const label patchi
591) const
592{
593 auto phasei = phases_.cbegin();
594
595 tmp<scalarField> tkappa
596 (
597 phasei().boundaryField()[patchi]*phasei().thermo().kappa(patchi)
598 );
599
600 for (++phasei; phasei != phases_.cend(); ++phasei)
601 {
602 tkappa.ref() +=
603 phasei().boundaryField()[patchi]*phasei().thermo().kappa(patchi);
604 }
605
606 return tkappa;
607}
608
609
610Foam::tmp<Foam::volScalarField> Foam::multiphaseMixtureThermo::alphahe() const
611{
613
614 tmp<volScalarField> talphaEff(phasei()*phasei().thermo().alphahe());
615
616 for (++phasei; phasei != phases_.end(); ++phasei)
617 {
618 talphaEff.ref() += phasei()*phasei().thermo().alphahe();
619 }
620
621 return talphaEff;
622}
623
624
625Foam::tmp<Foam::scalarField> Foam::multiphaseMixtureThermo::alphahe
626(
627 const label patchi
628) const
629{
631
632 tmp<scalarField> talphaEff
633 (
634 phasei().boundaryField()[patchi]
635 *phasei().thermo().alphahe(patchi)
636 );
637
638 for (++phasei; phasei != phases_.end(); ++phasei)
639 {
640 talphaEff.ref() +=
641 phasei().boundaryField()[patchi]
642 *phasei().thermo().alphahe(patchi);
643 }
644
645 return talphaEff;
646}
647
648
649Foam::tmp<Foam::volScalarField> Foam::multiphaseMixtureThermo::kappaEff
650(
651 const volScalarField& alphat
652) const
653{
654 auto phasei = phases_.cbegin();
655
656 tmp<volScalarField> tkappaEff(phasei()*phasei().thermo().kappaEff(alphat));
657
658 for (++phasei; phasei != phases_.cend(); ++phasei)
659 {
660 tkappaEff.ref() += phasei()*phasei().thermo().kappaEff(alphat);
661 }
662
663 return tkappaEff;
664}
665
666
667Foam::tmp<Foam::scalarField> Foam::multiphaseMixtureThermo::kappaEff
668(
669 const scalarField& alphat,
670 const label patchi
671) const
672{
673 auto phasei = phases_.cbegin();
674
675 tmp<scalarField> tkappaEff
676 (
677 phasei().boundaryField()[patchi]
678 *phasei().thermo().kappaEff(alphat, patchi)
679 );
680
681 for (++phasei; phasei != phases_.cend(); ++phasei)
682 {
683 tkappaEff.ref() +=
684 phasei().boundaryField()[patchi]
685 *phasei().thermo().kappaEff(alphat, patchi);
686 }
687
688 return tkappaEff;
689}
690
691
692Foam::tmp<Foam::volScalarField> Foam::multiphaseMixtureThermo::alphaEff
693(
694 const volScalarField& alphat
695) const
696{
697 auto phasei = phases_.cbegin();
698
699 tmp<volScalarField> talphaEff(phasei()*phasei().thermo().alphaEff(alphat));
700
701 for (++phasei; phasei != phases_.cend(); ++phasei)
702 {
703 talphaEff.ref() += phasei()*phasei().thermo().alphaEff(alphat);
704 }
705
706 return talphaEff;
707}
708
709
710Foam::tmp<Foam::scalarField> Foam::multiphaseMixtureThermo::alphaEff
711(
712 const scalarField& alphat,
713 const label patchi
714) const
715{
716 auto phasei = phases_.cbegin();
717
718 tmp<scalarField> talphaEff
719 (
720 phasei().boundaryField()[patchi]
721 *phasei().thermo().alphaEff(alphat, patchi)
722 );
723
724 for (++phasei; phasei != phases_.cend(); ++phasei)
725 {
726 talphaEff.ref() +=
727 phasei().boundaryField()[patchi]
728 *phasei().thermo().alphaEff(alphat, patchi);
729 }
730
731 return talphaEff;
732}
733
734
735Foam::tmp<Foam::volScalarField> Foam::multiphaseMixtureThermo::rCv() const
736{
737 auto phasei = phases_.cbegin();
738
740
741 for (++phasei; phasei != phases_.cend(); ++phasei)
742 {
743 trCv.ref() += phasei()/phasei().thermo().Cv();
744 }
745
746 return trCv;
747}
748
749
750Foam::tmp<Foam::surfaceScalarField>
752{
754 (
756 (
758 (
759 "surfaceTensionForce",
760 mesh_.time().timeName(),
761 mesh_
762 ),
763 mesh_,
764 dimensionedScalar(dimensionSet(1, -2, -2, 0, 0), Zero)
765 )
766 );
767
768 surfaceScalarField& stf = tstf.ref();
769 stf.setOriented();
770
771 forAllConstIters(phases_, phase1)
772 {
773 const phaseModel& alpha1 = *phase1;
774
775 auto phase2 = phase1;
776
777 for (++phase2; phase2 != phases_.cend(); ++phase2)
778 {
779 const phaseModel& alpha2 = *phase2;
780
781 auto sigma = sigmas_.cfind(interfacePair(alpha1, alpha2));
782
783 if (!sigma.good())
784 {
786 << "Cannot find interface " << interfacePair(alpha1, alpha2)
787 << " in list of sigma values"
788 << exit(FatalError);
789 }
790
791 stf += dimensionedScalar("sigma", dimSigma_, *sigma)
793 (
796 );
797 }
798 }
799
800 return tstf;
801}
802
803
805{
806 const Time& runTime = mesh_.time();
807
808 const dictionary& alphaControls = mesh_.solverDict("alpha");
809 label nAlphaSubCycles(alphaControls.get<label>("nAlphaSubCycles"));
810 scalar cAlpha(alphaControls.get<scalar>("cAlpha"));
811
812 volScalarField& alpha = phases_.first();
813
814 if (nAlphaSubCycles > 1)
815 {
816 surfaceScalarField rhoPhiSum(0.0*rhoPhi_);
817 dimensionedScalar totalDeltaT = runTime.deltaT();
818
819 for
820 (
822 !(++alphaSubCycle).end();
823 )
824 {
825 solveAlphas(cAlpha);
826 rhoPhiSum += (runTime.deltaT()/totalDeltaT)*rhoPhi_;
827 }
828
829 rhoPhi_ = rhoPhiSum;
830 }
831 else
832 {
833 solveAlphas(cAlpha);
834 }
835}
836
837
838Foam::tmp<Foam::surfaceVectorField> Foam::multiphaseMixtureThermo::nHatfv
839(
840 const volScalarField& alpha1,
842) const
843{
844 /*
845 // Cell gradient of alpha
846 volVectorField gradAlpha =
847 alpha2*fvc::grad(alpha1) - alpha1*fvc::grad(alpha2);
848
849 // Interpolated face-gradient of alpha
850 surfaceVectorField gradAlphaf = fvc::interpolate(gradAlpha);
851 */
852
853 surfaceVectorField gradAlphaf
854 (
857 );
858
859 // Face unit interface normal
860 return gradAlphaf/(mag(gradAlphaf) + deltaN_);
861}
862
863
864Foam::tmp<Foam::surfaceScalarField> Foam::multiphaseMixtureThermo::nHatf
865(
866 const volScalarField& alpha1,
868) const
869{
870 // Face unit interface normal flux
871 return nHatfv(alpha1, alpha2) & mesh_.Sf();
872}
873
874
875// Correction for the boundary condition on the unit normal nHat on
876// walls to produce the correct contact angle.
877
878// The dynamic contact angle is calculated from the component of the
879// velocity on the direction of the interface, parallel to the wall.
880
881void Foam::multiphaseMixtureThermo::correctContactAngle
882(
883 const phaseModel& alpha1,
884 const phaseModel& alpha2,
886) const
887{
888 const volScalarField::Boundary& gbf
890
891 const fvBoundaryMesh& boundary = mesh_.boundary();
892
893 forAll(boundary, patchi)
894 {
896 {
899
900 vectorField& nHatPatch = nHatb[patchi];
901
902 vectorField AfHatPatch
903 (
904 mesh_.Sf().boundaryField()[patchi]
905 /mesh_.magSf().boundaryField()[patchi]
906 );
907
908 const auto tp =
909 acap.thetaProps().cfind(interfacePair(alpha1, alpha2));
910
911 if (!tp.good())
912 {
914 << "Cannot find interface " << interfacePair(alpha1, alpha2)
915 << "\n in table of theta properties for patch "
916 << acap.patch().name()
917 << exit(FatalError);
918 }
919
920 const bool matched = (tp.key().first() == alpha1.name());
921
922 const scalar theta0 = degToRad(tp().theta0(matched));
923 scalarField theta(boundary[patchi].size(), theta0);
924
925 const scalar uTheta = tp().uTheta();
926
927 // Calculate the dynamic contact angle if required
928 if (uTheta > SMALL)
929 {
930 const scalar thetaA = degToRad(tp().thetaA(matched));
931 const scalar thetaR = degToRad(tp().thetaR(matched));
932
933 // Calculated the component of the velocity parallel to the wall
934 vectorField Uwall
935 (
936 U_.boundaryField()[patchi].patchInternalField()
937 - U_.boundaryField()[patchi]
938 );
939 Uwall -= (AfHatPatch & Uwall)*AfHatPatch;
940
941 // Find the direction of the interface parallel to the wall
942 vectorField nWall
943 (
944 nHatPatch - (AfHatPatch & nHatPatch)*AfHatPatch
945 );
946
947 // Normalise nWall
948 nWall /= (mag(nWall) + SMALL);
949
950 // Calculate Uwall resolved normal to the interface parallel to
951 // the interface
952 scalarField uwall(nWall & Uwall);
953
954 theta += (thetaA - thetaR)*tanh(uwall/uTheta);
955 }
956
957
958 // Reset nHatPatch to correspond to the contact angle
959
960 scalarField a12(nHatPatch & AfHatPatch);
961
962 scalarField b1(cos(theta));
963
964 scalarField b2(nHatPatch.size());
965
966 forAll(b2, facei)
967 {
968 b2[facei] = cos(acos(a12[facei]) - theta[facei]);
969 }
970
971 scalarField det(1.0 - a12*a12);
972
973 scalarField a((b1 - a12*b2)/det);
974 scalarField b((b2 - a12*b1)/det);
975
976 nHatPatch = a*AfHatPatch + b*nHatPatch;
977
978 nHatPatch /= (mag(nHatPatch) + deltaN_.value());
979 }
980 }
981}
982
983
984Foam::tmp<Foam::volScalarField> Foam::multiphaseMixtureThermo::K
985(
986 const phaseModel& alpha1,
987 const phaseModel& alpha2
988) const
989{
990 tmp<surfaceVectorField> tnHatfv = nHatfv(alpha1, alpha2);
991
992 correctContactAngle(alpha1, alpha2, tnHatfv.ref().boundaryFieldRef());
993
994 // Simple expression for curvature
995 return -fvc::div(tnHatfv & mesh_.Sf());
996}
997
998
999Foam::tmp<Foam::volScalarField>
1001{
1002 auto tnearInt = volScalarField::New
1003 (
1004 "nearInterface",
1006 mesh_,
1008 );
1009
1010 for (const phaseModel& phase : phases_)
1011 {
1012 tnearInt.ref() =
1013 max(tnearInt(), pos0(phase - 0.01)*pos0(0.99 - phase));
1014 }
1015
1016 return tnearInt;
1017}
1018
1019
1020void Foam::multiphaseMixtureThermo::solveAlphas
1021(
1022 const scalar cAlpha
1023)
1024{
1025 static label nSolves(-1);
1026 ++nSolves;
1027
1028 const word alphaScheme("div(phi,alpha)");
1029 const word alpharScheme("div(phirb,alpha)");
1030
1031 surfaceScalarField phic(mag(phi_/mesh_.magSf()));
1032 phic = min(cAlpha*phic, max(phic));
1033
1034 PtrList<surfaceScalarField> alphaPhiCorrs(phases_.size());
1035
1036 int phasei = 0;
1037 for (phaseModel& alpha : phases_)
1038 {
1039 alphaPhiCorrs.set
1040 (
1041 phasei,
1043 (
1044 phi_.name() + alpha.name(),
1045 fvc::flux
1046 (
1047 phi_,
1048 alpha,
1049 alphaScheme
1050 )
1051 )
1052 );
1053
1054 surfaceScalarField& alphaPhiCorr = alphaPhiCorrs[phasei];
1055
1056 for (phaseModel& alpha2 : phases_)
1057 {
1058 if (&alpha2 == &alpha) continue;
1059
1061
1062 alphaPhiCorr += fvc::flux
1063 (
1065 alpha,
1067 );
1068 }
1069
1071 (
1072 1.0/mesh_.time().deltaTValue(),
1074 alpha,
1075 phi_,
1076 alphaPhiCorr,
1077 zeroField(),
1078 zeroField(),
1079 oneField(),
1080 zeroField(),
1081 true
1082 );
1083
1084 ++phasei;
1085 }
1086
1087 MULES::limitSum(alphaPhiCorrs);
1088
1089 rhoPhi_ = Zero;
1090
1091 volScalarField sumAlpha
1092 (
1093 IOobject
1094 (
1095 "sumAlpha",
1096 mesh_.time().timeName(),
1097 mesh_
1098 ),
1099 mesh_,
1101 );
1102
1103
1105
1106
1107 phasei = 0;
1108
1109 for (phaseModel& alpha : phases_)
1110 {
1111 surfaceScalarField& alphaPhi = alphaPhiCorrs[phasei];
1112 alphaPhi += upwind<scalar>(mesh_, phi_).flux(alpha);
1113
1115 (
1116 IOobject
1117 (
1118 "Sp",
1119 mesh_.time().timeName(),
1120 mesh_
1121 ),
1122 mesh_,
1124 );
1125
1127 (
1128 IOobject
1129 (
1130 "Su",
1131 mesh_.time().timeName(),
1132 mesh_
1133 ),
1134 // Divergence term is handled explicitly to be
1135 // consistent with the explicit transport solution
1136 divU*min(alpha, scalar(1))
1137 );
1138
1139 {
1140 const scalarField& dgdt = alpha.dgdt();
1141
1142 forAll(dgdt, celli)
1143 {
1144 if (dgdt[celli] < 0.0 && alpha[celli] > 0.0)
1145 {
1146 Sp[celli] += dgdt[celli]*alpha[celli];
1147 Su[celli] -= dgdt[celli]*alpha[celli];
1148 }
1149 else if (dgdt[celli] > 0.0 && alpha[celli] < 1.0)
1150 {
1151 Sp[celli] -= dgdt[celli]*(1.0 - alpha[celli]);
1152 }
1153 }
1154 }
1155
1156 for (const phaseModel& alpha2 : phases_)
1157 {
1158 if (&alpha2 == &alpha) continue;
1159
1160 const scalarField& dgdt2 = alpha2.dgdt();
1161
1162 forAll(dgdt2, celli)
1163 {
1164 if (dgdt2[celli] > 0.0 && alpha2[celli] < 1.0)
1165 {
1166 Sp[celli] -= dgdt2[celli]*(1.0 - alpha2[celli]);
1167 Su[celli] += dgdt2[celli]*alpha[celli];
1168 }
1169 else if (dgdt2[celli] < 0.0 && alpha2[celli] > 0.0)
1170 {
1171 Sp[celli] += dgdt2[celli]*alpha2[celli];
1172 }
1173 }
1174 }
1175
1177 (
1179 alpha,
1180 alphaPhi,
1181 Sp,
1182 Su
1183 );
1184
1185 rhoPhi_ += fvc::interpolate(alpha.thermo().rho())*alphaPhi;
1186
1187 Info<< alpha.name() << " volume fraction, min, max = "
1188 << alpha.weightedAverage(mesh_.V()).value()
1189 << ' ' << min(alpha).value()
1190 << ' ' << max(alpha).value()
1191 << endl;
1192
1193 sumAlpha += alpha;
1194
1195 ++phasei;
1196 }
1197
1198 Info<< "Phase-sum volume fraction, min, max = "
1199 << sumAlpha.weightedAverage(mesh_.V()).value()
1200 << ' ' << min(sumAlpha).value()
1201 << ' ' << max(sumAlpha).value()
1202 << endl;
1203
1204 calcAlphas();
1205}
1206
1207
1208// ************************************************************************* //
CGAL::Exact_predicates_exact_constructions_kernel K
surfaceScalarField phic(mixture.cAlpha() *mag(alphaPhic/mesh.magSf()))
MULES: Multidimensional universal limiter for explicit solution.
volScalarField & he
Definition YEEqn.H:52
const volScalarField & alpha1
phaseModel & phase1
const volScalarField & alpha2
phaseModel & phase2
const dimensionSet & dimensions() const noexcept
Return dimensions.
dimensioned< Type > weightedAverage(const DimensionedField< scalar, GeoMesh > &weights, const label comm=UPstream::worldComm) const
Return the global weighted average.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
DimensionedField< scalar, volMesh > Internal
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef().
GeometricBoundaryField< vector, fvsPatchField, surfaceMesh > Boundary
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
@ NO_REGISTER
Do not request registration (bool: false).
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
An STL-conforming const_iterator.
Definition LList.H:461
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
dimensionedScalar deltaT() const
Return time step.
Definition TimeStateI.H:61
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
T & first()
Access first element of the list, position [0].
Definition UList.H:957
const_iterator cend() const noexcept
Return const_iterator to end traversing the constant UList.
Definition UListI.H:468
const_iterator cbegin() const noexcept
Return const_iterator to begin traversing the constant UList.
Definition UListI.H:424
Contact-angle boundary condition for multi-phase interface-capturing simulations. Used in conjunction...
volScalarField alpha_
Laminar thermal diffusivity [kg/m/s].
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
A fvBoundaryMesh is a fvPatch list with a reference to the associated fvMesh, with additional search ...
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the interpolation weighting factors.
tmp< volScalarField > nearInterface() const
Indicator of the proximity of the interface.
virtual word thermoName() const
Return the name of the thermo physics.
tmp< surfaceScalarField > surfaceTensionForce() const
virtual tmp< volScalarField > Cv() const
Heat capacity at constant volume [J/kg/K].
virtual tmp< volScalarField > rho() const
Density [kg/m^3].
virtual tmp< volScalarField > W() const
Molecular weight [kg/kmol].
virtual tmp< volScalarField > gamma() const
Gamma = Cp/Cv [].
virtual bool incompressible() const
Return true if the equation of state is incompressible.
virtual tmp< volScalarField > Cpv() const
Heat capacity at constant pressure/volume [J/kg/K].
virtual tmp< volScalarField > kappaEff(const volScalarField &alphat) const
Effective thermal diffusivity of mixture [J/m/s/K].
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.
virtual tmp< volScalarField > kappa() const
Thermal diffusivity for temperature of mixture [J/m/s/K].
tmp< volScalarField > rCv() const
Return the phase-averaged reciprocal Cv.
virtual tmp< volScalarField > Cp() const
Heat capacity at constant pressure [J/kg/K].
virtual void correct()
Update properties.
multiphaseMixtureThermo(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
virtual volScalarField & he()
Enthalpy/Internal energy [J/kg].
virtual tmp< volScalarField > hc() const
Chemical enthalpy [J/kg].
virtual tmp< volScalarField > alphaEff(const volScalarField &alphat) const
Effective thermal diffusivity of mixture [J/m/s/K].
virtual tmp< volScalarField > nu() const
Kinematic viscosity of mixture [m^2/s].
void correctRho(const volScalarField &dp)
Update densities for given pressure change.
virtual tmp< volScalarField > alphahe() const
Thermal diffusivity for energy of mixture [kg/m/s].
virtual tmp< volScalarField > CpByCpv() const
Heat capacity ratio [].
void solve()
Solve for the mixture phase-fractions.
virtual bool isochoric() const
Return true if the equation of state is isochoric.
const Time & time() const noexcept
Return time registry.
A class representing the concept of a field of 1 used to avoid unnecessary manipulations for objects ...
Definition oneField.H:50
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phaseModel.H:56
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phase.H:53
const dimensionedScalar & rho() const
Return const-access to phase1 density.
Definition phase.H:151
Basic thermodynamic properties based on compressibility.
Definition psiThermo.H:54
volScalarField mu_
Dynamic viscosity [kg/m/s].
Definition psiThermo.H:68
volScalarField psi_
Compressibility [s^2/m^2].
Definition psiThermo.H:63
Lookup type of boundary radiation properties.
Definition lookup.H:60
bool end() const
Return true if the number of sub-cycles has been reached.
Perform a subCycleTime on a field.
Definition subCycle.H:137
A class for managing temporary objects.
Definition tmp.H:75
Upwind differencing scheme class.
Definition upwind.H:55
A class for handling words, derived from Foam::string.
Definition word.H:66
A class representing the concept of a field of 0 used to avoid unnecessary manipulations for objects ...
Definition zeroField.H:50
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
volScalarField & p
thermo correct()
const tmp< volScalarField > & tCv
Definition EEqn.H:5
const volScalarField & Cv
Definition EEqn.H:8
const scalar gamma
Definition EEqn.H:9
const tmp< volScalarField > & tCp
Definition EEqn.H:4
const volScalarField & Cp
Definition EEqn.H:7
faceListList boundary
dynamicFvMesh & mesh
engineTime & runTime
surfaceScalarField phir(fvc::flux(UdmModel.Udm()))
word alpharScheme("div(phirb,alpha)")
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition error.H:688
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
Calculate the divergence of the given field.
Calculate the face-flux of the given field.
Calculate the gradient of the given field.
Calculate the mesh motion flux and convert fluxes from absolute to relative and back.
Calculate the snGrad of the given volField.
const cellShapeList & cells
volScalarField alphaEff("alphaEff", turbulence->nu()/Pr+alphat)
zeroField divU
Definition alphaSuSp.H:3
zeroField Su
Definition alphaSuSp.H:1
zeroField Sp
Definition alphaSuSp.H:2
word timeName
Definition getTimeIndex.H:3
label phasei
Definition pEqn.H:27
kappaEff
Definition TEqn.H:10
const expr V(m.psi().mesh().V())
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
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)
void limitSum(UPtrList< scalarField > &phiPsiCorrs)
Definition MULES.C:27
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
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.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition fvcSnGrad.C:40
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcDiv.C:42
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Definition fvcFlux.C:27
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition fvcMeshPhi.C:183
Namespace for OpenFOAM.
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
dimensionedScalar det(const dimensionedSphericalTensor &dt)
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
dimensionedScalar pos0(const dimensionedScalar &ds)
List< label > labelList
A List of labels.
Definition List.H:62
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar tanh(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
void Sp(fvMatrix< typename Expr::value_type > &m, const Expr2 &mult, const Expr &expression)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensionedScalar cbrt(const dimensionedScalar &ds)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &f1, const label comm)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
void Su(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar acos(const dimensionedScalar &ds)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
volScalarField & alpha
tmp< volScalarField > trho
volScalarField & h
volScalarField & b
psiReactionThermo & thermo
scalar T0
volScalarField & e
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
const dictionary & alphaControls
label nAlphaSubCycles(alphaControls.get< label >("nAlphaSubCycles"))
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition stdFoam.H:235
surfaceScalarField alphaPhi(phi.name()+alpha1.name(), fvc::flux(phi, alpha1, alphaScheme))
Unit conversion functions.