Loading...
Searching...
No Matches
pEqn.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ██╗████████╗██╗ ██╗ █████╗ ██████╗ █████╗ ███████╗██╗ ██╗
3 ██║╚══██╔══╝██║ ██║██╔══██╗██╔════╝██╔══██╗ ██╔════╝██║ ██║
4 ██║ ██║ ███████║███████║██║ ███████║█████╗█████╗ ██║ ██║
5 ██║ ██║ ██╔══██║██╔══██║██║ ██╔══██║╚════╝██╔══╝ ╚██╗ ██╔╝
6 ██║ ██║ ██║ ██║██║ ██║╚██████╗██║ ██║ ██║ ╚████╔╝
7 ╚═╝ ╚═╝ ╚═╝ ╚═╝╚═╝ ╚═╝ ╚═════╝╚═╝ ╚═╝ ╚═╝ ╚═══╝
8
9 * In real Time Highly Advanced Computational Applications for Finite Volumes
10 * Copyright (C) 2017 by the ITHACA-FV authors
11-------------------------------------------------------------------------------
12 License
13 This file is part of ITHACA-FV
14 ITHACA-FV is free software: you can redistribute it and/or modify
15 it under the terms of the GNU Lesser General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18 ITHACA-FV is distributed in the hope that it will be useful,
19 but WITHOUT ANY WARRANTY; without even the implied warranty of
20 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
21 GNU Lesser General Public License for more details.
22 You should have received a copy of the GNU Lesser General Public License
23 along with ITHACA-FV. If not, see <http://www.gnu.org/licenses/>.
24\*---------------------------------------------------------------------------*/
25// Solve pressure equation
26volScalarField rAU(1.0 / UEqn.A());
27
28volVectorField HbyA(constrainHbyA(rAU * UEqn.H(), U, p));
29surfaceScalarField phiHbyA
30(
31 "phiHbyA",
32 fvc::flux(HbyA)
33 + MRF.zeroFilter(fvc::interpolate(rAU)*fvc::ddtCorr(U, phi, _Uf))
34);
35
36MRF.makeRelative(phiHbyA);
37
38if (p.needReference())
39{
40 fvc::makeRelative(phiHbyA, U);
42 fvc::makeAbsolute(phiHbyA, U);
43}
44
45tmp<volScalarField> rAtU(rAU);
46
47if (pimple.consistent())
48{
49 rAtU = 1.0 / max(1.0 / rAU - UEqn.H1(), 0.1 / rAU);
50 phiHbyA +=
51 fvc::interpolate(rAtU() - rAU) * fvc::snGrad(p) * mesh.magSf();
52 HbyA -= (rAU - rAtU()) * fvc::grad(p);
53}
54
55if (pimple.nCorrPISO() <= 1)
56{
57 tUEqn.clear();
58}
59
60// Update the pressure BCs to ensure flux consistency
62
63// Non-orthogonal pressure corrector loop
64while (pimple.correctNonOrthogonal())
65{
66 // Pressure corrector
67 fvScalarMatrix pEqn
68 (
69 fvm::laplacian(rAtU(), p) == fvc::div(phiHbyA)
70 );
71 pEqn.setReference(pRefCell, pRefValue);
72 pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter())));
73
74 if (pimple.finalNonOrthogonalIter())
75 {
76 phi = phiHbyA - pEqn.flux();
77 }
78}
79
80#include "continuityErrs.H"
81
82// Explicitly relax pressure for momentum corrector
83p.relax();
84
85U = HbyA - rAtU() * fvc::grad(p);
86U.correctBoundaryConditions();
87fvOptions.correct(U);
88
89// Correct Uf if the mesh is moving
90fvc::correctUf(_Uf, U, phi);
91
92// Make the fluxes relative to the mesh motion
93fvc::makeRelative(phi, U);
fvScalarMatrix pEqn
Definition CFM.H:31
Foam::fvMesh & mesh
Definition createMesh.H:47
surfaceScalarField & phi
pimpleControl & pimple
tmp< fvVectorMatrix > tUEqn(fvm::div(phi, U)+MRF.DDt(U)+turbulence->divDevReff(U)==fvOptions(U))
fvVectorMatrix & UEqn
Definition UEqn.H:37
volScalarField & p
scalar pRefValue
label pRefCell
constrainPressure(p, U, phiHbyA, rAtU(), MRF)
adjustPhi(phiHbyA, U, p)
tmp< volScalarField > rAtU(rAU)
surfaceScalarField phiHbyA("phiHbyA", fvc::flux(HbyA))
volVectorField HbyA(constrainHbyA(rAU *UEqn.H(), U, p))
U
Definition pEqn.H:79
volScalarField rAU(1.0/UEqn.A())