[petsc-users] compare snes_mf_operator and snes_fd

Zou (Non-US), Ling ling.zou at inl.gov
Thu Jan 31 14:36:52 CST 2013


Thank you Jed. Yes I am using MOOSE. I've also tried snes_type test options
before and it worked like a charm to help me fix the Jacobian problem (as
you may remember I asked quite a lot snes_type test questions here before).

'-mat_mffd_type ds' is quite helpful as far as I can tell now. Here is a
table I made to show both '-mat_mffd_type ds' and
'-ksp_gmres_modifiedgramschmidt', both seems helpful.

'-mat_mffd_type ds'  |   '-ksp_gmres_modifiedgramschmidt'  |   # of
function call
No                                       No
                       36709
No                                       Yes
                        36014
Yes                                      No
                       22470
Yes                                     Yes
                       17862


Ultimately, I guess I need improve my Jacobian anyway.

Best,

Ling


On Thu, Jan 31, 2013 at 12:46 PM, Jed Brown <jedbrown at mcs.anl.gov> wrote:

> On Thu, Jan 31, 2013 at 1:29 PM, Matthew Knepley <knepley at gmail.com>wrote:
>
>> On Thu, Jan 31, 2013 at 2:25 PM, Zou (Non-US), Ling <ling.zou at inl.gov>wrote:
>>
>>>
>>>
>>> On Thu, Jan 31, 2013 at 11:28 AM, Matthew Knepley <knepley at gmail.com>wrote:
>>>
>>>> On Thu, Jan 31, 2013 at 1:16 PM, Zou (Non-US), Ling <ling.zou at inl.gov>wrote:
>>>>
>>>>> Thank you Matt and Barry. I didn't get a chance to reply you yesterday.
>>>>> Here are the new output files with -snes_view on.
>>>>>
>>>>
>>>> It seems clear that the matrix you are providing to snes_mf_operator is
>>>> not a good
>>>> preconditioner for the actual matrix obtained with snes_fd. Maybe you
>>>> have a bug in
>>>> your evaluation. Maybe you could try -snes_check_jacobian to see.
>>>>
>>>>     Matt
>>>>
>>>
>>> Thank you Matt. -snes_check_jacobian options seems not working (I am
>>> using PETSc 3.3p4).
>>>
>>
> That option is in petsc-dev, use -snes_type test or -snes_compare_explicit
> in petsc-3.3. If you are using MOOSE, then chances are you have not
> assembled an exact Jacobian thus these will show a difference even if
> everything is "working fine". Checking convergence with -snes_mf_operator
> -pc_type lu as you have done is a good test for whether an inexact Jacobian
> is still a good approximation. You might have to assembly an off-diagonal
> block, for example.
>
>
>>  However, now I got a clue what I need to improve. By the way, as ksp
>>> needs the Pmat as the matrix for preconditioning procedure, is there any
>>> way let ksp use the finite difference matrix provided by snes? or this is
>>> exactly what snes_fd is doing.
>>>
>>
>> That is what snes_fd is doing.
>>
>>
>>> Also, could you explain a bit more about the wired 'true resid norm'
>>> drops and increases behavior?
>>>
>>>     0 KSP unpreconditioned resid norm 7.527570931028e-02 true resid norm
>>> 7.527570931028e-02 ||r(i)||/||b|| 1.000000000000e+00
>>>     1 KSP unpreconditioned resid norm 7.217693018525e-06 true resid norm
>>> 7.217693018525e-06 ||r(i)||/||b|| 9.588342753138e-05
>>>     2 KSP unpreconditioned resid norm 1.052214184181e-07 true resid norm
>>> 1.410618177438e-02 ||r(i)||/||b|| 1.873935417365e-01
>>>     3 KSP unpreconditioned resid norm 1.023527631618e-07 true resid norm
>>> 1.410612986979e-02 ||r(i)||/||b|| 1.873928522101e-01
>>>     4 KSP unpreconditioned resid norm 1.930893544395e-08 true resid norm
>>> 1.408238386773e-02 ||r(i)||/||b|| 1.870773984964e-01
>>>
>>
>> This looks like you are losing orthogonality in the GMRES basis after
>> step 1. Maybe try *-ksp_gmres_modifiedgramschmidt*
>> You have an error in your Jacobian.
>>
>
> I would also be concerned about finite differencing error. Does the
> convergence behavior change at all if you use -mat_mffd_type ds?
>
>
>> *   *Matt
>>
>>
>>> Ling
>>>
>>>
>>>
>>>>
>>>>> Ling
>>>>>
>>>>>
>>>>> On Wed, Jan 30, 2013 at 6:40 PM, Matthew Knepley <knepley at gmail.com>wrote:
>>>>>
>>>>>> On Wed, Jan 30, 2013 at 6:30 PM, Zou (Non-US), Ling <ling.zou at inl.gov
>>>>>> > wrote:
>>>>>>
>>>>>>> Hi, All
>>>>>>>
>>>>>>> I am testing the performance of snes_mf_operator against snes_fd.
>>>>>>>
>>>>>>
>>>>>> You need to give -snes_view so we can see what solver is begin used.
>>>>>>
>>>>>>   Matt
>>>>>>
>>>>>>> I know snes_fd is for test/debugging and extremely slow, which is ok
>>>>>>> for my testing purpose. I then compared the code performance using
>>>>>>> snes_mf_operator against snes_fd. Of course, snes_mf_operator uses way less
>>>>>>> computing time then snes_fd, however, the snes_mf_operator non-linear
>>>>>>> solver performance is worse than snes_fd, in terms of non linear iteration
>>>>>>> in each time steps.
>>>>>>>
>>>>>>> Here is the PETSc Options Table entries taken from the log_summary
>>>>>>> when using snes_mf_operator
>>>>>>> #PETSc Option Table entries:
>>>>>>> -ksp_converged_reason
>>>>>>> -ksp_gmres_restart 300
>>>>>>> -ksp_monitor_true_residual
>>>>>>> -log_summary
>>>>>>> -m pipe_7eqn_2phase_step7_ps.i
>>>>>>> -mat_fd_type ds
>>>>>>> -pc_type lu
>>>>>>> -snes_mf_operator
>>>>>>> -snes_monitor
>>>>>>> #End of PETSc Option Table entries
>>>>>>>
>>>>>>> Here is the PETSc Options Table entries taken from the log_summary
>>>>>>> when using snes_fd
>>>>>>> #PETSc Option Table entries:
>>>>>>> -ksp_converged_reason
>>>>>>> -ksp_gmres_restart 300
>>>>>>> -ksp_monitor_true_residual
>>>>>>> -log_summary
>>>>>>> -m pipe_7eqn_2phase_step7_ps.i
>>>>>>> -mat_fd_type ds
>>>>>>> -pc_type lu
>>>>>>> -snes_fd
>>>>>>> -snes_monitor
>>>>>>> #End of PETSc Option Table entries
>>>>>>>
>>>>>>> The full code output along with log_summary are attached.
>>>>>>>
>>>>>>> I've noticed that when using snes_fd, the non-linear convergence is
>>>>>>> always good in each time step, around 3-4 non-linear steps with almost
>>>>>>> quadratic convergence rate. In each non-linear step, it uses only 1 linear
>>>>>>> step to converge as I used '-pc_type lu' and only 1 linear step is
>>>>>>> expected. Here is a piece of output I pulled out from the code output (very
>>>>>>> nice non-linear, linear performance but of course very expensive):
>>>>>>>
>>>>>>> DT: 1.234568e-05
>>>>>>>  Solving time step  7, time=4.34568e-05...
>>>>>>>   Initial |residual|_2 = 3.547156e+00
>>>>>>>   NL step  0, |residual|_2 = 3.547156e+00
>>>>>>>   0 SNES Function norm 3.547155872103e+00
>>>>>>>     0 KSP unpreconditioned resid norm 3.547155872103e+00 true resid
>>>>>>> norm 3.547155872103e+00 ||r(i)||/||b|| 1.000000000000e+00
>>>>>>>     1 KSP unpreconditioned resid norm 3.128472759493e-15 true resid
>>>>>>> norm 2.343197746412e-15 ||r(i)||/||b|| 6.605849392864e-16
>>>>>>>   Linear solve converged due to CONVERGED_RTOL iterations 1
>>>>>>>   NL step  1, |residual|_2 = 4.900005e-04
>>>>>>>   1 SNES Function norm 4.900004596844e-04
>>>>>>>     0 KSP unpreconditioned resid norm 4.900004596844e-04 true resid
>>>>>>> norm 4.900004596844e-04 ||r(i)||/||b|| 1.000000000000e+00
>>>>>>>     1 KSP unpreconditioned resid norm 5.026229113909e-18 true resid
>>>>>>> norm 1.400595243895e-17 ||r(i)||/||b|| 2.858354959089e-14
>>>>>>>   Linear solve converged due to CONVERGED_RTOL iterations 1
>>>>>>>   NL step  2, |residual|_2 = 1.171419e-06
>>>>>>>   2 SNES Function norm 1.171419468770e-06
>>>>>>>     0 KSP unpreconditioned resid norm 1.171419468770e-06 true resid
>>>>>>> norm 1.171419468770e-06 ||r(i)||/||b|| 1.000000000000e+00
>>>>>>>     1 KSP unpreconditioned resid norm 5.679448617332e-21 true resid
>>>>>>> norm 4.763172202015e-21 ||r(i)||/||b|| 4.066154207782e-15
>>>>>>>   Linear solve converged due to CONVERGED_RTOL iterations 1
>>>>>>>   NL step  3, |residual|_2 = 1.860041e-08
>>>>>>>   3 SNES Function norm 1.860041398803e-08
>>>>>>> Converged:1
>>>>>>>
>>>>>>> Back to the snes_mf_operator option, it behaviors differently. It
>>>>>>> generally takes more non-linear and linear steps. The 'KSP unpreconditioned
>>>>>>> resid norm' drops nicely however the 'true resid norm' seems to be a bit
>>>>>>> wired to me, drops then increases.
>>>>>>>
>>>>>>> DT: 1.524158e-05
>>>>>>>  Solving time step  9, time=7.24158e-05...
>>>>>>>   Initial |residual|_2 = 3.601003e+00
>>>>>>>   NL step  0, |residual|_2 = 3.601003e+00
>>>>>>>   0 SNES Function norm 3.601003423006e+00
>>>>>>>     0 KSP unpreconditioned resid norm 3.601003423006e+00 true resid
>>>>>>> norm 3.601003423006e+00 ||r(i)||/||b|| 1.000000000000e+00
>>>>>>>     1 KSP unpreconditioned resid norm 5.931429724028e-02 true resid
>>>>>>> norm 5.931429724028e-02 ||r(i)||/||b|| 1.647160257092e-02
>>>>>>>     2 KSP unpreconditioned resid norm 1.379343811770e-05 true resid
>>>>>>> norm 5.203950797327e+00 ||r(i)||/||b|| 1.445139086534e+00
>>>>>>>     3 KSP unpreconditioned resid norm 4.432805478482e-08 true resid
>>>>>>> norm 5.203984109211e+00 ||r(i)||/||b|| 1.445148337256e+00
>>>>>>>   Linear solve converged due to CONVERGED_RTOL iterations 3
>>>>>>>   NL step  1, |residual|_2 = 5.928815e-02
>>>>>>>   1 SNES Function norm 5.928815267199e-02
>>>>>>>     0 KSP unpreconditioned resid norm 5.928815267199e-02 true resid
>>>>>>> norm 5.928815267199e-02 ||r(i)||/||b|| 1.000000000000e+00
>>>>>>>     1 KSP unpreconditioned resid norm 3.276993782949e-06 true resid
>>>>>>> norm 3.276993782949e-06 ||r(i)||/||b|| 5.527232061148e-05
>>>>>>>     2 KSP unpreconditioned resid norm 2.082083269186e-08 true resid
>>>>>>> norm 1.551766076370e-05 ||r(i)||/||b|| 2.617329106129e-04
>>>>>>>   Linear solve converged due to CONVERGED_RTOL iterations 2
>>>>>>>   NL step  2, |residual|_2 = 3.340603e-05
>>>>>>>   2 SNES Function norm 3.340603450829e-05
>>>>>>>     0 KSP unpreconditioned resid norm 3.340603450829e-05 true resid
>>>>>>> norm 3.340603450829e-05 ||r(i)||/||b|| 1.000000000000e+00
>>>>>>>     1 KSP unpreconditioned resid norm 6.659426858789e-07 true resid
>>>>>>> norm 6.659426858789e-07 ||r(i)||/||b|| 1.993480207037e-02
>>>>>>>     2 KSP unpreconditioned resid norm 6.115119674466e-07 true resid
>>>>>>> norm 2.887921320245e-06 ||r(i)||/||b|| 8.644909109246e-02
>>>>>>>     3 KSP unpreconditioned resid norm 1.907116539439e-09 true resid
>>>>>>> norm 1.000874623281e-06 ||r(i)||/||b|| 2.996089293486e-02
>>>>>>>     4 KSP unpreconditioned resid norm 3.383211446515e-12 true resid
>>>>>>> norm 1.005586686459e-06 ||r(i)||/||b|| 3.010194718591e-02
>>>>>>>   Linear solve converged due to CONVERGED_RTOL iterations 4
>>>>>>>   NL step  3, |residual|_2 = 2.126180e-05
>>>>>>>   3 SNES Function norm 2.126179867301e-05
>>>>>>>     0 KSP unpreconditioned resid norm 2.126179867301e-05 true resid
>>>>>>> norm 2.126179867301e-05 ||r(i)||/||b|| 1.000000000000e+00
>>>>>>>     1 KSP unpreconditioned resid norm 2.724944027954e-06 true resid
>>>>>>> norm 2.724944027954e-06 ||r(i)||/||b|| 1.281615008147e-01
>>>>>>>     2 KSP unpreconditioned resid norm 7.933800605616e-10 true resid
>>>>>>> norm 2.776823963042e-06 ||r(i)||/||b|| 1.306015547295e-01
>>>>>>>     3 KSP unpreconditioned resid norm 6.130449965920e-11 true resid
>>>>>>> norm 2.777694372634e-06 ||r(i)||/||b|| 1.306424924510e-01
>>>>>>>     4 KSP unpreconditioned resid norm 2.090637685604e-13 true resid
>>>>>>> norm 2.777696567814e-06 ||r(i)||/||b|| 1.306425956963e-01
>>>>>>>   Linear solve converged due to CONVERGED_RTOL iterations 4
>>>>>>>   NL step  4, |residual|_2 = 2.863517e-06
>>>>>>>   4 SNES Function norm 2.863517221239e-06
>>>>>>>     0 KSP unpreconditioned resid norm 2.863517221239e-06 true resid
>>>>>>> norm 2.863517221239e-06 ||r(i)||/||b|| 1.000000000000e+00
>>>>>>>     1 KSP unpreconditioned resid norm 2.518692933040e-10 true resid
>>>>>>> norm 2.518692933039e-10 ||r(i)||/||b|| 8.795801590987e-05
>>>>>>>     2 KSP unpreconditioned resid norm 2.165272180327e-12 true resid
>>>>>>> norm 1.136392813468e-09 ||r(i)||/||b|| 3.968520967987e-04
>>>>>>>   Linear solve converged due to CONVERGED_RTOL iterations 2
>>>>>>>   NL step  5, |residual|_2 = 9.132390e-08
>>>>>>>   5 SNES Function norm 9.132390063388e-08
>>>>>>> Converged:1
>>>>>>>
>>>>>>>
>>>>>>> My questions:
>>>>>>> 1, Is it true? when using snes_fd, the real Jacobian matrix, say J,
>>>>>>> is explicitly constructed. when combined with -pc_type lu, the problem
>>>>>>> J (du) = -R
>>>>>>> is directly solved as (du) = J^{-1} * (-R)
>>>>>>> where J^{-1} is calculated from this explicitly constructed matrix
>>>>>>> J, using LU factorization.
>>>>>>>
>>>>>>> 2, what's the difference between snes_mf_operator and snes_fd?
>>>>>>> What I understand (might be wrong) is snes_mf_operator does not
>>>>>>> *explicitly construct* the matrix J, as it is a matrix free method. Is the
>>>>>>> finite differencing methods behind the matrix free operator
>>>>>>> in snes_mf_operator and the matrix construction in snes_fd are the same?
>>>>>>>
>>>>>>> 3, It seems that snes_mf_operator is preconditioned, while snes_fd
>>>>>>> is not. Why it says ' KSP unpreconditioned resid norm ' but I am expecting
>>>>>>> 'KSP preconditioned resid norm'. Also if it is 'unpreconditioned',
>>>>>>> should it be identical to the 'true resid norm'? Is it my fault, for
>>>>>>> example, giving a bad preconditioning matrix, makes the KSP not working
>>>>>>> well?
>>>>>>>
>>>>>>> I'd appreciate your help...there are too many (maybe bad) questions
>>>>>>> today. And please let me know if you may need more information.
>>>>>>>
>>>>>>> Best,
>>>>>>>
>>>>>>> Ling
>>>>>>>
>>>>>>
>>>>>>
>>>>>>
>>>>>> --
>>>>>> What most experimenters take for granted before they begin their
>>>>>> experiments is infinitely more interesting than any results to which their
>>>>>> experiments lead.
>>>>>> -- Norbert Wiener
>>>>>>
>>>>>
>>>>>
>>>>
>>>>
>>>> --
>>>> What most experimenters take for granted before they begin their
>>>> experiments is infinitely more interesting than any results to which their
>>>> experiments lead.
>>>> -- Norbert Wiener
>>>>
>>>
>>>
>>
>>
>> --
>> What most experimenters take for granted before they begin their
>> experiments is infinitely more interesting than any results to which their
>> experiments lead.
>> -- Norbert Wiener
>>
>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20130131/bdf15acc/attachment.html>


More information about the petsc-users mailing list