correctBuoyantPressure.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 "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 "fvcSnGrad.H"
34 #include "fvcReconstruct.H"
35 #include "fvcVolumeIntegrate.H"
36 #include "fvmDiv.H"
37 #include "fvmLaplacian.H"
38 
39 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
40 
41 void Foam::solvers::isothermalFluid::correctBuoyantPressure()
42 {
48 
49  // Local references to the buoyancy parameters
50  const volScalarField& gh = buoyancy->gh;
51  const surfaceScalarField& ghf = buoyancy->ghf;
52  const uniformDimensionedScalarField pRef = buoyancy->pRef;
53 
54  const volScalarField& psi = thermo.psi();
55  rho = thermo.rho();
56  rho.relax();
57 
58  fvVectorMatrix& UEqn = tUEqn.ref();
59 
60  // Thermodynamic density needs to be updated by psi*d(p) after the
61  // pressure solution
62  const volScalarField psip0(psi*p);
63 
65 
66  const volScalarField rAU("rAU", 1.0/UEqn.A());
67  const surfaceScalarField rhorAUf("rhorAUf", fvc::interpolate(rho*rAU));
68 
69  tmp<volScalarField> rAtU
70  (
72  ? volScalarField::New("rAtU", 1.0/(1.0/rAU - UEqn.H1()))
73  : tmp<volScalarField>(nullptr)
74  );
75 
76  tmp<surfaceScalarField> rhorAtUf
77  (
79  ? surfaceScalarField::New("rhoRAtUf", fvc::interpolate(rho*rAtU()))
80  : tmp<surfaceScalarField>(nullptr)
81  );
82 
83  const volScalarField& rAAtU = pimple.consistent() ? rAtU() : rAU;
84  const surfaceScalarField& rhorAAtUf =
85  pimple.consistent() ? rhorAtUf() : rhorAUf;
86 
88 
89  if (pimple.nCorrPiso() <= 1)
90  {
91  tUEqn.clear();
92  }
93 
95  (
96  "phiHbyA",
97  rhof*fvc::flux(HbyA)
98  + rhorAUf*fvc::ddtCorr(rho, U, phi, rhoUf)
99  );
100 
101  MRF.makeRelative(rhof, phiHbyA);
102 
103  const bool adjustMass =
105 
106  const surfaceScalarField ghGradRhof(-ghf*fvc::snGrad(rho)*mesh.magSf());
107 
108  phiHbyA += rhorAUf*ghGradRhof;
109 
110  tmp<fvScalarMatrix> tp_rghEqn;
111 
112  if (pimple.transonic())
113  {
114  const surfaceScalarField phidByPsi
115  (
117  (
118  fvc::relative(phiHbyA, rho, U)/rhof,
119  p_rgh
120  )
121  );
122 
123  const surfaceScalarField phid("phid", fvc::interpolate(psi)*phidByPsi);
124 
125  // Subtract the compressible part
126  // The resulting flux will be zero for a perfect gas
127  phiHbyA -= fvc::interpolate(psi*p_rgh)*phidByPsi;
128 
129  if (pimple.consistent())
130  {
131  const surfaceScalarField gradpf(fvc::snGrad(p_rgh)*mesh.magSf());
132  phiHbyA += (rhorAAtUf - rhorAUf)*gradpf;
133  HbyA += (rAAtU - rAU)*fvc::reconstruct(gradpf - ghGradRhof);
134  }
135 
136  // Update the pressure BCs to ensure flux consistency
137  constrainPressure(p_rgh, rho, U, phiHbyA, rhorAAtUf, MRF);
138 
140 
141  fvScalarMatrix p_rghDDtEqn
142  (
144  + fvc::div(phiHbyA) + fvm::div(phid, p_rgh)
145  ==
146  fvModels().sourceProxy(rho, p_rgh)
147  );
148 
149  while (pimple.correctNonOrthogonal())
150  {
151  tp_rghEqn = p_rghDDtEqn - fvm::laplacian(rhorAAtUf, p);
152  fvScalarMatrix& p_rghEqn = tp_rghEqn.ref();
153 
154  // Relax the pressure equation to ensure diagonal-dominance
155  p_rghEqn.relax();
156 
157  p_rghEqn.setReference
158  (
161  );
162 
163  fvConstraints().constrain(p_rghEqn);
164 
165  p_rghEqn.solve();
166  }
167  }
168  else
169  {
170  if (pimple.consistent())
171  {
172  const surfaceScalarField gradpf(fvc::snGrad(p_rgh)*mesh.magSf());
173  phiHbyA += (rhorAAtUf - rhorAUf)*gradpf;
174  HbyA += (rAAtU - rAU)*fvc::reconstruct(gradpf - ghGradRhof);
175  }
176 
177  // Update the pressure BCs to ensure flux consistency
178  constrainPressure(p_rgh, rho, U, phiHbyA, rhorAAtUf, MRF);
179 
181 
182  fvScalarMatrix p_rghDDtEqn
183  (
185  + fvc::div(phiHbyA)
186  ==
187  fvModels().sourceProxy(rho, p_rgh)
188  );
189 
190  while (pimple.correctNonOrthogonal())
191  {
192  tp_rghEqn = p_rghDDtEqn - fvm::laplacian(rhorAAtUf, p_rgh);
193  fvScalarMatrix& p_rghEqn = tp_rghEqn.ref();
194 
195  p_rghEqn.setReference
196  (
199  );
200 
201  fvConstraints().constrain(p_rghEqn);
202 
203  p_rghEqn.solve();
204  }
205  }
206 
207  const fvScalarMatrix& p_rghEqn = tp_rghEqn();
208 
209  phi = phiHbyA + p_rghEqn.flux();
210 
211  // Calculate and relax the net pressure-buoyancy force
212  netForce.ref().relax
213  (
214  fvc::reconstruct((ghGradRhof + p_rghEqn.flux()/rhorAAtUf)),
216  );
217 
218  // Correct the momentum source with the pressure gradient flux
219  // calculated from the relaxed pressure
220  U = HbyA + rAAtU*netForce();
223 
224  K = 0.5*magSqr(U);
225 
226  if (!mesh.schemes().steady())
227  {
228  p = p_rgh + rho*gh + pRef;
229 
230  const bool constrained = fvConstraints().constrain(p);
231 
232  // Thermodynamic density update
233  thermo_.correctRho(psi*p - psip0);
234 
235  if (constrained)
236  {
237  rho = thermo.rho();
238  }
239 
240  correctDensity();
241  }
242 
243  continuityErrors();
244 
245  p = p_rgh + rho*gh + pRef;
246 
247  if (mesh.schemes().steady())
248  {
249  if (fvConstraints().constrain(p))
250  {
251  p_rgh = p - rho*gh - pRef;
253  }
254  }
255 
256  // For steady compressible closed-volume cases adjust the pressure level
257  // to obey overall mass continuity
258  if (adjustMass && !thermo.incompressible())
259  {
262  p_rgh = p - rho*gh - pRef;
264  }
265 
266  // Optionally relax pressure for the thermophysics
267  p.relax();
268 
269  if (mesh.schemes().steady() || pimple.simpleRho() || adjustMass)
270  {
271  rho = thermo.rho();
272  }
273 
274  if (mesh.schemes().steady() || pimple.simpleRho())
275  {
276  rho.relax();
277  }
278 
279  // Correct rhoUf if the mesh is moving
281 
282  if (thermo.dpdt())
283  {
284  dpdt = fvc::ddt(p);
285  }
286 }
287 
288 
289 // ************************************************************************* //
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).
scalar relaxationFactor() const
Return the field relaxation factor read from fvSolution.
void correctBoundaryConditions()
Correct boundary field.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
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:1799
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 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: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
volScalarField rho_
The continuity density field.
const surfaceScalarField & phi
Mass-flux field.
volVectorField U_
Velocity field.
const volScalarField & p_rgh
Reference to the buoyant pressure for buoyant cases.
tmp< fvVectorMatrix > tUEqn
Cached momentum matrix.
tmp< volVectorField > netForce
Momentum equation net force source term.
autoPtr< surfaceVectorField > rhoUf
Pointer to the surface momentum field.
fluidThermo & thermo_
Reference to the fluid thermophysical properties.
volScalarField K
Kinetic energy field.
autoPtr< solvers::buoyancy > buoyancy
Pointer to the optional buoyancy force.
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.
volScalarField & p_rgh_
Reference to the buoyant pressure for buoyant cases.
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 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
tmp< VolField< typename outerProduct< vector, Type >::type > > reconstruct(const SurfaceField< Type > &ssf)
dimensioned< Type > domainIntegrate(const VolField< Type > &vf)
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:63
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
UniformDimensionedField< scalar > uniformDimensionedScalarField
SurfaceField< scalar > surfaceScalarField
void constrainHbyA(volVectorField &HbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:34
tmp< surfaceScalarField > constrainPhid(const tmp< surfaceScalarField > &tphiHbyA, const volScalarField &p)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:62
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:45
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
void magSqr(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
fvConstraints constrain(rhoEqn)