1 /*
2  *    This file is part of CasADi.
3  *
4  *    CasADi -- A symbolic framework for dynamic optimization.
5  *    Copyright (C) 2010-2014 Joel Andersson, Joris Gillis, Moritz Diehl,
6  *                            K.U. Leuven. All rights reserved.
7  *    Copyright (C) 2011-2014 Greg Horn
8  *
9  *    CasADi is free software; you can redistribute it and/or
10  *    modify it under the terms of the GNU Lesser General Public
11  *    License as published by the Free Software Foundation; either
12  *    version 3 of the License, or (at your option) any later version.
13  *
14  *    CasADi is distributed in the hope that it will be useful,
15  *    but WITHOUT ANY WARRANTY; without even the implied warranty of
16  *    MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  *    Lesser General Public License for more details.
18  *
19  *    You should have received a copy of the GNU Lesser General Public
20  *    License along with CasADi; if not, write to the Free Software
21  *    Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA
22  *
23  */
24 
25 
26 #include "rootfinder_impl.hpp"
27 #include "mx_node.hpp"
28 #include <iterator>
29 #include "linsol.hpp"
30 
31 #include "global_options.hpp"
32 
33 using namespace std;
34 namespace casadi {
35 
rootfinder_in()36   vector<string> rootfinder_in() {
37     vector<string> ret(rootfinder_n_in());
38     for (size_t i=0; i<ret.size(); ++i) ret[i]=rootfinder_in(i);
39     return ret;
40   }
41 
rootfinder_out()42   vector<string> rootfinder_out() {
43     vector<string> ret(rootfinder_n_out());
44     for (size_t i=0; i<ret.size(); ++i) ret[i]=rootfinder_out(i);
45     return ret;
46   }
47 
rootfinder_in(casadi_int ind)48   string rootfinder_in(casadi_int ind) {
49     switch (static_cast<RootfinderInput>(ind)) {
50     case ROOTFINDER_X0:  return "x0";
51     case ROOTFINDER_P:   return "p";
52     case ROOTFINDER_NUM_IN: break;
53     }
54     return string();
55   }
56 
rootfinder_out(casadi_int ind)57   string rootfinder_out(casadi_int ind) {
58     switch (static_cast<RootfinderOutput>(ind)) {
59     case ROOTFINDER_X:  return "x";
60     case ROOTFINDER_NUM_OUT: break;
61     }
62     return string();
63   }
64 
rootfinder_n_in()65   casadi_int rootfinder_n_in() {
66     return ROOTFINDER_NUM_IN;
67   }
68 
rootfinder_n_out()69   casadi_int rootfinder_n_out() {
70     return ROOTFINDER_NUM_OUT;
71   }
72 
rootfinder_options(const std::string & name)73   std::vector<std::string> rootfinder_options(const std::string& name) {
74     return Rootfinder::plugin_options(name).all();
75   }
76 
rootfinder_option_type(const std::string & name,const std::string & op)77   std::string rootfinder_option_type(const std::string& name, const std::string& op) {
78     return Rootfinder::plugin_options(name).type(op);
79   }
80 
rootfinder_option_info(const std::string & name,const std::string & op)81   std::string rootfinder_option_info(const std::string& name, const std::string& op) {
82     return Rootfinder::plugin_options(name).info(op);
83   }
84 
has_rootfinder(const string & name)85   bool has_rootfinder(const string& name) {
86     return Rootfinder::has_plugin(name);
87   }
88 
load_rootfinder(const string & name)89   void load_rootfinder(const string& name) {
90     Rootfinder::load_plugin(name);
91   }
92 
doc_rootfinder(const string & name)93   string doc_rootfinder(const string& name) {
94     return Rootfinder::getPlugin(name).doc;
95   }
96 
rootfinder(const string & name,const string & solver,const SXDict & rfp,const Dict & opts)97   Function rootfinder(const string& name, const string& solver,
98                       const SXDict& rfp, const Dict& opts) {
99     return rootfinder(name, solver, Rootfinder::create_oracle(rfp, opts), opts);
100   }
101 
rootfinder(const string & name,const string & solver,const MXDict & rfp,const Dict & opts)102   Function rootfinder(const string& name, const string& solver,
103                       const MXDict& rfp, const Dict& opts) {
104     return rootfinder(name, solver, Rootfinder::create_oracle(rfp, opts), opts);
105   }
106 
107   template<typename XType>
create_oracle(const std::map<std::string,XType> & d,const Dict & opts)108   Function Rootfinder::create_oracle(const std::map<std::string, XType>& d,
109                                  const Dict& opts) {
110     std::vector<XType> rfp_in(RFP_NUM_IN), rfp_out(RFP_NUM_OUT);
111     for (auto&& i : d) {
112       if (i.first=="x") {
113         rfp_in[RFP_X]=i.second;
114       } else if (i.first=="p") {
115         rfp_in[RFP_P]=i.second;
116       } else if (i.first=="g") {
117         rfp_out[RFP_G]=i.second;
118       } else {
119         casadi_error("No such field: " + i.first);
120       }
121     }
122 
123     // Options for the oracle
124     Dict oracle_options;
125     Dict::const_iterator it = opts.find("oracle_options");
126     if (it!=opts.end()) {
127       // "oracle_options" has been set
128       oracle_options = it->second;
129     } else if ((it=opts.find("verbose")) != opts.end()) {
130       // "oracle_options" has not been set, but "verbose" has
131       oracle_options["verbose"] = it->second;
132     }
133 
134     // Create oracle
135     return Function("rfp", rfp_in, rfp_out, {"x0", "p"}, {"x"}, oracle_options);
136   }
137 
rootfinder(const std::string & name,const std::string & solver,const Function & f,const Dict & opts)138   Function rootfinder(const std::string& name, const std::string& solver,
139                    const Function& f, const Dict& opts) {
140     // Make sure that residual function is sound
141     if (f.has_free()) {
142       casadi_error("Cannot create '" + name + "' since " + str(f.get_free()) + " are free.");
143     }
144     return Function::create(Rootfinder::instantiate(name, solver, f), opts);
145   }
146 
Rootfinder(const std::string & name,const Function & oracle)147   Rootfinder::Rootfinder(const std::string& name, const Function& oracle)
148     : OracleFunction(name, oracle) {
149 
150     // Default options
151     iin_ = 0;
152     iout_ = 0;
153     error_on_fail_ = true;
154   }
155 
~Rootfinder()156   Rootfinder::~Rootfinder() {
157   }
158 
159   const Options Rootfinder::options_
160   = {{&OracleFunction::options_},
161      {{"linear_solver",
162        {OT_STRING,
163         "User-defined linear solver class. Needed for sensitivities."}},
164       {"linear_solver_options",
165        {OT_DICT,
166         "Options to be passed to the linear solver."}},
167       {"constraints",
168        {OT_INTVECTOR,
169         "Constrain the unknowns. 0 (default): no constraint on ui, "
170         "1: ui >= 0.0, -1: ui <= 0.0, 2: ui > 0.0, -2: ui < 0.0."}},
171       {"implicit_input",
172        {OT_INT,
173         "Index of the input that corresponds to the actual root-finding"}},
174       {"implicit_output",
175        {OT_INT,
176         "Index of the output that corresponds to the actual root-finding"}},
177       {"jacobian_function",
178        {OT_FUNCTION,
179         "Function object for calculating the Jacobian (autogenerated by default)"}},
180       {"error_on_fail",
181        {OT_BOOL,
182         "When the numerical process returns unsuccessfully, raise an error (default false)."}}
183      }
184   };
185 
init(const Dict & opts)186   void Rootfinder::init(const Dict& opts) {
187 
188     // Default (temporary) options
189     Dict linear_solver_options;
190     string linear_solver = "qr";
191     Function jac; // Jacobian of f with respect to z
192 
193     // Read options
194     for (auto&& op : opts) {
195       if (op.first=="implicit_input") {
196         iin_ = op.second;
197       } else if (op.first=="implicit_output") {
198         iout_ = op.second;
199       } else if (op.first=="jacobian_function") {
200         jac = op.second;
201       } else if (op.first=="linear_solver_options") {
202         linear_solver_options = op.second;
203       } else if (op.first=="linear_solver") {
204         linear_solver = op.second.to_string();
205       } else if (op.first=="constraints") {
206         u_c_ = op.second;
207       } else if (op.first=="error_on_fail") {
208         error_on_fail_ = op.second;
209       }
210     }
211 
212     // Get the number of equations and check consistency
213     casadi_assert(iin_>=0 && iin_<oracle_.n_in() && oracle_.n_in()>0,
214                           "Implicit input not in range");
215     casadi_assert(iout_>=0 && iout_<oracle_.n_out() && oracle_.n_out()>0,
216                           "Implicit output not in range");
217     casadi_assert(oracle_.sparsity_out(iout_).is_dense()
218                           && oracle_.sparsity_out(iout_).is_column(),
219                           "Residual must be a dense vector");
220     casadi_assert(oracle_.sparsity_in(iin_).is_dense()
221                           && oracle_.sparsity_in(iin_).is_column(),
222                           "Unknown must be a dense vector");
223     n_ = oracle_.nnz_out(iout_);
224     casadi_assert(n_ == oracle_.nnz_in(iin_),
225       "Dimension mismatch. Input size is " + str(oracle_.nnz_in(iin_)) + ", "
226       "while output size is " + str(oracle_.nnz_out(iout_)));
227 
228     // Call the base class initializer
229     OracleFunction::init(opts);
230 
231     // Generate Jacobian if not provided
232     if (jac.is_null()) jac = oracle_.jacobian_old(iin_, iout_);
233     set_function(jac, "jac_f_z");
234     sp_jac_ = jac.sparsity_out(0);
235 
236     // Check for structural singularity in the Jacobian
237     casadi_assert(!sp_jac_.is_singular(),
238       "Rootfinder::init: singularity - the jacobian is structurally rank-deficient. "
239       "sprank(J)=" + str(sprank(sp_jac_)) + " (instead of " + str(sp_jac_.size1()) + ")");
240 
241     // Get the linear solver creator function
242     linsol_ = Linsol("linsol", linear_solver, sp_jac_, linear_solver_options);
243 
244     // Constraints
245     casadi_assert(u_c_.size()==n_ || u_c_.empty(),
246       "Constraint vector if supplied, must be of length n, but got "
247       + str(u_c_.size()) + " and n = " + str(n_));
248 
249     // Allocate sufficiently large work vectors
250     alloc(oracle_);
251     size_t sz_w = oracle_.sz_w();
252     if (!jac.is_null()) {
253       sz_w = max(sz_w, jac.sz_w());
254     }
255     alloc_w(sz_w + 2*static_cast<size_t>(n_));
256   }
257 
init_mem(void * mem) const258   int Rootfinder::init_mem(void* mem) const {
259     if (OracleFunction::init_mem(mem)) return 1;
260 
261     auto m = static_cast<RootfinderMemory*>(mem);
262 
263     // Problem has not been solved at this point
264     m->success = false;
265     m->unified_return_status = SOLVER_RET_UNKNOWN;
266 
267     return 0;
268   }
269 
eval(const double ** arg,double ** res,casadi_int * iw,double * w,void * mem) const270   int Rootfinder::eval(const double** arg, double** res,
271       casadi_int* iw, double* w, void* mem) const {
272     // Reset the solver, prepare for solution
273     setup(mem, arg, res, iw, w);
274 
275     // Solve the NLP
276     int ret = solve(mem);
277     auto m = static_cast<RootfinderMemory*>(mem);
278     if (error_on_fail_ && !m->success)
279       casadi_error("rootfinder process failed. "
280                    "Set 'error_on_fail' option to false to ignore this error.");
281 
282     return ret;
283   }
284 
set_work(void * mem,const double ** & arg,double ** & res,casadi_int * & iw,double * & w) const285   void Rootfinder::set_work(void* mem, const double**& arg, double**& res,
286                         casadi_int*& iw, double*& w) const {
287     auto m = static_cast<RootfinderMemory*>(mem);
288 
289     // Problem has not been solved at this point
290     m->success = false;
291     m->unified_return_status = SOLVER_RET_UNKNOWN;
292 
293     // Get input pointers
294     m->iarg = arg;
295     arg += n_in_;
296 
297     // Get output pointers
298     m->ires = res;
299     res += n_out_;
300   }
301 
302   Function Rootfinder
get_forward(casadi_int nfwd,const std::string & name,const std::vector<std::string> & inames,const std::vector<std::string> & onames,const Dict & opts) const303   ::get_forward(casadi_int nfwd, const std::string& name,
304                 const std::vector<std::string>& inames,
305                 const std::vector<std::string>& onames,
306                 const Dict& opts) const {
307     // Symbolic expression for the input
308     vector<MX> arg = mx_in(), res = mx_out();
309     vector<vector<MX>> fseed = fwd_seed<MX>(nfwd), fsens;
310     arg[iin_] = MX::sym(arg[iin_].name(), Sparsity(arg[iin_].size()));
311     for (auto&& e : fseed) e[iin_] = MX::sym(e[iin_].name(), e[iin_].size());
312     ad_forward(arg, res, fseed, fsens, false, false);
313 
314     // Construct return function
315     arg.insert(arg.end(), res.begin(), res.end());
316     vector<MX> v(nfwd);
317     for (casadi_int i=0; i<n_in_; ++i) {
318       for (casadi_int d=0; d<nfwd; ++d) v[d] = fseed[d][i];
319       arg.push_back(horzcat(v));
320     }
321     res.clear();
322     for (casadi_int i=0; i<n_out_; ++i) {
323       for (casadi_int d=0; d<nfwd; ++d) v[d] = fsens[d][i];
324       res.push_back(horzcat(v));
325     }
326     return Function(name, arg, res, inames, onames, opts);
327   }
328 
329   Function Rootfinder
get_reverse(casadi_int nadj,const std::string & name,const std::vector<std::string> & inames,const std::vector<std::string> & onames,const Dict & opts) const330   ::get_reverse(casadi_int nadj, const std::string& name,
331                 const std::vector<std::string>& inames,
332                 const std::vector<std::string>& onames,
333                 const Dict& opts) const {
334     // Symbolic expression for the input
335     vector<MX> arg = mx_in();
336     arg[iin_] = MX::sym(arg[iin_].name() + "_guess",
337                         Sparsity(arg[iin_].size()));
338     vector<MX> res = mx_out();
339     vector<vector<MX> > aseed = symbolicAdjSeed(nadj, res), asens;
340     ad_reverse(arg, res, aseed, asens, false, false);
341 
342     // Construct return function
343     arg.insert(arg.end(), res.begin(), res.end());
344     vector<MX> v(nadj);
345     for (casadi_int i=0; i<n_out_; ++i) {
346       for (casadi_int d=0; d<nadj; ++d) v[d] = aseed[d][i];
347       arg.push_back(horzcat(v));
348     }
349     res.clear();
350     for (casadi_int i=0; i<n_in_; ++i) {
351       for (casadi_int d=0; d<nadj; ++d) v[d] = asens[d][i];
352       res.push_back(horzcat(v));
353     }
354     return Function(name, arg, res, inames, onames, opts);
355   }
356 
357   int Rootfinder::
sp_forward(const bvec_t ** arg,bvec_t ** res,casadi_int * iw,bvec_t * w,void * mem) const358   sp_forward(const bvec_t** arg, bvec_t** res, casadi_int* iw, bvec_t* w, void* mem) const {
359     bvec_t* tmp1 = w; w += n_;
360     bvec_t* tmp2 = w; w += n_;
361 
362     // Propagate dependencies through the function
363     const bvec_t** arg1 = arg+n_in_;
364     copy(arg, arg+n_in_, arg1);
365     arg1[iin_] = nullptr;
366     bvec_t** res1 = res+n_out_;
367     fill_n(res1, n_out_, static_cast<bvec_t*>(nullptr));
368     res1[iout_] = tmp1;
369     oracle_(arg1, res1, iw, w, 0);
370 
371     // "Solve" in order to propagate to z
372     fill_n(tmp2, n_, 0);
373     sp_jac_.spsolve(tmp2, tmp1, false);
374     if (res[iout_]) copy(tmp2, tmp2+n_, res[iout_]);
375 
376     // Propagate to auxiliary outputs
377     if (n_out_>1) {
378       arg1[iin_] = tmp2;
379       copy(res, res+n_out_, res1);
380       res1[iout_] = nullptr;
381       oracle_(arg1, res1, iw, w, 0);
382     }
383     return 0;
384   }
385 
sp_reverse(bvec_t ** arg,bvec_t ** res,casadi_int * iw,bvec_t * w,void * mem) const386   int Rootfinder::sp_reverse(bvec_t** arg, bvec_t** res,
387       casadi_int* iw, bvec_t* w, void* mem) const {
388     bvec_t* tmp1 = w; w += n_;
389     bvec_t* tmp2 = w; w += n_;
390 
391     // Get & clear seed corresponding to implicitly defined variable
392     if (res[iout_]) {
393       copy(res[iout_], res[iout_]+n_, tmp1);
394       fill_n(res[iout_], n_, 0);
395     } else {
396       fill_n(tmp1, n_, 0);
397     }
398 
399     // Propagate dependencies from auxiliary outputs to z
400     bvec_t** res1 = res+n_out_;
401     copy(res, res+n_out_, res1);
402     res1[iout_] = nullptr;
403     bvec_t** arg1 = arg+n_in_;
404     copy(arg, arg+n_in_, arg1);
405     arg1[iin_] = tmp1;
406     if (n_out_>1) {
407       if (oracle_.rev(arg1, res1, iw, w, 0)) return 1;
408     }
409 
410     // "Solve" in order to get seed
411     fill_n(tmp2, n_, 0);
412     sp_jac_.spsolve(tmp2, tmp1, true);
413 
414     // Propagate dependencies through the function
415     for (casadi_int i=0; i<n_out_; ++i) res1[i] = nullptr;
416     res1[iout_] = tmp2;
417     arg1[iin_] = nullptr; // just a guess
418     if (oracle_.rev(arg1, res1, iw, w, 0)) return 1;
419     return 0;
420   }
421 
422   std::map<std::string, Rootfinder::Plugin> Rootfinder::solvers_;
423 
424   const std::string Rootfinder::infix_ = "rootfinder";
425 
426   void Rootfinder::
ad_forward(const std::vector<MX> & arg,const std::vector<MX> & res,const std::vector<std::vector<MX>> & fseed,std::vector<std::vector<MX>> & fsens,bool always_inline,bool never_inline) const427   ad_forward(const std::vector<MX>& arg, const std::vector<MX>& res,
428           const std::vector<std::vector<MX> >& fseed,
429           std::vector<std::vector<MX> >& fsens,
430           bool always_inline, bool never_inline) const {
431     // Number of directional derivatives
432     casadi_int nfwd = fseed.size();
433     fsens.resize(nfwd);
434 
435     // Quick return if no seeds
436     if (nfwd==0) return;
437 
438     // Propagate through f_
439     vector<MX> f_arg(arg);
440     f_arg.at(iin_) = res.at(iout_);
441     vector<MX> f_res(res);
442     f_res.at(iout_) = MX(size_in(iin_)); // zero residual
443     std::vector<std::vector<MX> > f_fseed(fseed);
444     for (casadi_int d=0; d<nfwd; ++d) {
445       f_fseed[d].at(iin_) = MX(size_in(iin_)); // ignore seeds for guess
446     }
447     oracle_->call_forward(f_arg, f_res, f_fseed, fsens,
448                           always_inline, never_inline);
449 
450     // Get expression of Jacobian
451     Function jac = get_function("jac_f_z");
452     MX J = jac(f_arg).front();
453 
454     // Solve for all the forward derivatives at once
455     vector<MX> rhs(nfwd);
456     for (casadi_int d=0; d<nfwd; ++d) rhs[d] = vec(fsens[d][iout_]);
457     rhs = horzsplit(J->get_solve(-horzcat(rhs), false, linsol_));
458     for (casadi_int d=0; d<nfwd; ++d) fsens[d][iout_] = reshape(rhs[d], size_in(iin_));
459 
460     // Propagate to auxiliary outputs
461     if (n_out_>1) {
462       for (casadi_int d=0; d<nfwd; ++d) f_fseed[d][iin_] = fsens[d][iout_];
463       oracle_->call_forward(f_arg, f_res, f_fseed, fsens,
464                             always_inline, never_inline);
465       for (casadi_int d=0; d<nfwd; ++d) fsens[d][iout_] = f_fseed[d][iin_]; // Otherwise overwritten
466     }
467   }
468 
469   void Rootfinder::
ad_reverse(const std::vector<MX> & arg,const std::vector<MX> & res,const std::vector<std::vector<MX>> & aseed,std::vector<std::vector<MX>> & asens,bool always_inline,bool never_inline) const470   ad_reverse(const std::vector<MX>& arg, const std::vector<MX>& res,
471           const std::vector<std::vector<MX> >& aseed,
472           std::vector<std::vector<MX> >& asens,
473           bool always_inline, bool never_inline) const {
474 
475     // Number of directional derivatives
476     casadi_int nadj = aseed.size();
477     asens.resize(nadj);
478 
479     // Quick return if no seeds
480     if (nadj==0) return;
481 
482     // Get expression of Jacobian
483     vector<MX> f_arg(arg);
484     f_arg[iin_] = res.at(iout_);
485     Function jac = get_function("jac_f_z");
486     MX J = jac(f_arg).front();
487 
488     // Get adjoint seeds for calling f
489     vector<MX> f_res(res);
490     f_res[iout_] = MX(size_in(iin_)); // zero residual
491     vector<vector<MX> > f_aseed(nadj);
492     for (casadi_int d=0; d<nadj; ++d) {
493       f_aseed[d].resize(n_out_);
494       for (casadi_int i=0; i<n_out_; ++i) f_aseed[d][i] = i==iout_ ? f_res[iout_] : aseed[d][i];
495     }
496 
497     // Propagate dependencies from auxiliary outputs
498     vector<MX> rhs(nadj);
499     vector<vector<MX> > asens_aux;
500     if (n_out_>1) {
501       oracle_->call_reverse(f_arg, f_res, f_aseed, asens_aux, always_inline, never_inline);
502       for (casadi_int d=0; d<nadj; ++d) rhs[d] = vec(asens_aux[d][iin_] + aseed[d][iout_]);
503     } else {
504       for (casadi_int d=0; d<nadj; ++d) rhs[d] = vec(aseed[d][iout_]);
505     }
506 
507     // Solve for all the adjoint seeds at once
508     rhs = horzsplit(J->get_solve(-horzcat(rhs), true, linsol_));
509     for (casadi_int d=0; d<nadj; ++d) {
510       for (casadi_int i=0; i<n_out_; ++i) {
511         if (i==iout_) {
512           f_aseed[d][i] = reshape(rhs[d], size_out(i));
513         } else {
514           // Avoid counting the auxiliary seeds twice
515           f_aseed[d][i] = MX(size_out(i));
516         }
517       }
518     }
519 
520     // No dependency on guess (1)
521     vector<MX> tmp(nadj);
522     for (casadi_int d=0; d<nadj; ++d) {
523       asens[d].resize(n_in_);
524       tmp[d] = asens[d][iin_].is_empty(true) ? MX(size_in(iin_)) : asens[d][iin_];
525     }
526 
527     // Propagate through f_
528     oracle_->call_reverse(f_arg, f_res, f_aseed, asens, always_inline, never_inline);
529 
530     // No dependency on guess (2)
531     for (casadi_int d=0; d<nadj; ++d) {
532       asens[d][iin_] = tmp[d];
533     }
534 
535     // Add contribution from auxiliary outputs
536     if (n_out_>1) {
537       for (casadi_int d=0; d<nadj; ++d) {
538         for (casadi_int i=0; i<n_in_; ++i) if (i!=iin_) asens[d][i] += asens_aux[d][i];
539       }
540     }
541   }
542 
get_stats(void * mem) const543   Dict Rootfinder::get_stats(void* mem) const {
544     Dict stats = OracleFunction::get_stats(mem);
545     auto m = static_cast<RootfinderMemory*>(mem);
546     stats["success"] = m->success;
547     stats["unified_return_status"] = string_from_UnifiedReturnStatus(m->unified_return_status);
548     return stats;
549   }
550 
serialize_body(SerializingStream & s) const551   void Rootfinder::serialize_body(SerializingStream &s) const {
552     OracleFunction::serialize_body(s);
553 
554     s.version("Rootfinder", 1);
555     s.pack("Rootfinder::n", n_);
556     s.pack("Rootfinder::linsol", linsol_);
557     s.pack("Rootfinder::sp_jac", sp_jac_);
558     s.pack("Rootfinder::u_c", u_c_);
559     s.pack("Rootfinder::iin", iin_);
560     s.pack("Rootfinder::iout", iout_);
561     s.pack("Rootfinder::error_on_fail", error_on_fail_);
562   }
563 
serialize_type(SerializingStream & s) const564   void Rootfinder::serialize_type(SerializingStream &s) const {
565     OracleFunction::serialize_type(s);
566     PluginInterface<Rootfinder>::serialize_type(s);
567   }
568 
deserialize(DeserializingStream & s)569   ProtoFunction* Rootfinder::deserialize(DeserializingStream& s) {
570     return PluginInterface<Rootfinder>::deserialize(s);
571   }
572 
Rootfinder(DeserializingStream & s)573   Rootfinder::Rootfinder(DeserializingStream & s) : OracleFunction(s) {
574     s.version("Rootfinder", 1);
575     s.unpack("Rootfinder::n", n_);
576     s.unpack("Rootfinder::linsol", linsol_);
577     s.unpack("Rootfinder::sp_jac", sp_jac_);
578     s.unpack("Rootfinder::u_c", u_c_);
579     s.unpack("Rootfinder::iin", iin_);
580     s.unpack("Rootfinder::iout", iout_);
581     s.unpack("Rootfinder::error_on_fail", error_on_fail_);
582   }
583 
584 } // namespace casadi
585