Skip to content

pass user jacobian to cvode #79

New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Closed
gasagna opened this issue Sep 6, 2016 · 11 comments
Closed

pass user jacobian to cvode #79

gasagna opened this issue Sep 6, 2016 · 11 comments

Comments

@gasagna
Copy link

gasagna commented Sep 6, 2016

Hi,

I am looking to implement the required code to pass the rhs jacobian to sundials to avoid expensive finite-difference approximations. According to the example cvde_Roberts_dns.jl, at line 26 (this link) this needs implementing a wrapper from Sundials._DlsMat to Matrix and the user jacobian wrapper.

My strategy for the latter is as follows:

# Type to store user rhs and jacobian, no data are provided, 
# as these can be defined using efficient closures in v0.5.
type FunJac{F, J}
    fun::F
    jac::J
end
funjac(f, J) = FunJac(f, J)

function cvodefunjac(t::Float64, 
                     x::N_Vector, 
                     ẋ::N_Vector, 
                     funjac::FunJac)
    funjac.fun(t, convert(Vector, x), convert(Vector, ẋ))
    return CV_SUCCESS
end

function cvodefunjac(N::Clong, 
                     t::realtype, 
                     x::N_Vector, 
                     ẋ::N_Vector, 
                     J::DlsMat,
                     userjac::FunJac,
                     tmp1::N_Vector, 
                     tmp2::N_Vector, 
                     tmp3::N_Vector)
    funjac.jac(t, convert(Matrix, J), convert(Vector, x))
    return CV_SUCCESS
end

i.e. the same FunJac type is used to store the rhs and the jacobian. Appropriate julia callbacks are then defined for the rhs and jacobian using cfunction and dispatch, along the lines of:

fun = cfunction(cvodefunjac, 
                Cint, 
                (realtype, 
                 N_Vector, 
                 N_Vector, 
                 Ref{typeof(funjac)}))

jac = cfunction(cvodefunjac, 
                Cint, 
                (Clong, 
                 realtype, 
                 N_Vector, 
                 N_Vector, 
                 DlsMat,
                 Ref{typeof(funjac)}, 
                 N_Vector, 
                 N_Vector, 
                 N_Vector))

ynv = NVector(copy(x₀))
CVodeInit(mem, fun, 0.0, convert(N_Vector, ynv))
CVodeSetUserData(mem, funjac)
CVDlsSetDenseJacFn(mem, jac)

Now the missing piece is how to define the convert(::Type{Matrix}, J::DlsMat) method called in the jacobian. Any pointers (no pun intended) are welcome. Is the code in cvode_Roberts_dns.jl a valid approach? It looks beyond my CJulia-fu skills...

I am actually writing this as part of a personal package but could contribute it to Sundials if requested.

@ChrisRackauckas
Copy link
Member

I think this is something that should be added to the simplified interface, which is being overhauled (see #75). We would love the contribution.

@gasagna
Copy link
Author

gasagna commented Sep 6, 2016

But I need a solution first 😅👍🏻

@stevengj
Copy link
Contributor

stevengj commented Sep 7, 2016

What exactly is the DlsMat data structure in C?

@acroy
Copy link
Contributor

acroy commented Sep 7, 2016

It is defined in include/sundials/sundials_direct.h (see here)

typedef struct _DlsMat {
  int type;
  long int M;
  long int N;
  long int ldim;
  long int mu;
  long int ml;
  long int s_mu;
  realtype *data;
  long int ldata;
  realtype **cols;
} *DlsMat;

@gasagna
Copy link
Author

gasagna commented Sep 7, 2016

Hi Steven, acroy,

In types_and_consts.jl we have

type _DlsMat
    _type::Cint
    M::Clong
    N::Clong
    ldim::Clong
    mu::Clong
    ml::Clong
    s_mu::Clong
    data::Ptr{realtype}
    ldata::Clong
    cols::Ptr{Ptr{realtype}}
end

typealias DlsMat Ptr{_DlsMat}

I think I somewhat figured it out. The convert method is

function convert(::Type{Matrix}, J::DlsMat)
    _dlsmat = unsafe_load(J)
    # own is false as memory is allocated by sundials
    pointer_to_array(_dlsmat.data, (_dlsmat.M, _dlsmat.N), false)
end

I have tests passing in a custom problem where I provide the rhs and the jacobian. I suspected I could have had troubles with the row/column-order difference between C and julia, but this seems to work.

@acroy
Copy link
Contributor

acroy commented Sep 7, 2016

While C uses row-major ordering, apparently the matrix in DlsMat.data is stored columnwise (according to the docs). So this might indeed have been the problem ...

@gasagna
Copy link
Author

gasagna commented Sep 7, 2016

Not really a problem! Much better 😏

@smldis
Copy link

smldis commented May 29, 2017

As a step forward I report here @ChrisRackauckas 's hints on solving this:

We need to grab has_jac from DiffEq and use that to define the Jacobian.
The same setup can be done to IDASpilsPrecSolveFn and IDASpilsPrecSetupFn.

@ChrisRackauckas
Copy link
Member

ChrisRackauckas commented May 29, 2017

Well, @gasagna's solution is great for the simplified interface. The common interface can just check has_jac(f) and if it does, put a closure over f(Val{:jac},t,u,du) and all of the other steps are the same. It looks like @gasagna has a full solution to this, so I'll look into this.

It looks like setting up the preconditioners has the same difficulties, so I think this solution can be applied there as well (for KINSOL, CVODE, and IDA). I'll probably take a deeper look into this in a week or so since this is probably the biggest thing the interface is missing for Sundials.

@ChrisRackauckas
Copy link
Member

The common interface now passes the Jacobian to CVODE. IDA is next and I'll close the issue.

@ChrisRackauckas
Copy link
Member

Thanks for your help @gasagna. I was able to use this and now the Jacobians pass to CVODE and IDA on the common interface bindings.

Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment
Projects
None yet
Development

No branches or pull requests

5 participants