Actual source code: sftype.c
petsc-3.14.3 2021-01-09
1: #include <petsc/private/sfimpl.h>
3: #if !defined(PETSC_HAVE_MPI_TYPE_GET_ENVELOPE)
4: #define MPI_Type_get_envelope(datatype,num_ints,num_addrs,num_dtypes,combiner) (*(num_ints)=0,*(num_addrs)=0,*(num_dtypes)=0,*(combiner)=0,1);SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP_SYS,"Need an MPI-2 implementation")
5: #define MPI_Type_get_contents(datatype,num_ints,num_addrs,num_dtypes,ints,addrs,dtypes) (*(ints)=0,*(addrs)=0,*(dtypes)=0,1);SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP_SYS,"Need an MPI-2 implementation")
6: #endif
7: #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. */
8: # define MPI_COMBINER_DUP 0
9: #endif
10: #if !defined(PETSC_HAVE_MPI_COMBINER_NAMED) && !defined(MPI_COMBINER_NAMED)
11: #define MPI_COMBINER_NAMED -2
12: #endif
13: #if !defined(PETSC_HAVE_MPI_COMBINER_CONTIGUOUS) && !defined(MPI_COMBINER_CONTIGUOUS) && MPI_VERSION < 2
14: # define MPI_COMBINER_CONTIGUOUS -1
15: #endif
17: static PetscErrorCode MPIPetsc_Type_free(MPI_Datatype *a)
18: {
19: PetscMPIInt nints,naddrs,ntypes,combiner;
23: MPI_Type_get_envelope(*a,&nints,&naddrs,&ntypes,&combiner);
25: if (combiner != MPI_COMBINER_NAMED) {
26: MPI_Type_free(a);
27: }
29: *a = MPI_DATATYPE_NULL;
30: return(0);
31: }
33: /*
34: Unwrap an MPI datatype recursively in case it is dupped or MPI_Type_contiguous(1,...)'ed from another type.
36: Input Arguments:
37: + a - the datatype to be unwrapped
39: Output Arguments:
40: + atype - the unwrapped datatype, which is either equal(=) to a or equivalent to a.
41: - flg - true if atype != a, which implies caller should MPIPetsc_Type_free(atype) after use. Note atype might be MPI builtin.
42: */
43: PetscErrorCode MPIPetsc_Type_unwrap(MPI_Datatype a,MPI_Datatype *atype,PetscBool *flg)
44: {
46: PetscMPIInt nints,naddrs,ntypes,combiner,ints[1];
47: MPI_Datatype types[1];
50: *flg = PETSC_FALSE;
51: *atype = a;
52: if (a == MPIU_INT || a == MPIU_REAL || a == MPIU_SCALAR) return(0);
53: MPI_Type_get_envelope(a,&nints,&naddrs,&ntypes,&combiner);
54: if (combiner == MPI_COMBINER_DUP) {
55: if (nints != 0 || naddrs != 0 || ntypes != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Unexpected returns from MPI_Type_get_envelope()");
56: MPI_Type_get_contents(a,0,0,1,ints,NULL,types);
57: /* Recursively unwrap dupped types. */
58: MPIPetsc_Type_unwrap(types[0],atype,flg);
59: if (*flg) {
60: /* If the recursive call returns a new type, then that means that atype[0] != types[0] and we're on the hook to
61: * free types[0]. Note that this case occurs if combiner(types[0]) is MPI_COMBINER_DUP, so we're safe to
62: * directly call MPI_Type_free rather than MPIPetsc_Type_free here. */
63: MPI_Type_free(&(types[0]));
64: }
65: /* In any case, it's up to the caller to free the returned type in this case. */
66: *flg = PETSC_TRUE;
67: } else if (combiner == MPI_COMBINER_CONTIGUOUS) {
68: if (nints != 1 || naddrs != 0 || ntypes != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Unexpected returns from MPI_Type_get_envelope()");
69: MPI_Type_get_contents(a,1,0,1,ints,NULL,types);
70: if (ints[0] == 1) { /* If a is created by MPI_Type_contiguous(1,..) */
71: MPIPetsc_Type_unwrap(types[0],atype,flg);
72: if (*flg) {MPIPetsc_Type_free(&(types[0]));}
73: *flg = PETSC_TRUE;
74: } else {
75: MPIPetsc_Type_free(&(types[0]));
76: }
77: }
78: return(0);
79: }
81: PetscErrorCode MPIPetsc_Type_compare(MPI_Datatype a,MPI_Datatype b,PetscBool *match)
82: {
84: MPI_Datatype atype,btype;
85: PetscMPIInt aintcount,aaddrcount,atypecount,acombiner;
86: PetscMPIInt bintcount,baddrcount,btypecount,bcombiner;
87: PetscBool freeatype, freebtype;
90: if (a == b) { /* this is common when using MPI builtin datatypes */
91: *match = PETSC_TRUE;
92: return(0);
93: }
94: MPIPetsc_Type_unwrap(a,&atype,&freeatype);
95: MPIPetsc_Type_unwrap(b,&btype,&freebtype);
96: *match = PETSC_FALSE;
97: if (atype == btype) {
98: *match = PETSC_TRUE;
99: goto free_types;
100: }
101: MPI_Type_get_envelope(atype,&aintcount,&aaddrcount,&atypecount,&acombiner);
102: MPI_Type_get_envelope(btype,&bintcount,&baddrcount,&btypecount,&bcombiner);
103: if (acombiner == bcombiner && aintcount == bintcount && aaddrcount == baddrcount && atypecount == btypecount && (aintcount > 0 || aaddrcount > 0 || atypecount > 0)) {
104: PetscMPIInt *aints,*bints;
105: MPI_Aint *aaddrs,*baddrs;
106: MPI_Datatype *atypes,*btypes;
107: PetscInt i;
108: PetscBool same;
109: PetscMalloc6(aintcount,&aints,bintcount,&bints,aaddrcount,&aaddrs,baddrcount,&baddrs,atypecount,&atypes,btypecount,&btypes);
110: MPI_Type_get_contents(atype,aintcount,aaddrcount,atypecount,aints,aaddrs,atypes);
111: MPI_Type_get_contents(btype,bintcount,baddrcount,btypecount,bints,baddrs,btypes);
112: PetscArraycmp(aints,bints,aintcount,&same);
113: if (same) {
114: PetscArraycmp(aaddrs,baddrs,aaddrcount,&same);
115: if (same) {
116: /* Check for identity first */
117: PetscArraycmp(atypes,btypes,atypecount,&same);
118: if (!same) {
119: /* If the atype or btype were not predefined data types, then the types returned from MPI_Type_get_contents
120: * will merely be equivalent to the types used in the construction, so we must recursively compare. */
121: for (i=0; i<atypecount; i++) {
122: MPIPetsc_Type_compare(atypes[i],btypes[i],&same);
123: if (!same) break;
124: }
125: }
126: }
127: }
128: for (i=0; i<atypecount; i++) {
129: MPIPetsc_Type_free(&(atypes[i]));
130: MPIPetsc_Type_free(&(btypes[i]));
131: }
132: PetscFree6(aints,bints,aaddrs,baddrs,atypes,btypes);
133: if (same) *match = PETSC_TRUE;
134: }
135: free_types:
136: if (freeatype) {
137: MPIPetsc_Type_free(&atype);
138: }
139: if (freebtype) {
140: MPIPetsc_Type_free(&btype);
141: }
142: return(0);
143: }
145: /* Check whether a was created via MPI_Type_contiguous from b
146: *
147: */
148: PetscErrorCode MPIPetsc_Type_compare_contig(MPI_Datatype a,MPI_Datatype b,PetscInt *n)
149: {
151: MPI_Datatype atype,btype;
152: PetscMPIInt aintcount,aaddrcount,atypecount,acombiner;
153: PetscBool freeatype,freebtype;
156: if (a == b) {*n = 1; return(0);}
157: *n = 0;
158: MPIPetsc_Type_unwrap(a,&atype,&freeatype);
159: MPIPetsc_Type_unwrap(b,&btype,&freebtype);
160: MPI_Type_get_envelope(atype,&aintcount,&aaddrcount,&atypecount,&acombiner);
161: if (acombiner == MPI_COMBINER_CONTIGUOUS && aintcount >= 1) {
162: PetscMPIInt *aints;
163: MPI_Aint *aaddrs;
164: MPI_Datatype *atypes;
165: PetscInt i;
166: PetscBool same;
167: PetscMalloc3(aintcount,&aints,aaddrcount,&aaddrs,atypecount,&atypes);
168: MPI_Type_get_contents(atype,aintcount,aaddrcount,atypecount,aints,aaddrs,atypes);
169: /* Check for identity first. */
170: if (atypes[0] == btype) {
171: *n = aints[0];
172: } else {
173: /* atypes[0] merely has to be equivalent to the type used to create atype. */
174: MPIPetsc_Type_compare(atypes[0],btype,&same);
175: if (same) *n = aints[0];
176: }
177: for (i=0; i<atypecount; i++) {
178: MPIPetsc_Type_free(&(atypes[i]));
179: }
180: PetscFree3(aints,aaddrs,atypes);
181: }
183: if (freeatype) {
184: MPIPetsc_Type_free(&atype);
185: }
186: if (freebtype) {
187: MPIPetsc_Type_free(&btype);
188: }
189: return(0);
190: }