[petsc-users] Using -ts_type sundials with -snes-fd

Jed Brown jedbrown at mcs.anl.gov
Fri Jun 22 02:30:09 CDT 2012


Geoff, just add -pc_type lu to the options below, it will use a (sparse)
direct solver as a preconditioner. Sundials is still using matrix-free
finite differencing for Jacobian-free Newton-Krylov, and we can't control
how often the Jacobian is recomputed, but we can choose whatever we want as
a preconditioner.

The Rosenbrock methods are pretty slick in that they only need a linear
solve per stage and that the Jacobian is only recomputed once per step. If
you have a problem where Jacobian evaluation is very expensive, we can make
methods with many stages (to further amortize setup cost; good unless it
causes more step rejections).

[...]
  PC Object:   1 MPI processes
    type: lu
      LU: out-of-place factorization
      tolerance for zero pivot 2.22045e-14
      matrix ordering: nd
      factor fill ratio given 5, needed 1
        Factored matrix follows:
          Matrix Object:           1 MPI processes
            type: seqaij
            rows=3, cols=3
            package used to perform factorization: petsc
            total: nonzeros=9, allocated nonzeros=9
            total number of mallocs used during MatSetValues calls =0
              using I-node routines: found 1 nodes, limit used is 5
    linear system matrix = precond matrix:
    Matrix Object:     1 MPI processes
      type: seqaij
      rows=3, cols=3
      total: nonzeros=9, allocated nonzeros=15
      total number of mallocs used during MatSetValues calls =0
        using I-node routines: found 1 nodes, limit used is 5
  Sundials no. of preconditioner evaluations 330
  Sundials no. of preconditioner solves 3497
  Sundials no. of Jacobian-vector product evaluations 1719
  Sundials no. of rhs calls for finite diff. Jacobian-vector evals 1719
steps 1000 (0 rejected, 0 SNES fails), ftime 1.97026e+10, nonlinits 1719,
linits 1719

On Thu, Jun 21, 2012 at 11:46 AM, Geoff Oxberry <goxberry at gmail.com> wrote:

> Peter,
>
> You are correct. Here is the output (now with petsc-dev changeset:
> 23668:a876ee8e66fd):
>
> goxberry$ ./ex8 -problem_type rober -ts_type sundials -ts_view
> TS Object: 1 MPI processes
>   type: sundials
>   maximum steps=1000
>   maximum time=1e+11
>   total number of nonlinear solver iterations=3739
>   total number of nonlinear solve failures=0
>   total number of linear solver iterations=3739
>   total number of rejected steps=0
>   Sundials integrater does not use SNES!
>   Sundials integrater type BDF: backward differentiation formula
>   Sundials abs tol 1e-06 rel tol 1e-06
>   Sundials linear solver tolerance factor 0.05
>   Sundials max dimension of Krylov subspace 5
>   Sundials using unmodified (classical) Gram-Schmidt for orthogonalization
> in GMRES
>   Sundials suggested factor for tolerance scaling 1
>   Sundials cumulative number of internal steps 1000
>   Sundials no. of calls to rhs function 1931
>   Sundials no. of calls to linear solver setup function 0
>   Sundials no. of error test failures 16
>   Sundials no. of nonlinear solver iterations 1930
>   Sundials no. of nonlinear convergence failure 204
>   Sundials no. of linear iterations 3739
>   Sundials no. of linear convergence failures 0
>   PC Object:   1 MPI processes
>     type: none
>   Sundials no. of preconditioner evaluations 0
>   Sundials no. of preconditioner solves 0
>   Sundials no. of Jacobian-vector product evaluations 3739
>   Sundials no. of rhs calls for finite diff. Jacobian-vector evals 3739
> steps 1000 (0 rejected, 0 SNES fails), ftime 744.845, nonlinits 3739,
> linits 3739
>
> I use Sundials mostly with finite difference Jacobians, but I also provide
> Jacobian matrices where I can (which would be another helpful feature).
>
> Cheers,
>
> Geoff
>
> On Jun 21, 2012, at 10:52 AM, Peter Brune wrote:
>
> I asked my question because it might be happening automatically.  His
> output looked like the equation was solved.  A brief google search on this
> makes it look like they do this anyways.  In addition:
>
>
> http://www.mcs.anl.gov/petsc/petsc-dev/src/ts/impls/implicit/sundials/sundials.c.html#TSView_Sundials
>
> directly refers to "Sundials no. of rhs calls for finite diff.
> Jacobian-vector evals."
>
> - Peter
>
> On Thu, Jun 21, 2012 at 9:48 AM, Matthew Knepley <knepley at gmail.com>wrote:
>
>> On Thu, Jun 21, 2012 at 8:45 AM, Peter Brune <prbrune at gmail.com> wrote:
>>
>>> What do you see when you run with -ts_view?
>>>
>>> - Peter
>>>
>>>
>>> On Thu, Jun 21, 2012 at 9:40 AM, Geoff Oxberry <goxberry at gmail.com>wrote:
>>>
>>>> Peter,
>>>>
>>>> Just wanted to make sure there wasn't some Sundials-specific option for
>>>> finite difference Jacobians that I was missing; despite reading the manual,
>>>> it's a large package, and it's easy to miss things. If that's the case, I'd
>>>> like to make a feature request for such an option.
>>>>
>>>
>> If I understand correctly, you want a MF Jacobian with Sundials. We can't
>> do that because Sundials is completely
>> closed package, which we cannot pry apart to insert something like this.
>> The alternative is to use the stuff solvers
>> we currently have in TS. I thought that you had used the Rosenbrock-W
>> stuff. Is this sufficient?
>>
>>   Thanks,
>>
>>       Matt
>>
>>
>>>  Geoff
>>>>
>>>> On Jun 21, 2012, at 9:53 AM, Peter Brune wrote:
>>>>
>>>> Note that in the code in ts/impls/implicit/sundials it says:
>>>>
>>>> This uses its own nonlinear solver and krylov method so PETSc SNES and
>>>> KSP options do not apply...
>>>>
>>>> - Peter
>>>> On Jun 21, 2012 7:59 AM, "Geoff Oxberry" <goxberry at gmail.com> wrote:
>>>>
>>>>> Running the following example from PETSC 3.3.0-dev (changeset:
>>>>> 23631:0e86ac5e4170)
>>>>>
>>>>> /path/to/petsc-dev/src/ts/examples/tutorials/ex8 -problem_type rober
>>>>> -snes_fd -ts_type sundials
>>>>>
>>>>> gives the following output
>>>>>
>>>>> steps 1000 (0 rejected, 0 SNES fails), ftime 744.845, nonlinits 3739,
>>>>> linits 3739
>>>>> WARNING! There are options you set that were not used!
>>>>> WARNING! could be spelling mistake, etc!
>>>>> Option left: name:-snes_fd no value
>>>>>
>>>>> Just to confirm, is it currently impossible to use a finite difference
>>>>> Jacobian matrix in concert with CVODE? If so, could this feature be
>>>>> implemented in a future release? I currently rely on Sundials to integrate
>>>>> stiff systems of ODEs, and for my application, it is impractical to derive
>>>>> an analytical Jacobian matrix. (It is an issue I've discussed both with Jed
>>>>> and Matt on another forum.)
>>>>>
>>>>> Cheers,
>>>>>
>>>>> Geoff
>>>>
>>>>
>>>>
>>>
>>
>>
>> --
>> 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/20120621/ea5153ff/attachment-0001.html>


More information about the petsc-users mailing list