[petsc-dev] Jacobian checking options

Dmitry Karpeyev karpeev at mcs.anl.gov
Mon Jan 20 12:39:26 CST 2014


Forwarding this to petsc-dev because of (*) below.

The behavior of -snes_check_jacobian, is, unfortunately, controlled by
several options with different prefixes.
This is because a mix of PETSc objects are involved.  In particular, we
have -snes_test_xxx and -mat_fd_xxx determining the choice of the
differencing parameter h in J(u)v ~ (F(u+hv) - F(u))/h.
In particular, -mat_fd_type ds|wp controls the order of magnitude of h:
(1) wp uses h = (1+||u||)*epsilon,
(2) ds (default) uses h = u_j*epsilon, where j is the column of the
Jacobian being computed (i.e., v = e_j).

Furthermore, epsilon is controlled by -snes_test_err (and no programmatic
option?).
By default epsilon is the square root of the machine epsilon. However, if h
ends up being too small (< 1e-16, which
generally should only possible with -mat_fd_type ds), then h =
0.1*sign(u_j)*epsilon is used.

(*) Maybe control of -snes_check_jacobian, -snes_test (do we need this to
be separate from -snes_check_jacobian?  maybe this can be a suboption that
causes an early exit from the solve?), and -mat_fd can be made more
consistent?

David, are you computing the Jacobian at u = 0? Then h = 0.1*epsilon, and
the differenced residual should be h^2/h =
(0.01*machine_epsilon)/(0.1*sqrt_machine_epsilon) and should evaluate to
zero, it seems to me.

Dmitry.



On Mon, Jan 20, 2014 at 10:48 AM, Andrs, David <david.andrs at inl.gov> wrote:

>  Is there a PETSc option to control the finite differencing parameter for
> -snes_check_jacobian?
>
>  The reason why I ask is that if I have a term like: u^2 (u is the
> velocity), then the FD value for this entry is something like 1e-5, but my
> exact one is correctly zero. Then the difference from PETSc is rather
> larger then 1e-8 (which PETSc suggest is the order of the difference it
> should be).
>
>  --
> David
>
>
> On Sat, Jan 18, 2014 at 10:37 AM, Dmitry Karpeyev <karpeev at mcs.anl.gov>wrote:
>
>> You might know this already:
>> Run the smallest (fewest dofs) case possible. Identify incorrect Jacobian
>> entries. Trace them back to the faulty mesh nodes/variables. Debug the
>> corresponding Jacobian code. -start_in_debugger (-start_in_debugger lldb on
>> macos) can be helpful there.
>>
>> Let me know, if you need further help.
>> Also let me know what test cases are good for testing the element-based
>> FD Jacobian code that I am writing.
>>
>> Cheers,
>> Dmitry
>> On Jan 18, 2014 11:26 AM, "Zou (Non-US), Ling" <ling.zou at inl.gov> wrote:
>>
>>>  Dmitry, this tool is VERY helpful.
>>> I showed you last time, the first several iterations were good as I
>>> could check the Jacobians for the very first iteration. It then became
>>> worse but I don't know why since I was not able to check the Jacobians.
>>>
>>>  With this new feature, I was able to find out in later steps the
>>> hand-coded Jacobians are far away from the FD Jacobians for whatever reason
>>> (bugs I bet). I am investigating my codes.
>>>
>>>  Best,
>>>
>>>  Ling
>>>
>>>
>>>
>>>
>>>
>>> On Tue, Jan 14, 2014 at 5:10 PM, Dmitry Karpeyev <karpeev at mcs.anl.gov>wrote:
>>>
>>>> Sounds good.  Let me know, if you need help.
>>>>
>>>>
>>>>  On Tue, Jan 14, 2014 at 5:02 PM, Zou (Non-US), Ling <ling.zou at inl.gov>wrote:
>>>>
>>>>>  Thank you Dmitry.
>>>>> I will at first work to update my petsc version, which may need some
>>>>> help from MOOSE team though.
>>>>> I will then update you the status of the Jacobian status.
>>>>>
>>>>>  Best,
>>>>>
>>>>>  Ling
>>>>>
>>>>>
>>>>>  On Tue, Jan 14, 2014 at 4:59 PM, Dmitry Karpeyev <karpeev at mcs.anl.gov
>>>>> > wrote:
>>>>>
>>>>>> Hi Ling,
>>>>>>
>>>>>>  You can check Jacobian at EVERY SNES solve in petsc-3.4 using
>>>>>> -snes_check_jacobian
>>>>>> and view the difference between two Jacobians using
>>>>>> -snes_check_jacobian_view.
>>>>>>
>>>>>>  Naturally, this is best done on the smallest problem possible, so
>>>>>> that
>>>>>> (a) the FD Jacobian calculation doesn't take too much time,
>>>>>> (b) matrix viewing doesn't overwhelm the output.
>>>>>>
>>>>>>  MOOSE works with petsc-3.4.
>>>>>>
>>>>>>  Cheers,
>>>>>> Dmitry.
>>>>>>
>>>>>
>>>>>
>>>>
>>>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-dev/attachments/20140120/ad2d1bcf/attachment.html>


More information about the petsc-dev mailing list