linsol_internal.cpp
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  * KU 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 "linsol_internal.hpp"
27 
28 namespace casadi {
29 
30  LinsolInternal::LinsolInternal(const std::string& name, const Sparsity& sp)
31  : ProtoFunction(name), sp_(sp) {
32  }
33 
35  }
36 
37  void LinsolInternal::init(const Dict& opts) {
38  // Call the base class initializer
39  ProtoFunction::init(opts);
40 
41  }
42 
43  void LinsolInternal::disp(std::ostream &stream, bool more) const {
44  stream << "Linear solver " << class_name();
45  if (more) {
46  stream << std::endl;
47  disp_more(stream);
48  }
49  }
50 
51  int LinsolInternal::init_mem(void* mem) const {
52  if (!mem) return 1;
53  if (ProtoFunction::init_mem(mem)) return 1;
54  auto m = static_cast<LinsolMemory*>(mem);
55  if (record_time_) {
56  m->add_stat("nfact");
57  m->add_stat("sfact");
58  m->add_stat("solve");
59  }
60  return 0;
61  }
62 
63  void LinsolInternal::linsol_eval_sx(const SXElem** arg, SXElem** res, casadi_int* iw, SXElem* w,
64  void* mem, bool tr, casadi_int nrhs) const {
65  casadi_error("eval_sx not defined for " + class_name());
66  }
67 
68  int LinsolInternal::solve(void* mem, const double* A, double* x, casadi_int nrhs, bool tr) const {
69  casadi_error("'solve' not defined for " + class_name());
70  }
71 
72 #if 0
73  casadi_int LinsolInternal::factorize(void* mem, const double* A) const {
74  // Symbolic factorization, if needed
75  if (needs_sfact(mem, A)) {
76  if (sfact(mem, A)) return 1;
77  }
78 
79  // Numeric factorization, if needed
80  if (needs_nfact(mem, A)) {
81  if (nfact(mem, A)) return 1;
82  }
83 
84  return 0;
85  }
86 
87  bool LinsolInternal::needs_sfact(void* mem, const double* A) const {
88  auto m = static_cast<LinsolMemory*>(mem);
89  return m->is_sfact;
90  }
91 
92  bool LinsolInternal::needs_nfact(void* mem, const double* A) const {
93  auto m = static_cast<LinsolMemory*>(mem);
94  return m->is_nfact;
95  }
96 #endif
97 
98  int LinsolInternal::nfact(void* mem, const double* A) const {
99  casadi_error("'nfact' not defined for " + class_name());
100  }
101 
102  casadi_int LinsolInternal::neig(void* mem, const double* A) const {
103  casadi_error("'neig' not defined for " + class_name());
104  }
105 
106  casadi_int LinsolInternal::rank(void* mem, const double* A) const {
107  casadi_error("'rank' not defined for " + class_name());
108  }
109 
110  void LinsolInternal::generate(CodeGenerator& g, const std::string& A, const std::string& x,
111  casadi_int nrhs, bool tr) const {
112  g << "#error " << class_name() << " does not support code generation\n";
113  }
114 
115  std::map<std::string, LinsolInternal::Plugin> LinsolInternal::solvers_;
116 
117 #ifdef CASADI_WITH_THREADSAFE_SYMBOLICS
118  std::mutex LinsolInternal::mutex_solvers_;
119 #endif // CASADI_WITH_THREADSAFE_SYMBOLICS
120 
121  const std::string LinsolInternal::infix_ = "linsol";
122 
123 
127  }
128 
131  s.pack("LinsolInternal::sp", sp_);
132  }
133 
135  s.unpack("LinsolInternal::sp", sp_);
136  }
137 
140  }
141 
142 } // namespace casadi
Helper class for C code generation.
Helper class for Serialization.
void unpack(Sparsity &e)
Reconstruct an object from the input stream.
virtual int solve(void *mem, const double *A, double *x, casadi_int nrhs, bool tr) const
void init(const Dict &opts) override
Initialize.
void disp(std::ostream &stream, bool more) const override
Display object.
virtual void linsol_eval_sx(const SXElem **arg, SXElem **res, casadi_int *iw, SXElem *w, void *mem, bool tr, casadi_int nrhs) const
Evaluate SX, possibly transposed.
static std::map< std::string, Plugin > solvers_
Collection of solvers.
virtual int nfact(void *mem, const double *A) const
Numeric factorization.
void serialize_body(SerializingStream &s) const override
Serialize an object without type information.
virtual int sfact(void *mem, const double *A) const
void serialize_type(SerializingStream &s) const override
Serialize type information.
virtual casadi_int rank(void *mem, const double *A) const
Matrix rank.
static ProtoFunction * deserialize(DeserializingStream &s)
Deserialize with type disambiguation.
virtual void disp_more(std::ostream &stream) const
Print more.
LinsolInternal(const std::string &name, const Sparsity &sp)
Constructor.
static const std::string infix_
Infix.
virtual void generate(CodeGenerator &g, const std::string &A, const std::string &x, casadi_int nrhs, bool tr) const
Generate C code.
~LinsolInternal() override
Destructor.
int init_mem(void *mem) const override
Initalize memory block.
virtual casadi_int neig(void *mem, const double *A) const
Number of negative eigenvalues.
void serialize_type(SerializingStream &s) const
Serialize type information.
static ProtoFunction * deserialize(DeserializingStream &s)
Deserialize with type disambiguation.
Base class for FunctionInternal and LinsolInternal.
virtual int init_mem(void *mem) const
Initalize memory block.
virtual void serialize_type(SerializingStream &s) const
Serialize type information.
virtual void serialize_body(SerializingStream &s) const
Serialize an object without type information.
virtual void init(const Dict &opts)
Initialize.
The basic scalar symbolic class of CasADi.
Definition: sx_elem.hpp:75
Helper class for Serialization.
void pack(const Sparsity &e)
Serializes an object to the output stream.
virtual std::string class_name() const =0
Readable name of the internal class.
General sparsity class.
Definition: sparsity.hpp:106
The casadi namespace.
Definition: archiver.cpp:28
GenericType::Dict Dict
C++ equivalent of Python's dict or MATLAB's struct.
void add_stat(const std::string &s)