alphaEqnSubCycle.H
Go to the documentation of this file.
1 #include "alphaControls.H"
2 
3 {
5  (
6  IOobject
7  (
8  "alphaPhi",
9  runTime.timeName(),
10  mesh
11  ),
12  mesh,
13  dimensionedScalar(phi.dimensions(), 0)
14  );
15 
16  surfaceScalarField phir(fvc::flux(UdmModel.Udm()));
17 
18  if (nAlphaSubCycles > 1)
19  {
20  dimensionedScalar totalDeltaT = runTime.deltaT();
21  surfaceScalarField alphaPhiSum
22  (
23  IOobject
24  (
25  "alphaPhiSum",
26  runTime.timeName(),
27  mesh
28  ),
29  mesh,
30  dimensionedScalar(phi.dimensions(), 0)
31  );
32 
33  for
34  (
35  subCycle<volScalarField> alphaSubCycle(alpha1, nAlphaSubCycles);
36  !(++alphaSubCycle).end();
37  )
38  {
39  #include "alphaEqn.H"
40  alphaPhiSum += (runTime.deltaT()/totalDeltaT)*alphaPhi;
41  }
42 
43  alphaPhi = alphaPhiSum;
44  }
45  else
46  {
47  #include "alphaEqn.H"
48  }
49 
50  // Apply the diffusion term separately to allow implicit solution
51  // and boundedness of the explicit advection
52  {
53  fvScalarMatrix alpha1Eqn
54  (
57  );
58 
59  alpha1Eqn.solve("alpha1Diffusion");
60 
61  alphaPhi += alpha1Eqn.flux();
62  alpha2 = 1.0 - alpha1;
63 
64  Info<< "Phase-1 volume fraction = "
65  << alpha1.weightedAverage(mesh.Vsc()).value()
66  << " Min(" << alpha1.name() << ") = " << min(alpha1).value()
67  << " Max(" << alpha1.name() << ") = " << max(alpha1).value()
68  << endl;
69  }
70 
71  rhoPhi = alphaPhi*(rho1 - rho2) + phi*rho2;
72  rho = mixture.rho();
73 }
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
surfaceScalarField phir(fvc::flux(UdmModel.Udm()))
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
alpha2
Definition: alphaEqn.H:115
phi
Definition: pEqn.H:104
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
dynamicFvMesh & mesh
volScalarField & alpha1(mixture.alpha1())
surfaceScalarField alphaPhi(phi.name()+alpha1.name(), fvc::flux(phi, alpha1, alphaScheme))
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar(dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar(dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::momentumTransportModel > turbulence(compressible::momentumTransportModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
const label nAlphaSubCycles(alphaControls.lookup< label >("nAlphaSubCycles"))
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
rhoPhi
Definition: rhoEqn.H:10
messageStream Info
const dimensionedScalar & rho2
Definition: createFields.H:44
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
phaseChangeTwoPhaseMixture & mixture
Definition: createFields.H:38
const dimensionedScalar & rho1
Definition: createFields.H:43