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-2022 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 namespace multiphaseInter
52 {
53  defineTypeNameAndDebug(multiphaseSystem, 0);
55 }
56 }
57 
58 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
59 
61 (
62  const fvMesh& mesh
63 )
64 :
66  cAlphas_(),
67  ddtAlphaMax_(0.0),
68  limitedPhiAlphas_(phaseModels_.size()),
69  Su_(phaseModels_.size()),
70  Sp_(phaseModels_.size())
71 {
72  label phasei = 0;
75  {
76  multiphaseInter::phaseModel& pm = iter()();
77  phases_.set(phasei++, &pm);
78  }
79 
80  mesh.solverDict("alpha").readEntry("cAlphas", cAlphas_);
81 
82  // Initiate Su and Sp
84  {
85  const multiphaseInter::phaseModel& pm = iter()();
86 
87  Su_.insert
88  (
89  pm.name(),
91  (
92  IOobject
93  (
94  "Su" + pm.name(),
95  mesh_.time().timeName(),
96  mesh_
97  ),
98  mesh_,
100  )
101  );
102 
103  Sp_.insert
104  (
105  pm.name(),
107  (
108  IOobject
109  (
110  "Sp" + pm.name(),
111  mesh_.time().timeName(),
112  mesh_
113  ),
114  mesh_,
116  )
117  );
118  }
119 }
120 
121 
122 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
125 {
126  this->alphaTransfer(Su_, Sp_);
127 }
128 
129 
131 {
132  const dictionary& alphaControls = mesh_.solverDict("alpha");
133  label nAlphaSubCycles(alphaControls.get<label>("nAlphaSubCycles"));
134 
135  volScalarField& alpha = phases_.first();
136 
137  if (nAlphaSubCycles > 1)
138  {
139  surfaceScalarField rhoPhiSum
140  (
141  IOobject
142  (
143  "rhoPhiSum",
144  mesh_.time().timeName(),
145  mesh_
146  ),
147  mesh_,
148  dimensionedScalar(rhoPhi_.dimensions(), Zero)
149  );
150 
151  dimensionedScalar totalDeltaT = mesh_.time().deltaT();
152 
153  for
154  (
156  !(++alphaSubCycle).end();
157  )
158  {
159  solveAlphas();
160  rhoPhiSum += (mesh_.time().deltaT()/totalDeltaT)*rhoPhi_;
161  }
162 
163  rhoPhi_ = rhoPhiSum;
164  }
165  else
166  {
167  solveAlphas();
168  }
169 
170 }
171 
173 {
174 
175  const dictionary& alphaControls = mesh_.solverDict("alpha");
176  alphaControls.readEntry("cAlphas", cAlphas_);
177  label nAlphaCorr(alphaControls.get<label>("nAlphaCorr"));
178 
179  PtrList<surfaceScalarField> phiAlphaCorrs(phases_.size());
180 
181  const surfaceScalarField& phi = this->phi();
182 
183  surfaceScalarField phic(mag((phi)/mesh_.magSf()));
184 
185  // Do not compress interface at non-coupled boundary faces
186  // (inlets, outlets etc.)
187  surfaceScalarField::Boundary& phicBf = phic.boundaryFieldRef();
188  forAll(phic.boundaryField(), patchi)
189  {
190  fvsPatchScalarField& phicp = phicBf[patchi];
191 
192  if (!phicp.coupled())
193  {
194  phicp == 0;
195  }
196  }
197 
198  for (int acorr=0; acorr<nAlphaCorr; acorr++)
199  {
200  label phasei = 0;
201  for (multiphaseInter::phaseModel& phase1 : phases_)
202  {
203  const volScalarField& alpha1 = phase1;
204 
205  phiAlphaCorrs.set
206  (
207  phasei,
209  (
210  "phi" + alpha1.name() + "Corr",
211  fvc::flux
212  (
213  phi,
214  alpha1,
215  "div(phi," + alpha1.name() + ')'
216  )
217  )
218  );
219 
220  surfaceScalarField& phiAlphaCorr = phiAlphaCorrs[phasei];
221 
222  for (multiphaseInter::phaseModel& phase2 : phases_)
223  {
224  const volScalarField& alpha2 = phase2;
225 
226  if (&phase2 == &phase1) continue;
227 
228  const phasePairKey key12(phase1.name(), phase2.name());
229 
230  if (!cAlphas_.found(key12))
231  {
233  << "Phase compression factor (cAlpha) not found for : "
234  << key12
235  << exit(FatalError);
236  }
237  scalar cAlpha = cAlphas_.find(key12)();
238 
239  phic = min(cAlpha*phic, max(phic));
240 
242 
243  word phirScheme
244  (
245  "div(phir," + alpha2.name() + ',' + alpha1.name() + ')'
246  );
247 
248  phiAlphaCorr += fvc::flux
249  (
250  -fvc::flux(-phir, alpha2, phirScheme),
251  alpha1,
252  phirScheme
253  );
254  }
255 
256  // Ensure that the flux at inflow BCs is preserved
257  forAll(phiAlphaCorr.boundaryField(), patchi)
258  {
259  fvsPatchScalarField& phiAlphaCorrp =
260  phiAlphaCorr.boundaryFieldRef()[patchi];
261 
262  if (!phiAlphaCorrp.coupled())
263  {
264  const scalarField& phi1p = phi.boundaryField()[patchi];
265  const scalarField& alpha1p =
266  alpha1.boundaryField()[patchi];
267 
268  forAll(phiAlphaCorrp, facei)
269  {
270  if (phi1p[facei] < 0)
271  {
272  phiAlphaCorrp[facei] = alpha1p[facei]*phi1p[facei];
273  }
274  }
275  }
276  }
277 
278  ++phasei;
279  }
280 
281  // Set Su and Sp to zero
282  for (const multiphaseInter::phaseModel& phase : phases_)
283  {
284  Su_[phase.name()] = dimensionedScalar("Su", dimless/dimTime, Zero);
285  Sp_[phase.name()] = dimensionedScalar("Sp", dimless/dimTime, Zero);
286 
287  // Add alpha*div(U)
288  //const volScalarField& alpha = phase;
289  //Su_[phase.name()] +=
290  // fvc::div(phi)*min(max(alpha, scalar(0)), scalar(1));
291  }
292 
293  // Fill Su and Sp
294  calculateSuSp();
295 
296  // Limit phiAlphaCorr on each phase
297  phasei = 0;
298  for (multiphaseInter::phaseModel& phase : phases_)
299  {
300  volScalarField& alpha1 = phase;
301 
302  surfaceScalarField& phiAlphaCorr = phiAlphaCorrs[phasei];
303 
304  volScalarField::Internal& Su = Su_[phase.name()];
305  volScalarField::Internal& Sp = Sp_[phase.name()];
306 
308  (
309  1.0/mesh_.time().deltaT().value(),
310  geometricOneField(),
311  alpha1,
312  phi,
313  phiAlphaCorr,
314  Sp,
315  Su,
316  oneField(),
317  zeroField(),
318  true
319  );
320  ++phasei;
321  }
322 
323  MULES::limitSum(phiAlphaCorrs);
324 
325  volScalarField sumAlpha
326  (
327  IOobject
328  (
329  "sumAlpha",
330  mesh_.time().timeName(),
331  mesh_
332  ),
333  mesh_,
335  );
336 
337  phasei = 0;
338  for (multiphaseInter::phaseModel& phase : phases_)
339  {
340  volScalarField& alpha1 = phase;
341 
342  const volScalarField::Internal& Su = Su_[phase.name()];
343 
344  const volScalarField::Internal& Sp = Sp_[phase.name()];
345 
346  surfaceScalarField& phiAlpha = phiAlphaCorrs[phasei];
347 
348  // Add a bounded upwind U-mean flux
349  //phiAlpha += upwind<scalar>(mesh_, phi).flux(alpha1);
350  fvScalarMatrix alpha1Eqn
351  (
352  fv::EulerDdtScheme<scalar>(mesh_).fvmDdt(alpha1)
353  + fv::gaussConvectionScheme<scalar>
354  (
355  mesh_,
356  phi,
357  upwind<scalar>(mesh_, phi)
358  ).fvmDiv(phi, alpha1)
359  ==
360  Su + fvm::Sp(Sp, alpha1)
361  );
362 
363  alpha1Eqn.boundaryManipulate(alpha1.boundaryFieldRef());
364 
365  alpha1Eqn.solve();
366 
367  phiAlpha += alpha1Eqn.flux();
368 
370  (
371  geometricOneField(),
372  alpha1,
373  phi,
374  phiAlpha,
375  Sp,
376  Su,
377  oneField(),
378  zeroField()
379  );
380 
381  phase.alphaPhi() = phiAlpha;
382 
383  ++phasei;
384  }
385 
386  if (acorr == nAlphaCorr - 1)
387  {
388  volScalarField sumAlpha
389  (
390  IOobject
391  (
392  "sumAlpha",
393  mesh_.time().timeName(),
394  mesh_
395  ),
396  mesh_,
398  );
399 
400  // Reset rhoPhi
401  rhoPhi_ = dimensionedScalar("rhoPhi", dimMass/dimTime, Zero);
402 
403  for (multiphaseInter::phaseModel& phase : phases_)
404  {
405  volScalarField& alpha1 = phase;
406  sumAlpha += alpha1;
407 
408  // Update rhoPhi
409  rhoPhi_ += fvc::interpolate(phase.rho()) * phase.alphaPhi();
410  }
411 
412  Info<< "Phase-sum volume fraction, min, max = "
413  << sumAlpha.weightedAverage(mesh_.V()).value()
414  << ' ' << min(sumAlpha).value()
415  << ' ' << max(sumAlpha).value()
416  << endl;
417 
418  volScalarField sumCorr(1.0 - sumAlpha);
419 
420  for (multiphaseInter::phaseModel& phase : phases_)
421  {
422  volScalarField& alpha = phase;
423  alpha += alpha*sumCorr;
424 
425  Info<< alpha.name() << " volume fraction = "
426  << alpha.weightedAverage(mesh_.V()).value()
427  << " Min(alpha) = " << min(alpha).value()
428  << " Max(alpha) = " << max(alpha).value()
429  << endl;
430  }
431  }
432  }
433 }
434 
435 
438 {
439  return phases_;
440 }
441 
442 
445 {
446  return phases_;
447 }
448 
449 
452 {
453  return phases_[i];
454 }
455 
456 
459 {
460  return phases_[i];
461 }
462 
463 
466 {
467  return ddtAlphaMax_;
468 }
469 
470 
471 Foam::scalar
473 {
474  auto iter = phaseModels_.cbegin();
475 
476  scalar maxVal = max(iter()->diffNo()).value();
477 
478  for (++iter; iter != phaseModels_.cend(); ++iter)
479  {
480  maxVal = max(maxVal, max(iter()->diffNo()).value());
481  }
482 
483  return maxVal * mesh_.time().deltaT().value();
484 }
485 
486 
489 {
490  return limitedPhiAlphas_;
491 }
492 
493 
496 {
497  return Su_;
498 }
499 
500 
503 {
504  return Sp_;
505 }
506 
507 
509 {
510  return true;
511 }
512 
513 
514 // ************************************************************************* //
const fvMesh & mesh() const
Return mesh.
fvsPatchField< scalar > fvsPatchScalarField
UPtrList< phaseModel > phases_
Unallocated phase list.
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:37
zeroField Su
Definition: alphaSuSp.H:1
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
label phasei
Definition: pEqn.H:27
virtual bool read()
Read thermophysical properties dictionary.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:578
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition: List.H:463
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
vtkSmartPointer< vtkFloatArray > zeroField(const word &name, const label size)
Create named field initialized to zero.
Calculate the snGrad of the given volField.
label size() const noexcept
The number of elements in table.
Definition: HashTableI.H:45
const word & name() const
The name of this phase.
Definition: phaseModel.H:122
const volScalarField & alpha2
const dimensionSet dimless
Dimensionless.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:361
const dictionary & solverDict(const word &name) const
Return the solver controls dictionary for the given field.
Definition: solution.C:363
phaseModel & phase2
DimensionedField< scalar, volMesh > Internal
The internal field type from which this GeometricField is derived.
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie...
Definition: UPtrListI.H:134
Area-weighted average a surfaceField creating a volField.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:84
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition: HashTableI.H:173
Calculate the first temporal derivative.
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect...
const fvMesh & mesh_
Reference to the mesh.
virtual void solve()
Solve for the phase fractions.
dynamicFvMesh & mesh
const dictionary & alphaControls
Definition: alphaControls.H:1
void limitSum(UPtrList< scalarField > &phiPsiCorrs)
Definition: MULES.C:27
multiphaseSystem(const fvMesh &)
Construct from fvMesh.
const phaseModel & phase(const label i) const
Constant access phase model i.
Calculate the gradient of the given field.
surfaceScalarField phir(fvc::flux(UdmModel.Udm()))
virtual bool coupled() const
True if this patch field is coupled.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Calculate the face-flux of the given field.
Calculate the matrix for the first temporal derivative.
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:328
dimensionedScalar ddtAlphaMax() const
Access to ddtAlphaMax.
Perform a subCycleTime on a field.
Definition: subCycle.H:133
const compressionFluxTable & limitedPhiAlphas() const
Access to compression fluxes for phaes.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
Calculate the divergence of the given field.
defineTypeNameAndDebug(interfaceCompositionModel, 0)
phaseModel & phase1
static word timeName(const scalar t, const int precision=precision_)
Return time name of given scalar time formatted with the given precision.
Definition: Time.C:760
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:193
scalarTable cAlphas_
Table for compression factors between phases.
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.
surfaceScalarField::Boundary & phicBf
Definition: alphaEqn.H:75
const word & name() const noexcept
Return const reference to name.
surfaceScalarField phic(mixture.cAlpha() *mag(alphaPhic/mesh.magSf()))
SuSpTable Su_
Su phase source terms.
phaseModelTable phaseModels_
Phase models.
label nAlphaSubCycles(alphaControls.get< label >("nAlphaSubCycles"))
SuSpTable Sp_
Sp phase source terms.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:79
IOobject(const IOobject &)=default
Copy construct.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
scalar maxDiffNo() const
Maximum diffusion number.
messageStream Info
Information stream (stdout output on master, null elsewhere)
defineRunTimeSelectionTable(interfaceCompositionModel, dictionary)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
void setSize(const label n)
Alias for resize()
Definition: UPtrList.H:337
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:166
const UPtrList< phaseModel > & phases() const
Return phases.
label nAlphaCorr(alphaControls.get< label >("nAlphaCorr"))
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
void calculateSuSp()
Calculate Sp and Su.
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
Calculate the finiteVolume matrix for implicit and explicit sources.
Namespace for OpenFOAM.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
zeroField Sp
Definition: alphaSuSp.H:2
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157
const volScalarField & alpha1