GNU Octave 11.1.0
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
 
Loading...
Searching...
No Matches
sparse-chol.cc
Go to the documentation of this file.
1////////////////////////////////////////////////////////////////////////
2//
3// Copyright (C) 1998-2026 The Octave Project Developers
4//
5// See the file COPYRIGHT.md in the top-level directory of this
6// distribution or <https://octave.org/copyright/>.
7//
8// This file is part of Octave.
9//
10// Octave is free software: you can redistribute it and/or modify it
11// under the terms of the GNU General Public License as published by
12// the Free Software Foundation, either version 3 of the License, or
13// (at your option) any later version.
14//
15// Octave is distributed in the hope that it will be useful, but
16// WITHOUT ANY WARRANTY; without even the implied warranty of
17// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18// GNU General Public License for more details.
19//
20// You should have received a copy of the GNU General Public License
21// along with Octave; see the file COPYING. If not, see
22// <https://www.gnu.org/licenses/>.
23//
24////////////////////////////////////////////////////////////////////////
25
26#if defined (HAVE_CONFIG_H)
27# include "config.h"
28#endif
29
30#include <cstddef>
31
32#include "CSparse.h"
33#include "MatrixType.h"
34#include "dRowVector.h"
35#include "dSparse.h"
36#include "oct-cmplx.h"
37#include "oct-error.h"
38#include "oct-sparse.h"
39#include "oct-spparms.h"
40#include "quit.h"
41#include "sparse-chol.h"
42#include "sparse-util.h"
43
46
47template <typename chol_type>
48class sparse_chol<chol_type>::sparse_chol_rep
49{
50public:
51
52 sparse_chol_rep ()
53 : m_is_pd (false), m_minor_p (0), m_perm (), m_rcond (0)
54#if defined (HAVE_CHOLMOD)
55 , m_L (nullptr), m_common ()
56#endif
57 { }
58
59 sparse_chol_rep (const chol_type& a, bool natural, bool force)
60 : m_is_pd (false), m_minor_p (0), m_perm (), m_rcond (0)
61#if defined (HAVE_CHOLMOD)
62 , m_L (nullptr), m_common ()
63#endif
64 {
65 init (a, natural, force);
66 }
67
68 sparse_chol_rep (const chol_type& a, octave_idx_type& info,
69 bool natural, bool force)
70 : m_is_pd (false), m_minor_p (0), m_perm (), m_rcond (0)
71#if defined (HAVE_CHOLMOD)
72 , m_L (nullptr), m_common ()
73#endif
74 {
75 info = init (a, natural, force);
76 }
77
78 OCTAVE_DISABLE_COPY_MOVE (sparse_chol_rep)
79
80 ~sparse_chol_rep ()
81 {
82#if defined (HAVE_CHOLMOD)
83 if (m_L)
84 CHOLMOD_NAME (free_sparse) (&m_L, &m_common);
85
86 CHOLMOD_NAME(finish) (&m_common);
87#endif
88 }
89
90#if defined (HAVE_CHOLMOD)
91 cholmod_sparse * L () const
92 {
93 return m_L;
94 }
95
96 cholmod_common * cc () const
97 {
98 cholmod_common *ret = const_cast<cholmod_common *> (&m_common);
99 return ret;
100 }
101
102#endif
103
104 octave_idx_type P () const
105 {
106#if defined (HAVE_CHOLMOD)
107 return (m_minor_p == from_size_t (m_L->ncol) ? 0 : m_minor_p + 1);
108#else
109 return 0;
110#endif
111 }
112
113 RowVector perm () const { return m_perm + 1; }
114
115 SparseMatrix Q () const;
116
117 bool is_positive_definite () const { return m_is_pd; }
118
119 double rcond () const { return m_rcond; }
120
121private:
122
123 bool m_is_pd;
124
125 octave_idx_type m_minor_p;
126
127 RowVector m_perm;
128
129 double m_rcond;
130
131#if defined (HAVE_CHOLMOD)
132 cholmod_sparse *m_L;
133
134 cholmod_common m_common;
135
136 void drop_zeros (const cholmod_sparse *S);
137#endif
138
139 octave_idx_type init (const chol_type& a, bool natural, bool force);
140};
141
142#if defined (HAVE_CHOLMOD)
143
144// Can't use CHOLMOD_NAME(drop)(0.0, S, cm) because it doesn't treat
145// complex matrices.
146
147template <typename chol_type>
148void
150{
151 if (! S)
152 return;
153
154 octave_idx_type *Sp = static_cast<octave_idx_type *> (S->p);
155 octave_idx_type *Si = static_cast<octave_idx_type *> (S->i);
156 chol_elt *Sx = static_cast<chol_elt *> (S->x);
157
158 octave_idx_type pdest = 0;
159 octave_idx_type ncol = S->ncol;
160
161 for (octave_idx_type k = 0; k < ncol; k++)
162 {
163 octave_idx_type p = Sp[k];
164 octave_idx_type pend = Sp[k+1];
165 Sp[k] = pdest;
166
167 for (; p < pend; p++)
168 {
169 chol_elt sik = Sx[p];
170
171 if (CHOLMOD_IS_NONZERO (sik))
172 {
173 if (p != pdest)
174 {
175 Si[pdest] = Si[p];
176 Sx[pdest] = sik;
177 }
178
179 pdest++;
180 }
181 }
182 }
183
184 Sp[ncol] = pdest;
185}
186
187// Must provide a specialization for this function.
188template <typename T>
189int
191
192template <>
193inline int
195{
196 return CHOLMOD_REAL;
197}
198
199template <>
200inline int
202{
203 return CHOLMOD_COMPLEX;
204}
205
206#endif
207
208template <typename chol_type>
211 bool natural, bool force)
212{
213 octave_idx_type info = 0;
214
215#if defined (HAVE_CHOLMOD)
216
217 octave_idx_type a_nr = a.rows ();
218 octave_idx_type a_nc = a.cols ();
219
220 if (a_nr != a_nc)
221 (*current_liboctave_error_handler)
222 ("sparse_chol requires square matrix");
223
224 cholmod_common *cm = &m_common;
225
226 // Setup initial parameters
227
228 CHOLMOD_NAME(start) (cm);
229 cm->prefer_zomplex = false;
230
231 double spu = sparse_params::get_key ("spumoni");
232
233 if (spu == 0.)
234 {
235 cm->print = -1;
236 SUITESPARSE_ASSIGN_FPTR (printf_func, cm->print_function, nullptr);
237 }
238 else
239 {
240 cm->print = static_cast<int> (spu) + 2;
241 SUITESPARSE_ASSIGN_FPTR (printf_func, cm->print_function,
243 }
244
245 cm->error_handler = &SparseCholError;
246
247 SUITESPARSE_ASSIGN_FPTR2 (divcomplex_func, cm->complex_divide,
248 divcomplex);
249
250 SUITESPARSE_ASSIGN_FPTR2 (hypot_func, cm->hypotenuse, hypot);
251
252 cm->final_asis = false;
253 cm->final_super = false;
254 cm->final_ll = true;
255 cm->final_pack = true;
256 cm->final_monotonic = true;
257 cm->final_resymbol = false;
258
259 cholmod_sparse A;
260 cholmod_sparse *ac = &A;
261 double dummy;
262
263 ac->nrow = a_nr;
264 ac->ncol = a_nc;
265
266 ac->p = a.cidx ();
267 ac->i = a.ridx ();
268 ac->nzmax = a.nnz ();
269 ac->packed = true;
270 ac->sorted = true;
271 ac->nz = nullptr;
272#if defined (OCTAVE_ENABLE_64)
273 ac->itype = CHOLMOD_LONG;
274#else
275 ac->itype = CHOLMOD_INT;
276#endif
277 ac->dtype = CHOLMOD_DOUBLE;
278 ac->stype = 1;
279 ac->xtype = get_xtype<chol_elt> ();
280
281 if (a_nr < 1)
282 ac->x = &dummy;
283 else
284 ac->x = a.data ();
285
286 // use natural ordering if no q output parameter
287 if (natural)
288 {
289 cm->nmethods = 1;
290 cm->method[0].ordering = CHOLMOD_NATURAL;
291 cm->postorder = false;
292 }
293
294 cholmod_factor *Lfactor = CHOLMOD_NAME(analyze) (ac, cm);
295 CHOLMOD_NAME(factorize) (ac, Lfactor, cm);
296
297 m_is_pd = cm->status == CHOLMOD_OK;
298 info = (m_is_pd ? 0 : cm->status);
299
300 if (m_is_pd || force)
301 {
302 m_rcond = CHOLMOD_NAME(rcond) (Lfactor, cm);
303
304 m_minor_p = Lfactor->minor;
305
306 m_L = CHOLMOD_NAME(factor_to_sparse) (Lfactor, cm);
307
308 if (m_minor_p > 0 && m_minor_p < a_nr)
309 {
310 std::size_t n1 = a_nr + 1;
311 m_L->p = CHOLMOD_NAME(realloc) (m_minor_p+1,
312 sizeof(octave_idx_type),
313 m_L->p, &n1, cm);
314
315 CHOLMOD_NAME(reallocate_sparse)
316 (static_cast<octave_idx_type *> (m_L->p)[m_minor_p],
317 m_L, cm);
318
319 m_L->ncol = m_minor_p;
320 }
321
322 drop_zeros (m_L);
323
324 if (! natural)
325 {
326 m_perm.resize (a_nr);
327 for (octave_idx_type i = 0; i < a_nr; i++)
328 m_perm(i) = static_cast<octave_idx_type *> (Lfactor->Perm)[i];
329 }
330 }
331
332 // NAME used to prefix statistics report from print_common
333 static char blank_name[] = " ";
334
335 CHOLMOD_NAME(print_common) (blank_name, cm);
336 CHOLMOD_NAME(free_factor) (&Lfactor, cm);
337
338 return info;
339
340#else
341
342 octave_unused_parameter (a);
343 octave_unused_parameter (natural);
344 octave_unused_parameter (force);
345
346 (*current_liboctave_error_handler)
347 ("support for CHOLMOD was unavailable or disabled when liboctave was built");
348
349 return info;
350
351#endif
352}
353
354template <typename chol_type>
357{
358#if defined (HAVE_CHOLMOD)
359
360 octave_idx_type n = m_L->nrow;
361 SparseMatrix p (n, n, n);
362
363 for (octave_idx_type i = 0; i < n; i++)
364 {
365 p.xcidx (i) = i;
366 p.xridx (i) = static_cast<octave_idx_type> (m_perm (i));
367 p.xdata (i) = 1;
368 }
369
370 p.xcidx (n) = n;
371
372 return p;
373
374#else
375
376 return SparseMatrix ();
377
378#endif
379}
380
381template <typename chol_type>
383 : m_rep (new typename sparse_chol<chol_type>::sparse_chol_rep ())
384{ }
385
386template <typename chol_type>
387sparse_chol<chol_type>::sparse_chol (const chol_type& a, bool natural,
388 bool force)
389 : m_rep (new typename
390 sparse_chol<chol_type>::sparse_chol_rep (a, natural, force))
391{ }
392
393template <typename chol_type>
395 octave_idx_type& info,
396 bool natural, bool force)
397 : m_rep (new typename
398 sparse_chol<chol_type>::sparse_chol_rep (a, info, natural, force))
399{ }
400
401template <typename chol_type>
403 octave_idx_type& info,
404 bool natural)
405 : m_rep (new typename
406 sparse_chol<chol_type>::sparse_chol_rep (a, info, natural, false))
407{ }
408
409template <typename chol_type>
411 octave_idx_type& info)
412 : m_rep (new typename
413 sparse_chol<chol_type>::sparse_chol_rep (a, info, false, false))
414{ }
415
416template <typename chol_type>
417chol_type
419{
420#if defined (HAVE_CHOLMOD)
421
422 cholmod_sparse *m = m_rep->L ();
423
424 octave_idx_type nc = m->ncol;
426 = from_suitesparse_long (CHOLMOD_NAME(nnz) (m, m_rep->cc ()));
427
428 chol_type ret (m->nrow, nc, nnz);
429
430 for (octave_idx_type j = 0; j < nc+1; j++)
431 ret.xcidx (j) = static_cast<octave_idx_type *> (m->p)[j];
432
433 for (octave_idx_type i = 0; i < nnz; i++)
434 {
435 ret.xridx (i) = static_cast<octave_idx_type *> (m->i)[i];
436 ret.xdata (i) = static_cast<chol_elt *> (m->x)[i];
437 }
438
439 return ret;
440
441#else
442
443 return chol_type ();
444
445#endif
446}
447
448template <typename chol_type>
451{
452 return m_rep->P ();
453}
454
455template <typename chol_type>
458{
459 return m_rep->perm ();
460}
461
462template <typename chol_type>
465{
466 return m_rep->Q ();
467}
468
469template <typename chol_type>
470bool
472{
473 return m_rep->is_positive_definite ();
474}
475
476template <typename chol_type>
477double
479{
480 return m_rep->rcond ();
481}
482
483template <typename chol_type>
484chol_type
486{
487 chol_type retval;
488
489#if defined (HAVE_CHOLMOD)
490
491 cholmod_sparse *m = m_rep->L ();
492 octave_idx_type n = m->ncol;
493 RowVector m_perm = m_rep->perm ();
494 double rcond2;
495 octave_idx_type info;
497 chol_type linv = L ().hermitian ().inverse (mattype, info, rcond2, 1, 0);
498
499 if (m_perm.numel () == n)
500 {
501 SparseMatrix Qc = Q ();
502
503 retval = Qc * linv * linv.hermitian () * Qc.transpose ();
504 }
505 else
506 retval = linv * linv.hermitian ();
507
508#endif
509
510 return retval;
511}
512
513template <typename chol_type>
514chol_type
515chol2inv (const chol_type& r)
516{
517 octave_idx_type r_nr = r.rows ();
518 octave_idx_type r_nc = r.cols ();
519 chol_type retval;
520
521 if (r_nr != r_nc)
522 (*current_liboctave_error_handler) ("U must be a square matrix");
523
524 MatrixType mattype (r);
525 int typ = mattype.type (false);
526 double rcond;
527 octave_idx_type info;
528 chol_type rtra, multip;
529
530 if (typ == MatrixType::Upper)
531 {
532 rtra = r.transpose ();
533 multip = (rtra*r);
534 }
535 else if (typ == MatrixType::Lower)
536 {
537 rtra = r.transpose ();
538 multip = (r*rtra);
539 }
540 else
541 (*current_liboctave_error_handler) ("U must be a triangular matrix");
542
543 MatrixType mattypenew (multip);
544 retval = multip.inverse (mattypenew, info, rcond, true, false);
545 return retval;
546}
547
548// SparseComplexMatrix specialization (the value for the NATURAL
549// parameter in the sparse_chol<T>::sparse_chol_rep constructor is
550// different from the default).
551
552template <>
555 octave_idx_type& info)
556 : m_rep (new sparse_chol<SparseComplexMatrix>::sparse_chol_rep (a, info,
557 true,
558 false))
559{ }
560
561// Instantiations we need.
562
564
566
569
572
573OCTAVE_END_NAMESPACE(math)
574OCTAVE_END_NAMESPACE(octave)
octave_idx_type numel() const
Number of elements in the array.
Definition Array-base.h:440
SparseMatrix transpose() const
Definition dSparse.h:131
SparseMatrix hermitian() const
Definition dSparse.h:135
octave_idx_type P() const
double rcond() const
bool is_positive_definite() const
SparseMatrix Q() const
RowVector perm() const
chol_type inverse() const
chol_type::element_type chol_elt
Definition sparse-chol.h:86
chol_type L() const
static double get_key(const std::string &key)
OCTAVE_BEGIN_NAMESPACE(octave) static octave_value daspk_fcn
F77_RET_T const F77_INT const F77_INT const F77_INT F77_DBLE const F77_INT F77_DBLE const F77_INT F77_DBLE * Q
F77_RET_T const F77_INT F77_CMPLX * A
OCTAVE_NORETURN liboctave_error_handler current_liboctave_error_handler
Definition lo-error.c:41
#define OCTAVE_API
Definition main.in.cc:55
#define CHOLMOD_NAME(name)
Definition oct-sparse.h:140
octave_idx_type from_suitesparse_long(SuiteSparse_long x)
Definition oct-sparse.h:200
octave_idx_type from_size_t(std::size_t x)
Definition oct-sparse.h:211
int get_xtype< Complex >()
template SparseMatrix chol2inv< SparseMatrix >(const SparseMatrix &r)
int get_xtype< double >()
chol_type chol2inv(const chol_type &r)
int get_xtype()
template SparseComplexMatrix chol2inv< SparseComplexMatrix >(const SparseComplexMatrix &r)
void SparseCholError(int status, char *file, int line, char *message)
int SparseCholPrint(const char *fmt,...)