181b34b04SJunchao Zhang static const char help[] = "Test PetscSF with MPI large count (more than 2 billion elements in messages)\n\n"; 281b34b04SJunchao Zhang 381b34b04SJunchao Zhang #include <petscsys.h> 481b34b04SJunchao Zhang #include <petscsf.h> 581b34b04SJunchao Zhang 681b34b04SJunchao Zhang int main(int argc,char **argv) 781b34b04SJunchao Zhang { 881b34b04SJunchao Zhang PetscErrorCode ierr; 981b34b04SJunchao Zhang PetscSF sf; 1081b34b04SJunchao Zhang PetscInt i,nroots,nleaves; 1181b34b04SJunchao Zhang PetscInt n = (1ULL<<31) + 1024; /* a little over 2G elements */ 1281b34b04SJunchao Zhang PetscSFNode *iremote = NULL; 1381b34b04SJunchao Zhang PetscMPIInt rank,size; 1481b34b04SJunchao Zhang char *rootdata=NULL,*leafdata=NULL; 1581b34b04SJunchao Zhang Vec x,y; 1681b34b04SJunchao Zhang VecScatter vscat; 1781b34b04SJunchao Zhang PetscInt rstart,rend; 1881b34b04SJunchao Zhang IS ix; 1981b34b04SJunchao Zhang const PetscScalar *xv; 2081b34b04SJunchao Zhang 2181b34b04SJunchao Zhang ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 2281b34b04SJunchao Zhang ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 2381b34b04SJunchao Zhang ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 24*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(size != 2,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"The test can only run with two MPI ranks"); 2581b34b04SJunchao Zhang 2681b34b04SJunchao Zhang /* Test PetscSF */ 2781b34b04SJunchao Zhang ierr = PetscSFCreate(PETSC_COMM_WORLD,&sf);CHKERRQ(ierr); 2881b34b04SJunchao Zhang ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr); 2981b34b04SJunchao Zhang 3081b34b04SJunchao Zhang if (!rank) { 3181b34b04SJunchao Zhang nroots = n; 3281b34b04SJunchao Zhang nleaves = 0; 3381b34b04SJunchao Zhang } else { 3481b34b04SJunchao Zhang nroots = 0; 3581b34b04SJunchao Zhang nleaves = n; 3681b34b04SJunchao Zhang ierr = PetscMalloc1(nleaves,&iremote);CHKERRQ(ierr); 3781b34b04SJunchao Zhang for (i=0; i<nleaves; i++) { 3881b34b04SJunchao Zhang iremote[i].rank = 0; 3981b34b04SJunchao Zhang iremote[i].index = i; 4081b34b04SJunchao Zhang } 4181b34b04SJunchao Zhang } 4281b34b04SJunchao Zhang ierr = PetscSFSetGraph(sf,nroots,nleaves,NULL,PETSC_COPY_VALUES,iremote,PETSC_COPY_VALUES);CHKERRQ(ierr); 4381b34b04SJunchao Zhang ierr = PetscMalloc2(nroots,&rootdata,nleaves,&leafdata);CHKERRQ(ierr); 4481b34b04SJunchao Zhang if (!rank) { 4581b34b04SJunchao Zhang memset(rootdata,11,nroots); 4681b34b04SJunchao Zhang rootdata[nroots-1] = 12; /* Use a different value at the end */ 4781b34b04SJunchao Zhang } 4881b34b04SJunchao Zhang 4981b34b04SJunchao Zhang ierr = PetscSFBcastBegin(sf,MPI_SIGNED_CHAR,rootdata,leafdata,MPI_REPLACE);CHKERRQ(ierr); /* rank 0->1, bcast rootdata to leafdata */ 5081b34b04SJunchao Zhang ierr = PetscSFBcastEnd(sf,MPI_SIGNED_CHAR,rootdata,leafdata,MPI_REPLACE);CHKERRQ(ierr); 5181b34b04SJunchao Zhang ierr = PetscSFReduceBegin(sf,MPI_SIGNED_CHAR,leafdata,rootdata,MPI_SUM);CHKERRQ(ierr); /* rank 1->0, add leafdata to rootdata */ 5281b34b04SJunchao Zhang ierr = PetscSFReduceEnd(sf,MPI_SIGNED_CHAR,leafdata,rootdata,MPI_SUM);CHKERRQ(ierr); 5381b34b04SJunchao Zhang if (!rank) { 54*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(rootdata[0] != 22 || rootdata[nroots-1] != 24,PETSC_COMM_SELF,PETSC_ERR_PLIB,"SF: wrong results"); 5581b34b04SJunchao Zhang } 5681b34b04SJunchao Zhang 5781b34b04SJunchao Zhang ierr = PetscFree2(rootdata,leafdata);CHKERRQ(ierr); 5881b34b04SJunchao Zhang ierr = PetscFree(iremote);CHKERRQ(ierr); 5981b34b04SJunchao Zhang ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 6081b34b04SJunchao Zhang 6181b34b04SJunchao Zhang /* Test VecScatter */ 6281b34b04SJunchao Zhang ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr); 6381b34b04SJunchao Zhang ierr = VecCreate(PETSC_COMM_WORLD,&y);CHKERRQ(ierr); 6481b34b04SJunchao Zhang ierr = VecSetSizes(x,rank==0? n : 64,PETSC_DECIDE);CHKERRQ(ierr); 6581b34b04SJunchao Zhang ierr = VecSetSizes(y,rank==0? 64 : n,PETSC_DECIDE);CHKERRQ(ierr); 6681b34b04SJunchao Zhang ierr = VecSetFromOptions(x);CHKERRQ(ierr); 6781b34b04SJunchao Zhang ierr = VecSetFromOptions(y);CHKERRQ(ierr); 6881b34b04SJunchao Zhang 6981b34b04SJunchao Zhang ierr = VecGetOwnershipRange(x,&rstart,&rend);CHKERRQ(ierr); 7081b34b04SJunchao Zhang ierr = ISCreateStride(PETSC_COMM_SELF,rend-rstart,rstart,1,&ix);CHKERRQ(ierr); 7181b34b04SJunchao Zhang ierr = VecScatterCreate(x,ix,y,ix,&vscat);CHKERRQ(ierr); 7281b34b04SJunchao Zhang 7381b34b04SJunchao Zhang ierr = VecSet(x,3.0);CHKERRQ(ierr); 7481b34b04SJunchao Zhang ierr = VecScatterBegin(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7581b34b04SJunchao Zhang ierr = VecScatterEnd(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7681b34b04SJunchao Zhang 7781b34b04SJunchao Zhang ierr = VecScatterBegin(vscat,y,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 7881b34b04SJunchao Zhang ierr = VecScatterEnd(vscat,y,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 7981b34b04SJunchao Zhang 8081b34b04SJunchao Zhang ierr = VecGetArrayRead(x,&xv);CHKERRQ(ierr); 81*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(xv[0] != 6.0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"VecScatter: wrong results"); 8281b34b04SJunchao Zhang ierr = VecRestoreArrayRead(x,&xv);CHKERRQ(ierr); 8381b34b04SJunchao Zhang ierr = VecDestroy(&x);CHKERRQ(ierr); 8481b34b04SJunchao Zhang ierr = VecDestroy(&y);CHKERRQ(ierr); 8581b34b04SJunchao Zhang ierr = VecScatterDestroy(&vscat);CHKERRQ(ierr); 8681b34b04SJunchao Zhang ierr = ISDestroy(&ix);CHKERRQ(ierr); 8781b34b04SJunchao Zhang 8881b34b04SJunchao Zhang ierr = PetscFinalize(); 8981b34b04SJunchao Zhang return ierr; 9081b34b04SJunchao Zhang } 9181b34b04SJunchao Zhang 9281b34b04SJunchao Zhang /**TEST 9381b34b04SJunchao Zhang test: 9481b34b04SJunchao Zhang requires: defined(PETSC_HAVE_MPI_LARGE_COUNT) defined(PETSC_USE_64BIT_INDICES) 9581b34b04SJunchao Zhang TODO: need a machine with big memory (~150GB) to run the test 9681b34b04SJunchao Zhang nsize: 2 9781b34b04SJunchao Zhang args: -sf_type {{basic neighbor}} 9881b34b04SJunchao Zhang 9981b34b04SJunchao Zhang TEST**/ 10081b34b04SJunchao Zhang 101