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 "scalarMatrices.H"
35 
36 #include "fvmDdt.H"
37 #include "fvmDiv.H"
38 #include "fvmSup.H"
39 
40 #include "fvcDdt.H"
41 #include "fvcDiv.H"
42 #include "fvcFlux.H"
43 #include "fvcSnGrad.H"
44 #include "fvcMeshPhi.H"
45 #include "fvcReconstruct.H"
46 
47 #include "pimpleNoLoopControl.H"
48 
49 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
50 
51 template<class BasePhaseSystem>
53 (
54  const phaseSystem::dmdtfTable& dmdtfs,
56 )
57 {
58  forAllConstIter(phaseSystem::dmdtfTable, dmdtfs, dmdtfIter)
59  {
60  const phaseInterface interface(*this, dmdtfIter.key());
61 
62  const volScalarField& dmdtf = *dmdtfIter();
63  const volScalarField dmdtf21(posPart(dmdtf));
64  const volScalarField dmdtf12(negPart(dmdtf));
65 
66  phaseModel& phase1 = this->phases()[interface.phase1().name()];
67  phaseModel& phase2 = this->phases()[interface.phase2().name()];
68 
69  if (!phase1.stationary())
70  {
71  *eqns[phase1.name()] +=
72  dmdtf21*phase2.U() + fvm::Sp(dmdtf12, phase1.URef());
73  }
74 
75  if (!phase2.stationary())
76  {
77  *eqns[phase2.name()] -=
78  dmdtf12*phase1.U() + fvm::Sp(dmdtf21, phase2.URef());
79  }
80  }
81 }
82 
83 
84 template<class BasePhaseSystem>
86 (
88  const tmp<surfaceScalarField>& field
89 ) const
90 {
91  if (result.valid())
92  {
93  result.ref() += field;
94  }
95  else
96  {
97  if (field.isTmp())
98  {
99  result = field;
100  }
101  else
102  {
103  result = field().clone();
104  }
105  }
106 }
107 
108 
109 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
110 
111 template<class BasePhaseSystem>
114 (
115  const fvMesh& mesh
116 )
117 :
118  BasePhaseSystem(mesh)
119 {
120  this->generateInterfacialModels(dragModels_);
121  this->generateInterfacialModels(virtualMassModels_);
122  this->generateInterfacialModels(liftModels_);
123  this->generateInterfacialModels(wallLubricationModels_);
124  this->generateInterfacialModels(turbulentDispersionModels_);
125 }
126 
127 
128 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
129 
130 template<class BasePhaseSystem>
133 {}
134 
135 
136 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
137 
138 template<class BasePhaseSystem>
141 {
142  // Create a momentum transfer matrix for each phase
144  (
146  );
147 
148  phaseSystem::momentumTransferTable& eqns = eqnsPtr();
149 
150  forAll(this->movingPhases(), movingPhasei)
151  {
152  const phaseModel& phase = this->movingPhases()[movingPhasei];
153 
154  eqns.insert
155  (
156  phase.name(),
158  );
159  }
160 
161  return eqnsPtr;
162 }
163 
164 
165 template<class BasePhaseSystem>
168 {
169  // Create a momentum transfer matrix for each phase
171  (
173  );
174 
175  phaseSystem::momentumTransferTable& eqns = eqnsPtr();
176 
177  forAll(this->movingPhases(), movingPhasei)
178  {
179  const phaseModel& phase = this->movingPhases()[movingPhasei];
180 
181  eqns.insert
182  (
183  phase.name(),
185  );
186  }
187 
188  return eqnsPtr;
189 }
190 
191 
192 template<class BasePhaseSystem>
195 {
196  PtrList<surfaceScalarField> Fs(this->phaseModels_.size());
197 
198  // Add the lift force
200  (
202  liftModels_,
203  liftModelIter
204  )
205  {
206  const phaseInterface& interface = liftModelIter()->interface();
207 
208  const volVectorField F(liftModelIter()->F());
209 
210  addField
211  (
212  interface.phase1(),
213  "F",
214  fvc::flux(F),
215  Fs
216  );
217  addField
218  (
219  interface.phase2(),
220  "F",
221  -fvc::flux(F),
222  Fs
223  );
224  }
225 
226  // Add the wall lubrication force
228  (
230  wallLubricationModels_,
231  wallLubricationModelIter
232  )
233  {
234  const phaseInterface& interface =
235  wallLubricationModelIter()->interface();
236 
237  const volVectorField F(wallLubricationModelIter()->F());
238 
239  addField
240  (
241  interface.phase1(),
242  "F",
243  fvc::flux(F),
244  Fs
245  );
246  addField
247  (
248  interface.phase2(),
249  "F",
250  -fvc::flux(F),
251  Fs
252  );
253  }
254 
255  // Add the phase pressure
256  forAll(this->movingPhases(), movingPhasei)
257  {
258  const phaseModel& phase = this->movingPhases()[movingPhasei];
259 
260  addField
261  (
262  phase,
263  "F",
264  phase.pPrimef()
265  *fvc::snGrad(phase)*this->mesh_.magSf(),
266  Fs
267  );
268  }
269 
270  // Add the turbulent dispersion force
272  (
274  turbulentDispersionModels_,
275  turbulentDispersionModelIter
276  )
277  {
278  const phaseInterface& interface =
279  turbulentDispersionModelIter()->interface();
280 
281  const surfaceScalarField Df
282  (
283  fvc::interpolate(turbulentDispersionModelIter()->D())
284  );
285 
286  const volScalarField alpha12(interface.phase1() + interface.phase2());
287  const surfaceScalarField snGradAlpha1By12
288  (
290  (
291  interface.phase1()
292  /max(alpha12, interface.phase1().residualAlpha())
293  )*this->mesh_.magSf()
294  );
295  const surfaceScalarField snGradAlpha2By12
296  (
298  (
299  interface.phase2()
300  /max(alpha12, interface.phase2().residualAlpha())
301  )*this->mesh_.magSf()
302  );
303 
304  addField(interface.phase1(), "F", Df*snGradAlpha1By12, Fs);
305  addField(interface.phase2(), "F", Df*snGradAlpha2By12, Fs);
306  }
307 
308  return Fs;
309 }
310 
311 
312 template<class BasePhaseSystem>
315 {
316  PtrList<surfaceScalarField> Ffs(this->phaseModels_.size());
317 
318  // Add the lift force
320  (
322  liftModels_,
323  liftModelIter
324  )
325  {
326  const phaseInterface& interface = liftModelIter()->interface();
327 
328  const surfaceScalarField Ff(liftModelIter()->Ff());
329 
330  addField
331  (
332  interface.phase1(),
333  "Ff",
334  Ff,
335  Ffs
336  );
337  addField
338  (
339  interface.phase2(),
340  "Ff",
341  -Ff,
342  Ffs
343  );
344  }
345 
346  // Add the wall lubrication force
348  (
350  wallLubricationModels_,
351  wallLubricationModelIter
352  )
353  {
354  const phaseInterface& interface =
355  wallLubricationModelIter()->interface();
356 
357  const surfaceScalarField Ff(wallLubricationModelIter()->Ff());
358 
359  addField
360  (
361  interface.phase1(),
362  "Ff",
363  Ff,
364  Ffs
365  );
366  addField
367  (
368  interface.phase2(),
369  "Ff",
370  -Ff,
371  Ffs
372  );
373  }
374 
375  // Add the phase pressure
376  forAll(this->movingPhases(), movingPhasei)
377  {
378  const phaseModel& phase = this->movingPhases()[movingPhasei];
379 
380  addField
381  (
382  phase,
383  "Ff",
384  phase.pPrimef()
385  *fvc::snGrad(phase)*this->mesh_.magSf(),
386  Ffs
387  );
388  }
389 
390  // Add the turbulent dispersion force
392  (
394  turbulentDispersionModels_,
395  turbulentDispersionModelIter
396  )
397  {
398  const phaseInterface& interface =
399  turbulentDispersionModelIter()->interface();
400 
401  const surfaceScalarField Df
402  (
403  fvc::interpolate(turbulentDispersionModelIter()->D())
404  );
405 
406  const volScalarField alpha12(interface.phase1() + interface.phase2());
407  const surfaceScalarField snGradAlpha1By12
408  (
410  (
411  interface.phase1()
412  /max(alpha12, interface.phase1().residualAlpha())
413  )*this->mesh_.magSf()
414  );
415  const surfaceScalarField snGradAlpha2By12
416  (
418  (
419  interface.phase2()
420  /max(alpha12, interface.phase2().residualAlpha())
421  )*this->mesh_.magSf()
422  );
423 
424  addField(interface.phase1(), "F", Df*snGradAlpha1By12, Ffs);
425  addField(interface.phase2(), "F", Df*snGradAlpha2By12, Ffs);
426  }
427 
428  return Ffs;
429 }
430 
431 
432 template<class BasePhaseSystem>
434 (
436 ) const
437 {
438  const label n = ADVs.size();
439 
440  scalarSquareMatrix AD(n);
441  scalarField source(n);
442  labelList pivotIndices(n);
443 
444  forAll(ADVs[0][0], ci)
445  {
446  for (label i=0; i<n; i++)
447  {
448  for (label j=0; j<n; j++)
449  {
450  AD(i, j) = ADVs[i][j][ci];
451  }
452  }
453 
454  // Calculate the inverse of AD using LD decomposition
455  // and back-substitution
456  LUDecompose(AD, pivotIndices);
457 
458  for (label j=0; j<n; j++)
459  {
460  source = Zero;
461  source[j] = 1;
462 
463  LUBacksubstitute(AD, pivotIndices, source);
464 
465  for (label i=0; i<n; i++)
466  {
467  ADVs[i][j][ci] = source[i];
468  }
469  }
470  }
471 }
472 
473 
474 template<class BasePhaseSystem>
475 template<template<class> class PatchField, class GeoMesh>
477 (
479 ) const
480 {
481  const label n = ADVs.size();
482 
484 
485  for (label i=0; i<n; i++)
486  {
487  ADps[i].setSize(n);
488 
489  for (label j=0; j<n; j++)
490  {
491  ADps[i].set(j, &ADVs[i][j]);
492 
493  ADVs[i][j].dimensions().reset(dimless/ADVs[i][j].dimensions());
494  }
495  }
496 
497  invADVs(ADps);
498 
499  forAll(ADVs[0][0].boundaryField(), patchi)
500  {
501  for (label i=0; i<n; i++)
502  {
503  for (label j=0; j<n; j++)
504  {
505  ADps[i].set(j, &ADVs[i][j].boundaryFieldRef()[patchi]);
506  }
507  }
508 
509  invADVs(ADps);
510  }
511 }
512 
513 
514 template<class BasePhaseSystem>
516 (
517  const PtrList<volScalarField>& As,
521 ) const
522 {
523  const label n = As.size();
524 
525  invADVs.setSize(n);
526  invADVfs.setSize(n);
527 
528  forAll(invADVs, i)
529  {
530  invADVs.set(i, new PtrList<volScalarField>(n));
531  invADVs[i].set(i, As[i].clone());
532 
533  invADVfs.set(i, new PtrList<surfaceScalarField>(n));
534  invADVfs[i].set(i, fvc::interpolate(As[i]));
535  }
536 
537  labelList movingPhases(this->phases().size(), -1);
538 
539  forAll(this->movingPhases(), movingPhasei)
540  {
541  movingPhases[this->movingPhases()[movingPhasei].index()] = movingPhasei;
542  }
543 
544  Kds_.clear();
545 
546  // Update the drag coefficients
548  (
550  dragModels_,
551  dragModelIter
552  )
553  {
554  const phaseInterface& interface = dragModelIter()->interface();
555 
556  tmp<volScalarField> tKd(dragModelIter()->K());
557  volScalarField& Kd = tKd.ref();
558  Kds_.insert(dragModelIter.key(), tKd.ptr());
559 
560  // Zero-gradient the drag coefficient to boundaries with fixed velocity
562  {
563  if
564  (
565  (
566  !interface.phase1().stationary()
567  && interface.phase1().U()()
568  .boundaryField()[patchi].fixesValue()
569  )
570  && (
571  !interface.phase2().stationary()
572  && interface.phase2().U()()
573  .boundaryField()[patchi].fixesValue()
574  )
575  )
576  {
577  Kd.boundaryFieldRef()[patchi] =
578  Kd.boundaryField()[patchi].patchInternalField();
579  }
580  }
581 
582  forAllConstIter(phaseInterface, interface, iter)
583  {
584  const phaseModel& phase = iter();
585  const phaseModel& otherPhase = iter.otherPhase();
586 
587  const label i = movingPhases[phase.index()];
588 
589  if (i != -1)
590  {
591  const volScalarField Kdij
592  (
593  (otherPhase/max(otherPhase, otherPhase.residualAlpha()))*Kd
594  );
595 
596  const surfaceScalarField Kdijf(fvc::interpolate(Kdij));
597 
598  invADVs[i][i] += Kdij;
599  invADVfs[i][i] += Kdijf;
600 
601  const label j = movingPhases[otherPhase.index()];
602 
603  if (j != -1)
604  {
605  invADVs[i].set(j, -Kdij);
606  invADVfs[i].set(j, -Kdijf);
607  }
608  }
609  }
610  }
611 
612  // Clear the Kds_ if they are not needed for the optional dragCorrection
613  const pimpleNoLoopControl& pimple = this->pimple();
614  if (!pimple.dict().lookupOrDefault<Switch>("dragCorrection", false))
615  {
616  Kds_.clear();
617  }
618 
619  for (label i=0; i<n; i++)
620  {
621  for (label j=0; j<n; j++)
622  {
623  if (!invADVs[i].set(j))
624  {
625  invADVs[i].set
626  (
627  j,
629  (
630  "0",
631  this->mesh(),
632  dimensionedScalar(As[0].dimensions(), 0)
633  )
634  );
635 
636  invADVfs[i].set
637  (
638  j,
640  (
641  "0",
642  this->mesh(),
643  dimensionedScalar(As[0].dimensions(), 0)
644  )
645  );
646  }
647  }
648  }
649 
650  // Cache the phase acceleration As and Hs
651  PtrList<volScalarField> ADUDts(movingPhases.size());
652  PtrList<volVectorField> HDUDts(movingPhases.size());
653 
655  (
657  virtualMassModels_,
658  VmIter
659  )
660  {
661  const phaseInterface& interface = VmIter()->interface();
662 
663  forAllConstIter(phaseInterface, interface, iter)
664  {
665  const phaseModel& phase = iter();
666  const label i = movingPhases[phase.index()];
667 
668  if (i != -1 && !ADUDts.set(i))
669  {
670  const fvVectorMatrix DUDt(phase.DUDt());
671  ADUDts.set(i, DUDt.A());
672  HDUDts.set(i, DUDt.H());
673  }
674  }
675  }
676 
677  // Add the virtual mass contributions
679  (
681  virtualMassModels_,
682  VmIter
683  )
684  {
685  const phaseInterface& interface = VmIter()->interface();
686  const volScalarField Vm(VmIter()->K());
687 
688  forAllConstIter(phaseInterface, interface, iter)
689  {
690  const phaseModel& phase = iter();
691  const phaseModel& otherPhase = iter.otherPhase();
692 
693  const label i = movingPhases[phase.index()];
694 
695  if (i != -1)
696  {
697  const volScalarField VmPhase
698  (
699  (otherPhase/max(otherPhase, otherPhase.residualAlpha()))*Vm
700  );
701 
702  {
703  const volScalarField AVm(VmPhase*ADUDts[i]);
704 
705  invADVs[i][i] += AVm;
706  invADVfs[i][i] += fvc::interpolate(AVm);
707 
708  addField
709  (
710  i,
711  "HVm",
712  VmPhase*HDUDts[i],
713  HVms
714  );
715  }
716 
717  const label j = movingPhases[otherPhase.index()];
718 
719  if (j != -1)
720  {
721  const volScalarField AVm(VmPhase*ADUDts[j]);
722 
723  invADVs[i][j] -= AVm;
724  invADVfs[i][j] -= fvc::interpolate(AVm);
725 
726  addField
727  (
728  i,
729  "HVm",
730  -VmPhase*HDUDts[j],
731  HVms
732  );
733  }
734  }
735  }
736  }
737 
740 }
741 
742 
743 template<class BasePhaseSystem>
746 (
747  const PtrList<surfaceScalarField>& Afs,
749 ) const
750 {
751  const label n = Afs.size();
752 
754 
755  forAll(invADVfs, i)
756  {
757  invADVfs.set(i, new PtrList<surfaceScalarField>(n));
758  invADVfs[i].set(i, Afs[i].clone());
759  }
760 
761  labelList movingPhases(this->phases().size(), -1);
762 
763  forAll(this->movingPhases(), movingPhasei)
764  {
765  movingPhases[this->movingPhases()[movingPhasei].index()] = movingPhasei;
766  }
767 
769  (
771  dragModels_,
772  dragModelIter
773  )
774  {
775  const phaseInterface& interface = dragModelIter()->interface();
776  const surfaceScalarField Kdf(dragModelIter()->Kf());
777 
778  forAllConstIter(phaseInterface, interface, iter)
779  {
780  const phaseModel& phase = iter();
781  const phaseModel& otherPhase = iter.otherPhase();
782 
783  const label i = movingPhases[phase.index()];
784 
785  if (i != -1)
786  {
787  const surfaceScalarField alphaf
788  (
789  fvc::interpolate(otherPhase)
790  );
791 
792  const surfaceScalarField Kdfij
793  (
794  (alphaf/max(alphaf, otherPhase.residualAlpha()))*Kdf
795  );
796 
797  invADVfs[i][i] += Kdfij;
798 
799  const label j = movingPhases[otherPhase.index()];
800 
801  if (j != -1)
802  {
803  invADVfs[i].set(j, -Kdfij);
804  }
805  }
806  }
807  }
808 
809  for (label i=0; i<n; i++)
810  {
811  for (label j=0; j<n; j++)
812  {
813  if (!invADVfs[i].set(j))
814  {
815  invADVfs[i].set
816  (
817  j,
819  (
820  "0",
821  this->mesh(),
822  dimensionedScalar(Afs[0].dimensions(), 0)
823  )
824  );
825  }
826  }
827  }
828 
829  // Cache the phase acceleration Afs and Hs
830  PtrList<volScalarField> AUgradUs(movingPhases.size());
831  PtrList<volVectorField> HUgradUs(movingPhases.size());
832 
834  (
836  virtualMassModels_,
837  VmIter
838  )
839  {
840  const phaseInterface& interface = VmIter()->interface();
841 
842  forAllConstIter(phaseInterface, interface, iter)
843  {
844  const phaseModel& phase = iter();
845  const label i = movingPhases[phase.index()];
846 
847  if (i != -1 && !AUgradUs.set(i))
848  {
849  const fvVectorMatrix UgradU(phase.UgradU());
850  AUgradUs.set(i, UgradU.A());
851  HUgradUs.set(i, UgradU.H());
852  }
853  }
854  }
855 
856  // Add the virtual mass contributions
858  (
860  virtualMassModels_,
861  VmIter
862  )
863  {
864  const phaseInterface& interface = VmIter()->interface();
865  const volScalarField Vm(VmIter()->K());
866 
867  forAllConstIter(phaseInterface, interface, iter)
868  {
869  const phaseModel& phase = iter();
870  const phaseModel& otherPhase = iter.otherPhase();
871 
872  const label i = movingPhases[phase.index()];
873 
874  if (i != -1)
875  {
876  const volScalarField VmPhase
877  (
878  (otherPhase/max(otherPhase, otherPhase.residualAlpha()))
879  *Vm
880  );
881 
882  const surfaceScalarField VmPhasef(fvc::interpolate(VmPhase));
883 
884  invADVfs[i][i] +=
885  byDt(VmPhasef) + fvc::interpolate(VmPhase*AUgradUs[i]);
886 
887  addField
888  (
889  i,
890  "HVmf",
891  VmPhasef
892  *byDt
893  (
895  (
896  this->MRF().absolute(phase.phi()().oldTime()),
897  phase.U()
898  )
899  )
900  + fvc::flux(VmPhase*HUgradUs[i]),
901  HVmfs
902  );
903 
904  const label j = movingPhases[otherPhase.index()];
905 
906  if (j != -1)
907  {
908  invADVfs[i][j] -=
909  byDt(VmPhasef) + fvc::interpolate(VmPhase*AUgradUs[j]);
910 
911  addField
912  (
913  i,
914  "HVmf",
915  -VmPhasef
916  *byDt
917  (
919  (
920  this->MRF().absolute
921  (
922  otherPhase.phi()().oldTime()
923  ),
924  otherPhase.U()
925  )
926  )
927  - fvc::flux(VmPhase*HUgradUs[j]),
928  HVmfs
929  );
930  }
931  }
932  }
933  }
934 
935  invADVs(invADVfs);
936 
937  return invADVfs;
938 }
939 
940 
941 template<class BasePhaseSystem>
943 implicitPhasePressure(const phaseModel& phase) const
944 {
945  return
946  this->mesh_.solution().solverDict(phase.volScalarField::name()).
947  template lookupOrDefault<Switch>
948  (
949  "implicitPhasePressure",
950  false
951  );
952 }
953 
954 
955 template<class BasePhaseSystem>
957 implicitPhasePressure() const
958 {
959  bool implicitPressure = false;
960 
961  forAll(this->phaseModels_, phasei)
962  {
963  const phaseModel& phase = this->phaseModels_[phasei];
964 
965  implicitPressure = implicitPressure || implicitPhasePressure(phase);
966  }
967 
968  return implicitPressure;
969 }
970 
971 
972 template<class BasePhaseSystem>
975 (
976  const PtrList<volScalarField>& rAs
977 ) const
978 {
979  tmp<surfaceScalarField> alphaDByAf;
980 
981  // Add the phase pressure
982  forAll(this->movingPhases(), movingPhasei)
983  {
984  const phaseModel& phase = this->movingPhases()[movingPhasei];
985 
986  addTmpField
987  (
988  alphaDByAf,
989  fvc::interpolate(max(phase, scalar(0)))
990  *fvc::interpolate(rAs[phase.index()])
991  *phase.pPrimef()
992  );
993  }
994 
995  // Add the turbulent dispersion
997  (
999  turbulentDispersionModels_,
1000  turbulentDispersionModelIter
1001  )
1002  {
1003  const phaseInterface& interface =
1004  turbulentDispersionModelIter()->interface();
1005 
1006  const surfaceScalarField alpha1f
1007  (
1008  fvc::interpolate(max(interface.phase1(), scalar(0)))
1009  );
1010 
1011  const surfaceScalarField alpha2f
1012  (
1013  fvc::interpolate(max(interface.phase2(), scalar(0)))
1014  );
1015 
1016  addTmpField
1017  (
1018  alphaDByAf,
1019  alpha1f*alpha2f
1020  /max(alpha1f + alpha2f, interface.phase1().residualAlpha())
1022  (
1023  max
1024  (
1025  rAs[interface.phase1().index()],
1026  rAs[interface.phase2().index()]
1027  )
1028  *turbulentDispersionModelIter()->D()
1029  )
1030  );
1031  }
1032 
1033  return alphaDByAf;
1034 }
1035 
1036 
1037 template<class BasePhaseSystem>
1040 {
1041  PtrList<surfaceScalarField> ddtCorrs(this->phaseModels_.size());
1042 
1043  // Add correction
1044  forAll(this->movingPhases(), movingPhasei)
1045  {
1046  const phaseModel& phase = this->movingPhases()[movingPhasei];
1047 
1048  addField
1049  (
1050  phase,
1051  "ddtCorr",
1052  fvc::ddtCorr
1053  (
1054  phase,
1055  phase.rho(),
1056  phase.U()(),
1057  phase.phi()(),
1058  phase.Uf()
1059  ),
1060  ddtCorrs
1061  );
1062  }
1063 
1064 
1065  const pimpleNoLoopControl& pimple = this->pimple();
1066  const Switch VmDdtCorr
1067  (
1068  pimple.dict().lookupOrDefault<Switch>("VmDdtCorrection", false)
1069  );
1070 
1071  // Optionally ddd virtual mass correction
1072  if (VmDdtCorr)
1073  {
1074  PtrList<volScalarField> VmDdtCoeffs(this->phaseModels_.size());
1075  PtrList<surfaceScalarField> VmDdtCorrs(this->phaseModels_.size());
1076 
1077  forAll(this->movingPhases(), movingPhasei)
1078  {
1079  const phaseModel& phase = this->movingPhases()[movingPhasei];
1080  const label phasei = phase.index();
1081 
1082  VmDdtCorrs.set
1083  (
1084  phasei,
1085  fvc::ddtCorr
1086  (
1087  phase.U()(),
1088  phase.phi()(),
1089  phase.Uf()
1090  )
1091  );
1092  }
1093 
1095  (
1097  virtualMassModels_,
1098  VmIter
1099  )
1100  {
1101  const phaseInterface& interface = VmIter()->interface();
1102  const volScalarField Vm(VmIter()->K());
1103 
1104  forAllConstIter(phaseInterface, interface, iter)
1105  {
1106  const phaseModel& phase = iter();
1107  const label phasei = phase.index();
1108  const phaseModel& otherPhase = iter.otherPhase();
1109  const label otherPhasei = otherPhase.index();
1110 
1111  const volScalarField VmPhase
1112  (
1113  (otherPhase/max(otherPhase, otherPhase.residualAlpha()))
1114  *Vm
1115  );
1116 
1117  addField
1118  (
1119  phase,
1120  "ddtCorr",
1121  fvc::interpolate(VmPhase)
1122  *(VmDdtCorrs[phasei] - VmDdtCorrs[otherPhasei]),
1123  ddtCorrs
1124  );
1125  }
1126  }
1127  }
1128 
1129  return ddtCorrs;
1130 }
1131 
1132 
1133 template<class BasePhaseSystem>
1135 (
1136  PtrList<volVectorField>& dragCorrs,
1137  PtrList<surfaceScalarField>& dragCorrfs
1138 ) const
1139 {
1140  labelList movingPhases(this->phases().size(), -1);
1141  PtrList<volVectorField> Uphis(this->movingPhases().size());
1142 
1143  forAll(this->movingPhases(), movingPhasei)
1144  {
1145  movingPhases[this->movingPhases()[movingPhasei].index()] = movingPhasei;
1146 
1147  Uphis.set
1148  (
1149  movingPhasei,
1150  fvc::reconstruct(this->movingPhases()[movingPhasei].phi())
1151  );
1152  }
1153 
1154  forAllConstIter(KdTable, Kds_, KdIter)
1155  {
1156  const volScalarField& K(*KdIter());
1157  const phaseInterface interface(*this, KdIter.key());
1158 
1159  forAllConstIter(phaseInterface, interface, iter)
1160  {
1161  const phaseModel& phase = iter();
1162  const phaseModel& otherPhase = iter.otherPhase();
1163 
1164  const label i = movingPhases[phase.index()];
1165  const label j = movingPhases[otherPhase.index()];
1166 
1167  if (i != -1)
1168  {
1169  const volScalarField K1
1170  (
1171  (otherPhase/max(otherPhase, otherPhase.residualAlpha()))*K
1172  );
1173 
1174  addField
1175  (
1176  i,
1177  "dragCorr",
1178  K1*(j == -1 ? -Uphis[i] : (Uphis[j] - Uphis[i])),
1179  dragCorrs
1180  );
1181 
1182  addField
1183  (
1184  i,
1185  "dragCorrf",
1187  *(j == -1 ? -phase.phi() : (otherPhase.phi() - phase.phi())),
1188  dragCorrfs
1189  );
1190  }
1191  }
1192  }
1193 }
1194 
1195 
1196 template<class BasePhaseSystem>
1198 {
1199  if (BasePhaseSystem::read())
1200  {
1201  bool readOK = true;
1202 
1203  // Read models ...
1204 
1205  return readOK;
1206  }
1207  else
1208  {
1209  return false;
1210  }
1211 }
1212 
1213 
1214 // ************************************************************************* //
static const Foam::dimensionedScalar D("D", Foam::dimTemperature, 257.14)
#define K1
Definition: SHA1.C:167
label n
#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 mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:47
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, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
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
void setSize(const label)
Reset size of List.
Definition: List.C:281
void invADVs(List< UPtrList< scalarField >> &ADVs) const
Invert the ADVs matrix inplace.
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.
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.
MomentumTransferPhaseSystem(const fvMesh &)
Construct from fvMesh.
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransferf()
As momentumTransfer, but for the face-based algorithm.
virtual 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.
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: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
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, 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
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
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:169
virtual tmp< fvVectorMatrix > DUDt() const =0
Return the substantive acceleration matrix.
virtual bool stationary() const =0
Return whether the phase is stationary.
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 volVectorField & URef()=0
Access 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.
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
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
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
const dimensionedScalar F
Faraday constant: default SI units: [C/kmol].
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< 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 zero Zero
Definition: zero.H:97
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
void addField(const label phasei, const word &name, tmp< GeoField > field, PtrList< GeoField > &fieldList)
Definition: phaseSystemI.H:263
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:855
const dimensionSet dimTime
dimensionedScalar negPart(const dimensionedScalar &ds)
T clone(const T &t)
Definition: List.H:55
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimMass
const dimensionSet dimVelocity
dimensionedScalar posPart(const dimensionedScalar &ds)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.