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