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-2022 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 #include "IOMRFZoneList.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 
49 /*---------------------------------------------------------------------------*\
50  Namespace fvc functions Declaration
51 \*---------------------------------------------------------------------------*/
52 
53 namespace fvc
54 {
55  tmp<surfaceScalarField> meshPhi
56  (
57  const volVectorField& U
58  );
59 
60  tmp<scalarField> meshPhi
61  (
62  const volVectorField& U,
63  const label patchi
64  );
65 
66  tmp<surfaceScalarField> meshPhi
67  (
68  const dimensionedScalar& rho,
69  const volVectorField& U
70  );
71 
72  tmp<surfaceScalarField> meshPhi
73  (
74  const volScalarField& rho,
75  const volVectorField& U
76  );
77 
78 
79  //- Make the given flux relative
80  void makeRelative
81  (
82  surfaceScalarField& phi,
83  const volVectorField& U
84  );
85 
86  //- Make the given flux relative
87  void makeRelative
88  (
89  surfaceScalarField& phi,
90  const dimensionedScalar& rho,
91  const volVectorField& U
92  );
93 
94  //- Make the given flux relative
95  void makeRelative
96  (
97  surfaceScalarField& phi,
98  const volScalarField& rho,
99  const volVectorField& U
100  );
101 
102 
103  //- Make the given flux absolute
104  void makeAbsolute
105  (
106  surfaceScalarField& phi,
107  const volVectorField& U
108  );
109 
110  //- Make the given flux absolute
111  void makeAbsolute
112  (
113  surfaceScalarField& phi,
114  const dimensionedScalar& rho,
115  const volVectorField& U
116  );
117 
118  //- Make the given flux absolute
119  void makeAbsolute
120  (
121  surfaceScalarField& phi,
122  const volScalarField& rho,
123  const volVectorField& U
124  );
125 
126 
127  //- Return the given absolute flux in relative form
128  tmp<surfaceScalarField> relative
129  (
130  const tmp<surfaceScalarField>& tphi,
131  const volVectorField& U
132  );
133 
134  //- Return the given absolute flux in relative form
135  tmp<surfaceScalarField> relative
136  (
137  const tmp<surfaceScalarField>& tphi,
138  const volScalarField& rho,
139  const volVectorField& U
140  );
141 
142 
143  //- Return the given relative flux in absolute form
144  tmp<surfaceScalarField> absolute
145  (
146  const tmp<surfaceScalarField>& tphi,
147  const volVectorField& U
148  );
149 
150  //- Return the given relative flux in absolute form
151  tmp<surfaceScalarField> absolute
152  (
153  const tmp<surfaceScalarField>& tphi,
154  const volScalarField& rho,
155  const volVectorField& U
156  );
157 
158  //- Return the given relative flux in absolute form
159  tmp<surfaceScalarField> absolute
160  (
161  const tmp<surfaceScalarField>& tphi,
162  const volScalarField& alpha,
163  const volScalarField& rho,
164  const volVectorField& U
165  );
166 
167  template<class MRFType>
168  void correctUf
169  (
170  autoPtr<surfaceVectorField>& Uf,
171  const volVectorField& U,
172  const surfaceScalarField& phi,
173  const MRFType& MRF
174  );
175 
176  template<class MRFType>
177  void correctRhoUf
178  (
179  autoPtr<surfaceVectorField>& rhoUf,
180  const volScalarField& rho,
181  const volVectorField& U,
182  const surfaceScalarField& phi,
183  const MRFType& MRF
184  );
185 
186  void correctUf
187  (
188  autoPtr<surfaceVectorField>& Uf,
189  const volVectorField& U,
190  const surfaceScalarField& phi
191  );
192 
193  void correctRhoUf
194  (
195  autoPtr<surfaceVectorField>& rhoUf,
196  const volScalarField& rho,
197  const volVectorField& U,
198  const surfaceScalarField& phi
199  );
200 }
201 
202 
203 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
204 
205 } // End namespace Foam
206 
207 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
208 
209 #ifdef NoRepository
210  #include "fvcMeshPhiTemplates.C"
211 #endif
212 
213 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
214 
215 #endif
216 
217 // ************************************************************************* //
IOMRFZoneList MRF(mesh)
label patchi
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
void correctRhoUf(autoPtr< surfaceVectorField > &rhoUf, const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const MRFType &MRF)
tmp< surfaceScalarField > meshPhi(const volVectorField &U)
Definition: fvcMeshPhi.C:34
tmp< surfaceScalarField > relative(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given absolute flux in relative form.
Definition: fvcMeshPhi.C:167
void correctUf(autoPtr< surfaceVectorField > &Uf, const volVectorField &U, const surfaceScalarField &phi, const MRFType &MRF)
void makeAbsolute(surfaceScalarField &phi, const volVectorField &U)
Make the given flux absolute.
Definition: fvcMeshPhi.C:128
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:202
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:89
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
VolField< vector > volVectorField
Definition: volFieldsFwd.H:62
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
SurfaceField< scalar > surfaceScalarField
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61