fvcMeshPhi.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) 2011-2017 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
26InNamespace
27 Foam::fvc
28
29Description
30 Calculate the mesh motion flux and convert fluxes from absolute to relative
31 and back.
32
33SourceFiles
34 fvcMeshPhi.C
35
36\*---------------------------------------------------------------------------*/
37
38#ifndef fvcMeshPhi_H
39#define fvcMeshPhi_H
40
41#include "volFieldsFwd.H"
42#include "surfaceFieldsFwd.H"
43#include "dimensionedTypes.H"
44
45// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46
47namespace Foam
48{
49
50/*---------------------------------------------------------------------------*\
51 Namespace fvc functions Declaration
52\*---------------------------------------------------------------------------*/
53
54namespace fvc
55{
56 tmp<surfaceScalarField> meshPhi
57 (
58 const volVectorField& U
59 );
60
61 tmp<surfaceScalarField> meshPhi
62 (
64 const volVectorField& U
65 );
66
67 tmp<surfaceScalarField> meshPhi
68 (
69 const volScalarField& rho,
70 const volVectorField& U
71 );
72
73
74 //- Make the given flux relative
75 void makeRelative
76 (
78 const volVectorField& U
79 );
80
81 //- Make the given flux relative
82 void makeRelative
83 (
86 const volVectorField& U
87 );
88
89 //- Make the given flux relative
90 void makeRelative
91 (
93 const volScalarField& rho,
94 const volVectorField& U
95 );
96
97
98 //- Make the given flux absolute
99 void makeAbsolute
100 (
102 const volVectorField& U
103 );
104
105 //- Make the given flux absolute
106 void makeAbsolute
107 (
109 const dimensionedScalar& rho,
110 const volVectorField& U
111 );
112
113 //- Make the given flux absolute
114 void makeAbsolute
115 (
117 const volScalarField& rho,
118 const volVectorField& U
119 );
120
121
122 //- Return the given absolute flux in relative form
123 tmp<surfaceScalarField> relative
124 (
125 const tmp<surfaceScalarField>& tphi,
126 const volVectorField& U
127 );
128
129 //- Return the given absolute flux in relative form
130 tmp<surfaceScalarField> relative
131 (
132 const tmp<surfaceScalarField>& tphi,
133 const volScalarField& rho,
134 const volVectorField& U
135 );
136
137
138 //- Return the given relative flux in absolute form
139 tmp<surfaceScalarField> absolute
140 (
141 const tmp<surfaceScalarField>& tphi,
142 const volVectorField& U
143 );
144
145 //- Return the given relative flux in absolute form
146 tmp<surfaceScalarField> absolute
147 (
148 const tmp<surfaceScalarField>& tphi,
149 const volScalarField& rho,
150 const volVectorField& U
151 );
152
153 void correctUf
154 (
155 autoPtr<surfaceVectorField>& Uf,
156 const volVectorField& U,
158 );
159
160 void correctRhoUf
161 (
162 autoPtr<surfaceVectorField>& rhoUf,
163 const volScalarField& rho,
164 const volVectorField& U,
166 );
167}
168
169
170// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
171
172} // End namespace Foam
173
174// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
175
176#endif
177
178// ************************************************************************* //
surfaceScalarField & phi
U
Definition: pEqn.H:72
autoPtr< surfaceVectorField > rhoUf
autoPtr< surfaceVectorField > Uf
void correctRhoUf(autoPtr< surfaceVectorField > &rhoUf, const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi)
Definition: fvcMeshPhi.C:243
tmp< surfaceScalarField > meshPhi(const volVectorField &U)
Definition: fvcMeshPhi.C:36
tmp< surfaceScalarField > relative(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given absolute flux in relative form.
Definition: fvcMeshPhi.C:155
void correctUf(autoPtr< surfaceVectorField > &Uf, const volVectorField &U, const surfaceScalarField &phi)
Definition: fvcMeshPhi.C:225
void makeAbsolute(surfaceScalarField &phi, const volVectorField &U)
Make the given flux absolute.
Definition: fvcMeshPhi.C:116
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:190
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:77
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:83
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField