/* --------------------------------------------------------------- purpose -- fill dvec[J] with the stack storage to solve for J in a forward solve created -- 97nov30, cca --------------------------------------------------------------- */ void ETree_forwSolveProfile ( ETree *etree, double dvec[] ) { int I, J, maxstack, nDJ, nfront, nUJ, stack ; int *bndwghts, *fch, *nodwghts, *sib ; Tree *tree ; /* --------------- check the input --------------- */ if ( etree == NULL || dvec == NULL ) { fprintf(stderr, "\n fatal error in ETree_forwSolveProfile(%p,%p)" "\n bad input\n", etree, dvec) ; exit(-1) ; } tree = ETree_tree(etree) ; nodwghts = ETree_nodwghts(etree) ; bndwghts = ETree_bndwghts(etree) ; nfront = ETree_nfront(etree) ; fch = ETree_fch(etree) ; sib = ETree_sib(etree) ; /* --------------------------------------------- loop over the nodes in a post-order traversal --------------------------------------------- */ maxstack = stack = 0 ; for ( J = Tree_postOTfirst(tree) ; J != -1 ; J = Tree_postOTnext(tree, J) ) { nDJ = nodwghts[J] ; nUJ = bndwghts[J] ; stack += nDJ + nUJ ; dvec[J] = stack ; if ( maxstack < stack ) { maxstack = stack ; } #if MYDEBUG > 0 fprintf(stdout, "\n working on front %d, nD = %d, nU = %d, stack = %d", J, nDJ, nUJ, stack) ; #endif for ( I = fch[J] ; I != -1 ; I = sib[I] ) { stack -= bndwghts[I] ; } stack -= nDJ ; } #if MYDEBUG >= 0 fprintf(stdout, "\n forward solve : final stack = %d, max stack = %d", stack, maxstack) ; #endif return ; }
/*--------------------------------------------------------------------*/ int main ( int argc, char *argv[] ) /* ------------------------------------------------------ (1) read in an ETree object. (2) read in an Graph object. (3) find the optimal domain/schur complement partition for a semi-implicit factorization created -- 96oct03, cca ------------------------------------------------------ */ { char *inETreeFileName, *inGraphFileName, *outIVfileName ; double alpha, nA21, nfent1, nfops1, nL11, nL22, nPhi, nV, t1, t2 ; Graph *graph ; int ii, inside, J, K, msglvl, nfind1, nfront, nJ, nleaves1, nnode1, nvtx, rc, sizeJ, totalgain, vsize, v, w ; int *adjJ, *compids, *nodwghts, *vadj, *vtxToFront, *vwghts ; IV *compidsIV ; IVL *symbfacIVL ; ETree *etree ; FILE *msgFile ; Tree *tree ; if ( argc != 7 ) { fprintf(stdout, "\n\n usage : %s msglvl msgFile inETreeFile inGraphFile alpha" "\n outIVfile " "\n msglvl -- message level" "\n msgFile -- message file" "\n inETreeFile -- input file, must be *.etreef or *.etreeb" "\n inGraphFile -- input file, must be *.graphf or *.graphb" "\n alpha -- weight parameter" "\n alpha = 0 --> minimize storage" "\n alpha = 1 --> minimize solve ops" "\n outIVfile -- output file for oldToNew vector," "\n must be *.ivf or *.ivb" "\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) ; } inETreeFileName = argv[3] ; inGraphFileName = argv[4] ; alpha = atof(argv[5]) ; outIVfileName = argv[6] ; fprintf(msgFile, "\n %s " "\n msglvl -- %d" "\n msgFile -- %s" "\n inETreeFile -- %s" "\n inGraphFile -- %s" "\n alpha -- %f" "\n outIVfile -- %s" "\n", argv[0], msglvl, argv[2], inETreeFileName, inGraphFileName, alpha, outIVfileName) ; fflush(msgFile) ; /* ------------------------ read in the ETree object ------------------------ */ if ( strcmp(inETreeFileName, "none") == 0 ) { fprintf(msgFile, "\n no file to read from") ; spoolesFatal(); } etree = ETree_new() ; MARKTIME(t1) ; rc = ETree_readFromFile(etree, inETreeFileName) ; MARKTIME(t2) ; fprintf(msgFile, "\n CPU %9.5f : read in etree from file %s", t2 - t1, inETreeFileName) ; if ( rc != 1 ) { fprintf(msgFile, "\n return value %d from ETree_readFromFile(%p,%s)", rc, etree, inETreeFileName) ; spoolesFatal(); } ETree_leftJustify(etree) ; fprintf(msgFile, "\n\n after reading ETree object from file %s", inETreeFileName) ; if ( msglvl > 2 ) { ETree_writeForHumanEye(etree, msgFile) ; } else { ETree_writeStats(etree, msgFile) ; } fflush(msgFile) ; nfront = ETree_nfront(etree) ; tree = ETree_tree(etree) ; nodwghts = ETree_nodwghts(etree) ; vtxToFront = ETree_vtxToFront(etree) ; /* ------------------------ read in the Graph object ------------------------ */ if ( strcmp(inGraphFileName, "none") == 0 ) { fprintf(msgFile, "\n no file to read from") ; spoolesFatal(); } graph = Graph_new() ; MARKTIME(t1) ; rc = Graph_readFromFile(graph, inGraphFileName) ; nvtx = graph->nvtx ; vwghts = graph->vwghts ; MARKTIME(t2) ; fprintf(msgFile, "\n CPU %9.5f : read in graph from file %s", t2 - t1, inGraphFileName) ; if ( rc != 1 ) { fprintf(msgFile, "\n return value %d from Graph_readFromFile(%p,%s)", rc, graph, inGraphFileName) ; spoolesFatal(); } fprintf(msgFile, "\n\n after reading Graph object from file %s", inGraphFileName) ; if ( msglvl > 2 ) { Graph_writeForHumanEye(graph, msgFile) ; } else { Graph_writeStats(graph, msgFile) ; } fflush(msgFile) ; /* ---------------------- compute the statistics ---------------------- */ nnode1 = etree->tree->n ; nfind1 = ETree_nFactorIndices(etree) ; nfent1 = ETree_nFactorEntries(etree, SPOOLES_SYMMETRIC) ; nfops1 = ETree_nFactorOps(etree, SPOOLES_REAL, SPOOLES_SYMMETRIC) ; nleaves1 = Tree_nleaves(etree->tree) ; fprintf(stdout, "\n root front %d has %d vertices", etree->tree->root, etree->nodwghtsIV->vec[etree->tree->root]) ; /* --------------------------------- create the symbolic factorization --------------------------------- */ symbfacIVL = SymbFac_initFromGraph(etree, graph) ; if ( msglvl > 2 ) { IVL_writeForHumanEye(symbfacIVL, msgFile) ; } else { IVL_writeStats(symbfacIVL, msgFile) ; } fflush(msgFile) ; /* -------------------------- find the optimal partition -------------------------- */ compidsIV = ETree_optPart(etree, graph, symbfacIVL, alpha, &totalgain, msglvl, msgFile) ; if ( msglvl > 2 ) { IV_writeForHumanEye(compidsIV, msgFile) ; } else { IV_writeStats(compidsIV, msgFile) ; } fflush(msgFile) ; compids = IV_entries(compidsIV) ; /* ------------------------------------------------------ compute the number of vertices in the schur complement ------------------------------------------------------ */ for ( J = 0, nPhi = nV = 0. ; J < nfront ; J++ ) { if ( compids[J] == 0 ) { nPhi += nodwghts[J] ; } nV += nodwghts[J] ; } /* -------------------------------------------- compute the number of entries in L11 and L22 -------------------------------------------- */ nL11 = nL22 = 0 ; for ( J = Tree_postOTfirst(tree) ; J != -1 ; J = Tree_postOTnext(tree, J) ) { nJ = nodwghts[J] ; if ( msglvl > 3 ) { fprintf(msgFile, "\n\n front %d, nJ = %d", J, nJ) ; } IVL_listAndSize(symbfacIVL, J, &sizeJ, &adjJ) ; for ( ii = 0, inside = 0 ; ii < sizeJ ; ii++ ) { w = adjJ[ii] ; K = vtxToFront[w] ; if ( msglvl > 3 ) { fprintf(msgFile, "\n w = %d, K = %d", w, K) ; } if ( K > J && compids[K] == compids[J] ) { inside += (vwghts == NULL) ? 1 : vwghts[w] ; if ( msglvl > 3 ) { fprintf(msgFile, ", inside") ; } } } if ( compids[J] != 0 ) { if ( msglvl > 3 ) { fprintf(msgFile, "\n inside = %d, adding %d to L11", inside, nJ*nJ + 2*nJ*inside) ; } nL11 += (nJ*(nJ+1))/2 + nJ*inside ; } else { if ( msglvl > 3 ) { fprintf(msgFile, "\n inside = %d, adding %d to L22", inside, (nJ*(nJ+1))/2 + nJ*inside) ; } nL22 += (nJ*(nJ+1))/2 + nJ*inside ; } } if ( msglvl > 0 ) { fprintf(msgFile, "\n |L| = %.0f, |L11| = %.0f, |L22| = %.0f", nfent1, nL11, nL22) ; } /* ------------------------------------ compute the number of entries in A21 ------------------------------------ */ nA21 = 0 ; if ( vwghts != NULL ) { for ( v = 0 ; v < nvtx ; v++ ) { J = vtxToFront[v] ; if ( compids[J] != 0 ) { Graph_adjAndSize(graph, v, &vsize, &vadj) ; for ( ii = 0 ; ii < vsize ; ii++ ) { w = vadj[ii] ; K = vtxToFront[w] ; if ( compids[K] == 0 ) { if ( msglvl > 3 ) { fprintf(msgFile, "\n A21 : v = %d, w = %d", v, w) ; } nA21 += vwghts[v] * vwghts[w] ; } } } } } else { for ( v = 0 ; v < nvtx ; v++ ) { J = vtxToFront[v] ; if ( compids[J] != 0 ) { Graph_adjAndSize(graph, v, &vsize, &vadj) ; for ( ii = 0 ; ii < vsize ; ii++ ) { w = vadj[ii] ; K = vtxToFront[w] ; if ( compids[K] == 0 ) { if ( msglvl > 3 ) { fprintf(msgFile, "\n A21 : v = %d, w = %d", v, w) ; } nA21++ ; } } } } } if ( msglvl > 0 ) { fprintf(msgFile, "\n |L| = %.0f, |L11| = %.0f, |L22| = %.0f, |A21| = %.0f", nfent1, nL11, nL22, nA21) ; fprintf(msgFile, "\n storage: explicit = %.0f, semi-implicit = %.0f, ratio = %.3f" "\n opcount: explicit = %.0f, semi-implicit = %.0f, ratio = %.3f", nfent1, nL11 + nA21 + nL22, nfent1/(nL11 + nA21 + nL22), 2*nfent1, 4*nL11 + 2*nA21 + 2*nL22, 2*nfent1/(4*nL11 + 2*nA21 + 2*nL22)) ; fprintf(msgFile, "\n ratios %8.3f %8.3f %8.3f", nPhi/nV, nfent1/(nL11 + nA21 + nL22), 2*nfent1/(4*nL11 + 2*nA21 + 2*nL22)) ; } /* ---------------- free the objects ---------------- */ ETree_free(etree) ; Graph_free(graph) ; IVL_free(symbfacIVL) ; fprintf(msgFile, "\n") ; fclose(msgFile) ; return(1) ; }
/* --------------------------------------------------------------- purpose -- fill dvec[J] with the active storage to eliminate J using the right-looking general sparse method symflag -- symmetry flag, one of SPOOLES_SYMMETRIC, SPOOLES_HERMITIAN or SPOOLES_NONSYMMETRIC created -- 98dec19, cca --------------------------------------------------------------- */ void ETree_FSstorageProfile ( ETree *etree, int symflag, IVL *symbfacIVL, double dvec[] ) { char *incore ; int ii, J, K, nDJ, nfront, nUJ, sizeJ, storage ; int *bndwghts, *indJ, *mark, *nodwghts, *stor, *vtxToFront ; Tree *tree ; /* --------------- check the input --------------- */ if ( etree == NULL || symbfacIVL == NULL || dvec == NULL ) { fprintf(stderr, "\n fatal error in ETree_FSstorageProfile(%p,%p,%p)" "\n bad input\n", etree, symbfacIVL, dvec) ; exit(-1) ; } tree = ETree_tree(etree) ; nodwghts = ETree_nodwghts(etree) ; bndwghts = ETree_bndwghts(etree) ; vtxToFront = ETree_vtxToFront(etree) ; nfront = ETree_nfront(etree) ; incore = CVinit(nfront, 'F') ; stor = IVinit(nfront, 0) ; mark = IVinit(nfront, -1) ; /* -------------------------------------------- compute the storage for each front's chevron -------------------------------------------- */ if ( symflag == SPOOLES_SYMMETRIC || symflag == SPOOLES_HERMITIAN ) { for ( J = 0 ; J < nfront ; J++ ) { nDJ = nodwghts[J] ; nUJ = bndwghts[J] ; stor[J] = (nDJ*(nDJ+1))/2 + nDJ*nUJ ; } } else { for ( J = 0 ; J < nfront ; J++ ) { nDJ = nodwghts[J] ; nUJ = bndwghts[J] ; stor[J] = nDJ*nDJ + 2*nDJ*nUJ ; } } /* --------------------------------------------- loop over the nodes in a post-order traversal --------------------------------------------- */ storage = 0 ; for ( J = Tree_postOTfirst(tree) ; J != -1 ; J = Tree_postOTnext(tree, J) ) { if ( incore[J] == 'F' ) { storage += stor[J] ; incore[J] = 'T' ; } IVL_listAndSize(symbfacIVL, J, &sizeJ, &indJ) ; mark[J] = J ; for ( ii = 0 ; ii < sizeJ ; ii++ ) { K = vtxToFront[indJ[ii]] ; if ( mark[K] != J ) { mark[K] = J ; if ( incore[K] == 'F' ) { storage += stor[K] ; incore[K] = 'T' ; } } } dvec[J] = storage ; storage -= stor[J] ; } IVfree(mark) ; IVfree(stor) ; CVfree(incore) ; return ; }
/* --------------------------------------------------------------- purpose -- fill dvec[J] with the active storage to eliminate J using the multifrontal method symflag -- symmetry flag, one of SPOOLES_SYMMETRIC, SPOOLES_HERMITIAN or SPOOLES_NONSYMMETRIC created -- 97may21, cca --------------------------------------------------------------- */ void ETree_MFstackProfile ( ETree *etree, int symflag, double dvec[] ) { int I, J, nDJ, nUI, nUJ, stack ; int *bndwghts, *fch, *nodwghts, *sib ; Tree *tree ; /* --------------- check the input --------------- */ if ( etree == NULL || dvec == NULL ) { fprintf(stderr, "\n fatal error in ETree_MFstackProfile(%p,%p)" "\n bad input\n", etree, dvec) ; exit(-1) ; } tree = ETree_tree(etree) ; nodwghts = ETree_nodwghts(etree) ; bndwghts = ETree_bndwghts(etree) ; fch = ETree_fch(etree) ; sib = ETree_sib(etree) ; /* --------------------------------------------- loop over the nodes in a post-order traversal --------------------------------------------- */ stack = 0 ; for ( J = Tree_postOTfirst(tree) ; J != -1 ; J = Tree_postOTnext(tree, J) ) { nDJ = nodwghts[J] ; nUJ = bndwghts[J] ; #if MYDEBUG > 0 fprintf(stdout, "\n working on front %d, nD = %d, nU = %d, stack = %d", J, nDJ, nUJ, stack) ; #endif if ( (I = fch[J]) != -1 ) { /* -------------------------------- remove last child from the stack -------------------------------- */ while ( sib[I] != -1 ) { I = sib[I] ; } nUI = bndwghts[I] ; if ( symflag == SPOOLES_SYMMETRIC || symflag == SPOOLES_HERMITIAN ) { stack -= (nUI*(nUI+1))/2 ; } else if ( symflag == SPOOLES_NONSYMMETRIC ) { stack -= nUI*nUI ; } #if MYDEBUG > 0 fprintf(stdout, "\n removing last child %d, stack = %d", I, stack) ; #endif } /* --------------------------------- add frontal matrix for J to stack and store as the storage for J --------------------------------- */ dvec[J] = stack + (nDJ + nUJ)*(nDJ + nUJ) ; if ( symflag == SPOOLES_SYMMETRIC || symflag == SPOOLES_HERMITIAN ) { dvec[J] = stack + ((nDJ + nUJ)*(nDJ + nUJ+1))/2 ; } else if ( symflag == SPOOLES_NONSYMMETRIC ) { dvec[J] = stack + (nDJ + nUJ)*(nDJ + nUJ) ; } #if MYDEBUG > 0 fprintf(stdout, "\n working storage = %.0f", dvec[J]) ; #endif if ( (I = fch[J]) != -1 ) { /* ------------------------------------ remove other children from the stack ------------------------------------ */ while ( sib[I] != -1 ) { nUI = bndwghts[I] ; if ( symflag == SPOOLES_SYMMETRIC || symflag == SPOOLES_HERMITIAN ) { stack -= (nUI*(nUI+1))/2 ; } else if ( symflag == SPOOLES_NONSYMMETRIC ) { stack -= nUI*nUI ; } #if MYDEBUG > 0 fprintf(stdout, "\n removing child %d, stack = %d", I, stack) ; #endif I = sib[I] ; } } /* -------------------------------- add update matrix for J to stack -------------------------------- */ if ( symflag == SPOOLES_SYMMETRIC || symflag == SPOOLES_HERMITIAN ) { stack += (nUJ*(nUJ+1))/2 ; } else if ( symflag == SPOOLES_NONSYMMETRIC ) { stack += nUJ*nUJ ; } #if MYDEBUG > 0 fprintf(stdout, "\n adding update matrix, stack = %d", stack) ; #endif } #if MYDEBUG >= 0 fprintf(stdout, "\n MF: final stack = %d", stack) ; #endif return ; }
/* --------------------------------------------------------------- purpose -- fill dvec[J] with the active storage to eliminate J using the left-looking general sparse method symflag -- symmetry flag, one of SPOOLES_SYMMETRIC, SPOOLES_HERMITIAN or SPOOLES_NONSYMMETRIC created -- 97may21, cca --------------------------------------------------------------- */ void ETree_GSstorageProfile ( ETree *etree, int symflag, IVL *symbfacIVL, int *vwghts, double dvec[] ) { int count, ii, I, J, K, nDJ, nfront, nUJ, sizeI, sizeJ, storage, v ; int *bndwghts, *head, *indI, *indJ, *link, *nodwghts, *offsets, *vtxToFront ; Tree *tree ; /* --------------- check the input --------------- */ if ( etree == NULL || symbfacIVL == NULL || dvec == NULL ) { fprintf(stderr, "\n fatal error in ETree_GSstorageProfile(%p,%p,%p,%p)" "\n bad input\n", etree, symbfacIVL, vwghts, dvec) ; exit(-1) ; } tree = ETree_tree(etree) ; nodwghts = ETree_nodwghts(etree) ; bndwghts = ETree_bndwghts(etree) ; vtxToFront = ETree_vtxToFront(etree) ; nfront = ETree_nfront(etree) ; head = IVinit(nfront, -1) ; link = IVinit(nfront, -1) ; offsets = IVinit(nfront, 0) ; /* --------------------------------------------- loop over the nodes in a post-order traversal --------------------------------------------- */ storage = 0 ; for ( J = Tree_postOTfirst(tree) ; J != -1 ; J = Tree_postOTnext(tree, J) ) { nDJ = nodwghts[J] ; nUJ = bndwghts[J] ; if ( symflag == SPOOLES_SYMMETRIC || symflag == SPOOLES_HERMITIAN ) { storage += (nDJ*(nDJ + 1))/2 + nDJ*nUJ ; } else if ( symflag == SPOOLES_NONSYMMETRIC ) { storage += nDJ*nDJ + 2*nDJ*nUJ ; } dvec[J] = storage ; #if MYDEBUG > 0 fprintf(stdout, "\n working on front %d, nD = %d, nU = %d, storage = %d", J, nDJ, nUJ, storage) ; #endif /* ----------------------------- loop over the updating fronts ----------------------------- */ while ( (I = head[J]) != -1 ) { head[J] = link[I] ; IVL_listAndSize(symbfacIVL, I, &sizeI, &indI) ; #if MYDEBUG > 0 fprintf(stdout, "\n updating front %d, offset = %d, sizeI = %d", I, offsets[I], sizeI) ; IVfprintf(stdout, sizeI, indI) ; #endif for ( ii = offsets[I], count = 0, K = -1 ; ii < sizeI ; ii++ ) { v = indI[ii] ; #if MYDEBUG > 0 fprintf(stdout, "\n ii = %d, v = %d, K = %d", ii, v, vtxToFront[v]) ; fflush(stdout) ; #endif K = vtxToFront[v] ; if ( K < 0 || K >= nfront ) { fprintf(stderr, "\n\n fatal error" "\n ii = %d, v = %d, K = %d", ii, v, K) ; exit(-1) ; } if ( (K = vtxToFront[v]) != J ) { #if MYDEBUG > 0 fprintf(stdout, "\n linking to next ancestor %d", K) ; #endif link[I] = head[K] ; head[K] = I ; offsets[I] = ii ; break ; } count += (vwghts == NULL) ? 1 : vwghts[v] ; #if MYDEBUG > 0 fprintf(stdout, "\n count = %d", count) ; fflush(stdout) ; #endif } if ( symflag == SPOOLES_SYMMETRIC || symflag == SPOOLES_HERMITIAN ) { storage -= count*nodwghts[I] ; } else if ( symflag == SPOOLES_NONSYMMETRIC ) { storage -= 2*count*nodwghts[I] ; } } if ( symflag == SPOOLES_SYMMETRIC || symflag == SPOOLES_HERMITIAN ) { storage -= (nDJ*(nDJ+1))/2 ; } else if ( symflag == SPOOLES_NONSYMMETRIC ) { storage -= nDJ*nDJ ; } if ( nUJ > 0 ) { IVL_listAndSize(symbfacIVL, J, &sizeJ, &indJ) ; for ( ii = 0 ; ii < sizeJ ; ii++ ) { v = indJ[ii] ; if ( (K = vtxToFront[v]) != J ) { break ; } } offsets[J] = ii ; #if MYDEBUG > 0 fprintf(stdout, "\n linking to next ancestor %d", K) ; #endif IVL_listAndSize(symbfacIVL, J, &sizeJ, &indJ) ; link[J] = head[K] ; head[K] = J ; } #if MYDEBUG > 0 fprintf(stdout, "\n at end of step %d, storage = %d", J, storage) ; #endif } #if MYDEBUG >= 0 fprintf(stdout, "\n GS: final storage = %d", storage) ; #endif IVfree(head) ; IVfree(link) ; IVfree(offsets) ; return ; }
/*--------------------------------------------------------------------*/ int main ( int argc, char *argv[] ) /* --------------------------------------------------------------- read in a ETree object, create an IV object with the same size, mark the vertices in the top level separator(s), write the IV object to a file created -- 96may02, cca --------------------------------------------------------------- */ { char *inETreeFileName, *outIVfileName ; double t1, t2 ; int msglvl, rc, J, K, ncomp, nfront, nvtx, v ; int *bndwghts, *compids, *fch, *map, *nodwghts, *par, *sib, *vtxToFront ; IV *compidsIV, *mapIV ; ETree *etree ; FILE *msgFile ; Tree *tree ; if ( argc != 5 ) { fprintf(stdout, "\n\n usage : %s msglvl msgFile inETreeFile outIVfile" "\n msglvl -- message level" "\n msgFile -- message file" "\n inETreeFile -- input file, must be *.etreef or *.etreeb" "\n outIVfile -- output file, must be *.ivf or *.ivb" "\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) ; } inETreeFileName = argv[3] ; outIVfileName = argv[4] ; fprintf(msgFile, "\n %s " "\n msglvl -- %d" "\n msgFile -- %s" "\n inETreeFile -- %s" "\n outIVfile -- %s" "\n", argv[0], msglvl, argv[2], inETreeFileName, outIVfileName) ; fflush(msgFile) ; /* ------------------------ read in the ETree object ------------------------ */ if ( strcmp(inETreeFileName, "none") == 0 ) { fprintf(msgFile, "\n no file to read from") ; exit(0) ; } etree = ETree_new() ; MARKTIME(t1) ; rc = ETree_readFromFile(etree, inETreeFileName) ; MARKTIME(t2) ; fprintf(msgFile, "\n CPU %9.5f : read in etree from file %s", t2 - t1, inETreeFileName) ; if ( rc != 1 ) { fprintf(msgFile, "\n return value %d from ETree_readFromFile(%p,%s)", rc, etree, inETreeFileName) ; exit(-1) ; } fprintf(msgFile, "\n\n after reading ETree object from file %s", inETreeFileName) ; if ( msglvl > 2 ) { ETree_writeForHumanEye(etree, msgFile) ; } else { ETree_writeStats(etree, msgFile) ; } fflush(msgFile) ; nfront = ETree_nfront(etree) ; nvtx = ETree_nvtx(etree) ; bndwghts = ETree_bndwghts(etree) ; vtxToFront = ETree_vtxToFront(etree) ; nodwghts = ETree_nodwghts(etree) ; par = ETree_par(etree) ; fch = ETree_fch(etree) ; sib = ETree_sib(etree) ; tree = ETree_tree(etree) ; /* ----------------------------------------- create the map from fronts to components, top level separator(s) are component zero ----------------------------------------- */ mapIV = IV_new() ; IV_init(mapIV, nfront, NULL) ; map = IV_entries(mapIV) ; ncomp = 0 ; for ( J = Tree_preOTfirst(tree) ; J != -1 ; J = Tree_preOTnext(tree, J) ) { if ( (K = par[J]) == -1 ) { map[J] = 0 ; } else if ( map[K] != 0 ) { map[J] = map[K] ; } else if ( J == fch[K] && sib[J] == -1 && bndwghts[J] == nodwghts[K] + bndwghts[K] ) { map[J] = 0 ; } else { map[J] = ++ncomp ; } } fprintf(msgFile, "\n\n mapIV object") ; if ( msglvl > 2 ) { IV_writeForHumanEye(mapIV, msgFile) ; } else { IV_writeStats(mapIV, msgFile) ; } /* ---------------------------------------- fill the map from vertices to components ---------------------------------------- */ compidsIV = IV_new() ; IV_init(compidsIV, nvtx, NULL) ; compids = IV_entries(compidsIV) ; for ( v = 0 ; v < nvtx ; v++ ) { compids[v] = map[vtxToFront[v]] ; } fprintf(msgFile, "\n\n compidsIV object") ; if ( msglvl > 2 ) { IV_writeForHumanEye(compidsIV, msgFile) ; } else { IV_writeStats(compidsIV, msgFile) ; } fflush(msgFile) ; /* ----------------------- write out the IV object ----------------------- */ if ( strcmp(outIVfileName, "none") != 0 ) { MARKTIME(t1) ; rc = IV_writeToFile(compidsIV, outIVfileName) ; MARKTIME(t2) ; fprintf(msgFile, "\n CPU %9.5f : write etree to file %s", t2 - t1, outIVfileName) ; } if ( rc != 1 ) { fprintf(msgFile, "\n return value %d from IV_writeToFile(%p,%s)", rc, compidsIV, outIVfileName) ; } /* ---------------- free the objects ---------------- */ ETree_free(etree) ; IV_free(mapIV) ; IV_free(compidsIV) ; fprintf(msgFile, "\n") ; fclose(msgFile) ; return(1) ; }
/*--------------------------------------------------------------------*/ int main ( int argc, char *argv[] ) /* ---------------------------------------- get statistics for a semi-implicit solve created -- 97dec11, cca ---------------------------------------- */ { char *inGraphFileName, *inETreeFileName, *inMapFileName ; double nA21, nL, nL11, nL22, nPhi, nV, t1, t2 ; ETree *etree ; int ii, inside, J, K, msglvl, nfront, nJ, nvtx, rc, sizeJ, v, vsize, w ; int *adjJ, *frontmap, *map, *nodwghts, *vadj, *vtxToFront, *vwghts ; IV *mapIV ; IVL *symbfacIVL ; Graph *graph ; FILE *msgFile ; Tree *tree ; if ( argc != 6 ) { fprintf(stdout, "\n\n usage : %s msglvl msgFile GraphFile ETreeFile mapFile " "\n msglvl -- message level" "\n msgFile -- message file" "\n GraphFile -- input graph file, must be *.graphf or *.graphb" "\n ETreeFile -- input ETree file, must be *.etreef or *.etreeb" "\n mapFile -- input map IV file, must be *.ivf or *.ivb" "\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) ; } inGraphFileName = argv[3] ; inETreeFileName = argv[4] ; inMapFileName = argv[5] ; fprintf(msgFile, "\n %s " "\n msglvl -- %d" "\n msgFile -- %s" "\n GraphFile -- %s" "\n ETreeFile -- %s" "\n mapFile -- %s" "\n", argv[0], msglvl, argv[2], inGraphFileName, inETreeFileName, inMapFileName) ; fflush(msgFile) ; /* ------------------------ read in the Graph object ------------------------ */ graph = Graph_new() ; if ( strcmp(inGraphFileName, "none") == 0 ) { fprintf(msgFile, "\n no file to read from") ; exit(0) ; } MARKTIME(t1) ; rc = Graph_readFromFile(graph, inGraphFileName) ; MARKTIME(t2) ; fprintf(msgFile, "\n CPU %9.5f : read in graph from file %s", t2 - t1, inGraphFileName) ; nvtx = graph->nvtx ; vwghts = graph->vwghts ; if ( rc != 1 ) { fprintf(msgFile, "\n return value %d from Graph_readFromFile(%p,%s)", rc, graph, inGraphFileName) ; exit(-1) ; } if ( msglvl > 2 ) { fprintf(msgFile, "\n\n after reading Graph object from file %s", inGraphFileName) ; Graph_writeForHumanEye(graph, msgFile) ; fflush(msgFile) ; } /* ------------------------ read in the ETree object ------------------------ */ etree = ETree_new() ; if ( strcmp(inETreeFileName, "none") == 0 ) { fprintf(msgFile, "\n no file to read from") ; exit(0) ; } MARKTIME(t1) ; rc = ETree_readFromFile(etree, inETreeFileName) ; MARKTIME(t2) ; fprintf(msgFile, "\n CPU %9.5f : read in etree from file %s", t2 - t1, inETreeFileName) ; nfront = ETree_nfront(etree) ; tree = ETree_tree(etree) ; vtxToFront = ETree_vtxToFront(etree) ; nodwghts = ETree_nodwghts(etree) ; nL = ETree_nFactorEntries(etree, 2) ; if ( rc != 1 ) { fprintf(msgFile, "\n return value %d from ETree_readFromFile(%p,%s)", rc, etree, inETreeFileName) ; exit(-1) ; } if ( msglvl > 2 ) { fprintf(msgFile, "\n\n after reading ETree object from file %s", inETreeFileName) ; ETree_writeForHumanEye(etree, msgFile) ; fflush(msgFile) ; } /* ------------------------- read in the map IV object ------------------------- */ mapIV = IV_new() ; if ( strcmp(inMapFileName, "none") == 0 ) { fprintf(msgFile, "\n no file to read from") ; exit(0) ; } MARKTIME(t1) ; rc = IV_readFromFile(mapIV, inMapFileName) ; MARKTIME(t2) ; fprintf(msgFile, "\n CPU %9.5f : read in mapIV from file %s", t2 - t1, inMapFileName) ; map = IV_entries(mapIV) ; if ( rc != 1 ) { fprintf(msgFile, "\n return value %d from IV_readFromFile(%p,%s)", rc, mapIV, inMapFileName) ; exit(-1) ; } if ( msglvl > 2 ) { fprintf(msgFile, "\n\n after reading IV object from file %s", inMapFileName) ; IV_writeForHumanEye(mapIV, msgFile) ; fflush(msgFile) ; } nV = nPhi = 0 ; if ( vwghts == NULL ) { for ( v = 0 ; v < nvtx ; v++ ) { nV++ ; if ( map[v] == 0 ) { nPhi++ ; } } } else { for ( v = 0 ; v < nvtx ; v++ ) { nV += vwghts[v] ; if ( map[v] == 0 ) { nPhi += vwghts[v] ; } } } fprintf(msgFile, "\n nPhi = %.0f, nV = %.0f", nPhi, nV) ; /* ------------------------- get the frontmap[] vector ------------------------- */ frontmap = IVinit(nfront, -1) ; for ( v = 0 ; v < nvtx ; v++ ) { J = vtxToFront[v] ; if ( frontmap[J] == -1 ) { frontmap[J] = map[v] ; } else if ( frontmap[J] != map[v] ) { fprintf(msgFile, "\n\n error, frontmap[%d] = %d, map[%d] = %d", J, frontmap[J], v, map[v]) ; } } /* ---------------------------------- compute the symbolic factorization ---------------------------------- */ symbfacIVL = SymbFac_initFromGraph(etree, graph) ; if ( msglvl > 2 ) { fprintf(msgFile, "\n\n symbolic factorization") ; IVL_writeForHumanEye(symbfacIVL, msgFile) ; fflush(msgFile) ; } /* -------------------------------------------- compute the number of entries in L11 and L22 -------------------------------------------- */ nL11 = nL22 = 0 ; for ( J = Tree_postOTfirst(tree) ; J != -1 ; J = Tree_postOTnext(tree, J) ) { nJ = nodwghts[J] ; if ( msglvl > 3 ) { fprintf(msgFile, "\n\n front %d, nJ = %d", J, nJ) ; } IVL_listAndSize(symbfacIVL, J, &sizeJ, &adjJ) ; for ( ii = 0, inside = 0 ; ii < sizeJ ; ii++ ) { w = adjJ[ii] ; K = vtxToFront[w] ; if ( msglvl > 3 ) { fprintf(msgFile, "\n w = %d, K = %d", w, K) ; } if ( K > J && frontmap[K] == frontmap[J] ) { inside += (vwghts == NULL) ? 1 : vwghts[w] ; if ( msglvl > 3 ) { fprintf(msgFile, ", inside") ; } } } if ( frontmap[J] != 0 ) { if ( msglvl > 3 ) { fprintf(msgFile, "\n inside = %d, adding %d to L11", inside, nJ*nJ + 2*nJ*inside) ; } nL11 += nJ*nJ + 2*nJ*inside ; } else { if ( msglvl > 3 ) { fprintf(msgFile, "\n inside = %d, adding %d to L22", inside, nJ*nJ + 2*nJ*inside) ; } nL22 += nJ*nJ + 2*nJ*inside ; } } if ( msglvl > 0 ) { fprintf(msgFile, "\n |L| = %.0f, |L11| = %.0f, |L22| = %.0f", nL, nL11, nL22) ; } /* ------------------------------------ compute the number of entries in A21 ------------------------------------ */ nA21 = 0 ; if ( vwghts != NULL ) { for ( v = 0 ; v < nvtx ; v++ ) { if ( map[v] == 0 ) { Graph_adjAndSize(graph, v, &vsize, &vadj) ; for ( ii = 0 ; ii < vsize ; ii++ ) { w = vadj[ii] ; if ( map[v] != map[w] ) { if ( msglvl > 3 ) { fprintf(msgFile, "\n A21 : v = %d, w = %d", v, w) ; } nA21 += vwghts[v] * vwghts[w] ; } } } } } else { for ( v = 0 ; v < nvtx ; v++ ) { if ( map[v] == 0 ) { Graph_adjAndSize(graph, v, &vsize, &vadj) ; for ( ii = 0 ; ii < vsize ; ii++ ) { w = vadj[ii] ; if ( map[v] != map[w] ) { if ( msglvl > 3 ) { fprintf(msgFile, "\n A21 : v = %d, w = %d", v, w) ; } nA21++ ; } } } } } if ( msglvl > 0 ) { fprintf(msgFile, "\n |L| = %.0f, |L11| = %.0f, |L22| = %.0f, |A21| = %.0f", nL, nL11, nL22, nA21) ; fprintf(msgFile, "\n storage: explicit = %.0f, semi-implicit = %.0f, ratio = %.3f" "\n opcount: explicit = %.0f, semi-implicit = %.0f, ratio = %.3f", nL, nL11 + nA21 + nL22, nL/(nL11 + nA21 + nL22), 2*nL, 4*nL11 + 2*nA21 + 2*nL22, 2*nL/(4*nL11 + 2*nA21 + 2*nL22)) ; fprintf(msgFile, "\n ratios %8.3f %8.3f %8.3f", nPhi/nV, nL/(nL11 + nA21 + nL22), 2*nL/(4*nL11 + 2*nA21 + 2*nL22)) ; } /* ------------------------ free the working storage ------------------------ */ Graph_free(graph) ; ETree_free(etree) ; IV_free(mapIV) ; IVL_free(symbfacIVL) ; IVfree(frontmap) ; fprintf(msgFile, "\n") ; fclose(msgFile) ; return(1) ; }