00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026 #include "fvcDdt.H"
00027 #include <finiteVolume/fvMesh.H>
00028 #include <finiteVolume/ddtScheme.H>
00029
00030
00031
00032 namespace Foam
00033 {
00034
00035
00036
00037 namespace fvc
00038 {
00039
00040
00041
00042 template<class Type>
00043 tmp<GeometricField<Type, fvPatchField, volMesh> >
00044 ddt
00045 (
00046 const dimensioned<Type> dt,
00047 const fvMesh& mesh
00048 )
00049 {
00050 return fv::ddtScheme<Type>::New
00051 (
00052 mesh,
00053 mesh.ddtScheme("ddt(" + dt.name() + ')')
00054 )().fvcDdt(dt);
00055 }
00056
00057
00058 template<class Type>
00059 tmp<GeometricField<Type, fvPatchField, volMesh> >
00060 ddt
00061 (
00062 const GeometricField<Type, fvPatchField, volMesh>& vf
00063 )
00064 {
00065 return fv::ddtScheme<Type>::New
00066 (
00067 vf.mesh(),
00068 vf.mesh().ddtScheme("ddt(" + vf.name() + ')')
00069 )().fvcDdt(vf);
00070 }
00071
00072
00073 template<class Type>
00074 tmp<GeometricField<Type, fvPatchField, volMesh> >
00075 ddt
00076 (
00077 const dimensionedScalar& rho,
00078 const GeometricField<Type, fvPatchField, volMesh>& vf
00079 )
00080 {
00081 return fv::ddtScheme<Type>::New
00082 (
00083 vf.mesh(),
00084 vf.mesh().ddtScheme("ddt(" + rho.name() + ',' + vf.name() + ')')
00085 )().fvcDdt(rho, vf);
00086 }
00087
00088
00089 template<class Type>
00090 tmp<GeometricField<Type, fvPatchField, volMesh> >
00091 ddt
00092 (
00093 const volScalarField& rho,
00094 const GeometricField<Type, fvPatchField, volMesh>& vf
00095 )
00096 {
00097 return fv::ddtScheme<Type>::New
00098 (
00099 vf.mesh(),
00100 vf.mesh().ddtScheme("ddt(" + rho.name() + ',' + vf.name() + ')')
00101 )().fvcDdt(rho, vf);
00102 }
00103
00104
00105 template<class Type>
00106 tmp<GeometricField<typename flux<Type>::type, fvsPatchField, surfaceMesh> >
00107 ddtPhiCorr
00108 (
00109 const volScalarField& rA,
00110 const GeometricField<Type, fvPatchField, volMesh>& U,
00111 const GeometricField
00112 <
00113 typename flux<Type>::type,
00114 fvsPatchField,
00115 surfaceMesh
00116 >& phi
00117 )
00118 {
00119 return fv::ddtScheme<Type>::New
00120 (
00121 U.mesh(),
00122 U.mesh().ddtScheme("ddt(" + U.name() + ')')
00123 )().fvcDdtPhiCorr(rA, U, phi);
00124 }
00125
00126
00127 template<class Type>
00128 tmp<GeometricField<typename flux<Type>::type, fvsPatchField, surfaceMesh> >
00129 ddtPhiCorr
00130 (
00131 const volScalarField& rA,
00132 const volScalarField& rho,
00133 const GeometricField<Type, fvPatchField, volMesh>& U,
00134 const GeometricField
00135 <
00136 typename flux<Type>::type,
00137 fvsPatchField,
00138 surfaceMesh
00139 >& phi
00140 )
00141 {
00142 return fv::ddtScheme<Type>::New
00143 (
00144 U.mesh(),
00145 U.mesh().ddtScheme("ddt(" + rho.name() + ',' + U.name() + ')')
00146 )().fvcDdtPhiCorr(rA, rho, U, phi);
00147 }
00148
00149
00150
00151
00152 }
00153
00154
00155
00156 }
00157
00158