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-2020 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 "fvMatrix.H"
46 
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 
375  eqn -=
376  Vm
377  *(
378  fvm::ddt(U)
379  + fvm::div(phi, U)
380  - fvm::Sp(fvc::div(phi), U)
381  - otherPhase.DUDt()
382  )
383  + this->MRF_.DDt(Vm, U - otherPhase.U());
384  }
385  }
386  }
387 
388  return eqnsPtr;
389 }
390 
391 
392 template<class BasePhaseSystem>
395 {
396  // Create a momentum transfer matrix for each phase
397  autoPtr<phaseSystem::momentumTransferTable> eqnsPtr
398  (
400  );
401 
402  phaseSystem::momentumTransferTable& eqns = eqnsPtr();
403 
404  forAll(this->movingPhases(), movingPhasei)
405  {
406  const phaseModel& phase = this->movingPhases()[movingPhasei];
407 
408  eqns.insert
409  (
410  phase.name(),
411  new fvVectorMatrix(phase.U(), dimMass*dimVelocity/dimTime)
412  );
413  }
414 
415  // Create U & grad(U) fields
416  PtrList<fvVectorMatrix> UgradUs(this->phaseModels_.size());
417  forAll(this->phaseModels_, phasei)
418  {
419  const phaseModel& phase = this->phaseModels_[phasei];
420 
421  if (!phase.stationary())
422  {
423  const volVectorField& U = phase.U();
424 
425  UgradUs.set
426  (
427  phasei,
428  new fvVectorMatrix
429  (
430  fvm::div(phase.phi(), U)
431  - fvm::Sp(fvc::div(phase.phi()), U)
432  + this->MRF().DDt(U)
433  )
434  );
435  }
436  }
437 
438  // Add the virtual mass force
439  forAllConstIter(VmTable, Vms_, VmIter)
440  {
441  const volScalarField& Vm(*VmIter());
442  const phasePair& pair(this->phasePairs_[VmIter.key()]);
443 
444  forAllConstIter(phasePair, pair, iter)
445  {
446  const phaseModel& phase = iter();
447  const phaseModel& otherPhase = iter.otherPhase();
448 
449  if (!phase.stationary())
450  {
451  *eqns[phase.name()] -=
452  Vm
453  *(
454  UgradUs[phase.index()]
455  - (UgradUs[otherPhase.index()] & otherPhase.U())
456  );
457  }
458  }
459  }
460 
461  return eqnsPtr;
462 }
463 
464 
465 template<class BasePhaseSystem>
468 {
469  PtrList<surfaceScalarField> AFfs(this->phaseModels_.size());
470 
471  // Add the implicit part of the drag force
472  forAllConstIter(KdfTable, Kdfs_, KdfIter)
473  {
474  const surfaceScalarField& Kf(*KdfIter());
475  const phasePair& pair(this->phasePairs_[KdfIter.key()]);
476 
477  forAllConstIter(phasePair, pair, iter)
478  {
479  addField(iter(), "AFf", Kf, AFfs);
480  }
481  }
482 
483  // Add the implicit part of the virtual mass force
484  forAllConstIter(VmfTable, Vmfs_, VmfIter)
485  {
486  const surfaceScalarField& Vmf(*VmfIter());
487  const phasePair& pair(this->phasePairs_[VmfIter.key()]);
488 
489  forAllConstIter(phasePair, pair, iter)
490  {
491  addField(iter(), "AFf", byDt(Vmf), AFfs);
492  }
493  }
494 
495  this->fillFields("AFf", dimDensity/dimTime, AFfs);
496 
497  return AFfs;
498 }
499 
500 
501 template<class BasePhaseSystem>
504 (
505  const PtrList<volScalarField>& rAUs
506 )
507 {
508  PtrList<surfaceScalarField> phiFs(this->phaseModels_.size());
509 
510  // Add the lift force
512  (
513  liftModelTable,
514  liftModels_,
515  liftModelIter
516  )
517  {
518  const volVectorField F(liftModelIter()->F<vector>());
519  const phasePair& pair(this->phasePairs_[liftModelIter.key()]);
520 
521  addField
522  (
523  pair.phase1(),
524  "phiF",
525  fvc::flux(rAUs[pair.phase1().index()]*F),
526  phiFs
527  );
528  addField
529  (
530  pair.phase2(),
531  "phiF",
532  -fvc::flux(rAUs[pair.phase2().index()]*F),
533  phiFs
534  );
535  }
536 
537  // Add the wall lubrication force
539  (
540  wallLubricationModelTable,
541  wallLubricationModels_,
542  wallLubricationModelIter
543  )
544  {
545  const volVectorField F(wallLubricationModelIter()->F<vector>());
546  const phasePair&
547  pair(this->phasePairs_[wallLubricationModelIter.key()]);
548 
549  addField
550  (
551  pair.phase1(),
552  "phiF",
553  fvc::flux(rAUs[pair.phase1().index()]*F),
554  phiFs
555  );
556  addField
557  (
558  pair.phase2(),
559  "phiF",
560  -fvc::flux(rAUs[pair.phase2().index()]*F),
561  phiFs
562  );
563  }
564 
565  // Add the phase pressure
566  forAll(this->phaseModels_, phasei)
567  {
568  const phaseModel& phase = this->phaseModels_[phasei];
569 
570  const surfaceScalarField pPrimeByAf
571  (
572  fvc::interpolate(rAUs[phasei]*phase.pPrime())
573  );
574 
575  const surfaceScalarField snGradAlpha1
576  (
577  fvc::snGrad(phase)*this->mesh_.magSf()
578  );
579 
580  addField(phase, "phiF", pPrimeByAf*snGradAlpha1, phiFs);
581  }
582 
583  // Add the turbulent dispersion force
585  (
586  turbulentDispersionModelTable,
587  turbulentDispersionModels_,
588  turbulentDispersionModelIter
589  )
590  {
591  const phasePair&
592  pair(this->phasePairs_[turbulentDispersionModelIter.key()]);
593 
594  const volScalarField D(turbulentDispersionModelIter()->D());
595 
596  const surfaceScalarField DByA1f
597  (
598  fvc::interpolate(rAUs[pair.phase1().index()]*D)
599  );
600  const surfaceScalarField DByA2f
601  (
602  fvc::interpolate(rAUs[pair.phase2().index()]*D)
603  );
604 
605  const volScalarField alpha12(pair.phase1() + pair.phase2());
606  const surfaceScalarField snGradAlpha1By12
607  (
609  (
610  pair.phase1()/max(alpha12, pair.phase1().residualAlpha())
611  )*this->mesh_.magSf()
612  );
613  const surfaceScalarField snGradAlpha2By12
614  (
616  (
617  pair.phase2()/max(alpha12, pair.phase2().residualAlpha())
618  )*this->mesh_.magSf()
619  );
620 
621  addField(pair.phase1(), "phiF", DByA1f*snGradAlpha1By12, phiFs);
622  addField(pair.phase2(), "phiF", DByA2f*snGradAlpha2By12, phiFs);
623  }
624 
625  if (this->fillFields_)
626  {
628  }
629 
630  return phiFs;
631 }
632 
633 
634 template<class BasePhaseSystem>
637 (
638  const PtrList<surfaceScalarField>& rAUfs
639 )
640 {
641  PtrList<surfaceScalarField> phiFfs(this->phaseModels_.size());
642 
643  // Add the explicit part of the virtual mass force
644  forAllConstIter(VmfTable, Vmfs_, VmfIter)
645  {
646  const surfaceScalarField& Vmf(*VmfIter());
647  const phasePair& pair(this->phasePairs_[VmfIter.key()]);
648 
649  forAllConstIter(phasePair, pair, iter)
650  {
651  addField
652  (
653  iter(),
654  "phiFf",
655  -rAUfs[iter().index()]*Vmf
656  *(
657  byDt(this->MRF().absolute(iter().phi()().oldTime()))
658  + iter.otherPhase().DUDtf()
659  ),
660  phiFfs
661  );
662  }
663  }
664 
665  // Add the lift force
667  (
668  liftModelTable,
669  liftModels_,
670  liftModelIter
671  )
672  {
673  const surfaceScalarField Ff(liftModelIter()->Ff());
674  const phasePair& pair(this->phasePairs_[liftModelIter.key()]);
675 
676  addField
677  (
678  pair.phase1(),
679  "phiFf",
680  rAUfs[pair.phase1().index()]*Ff,
681  phiFfs
682  );
683  addField
684  (
685  pair.phase2(),
686  "phiFf",
687  -rAUfs[pair.phase2().index()]*Ff,
688  phiFfs
689  );
690  }
691 
692  // Add the wall lubrication force
694  (
695  wallLubricationModelTable,
696  wallLubricationModels_,
697  wallLubricationModelIter
698  )
699  {
700  const surfaceScalarField Ff(wallLubricationModelIter()->Ff());
701  const phasePair&
702  pair(this->phasePairs_[wallLubricationModelIter.key()]);
703 
704  addField
705  (
706  pair.phase1(),
707  "phiFf",
708  rAUfs[pair.phase1().index()]*Ff,
709  phiFfs
710  );
711  addField
712  (
713  pair.phase2(),
714  "phiFf",
715  -rAUfs[pair.phase2().index()]*Ff,
716  phiFfs
717  );
718  }
719 
720  // Add the phase pressure
721  forAll(this->phaseModels_, phasei)
722  {
723  const phaseModel& phase = this->phaseModels_[phasei];
724 
725  const surfaceScalarField pPrimeByAf
726  (
727  rAUfs[phasei]*fvc::interpolate(phase.pPrime())
728  );
729 
730  const surfaceScalarField snGradAlpha1
731  (
732  fvc::snGrad(phase)*this->mesh_.magSf()
733  );
734 
735  addField(phase, "phiFf", pPrimeByAf*snGradAlpha1, phiFfs);
736  }
737 
738  // Add the turbulent dispersion force
740  (
741  turbulentDispersionModelTable,
742  turbulentDispersionModels_,
743  turbulentDispersionModelIter
744  )
745  {
746  const surfaceScalarField Ff(turbulentDispersionModelIter()->Ff());
747  const phasePair&
748  pair(this->phasePairs_[turbulentDispersionModelIter.key()]);
749 
750  addField
751  (
752  pair.phase1(),
753  "phiFf",
754  rAUfs[pair.phase1().index()]*Ff,
755  phiFfs
756  );
757  addField
758  (
759  pair.phase2(),
760  "phiFf",
761  -rAUfs[pair.phase2().index()]*Ff,
762  phiFfs
763  );
764  }
765 
766  if (this->fillFields_)
767  {
769  }
770 
771  return phiFfs;
772 }
773 
774 
775 template<class BasePhaseSystem>
778 (
779  const PtrList<volScalarField>& rAUs
780 ) const
781 {
782  PtrList<surfaceScalarField> phiKdPhis(this->phaseModels_.size());
783 
784  // Add the explicit part of the drag force
785  forAllConstIter(KdTable, Kds_, KdIter)
786  {
787  const volScalarField& K(*KdIter());
788  const phasePair& pair(this->phasePairs_[KdIter.key()]);
789 
790  forAllConstIter(phasePair, pair, iter)
791  {
792  addField
793  (
794  iter(),
795  "phiKdPhi",
796  -fvc::interpolate(rAUs[iter().index()]*K)
797  *this->MRF().absolute(iter.otherPhase().phi()),
798  phiKdPhis
799  );
800  }
801  }
802 
803  if (this->fillFields_)
804  {
805  this->fillFields
806  (
807  "phiKdPhi",
809  phiKdPhis
810  );
811  }
812 
813  return phiKdPhis;
814 }
815 
816 
817 template<class BasePhaseSystem>
820 (
821  const PtrList<surfaceScalarField>& rAUfs
822 ) const
823 {
824  PtrList<surfaceScalarField> phiKdPhifs(this->phaseModels_.size());
825 
826  // Add the explicit part of the drag force
827  forAllConstIter(KdfTable, Kdfs_, KdfIter)
828  {
829  const surfaceScalarField& Kf(*KdfIter());
830  const phasePair& pair(this->phasePairs_[KdfIter.key()]);
831 
832  forAllConstIter(phasePair, pair, iter)
833  {
834  addField
835  (
836  iter(),
837  "phiKdPhif",
838  -rAUfs[iter().index()]*Kf
839  *this->MRF().absolute(iter.otherPhase().phi()),
840  phiKdPhifs
841  );
842  }
843  }
844 
845  if (this->fillFields_)
846  {
847  this->fillFields
848  (
849  "phiKdPhif",
851  phiKdPhifs
852  );
853  }
854 
855  return phiKdPhifs;
856 }
857 
858 
859 template<class BasePhaseSystem>
862 (
863  const PtrList<volScalarField>& rAUs
864 ) const
865 {
866  PtrList<volVectorField> KdUByAs(this->phaseModels_.size());
867 
868  // Add the explicit part of the drag force
869  forAllConstIter(KdTable, Kds_, KdIter)
870  {
871  const volScalarField& K(*KdIter());
872  const phasePair& pair(this->phasePairs_[KdIter.key()]);
873 
874  forAllConstIter(phasePair, pair, iter)
875  {
876  addField
877  (
878  iter(),
879  "KdUByA",
880  -rAUs[iter().index()]*K*iter.otherPhase().U(),
881  KdUByAs
882  );
883  }
884  }
885 
886  if (this->fillFields_)
887  {
888  this->fillFields("KdUByA", dimVelocity, KdUByAs);
889  }
890 
891  return KdUByAs;
892 }
893 
894 
895 template<class BasePhaseSystem>
897 implicitPhasePressure(const phaseModel& phase) const
898 {
899  return
900  this->mesh_.solverDict(phase.volScalarField::name()).
901  template lookupOrDefault<Switch>
902  (
903  "implicitPhasePressure",
904  false
905  );
906 }
907 
908 
909 template<class BasePhaseSystem>
911 implicitPhasePressure() const
912 {
913  bool implicitPressure = false;
914 
915  forAll(this->phaseModels_, phasei)
916  {
917  const phaseModel& phase = this->phaseModels_[phasei];
918 
919  implicitPressure = implicitPressure || implicitPhasePressure(phase);
920  }
921 
922  return implicitPressure;
923 }
924 
925 
926 template<class BasePhaseSystem>
929 (
930  const PtrList<volScalarField>& rAUs,
931  const PtrList<surfaceScalarField>& rAUfs
932 ) const
933 {
934  PtrList<surfaceScalarField> DByAfs(this->phaseModels_.size());
935 
936  if (rAUfs.size())
937  {
938  // Add the phase pressure
939  forAll(this->phaseModels_, phasei)
940  {
941  const phaseModel& phase = this->phaseModels_[phasei];
942 
943  const surfaceScalarField pPrimeByAf
944  (
945  rAUfs[phasei]*fvc::interpolate(phase.pPrime())
946  );
947 
948  addField(phase, "DByAf", pPrimeByAf, DByAfs);
949 
950  forAll(this->phaseModels_, phasej)
951  {
952  if (phasej != phasei)
953  {
954  const phaseModel& phase2 = this->phaseModels_[phasej];
955 
956  addField
957  (
958  phase2,
959  "DByAf",
961  (
962  phase2
963  /max(1 - phase, phase2.residualAlpha())
964  )*pPrimeByAf,
965  DByAfs
966  );
967  }
968  }
969  }
970 
971  // Add the turbulent dispersion
973  (
974  turbulentDispersionModelTable,
975  turbulentDispersionModels_,
976  turbulentDispersionModelIter
977  )
978  {
979  const phasePair&
980  pair(this->phasePairs_[turbulentDispersionModelIter.key()]);
981 
982  const surfaceScalarField Df
983  (
984  fvc::interpolate(turbulentDispersionModelIter()->D())
985  );
986 
987  const surfaceScalarField alpha12f
988  (
989  fvc::interpolate(pair.phase1() + pair.phase2())
990  );
991 
992  addField
993  (
994  pair.phase1(),
995  "DByAf",
996  rAUfs[pair.phase1().index()]*Df
997  /max(alpha12f, pair.phase1().residualAlpha()),
998  DByAfs
999  );
1000  addField
1001  (
1002  pair.phase2(),
1003  "DByAf",
1004  rAUfs[pair.phase2().index()]*Df
1005  /max(alpha12f, pair.phase2().residualAlpha()),
1006  DByAfs
1007  );
1008  }
1009  }
1010  else
1011  {
1012  // Add the phase pressure
1013  forAll(this->phaseModels_, phasei)
1014  {
1015  const phaseModel& phase = this->phaseModels_[phasei];
1016 
1017  const surfaceScalarField pPrimeByAf
1018  (
1019  fvc::interpolate(rAUs[phasei]*phase.pPrime())
1020  );
1021 
1022  addField(phase, "DByAf", pPrimeByAf, DByAfs);
1023 
1024  forAll(this->phaseModels_, phasej)
1025  {
1026  if (phasej != phasei)
1027  {
1028  const phaseModel& phase2 = this->phaseModels_[phasej];
1029 
1030  addField
1031  (
1032  phase2,
1033  "DByAf",
1035  (
1036  phase2
1037  /max(1 - phase, phase2.residualAlpha())
1038  )*pPrimeByAf,
1039  DByAfs
1040  );
1041  }
1042  }
1043  }
1044 
1045  // Add the turbulent dispersion
1047  (
1048  turbulentDispersionModelTable,
1049  turbulentDispersionModels_,
1050  turbulentDispersionModelIter
1051  )
1052  {
1053  const phasePair&
1054  pair(this->phasePairs_[turbulentDispersionModelIter.key()]);
1055 
1056  const volScalarField D(turbulentDispersionModelIter()->D());
1057 
1058  const volScalarField alpha12(pair.phase1() + pair.phase2());
1059 
1060  addField
1061  (
1062  pair.phase1(),
1063  "DByAf",
1065  (
1066  rAUs[pair.phase1().index()]*D
1067  /max(alpha12, pair.phase1().residualAlpha())
1068  ),
1069  DByAfs
1070  );
1071  addField
1072  (
1073  pair.phase2(),
1074  "DByAf",
1076  (
1077  rAUs[pair.phase2().index()]*D
1078  /max(alpha12, pair.phase2().residualAlpha())
1079  ),
1080  DByAfs
1081  );
1082  }
1083  }
1084 
1085  return DByAfs;
1086 }
1087 
1088 
1089 template<class BasePhaseSystem>
1092 (
1093  const PtrList<volScalarField>& rAUs,
1094  const bool includeVirtualMass
1095 ) const
1096 {
1097  PtrList<surfaceScalarField> ddtCorrByAs(this->phaseModels_.size());
1098 
1099  // Construct phi differences
1100  PtrList<surfaceScalarField> phiCorrs(this->phaseModels_.size());
1101  forAll(this->phaseModels_, phasei)
1102  {
1103  const phaseModel& phase = this->phaseModels_[phasei];
1104 
1105  phiCorrs.set
1106  (
1107  phasei,
1108  this->MRF().absolute(phase.phi()().oldTime())
1109  -fvc::flux(phase.U()().oldTime())
1110  );
1111  }
1112 
1113  // Add correction
1114  forAll(this->phaseModels_, phasei)
1115  {
1116  const phaseModel& phase = this->phaseModels_[phasei];
1117  const volScalarField& alpha = phase;
1118 
1119  // Apply ddtPhiCorr filter in pure(ish) phases
1120  surfaceScalarField alphafBar
1121  (
1123  );
1124 
1125  tmp<surfaceScalarField> phiCorrCoeff = pos0(alphafBar - 0.99);
1126 
1127  surfaceScalarField::Boundary& phiCorrCoeffBf =
1128  phiCorrCoeff.ref().boundaryFieldRef();
1129 
1130  forAll(this->mesh_.boundary(), patchi)
1131  {
1132  // Set ddtPhiCorr to 0 on non-coupled boundaries
1133  if
1134  (
1135  !this->mesh_.boundary()[patchi].coupled()
1136  || isA<cyclicAMIFvPatch>(this->mesh_.boundary()[patchi])
1137  )
1138  {
1139  phiCorrCoeffBf[patchi] = 0;
1140  }
1141  }
1142 
1143  addField
1144  (
1145  phase,
1146  "ddtCorrByA",
1147  -phiCorrCoeff*phiCorrs[phasei]*fvc::interpolate
1148  (
1149  byDt(alpha.oldTime()*phase.rho()().oldTime()*rAUs[phasei])
1150  ),
1151  ddtCorrByAs
1152  );
1153  }
1154 
1155  // Add virtual mass correction
1156  if (includeVirtualMass)
1157  {
1158  forAllConstIter(VmTable, Vms_, VmIter)
1159  {
1160  const volScalarField& Vm(*VmIter());
1161  const phasePair& pair(this->phasePairs_[VmIter.key()]);
1162 
1163  forAllConstIter(phasePair, pair, iter)
1164  {
1165  const phaseModel& phase = iter();
1166  const phaseModel& otherPhase = iter.otherPhase();
1167 
1168  addField
1169  (
1170  iter(),
1171  "ddtCorrByA",
1172  -fvc::interpolate(Vm*byDt(rAUs[phase.index()]))
1173  *(
1174  phiCorrs[phase.index()]
1175  + this->MRF().absolute(otherPhase.phi())
1176  - fvc::flux(otherPhase.U())
1177  - phiCorrs[otherPhase.index()]
1178  ),
1179  ddtCorrByAs
1180  );
1181  }
1182  }
1183  }
1184 
1185  return ddtCorrByAs;
1186 }
1187 
1188 
1189 template<class BasePhaseSystem>
1191 (
1192  const PtrList<volScalarField>& rAUs,
1193  const PtrList<volVectorField>& KdUByAs,
1194  const PtrList<surfaceScalarField>& alphafs,
1195  const PtrList<surfaceScalarField>& phiKdPhis
1196 )
1197 {
1198  Info<< "Inverting drag systems: ";
1199 
1200  phaseSystem::phaseModelList& phases = this->phaseModels_;
1201 
1202  // Calculate the mean velocity from the current velocity
1203  // of the moving phases
1204  volVectorField Um(this->movingPhases()[0]*this->movingPhases()[0].U());
1205 
1206  for
1207  (
1208  label movingPhasei=1;
1209  movingPhasei<this->movingPhases().size();
1210  movingPhasei++
1211  )
1212  {
1213  Um +=
1214  this->movingPhases()[movingPhasei]
1215  *this->movingPhases()[movingPhasei].U();
1216  }
1217 
1218  // Remove the drag contributions from the velocity and flux of the phases
1219  // in preparation for the partial elimination of these terms
1220  forAll(phases, i)
1221  {
1222  if (!phases[i].stationary())
1223  {
1224  phases[i].URef() += KdUByAs[i];
1225  phases[i].phiRef() += phiKdPhis[i];
1226  }
1227  }
1228 
1229  {
1230  // Create drag coefficient matrices
1231  PtrList<PtrList<volScalarField>> KdByAs(phases.size());
1232  PtrList<PtrList<surfaceScalarField>> phiKds(phases.size());
1233 
1234  forAll(phases, i)
1235  {
1236  KdByAs.set
1237  (
1238  i,
1239  new PtrList<volScalarField>(phases.size())
1240  );
1241 
1242  phiKds.set
1243  (
1244  i,
1245  new PtrList<surfaceScalarField>(phases.size())
1246  );
1247  }
1248 
1249  forAllConstIter(KdTable, Kds_, KdIter)
1250  {
1251  const volScalarField& K(*KdIter());
1252  const phasePair& pair(this->phasePairs_[KdIter.key()]);
1253 
1254  const label phase1i = pair.phase1().index();
1255  const label phase2i = pair.phase2().index();
1256 
1257  addField
1258  (
1259  pair.phase2(),
1260  "KdByA",
1261  -rAUs[phase1i]*K,
1262  KdByAs[phase1i]
1263  );
1264  addField
1265  (
1266  pair.phase1(),
1267  "KdByA",
1268  -rAUs[phase2i]*K,
1269  KdByAs[phase2i]
1270  );
1271 
1272  addField
1273  (
1274  pair.phase2(),
1275  "phiKd",
1276  fvc::interpolate(KdByAs[phase1i][phase2i]),
1277  phiKds[phase1i]
1278  );
1279  addField
1280  (
1281  pair.phase1(),
1282  "phiKd",
1283  fvc::interpolate(KdByAs[phase2i][phase1i]),
1284  phiKds[phase2i]
1285  );
1286  }
1287 
1288  forAll(phases, i)
1289  {
1290  this->fillFields("KdByAs", dimless, KdByAs[i]);
1291  this->fillFields("phiKds", dimless, phiKds[i]);
1292 
1293  KdByAs[i][i] = 1;
1294  phiKds[i][i] = 1;
1295  }
1296 
1297  // Decompose
1298  for (label i = 0; i < phases.size(); i++)
1299  {
1300  for (label j = i + 1; j < phases.size(); j++)
1301  {
1302  KdByAs[i][j] /= KdByAs[i][i];
1303  phiKds[i][j] /= phiKds[i][i];
1304  for (label k = i + 1; k < phases.size(); ++ k)
1305  {
1306  KdByAs[j][k] -= KdByAs[j][i]*KdByAs[i][k];
1307  phiKds[j][k] -= phiKds[j][i]*phiKds[i][k];
1308  }
1309  }
1310  }
1311 
1312  {
1313  volScalarField detKdByAs(KdByAs[0][0]);
1314  surfaceScalarField detPhiKdfs(phiKds[0][0]);
1315 
1316  for (label i = 1; i < phases.size(); i++)
1317  {
1318  detKdByAs *= KdByAs[i][i];
1319  detPhiKdfs *= phiKds[i][i];
1320  }
1321 
1322  Info<< "Min cell/face det = " << gMin(detKdByAs.primitiveField())
1323  << "/" << gMin(detPhiKdfs.primitiveField()) << endl;
1324  }
1325 
1326  // Solve for the velocities and fluxes
1327  for (label i = 1; i < phases.size(); i++)
1328  {
1329  if (!phases[i].stationary())
1330  {
1331  for (label j = 0; j < i; j ++)
1332  {
1333  phases[i].URef() -= KdByAs[i][j]*phases[j].U();
1334  phases[i].phiRef() -= phiKds[i][j]*phases[j].phi();
1335  }
1336  }
1337  }
1338  for (label i = phases.size() - 1; i >= 0; i--)
1339  {
1340  if (!phases[i].stationary())
1341  {
1342  for (label j = phases.size() - 1; j > i; j--)
1343  {
1344  phases[i].URef() -= KdByAs[i][j]*phases[j].U();
1345  phases[i].phiRef() -= phiKds[i][j]*phases[j].phi();
1346  }
1347  phases[i].URef() /= KdByAs[i][i];
1348  phases[i].phiRef() /= phiKds[i][i];
1349  }
1350  }
1351  }
1352 
1353  this->setMixtureU(Um);
1354  this->setMixturePhi(alphafs, this->phi());
1355 }
1356 
1357 
1358 template<class BasePhaseSystem>
1360 (
1361  const PtrList<surfaceScalarField>& rAUfs,
1362  const PtrList<surfaceScalarField>& alphafs,
1363  const PtrList<surfaceScalarField>& phiKdPhifs
1364 )
1365 {
1366  Info<< "Inverting drag system: ";
1367 
1368  phaseSystem::phaseModelList& phases = this->phaseModels_;
1369 
1370  // Remove the drag contributions from the flux of the phases
1371  // in preparation for the partial elimination of these terms
1372  forAll(phases, i)
1373  {
1374  if (!phases[i].stationary())
1375  {
1376  phases[i].phiRef() += phiKdPhifs[i];
1377  }
1378  }
1379 
1380  {
1381  // Create drag coefficient matrix
1382  PtrList<PtrList<surfaceScalarField>> phiKdfs(phases.size());
1383 
1384  forAll(phases, phasei)
1385  {
1386  phiKdfs.set
1387  (
1388  phasei,
1389  new PtrList<surfaceScalarField>(phases.size())
1390  );
1391  }
1392 
1393  forAllConstIter(KdfTable, Kdfs_, KdfIter)
1394  {
1395  const surfaceScalarField& K(*KdfIter());
1396  const phasePair& pair(this->phasePairs_[KdfIter.key()]);
1397 
1398  const label phase1i = pair.phase1().index();
1399  const label phase2i = pair.phase2().index();
1400 
1401  addField
1402  (
1403  pair.phase2(),
1404  "phiKdf",
1405  -rAUfs[phase1i]*K,
1406  phiKdfs[phase1i]
1407  );
1408  addField
1409  (
1410  pair.phase1(),
1411  "phiKdf",
1412  -rAUfs[phase2i]*K,
1413  phiKdfs[phase2i]
1414  );
1415  }
1416 
1417  forAll(phases, phasei)
1418  {
1419  this->fillFields("phiKdf", dimless, phiKdfs[phasei]);
1420 
1421  phiKdfs[phasei][phasei] = 1;
1422  }
1423 
1424  // Decompose
1425  for (label i = 0; i < phases.size(); i++)
1426  {
1427  for (label j = i + 1; j < phases.size(); j++)
1428  {
1429  phiKdfs[i][j] /= phiKdfs[i][i];
1430  for (label k = i + 1; k < phases.size(); ++ k)
1431  {
1432  phiKdfs[j][k] -= phiKdfs[j][i]*phiKdfs[i][k];
1433  }
1434  }
1435  }
1436 
1437  {
1438  surfaceScalarField detPhiKdfs(phiKdfs[0][0]);
1439 
1440  for (label i = 1; i < phases.size(); i++)
1441  {
1442  detPhiKdfs *= phiKdfs[i][i];
1443  }
1444 
1445  Info<< "Min face det = "
1446  << gMin(detPhiKdfs.primitiveField()) << endl;
1447  }
1448 
1449  // Solve for the fluxes
1450  for (label i = 1; i < phases.size(); i++)
1451  {
1452  if (!phases[i].stationary())
1453  {
1454  for (label j = 0; j < i; j ++)
1455  {
1456  phases[i].phiRef() -= phiKdfs[i][j]*phases[j].phi();
1457  }
1458  }
1459  }
1460  for (label i = phases.size() - 1; i >= 0; i--)
1461  {
1462  if (!phases[i].stationary())
1463  {
1464  for (label j = phases.size() - 1; j > i; j--)
1465  {
1466  phases[i].phiRef() -= phiKdfs[i][j]*phases[j].phi();
1467  }
1468  phases[i].phiRef() /= phiKdfs[i][i];
1469  }
1470  }
1471  }
1472 
1473  this->setMixturePhi(alphafs, this->phi());
1474 }
1475 
1476 
1477 template<class BasePhaseSystem>
1479 {
1480  if (BasePhaseSystem::read())
1481  {
1482  bool readOK = true;
1483 
1484  // Read models ...
1485 
1486  return readOK;
1487  }
1488  else
1489  {
1490  return false;
1491  }
1492 }
1493 
1494 
1495 // ************************************************************************* //
const dimensionedScalar & F
Faraday constant: default SI units: [C/mol].
rho oldTime()
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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.
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:75
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:142
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.
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.
phi
Definition: pEqn.H:104
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.
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.
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:81
virtual PtrList< volVectorField > KdUByAs(const PtrList< volScalarField > &rAUs) const
Return the explicit part of the drag force for the cell-based.
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)
const dimensionSet dimForce
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
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 dimensionSet dimDensity
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
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
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
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
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
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
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)
zeroField Sp
Definition: alphaSuSp.H:2
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.
const dimensionSet dimVelocity