pEqn.H
Go to the documentation of this file.
1 {
2  if (rAU.valid())
3  {
4  rAU.ref() = 1.0/UEqn.A();
5  }
6  else
7  {
8  rAU = 1.0/UEqn.A();
9  }
10 
14  (
15  "phiHbyA",
17  + MRF.zeroFilter(fvc::interpolate(rho*rAU())*fvc::ddtCorr(U, phi, Uf))
18  );
19 
20  MRF.makeRelative(phiHbyA);
21 
23  (
24  (
25  mixture.surfaceTensionForce()
27  )*rAUf*mesh.magSf()
28  );
29 
31 
32  // Update the pressure BCs to ensure flux consistency
34 
35  // Cache the phase change pressure source
37 
38  // Make the fluxes relative to the mesh motion
40 
41  tmp<fvScalarMatrix> p_rghEqnComp1;
42  tmp<fvScalarMatrix> p_rghEqnComp2;
43 
44  if (pimple.transonic())
45  {
46  #include "rhofs.H"
47 
49  surfaceScalarField phid2("phid2", fvc::interpolate(psi2)*phi);
50 
51  p_rghEqnComp1 =
52  (
55  + (alpha1/rho1)
56  *correction
57  (
59  + fvm::div(phid1, p_rgh) - fvm::Sp(fvc::div(phid1), p_rgh)
60  )
61  );
62 
63  p_rghEqnComp2 =
64  (
67  + (alpha2/rho2)
68  *correction
69  (
71  + fvm::div(phid2, p_rgh) - fvm::Sp(fvc::div(phid2), p_rgh)
72  )
73  );
74  }
75  else
76  {
77  #include "rhofs.H"
78 
79  p_rghEqnComp1 =
80  (
84  );
85 
86  p_rghEqnComp2 =
87  (
91  );
92  }
93 
94  if (mesh.moving())
95  {
96  p_rghEqnComp1.ref() += fvc::div(mesh.phi())*alpha1;
97  p_rghEqnComp2.ref() += fvc::div(mesh.phi())*alpha2;
98  }
99 
100  p_rghEqnComp1.ref() *= pos(alpha1);
101  p_rghEqnComp2.ref() *= pos(alpha2);
102 
103  p_rghEqnComp1.ref() -=
104  (fvModels.source(alpha1, mixture.thermo1().rho())&rho1)/rho1;
105  p_rghEqnComp2.ref() -=
106  (fvModels.source(alpha2, mixture.thermo2().rho())&rho2)/rho2;
107 
108  if (pimple.transonic())
109  {
110  p_rghEqnComp1.ref().relax();
111  p_rghEqnComp2.ref().relax();
112  }
113 
114  // Cache p_rgh prior to solve for density update
116 
117  while (pimple.correctNonOrthogonal())
118  {
119  fvScalarMatrix p_rghEqnIncomp
120  (
122  == Sp_rgh
123  );
124 
125  solve
126  (
127  p_rghEqnComp1() + p_rghEqnComp2() + p_rghEqnIncomp
128  );
129 
130  if (pimple.finalNonOrthogonalIter())
131  {
132  dgdt =
133  (
134  alpha1*(p_rghEqnComp2 & p_rgh)
135  - alpha2*(p_rghEqnComp1 & p_rgh)
136  );
137 
138  phi = phiHbyA + p_rghEqnIncomp.flux();
139 
140  p = max(p_rgh + (alpha1*rho1 + alpha2*rho2)*gh, pMin);
141  p_rgh = p - (alpha1*rho1 + alpha2*rho2)*gh;
142  p_rgh.correctBoundaryConditions();
143 
144  U = HbyA
145  + rAU()*fvc::reconstruct((phig + p_rghEqnIncomp.flux())/rAUf);
146  U.correctBoundaryConditions();
148  }
149  }
150 
151  // Correct Uf if the mesh is moving
153 
154  // Update densities from change in p_rgh
155  mixture.thermo1().correctRho(psi1*(p_rgh - p_rgh_0));
156  mixture.thermo2().correctRho(psi2*(p_rgh - p_rgh_0));
157  mixture.correct();
158 
159  // Correct p_rgh for consistency with p and the updated densities
160  p_rgh = p - rho*gh;
161  p_rgh.correctBoundaryConditions();
162 
163  K = 0.5*magSqr(U);
164 }
surfaceScalarField alphaPhi1(alphaPhi1Header, phi *fvc::interpolate(alpha1))
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
rAU
Definition: pEqn.H:1
const volScalarField & rho1
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
tmp< fvScalarMatrix > p_rghEqnComp1
Definition: pEqn.H:41
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
pimpleNoLoopControl & pimple
p_rgh
Definition: pEqn.H:160
U
Definition: pEqn.H:72
surfaceScalarField rho1f(fvc::interpolate(rho1))
mixture MRF().makeRelative(phiHbyA)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
dimensionedScalar pMin("pMin", dimPressure, mixture)
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
MRF makeRelative(fvc::interpolate(rho), phiHbyA)
volScalarField & alpha1(mixture.alpha1())
twoPhaseChangeModel & phaseChange
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:59
fvMesh & mesh
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
p
Definition: pEqn.H:125
alpha2
Definition: alphaEqn.H:115
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:58
dimensionedScalar pos(const dimensionedScalar &ds)
tmp< fvScalarMatrix > p_rghEqnComp2
Definition: pEqn.H:42
surfaceScalarField rho2f(fvc::interpolate(rho2))
phiHbyA
Definition: pEqn.H:30
Foam::fvConstraints & fvConstraints
volScalarField p_rgh_0(p_rgh)
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:34
autoPtr< surfaceVectorField > Uf
dimensioned< scalar > magSqr(const dimensioned< Type > &)
const surfaceScalarField & ghf
const volScalarField & psi2
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.
const tmp< volScalarField::Internal > & Sp
Definition: alphaSuSp.H:7
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
Foam::fvModels & fvModels
phi
Definition: pEqn.H:98
volScalarField::Internal dgdt(IOobject("dgdt", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), alpha1 *fvc::div(phi))
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
rhoEqn solve()
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:202
const volScalarField & rho2
volVectorField & HbyA
Definition: pEqn.H:13
const surfaceScalarField phig(-rhorAUf *ghf *fvc::snGrad(rho) *mesh.magSf())
fvScalarMatrix Sp_rgh(phaseChange.Sp_rgh(rho, gh, p_rgh))
surfaceScalarField alphaPhi2("alphaPhi2", phi - alphaPhi1)
void correctUf(autoPtr< surfaceVectorField > &Uf, const volVectorField &U, const surfaceScalarField &phi, const MRFType &MRF)
const volScalarField & gh
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU))
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 volScalarField & psi1
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
constrainPressure(p, rho, U, phiHbyA, rhorAUf, MRF)
rho
Definition: pEqn.H:1
fvVectorMatrix & UEqn
Definition: UEqn.H:13