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