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 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 {
44  volVectorField& U = U_;
46 
47  fvVectorMatrix& UEqn = tUEqn.ref();
48  setrAU(UEqn);
49 
50  const surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU()));
51 
52  while (pimple.correct())
53  {
56  (
57  "phiHbyA",
60  );
61 
63 
64  const surfaceScalarField phig
65  (
66  (
69  )*rAUf*mesh.magSf()
70  );
71 
72  phiHbyA += phig;
73 
74  // Update the pressure BCs to ensure flux consistency
76 
77  PtrList<fvScalarMatrix> p_rghEqnComps(phases.size());
78 
79  forAll(phases, phasei)
80  {
81  const compressibleVoFphase& phase = phases[phasei];
82  const rhoThermo& thermo = phase.thermo();
83  const volScalarField& rho = phases[phasei].thermo().rho();
84 
85  p_rghEqnComps.set
86  (
87  phasei,
88  (
91  - (fvModels().source(phase, rho)&rho)
92  ).ptr()
93  );
94  }
95 
96  // Cache p_rgh prior to solve for density update
97  const volScalarField p_rgh_0(p_rgh);
98 
100  {
101  fvScalarMatrix p_rghEqnIncomp
102  (
104  - fvm::laplacian(rAUf, p_rgh)
105  );
106 
107  tmp<fvScalarMatrix> p_rghEqnComp;
108 
109  forAll(phases, phasei)
110  {
111  const compressibleVoFphase& phase = phases[phasei];
112 
113  tmp<fvScalarMatrix> p_rghEqnCompi
114  (
115  (max(phase, scalar(0))/phase.thermo().rho())
116  *p_rghEqnComps[phasei]
117  );
118 
119  if (phasei == 0)
120  {
121  p_rghEqnComp = p_rghEqnCompi;
122  }
123  else
124  {
125  p_rghEqnComp.ref() += p_rghEqnCompi;
126  }
127  }
128 
129  {
130  fvScalarMatrix p_rghEqn(p_rghEqnComp + p_rghEqnIncomp);
131 
132  fvConstraints().constrain(p_rghEqn);
133 
134  p_rghEqn.solve();
135  }
136 
138  {
139  forAll(phases, phasei)
140  {
141  compressibleVoFphase& phase = phases[phasei];
142 
143  phase.dgdt() =
144  pos0(phase)
145  *(p_rghEqnComps[phasei] & p_rgh)/phase.thermo().rho();
146  }
147 
148  phi = phiHbyA + p_rghEqnIncomp.flux();
149 
150  p = p_rgh + rho*buoyancy.gh;
152  p_rgh = p - rho*buoyancy.gh;
154 
155  U = HbyA
156  + rAU()*fvc::reconstruct((phig + p_rghEqnIncomp.flux())/rAUf);
159  }
160  }
161 
162  // Update densities from change in p_rgh
163  mixture.correctRho(p_rgh - p_rgh_0);
164  mixture.correct();
165 
166  // Correct p_rgh for consistency with p and the updated densities
167  p_rgh = p - rho*buoyancy.gh;
169  }
170 
171  // Correct Uf if the mesh is moving
173 
174  K = 0.5*magSqr(U);
175 
176  clearrAU();
177  tUEqn.clear();
178 }
179 
180 
181 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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:65
Single compressible phase derived from the VoFphase.
const volScalarField::Internal & dgdt() const
Return const-access to phase divergence.
const rhoThermo & thermo() const
Return const-access to phase rhoThermo.
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:963
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.
Definition: rhoThermo.H:55
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:85
pimpleNoLoopControl pimple
PIMPLE inner-loop controls.
Definition: solver.H:100
Foam::fvConstraints & fvConstraints() const
Return the fvConstraints that are created on demand.
Definition: solver.C:96
const fvMesh & mesh
Region mesh.
Definition: solver.H:94
volScalarField & p_rgh
Reference to the buoyant pressure for buoyant cases.
Definition: VoFSolver.H:107
const surfaceScalarField & phi
Reference to the mass-flux field.
Definition: VoFSolver.H:212
volVectorField U_
Velocity field.
Definition: VoFSolver.H:94
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:209
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
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
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:181
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)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:34
dimensioned< scalar > magSqr(const dimensioned< Type > &)
fluidMulticomponentThermo & thermo
Definition: createFields.H:31