xref: /petsc/src/vec/is/sf/tests/ex21.c (revision 58d68138c660dfb4e9f5b03334792cd4f2ffd7cc)
1 static const char help[] = "Test VecScatterCopy() on an SF with duplicated leaves \n\n";
2 
3 #include <petscvec.h>
4 #include <petscsf.h>
5 
6 /*
7   Contributed-by: "Hammond, Glenn E" <glenn.hammond@pnnl.gov>
8 */
9 int main(int argc, char *argv[]) {
10   PetscMPIInt size;
11   PetscInt    n;
12   PetscInt   *indices;
13   Vec         vec;
14   Vec         vec2;
15   IS          is;
16   IS          is2;
17   VecScatter  scatter;
18   VecScatter  scatter2;
19 
20   PetscFunctionBeginUser;
21   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
22   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
23   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This example requires 1 process");
24 
25   n = 4;
26   PetscCall(PetscMalloc1(n, &indices));
27   indices[0] = 0;
28   indices[1] = 1;
29   indices[2] = 2;
30   indices[3] = 3;
31   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, n, indices, PETSC_COPY_VALUES, &is));
32   PetscCall(PetscFree(indices));
33   PetscCall(VecCreateMPI(PETSC_COMM_WORLD, n, n, &vec));
34 
35   n = 4;
36   PetscCall(PetscMalloc1(n, &indices));
37   indices[0] = 0;
38   indices[1] = 0;
39   indices[2] = 1;
40   indices[3] = 1;
41   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, n, indices, PETSC_COPY_VALUES, &is2));
42   PetscCall(PetscFree(indices));
43   PetscCall(VecCreateMPI(PETSC_COMM_WORLD, n / 2, n / 2, &vec2));
44 
45   PetscCall(VecScatterCreate(vec, is, vec2, is2, &scatter));
46   PetscCall(ISDestroy(&is));
47   PetscCall(ISDestroy(&is2));
48 
49   PetscCall(VecScatterCopy(scatter, &scatter2));
50 
51   PetscCall(VecDestroy(&vec));
52   PetscCall(VecDestroy(&vec2));
53   PetscCall(VecScatterDestroy(&scatter));
54   PetscCall(VecScatterDestroy(&scatter2));
55   PetscCall(PetscFinalize());
56 }
57 
58 /*TEST
59   test:
60     nsize: 1
61     output_file: output/empty.out
62 TEST*/
63