EEqns.H
Go to the documentation of this file.
1 {
2  volScalarField& he1 = thermo1.he();
3  volScalarField& he2 = thermo2.he();
4 
5  volScalarField Cpv1("Cpv1", thermo1.Cpv());
6  volScalarField Cpv2("Cpv2", thermo2.Cpv());
7 
8  volScalarField Kh(fluid.Kh());
9 
11  (
13  - fvm::Sp(contErr1, he1)
14 
15  + fvc::ddt(alpha1, rho1, K1) + fvc::div(alphaRhoPhi1, K1)
16  - contErr1*K1
17  + (
18  he1.name() == thermo1.phasePropertyName("e")
20  : -alpha1*dpdt
21  )
22 
23  - fvm::laplacian
24  (
25  fvc::interpolate(alpha1)
27  he1
28  )
29  );
30 
31  E1Eqn.relax();
32 
33  E1Eqn -=
34  (
35  Kh*(thermo2.T() - thermo1.T())
36  + Kh*he1/Cpv1
37  - fvm::Sp(Kh/Cpv1, he1)
38  + alpha1*rho1*(U1&g)
39  + fvOptions(alpha1, rho1, he1)
40  );
41 
43  (
45  - fvm::Sp(contErr2, he2)
46 
47  + fvc::ddt(alpha2, rho2, K2) + fvc::div(alphaRhoPhi2, K2)
48  - contErr2*K2
49  + (
50  he2.name() == thermo2.phasePropertyName("e")
52  : -alpha2*dpdt
53  )
54 
55  - fvm::laplacian
56  (
57  fvc::interpolate(alpha2)
59  he2
60  )
61  );
62 
63  E2Eqn.relax();
64 
65  E2Eqn -=
66  (
67  Kh*(thermo1.T() - thermo2.T())
68  + Kh*he2/Cpv2
69  - fvm::Sp(Kh/Cpv2, he2)
70  + alpha2*rho2*(U2&g)
71  + fvOptions(alpha2, rho2, he2)
72  );
73 
74  fvOptions.constrain(E1Eqn);
75  E1Eqn.solve();
76 
77  fvOptions.constrain(E2Eqn);
78  E2Eqn.solve();
79 
80  thermo1.correct();
81  Info<< "min " << thermo1.T().name()
82  << " " << min(thermo1.T()).value() << endl;
83 
84  thermo2.correct();
85  Info<< "min " << thermo2.T().name()
86  << " " << min(thermo2.T()).value() << endl;
87 }
rhoThermo & thermo2
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
contErr1
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
multiphaseSystem & fluid
Definition: createFields.H:10
#define K1
Definition: SHA1.C:167
E1Eqn
Definition: EEqns.H:33
contErr2
#define K2
Definition: SHA1.C:168
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
surfaceScalarField & alphaPhi2
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
E2Eqn
Definition: EEqns.H:65
tmp< surfaceScalarField > interpolate(const RhoType &rho)
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
surfaceScalarField & alphaRhoPhi1
fv::options & fvOptions
rho1
Definition: pEqn.H:114
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:52
rhoThermo & thermo1
volScalarField Kh(fluid.Kh())
surfaceScalarField & alphaRhoPhi2
rho2
Definition: pEqn.H:115
volScalarField & dpdt
volVectorField & U1
alpha2
Definition: alphaEqn.H:112
phaseModel & phase1
const dimensionedVector & g
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
volScalarField & alpha1
volScalarField Cpv2("Cpv2", thermo2.Cpv())
volVectorField & U2
messageStream Info
volScalarField alphaEff("alphaEff", turbulence->nu()/Pr+alphat)
phaseModel & phase2
volScalarField & p
volScalarField Cpv1("Cpv1", thermo1.Cpv())
surfaceScalarField & alphaPhi1