xref: /petsc/src/vec/is/sf/tests/ex8.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
125f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
135f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
1497929ea7SJunchao Zhang 
155f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&x));
165f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(x));
175f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(x,PETSC_DECIDE,N));
185f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetOwnershipRange(x,&low,&high));
195f80ce2aSJacob Faibussowitsch   CHKERRQ(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] */
265f80ce2aSJacob Faibussowitsch   for (i=low; i<high; i++) CHKERRQ(VecSetValue(x,i,(PetscScalar)i,INSERT_VALUES));
275f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(x));
285f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyEnd(x));
2997929ea7SJunchao Zhang 
305f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\nTesting VecScatterCreateToZero\n"));
315f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterCreateToZero(x,&vscat,&y));
325f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)y,"y"));
3397929ea7SJunchao Zhang 
3497929ea7SJunchao Zhang   /* Test PetscSFBcastAndOp with op = MPI_REPLACE, which does y = x on rank 0 */
355f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD));
365f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD));
375f80ce2aSJacob Faibussowitsch   if (rank == 0) CHKERRQ(VecView(y,PETSC_VIEWER_STDOUT_SELF));
3897929ea7SJunchao Zhang 
3997929ea7SJunchao Zhang   /* Test PetscSFBcastAndOp with op = MPI_SUM, which does y += x */
405f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(vscat,x,y,ADD_VALUES,SCATTER_FORWARD));
415f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(vscat,x,y,ADD_VALUES,SCATTER_FORWARD));
425f80ce2aSJacob Faibussowitsch   if (rank == 0) CHKERRQ(VecView(y,PETSC_VIEWER_STDOUT_SELF));
4397929ea7SJunchao Zhang 
4497929ea7SJunchao Zhang   /* Test PetscSFReduce with op = MPI_REPLACE, which does x = y */
455f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(vscat,y,x,INSERT_VALUES,SCATTER_REVERSE));
465f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(vscat,y,x,INSERT_VALUES,SCATTER_REVERSE));
475f80ce2aSJacob Faibussowitsch   CHKERRQ(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*/
505f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(vscat,y,x,ADD_VALUES,SCATTER_REVERSE_LOCAL));
515f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(vscat,y,x,ADD_VALUES,SCATTER_REVERSE_LOCAL));
525f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(x,PETSC_VIEWER_STDOUT_WORLD));
5397929ea7SJunchao Zhang 
545f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&y));
555f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterDestroy(&vscat));
5697929ea7SJunchao Zhang 
5797929ea7SJunchao Zhang   /*-------------------------------------*/
5897929ea7SJunchao Zhang   /*       VecScatterCreateToAll         */
5997929ea7SJunchao Zhang   /*-------------------------------------*/
605f80ce2aSJacob Faibussowitsch   for (i=low; i<high; i++) CHKERRQ(VecSetValue(x,i,(PetscScalar)i,INSERT_VALUES));
615f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(x));
625f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyEnd(x));
6397929ea7SJunchao Zhang 
645f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\nTesting VecScatterCreateToAll\n"));
6597929ea7SJunchao Zhang 
665f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterCreateToAll(x,&vscat,&y));
675f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)y,"y"));
6897929ea7SJunchao Zhang 
6997929ea7SJunchao Zhang   /* Test PetscSFBcastAndOp with op = MPI_REPLACE, which does y = x on all ranks */
705f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD));
715f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD));
725f80ce2aSJacob Faibussowitsch   if (rank == 0) CHKERRQ(VecView(y,PETSC_VIEWER_STDOUT_SELF));
7397929ea7SJunchao Zhang 
7497929ea7SJunchao Zhang   /* Test PetscSFBcastAndOp with op = MPI_SUM, which does y += x */
755f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(vscat,x,y,ADD_VALUES,SCATTER_FORWARD));
765f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(vscat,x,y,ADD_VALUES,SCATTER_FORWARD));
775f80ce2aSJacob Faibussowitsch   if (rank == 0) CHKERRQ(VecView(y,PETSC_VIEWER_STDOUT_SELF));
7897929ea7SJunchao Zhang 
7997929ea7SJunchao Zhang   /* Test PetscSFReduce with op = MPI_REPLACE, which does x = y */
805f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(vscat,y,x,INSERT_VALUES,SCATTER_REVERSE));
815f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(vscat,y,x,INSERT_VALUES,SCATTER_REVERSE));
825f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(x,PETSC_VIEWER_STDOUT_WORLD));
8397929ea7SJunchao Zhang 
8497929ea7SJunchao Zhang   /* Test PetscSFReduce with op = MPI_SUM, which does x += size*y */
855f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(vscat,y,x,ADD_VALUES,SCATTER_REVERSE));
865f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(vscat,y,x,ADD_VALUES,SCATTER_REVERSE));
875f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(x,PETSC_VIEWER_STDOUT_WORLD));
885f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
895f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&y));
905f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterDestroy(&vscat));
9197929ea7SJunchao Zhang 
92*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
93*b122ec5aSJacob 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