Actual source code: ex152.c
petsc-3.12.5 2020-03-29
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: * patches in the PETSc distribution of ParMetis. See parmetis.py
8: *
9: *
10: * The bug was reported upstream, but has received no action so far.
11: *
12: * http://glaros.dtc.umn.edu/gkhome/node/837
13: *
14: */
16: #include <petscsys.h>
17: #include <parmetis.h>
19: #define CHKERRQPARMETIS(n) \
20: if (n == METIS_ERROR_INPUT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"ParMETIS error due to wrong inputs and/or options"); \
21: else if (n == METIS_ERROR_MEMORY) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"ParMETIS error due to insufficient memory"); \
22: else if (n == METIS_ERROR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"ParMETIS general error"); \
24: int main(int argc, char *argv[])
25: {
27: PetscBool flg;
28: PetscMPIInt rank, size;
29: int i, status;
30: idx_t ni,isize,*vtxdist, *xadj, *adjncy, *vwgt, *part;
31: idx_t wgtflag=0, numflag=0, ncon=1, ndims=3, edgecut=0;
32: idx_t options[5];
33: PetscReal *xyz;
34: real_t *sxyz, *tpwgts, ubvec[1];
35: MPI_Comm comm;
36: FILE *fp;
37: char fname[PETSC_MAX_PATH_LEN],prefix[PETSC_MAX_PATH_LEN] = "";
38: size_t red;
40: PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
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 ierr;
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: red = fread(vtxdist, sizeof(idx_t), size+1, fp);if (red != (size_t) (size+1)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"Unable to read from data file");
62: ni = vtxdist[rank+1]-vtxdist[rank];
64: PetscMalloc1(ni+1,&xadj);
66: red = fread(xadj, sizeof(idx_t), ni+1, fp);if (red != (size_t) (ni+1)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"Unable to read from data file");
68: PetscMalloc1(xadj[ni],&adjncy);
70: for (i=0; i<ni; i++) {
71: red = fread(&adjncy[xadj[i]], sizeof(idx_t), xadj[i+1]-xadj[i], fp);if (red != (size_t) (xadj[i+1]-xadj[i])) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"Unable to read from data file");
72: }
74: PetscFClose(PETSC_COMM_SELF,fp);
76: PetscSNPrintf(fname,sizeof(fname),"%s.%d.graph.xyz",prefix,rank);
77: PetscFOpen(PETSC_COMM_SELF,fname,"r",&fp);
79: PetscMalloc3(ni*ndims,&xyz,ni,&part,size,&tpwgts);
80: PetscMalloc1(ni*ndims,&sxyz);
82: red = fread(xyz, sizeof(PetscReal), ndims*ni, fp);if (red != (size_t) (ndims*ni)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"Unable to read from data file");
83: for (i=0; i<ni*ndims; i++) sxyz[i] = (size_t) xyz[i];
85: PetscFClose(PETSC_COMM_SELF,fp);
87: vwgt = NULL;
89: for (i = 0; i < size; i++) tpwgts[i] = 1. / size;
90: isize = size;
92: ubvec[0] = 1.05;
93: options[0] = 0;
94: options[1] = 2;
95: options[2] = 15;
96: options[3] = 0;
97: options[4] = 0;
99: MPI_Comm_dup(MPI_COMM_WORLD, &comm);
100: status = ParMETIS_V3_PartGeomKway(vtxdist, xadj, adjncy, vwgt, NULL, &wgtflag, &numflag, &ndims, sxyz, &ncon, &isize, tpwgts, ubvec,options, &edgecut, part, &comm);CHKERRQPARMETIS(status);
101: MPI_Comm_free(&comm);
103: PetscFree(vtxdist);
104: PetscFree(xadj);
105: PetscFree(adjncy);
106: PetscFree3(xyz,part,tpwgts);
107: PetscFree(sxyz);
108: PetscFinalize();
109: return ierr;
110: }
113: /*TEST
115: build:
116: requires: parmetis
118: test:
119: nsize: 2
120: requires: parmetis datafilespath !complex double !define(PETSC_USE_64BIT_INDICES)
121: args: -prefix ${DATAFILESPATH}/parmetis-test/testnp2
123: test:
124: suffix: 2
125: nsize: 4
126: requires: parmetis datafilespath !complex double !define(PETSC_USE_64BIT_INDICES)
127: args: -prefix ${DATAFILESPATH}/parmetis-test/testnp4
129: TEST*/