00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014 #include <parmetislib.h>
00015
00016
00017
00018
00019
00020 void ParMETIS_SerialNodeND(idxtype *vtxdist, idxtype *xadj, idxtype *adjncy, int *numflag,
00021 int *options, idxtype *order, idxtype *sizes, MPI_Comm *comm)
00022 {
00023 int i, npes, mype, seroptions[10];
00024 CtrlType ctrl;
00025 GraphType *agraph;
00026 idxtype *perm=NULL, *iperm=NULL;
00027 int *sendcount, *displs;
00028
00029 MPI_Comm_size(*comm, &npes);
00030 MPI_Comm_rank(*comm, &mype);
00031
00032 if (!ispow2(npes)) {
00033 if (mype == 0)
00034 printf("Error: The number of processors must be a power of 2!\n");
00035 return;
00036 }
00037
00038 if (*numflag == 1)
00039 ChangeNumbering(vtxdist, xadj, adjncy, order, npes, mype, 1);
00040
00041 SetUpCtrl(&ctrl, npes, options[OPTION_DBGLVL], *comm);
00042
00043 IFSET(ctrl.dbglvl, DBG_TIME, InitTimers(&ctrl));
00044 IFSET(ctrl.dbglvl, DBG_TIME, MPI_Barrier(ctrl.gcomm));
00045 IFSET(ctrl.dbglvl, DBG_TIME, starttimer(ctrl.TotalTmr));
00046
00047 IFSET(ctrl.dbglvl, DBG_TIME, MPI_Barrier(ctrl.gcomm));
00048 IFSET(ctrl.dbglvl, DBG_TIME, starttimer(ctrl.MoveTmr));
00049
00050 agraph = AssembleEntireGraph(&ctrl, vtxdist, xadj, adjncy);
00051
00052 IFSET(ctrl.dbglvl, DBG_TIME, MPI_Barrier(ctrl.gcomm));
00053 IFSET(ctrl.dbglvl, DBG_TIME, stoptimer(ctrl.MoveTmr));
00054
00055
00056 if (mype == 0) {
00057 perm = idxmalloc(agraph->nvtxs, "PAROMETISS: perm");
00058 iperm = idxmalloc(agraph->nvtxs, "PAROMETISS: iperm");
00059
00060 seroptions[0] = 0;
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071 METIS_NodeNDP(agraph->nvtxs, agraph->xadj, agraph->adjncy, npes, seroptions, perm, iperm, sizes);
00072 }
00073
00074 IFSET(ctrl.dbglvl, DBG_TIME, MPI_Barrier(ctrl.gcomm));
00075 IFSET(ctrl.dbglvl, DBG_TIME, starttimer(ctrl.MoveTmr));
00076
00077
00078 MPI_Bcast((void *)sizes, 2*npes, IDX_DATATYPE, 0, ctrl.gcomm);
00079
00080
00081 sendcount = imalloc(npes, "PAROMETISS: sendcount");
00082 displs = imalloc(npes, "PAROMETISS: displs");
00083 for (i=0; i<npes; i++) {
00084 sendcount[i] = vtxdist[i+1]-vtxdist[i];
00085 displs[i] = vtxdist[i];
00086 }
00087
00088 MPI_Scatterv((void *)iperm, sendcount, displs, IDX_DATATYPE, (void *)order, vtxdist[mype+1]-vtxdist[mype], IDX_DATATYPE, 0, ctrl.gcomm);
00089
00090 IFSET(ctrl.dbglvl, DBG_TIME, MPI_Barrier(ctrl.gcomm));
00091 IFSET(ctrl.dbglvl, DBG_TIME, stoptimer(ctrl.MoveTmr));
00092
00093 IFSET(ctrl.dbglvl, DBG_TIME, MPI_Barrier(ctrl.gcomm));
00094 IFSET(ctrl.dbglvl, DBG_TIME, stoptimer(ctrl.TotalTmr));
00095 IFSET(ctrl.dbglvl, DBG_TIME, PrintTimingInfo(&ctrl));
00096 IFSET(ctrl.dbglvl, DBG_TIME, MPI_Barrier(ctrl.gcomm));
00097
00098 GKfree((void **)&agraph->xadj, (void **)&agraph->adjncy, (void **)&perm, (void **)&iperm, (void **)&sendcount, (void **)&displs, LTERM);
00099 free(agraph);
00100 FreeCtrl(&ctrl);
00101
00102 if (*numflag == 1)
00103 ChangeNumbering(vtxdist, xadj, adjncy, order, npes, mype, 0);
00104
00105 }
00106
00107
00108
00109
00110
00111
00112 GraphType *AssembleEntireGraph(CtrlType *ctrl, idxtype *vtxdist, idxtype *xadj, idxtype *adjncy)
00113 {
00114 int i, gnvtxs, nvtxs, gnedges, nedges;
00115 int npes = ctrl->npes, mype = ctrl->mype;
00116 idxtype *axadj, *aadjncy;
00117 int *recvcounts, *displs;
00118 GraphType *agraph;
00119
00120 gnvtxs = vtxdist[npes];
00121 nvtxs = vtxdist[mype+1]-vtxdist[mype];
00122 nedges = xadj[nvtxs];
00123
00124 recvcounts = imalloc(npes, "AssembleGraph: recvcounts");
00125 displs = imalloc(npes+1, "AssembleGraph: displs");
00126
00127
00128 for (i=0; i<nvtxs; i++)
00129 xadj[i] = xadj[i+1]-xadj[i];
00130
00131 axadj = idxmalloc(gnvtxs+1, "AssembleEntireGraph: axadj");
00132
00133 for (i=0; i<npes; i++) {
00134 recvcounts[i] = vtxdist[i+1]-vtxdist[i];
00135 displs[i] = vtxdist[i];
00136 }
00137
00138
00139 MPI_Gatherv((void *)xadj, nvtxs, IDX_DATATYPE, axadj, recvcounts, displs, IDX_DATATYPE, 0, ctrl->comm);
00140
00141 MAKECSR(i, nvtxs, xadj);
00142 MAKECSR(i, gnvtxs, axadj);
00143
00144
00145
00146 MPI_Allgather((void *)(&nedges), 1, MPI_INT, (void *)recvcounts, 1, MPI_INT, ctrl->comm);
00147
00148 displs[0] = 0;
00149 for (i=1; i<npes+1; i++)
00150 displs[i] = displs[i-1] + recvcounts[i-1];
00151 gnedges = displs[npes];
00152
00153 aadjncy = idxmalloc(gnedges, "AssembleEntireGraph: aadjncy");
00154
00155
00156 MPI_Gatherv((void *)adjncy, nedges, IDX_DATATYPE, aadjncy, recvcounts, displs, IDX_DATATYPE, 0, ctrl->comm);
00157
00158
00159
00160 agraph = CreateGraph();
00161 agraph->nvtxs = gnvtxs;
00162 agraph->nedges = gnedges;
00163 agraph->xadj = axadj;
00164 agraph->adjncy = aadjncy;
00165
00166 return agraph;
00167 }