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 "steadyStateD2dt2Scheme.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 steadyStateD2dt2Scheme<Type>::fvcD2dt2
00045 (
00046 const GeometricField<Type, fvPatchField, volMesh>& vf
00047 )
00048 {
00049 return tmp<GeometricField<Type, fvPatchField, volMesh> >
00050 (
00051 new GeometricField<Type, fvPatchField, volMesh>
00052 (
00053 IOobject
00054 (
00055 "d2dt2("+vf.name()+')',
00056 mesh().time().timeName(),
00057 mesh(),
00058 IOobject::NO_READ,
00059 IOobject::NO_WRITE
00060 ),
00061 mesh(),
00062 dimensioned<Type>
00063 (
00064 "0",
00065 vf.dimensions()/dimTime/dimTime,
00066 pTraits<Type>::zero
00067 )
00068 )
00069 );
00070 }
00071
00072
00073 template<class Type>
00074 tmp<GeometricField<Type, fvPatchField, volMesh> >
00075 steadyStateD2dt2Scheme<Type>::fvcD2dt2
00076 (
00077 const volScalarField& rho,
00078 const GeometricField<Type, fvPatchField, volMesh>& vf
00079 )
00080 {
00081 return tmp<GeometricField<Type, fvPatchField, volMesh> >
00082 (
00083 new GeometricField<Type, fvPatchField, volMesh>
00084 (
00085 IOobject
00086 (
00087 "d2dt2("+rho.name()+','+vf.name()+')',
00088 mesh().time().timeName(),
00089 mesh(),
00090 IOobject::NO_READ,
00091 IOobject::NO_WRITE
00092 ),
00093 mesh(),
00094 dimensioned<Type>
00095 (
00096 "0",
00097 rho.dimensions()*vf.dimensions()/dimTime/dimTime,
00098 pTraits<Type>::zero
00099 )
00100 )
00101 );
00102 }
00103
00104
00105 template<class Type>
00106 tmp<fvMatrix<Type> >
00107 steadyStateD2dt2Scheme<Type>::fvmD2dt2
00108 (
00109 GeometricField<Type, fvPatchField, volMesh>& vf
00110 )
00111 {
00112 tmp<fvMatrix<Type> > tfvm
00113 (
00114 new fvMatrix<Type>
00115 (
00116 vf,
00117 vf.dimensions()*dimVol/dimTime/dimTime
00118 )
00119 );
00120
00121 return tfvm;
00122 }
00123
00124
00125 template<class Type>
00126 tmp<fvMatrix<Type> >
00127 steadyStateD2dt2Scheme<Type>::fvmD2dt2
00128 (
00129 const dimensionedScalar& rho,
00130 GeometricField<Type, fvPatchField, volMesh>& vf
00131 )
00132 {
00133 tmp<fvMatrix<Type> > tfvm
00134 (
00135 new fvMatrix<Type>
00136 (
00137 vf,
00138 rho.dimensions()*vf.dimensions()*dimVol/dimTime/dimTime
00139 )
00140 );
00141
00142 return tfvm;
00143 }
00144
00145
00146 template<class Type>
00147 tmp<fvMatrix<Type> >
00148 steadyStateD2dt2Scheme<Type>::fvmD2dt2
00149 (
00150 const volScalarField& rho,
00151 GeometricField<Type, fvPatchField, volMesh>& vf
00152 )
00153 {
00154 tmp<fvMatrix<Type> > tfvm
00155 (
00156 new fvMatrix<Type>
00157 (
00158 vf,
00159 rho.dimensions()*vf.dimensions()*dimVol/dimTime/dimTime
00160 )
00161 );
00162
00163 return tfvm;
00164 }
00165
00166
00167
00168
00169 }
00170
00171
00172
00173 }
00174
00175