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 "incompressibleFluid.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 "fvmLaplacian.H"
36 
37 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
38 
40 {
44 
45  fvVectorMatrix& UEqn = tUEqn.ref();
46 
47  volScalarField rAU(1.0/UEqn.A());
50  (
51  "phiHbyA",
54  );
55 
57 
58  if (p.needReference())
59  {
61  adjustPhi(phiHbyA, U, p);
63  }
64 
66 
67  if (pimple.consistent())
68  {
69  rAtU = 1.0/max(1.0/rAU - UEqn.H1(), 0.1/rAU);
70  phiHbyA +=
72  HbyA -= (rAU - rAtU())*fvc::grad(p);
73  }
74 
75  if (pimple.nCorrPiso() <= 1)
76  {
77  tUEqn.clear();
78  }
79 
80  // Update the pressure BCs to ensure flux consistency
81  constrainPressure(p, U, phiHbyA, rAtU(), MRF);
82 
83  // Non-orthogonal pressure corrector loop
85  {
86  fvScalarMatrix pEqn
87  (
88  fvm::laplacian(rAtU(), p) == fvc::div(phiHbyA)
89  );
90 
91  pEqn.setReference
92  (
95  );
96 
97  pEqn.solve();
98 
100  {
101  phi = phiHbyA - pEqn.flux();
102  }
103  }
104 
106 
107  // Explicitly relax pressure for momentum corrector
108  p.relax();
109 
110  U = HbyA - rAtU*fvc::grad(p);
113 
114  // Correct Uf if the mesh is moving
115  fvc::correctUf(Uf, U, phi, MRF);
116 
117  // Make the fluxes relative to the mesh motion
119 }
120 
121 
122 // ************************************************************************* //
For cases which do no have a pressure boundary adjust the balance of fluxes to obey continuity....
Generic GeometricField class.
void relax(const scalar alpha)
Relax field (for steady-state solution).
bool needReference() const
Does the field need a reference level for solution.
void correctBoundaryConditions()
Correct boundary field.
void makeRelative(volVectorField &U) const
Make the given absolute velocity relative within the MRF region.
Definition: MRFZoneList.C:157
bool consistent() const
Flag to indicate to relax pressure using the "consistent".
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
void setReference(const label celli, const Type &value, const bool forceReference=false)
Set reference level for solution.
Definition: fvMatrix.C:563
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
bool finalNonOrthogonalIter() const
Flag to indicate the last non-orthogonal iteration.
bool correctNonOrthogonal()
Non-orthogonal corrector loop.
label nCorrPiso() const
Maximum number of piso correctors.
Definition: pisoControlI.H:28
Provides controls for the pressure reference in closed-volume simulations.
scalar refValue() const
Return the pressure reference level.
label refCell() const
Return the cell in which the reference pressure is set.
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_
Pressure field.
const surfaceScalarField & phi
Reference to the volumetric-flux field.
volVectorField U_
Velocity field.
tmp< fvVectorMatrix > tUEqn
Cached momentum matrix.
autoPtr< surfaceVectorField > Uf
Pointer to the surface momentum field.
const volVectorField & U
Reference to the velocity field.
IOMRFZoneList MRF
MRF zone list.
virtual void correctPressure()
Construct the pressure equation.
surfaceScalarField phi_
Volumetric-flux field.
void continuityErrors()
Calculate and print the continuity errors.
const volScalarField & p
Reference to the pressure field.
A class for managing temporary objects.
Definition: tmp.H:55
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.
Calculate the snGrad of the given volField.
Calculate the matrix for the laplacian of the field.
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.
void correctUf(autoPtr< surfaceVectorField > &Uf, const volVectorField &U, const surfaceScalarField &phi, const MRFType &MRF)
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 makeAbsolute(surfaceScalarField &phi, const volVectorField &U)
Make the given flux absolute.
Definition: fvcMeshPhi.C:128
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
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
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:34