pEqn.H
Go to the documentation of this file.
1 {
2  volScalarField rAU("rAU", 1.0/UEqn.A());
6  (
7  "phiHbyA",
10  );
11  MRF.makeRelative(phiHbyA);
12 
14  (
15  (
16  mixture.surfaceTensionForce()
18  )*rAUf*mesh.magSf()
19  );
20 
22 
23  // Update the pressure BCs to ensure flux consistency
25 
26  tmp<fvScalarMatrix> p_rghEqnComp1;
27  tmp<fvScalarMatrix> p_rghEqnComp2;
28 
29  if (pimple.transonic())
30  {
32  surfaceScalarField phid2("phid2", fvc::interpolate(psi2)*phi);
33 
34  p_rghEqnComp1 =
35  fvc::ddt(rho1) + fvc::div(phi, rho1) - fvc::Sp(fvc::div(phi), rho1)
36  + correction
37  (
39  + fvm::div(phid1, p_rgh) - fvm::Sp(fvc::div(phid1), p_rgh)
40  );
41  deleteDemandDrivenData(p_rghEqnComp1.ref().faceFluxCorrectionPtr());
42  p_rghEqnComp1.ref().relax();
43 
44  p_rghEqnComp2 =
45  fvc::ddt(rho2) + fvc::div(phi, rho2) - fvc::Sp(fvc::div(phi), rho2)
46  + correction
47  (
49  + fvm::div(phid2, p_rgh) - fvm::Sp(fvc::div(phid2), p_rgh)
50  );
51  deleteDemandDrivenData(p_rghEqnComp2.ref().faceFluxCorrectionPtr());
52  p_rghEqnComp2.ref().relax();
53  }
54  else
55  {
56  #include "rhofs.H"
57 
58  p_rghEqnComp1 =
59  pos(alpha1)
60  *(
61  (
63  - (fvOptions(alpha1, mixture.thermo1().rho())&rho1)
64  )/rho1
67  );
68 
69  p_rghEqnComp2 =
70  pos(alpha2)
71  *(
72  (
74  - (fvOptions(alpha2, mixture.thermo2().rho())&rho2)
75  )/rho2
78  );
79  }
80 
81  // Cache p_rgh prior to solve for density update
83 
84  while (pimple.correctNonOrthogonal())
85  {
86  fvScalarMatrix p_rghEqnIncomp
87  (
90  );
91 
92  solve
93  (
94  p_rghEqnComp1() + p_rghEqnComp2() + p_rghEqnIncomp,
95  mesh.solver(p_rgh.select(pimple.finalInnerIter()))
96  );
97 
98  if (pimple.finalNonOrthogonalIter())
99  {
100  p = max(p_rgh + (alpha1*rho1 + alpha2*rho2)*gh, pMin);
101  p_rgh = p - (alpha1*rho1 + alpha2*rho2)*gh;
102 
103  dgdt =
104  (
105  alpha1*(p_rghEqnComp2 & p_rgh)
106  - alpha2*(p_rghEqnComp1 & p_rgh)
107  );
108 
109  phi = phiHbyA + p_rghEqnIncomp.flux();
110 
111  U = HbyA
112  + rAU*fvc::reconstruct((phig + p_rghEqnIncomp.flux())/rAUf);
113  U.correctBoundaryConditions();
114  fvOptions.correct(U);
115  }
116  }
117 
118  // Update densities from change in p_rgh
119  mixture.thermo1().correctRho(psi1*(p_rgh - p_rgh_0));
120  mixture.thermo2().correctRho(psi2*(p_rgh - p_rgh_0));
121 
123 
124  // Correct p_rgh for consistency with p and the updated densities
125  p_rgh = p - rho*gh;
126  p_rgh.correctBoundaryConditions();
127 
128  K = 0.5*magSqr(U);
129 }
tmp< fvScalarMatrix > p_rghEqnComp2
Definition: pEqn.H:30
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
psi2
Definition: TEqns.H:35
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
U
Definition: pEqn.H:83
p
Definition: pEqn.H:50
surfaceScalarField rho1f(fvc::interpolate(rho1))
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
phiHbyA
Definition: pEqn.H:20
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
surfaceScalarField & alphaPhi2
surfaceScalarField phig("phig", -rhorAUf *ghf *fvc::snGrad(rho) *mesh.magSf())
tmp< GeometricField< typename flux< Type >::type, fvsPatchField, surfaceMesh > > ddtCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
Definition: fvcDdt.C:170
const dictionary & pimple
const surfaceScalarField & alphaPhi1
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
CGAL::Exact_predicates_exact_constructions_kernel K
tmp< fvScalarMatrix > p_rghEqnComp1
Definition: pEqn.H:29
fv::options & fvOptions
const surfaceScalarField & ghf
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
dimensionedScalar pos(const dimensionedScalar &ds)
surfaceScalarField rho2f(fvc::interpolate(rho2))
dynamicFvMesh & mesh
IOMRFZoneList & MRF
rhoEqn solve()
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:33
alpha2
Definition: alphaEqn.H:112
p_rgh
Definition: pEqn.H:134
dimensioned< scalar > magSqr(const dimensioned< Type > &)
constrainPressure(p, rho, U, phiHbyA, rhorAUf, MRF)
Info<< "Reading field p_rgh\"<< endl;volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Creating phaseChangeTwoPhaseMixture\"<< endl;autoPtr< phaseChangeTwoPhaseMixture > mixture
Definition: createFields.H:33
volScalarField & alpha1
dimensionedScalar pMin("pMin", dimPressure, fluid)
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.
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU))
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
psi1
Definition: TEqns.H:34
volVectorField & HbyA
Definition: pEqn.H:13
fvVectorMatrix & UEqn
Definition: UEqn.H:13
const volScalarField & gh
volScalarField p_rgh_0(p_rgh)
const dimensionedScalar & rho2
Definition: createFields.H:40
phi
Definition: pEqn.H:18
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
const dimensionedScalar & rho1
Definition: createFields.H:39
rho
Definition: pEqn.H:1
void deleteDemandDrivenData(DataPtr &dataPtr)
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
volScalarField rAU(1.0/UEqn.A())
zeroField Sp
Definition: alphaSuSp.H:2