197929ea7SJunchao Zhang static char help[]= "Test VecScatterCreateToZero, VecScatterCreateToAll\n\n"; 297929ea7SJunchao Zhang 397929ea7SJunchao Zhang #include <petscvec.h> 497929ea7SJunchao Zhang int main(int argc,char **argv) 597929ea7SJunchao Zhang { 697929ea7SJunchao Zhang PetscInt i,N=10,low,high; 797929ea7SJunchao Zhang PetscMPIInt size,rank; 897929ea7SJunchao Zhang Vec x,y; 997929ea7SJunchao Zhang VecScatter vscat; 1097929ea7SJunchao Zhang 11*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 12*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 13*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 1497929ea7SJunchao Zhang 15*9566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&x)); 16*9566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(x)); 17*9566063dSJacob Faibussowitsch PetscCall(VecSetSizes(x,PETSC_DECIDE,N)); 18*9566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(x,&low,&high)); 19*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)x,"x")); 2097929ea7SJunchao Zhang 2197929ea7SJunchao Zhang /*-------------------------------------*/ 2297929ea7SJunchao Zhang /* VecScatterCreateToZero */ 2397929ea7SJunchao Zhang /*-------------------------------------*/ 2497929ea7SJunchao Zhang 2597929ea7SJunchao Zhang /* MPI vec x = [0, 1, 2, .., N-1] */ 26*9566063dSJacob Faibussowitsch for (i=low; i<high; i++) PetscCall(VecSetValue(x,i,(PetscScalar)i,INSERT_VALUES)); 27*9566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(x)); 28*9566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(x)); 2997929ea7SJunchao Zhang 30*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nTesting VecScatterCreateToZero\n")); 31*9566063dSJacob Faibussowitsch PetscCall(VecScatterCreateToZero(x,&vscat,&y)); 32*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)y,"y")); 3397929ea7SJunchao Zhang 3497929ea7SJunchao Zhang /* Test PetscSFBcastAndOp with op = MPI_REPLACE, which does y = x on rank 0 */ 35*9566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD)); 36*9566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD)); 37*9566063dSJacob Faibussowitsch if (rank == 0) PetscCall(VecView(y,PETSC_VIEWER_STDOUT_SELF)); 3897929ea7SJunchao Zhang 3997929ea7SJunchao Zhang /* Test PetscSFBcastAndOp with op = MPI_SUM, which does y += x */ 40*9566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(vscat,x,y,ADD_VALUES,SCATTER_FORWARD)); 41*9566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(vscat,x,y,ADD_VALUES,SCATTER_FORWARD)); 42*9566063dSJacob Faibussowitsch if (rank == 0) PetscCall(VecView(y,PETSC_VIEWER_STDOUT_SELF)); 4397929ea7SJunchao Zhang 4497929ea7SJunchao Zhang /* Test PetscSFReduce with op = MPI_REPLACE, which does x = y */ 45*9566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(vscat,y,x,INSERT_VALUES,SCATTER_REVERSE)); 46*9566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(vscat,y,x,INSERT_VALUES,SCATTER_REVERSE)); 47*9566063dSJacob Faibussowitsch PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 4897929ea7SJunchao Zhang 4997929ea7SJunchao Zhang /* Test PetscSFReduce with op = MPI_SUM, which does x += y on x's local part on rank 0*/ 50*9566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(vscat,y,x,ADD_VALUES,SCATTER_REVERSE_LOCAL)); 51*9566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(vscat,y,x,ADD_VALUES,SCATTER_REVERSE_LOCAL)); 52*9566063dSJacob Faibussowitsch PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 5397929ea7SJunchao Zhang 54*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 55*9566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&vscat)); 5697929ea7SJunchao Zhang 5797929ea7SJunchao Zhang /*-------------------------------------*/ 5897929ea7SJunchao Zhang /* VecScatterCreateToAll */ 5997929ea7SJunchao Zhang /*-------------------------------------*/ 60*9566063dSJacob Faibussowitsch for (i=low; i<high; i++) PetscCall(VecSetValue(x,i,(PetscScalar)i,INSERT_VALUES)); 61*9566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(x)); 62*9566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(x)); 6397929ea7SJunchao Zhang 64*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nTesting VecScatterCreateToAll\n")); 6597929ea7SJunchao Zhang 66*9566063dSJacob Faibussowitsch PetscCall(VecScatterCreateToAll(x,&vscat,&y)); 67*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)y,"y")); 6897929ea7SJunchao Zhang 6997929ea7SJunchao Zhang /* Test PetscSFBcastAndOp with op = MPI_REPLACE, which does y = x on all ranks */ 70*9566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD)); 71*9566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD)); 72*9566063dSJacob Faibussowitsch if (rank == 0) PetscCall(VecView(y,PETSC_VIEWER_STDOUT_SELF)); 7397929ea7SJunchao Zhang 7497929ea7SJunchao Zhang /* Test PetscSFBcastAndOp with op = MPI_SUM, which does y += x */ 75*9566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(vscat,x,y,ADD_VALUES,SCATTER_FORWARD)); 76*9566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(vscat,x,y,ADD_VALUES,SCATTER_FORWARD)); 77*9566063dSJacob Faibussowitsch if (rank == 0) PetscCall(VecView(y,PETSC_VIEWER_STDOUT_SELF)); 7897929ea7SJunchao Zhang 7997929ea7SJunchao Zhang /* Test PetscSFReduce with op = MPI_REPLACE, which does x = y */ 80*9566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(vscat,y,x,INSERT_VALUES,SCATTER_REVERSE)); 81*9566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(vscat,y,x,INSERT_VALUES,SCATTER_REVERSE)); 82*9566063dSJacob Faibussowitsch PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 8397929ea7SJunchao Zhang 8497929ea7SJunchao Zhang /* Test PetscSFReduce with op = MPI_SUM, which does x += size*y */ 85*9566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(vscat,y,x,ADD_VALUES,SCATTER_REVERSE)); 86*9566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(vscat,y,x,ADD_VALUES,SCATTER_REVERSE)); 87*9566063dSJacob Faibussowitsch PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 88*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 89*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 90*9566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&vscat)); 9197929ea7SJunchao Zhang 92*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 93b122ec5aSJacob Faibussowitsch return 0; 9497929ea7SJunchao Zhang } 9597929ea7SJunchao Zhang 9697929ea7SJunchao Zhang /*TEST 9797929ea7SJunchao Zhang 9897929ea7SJunchao Zhang testset: 9997929ea7SJunchao Zhang # N=10 is divisible by nsize, to trigger Allgather/Gather in SF 10097929ea7SJunchao Zhang nsize: 2 10126e8e884SScott Kruger # Exact numbers really matter here 10226e8e884SScott Kruger diff_args: -j 10397929ea7SJunchao Zhang filter: grep -v "type" 10497929ea7SJunchao Zhang output_file: output/ex8_1.out 10597929ea7SJunchao Zhang 10697929ea7SJunchao Zhang test: 10797929ea7SJunchao Zhang suffix: 1_standard 10897929ea7SJunchao Zhang 10997929ea7SJunchao Zhang test: 11097929ea7SJunchao Zhang suffix: 1_cuda 11126e8e884SScott Kruger # sf_backend cuda is not needed if compiling only with cuda 112bd46da1dSJunchao Zhang args: -vec_type cuda -sf_backend cuda 11397929ea7SJunchao Zhang requires: cuda 11497929ea7SJunchao Zhang 11597929ea7SJunchao Zhang test: 11626e8e884SScott Kruger suffix: 1_hip 117bd46da1dSJunchao Zhang args: -vec_type hip -sf_backend hip 11826e8e884SScott Kruger requires: hip 11926e8e884SScott Kruger 12026e8e884SScott Kruger test: 12197929ea7SJunchao Zhang suffix: 1_cuda_aware_mpi 12226e8e884SScott Kruger # sf_backend cuda is not needed if compiling only with cuda 123bd46da1dSJunchao Zhang args: -vec_type cuda -sf_backend cuda 124dfd57a17SPierre Jolivet requires: cuda defined(PETSC_HAVE_MPI_GPU_AWARE) 12597929ea7SJunchao Zhang 12697929ea7SJunchao Zhang testset: 12797929ea7SJunchao Zhang # N=10 is not divisible by nsize, to trigger Allgatherv/Gatherv in SF 12897929ea7SJunchao Zhang nsize: 3 12926e8e884SScott Kruger # Exact numbers really matter here 13026e8e884SScott Kruger diff_args: -j 13197929ea7SJunchao Zhang filter: grep -v "type" 13297929ea7SJunchao Zhang output_file: output/ex8_2.out 13397929ea7SJunchao Zhang 13497929ea7SJunchao Zhang test: 13597929ea7SJunchao Zhang suffix: 2_standard 13697929ea7SJunchao Zhang 13797929ea7SJunchao Zhang test: 13897929ea7SJunchao Zhang suffix: 2_cuda 13926e8e884SScott Kruger # sf_backend cuda is not needed if compiling only with cuda 140bd46da1dSJunchao Zhang args: -vec_type cuda -sf_backend cuda 14197929ea7SJunchao Zhang requires: cuda 14297929ea7SJunchao Zhang 14397929ea7SJunchao Zhang test: 14426e8e884SScott Kruger suffix: 2_hip 14526e8e884SScott Kruger # sf_backend hip is not needed if compiling only with hip 146bd46da1dSJunchao Zhang args: -vec_type hip -sf_backend hip 14726e8e884SScott Kruger requires: hip 14826e8e884SScott Kruger 14926e8e884SScott Kruger test: 15097929ea7SJunchao Zhang suffix: 2_cuda_aware_mpi 151bd46da1dSJunchao Zhang args: -vec_type cuda 152dfd57a17SPierre Jolivet requires: cuda defined(PETSC_HAVE_MPI_GPU_AWARE) 15397929ea7SJunchao Zhang 15497929ea7SJunchao Zhang TEST*/ 155