/* -------------------------------------------------------------------- inputComplex a number of (row,column, entry) triples into the matrix created -- 98jan28, cca -------------------------------------------------------------------- */ static void inputTriples ( InpMtx *inpmtx, int ntriples, int rowids[], int colids[], double entries[] ) { int nent ; int *ivec1, *ivec2 ; prepareToAddNewEntries(inpmtx, ntriples) ; nent = inpmtx->nent ; ivec1 = IV_entries(&inpmtx->ivec1IV) ; ivec2 = IV_entries(&inpmtx->ivec2IV) ; IVcopy(ntriples, ivec1 + nent, rowids) ; IVcopy(ntriples, ivec2 + nent, colids) ; IV_setSize(&inpmtx->ivec1IV, nent + ntriples) ; IV_setSize(&inpmtx->ivec2IV, nent + ntriples) ; if ( INPMTX_IS_REAL_ENTRIES(inpmtx) ) { double *dvec = DV_entries(&inpmtx->dvecDV) ; DVcopy(ntriples, dvec + nent, entries) ; DV_setSize(&inpmtx->dvecDV, nent + ntriples) ; } else if ( INPMTX_IS_COMPLEX_ENTRIES(inpmtx) ) { double *dvec = DV_entries(&inpmtx->dvecDV) ; ZVcopy(ntriples, dvec + 2*nent, entries) ; DV_setSize(&inpmtx->dvecDV, 2*(nent + ntriples)) ; } inpmtx->nent += ntriples ; inpmtx->storageMode = INPMTX_RAW_DATA ; return ; }
/* ----------------------------- input a chevron in the matrix created -- 98jan28, cca ----------------------------- */ static void inputChevron ( InpMtx *inpmtx, int chv, int chvsize, int chvind[], double chvent[] ) { int col, ii, jj, nent, offset, row ; int *ivec1, *ivec2 ; prepareToAddNewEntries(inpmtx, chvsize) ; nent = inpmtx->nent ; ivec1 = IV_entries(&inpmtx->ivec1IV) ; ivec2 = IV_entries(&inpmtx->ivec2IV) ; if ( INPMTX_IS_BY_ROWS(inpmtx) ) { for ( ii = 0, jj = nent ; ii < chvsize ; ii++, jj++ ) { if ( (offset = chvind[ii]) >= 0 ) { row = chv ; col = chv + offset ; } else { col = chv ; row = chv - offset ; } ivec1[jj] = row ; ivec2[jj] = col ; } } else if ( INPMTX_IS_BY_COLUMNS(inpmtx) ) { for ( ii = 0, jj = nent ; ii < chvsize ; ii++, jj++ ) { if ( (offset = chvind[ii]) >= 0 ) { row = chv ; col = chv + offset ; } else { col = chv ; row = chv - offset ; } ivec1[jj] = col ; ivec2[jj] = row ; } } else if ( INPMTX_IS_BY_CHEVRONS(inpmtx) ) { IVfill(chvsize, ivec1 + nent, chv) ; IVcopy(chvsize, ivec2 + nent, chvind) ; } IV_setSize(&inpmtx->ivec1IV, nent + chvsize) ; IV_setSize(&inpmtx->ivec2IV, nent + chvsize) ; if ( INPMTX_IS_REAL_ENTRIES(inpmtx) ) { double *dvec = DV_entries(&inpmtx->dvecDV) + nent ; DVcopy(chvsize, dvec, chvent) ; DV_setSize(&inpmtx->dvecDV, nent + chvsize) ; } else if ( INPMTX_IS_REAL_ENTRIES(inpmtx) ) { double *dvec = DV_entries(&inpmtx->dvecDV) + 2*nent ; ZVcopy(chvsize, dvec, chvent) ; DV_setSize(&inpmtx->dvecDV, 2*(nent + chvsize)) ; } inpmtx->nent += chvsize ; inpmtx->storageMode = INPMTX_RAW_DATA ; return ; }
/* --------------------------------- input a row in the matrix created -- 98jan28, cca --------------------------------- */ static void inputRow ( InpMtx *inpmtx, int row, int rowsize, int rowind[], double rowent[] ) { int col, ii, jj, nent ; int *ivec1, *ivec2 ; prepareToAddNewEntries(inpmtx, rowsize) ; nent = inpmtx->nent ; ivec1 = IV_entries(&inpmtx->ivec1IV) ; ivec2 = IV_entries(&inpmtx->ivec2IV) ; if ( INPMTX_IS_BY_ROWS(inpmtx) ) { /* row coordinates */ IVfill(rowsize, ivec1 + nent, row) ; IVcopy(rowsize, ivec2 + nent, rowind) ; } else if ( INPMTX_IS_BY_COLUMNS(inpmtx) ) { /* column coordinates */ IVfill(rowsize, ivec2 + nent, row) ; IVcopy(rowsize, ivec1 + nent, rowind) ; } else if ( INPMTX_IS_BY_CHEVRONS(inpmtx) ) { /* chevron coordinates */ for ( ii = 0, jj = nent ; ii < rowsize ; ii++, jj++ ) { col = rowind[ii] ; ivec1[ii] = (row <= col) ? row : col ; ivec2[ii] = col - row ; } } IV_setSize(&inpmtx->ivec1IV, nent + rowsize) ; IV_setSize(&inpmtx->ivec2IV, nent + rowsize) ; /* ----------------- input the entries ----------------- */ if ( INPMTX_IS_REAL_ENTRIES(inpmtx) ) { double *dvec = DV_entries(&inpmtx->dvecDV) ; DVcopy(rowsize, dvec + nent, rowent) ; DV_setSize(&inpmtx->dvecDV, nent + rowsize) ; } else if ( INPMTX_IS_COMPLEX_ENTRIES(inpmtx) ) { double *dvec = DV_entries(&inpmtx->dvecDV) ; ZVcopy(rowsize, dvec + 2*nent, rowent) ; DV_setSize(&inpmtx->dvecDV, 2*(nent + rowsize)) ; } inpmtx->storageMode = INPMTX_RAW_DATA ; inpmtx->nent += rowsize ; return ; }
/* ------------------------------------ input a complex column in the matrix created -- 98jan28, cca ------------------------------------ */ static void inputColumn ( InpMtx *inpmtx, int col, int colsize, int colind[], double colent[] ) { int ii, jj, nent, row ; int *ivec1, *ivec2 ; prepareToAddNewEntries(inpmtx, colsize) ; nent = inpmtx->nent ; ivec1 = IV_entries(&inpmtx->ivec1IV) ; ivec2 = IV_entries(&inpmtx->ivec2IV) ; if ( INPMTX_IS_BY_ROWS(inpmtx) ) { IVcopy(colsize, ivec1 + nent, colind) ; IVfill(colsize, ivec2 + nent, col) ; } else if ( INPMTX_IS_BY_COLUMNS(inpmtx) ) { IVfill(colsize, ivec1 + nent, col) ; IVcopy(colsize, ivec2 + nent, colind) ; } else if ( INPMTX_IS_BY_CHEVRONS(inpmtx) ) { for ( ii = 0, jj = nent ; ii < colsize ; ii++, jj++ ) { row = colind[jj] ; ivec1[jj] = (row <= col) ? row : col ; ivec2[jj] = col - row ; } } IV_setSize(&inpmtx->ivec1IV, nent + colsize) ; IV_setSize(&inpmtx->ivec2IV, nent + colsize) ; if ( INPMTX_IS_REAL_ENTRIES(inpmtx) ) { double *dvec = DV_entries(&inpmtx->dvecDV) + nent ; DVcopy(colsize, dvec, colent) ; DV_setSize(&inpmtx->dvecDV, nent + colsize) ; } else if ( INPMTX_IS_COMPLEX_ENTRIES(inpmtx) ) { double *dvec = DV_entries(&inpmtx->dvecDV) + 2*nent ; ZVcopy(colsize, dvec, colent) ; DV_setSize(&inpmtx->dvecDV, 2*(nent + colsize)) ; } inpmtx->nent = nent + colsize ; inpmtx->storageMode = INPMTX_RAW_DATA ; return ; }
/* ---------------------------------------------- sort the rows of the matrix in ascending order of the rowids[] vector. on return, rowids is in asending order. return value is the number of row swaps made. created -- 98apr15, cca ---------------------------------------------- */ int A2_sortRowsUp ( A2 *mtx, int nrow, int rowids[] ) { int ii, minrow, minrowid, nswap, target ; /* --------------- check the input --------------- */ if ( mtx == NULL || mtx->n1 < nrow || nrow < 0 || rowids == NULL ) { fprintf(stderr, "\n fatal error in A2_sortRowsUp(%p,%d,%p)" "\n bad input\n", mtx, nrow, rowids) ; if ( mtx != NULL ) { A2_writeStats(mtx, stderr) ; } exit(-1) ; } if ( ! (A2_IS_REAL(mtx) || A2_IS_COMPLEX(mtx)) ) { fprintf(stderr, "\n fatal error in A2_sortRowsUp(%p,%d,%p)" "\n bad type %d, must be SPOOLES_REAL or SPOOLES_COMPLEX\n", mtx, nrow, rowids, mtx->type) ; exit(-1) ; } nswap = 0 ; if ( mtx->inc1 == 1 ) { double *dvtmp ; int jcol, ncol ; int *ivtmp ; /* --------------------------------------------------- matrix is stored by columns, so permute each column --------------------------------------------------- */ ivtmp = IVinit(nrow, -1) ; if ( A2_IS_REAL(mtx) ) { dvtmp = DVinit(nrow, 0.0) ; } else if ( A2_IS_COMPLEX(mtx) ) { dvtmp = DVinit(2*nrow, 0.0) ; } IVramp(nrow, ivtmp, 0, 1) ; IV2qsortUp(nrow, rowids, ivtmp) ; ncol = mtx->n2 ; for ( jcol = 0 ; jcol < ncol ; jcol++ ) { if ( A2_IS_REAL(mtx) ) { DVcopy(nrow, dvtmp, A2_column(mtx, jcol)) ; DVgather(nrow, A2_column(mtx, jcol), dvtmp, ivtmp) ; } else if ( A2_IS_COMPLEX(mtx) ) { ZVcopy(nrow, dvtmp, A2_column(mtx, jcol)) ; ZVgather(nrow, A2_column(mtx, jcol), dvtmp, ivtmp) ; } } IVfree(ivtmp) ; DVfree(dvtmp) ; } else { /* ---------------------------------------- use a simple insertion sort to swap rows ---------------------------------------- */ for ( target = 0 ; target < nrow ; target++ ) { minrow = target ; minrowid = rowids[target] ; for ( ii = target + 1 ; ii < nrow ; ii++ ) { if ( minrowid > rowids[ii] ) { minrow = ii ; minrowid = rowids[ii] ; } } if ( minrow != target ) { rowids[minrow] = rowids[target] ; rowids[target] = minrowid ; A2_swapRows(mtx, target, minrow) ; nswap++ ; } } } return(nswap) ; }
/* ------------------------------------------------- sort the columns of the matrix in ascending order of the colids[] vector. on return, colids is in asending order. return value is the number of column swaps made. created -- 98apr15, cca ------------------------------------------------- */ int A2_sortColumnsUp ( A2 *mtx, int ncol, int colids[] ) { int ii, mincol, mincolid, nswap, target ; /* --------------- check the input --------------- */ if ( mtx == NULL || mtx->n2 < ncol || ncol < 0 || colids == NULL ) { fprintf(stderr, "\n fatal error in A2_sortColumnsUp(%p,%d,%p)" "\n bad input\n", mtx, ncol, colids) ; if ( mtx != NULL ) { A2_writeStats(mtx, stderr) ; } exit(-1) ; } if ( ! (A2_IS_REAL(mtx) || A2_IS_COMPLEX(mtx)) ) { fprintf(stderr, "\n fatal error in A2_sortColumnsUp(%p,%d,%p)" "\n bad type %d, must be SPOOLES_REAL or SPOOLES_COMPLEX\n", mtx, ncol, colids, mtx->type) ; exit(-1) ; } nswap = 0 ; if ( mtx->inc2 == 1 ) { double *dvtmp ; int irow, nrow ; int *ivtmp ; /* --------------------------------------------------- matrix is stored by rows, so permute each row --------------------------------------------------- */ ivtmp = IVinit(ncol, -1) ; if ( A2_IS_REAL(mtx) ) { dvtmp = DVinit(ncol, 0.0) ; } else if ( A2_IS_COMPLEX(mtx) ) { dvtmp = DVinit(2*ncol, 0.0) ; } IVramp(ncol, ivtmp, 0, 1) ; IV2qsortUp(ncol, colids, ivtmp) ; nrow = mtx->n1 ; for ( irow = 0 ; irow < nrow ; irow++ ) { if ( A2_IS_REAL(mtx) ) { DVcopy(ncol, dvtmp, A2_row(mtx, irow)) ; DVgather(ncol, A2_row(mtx, irow), dvtmp, ivtmp) ; } else if ( A2_IS_COMPLEX(mtx) ) { ZVcopy(ncol, dvtmp, A2_row(mtx, irow)) ; ZVgather(ncol, A2_row(mtx, irow), dvtmp, ivtmp) ; } } IVfree(ivtmp) ; DVfree(dvtmp) ; } else { /* ---------------------------------------- use a simple insertion sort to swap cols ---------------------------------------- */ for ( target = 0 ; target < ncol ; target++ ) { mincol = target ; mincolid = colids[target] ; for ( ii = target + 1 ; ii < ncol ; ii++ ) { if ( mincolid > colids[ii] ) { mincol = ii ; mincolid = colids[ii] ; } } if ( mincol != target ) { colids[mincol] = colids[target] ; colids[target] = mincolid ; A2_swapColumns(mtx, target, mincol) ; nswap++ ; } } } return(nswap) ; }
/*--------------------------------------------------------------------*/ int main ( int argc, char *argv[] ) /* ----------------------------- test the SubMtx_solve() method. created -- 98apr15, cca ----------------------------- */ { SubMtx *mtxA, *mtxB, *mtxX ; double idot, rdot, t1, t2 ; double *entB, *entX ; Drand *drand ; FILE *msgFile ; int inc1, inc2, mode, msglvl, ncolA, nentA, nrowA, ncolB, nrowB, ncolX, nrowX, seed, type ; if ( argc != 9 ) { fprintf(stdout, "\n\n usage : %s msglvl msgFile type mode nrowA nentA ncolB seed" "\n msglvl -- message level" "\n msgFile -- message file" "\n type -- type of matrix A" "\n 1 -- real" "\n 2 -- complex" "\n mode -- mode of matrix A" "\n 2 -- sparse stored by rows" "\n 3 -- sparse stored by columns" "\n 5 -- sparse stored by subrows" "\n 6 -- sparse stored by subcolumns" "\n 7 -- diagonal" "\n 8 -- block diagonal symmetric" "\n 9 -- block diagonal hermitian" "\n nrowA -- # of rows in matrix A" "\n nentA -- # of entries in matrix A" "\n ncolB -- # of columns in matrix B" "\n seed -- random number seed" "\n", argv[0]) ; return(0) ; } if ( (msglvl = atoi(argv[1])) < 0 ) { fprintf(stderr, "\n message level must be positive\n") ; spoolesFatal(); } if ( strcmp(argv[2], "stdout") == 0 ) { msgFile = stdout ; } else if ( (msgFile = fopen(argv[2], "a")) == NULL ) { fprintf(stderr, "\n unable to open file %s\n", argv[2]) ; return(-1) ; } type = atoi(argv[3]) ; mode = atoi(argv[4]) ; nrowA = atoi(argv[5]) ; nentA = atoi(argv[6]) ; ncolB = atoi(argv[7]) ; seed = atoi(argv[8]) ; fprintf(msgFile, "\n %% %s:" "\n %% msglvl = %d" "\n %% msgFile = %s" "\n %% type = %d" "\n %% mode = %d" "\n %% nrowA = %d" "\n %% nentA = %d" "\n %% ncolB = %d" "\n %% seed = %d", argv[0], msglvl, argv[2], type, mode, nrowA, nentA, ncolB, seed) ; ncolA = nrowA ; nrowB = nrowA ; nrowX = nrowA ; ncolX = ncolB ; /* ----------------------------- check for errors in the input ----------------------------- */ if ( nrowA <= 0 || nentA <= 0 || ncolB <= 0 ) { fprintf(stderr, "\n invalid input\n") ; spoolesFatal(); } switch ( type ) { case SPOOLES_REAL : switch ( mode ) { case SUBMTX_DENSE_SUBROWS : case SUBMTX_SPARSE_ROWS : case SUBMTX_DENSE_SUBCOLUMNS : case SUBMTX_SPARSE_COLUMNS : case SUBMTX_DIAGONAL : case SUBMTX_BLOCK_DIAGONAL_SYM : break ; default : fprintf(stderr, "\n invalid mode %d\n", mode) ; spoolesFatal(); } break ; case SPOOLES_COMPLEX : switch ( mode ) { case SUBMTX_DENSE_SUBROWS : case SUBMTX_SPARSE_ROWS : case SUBMTX_DENSE_SUBCOLUMNS : case SUBMTX_SPARSE_COLUMNS : case SUBMTX_DIAGONAL : case SUBMTX_BLOCK_DIAGONAL_SYM : case SUBMTX_BLOCK_DIAGONAL_HERM : break ; default : fprintf(stderr, "\n invalid mode %d\n", mode) ; spoolesFatal(); } break ; default : fprintf(stderr, "\n invalid type %d\n", type) ; spoolesFatal(); break ; } /* -------------------------------------- initialize the random number generator -------------------------------------- */ drand = Drand_new() ; Drand_init(drand) ; Drand_setSeed(drand, seed) ; Drand_setNormal(drand, 0.0, 1.0) ; /* ------------------------------ initialize the X SubMtx object ------------------------------ */ MARKTIME(t1) ; mtxX = SubMtx_new() ; SubMtx_initRandom(mtxX, type, SUBMTX_DENSE_COLUMNS, 0, 0, nrowX, ncolX, nrowX*ncolX, ++seed) ; SubMtx_denseInfo(mtxX, &nrowX, &ncolX, &inc1, &inc2, &entX) ; MARKTIME(t2) ; fprintf(msgFile, "\n %% CPU : %.3f to initialize X SubMtx object", t2 - t1) ; if ( msglvl > 1 ) { fprintf(msgFile, "\n\n %% X SubMtx object") ; fprintf(msgFile, "\n X = zeros(%d,%d) ;", nrowX, ncolX) ; SubMtx_writeForMatlab(mtxX, "X", msgFile) ; fflush(msgFile) ; } /* ------------------------------ initialize the B SubMtx object ------------------------------ */ MARKTIME(t1) ; mtxB = SubMtx_new() ; SubMtx_init(mtxB, type, SUBMTX_DENSE_COLUMNS, 0, 0, nrowB, ncolB, nrowB*ncolB) ; SubMtx_denseInfo(mtxB, &nrowB, &ncolB, &inc1, &inc2, &entB) ; switch ( mode ) { case SUBMTX_DENSE_SUBROWS : case SUBMTX_SPARSE_ROWS : case SUBMTX_DENSE_SUBCOLUMNS : case SUBMTX_SPARSE_COLUMNS : if ( SUBMTX_IS_REAL(mtxX) ) { DVcopy(nrowB*ncolB, entB, entX) ; } else if ( SUBMTX_IS_COMPLEX(mtxX) ) { ZVcopy(nrowB*ncolB, entB, entX) ; } break ; default : if ( SUBMTX_IS_REAL(mtxX) ) { DVzero(nrowB*ncolB, entB) ; } else if ( SUBMTX_IS_COMPLEX(mtxX) ) { DVzero(2*nrowB*ncolB, entB) ; } break ; } MARKTIME(t2) ; fprintf(msgFile, "\n %% CPU : %.3f to initialize B SubMtx object", t2 - t1) ; if ( msglvl > 1 ) { fprintf(msgFile, "\n\n %% B SubMtx object") ; fprintf(msgFile, "\n B = zeros(%d,%d) ;", nrowB, ncolB) ; SubMtx_writeForMatlab(mtxB, "B", msgFile) ; fflush(msgFile) ; } /* ------------------------------------- initialize the A matrix SubMtx object ------------------------------------- */ seed++ ; mtxA = SubMtx_new() ; switch ( mode ) { case SUBMTX_DENSE_SUBROWS : case SUBMTX_SPARSE_ROWS : SubMtx_initRandomLowerTriangle(mtxA, type, mode, 0, 0, nrowA, ncolA, nentA, seed, 1) ; break ; case SUBMTX_DENSE_SUBCOLUMNS : case SUBMTX_SPARSE_COLUMNS : SubMtx_initRandomUpperTriangle(mtxA, type, mode, 0, 0, nrowA, ncolA, nentA, seed, 1) ; break ; case SUBMTX_DIAGONAL : case SUBMTX_BLOCK_DIAGONAL_SYM : case SUBMTX_BLOCK_DIAGONAL_HERM : SubMtx_initRandom(mtxA, type, mode, 0, 0, nrowA, ncolA, nentA, seed) ; break ; default : fprintf(stderr, "\n fatal error in test_solve" "\n invalid mode = %d", mode) ; spoolesFatal(); } if ( msglvl > 1 ) { fprintf(msgFile, "\n\n %% A SubMtx object") ; fprintf(msgFile, "\n A = zeros(%d,%d) ;", nrowA, ncolA) ; SubMtx_writeForMatlab(mtxA, "A", msgFile) ; fflush(msgFile) ; } /* -------------------------------------------------------- compute B = A * X (for diagonal and block diagonal) or B = (I + A) * X (for lower and upper triangular) -------------------------------------------------------- */ if ( SUBMTX_IS_REAL(mtxA) ) { DV *colDV, *rowDV ; double value, *colX, *rowA, *pBij, *pXij ; int irowA, jcolX ; colDV = DV_new() ; DV_init(colDV, nrowA, NULL) ; colX = DV_entries(colDV) ; rowDV = DV_new() ; DV_init(rowDV, nrowA, NULL) ; rowA = DV_entries(rowDV) ; for ( jcolX = 0 ; jcolX < ncolB ; jcolX++ ) { SubMtx_fillColumnDV(mtxX, jcolX, colDV) ; for ( irowA = 0 ; irowA < nrowA ; irowA++ ) { SubMtx_fillRowDV(mtxA, irowA, rowDV) ; SubMtx_locationOfRealEntry(mtxX, irowA, jcolX, &pXij) ; SubMtx_locationOfRealEntry(mtxB, irowA, jcolX, &pBij) ; value = DVdot(nrowA, rowA, colX) ; switch ( mode ) { case SUBMTX_DENSE_SUBROWS : case SUBMTX_SPARSE_ROWS : case SUBMTX_DENSE_SUBCOLUMNS : case SUBMTX_SPARSE_COLUMNS : *pBij = *pXij + value ; break ; case SUBMTX_DIAGONAL : case SUBMTX_BLOCK_DIAGONAL_SYM : *pBij = value ; break ; } } } DV_free(colDV) ; DV_free(rowDV) ; } else if ( SUBMTX_IS_COMPLEX(mtxA) ) { ZV *colZV, *rowZV ; double *colX, *rowA, *pBIij, *pBRij, *pXIij, *pXRij ; int irowA, jcolX ; colZV = ZV_new() ; ZV_init(colZV, nrowA, NULL) ; colX = ZV_entries(colZV) ; rowZV = ZV_new() ; ZV_init(rowZV, nrowA, NULL) ; rowA = ZV_entries(rowZV) ; for ( jcolX = 0 ; jcolX < ncolB ; jcolX++ ) { SubMtx_fillColumnZV(mtxX, jcolX, colZV) ; for ( irowA = 0 ; irowA < nrowA ; irowA++ ) { SubMtx_fillRowZV(mtxA, irowA, rowZV) ; SubMtx_locationOfComplexEntry(mtxX, irowA, jcolX, &pXRij, &pXIij) ; SubMtx_locationOfComplexEntry(mtxB, irowA, jcolX, &pBRij, &pBIij) ; ZVdotU(nrowA, rowA, colX, &rdot, &idot) ; switch ( mode ) { case SUBMTX_DENSE_SUBROWS : case SUBMTX_SPARSE_ROWS : case SUBMTX_DENSE_SUBCOLUMNS : case SUBMTX_SPARSE_COLUMNS : *pBRij = *pXRij + rdot ; *pBIij = *pXIij + idot ; break ; case SUBMTX_DIAGONAL : case SUBMTX_BLOCK_DIAGONAL_SYM : case SUBMTX_BLOCK_DIAGONAL_HERM : *pBRij = rdot ; *pBIij = idot ; break ; } } } ZV_free(colZV) ; ZV_free(rowZV) ; } /* ---------------------- print out the matrices ---------------------- */ if ( msglvl > 1 ) { fprintf(msgFile, "\n\n %% X SubMtx object") ; fprintf(msgFile, "\n X = zeros(%d,%d) ;", nrowX, ncolX) ; SubMtx_writeForMatlab(mtxX, "X", msgFile) ; fprintf(msgFile, "\n\n %% A SubMtx object") ; fprintf(msgFile, "\n A = zeros(%d,%d) ;", nrowA, ncolA) ; SubMtx_writeForMatlab(mtxA, "A", msgFile) ; fprintf(msgFile, "\n\n %% B SubMtx object") ; fprintf(msgFile, "\n B = zeros(%d,%d) ;", nrowB, ncolB) ; SubMtx_writeForMatlab(mtxB, "B", msgFile) ; fflush(msgFile) ; } /* ----------------- check with matlab ----------------- */ if ( msglvl > 1 ) { switch ( mode ) { case SUBMTX_DENSE_SUBROWS : case SUBMTX_SPARSE_ROWS : case SUBMTX_DENSE_SUBCOLUMNS : case SUBMTX_SPARSE_COLUMNS : fprintf(msgFile, "\n\n emtx = abs(B - X - A*X) ;" "\n\n condA = cond(eye(%d,%d) + A)" "\n\n maxabsZ = max(max(abs(emtx))) ", nrowA, nrowA) ; fflush(msgFile) ; break ; case SUBMTX_DIAGONAL : case SUBMTX_BLOCK_DIAGONAL_SYM : case SUBMTX_BLOCK_DIAGONAL_HERM : fprintf(msgFile, "\n\n emtx = abs(B - A*X) ;" "\n\n condA = cond(A)" "\n\n maxabsZ = max(max(abs(emtx))) ") ; fflush(msgFile) ; break ; } } /* ---------------------------------------- compute the solve DY = B or (I + A)Y = B ---------------------------------------- */ SubMtx_solve(mtxA, mtxB) ; if ( msglvl > 1 ) { fprintf(msgFile, "\n\n %% Y SubMtx object") ; fprintf(msgFile, "\n Y = zeros(%d,%d) ;", nrowB, ncolB) ; SubMtx_writeForMatlab(mtxB, "Y", msgFile) ; fprintf(msgFile, "\n\n %% solerror = abs(Y - X) ;" "\n\n solerror = abs(Y - X) ;" "\n\n maxabserror = max(max(solerror)) ") ; fflush(msgFile) ; } /* ------------------------ free the working storage ------------------------ */ SubMtx_free(mtxA) ; SubMtx_free(mtxX) ; SubMtx_free(mtxB) ; Drand_free(drand) ; fprintf(msgFile, "\n") ; return(1) ; }