fvcMeshPhi.C
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 \*---------------------------------------------------------------------------*/
25 
26 #include "fvcMeshPhi.H"
27 #include "fvMesh.H"
28 #include "ddtScheme.H"
29 #include "surfaceInterpolate.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
34 (
35  const volVectorField& vf
36 )
37 {
39  (
40  vf.mesh(),
41  vf.mesh().schemes().ddt("ddt(" + vf.name() + ')')
42  ).ref().meshPhi(vf);
43 }
44 
45 
47 (
48  const volVectorField& vf,
49  const label patchi
50 )
51 {
53  (
54  vf.mesh(),
55  vf.mesh().schemes().ddt("ddt(" + vf.name() + ')')
56  ).ref().meshPhi(vf, patchi);
57 }
58 
59 
61 (
62  const dimensionedScalar& rho,
63  const volVectorField& vf
64 )
65 {
67  (
68  vf.mesh(),
69  vf.mesh().schemes().ddt("ddt(" + rho.name() + ',' + vf.name() + ')')
70  ).ref().meshPhi(vf);
71 }
72 
73 
75 (
76  const volScalarField& rho,
77  const volVectorField& vf
78 )
79 {
81  (
82  vf.mesh(),
83  vf.mesh().schemes().ddt("ddt(" + rho.name() + ',' + vf.name() + ')')
84  ).ref().meshPhi(vf);
85 }
86 
87 
89 (
90  surfaceScalarField& phi,
91  const volVectorField& U
92 )
93 {
94  if (phi.mesh().moving())
95  {
96  phi -= fvc::meshPhi(U);
97  }
98 }
99 
101 (
102  surfaceScalarField& phi,
103  const dimensionedScalar& rho,
104  const volVectorField& U
105 )
106 {
107  if (phi.mesh().moving())
108  {
109  phi -= rho*fvc::meshPhi(rho, U);
110  }
111 }
112 
114 (
115  surfaceScalarField& phi,
116  const volScalarField& rho,
117  const volVectorField& U
118 )
119 {
120  if (phi.mesh().moving())
121  {
123  }
124 }
125 
126 
128 (
129  surfaceScalarField& phi,
130  const volVectorField& U
131 )
132 {
133  if (phi.mesh().moving())
134  {
135  phi += fvc::meshPhi(U);
136  }
137 }
138 
140 (
141  surfaceScalarField& phi,
142  const dimensionedScalar& rho,
143  const volVectorField& U
144 )
145 {
146  if (phi.mesh().moving())
147  {
148  phi += rho*fvc::meshPhi(rho, U);
149  }
150 }
151 
153 (
154  surfaceScalarField& phi,
155  const volScalarField& rho,
156  const volVectorField& U
157 )
158 {
159  if (phi.mesh().moving())
160  {
162  }
163 }
164 
165 
167 (
168  const tmp<surfaceScalarField>& tphi,
169  const volVectorField& U
170 )
171 {
172  if (tphi().mesh().moving())
173  {
174  return tphi - fvc::meshPhi(U);
175  }
176  else
177  {
178  return tmp<surfaceScalarField>(tphi, true);
179  }
180 }
181 
182 
184 (
185  const tmp<surfaceScalarField>& tphi,
186  const volScalarField& rho,
187  const volVectorField& U
188 )
189 {
190  if (tphi().mesh().moving())
191  {
192  return tphi - fvc::interpolate(rho)*fvc::meshPhi(rho, U);
193  }
194  else
195  {
196  return tmp<surfaceScalarField>(tphi, true);
197  }
198 }
199 
200 
202 (
203  const tmp<surfaceScalarField>& tphi,
204  const volVectorField& U
205 )
206 {
207  if (tphi().mesh().moving())
208  {
209  const word phiName(tphi().name());
210 
212  (
213  phiName,
214  tphi + fvc::meshPhi(U)
215  );
216  }
217  else
218  {
219  return tmp<surfaceScalarField>(tphi, true);
220  }
221 }
222 
223 
225 (
226  const tmp<surfaceScalarField>& tphi,
227  const volScalarField& rho,
228  const volVectorField& U
229 )
230 {
231  if (tphi().mesh().moving())
232  {
233  const word phiName(tphi().name());
234 
236  (
237  phiName,
239  );
240  }
241  else
242  {
243  return tmp<surfaceScalarField>(tphi, true);
244  }
245 }
246 
247 
249 (
250  const tmp<surfaceScalarField>& tphi,
251  const volScalarField& alpha,
252  const volScalarField& rho,
253  const volVectorField& U
254 )
255 {
256  if (tphi().mesh().moving())
257  {
258  const word phiName(tphi().name());
259 
261  (
262  phiName,
263  tphi
265  );
266  }
267  else
268  {
269  return tmp<surfaceScalarField>(tphi, true);
270  }
271 }
272 
273 
275 (
277  const volVectorField& U,
278  const surfaceScalarField& phi
279 )
280 {
281  correctUf(Uf, U, phi, NullMRF());
282 }
283 
284 
286 (
288  const volScalarField& rho,
289  const volVectorField& U,
290  const surfaceScalarField& phi
291 )
292 {
293  correctRhoUf(rhoUf, rho, U, phi, NullMRF());
294 }
295 
296 
297 // ************************************************************************* //
const Mesh & mesh() const
Return mesh.
Generic GeometricField class.
const word & name() const
Return name.
Definition: IOobject.H:310
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Calculate the mesh motion flux and convert fluxes from absolute to relative and back.
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))
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
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
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
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47