View Issue Details
ID | Project | Category | View Status | Date Submitted | Last Update |
---|---|---|---|---|---|
0002020 | OpenFOAM | Bug | public | 2016-03-08 12:26 | 2016-03-08 15:28 |
Reporter | Shorty | Assigned To | henry | ||
Priority | high | Severity | crash | Reproducibility | always |
Status | resolved | Resolution | fixed | ||
Platform | GNU/Linux | OS | Ubuntu | OS Version | 14.04 |
Summary | 0002020: Field access somehow wrong | ||||
Description | Dear all, I am sure this is a critical bug. Using fixedShearStress BC will fail always due to division by zero here: operator==(tauHat*(tauHat & (tau0_*(1.0/(ry*nuEff)) + Uc))); The problem is that nuEff is initialized wrong (field is 0()); const label patchI = patch().index(); const surfaceScalarField& phi = db().lookupObject<surfaceScalarField>(phiName_); scalarField nuEff; if (phi.dimensions() == dimVelocity*dimArea) { const incompressible::turbulenceModel& turbModel = db().lookupObject<incompressible::turbulenceModel> ( "turbulenceModel" ); //- This is somehow wrong nuEff = turbModel.nuEff()()[patchI]; //- Output it Info<< nuEff -> 0() } Interesting that Info<< nuEff << is different to Info<< turbModel.nuEff()()[patchI] Maybe its correct but in my opinion nuEff should have the values of nu of each face at the patch patchI. | ||||
Steps To Reproduce | For incompressible solvers, it does not matter. I always used simpleFoam... division by zero. | ||||
Tags | No tags attached. | ||||
|
|
|
It is not clear from your description how your case is setup such that both the laminar and turbulent components of the viscosity are both zero. To study the problem further we will need a test-case for OpenFOAM-3.0.x or OpenFOAM-dev which reproduces the problem. |
|
Dear Henry, I am sorry for bad explanation. I will prepare a FOAM3.0 case but at least it should be the same. Just for clearness... the output of the field turbModel.nuEff() will show me all boundaries with non-zero values (everything is correct there). For example: The patch name is wall1, therefore I got the output: wall1 { type calculated; value nonuniform List<scalar> 500 ( 0.285685 0.285685 0.285685 0.285685 Please correct me if I am wrong: turbModel.nuEff()()[patchI] should return the face values of the patch (patchI), or am I wrong? But as far as I understand operator()()[patchI], it will return a scalar and not a scalar field. |
|
It is not clear to me what is wrong with your case setup and I will look into it whe you provide a version for OpenFOAM-3.0.x or dev. |
|
The code seems to use the patch index (patchi) to index into the internalField of a (vol)Field. Instead index into the boundaryField(), something like: const scalarField& bla = turbModel.nuEff()().boundaryField()[patchI]; |
|
scalarField nuEff(turbModel.nuEff(patch().index())); is cleaner and more efficient. |
|
I have made this change to OpenFOAM-dev: commit 9ef051dd18a345f96ea5e9380f82c5f5f20e96ab OpenFOAM-3.0.x: commit 09fa0b46eb4394f632f046b6c41fe57e8a140c7c |
|
Thanks henry, I changed it locally and now it is working. I think there is no need for a test case anymore now, or? |
|
No need for the test-case, I will close this report. |
|
Resolved in OpenFOAM-dev by commit 9ef051dd18a345f96ea5e9380f82c5f5f20e96ab Resolved in OpenFOAM-3.0.x by commit 09fa0b46eb4394f632f046b6c41fe57e8a140c7c |
Date Modified | Username | Field | Change |
---|---|---|---|
2016-03-08 12:26 | Shorty | New Issue | |
2016-03-08 12:26 | Shorty | File Added: ABLz0FSS.zip | |
2016-03-08 12:48 | henry | Note Added: 0006008 | |
2016-03-08 13:26 | Shorty | Note Added: 0006009 | |
2016-03-08 13:27 | Shorty | Note Edited: 0006009 | |
2016-03-08 13:43 | Shorty | Note Edited: 0006009 | |
2016-03-08 13:48 | henry | Note Added: 0006010 | |
2016-03-08 13:55 | MattijsJ | Note Added: 0006013 | |
2016-03-08 14:02 | henry | Note Added: 0006014 | |
2016-03-08 14:19 | henry | Note Added: 0006017 | |
2016-03-08 14:57 | Shorty | Note Added: 0006021 | |
2016-03-08 15:27 | henry | Note Added: 0006023 | |
2016-03-08 15:28 | henry | Note Added: 0006024 | |
2016-03-08 15:28 | henry | Status | new => resolved |
2016-03-08 15:28 | henry | Fixed in Version | => 3.0.x |
2016-03-08 15:28 | henry | Resolution | open => fixed |
2016-03-08 15:28 | henry | Assigned To | => henry |