xref: /petsc/src/ksp/pc/impls/bddc/bddcgraph.c (revision 32fe681d04f1d7625d6fc236e17224268fcf5f36)
1af0996ceSBarry Smith #include <petsc/private/petscimpl.h>
25e5bbd0aSStefano Zampini #include <petsc/private/pcbddcprivateimpl.h>
35e5bbd0aSStefano Zampini #include <petsc/private/pcbddcstructsimpl.h>
4674ae819SStefano Zampini 
5*32fe681dSStefano Zampini PetscErrorCode PCBDDCDestroyGraphCandidatesIS(void *ctx)
6*32fe681dSStefano Zampini {
7*32fe681dSStefano Zampini   PCBDDCGraphCandidates cand = (PCBDDCGraphCandidates)ctx;
8*32fe681dSStefano Zampini 
9*32fe681dSStefano Zampini   PetscFunctionBegin;
10*32fe681dSStefano Zampini   for (PetscInt i=0;i<cand->nfc;i++) PetscCall(ISDestroy(&cand->Faces[i]));
11*32fe681dSStefano Zampini   for (PetscInt i=0;i<cand->nec;i++) PetscCall(ISDestroy(&cand->Edges[i]));
12*32fe681dSStefano Zampini   PetscCall(PetscFree(cand->Faces));
13*32fe681dSStefano Zampini   PetscCall(PetscFree(cand->Edges));
14*32fe681dSStefano Zampini   PetscCall(ISDestroy(&cand->Vertices));
15*32fe681dSStefano Zampini   PetscCall(PetscFree(cand));
16*32fe681dSStefano Zampini   PetscFunctionReturn(0);
17*32fe681dSStefano Zampini }
18*32fe681dSStefano Zampini 
19d62866d3SStefano Zampini PetscErrorCode PCBDDCGraphGetDirichletDofsB(PCBDDCGraph graph, IS* dirdofs)
20d62866d3SStefano Zampini {
21d62866d3SStefano Zampini   PetscFunctionBegin;
22d62866d3SStefano Zampini   if (graph->dirdofsB) {
239566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)graph->dirdofsB));
24d62866d3SStefano Zampini   } else if (graph->has_dirichlet) {
25d62866d3SStefano Zampini     PetscInt i,size;
26d62866d3SStefano Zampini     PetscInt *dirdofs_idxs;
27d62866d3SStefano Zampini 
28d62866d3SStefano Zampini     size = 0;
29d62866d3SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
30d62866d3SStefano Zampini       if (graph->count[i] && graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) size++;
31d62866d3SStefano Zampini     }
32d62866d3SStefano Zampini 
339566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(size,&dirdofs_idxs));
34d62866d3SStefano Zampini     size = 0;
35d62866d3SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
36d62866d3SStefano Zampini       if (graph->count[i] && graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) dirdofs_idxs[size++] = i;
37d62866d3SStefano Zampini     }
389566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF,size,dirdofs_idxs,PETSC_OWN_POINTER,&graph->dirdofsB));
399566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)graph->dirdofsB));
40d62866d3SStefano Zampini   }
41d62866d3SStefano Zampini   *dirdofs = graph->dirdofsB;
42d62866d3SStefano Zampini   PetscFunctionReturn(0);
43d62866d3SStefano Zampini }
44d62866d3SStefano Zampini 
451b968477SStefano Zampini PetscErrorCode PCBDDCGraphGetDirichletDofs(PCBDDCGraph graph, IS* dirdofs)
461b968477SStefano Zampini {
471b968477SStefano Zampini   PetscFunctionBegin;
48a07ea27aSStefano Zampini   if (graph->dirdofs) {
499566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)graph->dirdofs));
50a07ea27aSStefano Zampini   } else if (graph->has_dirichlet) {
51a07ea27aSStefano Zampini     PetscInt i,size;
52a07ea27aSStefano Zampini     PetscInt *dirdofs_idxs;
53a07ea27aSStefano Zampini 
541b968477SStefano Zampini     size = 0;
551b968477SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
561b968477SStefano Zampini       if (graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) size++;
571b968477SStefano Zampini     }
581b968477SStefano Zampini 
599566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(size,&dirdofs_idxs));
601b968477SStefano Zampini     size = 0;
611b968477SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
621b968477SStefano Zampini       if (graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) dirdofs_idxs[size++] = i;
631b968477SStefano Zampini     }
649566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)graph->l2gmap),size,dirdofs_idxs,PETSC_OWN_POINTER,&graph->dirdofs));
659566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)graph->dirdofs));
661b968477SStefano Zampini   }
67a07ea27aSStefano Zampini   *dirdofs = graph->dirdofs;
681b968477SStefano Zampini   PetscFunctionReturn(0);
691b968477SStefano Zampini }
701b968477SStefano Zampini 
71e49050b4SStefano Zampini PetscErrorCode PCBDDCGraphASCIIView(PCBDDCGraph graph, PetscInt verbosity_level, PetscViewer viewer)
72674ae819SStefano Zampini {
732b510759SStefano Zampini   PetscInt       i,j,tabs;
7493fb5fd3SStefano Zampini   PetscInt*      queue_in_global_numbering;
75674ae819SStefano Zampini 
76674ae819SStefano Zampini   PetscFunctionBegin;
779566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushSynchronized(viewer));
789566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIGetTab(viewer,&tabs));
799566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"--------------------------------------------------\n"));
809566063dSJacob Faibussowitsch   PetscCall(PetscViewerFlush(viewer));
819566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"Local BDDC graph for subdomain %04d\n",PetscGlobalRank));
8263a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"Number of vertices %" PetscInt_FMT "\n",graph->nvtxs));
8363a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"Number of local subdomains %" PetscInt_FMT "\n",graph->n_local_subs ? graph->n_local_subs : 1));
8463a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"Custom minimal size %" PetscInt_FMT "\n",graph->custom_minimal_size));
85aefa1729SStefano Zampini   if (graph->maxcount != PETSC_MAX_INT) {
8663a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"Max count %" PetscInt_FMT "\n",graph->maxcount));
87aefa1729SStefano Zampini   }
8863a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"Topological two dim? %s (set %s)\n",PetscBools[graph->twodim],PetscBools[graph->twodimset]));
89d4890cceSStefano Zampini   if (verbosity_level > 2) {
90674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
9163a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"%" PetscInt_FMT ":\n",i));
9263a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"   which_dof: %" PetscInt_FMT "\n",graph->which_dof[i]));
9363a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"   special_dof: %" PetscInt_FMT "\n",graph->special_dof[i]));
9463a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"   neighbours: %" PetscInt_FMT "\n",graph->count[i]));
959566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
96674ae819SStefano Zampini       if (graph->count[i]) {
979566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"     set of neighbours:"));
98674ae819SStefano Zampini         for (j=0;j<graph->count[i];j++) {
9963a3b9bcSJacob Faibussowitsch           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer," %" PetscInt_FMT,graph->neighbours_set[i][j]));
100674ae819SStefano Zampini         }
1019566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"\n"));
102674ae819SStefano Zampini       }
1039566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISetTab(viewer,tabs));
1049566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
10551b0f6cfSStefano Zampini       if (graph->mirrors) {
10663a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"   mirrors: %" PetscInt_FMT "\n",graph->mirrors[i]));
10751b0f6cfSStefano Zampini         if (graph->mirrors[i]) {
1089566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
1099566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"     set of mirrors:"));
11051b0f6cfSStefano Zampini           for (j=0;j<graph->mirrors[i];j++) {
11163a3b9bcSJacob Faibussowitsch             PetscCall(PetscViewerASCIISynchronizedPrintf(viewer," %" PetscInt_FMT,graph->mirrors_set[i][j]));
11251b0f6cfSStefano Zampini           }
1139566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"\n"));
1149566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIISetTab(viewer,tabs));
1159566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
11651b0f6cfSStefano Zampini         }
11751b0f6cfSStefano Zampini       }
118d4890cceSStefano Zampini       if (verbosity_level > 3) {
1191633d1f0SStefano Zampini         if (graph->xadj) {
1209566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"   local adj list:"));
1219566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
122674ae819SStefano Zampini           for (j=graph->xadj[i];j<graph->xadj[i+1];j++) {
12363a3b9bcSJacob Faibussowitsch             PetscCall(PetscViewerASCIISynchronizedPrintf(viewer," %" PetscInt_FMT,graph->adjncy[j]));
124674ae819SStefano Zampini           }
1259566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"\n"));
1269566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIISetTab(viewer,tabs));
1279566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
128ec1c413dSStefano Zampini         } else {
1299566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"   no adj info\n"));
130674ae819SStefano Zampini         }
131e49050b4SStefano Zampini       }
132531609faSStefano Zampini       if (graph->n_local_subs) {
13363a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"   local sub id: %" PetscInt_FMT "\n",graph->local_subs[i]));
134531609faSStefano Zampini       }
13563a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"   interface subset id: %" PetscInt_FMT "\n",graph->subset[i]));
136674ae819SStefano Zampini       if (graph->subset[i] && graph->subset_ncc) {
13763a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"   ncc for subset: %" PetscInt_FMT "\n",graph->subset_ncc[graph->subset[i]-1]));
138674ae819SStefano Zampini       }
139674ae819SStefano Zampini     }
140e49050b4SStefano Zampini   }
14163a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"Total number of connected components %" PetscInt_FMT "\n",graph->ncc));
1429566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(graph->cptr[graph->ncc],&queue_in_global_numbering));
1439566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_in_global_numbering));
144674ae819SStefano Zampini   for (i=0;i<graph->ncc;i++) {
1451ce3d396SStefano Zampini     PetscInt node_num=graph->queue[graph->cptr[i]];
1461ce3d396SStefano Zampini     PetscBool printcc = PETSC_FALSE;
14763a3b9bcSJacob Faibussowitsch     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]));
1489566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
149674ae819SStefano Zampini     for (j=0;j<graph->count[node_num];j++) {
15063a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer," %" PetscInt_FMT,graph->neighbours_set[node_num][j]));
151674ae819SStefano Zampini     }
152d4890cceSStefano Zampini     if (verbosity_level > 1) {
1539566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"):"));
154c33722ecSStefano Zampini       if (verbosity_level > 2 || graph->twodim || graph->count[node_num] > 1 || (graph->count[node_num] == 1 && graph->special_dof[node_num] == PCBDDCGRAPH_NEUMANN_MARK)) {
1551ce3d396SStefano Zampini         printcc = PETSC_TRUE;
1561ce3d396SStefano Zampini       }
1571ce3d396SStefano Zampini       if (printcc) {
158674ae819SStefano Zampini         for (j=graph->cptr[i];j<graph->cptr[i+1];j++) {
15963a3b9bcSJacob Faibussowitsch           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer," %" PetscInt_FMT " (%" PetscInt_FMT ")",graph->queue[j],queue_in_global_numbering[j]));
160674ae819SStefano Zampini         }
1611ce3d396SStefano Zampini       }
162d4890cceSStefano Zampini     } else {
1639566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,")"));
164d4890cceSStefano Zampini     }
1659566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"\n"));
1669566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISetTab(viewer,tabs));
1679566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
168674ae819SStefano Zampini   }
1699566063dSJacob Faibussowitsch   PetscCall(PetscFree(queue_in_global_numbering));
1709566063dSJacob Faibussowitsch   PetscCall(PetscViewerFlush(viewer));
171674ae819SStefano Zampini   PetscFunctionReturn(0);
172674ae819SStefano Zampini }
173674ae819SStefano Zampini 
174c8272957SStefano Zampini PetscErrorCode PCBDDCGraphRestoreCandidatesIS(PCBDDCGraph graph, PetscInt *n_faces, IS *FacesIS[], PetscInt *n_edges, IS *EdgesIS[], IS *VerticesIS)
175c8272957SStefano Zampini {
176c8272957SStefano Zampini   PetscInt       i;
177*32fe681dSStefano Zampini   PetscContainer gcand;
178c8272957SStefano Zampini 
179c8272957SStefano Zampini   PetscFunctionBegin;
180*32fe681dSStefano Zampini   PetscCall(PetscObjectQuery((PetscObject)graph->l2gmap,"_PCBDDCGraphCandidatesIS",(PetscObject*)&gcand));
181*32fe681dSStefano Zampini   if (gcand) {
182*32fe681dSStefano Zampini     if (n_faces)    *n_faces = 0;
183*32fe681dSStefano Zampini     if (n_edges)    *n_edges = 0;
184*32fe681dSStefano Zampini     if (FacesIS)    *FacesIS = NULL;
185*32fe681dSStefano Zampini     if (EdgesIS)    *EdgesIS = NULL;
186*32fe681dSStefano Zampini     if (VerticesIS) *VerticesIS = NULL;
187*32fe681dSStefano Zampini   }
188c8272957SStefano Zampini   if (n_faces) {
189c8272957SStefano Zampini     if (FacesIS) {
190c8272957SStefano Zampini       for (i=0;i<*n_faces;i++) {
1919566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&((*FacesIS)[i])));
192c8272957SStefano Zampini       }
1939566063dSJacob Faibussowitsch       PetscCall(PetscFree(*FacesIS));
194c8272957SStefano Zampini     }
195c8272957SStefano Zampini     *n_faces = 0;
196c8272957SStefano Zampini   }
197c8272957SStefano Zampini   if (n_edges) {
198c8272957SStefano Zampini     if (EdgesIS) {
199c8272957SStefano Zampini       for (i=0;i<*n_edges;i++) {
2009566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&((*EdgesIS)[i])));
201c8272957SStefano Zampini       }
2029566063dSJacob Faibussowitsch       PetscCall(PetscFree(*EdgesIS));
203c8272957SStefano Zampini     }
204c8272957SStefano Zampini     *n_edges = 0;
205c8272957SStefano Zampini   }
206c8272957SStefano Zampini   if (VerticesIS) {
2079566063dSJacob Faibussowitsch     PetscCall(ISDestroy(VerticesIS));
208c8272957SStefano Zampini   }
209c8272957SStefano Zampini   PetscFunctionReturn(0);
210c8272957SStefano Zampini }
211c8272957SStefano Zampini 
212a873d5d0SStefano Zampini PetscErrorCode PCBDDCGraphGetCandidatesIS(PCBDDCGraph graph, PetscInt *n_faces, IS *FacesIS[], PetscInt *n_edges, IS *EdgesIS[], IS *VerticesIS)
213674ae819SStefano Zampini {
214674ae819SStefano Zampini   IS             *ISForFaces,*ISForEdges,ISForVertices;
215be12c134Sstefano_zampini   PetscInt       i,nfc,nec,nvc,*idx,*mark;
216*32fe681dSStefano Zampini   PetscContainer gcand;
217674ae819SStefano Zampini 
218674ae819SStefano Zampini   PetscFunctionBegin;
219*32fe681dSStefano Zampini   PetscCall(PetscObjectQuery((PetscObject)graph->l2gmap,"_PCBDDCGraphCandidatesIS",(PetscObject*)&gcand));
220*32fe681dSStefano Zampini   if (gcand) {
221*32fe681dSStefano Zampini     PCBDDCGraphCandidates cand;
222*32fe681dSStefano Zampini 
223*32fe681dSStefano Zampini     PetscCall(PetscContainerGetPointer(gcand,(void**)&cand));
224*32fe681dSStefano Zampini     if (n_faces)       *n_faces = cand->nfc;
225*32fe681dSStefano Zampini     if (FacesIS)       *FacesIS = cand->Faces;
226*32fe681dSStefano Zampini     if (n_edges)       *n_edges = cand->nec;
227*32fe681dSStefano Zampini     if (EdgesIS)       *EdgesIS = cand->Edges;
228*32fe681dSStefano Zampini     if (VerticesIS) *VerticesIS = cand->Vertices;
229*32fe681dSStefano Zampini     PetscFunctionReturn(0);
230*32fe681dSStefano Zampini   }
2319566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(graph->ncc,&mark));
232674ae819SStefano Zampini   /* loop on ccs to evalute number of faces, edges and vertices */
233674ae819SStefano Zampini   nfc = 0;
234674ae819SStefano Zampini   nec = 0;
235674ae819SStefano Zampini   nvc = 0;
236674ae819SStefano Zampini   for (i=0;i<graph->ncc;i++) {
2371e1f2d93SStefano Zampini     PetscInt repdof = graph->queue[graph->cptr[i]];
238be12c134Sstefano_zampini     if (graph->cptr[i+1]-graph->cptr[i] > graph->custom_minimal_size && graph->count[repdof] < graph->maxcount) {
239be12c134Sstefano_zampini       if (!graph->twodim && graph->count[repdof] == 1 && graph->special_dof[repdof] != PCBDDCGRAPH_NEUMANN_MARK) {
240674ae819SStefano Zampini         nfc++;
241be12c134Sstefano_zampini         mark[i] = 2;
242be12c134Sstefano_zampini       } else {
243674ae819SStefano Zampini         nec++;
244be12c134Sstefano_zampini         mark[i] = 1;
245674ae819SStefano Zampini       }
246674ae819SStefano Zampini     } else {
247674ae819SStefano Zampini       nvc += graph->cptr[i+1]-graph->cptr[i];
248674ae819SStefano Zampini     }
249674ae819SStefano Zampini   }
250adfc4fb2SStefano Zampini 
251674ae819SStefano Zampini   /* allocate IS arrays for faces, edges. Vertices need a single index set. */
252cf5a6209SStefano Zampini   if (FacesIS) {
2539566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nfc,&ISForFaces));
254cf5a6209SStefano Zampini   }
255cf5a6209SStefano Zampini   if (EdgesIS) {
2569566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nec,&ISForEdges));
257cf5a6209SStefano Zampini   }
258cf5a6209SStefano Zampini   if (VerticesIS) {
2599566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nvc,&idx));
260cf5a6209SStefano Zampini   }
261cf5a6209SStefano Zampini 
262674ae819SStefano Zampini   /* loop on ccs to compute index sets for faces and edges */
263acc113dbSStefano Zampini   if (!graph->queue_sorted) {
264acc113dbSStefano Zampini     PetscInt *queue_global;
265acc113dbSStefano Zampini 
2669566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(graph->cptr[graph->ncc],&queue_global));
2679566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_global));
268acc113dbSStefano Zampini     for (i=0;i<graph->ncc;i++) {
2699566063dSJacob Faibussowitsch       PetscCall(PetscSortIntWithArray(graph->cptr[i+1]-graph->cptr[i],&queue_global[graph->cptr[i]],&graph->queue[graph->cptr[i]]));
270acc113dbSStefano Zampini     }
2719566063dSJacob Faibussowitsch     PetscCall(PetscFree(queue_global));
272acc113dbSStefano Zampini     graph->queue_sorted = PETSC_TRUE;
273acc113dbSStefano Zampini   }
274674ae819SStefano Zampini   nfc = 0;
275674ae819SStefano Zampini   nec = 0;
276674ae819SStefano Zampini   for (i=0;i<graph->ncc;i++) {
277be12c134Sstefano_zampini     if (mark[i] == 2) {
278cf5a6209SStefano Zampini       if (FacesIS) {
2799566063dSJacob Faibussowitsch         PetscCall(ISCreateGeneral(PETSC_COMM_SELF,graph->cptr[i+1]-graph->cptr[i],&graph->queue[graph->cptr[i]],PETSC_USE_POINTER,&ISForFaces[nfc]));
280cf5a6209SStefano Zampini       }
281674ae819SStefano Zampini       nfc++;
282be12c134Sstefano_zampini     } else if (mark[i] == 1) {
283cf5a6209SStefano Zampini       if (EdgesIS) {
2849566063dSJacob Faibussowitsch         PetscCall(ISCreateGeneral(PETSC_COMM_SELF,graph->cptr[i+1]-graph->cptr[i],&graph->queue[graph->cptr[i]],PETSC_USE_POINTER,&ISForEdges[nec]));
285cf5a6209SStefano Zampini       }
286674ae819SStefano Zampini       nec++;
287674ae819SStefano Zampini     }
288674ae819SStefano Zampini   }
289be12c134Sstefano_zampini 
290674ae819SStefano Zampini   /* index set for vertices */
291cf5a6209SStefano Zampini   if (VerticesIS) {
292b8ffe317SStefano Zampini     nvc = 0;
293674ae819SStefano Zampini     for (i=0;i<graph->ncc;i++) {
294be12c134Sstefano_zampini       if (!mark[i]) {
295adfc4fb2SStefano Zampini         PetscInt j;
296adfc4fb2SStefano Zampini 
297674ae819SStefano Zampini         for (j=graph->cptr[i];j<graph->cptr[i+1];j++) {
298674ae819SStefano Zampini           idx[nvc]=graph->queue[j];
299674ae819SStefano Zampini           nvc++;
300674ae819SStefano Zampini         }
301674ae819SStefano Zampini       }
302674ae819SStefano Zampini     }
303674ae819SStefano Zampini     /* sort vertex set (by local ordering) */
3049566063dSJacob Faibussowitsch     PetscCall(PetscSortInt(nvc,idx));
3059566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF,nvc,idx,PETSC_OWN_POINTER,&ISForVertices));
306674ae819SStefano Zampini   }
3079566063dSJacob Faibussowitsch   PetscCall(PetscFree(mark));
308be12c134Sstefano_zampini 
309674ae819SStefano Zampini   /* get back info */
310a873d5d0SStefano Zampini   if (n_faces)       *n_faces = nfc;
311c8272957SStefano Zampini   if (FacesIS)       *FacesIS = ISForFaces;
312a873d5d0SStefano Zampini   if (n_edges)       *n_edges = nec;
313c8272957SStefano Zampini   if (EdgesIS)       *EdgesIS = ISForEdges;
314c8272957SStefano Zampini   if (VerticesIS) *VerticesIS = ISForVertices;
315674ae819SStefano Zampini   PetscFunctionReturn(0);
316674ae819SStefano Zampini }
317674ae819SStefano Zampini 
318674ae819SStefano Zampini PetscErrorCode PCBDDCGraphComputeConnectedComponents(PCBDDCGraph graph)
319674ae819SStefano Zampini {
3204a060362SStefano Zampini   PetscBool      adapt_interface_reduced;
321674ae819SStefano Zampini   MPI_Comm       interface_comm;
3224a060362SStefano Zampini   PetscMPIInt    size;
3238e4af1ccSStefano Zampini   PetscInt       i;
3241c7a958bSStefano Zampini   PetscBT        cornerp;
325674ae819SStefano Zampini 
326674ae819SStefano Zampini   PetscFunctionBegin;
327674ae819SStefano Zampini   /* compute connected components locally */
3289566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)(graph->l2gmap),&interface_comm));
3299566063dSJacob Faibussowitsch   PetscCall(PCBDDCGraphComputeConnectedComponentsLocal(graph));
3301c7a958bSStefano Zampini 
3311c7a958bSStefano Zampini   cornerp = NULL;
3321c7a958bSStefano Zampini   if (graph->active_coords) { /* face based corner selection */
3334f819b78SStefano Zampini     PetscBT   excluded;
3341c7a958bSStefano Zampini     PetscReal *wdist;
3351c7a958bSStefano Zampini     PetscInt  n_neigh,*neigh,*n_shared,**shared;
3361c7a958bSStefano Zampini     PetscInt  maxc, ns;
3371c7a958bSStefano Zampini 
3389566063dSJacob Faibussowitsch     PetscCall(PetscBTCreate(graph->nvtxs,&cornerp));
3399566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetInfo(graph->l2gmap,&n_neigh,&neigh,&n_shared,&shared));
3401c7a958bSStefano Zampini     for (ns = 1, maxc = 0; ns < n_neigh; ns++) maxc = PetscMax(maxc,n_shared[ns]);
3419566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(maxc*graph->cdim,&wdist));
3429566063dSJacob Faibussowitsch     PetscCall(PetscBTCreate(maxc,&excluded));
3431c7a958bSStefano Zampini 
3441c7a958bSStefano Zampini     for (ns = 1; ns < n_neigh; ns++) { /* first proc is self */
3451c7a958bSStefano Zampini       PetscReal *anchor,mdist;
3464f819b78SStefano Zampini       PetscInt  fst,j,k,d,cdim = graph->cdim,n = n_shared[ns];
34751ab8ad6SStefano Zampini       PetscInt  point1,point2,point3,point4;
3481c7a958bSStefano Zampini 
3491c7a958bSStefano Zampini       /* import coordinates on shared interface */
3509566063dSJacob Faibussowitsch       PetscCall(PetscBTMemzero(n,excluded));
3514f819b78SStefano Zampini       for (j=0,fst=-1,k=0;j<n;j++) {
3524f819b78SStefano Zampini         PetscBool skip = PETSC_FALSE;
3534f819b78SStefano Zampini         for (d=0;d<cdim;d++) {
3544f819b78SStefano Zampini           PetscReal c = graph->coords[shared[ns][j]*cdim+d];
3554f819b78SStefano Zampini           skip = (PetscBool)(skip || c == PETSC_MAX_REAL);
3564f819b78SStefano Zampini           wdist[k++] = c;
3574f819b78SStefano Zampini         }
3584f819b78SStefano Zampini         if (skip) {
3599566063dSJacob Faibussowitsch           PetscCall(PetscBTSet(excluded,j));
3604f819b78SStefano Zampini         } else if (fst == -1) fst = j;
3614f819b78SStefano Zampini       }
3624f819b78SStefano Zampini       if (fst == -1) continue;
3631c7a958bSStefano Zampini 
36451ab8ad6SStefano Zampini       /* the dofs are sorted by global numbering, so each rank starts from the same id
36551ab8ad6SStefano Zampini          and it will detect the same corners from the given set */
3661c7a958bSStefano Zampini 
3671c7a958bSStefano Zampini       /* find the farthest point from the starting one */
36851ab8ad6SStefano Zampini       anchor = wdist + fst*cdim;
3691c7a958bSStefano Zampini       mdist  = -1.0;
3704f819b78SStefano Zampini       point1 = fst;
3714f819b78SStefano Zampini       for (j=fst;j<n;j++) {
3721c7a958bSStefano Zampini         PetscReal dist = 0.0;
3731c7a958bSStefano Zampini 
3744f819b78SStefano Zampini         if (PetscUnlikely(PetscBTLookup(excluded,j))) continue;
3751c7a958bSStefano Zampini         for (d=0;d<cdim;d++) dist += (wdist[j*cdim+d]-anchor[d])*(wdist[j*cdim+d]-anchor[d]);
3761c7a958bSStefano Zampini         if (dist > mdist) { mdist = dist; point1 = j; }
3771c7a958bSStefano Zampini       }
3781c7a958bSStefano Zampini 
3791c7a958bSStefano Zampini       /* find the farthest point from point1 */
3801c7a958bSStefano Zampini       anchor = wdist + point1*cdim;
3811c7a958bSStefano Zampini       mdist  = -1.0;
3824f819b78SStefano Zampini       point2 = point1;
3834f819b78SStefano Zampini       for (j=fst;j<n;j++) {
3841c7a958bSStefano Zampini         PetscReal dist = 0.0;
3851c7a958bSStefano Zampini 
3864f819b78SStefano Zampini         if (PetscUnlikely(PetscBTLookup(excluded,j))) continue;
3871c7a958bSStefano Zampini         for (d=0;d<cdim;d++) dist += (wdist[j*cdim+d]-anchor[d])*(wdist[j*cdim+d]-anchor[d]);
3881c7a958bSStefano Zampini         if (dist > mdist) { mdist = dist; point2 = j; }
3891c7a958bSStefano Zampini       }
3901c7a958bSStefano Zampini 
3911c7a958bSStefano Zampini       /* find the third point maximizing the triangle area */
3921c7a958bSStefano Zampini       point3 = point2;
3931c7a958bSStefano Zampini       if (cdim > 2) {
3941c7a958bSStefano Zampini         PetscReal a = 0.0;
3951c7a958bSStefano Zampini 
3961c7a958bSStefano Zampini         for (d=0;d<cdim;d++) a += (wdist[point1*cdim+d]-wdist[point2*cdim+d])*(wdist[point1*cdim+d]-wdist[point2*cdim+d]);
397fff73eb4SStefano Zampini         a = PetscSqrtReal(a);
3981c7a958bSStefano Zampini         mdist = -1.0;
3994f819b78SStefano Zampini         for (j=fst;j<n;j++) {
400fff73eb4SStefano Zampini           PetscReal area,b = 0.0, c = 0.0,s;
4011c7a958bSStefano Zampini 
4024f819b78SStefano Zampini           if (PetscUnlikely(PetscBTLookup(excluded,j))) continue;
4031c7a958bSStefano Zampini           for (d=0;d<cdim;d++) {
4041c7a958bSStefano Zampini             b += (wdist[point1*cdim+d]-wdist[j*cdim+d])*(wdist[point1*cdim+d]-wdist[j*cdim+d]);
4051c7a958bSStefano Zampini             c += (wdist[point2*cdim+d]-wdist[j*cdim+d])*(wdist[point2*cdim+d]-wdist[j*cdim+d]);
4061c7a958bSStefano Zampini           }
407fff73eb4SStefano Zampini           b = PetscSqrtReal(b);
408fff73eb4SStefano Zampini           c = PetscSqrtReal(c);
409fff73eb4SStefano Zampini           s = 0.5*(a+b+c);
4104f819b78SStefano Zampini 
4114f819b78SStefano Zampini           /* Heron's formula, area squared */
4124f819b78SStefano Zampini           area = s*(s-a)*(s-b)*(s-c);
4131c7a958bSStefano Zampini           if (area > mdist) { mdist = area; point3 = j; }
4141c7a958bSStefano Zampini         }
4151c7a958bSStefano Zampini       }
4161c7a958bSStefano Zampini 
41751ab8ad6SStefano Zampini       /* find the farthest point from point3 different from point1 and point2 */
41851ab8ad6SStefano Zampini       anchor = wdist + point3*cdim;
41951ab8ad6SStefano Zampini       mdist  = -1.0;
42051ab8ad6SStefano Zampini       point4 = point3;
42151ab8ad6SStefano Zampini       for (j=fst;j<n;j++) {
42251ab8ad6SStefano Zampini         PetscReal dist = 0.0;
42351ab8ad6SStefano Zampini 
42451ab8ad6SStefano Zampini         if (PetscUnlikely(PetscBTLookup(excluded,j)) || j == point1 || j == point2 || j == point3) continue;
42551ab8ad6SStefano Zampini         for (d=0;d<cdim;d++) dist += (wdist[j*cdim+d]-anchor[d])*(wdist[j*cdim+d]-anchor[d]);
42651ab8ad6SStefano Zampini         if (dist > mdist) { mdist = dist; point4 = j; }
42751ab8ad6SStefano Zampini       }
42851ab8ad6SStefano Zampini 
4299566063dSJacob Faibussowitsch       PetscCall(PetscBTSet(cornerp,shared[ns][point1]));
4309566063dSJacob Faibussowitsch       PetscCall(PetscBTSet(cornerp,shared[ns][point2]));
4319566063dSJacob Faibussowitsch       PetscCall(PetscBTSet(cornerp,shared[ns][point3]));
43251ab8ad6SStefano Zampini       PetscCall(PetscBTSet(cornerp,shared[ns][point4]));
4334f819b78SStefano Zampini 
4341c7a958bSStefano Zampini       /* all dofs having the same coordinates will be primal */
4354f819b78SStefano Zampini       for (j=fst;j<n;j++) {
43651ab8ad6SStefano Zampini         PetscBool same[] = {PETSC_TRUE,PETSC_TRUE,PETSC_TRUE,PETSC_TRUE};
4371c7a958bSStefano Zampini 
4384f819b78SStefano Zampini         if (PetscUnlikely(PetscBTLookup(excluded,j))) continue;
4391c7a958bSStefano Zampini         for (d=0;d<cdim;d++) {
4401c7a958bSStefano Zampini           same[0] = (PetscBool)(same[0] && (PetscAbsReal(wdist[j*cdim + d]-wdist[point1*cdim+d]) < PETSC_SMALL));
4411c7a958bSStefano Zampini           same[1] = (PetscBool)(same[1] && (PetscAbsReal(wdist[j*cdim + d]-wdist[point2*cdim+d]) < PETSC_SMALL));
4421c7a958bSStefano Zampini           same[2] = (PetscBool)(same[2] && (PetscAbsReal(wdist[j*cdim + d]-wdist[point3*cdim+d]) < PETSC_SMALL));
44351ab8ad6SStefano Zampini           same[3] = (PetscBool)(same[3] && (PetscAbsReal(wdist[j*cdim + d]-wdist[point4*cdim+d]) < PETSC_SMALL));
4441c7a958bSStefano Zampini         }
44551ab8ad6SStefano Zampini         if (same[0] || same[1] || same[2] || same[3]) {
4469566063dSJacob Faibussowitsch           PetscCall(PetscBTSet(cornerp,shared[ns][j]));
4471c7a958bSStefano Zampini         }
4481c7a958bSStefano Zampini       }
4491c7a958bSStefano Zampini     }
4509566063dSJacob Faibussowitsch     PetscCall(PetscBTDestroy(&excluded));
4519566063dSJacob Faibussowitsch     PetscCall(PetscFree(wdist));
4529566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingRestoreInfo(graph->l2gmap,&n_neigh,&neigh,&n_shared,&shared));
4531c7a958bSStefano Zampini   }
4541c7a958bSStefano Zampini 
455674ae819SStefano Zampini   /* check consistency of connected components among neighbouring subdomains -> it adapt them in case it is needed */
4569566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(interface_comm,&size));
4574a060362SStefano Zampini   adapt_interface_reduced = PETSC_FALSE;
4584a060362SStefano Zampini   if (size > 1) {
4594a060362SStefano Zampini     PetscInt i;
4601c7a958bSStefano Zampini     PetscBool adapt_interface = cornerp ? PETSC_TRUE : PETSC_FALSE;
4611c7a958bSStefano Zampini     for (i=0;i<graph->n_subsets && !adapt_interface;i++) {
462674ae819SStefano Zampini       /* We are not sure that on a given subset of the local interface,
463674ae819SStefano Zampini          with two connected components, the latters be the same among sharing subdomains */
4641c7a958bSStefano Zampini       if (graph->subset_ncc[i] > 1) adapt_interface = PETSC_TRUE;
465674ae819SStefano Zampini     }
4661c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce(&adapt_interface,&adapt_interface_reduced,1,MPIU_BOOL,MPI_LOR,interface_comm));
4674a060362SStefano Zampini   }
468674ae819SStefano Zampini 
469674ae819SStefano Zampini   if (graph->n_subsets && adapt_interface_reduced) {
470ec1c413dSStefano Zampini     PetscBT     subset_cc_adapt;
471ec1c413dSStefano Zampini     MPI_Request *send_requests,*recv_requests;
472ec1c413dSStefano Zampini     PetscInt    *send_buffer,*recv_buffer;
473ec1c413dSStefano Zampini     PetscInt    sum_requests,start_of_recv,start_of_send;
474ec1c413dSStefano Zampini     PetscInt    *cum_recv_counts;
475ec1c413dSStefano Zampini     PetscInt    *labels;
476b3cdcd63SStefano Zampini     PetscInt    ncc,cum_queue,mss,mns,j,k,s;
4778875e3e6SStefano Zampini     PetscInt    **refine_buffer=NULL,*private_labels = NULL;
47849e09c70SStefano Zampini     PetscBool   *subset_has_corn,*recv_buffer_bool,*send_buffer_bool;
4795b1b9aeaSStefano Zampini 
4809566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(graph->n_subsets,&subset_has_corn));
4811c7a958bSStefano Zampini     if (cornerp) {
4821c7a958bSStefano Zampini       for (i=0;i<graph->n_subsets;i++) {
4831c7a958bSStefano Zampini         for (j=0;j<graph->subset_size[i];j++) {
4841c7a958bSStefano Zampini           if (PetscBTLookup(cornerp,graph->subset_idxs[i][j])) {
4851c7a958bSStefano Zampini             subset_has_corn[i] = PETSC_TRUE;
4861c7a958bSStefano Zampini             break;
4871c7a958bSStefano Zampini           }
4881c7a958bSStefano Zampini         }
4891c7a958bSStefano Zampini       }
4901c7a958bSStefano Zampini     }
4919566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(graph->nvtxs,&labels));
4929566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(labels,graph->nvtxs));
4931c7a958bSStefano Zampini     for (i=0,k=0;i<graph->ncc;i++) {
4941c7a958bSStefano Zampini       PetscInt s = 1;
4951c7a958bSStefano Zampini       for (j=graph->cptr[i];j<graph->cptr[i+1];j++) {
4961c7a958bSStefano Zampini         if (cornerp && PetscBTLookup(cornerp,graph->queue[j])) {
4971c7a958bSStefano Zampini           labels[graph->queue[j]] = k+s;
4981c7a958bSStefano Zampini           s += 1;
4991c7a958bSStefano Zampini         } else {
5001c7a958bSStefano Zampini           labels[graph->queue[j]] = k;
5011c7a958bSStefano Zampini         }
5021c7a958bSStefano Zampini       }
5031c7a958bSStefano Zampini       k += s;
5041c7a958bSStefano Zampini     }
505b3cdcd63SStefano Zampini 
506674ae819SStefano Zampini     /* allocate some space */
5079566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(graph->n_subsets+1,&cum_recv_counts));
5089566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(cum_recv_counts,graph->n_subsets+1));
509b3cdcd63SStefano Zampini 
510674ae819SStefano Zampini     /* first count how many neighbours per connected component I will receive from */
511674ae819SStefano Zampini     cum_recv_counts[0] = 0;
512b3cdcd63SStefano Zampini     for (i=0;i<graph->n_subsets;i++) cum_recv_counts[i+1] = cum_recv_counts[i]+graph->count[graph->subset_idxs[i][0]];
5139566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(graph->n_subsets,&send_buffer_bool));
5149566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(cum_recv_counts[graph->n_subsets],&recv_buffer_bool));
5159566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(cum_recv_counts[graph->n_subsets],&send_requests,cum_recv_counts[graph->n_subsets],&recv_requests));
516674ae819SStefano Zampini     for (i=0;i<cum_recv_counts[graph->n_subsets];i++) {
517674ae819SStefano Zampini       send_requests[i] = MPI_REQUEST_NULL;
518674ae819SStefano Zampini       recv_requests[i] = MPI_REQUEST_NULL;
519674ae819SStefano Zampini     }
520b3cdcd63SStefano Zampini 
521ec1c413dSStefano Zampini     /* exchange with my neighbours the number of my connected components on the subset of interface */
522674ae819SStefano Zampini     sum_requests = 0;
523674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
52449e09c70SStefano Zampini       send_buffer_bool[i] = (PetscBool)(graph->subset_ncc[i] > 1 || subset_has_corn[i]);
52549e09c70SStefano Zampini     }
52649e09c70SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
527ec1c413dSStefano Zampini       PetscMPIInt neigh,tag;
528b3cdcd63SStefano Zampini       PetscInt    count,*neighs;
529ec1c413dSStefano Zampini 
530b3cdcd63SStefano Zampini       count  = graph->count[graph->subset_idxs[i][0]];
531b3cdcd63SStefano Zampini       neighs = graph->neighbours_set[graph->subset_idxs[i][0]];
5329566063dSJacob Faibussowitsch       PetscCall(PetscMPIIntCast(2*graph->subset_ref_node[i],&tag));
533b3cdcd63SStefano Zampini       for (k=0;k<count;k++) {
5341c7a958bSStefano Zampini 
5359566063dSJacob Faibussowitsch         PetscCall(PetscMPIIntCast(neighs[k],&neigh));
5369566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Isend(send_buffer_bool + i,           1,MPIU_BOOL,neigh,tag,interface_comm,&send_requests[sum_requests]));
5379566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Irecv(recv_buffer_bool + sum_requests,1,MPIU_BOOL,neigh,tag,interface_comm,&recv_requests[sum_requests]));
538674ae819SStefano Zampini         sum_requests++;
539674ae819SStefano Zampini       }
540674ae819SStefano Zampini     }
5419566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(sum_requests,recv_requests,MPI_STATUSES_IGNORE));
5429566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(sum_requests,send_requests,MPI_STATUSES_IGNORE));
543b3cdcd63SStefano Zampini 
544b3cdcd63SStefano Zampini     /* determine the subsets I have to adapt (those having more than 1 cc) */
5459566063dSJacob Faibussowitsch     PetscCall(PetscBTCreate(graph->n_subsets,&subset_cc_adapt));
5469566063dSJacob Faibussowitsch     PetscCall(PetscBTMemzero(graph->n_subsets,subset_cc_adapt));
547674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
5481c7a958bSStefano Zampini       if (graph->subset_ncc[i] > 1 || subset_has_corn[i]) {
5499566063dSJacob Faibussowitsch         PetscCall(PetscBTSet(subset_cc_adapt,i));
550b3cdcd63SStefano Zampini         continue;
551b3cdcd63SStefano Zampini       }
552674ae819SStefano Zampini       for (j=cum_recv_counts[i];j<cum_recv_counts[i+1];j++) {
5531c7a958bSStefano Zampini          if (recv_buffer_bool[j]) {
5549566063dSJacob Faibussowitsch           PetscCall(PetscBTSet(subset_cc_adapt,i));
555674ae819SStefano Zampini           break;
556674ae819SStefano Zampini         }
557674ae819SStefano Zampini       }
558674ae819SStefano Zampini     }
5599566063dSJacob Faibussowitsch     PetscCall(PetscFree(send_buffer_bool));
5609566063dSJacob Faibussowitsch     PetscCall(PetscFree(recv_buffer_bool));
5619566063dSJacob Faibussowitsch     PetscCall(PetscFree(subset_has_corn));
562b3cdcd63SStefano Zampini 
563ec1c413dSStefano Zampini     /* determine send/recv buffers sizes */
564ec1c413dSStefano Zampini     j = 0;
565b3cdcd63SStefano Zampini     mss = 0;
566674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
567ec1c413dSStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i)) {
568b3cdcd63SStefano Zampini         j  += graph->subset_size[i];
569b3cdcd63SStefano Zampini         mss = PetscMax(graph->subset_size[i],mss);
570674ae819SStefano Zampini       }
571674ae819SStefano Zampini     }
572ec1c413dSStefano Zampini     k = 0;
573b3cdcd63SStefano Zampini     mns = 0;
574674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
575ec1c413dSStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i)) {
576b3cdcd63SStefano Zampini         k  += (cum_recv_counts[i+1]-cum_recv_counts[i])*graph->subset_size[i];
577b3cdcd63SStefano Zampini         mns = PetscMax(cum_recv_counts[i+1]-cum_recv_counts[i],mns);
578674ae819SStefano Zampini       }
579674ae819SStefano Zampini     }
5809566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(j,&send_buffer,k,&recv_buffer));
581ec1c413dSStefano Zampini 
582b3cdcd63SStefano Zampini     /* fill send buffer (order matters: subset_idxs ordered by global ordering) */
583ec1c413dSStefano Zampini     j = 0;
584b3cdcd63SStefano Zampini     for (i=0;i<graph->n_subsets;i++)
585b3cdcd63SStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i))
586b3cdcd63SStefano Zampini         for (k=0;k<graph->subset_size[i];k++)
587b3cdcd63SStefano Zampini           send_buffer[j++] = labels[graph->subset_idxs[i][k]];
588ec1c413dSStefano Zampini 
589674ae819SStefano Zampini     /* now exchange the data */
590674ae819SStefano Zampini     start_of_recv = 0;
591674ae819SStefano Zampini     start_of_send = 0;
592674ae819SStefano Zampini     sum_requests  = 0;
593674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
594ec1c413dSStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i)) {
595ec1c413dSStefano Zampini         PetscMPIInt neigh,tag;
596b3cdcd63SStefano Zampini         PetscInt    size_of_send = graph->subset_size[i];
597ec1c413dSStefano Zampini 
598b3cdcd63SStefano Zampini         j    = graph->subset_idxs[i][0];
5999566063dSJacob Faibussowitsch         PetscCall(PetscMPIIntCast(2*graph->subset_ref_node[i]+1,&tag));
600674ae819SStefano Zampini         for (k=0;k<graph->count[j];k++) {
6019566063dSJacob Faibussowitsch           PetscCall(PetscMPIIntCast(graph->neighbours_set[j][k],&neigh));
6029566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Isend(&send_buffer[start_of_send],size_of_send,MPIU_INT,neigh,tag,interface_comm,&send_requests[sum_requests]));
6039566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Irecv(&recv_buffer[start_of_recv],size_of_send,MPIU_INT,neigh,tag,interface_comm,&recv_requests[sum_requests]));
604ec1c413dSStefano Zampini           start_of_recv += size_of_send;
605674ae819SStefano Zampini           sum_requests++;
606674ae819SStefano Zampini         }
607674ae819SStefano Zampini         start_of_send += size_of_send;
608674ae819SStefano Zampini       }
609674ae819SStefano Zampini     }
6109566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(sum_requests,recv_requests,MPI_STATUSES_IGNORE));
611d52c4852SStefano Zampini 
612b3cdcd63SStefano Zampini     /* refine connected components */
613674ae819SStefano Zampini     start_of_recv = 0;
614b3cdcd63SStefano Zampini     /* allocate some temporary space */
615b3cdcd63SStefano Zampini     if (mss) {
6169566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(mss,&refine_buffer));
6179566063dSJacob Faibussowitsch       PetscCall(PetscMalloc2(mss*(mns+1),&refine_buffer[0],mss,&private_labels));
618b3cdcd63SStefano Zampini     }
619b3cdcd63SStefano Zampini     ncc = 0;
620b3cdcd63SStefano Zampini     cum_queue = 0;
621b3cdcd63SStefano Zampini     graph->cptr[0] = 0;
622674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
623ec1c413dSStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i)) {
624b3cdcd63SStefano Zampini         PetscInt subset_counter = 0;
625b3cdcd63SStefano Zampini         PetscInt sharingprocs = cum_recv_counts[i+1]-cum_recv_counts[i]+1; /* count myself */
626b3cdcd63SStefano Zampini         PetscInt buffer_size = graph->subset_size[i];
627ec1c413dSStefano Zampini 
628b3cdcd63SStefano Zampini         /* compute pointers */
629b3cdcd63SStefano Zampini         for (j=1;j<buffer_size;j++) refine_buffer[j] = refine_buffer[j-1] + sharingprocs;
630b3cdcd63SStefano Zampini         /* analyze contributions from subdomains that share the i-th subset
6311c7a958bSStefano Zampini            The structure of refine_buffer is suitable to find intersections of ccs among sharingprocs.
632b3cdcd63SStefano Zampini            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)
633b3cdcd63SStefano Zampini            sharing procs connected components:
634ec1c413dSStefano Zampini              neigh 0: [0 1 4], [2 3], labels [4,7]  (2 connected components)
635ec1c413dSStefano Zampini              neigh 1: [0 1], [2 3 4], labels [3 2]  (2 connected components)
636ec1c413dSStefano Zampini              neigh 2: [0 4], [1], [2 3], labels [1 5 6] (3 connected components)
637b3cdcd63SStefano Zampini            refine_buffer will be filled as:
638ec1c413dSStefano Zampini              [ 4, 3, 1;
639ec1c413dSStefano Zampini                4, 2, 1;
640ec1c413dSStefano Zampini                7, 2, 6;
641ec1c413dSStefano Zampini                4, 3, 5;
642ec1c413dSStefano Zampini                7, 2, 6; ];
643b3cdcd63SStefano Zampini            The connected components in local ordering are [0], [1], [2 3], [4] */
644b3cdcd63SStefano Zampini         /* fill temp_buffer */
645b3cdcd63SStefano Zampini         for (k=0;k<buffer_size;k++) refine_buffer[k][0] = labels[graph->subset_idxs[i][k]];
646b3cdcd63SStefano Zampini         for (j=0;j<sharingprocs-1;j++) {
647b3cdcd63SStefano Zampini           for (k=0;k<buffer_size;k++) refine_buffer[k][j+1] = recv_buffer[start_of_recv+k];
648b3cdcd63SStefano Zampini           start_of_recv += buffer_size;
649674ae819SStefano Zampini         }
6509566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(private_labels,buffer_size));
651b3cdcd63SStefano Zampini         for (j=0;j<buffer_size;j++) {
652b3cdcd63SStefano Zampini           if (!private_labels[j]) { /* found a new cc  */
653ec1c413dSStefano Zampini             PetscBool same_set;
654ec1c413dSStefano Zampini 
655b3cdcd63SStefano Zampini             graph->cptr[ncc] = cum_queue;
656b3cdcd63SStefano Zampini             ncc++;
657b3cdcd63SStefano Zampini             subset_counter++;
658b3cdcd63SStefano Zampini             private_labels[j] = subset_counter;
659b3cdcd63SStefano Zampini             graph->queue[cum_queue++] = graph->subset_idxs[i][j];
660b3cdcd63SStefano Zampini             for (k=j+1;k<buffer_size;k++) { /* check for other nodes in new cc */
661674ae819SStefano Zampini               same_set = PETSC_TRUE;
662b3cdcd63SStefano Zampini               for (s=0;s<sharingprocs;s++) {
663b3cdcd63SStefano Zampini                 if (refine_buffer[j][s] != refine_buffer[k][s]) {
664674ae819SStefano Zampini                   same_set = PETSC_FALSE;
665674ae819SStefano Zampini                   break;
666674ae819SStefano Zampini                 }
667674ae819SStefano Zampini               }
668674ae819SStefano Zampini               if (same_set) {
669b3cdcd63SStefano Zampini                 private_labels[k] = subset_counter;
670b3cdcd63SStefano Zampini                 graph->queue[cum_queue++] = graph->subset_idxs[i][k];
671674ae819SStefano Zampini               }
672674ae819SStefano Zampini             }
673674ae819SStefano Zampini           }
674674ae819SStefano Zampini         }
675b3cdcd63SStefano Zampini         graph->cptr[ncc]     = cum_queue;
676b3cdcd63SStefano Zampini         graph->subset_ncc[i] = subset_counter;
677b3cdcd63SStefano Zampini         graph->queue_sorted  = PETSC_FALSE;
678b3cdcd63SStefano Zampini       } else { /* this subset does not need to be adapted */
6799566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(graph->queue+cum_queue,graph->subset_idxs[i],graph->subset_size[i]));
680b3cdcd63SStefano Zampini         ncc++;
681b3cdcd63SStefano Zampini         cum_queue += graph->subset_size[i];
682b3cdcd63SStefano Zampini         graph->cptr[ncc] = cum_queue;
683674ae819SStefano Zampini       }
684674ae819SStefano Zampini     }
685b3cdcd63SStefano Zampini     graph->cptr[ncc] = cum_queue;
686b3cdcd63SStefano Zampini     graph->ncc       = ncc;
687b3cdcd63SStefano Zampini     if (mss) {
6889566063dSJacob Faibussowitsch       PetscCall(PetscFree2(refine_buffer[0],private_labels));
6899566063dSJacob Faibussowitsch       PetscCall(PetscFree(refine_buffer));
690b3cdcd63SStefano Zampini     }
6919566063dSJacob Faibussowitsch     PetscCall(PetscFree(labels));
6929566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(sum_requests,send_requests,MPI_STATUSES_IGNORE));
6939566063dSJacob Faibussowitsch     PetscCall(PetscFree2(send_requests,recv_requests));
6949566063dSJacob Faibussowitsch     PetscCall(PetscFree2(send_buffer,recv_buffer));
6959566063dSJacob Faibussowitsch     PetscCall(PetscFree(cum_recv_counts));
6969566063dSJacob Faibussowitsch     PetscCall(PetscBTDestroy(&subset_cc_adapt));
697674ae819SStefano Zampini   }
6989566063dSJacob Faibussowitsch   PetscCall(PetscBTDestroy(&cornerp));
6998e4af1ccSStefano Zampini 
700cc0a5675Sstefano_zampini   /* Determine if we are in 2D or 3D */
701cc0a5675Sstefano_zampini   if (!graph->twodimset) {
702cc0a5675Sstefano_zampini     PetscBool twodim = PETSC_TRUE;
703cc0a5675Sstefano_zampini     for (i=0;i<graph->ncc;i++) {
704cc0a5675Sstefano_zampini       PetscInt repdof = graph->queue[graph->cptr[i]];
705cc0a5675Sstefano_zampini       PetscInt ccsize = graph->cptr[i+1]-graph->cptr[i];
706cc0a5675Sstefano_zampini       if (graph->count[repdof] > 1 && ccsize > graph->custom_minimal_size) {
707cc0a5675Sstefano_zampini         twodim = PETSC_FALSE;
708cc0a5675Sstefano_zampini         break;
709cc0a5675Sstefano_zampini       }
710cc0a5675Sstefano_zampini     }
7111c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce(&twodim,&graph->twodim,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)graph->l2gmap)));
712cc0a5675Sstefano_zampini     graph->twodimset = PETSC_TRUE;
713cc0a5675Sstefano_zampini   }
714674ae819SStefano Zampini   PetscFunctionReturn(0);
715674ae819SStefano Zampini }
716674ae819SStefano Zampini 
7179fbee547SJacob Faibussowitsch static inline PetscErrorCode PCBDDCGraphComputeCC_Private(PCBDDCGraph graph,PetscInt pid,PetscInt* queue_tip,PetscInt n_prev,PetscInt* n_added)
718b3cdcd63SStefano Zampini {
719b3cdcd63SStefano Zampini   PetscInt       i,j,n;
720b3cdcd63SStefano Zampini   PetscInt       *xadj = graph->xadj,*adjncy = graph->adjncy;
721b3cdcd63SStefano Zampini   PetscBT        touched = graph->touched;
722d895b288Sstefano_zampini   PetscBool      havecsr = (PetscBool)(!!xadj);
723b3cdcd63SStefano Zampini   PetscBool      havesubs = (PetscBool)(!!graph->n_local_subs);
724b3cdcd63SStefano Zampini 
725b3cdcd63SStefano Zampini   PetscFunctionBegin;
726b3cdcd63SStefano Zampini   n = 0;
727b3cdcd63SStefano Zampini   if (havecsr && !havesubs) {
728b3cdcd63SStefano Zampini     for (i=-n_prev;i<0;i++) {
729b3cdcd63SStefano Zampini       PetscInt start_dof = queue_tip[i];
73054fffbccSStefano Zampini       /* we assume that if a dof has a size 1 adjacency list and the corresponding entry is negative, it is connected to all dofs */
73154fffbccSStefano Zampini       if (xadj[start_dof+1]-xadj[start_dof] == 1 && adjncy[xadj[start_dof]] < 0) {
73254fffbccSStefano Zampini         for (j=0;j<graph->subset_size[pid-1];j++) { /* pid \in [1,graph->n_subsets] */
73354fffbccSStefano Zampini           PetscInt dof = graph->subset_idxs[pid-1][j];
73454fffbccSStefano Zampini           if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid) {
7359566063dSJacob Faibussowitsch             PetscCall(PetscBTSet(touched,dof));
73654fffbccSStefano Zampini             queue_tip[n] = dof;
73754fffbccSStefano Zampini             n++;
73854fffbccSStefano Zampini           }
73954fffbccSStefano Zampini         }
74054fffbccSStefano Zampini       } else {
741b3cdcd63SStefano Zampini         for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
742b3cdcd63SStefano Zampini           PetscInt dof = adjncy[j];
743b3cdcd63SStefano Zampini           if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid) {
7449566063dSJacob Faibussowitsch             PetscCall(PetscBTSet(touched,dof));
745b3cdcd63SStefano Zampini             queue_tip[n] = dof;
746b3cdcd63SStefano Zampini             n++;
747b3cdcd63SStefano Zampini           }
748b3cdcd63SStefano Zampini         }
749b3cdcd63SStefano Zampini       }
75054fffbccSStefano Zampini     }
751b3cdcd63SStefano Zampini   } else if (havecsr && havesubs) {
752b3cdcd63SStefano Zampini     PetscInt sid = graph->local_subs[queue_tip[-n_prev]];
753b3cdcd63SStefano Zampini     for (i=-n_prev;i<0;i++) {
754b3cdcd63SStefano Zampini       PetscInt start_dof = queue_tip[i];
75554fffbccSStefano Zampini       /* 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 */
75654fffbccSStefano Zampini       if (xadj[start_dof+1]-xadj[start_dof] == 1 && adjncy[xadj[start_dof]] < 0) {
75754fffbccSStefano Zampini         for (j=0;j<graph->subset_size[pid-1];j++) { /* pid \in [1,graph->n_subsets] */
75854fffbccSStefano Zampini           PetscInt dof = graph->subset_idxs[pid-1][j];
75954fffbccSStefano Zampini           if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid && graph->local_subs[dof] == sid) {
7609566063dSJacob Faibussowitsch             PetscCall(PetscBTSet(touched,dof));
76154fffbccSStefano Zampini             queue_tip[n] = dof;
76254fffbccSStefano Zampini             n++;
76354fffbccSStefano Zampini           }
76454fffbccSStefano Zampini         }
76554fffbccSStefano Zampini       } else {
766b3cdcd63SStefano Zampini         for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
767b3cdcd63SStefano Zampini           PetscInt dof = adjncy[j];
768b3cdcd63SStefano Zampini           if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid && graph->local_subs[dof] == sid) {
7699566063dSJacob Faibussowitsch             PetscCall(PetscBTSet(touched,dof));
770b3cdcd63SStefano Zampini             queue_tip[n] = dof;
771b3cdcd63SStefano Zampini             n++;
772b3cdcd63SStefano Zampini           }
773b3cdcd63SStefano Zampini         }
774b3cdcd63SStefano Zampini       }
77554fffbccSStefano Zampini     }
7761c7a958bSStefano Zampini   } else if (havesubs) { /* sub info only */
777b3cdcd63SStefano Zampini     PetscInt sid = graph->local_subs[queue_tip[-n_prev]];
778b3cdcd63SStefano Zampini     for (j=0;j<graph->subset_size[pid-1];j++) { /* pid \in [1,graph->n_subsets] */
779b3cdcd63SStefano Zampini       PetscInt dof = graph->subset_idxs[pid-1][j];
780b3cdcd63SStefano Zampini       if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid && graph->local_subs[dof] == sid) {
7819566063dSJacob Faibussowitsch         PetscCall(PetscBTSet(touched,dof));
782b3cdcd63SStefano Zampini         queue_tip[n] = dof;
783b3cdcd63SStefano Zampini         n++;
784b3cdcd63SStefano Zampini       }
785b3cdcd63SStefano Zampini     }
7861c7a958bSStefano Zampini   } else {
7871c7a958bSStefano Zampini     for (j=0;j<graph->subset_size[pid-1];j++) { /* pid \in [1,graph->n_subsets] */
7881c7a958bSStefano Zampini       PetscInt dof = graph->subset_idxs[pid-1][j];
7891c7a958bSStefano Zampini       if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid) {
7909566063dSJacob Faibussowitsch         PetscCall(PetscBTSet(touched,dof));
7911c7a958bSStefano Zampini         queue_tip[n] = dof;
7921c7a958bSStefano Zampini         n++;
7931c7a958bSStefano Zampini       }
7941c7a958bSStefano Zampini     }
795b3cdcd63SStefano Zampini   }
796b3cdcd63SStefano Zampini   *n_added = n;
797b3cdcd63SStefano Zampini   PetscFunctionReturn(0);
798b3cdcd63SStefano Zampini }
799b3cdcd63SStefano Zampini 
800674ae819SStefano Zampini PetscErrorCode PCBDDCGraphComputeConnectedComponentsLocal(PCBDDCGraph graph)
801674ae819SStefano Zampini {
802b3cdcd63SStefano Zampini   PetscInt       ncc,cum_queue,n;
803b3cdcd63SStefano Zampini   PetscMPIInt    commsize;
804674ae819SStefano Zampini 
805674ae819SStefano Zampini   PetscFunctionBegin;
80628b400f6SJacob Faibussowitsch   PetscCheck(graph->setupcalled,PetscObjectComm((PetscObject)graph->l2gmap),PETSC_ERR_ORDER,"PCBDDCGraphSetUp should be called first");
807b3cdcd63SStefano Zampini   /* quiet return if there isn't any local info */
8081633d1f0SStefano Zampini   if (!graph->xadj && !graph->n_local_subs) {
809674ae819SStefano Zampini     PetscFunctionReturn(0);
810674ae819SStefano Zampini   }
811674ae819SStefano Zampini 
812674ae819SStefano Zampini   /* reset any previous search of connected components */
8139566063dSJacob Faibussowitsch   PetscCall(PetscBTMemzero(graph->nvtxs,graph->touched));
8149566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)graph->l2gmap),&commsize));
8154b2aedd3SStefano Zampini   if (commsize > graph->commsizelimit) {
816b3cdcd63SStefano Zampini     PetscInt i;
817674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
8180cf82fd6SStefano Zampini       if (graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK || !graph->count[i]) {
8199566063dSJacob Faibussowitsch         PetscCall(PetscBTSet(graph->touched,i));
820674ae819SStefano Zampini       }
8214a060362SStefano Zampini     }
822b3cdcd63SStefano Zampini   }
823674ae819SStefano Zampini 
824674ae819SStefano Zampini   /* begin search for connected components */
825674ae819SStefano Zampini   cum_queue = 0;
826674ae819SStefano Zampini   ncc = 0;
827674ae819SStefano Zampini   for (n=0;n<graph->n_subsets;n++) {
828b3cdcd63SStefano Zampini     PetscInt pid = n+1;  /* partition labeled by 0 is discarded */
829b3cdcd63SStefano Zampini     PetscInt found = 0,prev = 0,first = 0,ncc_pid = 0;
830b3cdcd63SStefano Zampini     while (found != graph->subset_size[n]) {
831b3cdcd63SStefano Zampini       PetscInt added = 0;
832b3cdcd63SStefano Zampini       if (!prev) { /* search for new starting dof */
833b3cdcd63SStefano Zampini         while (PetscBTLookup(graph->touched,graph->subset_idxs[n][first])) first++;
8349566063dSJacob Faibussowitsch         PetscCall(PetscBTSet(graph->touched,graph->subset_idxs[n][first]));
835b3cdcd63SStefano Zampini         graph->queue[cum_queue] = graph->subset_idxs[n][first];
836674ae819SStefano Zampini         graph->cptr[ncc] = cum_queue;
837b3cdcd63SStefano Zampini         prev = 1;
838b3cdcd63SStefano Zampini         cum_queue++;
839b3cdcd63SStefano Zampini         found++;
840674ae819SStefano Zampini         ncc_pid++;
841b3cdcd63SStefano Zampini         ncc++;
842674ae819SStefano Zampini       }
8439566063dSJacob Faibussowitsch       PetscCall(PCBDDCGraphComputeCC_Private(graph,pid,graph->queue + cum_queue,prev,&added));
844b3cdcd63SStefano Zampini       if (!added) {
845674ae819SStefano Zampini         graph->subset_ncc[n] = ncc_pid;
846b3cdcd63SStefano Zampini         graph->cptr[ncc] = cum_queue;
847b3cdcd63SStefano Zampini       }
848b3cdcd63SStefano Zampini       prev = added;
849b3cdcd63SStefano Zampini       found += added;
850b3cdcd63SStefano Zampini       cum_queue += added;
851b3cdcd63SStefano Zampini       if (added && found == graph->subset_size[n]) {
852b3cdcd63SStefano Zampini         graph->subset_ncc[n] = ncc_pid;
853b3cdcd63SStefano Zampini         graph->cptr[ncc] = cum_queue;
854b3cdcd63SStefano Zampini       }
855b3cdcd63SStefano Zampini     }
856674ae819SStefano Zampini   }
857674ae819SStefano Zampini   graph->ncc = ncc;
858acc113dbSStefano Zampini   graph->queue_sorted = PETSC_FALSE;
859674ae819SStefano Zampini   PetscFunctionReturn(0);
860674ae819SStefano Zampini }
861674ae819SStefano Zampini 
862674ae819SStefano Zampini PetscErrorCode PCBDDCGraphSetUp(PCBDDCGraph graph, PetscInt custom_minimal_size, IS neumann_is, IS dirichlet_is, PetscInt n_ISForDofs, IS ISForDofs[], IS custom_primal_vertices)
863674ae819SStefano Zampini {
8647a0e7b2cSstefano_zampini   IS             subset,subset_n;
865674ae819SStefano Zampini   MPI_Comm       comm;
866674ae819SStefano Zampini   const PetscInt *is_indices;
867dc456d91SStefano Zampini   PetscInt       n_neigh,*neigh,*n_shared,**shared,*queue_global;
868674ae819SStefano Zampini   PetscInt       i,j,k,s,total_counts,nodes_touched,is_size;
869b3cdcd63SStefano Zampini   PetscMPIInt    commsize;
870cc0a5675Sstefano_zampini   PetscBool      same_set,mirrors_found;
871674ae819SStefano Zampini 
872674ae819SStefano Zampini   PetscFunctionBegin;
8737a0e7b2cSstefano_zampini   PetscValidLogicalCollectiveInt(graph->l2gmap,custom_minimal_size,2);
8747a0e7b2cSstefano_zampini   if (neumann_is) {
8757a0e7b2cSstefano_zampini     PetscValidHeaderSpecific(neumann_is,IS_CLASSID,3);
8767a0e7b2cSstefano_zampini     PetscCheckSameComm(graph->l2gmap,1,neumann_is,3);
8777a0e7b2cSstefano_zampini   }
878a07ea27aSStefano Zampini   graph->has_dirichlet = PETSC_FALSE;
879a07ea27aSStefano Zampini   if (dirichlet_is) {
8807a0e7b2cSstefano_zampini     PetscValidHeaderSpecific(dirichlet_is,IS_CLASSID,4);
881a07ea27aSStefano Zampini     PetscCheckSameComm(graph->l2gmap,1,dirichlet_is,4);
882a07ea27aSStefano Zampini     graph->has_dirichlet = PETSC_TRUE;
883a07ea27aSStefano Zampini   }
8847a0e7b2cSstefano_zampini   PetscValidLogicalCollectiveInt(graph->l2gmap,n_ISForDofs,5);
8857a0e7b2cSstefano_zampini   for (i=0;i<n_ISForDofs;i++) {
8867a0e7b2cSstefano_zampini     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,6);
8877a0e7b2cSstefano_zampini     PetscCheckSameComm(graph->l2gmap,1,ISForDofs[i],6);
8887a0e7b2cSstefano_zampini   }
8897a0e7b2cSstefano_zampini   if (custom_primal_vertices) {
890ab8c8b98SStefano Zampini     PetscValidHeaderSpecific(custom_primal_vertices,IS_CLASSID,7);
8917a0e7b2cSstefano_zampini     PetscCheckSameComm(graph->l2gmap,1,custom_primal_vertices,7);
8927a0e7b2cSstefano_zampini   }
8939566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)(graph->l2gmap),&comm));
8949566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm,&commsize));
895b3cdcd63SStefano Zampini 
896674ae819SStefano Zampini   /* custom_minimal_size */
89714f95afaSStefano Zampini   graph->custom_minimal_size = custom_minimal_size;
898674ae819SStefano Zampini   /* get info l2gmap and allocate work vectors  */
8999566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetInfo(graph->l2gmap,&n_neigh,&neigh,&n_shared,&shared));
9002635a1d4SStefano Zampini   /* check if we have any local periodic nodes (periodic BCs) */
9012635a1d4SStefano Zampini   mirrors_found = PETSC_FALSE;
9022635a1d4SStefano Zampini   if (graph->nvtxs && n_neigh) {
9032635a1d4SStefano Zampini     for (i=0; i<n_shared[0]; i++) graph->count[shared[0][i]] += 1;
9042635a1d4SStefano Zampini     for (i=0; i<n_shared[0]; i++) {
9052635a1d4SStefano Zampini       if (graph->count[shared[0][i]] > 1) {
9062635a1d4SStefano Zampini         mirrors_found = PETSC_TRUE;
9072635a1d4SStefano Zampini         break;
9082635a1d4SStefano Zampini       }
9092635a1d4SStefano Zampini     }
9102635a1d4SStefano Zampini   }
91151b0f6cfSStefano Zampini   /* compute local mirrors (if any) */
91251b0f6cfSStefano Zampini   if (mirrors_found) {
9137a0e7b2cSstefano_zampini     IS       to,from;
91451b0f6cfSStefano Zampini     PetscInt *local_indices,*global_indices;
9157a0e7b2cSstefano_zampini 
9169566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF,graph->nvtxs,0,1,&to));
9179566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingApplyIS(graph->l2gmap,to,&from));
91851b0f6cfSStefano Zampini     /* get arrays of local and global indices */
9199566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(graph->nvtxs,&local_indices));
9209566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(to,(const PetscInt**)&is_indices));
9219566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(local_indices,is_indices,graph->nvtxs));
9229566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(to,(const PetscInt**)&is_indices));
9239566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(graph->nvtxs,&global_indices));
9249566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(from,(const PetscInt**)&is_indices));
9259566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(global_indices,is_indices,graph->nvtxs));
9269566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(from,(const PetscInt**)&is_indices));
92751b0f6cfSStefano Zampini     /* allocate space for mirrors */
9289566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(graph->nvtxs,&graph->mirrors,graph->nvtxs,&graph->mirrors_set));
9299566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(graph->mirrors,graph->nvtxs));
9300a545947SLisandro Dalcin     graph->mirrors_set[0] = NULL;
93151b0f6cfSStefano Zampini 
93251b0f6cfSStefano Zampini     k=0;
93351b0f6cfSStefano Zampini     for (i=0;i<n_shared[0];i++) {
93451b0f6cfSStefano Zampini       j=shared[0][i];
93551b0f6cfSStefano Zampini       if (graph->count[j] > 1) {
93651b0f6cfSStefano Zampini         graph->mirrors[j]++;
93751b0f6cfSStefano Zampini         k++;
93851b0f6cfSStefano Zampini       }
93951b0f6cfSStefano Zampini     }
94051b0f6cfSStefano Zampini     /* allocate space for set of mirrors */
9419566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(k,&graph->mirrors_set[0]));
94251b0f6cfSStefano Zampini     for (i=1;i<graph->nvtxs;i++)
94351b0f6cfSStefano Zampini       graph->mirrors_set[i]=graph->mirrors_set[i-1]+graph->mirrors[i-1];
94451b0f6cfSStefano Zampini 
94551b0f6cfSStefano Zampini     /* fill arrays */
9469566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(graph->mirrors,graph->nvtxs));
94751b0f6cfSStefano Zampini     for (j=0;j<n_shared[0];j++) {
94851b0f6cfSStefano Zampini       i=shared[0][j];
94951b0f6cfSStefano Zampini       if (graph->count[i] > 1)
95051b0f6cfSStefano Zampini         graph->mirrors_set[i][graph->mirrors[i]++]=global_indices[i];
95151b0f6cfSStefano Zampini     }
9529566063dSJacob Faibussowitsch     PetscCall(PetscSortIntWithArray(graph->nvtxs,global_indices,local_indices));
95351b0f6cfSStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
95451b0f6cfSStefano Zampini       if (graph->mirrors[i] > 0) {
9559566063dSJacob Faibussowitsch         PetscCall(PetscFindInt(graph->mirrors_set[i][0],graph->nvtxs,global_indices,&k));
95651b0f6cfSStefano Zampini         j = global_indices[k];
95751b0f6cfSStefano Zampini         while (k > 0 && global_indices[k-1] == j) k--;
95851b0f6cfSStefano Zampini         for (j=0;j<graph->mirrors[i];j++) {
95951b0f6cfSStefano Zampini           graph->mirrors_set[i][j]=local_indices[k+j];
96051b0f6cfSStefano Zampini         }
9619566063dSJacob Faibussowitsch         PetscCall(PetscSortInt(graph->mirrors[i],graph->mirrors_set[i]));
96251b0f6cfSStefano Zampini       }
96351b0f6cfSStefano Zampini     }
9649566063dSJacob Faibussowitsch     PetscCall(PetscFree(local_indices));
9659566063dSJacob Faibussowitsch     PetscCall(PetscFree(global_indices));
9669566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&to));
9679566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&from));
9687a0e7b2cSstefano_zampini   }
9699566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(graph->count,graph->nvtxs));
97051b0f6cfSStefano Zampini 
971674ae819SStefano Zampini   /* Count total number of neigh per node */
972674ae819SStefano Zampini   k = 0;
973674ae819SStefano Zampini   for (i=1;i<n_neigh;i++) {
974674ae819SStefano Zampini     k += n_shared[i];
975674ae819SStefano Zampini     for (j=0;j<n_shared[i];j++) {
976674ae819SStefano Zampini       graph->count[shared[i][j]] += 1;
977674ae819SStefano Zampini     }
978674ae819SStefano Zampini   }
979674ae819SStefano Zampini   /* Allocate space for storing the set of neighbours for each node */
980674ae819SStefano Zampini   if (graph->nvtxs) {
9819566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(k,&graph->neighbours_set[0]));
982674ae819SStefano Zampini   }
983674ae819SStefano Zampini   for (i=1;i<graph->nvtxs;i++) { /* dont count myself */
984674ae819SStefano Zampini     graph->neighbours_set[i]=graph->neighbours_set[i-1]+graph->count[i-1];
985674ae819SStefano Zampini   }
986674ae819SStefano Zampini   /* Get information for sharing subdomains */
9879566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(graph->count,graph->nvtxs));
988674ae819SStefano Zampini   for (i=1;i<n_neigh;i++) { /* dont count myself */
989674ae819SStefano Zampini     s = n_shared[i];
990674ae819SStefano Zampini     for (j=0;j<s;j++) {
991674ae819SStefano Zampini       k = shared[i][j];
992674ae819SStefano Zampini       graph->neighbours_set[k][graph->count[k]] = neigh[i];
993674ae819SStefano Zampini       graph->count[k] += 1;
994674ae819SStefano Zampini     }
995674ae819SStefano Zampini   }
996674ae819SStefano Zampini   /* sort set of sharing subdomains */
997674ae819SStefano Zampini   for (i=0;i<graph->nvtxs;i++) {
9989566063dSJacob Faibussowitsch     PetscCall(PetscSortRemoveDupsInt(&graph->count[i],graph->neighbours_set[i]));
999674ae819SStefano Zampini   }
10007fb0e2dbSStefano Zampini   /* free memory allocated by ISLocalToGlobalMappingGetInfo */
10019566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreInfo(graph->l2gmap,&n_neigh,&neigh,&n_shared,&shared));
10027fb0e2dbSStefano Zampini 
100367c9da69SStefano Zampini   /*
100467c9da69SStefano Zampini      Get info for dofs splitting
10055777c63cSStefano Zampini      User can specify just a subset; an additional field is considered as a complementary field
100667c9da69SStefano Zampini   */
10070c85b387SStefano Zampini   for (i=0,k=0;i<n_ISForDofs;i++) {
10080c85b387SStefano Zampini     PetscInt bs;
10090c85b387SStefano Zampini 
10109566063dSJacob Faibussowitsch     PetscCall(ISGetBlockSize(ISForDofs[i],&bs));
10110c85b387SStefano Zampini     k   += bs;
10120c85b387SStefano Zampini   }
10130c85b387SStefano Zampini   for (i=0;i<graph->nvtxs;i++) graph->which_dof[i] = k; /* by default a dof belongs to the complement set */
10140c85b387SStefano Zampini   for (i=0,k=0;i<n_ISForDofs;i++) {
10150c85b387SStefano Zampini     PetscInt bs;
10160c85b387SStefano Zampini 
10179566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(ISForDofs[i],&is_size));
10189566063dSJacob Faibussowitsch     PetscCall(ISGetBlockSize(ISForDofs[i],&bs));
10199566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(ISForDofs[i],(const PetscInt**)&is_indices));
10200c85b387SStefano Zampini     for (j=0;j<is_size/bs;j++) {
10210c85b387SStefano Zampini       PetscInt b;
10220c85b387SStefano Zampini 
10230c85b387SStefano Zampini       for (b=0;b<bs;b++) {
10240c85b387SStefano Zampini         PetscInt jj = bs*j + b;
10250c85b387SStefano Zampini 
10260c85b387SStefano Zampini         if (is_indices[jj] > -1 && is_indices[jj] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
10270c85b387SStefano Zampini           graph->which_dof[is_indices[jj]] = k+b;
10280c85b387SStefano Zampini         }
1029674ae819SStefano Zampini       }
103067c9da69SStefano Zampini     }
10319566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(ISForDofs[i],(const PetscInt**)&is_indices));
10320c85b387SStefano Zampini     k   += bs;
10335777c63cSStefano Zampini   }
10345777c63cSStefano Zampini 
1035674ae819SStefano Zampini   /* Take into account Neumann nodes */
1036674ae819SStefano Zampini   if (neumann_is) {
10379566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(neumann_is,&is_size));
10389566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(neumann_is,(const PetscInt**)&is_indices));
1039674ae819SStefano Zampini     for (i=0;i<is_size;i++) {
1040d50ed60aSStefano Zampini       if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
10417a0e7b2cSstefano_zampini         graph->special_dof[is_indices[i]] = PCBDDCGRAPH_NEUMANN_MARK;
1042674ae819SStefano Zampini       }
10433c629590SStefano Zampini     }
10449566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(neumann_is,(const PetscInt**)&is_indices));
1045674ae819SStefano Zampini   }
10467a0e7b2cSstefano_zampini   /* Take into account Dirichlet nodes (they overwrite any neumann boundary mark previously set) */
1047674ae819SStefano Zampini   if (dirichlet_is) {
10489566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(dirichlet_is,&is_size));
10499566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(dirichlet_is,(const PetscInt**)&is_indices));
1050674ae819SStefano Zampini     for (i=0;i<is_size;i++) {
1051d50ed60aSStefano Zampini       if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
10527a0e7b2cSstefano_zampini         if (commsize > graph->commsizelimit) { /* dirichlet nodes treated as internal */
10539566063dSJacob Faibussowitsch           PetscCall(PetscBTSet(graph->touched,is_indices[i]));
10547a0e7b2cSstefano_zampini           graph->subset[is_indices[i]] = 0;
10557a0e7b2cSstefano_zampini         }
10567a0e7b2cSstefano_zampini         graph->special_dof[is_indices[i]] = PCBDDCGRAPH_DIRICHLET_MARK;
1057674ae819SStefano Zampini       }
10583c629590SStefano Zampini     }
10599566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(dirichlet_is,(const PetscInt**)&is_indices));
10607fb0e2dbSStefano Zampini   }
106108a5cf49SStefano Zampini   /* mark local periodic nodes (if any) and adapt CSR graph (if any) */
106251b0f6cfSStefano Zampini   if (graph->mirrors) {
106351b0f6cfSStefano Zampini     for (i=0;i<graph->nvtxs;i++)
106451b0f6cfSStefano Zampini       if (graph->mirrors[i])
10650cf82fd6SStefano Zampini         graph->special_dof[i] = PCBDDCGRAPH_LOCAL_PERIODIC_MARK;
106651b0f6cfSStefano Zampini 
10671633d1f0SStefano Zampini     if (graph->xadj) {
106808a5cf49SStefano Zampini       PetscInt *new_xadj,*new_adjncy;
106951b0f6cfSStefano Zampini       /* sort CSR graph */
107060d4fc61SSatish Balay       for (i=0;i<graph->nvtxs;i++) {
10719566063dSJacob Faibussowitsch         PetscCall(PetscSortInt(graph->xadj[i+1]-graph->xadj[i],&graph->adjncy[graph->xadj[i]]));
107260d4fc61SSatish Balay       }
107351b0f6cfSStefano Zampini       /* adapt local CSR graph in case of local periodicity */
107451b0f6cfSStefano Zampini       k = 0;
107551b0f6cfSStefano Zampini       for (i=0;i<graph->nvtxs;i++)
107651b0f6cfSStefano Zampini         for (j=graph->xadj[i];j<graph->xadj[i+1];j++)
107751b0f6cfSStefano Zampini           k += graph->mirrors[graph->adjncy[j]];
107851b0f6cfSStefano Zampini 
10799566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(graph->nvtxs+1,&new_xadj));
10809566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(k+graph->xadj[graph->nvtxs],&new_adjncy));
108151b0f6cfSStefano Zampini       new_xadj[0] = 0;
108251b0f6cfSStefano Zampini       for (i=0;i<graph->nvtxs;i++) {
108351b0f6cfSStefano Zampini         k = graph->xadj[i+1]-graph->xadj[i];
10849566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(&new_adjncy[new_xadj[i]],&graph->adjncy[graph->xadj[i]],k));
108551b0f6cfSStefano Zampini         new_xadj[i+1] = new_xadj[i]+k;
108651b0f6cfSStefano Zampini         for (j=graph->xadj[i];j<graph->xadj[i+1];j++) {
108751b0f6cfSStefano Zampini           k = graph->mirrors[graph->adjncy[j]];
10889566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(&new_adjncy[new_xadj[i+1]],graph->mirrors_set[graph->adjncy[j]],k));
108951b0f6cfSStefano Zampini           new_xadj[i+1] += k;
109051b0f6cfSStefano Zampini         }
109151b0f6cfSStefano Zampini         k = new_xadj[i+1]-new_xadj[i];
10929566063dSJacob Faibussowitsch         PetscCall(PetscSortRemoveDupsInt(&k,&new_adjncy[new_xadj[i]]));
109351b0f6cfSStefano Zampini         new_xadj[i+1] = new_xadj[i]+k;
109451b0f6cfSStefano Zampini       }
109551b0f6cfSStefano Zampini       /* set new CSR into graph */
10969566063dSJacob Faibussowitsch       PetscCall(PetscFree(graph->xadj));
10979566063dSJacob Faibussowitsch       PetscCall(PetscFree(graph->adjncy));
109851b0f6cfSStefano Zampini       graph->xadj = new_xadj;
109951b0f6cfSStefano Zampini       graph->adjncy = new_adjncy;
110051b0f6cfSStefano Zampini     }
110108a5cf49SStefano Zampini   }
110251b0f6cfSStefano Zampini 
110317eb1463SStefano Zampini   /* mark special nodes (if any) -> each will become a single node equivalence class */
1104674ae819SStefano Zampini   if (custom_primal_vertices) {
11059566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(custom_primal_vertices,&is_size));
11069566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(custom_primal_vertices,(const PetscInt**)&is_indices));
11077a0e7b2cSstefano_zampini     for (i=0,j=0;i<is_size;i++) {
11087a0e7b2cSstefano_zampini       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 */
11097a0e7b2cSstefano_zampini         graph->special_dof[is_indices[i]] = PCBDDCGRAPH_SPECIAL_MARK-j;
11109b70a373SStefano Zampini         j++;
11119b70a373SStefano Zampini       }
11129b70a373SStefano Zampini     }
11139566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(custom_primal_vertices,(const PetscInt**)&is_indices));
111417eb1463SStefano Zampini   }
11159b70a373SStefano Zampini 
11164b2aedd3SStefano Zampini   /* mark interior nodes (if commsize > graph->commsizelimit) as touched and belonging to partition number 0 */
11174b2aedd3SStefano Zampini   if (commsize > graph->commsizelimit) {
1118674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
1119674ae819SStefano Zampini       if (!graph->count[i]) {
11209566063dSJacob Faibussowitsch         PetscCall(PetscBTSet(graph->touched,i));
1121674ae819SStefano Zampini         graph->subset[i] = 0;
1122674ae819SStefano Zampini       }
1123674ae819SStefano Zampini     }
1124b3cdcd63SStefano Zampini   }
1125b3cdcd63SStefano Zampini 
1126674ae819SStefano Zampini   /* init graph structure and compute default subsets */
1127674ae819SStefano Zampini   nodes_touched = 0;
1128674ae819SStefano Zampini   for (i=0;i<graph->nvtxs;i++) {
1129df48933bSStefano Zampini     if (PetscBTLookup(graph->touched,i)) {
1130674ae819SStefano Zampini       nodes_touched++;
1131674ae819SStefano Zampini     }
1132674ae819SStefano Zampini   }
1133674ae819SStefano Zampini   i = 0;
1134674ae819SStefano Zampini   graph->ncc = 0;
1135674ae819SStefano Zampini   total_counts = 0;
11367babac9bSStefano Zampini 
11377babac9bSStefano Zampini   /* allocated space for queues */
11384b2aedd3SStefano Zampini   if (commsize == graph->commsizelimit) {
11399566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(graph->nvtxs+1,&graph->cptr,graph->nvtxs,&graph->queue));
11407babac9bSStefano Zampini   } else {
11417babac9bSStefano Zampini     PetscInt nused = graph->nvtxs - nodes_touched;
11429566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nused+1,&graph->cptr,nused,&graph->queue));
11437babac9bSStefano Zampini   }
11447babac9bSStefano Zampini 
1145674ae819SStefano Zampini   while (nodes_touched<graph->nvtxs) {
1146674ae819SStefano Zampini     /*  find first untouched node in local ordering */
1147b3cdcd63SStefano Zampini     while (PetscBTLookup(graph->touched,i)) i++;
11489566063dSJacob Faibussowitsch     PetscCall(PetscBTSet(graph->touched,i));
1149674ae819SStefano Zampini     graph->subset[i] = graph->ncc+1;
1150674ae819SStefano Zampini     graph->cptr[graph->ncc] = total_counts;
1151674ae819SStefano Zampini     graph->queue[total_counts] = i;
1152674ae819SStefano Zampini     total_counts++;
1153674ae819SStefano Zampini     nodes_touched++;
1154674ae819SStefano Zampini     /* now find all other nodes having the same set of sharing subdomains */
1155674ae819SStefano Zampini     for (j=i+1;j<graph->nvtxs;j++) {
115674e413f5SStefano Zampini       /* check for same number of sharing subdomains, dof number and same special mark */
1157df48933bSStefano Zampini       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]) {
1158674ae819SStefano Zampini         /* check for same set of sharing subdomains */
1159674ae819SStefano Zampini         same_set = PETSC_TRUE;
1160674ae819SStefano Zampini         for (k=0;k<graph->count[j];k++) {
1161674ae819SStefano Zampini           if (graph->neighbours_set[i][k] != graph->neighbours_set[j][k]) {
1162674ae819SStefano Zampini             same_set = PETSC_FALSE;
1163674ae819SStefano Zampini           }
1164674ae819SStefano Zampini         }
1165ab8c8b98SStefano Zampini         /* I have found a friend of mine */
1166674ae819SStefano Zampini         if (same_set) {
11679566063dSJacob Faibussowitsch           PetscCall(PetscBTSet(graph->touched,j));
1168b3cdcd63SStefano Zampini           graph->subset[j] = graph->ncc+1;
1169674ae819SStefano Zampini           nodes_touched++;
1170674ae819SStefano Zampini           graph->queue[total_counts] = j;
1171674ae819SStefano Zampini           total_counts++;
1172674ae819SStefano Zampini         }
1173674ae819SStefano Zampini       }
1174674ae819SStefano Zampini     }
1175674ae819SStefano Zampini     graph->ncc++;
1176674ae819SStefano Zampini   }
1177b3cdcd63SStefano Zampini   /* 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 */
1178674ae819SStefano Zampini   graph->n_subsets = graph->ncc;
11799566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(graph->n_subsets,&graph->subset_ncc));
1180674ae819SStefano Zampini   for (i=0;i<graph->n_subsets;i++) {
1181674ae819SStefano Zampini     graph->subset_ncc[i] = 1;
1182674ae819SStefano Zampini   }
1183674ae819SStefano Zampini   /* final pointer */
1184674ae819SStefano Zampini   graph->cptr[graph->ncc] = total_counts;
1185ec1c413dSStefano Zampini 
1186b3cdcd63SStefano Zampini   /* For consistency reasons (among neighbours), I need to sort (by global ordering) each connected component */
1187ec1c413dSStefano Zampini   /* Get a reference node (min index in global ordering) for each subset for tagging messages */
11889566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(graph->ncc,&graph->subset_ref_node));
11899566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(graph->cptr[graph->ncc],&queue_global));
11909566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_global));
1191ec1c413dSStefano Zampini   for (j=0;j<graph->ncc;j++) {
11929566063dSJacob Faibussowitsch     PetscCall(PetscSortIntWithArray(graph->cptr[j+1]-graph->cptr[j],&queue_global[graph->cptr[j]],&graph->queue[graph->cptr[j]]));
1193dc456d91SStefano Zampini     graph->subset_ref_node[j] = graph->queue[graph->cptr[j]];
1194f0321c90SStefano Zampini   }
11959566063dSJacob Faibussowitsch   PetscCall(PetscFree(queue_global));
1196acc113dbSStefano Zampini   graph->queue_sorted = PETSC_TRUE;
119748bebe3cSStefano Zampini 
1198b3cdcd63SStefano Zampini   /* save information on subsets (needed when analyzing the connected components) */
11992f566a09SStefano Zampini   if (graph->ncc) {
12009566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(graph->ncc,&graph->subset_size,graph->ncc,&graph->subset_idxs));
12019566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(graph->cptr[graph->ncc],&graph->subset_idxs[0]));
12029566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(graph->subset_idxs[0],graph->cptr[graph->ncc]));
1203ec1c413dSStefano Zampini     for (j=1;j<graph->ncc;j++) {
1204b3cdcd63SStefano Zampini       graph->subset_size[j-1] = graph->cptr[j] - graph->cptr[j-1];
1205b3cdcd63SStefano Zampini       graph->subset_idxs[j] = graph->subset_idxs[j-1] + graph->subset_size[j-1];
1206ec1c413dSStefano Zampini     }
1207b3cdcd63SStefano Zampini     graph->subset_size[graph->ncc-1] = graph->cptr[graph->ncc] - graph->cptr[graph->ncc-1];
12089566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(graph->subset_idxs[0],graph->queue,graph->cptr[graph->ncc]));
1209ec1c413dSStefano Zampini   }
1210b3cdcd63SStefano Zampini 
1211f0321c90SStefano Zampini   /* renumber reference nodes */
12129566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)(graph->l2gmap)),graph->ncc,graph->subset_ref_node,PETSC_COPY_VALUES,&subset_n));
12139566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApplyIS(graph->l2gmap,subset_n,&subset));
12149566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&subset_n));
12159566063dSJacob Faibussowitsch   PetscCall(ISRenumber(subset,NULL,NULL,&subset_n));
12169566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&subset));
12179566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(subset_n,&k));
121863a3b9bcSJacob Faibussowitsch   PetscCheck(k == graph->ncc,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid size of new subset! %" PetscInt_FMT " != %" PetscInt_FMT,k,graph->ncc);
12199566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(subset_n,&is_indices));
12209566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(graph->subset_ref_node,is_indices,graph->ncc));
12219566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(subset_n,&is_indices));
12229566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&subset_n));
1223ec1c413dSStefano Zampini 
1224ec1c413dSStefano Zampini   /* free workspace */
1225b3cdcd63SStefano Zampini   graph->setupcalled = PETSC_TRUE;
1226674ae819SStefano Zampini   PetscFunctionReturn(0);
1227674ae819SStefano Zampini }
1228674ae819SStefano Zampini 
1229ab8c8b98SStefano Zampini PetscErrorCode PCBDDCGraphResetCoords(PCBDDCGraph graph)
1230ab8c8b98SStefano Zampini {
1231ab8c8b98SStefano Zampini   PetscFunctionBegin;
1232ab8c8b98SStefano Zampini   if (!graph) PetscFunctionReturn(0);
12339566063dSJacob Faibussowitsch   PetscCall(PetscFree(graph->coords));
1234ab8c8b98SStefano Zampini   graph->cdim  = 0;
1235ab8c8b98SStefano Zampini   graph->cnloc = 0;
1236ab8c8b98SStefano Zampini   graph->cloc  = PETSC_FALSE;
1237ab8c8b98SStefano Zampini   PetscFunctionReturn(0);
1238ab8c8b98SStefano Zampini }
1239ab8c8b98SStefano Zampini 
1240674ae819SStefano Zampini PetscErrorCode PCBDDCGraphResetCSR(PCBDDCGraph graph)
1241674ae819SStefano Zampini {
1242674ae819SStefano Zampini   PetscFunctionBegin;
124307d44595Sstefano_zampini   if (!graph) PetscFunctionReturn(0);
1244a1dbd327SStefano Zampini   if (graph->freecsr) {
12459566063dSJacob Faibussowitsch     PetscCall(PetscFree(graph->xadj));
12469566063dSJacob Faibussowitsch     PetscCall(PetscFree(graph->adjncy));
1247a1dbd327SStefano Zampini   } else {
1248a1dbd327SStefano Zampini     graph->xadj = NULL;
1249a1dbd327SStefano Zampini     graph->adjncy = NULL;
1250a1dbd327SStefano Zampini   }
1251c8272957SStefano Zampini   graph->freecsr = PETSC_FALSE;
1252575ad6abSStefano Zampini   graph->nvtxs_csr = 0;
1253674ae819SStefano Zampini   PetscFunctionReturn(0);
1254674ae819SStefano Zampini }
1255674ae819SStefano Zampini 
1256674ae819SStefano Zampini PetscErrorCode PCBDDCGraphReset(PCBDDCGraph graph)
1257674ae819SStefano Zampini {
1258674ae819SStefano Zampini   PetscFunctionBegin;
125907d44595Sstefano_zampini   if (!graph) PetscFunctionReturn(0);
12609566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingDestroy(&graph->l2gmap));
12619566063dSJacob Faibussowitsch   PetscCall(PetscFree(graph->subset_ncc));
12629566063dSJacob Faibussowitsch   PetscCall(PetscFree(graph->subset_ref_node));
1263674ae819SStefano Zampini   if (graph->nvtxs) {
12649566063dSJacob Faibussowitsch     PetscCall(PetscFree(graph->neighbours_set[0]));
1265674ae819SStefano Zampini   }
12669566063dSJacob Faibussowitsch   PetscCall(PetscBTDestroy(&graph->touched));
1267d0609cedSBarry Smith   PetscCall(PetscFree5(graph->count,graph->neighbours_set,graph->subset,graph->which_dof,graph->special_dof));
12689566063dSJacob Faibussowitsch   PetscCall(PetscFree2(graph->cptr,graph->queue));
126951b0f6cfSStefano Zampini   if (graph->mirrors) {
12709566063dSJacob Faibussowitsch     PetscCall(PetscFree(graph->mirrors_set[0]));
127151b0f6cfSStefano Zampini   }
12729566063dSJacob Faibussowitsch   PetscCall(PetscFree2(graph->mirrors,graph->mirrors_set));
1273b3cdcd63SStefano Zampini   if (graph->subset_idxs) {
12749566063dSJacob Faibussowitsch     PetscCall(PetscFree(graph->subset_idxs[0]));
1275ec1c413dSStefano Zampini   }
12769566063dSJacob Faibussowitsch   PetscCall(PetscFree2(graph->subset_size,graph->subset_idxs));
12779566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&graph->dirdofs));
12789566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&graph->dirdofsB));
12794f1b2e48SStefano Zampini   if (graph->n_local_subs) {
12809566063dSJacob Faibussowitsch     PetscCall(PetscFree(graph->local_subs));
12814f1b2e48SStefano Zampini   }
1282a07ea27aSStefano Zampini   graph->has_dirichlet       = PETSC_FALSE;
1283cc0a5675Sstefano_zampini   graph->twodimset           = PETSC_FALSE;
1284cc0a5675Sstefano_zampini   graph->twodim              = PETSC_FALSE;
1285674ae819SStefano Zampini   graph->nvtxs               = 0;
12867babac9bSStefano Zampini   graph->nvtxs_global        = 0;
1287674ae819SStefano Zampini   graph->n_subsets           = 0;
1288674ae819SStefano Zampini   graph->custom_minimal_size = 1;
12894f1b2e48SStefano Zampini   graph->n_local_subs        = 0;
1290be12c134Sstefano_zampini   graph->maxcount            = PETSC_MAX_INT;
1291fb597685SStefano Zampini   graph->setupcalled         = PETSC_FALSE;
1292674ae819SStefano Zampini   PetscFunctionReturn(0);
1293674ae819SStefano Zampini }
1294674ae819SStefano Zampini 
1295be12c134Sstefano_zampini PetscErrorCode PCBDDCGraphInit(PCBDDCGraph graph, ISLocalToGlobalMapping l2gmap, PetscInt N, PetscInt maxcount)
1296674ae819SStefano Zampini {
1297674ae819SStefano Zampini   PetscInt       n;
1298674ae819SStefano Zampini 
1299674ae819SStefano Zampini   PetscFunctionBegin;
1300674ae819SStefano Zampini   PetscValidPointer(graph,1);
1301674ae819SStefano Zampini   PetscValidHeaderSpecific(l2gmap,IS_LTOGM_CLASSID,2);
13027babac9bSStefano Zampini   PetscValidLogicalCollectiveInt(l2gmap,N,3);
1303be12c134Sstefano_zampini   PetscValidLogicalCollectiveInt(l2gmap,maxcount,4);
13048e61c736SStefano Zampini   /* raise an error if already allocated */
130528b400f6SJacob Faibussowitsch   PetscCheck(!graph->nvtxs_global,PetscObjectComm((PetscObject)l2gmap),PETSC_ERR_PLIB,"BDDCGraph already initialized");
1306674ae819SStefano Zampini   /* set number of vertices */
13079566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)l2gmap));
1308674ae819SStefano Zampini   graph->l2gmap = l2gmap;
13099566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetSize(l2gmap,&n));
1310674ae819SStefano Zampini   graph->nvtxs = n;
13117fb0e2dbSStefano Zampini   graph->nvtxs_global = N;
1312674ae819SStefano Zampini   /* allocate used space */
13139566063dSJacob Faibussowitsch   PetscCall(PetscBTCreate(graph->nvtxs,&graph->touched));
13149566063dSJacob Faibussowitsch   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));
1315674ae819SStefano Zampini   /* zeroes memory */
13169566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(graph->count,graph->nvtxs));
13179566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(graph->subset,graph->nvtxs));
131863602bcaSStefano Zampini   /* use -1 as a default value for which_dof array */
131963602bcaSStefano Zampini   for (n=0;n<graph->nvtxs;n++) graph->which_dof[n] = -1;
13209566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(graph->special_dof,graph->nvtxs));
1321674ae819SStefano Zampini   /* zeroes first pointer to neighbour set */
1322674ae819SStefano Zampini   if (graph->nvtxs) {
13230a545947SLisandro Dalcin     graph->neighbours_set[0] = NULL;
1324674ae819SStefano Zampini   }
1325674ae819SStefano Zampini   /* zeroes workspace for values of ncc */
13260a545947SLisandro Dalcin   graph->subset_ncc = NULL;
13270a545947SLisandro Dalcin   graph->subset_ref_node = NULL;
1328be12c134Sstefano_zampini   /* maxcount for cc */
1329be12c134Sstefano_zampini   graph->maxcount = maxcount;
1330674ae819SStefano Zampini   PetscFunctionReturn(0);
1331674ae819SStefano Zampini }
1332674ae819SStefano Zampini 
1333674ae819SStefano Zampini PetscErrorCode PCBDDCGraphDestroy(PCBDDCGraph* graph)
1334674ae819SStefano Zampini {
1335674ae819SStefano Zampini   PetscFunctionBegin;
13369566063dSJacob Faibussowitsch   PetscCall(PCBDDCGraphResetCSR(*graph));
13379566063dSJacob Faibussowitsch   PetscCall(PCBDDCGraphResetCoords(*graph));
13389566063dSJacob Faibussowitsch   PetscCall(PCBDDCGraphReset(*graph));
13399566063dSJacob Faibussowitsch   PetscCall(PetscFree(*graph));
1340674ae819SStefano Zampini   PetscFunctionReturn(0);
1341674ae819SStefano Zampini }
1342674ae819SStefano Zampini 
1343674ae819SStefano Zampini PetscErrorCode PCBDDCGraphCreate(PCBDDCGraph *graph)
1344674ae819SStefano Zampini {
1345674ae819SStefano Zampini   PCBDDCGraph    new_graph;
1346674ae819SStefano Zampini 
1347674ae819SStefano Zampini   PetscFunctionBegin;
13489566063dSJacob Faibussowitsch   PetscCall(PetscNew(&new_graph));
1349674ae819SStefano Zampini   new_graph->custom_minimal_size = 1;
13504b2aedd3SStefano Zampini   new_graph->commsizelimit = 1;
1351674ae819SStefano Zampini   *graph = new_graph;
1352674ae819SStefano Zampini   PetscFunctionReturn(0);
1353674ae819SStefano Zampini }
1354