CorrectPhi.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2015-2016 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26Class
27 Foam::CorrectPhi
28
29Description
30 Flux correction functions to ensure continuity.
31
32 Required during start-up, restart, mesh-motion etc. when non-conservative
33 fluxes may adversely affect the prediction-part of the solution algorithm
34 (the part before the first pressure solution which would ensure continuity).
35 This is particularly important for VoF and other multi-phase solver in
36 which non-conservative fluxes cause unboundedness of the phase-fraction.
37
38SourceFiles
39 CorrectPhi.C
40
41\*---------------------------------------------------------------------------*/
42
43#ifndef CorrectPhi_H
44#define CorrectPhi_H
45
46#include "volFieldsFwd.H"
47#include "surfaceFieldsFwd.H"
48
49// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50
51namespace Foam
52{
53 class pimpleControl;
54
55 //- If the mesh is moving correct the velocity BCs on the moving walls to
56 // ensure the corrected fluxes and velocity are consistent
58 (
61 );
62
63 //- If the mesh is moving correct the velocity BCs on the moving walls to
64 // ensure the corrected fluxes and velocity are consistent
66 (
67 const volScalarField& rho,
70 );
71
72 template<class RAUfType, class DivUType>
73 void CorrectPhi
74 (
77 const volScalarField& p,
78 const RAUfType& rAUf,
79 const DivUType& divU,
80 pimpleControl& pimple
81 );
82
83 template<class RAUfType, class DivRhoUType>
84 void CorrectPhi
85 (
88 const volScalarField& p,
89 const volScalarField& rho,
90 const volScalarField& psi,
91 const RAUfType& rAUf,
92 const DivRhoUType& divRhoU,
93 pimpleControl& pimple
94 );
95}
96
97// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
98
99#ifdef NoRepository
100 #include "CorrectPhi.C"
101#endif
102
103// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
104
105#endif
106
107// ************************************************************************* //
surfaceScalarField & phi
pimpleControl & pimple
Flux correction functions to ensure continuity.
U
Definition: pEqn.H:72
volScalarField & p
const volScalarField & psi
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU))
zeroField divU
Definition: alphaSuSp.H:3
Namespace for OpenFOAM.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:83
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
void correctUphiBCs(volVectorField &U, surfaceScalarField &phi)
If the mesh is moving correct the velocity BCs on the moving walls to.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField