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