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 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86
|
// CXSparse/MATLAB/CSparse/cs_qr_mex: sparse QR factorization
// CXSparse, Copyright (c) 2006-2022, Timothy A. Davis. All Rights Reserved.
// SPDX-License-Identifier: LGPL-2.1+
#include "cs_mex.h"
/* cs_qr: sparse QR factorization */
void mexFunction
(
int nargout,
mxArray *pargout [ ],
int nargin,
const mxArray *pargin [ ]
)
{
int64_t m, n, order, *p ;
if (nargout > 5 || nargin != 1)
{
mexErrMsgTxt ("Usage: [V,beta,p,R,q] = cs_qr(A)") ;
}
order = (nargout == 5) ? 3 : 0 ; /* determine ordering */
m = mxGetM (pargin [0]) ;
n = mxGetN (pargin [0]) ;
if (m < n) mexErrMsgTxt ("A must have # rows >= # columns") ;
if (mxIsComplex (pargin [0]))
{
#ifndef NCOMPLEX
cs_cls *S ;
cs_cln *N ;
cs_cl Amatrix, *A, *D ;
A = cs_cl_mex_get_sparse (&Amatrix, 0, pargin [0]) ; /* get A */
S = cs_cl_sqr (order, A, 1) ; /* symbolic QR ordering & analysis*/
N = cs_cl_qr (A, S) ; /* numeric QR factorization */
cs_cl_free (A->x) ;
if (!N) mexErrMsgTxt ("qr failed") ;
cs_cl_dropzeros (N->L) ; /* drop zeros from V and sort */
D = cs_cl_transpose (N->L, 1) ;
cs_cl_spfree (N->L) ;
N->L = cs_cl_transpose (D, 1) ;
cs_cl_spfree (D) ;
cs_cl_dropzeros (N->U) ; /* drop zeros from R and sort */
D = cs_cl_transpose (N->U, 1) ;
cs_cl_spfree (N->U) ;
N->U = cs_cl_transpose (D, 1) ;
cs_cl_spfree (D) ;
m = N->L->m ; /* m may be larger now */
p = cs_cl_pinv (S->pinv, m) ; /* p = pinv' */
pargout [0] = cs_cl_mex_put_sparse (&(N->L)) ; /* return V */
cs_dl_mex_put_double (n, N->B, &(pargout [1])) ; /* return beta */
pargout [2] = cs_dl_mex_put_int (p, m, 1, 1) ; /* return p */
pargout [3] = cs_cl_mex_put_sparse (&(N->U)) ; /* return R */
pargout [4] = cs_dl_mex_put_int (S->q, n, 1, 0) ; /* return q */
cs_cl_nfree (N) ;
cs_cl_sfree (S) ;
#else
mexErrMsgTxt ("complex matrices not supported") ;
#endif
}
else
{
cs_dls *S ;
cs_dln *N ;
cs_dl Amatrix, *A, *D ;
A = cs_dl_mex_get_sparse (&Amatrix, 0, 1, pargin [0]) ; /* get A */
S = cs_dl_sqr (order, A, 1) ; /* symbolic QR ordering & analysis*/
N = cs_dl_qr (A, S) ; /* numeric QR factorization */
if (!N) mexErrMsgTxt ("qr failed") ;
cs_dl_dropzeros (N->L) ; /* drop zeros from V and sort */
D = cs_dl_transpose (N->L, 1) ;
cs_dl_spfree (N->L) ;
N->L = cs_dl_transpose (D, 1) ;
cs_dl_spfree (D) ;
cs_dl_dropzeros (N->U) ; /* drop zeros from R and sort */
D = cs_dl_transpose (N->U, 1) ;
cs_dl_spfree (N->U) ;
N->U = cs_dl_transpose (D, 1) ;
cs_dl_spfree (D) ;
m = N->L->m ; /* m may be larger now */
p = cs_dl_pinv (S->pinv, m) ; /* p = pinv' */
pargout [0] = cs_dl_mex_put_sparse (&(N->L)) ; /* return V */
cs_dl_mex_put_double (n, N->B, &(pargout [1])) ; /* return beta */
pargout [2] = cs_dl_mex_put_int (p, m, 1, 1) ; /* return p */
pargout [3] = cs_dl_mex_put_sparse (&(N->U)) ; /* return R */
pargout [4] = cs_dl_mex_put_int (S->q, n, 1, 0) ; /* return q */
cs_dl_nfree (N) ;
cs_dl_sfree (S) ;
}
}
|