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) 2023-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 
27 #include "constrainHbyA.H"
28 #include "constrainPressure.H"
29 #include "adjustPhi.H"
30 #include "findRefCell.H"
31 #include "fvcMeshPhi.H"
32 #include "fvcFlux.H"
33 #include "fvcDdt.H"
34 #include "fvcDiv.H"
35 #include "fvcSnGrad.H"
36 #include "fvcSup.H"
37 #include "fvcReconstruct.H"
38 #include "fvmLaplacian.H"
39 
40 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
41 
43 {
45  volVectorField& U = U_;
47 
48  fvVectorMatrix& UEqn = tUEqn.ref();
49  setrAU(UEqn);
50 
51  const surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU()));
52 
53  while (pimple.correct())
54  {
57  (
58  "phiHbyA",
61  );
62 
64 
65  const surfaceScalarField phig
66  (
67  (
70  )*rAUf*mesh.magSf()
71  );
72 
73  phiHbyA += phig;
74 
75  // Update the pressure BCs to ensure flux consistency
77 
78  PtrList<fvScalarMatrix> p_rghEqnComps(phases.size());
79 
80  forAll(phases, phasei)
81  {
82  const compressibleVoFphase& phase = phases[phasei];
83  const rhoFluidThermo& thermo = phase.thermo();
84  const volScalarField& rho = phases[phasei].thermo().rho();
85 
86  p_rghEqnComps.set
87  (
88  phasei,
89  (
92  - fvModels().sourceProxy(phase, rho, p_rgh)
93  ).ptr()
94  );
95  }
96 
97  // Cache p_rgh prior to solve for density update
98  const volScalarField p_rgh_0(p_rgh);
99 
100  while (pimple.correctNonOrthogonal())
101  {
102  fvScalarMatrix p_rghEqnIncomp
103  (
105  - fvm::laplacian(rAUf, p_rgh)
106  );
107 
108  tmp<fvScalarMatrix> p_rghEqnComp;
109 
110  forAll(phases, phasei)
111  {
112  const compressibleVoFphase& phase = phases[phasei];
113 
114  tmp<fvScalarMatrix> p_rghEqnCompi
115  (
116  (max(phase, scalar(0))/phase.thermo().rho())
117  *p_rghEqnComps[phasei]
118  );
119 
120  if (phasei == 0)
121  {
122  p_rghEqnComp = p_rghEqnCompi;
123  }
124  else
125  {
126  p_rghEqnComp.ref() += p_rghEqnCompi;
127  }
128  }
129 
130  {
131  fvScalarMatrix p_rghEqn(p_rghEqnComp + p_rghEqnIncomp);
132 
133  fvConstraints().constrain(p_rghEqn);
134 
135  p_rghEqn.solve();
136  }
137 
139  {
140  forAll(phases, phasei)
141  {
142  compressibleVoFphase& phase = phases[phasei];
143 
144  phase.vDot() =
145  pos0(phase)
146  *(p_rghEqnComps[phasei] & p_rgh)/phase.thermo().rho();
147  }
148 
149  phi = phiHbyA + p_rghEqnIncomp.flux();
150 
155 
156  U = HbyA
157  + rAU()*fvc::reconstruct((phig + p_rghEqnIncomp.flux())/rAUf);
160  }
161  }
162 
163  // Update densities from change in p_rgh
164  mixture.correctRho(p_rgh - p_rgh_0);
165  mixture.correct();
166 
167  // Correct p_rgh for consistency with p and the updated densities
170  }
171 
172  // Correct Uf if the mesh is moving
174 
175  K = 0.5*magSqr(U);
176 
177  clearrAU();
178  tUEqn.clear();
179 }
180 
181 
182 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
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
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
bool set(const label) const
Is element set.
Definition: PtrListI.H:62
Single compressible phase derived from the VoFphase.
const volScalarField::Internal & vDot() const
Return const-access to phase divergence.
const rhoFluidThermo & thermo() const
Return const-access to phase rhoFluidThermo.
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 & 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.
Base-class for fluid thermodynamic properties based on density.
virtual tmp< volScalarField > rho() const =0
Density [kg/m^3].
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
UPtrListDictionary< compressibleVoFphase > & phases
Reference to the phases.
volScalarField & p
Reference to the mixture static pressure field.
virtual void pressureCorrector()
Construct and solve the pressure equation in the PISO loop.
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
Find the reference cell nearest (in index) to the given cell but which is not on a cyclic,...
Calculate the first temporal derivative.
Calculate the divergence of the given field.
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 field for explicit evaluation of implicit and explicit sources.
Calculate the matrix for the laplacian of the field.
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 > > Sp(const volScalarField &sp, const VolField< Type > &vf)
Definition: fvcSup.C:67
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
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 > > 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 pos0(const dimensionedScalar &ds)
void constrainHbyA(volVectorField &HbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:34
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
void magSqr(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
fluidMulticomponentThermo & thermo
Definition: createFields.H:31