1*97929ea7SJunchao Zhang static char help[]= "Test VecScatterCreateToZero, VecScatterCreateToAll\n\n"; 2*97929ea7SJunchao Zhang 3*97929ea7SJunchao Zhang #include <petscvec.h> 4*97929ea7SJunchao Zhang int main(int argc,char **argv) 5*97929ea7SJunchao Zhang { 6*97929ea7SJunchao Zhang PetscErrorCode ierr; 7*97929ea7SJunchao Zhang PetscInt i,N=10,low,high; 8*97929ea7SJunchao Zhang PetscMPIInt size,rank; 9*97929ea7SJunchao Zhang Vec x,y; 10*97929ea7SJunchao Zhang VecScatter vscat; 11*97929ea7SJunchao Zhang 12*97929ea7SJunchao Zhang ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 13*97929ea7SJunchao Zhang ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 14*97929ea7SJunchao Zhang ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 15*97929ea7SJunchao Zhang 16*97929ea7SJunchao Zhang ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr); 17*97929ea7SJunchao Zhang ierr = VecSetFromOptions(x);CHKERRQ(ierr); 18*97929ea7SJunchao Zhang ierr = VecSetSizes(x,PETSC_DECIDE,N);CHKERRQ(ierr); 19*97929ea7SJunchao Zhang ierr = VecGetOwnershipRange(x,&low,&high);CHKERRQ(ierr); 20*97929ea7SJunchao Zhang ierr = PetscObjectSetName((PetscObject)x,"x");CHKERRQ(ierr); 21*97929ea7SJunchao Zhang 22*97929ea7SJunchao Zhang /*-------------------------------------*/ 23*97929ea7SJunchao Zhang /* VecScatterCreateToZero */ 24*97929ea7SJunchao Zhang /*-------------------------------------*/ 25*97929ea7SJunchao Zhang 26*97929ea7SJunchao Zhang /* MPI vec x = [0, 1, 2, .., N-1] */ 27*97929ea7SJunchao Zhang for (i=low; i<high; i++) {ierr = VecSetValue(x,i,(PetscScalar)i,INSERT_VALUES);CHKERRQ(ierr);} 28*97929ea7SJunchao Zhang ierr = VecAssemblyBegin(x);CHKERRQ(ierr); 29*97929ea7SJunchao Zhang ierr = VecAssemblyEnd(x);CHKERRQ(ierr); 30*97929ea7SJunchao Zhang 31*97929ea7SJunchao Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"\nTesting VecScatterCreateToZero\n");CHKERRQ(ierr); 32*97929ea7SJunchao Zhang ierr = VecScatterCreateToZero(x,&vscat,&y);CHKERRQ(ierr); 33*97929ea7SJunchao Zhang ierr = PetscObjectSetName((PetscObject)y,"y");CHKERRQ(ierr); 34*97929ea7SJunchao Zhang 35*97929ea7SJunchao Zhang /* Test PetscSFBcastAndOp with op = MPI_REPLACE, which does y = x on rank 0 */ 36*97929ea7SJunchao Zhang ierr = VecScatterBegin(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 37*97929ea7SJunchao Zhang ierr = VecScatterEnd(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 38*97929ea7SJunchao Zhang if (!rank) {ierr = VecView(y,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);} 39*97929ea7SJunchao Zhang 40*97929ea7SJunchao Zhang /* Test PetscSFBcastAndOp with op = MPI_SUM, which does y += x */ 41*97929ea7SJunchao Zhang ierr = VecScatterBegin(vscat,x,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 42*97929ea7SJunchao Zhang ierr = VecScatterEnd(vscat,x,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 43*97929ea7SJunchao Zhang if (!rank) {ierr = VecView(y,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);} 44*97929ea7SJunchao Zhang 45*97929ea7SJunchao Zhang /* Test PetscSFReduce with op = MPI_REPLACE, which does x = y */ 46*97929ea7SJunchao Zhang ierr = VecScatterBegin(vscat,y,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 47*97929ea7SJunchao Zhang ierr = VecScatterEnd(vscat,y,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 48*97929ea7SJunchao Zhang ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 49*97929ea7SJunchao Zhang 50*97929ea7SJunchao Zhang /* Test PetscSFReduce with op = MPI_SUM, which does x += y on x's local part on rank 0*/ 51*97929ea7SJunchao Zhang ierr = VecScatterBegin(vscat,y,x,ADD_VALUES,SCATTER_REVERSE_LOCAL);CHKERRQ(ierr); 52*97929ea7SJunchao Zhang ierr = VecScatterEnd(vscat,y,x,ADD_VALUES,SCATTER_REVERSE_LOCAL);CHKERRQ(ierr); 53*97929ea7SJunchao Zhang ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 54*97929ea7SJunchao Zhang 55*97929ea7SJunchao Zhang ierr = VecDestroy(&y);CHKERRQ(ierr); 56*97929ea7SJunchao Zhang ierr = VecScatterDestroy(&vscat);CHKERRQ(ierr); 57*97929ea7SJunchao Zhang 58*97929ea7SJunchao Zhang /*-------------------------------------*/ 59*97929ea7SJunchao Zhang /* VecScatterCreateToAll */ 60*97929ea7SJunchao Zhang /*-------------------------------------*/ 61*97929ea7SJunchao Zhang for (i=low; i<high; i++) {ierr = VecSetValue(x,i,(PetscScalar)i,INSERT_VALUES);CHKERRQ(ierr);} 62*97929ea7SJunchao Zhang ierr = VecAssemblyBegin(x);CHKERRQ(ierr); 63*97929ea7SJunchao Zhang ierr = VecAssemblyEnd(x);CHKERRQ(ierr); 64*97929ea7SJunchao Zhang 65*97929ea7SJunchao Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"\nTesting VecScatterCreateToAll\n");CHKERRQ(ierr); 66*97929ea7SJunchao Zhang 67*97929ea7SJunchao Zhang ierr = VecScatterCreateToAll(x,&vscat,&y);CHKERRQ(ierr); 68*97929ea7SJunchao Zhang ierr = PetscObjectSetName((PetscObject)y,"y");CHKERRQ(ierr); 69*97929ea7SJunchao Zhang 70*97929ea7SJunchao Zhang /* Test PetscSFBcastAndOp with op = MPI_REPLACE, which does y = x on all ranks */ 71*97929ea7SJunchao Zhang ierr = VecScatterBegin(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 72*97929ea7SJunchao Zhang ierr = VecScatterEnd(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 73*97929ea7SJunchao Zhang if (!rank) {ierr = VecView(y,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);} 74*97929ea7SJunchao Zhang 75*97929ea7SJunchao Zhang /* Test PetscSFBcastAndOp with op = MPI_SUM, which does y += x */ 76*97929ea7SJunchao Zhang ierr = VecScatterBegin(vscat,x,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 77*97929ea7SJunchao Zhang ierr = VecScatterEnd(vscat,x,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 78*97929ea7SJunchao Zhang if (!rank) {ierr = VecView(y,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);} 79*97929ea7SJunchao Zhang 80*97929ea7SJunchao Zhang /* Test PetscSFReduce with op = MPI_REPLACE, which does x = y */ 81*97929ea7SJunchao Zhang ierr = VecScatterBegin(vscat,y,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 82*97929ea7SJunchao Zhang ierr = VecScatterEnd(vscat,y,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 83*97929ea7SJunchao Zhang ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 84*97929ea7SJunchao Zhang 85*97929ea7SJunchao Zhang /* Test PetscSFReduce with op = MPI_SUM, which does x += size*y */ 86*97929ea7SJunchao Zhang ierr = VecScatterBegin(vscat,y,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 87*97929ea7SJunchao Zhang ierr = VecScatterEnd(vscat,y,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 88*97929ea7SJunchao Zhang ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 89*97929ea7SJunchao Zhang ierr = VecDestroy(&x);CHKERRQ(ierr); 90*97929ea7SJunchao Zhang ierr = VecDestroy(&y);CHKERRQ(ierr); 91*97929ea7SJunchao Zhang ierr = VecScatterDestroy(&vscat);CHKERRQ(ierr); 92*97929ea7SJunchao Zhang 93*97929ea7SJunchao Zhang ierr = PetscFinalize(); 94*97929ea7SJunchao Zhang return ierr; 95*97929ea7SJunchao Zhang } 96*97929ea7SJunchao Zhang 97*97929ea7SJunchao Zhang /*TEST 98*97929ea7SJunchao Zhang 99*97929ea7SJunchao Zhang testset: 100*97929ea7SJunchao Zhang # N=10 is divisible by nsize, to trigger Allgather/Gather in SF 101*97929ea7SJunchao Zhang nsize: 2 102*97929ea7SJunchao Zhang filter: grep -v "type" 103*97929ea7SJunchao Zhang output_file: output/ex8_1.out 104*97929ea7SJunchao Zhang 105*97929ea7SJunchao Zhang test: 106*97929ea7SJunchao Zhang suffix: 1_standard 107*97929ea7SJunchao Zhang 108*97929ea7SJunchao Zhang test: 109*97929ea7SJunchao Zhang suffix: 1_cuda 110*97929ea7SJunchao Zhang args: -vec_type cuda -vecscatter_packongpu true 111*97929ea7SJunchao Zhang requires: cuda 112*97929ea7SJunchao Zhang 113*97929ea7SJunchao Zhang test: 114*97929ea7SJunchao Zhang suffix: 1_cuda_aware_mpi 115*97929ea7SJunchao Zhang args: -vec_type cuda -vecscatter_packongpu false 116*97929ea7SJunchao Zhang requires: cuda define(PETSC_HAVE_MPI_GPU_AWARE) 117*97929ea7SJunchao Zhang 118*97929ea7SJunchao Zhang testset: 119*97929ea7SJunchao Zhang # N=10 is not divisible by nsize, to trigger Allgatherv/Gatherv in SF 120*97929ea7SJunchao Zhang suffix: 2 121*97929ea7SJunchao Zhang nsize: 3 122*97929ea7SJunchao Zhang filter: grep -v "type" 123*97929ea7SJunchao Zhang output_file: output/ex8_2.out 124*97929ea7SJunchao Zhang 125*97929ea7SJunchao Zhang test: 126*97929ea7SJunchao Zhang suffix: 2_standard 127*97929ea7SJunchao Zhang 128*97929ea7SJunchao Zhang test: 129*97929ea7SJunchao Zhang suffix: 2_cuda 130*97929ea7SJunchao Zhang args: -vec_type cuda -vecscatter_packongpu true 131*97929ea7SJunchao Zhang requires: cuda 132*97929ea7SJunchao Zhang 133*97929ea7SJunchao Zhang test: 134*97929ea7SJunchao Zhang suffix: 2_cuda_aware_mpi 135*97929ea7SJunchao Zhang args: -vec_type cuda -vecscatter_packongpu false 136*97929ea7SJunchao Zhang requires: cuda define(PETSC_HAVE_MPI_GPU_AWARE) 137*97929ea7SJunchao Zhang 138*97929ea7SJunchao Zhang TEST*/ 139*97929ea7SJunchao Zhang 140