csparse_interface.hpp
1 /*
2  * This file is part of CasADi.
3  *
4  * CasADi -- A symbolic framework for dynamic optimization.
5  * Copyright (C) 2010-2023 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 #ifndef CASADI_CSPARSE_INTERFACE_HPP
27 #define CASADI_CSPARSE_INTERFACE_HPP
28 
39 #include <cs.h>
40 #include "casadi/core/linsol_internal.hpp"
41 #include <casadi/interfaces/csparse/casadi_linsol_csparse_export.h>
42 
43 namespace casadi {
44  struct CsparseMemory : public LinsolMemory {
45  // Destructor
46  ~CsparseMemory();
47 
48  // The linear system CSparse form (CCS)
49  cs A;
50 
51  // The symbolic factorization
52  css *S;
53 
54  // The numeric factorization
55  csn *N;
56 
57  // Temporary
58  std::vector<double> temp_;
59 
60  std::vector<int> colind, row;
61  };
62 
67  class CsparseInterface : public LinsolInternal {
68  public:
69 
70  // Create a linear solver given a sparsity pattern and a number of right hand sides
71  CsparseInterface(const std::string& name, const Sparsity& sp);
72 
74  static LinsolInternal* creator(const std::string& name, const Sparsity& sp) {
75  return new CsparseInterface(name, sp);
76  }
77 
78  // Destructor
79  ~CsparseInterface() override;
80 
81  // Initialize the solver
82  void init(const Dict& opts) override;
83 
85  void* alloc_mem() const override { return new CsparseMemory();}
86 
88  int init_mem(void* mem) const override;
89 
91  void free_mem(void *mem) const override { delete static_cast<CsparseMemory*>(mem);}
92 
93  // Symbolic factorization
94  int sfact(void* mem, const double* A) const override;
95 
96  // Factorize the linear system
97  int nfact(void* mem, const double* A) const override;
98 
99  // Solve the linear system
100  int solve(void* mem, const double* A, double* x, casadi_int nrhs, bool tr) const override;
101 
103  static const std::string meta_doc;
104 
105  // Get name of the plugin
106  const char* plugin_name() const override { return "csparse";}
107 
108  // Get name of the class
109  std::string class_name() const override { return "CsparseInterface";}
110 
112  static ProtoFunction* deserialize(DeserializingStream& s) { return new CsparseInterface(s); }
113 
114  protected:
116  explicit CsparseInterface(DeserializingStream& s) : LinsolInternal(s) {}
117  };
118 
119 } // namespace casadi
120 
122 
123 #endif // CASADI_CSPARSE_INTERFACE_HPP
The casadi namespace.
GenericType::Dict Dict
C++ equivalent of Python's dict or MATLAB's struct.