initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / finiteVolume / finiteVolume / fvc / fvcDdt.C
blobd85c9d953cb766ad1a247d02b9b351feca38d707
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 "fvcDdt.H"
28 #include "fvMesh.H"
29 #include "ddtScheme.H"
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 namespace Foam
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 namespace fvc
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 template<class Type>
44 tmp<GeometricField<Type, fvPatchField, volMesh> >
45 ddt
47     const dimensioned<Type> dt,
48     const fvMesh& mesh
51     return fv::ddtScheme<Type>::New
52     (
53         mesh,
54         mesh.ddtScheme("ddt(" + dt.name() + ')')
55     )().fvcDdt(dt);
59 template<class Type>
60 tmp<GeometricField<Type, fvPatchField, volMesh> >
61 ddt
63     const GeometricField<Type, fvPatchField, volMesh>& vf
66     return fv::ddtScheme<Type>::New
67     (
68         vf.mesh(),
69         vf.mesh().ddtScheme("ddt(" + vf.name() + ')')
70     )().fvcDdt(vf);
74 template<class Type>
75 tmp<GeometricField<Type, fvPatchField, volMesh> >
76 ddt
78     const dimensionedScalar& rho,
79     const GeometricField<Type, fvPatchField, volMesh>& vf
82     return fv::ddtScheme<Type>::New
83     (
84         vf.mesh(),
85         vf.mesh().ddtScheme("ddt(" + rho.name() + ',' + vf.name() + ')')
86     )().fvcDdt(rho, vf);
90 template<class Type>
91 tmp<GeometricField<Type, fvPatchField, volMesh> >
92 ddt
94     const volScalarField& rho,
95     const GeometricField<Type, fvPatchField, volMesh>& vf
98     return fv::ddtScheme<Type>::New
99     (
100         vf.mesh(),
101         vf.mesh().ddtScheme("ddt(" + rho.name() + ',' + vf.name() + ')')
102     )().fvcDdt(rho, vf);
106 template<class Type>
107 tmp<GeometricField<typename flux<Type>::type, fvsPatchField, surfaceMesh> >
108 ddtPhiCorr
110     const volScalarField& rA,
111     const GeometricField<Type, fvPatchField, volMesh>& U,
112     const GeometricField
113     <
114         typename flux<Type>::type,
115         fvsPatchField,
116         surfaceMesh
117     >& phi
120     return fv::ddtScheme<Type>::New
121     (
122         U.mesh(),
123         U.mesh().ddtScheme("ddt(" + U.name() + ')')
124     )().fvcDdtPhiCorr(rA, U, phi);
128 template<class Type>
129 tmp<GeometricField<typename flux<Type>::type, fvsPatchField, surfaceMesh> >
130 ddtPhiCorr
132     const volScalarField& rA,
133     const volScalarField& rho,
134     const GeometricField<Type, fvPatchField, volMesh>& U,
135     const GeometricField
136     <
137         typename flux<Type>::type,
138         fvsPatchField,
139         surfaceMesh
140     >& phi
143     return fv::ddtScheme<Type>::New
144     (
145         U.mesh(),
146         U.mesh().ddtScheme("ddt(" + rho.name() + ',' + U.name() + ')')
147     )().fvcDdtPhiCorr(rA, rho, U, phi);
151 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
153 } // End namespace fvc
155 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
157 } // End namespace Foam
159 // ************************************************************************* //