pEqn.H
Go to the documentation of this file.
1 {
2  volScalarField rAU("rAU", 1.0/UEqn.A());
3 
6 
8  (
9  "phiHbyA",
11  );
12 
13  if (ddtCorr)
14  {
15  surfaceScalarField faceMaskOld
16  (
17  localMin<scalar>(mesh).interpolate(cellMask.oldTime())
18  );
19  phiHbyA +=
20  MRF.zeroFilter
21  (
22  fvc::interpolate(rho*rAU)*faceMaskOld*fvc::ddtCorr(U, Uf)
23  );
24  }
25 
26  MRF.makeRelative(phiHbyA);
27 
29  (
30  (
31  mixture.surfaceTensionForce()
33  )*faceMask*rAUf*mesh.magSf()
34  );
35 
37 
38  // Update the pressure BCs to ensure flux consistency
40 
41  tmp<fvScalarMatrix> p_rghEqnComp1;
42  tmp<fvScalarMatrix> p_rghEqnComp2;
43 
44  if (pimple.transonic())
45  {
46  #include "rhofs.H"
47 
50 
52  pos(alpha1)
53  *(
54  (
56  - (fvOptions(alpha1, mixture.thermo1().rho())&rho1)
57  )/rho1
59  + (alpha1/rho1)
60  *correction
61  (
63  + fvm::div(phid1, p_rgh) - fvm::Sp(fvc::div(phid1), p_rgh)
64  )
65  );
66  p_rghEqnComp1.ref().relax();
67 
69  pos(alpha2)
70  *(
71  (
73  - (fvOptions(alpha2, mixture.thermo2().rho())&rho2)
74  )/rho2
76  + (alpha2/rho2)
77  *correction
78  (
80  + fvm::div(phid2, p_rgh) - fvm::Sp(fvc::div(phid2), p_rgh)
81  )
82  );
83  p_rghEqnComp2.ref().relax();
84  }
85  else
86  {
87  #include "rhofs.H"
88 
90  pos(alpha1)
91  *(
92  (
94  - (fvOptions(alpha1, mixture.thermo1().rho())&rho1)
95  )/rho1
98  );
99 
100  p_rghEqnComp2 =
101  pos(alpha2)
102  *(
103  (
105  - (fvOptions(alpha2, mixture.thermo2().rho())&rho2)
106  )/rho2
109  );
110  }
111 
112  // Cache p_rgh prior to solve for density update
114 
115  while (pimple.correctNonOrthogonal())
116  {
117  fvScalarMatrix p_rghEqnIncomp
118  (
121  );
122 
123  solve
124  (
125  p_rghEqnComp1() + p_rghEqnComp2() + p_rghEqnIncomp,
126  mesh.solver(p_rgh.select(pimple.finalInnerIter()))
127  );
128 
129  if (pimple.finalNonOrthogonalIter())
130  {
131  p = max(p_rgh + (alpha1*rho1 + alpha2*rho2)*gh, pMin);
132  p_rgh = p - (alpha1*rho1 + alpha2*rho2)*gh;
133 
134  dgdt =
135  (
138  );
139 
140  phi = phiHbyA + p_rghEqnIncomp.flux();
141 
142  U =
143  cellMask*
144  (
145  HbyA
146  + rAU*fvc::reconstruct((phig + p_rghEqnIncomp.flux())/rAUf)
147  );
148 
149  U.correctBoundaryConditions();
150  fvOptions.correct(U);
151  }
152  }
153 
154  {
155  Uf = fvc::interpolate(U);
156  surfaceVectorField n(mesh.Sf()/mesh.magSf());
157  Uf += n*(fvc::absolute(phi, U)/mesh.magSf() - (n & Uf));
158  }
159 
160  // Make the fluxes relative to the mesh motion
162 
163  // Zero faces H-I for transport Eq after pEq
164  phi *= faceMask;
165 
166  // Update densities from change in p_rgh
167  mixture.thermo1().correctRho(psi1*(p_rgh - p_rgh_0));
168  mixture.thermo2().correctRho(psi2*(p_rgh - p_rgh_0));
169 
171 
172  // Correct p_rgh for consistency with p and the updated densities
173  p_rgh = p - rho*gh;
174  p_rgh.correctBoundaryConditions();
175 
176  K = 0.5*magSqr(U);
177 }
Foam::fvc::reconstruct
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcReconstruct.C:49
Foam::constrainHbyA
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:28
Foam::fvc::snGrad
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:40
ghf
const surfaceScalarField & ghf
Definition: setRegionFluidFields.H:18
pMin
const dimensionedScalar & pMin
Definition: setRegionFluidFields.H:58
alpha2
const volScalarField & alpha2
Definition: setRegionFluidFields.H:9
p_rghEqnComp1
tmp< fvScalarMatrix > p_rghEqnComp1
Definition: pEqn.H:29
phiHbyA
phiHbyA
Definition: pEqn.H:20
Sp
zeroField Sp
Definition: alphaSuSp.H:2
faceMask
surfaceScalarField faceMask(localMin< scalar >(mesh).interpolate(cellMask))
Foam::correction
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
alpha1
const volScalarField & alpha1
Definition: setRegionFluidFields.H:8
gh
const volScalarField & gh
Definition: setRegionFluidFields.H:17
Foam::fac::div
tmp< GeometricField< Type, faPatchField, areaMesh > > div(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facDiv.C:43
psi2
const volScalarField & psi2
Definition: setRegionFluidFields.H:31
rho2f
surfaceScalarField rho2f(fvc::interpolate(rho2))
solve
CEqn solve()
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
K
CGAL::Exact_predicates_exact_constructions_kernel K
Definition: CGALTriangulation3DKernel.H:52
MRF
IOMRFZoneList & MRF
Definition: setRegionFluidFields.H:22
pimple
pimpleControl & pimple
Definition: setRegionFluidFields.H:56
HbyA
HbyA
Definition: pEqn.H:4
Foam::fvScalarMatrix
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:38
fvOptions
fv::options & fvOptions
Definition: setRegionFluidFields.H:23
rho2
volScalarField & rho2
Definition: setRegionFluidFields.H:30
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:53
makeRelative
MRF makeRelative(fvc::interpolate(rho), phiHbyA)
ddtCorr
ddtCorr
Definition: readControls.H:9
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Definition: hashSets.C:40
rho1
volScalarField & rho1
Definition: setRegionFluidFields.H:27
p
p
Definition: pEqn.H:50
rho1f
surfaceScalarField rho1f(fvc::interpolate(rho1))
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
alphaPhi1
const surfaceScalarField & alphaPhi1
Definition: setRegionFluidFields.H:13
psi1
const volScalarField & psi1
Definition: setRegionFluidFields.H:28
rAU
volScalarField rAU(1.0/UEqn.A())
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:58
U
U
Definition: pEqn.H:72
phig
surfaceScalarField phig("phig", -rhorAUf *ghf *fvc::snGrad(rho) *mesh.magSf())
p_rgh_0
volScalarField p_rgh_0(p_rgh)
p_rghEqnComp2
tmp< fvScalarMatrix > p_rghEqnComp2
Definition: pEqn.H:30
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:50
Foam::fac::ddt
tmp< GeometricField< Type, faPatchField, areaMesh > > ddt(const dimensioned< Type > dt, const faMesh &mesh)
Definition: facDdt.C:40
phi
phi
Definition: pEqn.H:18
UEqn
fvVectorMatrix & UEqn
Definition: UEqn.H:13
Uf
Uf
Definition: pEqn.H:107
p_rgh
p_rgh
Definition: pEqn.H:139
Foam::fac::laplacian
tmp< GeometricField< Type, faPatchField, areaMesh > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf, const word &name)
Definition: facLaplacian.C:40
Foam::surfaceVectorField
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
Definition: surfaceFieldsFwd.H:55
Foam::fvc::flux
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Definition: fvcFlux.C:27
rho
rho
Definition: pEqn.H:1
mixture
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\n"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture
Definition: createFields.H:39
Foam::fvc::absolute
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Definition: fvcMeshPhi.C:183
interpolate
mesh interpolate(rAU)
n
surfaceVectorField n(mesh.Sf()/mesh.magSf())
alphaPhi2
const surfaceScalarField & alphaPhi2
Definition: setRegionFluidFields.H:17
constrainPressure
constrainPressure(p_rgh, rho, U, phiHbyA, rhorAUf, MRF)
rAUf
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU))
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:170