csparse_interface.cpp
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 #include "csparse_interface.hpp"
27 #include "casadi/core/global_options.hpp"
28 
29 namespace casadi {
30 
31  extern "C"
32  int CASADI_LINSOL_CSPARSE_EXPORT
33  casadi_register_linsol_csparse(LinsolInternal::Plugin* plugin) {
34  plugin->creator = CsparseInterface::creator;
35  plugin->name = "csparse";
36  plugin->doc = CsparseInterface::meta_doc.c_str();
37  plugin->version = CASADI_VERSION;
38  plugin->options = &CsparseInterface::options_;
39  plugin->deserialize = &CsparseInterface::deserialize;
40  return 0;
41  }
42 
43  extern "C"
44  void CASADI_LINSOL_CSPARSE_EXPORT casadi_load_linsol_csparse() {
46  }
47 
48  CsparseInterface::CsparseInterface(const std::string& name, const Sparsity& sp)
49  : LinsolInternal(name, sp) {
50  }
51 
53  clear_mem();
54  }
55 
57  if (this->S) cs_sfree(this->S);
58  if (this->N) cs_nfree(this->N);
59  }
60 
61  void CsparseInterface::init(const Dict& opts) {
62  // Call the init method of the base class
64  }
65 
66  int CsparseInterface::init_mem(void* mem) const {
67  if (LinsolInternal::init_mem(mem)) return 1;
68  auto m = static_cast<CsparseMemory*>(mem);
69 
70  m->N = nullptr;
71  m->S = nullptr;
72  m->A.nzmax = this->nnz(); // maximum number of entries
73  m->A.m = this->nrow(); // number of rows
74  m->A.n = this->ncol(); // number of columns
75  m->colind.resize(this->ncol()+1);
76  m->row.resize(this->nnz());
77  copy_vector(this->colind(), m->colind);
78  copy_vector(this->row(), m->row);
79  m->A.p = get_ptr(m->colind); // row pointers (size n+1)
80  m->A.i = get_ptr(m->row); // row pointers (size n+1)
81  m->A.x = nullptr; // numerical values, size nzmax
82  m->A.nz = -1; // of entries in triplet matrix, -1 for compressed-column
83 
84  // Temporary
85  m->temp_.resize(m->A.n);
86  return 0;
87  }
88 
89  int CsparseInterface::sfact(void* mem, const double* A) const {
90  auto m = static_cast<CsparseMemory*>(mem);
91 
92  // Set the nonzeros of the matrix
93  m->A.x = const_cast<double*>(A);
94 
95  // ordering and symbolic analysis
96  casadi_int order = 0; // ordering?
97  if (m->S) cs_sfree(m->S);
98  m->S = cs_sqr(order, &m->A, 0);
99  return 0;
100  }
101 
102  int CsparseInterface::nfact(void* mem, const double* A) const {
103  auto m = static_cast<CsparseMemory*>(mem);
104 
105  // Set the nonzeros of the matrix
106  m->A.x = const_cast<double*>(A);
107 
108  // Make sure that all entries of the linear system are valid
109  for (casadi_int k=0; k<this->nnz(); ++k) {
110  casadi_assert(!isnan(A[k]),
111  "Nonzero " + str(k) + " is not-a-number");
112  casadi_assert(!isinf(A[k]),
113  "Nonzero " + str(k) + " is infinite");
114  }
115 
116  if (verbose_) {
117  uout() << "CsparseInterface::prepare: numeric factorization" << std::endl;
118  uout() << "linear system to be factorized = " << std::endl;
119  DM(sp_, std::vector<double>(A, A+nnz())).print_sparse(uout());
120  }
121 
122  double tol = 1e-8;
123 
124  if (m->N) cs_nfree(m->N);
125  m->N = cs_lu(&m->A, m->S, tol) ; // numeric LU factorization
126  if (m->N==nullptr) {
127  DM temp(sp_, std::vector<double>(A, A+nnz()));
128  temp = sparsify(temp);
129  if (temp.sparsity().is_singular()) {
130  std::stringstream ss;
131  ss << "CsparseInterface::prepare: factorization failed due to matrix"
132  " being singular. Matrix contains numerical zeros which are "
133  "structurally non-zero. Promoting these zeros to be structural "
134  "zeros, the matrix was found to be structurally rank deficient."
135  " sprank: " << sprank(temp.sparsity()) << " <-> " << temp.size2() << std::endl;
136  if (verbose_) {
137  ss << "Sparsity of the linear system: " << std::endl;
138  sp_.disp(ss, true); // print detailed
139  }
140  throw CasadiException(ss.str());
141  } else {
142  std::stringstream ss;
143  ss << "CsparseInterface::prepare: factorization failed, check if Jacobian is singular"
144  << std::endl;
145  if (verbose_) {
146  ss << "Sparsity of the linear system: " << std::endl;
147  sp_.disp(ss, true); // print detailed
148  }
149  throw CasadiException(ss.str());
150  }
151  }
152  casadi_assert_dev(m->N!=nullptr);
153  return 0;
154  }
155 
156  int CsparseInterface::solve(void* mem, const double* A, double* x,
157  casadi_int nrhs, bool tr) const {
158  auto m = static_cast<CsparseMemory*>(mem);
159  casadi_assert_dev(m->N!=nullptr);
160 
161  double *t = &m->temp_.front();
162 
163  for (casadi_int k=0; k<nrhs; ++k) {
164  if (tr) {
165  cs_pvec(m->S->q, x, t, m->A.n) ; // t = P2*b
166  casadi_assert_dev(m->N->U!=nullptr);
167  cs_utsolve(m->N->U, t) ; // t = U'\t
168  cs_ltsolve(m->N->L, t) ; // t = L'\t
169  cs_pvec(m->N->pinv, t, x, m->A.n) ; // x = P1*t
170  } else {
171  cs_ipvec(m->N->pinv, x, t, m->A.n) ; // t = P1\b
172  cs_lsolve(m->N->L, t) ; // t = L\t
173  cs_usolve(m->N->U, t) ; // t = U\t
174  cs_ipvec(m->S->q, t, x, m->A.n) ; // x = P2\t
175  }
176  x += ncol();
177  }
178  return 0;
179  }
180 
181 } // namespace casadi
Casadi exception class.
Definition: exception.hpp:77
static const std::string meta_doc
A documentation string.
int init_mem(void *mem) const override
Initalize memory block.
int sfact(void *mem, const double *A) const override
int solve(void *mem, const double *A, double *x, casadi_int nrhs, bool tr) const override
static ProtoFunction * deserialize(DeserializingStream &s)
Deserialize with type disambiguation.
void init(const Dict &opts) override
Initialize.
CsparseInterface(const std::string &name, const Sparsity &sp)
static LinsolInternal * creator(const std::string &name, const Sparsity &sp)
Create a new LinsolInternal.
int nfact(void *mem, const double *A) const override
Numeric factorization.
casadi_int size2() const
Get the second dimension (i.e. number of columns)
const casadi_int * colind() const
void init(const Dict &opts) override
Initialize.
const casadi_int * row() const
casadi_int nnz() const
casadi_int nrow() const
Get sparsity pattern.
casadi_int ncol() const
int init_mem(void *mem) const override
Initalize memory block.
void print_sparse(std::ostream &stream, bool truncate=true) const
Print sparse matrix style.
const Sparsity & sparsity() const
Const access the sparsity - reference to data member.
static void registerPlugin(const Plugin &plugin, bool needs_lock=true)
Register an integrator in the factory.
bool verbose_
Verbose printout.
static const Options options_
Options.
void clear_mem()
Clear all memory (called from destructor)
void disp(std::ostream &stream, bool more=false) const
Print a description of the object.
General sparsity class.
Definition: sparsity.hpp:106
bool is_singular() const
Check whether the sparsity-pattern indicates structural singularity.
Definition: sparsity.cpp:1299
The casadi namespace.
Definition: archiver.cpp:28
void copy_vector(const std::vector< S > &s, std::vector< D > &d)
int CASADI_LINSOL_CSPARSE_EXPORT casadi_register_linsol_csparse(LinsolInternal::Plugin *plugin)
std::string str(const T &v)
String representation, any type.
GenericType::Dict Dict
C++ equivalent of Python's dict or MATLAB's struct.
T * get_ptr(std::vector< T > &v)
Get a pointer to the data contained in the vector.
Matrix< double > DM
Definition: dm_fwd.hpp:33
void CASADI_LINSOL_CSPARSE_EXPORT casadi_load_linsol_csparse()
std::ostream & uout()