int Zoltan_Verify_Graph(MPI_Comm comm, indextype *vtxdist, indextype *xadj, indextype *adjncy, weighttype *vwgt, weighttype *adjwgt, int vwgt_dim, int ewgt_dim, int graph_type, int check_graph, int output_level) { int ierr, flag, cross_edges = 0, mesg_size, sum; int nprocs, proc, *proclist, errors, global_errors; int *perm=NULL; int free_adjncy_sort=0; ZOLTAN_COMM_OBJ *comm_plan; static char *yo = "Zoltan_Verify_Graph"; char msg[256]; ZOLTAN_GNO_TYPE num_obj; int nrecv; indextype *ptr, *ptr1, *ptr2; indextype global_i, global_j; indextype *sendgno=NULL, *recvgno=NULL, *adjncy_sort=NULL; weighttype *sendwgt, *recvwgt; ZOLTAN_GNO_TYPE num_duplicates, num_singletons; ZOLTAN_GNO_TYPE num_selfs, nedges, global_sum, num_zeros; ZOLTAN_GNO_TYPE i, j, ii, k; MPI_Datatype zoltan_gno_mpi_type; ierr = ZOLTAN_OK; zoltan_gno_mpi_type = Zoltan_mpi_gno_type(); /* Make sure all procs have same value of check_graph. */ MPI_Allreduce(&check_graph, &i, 1, MPI_INT, MPI_MAX, comm); check_graph = i; if (check_graph == 0) /* perform no error checking at all */ return ierr; /* Get number of procs and my rank */ MPI_Comm_size(comm, &nprocs); MPI_Comm_rank(comm, &proc); /* Check number of vertices (objects) */ num_obj = (ZOLTAN_GNO_TYPE)(vtxdist[proc+1] - vtxdist[proc]); MPI_Reduce(&num_obj, &global_sum, 1, zoltan_gno_mpi_type, MPI_SUM, 0, comm); if ((proc==0) && (global_sum==0)){ if (ierr == ZOLTAN_OK) ierr = ZOLTAN_WARN; if (output_level>0) ZOLTAN_PRINT_WARN(proc, yo, "No vertices in graph."); } /* Verify that vertex weights are non-negative */ num_zeros = 0; if (vwgt_dim){ for (i=0; i<num_obj; i++){ sum = 0; for (k=0; k<vwgt_dim; k++){ if (vwgt[i*vwgt_dim+k] < 0) { sprintf(msg, "Negative object weight of " TPL_WGT_SPEC " for object " ZOLTAN_GNO_SPEC ".", vwgt[i*vwgt_dim+k], i); ZOLTAN_PRINT_ERROR(proc, yo, msg); ierr = ZOLTAN_FATAL; } sum += vwgt[i*vwgt_dim+k]; } if (sum == 0){ num_zeros++; if (output_level>1) { sprintf(msg, "Zero vertex (object) weights for object " ZOLTAN_GNO_SPEC ".", i); ZOLTAN_PRINT_WARN(proc, yo, msg); } if (ierr == ZOLTAN_OK) ierr = ZOLTAN_WARN; } } MPI_Reduce(&num_zeros, &global_sum, 1, zoltan_gno_mpi_type, MPI_SUM, 0, comm); if ((proc==0) && (global_sum>0)){ if (ierr == ZOLTAN_OK) ierr = ZOLTAN_WARN; if (output_level>0){ sprintf(msg, ZOLTAN_GNO_SPEC " objects have zero weights.", global_sum); ZOLTAN_PRINT_WARN(proc, yo, msg); } } } /* Check number of edges */ nedges = (ZOLTAN_GNO_TYPE)xadj[num_obj]; MPI_Reduce(&nedges, &global_sum, 1, zoltan_gno_mpi_type, MPI_SUM, 0, comm); if ((proc==0) && (global_sum==0)){ if (ierr == ZOLTAN_OK) ierr = ZOLTAN_WARN; if (output_level>0) ZOLTAN_PRINT_WARN(proc, yo, "No edges in graph."); } /* Verify that edge weights are non-negative */ num_zeros = 0; if (ewgt_dim){ for (j=0; j<nedges; j++){ sum = 0; for (k=0; k<ewgt_dim; k++){ if (adjwgt[j*ewgt_dim+k] < 0) { sprintf(msg, "Negative edge weight of " TPL_WGT_SPEC " in edge " ZOLTAN_GNO_SPEC ".", adjwgt[j*ewgt_dim+k], j); ZOLTAN_PRINT_ERROR(proc, yo, msg); ierr = ZOLTAN_FATAL; } sum += adjwgt[j*ewgt_dim+k]; } if (sum == 0){ num_zeros++; if (output_level>1) { sprintf(msg, "Zero edge (communication) weights for edge " ZOLTAN_GNO_SPEC ".", j); ZOLTAN_PRINT_WARN(proc, yo, msg); } if (ierr == ZOLTAN_OK) ierr = ZOLTAN_WARN; } } MPI_Reduce(&num_zeros, &global_sum, 1, zoltan_gno_mpi_type, MPI_SUM, 0, comm); if ((proc==0) && (global_sum>0)){ if (ierr == ZOLTAN_OK) ierr = ZOLTAN_WARN; if (output_level>0){ sprintf(msg, ZOLTAN_GNO_SPEC " edges have zero weights.", global_sum); ZOLTAN_PRINT_WARN(proc, yo, msg); } } } /* Verify that the graph is symmetric (edge weights, too) */ /* Also check for self-edges and multi-edges */ /* Pre-processing: Check if edge lists are sorted. If not, */ /* make a copy and sort so we can save time in the lookups. */ flag = 0; /* Assume sorted. */ for (i=0; (i<num_obj) && (flag==0); i++){ for (ii=xadj[i]; ii<xadj[i+1]-1; ii++){ if (adjncy[ii] > adjncy[ii+1]){ flag = 1; /* Not sorted. */ break; } } } if (flag){ /* Need to sort. */ adjncy_sort = (indextype *) ZOLTAN_MALLOC(nedges*sizeof(indextype)); perm = (int *) ZOLTAN_MALLOC(nedges*sizeof(int)); free_adjncy_sort = 1; if (nedges && (!adjncy_sort || !perm)){ /* Out of memory. */ ZOLTAN_PRINT_ERROR(proc, yo, "Out of memory."); ierr = ZOLTAN_MEMERR; } for (k=0; k<nedges; k++){ adjncy_sort[k] = adjncy[k]; perm[k] = k; } if (sizeof(indextype) == sizeof(short)){ for (i=0; i<num_obj; i++) Zoltan_quicksort_list_inc_short((short *)adjncy_sort, perm, (int)xadj[i], (int)xadj[i+1]-1); } else if (sizeof(indextype) == sizeof(int)){ for (i=0; i<num_obj; i++) Zoltan_quicksort_list_inc_int((int *)adjncy_sort, perm, (int)xadj[i], (int)xadj[i+1]-1); } else if (sizeof(indextype) == sizeof(long)){ for (i=0; i<num_obj; i++) Zoltan_quicksort_list_inc_long((long *)adjncy_sort, perm, (int)xadj[i], (int)xadj[i+1]-1); } else if (sizeof(indextype) == sizeof(int64_t)){ for (i=0; i<num_obj; i++) Zoltan_quicksort_list_inc_long_long((int64_t*)adjncy_sort, perm, (int)xadj[i], (int)xadj[i+1]-1); } else{ ZOLTAN_PRINT_ERROR(proc, yo, "Error in third party library data type support."); ierr = ZOLTAN_MEMERR; } } else { /* Already sorted. */ adjncy_sort = adjncy; } /* First pass: Check on-processor edges and count # off-proc edges */ cross_edges = 0; num_selfs = 0; num_duplicates = 0; num_singletons = 0; for (i=0; i<num_obj; i++){ if (IS_GLOBAL_GRAPH(graph_type)){ global_i = vtxdist[proc]+i; } else{ /* graph_type == LOCAL_GRAPH */ global_i = i; /* A bit confusingly, global_i = i for local graphs */ } /* Singleton? */ if (xadj[i] == xadj[i+1]){ num_singletons++; if (output_level>1){ sprintf(msg, "Vertex " TPL_IDX_SPEC " has no edges.", global_i); ZOLTAN_PRINT_WARN(proc, yo, msg); } } for (ii=xadj[i]; ii<xadj[i+1]; ii++){ global_j = adjncy_sort[ii]; /* Valid vertex number? */ if ((IS_GLOBAL_GRAPH(graph_type) && ((global_j < vtxdist[0]) || (global_j >= vtxdist[nprocs]))) || (IS_LOCAL_GRAPH(graph_type) && ((global_j < 0) || (global_j >= num_obj)))){ sprintf(msg, "Edge to invalid vertex " TPL_IDX_SPEC " detected.", global_j); ZOLTAN_PRINT_ERROR(proc, yo, msg); ierr = ZOLTAN_FATAL; } /* Self edge? */ if (global_j == global_i){ num_selfs++; if (output_level>1){ sprintf(msg, "Self edge for vertex " TPL_IDX_SPEC " detected.", global_i); ZOLTAN_PRINT_WARN(proc, yo, msg); } } /* Duplicate edge? */ if ((ii+1<xadj[i+1]) && (adjncy_sort[ii]==adjncy_sort[ii+1])){ num_duplicates++; if (output_level>1){ sprintf(msg, "Duplicate edge (" TPL_IDX_SPEC "," TPL_IDX_SPEC ") detected.", global_i, global_j); ZOLTAN_PRINT_WARN(proc, yo, msg); } } /* Is global_j a local vertex? */ if (IS_LOCAL_GRAPH(graph_type) || (IS_GLOBAL_GRAPH(graph_type) && (global_j >= vtxdist[proc]) && (global_j < vtxdist[proc+1]))){ /* Check if (global_j, global_i) is an edge */ if (IS_GLOBAL_GRAPH(graph_type)) j = global_j - vtxdist[proc]; else /* graph_type == LOCAL_GRAPH */ j = global_j; /* Binary search for edge (global_j, global_i) */ ptr = (indextype *)bsearch(&global_i, &adjncy_sort[xadj[j]], (int)(xadj[j+1]-xadj[j]), sizeof(indextype), Zoltan_Compare_Indextypes); if (ptr){ /* OK, found edge (global_j, global_i) */ if ((adjncy_sort==adjncy) && ewgt_dim){ /* Compare weights */ /* EBEB For now, don't compare weights if we sorted edge lists. */ flag = 0; for (k=0; k<ewgt_dim; k++){ if (adjwgt[(ptr-adjncy_sort)*ewgt_dim+k] != adjwgt[ii*ewgt_dim+k]){ /* Numerically nonsymmetric */ flag = -1; if (ierr == ZOLTAN_OK) ierr = ZOLTAN_WARN; } } if (flag<0 && output_level>0){ sprintf(msg, "Graph is numerically nonsymmetric " "in edge (" TPL_IDX_SPEC "," TPL_IDX_SPEC ")", global_i, global_j); ZOLTAN_PRINT_WARN(proc, yo, msg); } } } else { /* bsearch failed */ sprintf(msg, "Graph is not symmetric. " "Edge (" TPL_IDX_SPEC "," TPL_IDX_SPEC ") exists, but no edge (" TPL_IDX_SPEC "," TPL_IDX_SPEC ").", global_i, global_j, global_j, global_i); ZOLTAN_PRINT_ERROR(proc, yo, msg); ierr = ZOLTAN_FATAL; } } else { cross_edges++; } } } /* Sum up warnings so far. */ MPI_Reduce(&num_selfs, &global_sum, 1, zoltan_gno_mpi_type, MPI_SUM, 0, comm); if ((proc==0) && (global_sum>0)){ ierr = ZOLTAN_WARN; if (output_level>0){ sprintf(msg, ZOLTAN_GNO_SPEC " self-edges in graph.", global_sum); ZOLTAN_PRINT_WARN(proc, yo, msg); } } MPI_Reduce(&num_duplicates, &global_sum, 1, zoltan_gno_mpi_type, MPI_SUM, 0, comm); if ((proc==0) && (global_sum>0)){ ierr = ZOLTAN_WARN; if (output_level>0){ sprintf(msg, ZOLTAN_GNO_SPEC " duplicate edges in graph.", global_sum); ZOLTAN_PRINT_WARN(proc, yo, msg); } } MPI_Reduce(&num_singletons, &global_sum, 1, zoltan_gno_mpi_type, MPI_SUM, 0, comm); if ((proc==0) && (global_sum>0)){ ierr = ZOLTAN_WARN; if (output_level>0){ sprintf(msg, ZOLTAN_GNO_SPEC " vertices in the graph are singletons (have no edges).", global_sum); ZOLTAN_PRINT_WARN(proc, yo, msg); } } /* Check if any processor has encountered an error so far */ errors = 0; if (ierr == ZOLTAN_WARN) errors |= 1; if (ierr == ZOLTAN_MEMERR) errors |= 2; if (ierr == ZOLTAN_FATAL) errors |= 4; MPI_Allreduce(&errors, &global_errors, 1, MPI_INT, MPI_BOR, comm); if (global_errors & 4){ /* Fatal error: return now */ if (free_adjncy_sort) ZOLTAN_FREE(&adjncy_sort); if (free_adjncy_sort) ZOLTAN_FREE(&perm); return ZOLTAN_FATAL; } if ((IS_GLOBAL_GRAPH(graph_type)) && (check_graph >= 2)) { /* Test for consistency across processors. */ /* Allocate space for off-proc data */ mesg_size = (2*sizeof(indextype)) + (ewgt_dim * sizeof(weighttype)); sendgno = (indextype *) ZOLTAN_MALLOC(cross_edges*mesg_size); recvgno = (indextype *) ZOLTAN_MALLOC(cross_edges*mesg_size); proclist = (int *) ZOLTAN_MALLOC(cross_edges*sizeof(int)); if (cross_edges && !(sendgno && recvgno && proclist)){ ZOLTAN_PRINT_ERROR(proc, yo, "Out of memory."); ierr = ZOLTAN_MEMERR; } /* Second pass: Copy data to send buffer */ nedges = 0; ptr1 = (indextype *) sendgno; for (i=0; i<num_obj; i++){ global_i = vtxdist[proc]+i; for (ii=xadj[i]; ii<xadj[i+1]; ii++){ global_j = adjncy[ii]; /* Is global_j off-proc? */ if ((global_j < vtxdist[proc]) || (global_j >= vtxdist[proc+1])){ /* Add to list */ k=0; while (global_j >= vtxdist[k+1]) k++; proclist[nedges++] = k; /* Copy (global_i, global_j) and corresponding weights to sendgno */ *ptr1++ = global_i; *ptr1++ = global_j; for (k=0; k<ewgt_dim; k++){ *ptr1++ = adjwgt[ii*ewgt_dim+k]; } } } } /* Do the irregular communication */ ierr = Zoltan_Comm_Create(&comm_plan, cross_edges, proclist, comm, TAG1, &nrecv); if (ierr != ZOLTAN_OK && ierr != ZOLTAN_WARN) { sprintf(msg, "Error %s returned from Zoltan_Comm_Create.", (ierr == ZOLTAN_MEMERR ? "ZOLTAN_MEMERR" : "ZOLTAN_FATAL")); ZOLTAN_PRINT_ERROR(proc, yo, msg); } else { if (nrecv != cross_edges){ sprintf(msg, "Incorrect number of edges to/from proc %d.", proc); ZOLTAN_PRINT_ERROR(proc, yo, msg); ierr = ZOLTAN_FATAL; } ierr = Zoltan_Comm_Do(comm_plan, TAG2, (char *)sendgno, mesg_size, (char *)recvgno); Zoltan_Comm_Destroy(&comm_plan); if (ierr != ZOLTAN_OK && ierr != ZOLTAN_WARN) { sprintf(msg, "Error %s returned from Zoltan_Comm_Do.", (ierr == ZOLTAN_MEMERR ? "ZOLTAN_MEMERR" : "ZOLTAN_FATAL")); ZOLTAN_PRINT_ERROR(proc, yo, msg); } else { /* Third pass: Compare on-proc data to the off-proc data we received */ /* sendgno and recvgno should contain the same data except (i,j) is */ /* (j,i) */ for (i=0, ptr1=sendgno; i<cross_edges; i++){ flag = 0; sendwgt = (weighttype *)(ptr1 + 2); for (j=0, ptr2=recvgno; j<cross_edges; j++){ recvwgt = (weighttype *)(ptr2 + 2); if ((ptr2[0] == ptr1[1]) && (ptr2[1] == ptr1[0])){ /* Found matching edge */ flag = 1; /* Check weights */ for (k=0; k<ewgt_dim; k++){ if (sendwgt[k] != recvwgt[k]){ flag = -1; ierr = ZOLTAN_WARN; } } if (flag<0 && output_level>0){ sprintf(msg, "Edge weight (" TPL_IDX_SPEC "," TPL_IDX_SPEC ") is not symmetric", ptr1[0], ptr1[1]); ZOLTAN_PRINT_WARN(proc, yo, msg); } } ptr2 = (indextype *)(recvwgt + ewgt_dim); } if (!flag){ sprintf(msg, "Graph is not symmetric. " "Edge (" TPL_IDX_SPEC "," TPL_IDX_SPEC ") exists, but not (" TPL_IDX_SPEC "," TPL_IDX_SPEC ").", ptr1[0], ptr1[1], ptr1[1], ptr1[0]); ZOLTAN_PRINT_ERROR(proc, yo, msg); ierr = ZOLTAN_FATAL; } ptr1 = (indextype *)(sendwgt + ewgt_dim); } } } /* Free memory */ ZOLTAN_FREE(&sendgno); ZOLTAN_FREE(&recvgno); ZOLTAN_FREE(&proclist); } if (free_adjncy_sort) ZOLTAN_FREE(&adjncy_sort); if (free_adjncy_sort) ZOLTAN_FREE(&perm); /* Compute global error code */ errors = 0; if (ierr == ZOLTAN_WARN) errors |= 1; if (ierr == ZOLTAN_MEMERR) errors |= 2; if (ierr == ZOLTAN_FATAL) errors |= 4; MPI_Allreduce(&errors, &global_errors, 1, MPI_INT, MPI_BOR, comm); if (global_errors & 4){ return ZOLTAN_FATAL; } else if (global_errors & 2){ return ZOLTAN_MEMERR; } else if (global_errors & 1){ return ZOLTAN_WARN; } else { if (proc==0 && output_level>0){ printf("ZOLTAN %s: The graph is valid with check_graph = %1d\n", yo, check_graph); } return ZOLTAN_OK; } }
int Zoltan_ParMetis_Order( ZZ *zz, /* Zoltan structure */ int num_obj, /* Number of (local) objects to order. */ ZOLTAN_ID_PTR gids, /* List of global ids (local to this proc) */ /* The application must allocate enough space */ ZOLTAN_ID_PTR lids, /* List of local ids (local to this proc) */ /* The application must allocate enough space */ ZOLTAN_ID_PTR rank, /* rank[i] is the rank of gids[i] */ int *iperm, ZOOS *order_opt /* Ordering options, parsed by Zoltan_Order */ ) { static char *yo = "Zoltan_ParMetis_Order"; int i, n, ierr; ZOLTAN_Output_Order ord; ZOLTAN_Third_Graph gr; #ifdef ZOLTAN_PARMETIS MPI_Comm comm = zz->Communicator;/* don't want to risk letting external packages changing our communicator */ #endif indextype numflag = 0; int timer_p = 0; int get_times = 0; int use_timers = 0; double times[5]; ZOLTAN_ID_PTR l_gids = NULL; ZOLTAN_ID_PTR l_lids = NULL; indextype options[MAX_PARMETIS_OPTIONS]; char alg[MAX_PARAM_STRING_LEN]; ZOLTAN_TRACE_ENTER(zz, yo); #ifdef ZOLTAN_PARMETIS #if TPL_USE_DATATYPE != TPL_METIS_DATATYPES #ifdef TPL_FLOAT_WEIGHT i = 1; #else i = 0; #endif if ((sizeof(indextype) != sizeof(idxtype)) || (sizeof(weighttype) != sizeof(idxtype)) || i){ ZOLTAN_THIRD_ERROR(ZOLTAN_FATAL, "Not supported: Multiple 3rd party libraries with incompatible " "data types."); return ZOLTAN_FATAL; } #endif #endif memset(&gr, 0, sizeof(ZOLTAN_Third_Graph)); memset(&ord, 0, sizeof(ZOLTAN_Output_Order)); memset(times, 0, sizeof(times)); ord.order_opt = order_opt; if (!order_opt){ /* If for some reason order_opt is NULL, allocate a new ZOOS here. */ /* This should really never happen. */ order_opt = (ZOOS *) ZOLTAN_MALLOC(sizeof(ZOOS)); strcpy(order_opt->method,"PARMETIS"); } ierr = Zoltan_Parmetis_Parse(zz, options, alg, NULL, NULL, &ord); /* ParMetis only computes the rank vector */ order_opt->return_args = RETURN_RANK; /* Check that num_obj equals the number of objects on this proc. */ /* This constraint may be removed in the future. */ n = zz->Get_Num_Obj(zz->Get_Num_Obj_Data, &ierr); if ((ierr!= ZOLTAN_OK) && (ierr!= ZOLTAN_WARN)){ /* Return error code */ ZOLTAN_PRINT_ERROR(zz->Proc, yo, "Get_Num_Obj returned error."); return(ZOLTAN_FATAL); } if (n != num_obj){ /* Currently this is a fatal error. */ ZOLTAN_PRINT_ERROR(zz->Proc, yo, "Input num_obj does not equal the " "number of objects."); return(ZOLTAN_FATAL); } /* Do not use weights for ordering */ gr.obj_wgt_dim = -1; gr.edge_wgt_dim = -1; gr.num_obj = num_obj; /* Check what ordering type is requested */ if (order_opt){ SET_GLOBAL_GRAPH(&gr.graph_type); /* GLOBAL by default */ #ifdef ZOLTAN_PARMETIS if ((strcmp(order_opt->method, "METIS") == 0)) #endif /* ZOLTAN_PARMETIS */ SET_LOCAL_GRAPH(&gr.graph_type); } gr.get_data = 1; if (IS_LOCAL_GRAPH(gr.graph_type) && zz->Num_Proc > 1) { ZOLTAN_PRINT_ERROR(zz->Proc, yo, "Serial ordering on more than 1 process: " "set ParMetis instead."); return(ZOLTAN_FATAL); } timer_p = Zoltan_Preprocess_Timer(zz, &use_timers); /* Start timer */ get_times = (zz->Debug_Level >= ZOLTAN_DEBUG_ATIME); if (get_times){ MPI_Barrier(zz->Communicator); times[0] = Zoltan_Time(zz->Timer); } ierr = Zoltan_Preprocess_Graph(zz, &l_gids, &l_lids, &gr, NULL, NULL, NULL); if ((ierr != ZOLTAN_OK) && (ierr != ZOLTAN_WARN)) { Zoltan_Third_Exit(&gr, NULL, NULL, NULL, NULL, NULL); return (ierr); } /* Allocate space for separator sizes */ if (IS_GLOBAL_GRAPH(gr.graph_type)) { if (Zoltan_TPL_Order_Init_Tree(&zz->TPL_Order, 2*zz->Num_Proc, zz->Num_Proc) != ZOLTAN_OK) { /* Not enough memory */ Zoltan_Third_Exit(&gr, NULL, NULL, NULL, NULL, &ord); ZOLTAN_THIRD_ERROR(ZOLTAN_MEMERR, "Out of memory."); } ord.sep_sizes = (indextype*)ZOLTAN_MALLOC((2*zz->Num_Proc+1)*sizeof(indextype)); if (ord.sep_sizes == NULL) { Zoltan_Third_Exit(&gr, NULL, NULL, NULL, NULL, &ord); ZOLTAN_THIRD_ERROR(ZOLTAN_MEMERR, "Out of memory."); } memset(ord.sep_sizes, 0, (2*zz->Num_Proc+1)*sizeof(int)); /* It seems parmetis don't initialize correctly */ } /* Allocate space for direct perm */ ord.rank = (indextype *) ZOLTAN_MALLOC(gr.num_obj*sizeof(indextype)); if (!ord.rank){ /* Not enough memory */ Zoltan_Third_Exit(&gr, NULL, NULL, NULL, NULL, &ord); ZOLTAN_THIRD_ERROR(ZOLTAN_MEMERR, "Out of memory."); } if (IS_LOCAL_GRAPH(gr.graph_type)){ /* Allocate space for inverse perm */ ord.iperm = (indextype *) ZOLTAN_MALLOC(gr.num_obj*sizeof(indextype)); if (!ord.iperm){ /* Not enough memory */ Zoltan_Third_Exit(&gr, NULL, NULL, NULL, NULL, &ord); ZOLTAN_THIRD_ERROR(ZOLTAN_MEMERR, "Out of memory."); } } else ord.iperm = NULL; /* Get a time here */ if (get_times) times[1] = Zoltan_Time(zz->Timer); #ifdef ZOLTAN_PARMETIS if (IS_GLOBAL_GRAPH(gr.graph_type)){ ZOLTAN_TRACE_DETAIL(zz, yo, "Calling the ParMETIS library"); ParMETIS_V3_NodeND (gr.vtxdist, gr.xadj, gr.adjncy, &numflag, options, ord.rank, ord.sep_sizes, &comm); ZOLTAN_TRACE_DETAIL(zz, yo, "Returned from the ParMETIS library"); } else #endif /* ZOLTAN_PARMETIS */ #if defined(ZOLTAN_METIS) || defined(ZOLTAN_PARMETIS) if (IS_LOCAL_GRAPH(gr.graph_type)) { /* Be careful : permutation parameters are in the opposite order */ indextype numobj = gr.num_obj; ZOLTAN_TRACE_DETAIL(zz, yo, "Calling the METIS library"); order_opt->return_args = RETURN_RANK|RETURN_IPERM; /* We provide directly all the permutations */ #if !defined(METIS_VER_MAJOR) || METIS_VER_MAJOR < 5 options[0] = 0; /* Use default options for METIS. */ METIS_NodeND(&numobj, gr.xadj, gr.adjncy, &numflag, options, ord.iperm, ord.rank); #else METIS_SetDefaultOptions(options); METIS_NodeND(&numobj, gr.xadj, gr.adjncy, NULL, options, ord.iperm, ord.rank); /* NULL is vwgt -- new interface in v4 */ #endif ZOLTAN_TRACE_DETAIL(zz, yo, "Returned from the METIS library"); } #endif /* ZOLTAN_METIS */ /* Get a time here */ if (get_times) times[2] = Zoltan_Time(zz->Timer); if (IS_GLOBAL_GRAPH(gr.graph_type)){ /* Update Elimination tree */ int numbloc; int start; int leaf; int *converttab; int levelmax; levelmax = mylog2(zz->Num_Proc) + 1; converttab = (int*)ZOLTAN_MALLOC(zz->Num_Proc*2*sizeof(int)); memset(converttab, 0, zz->Num_Proc*2*sizeof(int)); /* Determine the first node in each separator, store it in zz->TPL_Order.start */ for (numbloc = 0, start=0, leaf=0; numbloc < zz->Num_Proc /2; numbloc++) { int father; father = zz->Num_Proc + numbloc; converttab[start] = 2*numbloc; zz->TPL_Order.leaves[leaf++]=start; zz->TPL_Order.ancestor[start] = start + 2; converttab[start+1] = 2*numbloc+1; zz->TPL_Order.leaves[leaf++]=start+1; zz->TPL_Order.ancestor[start+1] = start + 2; start+=2; do { converttab[start] = father; if (father %2 == 0) { int nextoffset; int level; level = mylog2(2*zz->Num_Proc - 1 - father); nextoffset = (1<<(levelmax-level)); zz->TPL_Order.ancestor[start] = start+nextoffset; start++; break; } else { zz->TPL_Order.ancestor[start] = start+1; start++; father = zz->Num_Proc + father/2; } } while (father < 2*zz->Num_Proc - 1); } zz->TPL_Order.start[0] = 0; zz->TPL_Order.ancestor [2*zz->Num_Proc - 2] = -1; for (numbloc = 1 ; numbloc < 2*zz->Num_Proc ; numbloc++) { int oldblock=converttab[numbloc-1]; zz->TPL_Order.start[numbloc] = zz->TPL_Order.start[numbloc-1] + ord.sep_sizes[oldblock]; } ZOLTAN_FREE(&converttab); ZOLTAN_FREE(&ord.sep_sizes); zz->TPL_Order.leaves[zz->Num_Proc] = -1; zz->TPL_Order.nbr_leaves = zz->Num_Proc; zz->TPL_Order.nbr_blocks = 2*zz->Num_Proc-1; } else { /* No tree */ zz->TPL_Order.nbr_blocks = 0; zz->TPL_Order.start = NULL; zz->TPL_Order.ancestor = NULL; zz->TPL_Order.leaves = NULL; } /* Correct because no redistribution */ memcpy(gids, l_gids, n*zz->Num_GID*sizeof(ZOLTAN_ID_TYPE)); memcpy(lids, l_lids, n*zz->Num_LID*sizeof(ZOLTAN_ID_TYPE)); ierr = Zoltan_Postprocess_Graph (zz, l_gids, l_lids, &gr, NULL, NULL, NULL, &ord, NULL); ZOLTAN_FREE(&l_gids); ZOLTAN_FREE(&l_lids); /* Get a time here */ if (get_times) times[3] = Zoltan_Time(zz->Timer); if (get_times) Zoltan_Third_DisplayTime(zz, times); if (use_timers) ZOLTAN_TIMER_STOP(zz->ZTime, timer_p, zz->Communicator); if (sizeof(indextype) == sizeof(ZOLTAN_ID_TYPE)){ memcpy(rank, ord.rank, gr.num_obj*sizeof(indextype)); } else{ for (i=0; i < gr.num_obj; i++){ rank[i] = (ZOLTAN_ID_TYPE)ord.rank[i]; } } if ((ord.iperm != NULL) && (iperm != NULL)){ if (sizeof(indextype) == sizeof(int)){ memcpy(iperm, ord.iperm, gr.num_obj*sizeof(indextype)); } else{ for (i=0; i < gr.num_obj; i++){ iperm[i] = (int)ord.iperm[i]; } } } if (ord.iperm != NULL) ZOLTAN_FREE(&ord.iperm); ZOLTAN_FREE(&ord.rank); /* Free all other "graph" stuff */ Zoltan_Third_Exit(&gr, NULL, NULL, NULL, NULL, NULL); ZOLTAN_TRACE_EXIT(zz, yo); return (ZOLTAN_OK); }