View Issue Details
ID | Project | Category | View Status | Date Submitted | Last Update |
---|---|---|---|---|---|
0000255 | OpenFOAM | Bug | public | 2011-07-18 16:50 | 2011-07-23 18:49 |
Reporter | Assigned To | henry | |||
Priority | normal | Severity | tweak | Reproducibility | N/A |
Status | resolved | Resolution | fixed | ||
Platform | Linux | OS | Ubuntu | OS Version | 10.04 |
Summary | 0000255: missing coefficient | ||||
Description | In homogeneousDynSmagorinsky.C, the definition of cD seems inconsistent with the derivation of Lilly (1992). In Lilly's paper, there is an extra '2' coefficient in the denominator of equation 11 which is not seen in the openfoam implementation. The homogeneousDynSmagorinsky.C implementation is consistent with Fureby et al. (1997), however, the '2' coefficient in front of the Mij term seems to mysteriously vanish in Fureby's derivation between equations 6 and 7. | ||||
Steps To Reproduce | N/A | ||||
Additional Information | Lilly, D.K., 1992, A proposed modification of the Germano subgrid-scale closure method. Phys. Fluids, 4, 633-635. C. Fureby, G. Tabor, H. G. Weller, and A. D. Gosman, 1997, A comparative study of subgrid scale models in homogeneous isotropic turbulence, Phys. Fluids, 9, 1416. | ||||
Tags | No tags attached. | ||||
|
The cD in the homogeneousDynSmagorinsky model corresponds to 2*C in Lilly's paper which is convenient because then the factor of 2 in the expression for nuSgs is not required, i.e. nuSgs = cD*delta^2*||D|| and cD=<L.M>/<M.M> are consistent. If you believe otherwise and think there is actually an error in the homogeneousDynSmagorinsky model please reopen this thread and provide details. |
|
I apologize for the ambiguity of my initial explanation. I will try to be more explicit. I will label my own equations here, some copied directly from the code, with '(#)' OpenFOAM accounts for the turbulent and molecular deviatoric stresses in the incompressible solvers using the divDevReff (RAS) or divDevBeff (LES) term: (1) tmp<fvVectorMatrix> GenEddyVisc::divDevBeff(volVectorField& U) const { return ( - fvm::laplacian(nuEff(), U) - fvc::div(nuEff()*dev(T(fvc::grad(U)))) ); } where the effective viscosity 'nuEff' is the sum of the molecular viscosity and 'eddy viscosity' (2) virtual tmp<volScalarField> nuEff() const { return tmp<volScalarField> ( new volScalarField("nuEff", nuSgs() + nu()) ); } the turbulent stresses are thus interpreted as: (3) tau = -2 * nuSgs * D (4) D = 0.5 * (grad(U) + grad(U).T()) (5) divDevBeff = div(tau) this is all good and fine with me so far (the trace of the gradient tensor is zero by continuity for incompressible flows) from the homogeneousDynSmagorinsky code: (6) const volSymmTensorField D(dev(symm(gradU))); // as defined above in eqn (4) (7) cD(D) = average(LL && MM)/MMMM; (8) const volSymmTensorField MM ( sqr(delta())*(filter_(mag(D)*(D)) - 4*mag(filter_(D))*filter_(D)) ); (9) nuSgs_ = cD(D)*sqr(delta())*sqrt(magSqr(D)); therefore, piecing it together from eqns (3) and (9) (10) tau = -2 * cD(D) * sqr(delta()) * sqrt(magSqr(D)) * D which is the same as Lilly et al (1992) equation 5, and the MM term is consistent with Lilly's equation 9 (to be exact, it differs by a factor of sqrt(2) from the sqrt(magSqr(D)) term, but this is accounted for because the MM equation also differs from Lilly's equation 9 by the same factor, hence cD(D) differs from Lilly's definition by the reciprocal) therefore, I am of the opinion that the definition of cD(D) in eqn (8) should include the leading coefficient of 1/2 as in Lilly's equation 11. I hope that I was more clear this time. Regards, Perry |
|
If the cD includes the 1/2 coefficient then the nuSgs expression would need the factor of 2: nuSgs = cD*2*delta^2*||D|| what would be the advantage of first dividing by 2 and then multiplying by 2? |
|
the '2' coefficient you suggest is not needed because it is already there (implicitly) in the divDevBeff equation, see equation (1) and (3) above |
|
Thanks for the clarification. I have now had time to study the relevant papers and I agree with your finding and applied the appropriate correction: commit 90da72aa43cdb7f795f9abb3a522f0455f07b73d Thanks for reporting. |
Date Modified | Username | Field | Change |
---|---|---|---|
2011-07-18 16:50 |
|
New Issue | |
2011-07-19 12:19 | henry | Note Added: 0000565 | |
2011-07-19 12:19 | henry | Status | new => closed |
2011-07-19 12:19 | henry | Assigned To | => henry |
2011-07-19 12:19 | henry | Resolution | open => fixed |
2011-07-19 14:32 |
|
Note Added: 0000568 | |
2011-07-19 14:32 |
|
Status | closed => feedback |
2011-07-19 14:32 |
|
Resolution | fixed => reopened |
2011-07-19 14:44 | henry | Note Added: 0000569 | |
2011-07-19 14:55 |
|
Note Added: 0000570 | |
2011-07-19 14:55 |
|
Status | feedback => assigned |
2011-07-23 18:49 | henry | Note Added: 0000587 | |
2011-07-23 18:49 | henry | Status | assigned => resolved |
2011-07-23 18:49 | henry | Resolution | reopened => fixed |