Actual source code: ex2.c
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>
21: int main(int argc,char **argv)
22: {
23: PetscInt i,n,first,step;
24: IS set;
25: const PetscInt *indices;
27: PetscInitialize(&argc,&argv,(char*)0,help);
29: n = 10;
30: first = 3;
31: step = 2;
33: /*
34: Create stride index set, starting at 3 with a stride of 2
35: Note each processor is generating its own index set
36: (in this case they are all identical)
37: */
38: ISCreateStride(PETSC_COMM_SELF,n,first,step,&set);
39: ISView(set,PETSC_VIEWER_STDOUT_SELF);
41: /*
42: Extract indices from set.
43: */
44: ISGetIndices(set,&indices);
45: PetscPrintf(PETSC_COMM_WORLD,"Printing indices directly\n");
46: for (i=0; i<n; i++) {
47: PetscPrintf(PETSC_COMM_WORLD,"%" PetscInt_FMT "\n",indices[i]);
48: }
50: ISRestoreIndices(set,&indices);
52: /*
53: Determine information on stride
54: */
55: ISStrideGetInfo(set,&first,&step);
57: ISDestroy(&set);
58: PetscFinalize();
59: return 0;
60: }
62: /*TEST
64: test:
66: TEST*/