momentumTransferSystem.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-2025 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 
26 #include "momentumTransferSystem.H"
27 
28 #include "dragModel.H"
29 #include "virtualMassModel.H"
30 #include "liftModel.H"
31 #include "wallLubricationModel.H"
34 
35 #include "fvmDdt.H"
36 #include "fvmDiv.H"
37 #include "fvmSup.H"
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 "scalarMatrices.H"
46 #include "pimpleNoLoopControl.H"
47 
48 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
49 
50 namespace Foam
51 {
53 }
54 
55 
57 (
58  "momentumTransfer"
59 );
60 
61 
62 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
63 
64 Foam::IOobject Foam::momentumTransferSystem::io(const phaseSystem& fluid) const
65 {
66  typeIOobject<IOdictionary> result
67  (
69  fluid.mesh().time().constant(),
70  fluid.mesh(),
73  );
74 
75  // Check for and warn about momentum transfer models being in the old
76  // location in constant/phaseProperties
77  const wordList modelEntries
78  ({
79  modelName<blendedDragModel>(),
80  modelName<blendedVirtualMassModel>(),
81  modelName<blendedLiftModel>(),
82  modelName<blendedWallLubricationModel>(),
83  modelName<blendedTurbulentDispersionModel>()
84  });
85 
86  OStringStream modelEntriesString;
87  forAll(modelEntries, i)
88  {
89  modelEntriesString<< modelEntries[i];
90  if (i < modelEntries.size() - 2) modelEntriesString<< ", ";
91  if (i == modelEntries.size() - 2) modelEntriesString<< " and ";
92  }
93 
94  if (!result.headerOk())
95  {
96  result.readOpt() = IOobject::NO_READ;
97 
99  << "Specifying momentum transfer model entries - "
100  << modelEntriesString.str().c_str() << " - in "
101  << fluid.relativeObjectPath() << " is deprecated. These "
102  << "entries should now be specified in "
103  << result.relativeObjectPath() << "." << endl;
104  }
105  else
106  {
107  forAll(modelEntries, i)
108  {
109  if (!fluid.found(modelEntries[i])) continue;
110 
112  << "Momentum transfer model entries - "
113  << modelEntriesString.str().c_str() << " - in "
114  << fluid.relativeObjectPath() << " are no longer used. These "
115  << "entries are now read from " << result.relativeObjectPath()
116  << "." << endl;
117 
118  break;
119  }
120  }
121 
122  return result;
123 }
124 
125 
126 template<class ModelType>
127 const Foam::dictionary& Foam::momentumTransferSystem::modelsDict() const
128 {
129  const word key = modelName<ModelType>();
130 
131  return
132  readOpt() == IOobject::NO_READ && fluid_.found(key)
133  ? fluid_.subDict(key)
134  : subDict(key);
135 }
136 
137 
138 void Foam::momentumTransferSystem::addTmpField
139 (
140  tmp<surfaceScalarField>& result,
141  const tmp<surfaceScalarField>& field
142 ) const
143 {
144  if (result.valid())
145  {
146  result.ref() += field;
147  }
148  else
149  {
150  if (field.isTmp())
151  {
152  result = field;
153  }
154  else
155  {
156  result = field().clone();
157  }
158  }
159 }
160 
161 
162 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
163 
165 :
166  IOdictionary(io(fluid)),
167  fluid_(fluid),
168  dragModels_
169  (
171  (
172  fluid_,
173  modelsDict<blendedDragModel>()
174  )
175  ),
176  virtualMassModels_
177  (
179  (
180  fluid_,
181  modelsDict<blendedVirtualMassModel>()
182  )
183  ),
184  liftModels_
185  (
187  (
188  fluid_,
189  modelsDict<blendedLiftModel>()
190  )
191  ),
192  wallLubricationModels_
193  (
195  (
196  fluid_,
197  modelsDict<blendedWallLubricationModel>()
198  )
199  ),
200  turbulentDispersionModels_
201  (
203  (
204  fluid_,
205  modelsDict<blendedTurbulentDispersionModel>()
206  )
207  )
208 {}
209 
210 
211 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
212 
214 {}
215 
216 
217 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
218 
221 {
222  PtrList<surfaceScalarField> Fs(fluid_.phases().size());
223 
224  // Add the lift force
226  (
228  liftModels_,
229  liftModelIter
230  )
231  {
232  const phaseInterface& interface = liftModelIter()->interface();
233 
234  const volVectorField F(liftModelIter()->F());
235 
236  addField
237  (
238  interface.phase1(),
239  "F",
240  fvc::flux(F),
241  Fs
242  );
243  addField
244  (
245  interface.phase2(),
246  "F",
247  -fvc::flux(F),
248  Fs
249  );
250  }
251 
252  // Add the wall lubrication force
254  (
256  wallLubricationModels_,
257  wallLubricationModelIter
258  )
259  {
260  const phaseInterface& interface =
261  wallLubricationModelIter()->interface();
262 
263  const volVectorField F(wallLubricationModelIter()->F());
264 
265  addField
266  (
267  interface.phase1(),
268  "F",
269  fvc::flux(F),
270  Fs
271  );
272  addField
273  (
274  interface.phase2(),
275  "F",
276  -fvc::flux(F),
277  Fs
278  );
279  }
280 
281  // Add the phase pressure
282  forAll(fluid_.movingPhases(), movingPhasei)
283  {
284  const phaseModel& phase = fluid_.movingPhases()[movingPhasei];
285 
286  addField
287  (
288  phase,
289  "F",
290  phase.pPrimef()*fvc::snGrad(phase)*fluid_.mesh().magSf(),
291  Fs
292  );
293  }
294 
295  // Add the turbulent dispersion force
297  (
299  turbulentDispersionModels_,
300  turbulentDispersionModelIter
301  )
302  {
303  const phaseInterface& interface =
304  turbulentDispersionModelIter()->interface();
305 
306  const surfaceScalarField Df
307  (
308  fvc::interpolate(turbulentDispersionModelIter()->D())
309  );
310 
311  const volScalarField alpha12(interface.phase1() + interface.phase2());
312  const surfaceScalarField snGradAlpha1By12
313  (
315  (
316  interface.phase1()
317  /max(alpha12, interface.phase1().residualAlpha())
318  )*fluid_.mesh().magSf()
319  );
320  const surfaceScalarField snGradAlpha2By12
321  (
323  (
324  interface.phase2()
325  /max(alpha12, interface.phase2().residualAlpha())
326  )*fluid_.mesh().magSf()
327  );
328 
329  addField(interface.phase1(), "F", Df*snGradAlpha1By12, Fs);
330  addField(interface.phase2(), "F", Df*snGradAlpha2By12, Fs);
331  }
332 
333  return Fs;
334 }
335 
336 
339 {
340  PtrList<surfaceScalarField> Ffs(fluid_.phases().size());
341 
342  // Add the lift force
344  (
346  liftModels_,
347  liftModelIter
348  )
349  {
350  const phaseInterface& interface = liftModelIter()->interface();
351 
352  const surfaceScalarField Ff(liftModelIter()->Ff());
353 
354  addField
355  (
356  interface.phase1(),
357  "Ff",
358  Ff,
359  Ffs
360  );
361  addField
362  (
363  interface.phase2(),
364  "Ff",
365  -Ff,
366  Ffs
367  );
368  }
369 
370  // Add the wall lubrication force
372  (
374  wallLubricationModels_,
375  wallLubricationModelIter
376  )
377  {
378  const phaseInterface& interface =
379  wallLubricationModelIter()->interface();
380 
381  const surfaceScalarField Ff(wallLubricationModelIter()->Ff());
382 
383  addField
384  (
385  interface.phase1(),
386  "Ff",
387  Ff,
388  Ffs
389  );
390  addField
391  (
392  interface.phase2(),
393  "Ff",
394  -Ff,
395  Ffs
396  );
397  }
398 
399  // Add the phase pressure
400  forAll(fluid_.movingPhases(), movingPhasei)
401  {
402  const phaseModel& phase = fluid_.movingPhases()[movingPhasei];
403 
404  addField
405  (
406  phase,
407  "Ff",
408  phase.pPrimef()*fvc::snGrad(phase)*fluid_.mesh().magSf(),
409  Ffs
410  );
411  }
412 
413  // Add the turbulent dispersion force
415  (
417  turbulentDispersionModels_,
418  turbulentDispersionModelIter
419  )
420  {
421  const phaseInterface& interface =
422  turbulentDispersionModelIter()->interface();
423 
424  const surfaceScalarField Df
425  (
426  fvc::interpolate(turbulentDispersionModelIter()->D())
427  );
428 
429  const volScalarField alpha12(interface.phase1() + interface.phase2());
430  const surfaceScalarField snGradAlpha1By12
431  (
433  (
434  interface.phase1()
435  /max(alpha12, interface.phase1().residualAlpha())
436  )*fluid_.mesh().magSf()
437  );
438  const surfaceScalarField snGradAlpha2By12
439  (
441  (
442  interface.phase2()
443  /max(alpha12, interface.phase2().residualAlpha())
444  )*fluid_.mesh().magSf()
445  );
446 
447  addField(interface.phase1(), "F", Df*snGradAlpha1By12, Ffs);
448  addField(interface.phase2(), "F", Df*snGradAlpha2By12, Ffs);
449  }
450 
451  return Ffs;
452 }
453 
454 
455 void Foam::momentumTransferSystem::invADVs
456 (
458 ) const
459 {
460  const label n = ADVs.size();
461 
462  scalarSquareMatrix AD(n);
463  scalarField source(n);
464  labelList pivotIndices(n);
465 
466  forAll(ADVs[0][0], ci)
467  {
468  for (label i=0; i<n; i++)
469  {
470  for (label j=0; j<n; j++)
471  {
472  AD(i, j) = ADVs[i][j][ci];
473  }
474  }
475 
476  // Calculate the inverse of AD using LD decomposition
477  // and back-substitution
478  LUDecompose(AD, pivotIndices);
479 
480  for (label j=0; j<n; j++)
481  {
482  source = Zero;
483  source[j] = 1;
484 
485  LUBacksubstitute(AD, pivotIndices, source);
486 
487  for (label i=0; i<n; i++)
488  {
489  ADVs[i][j][ci] = source[i];
490  }
491  }
492  }
493 }
494 
495 
496 template<class GeoMesh>
497 void Foam::momentumTransferSystem::invADVs
498 (
499  PtrList<PtrList<GeometricField<scalar, GeoMesh>>>& ADVs
500 ) const
501 {
502  const label n = ADVs.size();
503 
504  List<UPtrList<scalarField>> ADps(n);
505 
506  for (label i=0; i<n; i++)
507  {
508  ADps[i].setSize(n);
509 
510  for (label j=0; j<n; j++)
511  {
512  ADps[i].set(j, &ADVs[i][j]);
513 
514  ADVs[i][j].dimensions().reset(dimless/ADVs[i][j].dimensions());
515  }
516  }
517 
518  invADVs(ADps);
519 
520  forAll(ADVs[0][0].boundaryField(), patchi)
521  {
522  for (label i=0; i<n; i++)
523  {
524  for (label j=0; j<n; j++)
525  {
526  ADps[i].set(j, &ADVs[i][j].boundaryFieldRef()[patchi]);
527  }
528  }
529 
530  invADVs(ADps);
531  }
532 }
533 
534 
535 void Foam::momentumTransferSystem::invADVs
536 (
537  const PtrList<volScalarField>& As,
541 ) const
542 {
543  const label n = As.size();
544 
545  invADVs.setSize(n);
546  invADVfs.setSize(n);
547 
548  forAll(invADVs, i)
549  {
550  invADVs.set(i, new PtrList<volScalarField>(n));
551  invADVs[i].set(i, As[i].clone());
552 
553  invADVfs.set(i, new PtrList<surfaceScalarField>(n));
554  invADVfs[i].set(i, fvc::interpolate(As[i]));
555  }
556 
557  labelList movingPhases(fluid_.phases().size(), -1);
558  forAll(fluid_.movingPhases(), movingPhasei)
559  {
560  movingPhases[fluid_.movingPhases()[movingPhasei].index()] =
561  movingPhasei;
562  }
563 
564  Kds_.clear();
565 
566  // Update the drag coefficients
568  (
570  dragModels_,
571  dragModelIter
572  )
573  {
574  const phaseInterface& interface = dragModelIter()->interface();
575 
576  tmp<volScalarField> tKd(dragModelIter()->K());
577  volScalarField& Kd = tKd.ref();
578  Kds_.insert(dragModelIter.key(), tKd.ptr());
579 
580  // Zero-gradient the drag coefficient to boundaries with fixed velocity
582  {
583  if
584  (
585  (
586  !interface.phase1().stationary()
587  && interface.phase1().U()()
588  .boundaryField()[patchi].fixesValue()
589  )
590  && (
591  !interface.phase2().stationary()
592  && interface.phase2().U()()
593  .boundaryField()[patchi].fixesValue()
594  )
595  )
596  {
597  Kd.boundaryFieldRef()[patchi] =
598  Kd.boundaryField()[patchi].patchInternalField();
599  }
600  }
601 
602  forAllConstIter(phaseInterface, interface, iter)
603  {
604  const phaseModel& phase = iter();
605  const phaseModel& otherPhase = iter.otherPhase();
606 
607  const label i = movingPhases[phase.index()];
608 
609  if (i != -1)
610  {
611  const volScalarField Kdij
612  (
613  (otherPhase/max(otherPhase, otherPhase.residualAlpha()))*Kd
614  );
615 
616  const surfaceScalarField Kdijf(fvc::interpolate(Kdij));
617 
618  invADVs[i][i] += Kdij;
619  invADVfs[i][i] += Kdijf;
620 
621  const label j = movingPhases[otherPhase.index()];
622 
623  if (j != -1)
624  {
625  invADVs[i].set(j, -Kdij);
626  invADVfs[i].set(j, -Kdijf);
627  }
628  }
629  }
630  }
631 
632  // Clear the Kds_ if they are not needed for the optional dragCorrection
633  const pimpleNoLoopControl& pimple = fluid_.pimple();
634  if (!pimple.dict().lookupOrDefault<Switch>("dragCorrection", false))
635  {
636  Kds_.clear();
637  }
638 
639  for (label i=0; i<n; i++)
640  {
641  for (label j=0; j<n; j++)
642  {
643  if (!invADVs[i].set(j))
644  {
645  invADVs[i].set
646  (
647  j,
649  (
650  "0",
651  fluid_.mesh(),
652  dimensionedScalar(As[0].dimensions(), 0)
653  )
654  );
655 
656  invADVfs[i].set
657  (
658  j,
660  (
661  "0",
662  fluid_.mesh(),
663  dimensionedScalar(As[0].dimensions(), 0)
664  )
665  );
666  }
667  }
668  }
669 
670  // Cache the phase acceleration As and Hs
671  PtrList<volScalarField> ADUDts(movingPhases.size());
672  PtrList<volVectorField> HDUDts(movingPhases.size());
673 
675  (
677  virtualMassModels_,
678  VmIter
679  )
680  {
681  const phaseInterface& interface = VmIter()->interface();
682 
683  forAllConstIter(phaseInterface, interface, iter)
684  {
685  const phaseModel& phase = iter();
686  const label i = movingPhases[phase.index()];
687 
688  if (i != -1 && !ADUDts.set(i))
689  {
690  const fvVectorMatrix DUDt(phase.DUDt());
691  ADUDts.set(i, DUDt.A());
692  HDUDts.set(i, DUDt.H());
693  }
694  }
695  }
696 
697  // Add the virtual mass contributions
699  (
701  virtualMassModels_,
702  VmIter
703  )
704  {
705  const phaseInterface& interface = VmIter()->interface();
706  const volScalarField Vm(VmIter()->K());
707 
708  forAllConstIter(phaseInterface, interface, iter)
709  {
710  const phaseModel& phase = iter();
711  const phaseModel& otherPhase = iter.otherPhase();
712 
713  const label i = movingPhases[phase.index()];
714 
715  if (i != -1)
716  {
717  const volScalarField VmPhase
718  (
719  (otherPhase/max(otherPhase, otherPhase.residualAlpha()))*Vm
720  );
721 
722  {
723  const volScalarField AVm(VmPhase*ADUDts[i]);
724 
725  invADVs[i][i] += AVm;
726  invADVfs[i][i] += fvc::interpolate(AVm);
727 
728  addField
729  (
730  i,
731  IOobject::groupName("HVm", phase.name()),
732  VmPhase*HDUDts[i],
733  HVms
734  );
735  }
736 
737  const label j = movingPhases[otherPhase.index()];
738 
739  if (j != -1)
740  {
741  const volScalarField AVm(VmPhase*ADUDts[j]);
742 
743  invADVs[i][j] -= AVm;
744  invADVfs[i][j] -= fvc::interpolate(AVm);
745 
746  addField
747  (
748  i,
749  IOobject::groupName("HVm", phase.name()),
750  -VmPhase*HDUDts[j],
751  HVms
752  );
753  }
754  }
755  }
756  }
757 
758  momentumTransferSystem::invADVs(invADVs);
759  momentumTransferSystem::invADVs(invADVfs);
760 }
761 
762 
765 (
766  const PtrList<surfaceScalarField>& Afs,
768 ) const
769 {
770  const label n = Afs.size();
771 
773 
774  forAll(invADVfs, i)
775  {
776  invADVfs.set(i, new PtrList<surfaceScalarField>(n));
777  invADVfs[i].set(i, Afs[i].clone());
778  }
779 
780  labelList movingPhases(fluid_.phases().size(), -1);
781  forAll(fluid_.movingPhases(), movingPhasei)
782  {
783  movingPhases[fluid_.movingPhases()[movingPhasei].index()] =
784  movingPhasei;
785  }
786 
788  (
790  dragModels_,
791  dragModelIter
792  )
793  {
794  const phaseInterface& interface = dragModelIter()->interface();
795  const surfaceScalarField Kdf(dragModelIter()->Kf());
796 
797  forAllConstIter(phaseInterface, interface, iter)
798  {
799  const phaseModel& phase = iter();
800  const phaseModel& otherPhase = iter.otherPhase();
801 
802  const label i = movingPhases[phase.index()];
803 
804  if (i != -1)
805  {
806  const surfaceScalarField alphaf
807  (
808  fvc::interpolate(otherPhase)
809  );
810 
811  const surfaceScalarField Kdfij
812  (
813  (alphaf/max(alphaf, otherPhase.residualAlpha()))*Kdf
814  );
815 
816  invADVfs[i][i] += Kdfij;
817 
818  const label j = movingPhases[otherPhase.index()];
819 
820  if (j != -1)
821  {
822  invADVfs[i].set(j, -Kdfij);
823  }
824  }
825  }
826  }
827 
828  for (label i=0; i<n; i++)
829  {
830  for (label j=0; j<n; j++)
831  {
832  if (!invADVfs[i].set(j))
833  {
834  invADVfs[i].set
835  (
836  j,
838  (
839  "0",
840  fluid_.mesh(),
841  dimensionedScalar(Afs[0].dimensions(), 0)
842  )
843  );
844  }
845  }
846  }
847 
848  // Cache the phase acceleration Afs and Hs
849  PtrList<volScalarField> AUgradUs(movingPhases.size());
850  PtrList<volVectorField> HUgradUs(movingPhases.size());
851 
853  (
855  virtualMassModels_,
856  VmIter
857  )
858  {
859  const phaseInterface& interface = VmIter()->interface();
860 
861  forAllConstIter(phaseInterface, interface, iter)
862  {
863  const phaseModel& phase = iter();
864  const label i = movingPhases[phase.index()];
865 
866  if (i != -1 && !AUgradUs.set(i))
867  {
868  const fvVectorMatrix UgradU(phase.UgradU());
869  AUgradUs.set(i, UgradU.A());
870  HUgradUs.set(i, UgradU.H());
871  }
872  }
873  }
874 
875  // Add the virtual mass contributions
877  (
879  virtualMassModels_,
880  VmIter
881  )
882  {
883  const phaseInterface& interface = VmIter()->interface();
884  const volScalarField Vm(VmIter()->K());
885 
886  forAllConstIter(phaseInterface, interface, iter)
887  {
888  const phaseModel& phase = iter();
889  const phaseModel& otherPhase = iter.otherPhase();
890 
891  const label i = movingPhases[phase.index()];
892 
893  if (i != -1)
894  {
895  const volScalarField VmPhase
896  (
897  (otherPhase/max(otherPhase, otherPhase.residualAlpha()))
898  *Vm
899  );
900 
901  const surfaceScalarField VmPhasef(fvc::interpolate(VmPhase));
902 
903  invADVfs[i][i] +=
904  byDt(VmPhasef) + fvc::interpolate(VmPhase*AUgradUs[i]);
905 
906  addField
907  (
908  i,
909  IOobject::groupName("HVmf", phase.name()),
910  VmPhasef
911  *byDt
912  (
914  (
915  fluid_.MRF().absolute
916  (
917  phase.phi()().oldTime()
918  ),
919  phase.U()
920  )
921  )
922  + fvc::flux(VmPhase*HUgradUs[i]),
923  HVmfs
924  );
925 
926  const label j = movingPhases[otherPhase.index()];
927 
928  if (j != -1)
929  {
930  invADVfs[i][j] -=
931  byDt(VmPhasef) + fvc::interpolate(VmPhase*AUgradUs[j]);
932 
933  addField
934  (
935  i,
936  IOobject::groupName("HVmf", phase.name()),
937  -VmPhasef
938  *byDt
939  (
941  (
942  fluid_.MRF().absolute
943  (
944  otherPhase.phi()().oldTime()
945  ),
946  otherPhase.U()
947  )
948  )
949  - fvc::flux(VmPhase*HUgradUs[j]),
950  HVmfs
951  );
952  }
953  }
954  }
955  }
956 
957  invADVs(invADVfs);
958 
959  return invADVfs;
960 }
961 
962 
965 (
966  const PtrList<volScalarField>& rAs
967 ) const
968 {
969  tmp<surfaceScalarField> alphaDByAf;
970 
971  // Add the phase pressure
972  forAll(fluid_.movingPhases(), movingPhasei)
973  {
974  const phaseModel& phase = fluid_.movingPhases()[movingPhasei];
975 
976  addTmpField
977  (
978  alphaDByAf,
979  fvc::interpolate(max(phase, scalar(0)))
980  *fvc::interpolate(rAs[phase.index()])
981  *phase.pPrimef()
982  );
983  }
984 
985  // Add the turbulent dispersion
987  (
989  turbulentDispersionModels_,
990  turbulentDispersionModelIter
991  )
992  {
993  const phaseInterface& interface =
994  turbulentDispersionModelIter()->interface();
995 
996  const surfaceScalarField alpha1f
997  (
998  fvc::interpolate(max(interface.phase1(), scalar(0)))
999  );
1000 
1001  const surfaceScalarField alpha2f
1002  (
1003  fvc::interpolate(max(interface.phase2(), scalar(0)))
1004  );
1005 
1006  addTmpField
1007  (
1008  alphaDByAf,
1009  alpha1f*alpha2f
1010  /max(alpha1f + alpha2f, interface.phase1().residualAlpha())
1012  (
1013  max
1014  (
1015  rAs[interface.phase1().index()],
1016  rAs[interface.phase2().index()]
1017  )
1018  *turbulentDispersionModelIter()->D()
1019  )
1020  );
1021  }
1022 
1023  return alphaDByAf;
1024 }
1025 
1026 
1029 {
1030  PtrList<surfaceScalarField> ddtCorrs(fluid_.phases().size());
1031 
1032  // Add correction
1033  forAll(fluid_.movingPhases(), movingPhasei)
1034  {
1035  const phaseModel& phase = fluid_.movingPhases()[movingPhasei];
1036 
1037  addField
1038  (
1039  phase,
1040  "ddtCorr",
1041  fvc::ddtCorr
1042  (
1043  phase,
1044  phase.rho(),
1045  phase.U()(),
1046  phase.phi()(),
1047  phase.Uf()
1048  ),
1049  ddtCorrs
1050  );
1051  }
1052 
1053 
1054  const pimpleNoLoopControl& pimple = fluid_.pimple();
1055  const Switch VmDdtCorr
1056  (
1057  pimple.dict().lookupOrDefault<Switch>("VmDdtCorrection", false)
1058  );
1059 
1060  // Optionally ddd virtual mass correction
1061  if (VmDdtCorr)
1062  {
1063  PtrList<volScalarField> VmDdtCoeffs(fluid_.phases().size());
1064  PtrList<surfaceScalarField> VmDdtCorrs(fluid_.phases().size());
1065 
1066  forAll(fluid_.movingPhases(), movingPhasei)
1067  {
1068  const phaseModel& phase = fluid_.movingPhases()[movingPhasei];
1069  const label phasei = phase.index();
1070 
1071  VmDdtCorrs.set
1072  (
1073  phasei,
1074  fvc::ddtCorr
1075  (
1076  phase.U()(),
1077  phase.phi()(),
1078  phase.Uf()
1079  )
1080  );
1081  }
1082 
1084  (
1086  virtualMassModels_,
1087  VmIter
1088  )
1089  {
1090  const phaseInterface& interface = VmIter()->interface();
1091  const volScalarField Vm(VmIter()->K());
1092 
1093  forAllConstIter(phaseInterface, interface, iter)
1094  {
1095  const phaseModel& phase = iter();
1096  const label phasei = phase.index();
1097  const phaseModel& otherPhase = iter.otherPhase();
1098  const label otherPhasei = otherPhase.index();
1099 
1100  const volScalarField VmPhase
1101  (
1102  (otherPhase/max(otherPhase, otherPhase.residualAlpha()))
1103  *Vm
1104  );
1105 
1106  addField
1107  (
1108  phase,
1109  "ddtCorr",
1110  fvc::interpolate(VmPhase)
1111  *(VmDdtCorrs[phasei] - VmDdtCorrs[otherPhasei]),
1112  ddtCorrs
1113  );
1114  }
1115  }
1116  }
1117 
1118  return ddtCorrs;
1119 }
1120 
1121 
1123 (
1124  PtrList<volVectorField>& dragCorrs,
1125  PtrList<surfaceScalarField>& dragCorrfs
1126 ) const
1127 {
1128  labelList movingPhases(fluid_.phases().size(), -1);
1129  forAll(fluid_.movingPhases(), movingPhasei)
1130  {
1131  movingPhases[fluid_.movingPhases()[movingPhasei].index()] =
1132  movingPhasei;
1133  }
1134 
1135  PtrList<volVectorField> Uphis(fluid_.movingPhases().size());
1136  forAll(fluid_.movingPhases(), movingPhasei)
1137  {
1138  Uphis.set
1139  (
1140  movingPhasei,
1141  fvc::reconstruct(fluid_.movingPhases()[movingPhasei].phi())
1142  );
1143  }
1144 
1145  forAllConstIter(KdTable, Kds_, KdIter)
1146  {
1147  const volScalarField& K(*KdIter());
1148  const phaseInterface interface(fluid_, KdIter.key());
1149 
1150  forAllConstIter(phaseInterface, interface, iter)
1151  {
1152  const phaseModel& phase = iter();
1153  const phaseModel& otherPhase = iter.otherPhase();
1154 
1155  const label i = movingPhases[phase.index()];
1156  const label j = movingPhases[otherPhase.index()];
1157 
1158  if (i != -1)
1159  {
1160  const volScalarField K1
1161  (
1162  (otherPhase/max(otherPhase, otherPhase.residualAlpha()))*K
1163  );
1164 
1165  addField
1166  (
1167  i,
1168  IOobject::groupName("dragCorr", phase.name()),
1169  K1*(j == -1 ? -Uphis[i] : (Uphis[j] - Uphis[i])),
1170  dragCorrs
1171  );
1172 
1173  addField
1174  (
1175  i,
1176  IOobject::groupName("dragCorrf", phase.name()),
1178  *(j == -1 ? -phase.phi() : (otherPhase.phi() - phase.phi())),
1179  dragCorrfs
1180  );
1181  }
1182  }
1183  }
1184 }
1185 
1186 
1188 {
1189  if (regIOobject::read())
1190  {
1191  bool readOK = true;
1192 
1193  // models ...
1194 
1195  return readOK;
1196  }
1197  else
1198  {
1199  return false;
1200  }
1201 }
1202 
1203 
1204 // ************************************************************************* //
#define K1
Definition: SHA1.C:167
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
A HashTable specialisation for hashing pointers.
Definition: HashPtrTable.H:67
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
@ MUST_READ_IF_MODIFIED
Definition: IOobject.H:119
static word groupName(Name name, const word &group)
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:125
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:62
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: UPtrList.H:66
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
T lookupOrDefault(const word &, const T &, const bool writeDefault=writeOptionalEntries > 0) const
Find and return a T, if not found return the given default.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
tmp< volScalarField > A() const
Return the central coefficient.
Definition: fvMatrix.C:808
tmp< VolField< Type > > H() const
Return the H operation source.
Definition: fvMatrix.C:868
Class which provides interfacial momentum transfer between a number of phases. Drag,...
static const word propertiesName
Default name of the phase properties dictionary.
void dragCorrs(PtrList< volVectorField > &dragCorrs, PtrList< surfaceScalarField > &dragCorrf) const
Set the cell and faces drag correction fields.
PtrList< surfaceScalarField > ddtCorrs() const
Return the flux corrections for the cell-based algorithm. These.
momentumTransferSystem(const phaseSystem &)
Construct from a phase system.
PtrList< surfaceScalarField > Ffs() const
As Fs, but for the face-based algorithm.
PtrList< surfaceScalarField > Fs() const
Return the explicit force fluxes for the cell-based algorithm, that.
virtual ~momentumTransferSystem()
Destructor.
PtrList< PtrList< surfaceScalarField > > invADVfs(const PtrList< surfaceScalarField > &Afs, PtrList< surfaceScalarField > &HVmfs) const
Return the inverse of the central + drag + virtual mass.
virtual tmp< surfaceScalarField > alphaDByAf(const PtrList< volScalarField > &rAs) const
Return the phase diffusivity.
virtual bool read()
Read base phaseProperties dictionary.
Class to represent an interface between phases. Derivations can further specify the configuration of ...
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:169
virtual tmp< fvVectorMatrix > DUDt() const =0
Return the substantive acceleration matrix.
virtual tmp< surfaceScalarField > phi() const =0
Return the volumetric flux.
label index() const
Return the index of the phase.
Definition: phaseModel.C:157
virtual tmp< volVectorField > U() const =0
Return the velocity.
virtual tmp< fvVectorMatrix > UgradU() const =0
Return the velocity transport matrix.
const word & name() const
Return the name of this phase.
Definition: phaseModel.C:145
virtual tmp< surfaceScalarField > pPrimef() const =0
Return the face-phase-pressure'.
virtual const volScalarField & rho() const =0
Return the density field.
Class to represent a system of phases.
Definition: phaseSystem.H:74
Pimple no-loop control class. Implements various option flags, but leaves loop controls to the deriva...
virtual bool read()
Read object.
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:183
bool isTmp() const
Return true if this is really a temporary object.
Definition: tmpI.H:169
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
A class for handling words, derived from string.
Definition: word.H:62
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
#define WarningInFunction
Report a warning using Foam::Warning.
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< 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
static const coefficient D("D", dimTemperature, 257.14)
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
List< word > wordList
A List of words.
Definition: fileName.H:54
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
const HashTable< dimensionSet > & dimensions()
Get the table of dimension sets.
Definition: dimensionSets.C:96
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
void addField(const label phasei, const word &name, tmp< GeoField > field, PtrList< GeoField > &fieldList)
Definition: phaseSystemI.H:265
void LUBacksubstitute(const scalarSquareMatrix &luMmatrix, const labelList &pivotIndices, List< Type > &source)
LU back-substitution with given source, returning the solution.
const dimensionSet dimless
void LUDecompose(scalarSquareMatrix &matrix, labelList &pivotIndices)
LU decompose the matrix with pivoting.
tmp< volScalarField > byDt(const volScalarField &vf)
Definition: phaseSystem.C:927
Foam::HashPtrTable< ModelType, Foam::phaseInterfaceKey, Foam::phaseInterfaceKey::hash > generateBlendedInterfacialModels(const phaseSystem &fluid, const dictionary &dict, const wordHashSet &ignoreKeys, const bool ignoreNonModelPhaseInterfaceTypes, const Args &... args)
defineTypeNameAndDebug(combustionModel, 0)
T clone(const T &t)
Definition: List.H:55
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)