Actual source code: iscomp.c

  1: #define PETSCVEC_DLL

 3:  #include petscsys.h
 4:  #include petscis.h

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

 11:    Collective on IS

 13:    Input Parameters:
 14: .  is1, is2 - The index sets being compared

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

 21:    Level: intermediate

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

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

 34:     Concepts: index sets^equal
 35:     Concepts: IS^equal

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


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

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

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

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

 80:       PetscMemcpy(a1,ptr1,sz1*sizeof(PetscInt));
 81:       PetscMemcpy(a2,ptr2,sz2*sizeof(PetscInt));

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

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