pEqn.H
Go to the documentation of this file.
1 {
2  if (pimple.nCorrPIMPLE() == 1)
3  {
4  p =
5  (
6  rho
7  - alphal*rhol0
8  - ((alphav*psiv + alphal*psil) - psi)*pSat
9  )/psi;
10  }
11 
13 
14  volScalarField rAU(1.0/UEqn.A());
16 
17  volVectorField HbyA("HbyA", U);
18  HbyA = rAU*UEqn.H();
19 
22 
24 
25  phi -= phiGradp/rhof;
26 
27  while (pimple.correctNonOrthogonal())
28  {
29  fvScalarMatrix pEqn
30  (
31  fvm::ddt(psi, p)
32  - (rhol0 + (psil - psiv)*pSat)*fvc::ddt(alphav) - pSat*fvc::ddt(psi)
33  + fvc::div(phi, rho)
36  );
37 
38  pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter())));
39 
40  if (pimple.finalNonOrthogonalIter())
41  {
42  phi += (phiGradp + pEqn.flux())/rhof;
43  }
44  }
45 
46  Info<< "Predicted p max-min : " << max(p).value()
47  << " " << min(p).value() << endl;
48 
49  rho == max
50  (
51  psi*p
52  + alphal*rhol0
53  + ((alphav*psiv + alphal*psil) - psi)*pSat,
54  rhoMin
55  );
56 
57  #include "alphavPsi.H"
58 
59  p =
60  (
61  rho
62  - alphal*rhol0
63  - ((alphav*psiv + alphal*psil) - psi)*pSat
64  )/psi;
65 
66  p.correctBoundaryConditions();
67 
68  Info<< "Phase-change corrected p max-min : " << max(p).value()
69  << " " << min(p).value() << endl;
70 
71  // Correct velocity
72 
73  U = HbyA - rAU*fvc::grad(p);
74 
75  // Remove the swirl component of velocity for "wedge" cases
76  if (pimple.dict().found("removeSwirl"))
77  {
78  label swirlCmpt(readLabel(pimple.dict().lookup("removeSwirl")));
79 
80  Info<< "Removing swirl component-" << swirlCmpt << " of U" << endl;
81  U.field().replace(swirlCmpt, 0.0);
82  }
83 
84  U.correctBoundaryConditions();
85 
86  Info<< "max(U) " << max(mag(U)).value() << endl;
87 }
Foam::fvc::snGrad
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
phiGradp
surfaceScalarField phiGradp(rhorAUf *mesh.magSf() *fvc::snGrad(p))
p
p
Definition: pEqn.H:62
Foam::MULES::interpolate
tmp< surfaceScalarField > interpolate(const RhoType &rho)
Definition: IMULESTemplates.C:40
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:261
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
phi
phi
Definition: pEqn.H:20
alphavPsi.H
U
U
Definition: pEqn.H:46
Foam::fvc::laplacian
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
rhoMin
PtrList< dimensionedScalar > rhoMin(fluidRegions.size())
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::fvScalarMatrix
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
Foam::Info
messageStream Info
Foam::fvc::ddtCorr
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
UEqn
tmp< fvVectorMatrix > UEqn(fvm::div(phi, U)+MRF.DDt(U)+turbulence->divDevReff(U)==fvOptions(U))
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
pSat
const dimensionedScalar & pSat
Definition: createFields.H:41
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
rhorAUf
surfaceScalarField rhorAUf("rhorAUf", fvc::interpolate(rho *rAU))
rho
rho
Definition: pEqn.H:3
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
rAU
volScalarField rAU("rAU", 1.0/UEqn().A())
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
HbyA
HbyA
Definition: pEqn.H:4
psi
const volScalarField & psi
Definition: setRegionFluidFields.H:13
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:52
alphal
alphal
Definition: alphavPsi.H:12
Foam::readLabel
label readLabel(Istream &is)
Definition: label.H:64
pimple
const dictionary & pimple
Definition: readFluidMultiRegionPIMPLEControls.H:1
Foam::fvc::ddt
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
rhof
surfaceScalarField rhof(fvc::interpolate(rho, "div(phi,rho)"))
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)