Actual source code: ex2.c

petsc-3.3-p7 2013-05-11
  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;
  8:     
  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>

 22: int main(int argc,char **argv)
 23: {
 25:   PetscInt       i,n,first,step;
 26:   IS             set;
 27:   const PetscInt *indices;

 29:   PetscInitialize(&argc,&argv,(char*)0,help);
 30: 
 31:   n     = 10;
 32:   first = 3;
 33:   step  = 2;

 35:   /*
 36:     Create stride index set, starting at 3 with a stride of 2
 37:     Note each processor is generating its own index set 
 38:     (in this case they are all identical)
 39:   */
 40:   ISCreateStride(PETSC_COMM_SELF,n,first,step,&set);
 41:   ISView(set,PETSC_VIEWER_STDOUT_SELF);

 43:   /*
 44:     Extract indices from set.
 45:   */
 46:   ISGetIndices(set,&indices);
 47:   PetscPrintf(PETSC_COMM_WORLD,"Printing indices directly\n");
 48:   for (i=0; i<n; i++) {
 49:     PetscPrintf(PETSC_COMM_WORLD,"%D\n",indices[i]);
 50:   }

 52:   ISRestoreIndices(set,&indices);

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