Actual source code: sftype.c

  1: #include <petsc/private/sfimpl.h>

  3: #if !defined(PETSC_HAVE_MPI_COMBINER_DUP) && !defined(MPI_COMBINER_DUP) /* We have no way to interpret output of MPI_Type_get_envelope without this. */
  4:   #define MPI_COMBINER_DUP 0
  5: #endif
  6: #if !defined(PETSC_HAVE_MPI_COMBINER_NAMED) && !defined(MPI_COMBINER_NAMED)
  7:   #define MPI_COMBINER_NAMED -2
  8: #endif
  9: #if !defined(PETSC_HAVE_MPI_COMBINER_CONTIGUOUS) && !defined(MPI_COMBINER_CONTIGUOUS) && MPI_VERSION < 2
 10:   #define MPI_COMBINER_CONTIGUOUS -1
 11: #endif

 13: static PetscErrorCode MPIPetsc_Type_free(MPI_Datatype *a)
 14: {
 15:   PetscMPIInt nints, naddrs, ntypes, combiner;

 17:   PetscFunctionBegin;
 18:   PetscCallMPI(MPI_Type_get_envelope(*a, &nints, &naddrs, &ntypes, &combiner));

 20:   if (combiner != MPI_COMBINER_NAMED) PetscCallMPI(MPI_Type_free(a));

 22:   *a = MPI_DATATYPE_NULL;
 23:   PetscFunctionReturn(PETSC_SUCCESS);
 24: }

 26: /*
 27:   Unwrap an MPI datatype recursively in case it is dupped or MPI_Type_contiguous(1,...)'ed from another type.

 29:    Input Parameter:
 30: .  a  - the datatype to be unwrapped

 32:    Output Parameters:
 33: + atype - the unwrapped datatype, which is either equal(=) to a or equivalent to a.
 34: - flg   - true if atype != a, which implies caller should MPIPetsc_Type_free(atype) after use. Note atype might be MPI builtin.
 35: */
 36: PetscErrorCode MPIPetsc_Type_unwrap(MPI_Datatype a, MPI_Datatype *atype, PetscBool *flg)
 37: {
 38:   PetscMPIInt  nints, naddrs, ntypes, combiner, ints[1];
 39:   MPI_Aint     addrs[1];
 40:   MPI_Datatype types[1];

 42:   PetscFunctionBegin;
 43:   *flg   = PETSC_FALSE;
 44:   *atype = a;
 45:   if (a == MPIU_INT || a == MPIU_REAL || a == MPIU_SCALAR) PetscFunctionReturn(PETSC_SUCCESS);
 46:   PetscCallMPI(MPI_Type_get_envelope(a, &nints, &naddrs, &ntypes, &combiner));
 47:   if (combiner == MPI_COMBINER_DUP) {
 48:     PetscCheck(nints == 0 && naddrs == 0 && ntypes == 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "Unexpected returns from MPI_Type_get_envelope()");
 49:     PetscCallMPI(MPI_Type_get_contents(a, 0, 0, 1, ints, addrs, types));
 50:     /* Recursively unwrap dupped types. */
 51:     PetscCall(MPIPetsc_Type_unwrap(types[0], atype, flg));
 52:     if (*flg) {
 53:       /* If the recursive call returns a new type, then that means that atype[0] != types[0] and we're on the hook to
 54:        * free types[0].  Note that this case occurs if combiner(types[0]) is MPI_COMBINER_DUP, so we're safe to
 55:        * directly call MPI_Type_free rather than MPIPetsc_Type_free here. */
 56:       PetscCallMPI(MPI_Type_free(&types[0]));
 57:     }
 58:     /* In any case, it's up to the caller to free the returned type in this case. */
 59:     *flg = PETSC_TRUE;
 60:   } else if (combiner == MPI_COMBINER_CONTIGUOUS) {
 61:     PetscCheck(nints == 1 && naddrs == 0 && ntypes == 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "Unexpected returns from MPI_Type_get_envelope()");
 62:     PetscCallMPI(MPI_Type_get_contents(a, 1, 0, 1, ints, addrs, types));
 63:     if (ints[0] == 1) { /* If a is created by MPI_Type_contiguous(1,..) */
 64:       PetscCall(MPIPetsc_Type_unwrap(types[0], atype, flg));
 65:       if (*flg) PetscCall(MPIPetsc_Type_free(&types[0]));
 66:       *flg = PETSC_TRUE;
 67:     } else {
 68:       PetscCall(MPIPetsc_Type_free(&types[0]));
 69:     }
 70:   }
 71:   PetscFunctionReturn(PETSC_SUCCESS);
 72: }

 74: PetscErrorCode MPIPetsc_Type_compare(MPI_Datatype a, MPI_Datatype b, PetscBool *match)
 75: {
 76:   MPI_Datatype atype, btype;
 77:   PetscMPIInt  aintcount, aaddrcount, atypecount, acombiner;
 78:   PetscMPIInt  bintcount, baddrcount, btypecount, bcombiner;
 79:   PetscBool    freeatype, freebtype;

 81:   PetscFunctionBegin;
 82:   if (a == b) { /* this is common when using MPI builtin datatypes */
 83:     *match = PETSC_TRUE;
 84:     PetscFunctionReturn(PETSC_SUCCESS);
 85:   }
 86:   PetscCall(MPIPetsc_Type_unwrap(a, &atype, &freeatype));
 87:   PetscCall(MPIPetsc_Type_unwrap(b, &btype, &freebtype));
 88:   *match = PETSC_FALSE;
 89:   if (atype == btype) {
 90:     *match = PETSC_TRUE;
 91:     goto free_types;
 92:   }
 93:   PetscCallMPI(MPI_Type_get_envelope(atype, &aintcount, &aaddrcount, &atypecount, &acombiner));
 94:   PetscCallMPI(MPI_Type_get_envelope(btype, &bintcount, &baddrcount, &btypecount, &bcombiner));
 95:   if (acombiner == bcombiner && aintcount == bintcount && aaddrcount == baddrcount && atypecount == btypecount && (aintcount > 0 || aaddrcount > 0 || atypecount > 0)) {
 96:     PetscMPIInt  *aints, *bints;
 97:     MPI_Aint     *aaddrs, *baddrs;
 98:     MPI_Datatype *atypes, *btypes;
 99:     PetscInt      i;
100:     PetscBool     same;
101:     PetscCall(PetscMalloc6(aintcount, &aints, bintcount, &bints, aaddrcount, &aaddrs, baddrcount, &baddrs, atypecount, &atypes, btypecount, &btypes));
102:     PetscCallMPI(MPI_Type_get_contents(atype, aintcount, aaddrcount, atypecount, aints, aaddrs, atypes));
103:     PetscCallMPI(MPI_Type_get_contents(btype, bintcount, baddrcount, btypecount, bints, baddrs, btypes));
104:     PetscCall(PetscArraycmp(aints, bints, aintcount, &same));
105:     if (same) {
106:       PetscCall(PetscArraycmp(aaddrs, baddrs, aaddrcount, &same));
107:       if (same) {
108:         /* Check for identity first */
109:         PetscCall(PetscArraycmp(atypes, btypes, atypecount, &same));
110:         if (!same) {
111:           /* If the atype or btype were not predefined data types, then the types returned from MPI_Type_get_contents
112:            * will merely be equivalent to the types used in the construction, so we must recursively compare. */
113:           for (i = 0; i < atypecount; i++) {
114:             PetscCall(MPIPetsc_Type_compare(atypes[i], btypes[i], &same));
115:             if (!same) break;
116:           }
117:         }
118:       }
119:     }
120:     for (i = 0; i < atypecount; i++) {
121:       PetscCall(MPIPetsc_Type_free(&atypes[i]));
122:       PetscCall(MPIPetsc_Type_free(&btypes[i]));
123:     }
124:     PetscCall(PetscFree6(aints, bints, aaddrs, baddrs, atypes, btypes));
125:     if (same) *match = PETSC_TRUE;
126:   }
127: free_types:
128:   if (freeatype) PetscCall(MPIPetsc_Type_free(&atype));
129:   if (freebtype) PetscCall(MPIPetsc_Type_free(&btype));
130:   PetscFunctionReturn(PETSC_SUCCESS);
131: }

133: /* Check whether a was created via MPI_Type_contiguous from b
134:  *
135:  */
136: PetscErrorCode MPIPetsc_Type_compare_contig(MPI_Datatype a, MPI_Datatype b, PetscInt *n)
137: {
138:   MPI_Datatype atype, btype;
139:   PetscMPIInt  aintcount, aaddrcount, atypecount, acombiner;
140:   PetscBool    freeatype, freebtype;

142:   PetscFunctionBegin;
143:   if (a == b) {
144:     *n = 1;
145:     PetscFunctionReturn(PETSC_SUCCESS);
146:   }
147:   *n = 0;
148:   PetscCall(MPIPetsc_Type_unwrap(a, &atype, &freeatype));
149:   PetscCall(MPIPetsc_Type_unwrap(b, &btype, &freebtype));
150:   PetscCallMPI(MPI_Type_get_envelope(atype, &aintcount, &aaddrcount, &atypecount, &acombiner));
151:   if (acombiner == MPI_COMBINER_CONTIGUOUS && aintcount >= 1) {
152:     PetscMPIInt  *aints;
153:     MPI_Aint     *aaddrs;
154:     MPI_Datatype *atypes;
155:     PetscInt      i;
156:     PetscBool     same;
157:     PetscCall(PetscMalloc3(aintcount, &aints, aaddrcount, &aaddrs, atypecount, &atypes));
158:     PetscCallMPI(MPI_Type_get_contents(atype, aintcount, aaddrcount, atypecount, aints, aaddrs, atypes));
159:     /* Check for identity first. */
160:     if (atypes[0] == btype) {
161:       *n = aints[0];
162:     } else {
163:       /* atypes[0] merely has to be equivalent to the type used to create atype. */
164:       PetscCall(MPIPetsc_Type_compare(atypes[0], btype, &same));
165:       if (same) *n = aints[0];
166:     }
167:     for (i = 0; i < atypecount; i++) PetscCall(MPIPetsc_Type_free(&atypes[i]));
168:     PetscCall(PetscFree3(aints, aaddrs, atypes));
169:   }

171:   if (freeatype) PetscCall(MPIPetsc_Type_free(&atype));
172:   if (freebtype) PetscCall(MPIPetsc_Type_free(&btype));
173:   PetscFunctionReturn(PETSC_SUCCESS);
174: }