initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / finiteVolume / finiteVolume / fvm / fvmDdt.C
blobc2e3e976af080d122cdf5c6add9d44e2931bdc14
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
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.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "volFields.H"
28 #include "surfaceFields.H"
29 #include "fvMatrix.H"
30 #include "ddtScheme.H"
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 namespace Foam
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 namespace fvm
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 template<class Type>
45 tmp<fvMatrix<Type> >
46 ddt
48     GeometricField<Type, fvPatchField, volMesh>& vf
51     return fv::ddtScheme<Type>::New
52     (
53         vf.mesh(),
54         vf.mesh().ddtScheme("ddt(" + vf.name() + ')')
55     )().fvmDdt(vf);
59 template<class Type>
60 tmp<fvMatrix<Type> >
61 ddt
63     const oneField&,
64     GeometricField<Type, fvPatchField, volMesh>& vf
67     return ddt(vf);
71 template<class Type>
72 tmp<fvMatrix<Type> >
73 ddt
75     const dimensionedScalar& rho,
76     GeometricField<Type, fvPatchField, volMesh>& vf
79     return fv::ddtScheme<Type>::New
80     (
81         vf.mesh(),
82         vf.mesh().ddtScheme("ddt(" + rho.name() + ',' + vf.name() + ')')
83     )().fvmDdt(rho, vf);
87 template<class Type>
88 tmp<fvMatrix<Type> >
89 ddt
91     const volScalarField& rho,
92     GeometricField<Type, fvPatchField, volMesh>& vf
95     return fv::ddtScheme<Type>::New
96     (
97         vf.mesh(),
98         vf.mesh().ddtScheme("ddt(" + rho.name() + ',' + vf.name() + ')')
99     )().fvmDdt(rho, vf);
103 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
105 } // End namespace fvm
107 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
109 } // End namespace Foam
111 // ************************************************************************* //