fvcMeshPhi.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2014 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 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
33 (
34  const volVectorField& vf
35 )
36 {
38  (
39  vf.mesh(),
40  vf.mesh().ddtScheme("ddt(" + vf.name() + ')')
41  )().meshPhi(vf);
42 }
43 
44 
46 (
47  const dimensionedScalar& rho,
48  const volVectorField& vf
49 )
50 {
52  (
53  vf.mesh(),
54  vf.mesh().ddtScheme("ddt(" + rho.name() + ',' + vf.name() + ')')
55  )().meshPhi(vf);
56 }
57 
58 
60 (
61  const volScalarField& rho,
62  const volVectorField& vf
63 )
64 {
66  (
67  vf.mesh(),
68  vf.mesh().ddtScheme("ddt(" + rho.name() + ',' + vf.name() + ')')
69  )().meshPhi(vf);
70 }
71 
72 
74 (
76  const volVectorField& U
77 )
78 {
79  if (phi.mesh().moving())
80  {
81  phi -= fvc::meshPhi(U);
82  }
83 }
84 
86 (
88  const dimensionedScalar& rho,
89  const volVectorField& U
90 )
91 {
92  if (phi.mesh().moving())
93  {
94  phi -= rho*fvc::meshPhi(rho, U);
95  }
96 }
97 
99 (
101  const volScalarField& rho,
102  const volVectorField& U
103 )
104 {
105  if (phi.mesh().moving())
106  {
108  }
109 }
110 
111 
113 (
115  const volVectorField& U
116 )
117 {
118  if (phi.mesh().moving())
119  {
120  phi += fvc::meshPhi(U);
121  }
122 }
123 
125 (
127  const dimensionedScalar& rho,
128  const volVectorField& U
129 )
130 {
131  if (phi.mesh().moving())
132  {
133  phi += rho*fvc::meshPhi(rho, U);
134  }
135 }
136 
138 (
140  const volScalarField& rho,
141  const volVectorField& U
142 )
143 {
144  if (phi.mesh().moving())
145  {
147  }
148 }
149 
150 
152 (
153  const tmp<surfaceScalarField>& tphi,
154  const volVectorField& U
155 )
156 {
157  if (tphi().mesh().moving())
158  {
159  return tphi - fvc::meshPhi(U);
160  }
161  else
162  {
163  return tmp<surfaceScalarField>(tphi, true);
164  }
165 }
166 
167 
169 (
170  const tmp<surfaceScalarField>& tphi,
171  const volScalarField& rho,
172  const volVectorField& U
173 )
174 {
175  if (tphi().mesh().moving())
176  {
177  return tphi - fvc::interpolate(rho)*fvc::meshPhi(rho, U);
178  }
179  else
180  {
181  return tmp<surfaceScalarField>(tphi, true);
182  }
183 }
184 
185 
187 (
188  const tmp<surfaceScalarField>& tphi,
189  const volVectorField& U
190 )
191 {
192  if (tphi().mesh().moving())
193  {
194  return tphi + fvc::meshPhi(U);
195  }
196  else
197  {
198  return tmp<surfaceScalarField>(tphi, true);
199  }
200 }
201 
202 
204 (
205  const tmp<surfaceScalarField>& tphi,
206  const volScalarField& rho,
207  const volVectorField& U
208 )
209 {
210  if (tphi().mesh().moving())
211  {
212  return tphi + fvc::interpolate(rho)*fvc::meshPhi(rho, U);
213  }
214  else
215  {
216  return tmp<surfaceScalarField>(tphi, true);
217  }
218 }
219 
220 
221 // ************************************************************************* //
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
fvcMeshPhi.H
Calculate the mesh motion flux and convert fluxes from absolute to relative and back.
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::MULES::interpolate
tmp< surfaceScalarField > interpolate(const RhoType &rho)
Definition: IMULESTemplates.C:40
Foam::fvc::meshPhi
tmp< surfaceScalarField > meshPhi(const volVectorField &U)
Definition: fvcMeshPhi.C:33
U
U
Definition: pEqn.H:46
Foam::fvc::makeRelative
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:74
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
fvMesh.H
rho
rho
Definition: pEqn.H:3
ddtScheme.H
Foam::fvc::relative
tmp< surfaceScalarField > relative(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given absolute flux in relative form.
Definition: fvcMeshPhi.C:152
Foam::fvc::makeAbsolute
void makeAbsolute(surfaceScalarField &phi, const volVectorField &U)
Make the given flux absolute.
Definition: fvcMeshPhi.C:113
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::fvc::absolute
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:187