1c4762a1bSJed Brown static const char help[] = "Test ParMETIS handling of negative weights.\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown /* Test contributed by John Fettig */ 4c4762a1bSJed Brown 5c4762a1bSJed Brown /* 6c4762a1bSJed Brown * This file implements two tests for a bug reported in ParMETIS. These tests are not expected to pass without the 7c4762a1bSJed Brown * patches in the PETSc distribution of ParMetis. See parmetis.py 8c4762a1bSJed Brown * 9c4762a1bSJed Brown * 10c4762a1bSJed Brown * The bug was reported upstream, but has received no action so far. 11c4762a1bSJed Brown * 12c4762a1bSJed Brown * http://glaros.dtc.umn.edu/gkhome/node/837 13c4762a1bSJed Brown * 14c4762a1bSJed Brown */ 15c4762a1bSJed Brown 16c4762a1bSJed Brown #include <petscsys.h> 17c4762a1bSJed Brown #include <parmetis.h> 18c4762a1bSJed Brown 19c4762a1bSJed Brown #define CHKERRQPARMETIS(n) \ 20c4762a1bSJed Brown if (n == METIS_ERROR_INPUT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"ParMETIS error due to wrong inputs and/or options"); \ 21c4762a1bSJed Brown else if (n == METIS_ERROR_MEMORY) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"ParMETIS error due to insufficient memory"); \ 22c4762a1bSJed Brown else if (n == METIS_ERROR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"ParMETIS general error"); \ 23c4762a1bSJed Brown 24c4762a1bSJed Brown int main(int argc, char *argv[]) 25c4762a1bSJed Brown { 26c4762a1bSJed Brown PetscErrorCode ierr; 27c4762a1bSJed Brown PetscBool flg; 28c4762a1bSJed Brown PetscMPIInt rank, size; 29c4762a1bSJed Brown int i, status; 30c4762a1bSJed Brown idx_t ni,isize,*vtxdist, *xadj, *adjncy, *vwgt, *part; 31c4762a1bSJed Brown idx_t wgtflag=0, numflag=0, ncon=1, ndims=3, edgecut=0; 32c4762a1bSJed Brown idx_t options[5]; 33c4762a1bSJed Brown PetscReal *xyz; 34c4762a1bSJed Brown real_t *sxyz, *tpwgts, ubvec[1]; 35c4762a1bSJed Brown MPI_Comm comm; 36c4762a1bSJed Brown FILE *fp; 37c4762a1bSJed Brown char fname[PETSC_MAX_PATH_LEN],prefix[PETSC_MAX_PATH_LEN] = ""; 38c4762a1bSJed Brown size_t red; 39c4762a1bSJed Brown 40c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 41c4762a1bSJed Brown #if defined(PETSC_USE_64BIT_INDICES) 421e1ea65dSPierre Jolivet ierr = PetscPrintf(PETSC_COMM_WORLD,"This example only works with 32 bit indices\n");CHKERRQ(ierr); 43c4762a1bSJed Brown ierr = PetscFinalize(); 44c4762a1bSJed Brown return ierr; 45c4762a1bSJed Brown #endif 46c4762a1bSJed Brown MPI_Comm_rank(PETSC_COMM_WORLD,&rank); 47c4762a1bSJed Brown MPI_Comm_size(PETSC_COMM_WORLD,&size); 48c4762a1bSJed Brown 49c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Parmetis test options","");CHKERRQ(ierr); 50c4762a1bSJed Brown ierr = PetscOptionsString("-prefix","Path and prefix of test file","",prefix,prefix,sizeof(prefix),&flg);CHKERRQ(ierr); 51c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Must specify -prefix"); 52c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 53c4762a1bSJed Brown 54c4762a1bSJed Brown ierr = PetscMalloc1(size+1,&vtxdist);CHKERRQ(ierr); 55c4762a1bSJed Brown 56c4762a1bSJed Brown ierr = PetscSNPrintf(fname,sizeof(fname),"%s.%d.graph",prefix,rank);CHKERRQ(ierr); 57c4762a1bSJed Brown 58c4762a1bSJed Brown ierr = PetscFOpen(PETSC_COMM_SELF,fname,"r",&fp);CHKERRQ(ierr); 59c4762a1bSJed Brown 60c4762a1bSJed Brown 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"); 61c4762a1bSJed Brown 62c4762a1bSJed Brown ni = vtxdist[rank+1]-vtxdist[rank]; 63c4762a1bSJed Brown 64c4762a1bSJed Brown ierr = PetscMalloc1(ni+1,&xadj);CHKERRQ(ierr); 65c4762a1bSJed Brown 66c4762a1bSJed Brown 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"); 67c4762a1bSJed Brown 68c4762a1bSJed Brown ierr = PetscMalloc1(xadj[ni],&adjncy);CHKERRQ(ierr); 69c4762a1bSJed Brown 70c4762a1bSJed Brown for (i=0; i<ni; i++) { 71c4762a1bSJed Brown 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"); 72c4762a1bSJed Brown } 73c4762a1bSJed Brown 74c4762a1bSJed Brown ierr = PetscFClose(PETSC_COMM_SELF,fp);CHKERRQ(ierr); 75c4762a1bSJed Brown 76c4762a1bSJed Brown ierr = PetscSNPrintf(fname,sizeof(fname),"%s.%d.graph.xyz",prefix,rank);CHKERRQ(ierr); 77c4762a1bSJed Brown ierr = PetscFOpen(PETSC_COMM_SELF,fname,"r",&fp);CHKERRQ(ierr); 78c4762a1bSJed Brown 79c4762a1bSJed Brown ierr = PetscMalloc3(ni*ndims,&xyz,ni,&part,size,&tpwgts);CHKERRQ(ierr); 80c4762a1bSJed Brown ierr = PetscMalloc1(ni*ndims,&sxyz);CHKERRQ(ierr); 81c4762a1bSJed Brown 82c4762a1bSJed Brown 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"); 83c4762a1bSJed Brown for (i=0; i<ni*ndims; i++) sxyz[i] = (size_t) xyz[i]; 84c4762a1bSJed Brown 85c4762a1bSJed Brown ierr = PetscFClose(PETSC_COMM_SELF,fp);CHKERRQ(ierr); 86c4762a1bSJed Brown 87c4762a1bSJed Brown vwgt = NULL; 88c4762a1bSJed Brown 89c4762a1bSJed Brown for (i = 0; i < size; i++) tpwgts[i] = 1. / size; 90c4762a1bSJed Brown isize = size; 91c4762a1bSJed Brown 92c4762a1bSJed Brown ubvec[0] = 1.05; 93c4762a1bSJed Brown options[0] = 0; 94c4762a1bSJed Brown options[1] = 2; 95c4762a1bSJed Brown options[2] = 15; 96c4762a1bSJed Brown options[3] = 0; 97c4762a1bSJed Brown options[4] = 0; 98c4762a1bSJed Brown 9955b25c41SPierre Jolivet ierr = MPI_Comm_dup(MPI_COMM_WORLD, &comm);CHKERRMPI(ierr); 100c4762a1bSJed Brown status = ParMETIS_V3_PartGeomKway(vtxdist, xadj, adjncy, vwgt, NULL, &wgtflag, &numflag, &ndims, sxyz, &ncon, &isize, tpwgts, ubvec,options, &edgecut, part, &comm);CHKERRQPARMETIS(status); 101ffc4695bSBarry Smith ierr = MPI_Comm_free(&comm);CHKERRMPI(ierr); 102c4762a1bSJed Brown 103c4762a1bSJed Brown ierr = PetscFree(vtxdist);CHKERRQ(ierr); 104c4762a1bSJed Brown ierr = PetscFree(xadj);CHKERRQ(ierr); 105c4762a1bSJed Brown ierr = PetscFree(adjncy);CHKERRQ(ierr); 106c4762a1bSJed Brown ierr = PetscFree3(xyz,part,tpwgts);CHKERRQ(ierr); 107c4762a1bSJed Brown ierr = PetscFree(sxyz);CHKERRQ(ierr); 108c4762a1bSJed Brown ierr = PetscFinalize(); 109c4762a1bSJed Brown return ierr; 110c4762a1bSJed Brown } 111c4762a1bSJed Brown 112c4762a1bSJed Brown /*TEST 113c4762a1bSJed Brown 114c4762a1bSJed Brown build: 115c4762a1bSJed Brown requires: parmetis 116c4762a1bSJed Brown 117c4762a1bSJed Brown test: 118c4762a1bSJed Brown nsize: 2 119*dfd57a17SPierre Jolivet requires: parmetis datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 120c4762a1bSJed Brown args: -prefix ${DATAFILESPATH}/parmetis-test/testnp2 121c4762a1bSJed Brown 122c4762a1bSJed Brown test: 123c4762a1bSJed Brown suffix: 2 124c4762a1bSJed Brown nsize: 4 125*dfd57a17SPierre Jolivet requires: parmetis datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 126c4762a1bSJed Brown args: -prefix ${DATAFILESPATH}/parmetis-test/testnp4 127c4762a1bSJed Brown 128c4762a1bSJed Brown TEST*/ 129