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 
31 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
32 
33 void Foam::solvers::incompressibleMultiphaseVoF::alphaSolve
34 (
35  const dictionary& alphaControls
36 )
37 {
38  const scalar cAlpha(alphaControls.lookup<scalar>("cAlpha"));
39 
40  const word alphaScheme("div(phi,alpha)");
41  const word alpharScheme("div(phirb,alpha)");
42 
44  phic = min(cAlpha*phic, max(phic));
45 
46  UPtrList<const volScalarField> alphas(phases.size());
47  PtrList<surfaceScalarField> alphaPhis(phases.size());
48 
49  forAll(phases, phasei)
50  {
51  const incompressibleVoFphase& alpha = phases[phasei];
52 
53  alphas.set(phasei, &alpha);
54 
55  alphaPhis.set
56  (
57  phasei,
59  (
60  "phi" + alpha.name() + "Corr",
61  fvc::flux
62  (
63  phi,
64  alpha,
65  alphaScheme
66  )
67  )
68  );
69 
70  surfaceScalarField& alphaPhi = alphaPhis[phasei];
71 
72  forAll(phases, phasej)
73  {
74  incompressibleVoFphase& alpha2 = phases[phasej];
75 
76  if (&alpha2 == &alpha) continue;
77 
78  surfaceScalarField phir(phic*mixture.nHatf(alpha, alpha2));
79 
80  alphaPhi += fvc::flux
81  (
82  -fvc::flux(-phir, alpha2, alpharScheme),
83  alpha,
84  alpharScheme
85  );
86  }
87 
88  // Limit alphaPhi for each phase
90  (
91  1.0/mesh.time().deltaT().value(),
92  geometricOneField(),
93  alpha,
94  phi,
95  alphaPhi,
96  zeroField(),
97  zeroField(),
98  oneField(),
99  zeroField(),
100  false
101  );
102  }
103 
104  MULES::limitSum(alphas, alphaPhis, phi);
105 
106  rhoPhi = Zero;
107 
108  volScalarField sumAlpha
109  (
110  IOobject
111  (
112  "sumAlpha",
113  mesh.time().name(),
114  mesh
115  ),
116  mesh,
118  );
119 
120  forAll(phases, phasei)
121  {
122  incompressibleVoFphase& alpha = phases[phasei];
123  surfaceScalarField& alphaPhi = alphaPhis[phasei];
124 
126  (
127  geometricOneField(),
128  alpha,
129  alphaPhi
130  );
131 
132  rhoPhi += alphaPhi*alpha.rho();
133 
134  Info<< alpha.name() << " volume fraction, min, max = "
135  << alpha.weightedAverage(mesh.V()).value()
136  << ' ' << min(alpha).value()
137  << ' ' << max(alpha).value()
138  << endl;
139 
140  sumAlpha += alpha;
141  }
142 
143  Info<< "Phase-sum volume fraction, min, max = "
144  << sumAlpha.weightedAverage(mesh.V()).value()
145  << ' ' << min(sumAlpha).value()
146  << ' ' << max(sumAlpha).value()
147  << endl;
148 
149  // Correct the sum of the phase-fractions to avoid 'drift'
150  volScalarField sumCorr(1.0 - sumAlpha);
151  forAll(phases, phasei)
152  {
153  incompressibleVoFphase& alpha = phases[phasei];
154  alpha += alpha*sumCorr;
155  }
156 }
157 
158 
160 {
161  const dictionary& alphaControls = mesh.solution().solverDict("alpha");
162 
163  const label nAlphaSubCycles(alphaControls.lookup<label>("nAlphaSubCycles"));
164 
165  if (nAlphaSubCycles > 1)
166  {
167  surfaceScalarField rhoPhiSum
168  (
169  IOobject
170  (
171  "rhoPhiSum",
172  runTime.name(),
173  mesh
174  ),
175  mesh,
176  dimensionedScalar(rhoPhi.dimensions(), 0)
177  );
178 
179  const dimensionedScalar totalDeltaT = runTime.deltaT();
180 
181  List<volScalarField*> alphaPtrs(phases.size());
182  forAll(phases, phasei)
183  {
184  alphaPtrs[phasei] = &phases[phasei];
185  }
186 
187  for
188  (
190  (
191  alphaPtrs,
193  );
194  !(++alphaSubCycle).end();
195  )
196  {
197  alphaSolve(alphaControls);
198  rhoPhiSum += (runTime.deltaT()/totalDeltaT)*rhoPhi;
199  }
200 
201  rhoPhi = rhoPhiSum;
202  }
203  else
204  {
205  alphaSolve(alphaControls);
206  }
207 
208  mixture.correct();
209 }
210 
211 
212 // ************************************************************************* //
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.
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
virtual void alphaPredictor()
Solve for the phase-fractions.
UPtrListDictionary< incompressibleVoFphase > & phases
Reference to the phases.
incompressibleMultiphaseVoFMixture & 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.
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.
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
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)