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
35//MRF.makeRelative(phiHbyA);
36
38
39// Update the pressure BCs to ensure flux consistency
41
42// Non-orthogonal pressure corrector loop
43while (piso.correctNonOrthogonal())
44{
45 // Pressure corrector
46 fvScalarMatrix pEqn
47 (
48 fvm::laplacian(rAU, p) == fvc::div(phiHbyA)
49 );
50 pEqn.setReference(pRefCell, pRefValue);
51 pEqn.solve(mesh.solver(p.select(piso.finalInnerIter())));
52
53 if (piso.finalNonOrthogonalIter())
54 {
55 phi = phiHbyA - pEqn.flux();
56 }
57}
58
59#include "continuityErrs.H"
60
61U = HbyA - rAU * fvc::grad(p);
62U.correctBoundaryConditions();
63//fvOptions.correct(U);
fvScalarMatrix pEqn
Definition CFM.H:31
Foam::fvMesh & mesh
Definition createMesh.H:47
surfaceScalarField & phi
fvVectorMatrix & UEqn
Definition UEqn.H:37
volScalarField & p
scalar pRefValue
label pRefCell
constrainPressure(p, U, phiHbyA, rAtU(), MRF)
adjustPhi(phiHbyA, U, p)
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())