1 #include <petsc/private/petscimpl.h> 2 #include <petsc/private/pcbddcprivateimpl.h> 3 #include <petsc/private/pcbddcstructsimpl.h> 4 #include <petsc/private/hashmapi.h> 5 #include <petscsf.h> 6 7 PetscErrorCode PCBDDCDestroyGraphCandidatesIS(void *ctx) 8 { 9 PCBDDCGraphCandidates cand = (PCBDDCGraphCandidates)ctx; 10 11 PetscFunctionBegin; 12 for (PetscInt i = 0; i < cand->nfc; i++) PetscCall(ISDestroy(&cand->Faces[i])); 13 for (PetscInt i = 0; i < cand->nec; i++) PetscCall(ISDestroy(&cand->Edges[i])); 14 PetscCall(PetscFree(cand->Faces)); 15 PetscCall(PetscFree(cand->Edges)); 16 PetscCall(ISDestroy(&cand->Vertices)); 17 PetscCall(PetscFree(cand)); 18 PetscFunctionReturn(PETSC_SUCCESS); 19 } 20 21 PetscErrorCode PCBDDCGraphGetDirichletDofsB(PCBDDCGraph graph, IS *dirdofs) 22 { 23 PetscFunctionBegin; 24 if (graph->dirdofsB) { 25 PetscCall(PetscObjectReference((PetscObject)graph->dirdofsB)); 26 } else if (graph->has_dirichlet) { 27 PetscInt i, size; 28 PetscInt *dirdofs_idxs; 29 30 size = 0; 31 for (i = 0; i < graph->nvtxs; i++) { 32 if (graph->nodes[i].count > 1 && graph->nodes[i].special_dof == PCBDDCGRAPH_DIRICHLET_MARK) size++; 33 } 34 35 PetscCall(PetscMalloc1(size, &dirdofs_idxs)); 36 size = 0; 37 for (i = 0; i < graph->nvtxs; i++) { 38 if (graph->nodes[i].count > 1 && graph->nodes[i].special_dof == PCBDDCGRAPH_DIRICHLET_MARK) dirdofs_idxs[size++] = i; 39 } 40 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, size, dirdofs_idxs, PETSC_OWN_POINTER, &graph->dirdofsB)); 41 PetscCall(PetscObjectReference((PetscObject)graph->dirdofsB)); 42 } 43 *dirdofs = graph->dirdofsB; 44 PetscFunctionReturn(PETSC_SUCCESS); 45 } 46 47 PetscErrorCode PCBDDCGraphGetDirichletDofs(PCBDDCGraph graph, IS *dirdofs) 48 { 49 PetscFunctionBegin; 50 if (graph->dirdofs) { 51 PetscCall(PetscObjectReference((PetscObject)graph->dirdofs)); 52 } else if (graph->has_dirichlet) { 53 PetscInt i, size; 54 PetscInt *dirdofs_idxs; 55 56 size = 0; 57 for (i = 0; i < graph->nvtxs; i++) { 58 if (graph->nodes[i].special_dof == PCBDDCGRAPH_DIRICHLET_MARK) size++; 59 } 60 61 PetscCall(PetscMalloc1(size, &dirdofs_idxs)); 62 size = 0; 63 for (i = 0; i < graph->nvtxs; i++) { 64 if (graph->nodes[i].special_dof == PCBDDCGRAPH_DIRICHLET_MARK) dirdofs_idxs[size++] = i; 65 } 66 PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)graph->l2gmap), size, dirdofs_idxs, PETSC_OWN_POINTER, &graph->dirdofs)); 67 PetscCall(PetscObjectReference((PetscObject)graph->dirdofs)); 68 } 69 *dirdofs = graph->dirdofs; 70 PetscFunctionReturn(PETSC_SUCCESS); 71 } 72 73 PetscErrorCode PCBDDCGraphASCIIView(PCBDDCGraph graph, PetscInt verbosity_level, PetscViewer viewer) 74 { 75 PetscInt i, j, tabs; 76 PetscInt *queue_in_global_numbering; 77 78 PetscFunctionBegin; 79 PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 80 PetscCall(PetscViewerASCIIGetTab(viewer, &tabs)); 81 PetscCall(PetscViewerASCIIPrintf(viewer, "--------------------------------------------------\n")); 82 PetscCall(PetscViewerFlush(viewer)); 83 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "Local BDDC graph for subdomain %04d (seq %d)\n", PetscGlobalRank, graph->seq_graph)); 84 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "Number of vertices %" PetscInt_FMT "\n", graph->nvtxs)); 85 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "Number of local subdomains %" PetscInt_FMT "\n", graph->n_local_subs ? graph->n_local_subs : 1)); 86 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "Custom minimal size %" PetscInt_FMT "\n", graph->custom_minimal_size)); 87 if (graph->maxcount != PETSC_MAX_INT) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "Max count %" PetscInt_FMT "\n", graph->maxcount)); 88 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "Topological two dim? %s (set %s)\n", PetscBools[graph->twodim], PetscBools[graph->twodimset])); 89 if (verbosity_level > 2) { 90 for (i = 0; i < graph->nvtxs; i++) { 91 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%" PetscInt_FMT ":\n", i)); 92 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " which_dof: %" PetscInt_FMT "\n", graph->nodes[i].which_dof)); 93 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " special_dof: %" PetscInt_FMT "\n", graph->nodes[i].special_dof)); 94 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " shared by: %" PetscInt_FMT "\n", graph->nodes[i].count)); 95 PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE)); 96 if (graph->nodes[i].count) { 97 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " set of neighbours:")); 98 for (j = 0; j < graph->nodes[i].count; j++) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " %" PetscInt_FMT, graph->nodes[i].neighbours_set[j])); 99 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\n")); 100 } 101 PetscCall(PetscViewerASCIISetTab(viewer, tabs)); 102 PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE)); 103 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " number of local groups: %" PetscInt_FMT "\n", graph->nodes[i].local_groups_count)); 104 PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE)); 105 if (graph->nodes[i].local_groups_count) { 106 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " groups:")); 107 for (j = 0; j < graph->nodes[i].local_groups_count; j++) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " %" PetscInt_FMT, graph->nodes[i].local_groups[j])); 108 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\n")); 109 } 110 PetscCall(PetscViewerASCIISetTab(viewer, tabs)); 111 PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE)); 112 113 if (verbosity_level > 3) { 114 if (graph->xadj) { 115 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " local adj list:")); 116 PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE)); 117 for (j = graph->xadj[i]; j < graph->xadj[i + 1]; j++) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " %" PetscInt_FMT, graph->adjncy[j])); 118 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\n")); 119 PetscCall(PetscViewerASCIISetTab(viewer, tabs)); 120 PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE)); 121 } else { 122 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " no adj info\n")); 123 } 124 } 125 if (graph->n_local_subs) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " local sub id: %" PetscInt_FMT "\n", graph->local_subs[i])); 126 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " interface subset id: %" PetscInt_FMT "\n", graph->nodes[i].subset)); 127 if (graph->nodes[i].subset && graph->subset_ncc) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " ncc for subset: %" PetscInt_FMT "\n", graph->subset_ncc[graph->nodes[i].subset - 1])); 128 } 129 } 130 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "Total number of connected components %" PetscInt_FMT "\n", graph->ncc)); 131 PetscCall(PetscMalloc1(graph->cptr[graph->ncc], &queue_in_global_numbering)); 132 PetscCall(ISLocalToGlobalMappingApply(graph->l2gmap, graph->cptr[graph->ncc], graph->queue, queue_in_global_numbering)); 133 for (i = 0; i < graph->ncc; i++) { 134 PetscInt node_num = graph->queue[graph->cptr[i]]; 135 PetscBool printcc = PETSC_FALSE; 136 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " cc %" PetscInt_FMT " (size %" PetscInt_FMT ", fid %" PetscInt_FMT ", neighs:", i, graph->cptr[i + 1] - graph->cptr[i], graph->nodes[node_num].which_dof)); 137 PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE)); 138 for (j = 0; j < graph->nodes[node_num].count; j++) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " %" PetscInt_FMT, graph->nodes[node_num].neighbours_set[j])); 139 if (verbosity_level > 1) { 140 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "):")); 141 if (verbosity_level > 2 || graph->twodim || graph->nodes[node_num].count > 2 || (graph->nodes[node_num].count == 2 && graph->nodes[node_num].special_dof == PCBDDCGRAPH_NEUMANN_MARK)) { printcc = PETSC_TRUE; } 142 if (printcc) { 143 for (j = graph->cptr[i]; j < graph->cptr[i + 1]; j++) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " %" PetscInt_FMT " (%" PetscInt_FMT ")", graph->queue[j], queue_in_global_numbering[j])); 144 } 145 } else { 146 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, ")")); 147 } 148 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\n")); 149 PetscCall(PetscViewerASCIISetTab(viewer, tabs)); 150 PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE)); 151 } 152 PetscCall(PetscFree(queue_in_global_numbering)); 153 PetscCall(PetscViewerFlush(viewer)); 154 PetscFunctionReturn(PETSC_SUCCESS); 155 } 156 157 PetscErrorCode PCBDDCGraphRestoreCandidatesIS(PCBDDCGraph graph, PetscInt *n_faces, IS *FacesIS[], PetscInt *n_edges, IS *EdgesIS[], IS *VerticesIS) 158 { 159 PetscInt i; 160 PetscContainer gcand; 161 162 PetscFunctionBegin; 163 PetscCall(PetscObjectQuery((PetscObject)graph->l2gmap, "_PCBDDCGraphCandidatesIS", (PetscObject *)&gcand)); 164 if (gcand) { 165 if (n_faces) *n_faces = 0; 166 if (n_edges) *n_edges = 0; 167 if (FacesIS) *FacesIS = NULL; 168 if (EdgesIS) *EdgesIS = NULL; 169 if (VerticesIS) *VerticesIS = NULL; 170 } 171 if (n_faces) { 172 if (FacesIS) { 173 for (i = 0; i < *n_faces; i++) PetscCall(ISDestroy(&((*FacesIS)[i]))); 174 PetscCall(PetscFree(*FacesIS)); 175 } 176 *n_faces = 0; 177 } 178 if (n_edges) { 179 if (EdgesIS) { 180 for (i = 0; i < *n_edges; i++) PetscCall(ISDestroy(&((*EdgesIS)[i]))); 181 PetscCall(PetscFree(*EdgesIS)); 182 } 183 *n_edges = 0; 184 } 185 if (VerticesIS) PetscCall(ISDestroy(VerticesIS)); 186 PetscFunctionReturn(PETSC_SUCCESS); 187 } 188 189 PetscErrorCode PCBDDCGraphGetCandidatesIS(PCBDDCGraph graph, PetscInt *n_faces, IS *FacesIS[], PetscInt *n_edges, IS *EdgesIS[], IS *VerticesIS) 190 { 191 IS *ISForFaces, *ISForEdges, ISForVertices; 192 PetscInt i, nfc, nec, nvc, *idx, *mark; 193 PetscContainer gcand; 194 195 PetscFunctionBegin; 196 PetscCall(PetscObjectQuery((PetscObject)graph->l2gmap, "_PCBDDCGraphCandidatesIS", (PetscObject *)&gcand)); 197 if (gcand) { 198 PCBDDCGraphCandidates cand; 199 200 PetscCall(PetscContainerGetPointer(gcand, (void **)&cand)); 201 if (n_faces) *n_faces = cand->nfc; 202 if (FacesIS) *FacesIS = cand->Faces; 203 if (n_edges) *n_edges = cand->nec; 204 if (EdgesIS) *EdgesIS = cand->Edges; 205 if (VerticesIS) *VerticesIS = cand->Vertices; 206 PetscFunctionReturn(PETSC_SUCCESS); 207 } 208 PetscCall(PetscCalloc1(graph->ncc, &mark)); 209 /* loop on ccs to evaluate number of faces, edges and vertices */ 210 nfc = 0; 211 nec = 0; 212 nvc = 0; 213 for (i = 0; i < graph->ncc; i++) { 214 PetscInt repdof = graph->queue[graph->cptr[i]]; 215 if (graph->cptr[i + 1] - graph->cptr[i] > graph->custom_minimal_size && graph->nodes[repdof].count <= graph->maxcount) { 216 if (!graph->twodim && graph->nodes[repdof].count == 2 && graph->nodes[repdof].special_dof != PCBDDCGRAPH_NEUMANN_MARK) { 217 nfc++; 218 mark[i] = 2; 219 } else { 220 nec++; 221 mark[i] = 1; 222 } 223 } else { 224 nvc += graph->cptr[i + 1] - graph->cptr[i]; 225 } 226 } 227 228 /* allocate IS arrays for faces, edges. Vertices need a single index set. */ 229 if (FacesIS) PetscCall(PetscMalloc1(nfc, &ISForFaces)); 230 if (EdgesIS) PetscCall(PetscMalloc1(nec, &ISForEdges)); 231 if (VerticesIS) PetscCall(PetscMalloc1(nvc, &idx)); 232 233 /* loop on ccs to compute index sets for faces and edges */ 234 if (!graph->queue_sorted) { 235 PetscInt *queue_global; 236 237 PetscCall(PetscMalloc1(graph->cptr[graph->ncc], &queue_global)); 238 PetscCall(ISLocalToGlobalMappingApply(graph->l2gmap, graph->cptr[graph->ncc], graph->queue, queue_global)); 239 for (i = 0; i < graph->ncc; i++) PetscCall(PetscSortIntWithArray(graph->cptr[i + 1] - graph->cptr[i], &queue_global[graph->cptr[i]], &graph->queue[graph->cptr[i]])); 240 PetscCall(PetscFree(queue_global)); 241 graph->queue_sorted = PETSC_TRUE; 242 } 243 nfc = 0; 244 nec = 0; 245 for (i = 0; i < graph->ncc; i++) { 246 if (mark[i] == 2) { 247 if (FacesIS) PetscCall(ISCreateGeneral(PETSC_COMM_SELF, graph->cptr[i + 1] - graph->cptr[i], &graph->queue[graph->cptr[i]], PETSC_USE_POINTER, &ISForFaces[nfc])); 248 nfc++; 249 } else if (mark[i] == 1) { 250 if (EdgesIS) PetscCall(ISCreateGeneral(PETSC_COMM_SELF, graph->cptr[i + 1] - graph->cptr[i], &graph->queue[graph->cptr[i]], PETSC_USE_POINTER, &ISForEdges[nec])); 251 nec++; 252 } 253 } 254 255 /* index set for vertices */ 256 if (VerticesIS) { 257 nvc = 0; 258 for (i = 0; i < graph->ncc; i++) { 259 if (!mark[i]) { 260 PetscInt j; 261 262 for (j = graph->cptr[i]; j < graph->cptr[i + 1]; j++) { 263 idx[nvc] = graph->queue[j]; 264 nvc++; 265 } 266 } 267 } 268 /* sort vertex set (by local ordering) */ 269 PetscCall(PetscSortInt(nvc, idx)); 270 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nvc, idx, PETSC_OWN_POINTER, &ISForVertices)); 271 } 272 PetscCall(PetscFree(mark)); 273 274 /* get back info */ 275 if (n_faces) *n_faces = nfc; 276 if (FacesIS) *FacesIS = ISForFaces; 277 if (n_edges) *n_edges = nec; 278 if (EdgesIS) *EdgesIS = ISForEdges; 279 if (VerticesIS) *VerticesIS = ISForVertices; 280 PetscFunctionReturn(PETSC_SUCCESS); 281 } 282 283 PetscErrorCode PCBDDCGraphComputeConnectedComponents(PCBDDCGraph graph) 284 { 285 PetscBool adapt_interface; 286 MPI_Comm interface_comm; 287 PetscBT cornerp = NULL; 288 289 PetscFunctionBegin; 290 PetscCall(PetscObjectGetComm((PetscObject)graph->l2gmap, &interface_comm)); 291 /* compute connected components locally */ 292 PetscCall(PCBDDCGraphComputeConnectedComponentsLocal(graph)); 293 if (graph->seq_graph) PetscFunctionReturn(PETSC_SUCCESS); 294 295 if (graph->active_coords && !graph->multi_element) { /* face based corner selection XXX multi_element */ 296 PetscBT excluded; 297 PetscReal *wdist; 298 PetscInt n_neigh, *neigh, *n_shared, **shared; 299 PetscInt maxc, ns; 300 301 PetscCall(PetscBTCreate(graph->nvtxs, &cornerp)); 302 PetscCall(ISLocalToGlobalMappingGetInfo(graph->l2gmap, &n_neigh, &neigh, &n_shared, &shared)); 303 for (ns = 1, maxc = 0; ns < n_neigh; ns++) maxc = PetscMax(maxc, n_shared[ns]); 304 PetscCall(PetscMalloc1(maxc * graph->cdim, &wdist)); 305 PetscCall(PetscBTCreate(maxc, &excluded)); 306 307 for (ns = 1; ns < n_neigh; ns++) { /* first proc is self */ 308 PetscReal *anchor, mdist; 309 PetscInt fst, j, k, d, cdim = graph->cdim, n = n_shared[ns]; 310 PetscInt point1, point2, point3, point4; 311 312 /* import coordinates on shared interface */ 313 PetscCall(PetscBTMemzero(n, excluded)); 314 for (j = 0, fst = -1, k = 0; j < n; j++) { 315 PetscBool skip = PETSC_FALSE; 316 for (d = 0; d < cdim; d++) { 317 PetscReal c = graph->coords[shared[ns][j] * cdim + d]; 318 skip = (PetscBool)(skip || c == PETSC_MAX_REAL); 319 wdist[k++] = c; 320 } 321 if (skip) { 322 PetscCall(PetscBTSet(excluded, j)); 323 } else if (fst == -1) fst = j; 324 } 325 if (fst == -1) continue; 326 327 /* the dofs are sorted by global numbering, so each rank starts from the same id 328 and it will detect the same corners from the given set */ 329 330 /* find the farthest point from the starting one */ 331 anchor = wdist + fst * cdim; 332 mdist = -1.0; 333 point1 = fst; 334 for (j = fst; j < n; j++) { 335 PetscReal dist = 0.0; 336 337 if (PetscUnlikely(PetscBTLookup(excluded, j))) continue; 338 for (d = 0; d < cdim; d++) dist += (wdist[j * cdim + d] - anchor[d]) * (wdist[j * cdim + d] - anchor[d]); 339 if (dist > mdist) { 340 mdist = dist; 341 point1 = j; 342 } 343 } 344 345 /* find the farthest point from point1 */ 346 anchor = wdist + point1 * cdim; 347 mdist = -1.0; 348 point2 = point1; 349 for (j = fst; j < n; j++) { 350 PetscReal dist = 0.0; 351 352 if (PetscUnlikely(PetscBTLookup(excluded, j))) continue; 353 for (d = 0; d < cdim; d++) dist += (wdist[j * cdim + d] - anchor[d]) * (wdist[j * cdim + d] - anchor[d]); 354 if (dist > mdist) { 355 mdist = dist; 356 point2 = j; 357 } 358 } 359 360 /* find the third point maximizing the triangle area */ 361 point3 = point2; 362 if (cdim > 2) { 363 PetscReal a = 0.0; 364 365 for (d = 0; d < cdim; d++) a += (wdist[point1 * cdim + d] - wdist[point2 * cdim + d]) * (wdist[point1 * cdim + d] - wdist[point2 * cdim + d]); 366 a = PetscSqrtReal(a); 367 mdist = -1.0; 368 for (j = fst; j < n; j++) { 369 PetscReal area, b = 0.0, c = 0.0, s; 370 371 if (PetscUnlikely(PetscBTLookup(excluded, j))) continue; 372 for (d = 0; d < cdim; d++) { 373 b += (wdist[point1 * cdim + d] - wdist[j * cdim + d]) * (wdist[point1 * cdim + d] - wdist[j * cdim + d]); 374 c += (wdist[point2 * cdim + d] - wdist[j * cdim + d]) * (wdist[point2 * cdim + d] - wdist[j * cdim + d]); 375 } 376 b = PetscSqrtReal(b); 377 c = PetscSqrtReal(c); 378 s = 0.5 * (a + b + c); 379 380 /* Heron's formula, area squared */ 381 area = s * (s - a) * (s - b) * (s - c); 382 if (area > mdist) { 383 mdist = area; 384 point3 = j; 385 } 386 } 387 } 388 389 /* find the farthest point from point3 different from point1 and point2 */ 390 anchor = wdist + point3 * cdim; 391 mdist = -1.0; 392 point4 = point3; 393 for (j = fst; j < n; j++) { 394 PetscReal dist = 0.0; 395 396 if (PetscUnlikely(PetscBTLookup(excluded, j)) || j == point1 || j == point2 || j == point3) continue; 397 for (d = 0; d < cdim; d++) dist += (wdist[j * cdim + d] - anchor[d]) * (wdist[j * cdim + d] - anchor[d]); 398 if (dist > mdist) { 399 mdist = dist; 400 point4 = j; 401 } 402 } 403 404 PetscCall(PetscBTSet(cornerp, shared[ns][point1])); 405 PetscCall(PetscBTSet(cornerp, shared[ns][point2])); 406 PetscCall(PetscBTSet(cornerp, shared[ns][point3])); 407 PetscCall(PetscBTSet(cornerp, shared[ns][point4])); 408 409 /* all dofs having the same coordinates will be primal */ 410 for (j = fst; j < n; j++) { 411 PetscBool same[] = {PETSC_TRUE, PETSC_TRUE, PETSC_TRUE, PETSC_TRUE}; 412 413 if (PetscUnlikely(PetscBTLookup(excluded, j))) continue; 414 for (d = 0; d < cdim; d++) { 415 same[0] = (PetscBool)(same[0] && (PetscAbsReal(wdist[j * cdim + d] - wdist[point1 * cdim + d]) < PETSC_SMALL)); 416 same[1] = (PetscBool)(same[1] && (PetscAbsReal(wdist[j * cdim + d] - wdist[point2 * cdim + d]) < PETSC_SMALL)); 417 same[2] = (PetscBool)(same[2] && (PetscAbsReal(wdist[j * cdim + d] - wdist[point3 * cdim + d]) < PETSC_SMALL)); 418 same[3] = (PetscBool)(same[3] && (PetscAbsReal(wdist[j * cdim + d] - wdist[point4 * cdim + d]) < PETSC_SMALL)); 419 } 420 if (same[0] || same[1] || same[2] || same[3]) PetscCall(PetscBTSet(cornerp, shared[ns][j])); 421 } 422 } 423 PetscCall(PetscBTDestroy(&excluded)); 424 PetscCall(PetscFree(wdist)); 425 PetscCall(ISLocalToGlobalMappingRestoreInfo(graph->l2gmap, &n_neigh, &neigh, &n_shared, &shared)); 426 } 427 428 /* Adapt connected components if needed */ 429 adapt_interface = (cornerp || graph->multi_element) ? PETSC_TRUE : PETSC_FALSE; 430 for (PetscInt i = 0; i < graph->n_subsets && !adapt_interface; i++) { 431 if (graph->subset_ncc[i] > 1) adapt_interface = PETSC_TRUE; 432 } 433 PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &adapt_interface, 1, MPIU_BOOL, MPI_LOR, interface_comm)); 434 if (adapt_interface) { 435 PetscSF msf; 436 const PetscInt *n_ref_sharing; 437 PetscInt *labels, *rootlabels, *mrlabels; 438 PetscInt nr, nmr, nrs, ncc, cum_queue; 439 440 PetscCall(PetscMalloc1(graph->nvtxs, &labels)); 441 PetscCall(PetscArrayzero(labels, graph->nvtxs)); 442 for (PetscInt i = 0, k = 0; i < graph->ncc; i++) { 443 PetscInt s = 1; 444 for (PetscInt j = graph->cptr[i]; j < graph->cptr[i + 1]; j++) { 445 if (cornerp && PetscBTLookup(cornerp, graph->queue[j])) { 446 labels[graph->queue[j]] = -(k + s + 1); 447 s += 1; 448 } else { 449 labels[graph->queue[j]] = -(k + 1); 450 } 451 } 452 k += s; 453 } 454 PetscCall(PetscSFGetGraph(graph->interface_ref_sf, &nr, NULL, NULL, NULL)); 455 PetscCall(PetscSFGetGraph(graph->interface_subset_sf, &nrs, NULL, NULL, NULL)); 456 PetscCall(PetscSFGetMultiSF(graph->interface_subset_sf, &msf)); 457 PetscCall(PetscSFGetGraph(msf, &nmr, NULL, NULL, NULL)); 458 PetscCall(PetscCalloc2(nmr, &mrlabels, nrs, &rootlabels)); 459 460 PetscCall(PetscSFComputeDegreeBegin(graph->interface_subset_sf, &n_ref_sharing)); 461 PetscCall(PetscSFComputeDegreeEnd(graph->interface_subset_sf, &n_ref_sharing)); 462 PetscCall(PetscSFGatherBegin(graph->interface_subset_sf, MPIU_INT, labels, mrlabels)); 463 PetscCall(PetscSFGatherEnd(graph->interface_subset_sf, MPIU_INT, labels, mrlabels)); 464 465 /* analyze contributions from processes 466 The structure of mrlabels is suitable to find intersections of ccs. 467 supposing the root subset has dimension 5 and leaves with labels: 468 0: [4 4 7 4 7], (2 connected components) 469 1: [3 2 2 3 2], (2 connected components) 470 2: [1 1 6 5 6], (3 connected components) 471 the multiroot data and the new labels corresponding to intersected connected components will be (column major) 472 473 4 4 7 4 7 474 mrlabels 3 2 2 3 2 475 1 1 6 5 6 476 --------- 477 rootlabels 0 1 2 3 2 478 */ 479 for (PetscInt i = 0, rcumlabels = 0, mcumlabels = 0; i < nr; i++) { 480 const PetscInt subset_size = graph->interface_ref_rsize[i]; 481 const PetscInt *n_sharing = n_ref_sharing + rcumlabels; 482 const PetscInt *mrbuffer = mrlabels + mcumlabels; 483 PetscInt *rbuffer = rootlabels + rcumlabels; 484 PetscInt subset_counter = 0; 485 486 for (PetscInt j = 0; j < subset_size; j++) { 487 if (!rbuffer[j]) { /* found a new cc */ 488 const PetscInt *jlabels = mrbuffer + j * n_sharing[0]; 489 rbuffer[j] = ++subset_counter; 490 491 for (PetscInt k = j + 1; k < subset_size; k++) { /* check for other nodes in new cc */ 492 PetscBool same_set = PETSC_TRUE; 493 const PetscInt *klabels = mrbuffer + k * n_sharing[0]; 494 495 for (PetscInt s = 0; s < n_sharing[0]; s++) { 496 if (jlabels[s] != klabels[s]) { 497 same_set = PETSC_FALSE; 498 break; 499 } 500 } 501 if (same_set) rbuffer[k] = subset_counter; 502 } 503 } 504 } 505 if (subset_size) { 506 rcumlabels += subset_size; 507 mcumlabels += n_sharing[0] * subset_size; 508 } 509 } 510 511 /* Now communicate the intersected labels */ 512 PetscCall(PetscSFBcastBegin(graph->interface_subset_sf, MPIU_INT, rootlabels, labels, MPI_REPLACE)); 513 PetscCall(PetscSFBcastEnd(graph->interface_subset_sf, MPIU_INT, rootlabels, labels, MPI_REPLACE)); 514 PetscCall(PetscFree2(mrlabels, rootlabels)); 515 516 /* and adapt local connected components */ 517 PetscInt *ocptr, *oqueue; 518 PetscBool *touched; 519 520 PetscCall(PetscMalloc3(graph->ncc + 1, &ocptr, graph->cptr[graph->ncc], &oqueue, graph->cptr[graph->ncc], &touched)); 521 PetscCall(PetscArraycpy(ocptr, graph->cptr, graph->ncc + 1)); 522 PetscCall(PetscArraycpy(oqueue, graph->queue, graph->cptr[graph->ncc])); 523 PetscCall(PetscArrayzero(touched, graph->cptr[graph->ncc])); 524 525 ncc = 0; 526 cum_queue = 0; 527 for (PetscInt i = 0; i < graph->ncc; i++) { 528 for (PetscInt j = ocptr[i]; j < ocptr[i + 1]; j++) { 529 const PetscInt jlabel = labels[oqueue[j]]; 530 531 if (jlabel) { 532 graph->cptr[ncc] = cum_queue; 533 ncc++; 534 for (PetscInt k = j; k < ocptr[i + 1]; k++) { /* check for other nodes in new cc */ 535 if (labels[oqueue[k]] == jlabel) { 536 graph->queue[cum_queue++] = oqueue[k]; 537 labels[oqueue[k]] = 0; 538 } 539 } 540 } 541 } 542 } 543 PetscCall(PetscFree3(ocptr, oqueue, touched)); 544 PetscCall(PetscFree(labels)); 545 graph->cptr[ncc] = cum_queue; 546 graph->queue_sorted = PETSC_FALSE; 547 graph->ncc = ncc; 548 } 549 PetscCall(PetscBTDestroy(&cornerp)); 550 551 /* Determine if we are in 2D or 3D */ 552 if (!graph->twodimset) { 553 PetscBool twodim = PETSC_TRUE; 554 for (PetscInt i = 0; i < graph->ncc; i++) { 555 PetscInt repdof = graph->queue[graph->cptr[i]]; 556 PetscInt ccsize = graph->cptr[i + 1] - graph->cptr[i]; 557 if (graph->nodes[repdof].count > 2 && ccsize > graph->custom_minimal_size) { 558 twodim = PETSC_FALSE; 559 break; 560 } 561 } 562 PetscCall(MPIU_Allreduce(&twodim, &graph->twodim, 1, MPIU_BOOL, MPI_LAND, PetscObjectComm((PetscObject)graph->l2gmap))); 563 graph->twodimset = PETSC_TRUE; 564 } 565 PetscFunctionReturn(PETSC_SUCCESS); 566 } 567 568 static inline PetscErrorCode PCBDDCGraphComputeCC_Private(PCBDDCGraph graph, PetscInt pid, PetscInt *PETSC_RESTRICT queue_tip, PetscInt n_prev, PetscInt *n_added) 569 { 570 PetscInt i, j, n = 0; 571 572 const PetscInt *PETSC_RESTRICT xadj = graph->xadj; 573 const PetscInt *PETSC_RESTRICT adjncy = graph->adjncy; 574 const PetscInt *PETSC_RESTRICT subset_idxs = graph->subset_idxs[pid - 1]; 575 const PetscInt *PETSC_RESTRICT local_subs = graph->local_subs; 576 const PetscInt subset_size = graph->subset_size[pid - 1]; 577 578 PCBDDCGraphNode *PETSC_RESTRICT nodes = graph->nodes; 579 580 const PetscBool havecsr = (PetscBool)(!!xadj); 581 const PetscBool havesubs = (PetscBool)(!!graph->n_local_subs); 582 583 PetscFunctionBegin; 584 if (havecsr && !havesubs) { 585 for (i = -n_prev; i < 0; i++) { 586 const PetscInt start_dof = queue_tip[i]; 587 588 /* we assume that if a dof has a size 1 adjacency list and the corresponding entry is negative, it is connected to all dofs */ 589 if (xadj[start_dof + 1] - xadj[start_dof] == 1 && adjncy[xadj[start_dof]] < 0) { 590 for (j = 0; j < subset_size; j++) { /* pid \in [1,graph->n_subsets] */ 591 const PetscInt dof = subset_idxs[j]; 592 593 if (!nodes[dof].touched && nodes[dof].subset == pid) { 594 nodes[dof].touched = PETSC_TRUE; 595 queue_tip[n] = dof; 596 n++; 597 } 598 } 599 } else { 600 for (j = xadj[start_dof]; j < xadj[start_dof + 1]; j++) { 601 const PetscInt dof = adjncy[j]; 602 603 if (!nodes[dof].touched && nodes[dof].subset == pid) { 604 nodes[dof].touched = PETSC_TRUE; 605 queue_tip[n] = dof; 606 n++; 607 } 608 } 609 } 610 } 611 } else if (havecsr && havesubs) { 612 const PetscInt sid = local_subs[queue_tip[-n_prev]]; 613 614 for (i = -n_prev; i < 0; i++) { 615 const PetscInt start_dof = queue_tip[i]; 616 617 /* we assume that if a dof has a size 1 adjacency list and the corresponding entry is negative, it is connected to all dofs belonging to the local sub */ 618 if (xadj[start_dof + 1] - xadj[start_dof] == 1 && adjncy[xadj[start_dof]] < 0) { 619 for (j = 0; j < subset_size; j++) { /* pid \in [1,graph->n_subsets] */ 620 const PetscInt dof = subset_idxs[j]; 621 622 if (!nodes[dof].touched && nodes[dof].subset == pid && local_subs[dof] == sid) { 623 nodes[dof].touched = PETSC_TRUE; 624 queue_tip[n] = dof; 625 n++; 626 } 627 } 628 } else { 629 for (j = xadj[start_dof]; j < xadj[start_dof + 1]; j++) { 630 const PetscInt dof = adjncy[j]; 631 632 if (!nodes[dof].touched && nodes[dof].subset == pid && local_subs[dof] == sid) { 633 nodes[dof].touched = PETSC_TRUE; 634 queue_tip[n] = dof; 635 n++; 636 } 637 } 638 } 639 } 640 } else if (havesubs) { /* sub info only */ 641 const PetscInt sid = local_subs[queue_tip[-n_prev]]; 642 643 for (j = 0; j < subset_size; j++) { /* pid \in [1,graph->n_subsets] */ 644 const PetscInt dof = subset_idxs[j]; 645 646 if (!nodes[dof].touched && nodes[dof].subset == pid && local_subs[dof] == sid) { 647 nodes[dof].touched = PETSC_TRUE; 648 queue_tip[n] = dof; 649 n++; 650 } 651 } 652 } else { 653 for (j = 0; j < subset_size; j++) { /* pid \in [1,graph->n_subsets] */ 654 const PetscInt dof = subset_idxs[j]; 655 656 if (!nodes[dof].touched && nodes[dof].subset == pid) { 657 nodes[dof].touched = PETSC_TRUE; 658 queue_tip[n] = dof; 659 n++; 660 } 661 } 662 } 663 *n_added = n; 664 PetscFunctionReturn(PETSC_SUCCESS); 665 } 666 667 PetscErrorCode PCBDDCGraphComputeConnectedComponentsLocal(PCBDDCGraph graph) 668 { 669 PetscInt ncc, cum_queue; 670 671 PetscFunctionBegin; 672 PetscCheck(graph->setupcalled, PetscObjectComm((PetscObject)graph->l2gmap), PETSC_ERR_ORDER, "PCBDDCGraphSetUp should be called first"); 673 /* quiet return if there isn't any local info */ 674 if (!graph->xadj && !graph->n_local_subs) PetscFunctionReturn(PETSC_SUCCESS); 675 676 /* reset any previous search of connected components */ 677 for (PetscInt i = 0; i < graph->nvtxs; i++) graph->nodes[i].touched = PETSC_FALSE; 678 if (!graph->seq_graph) { 679 for (PetscInt i = 0; i < graph->nvtxs; i++) { 680 if (graph->nodes[i].special_dof == PCBDDCGRAPH_DIRICHLET_MARK || graph->nodes[i].count < 2) graph->nodes[i].touched = PETSC_TRUE; 681 } 682 } 683 684 /* begin search for connected components */ 685 cum_queue = 0; 686 ncc = 0; 687 for (PetscInt n = 0; n < graph->n_subsets; n++) { 688 const PetscInt *subset_idxs = graph->subset_idxs[n]; 689 const PetscInt pid = n + 1; /* partition labeled by 0 is discarded */ 690 691 PetscInt found = 0, prev = 0, first = 0, ncc_pid = 0; 692 693 while (found != graph->subset_size[n]) { 694 PetscInt added = 0; 695 696 if (!prev) { /* search for new starting dof */ 697 while (graph->nodes[subset_idxs[first]].touched) first++; 698 graph->nodes[subset_idxs[first]].touched = PETSC_TRUE; 699 graph->queue[cum_queue] = subset_idxs[first]; 700 graph->cptr[ncc] = cum_queue; 701 prev = 1; 702 cum_queue++; 703 found++; 704 ncc_pid++; 705 ncc++; 706 } 707 PetscCall(PCBDDCGraphComputeCC_Private(graph, pid, graph->queue + cum_queue, prev, &added)); 708 if (!added) { 709 graph->subset_ncc[n] = ncc_pid; 710 graph->cptr[ncc] = cum_queue; 711 } 712 prev = added; 713 found += added; 714 cum_queue += added; 715 if (added && found == graph->subset_size[n]) { 716 graph->subset_ncc[n] = ncc_pid; 717 graph->cptr[ncc] = cum_queue; 718 } 719 } 720 } 721 graph->ncc = ncc; 722 graph->queue_sorted = PETSC_FALSE; 723 PetscFunctionReturn(PETSC_SUCCESS); 724 } 725 726 PetscErrorCode PCBDDCGraphSetUp(PCBDDCGraph graph, PetscInt custom_minimal_size, IS neumann_is, IS dirichlet_is, PetscInt n_ISForDofs, IS ISForDofs[], IS custom_primal_vertices) 727 { 728 IS subset; 729 MPI_Comm comm; 730 const PetscInt *is_indices; 731 PetscInt *queue_global, *nodecount, **nodeneighs; 732 PetscInt i, j, k, total_counts, nodes_touched, is_size, nvtxs = graph->nvtxs; 733 PetscMPIInt size, rank; 734 PetscBool same_set; 735 736 PetscFunctionBegin; 737 PetscValidLogicalCollectiveInt(graph->l2gmap, custom_minimal_size, 2); 738 if (neumann_is) { 739 PetscValidHeaderSpecific(neumann_is, IS_CLASSID, 3); 740 PetscCheckSameComm(graph->l2gmap, 1, neumann_is, 3); 741 } 742 graph->has_dirichlet = PETSC_FALSE; 743 if (dirichlet_is) { 744 PetscValidHeaderSpecific(dirichlet_is, IS_CLASSID, 4); 745 PetscCheckSameComm(graph->l2gmap, 1, dirichlet_is, 4); 746 graph->has_dirichlet = PETSC_TRUE; 747 } 748 PetscValidLogicalCollectiveInt(graph->l2gmap, n_ISForDofs, 5); 749 for (i = 0; i < n_ISForDofs; i++) { 750 PetscValidHeaderSpecific(ISForDofs[i], IS_CLASSID, 6); 751 PetscCheckSameComm(graph->l2gmap, 1, ISForDofs[i], 6); 752 } 753 if (custom_primal_vertices) { 754 PetscValidHeaderSpecific(custom_primal_vertices, IS_CLASSID, 7); 755 PetscCheckSameComm(graph->l2gmap, 1, custom_primal_vertices, 7); 756 } 757 for (i = 0; i < nvtxs; i++) graph->nodes[i].touched = PETSC_FALSE; 758 759 PetscCall(PetscObjectGetComm((PetscObject)graph->l2gmap, &comm)); 760 PetscCallMPI(MPI_Comm_size(comm, &size)); 761 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 762 763 /* custom_minimal_size */ 764 graph->custom_minimal_size = custom_minimal_size; 765 766 /* get node info from l2gmap */ 767 PetscCall(ISLocalToGlobalMappingGetNodeInfo(graph->l2gmap, NULL, &nodecount, &nodeneighs)); 768 769 /* Allocate space for storing the set of neighbours for each node */ 770 graph->multi_element = PETSC_FALSE; 771 for (i = 0; i < nvtxs; i++) { 772 graph->nodes[i].count = nodecount[i]; 773 if (!graph->seq_graph) { 774 PetscCall(PetscMalloc1(nodecount[i], &graph->nodes[i].neighbours_set)); 775 PetscCall(PetscArraycpy(graph->nodes[i].neighbours_set, nodeneighs[i], nodecount[i])); 776 777 if (!graph->multi_element) { 778 PetscInt nself; 779 for (j = 0, nself = 0; j < graph->nodes[i].count; j++) 780 if (graph->nodes[i].neighbours_set[j] == rank) nself++; 781 if (nself > 1) graph->multi_element = PETSC_TRUE; 782 } 783 } else { 784 PetscCall(PetscCalloc1(nodecount[i], &graph->nodes[i].neighbours_set)); 785 } 786 } 787 PetscCall(ISLocalToGlobalMappingRestoreNodeInfo(graph->l2gmap, NULL, &nodecount, &nodeneighs)); 788 PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &graph->multi_element, 1, MPIU_BOOL, MPI_LOR, comm)); 789 790 /* compute local groups */ 791 if (graph->multi_element) { 792 const PetscInt *idxs, *indegree; 793 IS is, lis; 794 PetscLayout layout; 795 PetscSF sf, multisf; 796 PetscInt n, nmulti, c, *multi_root_subs, *start; 797 798 PetscCheck(graph->local_subs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing local subdomain information"); 799 800 PetscCall(ISLocalToGlobalMappingGetIndices(graph->l2gmap, &idxs)); 801 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nvtxs, idxs, PETSC_USE_POINTER, &is)); 802 PetscCall(ISRenumber(is, NULL, &n, &lis)); 803 PetscCall(ISDestroy(&is)); 804 805 PetscCall(ISLocalToGlobalMappingRestoreIndices(graph->l2gmap, &idxs)); 806 PetscCall(ISGetIndices(lis, &idxs)); 807 PetscCall(PetscLayoutCreate(PETSC_COMM_SELF, &layout)); 808 PetscCall(PetscLayoutSetSize(layout, n)); 809 PetscCall(PetscSFCreate(PETSC_COMM_SELF, &sf)); 810 PetscCall(PetscSFSetGraphLayout(sf, layout, nvtxs, NULL, PETSC_OWN_POINTER, idxs)); 811 PetscCall(PetscLayoutDestroy(&layout)); 812 PetscCall(PetscSFGetMultiSF(sf, &multisf)); 813 PetscCall(PetscSFComputeDegreeBegin(sf, &indegree)); 814 PetscCall(PetscSFComputeDegreeEnd(sf, &indegree)); 815 PetscCall(PetscSFGetGraph(multisf, &nmulti, NULL, NULL, NULL)); 816 PetscCall(PetscMalloc2(nmulti, &multi_root_subs, n + 1, &start)); 817 start[0] = 0; 818 for (i = 0; i < n; i++) start[i + 1] = start[i] + indegree[i]; 819 PetscCall(PetscSFGatherBegin(sf, MPIU_INT, graph->local_subs, multi_root_subs)); 820 PetscCall(PetscSFGatherEnd(sf, MPIU_INT, graph->local_subs, multi_root_subs)); 821 for (i = 0; i < nvtxs; i++) { 822 PetscInt gid = idxs[i]; 823 824 graph->nodes[i].local_sub = graph->local_subs[i]; 825 for (j = 0, c = 0; j < graph->nodes[i].count; j++) { 826 if (graph->nodes[i].neighbours_set[j] == rank) c++; 827 } 828 PetscCheck(c == indegree[idxs[i]], PETSC_COMM_SELF, PETSC_ERR_PLIB, "%" PetscInt_FMT " != %" PetscInt_FMT, c, indegree[idxs[i]]); 829 PetscCall(PetscMalloc1(c, &graph->nodes[i].local_groups)); 830 for (j = 0; j < c; j++) graph->nodes[i].local_groups[j] = multi_root_subs[start[gid] + j]; 831 PetscCall(PetscSortInt(c, graph->nodes[i].local_groups)); 832 graph->nodes[i].local_groups_count = c; 833 } 834 PetscCall(PetscFree2(multi_root_subs, start)); 835 PetscCall(ISRestoreIndices(lis, &idxs)); 836 PetscCall(ISDestroy(&lis)); 837 PetscCall(PetscSFDestroy(&sf)); 838 } 839 840 /* 841 Get info for dofs splitting 842 User can specify just a subset; an additional field is considered as a complementary field 843 */ 844 for (i = 0, k = 0; i < n_ISForDofs; i++) { 845 PetscInt bs; 846 847 PetscCall(ISGetBlockSize(ISForDofs[i], &bs)); 848 k += bs; 849 } 850 for (i = 0; i < nvtxs; i++) graph->nodes[i].which_dof = k; /* by default a dof belongs to the complement set */ 851 for (i = 0, k = 0; i < n_ISForDofs; i++) { 852 PetscInt bs; 853 854 PetscCall(ISGetLocalSize(ISForDofs[i], &is_size)); 855 PetscCall(ISGetBlockSize(ISForDofs[i], &bs)); 856 PetscCall(ISGetIndices(ISForDofs[i], (const PetscInt **)&is_indices)); 857 for (j = 0; j < is_size / bs; j++) { 858 PetscInt b; 859 860 for (b = 0; b < bs; b++) { 861 PetscInt jj = bs * j + b; 862 863 if (is_indices[jj] > -1 && is_indices[jj] < nvtxs) { /* out of bounds indices (if any) are skipped */ 864 graph->nodes[is_indices[jj]].which_dof = k + b; 865 } 866 } 867 } 868 PetscCall(ISRestoreIndices(ISForDofs[i], (const PetscInt **)&is_indices)); 869 k += bs; 870 } 871 872 /* Take into account Neumann nodes */ 873 if (neumann_is) { 874 PetscCall(ISGetLocalSize(neumann_is, &is_size)); 875 PetscCall(ISGetIndices(neumann_is, (const PetscInt **)&is_indices)); 876 for (i = 0; i < is_size; i++) { 877 if (is_indices[i] > -1 && is_indices[i] < nvtxs) { /* out of bounds indices (if any) are skipped */ 878 graph->nodes[is_indices[i]].special_dof = PCBDDCGRAPH_NEUMANN_MARK; 879 } 880 } 881 PetscCall(ISRestoreIndices(neumann_is, (const PetscInt **)&is_indices)); 882 } 883 884 /* Take into account Dirichlet nodes (they overwrite any mark previously set) */ 885 if (dirichlet_is) { 886 PetscCall(ISGetLocalSize(dirichlet_is, &is_size)); 887 PetscCall(ISGetIndices(dirichlet_is, (const PetscInt **)&is_indices)); 888 for (i = 0; i < is_size; i++) { 889 if (is_indices[i] > -1 && is_indices[i] < nvtxs) { /* out of bounds indices (if any) are skipped */ 890 if (!graph->seq_graph) { /* dirichlet nodes treated as internal */ 891 graph->nodes[is_indices[i]].touched = PETSC_TRUE; 892 graph->nodes[is_indices[i]].subset = 0; 893 } 894 graph->nodes[is_indices[i]].special_dof = PCBDDCGRAPH_DIRICHLET_MARK; 895 } 896 } 897 PetscCall(ISRestoreIndices(dirichlet_is, (const PetscInt **)&is_indices)); 898 } 899 900 /* mark special nodes (if any) -> each will become a single dof equivalence class (i.e. point constraint for BDDC) */ 901 if (custom_primal_vertices) { 902 PetscCall(ISGetLocalSize(custom_primal_vertices, &is_size)); 903 PetscCall(ISGetIndices(custom_primal_vertices, (const PetscInt **)&is_indices)); 904 for (i = 0, j = 0; i < is_size; i++) { 905 if (is_indices[i] > -1 && is_indices[i] < nvtxs && graph->nodes[is_indices[i]].special_dof != PCBDDCGRAPH_DIRICHLET_MARK) { /* out of bounds indices (if any) are skipped */ 906 graph->nodes[is_indices[i]].special_dof = PCBDDCGRAPH_SPECIAL_MARK - j; 907 j++; 908 } 909 } 910 PetscCall(ISRestoreIndices(custom_primal_vertices, (const PetscInt **)&is_indices)); 911 } 912 913 /* mark interior nodes as touched and belonging to partition number 0 */ 914 if (!graph->seq_graph) { 915 for (i = 0; i < nvtxs; i++) { 916 if (graph->nodes[i].count < 2) { 917 graph->nodes[i].touched = PETSC_TRUE; 918 graph->nodes[i].subset = 0; 919 } 920 } 921 } 922 923 /* init graph structure and compute default subsets */ 924 nodes_touched = 0; 925 for (i = 0; i < nvtxs; i++) 926 if (graph->nodes[i].touched) nodes_touched++; 927 928 i = 0; 929 graph->ncc = 0; 930 total_counts = 0; 931 932 /* allocated space for queues */ 933 if (graph->seq_graph) { 934 PetscCall(PetscMalloc2(nvtxs + 1, &graph->cptr, nvtxs, &graph->queue)); 935 } else { 936 PetscInt nused = nvtxs - nodes_touched; 937 PetscCall(PetscMalloc2(nused + 1, &graph->cptr, nused, &graph->queue)); 938 } 939 940 while (nodes_touched < nvtxs) { 941 /* find first untouched node in local ordering */ 942 while (graph->nodes[i].touched) i++; 943 graph->nodes[i].touched = PETSC_TRUE; 944 graph->nodes[i].subset = graph->ncc + 1; 945 graph->cptr[graph->ncc] = total_counts; 946 graph->queue[total_counts] = i; 947 total_counts++; 948 nodes_touched++; 949 950 /* now find all other nodes having the same set of sharing subdomains */ 951 const PCBDDCGraphNode *nodei = &graph->nodes[i]; 952 const PetscInt icount = nodei->count; 953 const PetscInt iwhich_dof = nodei->which_dof; 954 const PetscInt ispecial_dof = nodei->special_dof; 955 const PetscInt ilocal_groups_count = nodei->local_groups_count; 956 const PetscInt *PETSC_RESTRICT ineighbours_set = nodei->neighbours_set; 957 const PetscInt *PETSC_RESTRICT ilocal_groups = nodei->local_groups; 958 for (j = i + 1; j < nvtxs; j++) { 959 PCBDDCGraphNode *PETSC_RESTRICT nodej = &graph->nodes[j]; 960 961 if (nodej->touched) continue; 962 /* check for same number of sharing subdomains, dof number and same special mark */ 963 if (icount == nodej->count && iwhich_dof == nodej->which_dof && ispecial_dof == nodej->special_dof) { 964 PetscBool mpi_shared = PETSC_TRUE; 965 966 /* check for same set of sharing subdomains */ 967 same_set = PETSC_TRUE; 968 for (k = 0; k < icount; k++) { 969 if (ineighbours_set[k] != nodej->neighbours_set[k]) { 970 same_set = PETSC_FALSE; 971 break; 972 } 973 } 974 975 if (graph->multi_element) { 976 mpi_shared = PETSC_FALSE; 977 for (k = 0; k < icount; k++) 978 if (ineighbours_set[k] != rank) { 979 mpi_shared = PETSC_TRUE; 980 break; 981 } 982 } 983 984 /* check for same local groups 985 shared dofs at the process boundaries will be handled differently */ 986 if (same_set && !mpi_shared) { 987 if (ilocal_groups_count != nodej->local_groups_count) same_set = PETSC_FALSE; 988 else { 989 for (k = 0; k < ilocal_groups_count; k++) { 990 if (ilocal_groups[k] != nodej->local_groups[k]) { 991 same_set = PETSC_FALSE; 992 break; 993 } 994 } 995 } 996 } 997 998 /* Add to subset */ 999 if (same_set) { 1000 nodej->touched = PETSC_TRUE; 1001 nodej->subset = graph->ncc + 1; 1002 nodes_touched++; 1003 graph->queue[total_counts] = j; 1004 total_counts++; 1005 } 1006 } 1007 } 1008 graph->ncc++; 1009 } 1010 graph->cptr[graph->ncc] = total_counts; 1011 1012 /* set default number of subsets */ 1013 graph->n_subsets = graph->ncc; 1014 PetscCall(PetscMalloc1(graph->n_subsets, &graph->subset_ncc)); 1015 for (i = 0; i < graph->n_subsets; i++) graph->subset_ncc[i] = 1; 1016 1017 PetscCall(PetscMalloc1(graph->ncc, &graph->subset_ref_node)); 1018 PetscCall(PetscMalloc1(graph->cptr[graph->ncc], &queue_global)); 1019 PetscCall(PetscMalloc2(graph->ncc, &graph->subset_size, graph->ncc, &graph->subset_idxs)); 1020 if (graph->multi_element) PetscCall(PetscMalloc1(graph->ncc, &graph->gsubset_size)); 1021 else graph->gsubset_size = graph->subset_size; 1022 PetscCall(ISLocalToGlobalMappingApply(graph->l2gmap, graph->cptr[graph->ncc], graph->queue, queue_global)); 1023 1024 PetscHMapI cnt_unique; 1025 1026 PetscCall(PetscHMapICreate(&cnt_unique)); 1027 for (j = 0; j < graph->ncc; j++) { 1028 PetscInt c = 0, ref_node = PETSC_MAX_INT; 1029 1030 for (k = graph->cptr[j]; k < graph->cptr[j + 1]; k++) { 1031 ref_node = PetscMin(ref_node, queue_global[k]); 1032 if (graph->multi_element) { 1033 PetscBool missing; 1034 PetscHashIter iter; 1035 1036 PetscCall(PetscHMapIPut(cnt_unique, queue_global[k], &iter, &missing)); 1037 if (missing) c++; 1038 } 1039 } 1040 graph->gsubset_size[j] = c; 1041 graph->subset_size[j] = graph->cptr[j + 1] - graph->cptr[j]; 1042 graph->subset_ref_node[j] = ref_node; 1043 if (graph->multi_element) PetscCall(PetscHMapIClear(cnt_unique)); 1044 } 1045 PetscCall(PetscHMapIDestroy(&cnt_unique)); 1046 1047 /* save information on subsets (needed when analyzing the connected components) */ 1048 if (graph->ncc) { 1049 PetscCall(PetscMalloc1(graph->cptr[graph->ncc], &graph->subset_idxs[0])); 1050 PetscCall(PetscArrayzero(graph->subset_idxs[0], graph->cptr[graph->ncc])); 1051 for (j = 1; j < graph->ncc; j++) { graph->subset_idxs[j] = graph->subset_idxs[j - 1] + graph->subset_size[j - 1]; } 1052 PetscCall(PetscArraycpy(graph->subset_idxs[0], graph->queue, graph->cptr[graph->ncc])); 1053 } 1054 1055 /* check consistency and create SF to analyze components on the interface between subdomains */ 1056 if (!graph->seq_graph) { 1057 PetscSF msf; 1058 PetscLayout map; 1059 const PetscInt *degree; 1060 PetscInt nr, nmr, *rdata; 1061 PetscBool valid = PETSC_TRUE; 1062 PetscInt subset_N; 1063 IS subset_n; 1064 const PetscInt *idxs; 1065 1066 PetscCall(ISCreateGeneral(comm, graph->n_subsets, graph->subset_ref_node, PETSC_USE_POINTER, &subset)); 1067 PetscCall(ISRenumber(subset, NULL, &subset_N, &subset_n)); 1068 PetscCall(ISDestroy(&subset)); 1069 1070 PetscCall(PetscSFCreate(comm, &graph->interface_ref_sf)); 1071 PetscCall(PetscLayoutCreateFromSizes(comm, PETSC_DECIDE, subset_N, 1, &map)); 1072 PetscCall(ISGetIndices(subset_n, &idxs)); 1073 PetscCall(PetscSFSetGraphLayout(graph->interface_ref_sf, map, graph->n_subsets, NULL, PETSC_OWN_POINTER, idxs)); 1074 PetscCall(ISRestoreIndices(subset_n, &idxs)); 1075 PetscCall(ISDestroy(&subset_n)); 1076 PetscCall(PetscLayoutDestroy(&map)); 1077 1078 PetscCall(PetscSFComputeDegreeBegin(graph->interface_ref_sf, °ree)); 1079 PetscCall(PetscSFComputeDegreeEnd(graph->interface_ref_sf, °ree)); 1080 PetscCall(PetscSFGetMultiSF(graph->interface_ref_sf, &msf)); 1081 PetscCall(PetscSFGetGraph(graph->interface_ref_sf, &nr, NULL, NULL, NULL)); 1082 PetscCall(PetscSFGetGraph(msf, &nmr, NULL, NULL, NULL)); 1083 PetscCall(PetscCalloc1(nmr, &rdata)); 1084 PetscCall(PetscSFGatherBegin(graph->interface_ref_sf, MPIU_INT, graph->gsubset_size, rdata)); 1085 PetscCall(PetscSFGatherEnd(graph->interface_ref_sf, MPIU_INT, graph->gsubset_size, rdata)); 1086 for (PetscInt i = 0, c = 0; i < nr && valid; i++) { 1087 for (PetscInt j = 0; j < degree[i]; j++) { 1088 if (rdata[j + c] != rdata[c]) valid = PETSC_FALSE; 1089 } 1090 c += degree[i]; 1091 } 1092 PetscCall(PetscFree(rdata)); 1093 PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &valid, 1, MPIU_BOOL, MPI_LAND, comm)); 1094 PetscCheck(valid, comm, PETSC_ERR_PLIB, "Initial local subsets are not consistent"); 1095 1096 /* Now create SF with each root extended to gsubset_size roots */ 1097 PetscInt mss = 0; 1098 const PetscSFNode *subs_remote; 1099 1100 PetscCall(PetscSFGetGraph(graph->interface_ref_sf, NULL, NULL, NULL, &subs_remote)); 1101 for (PetscInt i = 0; i < graph->n_subsets; i++) mss = PetscMax(graph->subset_size[i], mss); 1102 1103 PetscInt nri, nli, *start_rsize, *cum_rsize; 1104 PetscCall(PetscCalloc1(graph->n_subsets + 1, &start_rsize)); 1105 PetscCall(PetscCalloc1(nr, &graph->interface_ref_rsize)); 1106 PetscCall(PetscMalloc1(nr + 1, &cum_rsize)); 1107 PetscCall(PetscSFReduceBegin(graph->interface_ref_sf, MPIU_INT, graph->gsubset_size, graph->interface_ref_rsize, MPI_REPLACE)); 1108 PetscCall(PetscSFReduceEnd(graph->interface_ref_sf, MPIU_INT, graph->gsubset_size, graph->interface_ref_rsize, MPI_REPLACE)); 1109 1110 nri = 0; 1111 cum_rsize[0] = 0; 1112 for (PetscInt i = 0; i < nr; i++) { 1113 nri += graph->interface_ref_rsize[i]; 1114 cum_rsize[i + 1] = cum_rsize[i] + graph->interface_ref_rsize[i]; 1115 } 1116 nli = graph->cptr[graph->ncc]; 1117 PetscCall(PetscSFBcastBegin(graph->interface_ref_sf, MPIU_INT, cum_rsize, start_rsize, MPI_REPLACE)); 1118 PetscCall(PetscSFBcastEnd(graph->interface_ref_sf, MPIU_INT, cum_rsize, start_rsize, MPI_REPLACE)); 1119 PetscCall(PetscFree(cum_rsize)); 1120 1121 PetscInt *ilocal, *queue_global_uniq; 1122 PetscSFNode *iremote; 1123 PetscBool *touched; 1124 1125 PetscCall(PetscSFCreate(comm, &graph->interface_subset_sf)); 1126 PetscCall(PetscMalloc1(nli, &ilocal)); 1127 PetscCall(PetscMalloc1(nli, &iremote)); 1128 PetscCall(PetscMalloc2(mss, &queue_global_uniq, mss, &touched)); 1129 for (PetscInt i = 0, nli = 0; i < graph->n_subsets; i++) { 1130 const PetscMPIInt rr = subs_remote[i].rank; 1131 const PetscInt start = start_rsize[i]; 1132 const PetscInt subset_size = graph->subset_size[i]; 1133 const PetscInt gsubset_size = graph->gsubset_size[i]; 1134 const PetscInt *subset_idxs = graph->subset_idxs[i]; 1135 const PetscInt *lsub_queue_global = queue_global + graph->cptr[i]; 1136 1137 k = subset_size; 1138 PetscCall(PetscArrayzero(touched, subset_size)); 1139 PetscCall(PetscArraycpy(queue_global_uniq, lsub_queue_global, subset_size)); 1140 PetscCall(PetscSortRemoveDupsInt(&k, queue_global_uniq)); 1141 PetscCheck(k == gsubset_size, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid local subset %" PetscInt_FMT " size %" PetscInt_FMT " != %" PetscInt_FMT, i, k, gsubset_size); 1142 1143 PetscInt t = 0, j = 0; 1144 while (t < subset_size) { 1145 while (j < subset_size && touched[j]) j++; 1146 PetscCheck(j < subset_size, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unexpected %" PetscInt_FMT " >= %" PetscInt_FMT, j, subset_size); 1147 const PetscInt ls = graph->nodes[subset_idxs[j]].local_sub; 1148 1149 for (k = j; k < subset_size; k++) { 1150 if (graph->nodes[subset_idxs[k]].local_sub == ls) { 1151 PetscInt ig; 1152 1153 PetscCall(PetscFindInt(lsub_queue_global[k], gsubset_size, queue_global_uniq, &ig)); 1154 ilocal[nli] = subset_idxs[k]; 1155 iremote[nli].rank = rr; 1156 iremote[nli].index = start + ig; 1157 touched[k] = PETSC_TRUE; 1158 nli++; 1159 t++; 1160 } 1161 } 1162 } 1163 } 1164 PetscCheck(nli == graph->cptr[graph->ncc], PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid ilocal size %" PetscInt_FMT " != %" PetscInt_FMT, nli, graph->cptr[graph->ncc]); 1165 PetscCall(PetscSFSetGraph(graph->interface_subset_sf, nri, nli, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER)); 1166 PetscCall(PetscFree(start_rsize)); 1167 PetscCall(PetscFree2(queue_global_uniq, touched)); 1168 } 1169 PetscCall(PetscFree(queue_global)); 1170 1171 /* free workspace */ 1172 graph->setupcalled = PETSC_TRUE; 1173 PetscFunctionReturn(PETSC_SUCCESS); 1174 } 1175 1176 PetscErrorCode PCBDDCGraphResetCoords(PCBDDCGraph graph) 1177 { 1178 PetscFunctionBegin; 1179 if (!graph) PetscFunctionReturn(PETSC_SUCCESS); 1180 PetscCall(PetscFree(graph->coords)); 1181 graph->cdim = 0; 1182 graph->cnloc = 0; 1183 graph->cloc = PETSC_FALSE; 1184 PetscFunctionReturn(PETSC_SUCCESS); 1185 } 1186 1187 PetscErrorCode PCBDDCGraphResetCSR(PCBDDCGraph graph) 1188 { 1189 PetscFunctionBegin; 1190 if (!graph) PetscFunctionReturn(PETSC_SUCCESS); 1191 if (graph->freecsr) { 1192 PetscCall(PetscFree(graph->xadj)); 1193 PetscCall(PetscFree(graph->adjncy)); 1194 } else { 1195 graph->xadj = NULL; 1196 graph->adjncy = NULL; 1197 } 1198 graph->freecsr = PETSC_FALSE; 1199 graph->nvtxs_csr = 0; 1200 PetscFunctionReturn(PETSC_SUCCESS); 1201 } 1202 1203 PetscErrorCode PCBDDCGraphReset(PCBDDCGraph graph) 1204 { 1205 PetscFunctionBegin; 1206 if (!graph) PetscFunctionReturn(PETSC_SUCCESS); 1207 PetscCall(ISLocalToGlobalMappingDestroy(&graph->l2gmap)); 1208 PetscCall(PetscFree(graph->subset_ncc)); 1209 PetscCall(PetscFree(graph->subset_ref_node)); 1210 for (PetscInt i = 0; i < graph->nvtxs; i++) { 1211 PetscCall(PetscFree(graph->nodes[i].neighbours_set)); 1212 PetscCall(PetscFree(graph->nodes[i].local_groups)); 1213 } 1214 PetscCall(PetscFree(graph->nodes)); 1215 PetscCall(PetscFree2(graph->cptr, graph->queue)); 1216 if (graph->subset_idxs) PetscCall(PetscFree(graph->subset_idxs[0])); 1217 PetscCall(PetscFree2(graph->subset_size, graph->subset_idxs)); 1218 if (graph->multi_element) PetscCall(PetscFree(graph->gsubset_size)); 1219 PetscCall(PetscFree(graph->interface_ref_rsize)); 1220 PetscCall(PetscSFDestroy(&graph->interface_subset_sf)); 1221 PetscCall(PetscSFDestroy(&graph->interface_ref_sf)); 1222 PetscCall(ISDestroy(&graph->dirdofs)); 1223 PetscCall(ISDestroy(&graph->dirdofsB)); 1224 if (graph->n_local_subs) PetscCall(PetscFree(graph->local_subs)); 1225 graph->multi_element = PETSC_FALSE; 1226 graph->has_dirichlet = PETSC_FALSE; 1227 graph->twodimset = PETSC_FALSE; 1228 graph->twodim = PETSC_FALSE; 1229 graph->nvtxs = 0; 1230 graph->nvtxs_global = 0; 1231 graph->n_subsets = 0; 1232 graph->custom_minimal_size = 1; 1233 graph->n_local_subs = 0; 1234 graph->maxcount = PETSC_MAX_INT; 1235 graph->seq_graph = PETSC_FALSE; 1236 graph->setupcalled = PETSC_FALSE; 1237 PetscFunctionReturn(PETSC_SUCCESS); 1238 } 1239 1240 PetscErrorCode PCBDDCGraphInit(PCBDDCGraph graph, ISLocalToGlobalMapping l2gmap, PetscInt N, PetscInt maxcount) 1241 { 1242 PetscInt n; 1243 1244 PetscFunctionBegin; 1245 PetscAssertPointer(graph, 1); 1246 PetscValidHeaderSpecific(l2gmap, IS_LTOGM_CLASSID, 2); 1247 PetscValidLogicalCollectiveInt(l2gmap, N, 3); 1248 PetscValidLogicalCollectiveInt(l2gmap, maxcount, 4); 1249 /* raise an error if already allocated */ 1250 PetscCheck(!graph->nvtxs_global, PetscObjectComm((PetscObject)l2gmap), PETSC_ERR_PLIB, "BDDCGraph already initialized"); 1251 /* set number of vertices */ 1252 PetscCall(PetscObjectReference((PetscObject)l2gmap)); 1253 graph->l2gmap = l2gmap; 1254 PetscCall(ISLocalToGlobalMappingGetSize(l2gmap, &n)); 1255 graph->nvtxs = n; 1256 graph->nvtxs_global = N; 1257 /* allocate used space */ 1258 PetscCall(PetscCalloc1(graph->nvtxs, &graph->nodes)); 1259 /* use -1 as a default value for which_dof array */ 1260 for (n = 0; n < graph->nvtxs; n++) graph->nodes[n].which_dof = -1; 1261 1262 /* zeroes workspace for values of ncc */ 1263 graph->subset_ncc = NULL; 1264 graph->subset_ref_node = NULL; 1265 /* maxcount for cc */ 1266 graph->maxcount = maxcount; 1267 PetscFunctionReturn(PETSC_SUCCESS); 1268 } 1269 1270 PetscErrorCode PCBDDCGraphDestroy(PCBDDCGraph *graph) 1271 { 1272 PetscFunctionBegin; 1273 PetscCall(PCBDDCGraphResetCSR(*graph)); 1274 PetscCall(PCBDDCGraphResetCoords(*graph)); 1275 PetscCall(PCBDDCGraphReset(*graph)); 1276 PetscCall(PetscFree(*graph)); 1277 PetscFunctionReturn(PETSC_SUCCESS); 1278 } 1279 1280 PetscErrorCode PCBDDCGraphCreate(PCBDDCGraph *graph) 1281 { 1282 PCBDDCGraph new_graph; 1283 1284 PetscFunctionBegin; 1285 PetscCall(PetscNew(&new_graph)); 1286 new_graph->custom_minimal_size = 1; 1287 *graph = new_graph; 1288 PetscFunctionReturn(PETSC_SUCCESS); 1289 } 1290