[petsc-users] Questions about TS

Jed Brown jedbrown at mcs.anl.gov
Wed Oct 19 11:14:10 CDT 2011


On Wed, Oct 19, 2011 at 10:54, Blaise Bourdin <bourdin at math.lsu.edu> wrote:

> Hi,
>
> I am trying to use TS to solve a simple transient problem in an
> unstructured finite element f90 code.
>
> 1. Section 6.1.1 of the manual refers to a TSSetMatrices function that can
> be used to set the RHS and LHS matrices, but I can;t find it. Is this
> section outdated?
>

Yes, I must have missed this section when updating the documentation.


> 2. Since we are using unstructured finite elements, the LHS matrix is not
> the identity. As far as I understand, we have two possible choices:
>   - Use a mass lumping approximation of the variational identity matrix
> (mass matrix), M, and use M^{-1}K for the RHS matrix instead of K.
>

You can also write this as a special case of the choice below where you use
-ksp_type preonly -pc_type jacobi.


>   - Use an IMEX method where the implicit matrix is the variational
> identity M.
>  Is this right? What is the recommended way?
>

I would just do this because it's the most flexible. See the user's manual
section on IMEX methods. If you are interested in adaptive error control,
then you should also check out -ts_type rosw in petsc-dev. In any case, you
can write your mass matrix as well as any stiff terms that you want to treat
implicitly into TSSetIFunction(), provide an (approximate) Jacobian with
TSSetIJacobian(), and put the rest in TSSetRHSFunction(). You can be even
more sloppy about it with TSROSW.

Look at src/ts/examples/tutorials/ex22.c (has a Fortran twin, ex22f.F) or
ex25.c in petsc-dev.
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20111019/c6d14d5c/attachment.htm>


More information about the petsc-users mailing list