pEqn.H
Go to the documentation of this file.
1 {
2  volScalarField rAU("rAU", 1.0/UEqn.A());
6  (
7  "phiHbyA",
10  );
11 
13  (
14  (
15  interface.surfaceTensionForce()
17  )*rAUf*mesh.magSf()
18  );
19 
21 
22  // Update the pressure BCs to ensure flux consistency
24 
25  // Make the fluxes relative to the mesh motion
27 
28  tmp<fvScalarMatrix> p_rghEqnComp1;
29  tmp<fvScalarMatrix> p_rghEqnComp2;
30 
31  if (pimple.transonic())
32  {
34  surfaceScalarField phid2("phid2", fvc::interpolate(psi2)*phi);
35 
36  p_rghEqnComp1 =
37  fvc::ddt(rho1) + fvc::div(phi, rho1) - fvc::Sp(fvc::div(phi), rho1)
38  + correction
39  (
41  + fvm::div(phid1, p_rgh) - fvm::Sp(fvc::div(phid1), p_rgh)
42  );
43  deleteDemandDrivenData(p_rghEqnComp1.ref().faceFluxCorrectionPtr());
44  p_rghEqnComp1.ref().relax();
45 
46  p_rghEqnComp2 =
47  fvc::ddt(rho2) + fvc::div(phi, rho2) - fvc::Sp(fvc::div(phi), rho2)
48  + correction
49  (
51  + fvm::div(phid2, p_rgh) - fvm::Sp(fvc::div(phid2), p_rgh)
52  );
53  deleteDemandDrivenData(p_rghEqnComp2.ref().faceFluxCorrectionPtr());
54  p_rghEqnComp2.ref().relax();
55  }
56  else
57  {
58  p_rghEqnComp1 =
61 
62  p_rghEqnComp2 =
64  + fvc::div(phi, rho2) - fvc::Sp(fvc::div(phi), rho2);
65  }
66 
67  // Cache p_rgh prior to solve for density update
69 
70  while (pimple.correctNonOrthogonal())
71  {
72  fvScalarMatrix p_rghEqnIncomp
73  (
76  );
77 
78  solve
79  (
80  (
81  (max(alpha1, scalar(0))/rho1)*p_rghEqnComp1()
82  + (max(alpha2, scalar(0))/rho2)*p_rghEqnComp2()
83  )
84  + p_rghEqnIncomp,
85  mesh.solver(p_rgh.select(pimple.finalInnerIter()))
86  );
87 
88  if (pimple.finalNonOrthogonalIter())
89  {
90  p = max(p_rgh + (alpha1*rho1 + alpha2*rho2)*gh, pMin);
91  p_rgh = p - (alpha1*rho1 + alpha2*rho2)*gh;
92 
93  dgdt =
94  (
95  pos(alpha2)*(p_rghEqnComp2 & p_rgh)/rho2
96  - pos(alpha1)*(p_rghEqnComp1 & p_rgh)/rho1
97  );
98 
99  phi = phiHbyA + p_rghEqnIncomp.flux();
100 
101  U = HbyA
102  + rAU*fvc::reconstruct((phig + p_rghEqnIncomp.flux())/rAUf);
103  U.correctBoundaryConditions();
104  }
105  }
106 
107  {
108  Uf = fvc::interpolate(U);
109  surfaceVectorField n(mesh.Sf()/mesh.magSf());
110  Uf += n*(fvc::absolute(phi, U)/mesh.magSf() - (n & Uf));
111  }
112 
113  // Update densities from change in p_rgh
115  rho2 += psi2*(p_rgh - p_rgh_0);
116 
118 
119  // Correct p_rgh for consistency with p and the updated densities
120  p_rgh = p - rho*gh;
121  p_rgh.correctBoundaryConditions();
122 
123  K = 0.5*magSqr(U);
124 
125  Info<< "max(U) " << max(mag(U)).value() << endl;
126  Info<< "min(p_rgh) " << min(p_rgh).value() << endl;
127 }
tmp< fvScalarMatrix > p_rghEqnComp2
Definition: pEqn.H:29
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
psi2
Definition: TEqns.H:35
surfaceVectorField n(mesh.Sf()/mesh.magSf())
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
U
Definition: pEqn.H:83
p
Definition: pEqn.H:50
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
interfaceProperties interface(alpha1, U, mixture())
phiHbyA
Definition: pEqn.H:20
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
tmp< GeometricField< typename flux< Type >::type, fvsPatchField, surfaceMesh > > ddtCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
Definition: fvcDdt.C:155
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
const dictionary & pimple
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
Uf
Definition: pEqn.H:67
tmp< GeometricField< Type, fvPatchField, volMesh > > Sp(const volScalarField &sp, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:67
tmp< surfaceScalarField > interpolate(const RhoType &rho)
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:28
rho1
Definition: pEqn.H:114
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)
dynamicFvMesh & mesh
rhoEqn solve()
rho2
Definition: pEqn.H:115
surfaceScalarField phig("phig",-rhorAUf *ghf *fvc::snGrad(rho)*mesh.magSf())
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:120
dimensioned< scalar > magSqr(const dimensioned< Type > &)
constrainPressure(p, rho, U, phiHbyA, rhorAUf, MRF)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
volScalarField & alpha1
dimensionedScalar pMin("pMin", dimPressure, fluid)
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
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:188
volVectorField & HbyA
Definition: pEqn.H:13
fvVectorMatrix & UEqn
Definition: UEqn.H:13
const volScalarField & gh
volScalarField p_rgh_0(p_rgh)
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
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
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())
MRF makeRelative(fvc::interpolate(rho), phiHbyA)