Loading...
Searching...
No Matches
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 | 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
26\*---------------------------------------------------------------------------*/
27
28#include "fvcMeshPhi.H"
29#include "fvMesh.H"
30#include "ddtScheme.H"
31#include "surfaceInterpolate.H"
32
33// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
36(
37 const volVectorField& vf
38)
39{
41 (
42 vf.mesh(),
43 vf.mesh().ddtScheme("ddt(" + vf.name() + ')')
44 ).ref().meshPhi(vf);
45}
46
47
49(
51 const volVectorField& vf
52)
53{
55 (
56 vf.mesh(),
57 vf.mesh().ddtScheme("ddt(" + rho.name() + ',' + vf.name() + ')')
58 ).ref().meshPhi(vf);
59}
60
61
63(
64 const volScalarField& rho,
65 const volVectorField& vf
66)
67{
69 (
70 vf.mesh(),
71 vf.mesh().ddtScheme("ddt(" + rho.name() + ',' + vf.name() + ')')
72 ).ref().meshPhi(vf);
73}
74
75
77(
79 const volVectorField& U
80)
82 if (phi.mesh().moving())
83 {
84 phi -= fvc::meshPhi(U);
85 }
86}
87
89(
92 const volVectorField& U
93)
95 if (phi.mesh().moving())
96 {
98 }
99}
100
102(
104 const volScalarField& rho,
105 const volVectorField& U
106)
107{
108 if (phi.mesh().moving())
109 {
111 }
112}
113
114
116(
117 surfaceScalarField& phi,
118 const volVectorField& U
119)
121 if (phi.mesh().moving())
122 {
123 phi += fvc::meshPhi(U);
124 }
125}
126
128(
130 const dimensionedScalar& rho,
131 const volVectorField& U
132)
134 if (phi.mesh().moving())
135 {
136 phi += rho*fvc::meshPhi(rho, U);
137 }
138}
139
141(
143 const volScalarField& rho,
144 const volVectorField& U
145)
146{
147 if (phi.mesh().moving())
148 {
150 }
151}
152
153
155(
156 const tmp<surfaceScalarField>& tphi,
157 const volVectorField& U
158)
159{
160 if (tphi().mesh().moving())
161 {
162 return tphi - fvc::meshPhi(U);
163 }
164 else
165 {
166 return tmp<surfaceScalarField>(tphi, true);
167 }
168}
169
170
172(
173 const tmp<surfaceScalarField>& tphi,
174 const volScalarField& rho,
175 const volVectorField& U
176)
177{
178 if (tphi().mesh().moving())
179 {
180 return tphi - fvc::interpolate(rho)*fvc::meshPhi(rho, U);
181 }
182 else
183 {
184 return tmp<surfaceScalarField>(tphi, true);
185 }
186}
187
188
190(
191 const tmp<surfaceScalarField>& tphi,
192 const volVectorField& U
193)
194{
195 if (tphi().mesh().moving())
196 {
197 return tphi + fvc::meshPhi(U);
198 }
199 else
200 {
201 return tmp<surfaceScalarField>(tphi, true);
202 }
203}
204
205
207(
208 const tmp<surfaceScalarField>& tphi,
209 const volScalarField& rho,
210 const volVectorField& U
211)
212{
213 if (tphi().mesh().moving())
214 {
215 return tphi + fvc::interpolate(rho)*fvc::meshPhi(rho, U);
216 }
217 else
218 {
219 return tmp<surfaceScalarField>(tphi, true);
220 }
221}
222
223
225(
226 autoPtr<surfaceVectorField>& Uf,
227 const volVectorField& U,
228 const surfaceScalarField& phi
229)
230{
231 const fvMesh& mesh = U.mesh();
232
233 if (mesh.dynamic())
234 {
236 surfaceVectorField n(mesh.Sf()/mesh.magSf());
237 Uf() += n*(phi/mesh.magSf() - (n & Uf()));
238 }
239}
240
241
243(
244 autoPtr<surfaceVectorField>& rhoUf,
245 const volScalarField& rho,
246 const volVectorField& U,
247 const surfaceScalarField& phi
248)
249{
250 const fvMesh& mesh = U.mesh();
251
252 if (mesh.dynamic())
253 {
254 rhoUf() = fvc::interpolate(rho*U);
255 surfaceVectorField n(mesh.Sf()/mesh.magSf());
256 rhoUf() += n*(fvc::absolute(phi, rho, U)/mesh.magSf() - (n & rhoUf()));
257 }
258}
259
260
261// ************************************************************************* //
label n
const Mesh & mesh() const noexcept
Return const reference to mesh.
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
static tmp< ddtScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return a pointer to a new ddtScheme created on freestore.
Definition ddtScheme.C:43
A class for managing temporary objects.
Definition tmp.H:75
U
Definition pEqn.H:72
dynamicFvMesh & mesh
autoPtr< surfaceVectorField > rhoUf
autoPtr< surfaceVectorField > Uf
Calculate the mesh motion flux and convert fluxes from absolute to relative and back.
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &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)
Definition fvcMeshPhi.C:236
tmp< surfaceScalarField > meshPhi(const volVectorField &U)
Definition fvcMeshPhi.C:29
tmp< surfaceScalarField > relative(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given absolute flux in relative form.
Definition fvcMeshPhi.C:148
void correctUf(autoPtr< surfaceVectorField > &Uf, const volVectorField &U, const surfaceScalarField &phi)
Definition fvcMeshPhi.C:218
void makeAbsolute(surfaceScalarField &phi, const volVectorField &U)
Make the given flux absolute.
Definition fvcMeshPhi.C:109
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition fvcMeshPhi.C:183
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition fvcMeshPhi.C:70
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.