[petsc-users] Matrix-free vs finite differenced Jacobian approximation

Alexander Lindsay alexlindsay239 at gmail.com
Tue Dec 12 15:56:19 CST 2017


So I decided to look at the condition number of our matrix, running with
`-pc_type svd -pc_svd_monitor` and it was atrocious, roughly on the order
of 1e9. After doing some scaling we are down to a condition number of 1e3,
and both MF and FD operators now converge, regardless of the differencing
types chosen. I would say the problem was definitely on our end!

On Tue, Dec 12, 2017 at 2:49 PM, Matthew Knepley <knepley at gmail.com> wrote:

> On Tue, Dec 12, 2017 at 3:19 PM, Alexander Lindsay <
> alexlindsay239 at gmail.com> wrote:
>
>> I'm helping debug the finite strain models in the TensorMechanics module
>> in MOOSE, so unfortunately I don't have a nice small PetSc code I can hand
>> you guys :-(
>>
>> Hmm, interesting, if I run with `-snes_mf_operator -snes_fd
>> -mat_mffd_type ds`, I get DIVERGED_BREAKDOWN during the initial linear
>> solve.
>>
>
> So the MF operator always converges. The FD operator does not always
> converge, and factorization also can fail (DIVERGED_BREAKDOWN)
> so it seems that the FD operator is incorrect. Usually we have bugs with
> coloring, but I do not think coloring is used by -snes_fd. What happens
> if you get the coloring version by just deleting the FormJacobian pointer?
>
>    Thanks,
>
>      Matt
>
>
>> If I run with `-snes_fd -mat_fd_type ds`, then the solve converges.
>>
>> So summary:
>>
>> - J = B = finite-differenced, differencing type = wp : Solve fails due to
>> DIVERGED_LINE_SEARCH
>>
>> - J = B = finite-differenced, differencing type = ds : Solve converges in
>> 3 non-linear iterations
>>  0 Nonlinear |R| = 2.259203e-02
>>       0 Linear |R| = 2.259203e-02
>>       1 Linear |R| = 6.084393e-11
>>  1 Nonlinear |R| = 4.780691e-03
>>       0 Linear |R| = 4.780691e-03
>>       1 Linear |R| = 8.580132e-19
>>  2 Nonlinear |R| = 4.806625e-09
>>       0 Linear |R| = 4.806625e-09
>>       1 Linear |R| = 1.650725e-24
>>  3 Nonlinear |R| = 9.603678e-12
>>
>> - J = matrix-free, B = finite-differenced, mat_mffd_type = mat_fd_type =
>> wp: Solve converges in 2 non-linear iterations
>>  0 Nonlinear |R| = 2.259203e-02
>>       0 Linear |R| = 2.259203e-02
>>       1 Linear |R| = 2.258733e-02
>>       2 Linear |R| = 3.103342e-06
>>       3 Linear |R| = 6.779865e-12
>>  1 Nonlinear |R| = 7.497740e-06
>>       0 Linear |R| = 7.497740e-06
>>       1 Linear |R| = 8.265413e-12
>>  2 Nonlinear |R| = 7.993729e-12
>>
>> - J = matrix-free, B = finite-differenced, mat_mffd_type = ds,
>> mat_fd_type = wp: DIVERGED_BREAKDOWN in linear solve
>>
>> - J = matrix-free, B = finite-differenced, mat_mffd_type = wp,
>> mat_fd_type = ds: Solve converges in 2 non-linear iterations
>>  0 Nonlinear |R| = 2.259203e-02
>>       0 Linear |R| = 2.259203e-02
>>       1 Linear |R| = 4.635397e-03
>>       2 Linear |R| = 5.413676e-11
>>  1 Nonlinear |R| = 1.068626e-05
>>       0 Linear |R| = 1.068626e-05
>>       1 Linear |R| = 7.942385e-12
>>  2 Nonlinear |R| = 5.444448e-11
>>
>> - J = matrix-free, B = finite-differenced, mat_mffd_type = mat_fd_type =
>> ds: Solves converges in 3 non-linear iterations:
>>  0 Nonlinear |R| = 2.259203e-02
>>       0 Linear |R| = 2.259203e-02
>>       1 Linear |R| = 1.312921e-06
>>       2 Linear |R| = 7.714018e-09
>>  1 Nonlinear |R| = 4.780690e-03
>>       0 Linear |R| = 4.780690e-03
>>       1 Linear |R| = 7.773053e-09
>>  2 Nonlinear |R| = 1.226836e-08
>>       0 Linear |R| = 1.226836e-08
>>       1 Linear |R| = 1.546288e-14
>>  3 Nonlinear |R| = 1.295982e-10
>>
>>
>>
>>
>> On Tue, Dec 12, 2017 at 12:33 PM, Smith, Barry F. <bsmith at mcs.anl.gov>
>> wrote:
>>
>>>
>>>
>>> > On Dec 12, 2017, at 11:26 AM, Alexander Lindsay <
>>> alexlindsay239 at gmail.com> wrote:
>>> >
>>> > Ok, I'm going to go back on my original statement...the physics being
>>> run here is a sub-set of a much larger set of physics; for the current set
>>> the hand-coded Jacobian actually appears to be quite good.
>>> >
>>> > With hand-coded Jacobian, -pc_type lu, the convergence is perfect:
>>> >
>>> >  0 Nonlinear |R| = 2.259203e-02
>>> >       0 Linear |R| = 2.259203e-02
>>> >       1 Linear |R| = 1.129089e-10
>>> >  1 Nonlinear |R| = 6.295583e-11
>>> >
>>> > So yea I guess at this point I'm just curious about the different
>>> behavior between `-snes_fd` and `-snes_fd -snes_mf_operator`.
>>>
>>>   Now that you have provided the exact options you are using, yes it is
>>> very unexpected behavior. Is there any chance you can send us the code that
>>> reproduces this?
>>>
>>>    The code that does the differencing in -snes_fd is similar to the
>>> code that does the differencing for -snes_mf_operator so normally one
>>> expects similar behavior but there are a couple of options you can try. Run
>>> with -snes_mf_operator and -help  | grep mat_mffd  and this will show
>>> options to control the differencing for the matrix free. For -snes_fd you
>>> have the option -mat_fd_type  wp or ds
>>>
>>>
>>> > Does the hand-coded result change your opinion Matt that the rules for
>>> FormFunction/Jacobian might be being violated?
>>> >
>>> > I understand that a finite difference approximation of the true
>>> Jacobian is an approximation. However, in the absence of possible
>>> complications like Matt suggested where an on-the-fly calculation might
>>> stand a better chance of capturing the behavior, I would expect both
>>> `-snes_mf_operator -snes_fd` and `-snes_fd` to suffer from the same
>>> approximations, right?
>>> >
>>> > On Tue, Dec 12, 2017 at 9:43 AM, Matthew Knepley <knepley at gmail.com>
>>> wrote:
>>> > On Tue, Dec 12, 2017 at 11:30 AM, Alexander Lindsay <
>>> alexlindsay239 at gmail.com> wrote:
>>> > I'm not using any hand-coded Jacobians.
>>> >
>>> > This looks to me like the rules for FormFunction/Jacobian() are being
>>> broken. If the residual function
>>> > depends on some third variable, and it changes between calls
>>> independent of the solution U, then
>>> > the stored Jacobian could look wrong, but one done every time on the
>>> fly might converge.
>>> >
>>> >    Matt
>>> >
>>> > Case 1 options: -snes_fd -pc_type lu
>>> >
>>> > 0 Nonlinear |R| = 2.259203e-02
>>> >       0 Linear |R| = 2.259203e-02
>>> >       1 Linear |R| = 7.821248e-11
>>> >  1 Nonlinear |R| = 2.258733e-02
>>> >       0 Linear |R| = 2.258733e-02
>>> >       1 Linear |R| = 5.277296e-11
>>> >  2 Nonlinear |R| = 2.258733e-02
>>> >       0 Linear |R| = 2.258733e-02
>>> >       1 Linear |R| = 5.993971e-11
>>> > Nonlinear solve did not converge due to DIVERGED_LINE_SEARCH
>>> iterations 2
>>> >
>>> > Case 2 options: -snes_fd -snes_mf_operator -pc_type lu
>>> >
>>> >  0 Nonlinear |R| = 2.259203e-02
>>> >       0 Linear |R| = 2.259203e-02
>>> >       1 Linear |R| = 2.258733e-02
>>> >       2 Linear |R| = 3.103342e-06
>>> >       3 Linear |R| = 6.779865e-12
>>> >  1 Nonlinear |R| = 7.497740e-06
>>> >       0 Linear |R| = 7.497740e-06
>>> >       1 Linear |R| = 8.265413e-12
>>> >  2 Nonlinear |R| = 7.993729e-12
>>> > Nonlinear solve converged due to CONVERGED_FNORM_RELATIVE iterations 2
>>> >
>>> >
>>> > On Tue, Dec 12, 2017 at 9:12 AM, zakaryah . <zakaryah at gmail.com>
>>> wrote:
>>> > When you say "Jacobians are bad" and "debugging the Jacobians", do you
>>> mean that the hand-coded Jacobian is wrong?  In that case, why would you be
>>> surprised that the finite difference Jacobians, which are "correct" to
>>> approximation error, perform better?  Otherwise, what does "Jacobians are
>>> bad" mean - ill-conditioned?  Singular?  Not symmetric?  Not positive
>>> definite?
>>> >
>>> >
>>> >
>>> >
>>> > --
>>> > 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
>>> >
>>> > https://www.cse.buffalo.edu/~knepley/
>>> >
>>>
>>>
>>
>
>
> --
> 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
>
> https://www.cse.buffalo.edu/~knepley/ <http://www.caam.rice.edu/~mk51/>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20171212/f16d47ad/attachment.html>


More information about the petsc-users mailing list