CorrectPhi.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2015-2021 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 pressureReference;
52  class nonOrthogonalSolutionControl;
53 
54  //- If the mesh is moving correct the velocity BCs on the moving walls to
55  // ensure the corrected fluxes and velocity are consistent
56  void correctUphiBCs
57  (
60  const bool evaluateUBCs
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
65  void correctUphiBCs
66  (
67  const volScalarField& rho,
70  const bool evaluateUBCs
71  );
72 
73  template<class RAUfType, class DivUType>
74  void CorrectPhi
75  (
77  const volVectorField& U,
78  const volScalarField& p,
79  const RAUfType& rAUf,
80  const DivUType& divU,
82  nonOrthogonalSolutionControl& pcorrControl
83  );
84 
85  template<class RAUfType, class DivRhoUType>
86  void CorrectPhi
87  (
89  const volScalarField& p,
90  const volScalarField& rho,
91  const volScalarField& psi,
92  const RAUfType& rAUf,
93  const DivRhoUType& divRhoU,
94  nonOrthogonalSolutionControl& pcorrControl
95  );
96 }
97 
98 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
99 
100 #ifdef NoRepository
101  #include "CorrectPhi.C"
102 #endif
103 
104 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
105 
106 #endif
107 
108 // ************************************************************************* //
pressureReference & pressureReference
void correctUphiBCs(volVectorField &U, surfaceScalarField &phi, const bool evaluateUBCs)
If the mesh is moving correct the velocity BCs on the moving walls to.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:58
void CorrectPhi(surfaceScalarField &phi, const volVectorField &U, const volScalarField &p, const RAUfType &rAUf, const DivUType &divU, const pressureReference &pressureReference, nonOrthogonalSolutionControl &pcorrControl)
dimensionedScalar rAUf("rAUf", dimTime, 1.0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
phi
Definition: correctPhi.H:3
tmp< volScalarField > divU
Definition: alphaSuSp.H:9
U
Definition: pEqn.H:72
const volScalarField & psi
volScalarField & p
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Namespace for OpenFOAM.