FreeFOAM The Cross-Platform CFD Toolkit
fvcDdt.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 "fvcDdt.H"
27 #include <finiteVolume/fvMesh.H>
28 #include <finiteVolume/ddtScheme.H>
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 namespace fvc
38 {
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 template<class Type>
43 tmp<GeometricField<Type, fvPatchField, volMesh> >
44 ddt
45 (
46  const dimensioned<Type> dt,
47  const fvMesh& mesh
48 )
49 {
51  (
52  mesh,
53  mesh.ddtScheme("ddt(" + dt.name() + ')')
54  )().fvcDdt(dt);
55 }
56 
57 
58 template<class Type>
60 ddt
61 (
63 )
64 {
66  (
67  vf.mesh(),
68  vf.mesh().ddtScheme("ddt(" + vf.name() + ')')
69  )().fvcDdt(vf);
70 }
71 
72 
73 template<class Type>
75 ddt
76 (
77  const dimensionedScalar& rho,
79 )
80 {
82  (
83  vf.mesh(),
84  vf.mesh().ddtScheme("ddt(" + rho.name() + ',' + vf.name() + ')')
85  )().fvcDdt(rho, vf);
86 }
87 
88 
89 template<class Type>
91 ddt
92 (
93  const volScalarField& rho,
95 )
96 {
98  (
99  vf.mesh(),
100  vf.mesh().ddtScheme("ddt(" + rho.name() + ',' + vf.name() + ')')
101  )().fvcDdt(rho, vf);
102 }
103 
104 
105 template<class Type>
108 (
109  const volScalarField& rA,
111  const GeometricField
112  <
113  typename flux<Type>::type,
116  >& phi
117 )
118 {
120  (
121  U.mesh(),
122  U.mesh().ddtScheme("ddt(" + U.name() + ')')
123  )().fvcDdtPhiCorr(rA, U, phi);
124 }
125 
126 
127 template<class Type>
130 (
131  const volScalarField& rA,
132  const volScalarField& rho,
134  const GeometricField
135  <
136  typename flux<Type>::type,
139  >& phi
140 )
141 {
143  (
144  U.mesh(),
145  U.mesh().ddtScheme("ddt(" + rho.name() + ',' + U.name() + ')')
146  )().fvcDdtPhiCorr(rA, rho, U, phi);
147 }
148 
149 
150 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
151 
152 } // End namespace fvc
153 
154 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
155 
156 } // End namespace Foam
157 
158 // ************************ vim: set sw=4 sts=4 et: ************************ //