/* --------------------------------------- purpose -- initialize the EGraph object created -- 96oct24, cca --------------------------------------- */ void EGraph_init ( EGraph *egraph, int type, int nelem, int nvtx, int IVL_type ) { /* --------------- check the input --------------- */ if ( egraph == NULL || type < 0 || type > 1 || nelem <= 0 || nvtx <= 0 ) { fprintf(stderr, "\n fatal error in EGraph_init(%p,%d,%d,%d,%d)" "\n bad input\n", egraph, type, nelem, nvtx, IVL_type) ; exit(-1) ; } /* ---------------------------- clear the data in the object ---------------------------- */ EGraph_clearData(egraph) ; /* --------------------- initialize the object --------------------- */ egraph->type = type ; egraph->nelem = nelem ; egraph->nvtx = nvtx ; egraph->adjIVL = IVL_new() ; IVL_init1(egraph->adjIVL, IVL_type, nelem) ; if ( type == 1 ) { egraph->vwghts = IVinit(nvtx, 0) ; } return ; }
/*--------------------------------------------------------------------*/ int main ( int argc, char *argv[] ) /* -------------------------------------------------------------------- this program tests the Graph_MPI_Bcast() method (1) process root generates a random Graph object and computes its checksum (2) process root broadcasts the Graph object to the other processors (3) each process computes the checksum of its Graph object (4) the checksums are compared on root created -- 98sep10, cca -------------------------------------------------------------------- */ { char *buffer ; double chksum, t1, t2 ; double *sums ; Drand drand ; int iproc, length, loc, msglvl, myid, nitem, nproc, nvtx, root, seed, size, type, v ; int *list ; FILE *msgFile ; Graph *graph ; /* --------------------------------------------------------------- find out the identity of this process and the number of process --------------------------------------------------------------- */ MPI_Init(&argc, &argv) ; MPI_Comm_rank(MPI_COMM_WORLD, &myid) ; MPI_Comm_size(MPI_COMM_WORLD, &nproc) ; fprintf(stdout, "\n process %d of %d, argc = %d", myid, nproc, argc) ; fflush(stdout) ; if ( argc != 8 ) { fprintf(stdout, "\n\n usage : %s msglvl msgFile type nvtx nitem root seed " "\n msglvl -- message level" "\n msgFile -- message file" "\n type -- type of graph" "\n nvtx -- # of vertices" "\n nitem -- # of items used to generate graph" "\n root -- root processor for broadcast" "\n seed -- random number seed" "\n", argv[0]) ; return(0) ; } msglvl = atoi(argv[1]) ; if ( strcmp(argv[2], "stdout") == 0 ) { msgFile = stdout ; } else { length = strlen(argv[2]) + 1 + 4 ; buffer = CVinit(length, '\0') ; sprintf(buffer, "%s.%d", argv[2], myid) ; if ( (msgFile = fopen(buffer, "w")) == NULL ) { fprintf(stderr, "\n fatal error in %s" "\n unable to open file %s\n", argv[0], argv[2]) ; return(-1) ; } CVfree(buffer) ; } type = atoi(argv[3]) ; nvtx = atoi(argv[4]) ; nitem = atoi(argv[5]) ; root = atoi(argv[6]) ; seed = atoi(argv[7]) ; fprintf(msgFile, "\n %s " "\n msglvl -- %d" "\n msgFile -- %s" "\n type -- %d" "\n nvtx -- %d" "\n nitem -- %d" "\n root -- %d" "\n seed -- %d" "\n", argv[0], msglvl, argv[2], type, nvtx, nitem, root, seed) ; fflush(msgFile) ; /* ----------------------- set up the Graph object ----------------------- */ MARKTIME(t1) ; graph = Graph_new() ; if ( myid == root ) { InpMtx *inpmtx ; int nedges, totewght, totvwght, v ; int *adj, *vwghts ; IVL *adjIVL, *ewghtIVL ; /* ----------------------- generate a random graph ----------------------- */ inpmtx = InpMtx_new() ; InpMtx_init(inpmtx, INPMTX_BY_ROWS, INPMTX_INDICES_ONLY, nitem, 0) ; Drand_setDefaultFields(&drand) ; Drand_setSeed(&drand, seed) ; Drand_setUniform(&drand, 0, nvtx) ; Drand_fillIvector(&drand, nitem, InpMtx_ivec1(inpmtx)) ; Drand_fillIvector(&drand, nitem, InpMtx_ivec2(inpmtx)) ; InpMtx_setNent(inpmtx, nitem) ; InpMtx_changeStorageMode(inpmtx, INPMTX_BY_VECTORS) ; if ( msglvl > 2 ) { fprintf(msgFile, "\n\n inpmtx mtx filled with raw entries") ; InpMtx_writeForHumanEye(inpmtx, msgFile) ; fflush(msgFile) ; } adjIVL = InpMtx_fullAdjacency(inpmtx) ; if ( msglvl > 2 ) { fprintf(msgFile, "\n\n full adjacency structure") ; IVL_writeForHumanEye(adjIVL, msgFile) ; fflush(msgFile) ; } nedges = adjIVL->tsize ; if ( type == 1 || type == 3 ) { Drand_setUniform(&drand, 1, 10) ; vwghts = IVinit(nvtx, 0) ; Drand_fillIvector(&drand, nvtx, vwghts) ; totvwght = IVsum(nvtx, vwghts) ; if ( msglvl > 2 ) { fprintf(msgFile, "\n\n vertex weights") ; IVfprintf(msgFile, nvtx, vwghts) ; fflush(msgFile) ; } } else { vwghts = NULL ; totvwght = nvtx ; } if ( msglvl > 2 ) { fprintf(msgFile, "\n\n totvwght %d", totvwght) ; fflush(msgFile) ; } if ( type == 2 || type == 3 ) { ewghtIVL = IVL_new() ; IVL_init1(ewghtIVL, IVL_CHUNKED, nvtx) ; Drand_setUniform(&drand, 1, 100) ; totewght = 0 ; for ( v = 0 ; v < nvtx ; v++ ) { IVL_listAndSize(adjIVL, v, &size, &adj) ; IVL_setList(ewghtIVL, v, size, NULL) ; IVL_listAndSize(ewghtIVL, v, &size, &adj) ; Drand_fillIvector(&drand, size, adj) ; totewght += IVsum(size, adj) ; } if ( msglvl > 2 ) { fprintf(msgFile, "\n\n ewghtIVL") ; IVL_writeForHumanEye(ewghtIVL, msgFile) ; fflush(msgFile) ; } } else { ewghtIVL = NULL ; totewght = nedges ; } if ( msglvl > 2 ) { fprintf(msgFile, "\n\n totewght %d", totewght) ; fflush(msgFile) ; } Graph_init2(graph, type, nvtx, 0, nedges, totvwght, totewght, adjIVL, vwghts, ewghtIVL) ; InpMtx_free(inpmtx) ; } MARKTIME(t2) ; fprintf(msgFile, "\n CPU %8.3f : initialize the Graph object", t2 - t1) ; fflush(msgFile) ; if ( msglvl > 2 ) { Graph_writeForHumanEye(graph, msgFile) ; } else { Graph_writeStats(graph, msgFile) ; } fflush(msgFile) ; if ( myid == root ) { /* ---------------------------------------- compute the checksum of the Graph object ---------------------------------------- */ chksum = graph->type + graph->nvtx + graph->nvbnd + graph->nedges + graph->totvwght + graph->totewght ; for ( v = 0 ; v < nvtx ; v++ ) { IVL_listAndSize(graph->adjIVL, v, &size, &list) ; chksum += 1 + v + size + IVsum(size, list) ; } if ( graph->vwghts != NULL ) { chksum += IVsum(nvtx, graph->vwghts) ; } if ( graph->ewghtIVL != NULL ) { for ( v = 0 ; v < nvtx ; v++ ) { IVL_listAndSize(graph->ewghtIVL, v, &size, &list) ; chksum += 1 + v + size + IVsum(size, list) ; } } fprintf(msgFile, "\n\n local chksum = %12.4e", chksum) ; fflush(msgFile) ; } /* -------------------------- broadcast the Graph object -------------------------- */ MARKTIME(t1) ; graph = Graph_MPI_Bcast(graph, root, msglvl, msgFile, MPI_COMM_WORLD) ; MARKTIME(t2) ; fprintf(msgFile, "\n CPU %8.3f : broadcast the Graph object", t2 - t1) ; if ( msglvl > 2 ) { Graph_writeForHumanEye(graph, msgFile) ; } else { Graph_writeStats(graph, msgFile) ; } /* ---------------------------------------- compute the checksum of the Graph object ---------------------------------------- */ chksum = graph->type + graph->nvtx + graph->nvbnd + graph->nedges + graph->totvwght + graph->totewght ; for ( v = 0 ; v < nvtx ; v++ ) { IVL_listAndSize(graph->adjIVL, v, &size, &list) ; chksum += 1 + v + size + IVsum(size, list) ; } if ( graph->vwghts != NULL ) { chksum += IVsum(nvtx, graph->vwghts) ; } if ( graph->ewghtIVL != NULL ) { for ( v = 0 ; v < nvtx ; v++ ) { IVL_listAndSize(graph->ewghtIVL, v, &size, &list) ; chksum += 1 + v + size + IVsum(size, list) ; } } fprintf(msgFile, "\n\n local chksum = %12.4e", chksum) ; fflush(msgFile) ; /* --------------------------------------- gather the checksums from the processes --------------------------------------- */ sums = DVinit(nproc, 0.0) ; MPI_Gather((void *) &chksum, 1, MPI_DOUBLE, (void *) sums, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD) ; if ( myid == 0 ) { fprintf(msgFile, "\n\n sums") ; DVfprintf(msgFile, nproc, sums) ; for ( iproc = 0 ; iproc < nproc ; iproc++ ) { sums[iproc] -= chksum ; } fprintf(msgFile, "\n\n errors") ; DVfprintf(msgFile, nproc, sums) ; fprintf(msgFile, "\n\n maxerror = %12.4e", DVmax(nproc, sums, &loc)); } /* ---------------- free the objects ---------------- */ DVfree(sums) ; Graph_free(graph) ; /* ------------------------ exit the MPI environment ------------------------ */ MPI_Finalize() ; fprintf(msgFile, "\n") ; fclose(msgFile) ; return(0) ; }
/* ---------------------------------------------------------------- purpose -- to create and return an IVL object that holds the submatrix nonzero pattern for the lower triangular factor. NOTE: this method supercedes calling IVL_mapEntries() on the row adjacency structure. that gave problems when pivoting forced some fronts to have no eliminated columns. in some cases, solve aggregates were expected when none were forthcoming. created -- 98aug20, cca ---------------------------------------------------------------- */ IVL * FrontMtx_makeLowerBlockIVL ( FrontMtx *frontmtx, IV *rowmapIV ) { int count, ii, J, K, nrow, nfront, nJ ; int *rowmap, *rowind, *list, *mark ; IVL *lowerblockIVL ; /* --------------- check the input --------------- */ if ( frontmtx == NULL || rowmapIV == NULL ) { fprintf(stderr, "\n fatal error in FrontMtx_makeLowerBlockIVL()" "\n bad input\n") ; exit(-1) ; } nfront = FrontMtx_nfront(frontmtx) ; rowmap = IV_entries(rowmapIV) ; /* ----------------------------- set up the working storage and initialize the IVL object ----------------------------- */ mark = IVinit(nfront, -1) ; list = IVinit(nfront, -1) ; lowerblockIVL = IVL_new() ; IVL_init1(lowerblockIVL, IVL_CHUNKED, nfront) ; /* ------------------- fill the IVL object ------------------- */ for ( J = 0 ; J < nfront ; J++ ) { if ( (nJ = FrontMtx_frontSize(frontmtx, J)) > 0 ) { FrontMtx_rowIndices(frontmtx, J, &nrow, &rowind) ; if ( nrow > 0 ) { mark[J] = J ; count = 0 ; list[count++] = J ; for ( ii = nJ ; ii < nrow ; ii++ ) { K = rowmap[rowind[ii]] ; if ( mark[K] != J ) { mark[K] = J ; list[count++] = K ; } } IVL_setList(lowerblockIVL, J, count, list) ; } } } /* ------------------------ free the working storage ------------------------ */ IVfree(mark) ; IVfree(list) ; return(lowerblockIVL) ; }
/*--------------------------------------------------------------------*/ int main ( int argc, char *argv[] ) /* ------------------------------------------------- this program tests the IVL_MPI_allgather() method (1) each process generates the same owners[n] map (2) each process creates an IVL object and fills its owned lists with random numbers (3) the processes gather-all's the lists of ivl created -- 98apr03, cca ------------------------------------------------- */ { char *buffer ; double chksum, globalsum, t1, t2 ; Drand drand ; int ilist, length, myid, msglvl, nlist, nproc, rc, seed, size, tag ; int *list, *owners, *vec ; int stats[4], tstats[4] ; IV *ownersIV ; IVL *ivl ; FILE *msgFile ; /* --------------------------------------------------------------- find out the identity of this process and the number of process --------------------------------------------------------------- */ MPI_Init(&argc, &argv) ; MPI_Comm_rank(MPI_COMM_WORLD, &myid) ; MPI_Comm_size(MPI_COMM_WORLD, &nproc) ; fprintf(stdout, "\n process %d of %d, argc = %d", myid, nproc, argc) ; fflush(stdout) ; if ( argc != 5 ) { fprintf(stdout, "\n\n usage : %s msglvl msgFile n seed " "\n msglvl -- message level" "\n msgFile -- message file" "\n nlist -- number of lists in the IVL object" "\n seed -- random number seed" "\n", argv[0]) ; return(0) ; } msglvl = atoi(argv[1]) ; if ( strcmp(argv[2], "stdout") == 0 ) { msgFile = stdout ; } else { length = strlen(argv[2]) + 1 + 4 ; buffer = CVinit(length, '\0') ; sprintf(buffer, "%s.%d", argv[2], myid) ; if ( (msgFile = fopen(buffer, "w")) == NULL ) { fprintf(stderr, "\n fatal error in %s" "\n unable to open file %s\n", argv[0], argv[2]) ; return(-1) ; } CVfree(buffer) ; } nlist = atoi(argv[3]) ; seed = atoi(argv[4]) ; fprintf(msgFile, "\n %s " "\n msglvl -- %d" "\n msgFile -- %s" "\n nlist -- %d" "\n seed -- %d" "\n", argv[0], msglvl, argv[2], nlist, seed) ; fflush(msgFile) ; /* ---------------------------- generate the ownersIV object ---------------------------- */ MARKTIME(t1) ; ownersIV = IV_new() ; IV_init(ownersIV, nlist, NULL) ; owners = IV_entries(ownersIV) ; Drand_setDefaultFields(&drand) ; Drand_setSeed(&drand, seed) ; Drand_setUniform(&drand, 0, nproc) ; Drand_fillIvector(&drand, nlist, owners) ; MARKTIME(t2) ; fprintf(msgFile, "\n CPU %8.3f : initialize the ownersIV object", t2 - t1) ; fflush(msgFile) ; fprintf(msgFile, "\n\n ownersIV generated") ; if ( msglvl > 2 ) { IV_writeForHumanEye(ownersIV, msgFile) ; } else { IV_writeStats(ownersIV, msgFile) ; } fflush(msgFile) ; /* -------------------------------------------- set up the IVL object and fill owned entries -------------------------------------------- */ MARKTIME(t1) ; ivl = IVL_new() ; IVL_init1(ivl, IVL_CHUNKED, nlist) ; vec = IVinit(nlist, -1) ; Drand_setSeed(&drand, seed + myid) ; Drand_setUniform(&drand, 0, nlist) ; for ( ilist = 0 ; ilist < nlist ; ilist++ ) { if ( owners[ilist] == myid ) { size = (int) Drand_value(&drand) ; Drand_fillIvector(&drand, size, vec) ; IVL_setList(ivl, ilist, size, vec) ; } } MARKTIME(t2) ; fprintf(msgFile, "\n CPU %8.3f : initialize the IVL object", t2 - t1) ; fflush(msgFile) ; if ( msglvl > 2 ) { IVL_writeForHumanEye(ivl, msgFile) ; } else { IVL_writeStats(ivl, msgFile) ; } fflush(msgFile) ; /* -------------------------------------------- compute the local checksum of the ivl object -------------------------------------------- */ for ( ilist = 0, chksum = 0.0 ; ilist < nlist ; ilist++ ) { if ( owners[ilist] == myid ) { IVL_listAndSize(ivl, ilist, &size, &list) ; chksum += 1 + ilist + size + IVsum(size, list) ; } } fprintf(msgFile, "\n\n local partial chksum = %12.4e", chksum) ; fflush(msgFile) ; /* ----------------------- get the global checksum ----------------------- */ rc = MPI_Allreduce((void *) &chksum, (void *) &globalsum, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD) ; /* -------------------------------- execute the all-gather operation -------------------------------- */ tag = 47 ; IVzero(4, stats) ; IVL_MPI_allgather(ivl, ownersIV, stats, msglvl, msgFile, tag, MPI_COMM_WORLD) ; if ( msglvl > 0 ) { fprintf(msgFile, "\n\n return from IVL_MPI_allgather()") ; fprintf(msgFile, "\n local send stats : %10d messages with %10d bytes" "\n local recv stats : %10d messages with %10d bytes", stats[0], stats[2], stats[1], stats[3]) ; fflush(msgFile) ; } MPI_Reduce((void *) stats, (void *) tstats, 4, MPI_INT, MPI_SUM, 0, MPI_COMM_WORLD) ; if ( myid == 0 ) { fprintf(msgFile, "\n total send stats : %10d messages with %10d bytes" "\n total recv stats : %10d messages with %10d bytes", tstats[0], tstats[2], tstats[1], tstats[3]) ; fflush(msgFile) ; } if ( msglvl > 2 ) { fprintf(msgFile, "\n\n ivl") ; IVL_writeForHumanEye(ivl, msgFile) ; fflush(msgFile) ; } /* ----------------------------------------- compute the checksum of the entire object ----------------------------------------- */ for ( ilist = 0, chksum = 0.0 ; ilist < nlist ; ilist++ ) { IVL_listAndSize(ivl, ilist, &size, &list) ; chksum += 1 + ilist + size + IVsum(size, list) ; } fprintf(msgFile, "\n globalsum = %12.4e, chksum = %12.4e, error = %12.4e", globalsum, chksum, fabs(globalsum - chksum)) ; fflush(msgFile) ; /* ---------------- free the objects ---------------- */ IV_free(ownersIV) ; IVL_free(ivl) ; /* ------------------------ exit the MPI environment ------------------------ */ MPI_Finalize() ; fprintf(msgFile, "\n") ; fclose(msgFile) ; return(0) ; }
/* ------------------------------------------------------------------ this method is used during the setup for matrix-vector multiplies. each processor has computed the vertices it needs from other processors, these lists are contained in sendIVL. on return, recvIVL contains the lists of vertices this processor must send to all others. sendIVL -- on input, list[q] contains the vertices needed by this processor that are owned by q recvIVL -- on output, list[q] contains the vertices owned by this processor that are needed by q. note, if NULL on input, a new IVL object is allocated stats[] -- statistics vector stats[0] -- contains # of sends stats[1] -- contains # of receives stats[2] -- contains # of bytes sent stats[3] -- contains # of bytes received firsttag -- first tag for messages, tags in range [firsttag, firsttag+nproc-1] are used return value -- recvIVL created -- 98jul26, cca ------------------------------------------------------------------ */ IVL * IVL_MPI_alltoall ( IVL *sendIVL, IVL *recvIVL, int stats[], int msglvl, FILE *msgFile, int firsttag, MPI_Comm comm ) { int count, destination, lasttag, left, myid, nproc, offset, q, recvcount, right, sendcount, source, tag, tagbound ; int *incounts, *outcounts, *recvvec, *sendvec ; MPI_Status status ; /* --------------- check the input --------------- */ if ( sendIVL == NULL || stats == NULL || (msglvl > 0 && msgFile == NULL) ) { fprintf(msgFile, "\n fatal error in IVL_MPI_alltoall()" "\n bad input\n") ; exit(-1) ; } /* --------------------------------------- get id of self and number of processors --------------------------------------- */ MPI_Comm_rank(comm, &myid) ; MPI_Comm_size(comm, &nproc) ; if ( sendIVL->nlist != nproc ) { fprintf(msgFile, "\n fatal error in IVL_MPI_alltoall()" "\n sendIVL: nproc = %d, nlist = %d\n", nproc, sendIVL->nlist) ; exit(-1) ; } lasttag = firsttag + nproc ; if ( lasttag > (tagbound = maxTagMPI(comm)) ) { fprintf(stderr, "\n fatal error in IVL_MPI_alltoall()" "\n lasttag = %d, tag_bound = %d", lasttag, tagbound) ; exit(-1) ; } if ( recvIVL == NULL ) { recvIVL = IVL_new() ; } else { IVL_clearData(recvIVL) ; } IVL_init1(recvIVL, IVL_CHUNKED, nproc) ; /* ------------------------------------------ outcounts[] is sendIVL->sizes[] incounts[] will be recvIVL->sizes[] fill incounts via a call to MPI_Alltoall() and then initialize the recvIVL lists. ------------------------------------------ */ outcounts = sendIVL->sizes ; incounts = IVinit(nproc, 0) ; MPI_Alltoall((void *) outcounts, 1, MPI_INT, (void *) incounts, 1, MPI_INT, comm) ; for ( q = 0 ; q < nproc ; q++ ) { IVL_setList(recvIVL, q, incounts[q], NULL) ; } IVfree(incounts) ; /* --------------------------------------------------- load list myid of sendIVL into list myid of recvIVL --------------------------------------------------- */ IVL_listAndSize(sendIVL, myid, &sendcount, &sendvec) ; IVL_setList(recvIVL, myid, sendcount, sendvec) ; /* --------------------------------------------------------- now loop over the processes, send and receive information --------------------------------------------------------- */ for ( offset = 1, tag = firsttag ; offset < nproc ; offset++, tag++ ) { right = (myid + offset) % nproc ; if ( offset <= myid ) { left = myid - offset ; } else { left = nproc + myid - offset ; } IVL_listAndSize(sendIVL, right, &sendcount, &sendvec) ; IVL_listAndSize(recvIVL, left, &recvcount, &recvvec) ; if ( sendcount > 0 ) { destination = right ; stats[0]++ ; stats[2] += sendcount*sizeof(int) ; } else { destination = MPI_PROC_NULL ; } if ( recvcount > 0 ) { source = left ; stats[1]++ ; stats[3] += recvcount*sizeof(int) ; } else { source = MPI_PROC_NULL ; } if ( msglvl > 2 ) { fprintf(msgFile, "\n offset %d, recvcount %d, source %d, sendcount %d, destination %d", offset, recvcount, source, sendcount, destination) ; fflush(msgFile) ; } /* ----------------- do a send/receive ----------------- */ MPI_Sendrecv((void *) sendvec, sendcount, MPI_INT, destination, tag, (void *) recvvec, recvcount, MPI_INT, source, tag, comm, &status) ; if ( source != MPI_PROC_NULL ) { MPI_Get_count(&status, MPI_INT, &count) ; if ( count != recvcount ) { fprintf(stderr, "\n fatal error in IVL_MPI_alltoall()" "\n proc %d : source %d, count %d, recvcount %d\n", myid, source, count, recvcount) ; exit(-1) ; } } if ( msglvl > 2 ) { fprintf(msgFile, "\n send/recv completed") ; fflush(msgFile) ; } } return(recvIVL) ; }
/* -------------------------------------------------------------------- fill *pndom with ndom, the number of domains. fill *pnseg with nseg, the number of segments. domains are numbered in [0, ndom), segments in [ndom,ndom+nseg). return -- an IV object that contains the map from vertices to segments created -- 99feb29, cca -------------------------------------------------------------------- */ IV * GPart_domSegMap ( GPart *gpart, int *pndom, int *pnseg ) { FILE *msgFile ; Graph *g ; int adjdom, count, d, first, ierr, ii, jj1, jj2, last, ndom, msglvl, nextphi, nPhi, nPsi, nV, phi, phi0, phi1, phi2, phi3, psi, sigma, size, size0, size1, size2, v, vsize, w ; int *adj, *adj0, *adj1, *adj2, *compids, *dmark, *dsmap, *head, *link, *list, *offsets, *PhiToPsi, *PhiToV, *PsiToSigma, *vadj, *VtoPhi ; IV *dsmapIV ; IVL *PhiByPhi, *PhiByPowD, *PsiByPowD ; /* -------------------- set the initial time -------------------- */ icputimes = 0 ; MARKTIME(cputimes[icputimes]) ; /* --------------- check the input --------------- */ if ( gpart == NULL || (g = gpart->g) == NULL || pndom == NULL || pnseg == NULL ) { fprintf(stderr, "\n fatal error in GPart_domSegMap(%p,%p,%p)" "\n bad input\n", gpart, pndom, pnseg) ; exit(-1) ; } compids = IV_entries(&gpart->compidsIV) ; msglvl = gpart->msglvl ; msgFile = gpart->msgFile ; /* ------------------------ create the map IV object ------------------------ */ nV = g->nvtx ; dsmapIV = IV_new() ; IV_init(dsmapIV, nV, NULL) ; dsmap = IV_entries(dsmapIV) ; /* ---------------------------------- check compids[] and get the number of domains and interface vertices ---------------------------------- */ icputimes++ ; MARKTIME(cputimes[icputimes]) ; ndom = nPhi = 0 ; for ( v = 0 ; v < nV ; v++ ) { if ( (d = compids[v]) < 0 ) { fprintf(stderr, "\n fatal error in GPart_domSegMap(%p,%p,%p)" "\n compids[%d] = %d\n", gpart, pndom, pnseg, v, compids[v]) ; exit(-1) ; } else if ( d == 0 ) { nPhi++ ; } else if ( ndom < d ) { ndom = d ; } } *pndom = ndom ; if ( msglvl > 1 ) { fprintf(msgFile, "\n\n Inside GPart_domSegMap") ; fprintf(msgFile, "\n %d domains, %d Phi vertices", ndom, nPhi) ; } if ( ndom == 1 ) { IVfill(nV, dsmap, 0) ; *pndom = 1 ; *pnseg = 0 ; return(dsmapIV) ; } /* -------------------------------- get the maps PhiToV : [0,nPhi) |---> [0,nV) VtoPhi : [0,nV) |---> [0,nPhi) -------------------------------- */ icputimes++ ; MARKTIME(cputimes[icputimes]) ; PhiToV = IVinit(nPhi, -1) ; VtoPhi = IVinit(nV, -1) ; for ( v = 0, phi = 0 ; v < nV ; v++ ) { if ( (d = compids[v]) == 0 ) { PhiToV[phi] = v ; VtoPhi[v] = phi++ ; } } if ( phi != nPhi ) { fprintf(stderr, "\n fatal error in GPart_domSegMap(%p,%p,%p)" "\n phi = %d != %d = nPhi\n", gpart, pndom, pnseg, phi, nPhi) ; exit(-1) ; } if ( msglvl > 2 ) { fprintf(msgFile, "\n PhiToV(%d) :", nPhi) ; IVfp80(msgFile, nPhi, PhiToV, 15, &ierr) ; fflush(msgFile) ; } if ( msglvl > 3 ) { fprintf(msgFile, "\n VtoPhi(%d) :", nV) ; IVfp80(msgFile, nV, VtoPhi, 15, &ierr) ; fflush(msgFile) ; } /* --------------------------------------------------- create an IVL object, PhiByPowD, to hold lists from the interface vertices to their adjacent domains --------------------------------------------------- */ icputimes++ ; MARKTIME(cputimes[icputimes]) ; dmark = IVinit(ndom+1, -1) ; if ( nPhi >= ndom ) { list = IVinit(nPhi, -1) ; } else { list = IVinit(ndom, -1) ; } PhiByPowD = IVL_new() ; IVL_init1(PhiByPowD, IVL_CHUNKED, nPhi) ; for ( phi = 0 ; phi < nPhi ; phi++ ) { v = PhiToV[phi] ; Graph_adjAndSize(g, v, &vsize, &vadj) ; /* if ( phi == 0 ) { int ierr ; fprintf(msgFile, "\n adj(%d,%d) = ", v, phi) ; IVfp80(msgFile, vsize, vadj, 15, &ierr) ; fflush(msgFile) ; } */ count = 0 ; for ( ii = 0 ; ii < vsize ; ii++ ) { if ( (w = vadj[ii]) < nV && (d = compids[w]) > 0 && dmark[d] != phi ) { dmark[d] = phi ; list[count++] = d ; } } if ( count > 0 ) { IVqsortUp(count, list) ; IVL_setList(PhiByPowD, phi, count, list) ; } } if ( msglvl > 2 ) { fprintf(msgFile, "\n PhiByPowD : interface x adjacent domains") ; IVL_writeForHumanEye(PhiByPowD, msgFile) ; fflush(msgFile) ; } /* ------------------------------------------------------- create an IVL object, PhiByPhi to hold lists from the interface vertices to interface vertices. (s,t) are in the list if (s,t) is an edge in the graph and s and t do not share an adjacent domain ------------------------------------------------------- */ icputimes++ ; MARKTIME(cputimes[icputimes]) ; PhiByPhi = IVL_new() ; IVL_init1(PhiByPhi, IVL_CHUNKED, nPhi) ; offsets = IVinit(nPhi, 0) ; head = IVinit(nPhi, -1) ; link = IVinit(nPhi, -1) ; for ( phi1 = 0 ; phi1 < nPhi ; phi1++ ) { count = 0 ; if ( msglvl > 2 ) { v = PhiToV[phi1] ; Graph_adjAndSize(g, v, &vsize, &vadj) ; fprintf(msgFile, "\n checking out phi = %d, v = %d", phi1, v) ; fprintf(msgFile, "\n adj(%d) : ", v) ; IVfp80(msgFile, vsize, vadj, 10, &ierr) ; } /* ------------------------------------------------------------- get (phi1, phi0) edges that were previously put into the list ------------------------------------------------------------- */ if ( msglvl > 3 ) { if ( head[phi1] == -1 ) { fprintf(msgFile, "\n no previous edges") ; } else { fprintf(msgFile, "\n previous edges :") ; } } for ( phi0 = head[phi1] ; phi0 != -1 ; phi0 = nextphi ) { if ( msglvl > 3 ) { fprintf(msgFile, " %d", phi0) ; } nextphi = link[phi0] ; list[count++] = phi0 ; IVL_listAndSize(PhiByPhi, phi0, &size0, &adj0) ; if ( (ii = ++offsets[phi0]) < size0 ) { /* ---------------------------- link phi0 into the next list ---------------------------- */ phi2 = adj0[ii] ; link[phi0] = head[phi2] ; head[phi2] = phi0 ; } } /* -------------------------- get new edges (phi1, phi2) -------------------------- */ IVL_listAndSize(PhiByPowD, phi1, &size1, &adj1) ; v = PhiToV[phi1] ; Graph_adjAndSize(g, v, &vsize, &vadj) ; for ( ii = 0 ; ii < vsize ; ii++ ) { if ( (w = vadj[ii]) < nV && compids[w] == 0 && (phi2 = VtoPhi[w]) > phi1 ) { if ( msglvl > 3 ) { fprintf(msgFile, "\n checking out phi2 = %d", phi2) ; } /* -------------------------------------------------- see if phi1 and phi2 have a common adjacent domain -------------------------------------------------- */ IVL_listAndSize(PhiByPowD, phi2, &size2, &adj2) ; adjdom = 0 ; jj1 = jj2 = 0 ; while ( jj1 < size1 && jj2 < size2 ) { if ( adj1[jj1] < adj2[jj2] ) { jj1++ ; } else if ( adj1[jj1] > adj2[jj2] ) { jj2++ ; } else { if ( msglvl > 3 ) { fprintf(msgFile, ", common adj domain %d", adj1[jj1]) ; } adjdom = 1 ; break ; } } if ( adjdom == 0 ) { if ( msglvl > 3 ) { fprintf(msgFile, ", no adjacent domain") ; } list[count++] = phi2 ; } } } if ( count > 0 ) { /* --------------------- set the list for phi1 --------------------- */ IVqsortUp(count, list) ; IVL_setList(PhiByPhi, phi1, count, list) ; if ( msglvl > 2 ) { fprintf(msgFile, "\n edge list for %d :", phi1) ; IVfp80(msgFile, count, list, 15, &ierr) ; } for ( ii = 0, phi2 = -1 ; ii < count ; ii++ ) { if ( list[ii] > phi1 ) { offsets[phi1] = ii ; phi2 = list[ii] ; break ; } } if ( phi2 != -1 ) { if ( msglvl > 2 ) { fprintf(msgFile, "\n linking %d into list for %d", phi1, phi2) ; } link[phi1] = head[phi2] ; head[phi2] = phi1 ; } /* phi2 = list[0] ; link[phi1] = head[phi2] ; head[phi2] = phi1 ; */ } } if ( msglvl > 2 ) { fprintf(msgFile, "\n PhiByPhi : interface x interface") ; IVL_writeForHumanEye(PhiByPhi, msgFile) ; fflush(msgFile) ; } /* -------------------- get the PhiToPsi map -------------------- */ icputimes++ ; MARKTIME(cputimes[icputimes]) ; PhiToPsi = IVinit(nPhi, -1) ; nPsi = 0 ; for ( phi = 0 ; phi < nPhi ; phi++ ) { if ( PhiToPsi[phi] == -1 ) { /* --------------------------- phi not yet mapped to a psi --------------------------- */ first = last = 0 ; list[0] = phi ; PhiToPsi[phi] = nPsi ; while ( first <= last ) { phi2 = list[first++] ; IVL_listAndSize(PhiByPhi, phi2, &size, &adj) ; for ( ii = 0 ; ii < size ; ii++ ) { phi3 = adj[ii] ; if ( PhiToPsi[phi3] == -1 ) { PhiToPsi[phi3] = nPsi ; list[++last] = phi3 ; } } } nPsi++ ; } } if ( msglvl > 1 ) { fprintf(msgFile, "\n nPsi = %d", nPsi) ; fflush(msgFile) ; } if ( msglvl > 2 ) { fprintf(msgFile, "\n PhiToPsi(%d) :", nPhi) ; IVfp80(msgFile, nPhi, PhiToPsi, 15, &ierr) ; fflush(msgFile) ; } /* --------------------------------- create an IVL object, Psi --> 2^D --------------------------------- */ icputimes++ ; MARKTIME(cputimes[icputimes]) ; IVfill(nPsi, head, -1) ; IVfill(nPhi, link, -1) ; for ( phi = 0 ; phi < nPhi ; phi++ ) { psi = PhiToPsi[phi] ; link[phi] = head[psi] ; head[psi] = phi ; } PsiByPowD = IVL_new() ; IVL_init1(PsiByPowD, IVL_CHUNKED, nPsi) ; IVfill(ndom+1, dmark, -1) ; for ( psi = 0 ; psi < nPsi ; psi++ ) { count = 0 ; for ( phi = head[psi] ; phi != -1 ; phi = link[phi] ) { v = PhiToV[phi] ; Graph_adjAndSize(g, v, &size, &adj) ; for ( ii = 0 ; ii < size ; ii++ ) { if ( (w = adj[ii]) < nV && (d = compids[w]) > 0 && dmark[d] != psi ) { dmark[d] = psi ; list[count++] = d ; } } } if ( count > 0 ) { IVqsortUp(count, list) ; IVL_setList(PsiByPowD, psi, count, list) ; } } if ( msglvl > 2 ) { fprintf(msgFile, "\n PsiByPowD(%d) :", nPhi) ; IVL_writeForHumanEye(PsiByPowD, msgFile) ; fflush(msgFile) ; } icputimes++ ; MARKTIME(cputimes[icputimes]) ; /* ------------------------------------- now get the map Psi |---> Sigma that is the equivalence map over PhiByPowD ------------------------------------- */ icputimes++ ; MARKTIME(cputimes[icputimes]) ; PsiToSigma = IVL_equivMap1(PsiByPowD) ; *pnseg = 1 + IVmax(nPsi, PsiToSigma, &ii) ; if ( msglvl > 2 ) { fprintf(msgFile, "\n nSigma = %d", *pnseg) ; fprintf(msgFile, "\n PsiToSigma(%d) :", nPhi) ; IVfp80(msgFile, nPsi, PsiToSigma, 15, &ierr) ; fflush(msgFile) ; } /* -------------------------------------------------------------- now fill the map from the vertices to the domains and segments -------------------------------------------------------------- */ icputimes++ ; MARKTIME(cputimes[icputimes]) ; for ( v = 0 ; v < nV ; v++ ) { if ( (d = compids[v]) > 0 ) { dsmap[v] = d - 1 ; } else { phi = VtoPhi[v] ; psi = PhiToPsi[phi] ; sigma = PsiToSigma[psi] ; dsmap[v] = ndom + sigma ; } } /* ------------------------ free the working storage ------------------------ */ icputimes++ ; MARKTIME(cputimes[icputimes]) ; IVL_free(PhiByPhi) ; IVL_free(PhiByPowD) ; IVL_free(PsiByPowD) ; IVfree(PhiToV) ; IVfree(VtoPhi) ; IVfree(dmark) ; IVfree(list) ; IVfree(PhiToPsi) ; IVfree(head) ; IVfree(link) ; IVfree(offsets) ; IVfree(PsiToSigma) ; icputimes++ ; MARKTIME(cputimes[icputimes]) ; if ( msglvl > 1 ) { fprintf(msgFile, "\n domain/segment map timings split") ; fprintf(msgFile, "\n %9.5f : create the DSmap object" "\n %9.5f : get numbers of domain and interface vertices" "\n %9.5f : generate PhiToV and VtoPhi" "\n %9.5f : generate PhiByPowD" "\n %9.5f : generate PhiByPhi" "\n %9.5f : generate PhiToPsi" "\n %9.5f : generate PsiByPowD" "\n %9.5f : generate PsiToSigma" "\n %9.5f : generate dsmap" "\n %9.5f : free working storage" "\n %9.5f : total time", cputimes[1] - cputimes[0], cputimes[2] - cputimes[1], cputimes[3] - cputimes[2], cputimes[4] - cputimes[3], cputimes[5] - cputimes[4], cputimes[6] - cputimes[5], cputimes[7] - cputimes[6], cputimes[8] - cputimes[7], cputimes[9] - cputimes[8], cputimes[10] - cputimes[9], cputimes[11] - cputimes[0]) ; } return(dsmapIV) ; }
/* ---------------------------------- return an IVL object that contains the adjacency structure of A^TA. created -- 98jan28, cca ---------------------------------- */ IVL * InpMtx_adjForATA ( InpMtx *inpmtxA ) { InpMtx *inpmtxATA ; int firstcol, firstrow, irow, jvtx, lastcol, lastrow, loc, ncol, nent, nrow, size ; int *ind, *ivec1, *ivec2 ; IVL *adjIVL ; /* --------------- check the input --------------- */ if ( inpmtxA == NULL ) { fprintf(stderr, "\n fatal error in InpMtx_adjForATA(%p)" "\n NULL input\n", inpmtxA) ; exit(-1) ; } /* ---------------------------------------------------------- change the coordinate type and storage mode to row vectors ---------------------------------------------------------- */ InpMtx_changeCoordType(inpmtxA, INPMTX_BY_ROWS) ; InpMtx_changeStorageMode(inpmtxA, INPMTX_BY_VECTORS) ; nent = InpMtx_nent(inpmtxA) ; ivec1 = InpMtx_ivec1(inpmtxA) ; ivec2 = InpMtx_ivec2(inpmtxA) ; firstrow = IVmin(nent, ivec1, &loc) ; lastrow = IVmax(nent, ivec1, &loc) ; firstcol = IVmin(nent, ivec2, &loc) ; lastcol = IVmax(nent, ivec2, &loc) ; if ( firstrow < 0 || firstcol < 0 ) { fprintf(stderr, "\n fatal error" "\n firstrow = %d, firstcol = %d" "\n lastrow = %d, lastcol = %d", firstrow, firstcol, lastrow, lastcol) ; exit(-1) ; } nrow = 1 + lastrow ; ncol = 1 + lastcol ; /* ----------------------------------------------------------- create the new InpMtx object to hold the structure of A^TA ----------------------------------------------------------- */ inpmtxATA = InpMtx_new() ; InpMtx_init(inpmtxATA, INPMTX_BY_ROWS, INPMTX_INDICES_ONLY, 0, 0) ; for ( irow = 0 ; irow < nrow ; irow++ ) { InpMtx_vector(inpmtxA, irow, &size, &ind) ; InpMtx_inputMatrix(inpmtxATA, size, size, 1, size, ind, ind) ; } for ( jvtx = 0 ; jvtx < nrow ; jvtx++ ) { InpMtx_inputEntry(inpmtxATA, jvtx, jvtx) ; } InpMtx_changeStorageMode(inpmtxATA, INPMTX_BY_VECTORS) ; /* ------------------- fill the IVL object ------------------- */ adjIVL = IVL_new() ; IVL_init1(adjIVL, IVL_CHUNKED, nrow) ; for ( jvtx = 0 ; jvtx < ncol ; jvtx++ ) { InpMtx_vector(inpmtxATA, jvtx, &size, &ind) ; IVL_setList(adjIVL, jvtx, size, ind) ; } /* ------------------------------ free the working InpMtx object ------------------------------ */ InpMtx_free(inpmtxATA) ; return(adjIVL) ; }
/* --------------------------------------------- basic initializer for the Graph object type -- graph type nvtx -- # of vertices nvbnd -- # of boundary vertices nedges -- # of edges adjType -- IVL type for adjacency object ewghtType -- IVL type for edge weights object created -- 95sep27, cca --------------------------------------------- */ void Graph_init1 ( Graph *g, int type, int nvtx, int nvbnd, int nedges, int adjType, int ewghtType ) { int nvtot ; #if MYTRACE > 0 fprintf(stdout, "\n just inside Graph_init1(%p,%d,%d,%d,%d,%d,%d)", g, type, nvtx, nvbnd, nedges, adjType, ewghtType) ; fflush(stdout) ; #endif /* --------------- check the input --------------- */ if ( g == NULL ) { fprintf(stderr, "\n fatal error in Graph_init1(%p,%d,%d,%d,%d,%d,%d)" "\n graph is NULL\n", g, type, nvtx, nvbnd, nedges, adjType, ewghtType) ; exit(-1) ; } if ( type < 0 || type >= 4 ) { fprintf(stderr, "\n fatal error in Graph_init1(%p,%d,%d,%d,%d,%d,%d)" "\n invalid type = %d, must be in [0,3]\n", g, type, nvtx, nvbnd, nedges, adjType, ewghtType, type) ; exit(-1) ; } if ( nvtx <= 0 ) { fprintf(stderr, "\n fatal error in Graph_init1(%p,%d,%d,%d,%d,%d,%d)" "\n nvtx = %d, must be positive\n", g, type, nvtx, nvbnd, nedges, adjType, ewghtType, nvtx) ; exit(-1) ; } if ( nvbnd < 0 ) { fprintf(stderr, "\n fatal error in Graph_init1(%p,%d,%d,%d,%d,%d,%d)" "\n nvbnd = %d, must be nonnegative\n", g, type, nvtx, nvbnd, nedges, adjType, ewghtType, nvbnd) ; exit(-1) ; } if ( nedges < 0 ) { fprintf(stderr, "\n fatal error in Graph_init1(%p,%d,%d,%d,%d,%d,%d)" "\n nedges = %d, must be nonnegative\n", g, type, nvtx, nvbnd, nedges, adjType, ewghtType, nedges) ; exit(-1) ; } #if MYTRACE > 0 fprintf(stdout, "\n input checks out") ; fflush(stdout) ; #endif switch ( adjType ) { case IVL_CHUNKED : case IVL_SOLO : case IVL_UNKNOWN : break ; default : fprintf(stderr, "\n fatal error in Graph_init1(%p,%d,%d,%d,%d,%d,%d)" "\n invalid adjType = %d\n", g, type, nvtx, nvbnd, nedges, adjType, ewghtType, adjType) ; exit(-1) ; } switch ( ewghtType ) { case IVL_CHUNKED : case IVL_SOLO : case IVL_UNKNOWN : break ; default : fprintf(stderr, "\n fatal error in Graph_init1(%p,%d,%d,%d,%d,%d,%d)" "\n invalid ewghtType = %d\n", g, type, nvtx, nvbnd, nedges, adjType, ewghtType, ewghtType) ; exit(-1) ; } /* ------------------------- clear the data structures ------------------------- */ Graph_clearData(g) ; /* --------------------- set the scalar fields --------------------- */ g->type = type ; g->nvtx = nvtx ; g->nvbnd = nvbnd ; g->nedges = nedges ; nvtot = nvtx + nvbnd ; /* ---------------------------- set the adjacency IVL object ---------------------------- */ g->adjIVL = IVL_new() ; if ( nedges == 0 || adjType == IVL_UNKNOWN ) { IVL_init1(g->adjIVL, adjType, nvtot) ; } else { IVL_init2(g->adjIVL, adjType, nvtot, nedges) ; } if ( type % 2 == 1 ) { /* ----------------------------- set the vertex weights vector ----------------------------- */ g->vwghts = IVinit(nvtot, 0) ; } if ( type >= 2 ) { /* ------------------------------- set the edge weights IVL object ------------------------------- */ g->ewghtIVL = IVL_new() ; if ( nedges == 0 || ewghtType == IVL_UNKNOWN ) { IVL_init1(g->ewghtIVL, ewghtType, nvtot) ; } else { IVL_init2(g->ewghtIVL, ewghtType, nvtot, nedges) ; } } #if MYTRACE > 0 fprintf(stdout, "\n leaving Graph_init1(%p,%d,%d,%d,%d,%d,%d)", g, type, nvtx, nvbnd, nedges, adjType, ewghtType) ; fflush(stdout) ; #endif return ; }