Loading...
Searching...
No Matches
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
97 auto tmpPhi = surfaceScalarField::New
98 (
99 "phi",
101 fvc::interpolate(iter()()) * iter()->phi()
102 );
103
104 for (++iter; iter != phaseModels.cend(); ++iter)
105 {
106 tmpPhi.ref() += fvc::interpolate(iter()()) * iter()->phi();
107 }
108
109 return tmpPhi;
110}
111
112
113void Foam::multiphaseInterSystem::generatePairs(const dictTable& modelDicts)
114{
115 forAllConstIters(modelDicts, iter)
116 {
117 const phasePairKey& key = iter.key();
118
119 // pair already exists
120 if (phasePairs_.found(key))
121 {
122 // do nothing ...
123 }
124 else if (key.ordered())
125 {
126 // New ordered pair
127 phasePairs_.insert
128 (
129 key,
131 (
133 (
134 phaseModels_[key.first()],
135 phaseModels_[key.second()]
136 )
137 )
138 );
139 }
140 else
141 {
142 // New unordered pair
143 phasePairs_.insert
144 (
145 key,
147 (
148 new phasePair
149 (
150 phaseModels_[key.first()],
151 phaseModels_[key.second()]
152 )
154 );
155 }
156 }
157}
158
159
161{
162 forAllConstIters(phaseModels_, phaseIter1)
163 {
164 forAllConstIters(phaseModels_, phaseIter2)
165 {
166 if (phaseIter1()->name() != phaseIter2()->name())
167 {
169 (
170 phaseIter1()->name(),
171 phaseIter2()->name(),
172 true
173 );
174
175 phasePairKey keyInverse
176 (
177 phaseIter2()->name(),
178 phaseIter1()->name(),
179 true
180 );
181
182 if
183 (
184 !totalPhasePairs_.found(key)
185 && !totalPhasePairs_.found(keyInverse)
186 )
187 {
188 totalPhasePairs_.set
189 (
190 key,
192 (
193 new phasePair
194 (
195 phaseModels_[key.first()],
196 phaseModels_[key.second()]
197 )
198 )
199 );
200 }
201 }
203 }
204}
205
206
207// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
208
210(
211 const fvMesh& mesh
212)
213:
214 basicThermo(mesh, word::null, phasePropertiesName),
215 mesh_(mesh),
216 mu_
217 (
219 (
220 "mu",
221 mesh_.time().timeName(),
222 mesh_
223 ),
224 mesh_,
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{
263 rhoPhi_.setOriented();
264 phi_.setOriented();
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 * * * * * * * * * * * * * * * //
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
428 tmp<volScalarField> tmpCp
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
469 tmp<volScalarField> tmpCv
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
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 }
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 }
803 tmpAlpha.ref() += alphat;
804
805 return tmpAlpha;
806}
807
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
863 tmp<volScalarField> tmpNu
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
906}
907
912}
913
918}
919
922{
923 return turb_->nuEff();
924}
925
926
928{
929 return
931 this->kappa() + this->Cp()*turb_->mut()/Prt_
932 );
933}
934
935
937Foam::multiphaseInterSystem::kappaEff(const label patchi) const
938{
939 const scalarField Cp(this->Cp()().boundaryField()[patchi]);
940 const scalarField kappaEffp
941 (
942 this->kappa(patchi) + Cp*turb_->mut(patchi)/Prt_.value()
943 );
944
945 return tmp<scalarField>::New(kappaEffp);
946}
947
948
953
954
956Foam::multiphaseInterSystem::alphaEff(const label patchi) const
957{
958 return (this->alpha(patchi) + turb_->mut(patchi))/Prt_.value();
959}
960
965}
966
971}
972
977}
978
988 forAllIters(phaseModels_, iter)
989 {
990 iter()->correct();
991 }
992
993 calcMu();
994}
995
996
998{
999 forAllIters(phaseModels_, iter)
1001 iter()->correctTurbulence();
1002 }
1003}
1004
1005
1008{
1009 return phaseModels_;
1010}
1011
1012
1018
1019
1025
1026
1029{
1030 return totalPhasePairs_;
1031}
1032
1033
1035{
1036 forAllConstIters(phaseModels_, iter)
1037 {
1038 if (!iter()->thermo().incompressible())
1039 {
1040 return false;
1042 }
1043
1044 return true;
1045}
1046
1048bool Foam::multiphaseInterSystem::incompressible(const word phaseName) const
1049{
1050 return phaseModels_[phaseName]->thermo().incompressible();
1051}
1052
1053
1055{
1056 forAllConstIters(phaseModels_, iter)
1057 {
1058 if (!iter()->thermo().isochoric())
1059 {
1060 return false;
1062 }
1063
1064 return true;
1065}
1066
1067
1069{
1070 return mesh_;
1071}
1072
1073
1076{
1077 auto tstf = surfaceScalarField::New
1078 (
1079 "surfaceTensionForce",
1081 mesh_,
1082 dimensionedScalar({1, -2, -2, 0, 0, 0}, Zero)
1083 );
1084
1085 auto& stf = tstf.ref();
1086 stf.setOriented();
1087
1088 if (surfaceTensionModels_.size())
1089 {
1090 forAllConstIters(phaseModels_, iter1)
1091 {
1092 const volScalarField& alpha1 = iter1()();
1093
1094 auto iter2 = iter1;
1095
1096 for (++iter2; iter2 != phaseModels_.cend(); ++iter2)
1097 {
1098 const volScalarField& alpha2 = iter2()();
1099
1100 stf +=
1102 (
1103 surfaceTensionCoeff
1104 (
1105 phasePairKey(iter1()->name(), iter2()->name())
1106 )
1107 )
1109 (
1112 );
1113 }
1115 }
1116
1117 return tstf;
1118}
1119
1120
1121Foam::tmp<Foam::volVectorField> Foam::multiphaseInterSystem::U() const
1122{
1123 auto tstf = volVectorField::New
1124 (
1125 "U",
1127 mesh_,
1129 );
1130
1131 auto& stf = tstf.ref();
1132
1133 forAllConstIters(phaseModels_, iter)
1134 {
1135 stf += iter()() * iter()->U();
1137
1138 return tstf;
1139}
1140
1141
1150(
1151 const word& key
1152) const
1153{
1154 return 1.0/(phaseModels_[key]->thermo().rho());
1155}
1156
1157
1159{
1160 const scalarField& Vc = mesh_.V();
1161 scalarField& Udiag = UEqn.diag();
1162
1163 forAllConstIters(phaseModels_, iteri)
1164 {
1165 const multiphaseInter::phaseModel& phasei = iteri()();
1166
1167 auto iterk = iteri;
1168
1169 for (++iterk; iterk != phaseModels_.cend(); ++iterk)
1170 {
1171 if (iteri()().name() != iterk()().name())
1172 {
1173 const multiphaseInter::phaseModel& phasek = iterk()();
1174
1175 // Phase i and k
1176 const phasePairKey keyik
1177 (
1178 phasei.name(),
1179 phasek.name(),
1180 false
1181 );
1182
1183 if (interfacePorousModelTable_.found(keyik))
1184 {
1185 autoPtr<porousModel>& interfacePtr =
1186 interfacePorousModelTable_[keyik];
1187
1188 Udiag += Vc*interfacePtr->S();
1190 }
1191 }
1192 }
1193}
1194
1195
1196Foam::tmp<Foam::volScalarField> Foam::multiphaseInterSystem::K
1197(
1198 const volScalarField& alpha1,
1199 const volScalarField& alpha2
1200) const
1201{
1203
1204 // Simple expression for curvature
1205 return -fvc::div(tnHatfv.ref() & mesh_.Sf());
1206}
1207
1208
1210(
1211 const volScalarField& alpha1,
1212 const volScalarField& alpha2
1213) const
1214{
1215 return
1216 (
1217 pos(alpha1 - 0.1)*pos(0.9 - alpha1)
1218 *pos(alpha2 - 0.1)*pos(0.9 - alpha2)
1219 );
1220}
1221
1222
1225{
1226 auto tnearInt = volScalarField::New
1227 (
1228 "nearInterface",
1230 mesh_,
1232 );
1233
1234 auto& nearInt = tnearInt.ref();
1235
1236 forAllConstIters(phaseModels_, iter1)
1237 {
1238 const volScalarField& alpha1 = iter1()();
1239
1240 auto iter2 = iter1;
1241
1242 for (++iter2; iter2 != phaseModels_.cend(); ++iter2)
1243 {
1244 const volScalarField& alpha2 = iter2()();
1245
1246 nearInt +=
1247 (
1248 pos(alpha1 - 0.1)*pos(0.9 - alpha1)
1249 *pos(alpha2 - 0.1)*pos(0.9 - alpha2)
1250 );
1252 }
1253
1254 return tnearInt;
1255}
1256
1257
1258Foam::tmp<Foam::volVectorField> Foam::multiphaseInterSystem::nVolHatfv
1259(
1260 const volScalarField& alpha1,
1261 const volScalarField& alpha2
1262) const
1263{
1264 const volScalarField alpha1m
1265 (
1267 );
1268
1269 const volScalarField alpha2m
1270 (
1272 );
1273
1274 const volVectorField gradAlphaf
1275 (
1276 alpha2m*(fvc::grad(alpha1m))
1277 - alpha1m*(fvc::grad(alpha2m))
1278 );
1279
1280 const dimensionedScalar deltaN
1281 (
1282 "deltaN",
1283 1e-8/cbrt(average(mesh_.V()))
1285
1286 // Face unit interface normal
1287 return gradAlphaf/(mag(gradAlphaf) + deltaN);
1288}
1289
1290
1292(
1293 const volScalarField& alpha1,
1294 const volScalarField& alpha2
1295) const
1296{
1297
1298 const volScalarField alpha1b
1299 (
1301 );
1302
1303 const volScalarField alpha2b
1304 (
1305 clamp(alpha2, zero_one{})
1306 );
1307
1308 surfaceVectorField gradAlphaf
1309 (
1311 - fvc::interpolate(alpha1b)*fvc::interpolate(fvc::grad(alpha2b))
1312 );
1313
1314 const dimensionedScalar deltaN
1315 (
1316 "deltaN",
1317 1e-8/cbrt(average(mesh_.V()))
1319
1320 // Face unit interface normal
1321 return gradAlphaf/(mag(gradAlphaf) + deltaN);
1322}
1323
1324
1326(
1327 const volScalarField& alpha1,
1328 const volScalarField& alpha2
1329) const
1330{
1331 // Face unit interface normal flux
1332 return nHatfv(alpha1, alpha2) & mesh_.Sf();
1333}
1334
1335
1337{
1338 if (regIOobject::read())
1339 {
1340 return true;
1341 }
1342
1343 return false;
1344}
1345
1346
1347// ************************************************************************* //
CGAL::Exact_predicates_exact_constructions_kernel K
bool found
const volScalarField & alpha1
const volScalarField & alpha2
static tmp< GeometricField< scalar, fvsPatchField, surfaceMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvsPatchField< scalar >::calculatedType())
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:152
constexpr const_iterator cend() const noexcept
const_iterator to signal the end (for any HashTable)
@ NO_REGISTER
Do not request registration (bool: false).
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
Abstract base-class for fluid and solid thermodynamic properties.
Definition basicThermo.H:62
virtual const volScalarField & T() const
Temperature [K].
virtual volScalarField & p()
Pressure [Pa].
basicThermo(const basicThermo &)=delete
No copy construct.
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.
static const dictionary null
An empty dictionary, which is also the parent for all dictionaries.
Definition dictionary.H:487
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
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.
compressibleTurbulenceModel * turb_
Turbulence model.
virtual tmp< volScalarField > W() const
Molecular weight [kg/kmol] of the mixture.
dimensionedScalar Prt_
Turbulent Prandt number.
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< volScalarField > K(const volScalarField &alpha1, const volScalarField &alpha2) const
Interface curvature.
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].
HashTable< autoPtr< phasePair >, phasePairKey, phasePairKey::hash > phasePairTable
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.
virtual tmp< volScalarField > coeffs(const word &key) const
Return coefficients (1/rho).
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.
phasePairTable phasePairs_
Phase pairs.
surfaceScalarField phi_
Mixture total volumetric flux.
HashTable< autoPtr< multiphaseInter::phaseModel > > phaseModelTable
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.
multiphaseInterSystem(const fvMesh &mesh)
Construct from fvMesh.
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.
tmp< surfaceScalarField > nHatf(const volScalarField &alpha1, const volScalarField &alpha2) const
Interface normal surface vector.
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 tmp< volScalarField > kappaEff(const volScalarField &kappat) const
Effective thermal diffusivity for temperature.
phaseModelTable phaseModels_
Phase models.
HashTable< dictionary, phasePairKey, phasePairKey::hash > dictTable
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.
phasePairTable totalPhasePairs_
Total ordered phase pairs in the system.
static autoPtr< phaseModel > New(const multiphaseInterSystem &fluid, const word &phaseName)
Definition phaseModel.C:65
const word & name() const
The name of this phase.
Definition phaseModel.H:122
An ordered or unorder pair of phase names. Typically specified as follows.
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition phasePair.H:52
virtual bool read()
Read object.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
A class for handling words, derived from Foam::string.
Definition word.H:66
Represents 0/1 range or concept. Used for tagged dispatch or clamping.
Definition pTraits.H:51
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
const volScalarField & T
fvVectorMatrix & UEqn
Definition UEqn.H:13
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:688
Calculate the divergence of the given field.
Calculate the gradient of the given field.
Calculate the snGrad of the given volField.
auto & name
const cellShapeList & cells
word timeName
Definition getTimeIndex.H:3
label phasei
Definition pEqn.H:27
kappaEff
Definition TEqn.H:10
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
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
constexpr auto key(const Type &t) noexcept
Helper function to return the enum value.
Namespace for OpenFOAM.
const dimensionSet dimViscosity
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
dimensionedScalar pos(const dimensionedScalar &ds)
List< word > wordList
List of word.
Definition fileName.H:60
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
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
const dimensionSet dimVelocity
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
fvMatrix< vector > fvVectorMatrix
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensionedScalar cbrt(const dimensionedScalar &ds)
const dimensionSet dimVolume(pow3(dimLength))
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimDensity
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &f1, const label comm)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
volScalarField & nu
volScalarField & alpha
psiReactionThermo & thermo
scalar T0
volScalarField & e
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition stdFoam.H:214
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition stdFoam.H:235