multiphaseInterSystem.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) 2017-2022 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
29#include "surfaceTensionModel.H"
30#include "porousModel.H"
31
32#include "HashPtrTable.H"
33
34#include "surfaceInterpolate.H"
35#include "fvcGrad.H"
36#include "fvcSnGrad.H"
37#include "fvcDiv.H"
38#include "fvMatrix.H"
39
44
45// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
46
47namespace Foam
48{
50}
51
52/* * * * * * * * * * * * * * * private static data * * * * * * * * * * * * * */
53
54const Foam::word
56
57// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
58
60{
61 mu_ = mu()();
62}
63
64
67(
68 const wordList& phaseNames
69) const
70{
71 phaseModelTable phaseModels;
72
73 for (const word& phaseName : phaseNames)
74 {
75 phaseModels.insert
76 (
77 phaseName,
79 (
80 *this,
81 phaseName
82 )
83 );
84 }
85
86 return phaseModels;
87}
88
89
91(
92 const phaseModelTable& phaseModels
93) const
94{
95 auto iter = phaseModels.cbegin();
96
98 (
99 "phi",
100 fvc::interpolate(iter()()) * iter()->phi()
101 );
102
103 for (++iter; iter != phaseModels.cend(); ++iter)
104 {
105 tmpPhi.ref() += fvc::interpolate(iter()()) * iter()->phi();
106 }
107
108 return tmpPhi;
109}
110
111
113{
114 forAllConstIters(modelDicts, iter)
115 {
116 const phasePairKey& key = iter.key();
117
118 // pair already exists
119 if (phasePairs_.found(key))
120 {
121 // do nothing ...
122 }
123 else if (key.ordered())
124 {
125 // New ordered pair
126 phasePairs_.insert
127 (
128 key,
130 (
132 (
133 phaseModels_[key.first()],
134 phaseModels_[key.second()]
135 )
136 )
137 );
138 }
139 else
140 {
141 // New unordered pair
142 phasePairs_.insert
143 (
144 key,
146 (
147 new phasePair
148 (
149 phaseModels_[key.first()],
150 phaseModels_[key.second()]
151 )
152 )
153 );
154 }
155 }
156}
157
158
160{
161 forAllConstIters(phaseModels_, phaseIter1)
162 {
163 forAllConstIters(phaseModels_, phaseIter2)
164 {
165 if (phaseIter1()->name() != phaseIter2()->name())
166 {
167 phasePairKey key
168 (
169 phaseIter1()->name(),
170 phaseIter2()->name(),
171 true
172 );
173
174 phasePairKey keyInverse
175 (
176 phaseIter2()->name(),
177 phaseIter1()->name(),
178 true
179 );
180
181 if
182 (
183 !totalPhasePairs_.found(key)
184 && !totalPhasePairs_.found(keyInverse)
185 )
186 {
187 totalPhasePairs_.set
188 (
189 key,
191 (
192 new phasePair
193 (
194 phaseModels_[key.first()],
195 phaseModels_[key.second()]
196 )
197 )
198 );
199 }
200 }
201 }
202 }
203}
204
205
206// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
207
209(
210 const fvMesh& mesh
211)
212:
213 basicThermo(mesh, word::null, phasePropertiesName),
214 mesh_(mesh),
215 mu_
216 (
218 (
219 "mu",
220 mesh_.time().timeName(),
221 mesh_
222 ),
223 mesh_,
225 calculatedFvPatchScalarField::typeName
226 ),
227 phaseNames_(get<wordList>("phases")),
228 phi_
229 (
231 (
232 "phi",
233 mesh_.time().timeName(),
234 mesh_,
235 IOobject::READ_IF_PRESENT,
236 IOobject::AUTO_WRITE
237 ),
238 mesh_,
240 ),
241 rhoPhi_
242 (
244 (
245 "rhoPhi",
246 mesh_.time().timeName(),
247 mesh_
248 ),
249 mesh_,
251 ),
252 phaseModels_(generatePhaseModels(phaseNames_)),
253 phasePairs_(),
254 totalPhasePairs_(),
255 Prt_
256 (
257 dimensionedScalar::getOrAddToDict
258 (
259 "Prt", *this, 1.0
260 )
261 )
262{
265
266 // sub models
267 if (found("surfaceTension"))
268 {
270 (
271 "surfaceTension",
273 );
274 }
275 if (found("interfacePorous"))
276 {
278 (
279 "interfacePorous",
280 mesh_,
282 );
283 }
284
285 // Total phase pair
287
288 // Update mu_
289 calcMu();
290}
291
292
293// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
294
296{}
297
298
299// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
300
302(
303 const volScalarField& p,
304 const volScalarField& T
305) const
306{
308 return nullptr;
309}
310
311
313(
314 const scalarField& p,
315 const scalarField& T,
316 const labelList& cells
317) const
318{
320 return nullptr;
321}
322
323
325(
326 const scalarField& p,
327 const scalarField& T,
328 const label patchI
329) const
330{
332 return nullptr;
333}
334
335
337{
338 auto iter = phaseModels_.cbegin();
339
340 tmp<volScalarField> tAlphaHc
341 (
342 iter()() * iter()->hc()
343 );
344
345 for (++iter; iter != phaseModels_.cend(); ++iter)
346 {
347 tAlphaHc.ref() += iter()() * iter()->hc();
348 }
349
350 return tAlphaHc;
351}
352
353
355(
356 const scalarField& e,
357 const scalarField& p,
358 const scalarField& T0,
359 const labelList& cells
360) const
361{
363 return nullptr;
364}
365
366
368(
369 const scalarField& e,
370 const scalarField& p,
371 const scalarField& T0,
372 const label patchI
373) const
374{
376 return nullptr;
377}
378
379
381{
382 auto iter = phaseModels_.cbegin();
383
385 (
386 iter()() * iter()->rho()
387 );
388
389 for (++iter; iter != phaseModels_.cend(); ++iter)
390 {
391 tmpRho.ref() += iter()() * iter()->rho();
392 }
393
394 return tmpRho;
395}
396
397
399(
400 const label patchI
401) const
402{
403 auto iter = phaseModels_.cbegin();
404
405 tmp<scalarField> tmpRho
406 (
407 iter()().boundaryField()[patchI]
408 * iter()->rho()().boundaryField()[patchI]
409 );
410
411 for (++iter; iter != phaseModels_.cend(); ++iter)
412 {
413 tmpRho.ref() +=
414 (
415 iter()().boundaryField()[patchI]
416 * iter()->rho()().boundaryField()[patchI]
417 );
418 }
419
420 return tmpRho;
421}
422
423
425{
426 auto iter = phaseModels_.cbegin();
427
429 (
430 iter()() * iter()->Cp()
431 );
432
433 for (++iter; iter != phaseModels_.cend(); ++iter)
434 {
435 tmpCp.ref() += iter()() * iter()->Cp();
436 }
437
438 return tmpCp;
439}
440
441
443(
444 const scalarField& p,
445 const scalarField& T,
446 const label patchI
447) const
448{
449 auto iter = phaseModels_.cbegin();
450
451 tmp<scalarField> tmpCp
452 (
453 iter()() * iter()->Cp(p, T, patchI)
454 );
455
456 for (++iter; iter != phaseModels_.cend(); ++iter)
457 {
458 tmpCp.ref() += iter()() * iter()->Cp(p, T, patchI);
459 }
460
461 return tmpCp;
462}
463
464
466{
467 auto iter = phaseModels_.cbegin();
468
470 (
471 iter()() * iter()->Cv()
472 );
473
474 for (++iter; iter != phaseModels_.cend(); ++iter)
475 {
476 tmpCv.ref() += iter()() * iter()->Cv();
477 }
478
479 return tmpCv;
480}
481
482
484(
485 const scalarField& p,
486 const scalarField& T,
487 const label patchI
488) const
489{
490 auto iter = phaseModels_.cbegin();
491
492 tmp<scalarField> tmpCv
493 (
494 iter()() * iter()->Cv(p, T, patchI)
495 );
496
497 for (++iter; iter != phaseModels_.cend(); ++iter)
498 {
499 tmpCv.ref() += iter()() * iter()->Cv(p, T, patchI);
500 }
501
502 return tmpCv;
503}
504
505
507(
508 const scalarField& p,
509 const scalarField& T,
510 const labelList& cells
511) const
512{
514 return nullptr;
515}
516
517
519{
520 auto iter = phaseModels_.cbegin();
521
523 (
524 iter()() * iter()->Cp()
525 );
526
528 (
529 iter()() * iter()->Cv()
530 );
531
532 for (++iter; iter != phaseModels_.cend(); ++iter)
533 {
534 tmpCp.ref() += iter()() * iter()->Cp();
535 tmpCv.ref() += iter()() * iter()->Cv();
536 }
537
538 return (tmpCp/tmpCv);
539}
540
541
543(
544 const scalarField& p,
545 const scalarField& T,
546 const label patchI
547) const
548{
549 return
550 (
551 gamma()().boundaryField()[patchI]
552 );
553}
554
555
557{
558 auto iter = phaseModels_.cbegin();
559
561 (
562 iter()() * iter()->Cpv()
563 );
564
565 for (++iter; iter != phaseModels_.cend(); ++iter)
566 {
567 tmpCpv.ref() += iter()() * iter()->Cpv();
568 }
569
570 return tmpCpv;
571}
572
573
575(
576 const scalarField& p,
577 const scalarField& T,
578 const label patchI
579) const
580{
581 auto iter = phaseModels_.cbegin();
582
583 tmp<scalarField> tmpCpv
584 (
585 iter()() * iter()->Cpv(p, T, patchI)
586 );
587
588 for (++iter; iter != phaseModels_.cend(); ++iter)
589 {
590 tmpCpv.ref() += iter()() * iter()->Cpv(p, T, patchI);
591 }
592
593 return tmpCpv;
594}
595
596
598{
599 auto iter = phaseModels_.cbegin();
600
601 tmp<volScalarField> tmpCpByCpv
602 (
603 iter()() * iter()->CpByCpv()
604 );
605
606 for (++iter; iter != phaseModels_.cend(); ++iter)
607 {
608 tmpCpByCpv.ref() += iter()() * iter()->CpByCpv();
609 }
610
611 return tmpCpByCpv;
612}
613
614
616(
617 const scalarField& p,
618 const scalarField& T,
619 const label patchI
620) const
621{
622 auto iter = phaseModels_.cbegin();
623
624 tmp<scalarField> tmpCpv
625 (
626 iter()().boundaryField()[patchI]
627 * iter()->CpByCpv(p, T, patchI)
628 );
629
630 for (++iter; iter != phaseModels_.cend(); ++iter)
631 {
632 tmpCpv.ref() +=
633 (
634 iter()().boundaryField()[patchI]
635 * iter()->CpByCpv(p, T, patchI)
636 );
637 }
638
639 return tmpCpv;
640}
641
642
644{
646 return nullptr;
647}
648
649
651{
652 auto iter = phaseModels_.cbegin();
653
654 tmp<volScalarField> tmpkappa
655 (
656 iter()() * iter()->kappa()
657 );
658
659 for (++iter; iter != phaseModels_.cend(); ++iter)
660 {
661 tmpkappa.ref() += iter()() * iter()->kappa();
662 }
663
664 return tmpkappa;
665}
666
667
669(
670 const label patchI
671) const
672{
673 auto iter = phaseModels_.cbegin();
674
675 tmp<scalarField> tmpKappa
676 (
677 iter()().boundaryField()[patchI]
678 * iter()->kappa(patchI)
679 );
680
681 for (++iter; iter != phaseModels_.cend(); ++iter)
682 {
683 tmpKappa.ref() +=
684 (
685 iter()().boundaryField()[patchI]
686 * iter()->kappa(patchI)
687 );
688 }
689
690 return tmpKappa;
691}
692
693
695{
696 phaseModelTable::const_iterator phaseModelIter = phaseModels_.begin();
697
698 tmp<volScalarField> talphaEff
699 (
700 phaseModelIter()()*phaseModelIter()->alphahe()
701 );
702
703 for (; phaseModelIter != phaseModels_.end(); ++phaseModelIter)
704 {
705 talphaEff.ref() += phaseModelIter()()*phaseModelIter()->alphahe();
706 }
707
708 return talphaEff;
709}
710
711
713(
714 const label patchi
715) const
716{
717 phaseModelTable::const_iterator phaseModelIter = phaseModels_.begin();
718
719 tmp<scalarField> talphaEff
720 (
721 phaseModelIter()().boundaryField()[patchi]
722 *phaseModelIter()->alphahe(patchi)
723 );
724
725 for (; phaseModelIter != phaseModels_.end(); ++phaseModelIter)
726 {
727 talphaEff.ref() +=
728 phaseModelIter()().boundaryField()[patchi]
729 *phaseModelIter()->alphahe(patchi);
730 }
731
732 return talphaEff;
733}
734
735
737(
738 const volScalarField& kappat
739) const
740{
741 tmp<volScalarField> kappaEff(kappa() + kappat);
742 kappaEff.ref().rename("kappaEff");
743 return kappaEff;
744}
745
746
748(
749 const scalarField& kappat,
750 const label patchI
751) const
752{
753 return kappa(patchI) + kappat;
754}
755
756
758(
759 const volScalarField& alphat
760) const
761{
762 auto iter = phaseModels_.cbegin();
763
764 tmp<volScalarField> tmpAlpha
765 (
766 iter()() * iter()->alpha()
767 );
768
769 for (++iter; iter != phaseModels_.cend(); ++iter)
770 {
771 tmpAlpha.ref() += iter()() * iter()->alpha();
772 }
773
774 tmpAlpha.ref() += alphat;
775
776 return tmpAlpha;
777}
778
779
781(
782 const scalarField& alphat,
783 const label patchI
784) const
785{
786 auto iter = phaseModels_.cbegin();
787
788 tmp<scalarField> tmpAlpha
789 (
790 iter()().boundaryField()[patchI]
791 * iter()->alpha(patchI)
792 );
793
794 for (++iter; iter != phaseModels_.cend(); ++iter)
795 {
796 tmpAlpha.ref() +=
797 (
798 iter()().boundaryField()[patchI]
799 * iter()->alpha(patchI)
800 );
801 }
802
803 tmpAlpha.ref() += alphat;
804
805 return tmpAlpha;
806}
807
808
810{
811 return Prt_;
812}
813
814
816{
817 auto iter = phaseModels_.cbegin();
818
820 (
821 iter()() * iter()->mu()
822 );
823
824 for (++iter; iter != phaseModels_.cend(); ++iter)
825 {
826 tmpMu.ref() += iter()() * iter()->mu();
827 }
828
829 return tmpMu;
830}
831
832
834(
835 const label patchI
836) const
837{
838 auto iter = phaseModels_.cbegin();
839
840 tmp<scalarField> tmpMu
841 (
842 iter()().boundaryField()[patchI]
843 * iter()->mu(patchI)
844 );
845
846 for (++iter; iter != phaseModels_.cend(); ++iter)
847 {
848 tmpMu.ref() +=
849 (
850 iter()().boundaryField()[patchI]
851 * iter()->mu(patchI)
852 );
853 }
854
855 return tmpMu;
856}
857
858
860{
861 auto iter = phaseModels_.cbegin();
862
864 (
865 iter()() * iter()->nu()
866 );
867
868 for (++iter; iter != phaseModels_.cend(); ++iter)
869 {
870 tmpNu.ref() += iter()() * iter()->nu();
871 }
872
873 return tmpNu;
874}
875
876
878(
879 const label patchI
880) const
881{
882 auto iter = phaseModels_.cbegin();
883
884 tmp<scalarField> tmpNu
885 (
886 iter()().boundaryField()[patchI]
887 * iter()->nu(patchI)
888 );
889
890 for (++iter; iter != phaseModels_.cend(); ++iter)
891 {
892 tmpNu.ref() +=
893 (
894 iter()().boundaryField()[patchI]
895 * iter()->nu(patchI)
896 );
897 }
898
899 return tmpNu;
900}
901
902
904{
905 return turb_->mut();
906}
907
908
910{
911 return turb_->muEff();
912}
913
914
916{
917 return turb_->nut();
918}
919
920
922{
923 return turb_->nuEff();
924}
925
926
928{
930 (
931 this->kappa() + this->Cp()*turb_->mut()/Prt_
932 );
933
935}
936
937
940{
941 const scalarField Cp(this->Cp()().boundaryField()[patchi]);
942 const scalarField kappaEffp
943 (
944 this->kappa(patchi) + Cp*turb_->mut(patchi)/Prt_.value()
945 );
946
947 return tmp<scalarField>::New(kappaEffp);
948}
949
950
952{
953 return this->alpha() + turb_->mut()/Prt_;
954}
955
956
959{
960 return (this->alpha(patchi) + turb_->mut(patchi))/Prt_.value();
961}
962
963
965{
966 return phi_;
967}
968
969
971{
972 return phi_;
973}
974
975
977{
978 return rhoPhi_;
979}
980
981
983{
984 return rhoPhi_;
985}
986
987
989{
990 forAllIters(phaseModels_, iter)
991 {
992 iter()->correct();
993 }
994
995 calcMu();
996}
997
998
1000{
1001 forAllIters(phaseModels_, iter)
1002 {
1003 iter()->correctTurbulence();
1004 }
1005}
1006
1007
1010{
1011 return phaseModels_;
1012}
1013
1014
1017{
1018 return phaseModels_;
1019}
1020
1021
1024{
1025 return totalPhasePairs_;
1026}
1027
1028
1031{
1032 return totalPhasePairs_;
1033}
1034
1035
1037{
1038 forAllConstIters(phaseModels_, iter)
1039 {
1040 if (!iter()->thermo().incompressible())
1041 {
1042 return false;
1043 }
1044 }
1045
1046 return true;
1047}
1048
1049
1051{
1052 return phaseModels_[phaseName]->thermo().incompressible();
1053}
1054
1055
1057{
1058 forAllConstIters(phaseModels_, iter)
1059 {
1060 if (!iter()->thermo().isochoric())
1061 {
1062 return false;
1063 }
1064 }
1065
1066 return true;
1067}
1068
1069
1071{
1072 return mesh_;
1073}
1074
1075
1078{
1080 (
1081 IOobject
1082 (
1083 "surfaceTensionForce",
1084 mesh_.time().timeName(),
1085 mesh_
1086 ),
1087 mesh_,
1088 dimensionedScalar({1, -2, -2, 0, 0, 0}, Zero)
1089 );
1090
1091 auto& stf = tstf.ref();
1092 stf.setOriented();
1093
1094 if (surfaceTensionModels_.size())
1095 {
1096 forAllConstIters(phaseModels_, iter1)
1097 {
1098 const volScalarField& alpha1 = iter1()();
1099
1100 auto iter2 = iter1;
1101
1102 for (++iter2; iter2 != phaseModels_.cend(); ++iter2)
1103 {
1104 const volScalarField& alpha2 = iter2()();
1105
1106 stf +=
1108 (
1109 surfaceTensionCoeff
1110 (
1111 phasePairKey(iter1()->name(), iter2()->name())
1112 )
1113 )
1115 (
1118 );
1119 }
1120 }
1121 }
1122
1123 return tstf;
1124}
1125
1126
1128{
1129 auto tstf = tmp<volVectorField>::New
1130 (
1131 IOobject
1132 (
1133 "U",
1134 mesh_.time().timeName(),
1135 mesh_
1136 ),
1137 mesh_,
1139 );
1140
1141 auto& stf = tstf.ref();
1142
1143 forAllConstIters(phaseModels_, iter)
1144 {
1145 stf += iter()() * iter()->U();
1146 }
1147
1148 return tstf;
1149}
1150
1151
1154{
1155 return surfaceTensionModels_[key]->sigma();
1156}
1157
1158
1160(
1161 const word& key
1162) const
1163{
1164 return 1.0/(phaseModels_[key]->thermo().rho());
1165}
1166
1167
1169{
1170 const scalarField& Vc = mesh_.V();
1171 scalarField& Udiag = UEqn.diag();
1172
1173 forAllConstIters(phaseModels_, iteri)
1174 {
1175 const multiphaseInter::phaseModel& phasei = iteri()();
1176
1177 auto iterk = iteri;
1178
1179 for (++iterk; iterk != phaseModels_.cend(); ++iterk)
1180 {
1181 if (iteri()().name() != iterk()().name())
1182 {
1183 const multiphaseInter::phaseModel& phasek = iterk()();
1184
1185 // Phase i and k
1186 const phasePairKey keyik
1187 (
1188 phasei.name(),
1189 phasek.name(),
1190 false
1191 );
1192
1193 if (interfacePorousModelTable_.found(keyik))
1194 {
1195 autoPtr<porousModel>& interfacePtr =
1196 interfacePorousModelTable_[keyik];
1197
1198 Udiag += Vc*interfacePtr->S();
1199 }
1200 }
1201 }
1202 }
1203}
1204
1205
1207(
1208 const volScalarField& alpha1,
1209 const volScalarField& alpha2
1210) const
1211{
1212 tmp<surfaceVectorField> tnHatfv = nHatfv(alpha1, alpha2);
1213
1214 // Simple expression for curvature
1215 return -fvc::div(tnHatfv.ref() & mesh_.Sf());
1216}
1217
1218
1220(
1221 const volScalarField& alpha1,
1222 const volScalarField& alpha2
1223) const
1224{
1225 return
1226 (
1227 pos(alpha1 - 0.1)*pos(0.9 - alpha1)
1228 *pos(alpha2 - 0.1)*pos(0.9 - alpha2)
1229 );
1230}
1231
1232
1235{
1236 auto tnearInt = tmp<volScalarField>::New
1237 (
1238 IOobject
1239 (
1240 "nearInterface",
1241 mesh_.time().timeName(),
1242 mesh_
1243 ),
1244 mesh_,
1246 );
1247
1248 auto& nearInt = tnearInt.ref();
1249
1250 forAllConstIters(phaseModels_, iter1)
1251 {
1252 const volScalarField& alpha1 = iter1()();
1253
1254 auto iter2 = iter1;
1255
1256 for (++iter2; iter2 != phaseModels_.cend(); ++iter2)
1257 {
1258 const volScalarField& alpha2 = iter2()();
1259
1260 nearInt +=
1261 (
1262 pos(alpha1 - 0.1)*pos(0.9 - alpha1)
1263 *pos(alpha2 - 0.1)*pos(0.9 - alpha2)
1264 );
1265 }
1266 }
1267
1268 return tnearInt;
1269}
1270
1271
1273(
1274 const volScalarField& alpha1,
1275 const volScalarField& alpha2
1276) const
1277{
1278 const volScalarField alpha1m
1279 (
1280 min(max(alpha1, scalar(0)), scalar(1))
1281 );
1282
1283 const volScalarField alpha2m
1284 (
1285 min(max(alpha2, scalar(0)), scalar(1))
1286 );
1287
1288 const volVectorField gradAlphaf
1289 (
1290 alpha2m*(fvc::grad(alpha1m))
1291 - alpha1m*(fvc::grad(alpha2m))
1292 );
1293
1294 const dimensionedScalar deltaN
1295 (
1296 "deltaN",
1297 1e-8/cbrt(average(mesh_.V()))
1298 );
1299
1300 // Face unit interface normal
1301 return gradAlphaf/(mag(gradAlphaf) + deltaN);
1302}
1303
1304
1306(
1307 const volScalarField& alpha1,
1308 const volScalarField& alpha2
1309) const
1310{
1311
1312 const volScalarField alpha1b
1313 (
1314 min(max(alpha1, scalar(0)), scalar(1))
1315 );
1316
1317 const volScalarField alpha2b
1318 (
1319 min(max(alpha2, scalar(0)), scalar(1))
1320 );
1321
1322 surfaceVectorField gradAlphaf
1323 (
1325 - fvc::interpolate(alpha1b)*fvc::interpolate(fvc::grad(alpha2b))
1326 );
1327
1328 const dimensionedScalar deltaN
1329 (
1330 "deltaN",
1331 1e-8/cbrt(average(mesh_.V()))
1332 );
1333
1334 // Face unit interface normal
1335 return gradAlphaf/(mag(gradAlphaf) + deltaN);
1336}
1337
1338
1340(
1341 const volScalarField& alpha1,
1342 const volScalarField& alpha2
1343) const
1344{
1345 // Face unit interface normal flux
1346 return nHatfv(alpha1, alpha2) & mesh_.Sf();
1347}
1348
1349
1351{
1352 if (regIOobject::read())
1353 {
1354 return true;
1355 }
1356
1357 return false;
1358}
1359
1360
1361// ************************************************************************* //
CGAL::Exact_predicates_exact_constructions_kernel K
bool found
surfaceScalarField & phi
const volScalarField & alpha1
const volScalarField & alpha2
void setOriented(const bool oriented=true) noexcept
Set the oriented flag.
const_iterator cbegin() const
const_iterator set to the beginning of the HashTable
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition: HashTableI.H:180
constexpr const_iterator cend() const noexcept
const_iterator to signal the end (for any HashTable)
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
Abstract base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:66
A const_iterator for iterating across on values.
Definition: bitSet.H:505
const dictionary & coeffs() const
Return const dictionary of the model.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const surfaceScalarField & nHatf() const
scalarField & diag()
Definition: lduMatrix.C:192
void generatePairsTable()
Generate pair table.
tmp< surfaceScalarField > surfaceTensionForce() const
Calculate surface tension of the mixture.
tmp< volScalarField > nearInterface() const
Near Interface of alpha'n.
virtual ~multiphaseInterSystem()
Destructor.
virtual tmp< volScalarField > W() const
Molecular weight [kg/kmol] of the mixture.
static const word phasePropertiesName
Default name of the phase properties dictionary.
const dimensionedScalar & Prt() const
Return Prandt number.
virtual tmp< volScalarField > Cp() const
Return Cp of the mixture.
void generatePairs(const dictTable &modelDicts)
Generate pairs.
tmp< volScalarField > muEff() const
Return the effective dynamic viscosity.
const fvMesh & mesh_
Reference to the mesh.
virtual tmp< volScalarField > surfaceTensionCoeff(const phasePairKey &key) const
Return the surface tension coefficient.
const surfaceScalarField & phi() const
Constant access to the total flux.
virtual void correct()
Correct the mixture thermos.
tmp< volVectorField > U() const
Mixture U.
const phaseModelTable & phases() const
Constant access the phases.
virtual tmp< volScalarField > alphahe() const
Thermal diffusivity for energy of mixture [kg/m/s].
virtual tmp< volScalarField > mu() const
Dynamic viscosity of mixture [kg/m/s].
tmp< volScalarField > mut() const
Return the turbulent dynamic viscosity.
interfacePorousModelTable interfacePorousModelTable_
Interface porous models.
virtual tmp< volScalarField > nu() const
Kinematic viscosity of mixture [m^2/s].
const surfaceScalarField & rhoPhi() const
Constant access to the mixture mass flux.
surfaceTensionModelTable surfaceTensionModels_
Surface tension models.
const phasePairTable & totalPhasePairs() const
Constant access the total phase pairs.
virtual tmp< volScalarField > CpByCpv() const
Heat capacity ratio [].
virtual void correctTurbulence()
Correct the turbulence.
tmp< volScalarField > alphaEff() const
Effective thermal turbulent diffusivity of mixture [kg/m/s].
void generatePairsAndSubModels(const word &modelName, HashTable< autoPtr< modelType >, phasePairKey, phasePairKey::hash > &models)
Generate pairs and sub-model tables.
HashTable< autoPtr< multiphaseInter::phaseModel > > generatePhaseModels(const wordList &names) const
Generate the phases.
surfaceScalarField phi_
Mixture total volumetric flux.
virtual tmp< volScalarField > hc() const
Chemical enthalpy of the mixture [J/kg].
tmp< volScalarField > nuEff() const
Return the effective kinematic viscosity.
void calcMu()
Calculate and return the laminar viscosity.
virtual tmp< volScalarField > gamma() const
Gamma = Cp/Cv [].
virtual volScalarField & he()
Return access to the internal energy field [J/Kg].
virtual tmp< volScalarField > kappa() const
Thermal diffusivity for temperature 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.
volScalarField mu_
Dynamic viscocity.
void addInterfacePorosity(fvVectorMatrix &UEqn)
Add interface porosity on phasePair.
const fvMesh & mesh() const
Return mesh.
virtual tmp< volScalarField > rho() const
Return the mixture density.
tmp< volScalarField > nut() const
Return the turbulent kinematic viscosity.
surfaceScalarField rhoPhi_
Mixture total mass flux.
tmp< surfaceVectorField > nHatfv(const volScalarField &alpha1, const volScalarField &alpha2) const
Interface normal surface vector.
virtual tmp< scalarField > rhoEoS(const scalarField &p, const scalarField &T, const labelList &cells) const
Density from pressure and temperature.
tmp< volScalarField > kappaEff() const
Effective thermal turbulent diffusivity for temperature.
virtual tmp< volScalarField > Cv() const
Return Cv of the mixture.
tmp< volVectorField > nVolHatfv(const volScalarField &alpha1, const volScalarField &alpha2) const
Interface normal volField vector.
virtual bool read()
Read base phaseProperties dictionary.
virtual bool isochoric() const
Return true if the equation of state is isochoric for all phasses.
virtual tmp< volScalarField > Cpv() const
Heat capacity at constant pressure/volume [J/kg/K].
virtual bool incompressible() const
Return true if the equation of state is incompressible for all.
tmp< surfaceScalarField > generatePhi(const HashTable< autoPtr< multiphaseInter::phaseModel > > &phaseModels) const
Generate the mixture flux.
const word & name() const
The name of this phase.
Definition: phaseModel.H:114
An ordered or unorder pair of phase names. Typically specified as follows.
Definition: phasePairKey.H:68
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition: phasePair.H:56
virtual bool read()
Read object.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
volScalarField & p
const volScalarField & T
fvVectorMatrix & UEqn
Definition: UEqn.H:13
const volScalarField & mu
const volScalarField & Cv
Definition: EEqn.H:8
const scalar gamma
Definition: EEqn.H:9
const volScalarField & Cp
Definition: EEqn.H:7
dynamicFvMesh & mesh
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:517
Calculate the divergence of the given field.
Calculate the gradient of the given field.
Calculate the snGrad of the given volField.
const cellShapeList & cells
word timeName
Definition: getTimeIndex.H:3
label phasei
Definition: pEqn.H:27
kappaEff
Definition: TEqn.H:10
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:54
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:47
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
Namespace for OpenFOAM.
const dimensionSet dimViscosity
dimensionedScalar pos(const dimensionedScalar &ds)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
const dimensionSet dimless
Dimensionless.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
const dimensionSet dimVelocity
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
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:33
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
dimensionedScalar cbrt(const dimensionedScalar &ds)
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:60
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
const dimensionSet dimDensity
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
volScalarField & nu
volScalarField & alpha
scalar T0
Definition: createFields.H:22
volScalarField & e
Definition: createFields.H:11
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:260
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition: stdFoam.H:278