PetscErrorCode MatFactorNumeric_SeqSpooles(Mat F,Mat A,const MatFactorInfo *info) { Mat_Spooles *lu = (Mat_Spooles*)(F)->spptr; ChvManager *chvmanager ; Chv *rootchv ; IVL *adjIVL; PetscErrorCode ierr; PetscInt nz,nrow=A->rmap->n,irow,nedges,neqns=A->cmap->n,*ai,*aj,i,*diag=0,fierr; PetscScalar *av; double cputotal,facops; #if defined(PETSC_USE_COMPLEX) PetscInt nz_row,*aj_tmp; PetscScalar *av_tmp; #else PetscInt *ivec1,*ivec2,j; double *dvec; #endif PetscBool isSeqAIJ,isMPIAIJ; PetscFunctionBegin; if (lu->flg == DIFFERENT_NONZERO_PATTERN) { /* first numeric factorization */ (F)->ops->solve = MatSolve_SeqSpooles; (F)->assembled = PETSC_TRUE; /* set Spooles options */ ierr = SetSpoolesOptions(A, &lu->options);CHKERRQ(ierr); lu->mtxA = InpMtx_new(); } /* copy A to Spooles' InpMtx object */ ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr); ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isMPIAIJ);CHKERRQ(ierr); if (isSeqAIJ){ Mat_SeqAIJ *mat = (Mat_SeqAIJ*)A->data; ai=mat->i; aj=mat->j; av=mat->a; if (lu->options.symflag == SPOOLES_NONSYMMETRIC) { nz=mat->nz; } else { /* SPOOLES_SYMMETRIC || SPOOLES_HERMITIAN */ nz=(mat->nz + A->rmap->n)/2; diag=mat->diag; } } else { /* A is SBAIJ */ Mat_SeqSBAIJ *mat = (Mat_SeqSBAIJ*)A->data; ai=mat->i; aj=mat->j; av=mat->a; nz=mat->nz; } InpMtx_init(lu->mtxA, INPMTX_BY_ROWS, lu->options.typeflag, nz, 0); #if defined(PETSC_USE_COMPLEX) for (irow=0; irow<nrow; irow++) { if ( lu->options.symflag == SPOOLES_NONSYMMETRIC || !(isSeqAIJ || isMPIAIJ)){ nz_row = ai[irow+1] - ai[irow]; aj_tmp = aj + ai[irow]; av_tmp = av + ai[irow]; } else { nz_row = ai[irow+1] - diag[irow]; aj_tmp = aj + diag[irow]; av_tmp = av + diag[irow]; } for (i=0; i<nz_row; i++){ InpMtx_inputComplexEntry(lu->mtxA, irow, *aj_tmp++,PetscRealPart(*av_tmp),PetscImaginaryPart(*av_tmp)); av_tmp++; } } #else ivec1 = InpMtx_ivec1(lu->mtxA); ivec2 = InpMtx_ivec2(lu->mtxA); dvec = InpMtx_dvec(lu->mtxA); if ( lu->options.symflag == SPOOLES_NONSYMMETRIC || !isSeqAIJ){ for (irow = 0; irow < nrow; irow++){ for (i = ai[irow]; i<ai[irow+1]; i++) ivec1[i] = irow; } IVcopy(nz, ivec2, aj); DVcopy(nz, dvec, av); } else { nz = 0; for (irow = 0; irow < nrow; irow++){ for (j = diag[irow]; j<ai[irow+1]; j++) { ivec1[nz] = irow; ivec2[nz] = aj[j]; dvec[nz] = av[j]; nz++; } } } InpMtx_inputRealTriples(lu->mtxA, nz, ivec1, ivec2, dvec); #endif InpMtx_changeStorageMode(lu->mtxA, INPMTX_BY_VECTORS); if ( lu->options.msglvl > 0 ) { int err; printf("\n\n input matrix"); ierr = PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n input matrix");CHKERRQ(ierr); InpMtx_writeForHumanEye(lu->mtxA, lu->options.msgFile); err = fflush(lu->options.msgFile); if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fflush() failed on file"); } if ( lu->flg == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization */ /*--------------------------------------------------- find a low-fill ordering (1) create the Graph object (2) order the graph -------------------------------------------------------*/ if (lu->options.useQR){ adjIVL = InpMtx_adjForATA(lu->mtxA); } else { adjIVL = InpMtx_fullAdjacency(lu->mtxA); } nedges = IVL_tsize(adjIVL); lu->graph = Graph_new(); Graph_init2(lu->graph, 0, neqns, 0, nedges, neqns, nedges, adjIVL, NULL, NULL); if ( lu->options.msglvl > 2 ) { int err; if (lu->options.useQR){ ierr = PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n graph of A^T A");CHKERRQ(ierr); } else { ierr = PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n graph of the input matrix");CHKERRQ(ierr); } Graph_writeForHumanEye(lu->graph, lu->options.msgFile); err = fflush(lu->options.msgFile); if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fflush() failed on file"); } switch (lu->options.ordering) { case 0: lu->frontETree = orderViaBestOfNDandMS(lu->graph, lu->options.maxdomainsize, lu->options.maxzeros, lu->options.maxsize, lu->options.seed, lu->options.msglvl, lu->options.msgFile); break; case 1: lu->frontETree = orderViaMMD(lu->graph,lu->options.seed,lu->options.msglvl,lu->options.msgFile); break; case 2: lu->frontETree = orderViaMS(lu->graph, lu->options.maxdomainsize, lu->options.seed,lu->options.msglvl,lu->options.msgFile); break; case 3: lu->frontETree = orderViaND(lu->graph, lu->options.maxdomainsize, lu->options.seed,lu->options.msglvl,lu->options.msgFile); break; default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown Spooles's ordering"); } if ( lu->options.msglvl > 0 ) { int err; ierr = PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n front tree from ordering");CHKERRQ(ierr); ETree_writeForHumanEye(lu->frontETree, lu->options.msgFile); err = fflush(lu->options.msgFile); if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fflush() failed on file"); } /* get the permutation, permute the front tree */ lu->oldToNewIV = ETree_oldToNewVtxPerm(lu->frontETree); lu->oldToNew = IV_entries(lu->oldToNewIV); lu->newToOldIV = ETree_newToOldVtxPerm(lu->frontETree); if (!lu->options.useQR) ETree_permuteVertices(lu->frontETree, lu->oldToNewIV); /* permute the matrix */ if (lu->options.useQR){ InpMtx_permute(lu->mtxA, NULL, lu->oldToNew); } else { InpMtx_permute(lu->mtxA, lu->oldToNew, lu->oldToNew); if ( lu->options.symflag == SPOOLES_SYMMETRIC) { InpMtx_mapToUpperTriangle(lu->mtxA); } #if defined(PETSC_USE_COMPLEX) if ( lu->options.symflag == SPOOLES_HERMITIAN ) { InpMtx_mapToUpperTriangleH(lu->mtxA); } #endif InpMtx_changeCoordType(lu->mtxA, INPMTX_BY_CHEVRONS); } InpMtx_changeStorageMode(lu->mtxA, INPMTX_BY_VECTORS); /* get symbolic factorization */ if (lu->options.useQR){ lu->symbfacIVL = SymbFac_initFromGraph(lu->frontETree, lu->graph); IVL_overwrite(lu->symbfacIVL, lu->oldToNewIV); IVL_sortUp(lu->symbfacIVL); ETree_permuteVertices(lu->frontETree, lu->oldToNewIV); } else { lu->symbfacIVL = SymbFac_initFromInpMtx(lu->frontETree, lu->mtxA); } if ( lu->options.msglvl > 2 ) { int err; ierr = PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n old-to-new permutation vector");CHKERRQ(ierr); IV_writeForHumanEye(lu->oldToNewIV, lu->options.msgFile); ierr = PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n new-to-old permutation vector");CHKERRQ(ierr); IV_writeForHumanEye(lu->newToOldIV, lu->options.msgFile); ierr = PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n front tree after permutation");CHKERRQ(ierr); ETree_writeForHumanEye(lu->frontETree, lu->options.msgFile); ierr = PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n input matrix after permutation");CHKERRQ(ierr); InpMtx_writeForHumanEye(lu->mtxA, lu->options.msgFile); ierr = PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n symbolic factorization");CHKERRQ(ierr); IVL_writeForHumanEye(lu->symbfacIVL, lu->options.msgFile); err = fflush(lu->options.msgFile); if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fflush() failed on file"); } lu->frontmtx = FrontMtx_new(); lu->mtxmanager = SubMtxManager_new(); SubMtxManager_init(lu->mtxmanager, NO_LOCK, 0); } else { /* new num factorization using previously computed symbolic factor */ if (lu->options.pivotingflag) { /* different FrontMtx is required */ FrontMtx_free(lu->frontmtx); lu->frontmtx = FrontMtx_new(); } else { FrontMtx_clearData (lu->frontmtx); } SubMtxManager_free(lu->mtxmanager); lu->mtxmanager = SubMtxManager_new(); SubMtxManager_init(lu->mtxmanager, NO_LOCK, 0); /* permute mtxA */ if (lu->options.useQR){ InpMtx_permute(lu->mtxA, NULL, lu->oldToNew); } else { InpMtx_permute(lu->mtxA, lu->oldToNew, lu->oldToNew); if ( lu->options.symflag == SPOOLES_SYMMETRIC ) { InpMtx_mapToUpperTriangle(lu->mtxA); } InpMtx_changeCoordType(lu->mtxA, INPMTX_BY_CHEVRONS); } InpMtx_changeStorageMode(lu->mtxA, INPMTX_BY_VECTORS); if ( lu->options.msglvl > 2 ) { ierr = PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n input matrix after permutation");CHKERRQ(ierr); InpMtx_writeForHumanEye(lu->mtxA, lu->options.msgFile); } } /* end of if( lu->flg == DIFFERENT_NONZERO_PATTERN) */ if (lu->options.useQR){ FrontMtx_init(lu->frontmtx, lu->frontETree, lu->symbfacIVL, lu->options.typeflag, SPOOLES_SYMMETRIC, FRONTMTX_DENSE_FRONTS, SPOOLES_NO_PIVOTING, NO_LOCK, 0, NULL, lu->mtxmanager, lu->options.msglvl, lu->options.msgFile); } else { FrontMtx_init(lu->frontmtx, lu->frontETree, lu->symbfacIVL, lu->options.typeflag, lu->options.symflag, FRONTMTX_DENSE_FRONTS, lu->options.pivotingflag, NO_LOCK, 0, NULL, lu->mtxmanager, lu->options.msglvl, lu->options.msgFile); } if ( lu->options.symflag == SPOOLES_SYMMETRIC ) { /* || SPOOLES_HERMITIAN ? */ if ( lu->options.patchAndGoFlag == 1 ) { lu->frontmtx->patchinfo = PatchAndGoInfo_new(); PatchAndGoInfo_init(lu->frontmtx->patchinfo, 1, lu->options.toosmall, lu->options.fudge, lu->options.storeids, lu->options.storevalues); } else if ( lu->options.patchAndGoFlag == 2 ) { lu->frontmtx->patchinfo = PatchAndGoInfo_new(); PatchAndGoInfo_init(lu->frontmtx->patchinfo, 2, lu->options.toosmall, lu->options.fudge, lu->options.storeids, lu->options.storevalues); } } /* numerical factorization */ chvmanager = ChvManager_new(); ChvManager_init(chvmanager, NO_LOCK, 1); DVfill(10, lu->cpus, 0.0); if (lu->options.useQR){ facops = 0.0 ; FrontMtx_QR_factor(lu->frontmtx, lu->mtxA, chvmanager, lu->cpus, &facops, lu->options.msglvl, lu->options.msgFile); if ( lu->options.msglvl > 1 ) { ierr = PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n factor matrix");CHKERRQ(ierr); ierr = PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n facops = %9.2f", facops);CHKERRQ(ierr); } } else { IVfill(20, lu->stats, 0); rootchv = FrontMtx_factorInpMtx(lu->frontmtx, lu->mtxA, lu->options.tau, 0.0, chvmanager, &fierr, lu->cpus,lu->stats,lu->options.msglvl,lu->options.msgFile); if (rootchv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"\n matrix found to be singular"); if (fierr >= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"\n error encountered at front %D", fierr); if(lu->options.FrontMtxInfo){ ierr = PetscPrintf(PETSC_COMM_SELF,"\n %8d pivots, %8d pivot tests, %8d delayed rows and columns\n",lu->stats[0], lu->stats[1], lu->stats[2]);CHKERRQ(ierr); cputotal = lu->cpus[8] ; if ( cputotal > 0.0 ) { ierr = PetscPrintf(PETSC_COMM_SELF, "\n cpus cpus/totaltime" "\n initialize fronts %8.3f %6.2f" "\n load original entries %8.3f %6.2f" "\n update fronts %8.3f %6.2f" "\n assemble postponed data %8.3f %6.2f" "\n factor fronts %8.3f %6.2f" "\n extract postponed data %8.3f %6.2f" "\n store factor entries %8.3f %6.2f" "\n miscellaneous %8.3f %6.2f" "\n total time %8.3f \n", lu->cpus[0], 100.*lu->cpus[0]/cputotal, lu->cpus[1], 100.*lu->cpus[1]/cputotal, lu->cpus[2], 100.*lu->cpus[2]/cputotal, lu->cpus[3], 100.*lu->cpus[3]/cputotal, lu->cpus[4], 100.*lu->cpus[4]/cputotal, lu->cpus[5], 100.*lu->cpus[5]/cputotal, lu->cpus[6], 100.*lu->cpus[6]/cputotal, lu->cpus[7], 100.*lu->cpus[7]/cputotal, cputotal);CHKERRQ(ierr); } } } ChvManager_free(chvmanager); if ( lu->options.msglvl > 0 ) { int err; ierr = PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n factor matrix");CHKERRQ(ierr); FrontMtx_writeForHumanEye(lu->frontmtx, lu->options.msgFile); err = fflush(lu->options.msgFile); if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fflush() failed on file"); } if ( lu->options.symflag == SPOOLES_SYMMETRIC ) { /* || SPOOLES_HERMITIAN ? */ if ( lu->options.patchAndGoFlag == 1 ) { if ( lu->frontmtx->patchinfo->fudgeIV != NULL ) { if (lu->options.msglvl > 0 ){ ierr = PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n small pivots found at these locations");CHKERRQ(ierr); IV_writeForHumanEye(lu->frontmtx->patchinfo->fudgeIV, lu->options.msgFile); } } PatchAndGoInfo_free(lu->frontmtx->patchinfo); } else if ( lu->options.patchAndGoFlag == 2 ) { if (lu->options.msglvl > 0 ){ if ( lu->frontmtx->patchinfo->fudgeIV != NULL ) { ierr = PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n small pivots found at these locations");CHKERRQ(ierr); IV_writeForHumanEye(lu->frontmtx->patchinfo->fudgeIV, lu->options.msgFile); } if ( lu->frontmtx->patchinfo->fudgeDV != NULL ) { ierr = PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n perturbations");CHKERRQ(ierr); DV_writeForHumanEye(lu->frontmtx->patchinfo->fudgeDV, lu->options.msgFile); } } PatchAndGoInfo_free(lu->frontmtx->patchinfo); } } /* post-process the factorization */ FrontMtx_postProcess(lu->frontmtx, lu->options.msglvl, lu->options.msgFile); if ( lu->options.msglvl > 2 ) { int err; ierr = PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n factor matrix after post-processing");CHKERRQ(ierr); FrontMtx_writeForHumanEye(lu->frontmtx, lu->options.msgFile); err = fflush(lu->options.msgFile); if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fflush() failed on file"); } lu->flg = SAME_NONZERO_PATTERN; lu->CleanUpSpooles = PETSC_TRUE; PetscFunctionReturn(0); }
/* ------------------------------------------------------------------ purpose -- return an ETree object for a nested dissection ordering graph -- graph to order maxdomainsize -- used to control the incomplete nested dissection process. any subgraph whose weight is less than maxdomainsize is not split further. seed -- random number seed msglvl -- message level, 0 --> no output, 1 --> timings msgFile -- message file created -- 97nov06, cca ------------------------------------------------------------------ */ ETree * orderViaND ( Graph *graph, int maxdomainsize, int seed, int msglvl, FILE *msgFile ) { double t1, t2 ; DSTree *dstree ; ETree *etree ; int nvtx, Nvtx ; IV *eqmapIV, *stagesIV ; /* --------------- check the input --------------- */ if ( graph == NULL || maxdomainsize <= 0 || (msglvl > 0 && msgFile == NULL) ) { fprintf(stderr, "\n fatal error in orderViaND(%p,%d,%d,%d,%p)" "\n bad input\n", graph, maxdomainsize, seed, msglvl, msgFile) ; exit(-1) ; } /* ------------------------------ compress the graph if worth it ------------------------------ */ nvtx = graph->nvtx ; MARKTIME(t1) ; eqmapIV = Graph_equivMap(graph) ; MARKTIME(t2) ; if ( msglvl > 0 ) { fprintf(msgFile, "\n CPU %8.3f : get equivalence map", t2 - t1) ; fflush(msgFile) ; } Nvtx = 1 + IV_max(eqmapIV) ; if ( Nvtx <= COMPRESS_FRACTION * nvtx ) { MARKTIME(t1) ; graph = Graph_compress2(graph, eqmapIV, 1) ; MARKTIME(t2) ; if ( msglvl > 0 ) { fprintf(msgFile, "\n CPU %8.3f : compress graph", t2 - t1) ; fflush(msgFile) ; } } else { IV_free(eqmapIV) ; eqmapIV = NULL ; } MARKTIME(t1) ; IVL_sortUp(graph->adjIVL) ; MARKTIME(t2) ; if ( msglvl > 0 ) { fprintf(msgFile, "\n CPU %8.3f : sort adjacency", t2 - t1) ; fflush(msgFile) ; } /* ----------------------------- get the domain separator tree ----------------------------- */ { GPart *gpart ; DDsepInfo *info ; info = DDsepInfo_new() ; info->seed = seed ; info->maxcompweight = maxdomainsize ; info->alpha = 0.1 ; gpart = GPart_new() ; GPart_init(gpart, graph) ; GPart_setMessageInfo(gpart, msglvl, msgFile) ; dstree = GPart_RBviaDDsep(gpart, info) ; DSTree_renumberViaPostOT(dstree) ; if ( msglvl > 0 ) { DDsepInfo_writeCpuTimes(info, msgFile) ; } DDsepInfo_free(info) ; GPart_free(gpart) ; } /* --------------------- get the stages vector --------------------- */ stagesIV = DSTree_NDstages(dstree) ; DSTree_free(dstree) ; /* --------------------------------------------- order the vertices and extract the front tree --------------------------------------------- */ { MSMDinfo *info ; MSMD *msmd ; info = MSMDinfo_new() ; info->seed = seed ; info->compressFlag = 2 ; info->msglvl = msglvl ; info->msgFile = msgFile ; msmd = MSMD_new() ; MSMD_order(msmd, graph, IV_entries(stagesIV), info) ; etree = MSMD_frontETree(msmd) ; if ( msglvl > 0 ) { MSMDinfo_print(info, msgFile) ; } MSMDinfo_free(info) ; MSMD_free(msmd) ; IV_free(stagesIV) ; } /* ------------------------------------------------- expand the front tree if the graph was compressed ------------------------------------------------- */ if ( eqmapIV != NULL ) { ETree *etree2 = ETree_expand(etree, eqmapIV) ; ETree_free(etree) ; etree = etree2 ; Graph_free(graph) ; IV_free(eqmapIV) ; } else { MARKTIME(t1) ; IVL_sortUp(graph->adjIVL) ; MARKTIME(t2) ; if ( msglvl > 0 ) { fprintf(msgFile, "\n CPU %8.3f : sort adjacency", t2 - t1) ; fflush(msgFile) ; } } return(etree) ; }
/* -------------------------------------------------------- purpose -- return an ETree object for a multiple minimum degree ordering graph -- graph to order seed -- random number seed msglvl -- message level, 0 --> no output, 1 --> timings msgFile -- message file created -- 97nov08, cca -------------------------------------------------------- */ ETree * orderViaMMD ( Graph *graph, int seed, int msglvl, FILE *msgFile ) { double t1, t2 ; ETree *etree ; int nvtx, Nvtx ; IV *eqmapIV ; /* --------------- check the input --------------- */ if ( graph == NULL || (msglvl > 0 && msgFile == NULL) ) { fprintf(stderr, "\n fatal error in orderViaMMD(%p,%d,%d,%p)" "\n bad input\n", graph, seed, msglvl, msgFile) ; exit(-1) ; } /* ------------------------------ compress the graph if worth it ------------------------------ */ nvtx = graph->nvtx ; MARKTIME(t1) ; eqmapIV = Graph_equivMap(graph) ; MARKTIME(t2) ; if ( msglvl > 1 ) { fprintf(msgFile, "\n CPU %8.3f : get equivalence map", t2 - t1) ; fflush(msgFile) ; } Nvtx = 1 + IV_max(eqmapIV) ; if ( Nvtx <= COMPRESS_FRACTION * nvtx ) { MARKTIME(t1) ; graph = Graph_compress2(graph, eqmapIV, 1) ; MARKTIME(t2) ; if ( msglvl > 1 ) { fprintf(msgFile, "\n CPU %8.3f : compress graph", t2 - t1) ; fflush(msgFile) ; } } else { IV_free(eqmapIV) ; eqmapIV = NULL ; } MARKTIME(t1) ; IVL_sortUp(graph->adjIVL) ; MARKTIME(t2) ; if ( msglvl > 1 ) { fprintf(msgFile, "\n CPU %8.3f : sort adjacency", t2 - t1) ; fflush(msgFile) ; } /* --------------------------------------------- order the vertices and extract the front tree --------------------------------------------- */ { MSMDinfo *info ; MSMD *msmd ; info = MSMDinfo_new() ; info->seed = seed ; info->compressFlag = 2 ; info->msglvl = msglvl ; info->msgFile = msgFile ; msmd = MSMD_new() ; MSMD_order(msmd, graph, NULL, info) ; etree = MSMD_frontETree(msmd) ; if ( msglvl > 1 ) { MSMDinfo_print(info, msgFile) ; } MSMDinfo_free(info) ; MSMD_free(msmd) ; } /* ------------------------------------------------- expand the front tree if the graph was compressed ------------------------------------------------- */ if ( eqmapIV != NULL ) { ETree *etree2 = ETree_expand(etree, eqmapIV) ; ETree_free(etree) ; etree = etree2 ; Graph_free(graph) ; IV_free(eqmapIV) ; } else { MARKTIME(t1) ; IVL_sortUp(graph->adjIVL) ; MARKTIME(t2) ; if ( msglvl > 1 ) { fprintf(msgFile, "\n CPU %8.3f : sort adjacency", t2 - t1) ; fflush(msgFile) ; } } return(etree) ; }
/*--------------------------------------------------------------------*/ int main ( int argc, char *argv[] ) { /* -------------------------------------------------- QR all-in-one program (1) read in matrix entries and form InpMtx object of A and A^TA (2) form Graph object of A^TA (3) order matrix and form front tree (4) get the permutation, permute the matrix and front tree and get the symbolic factorization (5) compute the numeric factorization (6) read in right hand side entries (7) compute the solution created -- 98jun11, cca -------------------------------------------------- */ /*--------------------------------------------------------------------*/ char *matrixFileName, *rhsFileName ; ChvManager *chvmanager ; DenseMtx *mtxB, *mtxX ; double facops, imag, real, value ; double cpus[10] ; ETree *frontETree ; FILE *inputFile, *msgFile ; FrontMtx *frontmtx ; Graph *graph ; int ient, irow, jcol, jrhs, jrow, msglvl, neqns, nedges, nent, nrhs, nrow, seed, type ; InpMtx *mtxA ; IV *newToOldIV, *oldToNewIV ; IVL *adjIVL, *symbfacIVL ; SubMtxManager *mtxmanager ; /*--------------------------------------------------------------------*/ /* -------------------- get input parameters -------------------- */ if ( argc != 7 ) { fprintf(stdout, "\n usage: %s msglvl msgFile type matrixFileName rhsFileName seed" "\n msglvl -- message level" "\n msgFile -- message file" "\n type -- type of entries" "\n 1 (SPOOLES_REAL) -- real entries" "\n 2 (SPOOLES_COMPLEX) -- complex entries" "\n matrixFileName -- matrix file name, format" "\n nrow ncol nent" "\n irow jcol entry" "\n ..." "\n note: indices are zero based" "\n rhsFileName -- right hand side file name, format" "\n nrow " "\n entry[0]" "\n ..." "\n entry[nrow-1]" "\n seed -- random number seed, used for ordering" "\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) ; } type = atoi(argv[3]) ; matrixFileName = argv[4] ; rhsFileName = argv[5] ; seed = atoi(argv[6]) ; /*--------------------------------------------------------------------*/ /* -------------------------------------------- STEP 1: read the entries from the input file and create the InpMtx object of A -------------------------------------------- */ inputFile = fopen(matrixFileName, "r") ; fscanf(inputFile, "%d %d %d", &nrow, &neqns, &nent) ; mtxA = InpMtx_new() ; InpMtx_init(mtxA, INPMTX_BY_ROWS, type, nent, 0) ; if ( type == SPOOLES_REAL ) { for ( ient = 0 ; ient < nent ; ient++ ) { fscanf(inputFile, "%d %d %le", &irow, &jcol, &value) ; InpMtx_inputRealEntry(mtxA, irow, jcol, value) ; } } else { for ( ient = 0 ; ient < nent ; ient++ ) { fscanf(inputFile, "%d %d %le %le", &irow, &jcol, &real, &imag) ; InpMtx_inputComplexEntry(mtxA, irow, jcol, real, imag) ; } } fclose(inputFile) ; if ( msglvl > 1 ) { fprintf(msgFile, "\n\n input matrix") ; InpMtx_writeForHumanEye(mtxA, msgFile) ; fflush(msgFile) ; } /*--------------------------------------------------------------------*/ /* ---------------------------------------- STEP 2: read the right hand side entries ---------------------------------------- */ inputFile = fopen(rhsFileName, "r") ; fscanf(inputFile, "%d %d", &nrow, &nrhs) ; mtxB = DenseMtx_new() ; DenseMtx_init(mtxB, type, 0, 0, nrow, nrhs, 1, nrow) ; DenseMtx_zero(mtxB) ; if ( type == SPOOLES_REAL ) { for ( irow = 0 ; irow < nrow ; irow++ ) { fscanf(inputFile, "%d", &jrow) ; for ( jrhs = 0 ; jrhs < nrhs ; jrhs++ ) { fscanf(inputFile, "%le", &value) ; DenseMtx_setRealEntry(mtxB, jrow, jrhs, value) ; } } } else { for ( irow = 0 ; irow < nrow ; irow++ ) { fscanf(inputFile, "%d", &jrow) ; for ( jrhs = 0 ; jrhs < nrhs ; jrhs++ ) { fscanf(inputFile, "%le %le", &real, &imag) ; DenseMtx_setComplexEntry(mtxB, jrow, jrhs, real, imag) ; } } } fclose(inputFile) ; if ( msglvl > 1 ) { fprintf(msgFile, "\n\n rhs matrix in original ordering") ; DenseMtx_writeForHumanEye(mtxB, msgFile) ; fflush(msgFile) ; } /*--------------------------------------------------------------------*/ /* ------------------------------------------------- STEP 3 : find a low-fill ordering (1) create the Graph object for A^TA or A^HA (2) order the graph using multiple minimum degree ------------------------------------------------- */ graph = Graph_new() ; adjIVL = InpMtx_adjForATA(mtxA) ; nedges = IVL_tsize(adjIVL) ; Graph_init2(graph, 0, neqns, 0, nedges, neqns, nedges, adjIVL, NULL, NULL) ; if ( msglvl > 1 ) { fprintf(msgFile, "\n\n graph of A^T A") ; Graph_writeForHumanEye(graph, msgFile) ; fflush(msgFile) ; } frontETree = orderViaMMD(graph, seed, msglvl, msgFile) ; if ( msglvl > 1 ) { fprintf(msgFile, "\n\n front tree from ordering") ; ETree_writeForHumanEye(frontETree, msgFile) ; fflush(msgFile) ; } /*--------------------------------------------------------------------*/ /* ----------------------------------------------------- STEP 4: get the permutation, permute the matrix and front tree and get the symbolic factorization ----------------------------------------------------- */ oldToNewIV = ETree_oldToNewVtxPerm(frontETree) ; newToOldIV = ETree_newToOldVtxPerm(frontETree) ; InpMtx_permute(mtxA, NULL, IV_entries(oldToNewIV)) ; InpMtx_changeStorageMode(mtxA, INPMTX_BY_VECTORS) ; symbfacIVL = SymbFac_initFromGraph(frontETree, graph) ; IVL_overwrite(symbfacIVL, oldToNewIV) ; IVL_sortUp(symbfacIVL) ; ETree_permuteVertices(frontETree, oldToNewIV) ; if ( msglvl > 1 ) { fprintf(msgFile, "\n\n old-to-new permutation vector") ; IV_writeForHumanEye(oldToNewIV, msgFile) ; fprintf(msgFile, "\n\n new-to-old permutation vector") ; IV_writeForHumanEye(newToOldIV, msgFile) ; fprintf(msgFile, "\n\n front tree after permutation") ; ETree_writeForHumanEye(frontETree, msgFile) ; fprintf(msgFile, "\n\n input matrix after permutation") ; InpMtx_writeForHumanEye(mtxA, msgFile) ; fprintf(msgFile, "\n\n symbolic factorization") ; IVL_writeForHumanEye(symbfacIVL, msgFile) ; fflush(msgFile) ; } /*--------------------------------------------------------------------*/ /* ------------------------------------------ STEP 5: initialize the front matrix object ------------------------------------------ */ frontmtx = FrontMtx_new() ; mtxmanager = SubMtxManager_new() ; SubMtxManager_init(mtxmanager, NO_LOCK, 0) ; if ( type == SPOOLES_REAL ) { FrontMtx_init(frontmtx, frontETree, symbfacIVL, type, SPOOLES_SYMMETRIC, FRONTMTX_DENSE_FRONTS, SPOOLES_NO_PIVOTING, NO_LOCK, 0, NULL, mtxmanager, msglvl, msgFile) ; } else { FrontMtx_init(frontmtx, frontETree, symbfacIVL, type, SPOOLES_HERMITIAN, FRONTMTX_DENSE_FRONTS, SPOOLES_NO_PIVOTING, NO_LOCK, 0, NULL, mtxmanager, msglvl, msgFile) ; } /*--------------------------------------------------------------------*/ /* ----------------------------------------- STEP 6: compute the numeric factorization ----------------------------------------- */ chvmanager = ChvManager_new() ; ChvManager_init(chvmanager, NO_LOCK, 1) ; DVzero(10, cpus) ; facops = 0.0 ; FrontMtx_QR_factor(frontmtx, mtxA, chvmanager, cpus, &facops, msglvl, msgFile) ; ChvManager_free(chvmanager) ; if ( msglvl > 1 ) { fprintf(msgFile, "\n\n factor matrix") ; fprintf(msgFile, "\n facops = %9.2f", facops) ; FrontMtx_writeForHumanEye(frontmtx, msgFile) ; fflush(msgFile) ; } /*--------------------------------------------------------------------*/ /* -------------------------------------- STEP 7: post-process the factorization -------------------------------------- */ FrontMtx_postProcess(frontmtx, msglvl, msgFile) ; if ( msglvl > 1 ) { fprintf(msgFile, "\n\n factor matrix after post-processing") ; FrontMtx_writeForHumanEye(frontmtx, msgFile) ; fflush(msgFile) ; } /*--------------------------------------------------------------------*/ /* ------------------------------- STEP 8: solve the linear system ------------------------------- */ mtxX = DenseMtx_new() ; DenseMtx_init(mtxX, type, 0, 0, neqns, nrhs, 1, neqns) ; FrontMtx_QR_solve(frontmtx, mtxA, mtxX, mtxB, mtxmanager, cpus, msglvl, msgFile) ; if ( msglvl > 1 ) { fprintf(msgFile, "\n\n solution matrix in new ordering") ; DenseMtx_writeForHumanEye(mtxX, msgFile) ; fflush(msgFile) ; } /*--------------------------------------------------------------------*/ /* ------------------------------------------------------- STEP 9: permute the solution into the original ordering ------------------------------------------------------- */ DenseMtx_permuteRows(mtxX, newToOldIV) ; if ( msglvl > 0 ) { fprintf(msgFile, "\n\n solution matrix in original ordering") ; DenseMtx_writeForHumanEye(mtxX, msgFile) ; fflush(msgFile) ; } /*--------------------------------------------------------------------*/ /* ------------------------ free the working storage ------------------------ */ InpMtx_free(mtxA) ; FrontMtx_free(frontmtx) ; Graph_free(graph) ; DenseMtx_free(mtxX) ; DenseMtx_free(mtxB) ; ETree_free(frontETree) ; IV_free(newToOldIV) ; IV_free(oldToNewIV) ; IVL_free(symbfacIVL) ; SubMtxManager_free(mtxmanager) ; /*--------------------------------------------------------------------*/ return(1) ; }
/* ------------------------------------------------------------------- purpose -- given an InpMtx object that contains the structure of A, initialize the bridge data structure for the serial factor's and solve's. note: all parameters are pointers to be compatible with fortran's call by reference. return value -- 1 -- normal return -1 -- bridge is NULL -2 -- mtxA is NULL created -- 98sep17, cca ------------------------------------------------------------------- */ int BridgeMT_setup ( BridgeMT *bridge, InpMtx *mtxA ) { double t0, t1, t2 ; ETree *frontETree ; FILE *msgFile ; Graph *graph ; int compressed, msglvl, nedges, neqns, Neqns ; IV *eqmapIV ; IVL *adjIVL, *symbfacIVL ; MARKTIME(t0) ; /* -------------------- check the input data -------------------- */ if ( bridge == NULL ) { fprintf(stderr, "\n fatal error in BridgeMT_setup()" "\n data is NULL\n") ; return(-1) ; } if ( mtxA == NULL ) { fprintf(stderr, "\n fatal error in BridgeMT_setup()" "\n A is NULL\n") ; return(-2) ; } msglvl = bridge->msglvl ; msgFile = bridge->msgFile ; neqns = bridge->neqns ; if ( ! (INPMTX_IS_BY_ROWS(mtxA) || INPMTX_IS_BY_COLUMNS(mtxA)) ) { /* ------------------------------ change coordinate type to rows ------------------------------ */ InpMtx_changeCoordType(mtxA, INPMTX_BY_ROWS) ; } if ( ! INPMTX_IS_BY_VECTORS(mtxA) ) { /* ------------------------------ change storage mode to vectors ------------------------------ */ InpMtx_changeStorageMode(mtxA, INPMTX_BY_VECTORS) ; } /* --------------------------- create a Graph object for A --------------------------- */ MARKTIME(t1) ; graph = Graph_new() ; adjIVL = InpMtx_fullAdjacency(mtxA); nedges = bridge->nedges = IVL_tsize(adjIVL), Graph_init2(graph, 0, neqns, 0, nedges, neqns, nedges, adjIVL, NULL, NULL) ; MARKTIME(t2) ; bridge->cpus[0] += t2 - t1 ; if ( msglvl > 1 ) { fprintf(msgFile, "\n CPU %8.3f : time to create Graph", t2 - t1) ; fflush(msgFile) ; } if ( msglvl > 3 ) { fprintf(msgFile, "\n\n graph of the input matrix") ; Graph_writeForHumanEye(graph, msgFile) ; fflush(msgFile) ; } /* ------------------ compress the graph ------------------ */ MARKTIME(t1) ; eqmapIV = Graph_equivMap(graph) ; Neqns = bridge->Neqns = 1 + IV_max(eqmapIV) ; if ( msglvl > 2 ) { fprintf(msgFile, "\n\n graph's equivalence map") ; IV_writeForHumanEye(eqmapIV, msgFile) ; fflush(msgFile) ; } if ( Neqns < bridge->compressCutoff * neqns ) { Graph *cgraph ; /* ------------------ compress the graph ------------------ */ cgraph = Graph_compress2(graph, eqmapIV, 1) ; Graph_free(graph) ; graph = cgraph ; compressed = 1 ; bridge->Nedges = graph->nedges ; } else { compressed = 0 ; } MARKTIME(t2) ; bridge->cpus[1] += t2 - t1 ; if ( msglvl > 1 ) { fprintf(msgFile, "\n CPU %8.3f : time to create compressed graph", t2 - t1) ; fflush(msgFile) ; } if ( msglvl > 3 ) { fprintf(msgFile, "\n\n graph to order") ; Graph_writeForHumanEye(graph, msgFile) ; fflush(msgFile) ; } /* --------------- order the graph --------------- */ MARKTIME(t1) ; if ( bridge->maxdomainsize <= 0 ) { bridge->maxdomainsize = neqns/32 ; } if ( bridge->maxdomainsize <= 0 ) { bridge->maxdomainsize = 1 ; } if ( bridge->maxnzeros < 0 ) { bridge->maxnzeros = 0.01*neqns ; } if ( bridge->maxsize < 0 ) { bridge->maxsize = neqns ; } frontETree = orderViaBestOfNDandMS(graph, bridge->maxdomainsize, bridge->maxnzeros, bridge->maxsize, bridge->seed, msglvl, msgFile) ; MARKTIME(t2) ; bridge->cpus[2] += t2 - t1 ; if ( msglvl > 1 ) { fprintf(msgFile, "\n CPU %8.3f : time to order graph", t2 - t1) ; fflush(msgFile) ; } if ( msglvl > 3 ) { fprintf(msgFile, "\n\n front tree from ordering") ; ETree_writeForHumanEye(frontETree, msgFile) ; fflush(msgFile) ; } MARKTIME(t1) ; if ( compressed == 1 ) { ETree *etree ; IVL *tempIVL ; /* ---------------------------------------------------------- compute the symbolic factorization of the compressed graph ---------------------------------------------------------- */ tempIVL = SymbFac_initFromGraph(frontETree, graph) ; /* ------------------------------------------------------- expand the symbolic factorization to the original graph ------------------------------------------------------- */ symbfacIVL = IVL_expand(tempIVL, eqmapIV) ; IVL_free(tempIVL) ; /* --------------------- expand the front tree --------------------- */ etree = ETree_expand(frontETree, eqmapIV) ; ETree_free(frontETree) ; frontETree = etree ; } else { /* -------------------------------------------------------- compute the symbolic factorization of the original graph -------------------------------------------------------- */ symbfacIVL = SymbFac_initFromGraph(frontETree, graph) ; } MARKTIME(t2) ; bridge->frontETree = frontETree ; bridge->symbfacIVL = symbfacIVL ; /* ---------------------------------------------- get the old-to-new and new-to-old permutations ---------------------------------------------- */ bridge->oldToNewIV = ETree_oldToNewVtxPerm(frontETree) ; bridge->newToOldIV = ETree_newToOldVtxPerm(frontETree) ; if ( msglvl > 2 ) { fprintf(msgFile, "\n\n old-to-new permutation") ; IV_writeForHumanEye(bridge->oldToNewIV, msgFile) ; fprintf(msgFile, "\n\n new-to-old permutation") ; IV_writeForHumanEye(bridge->newToOldIV, msgFile) ; fflush(msgFile) ; } /* ------------------------------------------------------ overwrite the symbolic factorization with the permuted indices and sort the lists into ascending order ------------------------------------------------------ */ IVL_overwrite(symbfacIVL, bridge->oldToNewIV) ; IVL_sortUp(symbfacIVL) ; if ( msglvl > 2 ) { fprintf(msgFile, "\n\n symbolic factorization") ; IVL_writeForHumanEye(symbfacIVL, msgFile) ; fflush(msgFile) ; } /* -------------------------------------- permute the vertices in the front tree -------------------------------------- */ ETree_permuteVertices(frontETree, bridge->oldToNewIV) ; if ( msglvl > 2 ) { fprintf(msgFile, "\n\n permuted front etree") ; ETree_writeForHumanEye(frontETree, msgFile) ; fflush(msgFile) ; } MARKTIME(t2) ; bridge->cpus[3] += t2 - t1 ; if ( msglvl > 1 ) { fprintf(msgFile, "\n CPU %8.3f : time for symbolic factorization", t2 - t1) ; fflush(msgFile) ; } /* ------------------------ free the working storage ------------------------ */ Graph_free(graph) ; IV_free(eqmapIV) ; MARKTIME(t2) ; bridge->cpus[4] += t2 - t0 ; return(1) ; }