multiphaseSystem.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2017-2020 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "multiphaseSystem.H"
29 
31 #include "Time.H"
32 #include "subCycle.H"
33 #include "fvcMeshPhi.H"
34 
35 #include "surfaceInterpolate.H"
36 #include "fvcGrad.H"
37 #include "fvcSnGrad.H"
38 #include "fvcDiv.H"
39 #include "fvcDdt.H"
40 #include "fvcFlux.H"
41 #include "fvmDdt.H"
42 #include "fvcAverage.H"
43 #include "fvMatrix.H"
44 #include "fvmSup.H"
45 #include "CMULES.H"
46 
47 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51  defineTypeNameAndDebug(multiphaseSystem, 0);
52  defineRunTimeSelectionTable(multiphaseSystem, dictionary);
53 }
54 
55 
56 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
57 
59 (
60  const fvMesh& mesh
61 )
62 :
64  cAlphas_(),
65  ddtAlphaMax_(0.0),
66  limitedPhiAlphas_(phaseModels_.size()),
67  Su_(phaseModels_.size()),
68  Sp_(phaseModels_.size())
69 {
70  label phasei = 0;
71  phases_.setSize(phaseModels_.size());
72  forAllIters(phaseModels_, iter)
73  {
74  phaseModel& pm = iter()();
75  phases_.set(phasei++, &pm);
76  }
77 
78  mesh.solverDict("alpha").readEntry("cAlphas", cAlphas_);
79 
80  // Initiate Su and Sp
81  forAllConstIters(phaseModels_, iter)
82  {
83  const phaseModel& pm = iter()();
84 
85  Su_.insert
86  (
87  pm.name(),
89  (
90  IOobject
91  (
92  "Su" + pm.name(),
93  mesh_.time().timeName(),
94  mesh_
95  ),
96  mesh_,
98  )
99  );
100 
101  Sp_.insert
102  (
103  pm.name(),
105  (
106  IOobject
107  (
108  "Sp" + pm.name(),
109  mesh_.time().timeName(),
110  mesh_
111  ),
112  mesh_,
114  )
115  );
116  }
117 }
118 
119 
120 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
121 
123 {
124  this->alphaTransfer(Su_, Sp_);
125 }
126 
127 
129 {
130  const dictionary& alphaControls = mesh_.solverDict("alpha");
131  label nAlphaSubCycles(alphaControls.get<label>("nAlphaSubCycles"));
132 
133  volScalarField& alpha = phases_.first();
134 
135  if (nAlphaSubCycles > 1)
136  {
137  surfaceScalarField rhoPhiSum
138  (
139  IOobject
140  (
141  "rhoPhiSum",
142  mesh_.time().timeName(),
143  mesh_
144  ),
145  mesh_,
146  dimensionedScalar(rhoPhi_.dimensions(), Zero)
147  );
148 
149  dimensionedScalar totalDeltaT = mesh_.time().deltaT();
150 
151  for
152  (
154  !(++alphaSubCycle).end();
155  )
156  {
157  solveAlphas();
158  rhoPhiSum += (mesh_.time().deltaT()/totalDeltaT)*rhoPhi_;
159  }
160 
161  rhoPhi_ = rhoPhiSum;
162  }
163  else
164  {
165  solveAlphas();
166  }
167 
168 }
169 
170 void Foam::multiphaseSystem::solveAlphas()
171 {
172 
173  const dictionary& alphaControls = mesh_.solverDict("alpha");
174  alphaControls.readEntry("cAlphas", cAlphas_);
175  label nAlphaCorr(alphaControls.get<label>("nAlphaCorr"));
176 
177  PtrList<surfaceScalarField> phiAlphaCorrs(phases_.size());
178 
179  const surfaceScalarField& phi = this->phi();
180 
181  surfaceScalarField phic(mag((phi)/mesh_.magSf()));
182 
183  // Do not compress interface at non-coupled boundary faces
184  // (inlets, outlets etc.)
185  surfaceScalarField::Boundary& phicBf = phic.boundaryFieldRef();
186  forAll(phic.boundaryField(), patchi)
187  {
188  fvsPatchScalarField& phicp = phicBf[patchi];
189 
190  if (!phicp.coupled())
191  {
192  phicp == 0;
193  }
194  }
195 
196  for (int acorr=0; acorr<nAlphaCorr; acorr++)
197  {
198  label phasei = 0;
199  for (phaseModel& phase1 : phases_)
200  {
201  const volScalarField& alpha1 = phase1;
202 
203  phiAlphaCorrs.set
204  (
205  phasei,
207  (
208  "phi" + alpha1.name() + "Corr",
209  fvc::flux
210  (
211  phi,
212  alpha1,
213  "div(phi," + alpha1.name() + ')'
214  )
215  )
216  );
217 
218  surfaceScalarField& phiAlphaCorr = phiAlphaCorrs[phasei];
219 
220  for (phaseModel& phase2 : phases_)
221  {
222  const volScalarField& alpha2 = phase2;
223 
224  if (&phase2 == &phase1) continue;
225 
226  const phasePairKey key12(phase1.name(), phase2.name());
227 
228  if (!cAlphas_.found(key12))
229  {
231  << "Phase compression factor (cAlpha) not found for : "
232  << key12
233  << exit(FatalError);
234  }
235  scalar cAlpha = cAlphas_.find(key12)();
236 
237  phic = min(cAlpha*phic, max(phic));
238 
240 
241  word phirScheme
242  (
243  "div(phir," + alpha2.name() + ',' + alpha1.name() + ')'
244  );
245 
246  phiAlphaCorr += fvc::flux
247  (
248  -fvc::flux(-phir, alpha2, phirScheme),
249  alpha1,
250  phirScheme
251  );
252  }
253 
254  // Ensure that the flux at inflow BCs is preserved
255  forAll(phiAlphaCorr.boundaryField(), patchi)
256  {
257  fvsPatchScalarField& phiAlphaCorrp =
258  phiAlphaCorr.boundaryFieldRef()[patchi];
259 
260  if (!phiAlphaCorrp.coupled())
261  {
262  const scalarField& phi1p = phi.boundaryField()[patchi];
263  const scalarField& alpha1p =
264  alpha1.boundaryField()[patchi];
265 
266  forAll(phiAlphaCorrp, facei)
267  {
268  if (phi1p[facei] < 0)
269  {
270  phiAlphaCorrp[facei] = alpha1p[facei]*phi1p[facei];
271  }
272  }
273  }
274  }
275 
276  ++phasei;
277  }
278 
279  // Set Su and Sp to zero
280  for (const phaseModel& phase : phases_)
281  {
282  Su_[phase.name()] = dimensionedScalar("Su", dimless/dimTime, Zero);
283  Sp_[phase.name()] = dimensionedScalar("Sp", dimless/dimTime, Zero);
284 
285  // Add alpha*div(U)
286  //const volScalarField& alpha = phase;
287  //Su_[phase.name()] +=
288  // fvc::div(phi)*min(max(alpha, scalar(0)), scalar(1));
289  }
290 
291  // Fill Su and Sp
292  calculateSuSp();
293 
294  // Limit phiAlphaCorr on each phase
295  phasei = 0;
296  for (phaseModel& phase : phases_)
297  {
298  volScalarField& alpha1 = phase;
299 
300  surfaceScalarField& phiAlphaCorr = phiAlphaCorrs[phasei];
301 
302  volScalarField::Internal& Su = Su_[phase.name()];
303  volScalarField::Internal& Sp = Sp_[phase.name()];
304 
306  (
307  1.0/mesh_.time().deltaT().value(),
308  geometricOneField(),
309  alpha1,
310  phi,
311  phiAlphaCorr,
312  Sp,
313  Su,
314  oneField(),
315  zeroField(),
316  true
317  );
318  ++phasei;
319  }
320 
321  MULES::limitSum(phiAlphaCorrs);
322 
323  volScalarField sumAlpha
324  (
325  IOobject
326  (
327  "sumAlpha",
328  mesh_.time().timeName(),
329  mesh_
330  ),
331  mesh_,
333  );
334 
335  phasei = 0;
336  for (phaseModel& phase : phases_)
337  {
338  volScalarField& alpha1 = phase;
339 
340  const volScalarField::Internal& Su = Su_[phase.name()];
341 
342  const volScalarField::Internal& Sp = Sp_[phase.name()];
343 
344  surfaceScalarField& phiAlpha = phiAlphaCorrs[phasei];
345 
346  // Add a bounded upwind U-mean flux
347  //phiAlpha += upwind<scalar>(mesh_, phi).flux(alpha1);
348  fvScalarMatrix alpha1Eqn
349  (
350  fv::EulerDdtScheme<scalar>(mesh_).fvmDdt(alpha1)
351  + fv::gaussConvectionScheme<scalar>
352  (
353  mesh_,
354  phi,
355  upwind<scalar>(mesh_, phi)
356  ).fvmDiv(phi, alpha1)
357  ==
358  Su + fvm::Sp(Sp, alpha1)
359  );
360 
361  alpha1Eqn.boundaryManipulate(alpha1.boundaryFieldRef());
362 
363  alpha1Eqn.solve();
364 
365  phiAlpha += alpha1Eqn.flux();
366 
368  (
369  geometricOneField(),
370  alpha1,
371  phi,
372  phiAlpha,
373  Sp,
374  Su,
375  oneField(),
376  zeroField()
377  );
378 
379  phase.alphaPhi() = phiAlpha;
380 
381  ++phasei;
382  }
383 
384  if (acorr == nAlphaCorr - 1)
385  {
386  volScalarField sumAlpha
387  (
388  IOobject
389  (
390  "sumAlpha",
391  mesh_.time().timeName(),
392  mesh_
393  ),
394  mesh_,
396  );
397 
398  // Reset rhoPhi
399  rhoPhi_ = dimensionedScalar("rhoPhi", dimMass/dimTime, Zero);
400 
401  for (phaseModel& phase : phases_)
402  {
403  volScalarField& alpha1 = phase;
404  sumAlpha += alpha1;
405 
406  // Update rhoPhi
407  rhoPhi_ += fvc::interpolate(phase.rho()) * phase.alphaPhi();
408  }
409 
410  Info<< "Phase-sum volume fraction, min, max = "
411  << sumAlpha.weightedAverage(mesh_.V()).value()
412  << ' ' << min(sumAlpha).value()
413  << ' ' << max(sumAlpha).value()
414  << endl;
415 
416  volScalarField sumCorr(1.0 - sumAlpha);
417 
418  for (phaseModel& phase : phases_)
419  {
420  volScalarField& alpha = phase;
421  alpha += alpha*sumCorr;
422 
423  Info<< alpha.name() << " volume fraction = "
424  << alpha.weightedAverage(mesh_.V()).value()
425  << " Min(alpha) = " << min(alpha).value()
426  << " Max(alpha) = " << max(alpha).value()
427  << endl;
428  }
429  }
430  }
431 }
432 
433 
435 {
436  return phases_;
437 }
438 
439 
441 {
442  return phases_;
443 }
444 
445 
447 {
448  return phases_[i];
449 }
450 
451 
453 {
454  return phases_[i];
455 }
456 
457 
459 {
460  return ddtAlphaMax_;
461 }
462 
463 
465 {
466  auto iter = phaseModels_.cbegin();
467 
468  scalar maxVal = max(iter()->diffNo()).value();
469 
470  for (++iter; iter != phaseModels_.cend(); ++iter)
471  {
472  maxVal = max(maxVal, max(iter()->diffNo()).value());
473  }
474 
475  return maxVal * mesh_.time().deltaT().value();
476 }
477 
478 
481 {
482  return limitedPhiAlphas_;
483 }
484 
485 
487 {
488  return Su_;
489 }
490 
491 
493 {
494  return Sp_;
495 }
496 
497 
499 {
500  return true;
501 }
502 
503 
504 // ************************************************************************* //
Foam::phaseModel
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:54
Foam::multiphaseSystem::limitedPhiAlphas
const compressionFluxTable & limitedPhiAlphas() const
Access to compression fluxes for phaes.
Definition: multiphaseSystem.C:480
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
phicBf
surfaceScalarField::Boundary & phicBf
Definition: alphaEqn.H:75
Foam::fvc::flux
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Foam::MULES::explicitSolve
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
Definition: MULESTemplates.C:41
CMULES.H
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
Foam::multiphaseSystem::calculateSuSp
void calculateSuSp()
Calculate Sp and Su.
Definition: multiphaseSystem.C:122
Foam::multiphaseSystem::solve
void solve()
Solve for the mixture phase-fractions.
Definition: multiphaseSystem.C:821
Foam::MULES::limit
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)
Definition: MULESTemplates.C:581
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
fvcMeshPhi.H
Calculate the mesh motion flux and convert fluxes from absolute to relative and back.
alpha2
const volScalarField & alpha2
Definition: setRegionFluidFields.H:9
subCycle.H
Foam::defineRunTimeSelectionTable
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
Foam::fvsPatchScalarField
fvsPatchField< scalar > fvsPatchScalarField
Definition: fvsPatchFieldsFwd.H:45
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
fvcSnGrad.H
Calculate the snGrad of the given volField.
Foam::vtk::Tools::zeroField
vtkSmartPointer< vtkFloatArray > zeroField(const word &name, const label size)
Create named field initialized to zero.
Definition: foamVtkToolsTemplates.C:327
Sp
zeroField Sp
Definition: alphaSuSp.H:2
fvcDiv.H
Calculate the divergence of the given field.
Foam::dimensioned::name
const word & name() const
Return const reference to name.
Definition: dimensionedType.C:406
phasei
label phasei
Definition: pEqn.H:27
fvMatrix.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::multiphaseSystem::phase
const phaseModel & phase(const label i) const
Constant access phase model i.
Definition: multiphaseSystem.C:446
alpha1
const volScalarField & alpha1
Definition: setRegionFluidFields.H:8
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:107
nAlphaSubCycles
label nAlphaSubCycles(alphaControls.get< label >("nAlphaSubCycles"))
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
phic
surfaceScalarField phic(mixture.cAlpha() *mag(alphaPhic/mesh.magSf()))
Foam::multiphaseSystem::Su_
SuSpTable Su_
Su phase source terms.
Definition: multiphaseSystem.H:81
Foam::multiphaseSystem::maxDiffNo
scalar maxDiffNo() const
Maximum diffusion number.
Definition: multiphaseSystem.C:464
Su
zeroField Su
Definition: alphaSuSp.H:1
nAlphaCorr
label nAlphaCorr(alphaControls.get< label >("nAlphaCorr"))
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
phir
surfaceScalarField phir(fvc::flux(UdmModel.Udm()))
Foam::fvScalarMatrix
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:44
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::GeometricField< scalar, fvPatchField, volMesh >::Internal
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
Definition: GeometricField.H:107
Foam::MULES::limitSum
void limitSum(UPtrList< scalarField > &phiPsiCorrs)
Definition: MULES.C:34
Foam::UPtrList< Foam::phaseModel >
Foam::dictionary::readEntry
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
Definition: dictionaryTemplates.C:302
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::multiphaseSystem::ddtAlphaMax
dimensionedScalar ddtAlphaMax() const
Access to ddtAlphaMax.
Definition: multiphaseSystem.C:458
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::multiphaseSystem::Sp
SuSpTable & Sp()
Access Sp.
Definition: multiphaseSystem.C:492
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::multiphaseSystem::Sp_
SuSpTable Sp_
Sp phase source terms.
Definition: multiphaseSystem.H:84
forAllIters
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:223
fixedValueFvsPatchFields.H
phase2
phaseModel & phase2
Definition: setRegionFluidFields.H:6
Foam::multiphaseSystem::read
bool read()
Read base transportProperties dictionary.
Definition: multiphaseSystem.C:916
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::dimensioned< scalar >
stdFoam::end
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:121
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::fvc::interpolate
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.
Foam::multiphaseSystem::phases
const PtrDictionary< phaseModel > & phases() const
Return the phases.
Definition: multiphaseSystem.H:215
fvmSup.H
Calculate the matrix for implicit and explicit sources.
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::HashTable< surfaceScalarField >
Foam::phaseModel::name
const word & name() const
Definition: phaseModel.H:140
Time.H
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:54
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Definition: GeometricField.C:783
fvcGrad.H
Calculate the gradient of the given field.
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
fvcFlux.H
Calculate the face-flux of the given field.
fvcDdt.H
Calculate the first temporal derivative.
Foam::phaseSystem
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:66
phase1
phaseModel & phase1
Definition: setRegionFluidFields.H:5
fvcAverage.H
Area-weighted average a surfaceField creating a volField.
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::multiphaseSystem::Su
SuSpTable & Su()
Access Su.
Definition: multiphaseSystem.C:486
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
multiphaseSystem.H
Foam::multiphaseSystem::multiphaseSystem
multiphaseSystem(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
Definition: multiphaseSystem.C:361
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
fvmDdt.H
Calculate the matrix for the first temporal derivative.
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189
Foam::subCycle
Perform a subCycleTime on a field.
Definition: subCycle.H:123
alphaControls
const dictionary & alphaControls
Definition: alphaControls.H:1