xref: /petsc/src/vec/is/is/tutorials/ex5.c (revision 013f9f9e983abb680ff0f504877dbdb0e1ec26c3)
1 
2 static char help[] = "Demonstrates using ISLocalToGlobalMappings with block size.\n\n";
3 
4 #include <petscis.h>
5 #include <petscviewer.h>
6 
7 int main(int argc, char **argv)
8 {
9   PetscInt               i, n = 4, indices[] = {0, 3, 9, 12}, m = 2, input[] = {0, 2};
10   PetscInt               output[2], inglobals[13], outlocals[13];
11   ISLocalToGlobalMapping mapping;
12 
13   PetscFunctionBeginUser;
14   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
15 
16   /*
17       Create a local to global mapping. Each processor independently
18      creates a mapping
19   */
20   PetscCall(PetscIntView(n, indices, PETSC_VIEWER_STDOUT_WORLD));
21   PetscCall(ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD, 2, n, indices, PETSC_COPY_VALUES, &mapping));
22 
23   /*
24      Map a set of local indices to their global values
25   */
26   PetscCall(PetscIntView(m, input, PETSC_VIEWER_STDOUT_WORLD));
27   PetscCall(ISLocalToGlobalMappingApply(mapping, m, input, output));
28   PetscCall(PetscIntView(m, output, PETSC_VIEWER_STDOUT_WORLD));
29 
30   /*
31      Map some global indices to local, retaining the ones without a local index by -1
32   */
33   for (i = 0; i < 13; i++) inglobals[i] = i;
34   PetscCall(PetscIntView(13, inglobals, PETSC_VIEWER_STDOUT_WORLD));
35   PetscCall(ISGlobalToLocalMappingApply(mapping, IS_GTOLM_MASK, 13, inglobals, NULL, outlocals));
36   PetscCall(PetscIntView(13, outlocals, PETSC_VIEWER_STDOUT_WORLD));
37 
38   /*
39      Map some block global indices to local, dropping the ones without a local index.
40   */
41   PetscCall(PetscIntView(13, inglobals, PETSC_VIEWER_STDOUT_WORLD));
42   PetscCall(ISGlobalToLocalMappingApplyBlock(mapping, IS_GTOLM_DROP, 13, inglobals, &m, outlocals));
43   PetscCall(PetscIntView(m, outlocals, PETSC_VIEWER_STDOUT_WORLD));
44 
45   /*
46      Free the space used by the local to global mapping
47   */
48   PetscCall(ISLocalToGlobalMappingDestroy(&mapping));
49 
50   PetscCall(PetscFinalize());
51   return 0;
52 }
53 
54 /*TEST
55 
56    test:
57 
58 TEST*/
59