26 #include "csparse_interface.hpp"
27 #include "casadi/core/global_options.hpp"
32 int CASADI_LINSOL_CSPARSE_EXPORT
35 plugin->name =
"csparse";
37 plugin->version = CASADI_VERSION;
57 if (this->
S) cs_sfree(this->
S);
58 if (this->
N) cs_nfree(this->
N);
72 m->A.nzmax = this->
nnz();
73 m->A.m = this->
nrow();
74 m->A.n = this->
ncol();
75 m->colind.resize(this->
ncol()+1);
76 m->row.resize(this->
nnz());
85 m->temp_.resize(m->A.n);
93 m->
A.x =
const_cast<double*
>(A);
97 if (m->S) cs_sfree(m->S);
98 m->S = cs_sqr(order, &m->A, 0);
106 m->
A.x =
const_cast<double*
>(A);
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");
117 uout() <<
"CsparseInterface::prepare: numeric factorization" << std::endl;
118 uout() <<
"linear system to be factorized = " << std::endl;
124 if (m->N) cs_nfree(m->N);
125 m->N = cs_lu(&m->A, m->S, tol) ;
127 DM temp(
sp_, std::vector<double>(A, A+
nnz()));
128 temp = sparsify(temp);
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;
137 ss <<
"Sparsity of the linear system: " << std::endl;
142 std::stringstream ss;
143 ss <<
"CsparseInterface::prepare: factorization failed, check if Jacobian is singular"
146 ss <<
"Sparsity of the linear system: " << std::endl;
152 casadi_assert_dev(m->N!=
nullptr);
157 casadi_int nrhs,
bool tr)
const {
159 casadi_assert_dev(m->N!=
nullptr);
161 double *t = &m->temp_.front();
163 for (casadi_int k=0; k<nrhs; ++k) {
165 cs_pvec(m->S->q, x, t, m->A.n) ;
166 casadi_assert_dev(m->N->U!=
nullptr);
167 cs_utsolve(m->N->U, t) ;
168 cs_ltsolve(m->N->L, t) ;
169 cs_pvec(m->N->pinv, t, x, m->A.n) ;
171 cs_ipvec(m->N->pinv, x, t, m->A.n) ;
172 cs_lsolve(m->N->L, t) ;
173 cs_usolve(m->N->U, t) ;
174 cs_ipvec(m->S->q, t, x, m->A.n) ;
~CsparseInterface() override
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 nrow() const
Get sparsity pattern.
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.
bool is_singular() const
Check whether the sparsity-pattern indicates structural singularity.
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.
void CASADI_LINSOL_CSPARSE_EXPORT casadi_load_linsol_csparse()