File: cs_chol_mex.c

package info (click to toggle)
suitesparse 1%3A7.10.1%2Bdfsg-1
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 254,920 kB
  • sloc: ansic: 1,134,743; cpp: 46,133; makefile: 4,875; fortran: 2,087; java: 1,826; sh: 996; ruby: 725; python: 495; asm: 371; sed: 166; awk: 44
file content (61 lines) | stat: -rw-r--r-- 2,387 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
// CXSparse/MATLAB/CSparse/cs_chol_mex: sparse Cholesky factorization
// CXSparse, Copyright (c) 2006-2022, Timothy A. Davis. All Rights Reserved.
// SPDX-License-Identifier: LGPL-2.1+
#include "cs_mex.h"
/* cs_chol: sparse Cholesky factorization */
void mexFunction (int nargout, mxArray *pargout [ ], int nargin,
    const mxArray *pargin [ ])
{
    int64_t order, n, drop, *p ;
    if (nargout > 2 || nargin < 1 || nargin > 2)
    {
        mexErrMsgTxt ("Usage: [L,p] = cs_chol(A,drop)") ;
    }
    drop = (nargin == 1) ? 1 : mxGetScalar (pargin [1]) ;
    order = (nargout > 1) ? 1 : 0 ;                 /* determine ordering */
    if (mxIsComplex (pargin [0]))
    {
#ifndef NCOMPLEX
        cs_cl Amatrix, *A ;
        cs_cls *S ;
        cs_cln *N ;
        A = cs_cl_mex_get_sparse (&Amatrix, 1, pargin [0]) ;    /* get A */
        n = A->n ;
        S = cs_cl_schol (order, A) ;                /* symbolic Cholesky */
        N = cs_cl_chol (A, S) ;                     /* numeric Cholesky */
        if (!N) mexErrMsgTxt ("cs_chol failed: not positive definite\n") ;
        cs_cl_free (A->x) ;
        if (drop) cs_cl_dropzeros (N->L) ;          /* drop zeros if requested*/
        pargout [0] = cs_cl_mex_put_sparse (&(N->L)) ;      /* return L */
        if (nargout > 1)
        {
            p = cs_cl_pinv (S->pinv, n) ;           /* p=pinv' */
            pargout [1] = cs_dl_mex_put_int (p, n, 1, 1) ; /* return p */
        }
        cs_cl_nfree (N) ;
        cs_cl_sfree (S) ;
#else
        mexErrMsgTxt ("complex matrices not supported") ;
#endif
    }
    else
    {
        cs_dl Amatrix, *A ;
        cs_dls *S ;
        cs_dln *N ;
        A = cs_dl_mex_get_sparse (&Amatrix, 1, 1, pargin [0]) ; /* get A */
        n = A->n ;
        S = cs_dl_schol (order, A) ;                /* symbolic Cholesky */
        N = cs_dl_chol (A, S) ;                     /* numeric Cholesky */
        if (!N) mexErrMsgTxt ("cs_chol failed: not positive definite\n") ;
        if (drop) cs_dl_dropzeros (N->L) ;          /* drop zeros if requested*/
        pargout [0] = cs_dl_mex_put_sparse (&(N->L)) ; /* return L */
        if (nargout > 1)
        {
            p = cs_dl_pinv (S->pinv, n) ;                   /* p=pinv' */
            pargout [1] = cs_dl_mex_put_int (p, n, 1, 1) ; /* return p */
        }
        cs_dl_nfree (N) ;
        cs_dl_sfree (S) ;
    }
}