alphaPredictor.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) 2023-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 
27 #include "subCycle.H"
28 #include "MULES.H"
29 #include "fvcFlux.H"
30 
31 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
32 
33 void Foam::solvers::incompressibleMultiphaseVoF::alphaSolve()
34 {
35  const word alphaScheme("div(phi,alpha)");
36  const word alpharScheme("div(phirb,alpha)");
37 
39  phic = min(cAlpha*phic, max(phic));
40 
41  UPtrList<const volScalarField> alphas(phases.size());
42  PtrList<surfaceScalarField> alphaPhis(phases.size());
43 
44  forAll(phases, phasei)
45  {
46  const incompressibleVoFphase& alpha = phases[phasei];
47 
48  alphas.set(phasei, &alpha);
49 
50  alphaPhis.set
51  (
52  phasei,
54  (
55  "phi" + alpha.name() + "Corr",
56  fvc::flux
57  (
58  phi,
59  alpha,
60  alphaScheme
61  )
62  )
63  );
64 
65  surfaceScalarField& alphaPhi = alphaPhis[phasei];
66 
67  forAll(phases, phasej)
68  {
69  incompressibleVoFphase& alpha2 = phases[phasej];
70 
71  if (&alpha2 == &alpha) continue;
72 
73  surfaceScalarField phir(phic*mixture.nHatf(alpha, alpha2));
74 
75  alphaPhi += fvc::flux
76  (
77  -fvc::flux(-phir, alpha2, alpharScheme),
78  alpha,
79  alpharScheme
80  );
81  }
82 
83  // Limit alphaPhi for each phase
85  (
87  1.0/mesh.time().deltaT().value(),
88  geometricOneField(),
89  alpha,
90  phi,
91  alphaPhi,
92  zeroField(),
93  zeroField(),
94  oneField(),
95  zeroField(),
96  false
97  );
98  }
99 
100  MULES::limitSum(alphas, alphaPhis, phi);
101 
102  rhoPhi = Zero;
103 
104  volScalarField sumAlpha
105  (
106  IOobject
107  (
108  "sumAlpha",
109  mesh.time().name(),
110  mesh
111  ),
112  mesh,
114  );
115 
116  forAll(phases, phasei)
117  {
118  incompressibleVoFphase& alpha = phases[phasei];
119  surfaceScalarField& alphaPhi = alphaPhis[phasei];
120 
122  (
123  geometricOneField(),
124  alpha,
125  alphaPhi
126  );
127 
128  rhoPhi += alphaPhi*alpha.rho();
129 
130  Info<< alpha.name() << " volume fraction, min, max = "
131  << alpha.weightedAverage(mesh.V()).value()
132  << ' ' << min(alpha).value()
133  << ' ' << max(alpha).value()
134  << endl;
135 
136  sumAlpha += alpha;
137  }
138 
139  Info<< "Phase-sum volume fraction, min, max = "
140  << sumAlpha.weightedAverage(mesh.V()).value()
141  << ' ' << min(sumAlpha).value()
142  << ' ' << max(sumAlpha).value()
143  << endl;
144 
145  // Correct the sum of the phase-fractions to avoid 'drift'
146  volScalarField sumCorr(1.0 - sumAlpha);
147  forAll(phases, phasei)
148  {
149  incompressibleVoFphase& alpha = phases[phasei];
150  alpha += alpha*sumCorr;
151  }
152 }
153 
154 
156 {
157  const label nAlphaSubCycles = ceil(nAlphaSubCyclesPtr->value(alphaCoNum));
158 
159  if (nAlphaSubCycles > 1)
160  {
161  surfaceScalarField rhoPhiSum
162  (
163  IOobject
164  (
165  "rhoPhiSum",
166  runTime.name(),
167  mesh
168  ),
169  mesh,
170  dimensionedScalar(rhoPhi.dimensions(), 0)
171  );
172 
173  const dimensionedScalar totalDeltaT = runTime.deltaT();
174 
175  UPtrList<volScalarField> alphas(phases.size());
176  forAll(phases, phasei)
177  {
178  alphas.set(phasei, &phases[phasei]);
179  }
180 
181  for
182  (
184  (
185  alphas,
186  nAlphaSubCycles
187  );
188  !(++alphaSubCycle).end();
189  )
190  {
191  alphaSolve();
192  rhoPhiSum += (runTime.deltaT()/totalDeltaT)*rhoPhi;
193  }
194 
195  rhoPhi = rhoPhiSum;
196  }
197  else
198  {
199  alphaSolve();
200  }
201 
202  mixture.correct();
203 }
204 
205 
206 // ************************************************************************* //
MULES: Multidimensional universal limiter for explicit solution.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
IOobject(const word &name, const fileName &instance, const objectRegistry &registry, readOption r=NO_READ, writeOption w=NO_WRITE, bool registerObject=true)
Construct from name, instance, registry, io options.
Definition: IOobject.C:164
dimensionedScalar deltaT() const
Return time step.
Definition: TimeStateI.H:46
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
bool set(const label) const
Is element set.
Definition: UPtrListI.H:87
const Type & value() const
Return const reference to value.
const word & name() const
Return const reference to name.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:420
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
tmp< surfaceScalarField > nHatf(const volScalarField &alpha1, const volScalarField &alpha2) const
const fvMesh & mesh
Region mesh.
Definition: solver.H:101
const surfaceScalarField & phi
Reference to the mass-flux field.
Definition: VoFSolver.H:216
surfaceScalarField rhoPhi
Mass flux field.
Definition: VoFSolver.H:116
virtual void alphaPredictor()
Solve for the phase-fractions.
UPtrListDictionary< incompressibleVoFphase > & phases
Reference to the phases.
incompressibleMultiphaseVoFMixture & mixture
The compressible two-phase mixture.
MULES::control MULEScontrols
MULES controls.
Perform a subCycleTime on a field or list of fields.
Definition: subCycle.H:235
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
Calculate the face-flux of the given field.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
void limitSum(const UPtrList< const volScalarField > &psis, const PtrList< surfaceScalarField > &alphaPhiBDs, UPtrList< surfaceScalarField > &psiPhis, const surfaceScalarField &phi)
Definition: MULES.C:272
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)
tmp< SurfaceField< typename innerProduct< vector, Type >::type > > flux(const VolField< Type > &vf)
Return the face-flux field obtained from the given volVectorField.
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
const dimensionSet dimless
SurfaceField< scalar > surfaceScalarField
messageStream Info
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:62
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.