Actual source code: ex1.c

petsc-3.6.4 2016-04-12
Report Typos and Errors
  1: /*
  2:        Formatted test for ISGeneral routines.
  3: */

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

  7: #include <petscis.h>
  8: #include <petscviewer.h>

 12: int main(int argc,char **argv)
 13: {
 14:   PetscMPIInt    rank,size;
 15:   PetscInt       i,n,*indices;
 16:   const PetscInt *ii;
 17:   IS             is,newis;
 18:   PetscBool      flg;

 21:   PetscInitialize(&argc,&argv,(char*)0,help);
 22:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 23:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

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

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

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

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

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

 73:   /*
 74:      Thinks it is a different type?
 75:   */
 76:   PetscObjectTypeCompare((PetscObject)is,ISSTRIDE,&flg);
 77:   if (flg) SETERRQ(PETSC_COMM_SELF,1,"ISStride");
 78:   PetscObjectTypeCompare((PetscObject)is,ISBLOCK,&flg);
 79:   if (flg) SETERRQ(PETSC_COMM_SELF,1,"ISBlock");

 81:   ISDestroy(&is);

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