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