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-2023 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 "dragModel.H"
29 #include "virtualMassModel.H"
30 #include "liftModel.H"
31 #include "wallLubricationModel.H"
33 
34 #include "fvmDdt.H"
35 #include "fvmDiv.H"
36 #include "fvmSup.H"
37 
38 #include "fvcDdt.H"
39 #include "fvcDiv.H"
40 #include "fvcFlux.H"
41 #include "fvcSnGrad.H"
42 #include "fvcMeshPhi.H"
43 #include "fvcReconstruct.H"
44 
45 #include "pimpleNoLoopControl.H"
46 
47 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
48 
49 template<class BasePhaseSystem>
51 (
52  const phaseSystem::dmdtfTable& dmdtfs,
54 )
55 {
56  forAllConstIter(phaseSystem::dmdtfTable, dmdtfs, dmdtfIter)
57  {
58  const phaseInterface interface(*this, dmdtfIter.key());
59 
60  const volScalarField& dmdtf = *dmdtfIter();
61  const volScalarField dmdtf21(posPart(dmdtf));
62  const volScalarField dmdtf12(negPart(dmdtf));
63 
64  phaseModel& phase1 = this->phases()[interface.phase1().name()];
65  phaseModel& phase2 = this->phases()[interface.phase2().name()];
66 
67  if (!phase1.stationary())
68  {
69  *eqns[phase1.name()] +=
70  dmdtf21*phase2.U() + fvm::Sp(dmdtf12, phase1.URef());
71  }
72 
73  if (!phase2.stationary())
74  {
75  *eqns[phase2.name()] -=
76  dmdtf12*phase1.U() + fvm::Sp(dmdtf21, phase2.URef());
77  }
78  }
79 }
80 
81 
82 template<class BasePhaseSystem>
84 (
86  const tmp<surfaceScalarField>& field
87 ) const
88 {
89  if (result.valid())
90  {
91  result.ref() += field;
92  }
93  else
94  {
95  if (field.isTmp())
96  {
97  result = field;
98  }
99  else
100  {
101  result = field().clone();
102  }
103  }
104 }
105 
106 
107 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
108 
109 template<class BasePhaseSystem>
112 (
113  const fvMesh& mesh
114 )
115 :
116  BasePhaseSystem(mesh)
117 {
118  this->generateInterfacialModels(dragModels_);
119  this->generateInterfacialModels(virtualMassModels_);
120  this->generateInterfacialModels(liftModels_);
121  this->generateInterfacialModels(wallLubricationModels_);
122  this->generateInterfacialModels(turbulentDispersionModels_);
123 }
124 
125 
126 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
127 
128 template<class BasePhaseSystem>
131 {}
132 
133 
134 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
135 
136 template<class BasePhaseSystem>
139 {
140  // Create a momentum transfer matrix for each phase
142  (
144  );
145 
146  phaseSystem::momentumTransferTable& eqns = eqnsPtr();
147 
148  forAll(this->movingPhases(), movingPhasei)
149  {
150  const phaseModel& phase = this->movingPhases()[movingPhasei];
151 
152  eqns.insert
153  (
154  phase.name(),
156  );
157  }
158 
159  // Initialise Kds table if required
160  if (!Kds_.size())
161  {
163  (
165  dragModels_,
166  dragModelIter
167  )
168  {
169  const phaseInterface& interface = dragModelIter()->interface();
170 
171  Kds_.insert
172  (
173  dragModelIter.key(),
174  new volScalarField
175  (
176  IOobject
177  (
178  IOobject::groupName("Kd", interface.name()),
179  this->mesh().time().name(),
180  this->mesh()
181  ),
182  this->mesh(),
184  )
185  );
186  }
187  }
188 
189  // Update the drag coefficients
191  (
193  dragModels_,
194  dragModelIter
195  )
196  {
197  *Kds_[dragModelIter.key()] = dragModelIter()->K();
198 
199  // Zero-gradient the drag coefficient to boundaries with fixed velocity
200  const phaseInterface& interface = dragModelIter()->interface();
201  volScalarField& K = *Kds_[dragModelIter.key()];
202 
203  forAll(K.boundaryField(), patchi)
204  {
205  if
206  (
207  (
208  !interface.phase1().stationary()
209  && interface.phase1().U()()
210  .boundaryField()[patchi].fixesValue()
211  )
212  && (
213  !interface.phase2().stationary()
214  && interface.phase2().U()()
215  .boundaryField()[patchi].fixesValue()
216  )
217  )
218  {
219  K.boundaryFieldRef()[patchi] =
220  K.boundaryField()[patchi].patchInternalField();
221  }
222  }
223  }
224 
225  // Initialise Vms table if required
226  if (!Vms_.size())
227  {
229  (
231  virtualMassModels_,
232  virtualMassModelIter
233  )
234  {
235  const phaseInterface& interface =
236  virtualMassModelIter()->interface();
237 
238  Vms_.insert
239  (
240  interface,
241  new volScalarField
242  (
243  IOobject
244  (
245  IOobject::groupName("Vm", interface.name()),
246  this->mesh().time().name(),
247  this->mesh()
248  ),
249  this->mesh(),
251  )
252  );
253  }
254  }
255 
256  // Update the virtual mass coefficients
258  (
260  virtualMassModels_,
261  virtualMassModelIter
262  )
263  {
264  *Vms_[virtualMassModelIter.key()] = virtualMassModelIter()->K();
265  }
266 
267  // Add the virtual mass force
268  forAllConstIter(VmTable, Vms_, VmIter)
269  {
270  const volScalarField& Vm(*VmIter());
271  const phaseInterface interface(*this, VmIter.key());
272 
273  forAllConstIter(phaseInterface, interface, iter)
274  {
275  const phaseModel& phase = iter();
276  const phaseModel& otherPhase = iter.otherPhase();
277 
278  if (!phase.stationary())
279  {
280  fvVectorMatrix& eqn = *eqns[phase.name()];
281 
282  const volVectorField& U = eqn.psi();
283 
284  const surfaceScalarField& phi = phase.phiRef();
285  const tmp<surfaceScalarField> taphi(fvc::absolute(phi, U));
286  const surfaceScalarField& aphi(taphi());
287 
288  const volScalarField VmPhase
289  (
290  (otherPhase/max(otherPhase, otherPhase.residualAlpha()))
291  *Vm
292  );
293 
294  eqn -=
295  VmPhase
296  *(
297  fvm::ddt(U)
298  + fvm::div(aphi, U) - fvm::Sp(fvc::div(aphi), U)
299  - otherPhase.DUDt()
300  )
301  + this->MRF_.DDt(VmPhase, U - otherPhase.U());
302  }
303  }
304  }
305 
306  return eqnsPtr;
307 }
308 
309 
310 template<class BasePhaseSystem>
313 {
314  // Create a momentum transfer matrix for each phase
316  (
318  );
319 
320  phaseSystem::momentumTransferTable& eqns = eqnsPtr();
321 
322  forAll(this->movingPhases(), movingPhasei)
323  {
324  const phaseModel& phase = this->movingPhases()[movingPhasei];
325 
326  eqns.insert
327  (
328  phase.name(),
330  );
331  }
332 
333  // Initialise Kdfs table if required
334  if (!Kdfs_.size())
335  {
337  (
339  dragModels_,
340  dragModelIter
341  )
342  {
343  const phaseInterface& interface = dragModelIter()->interface();
344 
345  Kdfs_.insert
346  (
347  dragModelIter.key(),
349  (
350  IOobject
351  (
352  IOobject::groupName("Kdf", interface.name()),
353  this->mesh().time().name(),
354  this->mesh()
355  ),
356  this->mesh(),
358  )
359  );
360  }
361  }
362 
363  // Update the drag coefficients
365  (
367  dragModels_,
368  dragModelIter
369  )
370  {
371  *Kdfs_[dragModelIter.key()] = dragModelIter()->Kf();
372  }
373 
374  // Initialise Vms table if required
375  if (!Vms_.size())
376  {
378  (
380  virtualMassModels_,
381  virtualMassModelIter
382  )
383  {
384  const phaseInterface& interface =
385  virtualMassModelIter()->interface();
386 
387  Vms_.insert
388  (
389  interface,
390  new volScalarField
391  (
392  IOobject
393  (
394  IOobject::groupName("Vm", interface.name()),
395  this->mesh().time().name(),
396  this->mesh()
397  ),
398  this->mesh(),
400  )
401  );
402  }
403  }
404 
405  // Update the virtual mass coefficients
407  (
409  virtualMassModels_,
410  virtualMassModelIter
411  )
412  {
413  *Vms_[virtualMassModelIter.key()] = virtualMassModelIter()->K();
414  }
415 
416  // Create U & grad(U) fields
417  PtrList<fvVectorMatrix> UgradUs(this->phaseModels_.size());
418  forAll(this->phaseModels_, phasei)
419  {
420  const phaseModel& phase = this->phaseModels_[phasei];
421 
422  if (!phase.stationary())
423  {
424  const volVectorField& U = phase.URef();
425  const surfaceScalarField& phi = phase.phiRef();
426 
427  const tmp<surfaceScalarField> taphi(fvc::absolute(phi, U));
428  const surfaceScalarField& aphi(taphi());
429 
430  UgradUs.set
431  (
432  phasei,
433  new fvVectorMatrix
434  (
435  fvm::div(aphi, U) - fvm::Sp(fvc::div(aphi), U)
436  + this->MRF().DDt(U)
437  )
438  );
439  }
440  }
441 
442  // Add the virtual mass force
443  forAllConstIter(VmTable, Vms_, VmIter)
444  {
445  const volScalarField& Vm(*VmIter());
446  const phaseInterface interface(*this, VmIter.key());
447 
448  forAllConstIter(phaseInterface, interface, iter)
449  {
450  const phaseModel& phase = iter();
451  const phaseModel& otherPhase = iter.otherPhase();
452 
453  if (!phase.stationary())
454  {
455  const volScalarField VmPhase
456  (
457  (otherPhase/max(otherPhase, otherPhase.residualAlpha()))
458  *Vm
459  );
460 
461  *eqns[phase.name()] -=
462  VmPhase
463  *(
464  UgradUs[phase.index()]
465  - (UgradUs[otherPhase.index()] & otherPhase.U())
466  );
467  }
468  }
469  }
470 
471  return eqnsPtr;
472 }
473 
474 
475 template<class BasePhaseSystem>
478 {
479  PtrList<surfaceScalarField> KdVmfs(this->phaseModels_.size());
480 
481  // Add the implicit part of the drag force
482  forAllConstIter(KdfTable, Kdfs_, KdfIter)
483  {
484  const surfaceScalarField& Kf(*KdfIter());
485  const phaseInterface interface(*this, KdfIter.key());
486 
487  forAllConstIter(phaseInterface, interface, iter)
488  {
489  const phaseModel& phase = iter();
490  const phaseModel& otherPhase = iter.otherPhase();
491 
492  const surfaceScalarField alphaf
493  (
494  fvc::interpolate(otherPhase)
495  );
496 
497  addField
498  (
499  phase,
500  "KdVmf",
501  (alphaf/max(alphaf, otherPhase.residualAlpha()))
502  *Kf,
503  KdVmfs
504  );
505  }
506  }
507 
508  // Add the implicit part of the virtual mass force
509  forAllConstIter(VmTable, Vms_, VmIter)
510  {
511  const volScalarField& Vm(*VmIter());
512  const phaseInterface interface(*this, VmIter.key());
513 
514  forAllConstIter(phaseInterface, interface, iter)
515  {
516  const phaseModel& phase = iter();
517  const phaseModel& otherPhase = iter.otherPhase();
518 
519  const volScalarField VmPhase
520  (
521  (otherPhase/max(otherPhase, otherPhase.residualAlpha()))
522  *Vm
523  );
524 
525  addField(phase, "KdVmf", byDt(fvc::interpolate(VmPhase)), KdVmfs);
526  }
527  }
528 
529  this->fillFields("KdVmf", dimDensity/dimTime, KdVmfs);
530 
531  return KdVmfs;
532 }
533 
534 
535 template<class BasePhaseSystem>
538 {
539  PtrList<surfaceScalarField> Fs(this->phaseModels_.size());
540 
541  // Add the lift force
543  (
545  liftModels_,
546  liftModelIter
547  )
548  {
549  const phaseInterface& interface = liftModelIter()->interface();
550 
551  const volVectorField F(liftModelIter()->F());
552 
553  addField
554  (
555  interface.phase1(),
556  "F",
557  fvc::flux(F),
558  Fs
559  );
560  addField
561  (
562  interface.phase2(),
563  "F",
564  -fvc::flux(F),
565  Fs
566  );
567  }
568 
569  // Add the wall lubrication force
571  (
573  wallLubricationModels_,
574  wallLubricationModelIter
575  )
576  {
577  const phaseInterface& interface =
578  wallLubricationModelIter()->interface();
579 
580  const volVectorField F(wallLubricationModelIter()->F());
581 
582  addField
583  (
584  interface.phase1(),
585  "F",
586  fvc::flux(F),
587  Fs
588  );
589  addField
590  (
591  interface.phase2(),
592  "F",
593  -fvc::flux(F),
594  Fs
595  );
596  }
597 
598  // Add the phase pressure
599  forAll(this->phaseModels_, phasei)
600  {
601  const phaseModel& phase = this->phaseModels_[phasei];
602 
603  const tmp<volScalarField> pPrime(phase.pPrime());
604 
605  addField
606  (
607  phase,
608  "F",
609  fvc::interpolate(pPrime(), pPrime().name())
610  *fvc::snGrad(phase)*this->mesh_.magSf(),
611  Fs
612  );
613  }
614 
615  // Add the turbulent dispersion force
617  (
619  turbulentDispersionModels_,
620  turbulentDispersionModelIter
621  )
622  {
623  const phaseInterface& interface =
624  turbulentDispersionModelIter()->interface();
625 
626  const surfaceScalarField Df
627  (
628  fvc::interpolate(turbulentDispersionModelIter()->D())
629  );
630 
631  const volScalarField alpha12(interface.phase1() + interface.phase2());
632  const surfaceScalarField snGradAlpha1By12
633  (
635  (
636  interface.phase1()
637  /max(alpha12, interface.phase1().residualAlpha())
638  )*this->mesh_.magSf()
639  );
640  const surfaceScalarField snGradAlpha2By12
641  (
643  (
644  interface.phase2()
645  /max(alpha12, interface.phase2().residualAlpha())
646  )*this->mesh_.magSf()
647  );
648 
649  addField(interface.phase1(), "F", Df*snGradAlpha1By12, Fs);
650  addField(interface.phase2(), "F", Df*snGradAlpha2By12, Fs);
651  }
652 
653  if (this->fillFields_)
654  {
655  this->fillFields("F", dimArea*dimDensity*dimAcceleration, Fs);
656  }
657 
658  return Fs;
659 }
660 
661 
662 template<class BasePhaseSystem>
665 {
666  PtrList<surfaceScalarField> Ffs(this->phaseModels_.size());
667 
668  // Add the explicit part of the virtual mass force
669  forAllConstIter(VmTable, Vms_, VmIter)
670  {
671  const volScalarField& Vm(*VmIter());
672  const phaseInterface interface(*this, VmIter.key());
673 
674  forAllConstIter(phaseInterface, interface, iter)
675  {
676  const phaseModel& phase = iter();
677  const phaseModel& otherPhase = iter.otherPhase();
678 
679  const volScalarField VmPhase
680  (
681  (otherPhase/max(otherPhase, otherPhase.residualAlpha()))
682  *Vm
683  );
684 
685  addField
686  (
687  phase,
688  "Ff",
689  -fvc::interpolate(VmPhase)
690  *(
691  byDt
692  (
694  (
695  this->MRF().absolute(iter().phi()().oldTime()),
696  iter().U()
697  )
698  )
699  + otherPhase.DUDtf()
700  ),
701  Ffs
702  );
703  }
704  }
705 
706  // Add the lift force
708  (
710  liftModels_,
711  liftModelIter
712  )
713  {
714  const phaseInterface& interface = liftModelIter()->interface();
715 
716  const surfaceScalarField Ff(liftModelIter()->Ff());
717 
718  addField
719  (
720  interface.phase1(),
721  "Ff",
722  Ff,
723  Ffs
724  );
725  addField
726  (
727  interface.phase2(),
728  "Ff",
729  -Ff,
730  Ffs
731  );
732  }
733 
734  // Add the wall lubrication force
736  (
738  wallLubricationModels_,
739  wallLubricationModelIter
740  )
741  {
742  const phaseInterface& interface =
743  wallLubricationModelIter()->interface();
744 
745  const surfaceScalarField Ff(wallLubricationModelIter()->Ff());
746 
747  addField
748  (
749  interface.phase1(),
750  "Ff",
751  Ff,
752  Ffs
753  );
754  addField
755  (
756  interface.phase2(),
757  "Ff",
758  -Ff,
759  Ffs
760  );
761  }
762 
763  // Add the phase pressure
764  forAll(this->phaseModels_, phasei)
765  {
766  const phaseModel& phase = this->phaseModels_[phasei];
767 
768  addField
769  (
770  phase,
771  "Ff",
772  fvc::interpolate(phase.pPrime())
773  *fvc::snGrad(phase)*this->mesh_.magSf(),
774  Ffs
775  );
776  }
777 
778  // Add the turbulent dispersion force
780  (
782  turbulentDispersionModels_,
783  turbulentDispersionModelIter
784  )
785  {
786  const phaseInterface& interface =
787  turbulentDispersionModelIter()->interface();
788 
789  const surfaceScalarField Df
790  (
791  fvc::interpolate(turbulentDispersionModelIter()->D())
792  );
793 
794  const volScalarField alpha12(interface.phase1() + interface.phase2());
795  const surfaceScalarField snGradAlpha1By12
796  (
798  (
799  interface.phase1()
800  /max(alpha12, interface.phase1().residualAlpha())
801  )*this->mesh_.magSf()
802  );
803  const surfaceScalarField snGradAlpha2By12
804  (
806  (
807  interface.phase2()
808  /max(alpha12, interface.phase2().residualAlpha())
809  )*this->mesh_.magSf()
810  );
811 
812  addField(interface.phase1(), "F", Df*snGradAlpha1By12, Ffs);
813  addField(interface.phase2(), "F", Df*snGradAlpha2By12, Ffs);
814  }
815 
816  if (this->fillFields_)
817  {
818  this->fillFields("Ff", dimArea*dimDensity*dimAcceleration, Ffs);
819  }
820 
821  return Ffs;
822 }
823 
824 
825 template<class BasePhaseSystem>
828 {
829  PtrList<surfaceScalarField> KdPhis(this->phaseModels_.size());
830 
831  // Add the explicit part of the drag force
832  forAllConstIter(KdTable, Kds_, KdIter)
833  {
834  const volScalarField& K(*KdIter());
835  const phaseInterface interface(*this, KdIter.key());
836 
837  forAllConstIter(phaseInterface, interface, iter)
838  {
839  const phaseModel& phase = iter();
840  const phaseModel& otherPhase = iter.otherPhase();
841 
842  addField
843  (
844  phase,
845  "KdPhi",
847  (
848  -(otherPhase/max(otherPhase, otherPhase.residualAlpha()))
849  *K
850  )
852  (
853  this->MRF().absolute(otherPhase.phi()),
854  otherPhase.U()
855  ),
856  KdPhis
857  );
858  }
859  }
860 
861  if (this->fillFields_)
862  {
863  this->fillFields
864  (
865  "KdPhi",
867  KdPhis
868  );
869  }
870 
871  return KdPhis;
872 }
873 
874 
875 template<class BasePhaseSystem>
878 {
879  PtrList<surfaceScalarField> KdPhifs(this->phaseModels_.size());
880 
881  // Add the explicit part of the drag force
882  forAllConstIter(KdfTable, Kdfs_, KdfIter)
883  {
884  const surfaceScalarField& Kf(*KdfIter());
885  const phaseInterface interface(*this, KdfIter.key());
886 
887  forAllConstIter(phaseInterface, interface, iter)
888  {
889  const phaseModel& phase = iter();
890  const phaseModel& otherPhase = iter.otherPhase();
891 
892  const surfaceScalarField alphaf
893  (
894  fvc::interpolate(otherPhase)
895  );
896 
897  addField
898  (
899  phase,
900  "KdPhif",
901  -(alphaf/max(alphaf, otherPhase.residualAlpha()))
902  *Kf
904  (
905  this->MRF().absolute(otherPhase.phi()),
906  otherPhase.U()
907  ),
908  KdPhifs
909  );
910  }
911  }
912 
913  if (this->fillFields_)
914  {
915  this->fillFields
916  (
917  "KdPhif",
919  KdPhifs
920  );
921  }
922 
923  return KdPhifs;
924 }
925 
926 
927 template<class BasePhaseSystem>
930 {
931  PtrList<volScalarField> Kds(this->phaseModels_.size());
932 
933  forAllConstIter(KdTable, Kds_, KdIter)
934  {
935  const volScalarField& K(*KdIter());
936  const phaseInterface interface(*this, KdIter.key());
937 
938  forAllConstIter(phaseInterface, interface, iter)
939  {
940  const phaseModel& phase = iter();
941  const phaseModel& otherPhase = iter.otherPhase();
942 
943  addField
944  (
945  phase,
946  "Kd",
947  (otherPhase/max(otherPhase, otherPhase.residualAlpha()))
948  *K,
949  Kds
950  );
951  }
952  }
953 
954  if (this->fillFields_)
955  {
956  this->fillFields("Kd", dimDensity/dimTime, Kds);
957  }
958 
959  return Kds;
960 }
961 
962 
963 template<class BasePhaseSystem>
966 {
967  PtrList<volVectorField> KdUs(this->phaseModels_.size());
968 
969  // Add the explicit part of the drag force
970  forAllConstIter(KdTable, Kds_, KdIter)
971  {
972  const volScalarField& K(*KdIter());
973  const phaseInterface interface(*this, KdIter.key());
974 
975  forAllConstIter(phaseInterface, interface, iter)
976  {
977  const phaseModel& phase = iter();
978  const phaseModel& otherPhase = iter.otherPhase();
979 
980  addField
981  (
982  phase,
983  "KdU",
984  -(otherPhase/max(otherPhase, otherPhase.residualAlpha()))
985  *K*otherPhase.U(),
986  KdUs
987  );
988  }
989  }
990 
991  if (this->fillFields_)
992  {
993  this->fillFields("KdU", dimDensity*dimAcceleration, KdUs);
994  }
995 
996  return KdUs;
997 }
998 
999 
1000 template<class BasePhaseSystem>
1002 implicitPhasePressure(const phaseModel& phase) const
1003 {
1004  return
1005  this->mesh_.solution().solverDict(phase.volScalarField::name()).
1006  template lookupOrDefault<Switch>
1007  (
1008  "implicitPhasePressure",
1009  false
1010  );
1011 }
1012 
1013 
1014 template<class BasePhaseSystem>
1016 implicitPhasePressure() const
1017 {
1018  bool implicitPressure = false;
1019 
1020  forAll(this->phaseModels_, phasei)
1021  {
1022  const phaseModel& phase = this->phaseModels_[phasei];
1023 
1024  implicitPressure = implicitPressure || implicitPhasePressure(phase);
1025  }
1026 
1027  return implicitPressure;
1028 }
1029 
1030 
1031 template<class BasePhaseSystem>
1034 (
1035  const PtrList<volScalarField>& rAUs,
1036  const PtrList<surfaceScalarField>& rAUfs
1037 ) const
1038 {
1039  tmp<surfaceScalarField> alphaDByAf;
1040 
1041  // Add the phase pressure
1042  forAll(this->phaseModels_, phasei)
1043  {
1044  const phaseModel& phase = this->phaseModels_[phasei];
1045 
1046  const tmp<volScalarField> pPrime(phase.pPrime());
1047 
1048  addTmpField
1049  (
1050  alphaDByAf,
1051  fvc::interpolate(max(phase, scalar(0)))
1052  *(
1053  rAUfs.size()
1054 
1055  // Face-momentum form
1056  ? rAUfs[phasei]*fvc::interpolate(pPrime())
1057 
1058  // Cell-momentum form
1059  : fvc::interpolate(rAUs[phasei]*pPrime(), pPrime().name())
1060  )
1061  );
1062  }
1063 
1064  // Add the turbulent dispersion
1066  (
1068  turbulentDispersionModels_,
1069  turbulentDispersionModelIter
1070  )
1071  {
1072  const phaseInterface& interface =
1073  turbulentDispersionModelIter()->interface();
1074 
1075  const surfaceScalarField alpha1f
1076  (
1077  fvc::interpolate(max(interface.phase1(), scalar(0)))
1078  );
1079 
1080  const surfaceScalarField alpha2f
1081  (
1082  fvc::interpolate(max(interface.phase2(), scalar(0)))
1083  );
1084 
1085  addTmpField
1086  (
1087  alphaDByAf,
1088  alpha1f*alpha2f
1089  /max(alpha1f + alpha2f, interface.phase1().residualAlpha())
1090  *(
1091  rAUfs.size()
1092 
1093  // Face-momentum form
1094  ? max
1095  (
1096  rAUfs[interface.phase1().index()],
1097  rAUfs[interface.phase2().index()]
1098  )
1099  *fvc::interpolate(turbulentDispersionModelIter()->D())
1100 
1101  // Cell-momentum form
1103  (
1104  max
1105  (
1106  rAUs[interface.phase1().index()],
1107  rAUs[interface.phase2().index()]
1108  )
1109  *turbulentDispersionModelIter()->D()
1110  )
1111  )
1112  );
1113  }
1114 
1115  return alphaDByAf;
1116 }
1117 
1118 
1119 template<class BasePhaseSystem>
1122 {
1123  PtrList<surfaceScalarField> ddtCorrs(this->phaseModels_.size());
1124 
1125  // Add correction
1126  forAll(this->movingPhases(), movingPhasei)
1127  {
1128  const phaseModel& phase = this->movingPhases()[movingPhasei];
1129 
1130  addField
1131  (
1132  phase,
1133  "ddtCorr",
1134  fvc::ddtCorr
1135  (
1136  phase,
1137  phase.rho(),
1138  phase.U()(),
1139  phase.phi()(),
1140  phase.Uf()
1141  ),
1142  ddtCorrs
1143  );
1144  }
1145 
1146 
1147  const pimpleNoLoopControl& pimple = this->pimple();
1148  const Switch VmDdtCorr
1149  (
1150  pimple.dict().lookupOrDefault<Switch>("VmDdtCorrection", false)
1151  );
1152 
1153  // Optionally ddd virtual mass correction
1154  if (VmDdtCorr)
1155  {
1156  PtrList<volScalarField> VmDdtCoeffs(this->phaseModels_.size());
1157  PtrList<surfaceScalarField> VmDdtCorrs(this->phaseModels_.size());
1158 
1159  forAll(this->movingPhases(), movingPhasei)
1160  {
1161  const phaseModel& phase = this->movingPhases()[movingPhasei];
1162  const label phasei = phase.index();
1163 
1164  VmDdtCoeffs.set
1165  (
1166  phasei,
1167  fvm::ddt(phase.U()())().A()
1168  );
1169 
1170  VmDdtCorrs.set
1171  (
1172  phasei,
1173  fvc::ddtCorr
1174  (
1175  phase.U()(),
1176  phase.phi()(),
1177  phase.Uf()
1178  )
1179  );
1180  }
1181 
1182  forAllConstIter(VmTable, Vms_, VmIter)
1183  {
1184  const volScalarField& Vm(*VmIter());
1185  const phaseInterface interface(*this, VmIter.key());
1186 
1187  forAllConstIter(phaseInterface, interface, iter)
1188  {
1189  const phaseModel& phase = iter();
1190  const label phasei = phase.index();
1191  const phaseModel& otherPhase = iter.otherPhase();
1192  const label otherPhasei = otherPhase.index();
1193 
1194  const volScalarField VmPhase
1195  (
1196  (otherPhase/max(otherPhase, otherPhase.residualAlpha()))
1197  *Vm
1198  );
1199 
1200  addField
1201  (
1202  phase,
1203  "ddtCorr",
1204  fvc::interpolate(VmPhase)
1205  *(
1206  VmDdtCorrs[phasei]
1207  + (
1208  fvc::interpolate(VmDdtCoeffs[otherPhasei])
1209  *(
1210  otherPhase.Uf().valid()
1211  ? (this->mesh_.Sf() & otherPhase.Uf()())()
1212  : otherPhase.phi()()
1213  )
1214  - fvc::flux(VmDdtCoeffs[otherPhasei]*otherPhase.U())
1215  )
1216  - VmDdtCorrs[otherPhase.index()]
1217  ),
1218  ddtCorrs
1219  );
1220  }
1221  }
1222  }
1223 
1224  return ddtCorrs;
1225 }
1226 
1227 
1228 template<class BasePhaseSystem>
1230 (
1231  PtrList<volVectorField>& dragCorrs,
1232  PtrList<surfaceScalarField>& dragCorrfs
1233 ) const
1234 {
1235  const phaseSystem::phaseModelList& phases = this->phaseModels_;
1236 
1237  PtrList<volVectorField> Uphis(phases.size());
1238 
1239  forAll(phases, i)
1240  {
1241  if (!phases[i].stationary())
1242  {
1243  Uphis.set
1244  (
1245  i,
1246  fvc::reconstruct(phases[i].phi())
1247  );
1248  }
1249  }
1250 
1251  forAllConstIter(KdTable, Kds_, KdIter)
1252  {
1253  const volScalarField& K(*KdIter());
1254  const phaseInterface interface(*this, KdIter.key());
1255 
1256  const phaseModel& phase1 = interface.phase1();
1257  const phaseModel& phase2 = interface.phase2();
1258 
1259  const label phase1i = phase1.index();
1260  const label phase2i = phase2.index();
1261 
1262  if (!phase1.stationary())
1263  {
1264  const volScalarField K1
1265  (
1266  (phase2/max(phase2, phase2.residualAlpha()))*K
1267  );
1268 
1269  addField
1270  (
1271  phase1,
1272  "dragCorr",
1273  K1
1274  *(
1275  phase2.stationary()
1276  ? -Uphis[phase1i]
1277  : (Uphis[phase2i] - Uphis[phase1i])
1278  ),
1279  dragCorrs
1280  );
1281 
1282  addField
1283  (
1284  phase1,
1285  "dragCorrf",
1287  *(
1288  phase2.stationary()
1289  ? -phase1.phi()
1290  : (phase2.phi() - phase1.phi())
1291  ),
1292  dragCorrfs
1293  );
1294  }
1295 
1296  if (!phase2.stationary())
1297  {
1298  const volScalarField K2
1299  (
1300  (phase1/max(phase1, phase1.residualAlpha()))*K
1301  );
1302 
1303  addField
1304  (
1305  phase2,
1306  "dragCorr",
1307  K2
1308  *(
1309  phase1.stationary()
1310  ? -Uphis[phase2i]
1311  : (Uphis[phase1i] - Uphis[phase2i])
1312  ),
1313  dragCorrs
1314  );
1315 
1316  addField
1317  (
1318  phase2,
1319  "dragCorrf",
1321  *(
1322  phase1.stationary()
1323  ? -phase2.phi()
1324  : (phase1.phi() - phase2.phi())
1325  ),
1326  dragCorrfs
1327  );
1328  }
1329  }
1330 }
1331 
1332 
1333 template<class BasePhaseSystem>
1335 (
1336  const PtrList<volScalarField>& rAUs,
1337  const PtrList<volVectorField>& KdUs,
1338  const PtrList<surfaceScalarField>& alphafs,
1339  const PtrList<surfaceScalarField>& rAUfs,
1340  const PtrList<surfaceScalarField>& KdPhis
1341 )
1342 {
1343  Info<< "Inverting drag systems: ";
1344 
1345  phaseSystem::phaseModelList& phases = this->phaseModels_;
1346 
1347  // Calculate the mean velocity from the current velocity
1348  // of the moving phases
1349  volVectorField Um(this->movingPhases()[0]*this->movingPhases()[0].U());
1350 
1351  for
1352  (
1353  label movingPhasei=1;
1354  movingPhasei<this->movingPhases().size();
1355  movingPhasei++
1356  )
1357  {
1358  Um +=
1359  this->movingPhases()[movingPhasei]
1360  *this->movingPhases()[movingPhasei].U();
1361  }
1362 
1363  // Remove the drag contributions from the velocity and flux of the phases
1364  // in preparation for the partial elimination of these terms
1365  forAll(phases, i)
1366  {
1367  if (!phases[i].stationary())
1368  {
1369  phases[i].URef() += rAUs[i]*KdUs[i];
1370  phases[i].phiRef() += rAUfs[i]*KdPhis[i];
1371  }
1372  }
1373 
1374  {
1375  // Create drag coefficient matrices
1376  PtrList<PtrList<volScalarField>> KdByAs(phases.size());
1377  PtrList<PtrList<surfaceScalarField>> KdByAfs(phases.size());
1378 
1379  forAll(phases, i)
1380  {
1381  KdByAs.set
1382  (
1383  i,
1384  new PtrList<volScalarField>(phases.size())
1385  );
1386 
1387  KdByAfs.set
1388  (
1389  i,
1390  new PtrList<surfaceScalarField>(phases.size())
1391  );
1392  }
1393 
1394  forAllConstIter(KdTable, Kds_, KdIter)
1395  {
1396  const volScalarField& K(*KdIter());
1397  const phaseInterface interface(*this, KdIter.key());
1398 
1399  const label phase1i = interface.phase1().index();
1400  const label phase2i = interface.phase2().index();
1401 
1402  const volScalarField K1
1403  (
1404  interface.phase2()
1405  /max(interface.phase2(), interface.phase2().residualAlpha())
1406  *K
1407  );
1408 
1409  const volScalarField K2
1410  (
1411  interface.phase1()
1412  /max(interface.phase1(), interface.phase1().residualAlpha())
1413  *K
1414  );
1415 
1416  addField
1417  (
1418  interface.phase2(),
1419  "KdByA",
1420  -rAUs[phase1i]*K1,
1421  KdByAs[phase1i]
1422  );
1423  addField
1424  (
1425  interface.phase1(),
1426  "KdByA",
1427  -rAUs[phase2i]*K2,
1428  KdByAs[phase2i]
1429  );
1430 
1431  addField
1432  (
1433  interface.phase2(),
1434  "KdByAf",
1435  -rAUfs[phase1i]*fvc::interpolate(K1),
1436  KdByAfs[phase1i]
1437  );
1438  addField
1439  (
1440  interface.phase1(),
1441  "KdByAf",
1442  -rAUfs[phase2i]*fvc::interpolate(K2),
1443  KdByAfs[phase2i]
1444  );
1445  }
1446 
1447  forAll(phases, i)
1448  {
1449  this->fillFields("KdByAs", dimless, KdByAs[i]);
1450  this->fillFields("KdByAfs", dimless, KdByAfs[i]);
1451 
1452  KdByAs[i][i] = 1;
1453  KdByAfs[i][i] = 1;
1454  }
1455 
1456  // Decompose
1457  for (label i = 0; i < phases.size(); i++)
1458  {
1459  for (label j = i + 1; j < phases.size(); j++)
1460  {
1461  KdByAs[i][j] /= KdByAs[i][i];
1462  KdByAfs[i][j] /= KdByAfs[i][i];
1463  for (label k = i + 1; k < phases.size(); ++ k)
1464  {
1465  KdByAs[j][k] -= KdByAs[j][i]*KdByAs[i][k];
1466  KdByAfs[j][k] -= KdByAfs[j][i]*KdByAfs[i][k];
1467  }
1468  }
1469  }
1470 
1471  {
1472  volScalarField detKdByAs(KdByAs[0][0]);
1473  surfaceScalarField detPhiKdfs(KdByAfs[0][0]);
1474 
1475  for (label i = 1; i < phases.size(); i++)
1476  {
1477  detKdByAs *= KdByAs[i][i];
1478  detPhiKdfs *= KdByAfs[i][i];
1479  }
1480 
1481  Info<< "Min cell/face det = " << gMin(detKdByAs.primitiveField())
1482  << "/" << gMin(detPhiKdfs.primitiveField()) << endl;
1483  }
1484 
1485  // Solve for the velocities and fluxes
1486  for (label i = 1; i < phases.size(); i++)
1487  {
1488  if (!phases[i].stationary())
1489  {
1490  for (label j = 0; j < i; j ++)
1491  {
1492  phases[i].URef() -= KdByAs[i][j]*phases[j].U();
1493  phases[i].phiRef() -= KdByAfs[i][j]*phases[j].phi();
1494  }
1495  }
1496  }
1497  for (label i = phases.size() - 1; i >= 0; i--)
1498  {
1499  if (!phases[i].stationary())
1500  {
1501  for (label j = phases.size() - 1; j > i; j--)
1502  {
1503  phases[i].URef() -= KdByAs[i][j]*phases[j].U();
1504  phases[i].phiRef() -= KdByAfs[i][j]*phases[j].phi();
1505  }
1506  phases[i].URef() /= KdByAs[i][i];
1507  phases[i].phiRef() /= KdByAfs[i][i];
1508  }
1509  }
1510  }
1511 
1512  this->setMixtureU(Um);
1513  this->setMixturePhi(alphafs, this->phi());
1514 }
1515 
1516 
1517 template<class BasePhaseSystem>
1519 (
1520  const PtrList<surfaceScalarField>& rAUfs,
1521  const PtrList<surfaceScalarField>& alphafs,
1522  const PtrList<surfaceScalarField>& KdPhifs
1523 )
1524 {
1525  Info<< "Inverting drag system: ";
1526 
1527  phaseSystem::phaseModelList& phases = this->phaseModels_;
1528 
1529  // Remove the drag contributions from the flux of the phases
1530  // in preparation for the partial elimination of these terms
1531  forAll(phases, i)
1532  {
1533  if (!phases[i].stationary())
1534  {
1535  phases[i].phiRef() += rAUfs[i]*KdPhifs[i];
1536  }
1537  }
1538 
1539  {
1540  // Create drag coefficient matrix
1541  PtrList<PtrList<surfaceScalarField>> phiKdfs(phases.size());
1542 
1543  forAll(phases, phasei)
1544  {
1545  phiKdfs.set
1546  (
1547  phasei,
1548  new PtrList<surfaceScalarField>(phases.size())
1549  );
1550  }
1551 
1552  forAllConstIter(KdfTable, Kdfs_, KdfIter)
1553  {
1554  const surfaceScalarField& Kf(*KdfIter());
1555  const phaseInterface interface(*this, KdfIter.key());
1556 
1557  const label phase1i = interface.phase1().index();
1558  const label phase2i = interface.phase2().index();
1559 
1560  const surfaceScalarField alpha1f
1561  (
1562  fvc::interpolate(interface.phase1())
1563  );
1564 
1565  const surfaceScalarField alpha2f
1566  (
1567  fvc::interpolate(interface.phase2())
1568  );
1569 
1570  addField
1571  (
1572  interface.phase2(),
1573  "phiKdf",
1574  -rAUfs[phase1i]
1575  *alpha2f/max(alpha2f, interface.phase2().residualAlpha())
1576  *Kf,
1577  phiKdfs[phase1i]
1578  );
1579  addField
1580  (
1581  interface.phase1(),
1582  "phiKdf",
1583  -rAUfs[phase2i]
1584  *alpha1f/max(alpha1f, interface.phase1().residualAlpha())
1585  *Kf,
1586  phiKdfs[phase2i]
1587  );
1588  }
1589 
1590  forAll(phases, phasei)
1591  {
1592  this->fillFields("phiKdf", dimless, phiKdfs[phasei]);
1593 
1594  phiKdfs[phasei][phasei] = 1;
1595  }
1596 
1597  // Decompose
1598  for (label i = 0; i < phases.size(); i++)
1599  {
1600  for (label j = i + 1; j < phases.size(); j++)
1601  {
1602  phiKdfs[i][j] /= phiKdfs[i][i];
1603  for (label k = i + 1; k < phases.size(); ++ k)
1604  {
1605  phiKdfs[j][k] -= phiKdfs[j][i]*phiKdfs[i][k];
1606  }
1607  }
1608  }
1609 
1610  {
1611  surfaceScalarField detPhiKdfs(phiKdfs[0][0]);
1612 
1613  for (label i = 1; i < phases.size(); i++)
1614  {
1615  detPhiKdfs *= phiKdfs[i][i];
1616  }
1617 
1618  Info<< "Min face det = "
1619  << gMin(detPhiKdfs.primitiveField()) << endl;
1620  }
1621 
1622  // Solve for the fluxes
1623  for (label i = 1; i < phases.size(); i++)
1624  {
1625  if (!phases[i].stationary())
1626  {
1627  for (label j = 0; j < i; j ++)
1628  {
1629  phases[i].phiRef() -= phiKdfs[i][j]*phases[j].phi();
1630  }
1631  }
1632  }
1633  for (label i = phases.size() - 1; i >= 0; i--)
1634  {
1635  if (!phases[i].stationary())
1636  {
1637  for (label j = phases.size() - 1; j > i; j--)
1638  {
1639  phases[i].phiRef() -= phiKdfs[i][j]*phases[j].phi();
1640  }
1641  phases[i].phiRef() /= phiKdfs[i][i];
1642  }
1643  }
1644  }
1645 
1646  this->setMixturePhi(alphafs, this->phi());
1647 }
1648 
1649 
1650 template<class BasePhaseSystem>
1652 {
1653  if (BasePhaseSystem::read())
1654  {
1655  bool readOK = true;
1656 
1657  // Read models ...
1658 
1659  return readOK;
1660  }
1661  else
1662  {
1663  return false;
1664  }
1665 }
1666 
1667 
1668 // ************************************************************************* //
static const Foam::dimensionedScalar D("D", Foam::dimTemperature, 257.14)
label k
#define K1
Definition: SHA1.C:167
#define K2
Definition: SHA1.C:168
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Generic GeometricField class.
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
static word groupName(Name name, const word &group)
virtual void partialElimination(const PtrList< volScalarField > &rAUs, const PtrList< volVectorField > &KdUs, const PtrList< surfaceScalarField > &alphafs, const PtrList< surfaceScalarField > &rAUfs, const PtrList< surfaceScalarField > &KdPhis)
Solve the drag system for the velocities and fluxes.
virtual void dragCorrs(PtrList< volVectorField > &dragCorrs, PtrList< surfaceScalarField > &dragCorrf) const
Set the cell and faces drag correction fields.
virtual PtrList< surfaceScalarField > ddtCorrs() const
Return the flux corrections for the cell-based algorithm. These.
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransfer()
Return the momentum transfer matrices for the cell-based algorithm.
virtual tmp< surfaceScalarField > alphaDByAf(const PtrList< volScalarField > &rAUs, const PtrList< surfaceScalarField > &rAUfs) const
Return the phase diffusivity.
virtual PtrList< volVectorField > KdUs() const
Return the explicit part of the drag force for the cell-based.
virtual void partialEliminationf(const PtrList< surfaceScalarField > &rAUfs, const PtrList< surfaceScalarField > &alphafs, const PtrList< surfaceScalarField > &KdPhifs)
As partialElimination, but for the face-based algorithm. Only solves.
virtual PtrList< volScalarField > Kds() const
Return the implicit part of the drag force.
void addDmdtUfs(const phaseSystem::dmdtfTable &dmdtfs, phaseSystem::momentumTransferTable &eqns)
Add momentum transfer terms which result from bulk mass transfers.
virtual PtrList< surfaceScalarField > Ffs() const
As Fs, but for the face-based algorithm.
virtual bool implicitPhasePressure() const
Returns true if the phase pressure is treated implicitly.
virtual PtrList< surfaceScalarField > Fs() const
Return the explicit force fluxes for the cell-based algorithm, that.
virtual PtrList< surfaceScalarField > KdVmfs() const
Return implicit force coefficients on the faces, for the face-based.
MomentumTransferPhaseSystem(const fvMesh &)
Construct from fvMesh.
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransferf()
As momentumTransfer, but for the face-based algorithm.
virtual PtrList< surfaceScalarField > KdPhifs() const
As KdPhis, but for the face-based algorithm.
virtual PtrList< surfaceScalarField > KdPhis() const
Return the explicit drag force fluxes for the cell-based algorithm.
virtual bool read()
Read base phaseProperties dictionary.
void addTmpField(tmp< surfaceScalarField > &result, const tmp< surfaceScalarField > &field) const
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
bool set(const label) const
Is element set.
Definition: PtrListI.H:65
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
static const dimensionSet dimK
Coefficient dimensions.
Definition: dragModel.H:81
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
VolField< Type > & psi()
Definition: fvMatrix.H:289
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
Class to represent an interface between phases. Derivations can further specify the configuration of ...
const phaseModel & phase1() const
Return phase 1.
const phaseModel & phase2() const
Return phase 2.
virtual const autoPtr< surfaceVectorField > & Uf() const =0
Return the face velocity.
const dimensionedScalar & residualAlpha() const
Return the residual phase-fraction for given phase.
Definition: phaseModel.C:133
virtual surfaceScalarField & phiRef()=0
Access the volumetric flux.
virtual bool stationary() const =0
Return whether the phase is stationary.
virtual tmp< volScalarField > pPrime() const =0
Return the phase-pressure'.
virtual tmp< surfaceScalarField > phi() const =0
Return the volumetric flux.
virtual tmp< volVectorField > DUDt() const =0
Return the substantive acceleration.
virtual tmp< surfaceScalarField > DUDtf() const =0
Return the substantive acceleration on the faces.
label index() const
Return the index of the phase.
Definition: phaseModel.C:121
virtual tmp< volVectorField > U() const =0
Return the velocity.
virtual volVectorField & URef()=0
Access the velocity.
const word & name() const
Return the name of this phase.
Definition: phaseModel.C:109
virtual const volScalarField & rho() const =0
Return the density field.
Pimple no-loop control class. Implements various option flags, but leaves loop controls to the deriva...
virtual const dictionary & dict() const
Return the solution dictionary.
A class for managing temporary objects.
Definition: tmp.H:55
bool valid() const
Is this temporary object valid,.
Definition: tmpI.H:167
bool isTmp() const
Return true if this is really a temporary object.
Definition: tmpI.H:153
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
static const dimensionSet dimK
Coefficient dimensions.
IOMRFZoneList MRF(mesh)
pimpleControl pimple(mesh)
Calculate the first temporal derivative.
Calculate the divergence of the given field.
Calculate the face-flux of the given field.
Calculate the mesh motion flux and convert fluxes from absolute to relative and back.
Reconstruct volField from a face flux field.
Calculate the snGrad of the given volField.
Calculate the matrix for the first temporal derivative.
Calculate the matrix for the divergence of the given field and flux.
Calculate the matrix for implicit and explicit sources.
label patchi
K
Definition: pEqn.H:75
U
Definition: pEqn.H:72
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
const dimensionedScalar F
Faraday constant: default SI units: [C/mol].
tmp< SurfaceField< typename innerProduct< vector, Type >::type > > flux(const VolField< Type > &vf)
Return the face-flux field obtained from the given volVectorField.
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< VolField< typename outerProduct< vector, Type >::type > > reconstruct(const SurfaceField< Type > &ssf)
tmp< VolField< Type > > Sp(const volScalarField &sp, const VolField< Type > &vf)
Definition: fvcSup.C:67
tmp< VolField< Type > > DDt(const surfaceScalarField &phi, const VolField< Type > &psi)
Definition: fvcDDt.C:45
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:202
tmp< SurfaceField< Type > > snGrad(const VolField< Type > &vf, const word &name)
Definition: fvcSnGrad.C:45
tmp< SurfaceField< typename Foam::flux< Type >::type > > ddtCorr(const VolField< Type > &U, const SurfaceField< Type > &Uf)
Definition: fvcDdt.C:196
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
void addField(const Group &group, const word &name, tmp< GeoField > field, PtrList< GeoField > &fieldList)
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const dimensionSet dimless
messageStream Info
const dimensionSet dimAcceleration
tmp< volScalarField > byDt(const volScalarField &vf)
Definition: phaseSystem.C:853
const dimensionSet dimTime
dimensionedScalar negPart(const dimensionedScalar &ds)
const dimensionSet dimDensity
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimMass
const dimensionSet dimVelocity
const dimensionSet dimArea
Type gMin(const FieldField< Field, Type > &f)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
dimensionedScalar posPart(const dimensionedScalar &ds)