Actual source code: iscomp.c

  1: /*$Id: iscomp.c,v 1.35 2001/03/23 23:21:16 balay Exp $*/

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

  6: #undef __FUNCT__  
  8: /*@C
  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: int ISEqual(IS is1,IS is2,PetscTruth *flg)
 39: {
 40:   int        sz1,sz2,ierr,*ptr1,*ptr2,*a1,*a2;
 41:   PetscTruth flag;
 42:   MPI_Comm   comm;

 49: 
 50:   ISGetSize(is1,&sz1);
 51:   ISGetSize(is2,&sz2);
 52:   if (sz1 != sz2) {
 53:     *flg = PETSC_FALSE;
 54:   } else {
 55:     ISGetLocalSize(is1,&sz1);
 56:     ISGetLocalSize(is2,&sz2);

 58:     if (sz1 != sz2) {
 59:       flag = PETSC_FALSE;
 60:     } else {
 61:       ISGetIndices(is1,&ptr1);
 62:       ISGetIndices(is2,&ptr2);
 63: 
 64:       PetscMalloc((sz1+1)*sizeof(int),&a1);
 65:       PetscMalloc((sz2+1)*sizeof(int),&a2);

 67:       PetscMemcpy(a1,ptr1,sz1*sizeof(int));
 68:       PetscMemcpy(a2,ptr2,sz2*sizeof(int));

 70:       PetscSortInt(sz1,a1);
 71:       PetscSortInt(sz2,a2);
 72:       PetscMemcmp(a1,a2,sz1*sizeof(int),&flag);

 74:       ISRestoreIndices(is1,&ptr1);
 75:       ISRestoreIndices(is2,&ptr2);
 76: 
 77:       PetscFree(a1);
 78:       PetscFree(a2);
 79:     }
 80:     PetscObjectGetComm((PetscObject)is1,&comm);
 81:     MPI_Allreduce(&flag,flg,1,MPI_INT,MPI_MIN,comm);
 82:   }
 83:   return(0);
 84: }
 85: