phaseSystemSolve.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) 2013-2020 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 "phaseSystem.H"
27 
28 #include "MULES.H"
29 #include "subCycle.H"
30 
31 #include "fvcDdt.H"
32 #include "fvcDiv.H"
33 #include "fvcSnGrad.H"
34 #include "fvcFlux.H"
35 #include "fvcMeshPhi.H"
36 #include "fvcSup.H"
37 
38 #include "fvmDdt.H"
39 #include "fvmLaplacian.H"
40 #include "fvmSup.H"
41 
42 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
43 
45 (
46  const PtrList<volScalarField>& rAUs,
47  const PtrList<surfaceScalarField>& rAUfs
48 )
49 {
50  const dictionary& alphaControls = mesh_.solverDict("alpha");
51 
52  const label nAlphaSubCycles(alphaControls.lookup<label>("nAlphaSubCycles"));
53  const label nAlphaCorr(alphaControls.lookup<label>("nAlphaCorr"));
54 
55  const bool LTS = fv::localEulerDdt::enabled(mesh_);
56 
57  // Optional reference phase which is not solved for
58  // but obtained from the sum of the other phases
59  phaseModel* referencePhasePtr = nullptr;
60 
61  // The phases which are solved
62  // i.e. the moving phases less the optional reference phase
63  phaseModelPartialList solvePhases;
64 
66  {
67  referencePhasePtr = &phases()[referencePhaseName_];
68 
69  solvePhases.setSize(movingPhases().size() - 1);
70  label solvePhasesi = 0;
71  forAll(movingPhases(), movingPhasei)
72  {
73  if (&movingPhases()[movingPhasei] != referencePhasePtr)
74  {
75  solvePhases.set(solvePhasesi++, &movingPhases()[movingPhasei]);
76  }
77  }
78  }
79  else
80  {
81  solvePhases = movingPhases();
82  }
83 
84  // The phases included in the flux sum limit
85  // which is all moving phases if the number of solved phases is > 1
86  // otherwise it is just the solved phases
87  // as the flux sum limit is not needed in this case
88  phaseModelPartialList fluxPhases;
89  if (solvePhases.size() == 1)
90  {
91  fluxPhases = solvePhases;
92  }
93  else
94  {
95  fluxPhases = movingPhases();
96  }
97 
98  forAll(phases(), phasei)
99  {
100  phases()[phasei].correctBoundaryConditions();
101  }
102 
103  PtrList<surfaceScalarField> alphaPhiDbyA0s(phases().size());
104  if (implicitPhasePressure() && (rAUs.size() || rAUfs.size()))
105  {
106  const PtrList<surfaceScalarField> DByAfs(this->DByAfs(rAUs, rAUfs));
107 
108  forAll(solvePhases, solvePhasei)
109  {
110  phaseModel& phase = solvePhases[solvePhasei];
111  volScalarField& alpha = phase;
112 
113  alphaPhiDbyA0s.set
114  (
115  phase.index(),
116  DByAfs[phase.index()]
117  *fvc::snGrad(alpha, "bounded")*mesh_.magSf()
118  );
119  }
120  }
121 
122  // Calculate the void fraction
123  volScalarField alphaVoid
124  (
125  IOobject
126  (
127  "alphaVoid",
128  mesh_.time().timeName(),
129  mesh_
130  ),
131  mesh_,
133  );
134  forAll(stationaryPhases(), stationaryPhasei)
135  {
136  alphaVoid -= stationaryPhases()[stationaryPhasei];
137  }
138 
139  bool dilatation = false;
140  forAll(fluxPhases, fluxPhasei)
141  {
142  if (fluxPhases[fluxPhasei].divU().valid())
143  {
144  dilatation = true;
145  break;
146  }
147  }
148 
149  for (int acorr=0; acorr<nAlphaCorr; acorr++)
150  {
151  PtrList<volScalarField::Internal> Sps(phases().size());
152  PtrList<volScalarField::Internal> Sus(phases().size());
153 
154  forAll(fluxPhases, fluxPhasei)
155  {
156  phaseModel& phase = fluxPhases[fluxPhasei];
157  volScalarField& alpha = phase;
158  const label phasei = phase.index();
159 
160  Sps.set
161  (
162  phasei,
164  (
165  IOobject
166  (
167  "Sp",
168  mesh_.time().timeName(),
169  mesh_
170  ),
171  mesh_,
173  )
174  );
175 
176  Sus.set
177  (
178  phasei,
180  (
181  "Su",
182  min(alpha, scalar(1))
183  *fvc::div(fvc::absolute(phi_, phase.U()))
184  )
185  );
186 
187  if (dilatation)
188  {
189  // Construct the dilatation rate source term
191  (
193  (
194  "dgdt",
195  mesh_,
197  )
198  );
199 
200  forAll(phases(), phasej)
201  {
202  const phaseModel& phase2 = phases()[phasej];
203  const volScalarField& alpha2 = phase2;
204 
205  if (&phase2 != &phase)
206  {
207  if (phase.divU().valid())
208  {
209  dgdt += alpha2()*phase.divU()()();
210  }
211 
212  if (phase2.divU().valid())
213  {
214  dgdt -= alpha()*phase2.divU()()();
215  }
216  }
217  }
218 
219  volScalarField::Internal& Sp = Sps[phasei];
220  volScalarField::Internal& Su = Sus[phasei];
221 
222  forAll(dgdt, celli)
223  {
224  if (dgdt[celli] > 0)
225  {
226  Sp[celli] -= dgdt[celli]/max(1 - alpha[celli], 1e-4);
227  Su[celli] += dgdt[celli]/max(1 - alpha[celli], 1e-4);
228  }
229  else if (dgdt[celli] < 0)
230  {
231  Sp[celli] += dgdt[celli]/max(alpha[celli], 1e-4);
232  }
233  }
234  }
235  }
236 
237  tmp<volScalarField> trSubDeltaT;
238 
239  if (LTS && nAlphaSubCycles > 1)
240  {
241  trSubDeltaT =
243  }
244 
245  List<volScalarField*> alphaPtrs(phases().size());
246  forAll(phases(), phasei)
247  {
248  alphaPtrs[phasei] = &phases()[phasei];
249  }
250 
251  for
252  (
253  subCycle<volScalarField, subCycleFields> alphaSubCycle
254  (
255  alphaPtrs,
257  );
258  !(++alphaSubCycle).end();
259  )
260  {
261  // Generate face-alphas
262  PtrList<surfaceScalarField> alphafs(phases().size());
263  if (solvePhases.size() > 1)
264  {
265  forAll(phases(), phasei)
266  {
267  phaseModel& phase = phases()[phasei];
268  alphafs.set
269  (
270  phasei,
272  (
273  IOobject::groupName("alphaf", phase.name()),
274  upwind<scalar>(mesh_, phi_).interpolate(phase)
275  )
276  );
277  }
278  }
279 
280  // Create correction fluxes
281  PtrList<surfaceScalarField> alphaPhiCorrs(phases().size());
282 
283  if (solvePhases.size() > 1)
284  {
285  forAll(stationaryPhases(), stationaryPhasei)
286  {
287  phaseModel& phase = stationaryPhases()[stationaryPhasei];
288 
289  alphaPhiCorrs.set
290  (
291  phase.index(),
293  (
294  IOobject::groupName("alphaPhiCorr", phase.name()),
295  - upwind<scalar>(mesh_, phi_).flux(phase)
296  )
297  );
298  }
299  }
300 
301  forAll(fluxPhases, fluxPhasei)
302  {
303  phaseModel& phase = fluxPhases[fluxPhasei];
304  volScalarField& alpha = phase;
305 
306  alphaPhiCorrs.set
307  (
308  phase.index(),
310  (
311  IOobject::groupName("alphaPhiCorr", phase.name()),
312  fvc::flux(phi_, alpha, "div(phi," + alpha.name() + ')')
313  )
314  );
315 
316  surfaceScalarField& alphaPhiCorr = alphaPhiCorrs[phase.index()];
317 
318  forAll(phases(), phasei)
319  {
320  phaseModel& phase2 = phases()[phasei];
321  volScalarField& alpha2 = phase2;
322 
323  if (&phase2 == &phase) continue;
324 
325  surfaceScalarField phir(phase.phi() - phase2.phi());
326 
328  (
329  cAlphas_.find(phasePairKey(phase.name(), phase2.name()))
330  );
331 
332  if (cAlpha != cAlphas_.end())
333  {
334  surfaceScalarField phic
335  (
336  (mag(phi_) + mag(phir))/mesh_.magSf()
337  );
338 
339  phir +=
340  min(cAlpha()*phic, max(phic))
341  *nHatf(alpha, alpha2);
342  }
343 
344  word phirScheme
345  (
346  "div(phir," + alpha2.name() + ',' + alpha.name() + ')'
347  );
348 
349  alphaPhiCorr += fvc::flux
350  (
351  -fvc::flux(-phir, alpha2, phirScheme),
352  alpha,
353  phirScheme
354  );
355  }
356 
357  if (alphaPhiDbyA0s.set(phase.index()))
358  {
359  alphaPhiCorr +=
360  fvc::interpolate(max(alpha, scalar(0)))
361  *fvc::interpolate(max(1 - alpha, scalar(0)))
362  *alphaPhiDbyA0s[phase.index()];
363  }
364 
365  phase.correctInflowOutflow(alphaPhiCorr);
366 
368  (
369  geometricOneField(),
370  alpha,
371  phi_,
372  alphaPhiCorr,
373  Sps[phase.index()],
374  Sus[phase.index()],
375  min(alphaVoid.primitiveField(), phase.alphaMax())(),
376  zeroField(),
377  true
378  );
379  }
380 
381  if (solvePhases.size() > 1)
382  {
383  // Limit the flux sums, fixing those of the stationary phases
384  labelHashSet fixedAlphaPhiCorrs;
385  forAll(stationaryPhases(), stationaryPhasei)
386  {
387  fixedAlphaPhiCorrs.insert
388  (
389  stationaryPhases()[stationaryPhasei].index()
390  );
391  }
392  MULES::limitSum(alphafs, alphaPhiCorrs, fixedAlphaPhiCorrs);
393  }
394 
395  forAll(solvePhases, solvePhasei)
396  {
397  phaseModel& phase = solvePhases[solvePhasei];
398  volScalarField& alpha = phase;
399 
400  surfaceScalarField& alphaPhi = alphaPhiCorrs[phase.index()];
401  alphaPhi += upwind<scalar>(mesh_, phi_).flux(phase);
402  phase.correctInflowOutflow(alphaPhi);
403 
405  (
406  geometricOneField(),
407  alpha,
408  alphaPhi,
409  Sps[phase.index()],
410  Sus[phase.index()]
411  );
412 
413  if (alphaSubCycle.index() == 1)
414  {
415  phase.alphaPhiRef() = alphaPhi;
416  }
417  else
418  {
419  phase.alphaPhiRef() += alphaPhi;
420  }
421  }
422 
423  if (implicitPhasePressure() && (rAUs.size() || rAUfs.size()))
424  {
425  const PtrList<surfaceScalarField> DByAfs
426  (
427  this->DByAfs(rAUs, rAUfs)
428  );
429 
430  forAll(solvePhases, solvePhasei)
431  {
432  phaseModel& phase = solvePhases[solvePhasei];
433  volScalarField& alpha = phase;
434 
435  const surfaceScalarField alphaDbyA
436  (
437  fvc::interpolate(max(alpha, scalar(0)))
438  *fvc::interpolate(max(1 - alpha, scalar(0)))
439  *DByAfs[phase.index()]
440  );
441 
442  fvScalarMatrix alphaEqn
443  (
444  fvm::ddt(alpha) - fvc::ddt(alpha)
445  - fvm::laplacian(alphaDbyA, alpha, "bounded")
446  );
447 
448  alphaEqn.solve();
449 
450  phase.alphaPhiRef() += alphaEqn.flux();
451  }
452  }
453 
454  // Report the phase fractions and the phase fraction sum
455  forAll(solvePhases, solvePhasei)
456  {
457  phaseModel& phase = solvePhases[solvePhasei];
458 
459  Info<< phase.name() << " fraction, min, max = "
460  << phase.weightedAverage(mesh_.V()).value()
461  << ' ' << min(phase).value()
462  << ' ' << max(phase).value()
463  << endl;
464  }
465 
466  if (!referencePhasePtr)
467  {
469  (
470  IOobject
471  (
472  "sumAlphaMoving",
473  mesh_.time().timeName(),
474  mesh_
475  ),
476  mesh_,
478  );
479  forAll(movingPhases(), movingPhasei)
480  {
481  sumAlphaMoving += movingPhases()[movingPhasei];
482  }
483 
484  Info<< "Phase-sum volume fraction, min, max = "
485  << (sumAlphaMoving + 1 - alphaVoid)()
486  .weightedAverage(mesh_.V()).value()
487  << ' ' << min(sumAlphaMoving + 1 - alphaVoid).value()
488  << ' ' << max(sumAlphaMoving + 1 - alphaVoid).value()
489  << endl;
490 
491  // Correct the sum of the phase fractions to avoid drift
492  forAll(movingPhases(), movingPhasei)
493  {
494  movingPhases()[movingPhasei] *= alphaVoid/sumAlphaMoving;
495  }
496  }
497  }
498 
499  if (nAlphaSubCycles > 1)
500  {
501  forAll(solvePhases, solvePhasei)
502  {
503  phaseModel& phase = solvePhases[solvePhasei];
504 
505  phase.alphaPhiRef() /= nAlphaSubCycles;
506  }
507  }
508 
509  forAll(solvePhases, solvePhasei)
510  {
511  phaseModel& phase = solvePhases[solvePhasei];
512 
513  phase.alphaRhoPhiRef() =
514  fvc::interpolate(phase.rho())*phase.alphaPhi();
515 
516  phase.maxMin(0, 1);
517  }
518 
519  if (referencePhasePtr)
520  {
521  phaseModel& referencePhase = *referencePhasePtr;
522 
523  referencePhase.alphaPhiRef() = phi_;
524 
525  forAll(solvePhases, solvePhasei)
526  {
527  phaseModel& phase = solvePhases[solvePhasei];
528  referencePhase.alphaPhiRef() -= phase.alphaPhi();
529  }
530 
531  referencePhase.correctInflowOutflow(referencePhase.alphaPhiRef());
532  referencePhase.alphaRhoPhiRef() =
533  fvc::interpolate(referencePhase.rho())
534  *referencePhase.alphaPhi();
535 
536  volScalarField& referenceAlpha = referencePhase;
537  referenceAlpha = alphaVoid;
538 
539  forAll(solvePhases, solvePhasei)
540  {
541  referenceAlpha -= solvePhases[solvePhasei];
542  }
543  }
544  }
545 }
546 
547 
548 // ************************************************************************* //
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)
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const phaseModelPartialList & movingPhases() const
Return the models for phases that are moving.
Definition: phaseSystemI.H:49
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
const word & name() const
Return name.
Definition: IOobject.H:303
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &)
Return a temporary field constructed from name, mesh.
label phasei
Definition: pEqn.H:27
tmp< surfaceScalarField > nHatf(const volScalarField &alpha1, const volScalarField &alpha2) const
Normal to interface between two phases dotted with face areas.
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:112
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
PtrList< surfaceScalarField > alphafs(phases.size())
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
static tmp< volScalarField > localRSubDeltaT(const fvMesh &mesh, const label nAlphaSubCycles)
Calculate and return the reciprocal of the local sub-cycling.
Definition: localEulerDdt.C:70
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:177
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Calculate the matrix for the laplacian of the field.
Calculate the snGrad of the given volField.
surfaceScalarField phir(fvc::flux(UdmModel.Udm()))
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
alpha2
Definition: alphaEqn.H:115
const phaseModelPartialList & stationaryPhases() const
Return the models for phases that are stationary.
Definition: phaseSystemI.H:63
virtual void solve(const PtrList< volScalarField > &rAUs, const PtrList< surfaceScalarField > &rAUfs)
Solve for the phase fractions.
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:142
dictionary()
Construct top-level dictionary null.
Definition: dictionary.C:464
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:211
Calculate the first temporal derivative.
const phaseModelList & phases() const
Return the phase models.
Definition: phaseSystemI.H:35
static bool enabled(const fvMesh &mesh)
Return true if LTS is enabled.
Definition: localEulerDdt.C:37
void limitSum(UPtrList< scalarField > &phiPsiCorrs)
Definition: MULES.C:30
Calculate the face-flux of the given field.
static word groupName(Name name, const word &group)
friend class const_iterator
Declare friendship with the const_iterator.
Definition: HashTable.H:197
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
Calculate the matrix for the first temporal derivative.
word referencePhaseName_
Name of optional reference phase which is not solved for.
Definition: phaseSystem.H:145
static const word null
An empty word.
Definition: word.H:77
UPtrList< phaseModel > phaseModelPartialList
Definition: phaseSystem.H:83
surfaceScalarField alphaPhi(phi.name()+alpha1.name(), fvc::flux(phi, alpha1, alphaScheme))
const fvMesh & mesh_
Reference to the mesh.
Definition: phaseSystem.H:141
const iterator & end()
Definition: UILList.H:223
Calculate the field for explicit evaluation of implicit and explicit sources.
cAlphaTable cAlphas_
Interface compression coefficients.
Definition: phaseSystem.H:178
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:131
Calculate the divergence of the given field.
bool set(const label) const
Is element set.
Definition: UPtrListI.H:78
tmp< volScalarField > sumAlphaMoving() const
Return the sum of the phase fractions of the moving phases.
const label nAlphaSubCycles(alphaControls.lookup< label >("nAlphaSubCycles"))
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Calculate the mesh motion flux and convert fluxes from absolute to relative and back.
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:188
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
virtual PtrList< surfaceScalarField > DByAfs(const PtrList< volScalarField > &rAUs, const PtrList< surfaceScalarField > &rAUfs) const =0
Return the phase diffusivity.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
const label nAlphaCorr(alphaControls.lookup< label >("nAlphaCorr"))
MULES: Multidimensional universal limiter for explicit solution.
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
surfaceScalarField phi_
Total volumetric flux.
Definition: phaseSystem.H:166
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
zeroField divU
Definition: alphaSuSp.H:3
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Definition: fvcFlux.C:32
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
virtual bool implicitPhasePressure() const
Returns true if the phase pressure is treated implicitly.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
Calculate the matrix for implicit and explicit sources.