FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

steadyStateDdtScheme.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023     
00024 \*---------------------------------------------------------------------------*/
00025 
00026 #include "steadyStateDdtScheme.H"
00027 #include <finiteVolume/fvcDiv.H>
00028 #include <finiteVolume/fvMatrices.H>
00029 
00030 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00031 
00032 namespace Foam
00033 {
00034 
00035 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00036 
00037 namespace fv
00038 {
00039 
00040 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00041 
00042 template<class Type>
00043 tmp<GeometricField<Type, fvPatchField, volMesh> >
00044 steadyStateDdtScheme<Type>::fvcDdt
00045 (
00046     const dimensioned<Type>& dt
00047 )
00048 {
00049     return tmp<GeometricField<Type, fvPatchField, volMesh> >
00050     (
00051         new GeometricField<Type, fvPatchField, volMesh>
00052         (
00053             IOobject
00054             (
00055                 "ddt("+dt.name()+')',
00056                 mesh().time().timeName(),
00057                 mesh()
00058             ),
00059             mesh(),
00060             dimensioned<Type>
00061             (
00062                 "0",
00063                 dt.dimensions()/dimTime,
00064                 pTraits<Type>::zero
00065             )
00066         )
00067     );
00068 }
00069 
00070 
00071 template<class Type>
00072 tmp<GeometricField<Type, fvPatchField, volMesh> >
00073 steadyStateDdtScheme<Type>::fvcDdt
00074 (
00075     const GeometricField<Type, fvPatchField, volMesh>& vf
00076 )
00077 {
00078     return tmp<GeometricField<Type, fvPatchField, volMesh> >
00079     (
00080         new GeometricField<Type, fvPatchField, volMesh>
00081         (
00082             IOobject
00083             (
00084                 "ddt("+vf.name()+')',
00085                 mesh().time().timeName(),
00086                 mesh()
00087             ),
00088             mesh(),
00089             dimensioned<Type>
00090             (
00091                 "0",
00092                 vf.dimensions()/dimTime,
00093                 pTraits<Type>::zero
00094             )
00095         )
00096     );
00097 }
00098 
00099 
00100 template<class Type>
00101 tmp<GeometricField<Type, fvPatchField, volMesh> >
00102 steadyStateDdtScheme<Type>::fvcDdt
00103 (
00104     const dimensionedScalar& rho,
00105     const GeometricField<Type, fvPatchField, volMesh>& vf
00106 )
00107 {
00108     return tmp<GeometricField<Type, fvPatchField, volMesh> >
00109     (
00110         new GeometricField<Type, fvPatchField, volMesh>
00111         (
00112             IOobject
00113             (
00114                 "ddt("+rho.name()+','+vf.name()+')',
00115                 mesh().time().timeName(),
00116                 mesh()
00117             ),
00118             mesh(),
00119             dimensioned<Type>
00120             (
00121                 "0",
00122                 rho.dimensions()*vf.dimensions()/dimTime,
00123                 pTraits<Type>::zero
00124             )
00125         )
00126     );
00127 }
00128 
00129 
00130 template<class Type>
00131 tmp<GeometricField<Type, fvPatchField, volMesh> >
00132 steadyStateDdtScheme<Type>::fvcDdt
00133 (
00134     const volScalarField& rho,
00135     const GeometricField<Type, fvPatchField, volMesh>& vf
00136 )
00137 {
00138     return tmp<GeometricField<Type, fvPatchField, volMesh> >
00139     (
00140         new GeometricField<Type, fvPatchField, volMesh>
00141         (
00142             IOobject
00143             (
00144                 "ddt("+rho.name()+','+vf.name()+')',
00145                 mesh().time().timeName(),
00146                 mesh()
00147             ),
00148             mesh(),
00149             dimensioned<Type>
00150             (
00151                 "0",
00152                 rho.dimensions()*vf.dimensions()/dimTime,
00153                 pTraits<Type>::zero
00154             )
00155         )
00156     );
00157 }
00158 
00159 
00160 template<class Type>
00161 tmp<fvMatrix<Type> >
00162 steadyStateDdtScheme<Type>::fvmDdt
00163 (
00164     GeometricField<Type, fvPatchField, volMesh>& vf
00165 )
00166 {
00167     tmp<fvMatrix<Type> > tfvm
00168     (
00169         new fvMatrix<Type>
00170         (
00171             vf,
00172             vf.dimensions()*dimVol/dimTime
00173         )
00174     );
00175 
00176     return tfvm;
00177 }
00178 
00179 
00180 template<class Type>
00181 tmp<fvMatrix<Type> >
00182 steadyStateDdtScheme<Type>::fvmDdt
00183 (
00184     const dimensionedScalar& rho,
00185     GeometricField<Type, fvPatchField, volMesh>& vf
00186 )
00187 {
00188     tmp<fvMatrix<Type> > tfvm
00189     (
00190         new fvMatrix<Type>
00191         (
00192             vf,
00193             rho.dimensions()*vf.dimensions()*dimVol/dimTime
00194         )
00195     );
00196 
00197     return tfvm;
00198 }
00199 
00200 
00201 template<class Type>
00202 tmp<fvMatrix<Type> >
00203 steadyStateDdtScheme<Type>::fvmDdt
00204 (
00205     const volScalarField& rho,
00206     GeometricField<Type, fvPatchField, volMesh>& vf
00207 )
00208 {
00209     tmp<fvMatrix<Type> > tfvm
00210     (
00211         new fvMatrix<Type>
00212         (
00213             vf,
00214             rho.dimensions()*vf.dimensions()*dimVol/dimTime
00215         )
00216     );
00217 
00218     return tfvm;
00219 }
00220 
00221 
00222 template<class Type>
00223 tmp<typename steadyStateDdtScheme<Type>::fluxFieldType>
00224 steadyStateDdtScheme<Type>::fvcDdtPhiCorr
00225 (
00226     const volScalarField& rA,
00227     const GeometricField<Type, fvPatchField, volMesh>& U,
00228     const fluxFieldType& phi
00229 )
00230 {
00231     return tmp<fluxFieldType>
00232     (
00233         new fluxFieldType
00234         (
00235             IOobject
00236             (
00237                 "ddtPhiCorr("
00238               + rA.name() + ',' + U.name() + ',' + phi.name() + ')',
00239                 mesh().time().timeName(),
00240                 mesh()
00241             ),
00242             mesh(),
00243             dimensioned<typename flux<Type>::type>
00244             (
00245                 "0",
00246                 rA.dimensions()*phi.dimensions()/dimTime,
00247                 pTraits<typename flux<Type>::type>::zero
00248             )
00249         )
00250     );
00251 }
00252 
00253 
00254 template<class Type>
00255 tmp<typename steadyStateDdtScheme<Type>::fluxFieldType>
00256 steadyStateDdtScheme<Type>::fvcDdtPhiCorr
00257 (
00258     const volScalarField& rA,
00259     const volScalarField& rho,
00260     const GeometricField<Type, fvPatchField, volMesh>& U,
00261     const fluxFieldType& phi
00262 )
00263 {
00264     return tmp<fluxFieldType>
00265     (
00266         new fluxFieldType
00267         (
00268             IOobject
00269             (
00270                 "ddtPhiCorr("
00271               + rA.name() + ',' + rho.name()
00272               + ',' + U.name() + ',' + phi.name() + ')',
00273                 mesh().time().timeName(),
00274                 mesh()
00275             ),
00276             mesh(),
00277             dimensioned<typename flux<Type>::type>
00278             (
00279                 "0",
00280                 rA.dimensions()*rho.dimensions()*phi.dimensions()/dimTime,
00281                 pTraits<typename flux<Type>::type>::zero
00282             )
00283         )
00284     );
00285 }
00286 
00287 
00288 template<class Type>
00289 tmp<surfaceScalarField> steadyStateDdtScheme<Type>::meshPhi
00290 (
00291     const GeometricField<Type, fvPatchField, volMesh>& vf
00292 )
00293 {
00294     return tmp<surfaceScalarField>
00295     (
00296         new surfaceScalarField
00297         (
00298             IOobject
00299             (
00300                 "meshPhi",
00301                 mesh().time().timeName(),
00302                 mesh()
00303             ),
00304             mesh(),
00305             dimensionedScalar("0", dimVolume/dimTime, 0.0)
00306         )
00307     );
00308 }
00309 
00310 
00311 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00312 
00313 } // End namespace fv
00314 
00315 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00316 
00317 } // End namespace Foam
00318 
00319 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines