xref: /petsc/src/vec/is/sf/tests/ex8.c (revision dfd57a172ac9fa6c7b5fe6de6ab5df85cefc2996)
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;
13ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
14ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr);
1597929ea7SJunchao Zhang 
1697929ea7SJunchao Zhang   ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr);
1797929ea7SJunchao Zhang   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
1897929ea7SJunchao Zhang   ierr = VecSetSizes(x,PETSC_DECIDE,N);CHKERRQ(ierr);
1997929ea7SJunchao Zhang   ierr = VecGetOwnershipRange(x,&low,&high);CHKERRQ(ierr);
2097929ea7SJunchao Zhang   ierr = PetscObjectSetName((PetscObject)x,"x");CHKERRQ(ierr);
2197929ea7SJunchao Zhang 
2297929ea7SJunchao Zhang   /*-------------------------------------*/
2397929ea7SJunchao Zhang   /*       VecScatterCreateToZero        */
2497929ea7SJunchao Zhang   /*-------------------------------------*/
2597929ea7SJunchao Zhang 
2697929ea7SJunchao Zhang   /* MPI vec x = [0, 1, 2, .., N-1] */
2797929ea7SJunchao Zhang   for (i=low; i<high; i++) {ierr = VecSetValue(x,i,(PetscScalar)i,INSERT_VALUES);CHKERRQ(ierr);}
2897929ea7SJunchao Zhang   ierr = VecAssemblyBegin(x);CHKERRQ(ierr);
2997929ea7SJunchao Zhang   ierr = VecAssemblyEnd(x);CHKERRQ(ierr);
3097929ea7SJunchao Zhang 
3197929ea7SJunchao Zhang   ierr = PetscPrintf(PETSC_COMM_WORLD,"\nTesting VecScatterCreateToZero\n");CHKERRQ(ierr);
3297929ea7SJunchao Zhang   ierr = VecScatterCreateToZero(x,&vscat,&y);CHKERRQ(ierr);
3397929ea7SJunchao Zhang   ierr = PetscObjectSetName((PetscObject)y,"y");CHKERRQ(ierr);
3497929ea7SJunchao Zhang 
3597929ea7SJunchao Zhang   /* Test PetscSFBcastAndOp with op = MPI_REPLACE, which does y = x on rank 0 */
3697929ea7SJunchao Zhang   ierr = VecScatterBegin(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
3797929ea7SJunchao Zhang   ierr = VecScatterEnd(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
3897929ea7SJunchao Zhang   if (!rank) {ierr = VecView(y,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);}
3997929ea7SJunchao Zhang 
4097929ea7SJunchao Zhang   /* Test PetscSFBcastAndOp with op = MPI_SUM, which does y += x */
4197929ea7SJunchao Zhang   ierr = VecScatterBegin(vscat,x,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4297929ea7SJunchao Zhang   ierr = VecScatterEnd(vscat,x,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4397929ea7SJunchao Zhang   if (!rank) {ierr = VecView(y,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);}
4497929ea7SJunchao Zhang 
4597929ea7SJunchao Zhang   /* Test PetscSFReduce with op = MPI_REPLACE, which does x = y */
4697929ea7SJunchao Zhang   ierr = VecScatterBegin(vscat,y,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
4797929ea7SJunchao Zhang   ierr = VecScatterEnd(vscat,y,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
4897929ea7SJunchao Zhang   ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
4997929ea7SJunchao Zhang 
5097929ea7SJunchao Zhang   /* Test PetscSFReduce with op = MPI_SUM, which does x += y on x's local part on rank 0*/
5197929ea7SJunchao Zhang   ierr = VecScatterBegin(vscat,y,x,ADD_VALUES,SCATTER_REVERSE_LOCAL);CHKERRQ(ierr);
5297929ea7SJunchao Zhang   ierr = VecScatterEnd(vscat,y,x,ADD_VALUES,SCATTER_REVERSE_LOCAL);CHKERRQ(ierr);
5397929ea7SJunchao Zhang   ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
5497929ea7SJunchao Zhang 
5597929ea7SJunchao Zhang   ierr = VecDestroy(&y);CHKERRQ(ierr);
5697929ea7SJunchao Zhang   ierr = VecScatterDestroy(&vscat);CHKERRQ(ierr);
5797929ea7SJunchao Zhang 
5897929ea7SJunchao Zhang   /*-------------------------------------*/
5997929ea7SJunchao Zhang   /*       VecScatterCreateToAll         */
6097929ea7SJunchao Zhang   /*-------------------------------------*/
6197929ea7SJunchao Zhang   for (i=low; i<high; i++) {ierr = VecSetValue(x,i,(PetscScalar)i,INSERT_VALUES);CHKERRQ(ierr);}
6297929ea7SJunchao Zhang   ierr = VecAssemblyBegin(x);CHKERRQ(ierr);
6397929ea7SJunchao Zhang   ierr = VecAssemblyEnd(x);CHKERRQ(ierr);
6497929ea7SJunchao Zhang 
6597929ea7SJunchao Zhang   ierr = PetscPrintf(PETSC_COMM_WORLD,"\nTesting VecScatterCreateToAll\n");CHKERRQ(ierr);
6697929ea7SJunchao Zhang 
6797929ea7SJunchao Zhang   ierr = VecScatterCreateToAll(x,&vscat,&y);CHKERRQ(ierr);
6897929ea7SJunchao Zhang   ierr = PetscObjectSetName((PetscObject)y,"y");CHKERRQ(ierr);
6997929ea7SJunchao Zhang 
7097929ea7SJunchao Zhang   /* Test PetscSFBcastAndOp with op = MPI_REPLACE, which does y = x on all ranks */
7197929ea7SJunchao Zhang   ierr = VecScatterBegin(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7297929ea7SJunchao Zhang   ierr = VecScatterEnd(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7397929ea7SJunchao Zhang   if (!rank) {ierr = VecView(y,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);}
7497929ea7SJunchao Zhang 
7597929ea7SJunchao Zhang   /* Test PetscSFBcastAndOp with op = MPI_SUM, which does y += x */
7697929ea7SJunchao Zhang   ierr = VecScatterBegin(vscat,x,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7797929ea7SJunchao Zhang   ierr = VecScatterEnd(vscat,x,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7897929ea7SJunchao Zhang   if (!rank) {ierr = VecView(y,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);}
7997929ea7SJunchao Zhang 
8097929ea7SJunchao Zhang   /* Test PetscSFReduce with op = MPI_REPLACE, which does x = y */
8197929ea7SJunchao Zhang   ierr = VecScatterBegin(vscat,y,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8297929ea7SJunchao Zhang   ierr = VecScatterEnd(vscat,y,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8397929ea7SJunchao Zhang   ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
8497929ea7SJunchao Zhang 
8597929ea7SJunchao Zhang   /* Test PetscSFReduce with op = MPI_SUM, which does x += size*y */
8697929ea7SJunchao Zhang   ierr = VecScatterBegin(vscat,y,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8797929ea7SJunchao Zhang   ierr = VecScatterEnd(vscat,y,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8897929ea7SJunchao Zhang   ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
8997929ea7SJunchao Zhang   ierr = VecDestroy(&x);CHKERRQ(ierr);
9097929ea7SJunchao Zhang   ierr = VecDestroy(&y);CHKERRQ(ierr);
9197929ea7SJunchao Zhang   ierr = VecScatterDestroy(&vscat);CHKERRQ(ierr);
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
11326e8e884SScott Kruger         args: -vec_type cuda -sf_backend cuda -vecscatter_packongpu true
11497929ea7SJunchao Zhang         requires: cuda
11597929ea7SJunchao Zhang 
11697929ea7SJunchao Zhang       test:
11726e8e884SScott Kruger         suffix: 1_hip
11826e8e884SScott Kruger         args: -vec_type hip -sf_backend hip -vecscatter_packongpu true
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
12426e8e884SScott Kruger         args: -vec_type cuda -sf_backend cuda -vecscatter_packongpu false
125*dfd57a17SPierre 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
14126e8e884SScott Kruger         args: -vec_type cuda -sf_backend cuda -vecscatter_packongpu true
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
14726e8e884SScott Kruger         args: -vec_type hip -sf_backend hip -vecscatter_packongpu true
14826e8e884SScott Kruger         requires: hip
14926e8e884SScott Kruger 
15026e8e884SScott Kruger       test:
15197929ea7SJunchao Zhang         suffix: 2_cuda_aware_mpi
15297929ea7SJunchao Zhang         args: -vec_type cuda -vecscatter_packongpu false
153*dfd57a17SPierre Jolivet         requires: cuda defined(PETSC_HAVE_MPI_GPU_AWARE)
15497929ea7SJunchao Zhang 
15597929ea7SJunchao Zhang TEST*/
15697929ea7SJunchao Zhang 
157