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