compressibleMultiphaseMixture.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-2022 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 
27 #include "alphaContactAngleFvPatchScalarField.H"
28 #include "correctContactAngle.H"
29 #include "Time.H"
30 #include "subCycle.H"
31 #include "MULES.H"
32 #include "fvcDiv.H"
33 #include "fvcGrad.H"
34 #include "fvcSnGrad.H"
35 #include "fvcFlux.H"
36 #include "fvcMeshPhi.H"
37 #include "surfaceInterpolate.H"
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43  defineTypeNameAndDebug(compressibleMultiphaseMixture, 0);
44 }
45 
46 
47 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48 
49 void Foam::compressibleMultiphaseMixture::calcAlphas()
50 {
51  scalar level = 0.0;
52  alphas_ == 0.0;
53 
54  forAllIter(PtrDictionary<phaseModel>, phases_, phase)
55  {
56  alphas_ += level*phase();
57  level += 1.0;
58  }
59 }
60 
61 
62 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
63 
65 (
66  const volVectorField& U,
67  const surfaceScalarField& phi
68 )
69 :
70  IOdictionary
71  (
72  IOobject
73  (
74  "phaseProperties",
75  U.time().constant(),
76  U.db(),
77  IOobject::MUST_READ_IF_MODIFIED,
78  IOobject::NO_WRITE
79  )
80  ),
81 
82  mesh_(U.mesh()),
83 
84  p_
85  (
86  IOobject
87  (
88  "p",
89  U.mesh().time().timeName(),
90  U.mesh(),
91  IOobject::MUST_READ,
92  IOobject::AUTO_WRITE
93  ),
94  U.mesh()
95  ),
96 
97  T_
98  (
99  IOobject
100  (
101  "T",
102  U.mesh().time().timeName(),
103  U.mesh(),
104  IOobject::MUST_READ,
105  IOobject::AUTO_WRITE
106  ),
107  U.mesh()
108  ),
109 
110  phases_(lookup("phases"), phaseModel::iNew(p_, T_)),
111 
112  rho_
113  (
114  IOobject
115  (
116  "thermo:rho",
117  U.mesh().time().timeName(),
118  U.mesh(),
119  IOobject::NO_READ,
120  IOobject::AUTO_WRITE
121  ),
122  U.mesh(),
123  dimensionedScalar("rho", dimDensity, 0)
124  ),
125 
126  U_(U),
127 
128  phi_(phi),
129 
130  rhoPhi_
131  (
132  IOobject
133  (
134  "rhoPhi",
135  mesh_.time().timeName(),
136  mesh_,
137  IOobject::NO_READ,
138  IOobject::NO_WRITE
139  ),
140  mesh_,
142  ),
143 
144  alphas_
145  (
146  IOobject
147  (
148  "alphas",
149  mesh_.time().timeName(),
150  mesh_,
151  IOobject::NO_READ,
152  IOobject::AUTO_WRITE
153  ),
154  mesh_,
156  ),
157 
158  sigmas_(lookup("sigmas")),
159  dimSigma_(1, 0, -2, 0, 0),
160  deltaN_
161  (
162  "deltaN",
163  1e-8/pow(average(mesh_.V()), 1.0/3.0)
164  )
165 {
166  calcAlphas();
167  alphas_.write();
168  correct();
169 }
170 
171 
172 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
173 
175 {
176  forAllIter(PtrDictionary<phaseModel>, phases_, phasei)
177  {
178  phasei().correct();
179  }
180 }
181 
182 
184 {
185  PtrDictionary<phaseModel>::iterator phasei = phases_.begin();
186 
187  rho_ = phasei()*phasei().thermo().rho();
188 
189  for (++phasei; phasei != phases_.end(); ++phasei)
190  {
191  rho_ += phasei()*phasei().thermo().rho();
192  }
193 
194  forAllIter(PtrDictionary<phaseModel>, phases_, phasei)
195  {
196  phasei().Alpha() = phasei()*phasei().thermo().rho()/rho_;
197  }
198 }
199 
200 
202 {
203  forAllIter(PtrDictionary<phaseModel>, phases_, phasei)
204  {
205  phasei().thermo().rho() += phasei().thermo().psi()*dp;
206  }
207 }
208 
209 
211 {
212  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
213 
214  volScalarField mu(phasei().Alpha()*phasei().thermo().mu());
215 
216  for (++phasei; phasei != phases_.end(); ++phasei)
217  {
218  mu += phasei().Alpha()*phasei().thermo().mu();
219  }
220 
221  return mu/rho_;
222 }
223 
224 
226 (
227  const label patchi
228 ) const
229 {
230  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
231 
233  (
234  phasei().Alpha().boundaryField()[patchi]*phasei().thermo().mu(patchi)
235  );
236 
237  for (++phasei; phasei != phases_.end(); ++phasei)
238  {
239  mu +=
240  phasei().Alpha().boundaryField()[patchi]
241  *phasei().thermo().mu(patchi);
242  }
243 
244  return mu/rho_.boundaryField()[patchi];
245 }
246 
247 
249 (
250  const volScalarField& alphat
251 ) const
252 {
253  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
254 
255  tmp<volScalarField> talphaEff(phasei()*phasei().thermo().alphaEff(alphat));
256 
257  for (++phasei; phasei != phases_.end(); ++phasei)
258  {
259  talphaEff.ref() += phasei()*phasei().thermo().alphaEff(alphat);
260  }
261 
262  return talphaEff;
263 }
264 
265 
267 {
268  PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
269 
270  tmp<volScalarField> trCv(phasei()/phasei().thermo().Cv());
271 
272  for (++phasei; phasei != phases_.end(); ++phasei)
273  {
274  trCv.ref() += phasei()/phasei().thermo().Cv();
275  }
276 
277  return trCv;
278 }
279 
280 
283 {
284  tmp<surfaceScalarField> tstf
285  (
287  (
288  "surfaceTensionForce",
289  mesh_,
290  dimensionedScalar(dimensionSet(1, -2, -2, 0, 0), 0)
291  )
292  );
293 
294  surfaceScalarField& stf = tstf.ref();
295 
296  forAllConstIter(PtrDictionary<phaseModel>, phases_, phase1)
297  {
298  const phaseModel& alpha1 = phase1();
299 
301  ++phase2;
302 
303  for (; phase2 != phases_.end(); ++phase2)
304  {
305  const phaseModel& alpha2 = phase2();
306 
307  sigmaTable::const_iterator sigma =
308  sigmas_.find(interfacePair(alpha1, alpha2));
309 
310  if (sigma == sigmas_.end())
311  {
313  << "Cannot find interface " << interfacePair(alpha1, alpha2)
314  << " in list of sigma values"
315  << exit(FatalError);
316  }
317 
318  stf += dimensionedScalar(dimSigma_, sigma())
319  *fvc::interpolate(K(alpha1, alpha2))*
320  (
321  fvc::interpolate(alpha2)*fvc::snGrad(alpha1)
322  - fvc::interpolate(alpha1)*fvc::snGrad(alpha2)
323  );
324  }
325  }
326 
327  return tstf;
328 }
329 
330 
332 {
333  const Time& runTime = mesh_.time();
334 
335  const dictionary& alphaControls = mesh_.solution().solverDict("alpha");
336  label nAlphaSubCycles(alphaControls.lookup<label>("nAlphaSubCycles"));
337  scalar cAlpha(alphaControls.lookup<scalar>("cAlpha"));
338 
339  volScalarField& alpha = phases_.first();
340 
341  if (nAlphaSubCycles > 1)
342  {
343  surfaceScalarField rhoPhiSum(0.0*rhoPhi_);
344  dimensionedScalar totalDeltaT = runTime.deltaT();
345 
346  for
347  (
348  subCycle<volScalarField> alphaSubCycle(alpha, nAlphaSubCycles);
349  !(++alphaSubCycle).end();
350  )
351  {
352  solveAlphas(cAlpha);
353  rhoPhiSum += (runTime.deltaT()/totalDeltaT)*rhoPhi_;
354  }
355 
356  rhoPhi_ = rhoPhiSum;
357  }
358  else
359  {
360  solveAlphas(cAlpha);
361  }
362 
363  correct();
364 }
365 
366 
367 Foam::tmp<Foam::surfaceVectorField> Foam::compressibleMultiphaseMixture::nHatfv
368 (
369  const volScalarField& alpha1,
370  const volScalarField& alpha2
371 ) const
372 {
373  /*
374  // Cell gradient of alpha
375  volVectorField gradAlpha =
376  alpha2*fvc::grad(alpha1) - alpha1*fvc::grad(alpha2);
377 
378  // Interpolated face-gradient of alpha
379  surfaceVectorField gradAlphaf = fvc::interpolate(gradAlpha);
380  */
381 
382  surfaceVectorField gradAlphaf
383  (
385  - fvc::interpolate(alpha1)*fvc::interpolate(fvc::grad(alpha2))
386  );
387 
388  // Face unit interface normal
389  return gradAlphaf/(mag(gradAlphaf) + deltaN_);
390 }
391 
392 
393 Foam::tmp<Foam::surfaceScalarField> Foam::compressibleMultiphaseMixture::nHatf
394 (
395  const volScalarField& alpha1,
396  const volScalarField& alpha2
397 ) const
398 {
399  // Face unit interface normal flux
400  return nHatfv(alpha1, alpha2) & mesh_.Sf();
401 }
402 
403 
404 Foam::tmp<Foam::volScalarField> Foam::compressibleMultiphaseMixture::K
405 (
406  const phaseModel& alpha1,
407  const phaseModel& alpha2
408 ) const
409 {
410  tmp<surfaceVectorField> tnHatfv = nHatfv(alpha1, alpha2);
411 
413  (
414  alpha1,
415  alpha2,
416  U_.boundaryField(),
417  deltaN_,
418  tnHatfv.ref().boundaryFieldRef()
419  );
420 
421  // Simple expression for curvature
422  return -fvc::div(tnHatfv & mesh_.Sf());
423 }
424 
425 
428 {
429  tmp<volScalarField> tnearInt
430  (
432  (
433  "nearInterface",
434  mesh_,
436  )
437  );
438 
439  forAllConstIter(PtrDictionary<phaseModel>, phases_, phase)
440  {
441  tnearInt.ref() =
442  max(tnearInt(), pos0(phase() - 0.01)*pos0(0.99 - phase()));
443  }
444 
445  return tnearInt;
446 }
447 
448 
449 void Foam::compressibleMultiphaseMixture::solveAlphas
450 (
451  const scalar cAlpha
452 )
453 {
454  static label nSolves=-1;
455  nSolves++;
456 
457  word alphaScheme("div(phi,alpha)");
458  word alpharScheme("div(phirb,alpha)");
459 
460  surfaceScalarField phic(mag(phi_/mesh_.magSf()));
461  phic = min(cAlpha*phic, max(phic));
462 
463  UPtrList<const volScalarField> alphas(phases_.size());
464  PtrList<surfaceScalarField> alphaPhis(phases_.size());
465  int phasei = 0;
466 
467  forAllIter(PtrDictionary<phaseModel>, phases_, phase)
468  {
469  const phaseModel& alpha = phase();
470 
471  alphas.set(phasei, &alpha);
472 
473  alphaPhis.set
474  (
475  phasei,
477  (
478  phi_.name() + alpha.name(),
479  fvc::flux
480  (
481  phi_,
482  alpha,
484  )
485  )
486  );
487 
488  surfaceScalarField& alphaPhi = alphaPhis[phasei];
489 
490  forAllIter(PtrDictionary<phaseModel>, phases_, phase2)
491  {
492  phaseModel& alpha2 = phase2();
493 
494  if (&alpha2 == &alpha) continue;
495 
496  surfaceScalarField phir(phic*nHatf(alpha, alpha2));
497 
498  alphaPhi += fvc::flux
499  (
500  -fvc::flux(-phir, alpha2, alpharScheme),
501  alpha,
503  );
504  }
505 
506  // Limit alphaPhi for each phase
508  (
509  1.0/mesh_.time().deltaT().value(),
510  geometricOneField(),
511  alpha,
512  phi_,
513  alphaPhi,
514  zeroField(),
515  zeroField(),
516  oneField(),
517  zeroField(),
518  false
519  );
520 
521  phasei++;
522  }
523 
524  MULES::limitSum(alphas, alphaPhis, phi_);
525 
526  rhoPhi_ = dimensionedScalar(dimensionSet(1, 0, -1, 0, 0), 0);
527 
528  volScalarField sumAlpha
529  (
530  IOobject
531  (
532  "sumAlpha",
533  mesh_.time().timeName(),
534  mesh_
535  ),
536  mesh_,
538  );
539 
541 
542  phasei = 0;
543 
544  forAllIter(PtrDictionary<phaseModel>, phases_, phase)
545  {
546  phaseModel& alpha = phase();
547 
548  surfaceScalarField& alphaPhi = alphaPhis[phasei];
549 
551  (
552  IOobject
553  (
554  "Sp",
555  mesh_.time().timeName(),
556  mesh_
557  ),
558  mesh_,
559  dimensionedScalar(alpha.dgdt().dimensions(), 0)
560  );
561 
563  (
564  IOobject
565  (
566  "Su",
567  mesh_.time().timeName(),
568  mesh_
569  ),
570  // Divergence term is handled explicitly to be
571  // consistent with the explicit transport solution
572  divU.v()*min(alpha.v(), scalar(1))
573  );
574 
575  {
576  const scalarField& dgdt = alpha.dgdt();
577 
578  forAll(dgdt, celli)
579  {
580  if (dgdt[celli] < 0.0 && alpha[celli] > 0.0)
581  {
582  Sp[celli] += dgdt[celli]*alpha[celli];
583  Su[celli] -= dgdt[celli]*alpha[celli];
584  }
585  else if (dgdt[celli] > 0.0 && alpha[celli] < 1.0)
586  {
587  Sp[celli] -= dgdt[celli]*(1.0 - alpha[celli]);
588  }
589  }
590  }
591 
592  forAllConstIter(PtrDictionary<phaseModel>, phases_, phase2)
593  {
594  const phaseModel& alpha2 = phase2();
595 
596  if (&alpha2 == &alpha) continue;
597 
598  const scalarField& dgdt2 = alpha2.dgdt();
599 
600  forAll(dgdt2, celli)
601  {
602  if (dgdt2[celli] > 0.0 && alpha2[celli] < 1.0)
603  {
604  Sp[celli] -= dgdt2[celli]*(1.0 - alpha2[celli]);
605  Su[celli] += dgdt2[celli]*alpha[celli];
606  }
607  else if (dgdt2[celli] < 0.0 && alpha2[celli] > 0.0)
608  {
609  Sp[celli] += dgdt2[celli]*alpha2[celli];
610  }
611  }
612  }
613 
615  (
616  geometricOneField(),
617  alpha,
618  alphaPhi,
619  Sp,
620  Su
621  );
622 
623  rhoPhi_ += fvc::interpolate(alpha.thermo().rho())*alphaPhi;
624 
625  Info<< alpha.name() << " volume fraction, min, max = "
626  << alpha.weightedAverage(mesh_.V()).value()
627  << ' ' << min(alpha).value()
628  << ' ' << max(alpha).value()
629  << endl;
630 
631  sumAlpha += alpha;
632 
633  phasei++;
634  }
635 
636  Info<< "Phase-sum volume fraction, min, max = "
637  << sumAlpha.weightedAverage(mesh_.V()).value()
638  << ' ' << min(sumAlpha).value()
639  << ' ' << max(sumAlpha).value()
640  << endl;
641 
642  calcAlphas();
643 }
644 
645 
646 // ************************************************************************* //
virtual tmp< volScalarField > nu() const
Kinematic viscosity of mixture [m^2/s].
scalar Cv(const scalar p, const scalar T) const
Definition: HtoEthermo.H:2
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)
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
fluidReactionThermo & thermo
Definition: createFields.H:28
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
label phasei
Definition: pEqn.H:27
compressibleMultiphaseMixture(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
void correctContactAngle(const volScalarField &alpha1, const volScalarField &alpha2, const volVectorField::Boundary &Ubf, const dimensionedScalar &deltaN, surfaceVectorField::Boundary &nHatbf)
Correct the contact angle for the two volume fraction fields.
const word alphaScheme(mesh.schemes().div(divAlphaName)[1].wordToken())
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
static tmp< GeometricField< scalar, fvsPatchField, surfaceMesh > > New(const word &name, const Internal &, const PtrList< fvsPatchField< scalar >> &)
Return a temporary field constructed from name,.
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
surfaceScalarField phir(fvc::flux(mixture.Udm()))
Calculate the snGrad of the given volField.
const dimensionSet dimless
friend class iterator
Definition: LPtrList.H:84
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:59
fvMesh & mesh
T & first()
Return the first element of the list.
Definition: UListI.H:114
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
CGAL::Exact_predicates_exact_constructions_kernel K
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
tmp< volScalarField > nearInterface() const
Indicator of the proximity of the interface.
const dimensionSet dimTime
stressControl lookup("compactNormalStress") >> compactNormalStress
void limitSum(UPtrList< scalarField > &phiPsiCorrs)
Definition: MULES.C:30
void solve()
Solve for the mixture phase-fractions.
tmp< volScalarField > rCv() const
Return the phase-averaged reciprocal Cv.
Calculate the gradient of the given field.
phic
Definition: correctPhic.H:2
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Calculate the face-flux of the given field.
const dimensionSet dimDensity
void correctRho(const volScalarField &dp)
Update densities for given pressure change.
virtual void correctThermo()
Correct the thermodynamics of each phase.
surfaceScalarField alphaPhi(phi.name()+alpha1.name(), fvc::flux(phi, alpha1, alphaScheme))
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
word timeName
Definition: getTimeIndex.H:3
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
Calculate the divergence of the given field.
dimensionedScalar pos0(const dimensionedScalar &ds)
const dimensionedScalar mu
Atomic mass unit.
const label nAlphaSubCycles(alphaControls.lookup< label >("nAlphaSubCycles"))
defineTypeNameAndDebug(combustionModel, 0)
tmp< volScalarField > divU
Definition: alphaSuSp.H:9
const dimensionSet dimMass
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.
const tmp< volScalarField::Internal > & Sp
Definition: alphaSuSp.H:7
Calculate the mesh motion flux and convert fluxes from absolute to relative and back.
Internal & ref()
Return a reference to the dimensioned internal field.
virtual tmp< volScalarField > alphaEff(const volScalarField &alphat) const
Effective thermal diffusivity of mixture [W/m/K].
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
volScalarField::Internal dgdt(IOobject("dgdt", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), alpha1 *fvc::div(phi))
label patchi
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:202
virtual void correct()
Update properties.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual void correct()
Correct the fvModel.
Definition: fvModel.C:176
word alpharScheme("div(phirb,alpha)")
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
MULES: Multidimensional universal limiter for explicit solution.
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
const tmp< volScalarField::Internal > & Su
Definition: alphaSuSp.H:6
A class for managing temporary objects.
Definition: PtrList.H:53
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)
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
friend class const_iterator
Definition: LPtrList.H:87
tmp< surfaceScalarField > surfaceTensionForce() const
Namespace for OpenFOAM.