Loading...
Searching...
No Matches
alphaEqnSubCycle.H
Go to the documentation of this file.
1{
3 (
4 IOobject
5 (
6 "alphaPhi",
7 runTime.timeName(),
8 mesh
9 ),
10 mesh,
11 dimensionedScalar(phi.dimensions(), Zero)
12 );
13
14 surfaceScalarField phir(fvc::flux(UdmModel.Udm()));
15
17 {
18 dimensionedScalar totalDeltaT = runTime.deltaT();
19 surfaceScalarField alphaPhiSum
20 (
21 mesh.newIOobject("alphaPhiSum"),
22 mesh,
23 dimensionedScalar(phi.dimensions(), Zero)
24 );
25
26 for
27 (
28 subCycle<volScalarField> alphaSubCycle(alpha1, nAlphaSubCycles);
29 !(++alphaSubCycle).end();
30 )
31 {
32 #include "alphaEqn.H"
33 alphaPhiSum += (runTime.deltaT()/totalDeltaT)*alphaPhi;
34 }
35
36 alphaPhi = alphaPhiSum;
37 }
38 else
39 {
40 #include "alphaEqn.H"
41 }
42
43 // Apply the diffusion term separately to allow implicit solution
44 // and boundedness of the explicit advection
45 {
46 fvScalarMatrix alpha1Eqn
47 (
48 fvm::ddt(alpha1) - fvc::ddt(alpha1)
49 - fvm::laplacian(turbulence->nut(), alpha1)
50 );
51
52 alpha1Eqn.solve("alpha1Diffusion");
53
54 alphaPhi += alpha1Eqn.flux();
55 alpha2 = 1.0 - alpha1;
56
57 Info<< "Phase-1 volume fraction = "
58 << alpha1.weightedAverage(mesh.Vsc()).value()
59 << " Min(" << alpha1.name() << ") = " << min(alpha1).value()
60 << " Max(" << alpha1.name() << ") = " << max(alpha1).value()
61 << endl;
62 }
63
65 rho = mixture.rho();
66}
rhoPhi
Definition rhoEqn.H:10
const volScalarField & alpha1
volScalarField & rho2
const volScalarField & alpha2
volScalarField & rho1
dynamicFvMesh & mesh
engineTime & runTime
surfaceScalarField phir(fvc::flux(UdmModel.Udm()))
compressible::turbulenceModel & turbulence
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
fvMatrix< scalar > fvScalarMatrix
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\n"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture
label nAlphaSubCycles(alphaControls.get< label >("nAlphaSubCycles"))
surfaceScalarField alphaPhi(phi.name()+alpha1.name(), fvc::flux(phi, alpha1, alphaScheme))