/* ----------------------------- 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 real row in the matrix created -- 98jan28, cca ------------------------------ */ void InpMtx_inputRealRow ( InpMtx *inpmtx, int row, int rowsize, int rowind[], double rowent[] ) { /* -------------- check the data -------------- */ if ( inpmtx == NULL || row < 0 || rowsize < 0 || rowind == NULL || rowent == NULL ) { fprintf(stderr, "\n fatal error in InpMtx_inputRealRow(%p,%d,%d,%p,%p)" "\n bad input\n", inpmtx, row, rowsize, rowind, rowent) ; spoolesFatal(); } if ( ! INPMTX_IS_REAL_ENTRIES(inpmtx) ) { fprintf(stderr, "\n fatal error in InpMtx_inputRealRow(%p,%d,%d,%p,%p)" "\n inputMode is not SPOOLES_REAL\n", inpmtx, row, rowsize, rowind, rowent) ; spoolesFatal(); } inputRow(inpmtx, row, rowsize, rowind, rowent) ; return ; }
/* -------------------------------------------------------------------- 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 ; }
/* --------------------------------- set the present number of entries created -- 98jan28, cca -------------------------------- */ void InpMtx_setNent ( InpMtx *inpmtx, int newnent ) { /* --------------- check the input --------------- */ if ( inpmtx == NULL || newnent < 0 ) { fprintf(stderr, "\n fatal error in InpMtx_setNent(%p,%d)" "\n bad input\n", inpmtx, newnent) ; spoolesFatal(); } if ( inpmtx->maxnent < newnent ) { /* ------------------------------------------------------- newnent requested is more than maxnent, set new maxnent ------------------------------------------------------- */ InpMtx_setMaxnent(inpmtx, newnent) ; } inpmtx->nent = newnent ; IV_setSize(&inpmtx->ivec1IV, newnent) ; IV_setSize(&inpmtx->ivec2IV, newnent) ; if ( INPMTX_IS_REAL_ENTRIES(inpmtx) ) { DV_setSize(&(inpmtx->dvecDV), newnent) ; } else if ( INPMTX_IS_COMPLEX_ENTRIES(inpmtx) ) { DV_setSize(&inpmtx->dvecDV, 2*newnent) ; } return ; }
/* ------------------------------------------------------------- input a number of (row,column, entry) triples into the matrix created -- 98jan28, cca ------------------------------------------------------------- */ void InpMtx_inputRealTriples ( InpMtx *inpmtx, int ntriples, int rowids[], int colids[], double entries[] ) { /* -------------- check the data -------------- */ if ( inpmtx == NULL || ntriples < 0 || rowids == NULL || colids == NULL || entries == NULL ) { fprintf(stderr, "\n fatal error in InpMtx_inputRealTriples(%p,%d,%p,%p,%p)" "\n bad input\n", inpmtx, ntriples, rowids, colids, entries) ; spoolesFatal(); } if ( ! INPMTX_IS_REAL_ENTRIES(inpmtx) ) { fprintf(stderr, "\n fatal error in InpMtx_inputRealEntry(%p,%d,%p,%p,%p)" "\n coordType must be COMPLEX_REAL_ENTRIES\n", inpmtx, ntriples, rowids, colids, entries) ; spoolesFatal(); } inputTriples(inpmtx, ntriples, rowids, colids, entries) ; return ; }
/* ----------------------------- input a chevron in the matrix created -- 98jan28, cca ----------------------------- */ void InpMtx_inputRealChevron ( InpMtx *inpmtx, int chv, int chvsize, int chvind[], double chvent[] ) { /* -------------- check the data -------------- */ if ( inpmtx == NULL || chv < 0 || chvsize < 0 || chvind == NULL || chvent == NULL ) { fprintf(stderr, "\n fatal error in InpMtx_inputRealChevron(%p,%d,%d,%p,%p)" "\n bad input\n", inpmtx, chv, chvsize, chvind, chvent) ; spoolesFatal(); } if ( ! INPMTX_IS_REAL_ENTRIES(inpmtx) ) { fprintf(stderr, "\n fatal error in InpMtx_inputRealChevron(%p,%d,%d,%p,%p)" "\n inputMode must be SPOOLES_REAL\n", inpmtx, chv, chvsize, chvind, chvent) ; spoolesFatal(); } inputChevron(inpmtx, chv, chvsize, chvind, chvent) ; return ; }
/* -------------------------------------------------------------- sets the maximum numnber of entries. this methods resizes the ivec1[], ivece2[] and dvec[] vectors if newmaxnent != maxnent created -- 98jan28, cca -------------------------------------------------------------- */ void InpMtx_setMaxnent ( InpMtx *inpmtx, int newmaxnent ) { /* --------------- check the input --------------- */ if ( inpmtx == NULL || newmaxnent < 0 ) { fprintf(stderr, "\n fatal error in InpMtx_setMaxnent(%p, %d)" "\n bad input\n", inpmtx, newmaxnent) ; spoolesFatal(); } if ( newmaxnent != inpmtx->maxnent ) { IV_setMaxsize(&(inpmtx->ivec1IV), newmaxnent) ; IV_setMaxsize(&(inpmtx->ivec2IV), newmaxnent) ; if ( INPMTX_IS_REAL_ENTRIES(inpmtx) ) { DV_setMaxsize(&(inpmtx->dvecDV), newmaxnent) ; } else if ( INPMTX_IS_COMPLEX_ENTRIES(inpmtx) ) { DV_setMaxsize(&(inpmtx->dvecDV), 2*newmaxnent) ; } } inpmtx->maxnent = newmaxnent ; 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 ; }
/* --------------------------------------- given the data is in raw triples, sort and compress the data created -- 98jan28, cca modified -- 98sep04, cca test to see if the sort is necessary --------------------------------------- */ void InpMtx_sortAndCompress ( InpMtx *inpmtx ) { int ient, nent, sortMustBeDone ; int *ivec1, *ivec2 ; /* --------------- check the input --------------- */ if ( inpmtx == NULL ) { fprintf(stderr, "\n fatal error in InpMtx_sortAndCompress(%p)" "\n bad input\n", inpmtx) ; exit(-1) ; } if ( INPMTX_IS_SORTED(inpmtx) || INPMTX_IS_BY_VECTORS(inpmtx) || (nent = inpmtx->nent) == 0 ) { inpmtx->storageMode = INPMTX_SORTED ; return ; } ivec1 = InpMtx_ivec1(inpmtx) ; ivec2 = InpMtx_ivec2(inpmtx) ; sortMustBeDone = 0 ; for ( ient = 1 ; ient < nent ; ient++ ) { if ( ivec1[ient-1] > ivec1[ient] || ( ivec1[ient-1] == ivec1[ient] && ivec2[ient-1] > ivec2[ient] ) ) { sortMustBeDone = 1 ; break ; } } if ( sortMustBeDone == 1 ) { if ( INPMTX_IS_INDICES_ONLY(inpmtx) ) { inpmtx->nent = IV2sortUpAndCompress(nent, ivec1, ivec2) ; } else if ( INPMTX_IS_REAL_ENTRIES(inpmtx) ) { double *dvec = InpMtx_dvec(inpmtx) ; inpmtx->nent = IV2DVsortUpAndCompress(nent, ivec1, ivec2, dvec) ; } else if ( INPMTX_IS_COMPLEX_ENTRIES(inpmtx) ) { double *dvec = InpMtx_dvec(inpmtx) ; inpmtx->nent = IV2ZVsortUpAndCompress(nent, ivec1, ivec2, dvec) ; } } inpmtx->storageMode = INPMTX_SORTED ; return ; }
/* ---------------------------------- input a single entry in the matrix created -- 98jan28, cca ---------------------------------- */ static void inputEntry ( InpMtx *inpmtx, int row, int col, double real, double imag ) { int nent ; int *ivec1, *ivec2 ; prepareToAddNewEntries(inpmtx, 1) ; nent = inpmtx->nent ; ivec1 = IV_entries(&inpmtx->ivec1IV) ; ivec2 = IV_entries(&inpmtx->ivec2IV) ; if ( INPMTX_IS_BY_ROWS(inpmtx) ) { ivec1[nent] = row ; ivec2[nent] = col ; } else if ( INPMTX_IS_BY_COLUMNS(inpmtx) ) { ivec1[nent] = col ; ivec2[nent] = row ; } else if ( INPMTX_IS_BY_CHEVRONS(inpmtx) ) { if ( row <= col ) { ivec1[nent] = row ; ivec2[nent] = col - row ; } else { ivec1[nent] = col ; ivec2[nent] = col - row ; } } IV_setSize(&inpmtx->ivec1IV, nent + 1) ; IV_setSize(&inpmtx->ivec2IV, nent + 1) ; if ( INPMTX_IS_REAL_ENTRIES(inpmtx) ) { double *dvec = DV_entries(&inpmtx->dvecDV) ; dvec[nent] = real ; DV_setSize(&inpmtx->dvecDV, nent + 1) ; } else if ( INPMTX_IS_COMPLEX_ENTRIES(inpmtx) ) { double *dvec = DV_entries(&inpmtx->dvecDV) ; dvec[2*nent] = real ; dvec[2*nent+1] = imag ; DV_setSize(&inpmtx->dvecDV, 2*(nent + 1)) ; } inpmtx->nent++ ; inpmtx->storageMode = INPMTX_RAW_DATA ; 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 ; }
/* ----------------------- input a matrix created -- 98jan28, cca ----------------------- */ void InpMtx_inputRealMatrix ( InpMtx *inpmtx, int nrow, int ncol, int rowstride, int colstride, int rowind[], int colind[], double mtxent[] ) { /* -------------- check the data -------------- */ if ( inpmtx == NULL || nrow < 0 || ncol < 0 || rowstride < 1 || colstride < 1 || rowind == NULL || colind == NULL || mtxent == NULL ) { fprintf(stderr, "\n fatal error in InpMtx_inputRealMatrix(%p,%d,%d,%d,%d,%p,%p,%p)" "\n bad input\n", inpmtx, nrow, ncol, rowstride, colstride, rowind, colind, mtxent) ; spoolesFatal(); } if ( ! INPMTX_IS_REAL_ENTRIES(inpmtx) ) { fprintf(stderr, "\n fatal error in InpMtx_inputRealMatrix(%p,%d,%d,%d,%d,%p,%p,%p)" "\n inputMode must be SPOOLES_REAL\n", inpmtx, nrow, ncol, rowstride, colstride, rowind, colind, mtxent) ; spoolesFatal(); } if ( nrow == 0 || ncol == 0 ) { return ; } inputMatrix(inpmtx, nrow, ncol, rowstride, colstride, rowind, colind, mtxent) ; return ; }
/* --------------------------------------- input a single real entry in the matrix created -- 98jan28, cca --------------------------------------- */ void InpMtx_inputRealEntry ( InpMtx *inpmtx, int row, int col, double value ) { /* -------------- check the data -------------- */ if ( inpmtx == NULL || row < 0 || col < 0 ) { fprintf(stderr, "\n fatal error in InpMtx_inputRealEntry(%p,%d,%d,%e)" "\n bad inputReal\n", inpmtx, row, col, value) ; spoolesFatal(); } if ( !( INPMTX_IS_BY_ROWS(inpmtx) || INPMTX_IS_BY_COLUMNS(inpmtx) || INPMTX_IS_BY_CHEVRONS(inpmtx) ) ) { fprintf(stderr, "\n fatal error in InpMtx_inputRealEntry(%p,%d,%d,%e)" "\n bad coordType = %d\n", inpmtx, row, col, value, inpmtx->coordType) ; spoolesFatal(); } if ( ! INPMTX_IS_REAL_ENTRIES(inpmtx) ) { fprintf(stderr, "\n fatal error in InpMtx_inputRealEntry(%p,%d,%d,%e)" "\n inputMode is not SPOOLES_REAL\n", inpmtx, row, col, value) ; spoolesFatal(); } inputEntry(inpmtx, row, col, value, 0.0) ; return ; }
/* ------------------------------------------------------------ load entries from sigma*A chv -- pointer to the Chv object that holds the front pencil -- pointer to a Pencil that holds the matrix entries msglvl -- message level msgFile -- message file created -- 97jul18, cca ------------------------------------------------------------ */ void FrontMtx_loadEntries ( Chv *chv, Pencil *pencil, int msglvl, FILE *msgFile ) { InpMtx *inpmtxA, *inpmtxB ; double one[2] = {1.0,0.0} ; double *sigma ; double *chvent ; int chvsize, ichv, ncol, nD, nL, nU ; int *chvind, *colind ; /* --------------- check the input --------------- */ if ( chv == NULL || (msglvl > 0 && msgFile == NULL) ) { fprintf(stderr, "\n fatal error in FrontMtx_loadEntries(%p,%p,%d,%p)" "\n bad input\n", chv, pencil, msglvl, msgFile) ; exit(-1) ; } if ( msglvl > 3 ) { fprintf(msgFile, "\n\n # inside loadEntries for chv %d" ", sigma = %12.4e + i*%12.4e", chv->id, pencil->sigma[0], pencil->sigma[1]) ; fflush(msgFile) ; } Chv_dimensions(chv, &nD, &nL, &nU) ; Chv_columnIndices(chv, &ncol, &colind) ; /* ---------------------------------------- load the original entries, A + sigma * B ---------------------------------------- */ inpmtxA = pencil->inpmtxA ; sigma = pencil->sigma ; inpmtxB = pencil->inpmtxB ; if ( inpmtxA != NULL ) { int ii ; /* ------------------- load entries from A ------------------- */ for ( ii = 0 ; ii < nD ; ii++ ) { ichv = colind[ii] ; if ( INPMTX_IS_REAL_ENTRIES(inpmtxA) ) { InpMtx_realVector(inpmtxA, ichv, &chvsize, &chvind, &chvent) ; } else if ( INPMTX_IS_COMPLEX_ENTRIES(inpmtxA) ) { InpMtx_complexVector(inpmtxA, ichv, &chvsize, &chvind, &chvent) ; } if ( chvsize > 0 ) { if ( msglvl > 3 ) { int ierr ; fprintf(msgFile, "\n inpmtxA chevron %d : chvsize = %d", ichv, chvsize) ; fprintf(msgFile, "\n chvind") ; IVfp80(msgFile, chvsize, chvind, 80, &ierr) ; fprintf(msgFile, "\n chvent") ; if ( INPMTX_IS_REAL_ENTRIES(inpmtxA) ) { DVfprintf(msgFile, chvsize, chvent) ; } else if ( INPMTX_IS_COMPLEX_ENTRIES(inpmtxA) ) { DVfprintf(msgFile, 2*chvsize, chvent) ; } fflush(msgFile) ; } Chv_addChevron(chv, one, ichv, chvsize, chvind, chvent) ; } } } else { double *entries ; int ii, off, stride ; /* ----------------- load the identity ----------------- */ entries = Chv_entries(chv) ; if ( CHV_IS_REAL(chv) ) { if ( CHV_IS_SYMMETRIC(chv) || CHV_IS_HERMITIAN(chv) ) { stride = nD + chv->nU ; off = 0 ; for ( ii = 0 ; ii < nD ; ii++ ) { entries[off] += 1.0 ; off += stride ; stride-- ; } } else if ( CHV_IS_NONSYMMETRIC(chv) ) { stride = 2*nD + chv->nL + chv->nU - 2 ; off = nD + chv->nL - 1 ; for ( ii = 0 ; ii < nD ; ii++ ) { entries[off] += 1.0 ; off += stride ; stride -= 2 ; } } } else if ( CHV_IS_COMPLEX(chv) ) { if ( CHV_IS_SYMMETRIC(chv) || CHV_IS_HERMITIAN(chv) ) { stride = nD + chv->nU ; off = 0 ; for ( ii = 0 ; ii < nD ; ii++ ) { entries[2*off] += 1.0 ; off += stride ; stride-- ; } } else if ( CHV_IS_NONSYMMETRIC(chv) ) { stride = 2*nD + chv->nL + chv->nU - 2 ; off = nD + chv->nL - 1 ; for ( ii = 0 ; ii < nD ; ii++ ) { entries[2*off] += 1.0 ; off += stride ; stride -= 2 ; } } } } if ( inpmtxB != NULL ) { int ii ; /* ------------------------- load entries from sigma*B ------------------------- */ for ( ii = 0 ; ii < nD ; ii++ ) { ichv = colind[ii] ; if ( INPMTX_IS_REAL_ENTRIES(inpmtxB) ) { InpMtx_realVector(inpmtxB, ichv, &chvsize, &chvind, &chvent) ; } else if ( INPMTX_IS_COMPLEX_ENTRIES(inpmtxA) ) { InpMtx_complexVector(inpmtxB, ichv, &chvsize, &chvind, &chvent) ; } if ( chvsize > 0 ) { if ( msglvl > 3 ) { int ierr ; fprintf(msgFile, "\n inpmtxB chevron %d : chvsize = %d", ichv, chvsize) ; fprintf(msgFile, "\n chvind") ; IVfp80(msgFile, chvsize, chvind, 80, &ierr) ; fprintf(msgFile, "\n chvent") ; if ( INPMTX_IS_REAL_ENTRIES(inpmtxA) ) { DVfprintf(msgFile, chvsize, chvent) ; } else if ( INPMTX_IS_COMPLEX_ENTRIES(inpmtxA) ) { DVfprintf(msgFile, 2*chvsize, chvent) ; } } Chv_addChevron(chv, sigma, ichv, chvsize, chvind, chvent) ; } } } else { double *entries ; int ii, off, stride ; /* -------------------------------------- load a scalar multiple of the identity -------------------------------------- */ entries = Chv_entries(chv) ; if ( CHV_IS_REAL(chv) ) { if ( CHV_IS_SYMMETRIC(chv) ) { stride = nD + chv->nU ; off = 0 ; for ( ii = 0 ; ii < nD ; ii++ ) { entries[off] += sigma[0] ; off += stride ; stride-- ; } } else if ( CHV_IS_NONSYMMETRIC(chv) ) { stride = 2*nD + chv->nL + chv->nU - 2 ; off = nD + chv->nL - 1 ; for ( ii = 0 ; ii < nD ; ii++ ) { entries[off] += sigma[0] ; off += stride ; stride -= 2 ; } } } else if ( CHV_IS_COMPLEX(chv) ) { if ( CHV_IS_SYMMETRIC(chv) || CHV_IS_HERMITIAN(chv) ) { if ( CHV_IS_HERMITIAN(chv) && sigma[1] != 0.0 ) { fprintf(stderr, "\n fatal error in FrontMtx_loadEntries()" "\n chevron is hermitian" "\n sigma = %12.4e + %12.4e*i\n", sigma[0], sigma[1]) ; exit(-1) ; } stride = nD + chv->nU ; off = 0 ; for ( ii = 0 ; ii < nD ; ii++ ) { entries[2*off] += sigma[0] ; entries[2*off+1] += sigma[1] ; off += stride ; stride-- ; } } else if ( CHV_IS_NONSYMMETRIC(chv) ) { stride = 2*nD + chv->nL + chv->nU - 2 ; off = nD + chv->nL - 1 ; for ( ii = 0 ; ii < nD ; ii++ ) { entries[2*off] += sigma[0] ; entries[2*off+1] += sigma[1] ; off += stride ; stride -= 2 ; } } } } return ; }
/*--------------------------------------------------------------------*/ int main ( int argc, char *argv[] ) /* ------------------------------------------------------------------ generate a random matrix and test a matrix-matrix multiply method. the output is a matlab file to test correctness. created -- 98jan29, cca -------------------------------------------------------------------- */ { DenseMtx *X, *Y, *Y2 ; double alpha[2] ; double alphaImag, alphaReal, t1, t2 ; double *zvec ; Drand *drand ; int col, dataType, ii, msglvl, ncolA, nitem, nops, nrhs, nrowA, nrowX, nrowY, nthread, row, seed, storageMode, symflag, transposeflag ; int *colids, *rowids ; InpMtx *A ; FILE *msgFile ; if ( argc != 15 ) { fprintf(stdout, "\n\n %% usage : %s msglvl msgFile symflag storageMode " "\n %% nrow ncol nent nrhs seed alphaReal alphaImag nthread" "\n %% msglvl -- message level" "\n %% msgFile -- message file" "\n %% dataType -- type of matrix entries" "\n %% 1 -- real" "\n %% 2 -- complex" "\n %% symflag -- symmetry flag" "\n %% 0 -- symmetric" "\n %% 1 -- hermitian" "\n %% 2 -- nonsymmetric" "\n %% storageMode -- storage mode" "\n %% 1 -- by rows" "\n %% 2 -- by columns" "\n %% 3 -- by chevrons, (requires nrow = ncol)" "\n %% transpose -- transpose flag" "\n %% 0 -- Y := Y + alpha * A * X" "\n %% 1 -- Y := Y + alpha * A^H * X, nonsymmetric only" "\n %% 2 -- Y := Y + alpha * A^T * X, nonsymmetric only" "\n %% nrowA -- number of rows in A" "\n %% ncolA -- number of columns in A" "\n %% nitem -- number of items" "\n %% nrhs -- number of right hand sides" "\n %% seed -- random number seed" "\n %% alphaReal -- y := y + alpha*A*x" "\n %% alphaImag -- y := y + alpha*A*x" "\n %% nthread -- # of threads" "\n", argv[0]) ; return(0) ; } msglvl = atoi(argv[1]) ; if ( strcmp(argv[2], "stdout") == 0 ) { msgFile = stdout ; } else if ( (msgFile = fopen(argv[2], "a")) == NULL ) { fprintf(stderr, "\n fatal error in %s" "\n unable to open file %s\n", argv[0], argv[2]) ; return(-1) ; } dataType = atoi(argv[3]) ; symflag = atoi(argv[4]) ; storageMode = atoi(argv[5]) ; transposeflag = atoi(argv[6]) ; nrowA = atoi(argv[7]) ; ncolA = atoi(argv[8]) ; nitem = atoi(argv[9]) ; nrhs = atoi(argv[10]) ; seed = atoi(argv[11]) ; alphaReal = atof(argv[12]) ; alphaImag = atof(argv[13]) ; nthread = atoi(argv[14]) ; fprintf(msgFile, "\n %% %s " "\n %% msglvl -- %d" "\n %% msgFile -- %s" "\n %% dataType -- %d" "\n %% symflag -- %d" "\n %% storageMode -- %d" "\n %% transposeflag -- %d" "\n %% nrowA -- %d" "\n %% ncolA -- %d" "\n %% nitem -- %d" "\n %% nrhs -- %d" "\n %% seed -- %d" "\n %% alphaReal -- %e" "\n %% alphaImag -- %e" "\n %% nthread -- %d" "\n", argv[0], msglvl, argv[2], dataType, symflag, storageMode, transposeflag, nrowA, ncolA, nitem, nrhs, seed, alphaReal, alphaImag, nthread) ; fflush(msgFile) ; if ( dataType != 1 && dataType != 2 ) { fprintf(stderr, "\n invalid value %d for dataType\n", dataType) ; spoolesFatal(); } if ( symflag != 0 && symflag != 1 && symflag != 2 ) { fprintf(stderr, "\n invalid value %d for symflag\n", symflag) ; spoolesFatal(); } if ( storageMode != 1 && storageMode != 2 && storageMode != 3 ) { fprintf(stderr, "\n invalid value %d for storageMode\n", storageMode) ; spoolesFatal(); } if ( transposeflag < 0 || transposeflag > 2 ) { fprintf(stderr, "\n error, transposeflag = %d, must be 0, 1 or 2", transposeflag) ; spoolesFatal(); } if ( (transposeflag == 1 && symflag != 2) || (transposeflag == 2 && symflag != 2) ) { fprintf(stderr, "\n error, transposeflag = %d, symflag = %d", transposeflag, symflag) ; spoolesFatal(); } if ( transposeflag == 1 && dataType != 2 ) { fprintf(stderr, "\n error, transposeflag = %d, dataType = %d", transposeflag, dataType) ; spoolesFatal(); } if ( symflag == 1 && dataType != 2 ) { fprintf(stderr, "\n symflag = 1 (hermitian), dataType != 2 (complex)") ; spoolesFatal(); } if ( nrowA <= 0 || ncolA <= 0 || nitem <= 0 ) { fprintf(stderr, "\n invalid value: nrow = %d, ncol = %d, nitem = %d", nrowA, ncolA, nitem) ; spoolesFatal(); } if ( symflag < 2 && nrowA != ncolA ) { fprintf(stderr, "\n invalid data: symflag = %d, nrow = %d, ncol = %d", symflag, nrowA, ncolA) ; spoolesFatal(); } alpha[0] = alphaReal ; alpha[1] = alphaImag ; /* ---------------------------- initialize the matrix object ---------------------------- */ A = InpMtx_new() ; InpMtx_init(A, storageMode, dataType, 0, 0) ; drand = Drand_new() ; /* ---------------------------------- generate a vector of nitem triples ---------------------------------- */ rowids = IVinit(nitem, -1) ; Drand_setUniform(drand, 0, nrowA) ; Drand_fillIvector(drand, nitem, rowids) ; colids = IVinit(nitem, -1) ; Drand_setUniform(drand, 0, ncolA) ; Drand_fillIvector(drand, nitem, colids) ; Drand_setUniform(drand, 0.0, 1.0) ; if ( INPMTX_IS_REAL_ENTRIES(A) ) { zvec = DVinit(nitem, 0.0) ; Drand_fillDvector(drand, nitem, zvec) ; } else if ( INPMTX_IS_COMPLEX_ENTRIES(A) ) { zvec = ZVinit(nitem, 0.0, 0.0) ; Drand_fillDvector(drand, 2*nitem, zvec) ; } /* ----------------------------------- assemble the entries entry by entry ----------------------------------- */ if ( msglvl > 1 ) { fprintf(msgFile, "\n\n A = zeros(%d,%d) ;", nrowA, ncolA) ; } if ( symflag == 1 ) { /* ---------------- hermitian matrix ---------------- */ for ( ii = 0 ; ii < nitem ; ii++ ) { if ( rowids[ii] == colids[ii] ) { zvec[2*ii+1] = 0.0 ; } if ( rowids[ii] <= colids[ii] ) { row = rowids[ii] ; col = colids[ii] ; } else { row = colids[ii] ; col = rowids[ii] ; } InpMtx_inputComplexEntry(A, row, col, zvec[2*ii], zvec[2*ii+1]) ; } } else if ( symflag == 0 ) { /* ---------------- symmetric matrix ---------------- */ if ( INPMTX_IS_REAL_ENTRIES(A) ) { for ( ii = 0 ; ii < nitem ; ii++ ) { if ( rowids[ii] <= colids[ii] ) { row = rowids[ii] ; col = colids[ii] ; } else { row = colids[ii] ; col = rowids[ii] ; } InpMtx_inputRealEntry(A, row, col, zvec[ii]) ; } } else if ( INPMTX_IS_COMPLEX_ENTRIES(A) ) { for ( ii = 0 ; ii < nitem ; ii++ ) { if ( rowids[ii] <= colids[ii] ) { row = rowids[ii] ; col = colids[ii] ; } else { row = colids[ii] ; col = rowids[ii] ; } InpMtx_inputComplexEntry(A, row, col, zvec[2*ii], zvec[2*ii+1]) ; } } } else { /* ------------------- nonsymmetric matrix ------------------- */ if ( INPMTX_IS_REAL_ENTRIES(A) ) { for ( ii = 0 ; ii < nitem ; ii++ ) { InpMtx_inputRealEntry(A, rowids[ii], colids[ii], zvec[ii]) ; } } else if ( INPMTX_IS_COMPLEX_ENTRIES(A) ) { for ( ii = 0 ; ii < nitem ; ii++ ) { InpMtx_inputComplexEntry(A, rowids[ii], colids[ii], zvec[2*ii], zvec[2*ii+1]) ; } } } InpMtx_changeStorageMode(A, INPMTX_BY_VECTORS) ; DVfree(zvec) ; if ( symflag == 0 || symflag == 1 ) { if ( INPMTX_IS_REAL_ENTRIES(A) ) { nops = 4*A->nent*nrhs ; } else if ( INPMTX_IS_COMPLEX_ENTRIES(A) ) { nops = 16*A->nent*nrhs ; } } else { if ( INPMTX_IS_REAL_ENTRIES(A) ) { nops = 2*A->nent*nrhs ; } else if ( INPMTX_IS_COMPLEX_ENTRIES(A) ) { nops = 8*A->nent*nrhs ; } } if ( msglvl > 1 ) { /* ------------------------------------------- write the assembled matrix to a matlab file ------------------------------------------- */ InpMtx_writeForMatlab(A, "A", msgFile) ; if ( symflag == 0 ) { fprintf(msgFile, "\n for k = 1:%d" "\n for j = k+1:%d" "\n A(j,k) = A(k,j) ;" "\n end" "\n end", nrowA, ncolA) ; } else if ( symflag == 1 ) { fprintf(msgFile, "\n for k = 1:%d" "\n for j = k+1:%d" "\n A(j,k) = ctranspose(A(k,j)) ;" "\n end" "\n end", nrowA, ncolA) ; } } /* ------------------------------- generate dense matrices X and Y ------------------------------- */ if ( transposeflag == 0 ) { nrowX = ncolA ; nrowY = nrowA ; } else { nrowX = nrowA ; nrowY = ncolA ; } X = DenseMtx_new() ; Y = DenseMtx_new() ; Y2 = DenseMtx_new() ; if ( INPMTX_IS_REAL_ENTRIES(A) ) { DenseMtx_init(X, SPOOLES_REAL, 0, 0, nrowX, nrhs, 1, nrowX) ; Drand_fillDvector(drand, nrowX*nrhs, DenseMtx_entries(X)) ; DenseMtx_init(Y, SPOOLES_REAL, 0, 0, nrowY, nrhs, 1, nrowY) ; Drand_fillDvector(drand, nrowY*nrhs, DenseMtx_entries(Y)) ; DenseMtx_init(Y2, SPOOLES_REAL, 0, 0, nrowY, nrhs, 1, nrowY) ; DVcopy(nrowY*nrhs, DenseMtx_entries(Y2), DenseMtx_entries(Y)) ; } else if ( INPMTX_IS_COMPLEX_ENTRIES(A) ) { DenseMtx_init(X, SPOOLES_COMPLEX, 0, 0, nrowX, nrhs, 1, nrowX) ; Drand_fillDvector(drand, 2*nrowX*nrhs, DenseMtx_entries(X)) ; DenseMtx_init(Y, SPOOLES_COMPLEX, 0, 0, nrowY, nrhs, 1, nrowY) ; Drand_fillDvector(drand, 2*nrowY*nrhs, DenseMtx_entries(Y)) ; DenseMtx_init(Y2, SPOOLES_COMPLEX, 0, 0, nrowY, nrhs, 1, nrowY) ; DVcopy(2*nrowY*nrhs, DenseMtx_entries(Y2), DenseMtx_entries(Y)) ; } if ( msglvl > 1 ) { fprintf(msgFile, "\n X = zeros(%d,%d) ;", nrowX, nrhs) ; DenseMtx_writeForMatlab(X, "X", msgFile) ; fprintf(msgFile, "\n Y = zeros(%d,%d) ;", nrowY, nrhs) ; DenseMtx_writeForMatlab(Y, "Y", msgFile) ; } /* -------------------------------------------- perform the matrix-matrix multiply in serial -------------------------------------------- */ if ( msglvl > 1 ) { fprintf(msgFile, "\n alpha = %20.12e + %20.2e*i;", alpha[0], alpha[1]); fprintf(msgFile, "\n Z = zeros(%d,1) ;", nrowY) ; } if ( transposeflag == 0 ) { MARKTIME(t1) ; if ( symflag == 0 ) { InpMtx_sym_mmm(A, Y, alpha, X) ; } else if ( symflag == 1 ) { InpMtx_herm_mmm(A, Y, alpha, X) ; } else if ( symflag == 2 ) { InpMtx_nonsym_mmm(A, Y, alpha, X) ; } MARKTIME(t2) ; if ( msglvl > 1 ) { DenseMtx_writeForMatlab(Y, "Z", msgFile) ; fprintf(msgFile, "\n maxerr = max(Z - Y - alpha*A*X) ") ; fprintf(msgFile, "\n") ; } } else if ( transposeflag == 1 ) { MARKTIME(t1) ; InpMtx_nonsym_mmm_H(A, Y, alpha, X) ; MARKTIME(t2) ; if ( msglvl > 1 ) { DenseMtx_writeForMatlab(Y, "Z", msgFile) ; fprintf(msgFile, "\n maxerr = max(Z - Y - alpha*ctranspose(A)*X) ") ; fprintf(msgFile, "\n") ; } } else if ( transposeflag == 2 ) { MARKTIME(t1) ; InpMtx_nonsym_mmm_T(A, Y, alpha, X) ; MARKTIME(t2) ; if ( msglvl > 1 ) { DenseMtx_writeForMatlab(Y, "Z", msgFile) ; fprintf(msgFile, "\n maxerr = max(Z - Y - alpha*transpose(A)*X) ") ; fprintf(msgFile, "\n") ; } } fprintf(msgFile, "\n %% %d ops, %.3f time, %.3f serial mflops", nops, t2 - t1, 1.e-6*nops/(t2 - t1)) ; /* -------------------------------------------------------- perform the matrix-matrix multiply in multithreaded mode -------------------------------------------------------- */ if ( msglvl > 1 ) { fprintf(msgFile, "\n alpha = %20.12e + %20.2e*i;", alpha[0], alpha[1]); fprintf(msgFile, "\n Z = zeros(%d,1) ;", nrowY) ; } if ( transposeflag == 0 ) { MARKTIME(t1) ; if ( symflag == 0 ) { InpMtx_MT_sym_mmm(A, Y2, alpha, X, nthread, msglvl, msgFile) ; } else if ( symflag == 1 ) { InpMtx_MT_herm_mmm(A, Y2, alpha, X, nthread, msglvl, msgFile) ; } else if ( symflag == 2 ) { InpMtx_MT_nonsym_mmm(A, Y2, alpha, X, nthread, msglvl, msgFile) ; } MARKTIME(t2) ; if ( msglvl > 1 ) { DenseMtx_writeForMatlab(Y2, "Z2", msgFile) ; fprintf(msgFile, "\n maxerr2 = max(Z2 - Y - alpha*A*X) ") ; fprintf(msgFile, "\n") ; } } else if ( transposeflag == 1 ) { MARKTIME(t1) ; InpMtx_MT_nonsym_mmm_H(A, Y2, alpha, X, nthread, msglvl, msgFile) ; MARKTIME(t2) ; if ( msglvl > 1 ) { DenseMtx_writeForMatlab(Y2, "Z2", msgFile) ; fprintf(msgFile, "\n maxerr2 = max(Z2 - Y - alpha*ctranspose(A)*X) ") ; fprintf(msgFile, "\n") ; } } else if ( transposeflag == 2 ) { MARKTIME(t1) ; InpMtx_MT_nonsym_mmm_T(A, Y2, alpha, X, nthread, msglvl, msgFile) ; MARKTIME(t2) ; if ( msglvl > 1 ) { DenseMtx_writeForMatlab(Y2, "Z2", msgFile) ; fprintf(msgFile, "\n maxerr2 = max(Z2 - Y - alpha*transpose(A)*X) ") ; fprintf(msgFile, "\n") ; } } fprintf(msgFile, "\n %% %d ops, %.3f time, %.3f MT mflops", nops, t2 - t1, 1.e-6*nops/(t2 - t1)) ; /* ------------------------ free the working storage ------------------------ */ InpMtx_free(A) ; DenseMtx_free(X) ; DenseMtx_free(Y) ; DenseMtx_free(Y2) ; IVfree(rowids) ; IVfree(colids) ; Drand_free(drand) ; fclose(msgFile) ; return(1) ; }
/* ------------------------------------------- set up the nthread MTmvmObj data structures ------------------------------------------- */ static MTmvmObj * setup ( InpMtx *A, DenseMtx *Y, double alpha[], DenseMtx *X, int nthread ) { double *dvec ; int ithread, nentA, nextra, nlocal, offset ; int *ivec1, *ivec2 ; MTmvmObj *MTmvmObjs, *obj ; /* --------------------------------- allocate nthread MTmvmObj objects --------------------------------- */ ALLOCATE(MTmvmObjs, struct _MTmvmObj, nthread) ; for ( ithread = 0, obj = MTmvmObjs ; ithread < nthread ; ithread++, obj++ ) { obj->A = InpMtx_new() ; if ( ithread == 0 ) { obj->Y = Y ; } else { obj->Y = DenseMtx_new() ; } obj->alpha[0] = alpha[0] ; obj->alpha[1] = alpha[1] ; obj->X = X ; } /* ---------------------------------------- set up and zero the replicated Y objects ---------------------------------------- */ for ( ithread = 0, obj = MTmvmObjs ; ithread < nthread ; ithread++, obj++ ) { if ( ithread > 0 ) { DenseMtx_init(obj->Y, Y->type, Y->rowid, Y->colid, Y->nrow, Y->ncol, Y->inc1, Y->inc2) ; DenseMtx_zero(obj->Y) ; } } /* ------------------------------------- set up the partitioned InpMtx objects ------------------------------------- */ nentA = InpMtx_nent(A) ; nlocal = nentA / nthread ; nextra = nentA % nthread ; ivec1 = InpMtx_ivec1(A) ; ivec2 = InpMtx_ivec2(A) ; if ( INPMTX_IS_REAL_ENTRIES(A) || INPMTX_IS_COMPLEX_ENTRIES(A) ) { dvec = InpMtx_dvec(A) ; } else { dvec = NULL ; } offset = 0 ; for ( ithread = 0, obj = MTmvmObjs ; ithread < nthread ; ithread++, obj++ ) { InpMtx_init(obj->A, A->coordType, A->inputMode, 0, 0) ; obj->A->storageMode = A->storageMode ; if ( ithread < nextra ) { obj->A->nent = nlocal + 1 ; } else { obj->A->nent = nlocal ; } IV_init(&(obj->A->ivec1IV), obj->A->nent, ivec1 + offset) ; IV_init(&(obj->A->ivec2IV), obj->A->nent, ivec2 + offset) ; if ( INPMTX_IS_REAL_ENTRIES(A) ) { DV_init(&(obj->A->dvecDV), obj->A->nent, dvec + offset) ; } else if ( INPMTX_IS_COMPLEX_ENTRIES(A) ) { DV_init(&(obj->A->dvecDV), obj->A->nent, dvec + 2*offset) ; } offset += obj->A->nent ; } return(MTmvmObjs) ; }
/*--------------------------------------------------------------------*/ static void InpMtx_MT_mmm ( int flag, InpMtx *A, DenseMtx *Y, double alpha[], DenseMtx *X, int nthread, int msglvl, FILE *msgFile ) { double t1, t2 ; int myid, nent, rc ; MTmvmObj *MTmvmObjs, *obj ; /* ------------------------------- set up the nthread data objects ------------------------------- */ MARKTIME(t1) ; MTmvmObjs = setup(A, Y, alpha, X, nthread) ; MARKTIME(t2) ; if ( msglvl > 0 ) { fprintf(msgFile, "\n %% CPU %8.3f : setup time", t2 - t1) ; } #if THREAD_TYPE == TT_POSIX { pthread_t *tids ; pthread_attr_t attr ; void *status ; /* ##### NOTE: for SGI machines, this command must be present ##### for the thread scheduling to be efficient. ##### this is NOT a POSIX call, but SGI needs it anyway pthread_setconcurrency(nthread) ; */ pthread_attr_init(&attr) ; /* pthread_attr_setscope(&attr, PTHREAD_SCOPE_SYSTEM) ; */ pthread_attr_setscope(&attr, PTHREAD_SCOPE_PROCESS) ; ALLOCATE(tids, pthread_t, nthread) ; MARKTIME(t1) ; for ( myid = 0, obj = MTmvmObjs ; myid < nthread ; myid++, obj++ ) { switch ( flag ) { case NONSYM : rc = pthread_create(&tids[myid], &attr, worker_nonsym_mmm, obj) ; break ; case SYM : rc = pthread_create(&tids[myid], &attr, worker_sym_mmm, obj) ; break ; case HERM : rc = pthread_create(&tids[myid], &attr, worker_herm_mmm, obj) ; break ; case NONSYM_T : rc = pthread_create(&tids[myid], &attr, worker_nonsym_mmm_T, obj); break ; case NONSYM_H : rc = pthread_create(&tids[myid], &attr, worker_nonsym_mmm_H, obj); break ; } if ( rc != 0 ) { fprintf(stderr, "\n fatal error, myid = %d, rc = %d from pthread_create", myid, rc) ; exit(-1) ; } else if ( msglvl > 2 ) { fprintf(stderr, "\n %% thread %d created", myid) ; } } MARKTIME(t2) ; if ( msglvl > 0 ) { fprintf(msgFile, "\n %% CPU %8.3f : thread creation time", t2 - t1) ; } MARKTIME(t1) ; for ( myid = 0 ; myid < nthread ; myid++ ) { pthread_join(tids[myid], &status) ; } MARKTIME(t2) ; if ( msglvl > 0 ) { fprintf(msgFile, "\n %% CPU %8.3f : thread join time", t2 - t1) ; } FREE(tids) ; pthread_attr_destroy(&attr) ; } #endif /* ------------------------------------- accumulate the rhs hand side matrices ------------------------------------- */ MARKTIME(t1) ; nent = Y->nrow * Y->ncol ; for ( myid = 1, obj = MTmvmObjs + 1 ; myid < nthread ; myid++, obj++ ) { if ( INPMTX_IS_REAL_ENTRIES(A) ) { DVadd(nent, DenseMtx_entries(Y), DenseMtx_entries(obj->Y)) ; } else if ( INPMTX_IS_COMPLEX_ENTRIES(A) ) { DVadd(2*nent, DenseMtx_entries(Y), DenseMtx_entries(obj->Y)) ; } } MARKTIME(t2) ; if ( msglvl > 0 ) { fprintf(msgFile, "\n %% CPU %8.3f : time to accumulate rhs", t2 - t1) ; } /* --------------------------- release the data structures --------------------------- */ MARKTIME(t1) ; for ( myid = 0, obj = MTmvmObjs ; myid < nthread ; myid++, obj++ ) { InpMtx_free(obj->A) ; if ( myid > 0 ) { DenseMtx_free(obj->Y) ; } } FREE(MTmvmObjs) ; MARKTIME(t2) ; if ( msglvl > 0 ) { fprintf(msgFile, "\n %% CPU %8.3f : time to release and free data", t2 - t1) ; } return ; }
/*--------------------------------------------------------------------*/ int main ( int argc, char *argv[] ) /* --------------------------------------------------- read in (i, j, a(i,j)) triples, construct a InpMtx object and write it out to a file created -- 97oct17, cca --------------------------------------------------- */ { char *inFileName, *outFileName ; InpMtx *inpmtx ; FILE *inputFile, *msgFile ; int dataType, flag, ient, msglvl, ncol, nent, nrow, rc ; int *ivec1, *ivec2 ; if ( argc != 7 ) { fprintf(stdout, "\n\n usage : readAIJ msglvl msgFile dataType inputFile outFile flag" "\n msglvl -- message level" "\n msgFile -- message file" "\n dataType -- 0 for indices only, 1 for double, 2 for complex" "\n inputFile -- input file for a(i,j) entries" "\n the first line must be \"nrow ncol nentries\"" "\n if dataType == 0 then" "\n next lines are \"irow jcol\"" "\n else if dataType == 1 then" "\n next lines are \"irow jcol entry\"" "\n else if dataType == 2 then" "\n next lines are \"irow jcol realEntry imagEntry\"" "\n endif" "\n outFile -- output file, must be *.inpmtxf or *.inpmtxb" "\n flag -- flag for 0-based or 1-based addressing" "\n") ; return(0) ; } msglvl = atoi(argv[1]) ; if ( strcmp(argv[2], "stdout") == 0 ) { msgFile = stdout ; } else if ( (msgFile = fopen(argv[2], "a")) == NULL ) { fprintf(stderr, "\n fatal error in %s" "\n unable to open file %s\n", argv[0], argv[2]) ; return(-1) ; } dataType = atoi(argv[3]) ; inFileName = argv[4] ; outFileName = argv[5] ; flag = atoi(argv[6]) ; fprintf(msgFile, "\n readAIJ " "\n msglvl -- %d" "\n msgFile -- %s" "\n dataType -- %d" "\n inputFile -- %s" "\n outFile -- %s" "\n flag -- %d" "\n", msglvl, argv[2], dataType, inFileName, outFileName, flag) ; fflush(msgFile) ; /* ---------------------------- open the input file and read #rows #columns #entries ---------------------------- */ if ( (inputFile = fopen(inFileName, "r")) == NULL ) { fprintf(stderr, "\n fatal error in %s" "\n unable to open file %s\n", argv[0], inFileName) ; return(-1) ; } rc = fscanf(inputFile, "%d %d %d", &nrow, &ncol, &nent) ; if ( rc != 3 ) { fprintf(stderr, "\n fatal error in %s" "\n %d of 3 fields read on first line of file %s", argv[0], rc, inFileName) ; return(-1) ; } if ( msglvl > 1 ) { fprintf(msgFile, "\n\n read in nrow = %d, ncol = %d, nent = %d", nrow, ncol, nent) ; fflush(msgFile) ; } /* -------------------------------------------------- initialize the object set coordType = INPMTX_BY_ROWS --> row coordinates set inputMode = dataType -------------------------------------------------- */ inpmtx = InpMtx_new() ; InpMtx_init(inpmtx, INPMTX_BY_ROWS, dataType, nent, 0) ; /* ------------------------------------------------- read in the entries and load them into the object ------------------------------------------------- */ ivec1 = InpMtx_ivec1(inpmtx) ; ivec2 = InpMtx_ivec2(inpmtx) ; if ( INPMTX_IS_INDICES_ONLY(inpmtx) ) { for ( ient = 0 ; ient < nent ; ient++ ) { rc = fscanf(inputFile, "%d %d", ivec1 + ient, ivec2 + ient) ; if ( rc != 2 ) { fprintf(stderr, "\n fatal error in %s" "\n %d of 2 fields read on entry %d of file %s", argv[0], rc, ient, inFileName) ; return(-1) ; } if ( msglvl > 1 ) { fprintf(msgFile, "\n entry %d, row %d, column %d", ient, ivec1[ient], ivec2[ient]) ; fflush(msgFile) ; } } } else if ( INPMTX_IS_REAL_ENTRIES(inpmtx) ) { double *dvec = InpMtx_dvec(inpmtx) ; for ( ient = 0 ; ient < nent ; ient++ ) { rc = fscanf(inputFile, "%d %d %le", ivec1 + ient, ivec2 + ient, dvec + ient) ; if ( rc != 3 ) { fprintf(stderr, "\n fatal error in %s" "\n %d of 3 fields read on entry %d of file %s", argv[0], rc, ient, argv[3]) ; return(-1) ; } if ( msglvl > 1 ) { fprintf(msgFile, "\n entry %d, row %d, column %d, value %e", ient, ivec1[ient], ivec2[ient], dvec[ient]) ; fflush(msgFile) ; } } } else if ( INPMTX_IS_COMPLEX_ENTRIES(inpmtx) ) { double *dvec = InpMtx_dvec(inpmtx) ; for ( ient = 0 ; ient < nent ; ient++ ) { rc = fscanf(inputFile, "%d %d %le %le", ivec1 + ient, ivec2 + ient, dvec + 2*ient, dvec + 2*ient+1) ; if ( rc != 4 ) { fprintf(stderr, "\n fatal error in %s" "\n %d of 4 fields read on entry %d of file %s", argv[0], rc, ient, argv[3]) ; return(-1) ; } if ( msglvl > 1 ) { fprintf(msgFile, "\n entry %d, row %d, column %d, value %12.4e + %12.4e*i", ient, ivec1[ient], ivec2[ient], dvec[2*ient], dvec[2*ient+1]) ; fflush(msgFile) ; } } } inpmtx->nent = nent ; if ( flag == 1 ) { /* -------------------------------------------------- indices were in FORTRAN mode, decrement for C mode -------------------------------------------------- */ for ( ient = 0 ; ient < nent ; ient++ ) { ivec1[ient]-- ; ivec2[ient]-- ; } } /* ----------------------------- sort and compress the entries ----------------------------- */ InpMtx_changeStorageMode(inpmtx, 3) ; if ( msglvl > 1 ) { fprintf(msgFile, "\n\n sorted, compressed and vector form") ; InpMtx_writeForHumanEye(inpmtx, msgFile) ; fflush(msgFile) ; } /* --------------------------- write out the InpMtx object --------------------------- */ if ( strcmp(outFileName, "none") != 0 ) { rc = InpMtx_writeToFile(inpmtx, outFileName) ; fprintf(msgFile, "\n return value %d from InpMtx_writeToFile(%p,%s)", rc, inpmtx, outFileName) ; } /* --------------------- free the working data --------------------- */ InpMtx_free(inpmtx) ; fprintf(msgFile, "\n") ; fclose(msgFile) ; return(1) ; }
/* ---------------------------------- drop entries in the upper triangle created -- 98jan28, cca ---------------------------------- */ void InpMtx_dropUpperTriangle ( InpMtx *inpmtx ) { double *dvec ; int count, ii, nent ; int *ivec1, *ivec2 ; /* --------------- check the input --------------- */ if ( inpmtx == NULL ) { fprintf(stderr, "\n fatal error in InpMtx_dropUpperTriangle(%p)" "\n bad input\n", inpmtx) ; exit(-1) ; } if ( !( INPMTX_IS_BY_ROWS(inpmtx) || INPMTX_IS_BY_COLUMNS(inpmtx) || INPMTX_IS_BY_CHEVRONS(inpmtx) ) ) { fprintf(stderr, "\n fatal error in InpMtx_dropUpperTriangle(%p)" "\n bad coordType \n", inpmtx) ; exit(-1) ; } nent = inpmtx->nent ; ivec1 = InpMtx_ivec1(inpmtx) ; ivec2 = InpMtx_ivec2(inpmtx) ; count = 0 ; if ( INPMTX_IS_REAL_ENTRIES(inpmtx) || INPMTX_IS_COMPLEX_ENTRIES(inpmtx) ) { dvec = InpMtx_dvec(inpmtx) ; } if ( INPMTX_IS_BY_ROWS(inpmtx) ) { for ( ii = 0 ; ii < nent ; ii++ ) { if ( ivec1[ii] >= ivec2[ii] ) { ivec1[count] = ivec1[ii] ; ivec2[count] = ivec2[ii] ; if ( INPMTX_IS_REAL_ENTRIES(inpmtx) ) { dvec[count] = dvec[ii] ; } else if ( INPMTX_IS_COMPLEX_ENTRIES(inpmtx) ) { dvec[2*count] = dvec[2*ii] ; dvec[2*count+1] = dvec[2*ii+1] ; } count++ ; } } } else if ( INPMTX_IS_BY_COLUMNS(inpmtx) ) { for ( ii = 0 ; ii < nent ; ii++ ) { if ( ivec1[ii] <= ivec2[ii] ) { ivec1[count] = ivec1[ii] ; ivec2[count] = ivec2[ii] ; if ( INPMTX_IS_REAL_ENTRIES(inpmtx) ) { dvec[count] = dvec[ii] ; } else if ( INPMTX_IS_COMPLEX_ENTRIES(inpmtx) ) { dvec[2*count] = dvec[2*ii] ; dvec[2*count+1] = dvec[2*ii+1] ; } count++ ; } } } else if ( INPMTX_IS_BY_CHEVRONS(inpmtx) ) { for ( ii = 0 ; ii < nent ; ii++ ) { if ( ivec2[ii] <= 0 ) { ivec1[count] = ivec1[ii] ; ivec2[count] = ivec2[ii] ; if ( INPMTX_IS_REAL_ENTRIES(inpmtx) ) { dvec[count] = dvec[ii] ; } else if ( INPMTX_IS_COMPLEX_ENTRIES(inpmtx) ) { dvec[2*count] = dvec[2*ii] ; dvec[2*count+1] = dvec[2*ii+1] ; } count++ ; } } } inpmtx->nent = count ; IV_setSize(&inpmtx->ivec1IV, count) ; IV_setSize(&inpmtx->ivec2IV, count) ; if ( INPMTX_IS_REAL_ENTRIES(inpmtx) || INPMTX_IS_COMPLEX_ENTRIES(inpmtx) ) { DV_setSize(&inpmtx->dvecDV, count) ; } return ; }
/* ----------------------- input a matrix created -- 98jan28, cca ----------------------- */ static void inputMatrix ( InpMtx *inpmtx, int nrow, int ncol, int rowstride, int colstride, int rowind[], int colind[], double mtxent[] ) { int col, ii, jj, kk, nent, row ; int *ivec1, *ivec2 ; prepareToAddNewEntries(inpmtx, nrow*ncol) ; nent = inpmtx->nent ; ivec1 = IV_entries(&inpmtx->ivec1IV) ; ivec2 = IV_entries(&inpmtx->ivec2IV) ; if ( INPMTX_IS_BY_ROWS(inpmtx) ) { for ( jj = 0, kk = nent ; jj < ncol ; jj++ ) { col = colind[jj] ; for ( ii = 0 ; ii < nrow ; ii++, kk++ ) { row = rowind[ii] ; ivec1[kk] = row ; ivec2[kk] = col ; } } } else if ( INPMTX_IS_BY_COLUMNS(inpmtx) ) { for ( jj = 0, kk = nent ; jj < ncol ; jj++ ) { col = colind[jj] ; for ( ii = 0 ; ii < nrow ; ii++, kk++ ) { row = rowind[ii] ; ivec1[kk] = col ; ivec2[kk] = row ; } } } else if ( INPMTX_IS_BY_CHEVRONS(inpmtx) ) { for ( jj = 0, kk = nent ; jj < ncol ; jj++ ) { col = colind[jj] ; for ( ii = 0 ; ii < nrow ; ii++, kk++ ) { row = rowind[ii] ; if ( row <= col ) { ivec1[kk] = row ; } else { ivec1[kk] = col ; } ivec2[kk] = col - row ; } } } IV_setSize(&inpmtx->ivec1IV, nent + nrow*ncol) ; IV_setSize(&inpmtx->ivec2IV, nent + nrow*ncol) ; if ( INPMTX_IS_REAL_ENTRIES(inpmtx) ) { double *dvec = DV_entries(&inpmtx->dvecDV) ; int ij ; for ( jj = 0, kk = nent ; jj < ncol ; jj++ ) { for ( ii = 0 ; ii < nrow ; ii++, kk++ ) { ij = ii*rowstride + jj*colstride ; dvec[kk] = mtxent[ij] ; } } DV_setSize(&inpmtx->dvecDV, nent + nrow*ncol) ; } if ( INPMTX_IS_COMPLEX_ENTRIES(inpmtx) ) { double *dvec = DV_entries(&inpmtx->dvecDV) ; int ij ; for ( jj = 0, kk = nent ; jj < ncol ; jj++ ) { for ( ii = 0 ; ii < nrow ; ii++, kk++ ) { ij = ii*rowstride + jj*colstride ; dvec[2*kk] = mtxent[2*ij] ; dvec[2*kk+1] = mtxent[2*ij+1] ; } } DV_setSize(&inpmtx->dvecDV, 2*(nent + nrow*ncol)) ; } inpmtx->nent += nrow*ncol ; inpmtx->storageMode = INPMTX_RAW_DATA ; return ; }