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*327415f7SBarry Smith PetscFunctionBeginUser; 129566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 139566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 149566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 1597929ea7SJunchao Zhang 169566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&x)); 179566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(x)); 189566063dSJacob Faibussowitsch PetscCall(VecSetSizes(x,PETSC_DECIDE,N)); 199566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(x,&low,&high)); 209566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)x,"x")); 2197929ea7SJunchao Zhang 2297929ea7SJunchao Zhang /*-------------------------------------*/ 2397929ea7SJunchao Zhang /* VecScatterCreateToZero */ 2497929ea7SJunchao Zhang /*-------------------------------------*/ 2597929ea7SJunchao Zhang 2697929ea7SJunchao Zhang /* MPI vec x = [0, 1, 2, .., N-1] */ 279566063dSJacob Faibussowitsch for (i=low; i<high; i++) PetscCall(VecSetValue(x,i,(PetscScalar)i,INSERT_VALUES)); 289566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(x)); 299566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(x)); 3097929ea7SJunchao Zhang 319566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nTesting VecScatterCreateToZero\n")); 329566063dSJacob Faibussowitsch PetscCall(VecScatterCreateToZero(x,&vscat,&y)); 339566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)y,"y")); 3497929ea7SJunchao Zhang 3597929ea7SJunchao Zhang /* Test PetscSFBcastAndOp with op = MPI_REPLACE, which does y = x on rank 0 */ 369566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD)); 379566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD)); 389566063dSJacob Faibussowitsch if (rank == 0) PetscCall(VecView(y,PETSC_VIEWER_STDOUT_SELF)); 3997929ea7SJunchao Zhang 4097929ea7SJunchao Zhang /* Test PetscSFBcastAndOp with op = MPI_SUM, which does y += x */ 419566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(vscat,x,y,ADD_VALUES,SCATTER_FORWARD)); 429566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(vscat,x,y,ADD_VALUES,SCATTER_FORWARD)); 439566063dSJacob Faibussowitsch if (rank == 0) PetscCall(VecView(y,PETSC_VIEWER_STDOUT_SELF)); 4497929ea7SJunchao Zhang 4597929ea7SJunchao Zhang /* Test PetscSFReduce with op = MPI_REPLACE, which does x = y */ 469566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(vscat,y,x,INSERT_VALUES,SCATTER_REVERSE)); 479566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(vscat,y,x,INSERT_VALUES,SCATTER_REVERSE)); 489566063dSJacob Faibussowitsch PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 4997929ea7SJunchao Zhang 5097929ea7SJunchao Zhang /* Test PetscSFReduce with op = MPI_SUM, which does x += y on x's local part on rank 0*/ 519566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(vscat,y,x,ADD_VALUES,SCATTER_REVERSE_LOCAL)); 529566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(vscat,y,x,ADD_VALUES,SCATTER_REVERSE_LOCAL)); 539566063dSJacob Faibussowitsch PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 5497929ea7SJunchao Zhang 559566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 569566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&vscat)); 5797929ea7SJunchao Zhang 5897929ea7SJunchao Zhang /*-------------------------------------*/ 5997929ea7SJunchao Zhang /* VecScatterCreateToAll */ 6097929ea7SJunchao Zhang /*-------------------------------------*/ 619566063dSJacob Faibussowitsch for (i=low; i<high; i++) PetscCall(VecSetValue(x,i,(PetscScalar)i,INSERT_VALUES)); 629566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(x)); 639566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(x)); 6497929ea7SJunchao Zhang 659566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nTesting VecScatterCreateToAll\n")); 6697929ea7SJunchao Zhang 679566063dSJacob Faibussowitsch PetscCall(VecScatterCreateToAll(x,&vscat,&y)); 689566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)y,"y")); 6997929ea7SJunchao Zhang 7097929ea7SJunchao Zhang /* Test PetscSFBcastAndOp with op = MPI_REPLACE, which does y = x on all ranks */ 719566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD)); 729566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD)); 739566063dSJacob Faibussowitsch if (rank == 0) PetscCall(VecView(y,PETSC_VIEWER_STDOUT_SELF)); 7497929ea7SJunchao Zhang 7597929ea7SJunchao Zhang /* Test PetscSFBcastAndOp with op = MPI_SUM, which does y += x */ 769566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(vscat,x,y,ADD_VALUES,SCATTER_FORWARD)); 779566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(vscat,x,y,ADD_VALUES,SCATTER_FORWARD)); 789566063dSJacob Faibussowitsch if (rank == 0) PetscCall(VecView(y,PETSC_VIEWER_STDOUT_SELF)); 7997929ea7SJunchao Zhang 8097929ea7SJunchao Zhang /* Test PetscSFReduce with op = MPI_REPLACE, which does x = y */ 819566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(vscat,y,x,INSERT_VALUES,SCATTER_REVERSE)); 829566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(vscat,y,x,INSERT_VALUES,SCATTER_REVERSE)); 839566063dSJacob Faibussowitsch PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 8497929ea7SJunchao Zhang 8597929ea7SJunchao Zhang /* Test PetscSFReduce with op = MPI_SUM, which does x += size*y */ 869566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(vscat,y,x,ADD_VALUES,SCATTER_REVERSE)); 879566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(vscat,y,x,ADD_VALUES,SCATTER_REVERSE)); 889566063dSJacob Faibussowitsch PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 899566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 909566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 919566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&vscat)); 9297929ea7SJunchao Zhang 939566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 94b122ec5aSJacob Faibussowitsch return 0; 9597929ea7SJunchao Zhang } 9697929ea7SJunchao Zhang 9797929ea7SJunchao Zhang /*TEST 9897929ea7SJunchao Zhang 9997929ea7SJunchao Zhang testset: 10097929ea7SJunchao Zhang # N=10 is divisible by nsize, to trigger Allgather/Gather in SF 10197929ea7SJunchao Zhang nsize: 2 10226e8e884SScott Kruger # Exact numbers really matter here 10326e8e884SScott Kruger diff_args: -j 10497929ea7SJunchao Zhang filter: grep -v "type" 10597929ea7SJunchao Zhang output_file: output/ex8_1.out 10697929ea7SJunchao Zhang 10797929ea7SJunchao Zhang test: 10897929ea7SJunchao Zhang suffix: 1_standard 10997929ea7SJunchao Zhang 11097929ea7SJunchao Zhang test: 11197929ea7SJunchao Zhang suffix: 1_cuda 11226e8e884SScott Kruger # sf_backend cuda is not needed if compiling only with cuda 113bd46da1dSJunchao Zhang args: -vec_type cuda -sf_backend cuda 11497929ea7SJunchao Zhang requires: cuda 11597929ea7SJunchao Zhang 11697929ea7SJunchao Zhang test: 11726e8e884SScott Kruger suffix: 1_hip 118bd46da1dSJunchao Zhang args: -vec_type hip -sf_backend hip 11926e8e884SScott Kruger requires: hip 12026e8e884SScott Kruger 12126e8e884SScott Kruger test: 12297929ea7SJunchao Zhang suffix: 1_cuda_aware_mpi 12326e8e884SScott Kruger # sf_backend cuda is not needed if compiling only with cuda 124bd46da1dSJunchao Zhang args: -vec_type cuda -sf_backend cuda 125dfd57a17SPierre Jolivet requires: cuda defined(PETSC_HAVE_MPI_GPU_AWARE) 12697929ea7SJunchao Zhang 12797929ea7SJunchao Zhang testset: 12897929ea7SJunchao Zhang # N=10 is not divisible by nsize, to trigger Allgatherv/Gatherv in SF 12997929ea7SJunchao Zhang nsize: 3 13026e8e884SScott Kruger # Exact numbers really matter here 13126e8e884SScott Kruger diff_args: -j 13297929ea7SJunchao Zhang filter: grep -v "type" 13397929ea7SJunchao Zhang output_file: output/ex8_2.out 13497929ea7SJunchao Zhang 13597929ea7SJunchao Zhang test: 13697929ea7SJunchao Zhang suffix: 2_standard 13797929ea7SJunchao Zhang 13897929ea7SJunchao Zhang test: 13997929ea7SJunchao Zhang suffix: 2_cuda 14026e8e884SScott Kruger # sf_backend cuda is not needed if compiling only with cuda 141bd46da1dSJunchao Zhang args: -vec_type cuda -sf_backend cuda 14297929ea7SJunchao Zhang requires: cuda 14397929ea7SJunchao Zhang 14497929ea7SJunchao Zhang test: 14526e8e884SScott Kruger suffix: 2_hip 14626e8e884SScott Kruger # sf_backend hip is not needed if compiling only with hip 147bd46da1dSJunchao Zhang args: -vec_type hip -sf_backend hip 14826e8e884SScott Kruger requires: hip 14926e8e884SScott Kruger 15026e8e884SScott Kruger test: 15197929ea7SJunchao Zhang suffix: 2_cuda_aware_mpi 152bd46da1dSJunchao Zhang args: -vec_type cuda 153dfd57a17SPierre Jolivet requires: cuda defined(PETSC_HAVE_MPI_GPU_AWARE) 15497929ea7SJunchao Zhang 15597929ea7SJunchao Zhang TEST*/ 156