View Issue Details
ID | Project | Category | View Status | Date Submitted | Last Update |
---|---|---|---|---|---|
0000267 | OpenFOAM | Bug | public | 2011-08-03 11:57 | 2011-09-08 12:43 |
Reporter | Assigned To | henry | |||
Priority | normal | Severity | minor | Reproducibility | always |
Status | resolved | Resolution | fixed | ||
Platform | x84_64 | OS | Ubuntu | OS Version | 11.04 |
Summary | 0000267: zeta equation and fMu function in qZeta turbulence model wrong | ||||
Description | The implementation of the zeta equation is wrong in my opinion. It can be only a typo or by intention, but in the zeta equation tmp<fvScalarMatrix> zetaEqn ( fvm::ddt(zeta_) + fvm::div(phi_, zeta_) - fvm::laplacian(DzetaEff(), zeta_) == (2.0*C1_ - 1)*G*zeta_/q_ - fvm::Sp((2.0*C2_ - dimensionedScalar(1.0))*f2()*zeta_/q_, zeta_) + E ); the second to last term - fvm::Sp((2.0*C2_ - dimensionedScalar(1.0))*f2()*zeta_/q_, zeta_) should be - fvm::Sp((2.0*C2_*f2()_ - dimensionedScalar(1.0))*zeta_/q_, zeta_) in order to reproduce the relation C_zeta f_zeta = 2 C_eps f_eps - 1 as described in the paper just above equation (16) The other thing is an extra turbulence Reynolds number dependent factor in the model, that appears in the the fMu() function exp(-6.0/sqr(scalar(1) + Rt/50.0)) *(scalar(1) + 3.0*exp(-Rt/10.0)); The second factor scalar(1) + 3.0*exp(-Rt/10.0)) does not exist in the relation given in equation (18) of the paper. | ||||
Additional Information | Paper referred to is Gibson and Dafa'Alla, Two-Equation Model for Turbulent Wall Flow, AIAA Journal, vol. 33, No. 8, 1994 | ||||
Tags | Modelling, turbulence models | ||||
|
The implementation of the q-zeta in OpenFOAM is based on the final technical report on the work by Dafa'Alla, Juntasaro and Gibson which was subsequently presented and published: "Calculation of oscillating boundary layers with the q-[zeta] turbulence model" A.A. Dafa'Alla, E. Juntasaro and M.M. Gibson Engineering Turbulence Modelling and Experiments 3: Proceedings of the Third International Symposium, Crete, Greece, May 27-29, 1996, p. 141. Editors: Wolfgang Rodi and G. Bergeles in which alternative forms of the damping functions are presented. However - fvm::Sp((2.0*C2_ - dimensionedScalar(1.0))*f2()*zeta_/q_, zeta_) should be - fvm::SuSp((2.0*C2_*f2() - dimensionedScalar(1.0))*zeta_/q_, zeta_) and this is resolved by commit ebe324f10a942fd92f0fbb0aa9dade8765a18a90 |
Date Modified | Username | Field | Change |
---|---|---|---|
2011-08-03 11:57 |
|
New Issue | |
2011-08-03 11:59 |
|
Tag Attached: Modelling | |
2011-08-03 11:59 |
|
Tag Attached: turbulence models | |
2011-09-08 12:43 | henry | Note Added: 0000645 | |
2011-09-08 12:43 | henry | Status | new => resolved |
2011-09-08 12:43 | henry | Resolution | open => fixed |
2011-09-08 12:43 | henry | Assigned To | => henry |