FreeFOAM The Cross-Platform CFD Toolkit
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) 1991-2010 OpenCFD Ltd.
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 <finiteVolume/fvMesh.H>
28 #include <finiteVolume/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 (
87  surfaceScalarField& phi,
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 (
100  surfaceScalarField& phi,
101  const volScalarField& rho,
102  const volVectorField& U
103 )
104 {
105  if (phi.mesh().moving())
106  {
107  phi -= fvc::interpolate(rho)*fvc::meshPhi(rho, U);
108  }
109 }
110 
111 
113 (
114  surfaceScalarField& phi,
115  const volVectorField& U
116 )
117 {
118  if (phi.mesh().moving())
119  {
120  phi += fvc::meshPhi(U);
121  }
122 }
123 
125 (
126  surfaceScalarField& phi,
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 (
139  surfaceScalarField& phi,
140  const volScalarField& rho,
141  const volVectorField& U
142 )
143 {
144  if (phi.mesh().moving())
145  {
146  phi += fvc::interpolate(rho)*fvc::meshPhi(rho, U);
147  }
148 }
149 
150 
151 // ************************ vim: set sw=4 sts=4 et: ************************ //