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 "CoEulerDdtScheme.H"
00027 #include <finiteVolume/surfaceInterpolate.H>
00028 #include <finiteVolume/fvcDiv.H>
00029 #include <finiteVolume/fvMatrices.H>
00030
00031
00032
00033 namespace Foam
00034 {
00035
00036
00037
00038 namespace fv
00039 {
00040
00041
00042
00043 template<class Type>
00044 tmp<volScalarField> CoEulerDdtScheme<Type>::CorDeltaT() const
00045 {
00046 surfaceScalarField cofrDeltaT = CofrDeltaT();
00047
00048 tmp<volScalarField> tcorDeltaT
00049 (
00050 new volScalarField
00051 (
00052 IOobject
00053 (
00054 "CorDeltaT",
00055 cofrDeltaT.instance(),
00056 mesh()
00057 ),
00058 mesh(),
00059 dimensionedScalar("CorDeltaT", cofrDeltaT.dimensions(), 0.0),
00060 zeroGradientFvPatchScalarField::typeName
00061 )
00062 );
00063
00064 volScalarField& corDeltaT = tcorDeltaT();
00065
00066 const unallocLabelList& owner = mesh().owner();
00067 const unallocLabelList& neighbour = mesh().neighbour();
00068
00069 forAll(owner, faceI)
00070 {
00071 corDeltaT[owner[faceI]] =
00072 max(corDeltaT[owner[faceI]], cofrDeltaT[faceI]);
00073
00074 corDeltaT[neighbour[faceI]] =
00075 max(corDeltaT[neighbour[faceI]], cofrDeltaT[faceI]);
00076 }
00077
00078 forAll(corDeltaT.boundaryField(), patchi)
00079 {
00080 const fvsPatchScalarField& pcofrDeltaT =
00081 cofrDeltaT.boundaryField()[patchi];
00082
00083 const fvPatch& p = pcofrDeltaT.patch();
00084 const unallocLabelList& faceCells = p.patch().faceCells();
00085
00086 forAll(pcofrDeltaT, patchFacei)
00087 {
00088 corDeltaT[faceCells[patchFacei]] = max
00089 (
00090 corDeltaT[faceCells[patchFacei]],
00091 pcofrDeltaT[patchFacei]
00092 );
00093 }
00094 }
00095
00096 corDeltaT.correctBoundaryConditions();
00097
00098
00099
00100 return tcorDeltaT;
00101 }
00102
00103
00104 template<class Type>
00105 tmp<surfaceScalarField> CoEulerDdtScheme<Type>::CofrDeltaT() const
00106 {
00107 const dimensionedScalar& deltaT = mesh().time().deltaT();
00108
00109 const surfaceScalarField& phi =
00110 static_cast<const objectRegistry&>(mesh())
00111 .lookupObject<surfaceScalarField>(phiName_);
00112
00113 if (phi.dimensions() == dimensionSet(0, 3, -1, 0, 0))
00114 {
00115 surfaceScalarField Co
00116 (
00117 mesh().surfaceInterpolation::deltaCoeffs()
00118 *(mag(phi)/mesh().magSf())
00119 *deltaT
00120 );
00121
00122 return max(Co/maxCo_, scalar(1))/deltaT;
00123 }
00124 else if (phi.dimensions() == dimensionSet(1, 0, -1, 0, 0))
00125 {
00126 const volScalarField& rho =
00127 static_cast<const objectRegistry&>(mesh())
00128 .lookupObject<volScalarField>(rhoName_).oldTime();
00129
00130 surfaceScalarField Co
00131 (
00132 mesh().surfaceInterpolation::deltaCoeffs()
00133 *(mag(phi)/(fvc::interpolate(rho)*mesh().magSf()))
00134 *deltaT
00135 );
00136
00137 return max(Co/maxCo_, scalar(1))/deltaT;
00138 }
00139 else
00140 {
00141 FatalErrorIn("CoEulerDdtScheme<Type>::CofrDeltaT() const")
00142 << "Incorrect dimensions of phi: " << phi.dimensions()
00143 << abort(FatalError);
00144
00145 return tmp<surfaceScalarField>(NULL);
00146 }
00147 }
00148
00149
00150 template<class Type>
00151 tmp<GeometricField<Type, fvPatchField, volMesh> >
00152 CoEulerDdtScheme<Type>::fvcDdt
00153 (
00154 const dimensioned<Type>& dt
00155 )
00156 {
00157 volScalarField rDeltaT = CorDeltaT();
00158
00159 IOobject ddtIOobject
00160 (
00161 "ddt("+dt.name()+')',
00162 mesh().time().timeName(),
00163 mesh()
00164 );
00165
00166 if (mesh().moving())
00167 {
00168 tmp<GeometricField<Type, fvPatchField, volMesh> > tdtdt
00169 (
00170 new GeometricField<Type, fvPatchField, volMesh>
00171 (
00172 ddtIOobject,
00173 mesh(),
00174 dimensioned<Type>
00175 (
00176 "0",
00177 dt.dimensions()/dimTime,
00178 pTraits<Type>::zero
00179 )
00180 )
00181 );
00182
00183 tdtdt().internalField() =
00184 rDeltaT.internalField()*dt.value()*(1.0 - mesh().V0()/mesh().V());
00185
00186 return tdtdt;
00187 }
00188 else
00189 {
00190 return tmp<GeometricField<Type, fvPatchField, volMesh> >
00191 (
00192 new GeometricField<Type, fvPatchField, volMesh>
00193 (
00194 ddtIOobject,
00195 mesh(),
00196 dimensioned<Type>
00197 (
00198 "0",
00199 dt.dimensions()/dimTime,
00200 pTraits<Type>::zero
00201 ),
00202 calculatedFvPatchField<Type>::typeName
00203 )
00204 );
00205 }
00206 }
00207
00208
00209 template<class Type>
00210 tmp<GeometricField<Type, fvPatchField, volMesh> >
00211 CoEulerDdtScheme<Type>::fvcDdt
00212 (
00213 const GeometricField<Type, fvPatchField, volMesh>& vf
00214 )
00215 {
00216 volScalarField rDeltaT = CorDeltaT();
00217
00218 IOobject ddtIOobject
00219 (
00220 "ddt("+vf.name()+')',
00221 mesh().time().timeName(),
00222 mesh()
00223 );
00224
00225 if (mesh().moving())
00226 {
00227 return tmp<GeometricField<Type, fvPatchField, volMesh> >
00228 (
00229 new GeometricField<Type, fvPatchField, volMesh>
00230 (
00231 ddtIOobject,
00232 mesh(),
00233 rDeltaT.dimensions()*vf.dimensions(),
00234 rDeltaT.internalField()*
00235 (
00236 vf.internalField()
00237 - vf.oldTime().internalField()*mesh().V0()/mesh().V()
00238 ),
00239 rDeltaT.boundaryField()*
00240 (
00241 vf.boundaryField() - vf.oldTime().boundaryField()
00242 )
00243 )
00244 );
00245 }
00246 else
00247 {
00248 return tmp<GeometricField<Type, fvPatchField, volMesh> >
00249 (
00250 new GeometricField<Type, fvPatchField, volMesh>
00251 (
00252 ddtIOobject,
00253 rDeltaT*(vf - vf.oldTime())
00254 )
00255 );
00256 }
00257 }
00258
00259
00260 template<class Type>
00261 tmp<GeometricField<Type, fvPatchField, volMesh> >
00262 CoEulerDdtScheme<Type>::fvcDdt
00263 (
00264 const dimensionedScalar& rho,
00265 const GeometricField<Type, fvPatchField, volMesh>& vf
00266 )
00267 {
00268 volScalarField rDeltaT = CorDeltaT();
00269
00270 IOobject ddtIOobject
00271 (
00272 "ddt("+rho.name()+','+vf.name()+')',
00273 mesh().time().timeName(),
00274 mesh()
00275 );
00276
00277 if (mesh().moving())
00278 {
00279 return tmp<GeometricField<Type, fvPatchField, volMesh> >
00280 (
00281 new GeometricField<Type, fvPatchField, volMesh>
00282 (
00283 ddtIOobject,
00284 mesh(),
00285 rDeltaT.dimensions()*rho.dimensions()*vf.dimensions(),
00286 rDeltaT.internalField()*rho.value()*
00287 (
00288 vf.internalField()
00289 - vf.oldTime().internalField()*mesh().V0()/mesh().V()
00290 ),
00291 rDeltaT.boundaryField()*rho.value()*
00292 (
00293 vf.boundaryField() - vf.oldTime().boundaryField()
00294 )
00295 )
00296 );
00297 }
00298 else
00299 {
00300 return tmp<GeometricField<Type, fvPatchField, volMesh> >
00301 (
00302 new GeometricField<Type, fvPatchField, volMesh>
00303 (
00304 ddtIOobject,
00305 rDeltaT*rho*(vf - vf.oldTime())
00306 )
00307 );
00308 }
00309 }
00310
00311
00312 template<class Type>
00313 tmp<GeometricField<Type, fvPatchField, volMesh> >
00314 CoEulerDdtScheme<Type>::fvcDdt
00315 (
00316 const volScalarField& rho,
00317 const GeometricField<Type, fvPatchField, volMesh>& vf
00318 )
00319 {
00320 volScalarField rDeltaT = CorDeltaT();
00321
00322 IOobject ddtIOobject
00323 (
00324 "ddt("+rho.name()+','+vf.name()+')',
00325 mesh().time().timeName(),
00326 mesh()
00327 );
00328
00329 if (mesh().moving())
00330 {
00331 return tmp<GeometricField<Type, fvPatchField, volMesh> >
00332 (
00333 new GeometricField<Type, fvPatchField, volMesh>
00334 (
00335 ddtIOobject,
00336 mesh(),
00337 rDeltaT.dimensions()*rho.dimensions()*vf.dimensions(),
00338 rDeltaT.internalField()*
00339 (
00340 rho.internalField()*vf.internalField()
00341 - rho.oldTime().internalField()
00342 *vf.oldTime().internalField()*mesh().V0()/mesh().V()
00343 ),
00344 rDeltaT.boundaryField()*
00345 (
00346 rho.boundaryField()*vf.boundaryField()
00347 - rho.oldTime().boundaryField()
00348 *vf.oldTime().boundaryField()
00349 )
00350 )
00351 );
00352 }
00353 else
00354 {
00355 return tmp<GeometricField<Type, fvPatchField, volMesh> >
00356 (
00357 new GeometricField<Type, fvPatchField, volMesh>
00358 (
00359 ddtIOobject,
00360 rDeltaT*(rho*vf - rho.oldTime()*vf.oldTime())
00361 )
00362 );
00363 }
00364 }
00365
00366
00367 template<class Type>
00368 tmp<fvMatrix<Type> >
00369 CoEulerDdtScheme<Type>::fvmDdt
00370 (
00371 GeometricField<Type, fvPatchField, volMesh>& vf
00372 )
00373 {
00374 tmp<fvMatrix<Type> > tfvm
00375 (
00376 new fvMatrix<Type>
00377 (
00378 vf,
00379 vf.dimensions()*dimVol/dimTime
00380 )
00381 );
00382
00383 fvMatrix<Type>& fvm = tfvm();
00384
00385 scalarField rDeltaT = CorDeltaT()().internalField();
00386
00387 fvm.diag() = rDeltaT*mesh().V();
00388
00389 if (mesh().moving())
00390 {
00391 fvm.source() = rDeltaT*vf.oldTime().internalField()*mesh().V0();
00392 }
00393 else
00394 {
00395 fvm.source() = rDeltaT*vf.oldTime().internalField()*mesh().V();
00396 }
00397
00398 return tfvm;
00399 }
00400
00401
00402 template<class Type>
00403 tmp<fvMatrix<Type> >
00404 CoEulerDdtScheme<Type>::fvmDdt
00405 (
00406 const dimensionedScalar& rho,
00407 GeometricField<Type, fvPatchField, volMesh>& vf
00408 )
00409 {
00410 tmp<fvMatrix<Type> > tfvm
00411 (
00412 new fvMatrix<Type>
00413 (
00414 vf,
00415 rho.dimensions()*vf.dimensions()*dimVol/dimTime
00416 )
00417 );
00418 fvMatrix<Type>& fvm = tfvm();
00419
00420 scalarField rDeltaT = CorDeltaT()().internalField();
00421
00422 fvm.diag() = rDeltaT*rho.value()*mesh().V();
00423
00424 if (mesh().moving())
00425 {
00426 fvm.source() = rDeltaT
00427 *rho.value()*vf.oldTime().internalField()*mesh().V0();
00428 }
00429 else
00430 {
00431 fvm.source() = rDeltaT
00432 *rho.value()*vf.oldTime().internalField()*mesh().V();
00433 }
00434
00435 return tfvm;
00436 }
00437
00438
00439 template<class Type>
00440 tmp<fvMatrix<Type> >
00441 CoEulerDdtScheme<Type>::fvmDdt
00442 (
00443 const volScalarField& rho,
00444 GeometricField<Type, fvPatchField, volMesh>& vf
00445 )
00446 {
00447 tmp<fvMatrix<Type> > tfvm
00448 (
00449 new fvMatrix<Type>
00450 (
00451 vf,
00452 rho.dimensions()*vf.dimensions()*dimVol/dimTime
00453 )
00454 );
00455 fvMatrix<Type>& fvm = tfvm();
00456
00457 scalarField rDeltaT = CorDeltaT()().internalField();
00458
00459 fvm.diag() = rDeltaT*rho.internalField()*mesh().V();
00460
00461 if (mesh().moving())
00462 {
00463 fvm.source() = rDeltaT
00464 *rho.oldTime().internalField()
00465 *vf.oldTime().internalField()*mesh().V0();
00466 }
00467 else
00468 {
00469 fvm.source() = rDeltaT
00470 *rho.oldTime().internalField()
00471 *vf.oldTime().internalField()*mesh().V();
00472 }
00473
00474 return tfvm;
00475 }
00476
00477
00478 template<class Type>
00479 tmp<typename CoEulerDdtScheme<Type>::fluxFieldType>
00480 CoEulerDdtScheme<Type>::fvcDdtPhiCorr
00481 (
00482 const volScalarField& rA,
00483 const GeometricField<Type, fvPatchField, volMesh>& U,
00484 const fluxFieldType& phi
00485 )
00486 {
00487 IOobject ddtIOobject
00488 (
00489 "ddtPhiCorr(" + rA.name() + ',' + U.name() + ',' + phi.name() + ')',
00490 mesh().time().timeName(),
00491 mesh()
00492 );
00493
00494 if (mesh().moving())
00495 {
00496 return tmp<fluxFieldType>
00497 (
00498 new fluxFieldType
00499 (
00500 ddtIOobject,
00501 mesh(),
00502 dimensioned<typename flux<Type>::type>
00503 (
00504 "0",
00505 rA.dimensions()*phi.dimensions()/dimTime,
00506 pTraits<typename flux<Type>::type>::zero
00507 )
00508 )
00509 );
00510 }
00511 else
00512 {
00513 volScalarField rDeltaT = CorDeltaT();
00514
00515 return tmp<fluxFieldType>
00516 (
00517 new fluxFieldType
00518 (
00519 ddtIOobject,
00520 this->fvcDdtPhiCoeff(U.oldTime(), phi.oldTime())*
00521 (
00522 fvc::interpolate(rDeltaT*rA)*phi.oldTime()
00523 - (fvc::interpolate(rDeltaT*rA*U.oldTime()) & mesh().Sf())
00524 )
00525 )
00526 );
00527 }
00528 }
00529
00530
00531 template<class Type>
00532 tmp<typename CoEulerDdtScheme<Type>::fluxFieldType>
00533 CoEulerDdtScheme<Type>::fvcDdtPhiCorr
00534 (
00535 const volScalarField& rA,
00536 const volScalarField& rho,
00537 const GeometricField<Type, fvPatchField, volMesh>& U,
00538 const fluxFieldType& phi
00539 )
00540 {
00541 IOobject ddtIOobject
00542 (
00543 "ddtPhiCorr("
00544 + rA.name() + ',' + rho.name() + ',' + U.name() + ',' + phi.name() + ')',
00545 mesh().time().timeName(),
00546 mesh()
00547 );
00548
00549 if (mesh().moving())
00550 {
00551 return tmp<fluxFieldType>
00552 (
00553 new fluxFieldType
00554 (
00555 ddtIOobject,
00556 mesh(),
00557 dimensioned<typename flux<Type>::type>
00558 (
00559 "0",
00560 rA.dimensions()*rho.dimensions()*phi.dimensions()/dimTime,
00561 pTraits<typename flux<Type>::type>::zero
00562 )
00563 )
00564 );
00565 }
00566 else
00567 {
00568 volScalarField rDeltaT = CorDeltaT();
00569
00570 if
00571 (
00572 U.dimensions() == dimVelocity
00573 && phi.dimensions() == dimVelocity*dimArea
00574 )
00575 {
00576 return tmp<fluxFieldType>
00577 (
00578 new fluxFieldType
00579 (
00580 ddtIOobject,
00581 this->fvcDdtPhiCoeff(U.oldTime(), phi.oldTime())
00582 *(
00583 fvc::interpolate(rDeltaT*rA*rho.oldTime())*phi.oldTime()
00584 - (fvc::interpolate(rDeltaT*rA*rho.oldTime()*U.oldTime())
00585 & mesh().Sf())
00586 )
00587 )
00588 );
00589 }
00590 else if
00591 (
00592 U.dimensions() == dimVelocity
00593 && phi.dimensions() == dimDensity*dimVelocity*dimArea
00594 )
00595 {
00596 return tmp<fluxFieldType>
00597 (
00598 new fluxFieldType
00599 (
00600 ddtIOobject,
00601 this->fvcDdtPhiCoeff
00602 (
00603 U.oldTime(),
00604 phi.oldTime()/fvc::interpolate(rho.oldTime())
00605 )
00606 *(
00607 fvc::interpolate(rDeltaT*rA*rho.oldTime())
00608 *phi.oldTime()/fvc::interpolate(rho.oldTime())
00609 - (
00610 fvc::interpolate
00611 (
00612 rDeltaT*rA*rho.oldTime()*U.oldTime()
00613 ) & mesh().Sf()
00614 )
00615 )
00616 )
00617 );
00618 }
00619 else if
00620 (
00621 U.dimensions() == dimDensity*dimVelocity
00622 && phi.dimensions() == dimDensity*dimVelocity*dimArea
00623 )
00624 {
00625 return tmp<fluxFieldType>
00626 (
00627 new fluxFieldType
00628 (
00629 ddtIOobject,
00630 this->fvcDdtPhiCoeff
00631 (
00632 rho.oldTime(),
00633 U.oldTime(),
00634 phi.oldTime()
00635 )
00636 *(
00637 fvc::interpolate(rDeltaT*rA)*phi.oldTime()
00638 - (
00639 fvc::interpolate(rDeltaT*rA*U.oldTime())&mesh().Sf()
00640 )
00641 )
00642 )
00643 );
00644 }
00645 else
00646 {
00647 FatalErrorIn
00648 (
00649 "CoEulerDdtScheme<Type>::fvcDdtPhiCorr"
00650 ) << "dimensions of phi are not correct"
00651 << abort(FatalError);
00652
00653 return fluxFieldType::null();
00654 }
00655 }
00656 }
00657
00658
00659 template<class Type>
00660 tmp<surfaceScalarField> CoEulerDdtScheme<Type>::meshPhi
00661 (
00662 const GeometricField<Type, fvPatchField, volMesh>&
00663 )
00664 {
00665 return tmp<surfaceScalarField>
00666 (
00667 new surfaceScalarField
00668 (
00669 IOobject
00670 (
00671 "meshPhi",
00672 mesh().time().timeName(),
00673 mesh()
00674 ),
00675 mesh(),
00676 dimensionedScalar("0", dimVolume/dimTime, 0.0)
00677 )
00678 );
00679 }
00680
00681
00682
00683
00684 }
00685
00686
00687
00688 }
00689
00690