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