alphaEqn.H
Go to the documentation of this file.
1 {
2  word alphaScheme("div(phi,alpha)");
3  word alpharScheme("div(phirb,alpha)");
4 
6  (
7  IOobject
8  (
9  "phir",
10  runTime.timeName(),
11  mesh,
12  IOobject::NO_READ,
13  IOobject::NO_WRITE
14  ),
15  mixture.cAlpha()*mag(phi/mesh.magSf())*mixture.nHatf()
16  );
17 
18  for (int gCorr=0; gCorr<nAlphaCorr; gCorr++)
19  {
20  // Split the limiter into a surfaceScalarField
22  (
23  IOobject
24  (
25  "lambda",
26  mesh.time().timeName(),
27  mesh,
28  IOobject::NO_READ,
29  IOobject::NO_WRITE,
30  false
31  ),
32  mesh,
34  );
35 
36  // Create the complete convection flux for alpha1
38  (
39  fvc::flux
40  (
41  phi,
42  alpha1,
44  )
45  + fvc::flux
46  (
48  alpha1,
50  )
51  + fvc::flux
52  (
53  -fvc::flux(-phir, alpha3, alpharScheme),
54  alpha1,
56  )
57  );
58 
59  // Create the bounded (upwind) flux for alpha1
60  surfaceScalarField alphaPhi1BD
61  (
62  upwind<scalar>(mesh, phi).flux(alpha1)
63  );
64 
65  // Calculate the flux correction for alpha1
66  alphaPhi1 -= alphaPhi1BD;
67 
68  // Calculate the limiter for alpha1
69  if (LTS)
70  {
71  const volScalarField& rDeltaT =
72  fv::localEulerDdt::localRDeltaT(mesh);
73 
75  (
76  lambda,
77  rDeltaT,
78  geometricOneField(),
79  alpha1,
80  alphaPhi1BD,
81  alphaPhi1,
82  zeroField(),
83  zeroField(),
84  oneField(),
85  zeroField()
86  );
87  }
88  else
89  {
91  (
92  lambda,
93  1.0/runTime.deltaT().value(),
94  geometricOneField(),
95  alpha1,
96  alphaPhi1BD,
97  alphaPhi1,
98  zeroField(),
99  zeroField(),
100  oneField(),
101  zeroField()
102  );
103  }
104 
105  // Create the complete flux for alpha2
107  (
108  fvc::flux
109  (
110  phi,
111  alpha2,
113  )
114  + fvc::flux
115  (
117  alpha2,
119  )
120  );
121 
122  // Create the bounded (upwind) flux for alpha2
123  surfaceScalarField alphaPhi2BD
124  (
125  upwind<scalar>(mesh, phi).flux(alpha2)
126  );
127 
128  // Calculate the flux correction for alpha2
129  alphaPhi2 -= alphaPhi2BD;
130 
131  // Further limit the limiter for alpha2
132  if (LTS)
133  {
134  const volScalarField& rDeltaT =
135  fv::localEulerDdt::localRDeltaT(mesh);
136 
138  (
139  lambda,
140  rDeltaT,
141  geometricOneField(),
142  alpha2,
143  alphaPhi2BD,
144  alphaPhi2,
145  zeroField(),
146  zeroField(),
147  oneField(),
148  zeroField()
149  );
150  }
151  else
152  {
154  (
155  lambda,
156  1.0/runTime.deltaT().value(),
157  geometricOneField(),
158  alpha2,
159  alphaPhi2BD,
160  alphaPhi2,
161  zeroField(),
162  zeroField(),
163  oneField(),
164  zeroField()
165  );
166  }
167 
168  // Construct the limited fluxes
169  alphaPhi1 = alphaPhi1BD + lambda*alphaPhi1;
170  alphaPhi2 = alphaPhi2BD + lambda*alphaPhi2;
171 
172  // Solve for alpha1
173  solve(fvm::ddt(alpha1) + fvc::div(alphaPhi1));
174 
175  // Create the diffusion coefficients for alpha2<->alpha3
176  volScalarField Dc23(D23*max(alpha3, scalar(0))*pos0(alpha2));
177  volScalarField Dc32(D23*max(alpha2, scalar(0))*pos0(alpha3));
178 
179  // Add the diffusive flux for alpha3->alpha2
180  alphaPhi2 -= fvc::interpolate(Dc32)*mesh.magSf()*fvc::snGrad(alpha1);
181 
182  // Solve for alpha2
183  fvScalarMatrix alpha2Eqn
184  (
186  + fvc::div(alphaPhi2)
187  - fvm::laplacian(Dc23 + Dc32, alpha2)
188  );
189  alpha2Eqn.solve();
190 
191  // Construct the complete mass flux
192  rhoPhi =
193  alphaPhi1*(rho1 - rho3)
194  + (alphaPhi2 + alpha2Eqn.flux())*(rho2 - rho3)
195  + phi*rho3;
196 
197  alpha3 = 1.0 - alpha1 - alpha2;
198  }
199 
200  Info<< "Air phase volume fraction = "
201  << alpha1.weightedAverage(mesh.V()).value()
202  << " Min(" << alpha1.name() << ") = " << min(alpha1).value()
203  << " Max(" << alpha1.name() << ") = " << max(alpha1).value()
204  << endl;
205 
206  Info<< "Liquid phase volume fraction = "
207  << alpha2.weightedAverage(mesh.V()).value()
208  << " Min(" << alpha2.name() << ") = " << min(alpha2).value()
209  << " Max(" << alpha2.name() << ") = " << max(alpha2).value()
210  << endl;
211 }
surfaceScalarField alphaPhi1(alphaPhi1Header, phi *fvc::interpolate(alpha1))
const volScalarField & rho1
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
dimensionedScalar lambda(viscosity->lookup("lambda"))
volScalarField & alpha1(mixture.alpha1())
const dimensionSet dimless
fvMesh & mesh
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
rhoPhi
Definition: alphaEqn.H:25
alpha2
Definition: alphaEqn.H:115
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
phi
Definition: correctPhi.H:3
void limiter(surfaceScalarField &lambda, const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phiBD, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su, const PsiMaxType &psiMax, const PsiMinType &psiMin)
dimensionedScalar pos0(const dimensionedScalar &ds)
bool LTS
Definition: createRDeltaT.H:1
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.
rhoEqn solve()
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const volScalarField & rho2
surfaceScalarField alphaPhi2("alphaPhi2", phi - alphaPhi1)
word alpharScheme("div(phirb,alpha)")
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
const label nAlphaCorr(alphaControls.lookup< label >("nAlphaCorr"))
surfaceScalarField phir(IOobject("phir", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), mixture.cAlpha() *mag(phi/mesh.magSf()) *mixture.nHatf())
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, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45