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