Actual source code: iscomp.c

petsc-3.14.6 2021-03-30
Report Typos and Errors

  2: #include <petsc/private/isimpl.h>

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

  7:    Collective on IS

  9:    Input Parameters:
 10: .  is1, is2 - The index sets being compared

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

 17:    Level: intermediate

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

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


 31: .seealso: ISEqualUnsorted()
 32: @*/
 33: PetscErrorCode  ISEqual(IS is1,IS is2,PetscBool  *flg)
 34: {
 35:   PetscInt       sz1,sz2,*a1,*a2;
 36:   const PetscInt *ptr1,*ptr2;
 37:   PetscBool      flag;
 38:   MPI_Comm       comm;
 40:   PetscMPIInt    mflg;


 47:   if (is1 == is2) {
 48:     *flg = PETSC_TRUE;
 49:     return(0);
 50:   }

 52:   MPI_Comm_compare(PetscObjectComm((PetscObject)is1),PetscObjectComm((PetscObject)is2),&mflg);
 53:   if (mflg != MPI_CONGRUENT && mflg != MPI_IDENT) {
 54:     *flg = PETSC_FALSE;
 55:     return(0);
 56:   }

 58:   ISGetSize(is1,&sz1);
 59:   ISGetSize(is2,&sz2);
 60:   if (sz1 != sz2) *flg = PETSC_FALSE;
 61:   else {
 62:     ISGetLocalSize(is1,&sz1);
 63:     ISGetLocalSize(is2,&sz2);

 65:     if (sz1 != sz2) flag = PETSC_FALSE;
 66:     else {
 67:       ISGetIndices(is1,&ptr1);
 68:       ISGetIndices(is2,&ptr2);

 70:       PetscMalloc1(sz1,&a1);
 71:       PetscMalloc1(sz2,&a2);

 73:       PetscArraycpy(a1,ptr1,sz1);
 74:       PetscArraycpy(a2,ptr2,sz2);

 76:       PetscIntSortSemiOrdered(sz1,a1);
 77:       PetscIntSortSemiOrdered(sz2,a2);
 78:       PetscArraycmp(a1,a2,sz1,&flag);

 80:       ISRestoreIndices(is1,&ptr1);
 81:       ISRestoreIndices(is2,&ptr2);

 83:       PetscFree(a1);
 84:       PetscFree(a2);
 85:     }
 86:     PetscObjectGetComm((PetscObject)is1,&comm);
 87:     MPIU_Allreduce(&flag,flg,1,MPIU_BOOL,MPI_MIN,comm);
 88:   }
 89:   return(0);
 90: }

 92: /*@
 93:    ISEqualUnsorted  - Compares if two index sets have the same set of indices.

 95:    Collective on IS

 97:    Input Parameters:
 98: .  is1, is2 - The index sets being compared

100:    Output Parameters:
101: .  flg - output flag, either PETSC_TRUE (if both index sets have the
102:          same indices), or PETSC_FALSE if the index sets differ by size
103:          or by the set of indices)

105:    Level: intermediate

107:    Note:
108:    This routine does NOT sort the contents of the index sets before
109:    the comparision is made.


112: .seealso: ISEqual()
113: @*/
114: PetscErrorCode  ISEqualUnsorted(IS is1,IS is2,PetscBool  *flg)
115: {
116:   PetscInt       sz1,sz2;
117:   const PetscInt *ptr1,*ptr2;
118:   PetscBool      flag;
119:   MPI_Comm       comm;
121:   PetscMPIInt    mflg;


128:   if (is1 == is2) {
129:     *flg = PETSC_TRUE;
130:     return(0);
131:   }

133:   MPI_Comm_compare(PetscObjectComm((PetscObject)is1),PetscObjectComm((PetscObject)is2),&mflg);
134:   if (mflg != MPI_CONGRUENT && mflg != MPI_IDENT) {
135:     *flg = PETSC_FALSE;
136:     return(0);
137:   }

139:   ISGetSize(is1,&sz1);
140:   ISGetSize(is2,&sz2);
141:   if (sz1 != sz2) *flg = PETSC_FALSE;
142:   else {
143:     ISGetLocalSize(is1,&sz1);
144:     ISGetLocalSize(is2,&sz2);

146:     if (sz1 != sz2) flag = PETSC_FALSE;
147:     else {
148:       ISGetIndices(is1,&ptr1);
149:       ISGetIndices(is2,&ptr2);

151:       PetscArraycmp(ptr1,ptr2,sz1,&flag);

153:       ISRestoreIndices(is1,&ptr1);
154:       ISRestoreIndices(is2,&ptr2);
155:     }
156:     PetscObjectGetComm((PetscObject)is1,&comm);
157:     MPIU_Allreduce(&flag,flg,1,MPIU_BOOL,MPI_MIN,comm);
158:   }
159:   return(0);
160: }