Loading...
Searching...
No Matches
MomentumTransferPhaseSystem.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) 2015-2018 OpenFOAM Foundation
9 Copyright (C) 2020-2023 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
31#include "BlendedInterfacialModel.H"
32#include "dragModel.H"
33#include "virtualMassModel.H"
34#include "liftModel.H"
35#include "wallLubricationModel.H"
36#include "turbulentDispersionModel.H"
37
38#include "HashPtrTable.H"
39
40#include "fvmDdt.H"
41#include "fvmDiv.H"
42#include "fvmSup.H"
43#include "fvcAverage.H"
44#include "fvcDdt.H"
45#include "fvcDiv.H"
46#include "fvcFlux.H"
47#include "fvcSnGrad.H"
48#include "fvMatrix.H"
49
50
51// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
52
53template<class BasePhaseSystem>
55Foam::MomentumTransferPhaseSystem<BasePhaseSystem>::Kd
56(
57 const phasePairKey& key
58) const
59{
60 if (dragModels_.found(key))
61 {
62 return dragModels_[key]->K();
63 }
64
65 return volScalarField::New
66 (
67 IOobject::scopedName(dragModel::typeName, "K"),
68 IOobject::NO_REGISTER,
69 this->mesh_,
70 dimensionedScalar(dragModel::dimK)
71 );
72}
73
74
75template<class BasePhaseSystem>
77Foam::MomentumTransferPhaseSystem<BasePhaseSystem>::Kdf
78(
79 const phasePairKey& key
80) const
81{
82 if (dragModels_.found(key))
83 {
84 return dragModels_[key]->Kf();
85 }
86
87 return surfaceScalarField::New
88 (
89 IOobject::scopedName(dragModel::typeName, "K"),
90 IOobject::NO_REGISTER,
91 this->mesh_,
92 dimensionedScalar(dragModel::dimK)
93 );
94}
95
96
97template<class BasePhaseSystem>
99Foam::MomentumTransferPhaseSystem<BasePhaseSystem>::Vm
100(
101 const phasePairKey& key
102) const
103{
104 if (virtualMassModels_.found(key))
105 {
106 return virtualMassModels_[key]->K();
107 }
108
109 return volScalarField::New
110 (
111 IOobject::scopedName(virtualMassModel::typeName, "K"),
112 IOobject::NO_REGISTER,
113 this->mesh_,
114 dimensionedScalar(virtualMassModel::dimK)
115 );
116}
117
118
119template<class BasePhaseSystem>
121addMassTransferMomentumTransfer(phaseSystem::momentumTransferTable& eqns) const
122{
124 (
125 phaseSystem::phasePairTable,
126 this->phasePairs_,
127 phasePairIter
128 )
129 {
130 const phasePair& pair(phasePairIter());
131
132 if (pair.ordered())
133 {
134 continue;
135 }
136
137 // Note that the phase UEqn contains a continuity error term, which
138 // implicitly adds a mass transfer term of fvm::Sp(dmdt, U). These
139 // additions do not include this term.
140
141 const volScalarField dmdt(this->dmdt(pair));
142
143 if (!pair.phase1().stationary())
144 {
145 fvVectorMatrix& eqn = *eqns[pair.phase1().name()];
146 const volScalarField dmdt21(posPart(dmdt));
147
148 eqn += dmdt21*pair.phase2().U() - fvm::Sp(dmdt21, eqn.psi());
149 }
150
151 if (!pair.phase2().stationary())
152 {
153 fvVectorMatrix& eqn = *eqns[pair.phase2().name()];
154 const volScalarField dmdt12(negPart(dmdt));
155
156 eqn -= dmdt12*pair.phase1().U() - fvm::Sp(dmdt12, eqn.psi());
157 }
159}
160
161
162// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
163
164template<class BasePhaseSystem>
167(
168 const fvMesh& mesh
169)
170:
171 BasePhaseSystem(mesh)
172{
173 this->generatePairsAndSubModels
174 (
175 "drag",
176 dragModels_
177 );
178
179 this->generatePairsAndSubModels
180 (
181 "virtualMass",
182 virtualMassModels_
183 );
184
185 this->generatePairsAndSubModels
186 (
187 "lift",
188 liftModels_
189 );
190
191 this->generatePairsAndSubModels
192 (
193 "wallLubrication",
194 wallLubricationModels_
195 );
196
197 this->generatePairsAndSubModels
198 (
199 "turbulentDispersion",
200 turbulentDispersionModels_
201 );
202
204 (
206 dragModels_,
207 dragModelIter
208 )
209 {
210 const phasePair& pair(this->phasePairs_[dragModelIter.key()]);
211
212 Kds_.set
213 (
214 pair,
216 (
217 IOobject::groupName("Kd", pair.name()),
218 dragModelIter()->K()
219 )
220 );
221
222 Kdfs_.set
223 (
224 pair,
226 (
227 IOobject::groupName("Kdf", pair.name()),
228 dragModelIter()->Kf()
229 )
230 );
231 }
232
234 (
236 virtualMassModels_,
237 virtualMassModelIter
238 )
239 {
240 const phasePair& pair(this->phasePairs_[virtualMassModelIter.key()]);
241
242 Vms_.set
243 (
244 pair,
246 (
247 IOobject::groupName("Vm", pair.name()),
248 virtualMassModelIter()->K()
249 )
250 );
251
252 Vmfs_.set
253 (
254 pair,
256 (
257 IOobject::groupName("Vmf", pair.name()),
258 virtualMassModelIter()->Kf()
259 )
260 );
262}
263
264
265// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
266
267template<class BasePhaseSystem>
271
272
273// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
274
275template<class BasePhaseSystem>
278{
279 // Create a momentum transfer matrix for each phase
281 (
283 );
284
285 phaseSystem::momentumTransferTable& eqns = eqnsPtr();
286
287 forAll(this->movingPhases(), movingPhasei)
288 {
289 const phaseModel& phase = this->movingPhases()[movingPhasei];
290
291 eqns.set
292 (
293 phase.name(),
295 );
296 }
297
298 // Update the drag coefficients
300 (
301 dragModelTable,
302 dragModels_,
303 dragModelIter
304 )
305 {
306 *Kds_[dragModelIter.key()] = dragModelIter()->K();
307 *Kdfs_[dragModelIter.key()] = dragModelIter()->Kf();
308 }
309
310 // Add the implicit part of the drag force
311 forAllConstIter(KdTable, Kds_, KdIter)
312 {
313 const volScalarField& K(*KdIter());
314 const phasePair& pair(this->phasePairs_[KdIter.key()]);
315
316 forAllConstIter(phasePair, pair, iter)
317 {
318 if (!iter().stationary())
319 {
320 fvVectorMatrix& eqn = *eqns[iter().name()];
321
322 eqn -= fvm::Sp(K, eqn.psi());
323 }
324 }
325 }
326
327 // Update the virtual mass coefficients
329 (
330 virtualMassModelTable,
331 virtualMassModels_,
332 virtualMassModelIter
333 )
334 {
335 *Vms_[virtualMassModelIter.key()] = virtualMassModelIter()->K();
336 *Vmfs_[virtualMassModelIter.key()] = virtualMassModelIter()->Kf();
337 }
338
339 // Add the virtual mass force
340 forAllConstIter(VmTable, Vms_, VmIter)
341 {
342 const volScalarField& Vm(*VmIter());
343 const phasePair& pair(this->phasePairs_[VmIter.key()]);
344
345 forAllConstIter(phasePair, pair, iter)
346 {
347 const phaseModel& phase = iter();
348 const phaseModel& otherPhase = iter.otherPhase();
349
350 if (!phase.stationary())
351 {
352 fvVectorMatrix& eqn = *eqns[phase.name()];
353
354 const volVectorField& U = eqn.psi();
355 const tmp<surfaceScalarField> tphi(phase.phi());
356 const surfaceScalarField& phi = tphi();
357
358 eqn -=
359 Vm
360 *(
361 fvm::ddt(U)
362 + fvm::div(phi, U)
363 - fvm::Sp(fvc::div(phi), U)
364 - otherPhase.DUDt()
365 )
366 + this->MRF_.DDt(Vm, U - otherPhase.U());
367 }
368 }
369 }
370
371 // Add the source term due to mass transfer
372 addMassTransferMomentumTransfer(eqns);
374 return eqnsPtr;
375}
376
377
378template<class BasePhaseSystem>
381{
382 // Create a momentum transfer matrix for each phase
384 (
386 );
387
388 phaseSystem::momentumTransferTable& eqns = eqnsPtr();
389
390 forAll(this->movingPhases(), movingPhasei)
391 {
392 const phaseModel& phase = this->movingPhases()[movingPhasei];
393
394 eqns.set
395 (
396 phase.name(),
398 );
399 }
400
401 // Create U & grad(U) fields
402 PtrList<fvVectorMatrix> UgradUs(this->phaseModels_.size());
403 forAll(this->phaseModels_, phasei)
404 {
405 const phaseModel& phase = this->phaseModels_[phasei];
406
407 if (!phase.stationary())
408 {
409 const tmp<volVectorField> tU(phase.U());
410 const volVectorField& U = tU();
411
412 UgradUs.set
413 (
414 phasei,
416 (
417 fvm::div(phase.phi(), U)
418 - fvm::Sp(fvc::div(phase.phi()), U)
419 + this->MRF().DDt(U)
420 )
421 );
422 }
423 }
424
425 // Add the virtual mass force
426 forAllConstIter(VmTable, Vms_, VmIter)
427 {
428 const volScalarField& Vm(*VmIter());
429 const phasePair& pair(this->phasePairs_[VmIter.key()]);
430
431 forAllConstIter(phasePair, pair, iter)
432 {
433 const phaseModel& phase = iter();
434 const phaseModel& otherPhase = iter.otherPhase();
435
436 if (!phase.stationary())
437 {
438 *eqns[phase.name()] -=
439 Vm
440 *(
441 UgradUs[phase.index()]
442 - (UgradUs[otherPhase.index()] & otherPhase.U())
443 );
444 }
445 }
446 }
447
448 // Add the source term due to mass transfer
449 addMassTransferMomentumTransfer(eqns);
451 return eqnsPtr;
452}
453
454
455template<class BasePhaseSystem>
458{
459 PtrList<surfaceScalarField> AFfs(this->phaseModels_.size());
460
461 // Add the implicit part of the drag force
462 forAllConstIter(KdfTable, Kdfs_, KdfIter)
463 {
464 const surfaceScalarField& Kf(*KdfIter());
465 const phasePair& pair(this->phasePairs_[KdfIter.key()]);
466
467 forAllConstIter(phasePair, pair, iter)
468 {
469 this->addField(iter(), "AFf", Kf, AFfs);
470 }
471 }
472
473 // Add the implicit part of the virtual mass force
474 forAllConstIter(VmfTable, Vmfs_, VmfIter)
475 {
476 const surfaceScalarField& Vmf(*VmfIter());
477 const phasePair& pair(this->phasePairs_[VmfIter.key()]);
478
479 forAllConstIter(phasePair, pair, iter)
480 {
481 this->addField(iter(), "AFf", byDt(Vmf), AFfs);
482 }
483 }
484
485 if (this->fillFields_)
486 {
487 this->fillFields("AFf", dimDensity/dimTime, AFfs);
488 }
490 return AFfs;
491}
492
493
494template<class BasePhaseSystem>
497(
499)
500{
501 PtrList<surfaceScalarField> phiFs(this->phaseModels_.size());
502
503 // Add the lift force
505 (
506 liftModelTable,
507 liftModels_,
508 liftModelIter
509 )
510 {
511 const volVectorField F(liftModelIter()->F<vector>());
512 const phasePair& pair(this->phasePairs_[liftModelIter.key()]);
513
514 this->addField
515 (
516 pair.phase1(),
517 "phiF",
518 fvc::flux(rAUs[pair.phase1().index()]*F),
519 phiFs
520 );
521 this->addField
522 (
523 pair.phase2(),
524 "phiF",
525 - fvc::flux(rAUs[pair.phase2().index()]*F),
526 phiFs
527 );
528 }
529
530 // Add the wall lubrication force
532 (
533 wallLubricationModelTable,
534 wallLubricationModels_,
535 wallLubricationModelIter
536 )
537 {
538 const volVectorField F(wallLubricationModelIter()->F<vector>());
539 const phasePair&
540 pair(this->phasePairs_[wallLubricationModelIter.key()]);
541
542 this->addField
543 (
544 pair.phase1(),
545 "phiF",
546 fvc::flux(rAUs[pair.phase1().index()]*F),
547 phiFs
548 );
549 this->addField
550 (
551 pair.phase2(),
552 "phiF",
553 - fvc::flux(rAUs[pair.phase2().index()]*F),
554 phiFs
555 );
556 }
557
558 // Add the phase pressure
559 DByAfs_.clear();
560 forAll(this->phaseModels_, phasei)
561 {
562 const phaseModel& phase = this->phaseModels_[phasei];
563
564 const surfaceScalarField pPrimeByAf
565 (
567 );
568
570 (
571 fvc::snGrad(phase)*this->mesh_.magSf()
572 );
573
574 this->addField(phase, "phiF", pPrimeByAf*snGradAlpha1, phiFs);
575
576 const bool implicitPhasePressure =
577 this->mesh_.solverDict(phase.volScalarField::name()).
578 template getOrDefault<Switch>
579 (
580 "implicitPhasePressure",
581 false
582 );
583
584 if (implicitPhasePressure)
585 {
586 this->addField(phase, "DByAf", pPrimeByAf, DByAfs_);
587 }
588 }
589
590 // Add the turbulent dispersion force
592 (
593 turbulentDispersionModelTable,
594 turbulentDispersionModels_,
595 turbulentDispersionModelIter
596 )
597 {
598 const phasePair&
599 pair(this->phasePairs_[turbulentDispersionModelIter.key()]);
600
601 const volScalarField D(turbulentDispersionModelIter()->D());
602
603 const surfaceScalarField DByA1f
604 (
605 fvc::interpolate(rAUs[pair.phase1().index()]*D)
606 );
607 const surfaceScalarField DByA2f
608 (
609 fvc::interpolate(rAUs[pair.phase2().index()]*D)
610 );
611
613 (
614 fvc::snGrad(pair.phase1())*this->mesh_.magSf()
615 );
616
617 this->addField(pair.phase1(), "phiF", DByA1f*snGradAlpha1, phiFs);
618 this->addField(pair.phase2(), "phiF", - DByA2f*snGradAlpha1, phiFs);
619
620 if (DByAfs_.found(pair.phase1().name()))
621 {
622 this->addField(pair.phase1(), "DByAf", DByA1f, DByAfs_);
623 }
624 }
625
626 if (this->fillFields_)
627 {
629 }
631 return phiFs;
632}
633
634
635template<class BasePhaseSystem>
638(
640)
641{
642 PtrList<surfaceScalarField> phiFfs(this->phaseModels_.size());
643
644 // Add the explicit part of the virtual mass force
645 forAllConstIter(VmfTable, Vmfs_, VmfIter)
646 {
647 const surfaceScalarField& Vmf(*VmfIter());
648 const phasePair& pair(this->phasePairs_[VmfIter.key()]);
649
650 forAllConstIter(phasePair, pair, iter)
651 {
652 this->addField
653 (
654 iter(),
655 "phiFf",
656 - rAUfs[iter().index()]*Vmf
657 *(
658 byDt(this->MRF().absolute(iter().phi()().oldTime()))
659 + iter.otherPhase().DUDtf()
660 ),
661 phiFfs
662 );
663 }
664 }
665
666 // Add the lift force
668 (
669 liftModelTable,
670 liftModels_,
671 liftModelIter
672 )
673 {
674 const surfaceScalarField Ff(liftModelIter()->Ff());
675 const phasePair& pair(this->phasePairs_[liftModelIter.key()]);
676
677 this->addField
678 (
679 pair.phase1(),
680 "phiFs",
681 rAUfs[pair.phase1().index()]*Ff,
682 phiFfs
683 );
684 this->addField
685 (
686 pair.phase2(),
687 "phiFf",
688 - rAUfs[pair.phase2().index()]*Ff,
689 phiFfs
690 );
691 }
692
693 // Add the wall lubrication force
695 (
696 wallLubricationModelTable,
697 wallLubricationModels_,
698 wallLubricationModelIter
699 )
700 {
701 const surfaceScalarField Ff(wallLubricationModelIter()->Ff());
702 const phasePair&
703 pair(this->phasePairs_[wallLubricationModelIter.key()]);
704
705 this->addField
706 (
707 pair.phase1(),
708 "phiFf",
709 rAUfs[pair.phase1().index()]*Ff,
710 phiFfs
711 );
712 this->addField
713 (
714 pair.phase2(),
715 "phiFf",
716 - rAUfs[pair.phase2().index()]*Ff,
717 phiFfs
718 );
719 }
720
721 // Add the phase pressure
722 DByAfs_.clear();
723 forAll(this->phaseModels_, phasei)
724 {
725 const phaseModel& phase = this->phaseModels_[phasei];
726
727 const surfaceScalarField pPrimeByAf
728 (
730 );
731
733 (
734 fvc::snGrad(phase)*this->mesh_.magSf()
735 );
736
737 this->addField(phase, "phiFf", pPrimeByAf*snGradAlpha1, phiFfs);
738
739 const bool implicitPhasePressure =
740 this->mesh_.solverDict(phase.volScalarField::name()).
741 template getOrDefault<Switch>
742 (
743 "implicitPhasePressure",
744 false
745 );
746
747 if (implicitPhasePressure)
748 {
749 this->addField(phase, "DByAf", pPrimeByAf, DByAfs_);
750 }
751 }
752
753 // Add the turbulent dispersion force and phase pressure
755 (
756 turbulentDispersionModelTable,
757 turbulentDispersionModels_,
758 turbulentDispersionModelIter
759 )
760 {
761 const phasePair&
762 pair(this->phasePairs_[turbulentDispersionModelIter.key()]);
763
764 const volScalarField D(turbulentDispersionModelIter()->D());
765
766 const surfaceScalarField DByAf1
767 (
768 rAUfs[pair.phase1().index()]*fvc::interpolate(D)
769 );
770 const surfaceScalarField DByAf2
771 (
772 rAUfs[pair.phase2().index()]*fvc::interpolate(D)
773 );
774
776 (
777 fvc::snGrad(pair.phase1())*this->mesh_.magSf()
778 );
779
780 this->addField(pair.phase1(), "phiFf", DByAf1*snGradAlpha1, phiFfs);
781 this->addField(pair.phase2(), "phiFf", - DByAf2*snGradAlpha1, phiFfs);
782
783 if (DByAfs_.found(pair.phase1().name()))
784 {
785 this->addField(pair.phase1(), "DByAf", DByAf1, DByAfs_);
786 }
787 }
788
789 if (this->fillFields_)
790 {
792 }
794 return phiFfs;
795}
796
797
798template<class BasePhaseSystem>
801(
803) const
804{
805 PtrList<surfaceScalarField> phiKdPhis(this->phaseModels_.size());
806
807 // Add the explicit part of the drag force
808 forAllConstIter(KdTable, Kds_, KdIter)
809 {
810 const volScalarField& K(*KdIter());
811 const phasePair& pair(this->phasePairs_[KdIter.key()]);
812
813 forAllConstIter(phasePair, pair, iter)
814 {
815 this->addField
816 (
817 iter(),
818 "phiKdPhi",
819 - fvc::interpolate(rAUs[iter().index()]*K)
820 *this->MRF().absolute(iter.otherPhase().phi()),
821 phiKdPhis
822 );
823 }
824 }
825
826 if (this->fillFields_)
827 {
828 this->fillFields
829 (
830 "phiKdPhi",
832 phiKdPhis
833 );
834 }
836 return phiKdPhis;
837}
838
839
840template<class BasePhaseSystem>
843(
845) const
846{
847 PtrList<surfaceScalarField> phiKdPhifs(this->phaseModels_.size());
848
849 // Add the explicit part of the drag force
850 forAllConstIter(KdfTable, Kdfs_, KdfIter)
851 {
852 const surfaceScalarField& Kf(*KdfIter());
853 const phasePair& pair(this->phasePairs_[KdfIter.key()]);
854
855 forAllConstIter(phasePair, pair, iter)
856 {
857 this->addField
858 (
859 iter(),
860 "phiKdPhif",
861 - rAUfs[iter().index()]*Kf
862 *this->MRF().absolute(iter.otherPhase().phi()),
863 phiKdPhifs
864 );
865 }
866 }
867
868 if (this->fillFields_)
869 {
870 this->fillFields
871 (
872 "phiKdPhif",
874 phiKdPhifs
875 );
876 }
878 return phiKdPhifs;
879}
880
881
882template<class BasePhaseSystem>
885(
887) const
888{
889 PtrList<volVectorField> KdUByAs(this->phaseModels_.size());
890
891 // Add the explicit part of the drag force
892 forAllConstIter(KdTable, Kds_, KdIter)
893 {
894 const volScalarField& K(*KdIter());
895 const phasePair& pair(this->phasePairs_[KdIter.key()]);
896
897 forAllConstIter(phasePair, pair, iter)
898 {
899 this->addField
900 (
901 iter(),
902 "KdUByA",
903 - rAUs[iter().index()]*K*iter.otherPhase().U(),
904 KdUByAs
905 );
906 }
907 }
908
909 if (this->fillFields_)
910 {
911 this->fillFields("KdUByA", dimVelocity, KdUByAs);
912 }
914 return KdUByAs;
915}
916
917
918template<class BasePhaseSystem>
921(
923 const bool includeVirtualMass
924) const
925{
926 PtrList<surfaceScalarField> ddtCorrByAs(this->phaseModels_.size());
927
928 // Construct phi differences
929 PtrList<surfaceScalarField> phiCorrs(this->phaseModels_.size());
930 forAll(this->phaseModels_, phasei)
931 {
932 const phaseModel& phase = this->phaseModels_[phasei];
933
934 phiCorrs.set
935 (
936 phasei,
937 this->MRF().absolute(phase.phi()().oldTime())
938 - fvc::flux(phase.U()().oldTime())
939 );
940 }
941
942 // Add correction
943 forAll(this->phaseModels_, phasei)
944 {
945 const phaseModel& phase = this->phaseModels_[phasei];
946 const volScalarField& alpha = phase;
947
948 // Apply ddtPhiCorr filter in pure(ish) phases
949 surfaceScalarField alphafBar
950 (
952 );
953
954 tmp<surfaceScalarField> phiCorrCoeff = pos0(alphafBar - 0.99);
955
956 surfaceScalarField::Boundary& phiCorrCoeffBf =
957 phiCorrCoeff.ref().boundaryFieldRef();
958
959 forAll(this->mesh_.boundary(), patchi)
960 {
961 // Set ddtPhiCorr to 0 on non-coupled boundaries
962 if
963 (
964 !this->mesh_.boundary()[patchi].coupled()
965 || isA<cyclicAMIFvPatch>(this->mesh_.boundary()[patchi])
966 )
967 {
968 phiCorrCoeffBf[patchi] = 0;
969 }
970 }
971
972 this->addField
973 (
974 phase,
975 "ddtCorrByA",
976 - phiCorrCoeff*phiCorrs[phasei]*fvc::interpolate
977 (
978 byDt(alpha.oldTime()*phase.rho()().oldTime()*rAUs[phasei])
979 ),
980 ddtCorrByAs
981 );
982 }
983
984 // Add virtual mass correction
985 if (includeVirtualMass)
986 {
987 forAllConstIter(VmTable, Vms_, VmIter)
988 {
989 const volScalarField& Vm(*VmIter());
990 const phasePair& pair(this->phasePairs_[VmIter.key()]);
991
992 forAllConstIter(phasePair, pair, iter)
993 {
994 const phaseModel& phase = iter();
995 const phaseModel& otherPhase = iter.otherPhase();
996
997 this->addField
998 (
999 iter(),
1000 "ddtCorrByA",
1001 - fvc::interpolate(Vm*byDt(rAUs[phase.index()]))
1002 *(
1003 phiCorrs[phase.index()]
1004 + this->MRF().absolute(otherPhase.phi())
1005 - fvc::flux(otherPhase.U())
1006 - phiCorrs[otherPhase.index()]
1007 ),
1008 ddtCorrByAs
1009 );
1010 }
1011 }
1013
1014 return ddtCorrByAs;
1015}
1016
1017
1018template<class BasePhaseSystem>
1020(
1021 const PtrList<volScalarField>& rAUs
1022)
1023{
1024 Info<< "Inverting drag systems: ";
1025
1026 phaseSystem::phaseModelList& phases = this->phaseModels_;
1027
1028 // Create drag coefficient matrices
1029 PtrList<PtrList<volScalarField>> KdByAs(phases.size());
1030 PtrList<PtrList<surfaceScalarField>> phiKds(phases.size());
1031
1033 {
1034 KdByAs.set
1035 (
1036 phasei,
1037 new PtrList<volScalarField>(phases.size())
1038 );
1039
1040 phiKds.set
1041 (
1042 phasei,
1043 new PtrList<surfaceScalarField>(phases.size())
1044 );
1045 }
1046
1047 forAllConstIter(KdTable, Kds_, KdIter)
1048 {
1049 const volScalarField& K(*KdIter());
1050 const phasePair& pair(this->phasePairs_[KdIter.key()]);
1051
1052 const label phase1i = pair.phase1().index();
1053 const label phase2i = pair.phase2().index();
1054
1055 this->addField
1056 (
1057 pair.phase2(),
1058 "KdByA",
1059 - rAUs[phase1i]*K,
1060 KdByAs[phase1i]
1061 );
1062 this->addField
1063 (
1064 pair.phase1(),
1065 "KdByA",
1066 - rAUs[phase2i]*K,
1067 KdByAs[phase2i]
1068 );
1069
1070 this->addField
1071 (
1072 pair.phase2(),
1073 "phiKd",
1074 fvc::interpolate(KdByAs[phase1i][phase2i]),
1075 phiKds[phase1i]
1076 );
1077 this->addField
1078 (
1079 pair.phase1(),
1080 "phiKd",
1081 fvc::interpolate(KdByAs[phase2i][phase1i]),
1082 phiKds[phase2i]
1083 );
1084 }
1085
1087 {
1088 this->fillFields("KdByAs", dimless, KdByAs[phasei]);
1089 this->fillFields("phiKds", dimless, phiKds[phasei]);
1090
1091 KdByAs[phasei][phasei] = 1;
1092 phiKds[phasei][phasei] = 1;
1093 }
1094
1095 // Decompose
1096 for (label i = 0; i < phases.size(); ++ i)
1097 {
1098 for (label j = i + 1; j < phases.size(); ++ j)
1099 {
1100 KdByAs[i][j] /= KdByAs[i][i];
1101 phiKds[i][j] /= phiKds[i][i];
1102 for (label k = i + 1; k < phases.size(); ++ k)
1103 {
1104 KdByAs[j][k] -= KdByAs[j][i]*KdByAs[i][k];
1105 phiKds[j][k] -= phiKds[j][i]*phiKds[i][k];
1106 }
1107 }
1108 }
1109 {
1110 volScalarField detKdByAs(KdByAs[0][0]);
1111 surfaceScalarField detPhiKdfs(phiKds[0][0]);
1112 for (label i = 1; i < phases.size(); ++ i)
1113 {
1114 detKdByAs *= KdByAs[i][i];
1115 detPhiKdfs *= phiKds[i][i];
1116 }
1117 Info<< "Min cell/face det = " << gMin(detKdByAs.primitiveField())
1118 << "/" << gMin(detPhiKdfs.primitiveField()) << endl;
1119 }
1120
1121 // Solve for the velocities and fluxes
1122 for (label i = 1; i < phases.size(); ++ i)
1123 {
1124 if (!phases[i].stationary())
1125 {
1126 for (label j = 0; j < i; j ++)
1127 {
1128 phases[i].URef() -= KdByAs[i][j]*phases[j].U();
1129 phases[i].phiRef() -= phiKds[i][j]*phases[j].phi();
1130 }
1131 }
1132 }
1133 for (label i = phases.size() - 1; i >= 0; i --)
1134 {
1135 if (!phases[i].stationary())
1136 {
1137 for (label j = phases.size() - 1; j > i; j --)
1138 {
1139 phases[i].URef() -= KdByAs[i][j]*phases[j].U();
1140 phases[i].phiRef() -= phiKds[i][j]*phases[j].phi();
1141 }
1142 phases[i].URef() /= KdByAs[i][i];
1143 phases[i].phiRef() /= phiKds[i][i];
1144 }
1145 }
1146}
1147
1148
1149template<class BasePhaseSystem>
1151(
1153)
1154{
1155 Info<< "Inverting drag system: ";
1156
1157 phaseSystem::phaseModelList& phases = this->phaseModels_;
1158
1159 // Create drag coefficient matrix
1161
1163 {
1164 phiKdfs.set
1165 (
1166 phasei,
1168 );
1169 }
1170
1171 forAllConstIter(KdfTable, Kdfs_, KdfIter)
1172 {
1173 const surfaceScalarField& K(*KdfIter());
1174 const phasePair& pair(this->phasePairs_[KdfIter.key()]);
1175
1176 const label phase1i = pair.phase1().index();
1177 const label phase2i = pair.phase2().index();
1178
1179 this->addField
1180 (
1181 pair.phase2(),
1182 "phiKdf",
1183 - rAUfs[phase1i]*K,
1184 phiKdfs[phase1i]
1185 );
1186 this->addField
1187 (
1188 pair.phase1(),
1189 "phiKdf",
1190 - rAUfs[phase2i]*K,
1191 phiKdfs[phase2i]
1192 );
1193 }
1194
1196 {
1197 this->fillFields("phiKdf", dimless, phiKdfs[phasei]);
1198
1199 phiKdfs[phasei][phasei] = 1;
1200 }
1201
1202 // Decompose
1203 for (label i = 0; i < phases.size(); ++ i)
1204 {
1205 for (label j = i + 1; j < phases.size(); ++ j)
1206 {
1207 phiKdfs[i][j] /= phiKdfs[i][i];
1208 for (label k = i + 1; k < phases.size(); ++ k)
1209 {
1210 phiKdfs[j][k] -= phiKdfs[j][i]*phiKdfs[i][k];
1211 }
1212 }
1213 }
1214 {
1215 surfaceScalarField detPhiKdfs(phiKdfs[0][0]);
1216 for (label i = 1; i < phases.size(); ++ i)
1217 {
1218 detPhiKdfs *= phiKdfs[i][i];
1219 }
1220 Info<< "Min face det = " << gMin(detPhiKdfs.primitiveField()) << endl;
1221 }
1222
1223 // Solve for the fluxes
1224 for (label i = 1; i < phases.size(); ++ i)
1225 {
1226 if (!phases[i].stationary())
1227 {
1228 for (label j = 0; j < i; j ++)
1229 {
1230 phases[i].phiRef() -= phiKdfs[i][j]*phases[j].phi();
1231 }
1232 }
1233 }
1234 for (label i = phases.size() - 1; i >= 0; i --)
1235 {
1236 if (!phases[i].stationary())
1237 {
1238 for (label j = phases.size() - 1; j > i; j --)
1239 {
1240 phases[i].phiRef() -= phiKdfs[i][j]*phases[j].phi();
1241 }
1242 phases[i].phiRef() /= phiKdfs[i][i];
1244 }
1245}
1246
1247
1248template<class BasePhaseSystem>
1249const Foam::HashPtrTable<Foam::surfaceScalarField>&
1251{
1252 return DByAfs_;
1253}
1254
1255
1256template<class BasePhaseSystem>
1258{
1259 if (BasePhaseSystem::read())
1260 {
1261 bool readOK = true;
1262
1263 // Read models ...
1264
1265 return readOK;
1266 }
1267 else
1268 {
1269 return false;
1270 }
1271}
1272
1273
1274// ************************************************************************* //
CGAL::Exact_predicates_exact_constructions_kernel K
label k
IOMRFZoneList & MRF
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
GeometricBoundaryField< scalar, fvsPatchField, surfaceMesh > Boundary
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
A HashTable of pointers to objects of type <T>, with deallocation management of the pointers.
bool set(const Key &key, T *ptr)
Assign a new entry, overwrites existing.
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Class which models interfacial momenum transfer between a number of phases. Drag, virtual mass,...
virtual PtrList< surfaceScalarField > AFfs() const
Return implicit force coefficients on the faces, for the face-based.
HashTable< autoPtr< BlendedInterfacialModel< dragModel > >, phasePairKey, phasePairKey::hash > dragModelTable
HashPtrTable< volScalarField, phasePairKey, phasePairKey::hash > KdTable
virtual PtrList< surfaceScalarField > phiFs(const PtrList< volScalarField > &rAUs)
Return the explicit force fluxes for the cell-based algorithm, that.
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransfer()
Return the momentum transfer matrices for the cell-based algorithm.
virtual PtrList< surfaceScalarField > ddtCorrByAs(const PtrList< volScalarField > &rAUs, const bool includeVirtualMass=false) const
Return the flux corrections for the cell-based algorithm. These.
virtual void partialElimination(const PtrList< volScalarField > &rAUs)
Solve the drag system for the velocities and fluxes.
virtual PtrList< surfaceScalarField > phiKdPhifs(const PtrList< surfaceScalarField > &rAUfs) const
As phiKdPhis, but for the face-based algorithm.
HashTable< autoPtr< BlendedInterfacialModel< turbulentDispersionModel > >, phasePairKey, phasePairKey::hash > turbulentDispersionModelTable
HashPtrTable< surfaceScalarField, phasePairKey, phasePairKey::hash > VmfTable
HashPtrTable< surfaceScalarField, phasePairKey, phasePairKey::hash > KdfTable
HashTable< autoPtr< BlendedInterfacialModel< virtualMassModel > >, phasePairKey, phasePairKey::hash > virtualMassModelTable
HashPtrTable< volScalarField, phasePairKey, phasePairKey::hash > VmTable
virtual PtrList< volVectorField > KdUByAs(const PtrList< volScalarField > &rAUs) const
Return the explicit part of the drag force for the cell-based.
virtual const HashPtrTable< surfaceScalarField > & DByAfs() const
Return the phase diffusivities divided by the momentum coefficients.
virtual void partialEliminationf(const PtrList< surfaceScalarField > &rAUfs)
As partialElimination, but for the face-based algorithm. Only solves.
virtual PtrList< surfaceScalarField > phiFfs(const PtrList< surfaceScalarField > &rAUfs)
As phiFs, but for the face-based algorithm.
MomentumTransferPhaseSystem(const fvMesh &)
Construct from fvMesh.
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransferf()
As momentumTransfer, but for the face-based algorithm.
virtual PtrList< surfaceScalarField > phiKdPhis(const PtrList< volScalarField > &rAUs) const
Return the explicit drag force fluxes for the cell-based algorithm.
HashTable< autoPtr< BlendedInterfacialModel< wallLubricationModel > >, phasePairKey, phasePairKey::hash > wallLubricationModelTable
HashTable< autoPtr< BlendedInterfacialModel< liftModel > >, phasePairKey, phasePairKey::hash > liftModelTable
virtual bool read()
Read base phaseProperties dictionary.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
Definition PtrList.H:171
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition fvMatrix.H:487
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const word & name() const
The name of this phase.
Definition phaseModel.H:122
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phaseModel.H:56
const volVectorField & U() const
Definition phaseModel.H:198
const surfaceScalarField & phi() const
Definition phaseModel.H:218
const phaseModel & otherPhase() const
Return the other phase in this two-phase system.
Definition phaseModel.C:211
virtual tmp< volVectorField > DUDt() const =0
Return the substantive acceleration.
label index() const
Return the index of the phase.
Definition phaseModel.C:135
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition phasePair.H:52
virtual word name() const
Pair name.
Definition phasePair.C:62
const multiphaseInter::phaseModel & phase1() const
Definition phasePairI.H:23
const multiphaseInter::phaseModel & phase2() const
Definition phasePairI.H:29
PtrListDictionary< phaseModel > phaseModelList
Definition phaseSystem.H:83
HashPtrTable< fvVectorMatrix > momentumTransferTable
Definition phaseSystem.H:77
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phase.H:53
const word & name() const
Definition phase.H:114
const dimensionedScalar & rho() const
Return const-access to phase1 density.
Definition phase.H:151
A class for managing temporary objects.
Definition tmp.H:75
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
U
Definition pEqn.H:72
dynamicFvMesh & mesh
Info<< "Creating field kinetic energy K\n"<< endl;volScalarField K("K", 0.5 *magSqr(U));if(U.nOldTimes()){ volVectorField *Uold=&U.oldTime();volScalarField *Kold=&K.oldTime(); *Kold==0.5 *magSqr(*Uold);while(Uold->nOldTimes()) { Uold=&Uold-> oldTime()
Area-weighted average a surfaceField creating a volField.
Calculate the first temporal derivative.
Calculate the divergence of the given field.
Calculate the face-flux of the given field.
Calculate the snGrad of the given volField.
Calculate the matrix for the first temporal derivative.
Calculate the matrix for the divergence of the given field and flux.
Calculate the finiteVolume matrix for implicit and explicit sources.
label phasei
Definition pEqn.H:27
PtrList< surfaceScalarField > phiFs(fluid.phiFs(rAUs))
fluid fillFields("rAU", dimTime/dimDensity, rAUs)
PtrList< surfaceScalarField > phiFfs(fluid.phiFfs(rAUfs))
PtrList< volScalarField > rAUs
Definition pEqn.H:4
PtrList< surfaceScalarField > rAUfs
Definition pEqn.H:30
PtrList< surfaceScalarField > AFfs(fluid.AFfs())
volVectorField F(fluid.F())
surfaceScalarField snGradAlpha1(fvc::snGrad(alpha1) *mesh.magSf())
surfaceScalarField Vmf("Vmf", fluid.Vmf())
surfaceScalarField Ff(fluid.Ff())
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< Type, fvPatchField, volMesh > > average(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Area-weighted average a surfaceField creating a volField.
Definition fvcAverage.C:41
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< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition fvmDiv.C:41
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition fvmDdt.C:41
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
constexpr auto key(const Type &t) noexcept
Helper function to return the enum value.
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
dimensionedScalar pos0(const dimensionedScalar &ds)
tmp< volScalarField > byDt(const volScalarField &vf)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimVelocity
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensionedScalar negPart(const dimensionedScalar &ds)
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
fvMatrix< vector > fvVectorMatrix
const dimensionSet dimForce
Type gMin(const FieldField< Field, Type > &f)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimDensity
dimensionedScalar posPart(const dimensionedScalar &ds)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
volScalarField & alpha
const dimensionedScalar & D
multiphaseSystem::phaseModelList & phases
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object.
Definition stdFoam.H:356