pressureCorrector.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2022-2025 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "compressibleVoF.H"
27 #include "constrainHbyA.H"
28 #include "constrainPressure.H"
29 #include "adjustPhi.H"
30 #include "fvcMeshPhi.H"
31 #include "fvcFlux.H"
32 #include "fvcDdt.H"
33 #include "fvcSnGrad.H"
34 #include "fvcReconstruct.H"
35 #include "fvmDiv.H"
36 #include "fvmSup.H"
37 #include "fvmLaplacian.H"
38 
39 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
40 
42 {
44  volVectorField& U = U_;
46 
47  const volScalarField& rho1 = mixture.rho1();
48  const volScalarField& rho2 = mixture.rho2();
49 
50  const volScalarField& psi1 = mixture.thermo1().psi();
51  const volScalarField& psi2 = mixture.thermo2().psi();
52 
53  fvVectorMatrix& UEqn = tUEqn.ref();
54  setrAU(UEqn);
55 
56  const surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU()));
57 
58  while (pimple.correct())
59  {
62  (
63  "phiHbyA",
66  );
67 
69 
70  const surfaceScalarField phig
71  (
72  (
75  )*rAUf*mesh.magSf()
76  );
77 
78  phiHbyA += phig;
79 
80  // Update the pressure BCs to ensure flux consistency
82 
83  // Cache any sources
84  fvScalarMatrix p_rghEqnSource
85  (
86  fvModels().sourceProxy(alpha1, rho1, p_rgh)/rho1
87  + fvModels().sourceProxy(alpha2, rho2, p_rgh)/rho2
88  );
89 
90  // Make the fluxes relative to the mesh motion
92 
93  tmp<fvScalarMatrix> p_rghEqnComp1;
94  tmp<fvScalarMatrix> p_rghEqnComp2;
95 
96  if (pimple.transonic())
97  {
98  const surfaceScalarField rho1f(fvc::interpolate(rho1));
99  const surfaceScalarField rho2f(fvc::interpolate(rho2));
100 
101  surfaceScalarField phid1("phid1", fvc::interpolate(psi1)*phi);
102  surfaceScalarField phid2("phid2", fvc::interpolate(psi2)*phi);
103 
104  p_rghEqnComp1 =
105  (
106  (fvc::ddt(alpha1, rho1) + fvc::div(alphaPhi1*rho1f))/rho1
108  + (alpha1/rho1)
109  *correction
110  (
111  psi1*fvm::ddt(p_rgh)
112  + fvm::div(phid1, p_rgh) - fvm::Sp(fvc::div(phid1), p_rgh)
113  )
114  );
115 
116  p_rghEqnComp2 =
117  (
118  (fvc::ddt(alpha2, rho2) + fvc::div(alphaPhi2*rho2f))/rho2
120  + (alpha2/rho2)
121  *correction
122  (
123  psi2*fvm::ddt(p_rgh)
124  + fvm::div(phid2, p_rgh) - fvm::Sp(fvc::div(phid2), p_rgh)
125  )
126  );
127  }
128  else
129  {
130  const surfaceScalarField rho1f(fvc::interpolate(rho1));
131  const surfaceScalarField rho2f(fvc::interpolate(rho2));
132 
133  p_rghEqnComp1 =
134  (
135  (fvc::ddt(alpha1, rho1) + fvc::div(alphaPhi1*rho1f))/rho1
137  + (alpha1*psi1/rho1)*correction(fvm::ddt(p_rgh))
138  );
139 
140  p_rghEqnComp2 =
141  (
142  (fvc::ddt(alpha2, rho2) + fvc::div(alphaPhi2*rho2f))/rho2
144  + (alpha2*psi2/rho2)*correction(fvm::ddt(p_rgh))
145  );
146  }
147 
148  if (mesh.moving())
149  {
150  p_rghEqnComp1.ref() += fvc::div(mesh.phi())*alpha1;
151  p_rghEqnComp2.ref() += fvc::div(mesh.phi())*alpha2;
152  }
153 
154  p_rghEqnComp1.ref() *= pos(alpha1);
155  p_rghEqnComp2.ref() *= pos(alpha2);
156 
157  if (pimple.transonic())
158  {
159  p_rghEqnComp1.ref().relax();
160  p_rghEqnComp2.ref().relax();
161  }
162 
163  // Cache p_rgh prior to solve for density update
164  const volScalarField p_rgh_0(p_rgh);
165 
166  while (pimple.correctNonOrthogonal())
167  {
168  fvScalarMatrix p_rghEqnIncomp
169  (
171  == p_rghEqnSource
172  );
173 
174  {
175  fvScalarMatrix p_rghEqn
176  (
177  p_rghEqnComp1() + p_rghEqnComp2() + p_rghEqnIncomp
178  );
179 
180  fvConstraints().constrain(p_rghEqn);
181 
182  p_rghEqn.solve();
183  }
184 
186  {
187  vDot =
188  (
189  alpha1*(p_rghEqnComp2 & p_rgh)
190  - alpha2*(p_rghEqnComp1 & p_rgh)
191  );
192 
193  phi = phiHbyA + p_rghEqnIncomp.flux();
194 
199 
200  U = HbyA
201  + rAU()*fvc::reconstruct((phig + p_rghEqnIncomp.flux())/rAUf);
204  }
205  }
206 
207  // Update densities from change in p_rgh
208  mixture_.thermo1().correctRho(psi1*(p_rgh - p_rgh_0));
209  mixture_.thermo2().correctRho(psi2*(p_rgh - p_rgh_0));
210  mixture_.correct();
211 
212  // Correct p_rgh for consistency with p and the updated densities
215  }
216 
217  // Correct Uf if the mesh is moving
219 
220  K = 0.5*magSqr(U);
221 
222  clearrAU();
223  tUEqn.clear();
224 }
225 
226 
227 // ************************************************************************* //
For cases which do no have a pressure boundary adjust the balance of fluxes to obey continuity....
Generic GeometricField class.
void correctBoundaryConditions()
Correct boundary field.
void makeRelative(volVectorField &U) const
Make the given absolute velocity relative within the MRF region.
Definition: MRFZoneList.C:157
virtual void correct()
Update mixture properties.
const rhoFluidThermo & thermo1() const
Return the thermo for phase 1.
const rhoFluidThermo & thermo2() const
Return the thermo for phase 2.
bool transonic() const
Flag to indicate to solve using transonic algorithm.
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:58
tmp< SurfaceField< Type > > flux() const
Return the face-flux field from the matrix.
Definition: fvMatrix.C:964
const surfaceScalarField & phi() const
Return cell face motion fluxes.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
bool finalNonOrthogonalIter() const
Flag to indicate the last non-orthogonal iteration.
bool correct()
Piso loop within outer loop.
bool correctNonOrthogonal()
Non-orthogonal corrector loop.
bool moving() const
Is mesh moving.
Definition: polyMesh.H:473
virtual void correctRho(const volScalarField &deltaRho)
Add the given density correction to the density field.
Foam::fvModels & fvModels() const
Return the fvModels that are created on demand.
Definition: solver.C:96
pimpleNoLoopControl pimple
PIMPLE inner-loop controls.
Definition: solver.H:107
Foam::fvConstraints & fvConstraints() const
Return the fvConstraints that are created on demand.
Definition: solver.C:107
const fvMesh & mesh
Region mesh.
Definition: solver.H:101
const surfaceScalarField & phi
Reference to the mass-flux field.
Definition: VoFSolver.H:216
volVectorField U_
Velocity field.
Definition: VoFSolver.H:94
const volScalarField & p_rgh
Reference to the buoyant pressure for buoyant cases.
Definition: VoFSolver.H:210
tmp< fvVectorMatrix > tUEqn
Cached momentum matrix.
Definition: VoFSolver.H:145
autoPtr< surfaceVectorField > Uf
Pointer to the surface momentum field.
Definition: VoFSolver.H:133
autoPtr< volScalarField > rAU
Inverse momentum equation diagonal.
Definition: VoFSolver.H:129
const volVectorField & U
Reference to the velocity field.
Definition: VoFSolver.H:213
void clearrAU()
Clear the cached rAU is no longer needed.
Definition: VoFSolver.C:66
IOMRFZoneList MRF
MRF zone list.
Definition: VoFSolver.H:122
surfaceScalarField phi_
Volumetric flux field.
Definition: VoFSolver.H:97
const volScalarField & rho
Reference to the mixture continuity density field.
Definition: VoFSolver.H:110
volScalarField & p_rgh_
Reference to the buoyant pressure for buoyant cases.
Definition: VoFSolver.H:107
void setrAU(const fvVectorMatrix &UEqn)
Set or update the cached rAU.
Definition: VoFSolver.C:53
Buoyancy related data for the Foam::solvers::isothermalFluid solver module when solving buoyant cases...
Definition: buoyancy.H:70
volScalarField gh
(g & h) - ghRef
Definition: buoyancy.H:95
uniformDimensionedScalarField pRef
Reference pressure.
Definition: buoyancy.H:89
surfaceScalarField ghf
(g & hf) - ghRef
Definition: buoyancy.H:98
volScalarField & p
Reference to the mixture static pressure field.
volScalarField::Internal vDot
Compressibility source.
volScalarField K
Kinetic energy field.
virtual void pressureCorrector()
Construct and solve the pressure equation in the PISO loop.
compressibleTwoPhaseVoFMixture & mixture_
The compressible two-phase mixture.
surfaceScalarField alphaPhi1
volScalarField & alpha1
Reference to the phase1-fraction.
surfaceScalarField alphaPhi2
volScalarField & alpha2
Reference to the phase2-fraction.
virtual tmp< surfaceScalarField > surfaceTensionForce() const
Return the interface surface tension force for the momentum equation.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
fvVectorMatrix & UEqn
Definition: UEqn.H:11
Calculate the first temporal derivative.
Calculate the face-flux of the given field.
Calculate the mesh motion flux and convert fluxes from absolute to relative and back.
Reconstruct volField from a face flux field.
Calculate the snGrad of the given volField.
Calculate the matrix for the divergence of the given field and flux.
Calculate the matrix for the laplacian of the field.
Calculate the matrix for implicit and explicit sources.
volVectorField & HbyA
Definition: pEqn.H:13
surfaceScalarField phiHbyA("phiHbyA", fvc::interpolate(rho) *fvc::flux(HbyA))
tmp< SurfaceField< typename innerProduct< vector, Type >::type > > flux(const VolField< Type > &vf)
Return the face-flux field obtained from the given volVectorField.
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< VolField< Type > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
tmp< VolField< typename outerProduct< vector, Type >::type > > reconstruct(const SurfaceField< Type > &ssf)
void correctUf(autoPtr< surfaceVectorField > &Uf, const volVectorField &U, const surfaceScalarField &phi, const MRFType &MRF)
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:202
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:89
tmp< SurfaceField< Type > > snGrad(const VolField< Type > &vf, const word &name)
Definition: fvcSnGrad.C:45
tmp< SurfaceField< typename Foam::flux< Type >::type > > ddtCorr(const VolField< Type > &U, const SurfaceField< Type > &Uf)
Definition: fvcDdt.C:196
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvmLaplacian.C:47
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
void constrainPressure(volScalarField &p, const RhoType &rho, const volVectorField &U, const surfaceScalarField &phiHbyA, const RAUType &rhorAU, const MRFType &MRF)
dimensionedScalar pos(const dimensionedScalar &ds)
void constrainHbyA(volVectorField &HbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:34
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
void magSqr(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)