Actual source code: iscomp.c

petsc-3.3-p7 2013-05-11
  2: #include <petscis.h>    /*I "petscis.h"  I*/

  6: /*@
  7:    ISEqual  - Compares if two index sets have the same set of indices.

  9:    Collective on IS

 11:    Input Parameters:
 12: .  is1, is2 - The index sets being compared

 14:    Output Parameters:
 15: .  flg - output flag, either PETSC_TRUE (if both index sets have the
 16:          same indices), or PETSC_FALSE if the index sets differ by size 
 17:          or by the set of indices)

 19:    Level: intermediate

 21:    Note: 
 22:    This routine sorts the contents of the index sets before
 23:    the comparision is made, so the order of the indices on a processor is immaterial.

 25:    Each processor has to have the same indices in the two sets, for example,
 26: $           Processor 
 27: $             0      1
 28: $    is1 = {0, 1} {2, 3}
 29: $    is2 = {2, 3} {0, 1}
 30:    will return false.

 32:     Concepts: index sets^equal
 33:     Concepts: IS^equal

 35: @*/
 36: PetscErrorCode  ISEqual(IS is1,IS is2,PetscBool  *flg)
 37: {
 38:   PetscInt       sz1,sz2,*a1,*a2;
 39:   const PetscInt *ptr1,*ptr2;
 40:   PetscBool      flag;
 41:   MPI_Comm       comm;
 43:   PetscMPIInt    mflg;


 50:   if (is1 == is2) {
 51:     *flg = PETSC_TRUE;
 52:     return(0);
 53:   }

 55:   MPI_Comm_compare(((PetscObject)is1)->comm,((PetscObject)is2)->comm,&mflg);
 56:   if (mflg != MPI_CONGRUENT && mflg != MPI_IDENT) {
 57:     *flg = PETSC_FALSE;
 58:     return(0);
 59:   }

 61:   ISGetSize(is1,&sz1);
 62:   ISGetSize(is2,&sz2);
 63:   if (sz1 != sz2) {
 64:     *flg = PETSC_FALSE;
 65:   } else {
 66:     ISGetLocalSize(is1,&sz1);
 67:     ISGetLocalSize(is2,&sz2);

 69:     if (sz1 != sz2) {
 70:       flag = PETSC_FALSE;
 71:     } else {
 72:       ISGetIndices(is1,&ptr1);
 73:       ISGetIndices(is2,&ptr2);
 74: 
 75:       PetscMalloc(sz1*sizeof(PetscInt),&a1);
 76:       PetscMalloc(sz2*sizeof(PetscInt),&a2);

 78:       PetscMemcpy(a1,ptr1,sz1*sizeof(PetscInt));
 79:       PetscMemcpy(a2,ptr2,sz2*sizeof(PetscInt));

 81:       PetscSortInt(sz1,a1);
 82:       PetscSortInt(sz2,a2);
 83:       PetscMemcmp(a1,a2,sz1*sizeof(PetscInt),&flag);

 85:       ISRestoreIndices(is1,&ptr1);
 86:       ISRestoreIndices(is2,&ptr2);
 87: 
 88:       PetscFree(a1);
 89:       PetscFree(a2);
 90:     }
 91:     PetscObjectGetComm((PetscObject)is1,&comm);
 92:     MPI_Allreduce(&flag,flg,1,MPI_INT,MPI_MIN,comm);
 93:   }
 94:   return(0);
 95: }
 96: