MomentumTransferPhaseSystem.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2015-2021 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 
29 #include "dragModel.H"
30 #include "virtualMassModel.H"
31 #include "liftModel.H"
32 #include "wallLubricationModel.H"
34 
35 #include "HashPtrTable.H"
36 
37 #include "fvmDdt.H"
38 #include "fvmDiv.H"
39 #include "fvmSup.H"
40 #include "fvcAverage.H"
41 #include "fvcDdt.H"
42 #include "fvcDiv.H"
43 #include "fvcFlux.H"
44 #include "fvcSnGrad.H"
45 #include "fvcMeshPhi.H"
46 #include "fvMatrix.H"
47 
48 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
49 
50 template<class BasePhaseSystem>
53 (
54  const phasePairKey& key
55 ) const
56 {
57  if (dragModels_.found(key))
58  {
59  return dragModels_[key]->K();
60  }
61  else
62  {
63  return volScalarField::New
64  (
65  dragModel::typeName + ":K",
66  this->mesh_,
67  dimensionedScalar(dragModel::dimK, 0)
68  );
69  }
70 }
71 
72 
73 template<class BasePhaseSystem>
76 (
77  const phasePairKey& key
78 ) const
79 {
80  if (dragModels_.found(key))
81  {
82  return dragModels_[key]->Kf();
83  }
84  else
85  {
87  (
88  dragModel::typeName + ":K",
89  this->mesh_,
90  dimensionedScalar(dragModel::dimK, 0)
91  );
92  }
93 }
94 
95 
96 template<class BasePhaseSystem>
99 (
100  const phasePairKey& key
101 ) const
102 {
103  if (virtualMassModels_.found(key))
104  {
105  return virtualMassModels_[key]->K();
106  }
107  else
108  {
109  return volScalarField::New
110  (
111  virtualMassModel::typeName + ":K",
112  this->mesh_,
113  dimensionedScalar(virtualMassModel::dimK, 0)
114  );
115  }
116 }
117 
118 
119 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
120 
121 template<class BasePhaseSystem>
123 (
124  const phaseSystem::dmdtfTable& dmdtfs,
125  phaseSystem::momentumTransferTable& eqns
126 )
127 {
128  forAllConstIter(phaseSystem::dmdtfTable, dmdtfs, dmdtfIter)
129  {
130  const phasePairKey& key = dmdtfIter.key();
131  const phasePair& pair(this->phasePairs_[key]);
132 
133  const volScalarField dmdtf(Pair<word>::compare(pair, key)**dmdtfIter());
134 
135  phaseModel& phase1 = this->phases()[pair.phase1().name()];
136  phaseModel& phase2 = this->phases()[pair.phase2().name()];
137 
138  const volScalarField dmdtf21(posPart(dmdtf));
139  const volScalarField dmdtf12(negPart(dmdtf));
140 
141  if (!phase1.stationary())
142  {
143  *eqns[phase1.name()] +=
144  dmdtf21*phase2.U() + fvm::Sp(dmdtf12, phase1.URef());
145  }
146 
147  if (!phase2.stationary())
148  {
149  *eqns[phase2.name()] -=
150  dmdtf12*phase1.U() + fvm::Sp(dmdtf21, phase2.URef());
151  }
152  }
153 }
154 
155 
156 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
157 
158 template<class BasePhaseSystem>
161 (
162  const fvMesh& mesh
163 )
164 :
165  BasePhaseSystem(mesh)
166 {
167  this->generatePairsAndSubModels
168  (
169  "drag",
170  dragModels_
171  );
172 
173  this->generatePairsAndSubModels
174  (
175  "virtualMass",
176  virtualMassModels_
177  );
178 
179  this->generatePairsAndSubModels
180  (
181  "lift",
182  liftModels_
183  );
184 
185  this->generatePairsAndSubModels
186  (
187  "wallLubrication",
188  wallLubricationModels_
189  );
190 
191  this->generatePairsAndSubModels
192  (
193  "turbulentDispersion",
194  turbulentDispersionModels_
195  );
196 
198  (
199  dragModelTable,
200  dragModels_,
201  dragModelIter
202  )
203  {
204  const phasePair& pair(this->phasePairs_[dragModelIter.key()]);
205 
206  Kds_.insert
207  (
208  pair,
209  new volScalarField
210  (
211  IOobject
212  (
213  IOobject::groupName("Kd", pair.name()),
214  this->mesh().time().timeName(),
215  this->mesh()
216  ),
217  this->mesh(),
219  )
220  );
221 
222  Kdfs_.insert
223  (
224  pair,
226  (
227  IOobject
228  (
229  IOobject::groupName("Kdf", pair.name()),
230  this->mesh().time().timeName(),
231  this->mesh()
232  ),
233  this->mesh(),
235  )
236  );
237  }
238 
240  (
241  virtualMassModelTable,
242  virtualMassModels_,
243  virtualMassModelIter
244  )
245  {
246  const phasePair& pair(this->phasePairs_[virtualMassModelIter.key()]);
247 
248  Vms_.insert
249  (
250  pair,
251  new volScalarField
252  (
253  IOobject
254  (
255  IOobject::groupName("Vm", pair.name()),
256  this->mesh().time().timeName(),
257  this->mesh()
258  ),
259  this->mesh(),
261  )
262  );
263 
264  Vmfs_.insert
265  (
266  pair,
268  (
269  IOobject
270  (
271  IOobject::groupName("Vmf", pair.name()),
272  this->mesh().time().timeName(),
273  this->mesh()
274  ),
275  this->mesh(),
277  )
278  );
279  }
280 }
281 
282 
283 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
284 
285 template<class BasePhaseSystem>
288 {}
289 
290 
291 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
292 
293 template<class BasePhaseSystem>
296 {
297  // Create a momentum transfer matrix for each phase
298  autoPtr<phaseSystem::momentumTransferTable> eqnsPtr
299  (
301  );
302 
303  phaseSystem::momentumTransferTable& eqns = eqnsPtr();
304 
305  forAll(this->movingPhases(), movingPhasei)
306  {
307  const phaseModel& phase = this->movingPhases()[movingPhasei];
308 
309  eqns.insert
310  (
311  phase.name(),
312  new fvVectorMatrix(phase.U(), dimMass*dimVelocity/dimTime)
313  );
314  }
315 
316  // Update the drag coefficients
318  (
319  dragModelTable,
320  dragModels_,
321  dragModelIter
322  )
323  {
324  *Kds_[dragModelIter.key()] = dragModelIter()->K();
325  *Kdfs_[dragModelIter.key()] = dragModelIter()->Kf();
326  }
327 
328  // Add the implicit part of the drag force
329  forAllConstIter(KdTable, Kds_, KdIter)
330  {
331  const volScalarField& K(*KdIter());
332  const phasePair& pair(this->phasePairs_[KdIter.key()]);
333 
334  forAllConstIter(phasePair, pair, iter)
335  {
336  if (!iter().stationary())
337  {
338  fvVectorMatrix& eqn = *eqns[iter().name()];
339 
340  eqn -= fvm::Sp(K, eqn.psi());
341  }
342  }
343  }
344 
345  // Update the virtual mass coefficients
347  (
348  virtualMassModelTable,
349  virtualMassModels_,
350  virtualMassModelIter
351  )
352  {
353  *Vms_[virtualMassModelIter.key()] = virtualMassModelIter()->K();
354  *Vmfs_[virtualMassModelIter.key()] = virtualMassModelIter()->Kf();
355  }
356 
357  // Add the virtual mass force
358  forAllConstIter(VmTable, Vms_, VmIter)
359  {
360  const volScalarField& Vm(*VmIter());
361  const phasePair& pair(this->phasePairs_[VmIter.key()]);
362 
363  forAllConstIter(phasePair, pair, iter)
364  {
365  const phaseModel& phase = iter();
366  const phaseModel& otherPhase = iter.otherPhase();
367 
368  if (!phase.stationary())
369  {
370  fvVectorMatrix& eqn = *eqns[phase.name()];
371 
372  const volVectorField& U = eqn.psi();
373  const surfaceScalarField& phi = phase.phi();
374  const tmp<surfaceScalarField> taphi(fvc::absolute(phi, U));
375  const surfaceScalarField& aphi(taphi());
376 
377  eqn -=
378  Vm
379  *(
380  fvm::ddt(U)
381  + fvm::div(aphi, U) - fvm::Sp(fvc::div(aphi), U)
382  - otherPhase.DUDt()
383  )
384  + this->MRF_.DDt(Vm, U - otherPhase.U());
385  }
386  }
387  }
388 
389  return eqnsPtr;
390 }
391 
392 
393 template<class BasePhaseSystem>
396 {
397  // Create a momentum transfer matrix for each phase
398  autoPtr<phaseSystem::momentumTransferTable> eqnsPtr
399  (
401  );
402 
403  phaseSystem::momentumTransferTable& eqns = eqnsPtr();
404 
405  forAll(this->movingPhases(), movingPhasei)
406  {
407  const phaseModel& phase = this->movingPhases()[movingPhasei];
408 
409  eqns.insert
410  (
411  phase.name(),
412  new fvVectorMatrix(phase.U(), dimMass*dimVelocity/dimTime)
413  );
414  }
415 
416  // Create U & grad(U) fields
417  PtrList<fvVectorMatrix> UgradUs(this->phaseModels_.size());
418  forAll(this->phaseModels_, phasei)
419  {
420  const phaseModel& phase = this->phaseModels_[phasei];
421 
422  if (!phase.stationary())
423  {
424  const volVectorField& U = phase.U();
425  const surfaceScalarField& phi = phase.phi();
426  const tmp<surfaceScalarField> taphi(fvc::absolute(phi, U));
427  const surfaceScalarField& aphi(taphi());
428 
429  UgradUs.set
430  (
431  phasei,
432  new fvVectorMatrix
433  (
434  fvm::div(aphi, U) - fvm::Sp(fvc::div(aphi), U)
435  + this->MRF().DDt(U)
436  )
437  );
438  }
439  }
440 
441  // Add the virtual mass force
442  forAllConstIter(VmTable, Vms_, VmIter)
443  {
444  const volScalarField& Vm(*VmIter());
445  const phasePair& pair(this->phasePairs_[VmIter.key()]);
446 
447  forAllConstIter(phasePair, pair, iter)
448  {
449  const phaseModel& phase = iter();
450  const phaseModel& otherPhase = iter.otherPhase();
451 
452  if (!phase.stationary())
453  {
454  *eqns[phase.name()] -=
455  Vm
456  *(
457  UgradUs[phase.index()]
458  - (UgradUs[otherPhase.index()] & otherPhase.U())
459  );
460  }
461  }
462  }
463 
464  return eqnsPtr;
465 }
466 
467 
468 template<class BasePhaseSystem>
471 {
472  PtrList<surfaceScalarField> AFfs(this->phaseModels_.size());
473 
474  // Add the implicit part of the drag force
475  forAllConstIter(KdfTable, Kdfs_, KdfIter)
476  {
477  const surfaceScalarField& Kf(*KdfIter());
478  const phasePair& pair(this->phasePairs_[KdfIter.key()]);
479 
480  forAllConstIter(phasePair, pair, iter)
481  {
482  addField(iter(), "AFf", Kf, AFfs);
483  }
484  }
485 
486  // Add the implicit part of the virtual mass force
487  forAllConstIter(VmfTable, Vmfs_, VmfIter)
488  {
489  const surfaceScalarField& Vmf(*VmfIter());
490  const phasePair& pair(this->phasePairs_[VmfIter.key()]);
491 
492  forAllConstIter(phasePair, pair, iter)
493  {
494  addField(iter(), "AFf", byDt(Vmf), AFfs);
495  }
496  }
497 
498  this->fillFields("AFf", dimDensity/dimTime, AFfs);
499 
500  return AFfs;
501 }
502 
503 
504 template<class BasePhaseSystem>
507 (
508  const PtrList<volScalarField>& rAUs
509 )
510 {
511  PtrList<surfaceScalarField> phiFs(this->phaseModels_.size());
512 
513  // Add the lift force
515  (
516  liftModelTable,
517  liftModels_,
518  liftModelIter
519  )
520  {
521  const volVectorField F(liftModelIter()->F<vector>());
522  const phasePair& pair(this->phasePairs_[liftModelIter.key()]);
523 
524  addField
525  (
526  pair.phase1(),
527  "phiF",
528  fvc::flux(rAUs[pair.phase1().index()]*F),
529  phiFs
530  );
531  addField
532  (
533  pair.phase2(),
534  "phiF",
535  -fvc::flux(rAUs[pair.phase2().index()]*F),
536  phiFs
537  );
538  }
539 
540  // Add the wall lubrication force
542  (
543  wallLubricationModelTable,
544  wallLubricationModels_,
545  wallLubricationModelIter
546  )
547  {
548  const volVectorField F(wallLubricationModelIter()->F<vector>());
549  const phasePair&
550  pair(this->phasePairs_[wallLubricationModelIter.key()]);
551 
552  addField
553  (
554  pair.phase1(),
555  "phiF",
556  fvc::flux(rAUs[pair.phase1().index()]*F),
557  phiFs
558  );
559  addField
560  (
561  pair.phase2(),
562  "phiF",
563  -fvc::flux(rAUs[pair.phase2().index()]*F),
564  phiFs
565  );
566  }
567 
568  // Add the phase pressure
569  forAll(this->phaseModels_, phasei)
570  {
571  const phaseModel& phase = this->phaseModels_[phasei];
572 
573  const surfaceScalarField pPrimeByAf
574  (
575  fvc::interpolate(rAUs[phasei]*phase.pPrime())
576  );
577 
578  const surfaceScalarField snGradAlpha1
579  (
580  fvc::snGrad(phase)*this->mesh_.magSf()
581  );
582 
583  addField(phase, "phiF", pPrimeByAf*snGradAlpha1, phiFs);
584  }
585 
586  // Add the turbulent dispersion force
588  (
589  turbulentDispersionModelTable,
590  turbulentDispersionModels_,
591  turbulentDispersionModelIter
592  )
593  {
594  const phasePair&
595  pair(this->phasePairs_[turbulentDispersionModelIter.key()]);
596 
597  const volScalarField D(turbulentDispersionModelIter()->D());
598 
599  const surfaceScalarField DByA1f
600  (
601  fvc::interpolate(rAUs[pair.phase1().index()]*D)
602  );
603  const surfaceScalarField DByA2f
604  (
605  fvc::interpolate(rAUs[pair.phase2().index()]*D)
606  );
607 
608  const volScalarField alpha12(pair.phase1() + pair.phase2());
609  const surfaceScalarField snGradAlpha1By12
610  (
612  (
613  pair.phase1()/max(alpha12, pair.phase1().residualAlpha())
614  )*this->mesh_.magSf()
615  );
616  const surfaceScalarField snGradAlpha2By12
617  (
619  (
620  pair.phase2()/max(alpha12, pair.phase2().residualAlpha())
621  )*this->mesh_.magSf()
622  );
623 
624  addField(pair.phase1(), "phiF", DByA1f*snGradAlpha1By12, phiFs);
625  addField(pair.phase2(), "phiF", DByA2f*snGradAlpha2By12, phiFs);
626  }
627 
628  if (this->fillFields_)
629  {
631  }
632 
633  return phiFs;
634 }
635 
636 
637 template<class BasePhaseSystem>
640 (
641  const PtrList<surfaceScalarField>& rAUfs
642 )
643 {
644  PtrList<surfaceScalarField> phiFfs(this->phaseModels_.size());
645 
646  // Add the explicit part of the virtual mass force
647  forAllConstIter(VmfTable, Vmfs_, VmfIter)
648  {
649  const surfaceScalarField& Vmf(*VmfIter());
650  const phasePair& pair(this->phasePairs_[VmfIter.key()]);
651 
652  forAllConstIter(phasePair, pair, iter)
653  {
654  addField
655  (
656  iter(),
657  "phiFf",
658  -rAUfs[iter().index()]*Vmf
659  *(
660  byDt
661  (
663  (
664  this->MRF().absolute(iter().phi()().oldTime()),
665  iter().U()
666  )
667  )
668  + iter.otherPhase().DUDtf()
669  ),
670  phiFfs
671  );
672  }
673  }
674 
675  // Add the lift force
677  (
678  liftModelTable,
679  liftModels_,
680  liftModelIter
681  )
682  {
683  const surfaceScalarField Ff(liftModelIter()->Ff());
684  const phasePair& pair(this->phasePairs_[liftModelIter.key()]);
685 
686  addField
687  (
688  pair.phase1(),
689  "phiFf",
690  rAUfs[pair.phase1().index()]*Ff,
691  phiFfs
692  );
693  addField
694  (
695  pair.phase2(),
696  "phiFf",
697  -rAUfs[pair.phase2().index()]*Ff,
698  phiFfs
699  );
700  }
701 
702  // Add the wall lubrication force
704  (
705  wallLubricationModelTable,
706  wallLubricationModels_,
707  wallLubricationModelIter
708  )
709  {
710  const surfaceScalarField Ff(wallLubricationModelIter()->Ff());
711  const phasePair&
712  pair(this->phasePairs_[wallLubricationModelIter.key()]);
713 
714  addField
715  (
716  pair.phase1(),
717  "phiFf",
718  rAUfs[pair.phase1().index()]*Ff,
719  phiFfs
720  );
721  addField
722  (
723  pair.phase2(),
724  "phiFf",
725  -rAUfs[pair.phase2().index()]*Ff,
726  phiFfs
727  );
728  }
729 
730  // Add the phase pressure
731  forAll(this->phaseModels_, phasei)
732  {
733  const phaseModel& phase = this->phaseModels_[phasei];
734 
735  const surfaceScalarField pPrimeByAf
736  (
737  rAUfs[phasei]*fvc::interpolate(phase.pPrime())
738  );
739 
740  const surfaceScalarField snGradAlpha1
741  (
742  fvc::snGrad(phase)*this->mesh_.magSf()
743  );
744 
745  addField(phase, "phiFf", pPrimeByAf*snGradAlpha1, phiFfs);
746  }
747 
748  // Add the turbulent dispersion force
750  (
751  turbulentDispersionModelTable,
752  turbulentDispersionModels_,
753  turbulentDispersionModelIter
754  )
755  {
756  const surfaceScalarField Ff(turbulentDispersionModelIter()->Ff());
757  const phasePair&
758  pair(this->phasePairs_[turbulentDispersionModelIter.key()]);
759 
760  addField
761  (
762  pair.phase1(),
763  "phiFf",
764  rAUfs[pair.phase1().index()]*Ff,
765  phiFfs
766  );
767  addField
768  (
769  pair.phase2(),
770  "phiFf",
771  -rAUfs[pair.phase2().index()]*Ff,
772  phiFfs
773  );
774  }
775 
776  if (this->fillFields_)
777  {
779  }
780 
781  return phiFfs;
782 }
783 
784 
785 template<class BasePhaseSystem>
788 (
789  const PtrList<volScalarField>& rAUs
790 ) const
791 {
792  PtrList<surfaceScalarField> phiKdPhis(this->phaseModels_.size());
793 
794  // Add the explicit part of the drag force
795  forAllConstIter(KdTable, Kds_, KdIter)
796  {
797  const volScalarField& K(*KdIter());
798  const phasePair& pair(this->phasePairs_[KdIter.key()]);
799 
800  forAllConstIter(phasePair, pair, iter)
801  {
802  addField
803  (
804  iter(),
805  "phiKdPhi",
806  -fvc::interpolate(rAUs[iter().index()]*K)
808  (
809  this->MRF().absolute(iter.otherPhase().phi()),
810  iter.otherPhase().U()
811  ),
812  phiKdPhis
813  );
814  }
815  }
816 
817  if (this->fillFields_)
818  {
819  this->fillFields
820  (
821  "phiKdPhi",
823  phiKdPhis
824  );
825  }
826 
827  return phiKdPhis;
828 }
829 
830 
831 template<class BasePhaseSystem>
834 (
835  const PtrList<surfaceScalarField>& rAUfs
836 ) const
837 {
838  PtrList<surfaceScalarField> phiKdPhifs(this->phaseModels_.size());
839 
840  // Add the explicit part of the drag force
841  forAllConstIter(KdfTable, Kdfs_, KdfIter)
842  {
843  const surfaceScalarField& Kf(*KdfIter());
844  const phasePair& pair(this->phasePairs_[KdfIter.key()]);
845 
846  forAllConstIter(phasePair, pair, iter)
847  {
848  addField
849  (
850  iter(),
851  "phiKdPhif",
852  -rAUfs[iter().index()]*Kf
854  (
855  this->MRF().absolute(iter.otherPhase().phi()),
856  iter.otherPhase().U()
857  ),
858  phiKdPhifs
859  );
860  }
861  }
862 
863  if (this->fillFields_)
864  {
865  this->fillFields
866  (
867  "phiKdPhif",
869  phiKdPhifs
870  );
871  }
872 
873  return phiKdPhifs;
874 }
875 
876 
877 template<class BasePhaseSystem>
880 (
881  const PtrList<volScalarField>& rAUs
882 ) const
883 {
884  PtrList<volVectorField> KdUByAs(this->phaseModels_.size());
885 
886  // Add the explicit part of the drag force
887  forAllConstIter(KdTable, Kds_, KdIter)
888  {
889  const volScalarField& K(*KdIter());
890  const phasePair& pair(this->phasePairs_[KdIter.key()]);
891 
892  forAllConstIter(phasePair, pair, iter)
893  {
894  addField
895  (
896  iter(),
897  "KdUByA",
898  -rAUs[iter().index()]*K*iter.otherPhase().U(),
899  KdUByAs
900  );
901  }
902  }
903 
904  if (this->fillFields_)
905  {
906  this->fillFields("KdUByA", dimVelocity, KdUByAs);
907  }
908 
909  return KdUByAs;
910 }
911 
912 
913 template<class BasePhaseSystem>
915 implicitPhasePressure(const phaseModel& phase) const
916 {
917  return
918  this->mesh_.solverDict(phase.volScalarField::name()).
919  template lookupOrDefault<Switch>
920  (
921  "implicitPhasePressure",
922  false
923  );
924 }
925 
926 
927 template<class BasePhaseSystem>
929 implicitPhasePressure() const
930 {
931  bool implicitPressure = false;
932 
933  forAll(this->phaseModels_, phasei)
934  {
935  const phaseModel& phase = this->phaseModels_[phasei];
936 
937  implicitPressure = implicitPressure || implicitPhasePressure(phase);
938  }
939 
940  return implicitPressure;
941 }
942 
943 
944 template<class BasePhaseSystem>
947 (
948  const PtrList<volScalarField>& rAUs,
949  const PtrList<surfaceScalarField>& rAUfs
950 ) const
951 {
952  PtrList<surfaceScalarField> DByAfs(this->phaseModels_.size());
953 
954  if (rAUfs.size())
955  {
956  // Add the phase pressure
957  forAll(this->phaseModels_, phasei)
958  {
959  const phaseModel& phase = this->phaseModels_[phasei];
960 
961  const surfaceScalarField pPrimeByAf
962  (
963  rAUfs[phasei]*fvc::interpolate(phase.pPrime())
964  );
965 
966  addField(phase, "DByAf", pPrimeByAf, DByAfs);
967 
968  forAll(this->phaseModels_, phasej)
969  {
970  if (phasej != phasei)
971  {
972  const phaseModel& phase2 = this->phaseModels_[phasej];
973 
974  addField
975  (
976  phase2,
977  "DByAf",
979  (
980  phase2
981  /max(1 - phase, phase2.residualAlpha())
982  )*pPrimeByAf,
983  DByAfs
984  );
985  }
986  }
987  }
988 
989  // Add the turbulent dispersion
991  (
992  turbulentDispersionModelTable,
993  turbulentDispersionModels_,
994  turbulentDispersionModelIter
995  )
996  {
997  const phasePair&
998  pair(this->phasePairs_[turbulentDispersionModelIter.key()]);
999 
1000  const surfaceScalarField Df
1001  (
1002  fvc::interpolate(turbulentDispersionModelIter()->D())
1003  );
1004 
1005  const surfaceScalarField alpha12f
1006  (
1007  fvc::interpolate(pair.phase1() + pair.phase2())
1008  );
1009 
1010  addField
1011  (
1012  pair.phase1(),
1013  "DByAf",
1014  rAUfs[pair.phase1().index()]*Df
1015  /max(alpha12f, pair.phase1().residualAlpha()),
1016  DByAfs
1017  );
1018  addField
1019  (
1020  pair.phase2(),
1021  "DByAf",
1022  rAUfs[pair.phase2().index()]*Df
1023  /max(alpha12f, pair.phase2().residualAlpha()),
1024  DByAfs
1025  );
1026  }
1027  }
1028  else
1029  {
1030  // Add the phase pressure
1031  forAll(this->phaseModels_, phasei)
1032  {
1033  const phaseModel& phase = this->phaseModels_[phasei];
1034 
1035  const surfaceScalarField pPrimeByAf
1036  (
1037  fvc::interpolate(rAUs[phasei]*phase.pPrime())
1038  );
1039 
1040  addField(phase, "DByAf", pPrimeByAf, DByAfs);
1041 
1042  forAll(this->phaseModels_, phasej)
1043  {
1044  if (phasej != phasei)
1045  {
1046  const phaseModel& phase2 = this->phaseModels_[phasej];
1047 
1048  addField
1049  (
1050  phase2,
1051  "DByAf",
1053  (
1054  phase2
1055  /max(1 - phase, phase2.residualAlpha())
1056  )*pPrimeByAf,
1057  DByAfs
1058  );
1059  }
1060  }
1061  }
1062 
1063  // Add the turbulent dispersion
1065  (
1066  turbulentDispersionModelTable,
1067  turbulentDispersionModels_,
1068  turbulentDispersionModelIter
1069  )
1070  {
1071  const phasePair&
1072  pair(this->phasePairs_[turbulentDispersionModelIter.key()]);
1073 
1074  const volScalarField D(turbulentDispersionModelIter()->D());
1075 
1076  const volScalarField alpha12(pair.phase1() + pair.phase2());
1077 
1078  addField
1079  (
1080  pair.phase1(),
1081  "DByAf",
1083  (
1084  rAUs[pair.phase1().index()]*D
1085  /max(alpha12, pair.phase1().residualAlpha())
1086  ),
1087  DByAfs
1088  );
1089  addField
1090  (
1091  pair.phase2(),
1092  "DByAf",
1094  (
1095  rAUs[pair.phase2().index()]*D
1096  /max(alpha12, pair.phase2().residualAlpha())
1097  ),
1098  DByAfs
1099  );
1100  }
1101  }
1102 
1103  return DByAfs;
1104 }
1105 
1106 
1107 template<class BasePhaseSystem>
1110 (
1111  const PtrList<volScalarField>& rAUs,
1112  const bool includeVirtualMass
1113 ) const
1114 {
1115  PtrList<surfaceScalarField> ddtCorrByAs(this->phaseModels_.size());
1116 
1117  // Construct phi differences
1118  PtrList<surfaceScalarField> phiCorrs(this->phaseModels_.size());
1119 
1120  forAll(this->movingPhases(), movingPhasei)
1121  {
1122  const phaseModel& phase = this->movingPhases()[movingPhasei];
1123  const label phasei = phase.index();
1124 
1125  phiCorrs.set
1126  (
1127  phasei,
1128  this->MRF().zeroFilter
1129  (
1130  (
1131  phase.Uf().valid()
1132  ? (this->mesh_.Sf() & phase.Uf()().oldTime())()
1133  : phase.phi()().oldTime()
1134  )
1135  - fvc::flux(phase.U()().oldTime())
1136  )()
1137  );
1138  }
1139 
1140  // Add correction
1141  forAll(this->movingPhases(), movingPhasei)
1142  {
1143  const phaseModel& phase = this->movingPhases()[movingPhasei];
1144  const label phasei = phase.index();
1145  const volScalarField& alpha = phase;
1146 
1147  // Apply ddtPhiCorr filter in pure(ish) phases
1148  surfaceScalarField alphafBar
1149  (
1151  );
1152 
1153  tmp<surfaceScalarField> phiCorrCoeff = pos0(alphafBar - 0.99);
1154 
1155  surfaceScalarField::Boundary& phiCorrCoeffBf =
1156  phiCorrCoeff.ref().boundaryFieldRef();
1157 
1158  forAll(this->mesh_.boundary(), patchi)
1159  {
1160  // Set ddtPhiCorr to 0 on non-coupled boundaries
1161  if
1162  (
1163  !this->mesh_.boundary()[patchi].coupled()
1164  || isA<cyclicAMIFvPatch>(this->mesh_.boundary()[patchi])
1165  )
1166  {
1167  phiCorrCoeffBf[patchi] = 0;
1168  }
1169  }
1170 
1171  addField
1172  (
1173  phase,
1174  "ddtCorrByA",
1175  -phiCorrCoeff*phiCorrs[phasei]*fvc::interpolate
1176  (
1177  byDt(alpha.oldTime()*phase.rho()().oldTime()*rAUs[phasei])
1178  ),
1179  ddtCorrByAs
1180  );
1181  }
1182 
1183  // Add virtual mass correction
1184  if (includeVirtualMass)
1185  {
1186  forAllConstIter(VmTable, Vms_, VmIter)
1187  {
1188  const volScalarField& Vm(*VmIter());
1189  const phasePair& pair(this->phasePairs_[VmIter.key()]);
1190 
1191  forAllConstIter(phasePair, pair, iter)
1192  {
1193  const phaseModel& phase = iter();
1194  const phaseModel& otherPhase = iter.otherPhase();
1195 
1196  addField
1197  (
1198  iter(),
1199  "ddtCorrByA",
1200  -fvc::interpolate(Vm*byDt(rAUs[phase.index()]))
1201  *(
1202  phiCorrs[phase.index()]
1203  + fvc::absolute
1204  (
1205  this->MRF().absolute(otherPhase.phi()),
1206  otherPhase.U()
1207  )
1208  - fvc::flux(otherPhase.U())
1209  - phiCorrs[otherPhase.index()]
1210  ),
1211  ddtCorrByAs
1212  );
1213  }
1214  }
1215  }
1216 
1217  return ddtCorrByAs;
1218 }
1219 
1220 
1221 template<class BasePhaseSystem>
1223 (
1224  const PtrList<volScalarField>& rAUs,
1225  const PtrList<volVectorField>& KdUByAs,
1226  const PtrList<surfaceScalarField>& alphafs,
1227  const PtrList<surfaceScalarField>& phiKdPhis
1228 )
1229 {
1230  Info<< "Inverting drag systems: ";
1231 
1232  phaseSystem::phaseModelList& phases = this->phaseModels_;
1233 
1234  // Calculate the mean velocity from the current velocity
1235  // of the moving phases
1236  volVectorField Um(this->movingPhases()[0]*this->movingPhases()[0].U());
1237 
1238  for
1239  (
1240  label movingPhasei=1;
1241  movingPhasei<this->movingPhases().size();
1242  movingPhasei++
1243  )
1244  {
1245  Um +=
1246  this->movingPhases()[movingPhasei]
1247  *this->movingPhases()[movingPhasei].U();
1248  }
1249 
1250  // Remove the drag contributions from the velocity and flux of the phases
1251  // in preparation for the partial elimination of these terms
1252  forAll(phases, i)
1253  {
1254  if (!phases[i].stationary())
1255  {
1256  phases[i].URef() += KdUByAs[i];
1257  phases[i].phiRef() += phiKdPhis[i];
1258  }
1259  }
1260 
1261  {
1262  // Create drag coefficient matrices
1263  PtrList<PtrList<volScalarField>> KdByAs(phases.size());
1264  PtrList<PtrList<surfaceScalarField>> phiKds(phases.size());
1265 
1266  forAll(phases, i)
1267  {
1268  KdByAs.set
1269  (
1270  i,
1271  new PtrList<volScalarField>(phases.size())
1272  );
1273 
1274  phiKds.set
1275  (
1276  i,
1277  new PtrList<surfaceScalarField>(phases.size())
1278  );
1279  }
1280 
1281  forAllConstIter(KdTable, Kds_, KdIter)
1282  {
1283  const volScalarField& K(*KdIter());
1284  const phasePair& pair(this->phasePairs_[KdIter.key()]);
1285 
1286  const label phase1i = pair.phase1().index();
1287  const label phase2i = pair.phase2().index();
1288 
1289  addField
1290  (
1291  pair.phase2(),
1292  "KdByA",
1293  -rAUs[phase1i]*K,
1294  KdByAs[phase1i]
1295  );
1296  addField
1297  (
1298  pair.phase1(),
1299  "KdByA",
1300  -rAUs[phase2i]*K,
1301  KdByAs[phase2i]
1302  );
1303 
1304  addField
1305  (
1306  pair.phase2(),
1307  "phiKd",
1308  fvc::interpolate(KdByAs[phase1i][phase2i]),
1309  phiKds[phase1i]
1310  );
1311  addField
1312  (
1313  pair.phase1(),
1314  "phiKd",
1315  fvc::interpolate(KdByAs[phase2i][phase1i]),
1316  phiKds[phase2i]
1317  );
1318  }
1319 
1320  forAll(phases, i)
1321  {
1322  this->fillFields("KdByAs", dimless, KdByAs[i]);
1323  this->fillFields("phiKds", dimless, phiKds[i]);
1324 
1325  KdByAs[i][i] = 1;
1326  phiKds[i][i] = 1;
1327  }
1328 
1329  // Decompose
1330  for (label i = 0; i < phases.size(); i++)
1331  {
1332  for (label j = i + 1; j < phases.size(); j++)
1333  {
1334  KdByAs[i][j] /= KdByAs[i][i];
1335  phiKds[i][j] /= phiKds[i][i];
1336  for (label k = i + 1; k < phases.size(); ++ k)
1337  {
1338  KdByAs[j][k] -= KdByAs[j][i]*KdByAs[i][k];
1339  phiKds[j][k] -= phiKds[j][i]*phiKds[i][k];
1340  }
1341  }
1342  }
1343 
1344  {
1345  volScalarField detKdByAs(KdByAs[0][0]);
1346  surfaceScalarField detPhiKdfs(phiKds[0][0]);
1347 
1348  for (label i = 1; i < phases.size(); i++)
1349  {
1350  detKdByAs *= KdByAs[i][i];
1351  detPhiKdfs *= phiKds[i][i];
1352  }
1353 
1354  Info<< "Min cell/face det = " << gMin(detKdByAs.primitiveField())
1355  << "/" << gMin(detPhiKdfs.primitiveField()) << endl;
1356  }
1357 
1358  // Solve for the velocities and fluxes
1359  for (label i = 1; i < phases.size(); i++)
1360  {
1361  if (!phases[i].stationary())
1362  {
1363  for (label j = 0; j < i; j ++)
1364  {
1365  phases[i].URef() -= KdByAs[i][j]*phases[j].U();
1366  phases[i].phiRef() -= phiKds[i][j]*phases[j].phi();
1367  }
1368  }
1369  }
1370  for (label i = phases.size() - 1; i >= 0; i--)
1371  {
1372  if (!phases[i].stationary())
1373  {
1374  for (label j = phases.size() - 1; j > i; j--)
1375  {
1376  phases[i].URef() -= KdByAs[i][j]*phases[j].U();
1377  phases[i].phiRef() -= phiKds[i][j]*phases[j].phi();
1378  }
1379  phases[i].URef() /= KdByAs[i][i];
1380  phases[i].phiRef() /= phiKds[i][i];
1381  }
1382  }
1383  }
1384 
1385  this->setMixtureU(Um);
1386  this->setMixturePhi(alphafs, this->phi());
1387 }
1388 
1389 
1390 template<class BasePhaseSystem>
1392 (
1393  const PtrList<surfaceScalarField>& rAUfs,
1394  const PtrList<surfaceScalarField>& alphafs,
1395  const PtrList<surfaceScalarField>& phiKdPhifs
1396 )
1397 {
1398  Info<< "Inverting drag system: ";
1399 
1400  phaseSystem::phaseModelList& phases = this->phaseModels_;
1401 
1402  // Remove the drag contributions from the flux of the phases
1403  // in preparation for the partial elimination of these terms
1404  forAll(phases, i)
1405  {
1406  if (!phases[i].stationary())
1407  {
1408  phases[i].phiRef() += phiKdPhifs[i];
1409  }
1410  }
1411 
1412  {
1413  // Create drag coefficient matrix
1414  PtrList<PtrList<surfaceScalarField>> phiKdfs(phases.size());
1415 
1416  forAll(phases, phasei)
1417  {
1418  phiKdfs.set
1419  (
1420  phasei,
1421  new PtrList<surfaceScalarField>(phases.size())
1422  );
1423  }
1424 
1425  forAllConstIter(KdfTable, Kdfs_, KdfIter)
1426  {
1427  const surfaceScalarField& K(*KdfIter());
1428  const phasePair& pair(this->phasePairs_[KdfIter.key()]);
1429 
1430  const label phase1i = pair.phase1().index();
1431  const label phase2i = pair.phase2().index();
1432 
1433  addField
1434  (
1435  pair.phase2(),
1436  "phiKdf",
1437  -rAUfs[phase1i]*K,
1438  phiKdfs[phase1i]
1439  );
1440  addField
1441  (
1442  pair.phase1(),
1443  "phiKdf",
1444  -rAUfs[phase2i]*K,
1445  phiKdfs[phase2i]
1446  );
1447  }
1448 
1449  forAll(phases, phasei)
1450  {
1451  this->fillFields("phiKdf", dimless, phiKdfs[phasei]);
1452 
1453  phiKdfs[phasei][phasei] = 1;
1454  }
1455 
1456  // Decompose
1457  for (label i = 0; i < phases.size(); i++)
1458  {
1459  for (label j = i + 1; j < phases.size(); j++)
1460  {
1461  phiKdfs[i][j] /= phiKdfs[i][i];
1462  for (label k = i + 1; k < phases.size(); ++ k)
1463  {
1464  phiKdfs[j][k] -= phiKdfs[j][i]*phiKdfs[i][k];
1465  }
1466  }
1467  }
1468 
1469  {
1470  surfaceScalarField detPhiKdfs(phiKdfs[0][0]);
1471 
1472  for (label i = 1; i < phases.size(); i++)
1473  {
1474  detPhiKdfs *= phiKdfs[i][i];
1475  }
1476 
1477  Info<< "Min face det = "
1478  << gMin(detPhiKdfs.primitiveField()) << endl;
1479  }
1480 
1481  // Solve for the fluxes
1482  for (label i = 1; i < phases.size(); i++)
1483  {
1484  if (!phases[i].stationary())
1485  {
1486  for (label j = 0; j < i; j ++)
1487  {
1488  phases[i].phiRef() -= phiKdfs[i][j]*phases[j].phi();
1489  }
1490  }
1491  }
1492  for (label i = phases.size() - 1; i >= 0; i--)
1493  {
1494  if (!phases[i].stationary())
1495  {
1496  for (label j = phases.size() - 1; j > i; j--)
1497  {
1498  phases[i].phiRef() -= phiKdfs[i][j]*phases[j].phi();
1499  }
1500  phases[i].phiRef() /= phiKdfs[i][i];
1501  }
1502  }
1503  }
1504 
1505  this->setMixturePhi(alphafs, this->phi());
1506 }
1507 
1508 
1509 template<class BasePhaseSystem>
1511 {
1512  if (BasePhaseSystem::read())
1513  {
1514  bool readOK = true;
1515 
1516  // Read models ...
1517 
1518  return readOK;
1519  }
1520  else
1521  {
1522  return false;
1523  }
1524 }
1525 
1526 
1527 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
tmp< GeometricField< Type, fvPatchField, volMesh > > DDt(const surfaceScalarField &phi, const GeometricField< Type, fvPatchField, volMesh > &psi)
Definition: fvcDDt.C:45
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
virtual PtrList< surfaceScalarField > phiFs(const PtrList< volScalarField > &rAUs)
Return the explicit force fluxes for the cell-based algorithm, that.
virtual tmp< surfaceScalarField > Kdf(const phasePairKey &key) const
Return the face drag coefficient for the phase pair.
const dimensionedScalar F
Faraday constant: default SI units: [C/mol].
label phasei
Definition: pEqn.H:27
IOMRFZoneList & MRF
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Type gMin(const FieldField< Field, Type > &f)
HashPtrTable< fvVectorMatrix > momentumTransferTable
Definition: phaseSystem.H:78
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransferf()
As momentumTransfer, but for the face-based algorithm.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
static int compare(const Pair< word > &a, const Pair< word > &b)
Compare Pairs.
Definition: Pair.H:153
virtual bool read()
Read base phaseProperties dictionary.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
fluid fillFields("rAU", dimTime/dimDensity, rAUs)
MomentumTransferPhaseSystem(const fvMesh &)
Construct from fvMesh.
Calculate the snGrad of the given volField.
const dimensionSet dimless
label k
Boltzmann constant.
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:58
CGAL::Exact_predicates_exact_constructions_kernel K
static const dimensionSet dimK
Coefficient dimensions.
dimensionedScalar posPart(const dimensionedScalar &ds)
Area-weighted average a surfaceField creating a volField.
virtual bool implicitPhasePressure() const
Returns true if the phase pressure is treated implicitly.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Calculate the first temporal derivative.
const dimensionSet dimTime
PtrList< surfaceScalarField > phiFs(fluid.phiFs(rAUs))
dynamicFvMesh & mesh
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
Calculate the face-flux of the given field.
static word groupName(Name name, const word &group)
virtual void partialEliminationf(const PtrList< surfaceScalarField > &rAUfs, const PtrList< surfaceScalarField > &alphafs, const PtrList< surfaceScalarField > &phiKdPhifs)
As partialElimination, but for the face-based algorithm. Only solves.
virtual tmp< volScalarField > Vm(const phasePairKey &key) const
Return the virtual mass coefficient for the phase pair.
const dimensionSet dimDensity
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
Calculate the matrix for the first temporal derivative.
PtrListDictionary< phaseModel > phaseModelList
Definition: phaseSystem.H:84
phi
Definition: correctPhi.H:3
virtual PtrList< volVectorField > KdUByAs(const PtrList< volScalarField > &rAUs) const
Return the explicit part of the drag force for the cell-based.
const dimensionSet dimForce
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
virtual PtrList< surfaceScalarField > phiKdPhis(const PtrList< volScalarField > &rAUs) const
Return the explicit drag force fluxes for the cell-based algorithm.
Calculate the divergence of the given field.
dimensionedScalar pos0(const dimensionedScalar &ds)
PtrList< surfaceScalarField > phiFfs(fluid.phiFfs(rAUfs))
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
const dimensionSet dimVelocity
const dimensionSet dimMass
virtual PtrList< surfaceScalarField > ddtCorrByAs(const PtrList< volScalarField > &rAUs, const bool includeVirtualMass=false) const
Return the flux corrections for the cell-based algorithm. These.
tmp< volScalarField > byDt(const volScalarField &vf)
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.
const tmp< volScalarField::Internal > & Sp
Definition: alphaSuSp.H:7
Calculate the mesh motion flux and convert fluxes from absolute to relative and back.
rho oldTime()
Calculate the matrix for the divergence of the given field and flux.
virtual ~MomentumTransferPhaseSystem()
Destructor.
virtual PtrList< surfaceScalarField > phiFfs(const PtrList< surfaceScalarField > &rAUfs)
As phiFs, but for the face-based algorithm.
label patchi
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:188
U
Definition: pEqn.H:72
static const dimensionSet dimK
Coefficient dimensions.
Definition: dragModel.H:93
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
virtual void partialElimination(const PtrList< volScalarField > &rAUs, const PtrList< volVectorField > &KdUByAs, const PtrList< surfaceScalarField > &alphafs, const PtrList< surfaceScalarField > &phiKdPhis)
Solve the drag system for the velocities and fluxes.
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:45
phaseSystem::phaseModelList & phases
Definition: createFields.H:12
messageStream Info
tmp< GeometricField< Type, fvPatchField, volMesh > > average(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Area-weighted average a surfaceField creating a volField.
Definition: fvcAverage.C:46
virtual tmp< volScalarField > Kd(const phasePairKey &key) const
Return the drag coefficient for the phase pair.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
A class for managing temporary objects.
Definition: PtrList.H:53
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Definition: fvcFlux.C:32
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
virtual PtrList< surfaceScalarField > AFfs() const
Return implicit force coefficients on the faces, for the face-based.
virtual PtrList< surfaceScalarField > DByAfs(const PtrList< volScalarField > &rAUs, const PtrList< surfaceScalarField > &rAUfs) const
Return the phase diffusivity.
void addDmdtUfs(const phaseSystem::dmdtfTable &dmdtfs, phaseSystem::momentumTransferTable &eqns)
Add momentum transfer terms which result from bulk mass transfers.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
Calculate the matrix for implicit and explicit sources.
dimensionedScalar negPart(const dimensionedScalar &ds)
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransfer()
Return the momentum transfer matrices for the cell-based algorithm.
virtual PtrList< surfaceScalarField > phiKdPhifs(const PtrList< surfaceScalarField > &rAUfs) const
As phiKdPhis, but for the face-based algorithm.