xref: /petsc/src/vec/vec/tests/ex24.c (revision 732aec7a18f2199fb53bb9a2f3aef439a834ce31)
1 static char help[] = "Scatters from a parallel vector to a sequential vector.\n\
2 Tests where the local part of the scatter is a copy.\n\n";
3 
4 #include <petscvec.h>
5 
main(int argc,char ** argv)6 int main(int argc, char **argv)
7 {
8   PetscMPIInt size, rank;
9   PetscInt    n = 5, i, *blks, bs = 1, m = 2;
10   PetscScalar value;
11   Vec         x, y;
12   IS          is1, is2;
13   VecScatter  ctx = 0;
14   PetscViewer sviewer;
15 
16   PetscFunctionBeginUser;
17   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
18 
19   PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL));
20   PetscCall(PetscOptionsGetInt(NULL, NULL, "-bs", &bs, NULL));
21 
22   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
23   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
24 
25   /* create two vectors */
26   PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
27   PetscCall(VecSetSizes(x, PETSC_DECIDE, size * bs * n));
28   PetscCall(VecSetFromOptions(x));
29 
30   /* create two index sets */
31   if (rank < size - 1) m = n + 2;
32   else m = n;
33 
34   PetscCall(PetscMalloc1(m, &blks));
35   blks[0] = n * rank;
36   for (i = 1; i < m; i++) blks[i] = blks[i - 1] + 1;
37   PetscCall(ISCreateBlock(PETSC_COMM_SELF, bs, m, blks, PETSC_COPY_VALUES, &is1));
38   PetscCall(PetscFree(blks));
39 
40   PetscCall(VecCreateSeq(PETSC_COMM_SELF, bs * m, &y));
41   PetscCall(ISCreateStride(PETSC_COMM_SELF, bs * m, 0, 1, &is2));
42 
43   /* each processor inserts the entire vector */
44   /* this is redundant but tests assembly */
45   for (i = 0; i < bs * n * size; i++) {
46     value = (PetscScalar)i;
47     PetscCall(VecSetValues(x, 1, &i, &value, INSERT_VALUES));
48   }
49   PetscCall(VecAssemblyBegin(x));
50   PetscCall(VecAssemblyEnd(x));
51   PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
52 
53   PetscCall(VecScatterCreate(x, is1, y, is2, &ctx));
54   PetscCall(VecScatterBegin(ctx, x, y, INSERT_VALUES, SCATTER_FORWARD));
55   PetscCall(VecScatterEnd(ctx, x, y, INSERT_VALUES, SCATTER_FORWARD));
56 
57   PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
58   PetscCall(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &sviewer));
59   PetscCall(VecView(y, sviewer));
60   PetscCall(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &sviewer));
61   PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD));
62 
63   PetscCall(VecScatterDestroy(&ctx));
64 
65   PetscCall(VecDestroy(&x));
66   PetscCall(VecDestroy(&y));
67   PetscCall(ISDestroy(&is1));
68   PetscCall(ISDestroy(&is2));
69 
70   PetscCall(PetscFinalize());
71   return 0;
72 }
73 
74 /*TEST
75 
76    testset:
77       nsize: 3
78       output_file: output/ex24_1.out
79       filter: grep -v "  type:"
80       test:
81         suffix: standard
82         args: -vec_type standard
83       test:
84         requires: cuda
85         suffix: cuda
86         args: -vec_type cuda
87       test:
88         requires: viennacl
89         suffix: viennacl
90         args: -vec_type viennacl
91 
92 TEST*/
93