alphaEqns.H
Go to the documentation of this file.
1 {
2  word alphaScheme("div(phi,alpha)");
3  word alpharScheme("div(phirb,alpha)");
4 
6 
7  for (int gCorr=0; gCorr<nAlphaCorr; gCorr++)
8  {
9  volScalarField::Internal Sp
10  (
11  IOobject
12  (
13  "Sp",
14  runTime.timeName(),
15  mesh
16  ),
17  mesh,
18  dimensionedScalar("Sp", dgdt.dimensions(), 0.0)
19  );
20 
21  volScalarField::Internal Su
22  (
23  IOobject
24  (
25  "Su",
26  runTime.timeName(),
27  mesh
28  ),
29  // Divergence term is handled explicitly to be
30  // consistent with the explicit transport solution
31  divU*min(alpha1, scalar(1))
32  );
33 
34  forAll(dgdt, celli)
35  {
36  if (dgdt[celli] > 0.0 && alpha1[celli] > 0.0)
37  {
38  Sp[celli] -= dgdt[celli]*alpha1[celli];
39  Su[celli] += dgdt[celli]*alpha1[celli];
40  }
41  else if (dgdt[celli] < 0.0 && alpha1[celli] < 1.0)
42  {
43  Sp[celli] += dgdt[celli]*(1.0 - alpha1[celli]);
44  }
45  }
46 
47 
49  (
50  fvc::flux
51  (
52  phi,
53  alpha1,
54  alphaScheme
55  )
56  + fvc::flux
57  (
59  alpha1,
61  )
62  );
63 
65  (
66  geometricOneField(),
67  alpha1,
68  phi,
69  alphaPhi1,
70  Sp,
71  Su,
72  1,
73  0
74  );
75 
78  rhoPhi = alphaPhi1*(rho1f - rho2f) + phi*rho2f;
79 
80  alpha2 = scalar(1) - alpha1;
81  }
82 
83  Info<< "Liquid phase volume fraction = "
84  << alpha1.weightedAverage(mesh.V()).value()
85  << " Min(" << alpha1.name() << ") = " << min(alpha1).value()
86  << " Min(" << alpha2.name() << ") = " << min(alpha2).value()
87  << endl;
88 }
surfaceScalarField & phi
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
interfaceProperties interface(alpha1, U, mixture())
word alpharScheme("div(phirb,alpha)")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
tmp< GeometricField< Type, fvPatchField, volMesh > > Sp(const volScalarField &sp, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:67
tmp< surfaceScalarField > interpolate(const RhoType &rho)
label nAlphaCorr(readLabel(alphaControls.lookup("nAlphaCorr")))
surfaceScalarField phir(phic *interface.nHatf())
rho1
Definition: pEqn.H:114
dynamicFvMesh & mesh
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
rho2
Definition: pEqn.H:115
alpha2
Definition: alphaEqn.H:112
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
volScalarField & alpha1
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
rhoPhi
Definition: rhoEqn.H:10
messageStream Info
volScalarField divU(fvc::div(fvc::absolute(phi, U)))
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
tmp< GeometricField< Type, fvPatchField, volMesh > > Su(const GeometricField< Type, fvPatchField, volMesh > &su, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:44
surfaceScalarField & alphaPhi1