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  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);
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 = 1.0 - alpha1;
229 
230  correctInterface();
231  }
232 
233  for (int aCorr=0; aCorr<nAlphaCorr; aCorr++)
234  {
235  // Split operator
236  tmp<surfaceScalarField> talphaPhi1Un
237  (
238  alphaPhi
239  (
240  phiCN(),
241  (cnCoeff*alpha1 + (1.0 - cnCoeff)*alpha1.oldTime())(),
243  )
244  );
245 
246  if (MULESCorr)
247  {
248  tmp<surfaceScalarField> talphaPhi1Corr(talphaPhi1Un() - alphaPhi1);
249  volScalarField alpha10("alpha10", alpha1);
250 
251  if (divU.valid())
252  {
254  (
255  geometricOneField(),
256  alpha1,
257  talphaPhi1Un(),
258  talphaPhi1Corr.ref(),
259  Sp(),
260  (-Sp()*alpha1)(),
261  oneField(),
262  zeroField()
263  );
264  }
265  else
266  {
268  (
269  geometricOneField(),
270  alpha1,
271  talphaPhi1Un(),
272  talphaPhi1Corr.ref(),
273  oneField(),
274  zeroField()
275  );
276  }
277 
278  // Under-relax the correction for all but the 1st corrector
279  if (aCorr == 0)
280  {
281  alphaPhi1 += talphaPhi1Corr();
282  }
283  else
284  {
285  alpha1 = 0.5*alpha1 + 0.5*alpha10;
286  alphaPhi1 += 0.5*talphaPhi1Corr();
287  }
288  }
289  else
290  {
291  alphaPhi1 = talphaPhi1Un;
292 
293  if (divU.valid())
294  {
296  (
297  geometricOneField(),
298  alpha1,
299  phiCN,
300  alphaPhi1,
301  Sp(),
302  (Su() + divU()*min(alpha1(), scalar(1)))(),
303  oneField(),
304  zeroField()
305  );
306  }
307  else
308  {
310  (
311  geometricOneField(),
312  alpha1,
313  phiCN,
314  alphaPhi1,
315  oneField(),
316  zeroField()
317  );
318  }
319  }
320 
321  alpha2 = 1.0 - alpha1;
322 
323  // Correct only the mixture interface for the interface compression flux
324  correctInterface();
325  }
326 
328  {
329  talphaPhi1Corr0 = alphaPhi1 - talphaPhi1Corr0;
330 
331  // Register alphaPhiCorr0.<phase1> for redistribution
332  talphaPhi1Corr0.ref().rename
333  (
334  IOobject::groupName("alphaPhiCorr0", alpha1.group())
335  );
336  talphaPhi1Corr0.ref().checkIn();
337  }
338  else
339  {
340  talphaPhi1Corr0.clear();
341  }
342 
343  if
344  (
345  word(mesh.schemes().ddt("ddt(rho,U)"))
346  != fv::EulerDdtScheme<vector>::typeName
347  && word(mesh.schemes().ddt("ddt(rho,U)"))
348  != fv::localEulerDdtScheme<vector>::typeName
349  )
350  {
351  if (ocCoeff > 0)
352  {
353  // Calculate the end-of-time-step alpha flux
354  alphaPhi1 =
355  (alphaPhi1 - (1.0 - cnCoeff)*alphaPhi1.oldTime())/cnCoeff;
356  }
357  }
358 
359  Info<< "Phase-1 volume fraction = "
360  << alpha1.weightedAverage(mesh.Vsc()).value()
361  << " Min(" << alpha1.name() << ") = " << min(alpha1).value()
362  << " Max(" << alpha1.name() << ") = " << max(alpha1).value()
363  << endl;
364 }
365 
366 
368 {
369  const dictionary& alphaControls = mesh.solution().solverDict(alpha1.name());
370 
371  const label nAlphaSubCycles(alphaControls.lookup<label>("nAlphaSubCycles"));
372 
373 
374  if (nAlphaSubCycles > 1)
375  {
376  dimensionedScalar totalDeltaT = runTime.deltaT();
377  tmp<volScalarField> trSubDeltaT;
378 
379  if (LTS)
380  {
381  trSubDeltaT =
383  }
384 
385  // Create a temporary alphaPhi1 to accumulate the sub-cycled alphaPhi1
386  tmp<surfaceScalarField> talphaPhi1
387  (
389  (
390  "alphaPhi1",
391  mesh,
392  dimensionedScalar(alphaPhi1.dimensions(), 0)
393  )
394  );
395 
396  List<volScalarField*> alphaPtrs({&alpha1, &alpha2});
397 
398  for
399  (
401  (
402  alphaPtrs,
404  );
405  !(++alphaSubCycle).end();
406  )
407  {
408  alphaSolve(alphaControls);
409  talphaPhi1.ref() += (runTime.deltaT()/totalDeltaT)*alphaPhi1;
410  }
411 
412  alphaPhi1 = talphaPhi1();
413  }
414  else
415  {
416  alphaSolve(alphaControls);
417  }
418 }
419 
420 
421 // ************************************************************************* //
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 >> &)
Return a temporary field constructed from name,.
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:113
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:160
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1748
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:67
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:94
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:306
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
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
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
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
bool isType(const Type &t)
Check the typeid.
Definition: typeInfo.H:131
messageStream Info
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const wordHashSet compressionSchemes
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
error FatalError