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