correctPressure.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-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 
26 #include "isothermalFluid.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 "fvcGrad.H"
34 #include "fvcSnGrad.H"
35 #include "fvcReconstruct.H"
36 #include "fvcVolumeIntegrate.H"
37 #include "fvmDiv.H"
38 #include "fvmLaplacian.H"
39 
40 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
41 
42 void Foam::solvers::isothermalFluid::correctPressure()
43 {
48 
49  const volScalarField& psi = thermo.psi();
50  rho = thermo.rho();
51  rho.relax();
52 
53  fvVectorMatrix& UEqn = tUEqn.ref();
54 
55  // Thermodynamic density needs to be updated by psi*d(p) after the
56  // pressure solution
57  const volScalarField psip0(psi*p);
58 
60 
61  const volScalarField rAU("rAU", 1.0/UEqn.A());
62  const surfaceScalarField rhorAUf("rhorAUf", fvc::interpolate(rho*rAU));
63 
64  tmp<volScalarField> rAtU
65  (
67  ? volScalarField::New("rAtU", 1.0/(1.0/rAU - UEqn.H1()))
68  : tmp<volScalarField>(nullptr)
69  );
70 
71  tmp<surfaceScalarField> rhorAtUf
72  (
74  ? surfaceScalarField::New("rhoRAtUf", fvc::interpolate(rho*rAtU()))
75  : tmp<surfaceScalarField>(nullptr)
76  );
77 
78  const volScalarField& rAAtU = pimple.consistent() ? rAtU() : rAU;
79  const surfaceScalarField& rhorAAtUf =
80  pimple.consistent() ? rhorAtUf() : rhorAUf;
81 
83 
84  if (pimple.nCorrPiso() <= 1)
85  {
86  tUEqn.clear();
87  }
88 
90  (
91  "phiHbyA",
92  rhof*fvc::flux(HbyA)
93  + rhorAUf*fvc::ddtCorr(rho, U, phi, rhoUf)
94  );
95 
96  MRF.makeRelative(rhof, phiHbyA);
97 
98  bool adjustMass = false;
99 
100  if (pimple.transonic())
101  {
102  const surfaceScalarField phidByPsi
103  (
105  (
106  fvc::relative(phiHbyA, rho, U)/rhof,
107  p
108  )
109  );
110 
111  const surfaceScalarField phid("phid", fvc::interpolate(psi)*phidByPsi);
112 
113  // Subtract the compressible part
114  // The resulting flux will be zero for a perfect gas
115  phiHbyA -= fvc::interpolate(psi*p)*phidByPsi;
116 
117  if (pimple.consistent())
118  {
119  phiHbyA += (rhorAAtUf - rhorAUf)*fvc::snGrad(p)*mesh.magSf();
120  HbyA += (rAAtU - rAU)*fvc::grad(p);
121  }
122 
123  // Update the pressure BCs to ensure flux consistency
124  constrainPressure(p, rho, U, phiHbyA, rhorAAtUf, MRF);
125 
127 
128  fvScalarMatrix pDDtEqn
129  (
131  + fvc::div(phiHbyA) + fvm::div(phid, p)
132  ==
133  fvModels().source(psi, p, rho.name())
134  );
135 
136  while (pimple.correctNonOrthogonal())
137  {
138  fvScalarMatrix pEqn(pDDtEqn - fvm::laplacian(rhorAAtUf, p));
139 
140  // Relax the pressure equation to ensure diagonal-dominance
141  pEqn.relax();
142 
143  pEqn.setReference
144  (
147  );
148 
149  fvConstraints().constrain(pEqn);
150 
151  pEqn.solve();
152 
154  {
155  phi = phiHbyA + pEqn.flux();
156  }
157  }
158  }
159  else
160  {
161  if (pimple.consistent())
162  {
163  phiHbyA += (rhorAAtUf - rhorAUf)*fvc::snGrad(p)*mesh.magSf();
164  HbyA += (rAAtU - rAU)*fvc::grad(p);
165  }
166 
167  // Update the pressure BCs to ensure flux consistency
168  constrainPressure(p, rho, U, phiHbyA, rhorAAtUf, MRF);
169 
171 
172  if (mesh.schemes().steady())
173  {
174  adjustMass = adjustPhi(phiHbyA, U, p);
175  }
176 
177  fvScalarMatrix pDDtEqn
178  (
180  + fvc::div(phiHbyA)
181  ==
182  fvModels().source(psi, p, rho.name())
183  );
184 
185  while (pimple.correctNonOrthogonal())
186  {
187  fvScalarMatrix pEqn(pDDtEqn - fvm::laplacian(rhorAAtUf, p));
188 
189  pEqn.setReference
190  (
193  );
194 
195  fvConstraints().constrain(pEqn);
196 
197  pEqn.solve();
198 
200  {
201  phi = phiHbyA + pEqn.flux();
202  }
203  }
204  }
205 
206  if (!mesh.schemes().steady())
207  {
208  const bool constrained = fvConstraints().constrain(p);
209 
210  // Thermodynamic density update
211  thermo_.correctRho(psi*p - psip0);
212 
213  if (constrained)
214  {
215  rho = thermo.rho();
216  }
217 
218  correctDensity();
219  }
220 
221  continuityErrors();
222 
223  // Explicitly relax pressure for momentum corrector
224  p.relax();
225 
226  U = HbyA - rAAtU*fvc::grad(p);
229  K = 0.5*magSqr(U);
230 
231  if (mesh.schemes().steady())
232  {
234  }
235 
236  // For steady compressible closed-volume cases adjust the pressure level
237  // to obey overall mass continuity
238  if (adjustMass && !thermo.incompressible())
239  {
243  }
244 
245  if (mesh.schemes().steady() || pimple.simpleRho() || adjustMass)
246  {
247  rho = thermo.rho();
248  }
249 
250  if (mesh.schemes().steady() || pimple.simpleRho())
251  {
252  rho.relax();
253  }
254 
255  // Correct rhoUf if the mesh is moving
257 
258  if (thermo.dpdt())
259  {
260  dpdt = fvc::ddt(p);
261  }
262 }
263 
264 
265 // ************************************************************************* //
For cases which do no have a pressure boundary adjust the balance of fluxes to obey continuity....
void relax(const scalar alpha)
Relax field (for steady-state solution).
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
void correctBoundaryConditions()
Correct boundary field.
const word & name() const
Return name.
Definition: IOobject.H:310
void makeRelative(volVectorField &U) const
Make the given absolute velocity relative within the MRF region.
Definition: MRFZoneList.C:157
virtual Switch dpdt() const =0
Should the dpdt term be included in the enthalpy equation.
virtual bool incompressible() const =0
Return true if the equation of state is incompressible.
virtual tmp< volScalarField > rho() const =0
Density [kg/m^3].
bool transonic() const
Flag to indicate to solve using transonic algorithm.
bool consistent() const
Flag to indicate to relax pressure using the "consistent".
virtual void correctRho(const volScalarField &deltaRho)=0
Add the given density correction to the density field.
virtual const volScalarField & psi() const =0
Compressibility [s^2/m^2].
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1750
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
bool steady() const
Return true if the default ddtScheme is steadyState.
Definition: fvSchemes.H:147
bool finalNonOrthogonalIter() const
Flag to indicate the last non-orthogonal iteration.
bool simpleRho() const
Switch to indicate whether to update the density in simple mode.
bool correctNonOrthogonal()
Non-orthogonal corrector loop.
label nCorrPiso() const
Maximum number of piso correctors.
Definition: pisoControlI.H:28
scalar refValue() const
Return the pressure reference level.
label refCell() const
Return the cell in which the reference pressure is set.
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 rho_
The continuity density field.
const surfaceScalarField & phi
Mass-flux field.
volVectorField U_
Velocity field.
tmp< fvVectorMatrix > tUEqn
Cached momentum matrix.
autoPtr< surfaceVectorField > rhoUf
Pointer to the surface momentum field.
fluidThermo & thermo_
Reference to the fluid thermophysical properties.
volScalarField K
Kinetic energy field.
const volVectorField & U
Velocity field.
IOMRFZoneList MRF
MRF zone list.
volScalarField & p_
Reference to the pressure field.
surfaceScalarField phi_
Mass-flux field.
const volScalarField & rho
Reference to the continuity density field.
volScalarField::Internal dpdt
Rate of change of the pressure.
Foam::pressureReference pressureReference
Pressure reference.
const fluidThermo & thermo
Reference to the fluid thermophysical properties.
const volScalarField & p
Reference to the pressure field.
dimensionedScalar initialMass
Initial mass in the region.
fvVectorMatrix & UEqn
Definition: UEqn.H:11
Calculate the first temporal derivative.
Calculate the face-flux of the given field.
Calculate the gradient 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.
Volume integrate volField creating a volField.
Calculate the matrix for the divergence of the given field and flux.
Calculate the matrix for the laplacian of the field.
const volScalarField & psi
volScalarField rAU(1.0/UEqn.A())
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
void correctRhoUf(autoPtr< surfaceVectorField > &rhoUf, const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const MRFType &MRF)
tmp< surfaceScalarField > relative(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given absolute flux in relative form.
Definition: fvcMeshPhi.C:167
dimensioned< Type > domainIntegrate(const VolField< Type > &vf)
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
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 > > 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)
bool adjustPhi(surfaceScalarField &phi, const volVectorField &U, volScalarField &p)
Adjust the balance of fluxes to obey continuity.
Definition: adjustPhi.C:35
VolField< vector > volVectorField
Definition: volFieldsFwd.H:62
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
SurfaceField< scalar > surfaceScalarField
tmp< surfaceScalarField > constrainPhid(const tmp< surfaceScalarField > &tphiHbyA, const volScalarField &p)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:45
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 > &)