View Issue Details
ID | Project | Category | View Status | Date Submitted | Last Update |
---|---|---|---|---|---|
0000371 | OpenFOAM | Bug | public | 2011-12-22 22:59 | 2011-12-23 09:25 |
Reporter | albertop | Assigned To | |||
Priority | normal | Severity | major | Reproducibility | always |
Status | resolved | Resolution | fixed | ||
Platform | x86_64 | OS | openSUSE | OS Version | 12.1 |
Summary | 0000371: Incorrect behaviour of PIMPLE loop in unsteady simulations | ||||
Description | PIMPLE solvers seem to have an incorrect behaviour when the residualControl feature is used: - The solver correctly performs the number of PIMPLE correctors if the solution does *not* converge. - When the solution converges, the PIMPLE iteration is interrupted, the time is advanced, but the number of PIMPLE iterations is not re-initialized. | ||||
Steps To Reproduce | Run, for example, the tutorial "angledDuct" for rhoPimpleFoam. Observe what happens at the time steps from Time = 6 to Time = 10. - The solution converges in 45 iterations (50 is the maximum allowed) at Time = 6. The loop is interrupted, and time advanced (correct behavior). - The new time-step starts (Time = 7), but the code states "PIMPLE: iteration 46", then it performs one iteration, and it states the solution converged in 47 iterations. - The time is advanced again (Time = 8), and the code states "PIMPLE: iteration 48", it does one iteration, and states "PIMPLE: converged in 49 iterations". - Same behaviour at Time = 9, but the solution does not converge "in 50 iterations". At this point, Time is increased again, and PIMPLE iterations start from zero for the new time-step. See attached log for the output of the solver. | ||||
Additional Information | If my understanding is correct, the expected behavior should be that at the end of the time-step, if the solution is converged, the number of PIMPLE iterations is reset to zero, so that 50 iterations are allowed for the new time-step. The current behaviour leads to incorrect results in transient simulations. | ||||
Tags | No tags attached. | ||||
|
output.log (8,340 bytes)
PIMPLE: iteration 42 DILUPBiCG: Solving for Ux, Initial residual = 1.06776e-05, Final residual = 3.21993e-07, No Iterations 1 DILUPBiCG: Solving for Uy, Initial residual = 1.22935e-05, Final residual = 1.38374e-06, No Iterations 1 DILUPBiCG: Solving for Uz, Initial residual = 0.000132657, Final residual = 4.34948e-06, No Iterations 1 DILUPBiCG: Solving for h, Initial residual = 8.6855e-05, Final residual = 3.13538e-06, No Iterations 1 DICPCG: Solving for p, Initial residual = 4.89827e-05, Final residual = 8.20175e-07, No Iterations 72 diagonal: Solving for rho, Initial residual = 0, Final residual = 0, No Iterations 0 time step continuity errors : sum local = 0.000160787, global = 6.49057e-07, cumulative = 0.172084 rho max/min : 1.19331 1.1614 PIMPLE: iteration 43 DILUPBiCG: Solving for Ux, Initial residual = 8.8961e-06, Final residual = 8.8961e-06, No Iterations 0 DILUPBiCG: Solving for Uy, Initial residual = 1.02509e-05, Final residual = 1.07183e-06, No Iterations 1 DILUPBiCG: Solving for Uz, Initial residual = 0.000107678, Final residual = 3.72043e-06, No Iterations 1 DILUPBiCG: Solving for h, Initial residual = 7.87623e-05, Final residual = 2.90994e-06, No Iterations 1 DICPCG: Solving for p, Initial residual = 7.77476e-05, Final residual = 8.92265e-07, No Iterations 74 diagonal: Solving for rho, Initial residual = 0, Final residual = 0, No Iterations 0 time step continuity errors : sum local = 0.000174919, global = -1.6737e-06, cumulative = 0.172082 rho max/min : 1.19331 1.1614 PIMPLE: iteration 44 DILUPBiCG: Solving for Ux, Initial residual = 1.00314e-05, Final residual = 3.44023e-07, No Iterations 1 DILUPBiCG: Solving for Uy, Initial residual = 8.15748e-06, Final residual = 8.15748e-06, No Iterations 0 DILUPBiCG: Solving for Uz, Initial residual = 9.34829e-05, Final residual = 3.00462e-06, No Iterations 1 DILUPBiCG: Solving for h, Initial residual = 7.02682e-05, Final residual = 2.93709e-06, No Iterations 1 DICPCG: Solving for p, Initial residual = 6.18899e-05, Final residual = 9.70891e-07, No Iterations 45 diagonal: Solving for rho, Initial residual = 0, Final residual = 0, No Iterations 0 time step continuity errors : sum local = 0.000190324, global = 1.40623e-05, cumulative = 0.172096 rho max/min : 1.19331 1.1614 PIMPLE: converged in 45 iterations ExecutionTime = 31.91 s ClockTime = 32 s Courant Number mean: 7181.37 max: 16740 Time = 7 diagonal: Solving for rho, Initial residual = 0, Final residual = 0, No Iterations 0 PIMPLE: iteration 46 DILUPBiCG: Solving for Ux, Initial residual = 8.05191e-06, Final residual = 8.05191e-06, No Iterations 0 DILUPBiCG: Solving for Uy, Initial residual = 8.65461e-06, Final residual = 8.65461e-06, No Iterations 0 DILUPBiCG: Solving for Uz, Initial residual = 6.08186e-05, Final residual = 3.36791e-06, No Iterations 1 DILUPBiCG: Solving for h, Initial residual = 4.81271e-05, Final residual = 1.65391e-06, No Iterations 1 DICPCG: Solving for p, Initial residual = 4.70836e-05, Final residual = 9.72231e-07, No Iterations 25 diagonal: Solving for rho, Initial residual = 0, Final residual = 0, No Iterations 0 time step continuity errors : sum local = 0.000190587, global = 8.0717e-06, cumulative = 0.172104 rho max/min : 1.19331 1.1614 PIMPLE: converged in 47 iterations ExecutionTime = 32.01 s ClockTime = 32 s Courant Number mean: 7181.39 max: 16740 Time = 8 diagonal: Solving for rho, Initial residual = 0, Final residual = 0, No Iterations 0 PIMPLE: iteration 48 DILUPBiCG: Solving for Ux, Initial residual = 7.22976e-06, Final residual = 7.22976e-06, No Iterations 0 DILUPBiCG: Solving for Uy, Initial residual = 8.51094e-06, Final residual = 8.51094e-06, No Iterations 0 DILUPBiCG: Solving for Uz, Initial residual = 3.84221e-05, Final residual = 3.42924e-06, No Iterations 1 DILUPBiCG: Solving for h, Initial residual = 3.36357e-05, Final residual = 1.09803e-06, No Iterations 1 DICPCG: Solving for p, Initial residual = 2.84791e-05, Final residual = 9.55147e-07, No Iterations 22 diagonal: Solving for rho, Initial residual = 0, Final residual = 0, No Iterations 0 time step continuity errors : sum local = 0.000187239, global = 6.55787e-06, cumulative = 0.172111 rho max/min : 1.19331 1.1614 PIMPLE: converged in 49 iterations ExecutionTime = 32.11 s ClockTime = 32 s Courant Number mean: 7181.4 max: 16740 Time = 9 diagonal: Solving for rho, Initial residual = 0, Final residual = 0, No Iterations 0 PIMPLE: iteration 50 DILUPBiCG: Solving for Ux, Initial residual = 6.716e-06, Final residual = 6.716e-06, No Iterations 0 DILUPBiCG: Solving for Uy, Initial residual = 8.14389e-06, Final residual = 8.14389e-06, No Iterations 0 DILUPBiCG: Solving for Uz, Initial residual = 2.97565e-05, Final residual = 3.85919e-06, No Iterations 1 DILUPBiCG: Solving for h, Initial residual = 2.57404e-05, Final residual = 8.36795e-07, No Iterations 1 DICPCG: Solving for p, Initial residual = 2.53165e-05, Final residual = 9.34994e-07, No Iterations 16 diagonal: Solving for rho, Initial residual = 0, Final residual = 0, No Iterations 0 time step continuity errors : sum local = 0.000183288, global = 9.17045e-06, cumulative = 0.17212 rho max/min : 1.19331 1.1614 DILUPBiCG: Solving for epsilon, Initial residual = 0.0243377, Final residual = 3.0661e-06, No Iterations 4 DILUPBiCG: Solving for k, Initial residual = 0.104572, Final residual = 5.67177e-06, No Iterations 5 PIMPLE: not converged within 50 iterations ExecutionTime = 32.24 s ClockTime = 32 s Courant Number mean: 7181.41 max: 16740 Time = 10 diagonal: Solving for rho, Initial residual = 0, Final residual = 0, No Iterations 0 PIMPLE: iteration 1 DILUPBiCG: Solving for Ux, Initial residual = 0.000331171, Final residual = 3.22179e-05, No Iterations 1 DILUPBiCG: Solving for Uy, Initial residual = 0.000200863, Final residual = 1.7147e-05, No Iterations 1 DILUPBiCG: Solving for Uz, Initial residual = 0.00355924, Final residual = 0.000231287, No Iterations 1 DILUPBiCG: Solving for h, Initial residual = 0.000300133, Final residual = 2.30866e-05, No Iterations 1 DICPCG: Solving for p, Initial residual = 0.00494924, Final residual = 4.76067e-05, No Iterations 46 diagonal: Solving for rho, Initial residual = 0, Final residual = 0, No Iterations 0 time step continuity errors : sum local = 0.00932971, global = -5.98021e-05, cumulative = 0.17206 rho max/min : 1.19329 1.16137 PIMPLE: iteration 2 DILUPBiCG: Solving for Ux, Initial residual = 0.000353627, Final residual = 2.88376e-05, No Iterations 1 DILUPBiCG: Solving for Uy, Initial residual = 0.000191335, Final residual = 1.85767e-05, No Iterations 1 DILUPBiCG: Solving for Uz, Initial residual = 0.00270185, Final residual = 0.000255215, No Iterations 1 DILUPBiCG: Solving for h, Initial residual = 0.00164716, Final residual = 1.4852e-05, No Iterations 2 DICPCG: Solving for p, Initial residual = 0.00568365, Final residual = 4.74496e-05, No Iterations 44 diagonal: Solving for rho, Initial residual = 0, Final residual = 0, No Iterations 0 time step continuity errors : sum local = 0.0092973, global = 0.000964444, cumulative = 0.173025 rho max/min : 1.19328 1.16136 PIMPLE: iteration 3 DILUPBiCG: Solving for Ux, Initial residual = 0.000376318, Final residual = 4.13185e-06, No Iterations 2 DILUPBiCG: Solving for Uy, Initial residual = 0.000180439, Final residual = 1.75917e-05, No Iterations 1 DILUPBiCG: Solving for Uz, Initial residual = 0.0026057, Final residual = 3.92159e-05, No Iterations 2 DILUPBiCG: Solving for h, Initial residual = 0.0019651, Final residual = 0.000127436, No Iterations 1 DICPCG: Solving for p, Initial residual = 0.00533171, Final residual = 4.95756e-05, No Iterations 73 diagonal: Solving for rho, Initial residual = 0, Final residual = 0, No Iterations 0 time step continuity errors : sum local = 0.00972879, global = -0.000147045, cumulative = 0.172878 rho max/min : 1.19327 1.16137 PIMPLE: iteration 4 DILUPBiCG: Solving for Ux, Initial residual = 0.000336712, Final residual = 5.87149e-06, No Iterations 2 DILUPBiCG: Solving for Uy, Initial residual = 0.000161355, Final residual = 1.68774e-06, No Iterations 2 DILUPBiCG: Solving for Uz, Initial residual = 0.00207619, Final residual = 2.98622e-05, No Iterations 2 DILUPBiCG: Solving for h, Initial residual = 0.00194045, Final residual = 2.13209e-05, No Iterations 2 |
|
I attach a proposed patch, where corr_ is reset to zero in pimpleControl.C, if convergence criteria are satisfied: 00197 if (criteriaSatisfied()) 00198 { 00199 Info<< algorithmName_ << ": converged in " << corr_ << " iterations" 00200 << endl; 00201 completed = true; corr_ = 0; // Added code 00202 } |
|
pimpleControl.C.patch (343 bytes)
*** pimpleControlOld.C 2011-12-23 00:57:19.000000000 +0100 --- pimpleControl.C 2011-12-23 01:01:53.332168358 +0100 *************** *** 199,204 **** --- 199,205 ---- Info<< algorithmName_ << ": converged in " << corr_ << " iterations" << endl; completed = true; + corr_ = 0; } else { |
|
Thanks for the report - bug corrected in commit: a22b8fd6d4d3607b5d8a4f00776997ba59b0f4c7 |
Date Modified | Username | Field | Change |
---|---|---|---|
2011-12-22 22:59 | albertop | New Issue | |
2011-12-22 22:59 | albertop | File Added: output.log | |
2011-12-23 00:05 | albertop | Note Added: 0000865 | |
2011-12-23 00:06 | albertop | File Added: pimpleControl.C.patch | |
2011-12-23 09:25 |
|
Note Added: 0000866 | |
2011-12-23 09:25 |
|
Status | new => resolved |
2011-12-23 09:25 |
|
Fixed in Version | => 2.1.x |
2011-12-23 09:25 |
|
Resolution | open => fixed |
2011-12-23 09:25 |
|
Assigned To | => user2 |