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 "CMULES.H"
29 #include "fvcFlux.H"
30 #include "fvcMeshPhi.H"
31 
32 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
33 
34 void Foam::solvers::compressibleMultiphaseVoF::alphaSolve()
35 {
36  const word alphaScheme("div(phi,alpha)");
37  const word alpharScheme("div(phirb,alpha)");
38 
40  phic = min(cAlpha*phic, max(phic));
41 
42  UPtrList<const volScalarField> alphas(phases.size());
43  PtrList<surfaceScalarField> alphaPhis(phases.size());
44 
45  forAll(phases, phasei)
46  {
47  const compressibleVoFphase& alpha = phases[phasei];
48 
49  alphas.set(phasei, &alpha);
50 
51  alphaPhis.set
52  (
53  phasei,
55  (
56  "phi" + alpha.name() + "Corr",
57  fvc::flux
58  (
59  phi,
60  alpha,
61  alphaScheme
62  )
63  )
64  );
65 
66  surfaceScalarField& alphaPhi = alphaPhis[phasei];
67 
68  forAll(phases, phasej)
69  {
70  compressibleVoFphase& alpha2 = phases[phasej];
71 
72  if (&alpha2 == &alpha) continue;
73 
74  surfaceScalarField phir(phic*mixture.nHatf(alpha, alpha2));
75 
76  alphaPhi += fvc::flux
77  (
78  -fvc::flux(-phir, alpha2, alpharScheme),
79  alpha,
80  alpharScheme
81  );
82  }
83 
84  // Limit alphaPhi for each phase
86  (
88  1.0/mesh.time().deltaT().value(),
89  geometricOneField(),
90  alpha,
91  phi,
92  alphaPhi,
93  zeroField(),
94  zeroField(),
95  oneField(),
96  zeroField(),
97  false
98  );
99  }
100 
101  MULES::limitSum(alphas, alphaPhis, phi);
102 
103  rhoPhi = Zero;
104 
105  volScalarField sumAlpha
106  (
107  IOobject
108  (
109  "sumAlpha",
110  mesh.time().name(),
111  mesh
112  ),
113  mesh,
115  );
116 
118 
119  forAll(phases, phasei)
120  {
121  compressibleVoFphase& alpha = phases[phasei];
122 
123  surfaceScalarField& alphaPhi = alphaPhis[phasei];
124 
126  (
127  IOobject
128  (
129  "Sp",
130  mesh.time().name(),
131  mesh
132  ),
133  mesh,
134  dimensionedScalar(alpha.vDot().dimensions(), 0)
135  );
136 
138  (
139  IOobject
140  (
141  "Su",
142  mesh.time().name(),
143  mesh
144  ),
145  // Divergence term is handled explicitly to be
146  // consistent with the explicit transport solution
147  divU.v()*min(alpha.v(), scalar(1))
148  );
149 
150  {
151  const scalarField& vDot = alpha.vDot();
152 
153  forAll(vDot, celli)
154  {
155  if (vDot[celli] < 0.0 && alpha[celli] > 0.0)
156  {
157  Sp[celli] += vDot[celli]*alpha[celli];
158  Su[celli] -= vDot[celli]*alpha[celli];
159  }
160  else if (vDot[celli] > 0.0 && alpha[celli] < 1.0)
161  {
162  Sp[celli] -= vDot[celli]*(1.0 - alpha[celli]);
163  }
164  }
165  }
166 
167 
168  forAll(phases, phasej)
169  {
170  const compressibleVoFphase& alpha2 = phases[phasej];
171 
172  if (&alpha2 == &alpha) continue;
173 
174  const scalarField& vDot2 = alpha2.vDot();
175 
176  forAll(vDot2, celli)
177  {
178  if (vDot2[celli] > 0.0 && alpha2[celli] < 1.0)
179  {
180  Sp[celli] -= vDot2[celli]*(1.0 - alpha2[celli]);
181  Su[celli] += vDot2[celli]*alpha[celli];
182  }
183  else if (vDot2[celli] < 0.0 && alpha2[celli] > 0.0)
184  {
185  Sp[celli] += vDot2[celli]*alpha2[celli];
186  }
187  }
188  }
189 
191  (
192  geometricOneField(),
193  alpha,
194  alphaPhi,
195  Sp,
196  Su
197  );
198 
199  rhoPhi += fvc::interpolate(alpha.thermo().rho())*alphaPhi;
200 
201  Info<< alpha.name() << " volume fraction, min, max = "
202  << alpha.weightedAverage(mesh.V()).value()
203  << ' ' << min(alpha).value()
204  << ' ' << max(alpha).value()
205  << endl;
206 
207  sumAlpha += alpha;
208  }
209 
210  // Correct the sum of the phase-fractions to avoid 'drift'
211  const volScalarField sumCorr(1.0 - sumAlpha);
212  forAll(phases, phasei)
213  {
214  compressibleVoFphase& alpha = phases[phasei];
215  alpha += alpha*sumCorr;
216  }
217 }
218 
219 
221 {
222  const label nAlphaSubCycles = ceil(nAlphaSubCyclesPtr->value(alphaCoNum));
223 
224  if (nAlphaSubCycles > 1)
225  {
226  surfaceScalarField rhoPhiSum
227  (
228  IOobject
229  (
230  "rhoPhiSum",
231  runTime.name(),
232  mesh
233  ),
234  mesh,
235  dimensionedScalar(rhoPhi.dimensions(), 0)
236  );
237 
238  const dimensionedScalar totalDeltaT = runTime.deltaT();
239 
240  UPtrList<volScalarField> alphas(phases.size());
241  forAll(phases, phasei)
242  {
243  alphas.set(phasei, &phases[phasei]);
244  }
245 
246  for
247  (
249  (
250  alphas,
251  nAlphaSubCycles
252  );
253  !(++alphaSubCycle).end();
254  )
255  {
256  alphaSolve();
257  rhoPhiSum += (runTime.deltaT()/totalDeltaT)*rhoPhi;
258  }
259 
260  rhoPhi = rhoPhiSum;
261  }
262  else
263  {
264  alphaSolve();
265  }
266 
267  mixture.correct();
268 }
269 
270 
271 // ************************************************************************* //
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Generic GeometricField class.
DimensionedField< Type, GeoMesh, PrimitiveField > Internal
Type of the internal field from which this GeometricField is derived.
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
autoPtr< volScalarField > divU
Pointer to the momentum divergence field.
Definition: VoFSolver.H:138
const volVectorField & U
Reference to the velocity field.
Definition: VoFSolver.H:213
UPtrListDictionary< compressibleVoFphase > & phases
Reference to the phases.
virtual void alphaPredictor()
Solve for the phase-fractions.
compressibleMultiphaseVoFMixture & 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.
Calculate the mesh motion flux and convert fluxes from absolute to relative and back.
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 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 > > 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
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)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
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.