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 
26 #include "twoPhaseSolver.H"
27 #include "subCycle.H"
28 #include "interfaceCompression.H"
29 #include "CMULES.H"
30 #include "CrankNicolsonDdtScheme.H"
31 #include "fvcFlux.H"
32 #include "fvmSup.H"
33 
34 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
35 
37 (
38  const surfaceScalarField& phi,
39  const volScalarField& alpha,
41 )
42 {
43  const word alphaScheme(mesh.schemes().div(divAlphaName)[1].wordToken());
44 
45  ITstream compressionScheme
46  (
47  compressionSchemes.found(alphaScheme)
49  : ITstream
50  (
52  tokenList
53  {
54  word("Gauss"),
55  word("interfaceCompression"),
56  alphaScheme,
57  alphaControls.lookup<scalar>("cAlpha")
58  }
59  )
60  );
61 
62  return fvc::flux
63  (
64  phi,
65  alpha,
66  compressionScheme
67  );
68 }
69 
70 
71 void Foam::solvers::twoPhaseSolver::alphaSolve
72 (
74 )
75 {
76  const label nAlphaSubCycles(alphaControls.lookup<label>("nAlphaSubCycles"));
77 
78  const label nAlphaCorr(alphaControls.lookup<label>("nAlphaCorr"));
79 
80  const bool MULESCorr
81  (
82  alphaControls.lookupOrDefault<Switch>("MULESCorr", false)
83  );
84 
85  // Apply the compression correction from the previous iteration
86  // Improves efficiency for steady-simulations but can only be applied
87  // once the alpha field is reasonably steady, i.e. fully developed
88  const bool alphaApplyPrevCorr
89  (
90  alphaControls.lookupOrDefault<Switch>("alphaApplyPrevCorr", false)
91  );
92 
93 
94  // Set the off-centering coefficient according to ddt scheme
95  scalar ocCoeff = 0;
96  {
97  tmp<fv::ddtScheme<scalar>> tddtAlpha
98  (
100  (
101  mesh,
102  mesh.schemes().ddt("ddt(alpha)")
103  )
104  );
105  const fv::ddtScheme<scalar>& ddtAlpha = tddtAlpha();
106 
107  if
108  (
111  )
112  {
113  ocCoeff = 0;
114  }
115  else if (isType<fv::CrankNicolsonDdtScheme<scalar>>(ddtAlpha))
116  {
117  if (nAlphaSubCycles > 1)
118  {
120  << "Sub-cycling is not supported "
121  "with the CrankNicolson ddt scheme"
122  << exit(FatalError);
123  }
124 
125  if
126  (
127  alphaRestart
128  || mesh.time().timeIndex() > mesh.time().startTimeIndex() + 1
129  )
130  {
131  ocCoeff =
132  refCast<const fv::CrankNicolsonDdtScheme<scalar>>(ddtAlpha)
133  .ocCoeff();
134  }
135  }
136  else
137  {
139  << "Only Euler and CrankNicolson ddt schemes are supported"
140  << exit(FatalError);
141  }
142  }
143 
144  // Set the time blending factor, 1 for Euler
145  const scalar cnCoeff = 1.0/(1.0 + ocCoeff);
146 
147  tmp<surfaceScalarField> phiCN(phi);
148 
149  // Calculate the Crank-Nicolson off-centred volumetric flux
150  if (ocCoeff > 0)
151  {
153  (
154  "phiCN",
155  cnCoeff*phi + (1.0 - cnCoeff)*phi.oldTime()
156  );
157  }
158 
159  tmp<volScalarField> divU;
160 
161  if (divergent())
162  {
163  divU =
164  (
165  mesh.moving()
166  ? fvc::div(phiCN() + mesh.phi())
167  : fvc::div(phiCN())
168  );
169  }
170 
171  tmp<volScalarField::Internal> Su;
172  tmp<volScalarField::Internal> Sp;
173 
174  alphaSuSp(Su, Sp, alphaControls);
175 
176  if (MULESCorr)
177  {
178  fvScalarMatrix alpha1Eqn
179  (
180  (
181  LTS
182  ? fv::localEulerDdtScheme<scalar>(mesh).fvmDdt(alpha1)
183  : fv::EulerDdtScheme<scalar>(mesh).fvmDdt(alpha1)
184  )
185  + fv::gaussConvectionScheme<scalar>
186  (
187  mesh,
188  phiCN,
189  upwind<scalar>(mesh, phiCN)
190  ).fvmDiv(phiCN, alpha1)
191  );
192 
193  if (divU.valid())
194  {
195  alpha1Eqn -= Su() + fvm::Sp(Sp() + divU(), alpha1);
196  }
197 
198  alpha1Eqn.solve();
199 
200  Info<< "Phase-1 volume fraction = "
201  << alpha1.weightedAverage(mesh.Vsc()).value()
202  << " Min(" << alpha1.name() << ") = " << min(alpha1).value()
203  << " Max(" << alpha1.name() << ") = " << max(alpha1).value()
204  << endl;
205 
206  tmp<surfaceScalarField> talphaPhi1UD(alpha1Eqn.flux());
207  alphaPhi1 = talphaPhi1UD();
208 
209  if (alphaApplyPrevCorr && talphaPhi1Corr0.valid())
210  {
211  Info<< "Applying the previous iteration compression flux" << endl;
213  (
214  geometricOneField(),
215  alpha1,
216  alphaPhi1,
217  talphaPhi1Corr0.ref(),
218  oneField(),
219  zeroField()
220  );
221 
222  alphaPhi1 += talphaPhi1Corr0();
223  }
224 
225  // Cache the upwind-flux
226  talphaPhi1Corr0 = talphaPhi1UD;
227 
228  alpha2 = scalar(1) - alpha1;
229  alphaPhi2 = phi - alphaPhi1;
230 
231  correctInterface();
232  }
233 
234  for (int aCorr=0; aCorr<nAlphaCorr; aCorr++)
235  {
236  tmp<volScalarField> talpha1CN(alpha1);
237 
238  if (ocCoeff > 0)
239  {
240  // Preserve the BCs of alpha1 in alpha1CN for interpolation
241  talpha1CN = alpha1.clone();
242  talpha1CN.ref() ==
243  (cnCoeff*alpha1 + (1.0 - cnCoeff)*alpha1.oldTime())();
244  }
245 
246  // Split operator
247  tmp<surfaceScalarField> talphaPhi1Un
248  (
249  alphaPhi
250  (
251  phiCN(),
252  talpha1CN(),
254  )
255  );
256 
257  if (MULESCorr)
258  {
259  tmp<surfaceScalarField> talphaPhi1Corr(talphaPhi1Un() - alphaPhi1);
260  volScalarField alpha10("alpha10", alpha1);
261 
262  if (divU.valid())
263  {
265  (
266  geometricOneField(),
267  alpha1,
268  talphaPhi1Un(),
269  talphaPhi1Corr.ref(),
270  (Sp() + divU())(),
271  (-(Sp() + divU())*alpha1)(),
272  oneField(),
273  zeroField()
274  );
275  }
276  else
277  {
279  (
280  geometricOneField(),
281  alpha1,
282  talphaPhi1Un(),
283  talphaPhi1Corr.ref(),
284  oneField(),
285  zeroField()
286  );
287  }
288 
289  // Under-relax the correction for all but the 1st corrector
290  if (aCorr == 0)
291  {
292  alphaPhi1 += talphaPhi1Corr();
293  }
294  else
295  {
296  alpha1 = 0.5*alpha1 + 0.5*alpha10;
297  alphaPhi1 += 0.5*talphaPhi1Corr();
298  }
299  }
300  else
301  {
302  alphaPhi1 = talphaPhi1Un;
303 
304  if (divU.valid())
305  {
307  (
308  geometricOneField(),
309  alpha1,
310  phiCN,
311  alphaPhi1,
312  Sp(),
313  (Su() + divU()*min(alpha1(), scalar(1)))(),
314  oneField(),
315  zeroField()
316  );
317  }
318  else
319  {
321  (
322  geometricOneField(),
323  alpha1,
324  phiCN,
325  alphaPhi1,
326  oneField(),
327  zeroField()
328  );
329  }
330  }
331 
332  alpha2 = scalar(1) - alpha1;
333  alphaPhi2 = phi - alphaPhi1;
334 
335  // Correct only the mixture interface for the interface compression flux
336  correctInterface();
337  }
338 
340  {
341  talphaPhi1Corr0 = alphaPhi1 - talphaPhi1Corr0;
342 
343  // Register alphaPhiCorr0.<phase1> for redistribution
344  talphaPhi1Corr0.ref().rename
345  (
346  IOobject::groupName("alphaPhiCorr0", alpha1.group())
347  );
348  talphaPhi1Corr0.ref().checkIn();
349  }
350  else
351  {
352  talphaPhi1Corr0.clear();
353  }
354 
355  if
356  (
357  word(mesh.schemes().ddt("ddt(rho,U)"))
358  != fv::EulerDdtScheme<vector>::typeName
359  && word(mesh.schemes().ddt("ddt(rho,U)"))
360  != fv::localEulerDdtScheme<vector>::typeName
361  )
362  {
363  if (ocCoeff > 0)
364  {
365  // Calculate the end-of-time-step alpha flux
366  alphaPhi1 =
367  (alphaPhi1 - (1.0 - cnCoeff)*alphaPhi1.oldTime())/cnCoeff;
368  alphaPhi2 = phi - alphaPhi1;
369  }
370  }
371 
372  Info<< "Phase-1 volume fraction = "
373  << alpha1.weightedAverage(mesh.Vsc()).value()
374  << " Min(" << alpha1.name() << ") = " << min(alpha1).value()
375  << " Max(" << alpha1.name() << ") = " << max(alpha1).value()
376  << endl;
377 }
378 
379 
381 {
382  const dictionary& alphaControls = mesh.solution().solverDict(alpha1.name());
383 
384  const label nAlphaSubCycles(alphaControls.lookup<label>("nAlphaSubCycles"));
385 
386 
387  if (nAlphaSubCycles > 1)
388  {
389  dimensionedScalar totalDeltaT = runTime.deltaT();
390  tmp<volScalarField> trSubDeltaT;
391 
392  if (LTS)
393  {
394  trSubDeltaT =
396  }
397 
398  // Create a temporary alphaPhi1 to accumulate the sub-cycled alphaPhi1
399  tmp<surfaceScalarField> talphaPhi1
400  (
402  (
403  "alphaPhi1",
404  mesh,
405  dimensionedScalar(alphaPhi1.dimensions(), 0)
406  )
407  );
408 
409  List<volScalarField*> alphaPtrs({&alpha1, &alpha2});
410 
411  for
412  (
414  (
415  alphaPtrs,
417  );
418  !(++alphaSubCycle).end();
419  )
420  {
421  alphaSolve(alphaControls);
422  talphaPhi1.ref() += (runTime.deltaT()/totalDeltaT)*alphaPhi1;
423  }
424 
425  alphaPhi1 = talphaPhi1();
426  alphaPhi2 = phi - talphaPhi1();
427  }
428  else
429  {
430  alphaSolve(alphaControls);
431  }
432 }
433 
434 
435 // ************************************************************************* //
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
const label nAlphaSubCycles(alphaControls.lookup< label >("nAlphaSubCycles"))
const dictionary & alphaControls
Definition: alphaControls.H:1
const bool MULESCorr(alphaControls.lookupOrDefault< Switch >("MULESCorr", false))
const label nAlphaCorr(alphaControls.lookup< label >("nAlphaCorr"))
const bool alphaApplyPrevCorr(alphaControls.lookupOrDefault< Switch >("alphaApplyPrevCorr", false))
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:138
static word groupName(Name name, const word &group)
Input token stream.
Definition: ITstream.H:53
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
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1751
ITstream & div(const word &name) const
Definition: fvSchemes.C:479
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
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:212
void alphaPredictor()
Solve for the phase-fractions.
virtual tmp< surfaceScalarField > alphaPhi(const surfaceScalarField &phi, const volScalarField &alpha, const dictionary &alphaControls)
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:181
A class for handling words, derived from string.
Definition: word.H:62
#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 explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiCorr, 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:257
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
bool isType(const Type &t)
Check the typeid.
Definition: typeInfo.H:158
messageStream Info
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const wordHashSet compressionSchemes
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
error FatalError
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.