Actual source code: iscomp.c
petsc-3.8.4 2018-03-24
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.
30: Concepts: index sets^equal
31: Concepts: IS^equal
33: @*/
34: PetscErrorCode ISEqual(IS is1,IS is2,PetscBool *flg)
35: {
36: PetscInt sz1,sz2,*a1,*a2;
37: const PetscInt *ptr1,*ptr2;
38: PetscBool flag;
39: MPI_Comm comm;
41: PetscMPIInt mflg;
48: if (is1 == is2) {
49: *flg = PETSC_TRUE;
50: return(0);
51: }
53: MPI_Comm_compare(PetscObjectComm((PetscObject)is1),PetscObjectComm((PetscObject)is2),&mflg);
54: if (mflg != MPI_CONGRUENT && mflg != MPI_IDENT) {
55: *flg = PETSC_FALSE;
56: return(0);
57: }
59: ISGetSize(is1,&sz1);
60: ISGetSize(is2,&sz2);
61: if (sz1 != sz2) *flg = PETSC_FALSE;
62: else {
63: ISGetLocalSize(is1,&sz1);
64: ISGetLocalSize(is2,&sz2);
66: if (sz1 != sz2) flag = PETSC_FALSE;
67: else {
68: ISGetIndices(is1,&ptr1);
69: ISGetIndices(is2,&ptr2);
71: PetscMalloc1(sz1,&a1);
72: PetscMalloc1(sz2,&a2);
74: PetscMemcpy(a1,ptr1,sz1*sizeof(PetscInt));
75: PetscMemcpy(a2,ptr2,sz2*sizeof(PetscInt));
77: PetscSortInt(sz1,a1);
78: PetscSortInt(sz2,a2);
79: PetscMemcmp(a1,a2,sz1*sizeof(PetscInt),&flag);
81: ISRestoreIndices(is1,&ptr1);
82: ISRestoreIndices(is2,&ptr2);
84: PetscFree(a1);
85: PetscFree(a2);
86: }
87: PetscObjectGetComm((PetscObject)is1,&comm);
88: MPIU_Allreduce(&flag,flg,1,MPIU_BOOL,MPI_MIN,comm);
89: }
90: return(0);
91: }