Actual source code: garbage.c
  1: #include <petsc/private/garbagecollector.h>
  3: /* Fetches garbage hashmap from communicator */
  4: static PetscErrorCode GarbageGetHMap_Private(MPI_Comm comm, PetscGarbage *garbage)
  5: {
  6:   PetscMPIInt  flag;
  7:   PetscHMapObj garbage_map;
  9:   PetscFunctionBegin;
 10:   PetscCallMPI(MPI_Comm_get_attr(comm, Petsc_Garbage_HMap_keyval, garbage, &flag));
 11:   if (!flag) {
 12:     /* No garbage,create one */
 13:     PetscCall(PetscHMapObjCreate(&garbage_map));
 14:     garbage->map = garbage_map;
 15:     PetscCallMPI(MPI_Comm_set_attr(comm, Petsc_Garbage_HMap_keyval, garbage->ptr));
 16:   }
 17:   PetscFunctionReturn(PETSC_SUCCESS);
 18: }
 20: /*@C
 21:     PetscObjectDelayedDestroy - Adds an object to a data structure for
 22:     later destruction.
 24:     Not Collective
 26:     Input Parameter:
 27: .   obj - object to be destroyed
 29:     Level: developer
 31:     Notes:
 32:     Analogue to `PetscObjectDestroy()` for use in managed languages.
 34:     A PETSc object is given a creation index at initialisation based on
 35:     the communicator it was created on and the order in which it is
 36:     created. When this function is passed a PETSc object, a pointer to
 37:     the object is stashed on a garbage dictionary (`PetscHMapObj`) which is
 38:     keyed by its creation index.
 40:     Objects stashed on this garbage dictionary can later be destroyed
 41:     with a call to `PetscGarbageCleanup()`.
 43:     This function is intended for use with managed languages such as
 44:     Python or Julia, which may not destroy objects in a deterministic
 45:     order.
 47: .seealso: `PetscGarbageCleanup()`, `PetscObjectDestroy()`
 48: @*/
 49: PetscErrorCode PetscObjectDelayedDestroy(PetscObject *obj)
 50: {
 51:   MPI_Comm     petsc_comm;
 52:   PetscInt     count;
 53:   PetscGarbage garbage;
 55:   PetscFunctionBegin;
 57:   /* Don't stash NULL pointers */
 58:   if (*obj != NULL) {
 59:     /* Elaborate check for getting non-cyclic reference counts */
 60:     if (!(*obj)->non_cyclic_references) {
 61:       count = --(*obj)->refct;
 62:     } else {
 63:       PetscCall((*obj)->non_cyclic_references(*obj, &count));
 64:       --count;
 65:       --(*obj)->refct;
 66:     }
 67:     /* Only stash if the (non-cyclic) reference count hits 0 */
 68:     if (count == 0) {
 69:       (*obj)->refct = 1;
 70:       PetscCall(PetscObjectGetComm(*obj, &petsc_comm));
 71:       PetscCall(GarbageGetHMap_Private(petsc_comm, &garbage));
 72:       PetscCall(PetscHMapObjSet(garbage.map, (*obj)->cidx, *obj));
 73:     }
 74:   }
 75:   *obj = NULL;
 76:   PetscFunctionReturn(PETSC_SUCCESS);
 77: }
 79: /* Performs the intersection of 2 sorted arrays seta and setb of lengths
 80:    lena and lenb respectively,returning the result in seta and lena
 81:    This is an O(n) operation */
 82: static PetscErrorCode GarbageKeySortedIntersect_Private(PetscInt64 seta[], PetscInt *lena, PetscInt64 setb[], PetscInt lenb)
 83: {
 84:   /* The arrays seta and setb MUST be sorted! */
 85:   PetscInt ii, jj = 0, counter = 0;
 87:   PetscFunctionBegin;
 88:   if (PetscDefined(USE_DEBUG)) {
 89:     PetscBool sorted = PETSC_FALSE;
 90:     /* In debug mode check whether the array are sorted */
 91:     PetscCall(PetscSortedInt64(*lena, seta, &sorted));
 92:     PetscCheck(sorted, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Provided array in argument 1 is not sorted");
 93:     PetscCall(PetscSortedInt64(lenb, setb, &sorted));
 94:     PetscCheck(sorted, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Provided array in argument 3 is not sorted");
 95:   }
 96:   for (ii = 0; ii < *lena; ii++) {
 97:     while (jj < lenb && seta[ii] > setb[jj]) jj++;
 98:     if (jj >= lenb) break;
 99:     if (seta[ii] == setb[jj]) {
100:       seta[counter] = seta[ii];
101:       counter++;
102:     }
103:   }
105:   *lena = counter;
106:   PetscFunctionReturn(PETSC_SUCCESS);
107: }
109: /* Wrapper to create MPI reduce operator for set intersection */
110: void PetscGarbageKeySortedIntersect(void *inset, void *inoutset, PetscMPIInt *length, MPI_Datatype *dtype)
111: {
112:   PetscInt64 *seta, *setb;
114:   seta = (PetscInt64 *)inoutset;
115:   setb = (PetscInt64 *)inset;
117:   PetscCallAbort(PETSC_COMM_SELF, GarbageKeySortedIntersect_Private(&seta[1], (PetscInt *)&seta[0], &setb[1], (PetscInt)setb[0]));
118: }
120: /* Performs a collective allreduce intersection of one array per rank */
121: PetscErrorCode GarbageKeyAllReduceIntersect_Private(MPI_Comm comm, PetscInt64 *set, PetscInt *entries)
122: {
123:   PetscInt     ii, max_entries;
124:   PetscInt64  *sendset, *recvset;
125:   MPI_Datatype keyset_type;
127:   PetscFunctionBegin;
128:   /* Sort keys first for use with `GarbageKeySortedIntersect_Private()`*/
129:   PetscCall(PetscSortInt64(*entries, set));
131:   /* Get the maximum size of all key sets */
132:   PetscCall(MPIU_Allreduce(entries, &max_entries, 1, MPIU_INT, MPI_MAX, comm));
133:   PetscCall(PetscMalloc1(max_entries + 1, &sendset));
134:   PetscCall(PetscMalloc1(max_entries + 1, &recvset));
135:   sendset[0] = (PetscInt64)*entries;
136:   for (ii = 1; ii < *entries + 1; ii++) sendset[ii] = set[ii - 1];
138:   /* Create a custom data type to hold the set */
139:   PetscCallMPI(MPI_Type_contiguous(max_entries + 1, MPIU_INT64, &keyset_type));
140:   /* PetscCallMPI(MPI_Type_set_name(keyset_type,"PETSc garbage key set type")); */
141:   PetscCallMPI(MPI_Type_commit(&keyset_type));
143:   /* Perform custom intersect reduce operation over sets */
144:   PetscCallMPI(MPI_Allreduce(sendset, recvset, 1, keyset_type, Petsc_Garbage_SetIntersectOp, comm));
146:   PetscCallMPI(MPI_Type_free(&keyset_type));
148:   *entries = (PetscInt)recvset[0];
149:   for (ii = 0; ii < *entries; ii++) set[ii] = recvset[ii + 1];
151:   PetscCall(PetscFree(sendset));
152:   PetscCall(PetscFree(recvset));
153:   PetscFunctionReturn(PETSC_SUCCESS);
154: }
156: /*@C
157:     PetscGarbageCleanup - Destroys objects placed in the garbage by
158:     `PetscObjectDelayedDestroy()`.
160:     Collective
162:     Input Parameter:
163: .   comm      - MPI communicator over which to perform collective cleanup
165:     Level: developer
167:     Notes:
168:     Implements a collective garbage collection.
169:     A per- MPI communicator garbage dictionary is created to store
170:     references to objects destroyed using `PetscObjectDelayedDestroy()`.
171:     Objects that appear in this dictionary on all MPI processes can be destroyed
172:     by calling `PetscGarbageCleanup()`.
174:     This is done as follows:
175:     1.  Keys of the garbage dictionary, which correspond to the creation
176:         indices of the objects stashed, are sorted.
177:     2.  A collective intersection of dictionary keys is performed by all
178:         ranks in the communicator.
179:     3.  The intersection is broadcast back to all ranks in the
180:         communicator.
181:     4.  The objects on the dictionary are collectively destroyed in
182:         creation index order using a call to PetscObjectDestroy().
184:     This function is intended for use with managed languages such as
185:     Python or Julia, which may not destroy objects in a deterministic
186:     order.
188: .seealso: PetscObjectDelayedDestroy()
189: @*/
190: PetscErrorCode PetscGarbageCleanup(MPI_Comm comm)
191: {
192:   PetscInt     ii, entries, offset;
193:   PetscInt64  *keys;
194:   PetscObject  obj;
195:   PetscGarbage garbage;
197:   PetscFunctionBegin;
198:   /* Duplicate comm to prevent it being cleaned up by PetscObjectDestroy() */
199:   PetscCall(PetscCommDuplicate(comm, &comm, NULL));
201:   /* Grab garbage from comm and remove it
202:    this avoids calling PetscCommDestroy() and endlessly recursing */
203:   PetscCall(GarbageGetHMap_Private(comm, &garbage));
204:   PetscCallMPI(MPI_Comm_delete_attr(comm, Petsc_Garbage_HMap_keyval));
206:   /* Get keys from garbage hash map */
207:   PetscCall(PetscHMapObjGetSize(garbage.map, &entries));
208:   PetscCall(PetscMalloc1(entries, &keys));
209:   offset = 0;
210:   PetscCall(PetscHMapObjGetKeys(garbage.map, &offset, keys));
212:   /* Gather and intersect */
213:   PetscCall(GarbageKeyAllReduceIntersect_Private(comm, keys, &entries));
215:   /* Collectively destroy objects objects that appear in garbage in
216:      creation index order */
217:   for (ii = 0; ii < entries; ii++) {
218:     PetscCall(PetscHMapObjGet(garbage.map, keys[ii], &obj));
219:     PetscCall(PetscObjectDestroy(&obj));
220:     PetscCall(PetscFree(obj));
221:     PetscCall(PetscHMapObjDel(garbage.map, keys[ii]));
222:   }
223:   PetscCall(PetscFree(keys));
225:   /* Put garbage back */
226:   PetscCallMPI(MPI_Comm_set_attr(comm, Petsc_Garbage_HMap_keyval, garbage.ptr));
227:   PetscCall(PetscCommDestroy(&comm));
228:   PetscFunctionReturn(PETSC_SUCCESS);
229: }
231: /* Utility function for printing the contents of the garbage on a given comm */
232: PetscErrorCode PetscGarbageView(MPI_Comm comm, PetscViewer viewer)
233: {
234:   char         text[64];
235:   PetscInt     ii, entries, offset;
236:   PetscInt64  *keys;
237:   PetscObject  obj;
238:   PetscGarbage garbage;
239:   PetscMPIInt  rank;
241:   PetscFunctionBegin;
242:   PetscCall(PetscPrintf(comm, "PETSc garbage on "));
243:   if (comm == PETSC_COMM_WORLD) {
244:     PetscCall(PetscPrintf(comm, "PETSC_COMM_WORLD\n"));
245:   } else if (comm == PETSC_COMM_SELF) {
246:     PetscCall(PetscPrintf(comm, "PETSC_COMM_SELF\n"));
247:   } else {
248:     PetscCall(PetscPrintf(comm, "UNKNOWN_COMM\n"));
249:   }
250:   PetscCall(PetscCommDuplicate(comm, &comm, NULL));
251:   PetscCall(GarbageGetHMap_Private(comm, &garbage));
253:   /* Get keys from garbage hash map and sort */
254:   PetscCall(PetscHMapObjGetSize(garbage.map, &entries));
255:   PetscCall(PetscMalloc1(entries, &keys));
256:   offset = 0;
257:   PetscCall(PetscHMapObjGetKeys(garbage.map, &offset, keys));
259:   /* Pretty print entries in a table */
260:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
261:   PetscCall(PetscSynchronizedPrintf(comm, "Rank %i:: ", rank));
262:   PetscCall(PetscFormatConvert("Total entries: %D\n", text));
263:   PetscCall(PetscSynchronizedPrintf(comm, text, entries));
264:   if (entries) {
265:     PetscCall(PetscSynchronizedPrintf(comm, "| Key   | Type                   | Name                             | Object ID |\n"));
266:     PetscCall(PetscSynchronizedPrintf(comm, "|-------|------------------------|----------------------------------|-----------|\n"));
267:   }
268:   for (ii = 0; ii < entries; ii++) {
269:     PetscCall(PetscHMapObjGet(garbage.map, keys[ii], &obj));
270:     PetscCall(PetscFormatConvert("| %5" PetscInt64_FMT " | %-22s | %-32s | %6D    |\n", text));
271:     PetscCall(PetscSynchronizedPrintf(comm, text, keys[ii], obj->class_name, obj->description, obj->id));
272:   }
273:   PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT));
275:   PetscCall(PetscFree(keys));
276:   PetscCall(PetscCommDestroy(&comm));
277:   PetscFunctionReturn(PETSC_SUCCESS);
278: }