GNU Octave 10.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-2025 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 "lo-error.h"
37#include "oct-cmplx.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
45
47
48template <typename chol_type>
49class sparse_chol<chol_type>::sparse_chol_rep
50{
51public:
52
53 sparse_chol_rep ()
54 : m_is_pd (false), m_minor_p (0), m_perm (), m_rcond (0)
55#if defined (HAVE_CHOLMOD)
56 , m_L (nullptr), m_common ()
57#endif
58 { }
59
60 sparse_chol_rep (const chol_type& a, bool natural, bool force)
61 : m_is_pd (false), m_minor_p (0), m_perm (), m_rcond (0)
62#if defined (HAVE_CHOLMOD)
63 , m_L (nullptr), m_common ()
64#endif
65 {
66 init (a, natural, force);
67 }
68
69 sparse_chol_rep (const chol_type& a, octave_idx_type& info,
70 bool natural, bool force)
71 : m_is_pd (false), m_minor_p (0), m_perm (), m_rcond (0)
72#if defined (HAVE_CHOLMOD)
73 , m_L (nullptr), m_common ()
74#endif
75 {
76 info = init (a, natural, force);
77 }
78
79 OCTAVE_DISABLE_COPY_MOVE (sparse_chol_rep)
80
81 ~sparse_chol_rep ()
82 {
83#if defined (HAVE_CHOLMOD)
84 if (m_L)
85 CHOLMOD_NAME (free_sparse) (&m_L, &m_common);
86
87 CHOLMOD_NAME(finish) (&m_common);
88#endif
89 }
90
91#if defined (HAVE_CHOLMOD)
92 cholmod_sparse * L () const
93 {
94 return m_L;
95 }
96
97 cholmod_common * cc () const
98 {
99 cholmod_common *ret = const_cast<cholmod_common *> (&m_common);
100 return ret;
101 }
102
103#endif
104
105 octave_idx_type P () const
106 {
107#if defined (HAVE_CHOLMOD)
108 return (m_minor_p == from_size_t (m_L->ncol) ? 0 : m_minor_p + 1);
109#else
110 return 0;
111#endif
112 }
113
114 RowVector perm () const { return m_perm + 1; }
115
116 SparseMatrix Q () const;
117
118 bool is_positive_definite () const { return m_is_pd; }
119
120 double rcond () const { return m_rcond; }
121
122private:
123
124 bool m_is_pd;
125
126 octave_idx_type m_minor_p;
127
128 RowVector m_perm;
129
130 double m_rcond;
131
132#if defined (HAVE_CHOLMOD)
133 cholmod_sparse *m_L;
134
135 cholmod_common m_common;
136
137 void drop_zeros (const cholmod_sparse *S);
138#endif
139
140 octave_idx_type init (const chol_type& a, bool natural, bool force);
141};
142
143#if defined (HAVE_CHOLMOD)
144
145// Can't use CHOLMOD_NAME(drop)(0.0, S, cm) because it doesn't treat
146// complex matrices.
147
148template <typename chol_type>
149void
151{
152 if (! S)
153 return;
154
155 octave_idx_type *Sp = static_cast<octave_idx_type *> (S->p);
156 octave_idx_type *Si = static_cast<octave_idx_type *> (S->i);
157 chol_elt *Sx = static_cast<chol_elt *> (S->x);
158
159 octave_idx_type pdest = 0;
160 octave_idx_type ncol = S->ncol;
161
162 for (octave_idx_type k = 0; k < ncol; k++)
163 {
164 octave_idx_type p = Sp[k];
165 octave_idx_type pend = Sp[k+1];
166 Sp[k] = pdest;
167
168 for (; p < pend; p++)
169 {
170 chol_elt sik = Sx[p];
171
172 if (CHOLMOD_IS_NONZERO (sik))
173 {
174 if (p != pdest)
175 {
176 Si[pdest] = Si[p];
177 Sx[pdest] = sik;
178 }
179
180 pdest++;
181 }
182 }
183 }
184
185 Sp[ncol] = pdest;
186}
187
188// Must provide a specialization for this function.
189template <typename T>
190int
192
193template <>
194inline int
196{
197 return CHOLMOD_REAL;
198}
199
200template <>
201inline int
203{
204 return CHOLMOD_COMPLEX;
205}
206
207#endif
208
209template <typename chol_type>
212 bool natural, bool force)
213{
214 octave_idx_type info = 0;
215
216#if defined (HAVE_CHOLMOD)
217
218 octave_idx_type a_nr = a.rows ();
219 octave_idx_type a_nc = a.cols ();
220
221 if (a_nr != a_nc)
222 (*current_liboctave_error_handler)
223 ("sparse_chol requires square matrix");
224
225 cholmod_common *cm = &m_common;
226
227 // Setup initial parameters
228
229 CHOLMOD_NAME(start) (cm);
230 cm->prefer_zomplex = false;
231
232 double spu = sparse_params::get_key ("spumoni");
233
234 if (spu == 0.)
235 {
236 cm->print = -1;
237 SUITESPARSE_ASSIGN_FPTR (printf_func, cm->print_function, nullptr);
238 }
239 else
240 {
241 cm->print = static_cast<int> (spu) + 2;
242 SUITESPARSE_ASSIGN_FPTR (printf_func, cm->print_function,
244 }
245
246 cm->error_handler = &SparseCholError;
247
248 SUITESPARSE_ASSIGN_FPTR2 (divcomplex_func, cm->complex_divide,
249 divcomplex);
250
251 SUITESPARSE_ASSIGN_FPTR2 (hypot_func, cm->hypotenuse, hypot);
252
253 cm->final_asis = false;
254 cm->final_super = false;
255 cm->final_ll = true;
256 cm->final_pack = true;
257 cm->final_monotonic = true;
258 cm->final_resymbol = false;
259
260 cholmod_sparse A;
261 cholmod_sparse *ac = &A;
262 double dummy;
263
264 ac->nrow = a_nr;
265 ac->ncol = a_nc;
266
267 ac->p = a.cidx ();
268 ac->i = a.ridx ();
269 ac->nzmax = a.nnz ();
270 ac->packed = true;
271 ac->sorted = true;
272 ac->nz = nullptr;
273#if defined (OCTAVE_ENABLE_64)
274 ac->itype = CHOLMOD_LONG;
275#else
276 ac->itype = CHOLMOD_INT;
277#endif
278 ac->dtype = CHOLMOD_DOUBLE;
279 ac->stype = 1;
280 ac->xtype = get_xtype<chol_elt> ();
281
282 if (a_nr < 1)
283 ac->x = &dummy;
284 else
285 ac->x = a.data ();
286
287 // use natural ordering if no q output parameter
288 if (natural)
289 {
290 cm->nmethods = 1;
291 cm->method[0].ordering = CHOLMOD_NATURAL;
292 cm->postorder = false;
293 }
294
295 cholmod_factor *Lfactor = CHOLMOD_NAME(analyze) (ac, cm);
296 CHOLMOD_NAME(factorize) (ac, Lfactor, cm);
297
298 m_is_pd = cm->status == CHOLMOD_OK;
299 info = (m_is_pd ? 0 : cm->status);
300
301 if (m_is_pd || force)
302 {
303 m_rcond = CHOLMOD_NAME(rcond) (Lfactor, cm);
304
305 m_minor_p = Lfactor->minor;
306
307 m_L = CHOLMOD_NAME(factor_to_sparse) (Lfactor, cm);
308
309 if (m_minor_p > 0 && m_minor_p < a_nr)
310 {
311 std::size_t n1 = a_nr + 1;
312 m_L->p = CHOLMOD_NAME(realloc) (m_minor_p+1,
313 sizeof(octave_idx_type),
314 m_L->p, &n1, cm);
315
316 CHOLMOD_NAME(reallocate_sparse)
317 (static_cast<octave_idx_type *> (m_L->p)[m_minor_p],
318 m_L, cm);
319
320 m_L->ncol = m_minor_p;
321 }
322
323 drop_zeros (m_L);
324
325 if (! natural)
326 {
327 m_perm.resize (a_nr);
328 for (octave_idx_type i = 0; i < a_nr; i++)
329 m_perm(i) = static_cast<octave_idx_type *> (Lfactor->Perm)[i];
330 }
331 }
332
333 // NAME used to prefix statistics report from print_common
334 static char blank_name[] = " ";
335
336 CHOLMOD_NAME(print_common) (blank_name, cm);
337 CHOLMOD_NAME(free_factor) (&Lfactor, cm);
338
339 return info;
340
341#else
342
343 octave_unused_parameter (a);
344 octave_unused_parameter (natural);
345 octave_unused_parameter (force);
346
347 (*current_liboctave_error_handler)
348 ("support for CHOLMOD was unavailable or disabled when liboctave was built");
349
350 return info;
351
352#endif
353}
354
355template <typename chol_type>
358{
359#if defined (HAVE_CHOLMOD)
360
361 octave_idx_type n = m_L->nrow;
362 SparseMatrix p (n, n, n);
363
364 for (octave_idx_type i = 0; i < n; i++)
365 {
366 p.xcidx (i) = i;
367 p.xridx (i) = static_cast<octave_idx_type> (m_perm (i));
368 p.xdata (i) = 1;
369 }
370
371 p.xcidx (n) = n;
372
373 return p;
374
375#else
376
377 return SparseMatrix ();
378
379#endif
380}
381
382template <typename chol_type>
384 : m_rep (new typename sparse_chol<chol_type>::sparse_chol_rep ())
385{ }
386
387template <typename chol_type>
388sparse_chol<chol_type>::sparse_chol (const chol_type& a, bool natural,
389 bool force)
390 : m_rep (new typename
391 sparse_chol<chol_type>::sparse_chol_rep (a, natural, force))
392{ }
393
394template <typename chol_type>
396 octave_idx_type& info,
397 bool natural, bool force)
398 : m_rep (new typename
399 sparse_chol<chol_type>::sparse_chol_rep (a, info, natural, force))
400{ }
401
402template <typename chol_type>
404 octave_idx_type& info,
405 bool natural)
406 : m_rep (new typename
407 sparse_chol<chol_type>::sparse_chol_rep (a, info, natural, false))
408{ }
409
410template <typename chol_type>
412 octave_idx_type& info)
413 : m_rep (new typename
414 sparse_chol<chol_type>::sparse_chol_rep (a, info, false, false))
415{ }
416
417template <typename chol_type>
418chol_type
420{
421#if defined (HAVE_CHOLMOD)
422
423 cholmod_sparse *m = m_rep->L ();
424
425 octave_idx_type nc = m->ncol;
427 = from_suitesparse_long (CHOLMOD_NAME(nnz) (m, m_rep->cc ()));
428
429 chol_type ret (m->nrow, nc, nnz);
430
431 for (octave_idx_type j = 0; j < nc+1; j++)
432 ret.xcidx (j) = static_cast<octave_idx_type *> (m->p)[j];
433
434 for (octave_idx_type i = 0; i < nnz; i++)
435 {
436 ret.xridx (i) = static_cast<octave_idx_type *> (m->i)[i];
437 ret.xdata (i) = static_cast<chol_elt *> (m->x)[i];
438 }
439
440 return ret;
441
442#else
443
444 return chol_type ();
445
446#endif
447}
448
449template <typename chol_type>
452{
453 return m_rep->P ();
454}
455
456template <typename chol_type>
459{
460 return m_rep->perm ();
461}
462
463template <typename chol_type>
466{
467 return m_rep->Q ();
468}
469
470template <typename chol_type>
471bool
473{
474 return m_rep->is_positive_definite ();
475}
476
477template <typename chol_type>
478double
480{
481 return m_rep->rcond ();
482}
483
484template <typename chol_type>
485chol_type
487{
488 chol_type retval;
489
490#if defined (HAVE_CHOLMOD)
491
492 cholmod_sparse *m = m_rep->L ();
493 octave_idx_type n = m->ncol;
494 RowVector m_perm = m_rep->perm ();
495 double rcond2;
496 octave_idx_type info;
498 chol_type linv = L ().hermitian ().inverse (mattype, info, rcond2, 1, 0);
499
500 if (m_perm.numel () == n)
501 {
502 SparseMatrix Qc = Q ();
503
504 retval = Qc * linv * linv.hermitian () * Qc.transpose ();
505 }
506 else
507 retval = linv * linv.hermitian ();
508
509#endif
510
511 return retval;
512}
513
514template <typename chol_type>
515chol_type
516chol2inv (const chol_type& r)
517{
518 octave_idx_type r_nr = r.rows ();
519 octave_idx_type r_nc = r.cols ();
520 chol_type retval;
521
522 if (r_nr != r_nc)
523 (*current_liboctave_error_handler) ("U must be a square matrix");
524
525 MatrixType mattype (r);
526 int typ = mattype.type (false);
527 double rcond;
528 octave_idx_type info;
529 chol_type rtra, multip;
530
531 if (typ == MatrixType::Upper)
532 {
533 rtra = r.transpose ();
534 multip = (rtra*r);
535 }
536 else if (typ == MatrixType::Lower)
537 {
538 rtra = r.transpose ();
539 multip = (r*rtra);
540 }
541 else
542 (*current_liboctave_error_handler) ("U must be a triangular matrix");
543
544 MatrixType mattypenew (multip);
545 retval = multip.inverse (mattypenew, info, rcond, true, false);
546 return retval;
547}
548
549// SparseComplexMatrix specialization (the value for the NATURAL
550// parameter in the sparse_chol<T>::sparse_chol_rep constructor is
551// different from the default).
552
553template <>
556 octave_idx_type& info)
557 : m_rep (new sparse_chol<SparseComplexMatrix>::sparse_chol_rep (a, info,
558 true,
559 false))
560{ }
561
562// Instantiations we need.
563
565
567
570
573
574OCTAVE_END_NAMESPACE(math)
575OCTAVE_END_NAMESPACE(octave)
octave_idx_type numel() const
Number of elements in the array.
Definition Array.h:418
SparseMatrix transpose() const
Definition dSparse.h:125
SparseMatrix hermitian() const
Definition dSparse.h:129
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:87
chol_type L() const
static double get_key(const std::string &key)
OCTAVE_BEGIN_NAMESPACE(octave) static octave_value daspk_fcn
OCTAVE_NORETURN liboctave_error_handler current_liboctave_error_handler
Definition lo-error.c:41
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
#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,...)