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 
26 #include "twoPhaseSolver.H"
27 #include "subCycle.H"
28 #include "CMULES.H"
29 #include "CrankNicolsonDdtScheme.H"
30 #include "fvcFlux.H"
31 #include "fvmSup.H"
32 
33 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
34 
36 (
37  const surfaceScalarField& phi,
38  const volScalarField& alpha
39 )
40 {
41  return fvc::flux
42  (
43  phi,
44  alpha,
46  );
47 }
48 
49 
50 void Foam::solvers::twoPhaseSolver::alphaSolve(const label nAlphaSubCycles)
51 {
52  // Set the off-centering coefficient according to ddt scheme
53  scalar ocCoeff = 0;
54  {
55  tmp<fv::ddtScheme<scalar>> tddtAlpha
56  (
58  (
59  mesh,
60  mesh.schemes().ddt("ddt(alpha)")
61  )
62  );
63  const fv::ddtScheme<scalar>& ddtAlpha = tddtAlpha();
64 
65  if
66  (
69  )
70  {
71  ocCoeff = 0;
72  }
73  else if (isType<fv::CrankNicolsonDdtScheme<scalar>>(ddtAlpha))
74  {
75  if (nAlphaSubCycles > 1)
76  {
78  << "Sub-cycling is not supported "
79  "with the CrankNicolson ddt scheme"
80  << exit(FatalError);
81  }
82 
83  if
84  (
85  alphaRestart
86  || mesh.time().timeIndex() > mesh.time().startTimeIndex() + 1
87  )
88  {
89  ocCoeff =
90  refCast<const fv::CrankNicolsonDdtScheme<scalar>>(ddtAlpha)
91  .ocCoeff();
92  }
93  }
94  else
95  {
97  << "Only Euler and CrankNicolson ddt schemes are supported"
98  << exit(FatalError);
99  }
100  }
101 
102  // Set the time blending factor, 1 for Euler
103  const scalar cnCoeff = 1.0/(1.0 + ocCoeff);
104 
105  tmp<surfaceScalarField> phiCN(phi);
106 
107  // Calculate the Crank-Nicolson off-centred volumetric flux
108  if (ocCoeff > 0)
109  {
111  (
112  "phiCN",
113  cnCoeff*phi + (1.0 - cnCoeff)*phi.oldTime()
114  );
115  }
116 
117  tmp<volScalarField> divU;
118 
119  if (divergent())
120  {
121  divU =
122  (
123  mesh.moving()
124  ? fvc::div(phiCN() + mesh.phi())
125  : fvc::div(phiCN())
126  );
127  }
128 
129  tmp<volScalarField::Internal> Su;
130  tmp<volScalarField::Internal> Sp;
131 
132  alphaSuSp(Su, Sp);
133 
134  if (MULESCorr)
135  {
136  fvScalarMatrix alpha1Eqn
137  (
138  (
139  LTS
140  ? fv::localEulerDdtScheme<scalar>(mesh).fvmDdt(alpha1)
141  : fv::EulerDdtScheme<scalar>(mesh).fvmDdt(alpha1)
142  )
143  + fv::gaussConvectionScheme<scalar>
144  (
145  mesh,
146  phiCN,
147  upwind<scalar>(mesh, phiCN)
148  ).fvmDiv(phiCN, alpha1)
149  );
150 
151  if (divU.valid())
152  {
153  alpha1Eqn -= Su() + fvm::Sp(Sp() + divU(), alpha1);
154  }
155 
156  alpha1Eqn.solve();
157 
158  Info<< "Phase-1 volume fraction = "
159  << alpha1.weightedAverage(mesh.Vsc()).value()
160  << " Min(" << alpha1.name() << ") = " << min(alpha1).value()
161  << " Max(" << alpha1.name() << ") = " << max(alpha1).value()
162  << endl;
163 
164  tmp<surfaceScalarField> talphaPhi1UD(alpha1Eqn.flux());
165  alphaPhi1 = talphaPhi1UD();
166 
167  if (alphaApplyPrevCorr && talphaPhi1Corr0.valid())
168  {
169  Info<< "Applying the previous iteration compression flux" << endl;
171  (
172  MULEScontrols,
173  geometricOneField(),
174  alpha1,
175  alphaPhi1,
176  talphaPhi1Corr0.ref(),
177  oneField(),
178  zeroField()
179  );
180 
181  alphaPhi1 += talphaPhi1Corr0();
182  }
183 
184  // Cache the upwind-flux
185  talphaPhi1Corr0 = talphaPhi1UD;
186 
187  alpha2 = scalar(1) - alpha1;
188  alphaPhi2 = phi - alphaPhi1;
189 
190  correctInterface();
191  }
192 
193  for (int aCorr=0; aCorr<nAlphaCorr; aCorr++)
194  {
195  tmp<volScalarField> talpha1CN(alpha1);
196 
197  if (ocCoeff > 0)
198  {
199  // Preserve the BCs of alpha1 in alpha1CN for interpolation
200  talpha1CN = alpha1.clone();
201  talpha1CN.ref() ==
202  (cnCoeff*alpha1 + (1.0 - cnCoeff)*alpha1.oldTime())();
203  }
204 
205  // Split operator
206  tmp<surfaceScalarField> talphaPhi1Un(alphaPhi(phiCN(), talpha1CN()));
207 
208  if (MULESCorr)
209  {
210  tmp<surfaceScalarField> talphaPhi1Corr(talphaPhi1Un() - alphaPhi1);
211  volScalarField alpha10("alpha10", alpha1);
212 
213  if (divU.valid())
214  {
216  (
217  MULEScontrols,
218  geometricOneField(),
219  alpha1,
220  talphaPhi1Un(),
221  talphaPhi1Corr.ref(),
222  (Sp() + divU())(),
223  oneField(),
224  zeroField()
225  );
226  }
227  else
228  {
230  (
231  MULEScontrols,
232  geometricOneField(),
233  alpha1,
234  talphaPhi1Un(),
235  talphaPhi1Corr.ref(),
236  oneField(),
237  zeroField()
238  );
239  }
240 
241  // Under-relax the correction for all but the 1st corrector
242  if (aCorr == 0)
243  {
244  alphaPhi1 += talphaPhi1Corr();
245  }
246  else
247  {
248  alpha1 = 0.5*alpha1 + 0.5*alpha10;
249  alphaPhi1 += 0.5*talphaPhi1Corr();
250  }
251  }
252  else
253  {
254  alphaPhi1 = talphaPhi1Un;
255 
256  if (divU.valid())
257  {
259  (
260  MULEScontrols,
261  geometricOneField(),
262  alpha1,
263  phiCN,
264  alphaPhi1,
265  Sp(),
266  (Su() + divU()*min(alpha1(), scalar(1)))(),
267  oneField(),
268  zeroField()
269  );
270  }
271  else
272  {
274  (
275  MULEScontrols,
276  geometricOneField(),
277  alpha1,
278  phiCN,
279  alphaPhi1,
280  oneField(),
281  zeroField()
282  );
283  }
284  }
285 
286  alpha2 = scalar(1) - alpha1;
287  alphaPhi2 = phi - alphaPhi1;
288 
289  // Correct only the mixture interface for the interface compression flux
290  correctInterface();
291  }
292 
293  if (alphaApplyPrevCorr && MULESCorr)
294  {
295  talphaPhi1Corr0 = alphaPhi1 - talphaPhi1Corr0;
296 
297  // Register alphaPhiCorr0.<phase1> for redistribution
298  talphaPhi1Corr0.ref().rename
299  (
300  IOobject::groupName("alphaPhiCorr0", alpha1.group())
301  );
302  talphaPhi1Corr0.ref().checkIn();
303  }
304  else
305  {
306  talphaPhi1Corr0.clear();
307  }
308 
309  if
310  (
311  word(mesh.schemes().ddt("ddt(rho,U)"))
312  != fv::EulerDdtScheme<vector>::typeName
313  && word(mesh.schemes().ddt("ddt(rho,U)"))
314  != fv::localEulerDdtScheme<vector>::typeName
315  )
316  {
317  if (ocCoeff > 0)
318  {
319  // Calculate the end-of-time-step alpha flux
320  alphaPhi1 =
321  (alphaPhi1 - (1.0 - cnCoeff)*alphaPhi1.oldTime())/cnCoeff;
322  alphaPhi2 = phi - alphaPhi1;
323  }
324  }
325 
326  Info<< "Phase-1 volume fraction = "
327  << alpha1.weightedAverage(mesh.Vsc()).value()
328  << " Min(" << alpha1.name() << ") = " << min(alpha1).value()
329  << " Max(" << alpha1.name() << ") = " << max(alpha1).value()
330  << endl;
331 }
332 
333 
335 {
336  const label nAlphaSubCycles = ceil(nAlphaSubCyclesPtr->value(alphaCoNum));
337 
338  if (nAlphaSubCycles > 1)
339  {
340  dimensionedScalar totalDeltaT = runTime.deltaT();
341  tmp<volScalarField> trSubDeltaT;
342 
343  if (LTS)
344  {
345  trSubDeltaT =
346  fv::localEulerDdt::localRSubDeltaT(mesh, nAlphaSubCycles);
347  }
348 
349  // Create a temporary alphaPhi1 to accumulate the sub-cycled alphaPhi1
350  tmp<surfaceScalarField> talphaPhi1
351  (
353  (
354  "alphaPhi1",
355  mesh,
356  dimensionedScalar(alphaPhi1.dimensions(), 0)
357  )
358  );
359 
360  UPtrList<volScalarField> alphas({&alpha1, &alpha2});
361 
362  for
363  (
365  (
366  alphas,
367  nAlphaSubCycles
368  );
369  !(++alphaSubCycle).end();
370  )
371  {
372  alphaSolve(nAlphaSubCycles);
373  talphaPhi1.ref() += (runTime.deltaT()/totalDeltaT)*alphaPhi1;
374  }
375 
376  alphaPhi1 = talphaPhi1();
377  alphaPhi2 = phi - talphaPhi1();
378  }
379  else
380  {
381  alphaSolve(nAlphaSubCycles);
382  }
383 }
384 
385 
386 // ************************************************************************* //
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
Generic GeometricField class.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
static word groupName(Name name, const word &group)
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:28
virtual label startTimeIndex() const
Return start time index.
Definition: Time.C:774
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
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:420
const surfaceScalarField & phi() const
Return cell face motion fluxes.
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1799
tmp< DimensionedField< scalar, volMesh > > Vsc() const
Return sub-cycle cell volumes.
ITstream & ddt(const word &name) const
Definition: fvSchemes.C:400
ITstream & div(const word &name) const
Definition: fvSchemes.C:423
Second-oder Crank-Nicolson implicit ddt using the current and previous time-step fields as well as th...
Basic first-order Euler implicit/explicit ddt using only the current and previous time-step values.
Abstract base class for ddt schemes.
Definition: ddtScheme.H:68
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
bool moving() const
Is mesh moving.
Definition: polyMesh.H:473
const fvMesh & mesh
Region mesh.
Definition: solver.H:101
const word divAlphaName
Name of the alpha convection scheme.
Definition: VoFSolver.H:85
const surfaceScalarField & phi
Reference to the mass-flux field.
Definition: VoFSolver.H:216
void alphaPredictor()
Solve for the phase-fractions.
virtual tmp< surfaceScalarField > alphaPhi(const surfaceScalarField &phi, const volScalarField &alpha)
Perform a subCycleTime on a field or list of fields.
Definition: subCycle.H:235
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
Calculate the face-flux of the given field.
Calculate the matrix for implicit and explicit sources.
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 explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &psiPhi, const SpType &Sp, const SuType &Su)
tmp< SurfaceField< typename innerProduct< vector, Type >::type > > flux(const VolField< Type > &vf)
Return the face-flux field obtained from the given volVectorField.
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< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
bool isType(const Type &t)
Check the typeid.
Definition: typeInfo.H:163
messageStream Info
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:62
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
error FatalError
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.