void mexFunction ( int nargout, mxArray *pargout [ ], int nargin, const mxArray *pargin [ ] ) { Int *Ap, *Ai, *E, *Bp, *Bi, *HPinv ; double *Ax, *Bx, dummy, tol ; Int m, n, anz, bnz, is_complex, econ, A_complex, B_complex ; spqr_mx_options opts ; cholmod_sparse *A, Amatrix, *R, *Q, *Csparse, Bsmatrix, *Bsparse, *H ; cholmod_dense *Cdense, Bdmatrix, *Bdense, *HTau ; cholmod_common Common, *cc ; #ifdef TIMING double t0 = (nargout > 3) ? spqr_time ( ) : 0 ; #endif // ------------------------------------------------------------------------- // start CHOLMOD and set parameters // ------------------------------------------------------------------------- cc = &Common ; cholmod_l_start (cc) ; spqr_mx_config (SPUMONI, cc) ; // ------------------------------------------------------------------------- // check inputs // ------------------------------------------------------------------------- // nargin can be 1, 2, or 3 // nargout can be 0, 1, 2, 3, or 4 if (nargout > 4) { mexErrMsgIdAndTxt ("MATLAB:maxlhs", "Too many output arguments") ; } if (nargin < 1) { mexErrMsgIdAndTxt ("MATLAB:minrhs", "Not enough input arguments") ; } if (nargin > 3) { mexErrMsgIdAndTxt ("MATLAB:maxrhs", "Too many input arguments") ; } // ------------------------------------------------------------------------- // get the input matrix A // ------------------------------------------------------------------------- if (!mxIsSparse (pargin [0])) { mexErrMsgIdAndTxt ("QR:invalidInput", "A must be sparse") ; } A = spqr_mx_get_sparse (pargin [0], &Amatrix, &dummy) ; Ap = (Int *) A->p ; Ai = (Int *) A->i ; m = A->nrow ; n = A->ncol ; A_complex = mxIsComplex (pargin [0]) ; // ------------------------------------------------------------------------- // determine usage and parameters // ------------------------------------------------------------------------- if (nargin == 1) { // --------------------------------------------------------------------- // [ ] = qr (A) // --------------------------------------------------------------------- spqr_mx_get_options (NULL, &opts, m, nargout, cc) ; // R = qr (A) // [Q,R] = qr (A) // [Q,R,E] = qr (A) } else if (nargin == 2) { // --------------------------------------------------------------------- // [ ] = qr (A,0), [ ] = qr (A,opts), or [ ] = qr (A,B) // --------------------------------------------------------------------- if (is_zero (pargin [1])) { // ----------------------------------------------------------------- // [ ... ] = qr (A,0) // ----------------------------------------------------------------- spqr_mx_get_options (NULL, &opts, m, nargout, cc) ; opts.econ = n ; opts.permvector = TRUE ; // R = qr (A,0) // [Q,R] = qr (A,0) // [Q,R,E] = qr (A,0) } else if (mxIsEmpty (pargin [1]) || mxIsStruct (pargin [1])) { // ----------------------------------------------------------------- // [ ] = qr (A,opts) // ----------------------------------------------------------------- spqr_mx_get_options (pargin [1], &opts, m, nargout, cc) ; // R = qr (A,opts) // [Q,R] = qr (A,opts) // [Q,R,E] = qr (A,opts) } else { // ----------------------------------------------------------------- // [ ] = qr (A,B) // ----------------------------------------------------------------- spqr_mx_get_options (NULL, &opts, m, nargout, cc) ; opts.haveB = TRUE ; opts.Qformat = SPQR_Q_DISCARD ; if (nargout <= 1) { mexErrMsgIdAndTxt ("MATLAB:minlhs", "Not enough output arguments") ; } // [C,R] = qr (A,B) // [C,R,E] = qr (A,B) } } else // if (nargin == 3) { // --------------------------------------------------------------------- // [ ] = qr (A,B,opts) // --------------------------------------------------------------------- if (is_zero (pargin [2])) { // ----------------------------------------------------------------- // [ ] = qr (A,B,0) // ----------------------------------------------------------------- spqr_mx_get_options (NULL, &opts, m, nargout, cc) ; opts.econ = n ; opts.permvector = TRUE ; opts.haveB = TRUE ; opts.Qformat = SPQR_Q_DISCARD ; if (nargout <= 1) { mexErrMsgIdAndTxt ("MATLAB:minlhs", "Not enough output arguments") ; } // [C,R] = qr (A,B,0) // [C,R,E] = qr (A,B,0) } else if (mxIsEmpty (pargin [2]) || mxIsStruct (pargin [2])) { // ----------------------------------------------------------------- // [ ] = qr (A,B,opts) // ----------------------------------------------------------------- spqr_mx_get_options (pargin [2], &opts, m, nargout, cc) ; opts.haveB = TRUE ; opts.Qformat = SPQR_Q_DISCARD ; if (nargout <= 1) { mexErrMsgIdAndTxt ("MATLAB:minlhs", "Not enough output arguments") ; } // [C,R] = qr (A,B,opts) // [C,R,E] = qr (A,B,opts) } else { mexErrMsgIdAndTxt ("QR:invalidInput", "Invalid opts argument") ; } } int order = opts.ordering ; tol = opts.tol ; econ = opts.econ ; // ------------------------------------------------------------------------- // get A and convert to merged-complex if needed // ------------------------------------------------------------------------- if (opts.haveB) { B_complex = mxIsComplex (pargin [1]) ; } else { B_complex = FALSE ; } is_complex = (A_complex || B_complex) ; Ax = spqr_mx_merge_if_complex (pargin [0], is_complex, &anz, cc) ; if (is_complex) { // A has been converted from real or zomplex to complex A->x = Ax ; A->z = NULL ; A->xtype = CHOLMOD_COMPLEX ; } // ------------------------------------------------------------------------- // analyze, factorize, and get the results // ------------------------------------------------------------------------- if (opts.haveB) { // --------------------------------------------------------------------- // get B, and convert to complex if necessary // --------------------------------------------------------------------- if (!mxIsNumeric (pargin [1])) { mexErrMsgIdAndTxt ("QR:invalidInput", "invalid non-numeric B") ; } if (mxGetM (pargin [1]) != m) { mexErrMsgIdAndTxt ("QR:invalidInput", "A and B must have the same number of rows") ; } // convert from real or zomplex to complex Bx = spqr_mx_merge_if_complex (pargin [1], is_complex, &bnz, cc) ; int B_is_sparse = mxIsSparse (pargin [1]) ; if (B_is_sparse) { Bsparse = spqr_mx_get_sparse (pargin [1], &Bsmatrix, &dummy) ; Bdense = NULL ; if (is_complex) { // Bsparse has been converted from real or zomplex to complex Bsparse->x = Bx ; Bsparse->z = NULL ; Bsparse->xtype = CHOLMOD_COMPLEX ; } } else { Bsparse = NULL ; Bdense = spqr_mx_get_dense (pargin [1], &Bdmatrix, &dummy) ; if (is_complex) { // Bdense has been converted from real or zomplex to complex Bdense->x = Bx ; Bdense->z = NULL ; Bdense->xtype = CHOLMOD_COMPLEX ; } } // --------------------------------------------------------------------- // [C,R,E] = qr (A,B,...) or [C,R] = qr (A,B,...) // --------------------------------------------------------------------- if (is_complex) { // ----------------------------------------------------------------- // [C,R,E] = qr (A,B): complex case // ----------------------------------------------------------------- if (B_is_sparse) { // B and C are both sparse and complex SuiteSparseQR <Complex> (order, tol, econ, A, Bsparse, &Csparse, &R, &E, cc) ; pargout [0] = spqr_mx_put_sparse (&Csparse, cc) ; } else { // B and C are both dense and complex SuiteSparseQR <Complex> (order, tol, econ, A, Bdense, &Cdense, &R, &E, cc) ; pargout [0] = spqr_mx_put_dense (&Cdense, cc) ; } } else { // ----------------------------------------------------------------- // [C,R,E] = qr (A,B): real case // ----------------------------------------------------------------- if (B_is_sparse) { // B and C are both sparse and real SuiteSparseQR <double> (order, tol, econ, A, Bsparse, &Csparse, &R, &E, cc) ; pargout [0] = spqr_mx_put_sparse (&Csparse, cc) ; } else { // B and C are both dense and real SuiteSparseQR <double> (order, tol, econ, A, Bdense, &Cdense, &R, &E, cc) ; pargout [0] = spqr_mx_put_dense (&Cdense, cc) ; } } pargout [1] = spqr_mx_put_sparse (&R, cc) ; } else if (nargout <= 1) { // --------------------------------------------------------------------- // R = qr (A) or R = qr (A,opts) // --------------------------------------------------------------------- if (is_complex) { SuiteSparseQR <Complex> (0, tol, econ, A, &R, NULL, cc) ; } else { SuiteSparseQR <double> (0, tol, econ, A, &R, NULL, cc) ; } pargout [0] = spqr_mx_put_sparse (&R, cc) ; } else { // --------------------------------------------------------------------- // [Q,R,E] = qr (A,...) or [Q,R] = qr (A,...) // --------------------------------------------------------------------- if (opts.Qformat == SPQR_Q_DISCARD) { // ----------------------------------------------------------------- // Q is discarded, and Q = [ ] is returned as a placeholder // ----------------------------------------------------------------- if (is_complex) { SuiteSparseQR <Complex> (order, tol, econ, A, &R, &E, cc); } else { SuiteSparseQR <double> (order, tol, econ, A, &R, &E, cc) ; } pargout [0] = mxCreateDoubleMatrix (0, 0, mxREAL) ; } else if (opts.Qformat == SPQR_Q_MATRIX) { // ----------------------------------------------------------------- // Q is a sparse matrix // ----------------------------------------------------------------- if (is_complex) { SuiteSparseQR <Complex> (order, tol, econ, A, &Q, &R, &E, cc) ; } else { SuiteSparseQR <double> (order, tol, econ, A, &Q, &R, &E, cc) ; } pargout [0] = spqr_mx_put_sparse (&Q, cc) ; } else { // ----------------------------------------------------------------- // H is kept, and Q is a struct containing H, Tau, and P // ----------------------------------------------------------------- mxArray *Tau, *P, *Hmatlab ; if (is_complex) { SuiteSparseQR <Complex> (order, tol, econ, A, &R, &E, &H, &HPinv, &HTau, cc) ; } else { SuiteSparseQR <double> (order, tol, econ, A, &R, &E, &H, &HPinv, &HTau, cc) ; } Tau = spqr_mx_put_dense (&HTau, cc) ; Hmatlab = spqr_mx_put_sparse (&H, cc) ; // Q.P contains the inverse row permutation P = mxCreateDoubleMatrix (1, m, mxREAL) ; double *Tx = mxGetPr (P) ; for (Int i = 0 ; i < m ; i++) { Tx [i] = HPinv [i] + 1 ; } // return Q const char *Qstruct [ ] = { "H", "Tau", "P" } ; pargout [0] = mxCreateStructMatrix (1, 1, 3, Qstruct) ; mxSetFieldByNumber (pargout [0], 0, 0, Hmatlab) ; mxSetFieldByNumber (pargout [0], 0, 1, Tau) ; mxSetFieldByNumber (pargout [0], 0, 2, P) ; } pargout [1] = spqr_mx_put_sparse (&R, cc) ; } // ------------------------------------------------------------------------- // return E // ------------------------------------------------------------------------- if (nargout > 2) { pargout [2] = spqr_mx_put_permutation (E, n, opts.permvector, cc) ; } // ------------------------------------------------------------------------- // free copy of merged-complex, if needed // ------------------------------------------------------------------------- if (is_complex) { // this was allocated by merge_if_complex cholmod_l_free (anz, sizeof (Complex), Ax, cc) ; if (opts.haveB) { cholmod_l_free (bnz, sizeof (Complex), Bx, cc) ; } } // ------------------------------------------------------------------------- // info output // ------------------------------------------------------------------------- if (nargout > 3) { #ifdef TIMING double flops = cc->other1 [0] ; double t = spqr_time ( ) - t0 ; #else double flops = -1 ; double t = -1 ; #endif pargout [3] = spqr_mx_info (cc, t, flops) ; } cholmod_l_finish (cc) ; if (opts.spumoni > 0) spqr_mx_spumoni (&opts, is_complex, cc) ; }
void mexFunction ( int nargout, mxArray *pargout [ ], int nargin, const mxArray *pargin [ ] ) { Int *Bp, *Bi ; double *Ax, *Bx, dummy ; Int m, n, k, bncols, p, i, rank, A_complex, B_complex, is_complex, anz, bnz ; spqr_mx_options opts ; cholmod_sparse *A, Amatrix, *Xsparse ; cholmod_dense *Xdense ; cholmod_common Common, *cc ; char msg [LEN+1] ; #ifdef TIMING double t0 = (nargout > 1) ? spqr_time ( ) : 0 ; #endif // ------------------------------------------------------------------------- // start CHOLMOD and set parameters // ------------------------------------------------------------------------- cc = &Common ; cholmod_l_start (cc) ; spqr_mx_config (SPUMONI, cc) ; // ------------------------------------------------------------------------- // check inputs // ------------------------------------------------------------------------- if (nargout > 2) { mexErrMsgIdAndTxt ("MATLAB:maxlhs", "Too many output arguments") ; } if (nargin < 2) { mexErrMsgIdAndTxt ("MATLAB:minrhs", "Not enough input arguments") ; } if (nargin > 3) { mexErrMsgIdAndTxt ("MATLAB:maxrhs", "Too many input arguments") ; } // ------------------------------------------------------------------------- // get the input matrix A (must be sparse) // ------------------------------------------------------------------------- if (!mxIsSparse (pargin [0])) { mexErrMsgIdAndTxt ("QR:invalidInput", "A must be sparse") ; } A = spqr_mx_get_sparse (pargin [0], &Amatrix, &dummy) ; m = A->nrow ; n = A->ncol ; A_complex = mxIsComplex (pargin [0]) ; B_complex = mxIsComplex (pargin [1]) ; is_complex = (A_complex || B_complex) ; Ax = spqr_mx_merge_if_complex (pargin [0], is_complex, &anz, cc) ; if (is_complex) { // A has been converted from real or zomplex to complex A->x = Ax ; A->z = NULL ; A->xtype = CHOLMOD_COMPLEX ; } // ------------------------------------------------------------------------- // determine usage and parameters // ------------------------------------------------------------------------- spqr_mx_get_options ((nargin < 3) ? NULL : pargin [2], &opts, m, 3, cc) ; opts.Qformat = SPQR_Q_DISCARD ; opts.econ = 0 ; opts.permvector = TRUE ; opts.haveB = TRUE ; // ------------------------------------------------------------------------- // get the input matrix B (sparse or dense) // ------------------------------------------------------------------------- if (!mxIsNumeric (pargin [1])) { mexErrMsgIdAndTxt ("QR:invalidInput", "invalid non-numeric B") ; } if (mxGetM (pargin [1]) != m) { mexErrMsgIdAndTxt ("QR:invalidInput", "A and B must have the same number of rows") ; } cholmod_sparse Bsmatrix, *Bsparse ; cholmod_dense Bdmatrix, *Bdense ; // convert from real or zomplex to complex Bx = spqr_mx_merge_if_complex (pargin [1], is_complex, &bnz, cc) ; int B_is_sparse = mxIsSparse (pargin [1]) ; if (B_is_sparse) { Bsparse = spqr_mx_get_sparse (pargin [1], &Bsmatrix, &dummy) ; Bdense = NULL ; if (is_complex) { // Bsparse has been converted from real or zomplex to complex Bsparse->x = Bx ; Bsparse->z = NULL ; Bsparse->xtype = CHOLMOD_COMPLEX ; } } else { Bsparse = NULL ; Bdense = spqr_mx_get_dense (pargin [1], &Bdmatrix, &dummy) ; if (is_complex) { // Bdense has been converted from real or zomplex to complex Bdense->x = Bx ; Bdense->z = NULL ; Bdense->xtype = CHOLMOD_COMPLEX ; } } // ------------------------------------------------------------------------- // X = A\B // ------------------------------------------------------------------------- if (opts.min2norm && m < n) { #ifndef NEXPERT // This requires SuiteSparseQR_expert.cpp if (is_complex) { if (B_is_sparse) { // X and B are both sparse and complex Xsparse = SuiteSparseQR_min2norm <Complex> (opts.ordering, opts.tol, A, Bsparse, cc) ; pargout [0] = spqr_mx_put_sparse (&Xsparse, cc) ; } else { // X and B are both dense and complex Xdense = SuiteSparseQR_min2norm <Complex> (opts.ordering, opts.tol, A, Bdense, cc) ; pargout [0] = spqr_mx_put_dense (&Xdense, cc) ; } } else { if (B_is_sparse) { // X and B are both sparse and real Xsparse = SuiteSparseQR_min2norm <double> (opts.ordering, opts.tol, A, Bsparse, cc) ; pargout [0] = spqr_mx_put_sparse (&Xsparse, cc) ; } else { // X and B are both dense and real Xdense = SuiteSparseQR_min2norm <double> (opts.ordering, opts.tol, A, Bdense, cc) ; pargout [0] = spqr_mx_put_dense (&Xdense, cc) ; } } #else mexErrMsgIdAndTxt ("QR:notInstalled", "min2norm method not installed") ; #endif } else { if (is_complex) { if (B_is_sparse) { // X and B are both sparse and complex Xsparse = SuiteSparseQR <Complex> (opts.ordering, opts.tol, A, Bsparse, cc) ; pargout [0] = spqr_mx_put_sparse (&Xsparse, cc) ; } else { // X and B are both dense and complex Xdense = SuiteSparseQR <Complex> (opts.ordering, opts.tol, A, Bdense, cc) ; pargout [0] = spqr_mx_put_dense (&Xdense, cc) ; } } else { if (B_is_sparse) { // X and B are both sparse and real Xsparse = SuiteSparseQR <double> (opts.ordering, opts.tol, A, Bsparse, cc) ; pargout [0] = spqr_mx_put_sparse (&Xsparse, cc) ; } else { // X and B are both dense and real Xdense = SuiteSparseQR <double> (opts.ordering, opts.tol, A, Bdense, cc) ; pargout [0] = spqr_mx_put_dense (&Xdense, cc) ; } } } // ------------------------------------------------------------------------- // info output // ------------------------------------------------------------------------- if (nargout > 1) { #ifdef TIMING double flops = cc->other1 [0] ; double t = spqr_time ( ) - t0 ; #else double flops = -1 ; double t = -1 ; #endif pargout [1] = spqr_mx_info (cc, t, flops) ; } // ------------------------------------------------------------------------- // warn if rank deficient // ------------------------------------------------------------------------- rank = cc->SPQR_istat [4] ; if (rank < MIN (m,n)) { // snprintf would be safer, but Windows is oblivious to safety ... // (Visual Studio C++ 2008 does not recognize snprintf!) sprintf (msg, "rank deficient. rank = %ld tol = %g\n", rank, cc->SPQR_xstat [1]) ; mexWarnMsgIdAndTxt ("MATLAB:rankDeficientMatrix", msg) ; } if (is_complex) { // free the merged complex copies of A and B cholmod_l_free (anz, sizeof (Complex), Ax, cc) ; cholmod_l_free (bnz, sizeof (Complex), Bx, cc) ; } cholmod_l_finish (cc) ; if (opts.spumoni > 0) spqr_mx_spumoni (&opts, is_complex, cc) ; }
void mexFunction ( int nargout, mxArray *pargout [ ], int nargin, const mxArray *pargin [ ] ) { mxArray *Hmatlab, *Tau, *P ; Long *HPinv, *Yp, *Yi ; double *Hx, *Xx, *Tx, *Px, dummy ; Long m, n, k, nh, nb, p, i, method, mh, gotP, X_is_sparse, is_complex, hnz, tnz, xnz, inuse, count ; cholmod_sparse *Ysparse, *H, Hmatrix, *Xsparse, Xsmatrix ; cholmod_dense *Ydense, *Xdense, Xdmatrix, *HTau, HTau_matrix ; cholmod_common Common, *cc ; // ------------------------------------------------------------------------- // start CHOLMOD and set parameters // ------------------------------------------------------------------------- cc = &Common ; cholmod_l_start (cc) ; spqr_mx_config (SPUMONI, cc) ; // ------------------------------------------------------------------------- // check inputs // ------------------------------------------------------------------------- // nargin can be 2 or 3 // nargout can be 0 or 1 if (nargout > 1) { mexErrMsgIdAndTxt ("MATLAB:maxlhs", "Too many output arguments") ; } if (nargin < 2) { mexErrMsgIdAndTxt ("MATLAB:minrhs", "Not enough input arguments") ; } if (nargin > 3) { mexErrMsgIdAndTxt ("MATLAB:maxrhs", "Too many input arguments") ; } if (!mxIsStruct (pargin [0])) { mexErrMsgIdAndTxt ("QR:invalidInput", "invalid Q (must be a struct)") ; } // ------------------------------------------------------------------------- // get H, Tau, and P from the Q struct // ------------------------------------------------------------------------- i = mxGetFieldNumber (pargin [0], "H") ; if (i < 0) { mexErrMsgIdAndTxt ("QR:invalidInput", "invalid Q struct") ; } Hmatlab = mxGetFieldByNumber (pargin [0], 0, i) ; nh = mxGetN (Hmatlab) ; if (!mxIsSparse (Hmatlab)) { mexErrMsgIdAndTxt ("QR:invalidInput", "H must be sparse") ; } i = mxGetFieldNumber (pargin [0], "Tau") ; if (i < 0) { mexErrMsgIdAndTxt ("QR:invalidInput", "invalid Q struct") ; } Tau = mxGetFieldByNumber (pargin [0], 0, i) ; if (nh != mxGetNumberOfElements (Tau)) { mexErrMsgIdAndTxt ("QR:invalidInput", "H and Tau must have the same number of columns") ; } is_complex = mxIsComplex (Tau) || mxIsComplex (Hmatlab) || mxIsComplex (pargin [1]) ; // ------------------------------------------------------------------------- // get the Householder vectors // ------------------------------------------------------------------------- H = spqr_mx_get_sparse (Hmatlab, &Hmatrix, &dummy) ; mh = H->nrow ; Hx = spqr_mx_merge_if_complex (Hmatlab, is_complex, &hnz, cc) ; if (is_complex) { // H has been converted from real or zomplex to complex H->x = Hx ; H->z = NULL ; H->xtype = CHOLMOD_COMPLEX ; } // ------------------------------------------------------------------------- // get Tau // ------------------------------------------------------------------------- HTau = spqr_mx_get_dense (Tau, &HTau_matrix, &dummy) ; Tx = spqr_mx_merge_if_complex (Tau, is_complex, &tnz, cc) ; if (is_complex) { // HTau has been converted from real or zomplex to complex HTau->x = Tx ; HTau->z = NULL ; HTau->xtype = CHOLMOD_COMPLEX ; } // ------------------------------------------------------------------------- // get method // ------------------------------------------------------------------------- if (nargin < 3) { method = 0 ; } else { method = (Long) mxGetScalar (pargin [2]) ; if (method < 0 || method > 3) { mexErrMsgIdAndTxt ("QR:invalidInput", "invalid method") ; } } // ------------------------------------------------------------------------- // get X // ------------------------------------------------------------------------- m = mxGetM (pargin [1]) ; n = mxGetN (pargin [1]) ; X_is_sparse = mxIsSparse (pargin [1]) ; Xsparse = NULL ; if (X_is_sparse) { Xsparse = spqr_mx_get_sparse (pargin [1], &Xsmatrix, &dummy) ; } else { Xdense = spqr_mx_get_dense (pargin [1], &Xdmatrix, &dummy) ; } Xx = spqr_mx_merge_if_complex (pargin [1], is_complex, &xnz, cc) ; if (is_complex) { // X has been converted from real or zomplex to complex if (X_is_sparse) { Xsparse->x = Xx ; Xsparse->xtype = CHOLMOD_COMPLEX ; } else { Xdense->x = Xx ; Xdense->xtype = CHOLMOD_COMPLEX ; } } if (method == 0 || method == 1) { if (mh != m) { mexErrMsgIdAndTxt ("QR:invalidInput", "H and X must have same number of rows") ; } } else { if (mh != n) { mexErrMsgIdAndTxt ("QR:invalidInput", "# of cols of X must equal # of rows of H") ; } } // ------------------------------------------------------------------------- // get P // ------------------------------------------------------------------------- i = mxGetFieldNumber (pargin [0], "P") ; gotP = (i >= 0) ; HPinv = NULL ; if (gotP) { // get P from the H struct P = mxGetFieldByNumber (pargin [0], 0, i) ; if (mxGetNumberOfElements (P) != mh) { mexErrMsgIdAndTxt ("QR:invalidInput", "P must be a vector of length equal to # rows of H") ; } HPinv = (Long *) cholmod_l_malloc (mh, sizeof (Long), cc) ; Px = mxGetPr (P) ; for (i = 0 ; i < mh ; i++) { HPinv [i] = (Long) (Px [i] - 1) ; if (HPinv [i] < 0 || HPinv [i] >= mh) { mexErrMsgIdAndTxt ("QR:invalidInput", "invalid permutation") ; } } } // ------------------------------------------------------------------------- // Y = Q'*X, Q*X, X*Q or X*Q' // ------------------------------------------------------------------------- if (is_complex) { if (X_is_sparse) { Ysparse = SuiteSparseQR_qmult <Complex> (method, H, HTau, HPinv, Xsparse, cc) ; pargout [0] = spqr_mx_put_sparse (&Ysparse, cc) ; } else { Ydense = SuiteSparseQR_qmult <Complex> (method, H, HTau, HPinv, Xdense, cc) ; pargout [0] = spqr_mx_put_dense (&Ydense, cc) ; } } else { if (X_is_sparse) { Ysparse = SuiteSparseQR_qmult <double> (method, H, HTau, HPinv, Xsparse, cc) ; pargout [0] = spqr_mx_put_sparse (&Ysparse, cc) ; } else { Ydense = SuiteSparseQR_qmult <double> (method, H, HTau, HPinv, Xdense, cc) ; pargout [0] = spqr_mx_put_dense (&Ydense, cc) ; } } // ------------------------------------------------------------------------- // free workspace // ------------------------------------------------------------------------- cholmod_l_free (mh, sizeof (Long), HPinv, cc) ; if (is_complex) { // free the merged copies of the real parts of the H and Tau matrices cholmod_l_free (hnz, sizeof (Complex), Hx, cc) ; cholmod_l_free (tnz, sizeof (Complex), Tx, cc) ; cholmod_l_free (xnz, sizeof (Complex), Xx, cc) ; } cholmod_l_finish (cc) ; #if 0 // malloc count for testing only ... spqr_mx_get_usage (pargout [0], 1, &inuse, &count, cc) ; if (inuse != cc->memory_inuse || count != cc->malloc_count) { mexErrMsgIdAndTxt ("QR:internalError", "memory leak!") ; } #endif }