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 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 dictionary& alphaControls
37 )
38 {
39  const scalar cAlpha(alphaControls.lookup<scalar>("cAlpha"));
40 
41  const word alphaScheme("div(phi,alpha)");
42  const word alpharScheme("div(phirb,alpha)");
43 
45  phic = min(cAlpha*phic, max(phic));
46 
47  UPtrList<const volScalarField> alphas(phases.size());
48  PtrList<surfaceScalarField> alphaPhis(phases.size());
49 
50  forAll(phases, phasei)
51  {
52  const compressibleVoFphase& alpha = phases[phasei];
53 
54  alphas.set(phasei, &alpha);
55 
56  alphaPhis.set
57  (
58  phasei,
60  (
61  "phi" + alpha.name() + "Corr",
62  fvc::flux
63  (
64  phi,
65  alpha,
66  alphaScheme
67  )
68  )
69  );
70 
71  surfaceScalarField& alphaPhi = alphaPhis[phasei];
72 
73  forAll(phases, phasej)
74  {
75  compressibleVoFphase& alpha2 = phases[phasej];
76 
77  if (&alpha2 == &alpha) continue;
78 
79  surfaceScalarField phir(phic*mixture.nHatf(alpha, alpha2));
80 
81  alphaPhi += fvc::flux
82  (
83  -fvc::flux(-phir, alpha2, alpharScheme),
84  alpha,
85  alpharScheme
86  );
87  }
88 
89  // Limit alphaPhi for each phase
91  (
92  1.0/mesh.time().deltaT().value(),
93  geometricOneField(),
94  alpha,
95  phi,
96  alphaPhi,
97  zeroField(),
98  zeroField(),
99  oneField(),
100  zeroField(),
101  false
102  );
103  }
104 
105  MULES::limitSum(alphas, alphaPhis, phi);
106 
107  rhoPhi = Zero;
108 
109  volScalarField sumAlpha
110  (
111  IOobject
112  (
113  "sumAlpha",
114  mesh.time().name(),
115  mesh
116  ),
117  mesh,
119  );
120 
122 
123  forAll(phases, phasei)
124  {
125  compressibleVoFphase& alpha = phases[phasei];
126 
127  surfaceScalarField& alphaPhi = alphaPhis[phasei];
128 
130  (
131  IOobject
132  (
133  "Sp",
134  mesh.time().name(),
135  mesh
136  ),
137  mesh,
138  dimensionedScalar(alpha.dgdt().dimensions(), 0)
139  );
140 
142  (
143  IOobject
144  (
145  "Su",
146  mesh.time().name(),
147  mesh
148  ),
149  // Divergence term is handled explicitly to be
150  // consistent with the explicit transport solution
151  divU.v()*min(alpha.v(), scalar(1))
152  );
153 
154  {
155  const scalarField& dgdt = alpha.dgdt();
156 
157  forAll(dgdt, celli)
158  {
159  if (dgdt[celli] < 0.0 && alpha[celli] > 0.0)
160  {
161  Sp[celli] += dgdt[celli]*alpha[celli];
162  Su[celli] -= dgdt[celli]*alpha[celli];
163  }
164  else if (dgdt[celli] > 0.0 && alpha[celli] < 1.0)
165  {
166  Sp[celli] -= dgdt[celli]*(1.0 - alpha[celli]);
167  }
168  }
169  }
170 
171 
172  forAll(phases, phasej)
173  {
174  const compressibleVoFphase& alpha2 = phases[phasej];
175 
176  if (&alpha2 == &alpha) continue;
177 
178  const scalarField& dgdt2 = alpha2.dgdt();
179 
180  forAll(dgdt2, celli)
181  {
182  if (dgdt2[celli] > 0.0 && alpha2[celli] < 1.0)
183  {
184  Sp[celli] -= dgdt2[celli]*(1.0 - alpha2[celli]);
185  Su[celli] += dgdt2[celli]*alpha[celli];
186  }
187  else if (dgdt2[celli] < 0.0 && alpha2[celli] > 0.0)
188  {
189  Sp[celli] += dgdt2[celli]*alpha2[celli];
190  }
191  }
192  }
193 
195  (
196  geometricOneField(),
197  alpha,
198  alphaPhi,
199  Sp,
200  Su
201  );
202 
203  rhoPhi += fvc::interpolate(alpha.thermo().rho())*alphaPhi;
204 
205  Info<< alpha.name() << " volume fraction, min, max = "
206  << alpha.weightedAverage(mesh.V()).value()
207  << ' ' << min(alpha).value()
208  << ' ' << max(alpha).value()
209  << endl;
210 
211  sumAlpha += alpha;
212  }
213 
214  // Correct the sum of the phase-fractions to avoid 'drift'
215  const volScalarField sumCorr(1.0 - sumAlpha);
216  forAll(phases, phasei)
217  {
218  compressibleVoFphase& alpha = phases[phasei];
219  alpha += alpha*sumCorr;
220  }
221 }
222 
223 
225 {
226  const dictionary& alphaControls = mesh.solution().solverDict("alpha");
227 
228  const label nAlphaSubCycles(alphaControls.lookup<label>("nAlphaSubCycles"));
229 
230  if (nAlphaSubCycles > 1)
231  {
232  surfaceScalarField rhoPhiSum
233  (
234  IOobject
235  (
236  "rhoPhiSum",
237  runTime.name(),
238  mesh
239  ),
240  mesh,
241  dimensionedScalar(rhoPhi.dimensions(), 0)
242  );
243 
244  const dimensionedScalar totalDeltaT = runTime.deltaT();
245 
246  List<volScalarField*> alphaPtrs(phases.size());
247  forAll(phases, phasei)
248  {
249  alphaPtrs[phasei] = &phases[phasei];
250  }
251 
252  for
253  (
255  (
256  alphaPtrs,
258  );
259  !(++alphaSubCycle).end();
260  )
261  {
262  alphaSolve(alphaControls);
263  rhoPhiSum += (runTime.deltaT()/totalDeltaT)*rhoPhi;
264  }
265 
266  rhoPhi = rhoPhiSum;
267  }
268  else
269  {
270  alphaSolve(alphaControls);
271  }
272 
273  mixture.correct();
274 }
275 
276 
277 // ************************************************************************* //
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const label nAlphaSubCycles(alphaControls.lookup< label >("nAlphaSubCycles"))
const dictionary & alphaControls
Definition: alphaControls.H:1
Generic GeometricField class.
DimensionedField< Type, GeoMesh > 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:167
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
dimensionedScalar deltaT() const
Return time step.
Definition: TimeStateI.H:46
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
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:402
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:94
const surfaceScalarField & phi
Reference to the mass-flux field.
Definition: VoFSolver.H:212
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:209
UPtrListDictionary< compressibleVoFphase > & phases
Reference to the phases.
virtual void alphaPredictor()
Solve for the phase-fractions.
compressibleMultiphaseVoFMixture & mixture
The compressible two-phase mixture.
Perform a subCycleTime on a field or list of fields.
Definition: subCycle.H:235
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 explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
void limit(const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phi, surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su, const PsiMaxType &psiMax, const PsiMinType &psiMin, const bool returnCorr)
void limitSum(UPtrList< scalarField > &phiPsiCorrs)
Definition: MULES.C:30
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
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
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:251
const dimensionSet dimless
SurfaceField< scalar > surfaceScalarField
messageStream Info
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)