Most efficient way to evaluate a casadi::Function from C++

425 views
Skip to first unread message

Arturo Laurenzi

unread,
Aug 2, 2021, 6:25:16 AM8/2/21
to CasADi
Hi all, 
I am looking for a way to minimize copies to/from my linear algebra library of choice (in my case it's Eigen) and Casadi's DM type, and more in general to achieve best performance. I have:
 - inputs as Eigen::VectorXd (from which I can extract the double* buffer of course)
 - outputs as dense Eigen::MatrixXd 

If I understand correctly, Casadi's output will always be in a sparse format (is that true btw?). Therefore, I believe the best I can obtain is (i) to avoid any copying around as far as the input vectors are concerned by using raw buffers, and (ii) to make a single copy of all non-zeroes elements from the output DMs to "my" Eigen::MatrixXd.

To this end, I found this overload
```
/** \brief Evaluate memory-less, numerically */ 
int operator()(const double** arg, double** res, casadi_int* iw, double* w, int mem) const;
```
however, I don't understand the meaning of the last three arguments, and I'd need your assistance for that.

Is this the right path?

Thanks,
Arturo

Joris Gillis

unread,
Aug 2, 2021, 6:43:19 AM8/2/21
to CasADi
Dear Arturo,

You can always add densify() as last operation of you expression graph; casadi will then dump the dense equivalent column first into the indicated memory buffer.
sz_arg, sz_res, sz_iw and sz_w need to be queried from the Function and
 used to allocate enough space. Note that sz_arg and sz_res may well be bigger than n_in and n_out.
Simply leave the the trailing memory uninitialized.
For mem argument, use the output of checkout() and do not forget to release() afterwards.

This API is basically the same as the CasADi codegen API.

best,
  Joris

Arturo Laurenzi

unread,
Aug 5, 2021, 12:19:12 PM8/5/21
to CasADi
Thanks for the prompt reply!
Still I believe I'm missing something about arg and res. Before your reply, as in my application inputs are always dense column vectors, I was allocating arg as an array of n_in pointers to double, 
each of which was then allocated as an array of size f.size1_in(i). Same story for the results, except I have a number n_out of pointers to double, each of which points to a location f.sparsity_out(i).nnz() large
(outputs can be sparse in my application)

Roughly:
double** arg = new double*[f.n_in()]; 
for(int i = 0; i < f.n_in(); ++i) 
  arg[i] = new double[f.size1_in(i)];

double** res = new double*[f.n_out()]; 
for(int i = 0; i < f.n_out(); ++i) 
  res[i] = new double[f.sparsity_out(i).nnz()];

Now I understood that sz_arg and sz_res can be larger than what I provide. However, what should be the length of the trailing arrays, as f.size1_in(i) or f.sparsity_out(i) can't be called
for i >= n_in ?

Last thing, what is the purpose of the checkout/release mechanism? Should I call release every time I call the function (after consuming its output maybe)?

Joris Gillis

unread,
Aug 5, 2021, 2:45:16 PM8/5/21
to CasADi
what should be the length of the trailing arrays, as f.size1_in(i) or f.sparsity_out(i) can't be called
for i >= n_in ?
Nothing. No need to assign the trailing elements to anything; keep using your old stop condition.
You can do checkout once and keep using that indefinitely, and release at program termination.

 joris

--
Sent from CasADi's user forum at http://forum.casadi.org.
---
You received this message because you are subscribed to a topic in the Google Groups "CasADi" group.
To unsubscribe from this topic, visit https://groups.google.com/d/topic/casadi-users/uusT4p7zZDM/unsubscribe.
To unsubscribe from this group and all its topics, send an email to casadi-users...@googlegroups.com.
To view this discussion on the web visit https://groups.google.com/d/msgid/casadi-users/454303bd-5088-426c-888e-539bb6a91664n%40googlegroups.com.

Reply all
Reply to author
Forward
0 new messages