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 "EulerD2dt2Scheme.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 EulerD2dt2Scheme<Type>::fvcD2dt2
00045 (
00046 const GeometricField<Type, fvPatchField, volMesh>& vf
00047 )
00048 {
00049 dimensionedScalar rDeltaT2 =
00050 4.0/sqr(mesh().time().deltaT() + mesh().time().deltaT0());
00051
00052 IOobject d2dt2IOobject
00053 (
00054 "d2dt2("+vf.name()+')',
00055 mesh().time().timeName(),
00056 mesh(),
00057 IOobject::NO_READ,
00058 IOobject::NO_WRITE
00059 );
00060
00061 scalar deltaT = mesh().time().deltaT().value();
00062 scalar deltaT0 = mesh().time().deltaT0().value();
00063
00064 scalar coefft = (deltaT + deltaT0)/(2*deltaT);
00065 scalar coefft00 = (deltaT + deltaT0)/(2*deltaT0);
00066 scalar coefft0 = coefft + coefft00;
00067
00068 if (mesh().moving())
00069 {
00070 scalar halfRdeltaT2 = rDeltaT2.value()/2.0;
00071
00072 scalarField VV0 = mesh().V() + mesh().V0();
00073 scalarField V0V00 = mesh().V0() + mesh().V00();
00074
00075 return tmp<GeometricField<Type, fvPatchField, volMesh> >
00076 (
00077 new GeometricField<Type, fvPatchField, volMesh>
00078 (
00079 d2dt2IOobject,
00080 mesh(),
00081 rDeltaT2.dimensions()*vf.dimensions(),
00082 halfRdeltaT2*
00083 (
00084 coefft*VV0*vf.internalField()
00085
00086 - (coefft*VV0 + coefft00*V0V00)
00087 *vf.oldTime().internalField()
00088
00089 + (coefft00*V0V00)*vf.oldTime().oldTime().internalField()
00090 )/mesh().V(),
00091 rDeltaT2.value()*
00092 (
00093 coefft*vf.boundaryField()
00094 - coefft0*vf.oldTime().boundaryField()
00095 + coefft00*vf.oldTime().oldTime().boundaryField()
00096 )
00097 )
00098 );
00099 }
00100 else
00101 {
00102 return tmp<GeometricField<Type, fvPatchField, volMesh> >
00103 (
00104 new GeometricField<Type, fvPatchField, volMesh>
00105 (
00106 d2dt2IOobject,
00107 rDeltaT2*
00108 (
00109 coefft*vf
00110 - coefft0*vf.oldTime()
00111 + coefft00*vf.oldTime().oldTime()
00112 )
00113 )
00114 );
00115 }
00116 }
00117
00118
00119 template<class Type>
00120 tmp<GeometricField<Type, fvPatchField, volMesh> >
00121 EulerD2dt2Scheme<Type>::fvcD2dt2
00122 (
00123 const volScalarField& rho,
00124 const GeometricField<Type, fvPatchField, volMesh>& vf
00125 )
00126 {
00127 dimensionedScalar rDeltaT2 =
00128 4.0/sqr(mesh().time().deltaT() + mesh().time().deltaT0());
00129
00130 IOobject d2dt2IOobject
00131 (
00132 "d2dt2("+rho.name()+','+vf.name()+')',
00133 mesh().time().timeName(),
00134 mesh(),
00135 IOobject::NO_READ,
00136 IOobject::NO_WRITE
00137 );
00138
00139 scalar deltaT = mesh().time().deltaT().value();
00140 scalar deltaT0 = mesh().time().deltaT0().value();
00141
00142 scalar coefft = (deltaT + deltaT0)/(2*deltaT);
00143 scalar coefft00 = (deltaT + deltaT0)/(2*deltaT0);
00144
00145 if (mesh().moving())
00146 {
00147 scalar halfRdeltaT2 = 0.5*rDeltaT2.value();
00148 scalar quarterRdeltaT2 = 0.25*rDeltaT2.value();
00149
00150 scalarField VV0rhoRho0 =
00151 (mesh().V() + mesh().V0())
00152 *(rho.internalField() + rho.oldTime().internalField());
00153
00154 scalarField V0V00rho0Rho00 =
00155 (mesh().V0() + mesh().V00())
00156 *(
00157 rho.oldTime().internalField()
00158 + rho.oldTime().oldTime().internalField()
00159 );
00160
00161 return tmp<GeometricField<Type, fvPatchField, volMesh> >
00162 (
00163 new GeometricField<Type, fvPatchField, volMesh>
00164 (
00165 d2dt2IOobject,
00166 mesh(),
00167 rDeltaT2.dimensions()*rho.dimensions()*vf.dimensions(),
00168 quarterRdeltaT2*
00169 (
00170 coefft*VV0rhoRho0*vf.internalField()
00171
00172 - (coefft*VV0rhoRho0 + coefft00*V0V00rho0Rho00)
00173 *vf.oldTime().internalField()
00174
00175 + (coefft00*V0V00rho0Rho00)
00176 *vf.oldTime().oldTime().internalField()
00177 )/mesh().V(),
00178 halfRdeltaT2*
00179 (
00180 coefft
00181 *(rho.boundaryField() + rho.oldTime().boundaryField())
00182 *vf.boundaryField()
00183
00184 - (
00185 coefft
00186 *(
00187 rho.boundaryField()
00188 + rho.oldTime().boundaryField()
00189 )
00190 + coefft00
00191 *(
00192 rho.oldTime().boundaryField()
00193 + rho.oldTime().oldTime().boundaryField()
00194 )
00195 )*vf.oldTime().boundaryField()
00196
00197 + coefft00
00198 *(
00199 rho.oldTime().boundaryField()
00200 + rho.oldTime().oldTime().boundaryField()
00201 )*vf.oldTime().oldTime().boundaryField()
00202 )
00203 )
00204 );
00205 }
00206 else
00207 {
00208 dimensionedScalar halfRdeltaT2 = 0.5*rDeltaT2;
00209
00210 volScalarField rhoRho0 = rho + rho.oldTime();
00211 volScalarField rho0Rho00 = rho.oldTime() +rho.oldTime().oldTime();
00212
00213 return tmp<GeometricField<Type, fvPatchField, volMesh> >
00214 (
00215 new GeometricField<Type, fvPatchField, volMesh>
00216 (
00217 d2dt2IOobject,
00218 halfRdeltaT2*
00219 (
00220 coefft*rhoRho0*vf
00221 - (coefft*rhoRho0 + coefft00*rho0Rho00)*vf.oldTime()
00222 + coefft00*rho0Rho00*vf.oldTime().oldTime()
00223 )
00224 )
00225 );
00226 }
00227 }
00228
00229
00230 template<class Type>
00231 tmp<fvMatrix<Type> >
00232 EulerD2dt2Scheme<Type>::fvmD2dt2
00233 (
00234 GeometricField<Type, fvPatchField, volMesh>& vf
00235 )
00236 {
00237 tmp<fvMatrix<Type> > tfvm
00238 (
00239 new fvMatrix<Type>
00240 (
00241 vf,
00242 vf.dimensions()*dimVol/dimTime/dimTime
00243 )
00244 );
00245
00246 fvMatrix<Type>& fvm = tfvm();
00247
00248 scalar deltaT = mesh().time().deltaT().value();
00249 scalar deltaT0 = mesh().time().deltaT0().value();
00250
00251 scalar coefft = (deltaT + deltaT0)/(2*deltaT);
00252 scalar coefft00 = (deltaT + deltaT0)/(2*deltaT0);
00253 scalar coefft0 = coefft + coefft00;
00254
00255 scalar rDeltaT2 = 4.0/sqr(deltaT + deltaT0);
00256
00257 if (mesh().moving())
00258 {
00259 scalar halfRdeltaT2 = rDeltaT2/2.0;
00260
00261 scalarField VV0 = mesh().V() + mesh().V0();
00262 scalarField V0V00 = mesh().V0() + mesh().V00();
00263
00264 fvm.diag() = (coefft*halfRdeltaT2)*VV0;
00265
00266 fvm.source() = halfRdeltaT2*
00267 (
00268 (coefft*VV0 + coefft00*V0V00)
00269 *vf.oldTime().internalField()
00270
00271 - (coefft00*V0V00)*vf.oldTime().oldTime().internalField()
00272 );
00273 }
00274 else
00275 {
00276 fvm.diag() = (coefft*rDeltaT2)*mesh().V();
00277
00278 fvm.source() = rDeltaT2*mesh().V()*
00279 (
00280 coefft0*vf.oldTime().internalField()
00281 - coefft00*vf.oldTime().oldTime().internalField()
00282 );
00283 }
00284
00285 return tfvm;
00286 }
00287
00288
00289 template<class Type>
00290 tmp<fvMatrix<Type> >
00291 EulerD2dt2Scheme<Type>::fvmD2dt2
00292 (
00293 const dimensionedScalar& rho,
00294 GeometricField<Type, fvPatchField, volMesh>& vf
00295 )
00296 {
00297 tmp<fvMatrix<Type> > tfvm
00298 (
00299 new fvMatrix<Type>
00300 (
00301 vf,
00302 rho.dimensions()*vf.dimensions()*dimVol
00303 /dimTime/dimTime
00304 )
00305 );
00306
00307 fvMatrix<Type>& fvm = tfvm();
00308
00309 scalar deltaT = mesh().time().deltaT().value();
00310 scalar deltaT0 = mesh().time().deltaT0().value();
00311
00312 scalar coefft = (deltaT + deltaT0)/(2*deltaT);
00313 scalar coefft00 = (deltaT + deltaT0)/(2*deltaT0);
00314
00315 scalar rDeltaT2 = 4.0/sqr(deltaT + deltaT0);
00316
00317 if (mesh().moving())
00318 {
00319 scalar halfRdeltaT2 = 0.5*rDeltaT2;
00320
00321 scalarField VV0 = mesh().V() + mesh().V0();
00322
00323 scalarField V0V00 = mesh().V0() + mesh().V00();
00324
00325 fvm.diag() = rho.value()*(coefft*halfRdeltaT2)*VV0;
00326
00327 fvm.source() = halfRdeltaT2*rho.value()*
00328 (
00329 (coefft*VV0 + coefft00*V0V00)
00330 *vf.oldTime().internalField()
00331
00332 - (coefft00*V0V00)*vf.oldTime().oldTime().internalField()
00333 );
00334 }
00335 else
00336 {
00337 fvm.diag() = (coefft*rDeltaT2)*mesh().V()*rho.value();
00338
00339 fvm.source() = rDeltaT2*mesh().V()*rho.value()*
00340 (
00341 (coefft + coefft00)*vf.oldTime().internalField()
00342 - coefft00*vf.oldTime().oldTime().internalField()
00343 );
00344 }
00345
00346 return tfvm;
00347 }
00348
00349
00350 template<class Type>
00351 tmp<fvMatrix<Type> >
00352 EulerD2dt2Scheme<Type>::fvmD2dt2
00353 (
00354 const volScalarField& rho,
00355 GeometricField<Type, fvPatchField, volMesh>& vf
00356 )
00357 {
00358 tmp<fvMatrix<Type> > tfvm
00359 (
00360 new fvMatrix<Type>
00361 (
00362 vf,
00363 rho.dimensions()*vf.dimensions()*dimVol
00364 /dimTime/dimTime
00365 )
00366 );
00367
00368 fvMatrix<Type>& fvm = tfvm();
00369
00370 scalar deltaT = mesh().time().deltaT().value();
00371 scalar deltaT0 = mesh().time().deltaT0().value();
00372
00373 scalar coefft = (deltaT + deltaT0)/(2*deltaT);
00374 scalar coefft00 = (deltaT + deltaT0)/(2*deltaT0);
00375
00376 scalar rDeltaT2 = 4.0/sqr(deltaT + deltaT0);
00377
00378 if (mesh().moving())
00379 {
00380 scalar quarterRdeltaT2 = 0.25*rDeltaT2;
00381
00382 scalarField VV0rhoRho0 =
00383 (mesh().V() + mesh().V0())
00384 *(rho.internalField() + rho.oldTime().internalField());
00385
00386 scalarField V0V00rho0Rho00 =
00387 (mesh().V0() + mesh().V00())
00388 *(
00389 rho.oldTime().internalField()
00390 + rho.oldTime().oldTime().internalField()
00391 );
00392
00393 fvm.diag() = (coefft*quarterRdeltaT2)*VV0rhoRho0;
00394
00395 fvm.source() = quarterRdeltaT2*
00396 (
00397 (coefft*VV0rhoRho0 + coefft00*V0V00rho0Rho00)
00398 *vf.oldTime().internalField()
00399
00400 - (coefft00*V0V00rho0Rho00)
00401 *vf.oldTime().oldTime().internalField()
00402 );
00403 }
00404 else
00405 {
00406 scalar halfRdeltaT2 = 0.5*rDeltaT2;
00407
00408 scalarField rhoRho0 =
00409 (rho.internalField() + rho.oldTime().internalField());
00410
00411 scalarField rho0Rho00 =
00412 (
00413 rho.oldTime().internalField()
00414 + rho.oldTime().oldTime().internalField()
00415 );
00416
00417 fvm.diag() = (coefft*halfRdeltaT2)*mesh().V()*rhoRho0;
00418
00419 fvm.source() = halfRdeltaT2*mesh().V()*
00420 (
00421 (coefft*rhoRho0 + coefft00*rho0Rho00)
00422 *vf.oldTime().internalField()
00423
00424 - (coefft00*rho0Rho00)
00425 *vf.oldTime().oldTime().internalField()
00426 );
00427 }
00428
00429 return tfvm;
00430 }
00431
00432
00433
00434
00435 }
00436
00437
00438
00439 }
00440
00441