Loading...
Searching...
No Matches
pEqn.H
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/>.
24Description
25 Example of steady NS Reduction Problem solved by the use of the SIMPLE algorithm
26SourceFiles
27 fsiBasic.C
28\*---------------------------------------------------------------------------*/
29volScalarField rAU(1.0 / UEqn.A());
30volVectorField HbyA(constrainHbyA(rAU*UEqn.H(), U, p));
31surfaceScalarField phiHbyA("phiHbyA", fvc::flux(HbyA));
32
33if (pimple.ddtCorr())
34{
35 phiHbyA += MRF.zeroFilter(fvc::interpolate(rAU) * fvc::ddtCorr(U, phi, Uf));
36}
37else
38{
39 phiHbyA += MRF.zeroFilter(fvc::interpolate(rAU));
40}
41
42MRF.makeRelative(phiHbyA);
43
44if (p.needReference())
45{
46 fvc::makeRelative(phiHbyA, U);
47 adjustPhi(phiHbyA, U, p);
48 fvc::makeAbsolute(phiHbyA, U);
49}
50
51tmp<volScalarField> rAtU(rAU);
52
53if (pimple.consistent())
54{
55 rAtU = 1.0 / max(1.0 / rAU - UEqn.H1(), 0.1 / rAU);
56 phiHbyA +=
57 fvc::interpolate(rAtU() - rAU) * fvc::snGrad(p) * mesh.magSf();
58 HbyA -= (rAU - rAtU()) * fvc::grad(p);
59}
60
61if (pimple.nCorrPISO() <= 1)
62{
63 tUEqn.clear();
64}
65
66// Update the pressure BCs to ensure flux consistency
67constrainPressure(p, U, phiHbyA, rAtU(), MRF);
68
69// Non-orthogonal pressure corrector loop
70while (pimple.correctNonOrthogonal())
71{
72 fvScalarMatrix pEqn
73 (
74 fvm::laplacian(rAtU(), p) == fvc::div(phiHbyA)
75 );
76 pEqn.setReference(pRefCell, pRefValue);
77 pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter())));
78
79 if (pimple.finalNonOrthogonalIter())
80 {
81 phi = phiHbyA - pEqn.flux();
82 }
83}
84
85#include "continuityErrs.H"
86
87// Explicitly relax pressure for momentum corrector
88p.relax();
89
90U = HbyA - rAtU * fvc::grad(p);
91U.correctBoundaryConditions();
92// std::cout << "============" << "the divergence of the velocity is "<< "\n";
93// Info<< max(mag(fvc::div(U))) << "\n";
94fvOptions.correct(U);
95
96// Correct Uf if the mesh is moving
97fvc::correctUf(Uf, U, phi);
98// {
99// Uf = fvc::interpolate(U);
100// surfaceVectorField n(mesh.Sf()/mesh.magSf());
101// Uf += n*(phi/mesh.magSf() - (n & Uf));
102// }
103
104// Make the fluxes relative to the mesh motion
105fvc::makeRelative(phi, U);
T max(Eigen::SparseMatrix< T > &mat, label &ind_row, label &ind_col)
Find the maximum of a sparse Matrix (Useful for DEIM).