FreeFOAM The Cross-Platform CFD Toolkit
fvmDdt.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 <finiteVolume/volFields.H>
28 #include <finiteVolume/fvMatrix.H>
29 #include <finiteVolume/ddtScheme.H>
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace fvm
39 {
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 template<class Type>
44 tmp<fvMatrix<Type> >
45 ddt
46 (
48 )
49 {
51  (
52  vf.mesh(),
53  vf.mesh().ddtScheme("ddt(" + vf.name() + ')')
54  )().fvmDdt(vf);
55 }
56 
57 
58 template<class Type>
60 ddt
61 (
62  const geometricOneField&,
64 )
65 {
66  return ddt(vf);
67 }
68 
69 
70 template<class Type>
72 ddt
73 (
74  const dimensionedScalar& rho,
76 )
77 {
79  (
80  vf.mesh(),
81  vf.mesh().ddtScheme("ddt(" + rho.name() + ',' + vf.name() + ')')
82  )().fvmDdt(rho, vf);
83 }
84 
85 
86 template<class Type>
88 ddt
89 (
90  const volScalarField& rho,
92 )
93 {
95  (
96  vf.mesh(),
97  vf.mesh().ddtScheme("ddt(" + rho.name() + ',' + vf.name() + ')')
98  )().fvmDdt(rho, vf);
99 }
100 
101 
102 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
103 
104 } // End namespace fvm
105 
106 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
107 
108 } // End namespace Foam
109 
110 // ************************ vim: set sw=4 sts=4 et: ************************ //