phaseSystemSolve.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) 2013-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 "phaseSystem.H"
27 #include "momentumTransferSystem.H"
28 
29 #include "subCycle.H"
30 
31 #include "fvcDdt.H"
32 #include "fvcDiv.H"
33 #include "fvcSnGrad.H"
34 #include "fvcFlux.H"
35 #include "fvcMeshPhi.H"
36 #include "fvcSup.H"
37 
38 #include "fvmDdt.H"
39 #include "fvmLaplacian.H"
40 #include "fvmSup.H"
41 
42 #include "upwind.H"
43 
44 #include "CMULES.H"
45 
46 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
47 
49 {
50  forAll(phases(), phasei)
51  {
52  if
53  (
54  mesh()
55  .solution()
56  .solverDict(phases()[phasei].volScalarField::name())
57  .lookupOrDefault<Switch>("implicitPhasePressure", false)
58  )
59  {
60  return true;
61  }
62  }
63 
64  return false;
65 }
66 
67 
69 (
70  const alphaControl& alphaControls,
71  const PtrList<volScalarField>& rAs,
72  const momentumTransferSystem& mts
73 )
74 {
75  const bool boundedPredictor =
76  !alphaControls.MULES.globalBounds
77  && alphaControls.MULES.extremaCoeff == 0;
78 
79  MULES::control MULESBD(alphaControls.MULES);
80  MULESBD.globalBounds = true;
81 
82  const label nAlphaSubCycles = alphaControls.nAlphaSubCycles;
83  const bool LTS = fv::localEulerDdt::enabled(mesh_);
84 
85  // Optional reference phase which is not solved for
86  // but obtained from the sum of the other phases
87  phaseModel* referencePhasePtr = nullptr;
88 
89  // The phases which are solved
90  // i.e. the moving phases less the optional reference phase
91  phaseModelPartialList solvePhases;
92  labelList solveMovingPhaseIndices;
93  label referenceMovingPhaseIndex = -1;
94 
95  if (referencePhaseName_ != word::null)
96  {
97  referencePhasePtr = &phases()[referencePhaseName_];
98 
99  solvePhases.setSize(movingPhases().size() - 1);
100  solveMovingPhaseIndices.setSize(solvePhases.size());
101 
102  label solvePhasesi = 0;
103  forAll(movingPhases(), movingPhasei)
104  {
105  if (&movingPhases()[movingPhasei] != referencePhasePtr)
106  {
107  solveMovingPhaseIndices[solvePhasesi] = movingPhasei;
108  solvePhases.set(solvePhasesi++, &movingPhases()[movingPhasei]);
109  }
110  else
111  {
112  referenceMovingPhaseIndex = movingPhasei;
113  }
114  }
115  }
116  else
117  {
118  solvePhases = movingPhases();
119  solveMovingPhaseIndices = identityMap(solvePhases.size());
120  }
121 
122  forAll(phases(), phasei)
123  {
124  phases()[phasei].correctBoundaryConditions();
125  }
126 
127  // Create sub-list of alphas and phis for the moving phases
128  UPtrList<const volScalarField> alphasMoving(movingPhases().size());
129  UPtrList<const surfaceScalarField> phisMoving(movingPhases().size());
130  forAll(movingPhases(), movingPhasei)
131  {
132  alphasMoving.set(movingPhasei, &movingPhases()[movingPhasei]);
133  phisMoving.set(movingPhasei, &movingPhases()[movingPhasei].phi()());
134  }
135 
136  // Calculate the void fraction
137  volScalarField alphaVoid
138  (
139  IOobject
140  (
141  "alphaVoid",
142  mesh_.time().name(),
143  mesh_
144  ),
145  mesh_,
147  );
148  forAll(stationaryPhases(), stationaryPhasei)
149  {
150  alphaVoid -= stationaryPhases()[stationaryPhasei];
151  }
152 
153  // Calculate the effective flux of the moving phases
154  tmp<surfaceScalarField> tphiMoving(phi_);
155  if (stationaryPhases().size())
156  {
157  tphiMoving = phi_/upwind<scalar>(mesh_, phi_).interpolate(alphaVoid);
158  }
159  const surfaceScalarField& phiMoving = tphiMoving();
160 
161  bool dilatation = false;
162  forAll(movingPhases(), movingPhasei)
163  {
164  if (movingPhases()[movingPhasei].divU().valid())
165  {
166  dilatation = true;
167  break;
168  }
169  }
170 
171  PtrList<volScalarField::Internal> Sps(phases().size());
172  PtrList<volScalarField::Internal> Sus(phases().size());
173 
174  forAll(movingPhases(), movingPhasei)
175  {
176  const phaseModel& phase = movingPhases()[movingPhasei];
177  const volScalarField& alpha = phase;
178  const label phasei = phase.index();
179 
180  Sps.set
181  (
182  phasei,
184  (
185  IOobject
186  (
187  IOobject::groupName("Sp", phase.name()),
188  mesh_.time().name(),
189  mesh_
190  ),
191  mesh_,
193  )
194  );
195 
196  Sus.set
197  (
198  phasei,
200  (
201  IOobject::groupName("Su", phase.name()),
202  min(alpha.v(), scalar(1))
203  *fvc::div(fvc::absolute(phi_, phase.U()))->v()
204  )
205  );
206 
207  if (dilatation)
208  {
209  // Construct the dilatation rate source term
211  (
213  (
214  "vDot",
215  mesh_,
217  )
218  );
219 
220  forAll(phases(), phasej)
221  {
222  const phaseModel& phase2 = phases()[phasej];
223  const volScalarField& alpha2 = phase2;
224 
225  if (&phase2 != &phase)
226  {
227  if (!phase.stationary() && phase.divU().valid())
228  {
229  vDot += alpha2()*phase.divU()()();
230  }
231 
232  if (!phase2.stationary() && phase2.divU().valid())
233  {
234  vDot -= alpha()*phase2.divU()()();
235  }
236  }
237  }
238 
239  volScalarField::Internal& Sp = Sps[phasei];
240  volScalarField::Internal& Su = Sus[phasei];
241 
242  forAll(vDot, celli)
243  {
244  if (vDot[celli] > 0)
245  {
246  Sp[celli] -=
247  vDot[celli]
248  /max
249  (
250  1 - alpha[celli],
251  alphaControls.vDotResidualAlpha
252  );
253  Su[celli] +=
254  vDot[celli]
255  /max
256  (
257  1 - alpha[celli],
258  alphaControls.vDotResidualAlpha
259  );
260  }
261  else if (vDot[celli] < 0)
262  {
263  Sp[celli] +=
264  vDot[celli]
265  /max
266  (
267  alpha[celli],
268  alphaControls.vDotResidualAlpha
269  );
270  }
271  }
272  }
273  }
274 
275  tmp<volScalarField> trSubDeltaT;
276 
277  if (LTS && nAlphaSubCycles > 1)
278  {
279  trSubDeltaT =
280  fv::localEulerDdt::localRSubDeltaT(mesh_, nAlphaSubCycles);
281  }
282 
283  // Cache a list of phase-fractions
284  UPtrList<volScalarField> alphas(phases().size());
285  forAll(phases(), phasei)
286  {
287  alphas.set(phasei, &phases()[phasei]);
288  }
289 
290  // Optional current phase-pressure diffusion fluxes
291  PtrList<surfaceScalarField> alphaPhiDByA(movingPhases().size());
292 
293  if (implicitPhasePressure() && rAs.size())
294  {
295  // Cache the phase-pressure diffusion coefficient
296  const surfaceScalarField alphaDByAf(mts.alphaDByAf(rAs));
297 
298  forAll(movingPhases(), movingPhasei)
299  {
300  const phaseModel& phase = movingPhases()[movingPhasei];
301  const volScalarField& alpha = phase;
302 
303  alphaPhiDByA.set
304  (
305  movingPhasei,
306  alphaDByAf*fvc::snGrad(alpha, "bounded")*mesh_.magSf()
307  );
308  }
309  }
310 
311  for
312  (
314  (
315  alphas,
316  nAlphaSubCycles
317  );
318  !(++alphaSubCycle).end();
319  )
320  {
321  // Bounded implicit predictor phase-fractions
322  PtrList<volScalarField> alphaPreds(movingPhases().size());
323 
324  // Bounded implicit predictor phase-fluxes
325  PtrList<surfaceScalarField> alphaPhiPreds(movingPhases().size());
326 
327  if (alphaControls.MULESCorr)
328  {
329  forAll(solvePhases, solvePhasei)
330  {
331  phaseModel& phase = solvePhases[solvePhasei];
332  volScalarField& alpha = phase;
333 
334  // Bounded implicit predictor matrix using the mean-flux only
335  // to ensure boundedness and phase-consistency
336  fvScalarMatrix alphaEqn
337  (
338  (
339  LTS
340  ? fv::localEulerDdtScheme<scalar>(mesh_).fvmDdt(alpha)
341  : fv::EulerDdtScheme<scalar>(mesh_).fvmDdt(alpha)
342  )
344  (
345  mesh_,
346  phiMoving,
347  upwind<scalar>(mesh_, phiMoving)
348  ).fvmDiv(phiMoving, alpha)
349  ==
350  Sus[phase.index()]
351  + fvm::Sp(Sps[phase.index()], alpha)
352  );
353 
354  alphaEqn.solve();
355 
356  alphaPreds.set
357  (
358  solveMovingPhaseIndices[solvePhasei],
359  alpha
360  );
361 
362  alphaPhiPreds.set
363  (
364  solveMovingPhaseIndices[solvePhasei],
365  alphaEqn.flux()
366  );
367 
368  Info<< phase.name() << " predicted fraction, min, max = "
369  << phase.weightedAverage(mesh_.Vsc()).value()
370  << ' ' << min(phase).value()
371  << ' ' << max(phase).value()
372  << endl;
373  }
374 
375  // Calculate the optional reference phase fraction and flux
376  // from the sum of the other phases
377  if (referencePhasePtr)
378  {
379  volScalarField& referenceAlpha = *referencePhasePtr;
380  referenceAlpha = alphaVoid;
381 
382  forAll(solvePhases, solvePhasei)
383  {
384  referenceAlpha -= solvePhases[solvePhasei];
385  }
386 
387  alphaPreds.set
388  (
389  referenceMovingPhaseIndex,
390  referenceAlpha
391  );
392 
393  alphaPhiPreds.set(referenceMovingPhaseIndex, phi_);
394 
395  forAll(solvePhases, solvePhasei)
396  {
397  alphaPhiPreds[referenceMovingPhaseIndex]
398  -= alphaPhiPreds[solveMovingPhaseIndices[solvePhasei]];
399  }
400  }
401 
402  // Cache the initial predicted phase-fluxes in the phase
403  forAll(solvePhases, solvePhasei)
404  {
405  phaseModel& phase = solvePhases[solvePhasei];
406  const label movingPhasei = solveMovingPhaseIndices[solvePhasei];
407 
408  if (alphaSubCycle.index() == 1)
409  {
410  phase.alphaPhiRef() = alphaPhiPreds[movingPhasei];
411  }
412  else
413  {
414  phase.alphaPhiRef() += alphaPhiPreds[movingPhasei];
415  }
416  }
417  }
418 
419  // Corrector loop
420  for (int acorr=0; acorr<alphaControls.nAlphaCorr; acorr++)
421  {
422  // High-order transport and optional compression phase-fluxes
423  PtrList<surfaceScalarField> alphaPhis(movingPhases().size());
424 
425  // Bounded phase-fluxes
426  PtrList<surfaceScalarField> alphaPhiBDs(movingPhases().size());
427 
428  forAll(movingPhases(), movingPhasei)
429  {
430  const phaseModel& phase = movingPhases()[movingPhasei];
431  const volScalarField& alpha = phase;
432 
433  alphaPhis.set
434  (
435  movingPhasei,
437  (
438  IOobject::groupName("alphaPhi", phase.name()),
439  fvc::flux
440  (
441  phase.phi()(),
442  alpha,
443  "div(phi," + alpha.name() + ')'
444  )
445  )
446  );
447 
448  surfaceScalarField& alphaPhi = alphaPhis[movingPhasei];
449 
450  if (!cAlphas_.empty())
451  {
452  forAll(phases(), phasei)
453  {
454  const phaseModel& phase2 = phases()[phasei];
455  const volScalarField& alpha2 = phase2;
456 
457  if (&phase2 == &phase) continue;
458 
460  (
461  cAlphas_.find(phaseInterface(phase, phase2))
462  );
463 
464  if (cAlpha != cAlphas_.end())
465  {
466  const surfaceScalarField phir
467  (
468  phase.phi() - phase2.phi()
469  );
470 
471  const surfaceScalarField phic
472  (
473  (mag(phi_) + mag(phir))/mesh_.magSf()
474  );
475 
476  const surfaceScalarField phirc
477  (
478  min(cAlpha()*phic, max(phic))
479  *nHatf(alpha, alpha2)
480  );
481 
482  const word phirScheme
483  (
484  "div(phir,"
485  + alpha2.name() + ',' + alpha.name()
486  + ')'
487  );
488 
489  alphaPhi += fvc::flux
490  (
491  -fvc::flux(-phirc, alpha2, phirScheme),
492  alpha,
493  phirScheme
494  );
495  }
496  }
497  }
498 
499  // Add the optional phase-pressure diffusion flux
500  if (alphaPhiDByA.set(movingPhasei))
501  {
502  alphaPhi += alphaPhiDByA[movingPhasei];
503  }
504 
505  // Correct the inflow and outflow boundary conditions
506  // of the phase flux
507  phase.correctInflowOutflow(alphaPhi);
508 
509  // Construct the optional bounded phase-flux
510  if (boundedPredictor)
511  {
512  alphaPhiBDs.set
513  (
514  movingPhasei,
516  (
517  IOobject::groupName("alphaPhiBD", phase.name()),
518  upwind<scalar>(mesh_, phase.phi()()).flux(alpha)
519  )
520  );
521 
522  const surfaceScalarField::Boundary& alphaPhiBf =
523  alphaPhi.boundaryField();
524 
525  surfaceScalarField::Boundary& alphaPhiBDBf =
526  alphaPhiBDs[movingPhasei].boundaryFieldRef();
527 
528  // For non-coupled boundaries
529  // set the bounded flux to the high-order flux
530  forAll(alphaPhiBDBf, patchi)
531  {
532  fvsPatchScalarField& alphaPhiBDPf =
533  alphaPhiBDBf[patchi];
534 
535  if (!alphaPhiBDPf.coupled())
536  {
537  alphaPhiBDPf = alphaPhiBf[patchi];
538  }
539  }
540  }
541  }
542 
543  if (alphaControls.MULESCorr)
544  {
545  // Set local reference to the bounded or high-order phase-fluxes
546  PtrList<surfaceScalarField>& alphaPhiCorrs =
547  boundedPredictor ? alphaPhiBDs : alphaPhis;
548 
549  forAll(movingPhases(), movingPhasei)
550  {
551  const phaseModel& phase = movingPhases()[movingPhasei];
552  const volScalarField& alpha = phase;
553 
554  // Calculate the phase-flux correction
555  // with respect to the bounded implicit prediction
556  alphaPhiCorrs[movingPhasei] -= alphaPhiPreds[movingPhasei];
557 
558  // Limit the phase-flux correction
560  (
561  boundedPredictor
562  ? MULESBD
563  : alphaControls.MULES,
565  alpha,
566  alphaPhiPreds[movingPhasei],
567  alphaPhiCorrs[movingPhasei],
568  Sps[phase.index()],
569  min(alphaVoid.primitiveField(), phase.alphaMax())(),
570  zeroField()
571  );
572  }
573 
574  // Limit the flux corrections of the phases
575  // to ensure the phase fractions sum to 1
576  MULES::limitSumCorr(alphasMoving, alphaPhiCorrs, phiMoving);
577 
578  // Correct the phase-fractions from the phase-flux corrections
579  forAll(solvePhases, solvePhasei)
580  {
581  phaseModel& phase = solvePhases[solvePhasei];
582  volScalarField& alpha = phase;
583  const label movingPhasei =
584  solveMovingPhaseIndices[solvePhasei];
585 
586  const surfaceScalarField& alphaPhiCorr =
587  alphaPhiCorrs[movingPhasei];
588 
590  (
592  alpha,
593  alphaPhiCorr,
594  Sps[phase.index()]
595  );
596  }
597  }
598  else
599  {
600  forAll(movingPhases(), movingPhasei)
601  {
602  const phaseModel& phase = movingPhases()[movingPhasei];
603  const volScalarField& alpha = phase;
604 
606  (
607  boundedPredictor
608  ? MULESBD
609  : alphaControls.MULES,
611  alpha,
612  phiMoving,
613  boundedPredictor
614  ? alphaPhiBDs[movingPhasei]
615  : alphaPhis[movingPhasei],
616  Sps[phase.index()],
617  Sus[phase.index()],
618  min(alphaVoid.primitiveField(), phase.alphaMax())(),
619  zeroField(),
620  false
621  );
622  }
623 
624  // Limit the flux of the phases
625  // to ensure the phase fractions sum to 1
627  (
628  alphasMoving,
629  boundedPredictor ? alphaPhiBDs : alphaPhis,
630  phiMoving
631  );
632 
633  forAll(solvePhases, solvePhasei)
634  {
635  phaseModel& phase = solvePhases[solvePhasei];
636  volScalarField& alpha = phase;
637 
638  const surfaceScalarField& alphaPhi =
639  (boundedPredictor ? alphaPhiBDs : alphaPhis)
640  [solveMovingPhaseIndices[solvePhasei]];
641 
643  (
645  alpha,
646  alphaPhi,
647  Sps[phase.index()],
648  Sus[phase.index()]
649  );
650  }
651  }
652 
653  // Update the phase-fraction of the reference phase
654  if (referencePhasePtr)
655  {
656  volScalarField& referenceAlpha = *referencePhasePtr;
657  referenceAlpha = alphaVoid;
658 
659  forAll(solvePhases, solvePhasei)
660  {
661  referenceAlpha -= solvePhases[solvePhasei];
662  }
663  }
664 
665  if (boundedPredictor)
666  {
667  // Limit the high-order phase-fluxes
668  // to ensure the phase fractions sum to 1
669  MULES::limitSum(alphasMoving, alphaPhis, phiMoving);
670 
671  // Calculate the phase-flux corrections
672  // with respect to the current prediction
673  if (alphaControls.MULESCorr)
674  {
675  forAll(movingPhases(), movingPhasei)
676  {
677  alphaPhis[movingPhasei] -=
678  (
679  alphaPhiPreds[movingPhasei]
680  + alphaPhiBDs[movingPhasei]
681  );
682  }
683  }
684  else
685  {
686  forAll(movingPhases(), movingPhasei)
687  {
688  alphaPhis[movingPhasei] -= alphaPhiBDs[movingPhasei];
689  }
690  }
691 
692  // Limit the phase-flux corrections
693  forAll(movingPhases(), movingPhasei)
694  {
695  const phaseModel& phase = movingPhases()[movingPhasei];
696  const volScalarField& alpha = phase;
697 
699  (
700  alphaControls.MULES,
702  alpha,
703  alphaControls.MULESCorr
704  ? alphaPhiPreds[movingPhasei]
705  : alphaPhiBDs[movingPhasei],
706  alphaPhis[movingPhasei],
707  Sps[phase.index()],
708  min(alphaVoid.primitiveField(), phase.alphaMax())(),
709  zeroField()
710  );
711  }
712 
713  // Limit the flux corrections of the phases
714  // to ensure the phase fractions sum to 1
715  MULES::limitSumCorr(alphasMoving, alphaPhis, phiMoving);
716 
717  // Correct the phase-fractions from the phase-flux corrections
718  forAll(solvePhases, solvePhasei)
719  {
720  phaseModel& phase = solvePhases[solvePhasei];
721  volScalarField& alpha = phase;
722 
723  const surfaceScalarField& alphaPhi =
724  alphaPhis[solveMovingPhaseIndices[solvePhasei]];
725 
727  (
729  alpha,
730  alphaPhi,
731  Sps[phase.index()]
732  );
733  }
734  }
735 
736  // Add the optional implicit phase pressure contribution
737  if (implicitPhasePressure() && rAs.size())
738  {
739  // Cache the phase-pressure diffusion coefficient
740  const surfaceScalarField alphaDByAf(mts.alphaDByAf(rAs));
741 
742  // Add the implicit phase-pressure diffusion contribution
743  // to the phase-fractions and phase-fluxes
744  forAll(solvePhases, solvePhasei)
745  {
746  phaseModel& phase = solvePhases[solvePhasei];
747  volScalarField& alpha = phase;
748 
749  fvScalarMatrix alphaEqn
750  (
752  - fvm::laplacian(alphaDByAf, alpha, "bounded")
753  );
754 
755  alphaEqn.solve
756  (
757  mesh_.solution().solverDict("alpha")
758  .optionalSubDict("phasePressure")
759  );
760 
761  alphaPhis[solveMovingPhaseIndices[solvePhasei]] +=
762  alphaEqn.flux();
763  }
764  }
765 
766  // Report the phase fractions and the phase fraction sum
767  forAll(solvePhases, solvePhasei)
768  {
769  const phaseModel& phase = solvePhases[solvePhasei];
770 
771  Info<< phase.name() << " fraction, min, max = "
772  << phase.weightedAverage(mesh_.Vsc()).value()
773  << ' ' << min(phase).value()
774  << ' ' << max(phase).value()
775  << endl;
776  }
777 
778  // Calculate the reference phase-fraction from the other phases
779  // or re-scale the all the phase fractions to sum to 1 exactly
780  if (referencePhasePtr)
781  {
782  volScalarField& referenceAlpha = *referencePhasePtr;
783  referenceAlpha = alphaVoid;
784 
785  forAll(solvePhases, solvePhasei)
786  {
787  referenceAlpha -= solvePhases[solvePhasei];
788  }
789  }
790  else
791  {
792  volScalarField::Internal sumAlphaMoving
793  (
794  IOobject
795  (
796  "sumAlphaMoving",
797  mesh_.time().name(),
798  mesh_
799  ),
800  mesh_,
802  );
803 
804  forAll(movingPhases(), movingPhasei)
805  {
806  sumAlphaMoving += movingPhases()[movingPhasei];
807  }
808 
809  Info<< "Phase-sum volume fraction, min, max = "
810  << (sumAlphaMoving + 1 - alphaVoid())()
811  .weightedAverage(mesh_.Vsc()).value()
812  << ' ' << min(sumAlphaMoving + 1 - alphaVoid()).value()
813  << ' ' << max(sumAlphaMoving + 1 - alphaVoid()).value()
814  << endl;
815 
816  if (alphaControls.clip)
817  {
818  // Scale moving phase-fractions to sum to alphaVoid
819  forAll(movingPhases(), movingPhasei)
820  {
821  movingPhases()[movingPhasei].internalFieldRef() *=
822  alphaVoid()/sumAlphaMoving;
823  }
824  }
825  }
826 
827  if (alphaControls.MULESCorr)
828  {
829  // For the semi-implicit algorithm relax the phase-fractions
830  // and the phase-flux accumulation
831  if (boundedPredictor)
832  {
833  forAll(movingPhases(), movingPhasei)
834  {
835  phaseModel& phase = movingPhases()[movingPhasei];
836  volScalarField& alpha = phase;
837 
838  alpha +=
839  (1 - alphaControls.MULESCorrRelax)
840  *(alphaPreds[movingPhasei] - alpha);
841 
842  alphaPreds[movingPhasei] = alpha;
843 
844  const surfaceScalarField alphaPhiInc
845  (
846  alphaControls.MULESCorrRelax
847  *(
848  alphaPhiBDs[movingPhasei]
849  + alphaPhis[movingPhasei]
850  )
851  );
852  phase.alphaPhiRef() += alphaPhiInc;
853  alphaPhiPreds[movingPhasei] += alphaPhiInc;
854  }
855  }
856  else
857  {
858  forAll(movingPhases(), movingPhasei)
859  {
860  phaseModel& phase = movingPhases()[movingPhasei];
861  volScalarField& alpha = phase;
862 
863  alpha +=
864  alphaControls.MULESCorrRelax
865  *(alphaPreds[movingPhasei] - alpha);
866 
867  alphaPreds[movingPhasei] = alpha;
868 
869  const surfaceScalarField alphaPhiInc
870  (
871  alphaControls.MULESCorrRelax*alphaPhis[movingPhasei]
872  );
873  phase.alphaPhiRef() += alphaPhiInc;
874  alphaPhiPreds[movingPhasei] += alphaPhiInc;
875  }
876  }
877  }
878  else if (acorr == alphaControls.nAlphaCorr-1)
879  {
880  // Accumulate the phase-fluxes
881  if (boundedPredictor)
882  {
883  forAll(movingPhases(), movingPhasei)
884  {
885  phaseModel& phase = movingPhases()[movingPhasei];
886 
887  if (alphaSubCycle.index() == 1)
888  {
889  phase.alphaPhiRef() = alphaPhiBDs[movingPhasei];
890  phase.alphaPhiRef() += alphaPhis[movingPhasei];
891  }
892  else
893  {
894  phase.alphaPhiRef() += alphaPhiBDs[movingPhasei];
895  phase.alphaPhiRef() += alphaPhis[movingPhasei];
896  }
897  }
898  }
899  else
900  {
901  forAll(movingPhases(), movingPhasei)
902  {
903  phaseModel& phase = movingPhases()[movingPhasei];
904 
905  if (alphaSubCycle.index() == 1)
906  {
907  phase.alphaPhiRef() = alphaPhis[movingPhasei];
908  }
909  else
910  {
911  phase.alphaPhiRef() += alphaPhis[movingPhasei];
912  }
913  }
914  }
915  }
916  }
917  }
918 
919  if (nAlphaSubCycles > 1)
920  {
921  forAll(solvePhases, solvePhasei)
922  {
923  solvePhases[solvePhasei].alphaPhiRef() /= nAlphaSubCycles;
924  }
925  }
926 
927  forAll(solvePhases, solvePhasei)
928  {
929  phaseModel& phase = solvePhases[solvePhasei];
930  volScalarField& alpha = phase;
931 
932  phase.alphaRhoPhiRef() =
933  fvc::interpolate(phase.rho())*phase.alphaPhi();
934 
935  if (alphaControls.clip)
936  {
937  // Clip the phase-fractions between 0 and alphaMax
938  alpha.maxMin(0, phase.alphaMax());
939 
940  if (stationaryPhases().size())
941  {
942  alpha = min(alpha, alphaVoid);
943  }
944  }
945  }
946 
947  // Calculate the optional reference phase fraction and flux
948  // from the sum of the other phases
949  if (referencePhasePtr)
950  {
951  phaseModel& referencePhase = *referencePhasePtr;
952 
953  referencePhase.alphaPhiRef() = phi_;
954 
955  forAll(solvePhases, solvePhasei)
956  {
957  referencePhase.alphaPhiRef() -=
958  solvePhases[solvePhasei].alphaPhi();
959  }
960 
961  referencePhase.alphaRhoPhiRef() =
962  fvc::interpolate(referencePhase.rho())
963  *referencePhase.alphaPhi();
964 
965  volScalarField& referenceAlpha = referencePhase;
966  referenceAlpha = alphaVoid;
967 
968  forAll(solvePhases, solvePhasei)
969  {
970  referenceAlpha -= solvePhases[solvePhasei];
971  }
972  }
973 }
974 
975 
976 // ************************************************************************* //
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dimensioned< Type > weightedAverage(const DimensionedField< scalar, GeoMesh, PrimitiveField2 > &) const
Calculate and return weighted average.
static tmp< DimensionedField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Mesh &mesh, const dimensionSet &, const PrimitiveField< Type > &)
Return a temporary field constructed from name, mesh,.
Generic GeometricBoundaryField class.
Generic GeometricField class.
void maxMin(const dimensioned< Type > &minDt, const dimensioned< Type > &maxDt)
const Boundary & boundaryField() const
Return const-reference to the boundary field.
const Internal::FieldType & primitiveField() const
Return a const-reference to the primitive field.
const Internal & v() const
Return a const-reference to the dimensioned internal field.
Definition: volFieldsI.H:29
An STL-conforming const_iterator.
Definition: HashTable.H:498
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const word & name() const
Return name.
Definition: IOobject.H:307
static word groupName(Name name, const word &group)
void setSize(const label)
Reset size of List.
Definition: List.C:281
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
bool set(const label) const
Is element set.
Definition: UPtrListI.H:87
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
void setSize(const label)
Reset size of UPtrList. This can only be used to set the size.
Definition: UPtrList.C:61
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:58
tmp< SurfaceField< Type > > flux() const
Return the face-flux field from the matrix.
Definition: fvMatrix.C:964
Basic first-order Euler implicit/explicit ddt using only the current and previous time-step values.
Basic second-order convection using face-gradients and Gauss' theorem.
Local time-step first-order Euler implicit/explicit ddt.
static tmp< volScalarField > localRSubDeltaT(const fvMesh &mesh, const label nAlphaSubCycles)
Calculate and return the reciprocal of the local sub-cycling.
Definition: localEulerDdt.C:70
static bool enabled(const fvMesh &mesh)
Return true if LTS is enabled.
Definition: localEulerDdt.C:37
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:86
virtual bool coupled() const
Return true if this patch field is coupled.
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
virtual tmp< SurfaceField< Type > > flux(const VolField< Type > &) const
Return the interpolation weighting factors.
Class which provides interfacial momentum transfer between a number of phases. Drag,...
virtual tmp< surfaceScalarField > alphaDByAf(const PtrList< volScalarField > &rAs) const
Return the phase diffusivity.
Class to represent an interface between phases. Derivations can further specify the configuration of ...
scalar alphaMax() const
Return the maximum phase-fraction (e.g. packing limit)
Definition: phaseModel.C:175
virtual bool stationary() const =0
Return whether the phase is stationary.
virtual const autoPtr< volScalarField > & divU() const =0
Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
virtual tmp< surfaceScalarField > phi() const =0
Return the volumetric flux.
label index() const
Return the index of the phase.
Definition: phaseModel.C:157
virtual surfaceScalarField & alphaRhoPhiRef()=0
Access the mass flux of the phase.
virtual tmp< volVectorField > U() const =0
Return the velocity.
virtual surfaceScalarField & alphaPhiRef()=0
Access the volumetric flux of the phase.
void correctInflowOutflow(surfaceScalarField &alphaPhi) const
Ensure that the flux at inflow/outflow BCs is preserved.
Definition: phaseModel.C:251
const word & name() const
Return the name of this phase.
Definition: phaseModel.C:145
virtual const volScalarField & rho() const =0
Return the density field.
virtual tmp< surfaceScalarField > alphaPhi() const =0
Return the volumetric flux of the phase.
const phaseModelList & phases() const
Return the phase models.
Definition: phaseSystemI.H:104
bool implicitPhasePressure() const
Returns true if the phase pressure is treated implicitly.
const fvMesh & mesh() const
Return the mesh.
Definition: phaseSystemI.H:91
void solve(const alphaControl &alphaControls, const PtrList< volScalarField > &rAs, const momentumTransferSystem &mts)
Solve for the phase fractions.
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:51
Perform a subCycleTime on a field or list of fields.
Definition: subCycle.H:235
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &, const tmp< surfaceScalarField > &, const tmp< surfaceScalarField > &)
Return the face-interpolate of the given cell field.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
A class representing the concept of a field of 0 used to avoid unnecessary manipulations for objects ...
Definition: zeroField.H:53
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.
Calculate the snGrad of the given volField.
Calculate the field for explicit evaluation of implicit and explicit sources.
Calculate the matrix for the first temporal derivative.
Calculate the matrix for the laplacian of the field.
Calculate the matrix for implicit and explicit sources.
label patchi
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiCorr, const SpType &Sp)
void limitSum(const UPtrList< const volScalarField > &psis, const PtrList< surfaceScalarField > &alphaPhiBDs, UPtrList< surfaceScalarField > &psiPhis, const surfaceScalarField &phi)
Definition: MULES.C:272
void limitSumCorr(UPtrList< scalarField > &psiPhiCorrs)
Definition: MULES.C:78
void limitCorr(const control &controls, const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phiBD, surfaceScalarField &phiCorr, const SpType &Sp, const PsiMaxType &psiMax, const PsiMinType &psiMin)
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &psiPhi, const SpType &Sp, const SuType &Su)
void limit(const control &controls, const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phi, surfaceScalarField &psiPhi, const SpType &Sp, const SuType &Su, const PsiMaxType &psiMax, const PsiMinType &psiMin, const bool returnCorr)
bool valid(const PtrList< ModelType > &l)
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< Type > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
tmp< VolField< Type > > Su(const VolField< Type > &su, const VolField< Type > &vf)
Definition: fvcSup.C:44
tmp< VolField< Type > > Sp(const volScalarField &sp, const VolField< Type > &vf)
Definition: fvcSup.C:67
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< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvmLaplacian.C:47
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
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:258
const dimensionSet dimless
messageStream Info
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
labelList identityMap(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Switch globalBounds
Optional switch to select global bounds only.
Definition: MULES.H:80
scalar extremaCoeff
Optional coefficient to relax the local boundedness constraint.
Definition: MULES.H:90
alpha solution control structure
Definition: phaseSystem.H:81
Switch clip
Optionally clip the phase-fractions.
Definition: phaseSystem.H:109
label nAlphaCorr
Number of alpha correctors.
Definition: phaseSystem.H:92
scalar vDotResidualAlpha
Compressibility source stabilisation tolerance.
Definition: phaseSystem.H:101
MULES::control MULES
MULES controls.
Definition: phaseSystem.H:104
bool MULESCorr
Semi-implicit MULES switch.
Definition: phaseSystem.H:95
scalar MULESCorrRelax
Explicit correction relaxation factor (defaults to 0.5)
Definition: phaseSystem.H:98
label nAlphaSubCycles
Number of alpha equation sub-cycles.
Definition: phaseSystem.H:87