void mexFunction ( int nargout, mxArray *pargout [ ], int nargin, const mxArray *pargin [ ] ) { UF_long i, m, n, *Ap, *Ai, *P, nc, result, spumoni, full ; double *Pout, *InfoOut, Control [AMD_CONTROL], Info [AMD_INFO], *ControlIn ; mxArray *A ; /* --------------------------------------------------------------------- */ /* get control parameters */ /* --------------------------------------------------------------------- */ amd_malloc = mxMalloc ; amd_free = mxFree ; amd_calloc = mxCalloc ; amd_realloc = mxRealloc ; amd_printf = mexPrintf ; spumoni = 0 ; if (nargin == 0) { /* get the default control parameters, and return */ pargout [0] = mxCreateDoubleMatrix (AMD_CONTROL, 1, mxREAL) ; amd_l_defaults (mxGetPr (pargout [0])) ; if (nargout == 0) { amd_l_control (mxGetPr (pargout [0])) ; } return ; } amd_l_defaults (Control) ; if (nargin > 1) { ControlIn = mxGetPr (pargin [1]) ; nc = mxGetM (pargin [1]) * mxGetN (pargin [1]) ; Control [AMD_DENSE] = (nc > 0) ? ControlIn [AMD_DENSE] : AMD_DEFAULT_DENSE ; Control [AMD_AGGRESSIVE] = (nc > 1) ? ControlIn [AMD_AGGRESSIVE] : AMD_DEFAULT_AGGRESSIVE ; spumoni = (nc > 2) ? (ControlIn [2] != 0) : 0 ; } if (spumoni > 0) { amd_l_control (Control) ; } /* --------------------------------------------------------------------- */ /* get inputs */ /* --------------------------------------------------------------------- */ if (nargout > 2 || nargin > 2) { mexErrMsgTxt ("Usage: p = amd (A)\nor [p, Info] = amd (A, Control)") ; } A = (mxArray *) pargin [0] ; n = mxGetN (A) ; m = mxGetM (A) ; if (spumoni > 0) { mexPrintf (" input matrix A is %d-by-%d\n", m, n) ; } if (mxGetNumberOfDimensions (A) != 2) { mexErrMsgTxt ("amd: A must be 2-dimensional") ; } if (m != n) { mexErrMsgTxt ("amd: A must be square") ; } /* --------------------------------------------------------------------- */ /* allocate workspace for output permutation */ /* --------------------------------------------------------------------- */ P = mxMalloc ((n+1) * sizeof (UF_long)) ; /* --------------------------------------------------------------------- */ /* if A is full, convert to a sparse matrix */ /* --------------------------------------------------------------------- */ full = !mxIsSparse (A) ; if (full) { if (spumoni > 0) { mexPrintf ( " input matrix A is full (sparse copy of A will be created)\n"); } mexCallMATLAB (1, &A, 1, (mxArray **) pargin, "sparse") ; } Ap = (UF_long *) mxGetJc (A) ; Ai = (UF_long *) mxGetIr (A) ; if (spumoni > 0) { mexPrintf (" input matrix A has %d nonzero entries\n", Ap [n]) ; } /* --------------------------------------------------------------------- */ /* order the matrix */ /* --------------------------------------------------------------------- */ result = amd_l_order (n, Ap, Ai, P, Control, Info) ; /* --------------------------------------------------------------------- */ /* if A is full, free the sparse copy of A */ /* --------------------------------------------------------------------- */ if (full) { mxDestroyArray (A) ; } /* --------------------------------------------------------------------- */ /* print results (including return value) */ /* --------------------------------------------------------------------- */ if (spumoni > 0) { amd_l_info (Info) ; } /* --------------------------------------------------------------------- */ /* check error conditions */ /* --------------------------------------------------------------------- */ if (result == AMD_OUT_OF_MEMORY) { mexErrMsgTxt ("amd: out of memory") ; } else if (result == AMD_INVALID) { mexErrMsgTxt ("amd: input matrix A is corrupted") ; } /* --------------------------------------------------------------------- */ /* copy the outputs to MATLAB */ /* --------------------------------------------------------------------- */ /* output permutation, P */ pargout [0] = mxCreateDoubleMatrix (1, n, mxREAL) ; Pout = mxGetPr (pargout [0]) ; for (i = 0 ; i < n ; i++) { Pout [i] = P [i] + 1 ; /* change to 1-based indexing for MATLAB */ } mxFree (P) ; /* Info */ if (nargout > 1) { pargout [1] = mxCreateDoubleMatrix (AMD_INFO, 1, mxREAL) ; InfoOut = mxGetPr (pargout [1]) ; for (i = 0 ; i < AMD_INFO ; i++) { InfoOut [i] = Info [i] ; } } }
int main (void) { /* The symmetric can_24 Harwell/Boeing matrix, including upper and lower * triangular parts, and the diagonal entries. Note that this matrix is * 0-based, with row and column indices in the range 0 to n-1. */ Long n = 24, nz, Ap [ ] = { 0, 9, 15, 21, 27, 33, 39, 48, 57, 61, 70, 76, 82, 88, 94, 100, 106, 110, 119, 128, 137, 143, 152, 156, 160 }, Ai [ ] = { /* column 0: */ 0, 5, 6, 12, 13, 17, 18, 19, 21, /* column 1: */ 1, 8, 9, 13, 14, 17, /* column 2: */ 2, 6, 11, 20, 21, 22, /* column 3: */ 3, 7, 10, 15, 18, 19, /* column 4: */ 4, 7, 9, 14, 15, 16, /* column 5: */ 0, 5, 6, 12, 13, 17, /* column 6: */ 0, 2, 5, 6, 11, 12, 19, 21, 23, /* column 7: */ 3, 4, 7, 9, 14, 15, 16, 17, 18, /* column 8: */ 1, 8, 9, 14, /* column 9: */ 1, 4, 7, 8, 9, 13, 14, 17, 18, /* column 10: */ 3, 10, 18, 19, 20, 21, /* column 11: */ 2, 6, 11, 12, 21, 23, /* column 12: */ 0, 5, 6, 11, 12, 23, /* column 13: */ 0, 1, 5, 9, 13, 17, /* column 14: */ 1, 4, 7, 8, 9, 14, /* column 15: */ 3, 4, 7, 15, 16, 18, /* column 16: */ 4, 7, 15, 16, /* column 17: */ 0, 1, 5, 7, 9, 13, 17, 18, 19, /* column 18: */ 0, 3, 7, 9, 10, 15, 17, 18, 19, /* column 19: */ 0, 3, 6, 10, 17, 18, 19, 20, 21, /* column 20: */ 2, 10, 19, 20, 21, 22, /* column 21: */ 0, 2, 6, 10, 11, 19, 20, 21, 22, /* column 22: */ 2, 20, 21, 22, /* column 23: */ 6, 11, 12, 23 } ; Long P [24], Pinv [24], i, j, k, jnew, p, inew, result ; double Control [AMD_CONTROL], Info [AMD_INFO] ; char A [24][24] ; /* here is an example of how to use AMD_VERSION. This code will work in * any version of AMD. */ #if defined(AMD_VERSION) && (AMD_VERSION >= AMD_VERSION_CODE(1,2)) printf ("AMD version %d.%d.%d, date: %s\n", AMD_MAIN_VERSION, AMD_SUB_VERSION, AMD_SUBSUB_VERSION, AMD_DATE) ; #else printf ("AMD version: 1.1 or earlier\n") ; #endif printf ("AMD demo, with the 24-by-24 Harwell/Boeing matrix, can_24:\n") ; /* get the default parameters, and print them */ amd_l_defaults (Control) ; amd_l_control (Control) ; /* print the input matrix */ nz = Ap [n] ; printf ("\nInput matrix: %ld-by-%ld, with %ld entries.\n" " Note that for a symmetric matrix such as this one, only the\n" " strictly lower or upper triangular parts would need to be\n" " passed to AMD, since AMD computes the ordering of A+A'. The\n" " diagonal entries are also not needed, since AMD ignores them.\n" , n, n, nz) ; for (j = 0 ; j < n ; j++) { printf ("\nColumn: %ld, number of entries: %ld, with row indices in" " Ai [%ld ... %ld]:\n row indices:", j, Ap [j+1] - Ap [j], Ap [j], Ap [j+1]-1) ; for (p = Ap [j] ; p < Ap [j+1] ; p++) { i = Ai [p] ; printf (" %ld", i) ; } printf ("\n") ; } /* print a character plot of the input matrix. This is only reasonable * because the matrix is small. */ printf ("\nPlot of input matrix pattern:\n") ; for (j = 0 ; j < n ; j++) { for (i = 0 ; i < n ; i++) A [i][j] = '.' ; for (p = Ap [j] ; p < Ap [j+1] ; p++) { i = Ai [p] ; A [i][j] = 'X' ; } } printf (" ") ; for (j = 0 ; j < n ; j++) printf (" %1ld", j % 10) ; printf ("\n") ; for (i = 0 ; i < n ; i++) { printf ("%2ld: ", i) ; for (j = 0 ; j < n ; j++) { printf (" %c", A [i][j]) ; } printf ("\n") ; } /* order the matrix */ result = amd_l_order (n, Ap, Ai, P, Control, Info) ; printf ("return value from amd_l_order: %ld (should be %d)\n", result, AMD_OK) ; /* print the statistics */ amd_l_info (Info) ; if (result != AMD_OK) { printf ("AMD failed\n") ; exit (1) ; } /* print the permutation vector, P, and compute the inverse permutation */ printf ("Permutation vector:\n") ; for (k = 0 ; k < n ; k++) { /* row/column j is the kth row/column in the permuted matrix */ j = P [k] ; Pinv [j] = k ; printf (" %2ld", j) ; } printf ("\n\n") ; printf ("Inverse permutation vector:\n") ; for (j = 0 ; j < n ; j++) { k = Pinv [j] ; printf (" %2ld", k) ; } printf ("\n\n") ; /* print a character plot of the permuted matrix. */ printf ("\nPlot of permuted matrix pattern:\n") ; for (jnew = 0 ; jnew < n ; jnew++) { j = P [jnew] ; for (inew = 0 ; inew < n ; inew++) A [inew][jnew] = '.' ; for (p = Ap [j] ; p < Ap [j+1] ; p++) { inew = Pinv [Ai [p]] ; A [inew][jnew] = 'X' ; } } printf (" ") ; for (j = 0 ; j < n ; j++) printf (" %1ld", j % 10) ; printf ("\n") ; for (i = 0 ; i < n ; i++) { printf ("%2ld: ", i) ; for (j = 0 ; j < n ; j++) { printf (" %c", A [i][j]) ; } printf ("\n") ; } return (0) ; }