Actual source code: ex8.c
petsc-3.14.6 2021-03-30
1: static char help[]= "Test VecScatterCreateToZero, VecScatterCreateToAll\n\n";
3: #include <petscvec.h>
4: int main(int argc,char **argv)
5: {
6: PetscErrorCode ierr;
7: PetscInt i,N=10,low,high;
8: PetscMPIInt size,rank;
9: Vec x,y;
10: VecScatter vscat;
12: PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
13: MPI_Comm_size(PETSC_COMM_WORLD,&size);
14: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
16: VecCreate(PETSC_COMM_WORLD,&x);
17: VecSetFromOptions(x);
18: VecSetSizes(x,PETSC_DECIDE,N);
19: VecGetOwnershipRange(x,&low,&high);
20: PetscObjectSetName((PetscObject)x,"x");
22: /*-------------------------------------*/
23: /* VecScatterCreateToZero */
24: /*-------------------------------------*/
26: /* MPI vec x = [0, 1, 2, .., N-1] */
27: for (i=low; i<high; i++) {VecSetValue(x,i,(PetscScalar)i,INSERT_VALUES);}
28: VecAssemblyBegin(x);
29: VecAssemblyEnd(x);
31: PetscPrintf(PETSC_COMM_WORLD,"\nTesting VecScatterCreateToZero\n");
32: VecScatterCreateToZero(x,&vscat,&y);
33: PetscObjectSetName((PetscObject)y,"y");
35: /* Test PetscSFBcastAndOp with op = MPI_REPLACE, which does y = x on rank 0 */
36: VecScatterBegin(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD);
37: VecScatterEnd(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD);
38: if (!rank) {VecView(y,PETSC_VIEWER_STDOUT_SELF);}
40: /* Test PetscSFBcastAndOp with op = MPI_SUM, which does y += x */
41: VecScatterBegin(vscat,x,y,ADD_VALUES,SCATTER_FORWARD);
42: VecScatterEnd(vscat,x,y,ADD_VALUES,SCATTER_FORWARD);
43: if (!rank) {VecView(y,PETSC_VIEWER_STDOUT_SELF);}
45: /* Test PetscSFReduce with op = MPI_REPLACE, which does x = y */
46: VecScatterBegin(vscat,y,x,INSERT_VALUES,SCATTER_REVERSE);
47: VecScatterEnd(vscat,y,x,INSERT_VALUES,SCATTER_REVERSE);
48: VecView(x,PETSC_VIEWER_STDOUT_WORLD);
50: /* Test PetscSFReduce with op = MPI_SUM, which does x += y on x's local part on rank 0*/
51: VecScatterBegin(vscat,y,x,ADD_VALUES,SCATTER_REVERSE_LOCAL);
52: VecScatterEnd(vscat,y,x,ADD_VALUES,SCATTER_REVERSE_LOCAL);
53: VecView(x,PETSC_VIEWER_STDOUT_WORLD);
55: VecDestroy(&y);
56: VecScatterDestroy(&vscat);
58: /*-------------------------------------*/
59: /* VecScatterCreateToAll */
60: /*-------------------------------------*/
61: for (i=low; i<high; i++) {VecSetValue(x,i,(PetscScalar)i,INSERT_VALUES);}
62: VecAssemblyBegin(x);
63: VecAssemblyEnd(x);
65: PetscPrintf(PETSC_COMM_WORLD,"\nTesting VecScatterCreateToAll\n");
67: VecScatterCreateToAll(x,&vscat,&y);
68: PetscObjectSetName((PetscObject)y,"y");
70: /* Test PetscSFBcastAndOp with op = MPI_REPLACE, which does y = x on all ranks */
71: VecScatterBegin(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD);
72: VecScatterEnd(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD);
73: if (!rank) {VecView(y,PETSC_VIEWER_STDOUT_SELF);}
75: /* Test PetscSFBcastAndOp with op = MPI_SUM, which does y += x */
76: VecScatterBegin(vscat,x,y,ADD_VALUES,SCATTER_FORWARD);
77: VecScatterEnd(vscat,x,y,ADD_VALUES,SCATTER_FORWARD);
78: if (!rank) {VecView(y,PETSC_VIEWER_STDOUT_SELF);}
80: /* Test PetscSFReduce with op = MPI_REPLACE, which does x = y */
81: VecScatterBegin(vscat,y,x,INSERT_VALUES,SCATTER_REVERSE);
82: VecScatterEnd(vscat,y,x,INSERT_VALUES,SCATTER_REVERSE);
83: VecView(x,PETSC_VIEWER_STDOUT_WORLD);
85: /* Test PetscSFReduce with op = MPI_SUM, which does x += size*y */
86: VecScatterBegin(vscat,y,x,ADD_VALUES,SCATTER_REVERSE);
87: VecScatterEnd(vscat,y,x,ADD_VALUES,SCATTER_REVERSE);
88: VecView(x,PETSC_VIEWER_STDOUT_WORLD);
89: VecDestroy(&x);
90: VecDestroy(&y);
91: VecScatterDestroy(&vscat);
93: PetscFinalize();
94: return ierr;
95: }
97: /*TEST
99: testset:
100: # N=10 is divisible by nsize, to trigger Allgather/Gather in SF
101: nsize: 2
102: filter: grep -v "type"
103: output_file: output/ex8_1.out
105: test:
106: suffix: 1_standard
108: test:
109: suffix: 1_cuda
110: args: -vec_type cuda -vecscatter_packongpu true
111: requires: cuda
113: test:
114: suffix: 1_cuda_aware_mpi
115: args: -vec_type cuda -vecscatter_packongpu false
116: requires: cuda define(PETSC_HAVE_MPI_GPU_AWARE)
118: testset:
119: # N=10 is not divisible by nsize, to trigger Allgatherv/Gatherv in SF
120: suffix: 2
121: nsize: 3
122: filter: grep -v "type"
123: output_file: output/ex8_2.out
125: test:
126: suffix: 2_standard
128: test:
129: suffix: 2_cuda
130: args: -vec_type cuda -vecscatter_packongpu true
131: requires: cuda
133: test:
134: suffix: 2_cuda_aware_mpi
135: args: -vec_type cuda -vecscatter_packongpu false
136: requires: cuda define(PETSC_HAVE_MPI_GPU_AWARE)
138: TEST*/