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
26{
27 volScalarField rAU("rAU", 1.0 / UEqn.A());
28 surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU));
29 volVectorField HbyA(constrainHbyA(rAU * UEqn.H(), U, p_rgh));
30
31 surfaceScalarField phig(-rAUf* ghf * fvc::snGrad(rhok)*mesh.magSf());
32
33 surfaceScalarField phiHbyA
34 (
35 "phiHbyA",
36 fvc::flux(HbyA)
37 + rAUf * fvc::ddtCorr(U, phi)
38 + phig
39 );
40
41 MRF.makeRelative(phiHbyA);
42
43 // Update the pressure BCs to ensure flux consistency
45
46 while (pimple.correctNonOrthogonal())
47 {
48 fvScalarMatrix p_rghEqn
49 (
50 fvm::laplacian(rAUf, p_rgh) == fvc::div(phiHbyA)
51 );
52 p_rghEqn.setReference(pRefCell, getRefCellValue(p_rgh, pRefCell));
53 p_rghEqn.solve(mesh.solver(p_rgh.select(pimple.finalInnerIter())));
54
55 if (pimple.finalNonOrthogonalIter())
56 {
57 // Calculate the conservative fluxes
58 phi = phiHbyA - p_rghEqn.flux();
59 // Explicitly relax pressure for momentum corrector
60 p_rgh.relax();
61 // Correct the momentum source with the pressure gradient flux calculated from the relaxed pressure
62 U = HbyA + rAU * fvc::reconstruct((phig - p_rghEqn.flux()) / rAUf);
63 U.correctBoundaryConditions();
64 fvOptions.correct(U);
65 }
66 }
67
68#include "continuityErrs.H"
69
70 p = p_rgh + rhok * gh;
71
72
73 if (p_rgh.needReference())
74 {
75 p += dimensionedScalar
76 (
77 "p",
78 p.dimensions(),
79 pRefValue - getRefCellValue(p, pRefCell)
80 );
81 p_rgh = p - rhok * gh;
82 }
83}
84
85
Foam::fvMesh & mesh
Definition createMesh.H:47
surfaceScalarField & phi
volScalarField & gh
volScalarField & p_rgh
pimpleControl & pimple
volScalarField & rhok
surfaceScalarField & ghf
surfaceScalarField phig(-rAUf *ghf *fvc::snGrad(rhok) *mesh.magSf())
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU))
fvVectorMatrix & UEqn
Definition UEqn.H:37
volScalarField & p
scalar pRefValue
label pRefCell
constrainPressure(p, U, phiHbyA, rAtU(), MRF)
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())