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