Actual source code: ex2.c
petsc-3.7.5 2017-01-01
2: static char help[] = "Demonstrates creating a stride index set.\n\n";
4: /*T
5: Concepts: index sets^creating a stride index set;
6: Concepts: stride^creating a stride index set;
7: Concepts: IS^creating a stride index set;
9: Description: Creates an index set based on a stride. Views that index set
10: and then destroys it.
11: T*/
13: /*
14: Include petscis.h so we can use PETSc IS objects. Note that this automatically
15: includes petscsys.h.
16: */
18: #include <petscis.h>
19: #include <petscviewer.h>
23: int main(int argc,char **argv)
24: {
26: PetscInt i,n,first,step;
27: IS set;
28: const PetscInt *indices;
30: PetscInitialize(&argc,&argv,(char*)0,help);
32: n = 10;
33: first = 3;
34: step = 2;
36: /*
37: Create stride index set, starting at 3 with a stride of 2
38: Note each processor is generating its own index set
39: (in this case they are all identical)
40: */
41: ISCreateStride(PETSC_COMM_SELF,n,first,step,&set);
42: ISView(set,PETSC_VIEWER_STDOUT_SELF);
44: /*
45: Extract indices from set.
46: */
47: ISGetIndices(set,&indices);
48: PetscPrintf(PETSC_COMM_WORLD,"Printing indices directly\n");
49: for (i=0; i<n; i++) {
50: PetscPrintf(PETSC_COMM_WORLD,"%D\n",indices[i]);
51: }
53: ISRestoreIndices(set,&indices);
55: /*
56: Determine information on stride
57: */
58: ISStrideGetInfo(set,&first,&step);
59: if (first != 3 || step != 2) SETERRQ(PETSC_COMM_SELF,1,"Stride info not correct!\n");
60: ISDestroy(&set);
61: PetscFinalize();
62: return 0;
63: }