Actual source code: ex2.c

  2: /*      "$Id: ex2.c,v 1.20 2001/03/23 23:21:14 balay Exp $"; */

  4: static char help[] = "Demonstrates creating a stride index set.nn";

  6: /*T
  7:     Concepts: index sets^creating a stride index set;
  8:     Concepts: stride^creating a stride index set;
  9:     Concepts: IS^creating a stride index set;
 10:     
 11:     Description: Creates an index set based on a stride. Views that index set
 12:     and then destroys it.
 13: T*/

 15: /*
 16:   Include petscis.h so we can use PETSc IS objects. Note that this automatically 
 17:   includes petsc.h.
 18: */

 20:  #include petscis.h

 22: int main(int argc,char **argv)
 23: {
 24:   int i,n,ierr, *indices,first,step;
 25:   IS  set;

 27:   PetscInitialize(&argc,&argv,(char*)0,help);
 28: 
 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 directlyn");
 46:   for (i=0; i<n; i++) {
 47:     PetscPrintf(PETSC_COMM_WORLD,"%dn",indices[i]);
 48:   }

 50:   ISRestoreIndices(set,&indices);

 52:   /*
 53:       Determine information on stride
 54:   */
 55:   ISStrideGetInfo(set,&first,&step);
 56:   if (first != 3 || step != 2) SETERRQ(1,"Stride info not correct!n");
 57:   ISDestroy(set);
 58:   PetscFinalize();
 59:   return 0;
 60: }