xref: /petsc/src/vec/is/sf/tests/ex8.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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   PetscErrorCode     ierr;
797929ea7SJunchao Zhang   PetscInt           i,N=10,low,high;
897929ea7SJunchao Zhang   PetscMPIInt        size,rank;
997929ea7SJunchao Zhang   Vec                x,y;
1097929ea7SJunchao Zhang   VecScatter         vscat;
1197929ea7SJunchao Zhang 
1297929ea7SJunchao Zhang   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
13*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
14*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
1597929ea7SJunchao Zhang 
16*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&x));
17*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(x));
18*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(x,PETSC_DECIDE,N));
19*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetOwnershipRange(x,&low,&high));
20*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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] */
27*5f80ce2aSJacob Faibussowitsch   for (i=low; i<high; i++) CHKERRQ(VecSetValue(x,i,(PetscScalar)i,INSERT_VALUES));
28*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(x));
29*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyEnd(x));
3097929ea7SJunchao Zhang 
31*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\nTesting VecScatterCreateToZero\n"));
32*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterCreateToZero(x,&vscat,&y));
33*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)y,"y"));
3497929ea7SJunchao Zhang 
3597929ea7SJunchao Zhang   /* Test PetscSFBcastAndOp with op = MPI_REPLACE, which does y = x on rank 0 */
36*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD));
37*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD));
38*5f80ce2aSJacob Faibussowitsch   if (rank == 0) CHKERRQ(VecView(y,PETSC_VIEWER_STDOUT_SELF));
3997929ea7SJunchao Zhang 
4097929ea7SJunchao Zhang   /* Test PetscSFBcastAndOp with op = MPI_SUM, which does y += x */
41*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(vscat,x,y,ADD_VALUES,SCATTER_FORWARD));
42*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(vscat,x,y,ADD_VALUES,SCATTER_FORWARD));
43*5f80ce2aSJacob Faibussowitsch   if (rank == 0) CHKERRQ(VecView(y,PETSC_VIEWER_STDOUT_SELF));
4497929ea7SJunchao Zhang 
4597929ea7SJunchao Zhang   /* Test PetscSFReduce with op = MPI_REPLACE, which does x = y */
46*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(vscat,y,x,INSERT_VALUES,SCATTER_REVERSE));
47*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(vscat,y,x,INSERT_VALUES,SCATTER_REVERSE));
48*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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*/
51*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(vscat,y,x,ADD_VALUES,SCATTER_REVERSE_LOCAL));
52*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(vscat,y,x,ADD_VALUES,SCATTER_REVERSE_LOCAL));
53*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(x,PETSC_VIEWER_STDOUT_WORLD));
5497929ea7SJunchao Zhang 
55*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&y));
56*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterDestroy(&vscat));
5797929ea7SJunchao Zhang 
5897929ea7SJunchao Zhang   /*-------------------------------------*/
5997929ea7SJunchao Zhang   /*       VecScatterCreateToAll         */
6097929ea7SJunchao Zhang   /*-------------------------------------*/
61*5f80ce2aSJacob Faibussowitsch   for (i=low; i<high; i++) CHKERRQ(VecSetValue(x,i,(PetscScalar)i,INSERT_VALUES));
62*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(x));
63*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyEnd(x));
6497929ea7SJunchao Zhang 
65*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\nTesting VecScatterCreateToAll\n"));
6697929ea7SJunchao Zhang 
67*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterCreateToAll(x,&vscat,&y));
68*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)y,"y"));
6997929ea7SJunchao Zhang 
7097929ea7SJunchao Zhang   /* Test PetscSFBcastAndOp with op = MPI_REPLACE, which does y = x on all ranks */
71*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD));
72*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD));
73*5f80ce2aSJacob Faibussowitsch   if (rank == 0) CHKERRQ(VecView(y,PETSC_VIEWER_STDOUT_SELF));
7497929ea7SJunchao Zhang 
7597929ea7SJunchao Zhang   /* Test PetscSFBcastAndOp with op = MPI_SUM, which does y += x */
76*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(vscat,x,y,ADD_VALUES,SCATTER_FORWARD));
77*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(vscat,x,y,ADD_VALUES,SCATTER_FORWARD));
78*5f80ce2aSJacob Faibussowitsch   if (rank == 0) CHKERRQ(VecView(y,PETSC_VIEWER_STDOUT_SELF));
7997929ea7SJunchao Zhang 
8097929ea7SJunchao Zhang   /* Test PetscSFReduce with op = MPI_REPLACE, which does x = y */
81*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(vscat,y,x,INSERT_VALUES,SCATTER_REVERSE));
82*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(vscat,y,x,INSERT_VALUES,SCATTER_REVERSE));
83*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(x,PETSC_VIEWER_STDOUT_WORLD));
8497929ea7SJunchao Zhang 
8597929ea7SJunchao Zhang   /* Test PetscSFReduce with op = MPI_SUM, which does x += size*y */
86*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(vscat,y,x,ADD_VALUES,SCATTER_REVERSE));
87*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(vscat,y,x,ADD_VALUES,SCATTER_REVERSE));
88*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(x,PETSC_VIEWER_STDOUT_WORLD));
89*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
90*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&y));
91*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterDestroy(&vscat));
9297929ea7SJunchao Zhang 
9397929ea7SJunchao Zhang   ierr = PetscFinalize();
9497929ea7SJunchao Zhang   return ierr;
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