initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / finiteVolume / finiteVolume / ddtSchemes / steadyStateDdtScheme / steadyStateDdtScheme.C
blob7bd7b6f8a83f60424124f4171e2f456e0474070a
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
24     
25 \*---------------------------------------------------------------------------*/
27 #include "steadyStateDdtScheme.H"
28 #include "fvcDiv.H"
29 #include "fvMatrices.H"
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 namespace Foam
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 namespace fv
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 template<class Type>
44 tmp<GeometricField<Type, fvPatchField, volMesh> >
45 steadyStateDdtScheme<Type>::fvcDdt
47     const dimensioned<Type>& dt
50     return tmp<GeometricField<Type, fvPatchField, volMesh> >
51     (
52         new GeometricField<Type, fvPatchField, volMesh>
53         (
54             IOobject
55             (
56                 "ddt("+dt.name()+')',
57                 mesh().time().timeName(),
58                 mesh()
59             ),
60             mesh(),
61             dimensioned<Type>
62             (
63                 "0",
64                 dt.dimensions()/dimTime,
65                 pTraits<Type>::zero
66             )
67         )
68     );
72 template<class Type>
73 tmp<GeometricField<Type, fvPatchField, volMesh> >
74 steadyStateDdtScheme<Type>::fvcDdt
76     const GeometricField<Type, fvPatchField, volMesh>& vf
79     return tmp<GeometricField<Type, fvPatchField, volMesh> >
80     (
81         new GeometricField<Type, fvPatchField, volMesh>
82         (
83             IOobject
84             (
85                 "ddt("+vf.name()+')',
86                 mesh().time().timeName(),
87                 mesh()
88             ),
89             mesh(),
90             dimensioned<Type>
91             (
92                 "0",
93                 vf.dimensions()/dimTime,
94                 pTraits<Type>::zero
95             )
96         )
97     );
101 template<class Type>
102 tmp<GeometricField<Type, fvPatchField, volMesh> >
103 steadyStateDdtScheme<Type>::fvcDdt
105     const dimensionedScalar& rho,
106     const GeometricField<Type, fvPatchField, volMesh>& vf
109     return tmp<GeometricField<Type, fvPatchField, volMesh> >
110     (
111         new GeometricField<Type, fvPatchField, volMesh>
112         (
113             IOobject
114             (
115                 "ddt("+rho.name()+','+vf.name()+')',
116                 mesh().time().timeName(),
117                 mesh()
118             ),
119             mesh(),
120             dimensioned<Type>
121             (
122                 "0",
123                 rho.dimensions()*vf.dimensions()/dimTime,
124                 pTraits<Type>::zero
125             )
126         )
127     );
131 template<class Type>
132 tmp<GeometricField<Type, fvPatchField, volMesh> >
133 steadyStateDdtScheme<Type>::fvcDdt
135     const volScalarField& rho,
136     const GeometricField<Type, fvPatchField, volMesh>& vf
139     return tmp<GeometricField<Type, fvPatchField, volMesh> >
140     (
141         new GeometricField<Type, fvPatchField, volMesh>
142         (
143             IOobject
144             (
145                 "ddt("+rho.name()+','+vf.name()+')',
146                 mesh().time().timeName(),
147                 mesh()
148             ),
149             mesh(),
150             dimensioned<Type>
151             (
152                 "0",
153                 rho.dimensions()*vf.dimensions()/dimTime,
154                 pTraits<Type>::zero
155             )
156         )
157     );
161 template<class Type>
162 tmp<fvMatrix<Type> >
163 steadyStateDdtScheme<Type>::fvmDdt
165     GeometricField<Type, fvPatchField, volMesh>& vf
168     tmp<fvMatrix<Type> > tfvm
169     (
170         new fvMatrix<Type>
171         (
172             vf,
173             vf.dimensions()*dimVol/dimTime
174         )
175     );
177     return tfvm;
181 template<class Type>
182 tmp<fvMatrix<Type> >
183 steadyStateDdtScheme<Type>::fvmDdt
185     const dimensionedScalar& rho,
186     GeometricField<Type, fvPatchField, volMesh>& vf
189     tmp<fvMatrix<Type> > tfvm
190     (
191         new fvMatrix<Type>
192         (
193             vf,
194             rho.dimensions()*vf.dimensions()*dimVol/dimTime
195         )
196     );
198     return tfvm;
202 template<class Type>
203 tmp<fvMatrix<Type> >
204 steadyStateDdtScheme<Type>::fvmDdt
206     const volScalarField& rho,
207     GeometricField<Type, fvPatchField, volMesh>& vf
210     tmp<fvMatrix<Type> > tfvm
211     (
212         new fvMatrix<Type>
213         (
214             vf,
215             rho.dimensions()*vf.dimensions()*dimVol/dimTime
216         )
217     );
219     return tfvm;
223 template<class Type>
224 tmp<typename steadyStateDdtScheme<Type>::fluxFieldType>
225 steadyStateDdtScheme<Type>::fvcDdtPhiCorr
227     const volScalarField& rA,
228     const GeometricField<Type, fvPatchField, volMesh>& U,
229     const fluxFieldType& phi
232     return tmp<fluxFieldType>
233     (
234         new fluxFieldType
235         (
236             IOobject
237             (
238                 "ddtPhiCorr("
239               + rA.name() + ',' + U.name() + ',' + phi.name() + ')',
240                 mesh().time().timeName(),
241                 mesh()
242             ),
243             mesh(),
244             dimensioned<typename flux<Type>::type>
245             (
246                 "0",
247                 rA.dimensions()*phi.dimensions()/dimTime,
248                 pTraits<typename flux<Type>::type>::zero
249             )
250         )
251     );
255 template<class Type>
256 tmp<typename steadyStateDdtScheme<Type>::fluxFieldType>
257 steadyStateDdtScheme<Type>::fvcDdtPhiCorr
259     const volScalarField& rA,
260     const volScalarField& rho,
261     const GeometricField<Type, fvPatchField, volMesh>& U,
262     const fluxFieldType& phi
265     return tmp<fluxFieldType>
266     (
267         new fluxFieldType
268         (
269             IOobject
270             (
271                 "ddtPhiCorr("
272               + rA.name() + ',' + rho.name()
273               + ',' + U.name() + ',' + phi.name() + ')',
274                 mesh().time().timeName(),
275                 mesh()
276             ),
277             mesh(),
278             dimensioned<typename flux<Type>::type>
279             (
280                 "0",
281                 rA.dimensions()*rho.dimensions()*phi.dimensions()/dimTime,
282                 pTraits<typename flux<Type>::type>::zero
283             )
284         )
285     );
289 template<class Type>
290 tmp<surfaceScalarField> steadyStateDdtScheme<Type>::meshPhi
292     const GeometricField<Type, fvPatchField, volMesh>& vf
295     return tmp<surfaceScalarField>
296     (
297         new surfaceScalarField
298         (
299             IOobject
300             (
301                 "meshPhi",
302                 mesh().time().timeName(),
303                 mesh()
304             ),
305             mesh(),
306             dimensionedScalar("0", dimVolume/dimTime, 0.0)
307         )
308     );
312 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
314 } // End namespace fv
316 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
318 } // End namespace Foam
320 // ************************************************************************* //