View Issue Details
ID | Project | Category | View Status | Date Submitted | Last Update |
---|---|---|---|---|---|
0003068 | OpenFOAM | Bug | public | 2018-09-06 13:32 | 2018-09-08 01:36 |
Reporter | sose | Assigned To | henry | ||
Priority | normal | Severity | major | Reproducibility | always |
Status | closed | Resolution | suspended | ||
Summary | 0003068: wrong coefficients for meshPhi in backward ddtScheme | ||||
Description | For the backward ddtScheme the implementation of meshPhi is valid only for constant timestep size. To account correctly for the variable timestep size the coefficients should be as shown below. template<class Type> tmp<surfaceScalarField> backwardDdtScheme<Type>::meshPhi ( const GeometricField<Type, fvPatchField, volMesh>& vf ) { ... current implementation // Coefficient for t-3/2 (between times 0 and 00) scalar coefft0_00 = deltaT/(deltaT + deltaT0); // Coefficient for t-1/2 (between times n and 0) scalar coefftn_0 = 1 + coefft0_00; ... the coefficients should have been scalar coefft0_00 = deltaT*deltaT/(deltaT + deltaT0)/deltaT0; scalar coefftn_0 = 1 + deltaT/(deltaT + deltaT0); } | ||||
Steps To Reproduce | for backward scheme, we have sum(meshPhi) = (1/h) (c * V - c0 * V0 + c00 * V00) where c = 1 + h/(h+h0) ; c0 = c + c00 ; c00 = h*h/h0/(h+h0) this can be expanded to sum(meshPhi) = (1/h) (c*(V - V0) + (c-c0) *(V0 - V00) + (c00+c-c0)*V00) in terms of swept-volume, we have mesh.phi() defined as mesh.phi() := deltaV/h and the expression for each cell face becomes meshPhi = c*mesh.phi() + (c-c0) * mesh.phi().oldTime() + (c00+c-c0)*V00 where it is clear that c = 1 + h/(h+h0) this is coefftn_0 (c-c0) = c00 = h*h/h0/(h+h0) this is coefft0_00 (c00+c-c0) = 0 | ||||
Tags | No tags attached. | ||||
|
Can you provide a test-case which demonstrates the need for this change? |
|
thank you for your quick response, in my opinion there is no need for a test case. The mathematical proof shown in the description is more than enough. |
|
This is a significant change to a piece of code which has already been tested so we will need a case which reproduces the problem and demonstrates that your change is correct. Presumably you already have a case which demonstrates a need for this change? |
|
As already indicated this bug affects only variable time stepping. For constant time step, the new and old formulations produce the same value for coefft0_00. I understand your point, but unfortunately I don't have a test case. I have only a mathematical proof, and the math is very basic as shown in the description ... PS: the second line from bottom should be " - (c-c0) = c00 ", please note that the "minus" is already accounted for in the current implementation, except that the definition of c00 is missing "h/h0", and for a constant time step this ratio is 1 which does not change anything. For variable time step, however, the resulting coefficient will be different. |
|
If it not possible to demonstrate this change then I cannot make it. This part of the code has been developed and tested carefully and used successfully regularly, it would be insane to make changes here without any means to demonstrate either the need for or correctness of the change. I will suspend this report until you can provide a case which reproduces the problem and demonstrates that your change is correct. |
Date Modified | Username | Field | Change |
---|---|---|---|
2018-09-06 13:32 | sose | New Issue | |
2018-09-06 13:45 | henry | Note Added: 0010040 | |
2018-09-06 13:50 | sose | Note Added: 0010041 | |
2018-09-06 14:22 | henry | Note Added: 0010042 | |
2018-09-06 14:41 | sose | Note Added: 0010043 | |
2018-09-06 16:23 | henry | Note Added: 0010045 | |
2018-09-06 16:24 | henry | Assigned To | => henry |
2018-09-06 16:24 | henry | Status | new => closed |
2018-09-06 16:24 | henry | Resolution | open => suspended |