Actual source code: ex1.c

  1: /*$Id: ex1.c,v 1.41 2001/04/10 19:34:49 bsmith Exp $*/
  2: /*
  3:        Formatted test for ISGeneral routines.
  4: */

  6: static char help[] = "Tests IS general routines.\n\n";

 8:  #include petscis.h

 12: int main(int argc,char **argv)
 13: {
 14:   int        i,n,ierr,*indices,rank,size,*ii;
 15:   IS         is,newis;
 16:   PetscTruth flg;

 18:   PetscInitialize(&argc,&argv,(char*)0,help);
 19:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 20:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

 22:   /*
 23:      Test IS of size 0 
 24:   */
 25:   ISCreateGeneral(PETSC_COMM_SELF,0,&n,&is);
 26:   ISGetSize(is,&n);
 27:   if (n != 0) SETERRQ(1,"ISGetSize");
 28:   ISDestroy(is);

 30:   /*
 31:      Create large IS and test ISGetIndices()
 32:   */
 33:   n = 10000 + rank;
 34:   PetscMalloc(n*sizeof(int),&indices);
 35:   for (i=0; i<n; i++) {
 36:     indices[i] = rank + i;
 37:   }
 38:   ISCreateGeneral(PETSC_COMM_SELF,n,indices,&is);
 39:   ISGetIndices(is,&ii);
 40:   for (i=0; i<n; i++) {
 41:     if (ii[i] != indices[i]) SETERRQ(1,"ISGetIndices");
 42:   }
 43:   ISRestoreIndices(is,&ii);

 45:   /* 
 46:      Check identity and permutation 
 47:   */
 48:   ISPermutation(is,&flg);
 49:   if (flg == PETSC_TRUE) SETERRQ(1,"ISPermutation");
 50:   ISIdentity(is,&flg);
 51:   if (flg != PETSC_TRUE) SETERRQ(1,"ISIdentity");
 52:   ISSetPermutation(is);
 53:   ISSetIdentity(is);
 54:   ISPermutation(is,&flg);
 55:   if (flg != PETSC_TRUE) SETERRQ(1,"ISPermutation");
 56:   ISIdentity(is,&flg);
 57:   if (flg != PETSC_TRUE) SETERRQ(1,"ISIdentity");

 59:   /*
 60:      Check equality of index sets 
 61:   */
 62:   ISEqual(is,is,&flg);
 63:   if (flg != PETSC_TRUE) SETERRQ(1,"ISEqual");

 65:   /*
 66:      Sorting 
 67:   */
 68:   ISSort(is);
 69:   ISSorted(is,&flg);
 70:   if (flg != PETSC_TRUE) SETERRQ(1,"ISSort");

 72:   /*
 73:      Thinks it is a different type?
 74:   */
 75:   ISStride(is,&flg);
 76:   if (flg == PETSC_TRUE) SETERRQ(1,"ISStride");
 77:   ISBlock(is,&flg);
 78:   if (flg == PETSC_TRUE) SETERRQ(1,"ISBlock");

 80:   ISDestroy(is);

 82:   /*
 83:      Inverting permutation
 84:   */
 85:   for (i=0; i<n; i++) {
 86:     indices[i] = n - i - 1;
 87:   }
 88:   ISCreateGeneral(PETSC_COMM_SELF,n,indices,&is);
 89:   PetscFree(indices);
 90:   ISSetPermutation(is);
 91:   ISInvertPermutation(is,PETSC_DECIDE,&newis);
 92:   ISGetIndices(newis,&ii);
 93:   for (i=0; i<n; i++) {
 94:     if (ii[i] != n - i - 1) SETERRQ(1,"ISInvertPermutation");
 95:   }
 96:   ISRestoreIndices(newis,&ii);
 97:   ISDestroy(newis);
 98:   ISDestroy(is);
 99:   PetscFinalize();
100:   return 0;
101: }
102: