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-2019 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 
28 #include "BlendedInterfacialModel.H"
29 #include "dragModel.H"
30 #include "virtualMassModel.H"
31 #include "liftModel.H"
32 #include "wallLubricationModel.H"
33 #include "turbulentDispersionModel.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 template<class BasePhaseSystem>
121 addMassTransferMomentumTransfer(phaseSystem::momentumTransferTable& eqns) const
122 {
124  (
125  phaseSystem::phasePairTable,
126  this->phasePairs_,
127  phasePairIter
128  )
129  {
130  const phasePair& pair(phasePairIter());
131 
132  if (pair.ordered())
133  {
134  continue;
135  }
136 
137  // Note that the phase UEqn contains a continuity error term, which
138  // implicitly adds a mass transfer term of fvm::Sp(dmdt, U). These
139  // additions do not include this term.
140 
141  const volScalarField dmdt(this->dmdt(pair));
142 
143  if (!pair.phase1().stationary())
144  {
145  fvVectorMatrix& eqn = *eqns[pair.phase1().name()];
146  const volScalarField dmdt21(posPart(dmdt));
147 
148  eqn += dmdt21*pair.phase2().U() - fvm::Sp(dmdt21, eqn.psi());
149  }
150 
151  if (!pair.phase2().stationary())
152  {
153  fvVectorMatrix& eqn = *eqns[pair.phase2().name()];
154  const volScalarField dmdt12(negPart(dmdt));
155 
156  eqn -= dmdt12*pair.phase1().U() - fvm::Sp(dmdt12, eqn.psi());
157  }
158  }
159 }
160 
161 
162 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
163 
164 template<class BasePhaseSystem>
167 (
168  const fvMesh& mesh
169 )
170 :
171  BasePhaseSystem(mesh)
172 {
173  this->generatePairsAndSubModels
174  (
175  "drag",
176  dragModels_
177  );
178 
179  this->generatePairsAndSubModels
180  (
181  "virtualMass",
182  virtualMassModels_
183  );
184 
185  this->generatePairsAndSubModels
186  (
187  "lift",
188  liftModels_
189  );
190 
191  this->generatePairsAndSubModels
192  (
193  "wallLubrication",
194  wallLubricationModels_
195  );
196 
197  this->generatePairsAndSubModels
198  (
199  "turbulentDispersion",
200  turbulentDispersionModels_
201  );
202 
204  (
205  dragModelTable,
206  dragModels_,
207  dragModelIter
208  )
209  {
210  const phasePair& pair(this->phasePairs_[dragModelIter.key()]);
211 
212  Kds_.insert
213  (
214  pair,
215  new volScalarField
216  (
217  IOobject::groupName("Kd", pair.name()),
218  dragModelIter()->K()
219  )
220  );
221 
222  Kdfs_.insert
223  (
224  pair,
226  (
227  IOobject::groupName("Kdf", pair.name()),
228  dragModelIter()->Kf()
229  )
230  );
231  }
232 
234  (
235  virtualMassModelTable,
236  virtualMassModels_,
237  virtualMassModelIter
238  )
239  {
240  const phasePair& pair(this->phasePairs_[virtualMassModelIter.key()]);
241 
242  Vms_.insert
243  (
244  pair,
245  new volScalarField
246  (
247  IOobject::groupName("Vm", pair.name()),
248  virtualMassModelIter()->K()
249  )
250  );
251 
252  Vmfs_.insert
253  (
254  pair,
256  (
257  IOobject::groupName("Vmf", pair.name()),
258  virtualMassModelIter()->Kf()
259  )
260  );
261  }
262 }
263 
264 
265 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
266 
267 template<class BasePhaseSystem>
270 {}
271 
272 
273 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
274 
275 template<class BasePhaseSystem>
278 {
279  // Create a momentum transfer matrix for each phase
280  autoPtr<phaseSystem::momentumTransferTable> eqnsPtr
281  (
283  );
284 
285  phaseSystem::momentumTransferTable& eqns = eqnsPtr();
286 
287  forAll(this->movingPhases(), movingPhasei)
288  {
289  const phaseModel& phase = this->movingPhases()[movingPhasei];
290 
291  eqns.insert
292  (
293  phase.name(),
294  new fvVectorMatrix(phase.U(), dimMass*dimVelocity/dimTime)
295  );
296  }
297 
298  // Update the drag coefficients
300  (
301  dragModelTable,
302  dragModels_,
303  dragModelIter
304  )
305  {
306  *Kds_[dragModelIter.key()] = dragModelIter()->K();
307  *Kdfs_[dragModelIter.key()] = dragModelIter()->Kf();
308  }
309 
310  // Add the implicit part of the drag force
311  forAllConstIter(KdTable, Kds_, KdIter)
312  {
313  const volScalarField& K(*KdIter());
314  const phasePair& pair(this->phasePairs_[KdIter.key()]);
315 
316  forAllConstIter(phasePair, pair, iter)
317  {
318  if (!iter().stationary())
319  {
320  fvVectorMatrix& eqn = *eqns[iter().name()];
321 
322  eqn -= fvm::Sp(K, eqn.psi());
323  }
324  }
325  }
326 
327  // Update the virtual mass coefficients
329  (
330  virtualMassModelTable,
331  virtualMassModels_,
332  virtualMassModelIter
333  )
334  {
335  *Vms_[virtualMassModelIter.key()] = virtualMassModelIter()->K();
336  *Vmfs_[virtualMassModelIter.key()] = virtualMassModelIter()->Kf();
337  }
338 
339  // Add the virtual mass force
340  forAllConstIter(VmTable, Vms_, VmIter)
341  {
342  const volScalarField& Vm(*VmIter());
343  const phasePair& pair(this->phasePairs_[VmIter.key()]);
344 
345  forAllConstIter(phasePair, pair, iter)
346  {
347  const phaseModel& phase = iter();
348  const phaseModel& otherPhase = iter.otherPhase();
349 
350  if (!phase.stationary())
351  {
352  fvVectorMatrix& eqn = *eqns[phase.name()];
353 
354  const volVectorField& U = eqn.psi();
355  const surfaceScalarField& phi = phase.phi();
356 
357  eqn -=
358  Vm
359  *(
360  fvm::ddt(U)
361  + fvm::div(phi, U)
362  - fvm::Sp(fvc::div(phi), U)
363  - otherPhase.DUDt()
364  )
365  + this->MRF_.DDt(Vm, U - otherPhase.U());
366  }
367  }
368  }
369 
370  // Add the source term due to mass transfer
371  addMassTransferMomentumTransfer(eqns);
372 
373  return eqnsPtr;
374 }
375 
376 
377 template<class BasePhaseSystem>
380 {
381  // Create a momentum transfer matrix for each phase
382  autoPtr<phaseSystem::momentumTransferTable> eqnsPtr
383  (
385  );
386 
387  phaseSystem::momentumTransferTable& eqns = eqnsPtr();
388 
389  forAll(this->movingPhases(), movingPhasei)
390  {
391  const phaseModel& phase = this->movingPhases()[movingPhasei];
392 
393  eqns.insert
394  (
395  phase.name(),
396  new fvVectorMatrix(phase.U(), dimMass*dimVelocity/dimTime)
397  );
398  }
399 
400  // Create U & grad(U) fields
401  PtrList<fvVectorMatrix> UgradUs(this->phaseModels_.size());
402  forAll(this->phaseModels_, phasei)
403  {
404  const phaseModel& phase = this->phaseModels_[phasei];
405 
406  if (!phase.stationary())
407  {
408  const volVectorField& U = phase.U();
409 
410  UgradUs.set
411  (
412  phasei,
413  new fvVectorMatrix
414  (
415  fvm::div(phase.phi(), U)
416  - fvm::Sp(fvc::div(phase.phi()), U)
417  + this->MRF().DDt(U)
418  )
419  );
420  }
421  }
422 
423  // Add the virtual mass force
424  forAllConstIter(VmTable, Vms_, VmIter)
425  {
426  const volScalarField& Vm(*VmIter());
427  const phasePair& pair(this->phasePairs_[VmIter.key()]);
428 
429  forAllConstIter(phasePair, pair, iter)
430  {
431  const phaseModel& phase = iter();
432  const phaseModel& otherPhase = iter.otherPhase();
433 
434  if (!phase.stationary())
435  {
436  *eqns[phase.name()] -=
437  Vm
438  *(
439  UgradUs[phase.index()]
440  - (UgradUs[otherPhase.index()] & otherPhase.U())
441  );
442  }
443  }
444  }
445 
446  // Add the source term due to mass transfer
447  addMassTransferMomentumTransfer(eqns);
448 
449  return eqnsPtr;
450 }
451 
452 
453 template<class BasePhaseSystem>
456 {
457  PtrList<surfaceScalarField> AFfs(this->phaseModels_.size());
458 
459  // Add the implicit part of the drag force
460  forAllConstIter(KdfTable, Kdfs_, KdfIter)
461  {
462  const surfaceScalarField& Kf(*KdfIter());
463  const phasePair& pair(this->phasePairs_[KdfIter.key()]);
464 
465  forAllConstIter(phasePair, pair, iter)
466  {
467  this->addField(iter(), "AFf", Kf, AFfs);
468  }
469  }
470 
471  // Add the implicit part of the virtual mass force
472  forAllConstIter(VmfTable, Vmfs_, VmfIter)
473  {
474  const surfaceScalarField& Vmf(*VmfIter());
475  const phasePair& pair(this->phasePairs_[VmfIter.key()]);
476 
477  forAllConstIter(phasePair, pair, iter)
478  {
479  this->addField(iter(), "AFf", byDt(Vmf), AFfs);
480  }
481  }
482 
483  if (this->fillFields_)
484  {
485  this->fillFields("AFf", dimDensity/dimTime, AFfs);
486  }
487 
488  return AFfs;
489 }
490 
491 
492 template<class BasePhaseSystem>
495 (
496  const PtrList<volScalarField>& rAUs
497 )
498 {
499  PtrList<surfaceScalarField> phiFs(this->phaseModels_.size());
500 
501  // Add the lift force
503  (
504  liftModelTable,
505  liftModels_,
506  liftModelIter
507  )
508  {
509  const volVectorField F(liftModelIter()->F<vector>());
510  const phasePair& pair(this->phasePairs_[liftModelIter.key()]);
511 
512  this->addField
513  (
514  pair.phase1(),
515  "phiF",
516  fvc::flux(rAUs[pair.phase1().index()]*F),
517  phiFs
518  );
519  this->addField
520  (
521  pair.phase2(),
522  "phiF",
523  - fvc::flux(rAUs[pair.phase2().index()]*F),
524  phiFs
525  );
526  }
527 
528  // Add the wall lubrication force
530  (
531  wallLubricationModelTable,
532  wallLubricationModels_,
533  wallLubricationModelIter
534  )
535  {
536  const volVectorField F(wallLubricationModelIter()->F<vector>());
537  const phasePair&
538  pair(this->phasePairs_[wallLubricationModelIter.key()]);
539 
540  this->addField
541  (
542  pair.phase1(),
543  "phiF",
544  fvc::flux(rAUs[pair.phase1().index()]*F),
545  phiFs
546  );
547  this->addField
548  (
549  pair.phase2(),
550  "phiF",
551  - fvc::flux(rAUs[pair.phase2().index()]*F),
552  phiFs
553  );
554  }
555 
556  // Add the phase pressure
557  DByAfs_.clear();
558  forAll(this->phaseModels_, phasei)
559  {
560  const phaseModel& phase = this->phaseModels_[phasei];
561 
562  const surfaceScalarField pPrimeByAf
563  (
564  fvc::interpolate(rAUs[phasei]*phase.pPrime())
565  );
566 
568  (
569  fvc::snGrad(phase)*this->mesh_.magSf()
570  );
571 
572  this->addField(phase, "phiF", pPrimeByAf*snGradAlpha1, phiFs);
573 
574  const bool implicitPhasePressure =
575  this->mesh_.solverDict(phase.volScalarField::name()).
576  template lookupOrDefault<Switch>
577  (
578  "implicitPhasePressure",
579  false
580  );
581 
582  if (implicitPhasePressure)
583  {
584  this->addField(phase, "DByAf", pPrimeByAf, DByAfs_);
585  }
586  }
587 
588  // Add the turbulent dispersion force
590  (
591  turbulentDispersionModelTable,
592  turbulentDispersionModels_,
593  turbulentDispersionModelIter
594  )
595  {
596  const phasePair&
597  pair(this->phasePairs_[turbulentDispersionModelIter.key()]);
598 
599  const volScalarField D(turbulentDispersionModelIter()->D());
600 
601  const surfaceScalarField DByA1f
602  (
603  fvc::interpolate(rAUs[pair.phase1().index()]*D)
604  );
605  const surfaceScalarField DByA2f
606  (
607  fvc::interpolate(rAUs[pair.phase2().index()]*D)
608  );
609 
611  (
612  fvc::snGrad(pair.phase1())*this->mesh_.magSf()
613  );
614 
615  this->addField(pair.phase1(), "phiF", DByA1f*snGradAlpha1, phiFs);
616  this->addField(pair.phase2(), "phiF", - DByA2f*snGradAlpha1, phiFs);
617 
618  if (DByAfs_.found(pair.phase1().name()))
619  {
620  this->addField(pair.phase1(), "DByAf", DByA1f, DByAfs_);
621  }
622  }
623 
624  if (this->fillFields_)
625  {
627  }
628 
629  return phiFs;
630 }
631 
632 
633 template<class BasePhaseSystem>
636 (
637  const PtrList<surfaceScalarField>& rAUfs
638 )
639 {
640  PtrList<surfaceScalarField> phiFfs(this->phaseModels_.size());
641 
642  // Add the explicit part of the virtual mass force
643  forAllConstIter(VmfTable, Vmfs_, VmfIter)
644  {
645  const surfaceScalarField& Vmf(*VmfIter());
646  const phasePair& pair(this->phasePairs_[VmfIter.key()]);
647 
648  forAllConstIter(phasePair, pair, iter)
649  {
650  this->addField
651  (
652  iter(),
653  "phiFf",
654  - rAUfs[iter().index()]*Vmf
655  *(
656  byDt(this->MRF().absolute(iter().phi()().oldTime()))
657  + iter.otherPhase().DUDtf()
658  ),
659  phiFfs
660  );
661  }
662  }
663 
664  // Add the lift force
666  (
667  liftModelTable,
668  liftModels_,
669  liftModelIter
670  )
671  {
672  const surfaceScalarField Ff(liftModelIter()->Ff());
673  const phasePair& pair(this->phasePairs_[liftModelIter.key()]);
674 
675  this->addField
676  (
677  pair.phase1(),
678  "phiFs",
679  rAUfs[pair.phase1().index()]*Ff,
680  phiFfs
681  );
682  this->addField
683  (
684  pair.phase2(),
685  "phiFf",
686  - rAUfs[pair.phase2().index()]*Ff,
687  phiFfs
688  );
689  }
690 
691  // Add the wall lubrication force
693  (
694  wallLubricationModelTable,
695  wallLubricationModels_,
696  wallLubricationModelIter
697  )
698  {
699  const surfaceScalarField Ff(wallLubricationModelIter()->Ff());
700  const phasePair&
701  pair(this->phasePairs_[wallLubricationModelIter.key()]);
702 
703  this->addField
704  (
705  pair.phase1(),
706  "phiFf",
707  rAUfs[pair.phase1().index()]*Ff,
708  phiFfs
709  );
710  this->addField
711  (
712  pair.phase2(),
713  "phiFf",
714  - rAUfs[pair.phase2().index()]*Ff,
715  phiFfs
716  );
717  }
718 
719  // Add the phase pressure
720  DByAfs_.clear();
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 
731  (
732  fvc::snGrad(phase)*this->mesh_.magSf()
733  );
734 
735  this->addField(phase, "phiFf", pPrimeByAf*snGradAlpha1, phiFfs);
736 
737  const bool implicitPhasePressure =
738  this->mesh_.solverDict(phase.volScalarField::name()).
739  template lookupOrDefault<Switch>
740  (
741  "implicitPhasePressure",
742  false
743  );
744 
745  if (implicitPhasePressure)
746  {
747  this->addField(phase, "DByAf", pPrimeByAf, DByAfs_);
748  }
749  }
750 
751  // Add the turbulent dispersion force and phase pressure
753  (
754  turbulentDispersionModelTable,
755  turbulentDispersionModels_,
756  turbulentDispersionModelIter
757  )
758  {
759  const phasePair&
760  pair(this->phasePairs_[turbulentDispersionModelIter.key()]);
761 
762  const volScalarField D(turbulentDispersionModelIter()->D());
763 
764  const surfaceScalarField DByAf1
765  (
766  rAUfs[pair.phase1().index()]*fvc::interpolate(D)
767  );
768  const surfaceScalarField DByAf2
769  (
770  rAUfs[pair.phase2().index()]*fvc::interpolate(D)
771  );
772 
774  (
775  fvc::snGrad(pair.phase1())*this->mesh_.magSf()
776  );
777 
778  this->addField(pair.phase1(), "phiFf", DByAf1*snGradAlpha1, phiFfs);
779  this->addField(pair.phase2(), "phiFf", - DByAf2*snGradAlpha1, phiFfs);
780 
781  if (DByAfs_.found(pair.phase1().name()))
782  {
783  this->addField(pair.phase1(), "DByAf", DByAf1, DByAfs_);
784  }
785  }
786 
787  if (this->fillFields_)
788  {
790  }
791 
792  return phiFfs;
793 }
794 
795 
796 template<class BasePhaseSystem>
799 (
800  const PtrList<volScalarField>& rAUs
801 ) const
802 {
803  PtrList<surfaceScalarField> phiKdPhis(this->phaseModels_.size());
804 
805  // Add the explicit part of the drag force
806  forAllConstIter(KdTable, Kds_, KdIter)
807  {
808  const volScalarField& K(*KdIter());
809  const phasePair& pair(this->phasePairs_[KdIter.key()]);
810 
811  forAllConstIter(phasePair, pair, iter)
812  {
813  this->addField
814  (
815  iter(),
816  "phiKdPhi",
817  - fvc::interpolate(rAUs[iter().index()]*K)
818  *this->MRF().absolute(iter.otherPhase().phi()),
819  phiKdPhis
820  );
821  }
822  }
823 
824  if (this->fillFields_)
825  {
826  this->fillFields
827  (
828  "phiKdPhi",
830  phiKdPhis
831  );
832  }
833 
834  return phiKdPhis;
835 }
836 
837 
838 template<class BasePhaseSystem>
841 (
842  const PtrList<surfaceScalarField>& rAUfs
843 ) const
844 {
845  PtrList<surfaceScalarField> phiKdPhifs(this->phaseModels_.size());
846 
847  // Add the explicit part of the drag force
848  forAllConstIter(KdfTable, Kdfs_, KdfIter)
849  {
850  const surfaceScalarField& Kf(*KdfIter());
851  const phasePair& pair(this->phasePairs_[KdfIter.key()]);
852 
853  forAllConstIter(phasePair, pair, iter)
854  {
855  this->addField
856  (
857  iter(),
858  "phiKdPhif",
859  - rAUfs[iter().index()]*Kf
860  *this->MRF().absolute(iter.otherPhase().phi()),
861  phiKdPhifs
862  );
863  }
864  }
865 
866  if (this->fillFields_)
867  {
868  this->fillFields
869  (
870  "phiKdPhif",
872  phiKdPhifs
873  );
874  }
875 
876  return phiKdPhifs;
877 }
878 
879 
880 template<class BasePhaseSystem>
883 (
884  const PtrList<volScalarField>& rAUs
885 ) const
886 {
887  PtrList<volVectorField> KdUByAs(this->phaseModels_.size());
888 
889  // Add the explicit part of the drag force
890  forAllConstIter(KdTable, Kds_, KdIter)
891  {
892  const volScalarField& K(*KdIter());
893  const phasePair& pair(this->phasePairs_[KdIter.key()]);
894 
895  forAllConstIter(phasePair, pair, iter)
896  {
897  this->addField
898  (
899  iter(),
900  "KdUByA",
901  - rAUs[iter().index()]*K*iter.otherPhase().U(),
902  KdUByAs
903  );
904  }
905  }
906 
907  if (this->fillFields_)
908  {
909  this->fillFields("KdUByA", dimVelocity, KdUByAs);
910  }
911 
912  return KdUByAs;
913 }
914 
915 
916 template<class BasePhaseSystem>
919 (
920  const PtrList<volScalarField>& rAUs,
921  const bool includeVirtualMass
922 ) const
923 {
924  PtrList<surfaceScalarField> ddtCorrByAs(this->phaseModels_.size());
925 
926  // Construct phi differences
927  PtrList<surfaceScalarField> phiCorrs(this->phaseModels_.size());
928  forAll(this->phaseModels_, phasei)
929  {
930  const phaseModel& phase = this->phaseModels_[phasei];
931 
932  phiCorrs.set
933  (
934  phasei,
935  this->MRF().absolute(phase.phi()().oldTime())
936  - fvc::flux(phase.U()().oldTime())
937  );
938  }
939 
940  // Add correction
941  forAll(this->phaseModels_, phasei)
942  {
943  const phaseModel& phase = this->phaseModels_[phasei];
944  const volScalarField& alpha = phase;
945 
946  // Apply ddtPhiCorr filter in pure(ish) phases
947  surfaceScalarField alphafBar
948  (
950  );
951 
952  tmp<surfaceScalarField> phiCorrCoeff = pos0(alphafBar - 0.99);
953 
954  surfaceScalarField::Boundary& phiCorrCoeffBf =
955  phiCorrCoeff.ref().boundaryFieldRef();
956 
957  forAll(this->mesh_.boundary(), patchi)
958  {
959  // Set ddtPhiCorr to 0 on non-coupled boundaries
960  if
961  (
962  !this->mesh_.boundary()[patchi].coupled()
963  || isA<cyclicAMIFvPatch>(this->mesh_.boundary()[patchi])
964  )
965  {
966  phiCorrCoeffBf[patchi] = 0;
967  }
968  }
969 
970  this->addField
971  (
972  phase,
973  "ddtCorrByA",
974  - phiCorrCoeff*phiCorrs[phasei]*fvc::interpolate
975  (
976  byDt(alpha.oldTime()*phase.rho()().oldTime()*rAUs[phasei])
977  ),
978  ddtCorrByAs
979  );
980  }
981 
982  // Add virtual mass correction
983  if (includeVirtualMass)
984  {
985  forAllConstIter(VmTable, Vms_, VmIter)
986  {
987  const volScalarField& Vm(*VmIter());
988  const phasePair& pair(this->phasePairs_[VmIter.key()]);
989 
990  forAllConstIter(phasePair, pair, iter)
991  {
992  const phaseModel& phase = iter();
993  const phaseModel& otherPhase = iter.otherPhase();
994 
995  this->addField
996  (
997  iter(),
998  "ddtCorrByA",
999  - fvc::interpolate(Vm*byDt(rAUs[phase.index()]))
1000  *(
1001  phiCorrs[phase.index()]
1002  + this->MRF().absolute(otherPhase.phi())
1003  - fvc::flux(otherPhase.U())
1004  - phiCorrs[otherPhase.index()]
1005  ),
1006  ddtCorrByAs
1007  );
1008  }
1009  }
1010  }
1011 
1012  return ddtCorrByAs;
1013 }
1014 
1015 
1016 template<class BasePhaseSystem>
1018 (
1019  const PtrList<volScalarField>& rAUs
1020 )
1021 {
1022  Info<< "Inverting drag systems: ";
1023 
1024  phaseSystem::phaseModelList& phases = this->phaseModels_;
1025 
1026  // Create drag coefficient matrices
1027  PtrList<PtrList<volScalarField>> KdByAs(phases.size());
1028  PtrList<PtrList<surfaceScalarField>> phiKds(phases.size());
1029 
1030  forAll(phases, phasei)
1031  {
1032  KdByAs.set
1033  (
1034  phasei,
1035  new PtrList<volScalarField>(phases.size())
1036  );
1037 
1038  phiKds.set
1039  (
1040  phasei,
1041  new PtrList<surfaceScalarField>(phases.size())
1042  );
1043  }
1044 
1045  forAllConstIter(KdTable, Kds_, KdIter)
1046  {
1047  const volScalarField& K(*KdIter());
1048  const phasePair& pair(this->phasePairs_[KdIter.key()]);
1049 
1050  const label phase1i = pair.phase1().index();
1051  const label phase2i = pair.phase2().index();
1052 
1053  this->addField
1054  (
1055  pair.phase2(),
1056  "KdByA",
1057  - rAUs[phase1i]*K,
1058  KdByAs[phase1i]
1059  );
1060  this->addField
1061  (
1062  pair.phase1(),
1063  "KdByA",
1064  - rAUs[phase2i]*K,
1065  KdByAs[phase2i]
1066  );
1067 
1068  this->addField
1069  (
1070  pair.phase2(),
1071  "phiKd",
1072  fvc::interpolate(KdByAs[phase1i][phase2i]),
1073  phiKds[phase1i]
1074  );
1075  this->addField
1076  (
1077  pair.phase1(),
1078  "phiKd",
1079  fvc::interpolate(KdByAs[phase2i][phase1i]),
1080  phiKds[phase2i]
1081  );
1082  }
1083 
1084  forAll(phases, phasei)
1085  {
1086  this->fillFields("KdByAs", dimless, KdByAs[phasei]);
1087  this->fillFields("phiKds", dimless, phiKds[phasei]);
1088 
1089  KdByAs[phasei][phasei] = 1;
1090  phiKds[phasei][phasei] = 1;
1091  }
1092 
1093  // Decompose
1094  for (label i = 0; i < phases.size(); ++ i)
1095  {
1096  for (label j = i + 1; j < phases.size(); ++ j)
1097  {
1098  KdByAs[i][j] /= KdByAs[i][i];
1099  phiKds[i][j] /= phiKds[i][i];
1100  for (label k = i + 1; k < phases.size(); ++ k)
1101  {
1102  KdByAs[j][k] -= KdByAs[j][i]*KdByAs[i][k];
1103  phiKds[j][k] -= phiKds[j][i]*phiKds[i][k];
1104  }
1105  }
1106  }
1107  {
1108  volScalarField detKdByAs(KdByAs[0][0]);
1109  surfaceScalarField detPhiKdfs(phiKds[0][0]);
1110  for (label i = 1; i < phases.size(); ++ i)
1111  {
1112  detKdByAs *= KdByAs[i][i];
1113  detPhiKdfs *= phiKds[i][i];
1114  }
1115  Info<< "Min cell/face det = " << gMin(detKdByAs.primitiveField())
1116  << "/" << gMin(detPhiKdfs.primitiveField()) << endl;
1117  }
1118 
1119  // Solve for the velocities and fluxes
1120  for (label i = 1; i < phases.size(); ++ i)
1121  {
1122  if (!phases[i].stationary())
1123  {
1124  for (label j = 0; j < i; j ++)
1125  {
1126  phases[i].URef() -= KdByAs[i][j]*phases[j].U();
1127  phases[i].phiRef() -= phiKds[i][j]*phases[j].phi();
1128  }
1129  }
1130  }
1131  for (label i = phases.size() - 1; i >= 0; i --)
1132  {
1133  if (!phases[i].stationary())
1134  {
1135  for (label j = phases.size() - 1; j > i; j --)
1136  {
1137  phases[i].URef() -= KdByAs[i][j]*phases[j].U();
1138  phases[i].phiRef() -= phiKds[i][j]*phases[j].phi();
1139  }
1140  phases[i].URef() /= KdByAs[i][i];
1141  phases[i].phiRef() /= phiKds[i][i];
1142  }
1143  }
1144 }
1145 
1146 
1147 template<class BasePhaseSystem>
1149 (
1150  const PtrList<surfaceScalarField>& rAUfs
1151 )
1152 {
1153  Info<< "Inverting drag system: ";
1154 
1155  phaseSystem::phaseModelList& phases = this->phaseModels_;
1156 
1157  // Create drag coefficient matrix
1158  PtrList<PtrList<surfaceScalarField>> phiKdfs(phases.size());
1159 
1160  forAll(phases, phasei)
1161  {
1162  phiKdfs.set
1163  (
1164  phasei,
1165  new PtrList<surfaceScalarField>(phases.size())
1166  );
1167  }
1168 
1169  forAllConstIter(KdfTable, Kdfs_, KdfIter)
1170  {
1171  const surfaceScalarField& K(*KdfIter());
1172  const phasePair& pair(this->phasePairs_[KdfIter.key()]);
1173 
1174  const label phase1i = pair.phase1().index();
1175  const label phase2i = pair.phase2().index();
1176 
1177  this->addField
1178  (
1179  pair.phase2(),
1180  "phiKdf",
1181  - rAUfs[phase1i]*K,
1182  phiKdfs[phase1i]
1183  );
1184  this->addField
1185  (
1186  pair.phase1(),
1187  "phiKdf",
1188  - rAUfs[phase2i]*K,
1189  phiKdfs[phase2i]
1190  );
1191  }
1192 
1193  forAll(phases, phasei)
1194  {
1195  this->fillFields("phiKdf", dimless, phiKdfs[phasei]);
1196 
1197  phiKdfs[phasei][phasei] = 1;
1198  }
1199 
1200  // Decompose
1201  for (label i = 0; i < phases.size(); ++ i)
1202  {
1203  for (label j = i + 1; j < phases.size(); ++ j)
1204  {
1205  phiKdfs[i][j] /= phiKdfs[i][i];
1206  for (label k = i + 1; k < phases.size(); ++ k)
1207  {
1208  phiKdfs[j][k] -= phiKdfs[j][i]*phiKdfs[i][k];
1209  }
1210  }
1211  }
1212  {
1213  surfaceScalarField detPhiKdfs(phiKdfs[0][0]);
1214  for (label i = 1; i < phases.size(); ++ i)
1215  {
1216  detPhiKdfs *= phiKdfs[i][i];
1217  }
1218  Info<< "Min face det = " << gMin(detPhiKdfs.primitiveField()) << endl;
1219  }
1220 
1221  // Solve for the fluxes
1222  for (label i = 1; i < phases.size(); ++ i)
1223  {
1224  if (!phases[i].stationary())
1225  {
1226  for (label j = 0; j < i; j ++)
1227  {
1228  phases[i].phiRef() -= phiKdfs[i][j]*phases[j].phi();
1229  }
1230  }
1231  }
1232  for (label i = phases.size() - 1; i >= 0; i --)
1233  {
1234  if (!phases[i].stationary())
1235  {
1236  for (label j = phases.size() - 1; j > i; j --)
1237  {
1238  phases[i].phiRef() -= phiKdfs[i][j]*phases[j].phi();
1239  }
1240  phases[i].phiRef() /= phiKdfs[i][i];
1241  }
1242  }
1243 }
1244 
1245 
1246 template<class BasePhaseSystem>
1249 {
1250  return DByAfs_;
1251 }
1252 
1253 
1254 template<class BasePhaseSystem>
1256 {
1257  if (BasePhaseSystem::read())
1258  {
1259  bool readOK = true;
1260 
1261  // Read models ...
1262 
1263  return readOK;
1264  }
1265  else
1266  {
1267  return false;
1268  }
1269 }
1270 
1271 
1272 // ************************************************************************* //
Class which models interfacial momenum transfer between a number of phases. Drag, virtual mass...
surfaceScalarField Vmf("Vmf", fluid.Vmf())
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.
surfaceScalarField & phi
label phasei
Definition: pEqn.H:27
IOMRFZoneList & MRF
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
virtual void partialEliminationf(const PtrList< surfaceScalarField > &rAUfs)
As partialElimination, but for the face-based algorithm. Only solves.
virtual bool read()
Read base phaseProperties dictionary.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
MomentumTransferPhaseSystem(const fvMesh &)
Construct from fvMesh.
Calculate the snGrad of the given volField.
label k
Boltzmann constant.
A HashTable specialization for hashing pointers.
Definition: HashPtrTable.H:50
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
PtrList< surfaceScalarField > AFfs(fluid.AFfs())
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
CGAL::Exact_predicates_exact_constructions_kernel K
dimensionedScalar posPart(const dimensionedScalar &ds)
Area-weighted average a surfaceField creating a volField.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Calculate the first temporal derivative.
PtrList< surfaceScalarField > phiFs(fluid.phiFs(rAUs))
virtual const HashPtrTable< surfaceScalarField > & DByAfs() const
Return the phase diffusivities divided by the momentum coefficients.
Calculate the face-flux of the given field.
static word groupName(Name name, const word &group)
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.
volVectorField F(fluid.F())
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
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
surfaceScalarField snGradAlpha1(fvc::snGrad(alpha1) *mesh.magSf())
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
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:45
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
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
fluid fillFields("rAU", dimTime/dimDensity, rAUs)
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
surfaceScalarField Ff(fluid.Ff())
virtual PtrList< surfaceScalarField > AFfs() const
Return implicit force coefficients on the faces, for the face-based.
PtrList< surfaceScalarField > phiFfs(fluid.phiFfs(rAUfs))
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 void partialElimination(const PtrList< volScalarField > &rAUs)
Solve the drag system for the velocities and fluxes.
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