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-------------------------------------------------------------------------------
12License
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
26// Solve pressure equation
27volScalarField rAU(1.0 / UEqn.A());
28volVectorField HbyA(constrainHbyA(rAU* UEqn.H(), U, p));
29surfaceScalarField phiHbyA
30(
31 "phiHbyA",
32 fvc::flux(HbyA)
33 + fvc::interpolate(rAU)*fvc::ddtCorr(U, phi)
34 );
35
36MRF.makeRelative(phiHbyA);
37
39
40tmp<volScalarField> rAtU(rAU);
41
42if (pimple.consistent())
43{
44 rAtU = 1.0 / max(1.0 / rAU - UEqn.H1(), 0.1 / rAU);
45 phiHbyA +=
46 fvc::interpolate(rAtU() - rAU) * fvc::snGrad(p) * mesh.magSf();
47 HbyA -= (rAU - rAtU()) * fvc::grad(p);
48}
49
50if (pimple.nCorrPISO() <= 1)
51{
52 tUEqn.clear();
53}
54
55// Update the pressure BCs to ensure flux consistency
57
58// Non-orthogonal pressure corrector loop
59while (pimple.correctNonOrthogonal())
60{
61 // Pressure corrector
62 fvScalarMatrix pEqn
63 (
64 fvm::laplacian(rAtU(), p) == fvc::div(phiHbyA)
65 );
66 pEqn.setReference(pRefCell, pRefValue);
67 pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter())));
68
69 if (pimple.finalNonOrthogonalIter())
70 {
71 phi = phiHbyA - pEqn.flux();
72 }
73}
74
75#include "continuityErrs.H"
76
77// Explicitly relax pressure for momentum corrector
78p.relax();
79
80U = HbyA - rAtU() * fvc::grad(p);
81U.correctBoundaryConditions();
82fvOptions.correct(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())