Actual source code: ex152.c
petsc-3.6.1 2015-08-06
1: static const char help[] = "Test ParMETIS handling of negative weights.\n\n";
3: /* Test contributed by John Fettig */
5: /*
6: * This file implements two tests for a bug reported in ParMETIS. These tests are not expected to pass without the
7: * following two patches.
8: *
9: * http://petsc.cs.iit.edu/petsc/externalpackages/parmetis-4.0.2/rev/2dd2eae596ac
10: * http://petsc.cs.iit.edu/petsc/externalpackages/parmetis-4.0.2/rev/1c2b9fe39201
11: *
12: * The bug was reported upstream, but has received no action so far.
13: *
14: * http://glaros.dtc.umn.edu/gkhome/node/837
15: *
16: */
18: #include <petscsys.h>
19: #include <parmetis.h>
21: #define CHKERRQPARMETIS(n) \
22: if (n == METIS_ERROR_INPUT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"ParMETIS error due to wrong inputs and/or options"); \
23: else if (n == METIS_ERROR_MEMORY) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"ParMETIS error due to insufficient memory"); \
24: else if (n == METIS_ERROR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"ParMETIS general error"); \
26: int main(int argc, char *argv[])
27: {
29: PetscBool flg;
30: PetscMPIInt rank, size;
31: int i, status;
32: idx_t ni,isize,*vtxdist, *xadj, *adjncy, *vwgt, *part;
33: idx_t wgtflag=0, numflag=0, ncon=1, ndims=3, edgecut=0;
34: idx_t options[5];
35: real_t *xyz, *tpwgts, ubvec[1];
36: MPI_Comm comm;
37: FILE *fp;
38: char fname[PETSC_MAX_PATH_LEN],prefix[PETSC_MAX_PATH_LEN] = "";
40: PetscInitialize(&argc,&argv,NULL,help);
41: #if defined(PETSC_USE_64BIT_INDICES)
42: PetscPrintf(PETSC_COMM_WORLD,"This example only works with 32 bit indices\n");
43: PetscFinalize();
44: return 0;
45: #endif
46: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
47: MPI_Comm_size(PETSC_COMM_WORLD,&size);
49: PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Parmetis test options","");
50: PetscOptionsString("-prefix","Path and prefix of test file","",prefix,prefix,sizeof(prefix),&flg);
51: if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Must specify -prefix");
52: PetscOptionsEnd();
54: PetscMalloc1(size+1,&vtxdist);
56: PetscSNPrintf(fname,sizeof(fname),"%s.%d.graph",prefix,rank);
58: PetscFOpen(PETSC_COMM_SELF,fname,"r",&fp);
60: fread(vtxdist, sizeof(idx_t), size+1, fp);
62: ni = vtxdist[rank+1]-vtxdist[rank];
64: PetscMalloc1(ni+1,&xadj);
66: fread(xadj, sizeof(idx_t), ni+1, fp);
68: PetscMalloc1(xadj[ni],&adjncy);
70: for (i=0; i<ni; i++) fread(&adjncy[xadj[i]], sizeof(idx_t), xadj[i+1]-xadj[i], fp);
72: PetscFClose(PETSC_COMM_SELF,fp);
74: PetscSNPrintf(fname,sizeof(fname),"%s.%d.graph.xyz",prefix,rank);
75: PetscFOpen(PETSC_COMM_SELF,fname,"r",&fp);
77: PetscMalloc3(ni*ndims,&xyz,ni,&part,size,&tpwgts);
79: fread(xyz, sizeof(real_t), ndims*ni, fp);
80: PetscFClose(PETSC_COMM_SELF,fp);
82: vwgt = NULL;
84: for (i = 0; i < size; i++) tpwgts[i] = 1. / size;
85: isize = size;
87: ubvec[0] = 1.05;
88: options[0] = 1;
89: options[1] = 2;
90: options[2] = 15;
91: options[3] = 0;
92: options[4] = 0;
94: MPI_Comm_dup(MPI_COMM_WORLD, &comm);
95: status = ParMETIS_V3_PartGeomKway(vtxdist, xadj, adjncy, vwgt, NULL, &wgtflag, &numflag, &ndims, xyz, &ncon, &isize, tpwgts, ubvec,options, &edgecut, part, &comm);CHKERRQPARMETIS(status);
96: MPI_Comm_free(&comm);
98: PetscFree(vtxdist);
99: PetscFree(xadj);
100: PetscFree(adjncy);
101: PetscFree3(xyz,part,tpwgts);
102: PetscFinalize();
103: return 0;
104: }