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