xref: /petsc/src/sys/objects/garbage.c (revision 51b144c619aff302b570817d6f78637b8418d403) !
162e5d2d2SJDBetteridge #include <petsc/private/garbagecollector.h>
262e5d2d2SJDBetteridge 
362e5d2d2SJDBetteridge /* Fetches garbage hashmap from communicator */
GarbageGetHMap_Private(MPI_Comm comm,PetscGarbage * garbage)462e5d2d2SJDBetteridge static PetscErrorCode GarbageGetHMap_Private(MPI_Comm comm, PetscGarbage *garbage)
562e5d2d2SJDBetteridge {
6*b8b5be36SMartin Diehl   PetscMPIInt  iflg;
762e5d2d2SJDBetteridge   PetscHMapObj garbage_map;
862e5d2d2SJDBetteridge 
962e5d2d2SJDBetteridge   PetscFunctionBegin;
10*b8b5be36SMartin Diehl   PetscCallMPI(MPI_Comm_get_attr(comm, Petsc_Garbage_HMap_keyval, garbage, &iflg));
11*b8b5be36SMartin Diehl   if (!iflg) {
1262e5d2d2SJDBetteridge     /* No garbage,create one */
1362e5d2d2SJDBetteridge     PetscCall(PetscHMapObjCreate(&garbage_map));
1462e5d2d2SJDBetteridge     garbage->map = garbage_map;
1562e5d2d2SJDBetteridge     PetscCallMPI(MPI_Comm_set_attr(comm, Petsc_Garbage_HMap_keyval, garbage->ptr));
1662e5d2d2SJDBetteridge   }
173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1862e5d2d2SJDBetteridge }
1962e5d2d2SJDBetteridge 
2062e5d2d2SJDBetteridge /*@C
2162e5d2d2SJDBetteridge   PetscObjectDelayedDestroy - Adds an object to a data structure for
2262e5d2d2SJDBetteridge   later destruction.
2362e5d2d2SJDBetteridge 
2462e5d2d2SJDBetteridge   Not Collective
2562e5d2d2SJDBetteridge 
262fe279fdSBarry Smith   Input Parameter:
2762e5d2d2SJDBetteridge . obj - object to be destroyed
2862e5d2d2SJDBetteridge 
2921532e8aSBarry Smith   Level: developer
3021532e8aSBarry Smith 
3162e5d2d2SJDBetteridge   Notes:
3262e5d2d2SJDBetteridge   Analogue to `PetscObjectDestroy()` for use in managed languages.
3362e5d2d2SJDBetteridge 
3462e5d2d2SJDBetteridge   A PETSc object is given a creation index at initialisation based on
3562e5d2d2SJDBetteridge   the communicator it was created on and the order in which it is
3662e5d2d2SJDBetteridge   created. When this function is passed a PETSc object, a pointer to
3721532e8aSBarry Smith   the object is stashed on a garbage dictionary (`PetscHMapObj`) which is
3862e5d2d2SJDBetteridge   keyed by its creation index.
3962e5d2d2SJDBetteridge 
4062e5d2d2SJDBetteridge   Objects stashed on this garbage dictionary can later be destroyed
4162e5d2d2SJDBetteridge   with a call to `PetscGarbageCleanup()`.
4262e5d2d2SJDBetteridge 
4362e5d2d2SJDBetteridge   This function is intended for use with managed languages such as
4462e5d2d2SJDBetteridge   Python or Julia, which may not destroy objects in a deterministic
4562e5d2d2SJDBetteridge   order.
4662e5d2d2SJDBetteridge 
47b0ed7cdaSConnor Ward   Serial objects (that have a communicator with size 1) are destroyed
48b0ed7cdaSConnor Ward   eagerly since deadlocks cannot occur.
49b0ed7cdaSConnor Ward 
5021532e8aSBarry Smith .seealso: `PetscGarbageCleanup()`, `PetscObjectDestroy()`
5162e5d2d2SJDBetteridge @*/
PetscObjectDelayedDestroy(PetscObject * obj)5262e5d2d2SJDBetteridge PetscErrorCode PetscObjectDelayedDestroy(PetscObject *obj)
5362e5d2d2SJDBetteridge {
54b0ed7cdaSConnor Ward   MPI_Comm     comm;
55b0ed7cdaSConnor Ward   PetscMPIInt  size;
5662e5d2d2SJDBetteridge   PetscInt     count;
5762e5d2d2SJDBetteridge   PetscGarbage garbage;
5862e5d2d2SJDBetteridge 
5962e5d2d2SJDBetteridge   PetscFunctionBegin;
604f572ea9SToby Isaac   PetscAssertPointer(obj, 1);
6162e5d2d2SJDBetteridge   /* Don't stash NULL pointers */
6262e5d2d2SJDBetteridge   if (*obj != NULL) {
6362e5d2d2SJDBetteridge     /* Elaborate check for getting non-cyclic reference counts */
6462e5d2d2SJDBetteridge     if (!(*obj)->non_cyclic_references) {
6562e5d2d2SJDBetteridge       count = --(*obj)->refct;
6662e5d2d2SJDBetteridge     } else {
6762e5d2d2SJDBetteridge       PetscCall((*obj)->non_cyclic_references(*obj, &count));
6862e5d2d2SJDBetteridge       --count;
6962e5d2d2SJDBetteridge       --(*obj)->refct;
7062e5d2d2SJDBetteridge     }
7162e5d2d2SJDBetteridge     /* Only stash if the (non-cyclic) reference count hits 0 */
7262e5d2d2SJDBetteridge     if (count == 0) {
7362e5d2d2SJDBetteridge       (*obj)->refct = 1;
74b0ed7cdaSConnor Ward       PetscCall(PetscObjectGetComm(*obj, &comm));
75b0ed7cdaSConnor Ward       PetscCallMPI(MPI_Comm_size(comm, &size));
76b0ed7cdaSConnor Ward       /* Eagerly destroy serial objects */
77b0ed7cdaSConnor Ward       if (size == 1) {
78b0ed7cdaSConnor Ward         PetscCall(PetscObjectDestroy(obj));
79b0ed7cdaSConnor Ward       } else {
80b0ed7cdaSConnor Ward         PetscCall(GarbageGetHMap_Private(comm, &garbage));
8162e5d2d2SJDBetteridge         PetscCall(PetscHMapObjSet(garbage.map, (*obj)->cidx, *obj));
8262e5d2d2SJDBetteridge       }
8362e5d2d2SJDBetteridge     }
84b0ed7cdaSConnor Ward   }
8562e5d2d2SJDBetteridge   *obj = NULL;
863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8762e5d2d2SJDBetteridge }
8862e5d2d2SJDBetteridge 
8962e5d2d2SJDBetteridge /* Performs the intersection of 2 sorted arrays seta and setb of lengths
9062e5d2d2SJDBetteridge    lena and lenb respectively,returning the result in seta and lena
9162e5d2d2SJDBetteridge    This is an O(n) operation */
GarbageKeySortedIntersect_Private(PetscInt64 seta[],PetscInt * lena,PetscInt64 setb[],PetscInt lenb)9262e5d2d2SJDBetteridge static PetscErrorCode GarbageKeySortedIntersect_Private(PetscInt64 seta[], PetscInt *lena, PetscInt64 setb[], PetscInt lenb)
9362e5d2d2SJDBetteridge {
9462e5d2d2SJDBetteridge   /* The arrays seta and setb MUST be sorted! */
95b0271995SJDBetteridge   PetscInt ii, jj = 0, counter = 0;
9662e5d2d2SJDBetteridge 
9762e5d2d2SJDBetteridge   PetscFunctionBegin;
9862e5d2d2SJDBetteridge   if (PetscDefined(USE_DEBUG)) {
9962e5d2d2SJDBetteridge     PetscBool sorted = PETSC_FALSE;
10062e5d2d2SJDBetteridge     /* In debug mode check whether the array are sorted */
10162e5d2d2SJDBetteridge     PetscCall(PetscSortedInt64(*lena, seta, &sorted));
10279528a3fSJDBetteridge     PetscCheck(sorted, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Provided array in argument 1 is not sorted");
10362e5d2d2SJDBetteridge     PetscCall(PetscSortedInt64(lenb, setb, &sorted));
10479528a3fSJDBetteridge     PetscCheck(sorted, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Provided array in argument 3 is not sorted");
10562e5d2d2SJDBetteridge   }
10662e5d2d2SJDBetteridge   for (ii = 0; ii < *lena; ii++) {
107aa624791SPierre Jolivet     while (jj < lenb && seta[ii] > setb[jj]) jj++;
108b0271995SJDBetteridge     if (jj >= lenb) break;
109b0271995SJDBetteridge     if (seta[ii] == setb[jj]) {
11062e5d2d2SJDBetteridge       seta[counter] = seta[ii];
11162e5d2d2SJDBetteridge       counter++;
11262e5d2d2SJDBetteridge     }
11362e5d2d2SJDBetteridge   }
114b0271995SJDBetteridge 
11562e5d2d2SJDBetteridge   *lena = counter;
1163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11762e5d2d2SJDBetteridge }
11862e5d2d2SJDBetteridge 
11962e5d2d2SJDBetteridge /* Wrapper to create MPI reduce operator for set intersection */
PetscGarbageKeySortedIntersect(void * inset,void * inoutset,PetscMPIInt * length,MPI_Datatype * dtype)12062e5d2d2SJDBetteridge void PetscGarbageKeySortedIntersect(void *inset, void *inoutset, PetscMPIInt *length, MPI_Datatype *dtype)
12162e5d2d2SJDBetteridge {
12262e5d2d2SJDBetteridge   PetscInt64 *seta, *setb;
1234f2d7745SStefano Zampini   PetscInt    lena = 0, lenb = 0;
12462e5d2d2SJDBetteridge 
12562e5d2d2SJDBetteridge   seta = (PetscInt64 *)inoutset;
12662e5d2d2SJDBetteridge   setb = (PetscInt64 *)inset;
12762e5d2d2SJDBetteridge 
1284f2d7745SStefano Zampini   PetscCallAbort(PETSC_COMM_SELF, PetscIntCast(seta[0], &lena));
1294f2d7745SStefano Zampini   PetscCallAbort(PETSC_COMM_SELF, PetscIntCast(setb[0], &lenb));
1304f2d7745SStefano Zampini   PetscCallAbort(PETSC_COMM_SELF, GarbageKeySortedIntersect_Private(seta + 1, &lena, setb + 1, lenb));
1314f2d7745SStefano Zampini   seta[0] = lena;
13262e5d2d2SJDBetteridge }
13362e5d2d2SJDBetteridge 
13462e5d2d2SJDBetteridge /* Performs a collective allreduce intersection of one array per rank */
GarbageKeyAllReduceIntersect_Private(MPI_Comm comm,PetscInt64 * set,PetscInt * entries)1358e89d99cSJDBetteridge PetscErrorCode GarbageKeyAllReduceIntersect_Private(MPI_Comm comm, PetscInt64 *set, PetscInt *entries)
13662e5d2d2SJDBetteridge {
13762e5d2d2SJDBetteridge   PetscInt     ii, max_entries;
13862e5d2d2SJDBetteridge   PetscInt64  *sendset, *recvset;
13962e5d2d2SJDBetteridge   MPI_Datatype keyset_type;
1406497c311SBarry Smith   PetscMPIInt  imax_entries;
14162e5d2d2SJDBetteridge 
14262e5d2d2SJDBetteridge   PetscFunctionBegin;
14362e5d2d2SJDBetteridge   /* Sort keys first for use with `GarbageKeySortedIntersect_Private()`*/
14462e5d2d2SJDBetteridge   PetscCall(PetscSortInt64(*entries, set));
14562e5d2d2SJDBetteridge 
14662e5d2d2SJDBetteridge   /* Get the maximum size of all key sets */
147462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(entries, &max_entries, 1, MPIU_INT, MPI_MAX, comm));
14862e5d2d2SJDBetteridge   PetscCall(PetscMalloc1(max_entries + 1, &sendset));
14962e5d2d2SJDBetteridge   PetscCall(PetscMalloc1(max_entries + 1, &recvset));
1504f2d7745SStefano Zampini   sendset[0] = *entries;
15162e5d2d2SJDBetteridge   for (ii = 1; ii < *entries + 1; ii++) sendset[ii] = set[ii - 1];
15262e5d2d2SJDBetteridge 
15362e5d2d2SJDBetteridge   /* Create a custom data type to hold the set */
1546497c311SBarry Smith   PetscCall(PetscMPIIntCast(max_entries, &imax_entries));
1556497c311SBarry Smith   PetscCallMPI(MPI_Type_contiguous(imax_entries + 1, MPIU_INT64, &keyset_type));
15662e5d2d2SJDBetteridge   /* PetscCallMPI(MPI_Type_set_name(keyset_type,"PETSc garbage key set type")); */
15762e5d2d2SJDBetteridge   PetscCallMPI(MPI_Type_commit(&keyset_type));
15862e5d2d2SJDBetteridge 
15962e5d2d2SJDBetteridge   /* Perform custom intersect reduce operation over sets */
1606a210b70SBarry Smith   PetscCallMPI(MPIU_Allreduce(sendset, recvset, 1, keyset_type, Petsc_Garbage_SetIntersectOp, comm));
16162e5d2d2SJDBetteridge 
16262e5d2d2SJDBetteridge   PetscCallMPI(MPI_Type_free(&keyset_type));
16362e5d2d2SJDBetteridge 
1644f2d7745SStefano Zampini   PetscCall(PetscIntCast(recvset[0], entries));
16562e5d2d2SJDBetteridge   for (ii = 0; ii < *entries; ii++) set[ii] = recvset[ii + 1];
16662e5d2d2SJDBetteridge 
16762e5d2d2SJDBetteridge   PetscCall(PetscFree(sendset));
16862e5d2d2SJDBetteridge   PetscCall(PetscFree(recvset));
1693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17062e5d2d2SJDBetteridge }
17162e5d2d2SJDBetteridge 
17262e5d2d2SJDBetteridge /*@C
17362e5d2d2SJDBetteridge   PetscGarbageCleanup - Destroys objects placed in the garbage by
17421532e8aSBarry Smith   `PetscObjectDelayedDestroy()`.
17562e5d2d2SJDBetteridge 
17662e5d2d2SJDBetteridge   Collective
17762e5d2d2SJDBetteridge 
1782fe279fdSBarry Smith   Input Parameter:
17921532e8aSBarry Smith . comm - MPI communicator over which to perform collective cleanup
18021532e8aSBarry Smith 
18121532e8aSBarry Smith   Level: developer
18262e5d2d2SJDBetteridge 
18362e5d2d2SJDBetteridge   Notes:
18462e5d2d2SJDBetteridge   Implements a collective garbage collection.
18562e5d2d2SJDBetteridge   A per- MPI communicator garbage dictionary is created to store
18621532e8aSBarry Smith   references to objects destroyed using `PetscObjectDelayedDestroy()`.
18721532e8aSBarry Smith   Objects that appear in this dictionary on all MPI processes can be destroyed
18821532e8aSBarry Smith   by calling `PetscGarbageCleanup()`.
18962e5d2d2SJDBetteridge 
19010450e9eSJacob Faibussowitsch   This is done as follows\:
19162e5d2d2SJDBetteridge   1.  Keys of the garbage dictionary, which correspond to the creation
19262e5d2d2SJDBetteridge   indices of the objects stashed, are sorted.
19362e5d2d2SJDBetteridge   2.  A collective intersection of dictionary keys is performed by all
19462e5d2d2SJDBetteridge   ranks in the communicator.
19562e5d2d2SJDBetteridge   3.  The intersection is broadcast back to all ranks in the
19662e5d2d2SJDBetteridge   communicator.
19762e5d2d2SJDBetteridge   4.  The objects on the dictionary are collectively destroyed in
19862e5d2d2SJDBetteridge   creation index order using a call to PetscObjectDestroy().
19962e5d2d2SJDBetteridge 
20062e5d2d2SJDBetteridge   This function is intended for use with managed languages such as
20162e5d2d2SJDBetteridge   Python or Julia, which may not destroy objects in a deterministic
20262e5d2d2SJDBetteridge   order.
20362e5d2d2SJDBetteridge 
204aec76313SJacob Faibussowitsch .seealso: `PetscObjectDelayedDestroy()`
20562e5d2d2SJDBetteridge @*/
PetscGarbageCleanup(MPI_Comm comm)20662e5d2d2SJDBetteridge PetscErrorCode PetscGarbageCleanup(MPI_Comm comm)
20762e5d2d2SJDBetteridge {
20862e5d2d2SJDBetteridge   PetscInt     ii, entries, offset;
20962e5d2d2SJDBetteridge   PetscInt64  *keys;
21062e5d2d2SJDBetteridge   PetscObject  obj;
21162e5d2d2SJDBetteridge   PetscGarbage garbage;
21262e5d2d2SJDBetteridge 
21362e5d2d2SJDBetteridge   PetscFunctionBegin;
21462e5d2d2SJDBetteridge   /* Duplicate comm to prevent it being cleaned up by PetscObjectDestroy() */
21562e5d2d2SJDBetteridge   PetscCall(PetscCommDuplicate(comm, &comm, NULL));
21662e5d2d2SJDBetteridge 
21762e5d2d2SJDBetteridge   /* Grab garbage from comm and remove it
21862e5d2d2SJDBetteridge    this avoids calling PetscCommDestroy() and endlessly recursing */
21962e5d2d2SJDBetteridge   PetscCall(GarbageGetHMap_Private(comm, &garbage));
22062e5d2d2SJDBetteridge   PetscCallMPI(MPI_Comm_delete_attr(comm, Petsc_Garbage_HMap_keyval));
22162e5d2d2SJDBetteridge 
22262e5d2d2SJDBetteridge   /* Get keys from garbage hash map */
22362e5d2d2SJDBetteridge   PetscCall(PetscHMapObjGetSize(garbage.map, &entries));
22462e5d2d2SJDBetteridge   PetscCall(PetscMalloc1(entries, &keys));
22562e5d2d2SJDBetteridge   offset = 0;
22662e5d2d2SJDBetteridge   PetscCall(PetscHMapObjGetKeys(garbage.map, &offset, keys));
22762e5d2d2SJDBetteridge 
22862e5d2d2SJDBetteridge   /* Gather and intersect */
22962e5d2d2SJDBetteridge   PetscCall(GarbageKeyAllReduceIntersect_Private(comm, keys, &entries));
23062e5d2d2SJDBetteridge 
23115229ffcSPierre Jolivet   /* Collectively destroy objects that appear in garbage in
23262e5d2d2SJDBetteridge      creation index order */
23362e5d2d2SJDBetteridge   for (ii = 0; ii < entries; ii++) {
23462e5d2d2SJDBetteridge     PetscCall(PetscHMapObjGet(garbage.map, keys[ii], &obj));
23562e5d2d2SJDBetteridge     PetscCall(PetscObjectDestroy(&obj));
23662e5d2d2SJDBetteridge     PetscCall(PetscFree(obj));
23762e5d2d2SJDBetteridge     PetscCall(PetscHMapObjDel(garbage.map, keys[ii]));
23862e5d2d2SJDBetteridge   }
23962e5d2d2SJDBetteridge   PetscCall(PetscFree(keys));
24062e5d2d2SJDBetteridge 
24162e5d2d2SJDBetteridge   /* Put garbage back */
24262e5d2d2SJDBetteridge   PetscCallMPI(MPI_Comm_set_attr(comm, Petsc_Garbage_HMap_keyval, garbage.ptr));
24362e5d2d2SJDBetteridge   PetscCall(PetscCommDestroy(&comm));
2443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
24562e5d2d2SJDBetteridge }
24662e5d2d2SJDBetteridge 
24762e5d2d2SJDBetteridge /* Utility function for printing the contents of the garbage on a given comm */
PetscGarbageView(MPI_Comm comm,PetscViewer viewer)24862e5d2d2SJDBetteridge PetscErrorCode PetscGarbageView(MPI_Comm comm, PetscViewer viewer)
24962e5d2d2SJDBetteridge {
25062e5d2d2SJDBetteridge   char         text[64];
25162e5d2d2SJDBetteridge   PetscInt     ii, entries, offset;
25262e5d2d2SJDBetteridge   PetscInt64  *keys;
25362e5d2d2SJDBetteridge   PetscObject  obj;
25462e5d2d2SJDBetteridge   PetscGarbage garbage;
25562e5d2d2SJDBetteridge   PetscMPIInt  rank;
25662e5d2d2SJDBetteridge 
25762e5d2d2SJDBetteridge   PetscFunctionBegin;
25862e5d2d2SJDBetteridge   PetscCall(PetscPrintf(comm, "PETSc garbage on "));
25962e5d2d2SJDBetteridge   if (comm == PETSC_COMM_WORLD) {
26062e5d2d2SJDBetteridge     PetscCall(PetscPrintf(comm, "PETSC_COMM_WORLD\n"));
26162e5d2d2SJDBetteridge   } else if (comm == PETSC_COMM_SELF) {
26262e5d2d2SJDBetteridge     PetscCall(PetscPrintf(comm, "PETSC_COMM_SELF\n"));
26362e5d2d2SJDBetteridge   } else {
26462e5d2d2SJDBetteridge     PetscCall(PetscPrintf(comm, "UNKNOWN_COMM\n"));
26562e5d2d2SJDBetteridge   }
26662e5d2d2SJDBetteridge   PetscCall(PetscCommDuplicate(comm, &comm, NULL));
26762e5d2d2SJDBetteridge   PetscCall(GarbageGetHMap_Private(comm, &garbage));
26862e5d2d2SJDBetteridge 
26962e5d2d2SJDBetteridge   /* Get keys from garbage hash map and sort */
27062e5d2d2SJDBetteridge   PetscCall(PetscHMapObjGetSize(garbage.map, &entries));
27162e5d2d2SJDBetteridge   PetscCall(PetscMalloc1(entries, &keys));
27262e5d2d2SJDBetteridge   offset = 0;
27362e5d2d2SJDBetteridge   PetscCall(PetscHMapObjGetKeys(garbage.map, &offset, keys));
27462e5d2d2SJDBetteridge 
27562e5d2d2SJDBetteridge   /* Pretty print entries in a table */
27662e5d2d2SJDBetteridge   PetscCallMPI(MPI_Comm_rank(comm, &rank));
27762e5d2d2SJDBetteridge   PetscCall(PetscSynchronizedPrintf(comm, "Rank %i:: ", rank));
278da0802e2SStefano Zampini   PetscCall(PetscFormatConvert("Total entries: %" PetscInt_FMT "\n", text));
27962e5d2d2SJDBetteridge   PetscCall(PetscSynchronizedPrintf(comm, text, entries));
28062e5d2d2SJDBetteridge   if (entries) {
28162e5d2d2SJDBetteridge     PetscCall(PetscSynchronizedPrintf(comm, "| Key   | Type                   | Name                             | Object ID |\n"));
28262e5d2d2SJDBetteridge     PetscCall(PetscSynchronizedPrintf(comm, "|-------|------------------------|----------------------------------|-----------|\n"));
28362e5d2d2SJDBetteridge   }
28462e5d2d2SJDBetteridge   for (ii = 0; ii < entries; ii++) {
28562e5d2d2SJDBetteridge     PetscCall(PetscHMapObjGet(garbage.map, keys[ii], &obj));
286da0802e2SStefano Zampini     PetscCall(PetscFormatConvert("| %5" PetscInt64_FMT " | %-22s | %-32s | %6" PetscInt_FMT "    |\n", text));
28762e5d2d2SJDBetteridge     PetscCall(PetscSynchronizedPrintf(comm, text, keys[ii], obj->class_name, obj->description, obj->id));
28862e5d2d2SJDBetteridge   }
28962e5d2d2SJDBetteridge   PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT));
29062e5d2d2SJDBetteridge 
29162e5d2d2SJDBetteridge   PetscCall(PetscFree(keys));
29262e5d2d2SJDBetteridge   PetscCall(PetscCommDestroy(&comm));
2933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
29462e5d2d2SJDBetteridge }
295