UEqns.H
Go to the documentation of this file.
1 Info<< "Constructing face momentum equations" << endl;
2 
3 fvVectorMatrix U1Eqn(U1, rho1.dimensions()*U1.dimensions()*dimVol/dimTime);
4 fvVectorMatrix U2Eqn(U2, rho2.dimensions()*U2.dimensions()*dimVol/dimTime);
5 
6 {
7  volScalarField Vm(fluid.Vm());
8 
9  fvVectorMatrix UgradU1
10  (
12  + MRF.DDt(U1)
13  );
14 
15  fvVectorMatrix UgradU2
16  (
18  + MRF.DDt(U2)
19  );
20 
21  const volScalarField dmdt21(posPart(fluid.dmdt()));
22  const volScalarField dmdt12(negPart(fluid.dmdt()));
23 
24  {
25  U1Eqn =
26  (
28  + fvm::Sp(dmdt21, U1) - dmdt21*U2
29  + MRF.DDt(alpha1*rho1, U1)
30  + phase1.turbulence().divDevRhoReff(U1)
31  + Vm*(UgradU1 - (UgradU2 & U2))
32  - fvOptions(alpha1, rho1, U1)
33  );
34  U1Eqn.relax();
35  fvOptions.constrain(U1Eqn);
36  U1.correctBoundaryConditions();
37  fvOptions.correct(U1);
38  }
39 
40  {
41  U2Eqn =
42  (
44  - fvm::Sp(dmdt12, U2) + dmdt12*U1
45  + MRF.DDt(alpha2*rho2, U2)
46  + phase2.turbulence().divDevRhoReff(U2)
47  + Vm*(UgradU2 - (UgradU1 & U1))
48  - fvOptions(alpha2, rho2, U2)
49  );
50  U2Eqn.relax();
51  fvOptions.constrain(U2Eqn);
52  U2.correctBoundaryConditions();
53  fvOptions.correct(U2);
54  }
55 }
surfaceScalarField & phi2
multiphaseSystem & fluid
Definition: createFields.H:10
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:253
const dimensionSet dimVol(dimVolume)
Definition: dimensionSets.H:59
tmp< GeometricField< Type, fvPatchField, volMesh > > Sp(const volScalarField &sp, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:67
surfaceScalarField & phi1
surfaceScalarField & alphaRhoPhi1
dimensionedScalar posPart(const dimensionedScalar &ds)
fv::options & fvOptions
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
IOMRFZoneList & MRF
surfaceScalarField & alphaRhoPhi2
volVectorField & U1
alpha2
Definition: alphaEqn.H:112
phaseModel & phase1
volScalarField & alpha1
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:45
volVectorField & U2
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
phaseModel & phase2
dimensionedScalar negPart(const dimensionedScalar &ds)