xref: /petsc/src/ksp/pc/impls/bddc/bddcgraph.c (revision 5e5bbd0a99ea6ef6b187abb9b4a6405ccdf6db1e)
1af0996ceSBarry Smith #include <petsc/private/petscimpl.h>
2*5e5bbd0aSStefano Zampini #include <petsc/private/pcbddcprivateimpl.h>
3*5e5bbd0aSStefano Zampini #include <petsc/private/pcbddcstructsimpl.h>
4674ae819SStefano Zampini 
5d62866d3SStefano Zampini PetscErrorCode PCBDDCGraphGetDirichletDofsB(PCBDDCGraph graph, IS* dirdofs)
6d62866d3SStefano Zampini {
7d62866d3SStefano Zampini   PetscFunctionBegin;
8d62866d3SStefano Zampini   if (graph->dirdofsB) {
99566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)graph->dirdofsB));
10d62866d3SStefano Zampini   } else if (graph->has_dirichlet) {
11d62866d3SStefano Zampini     PetscInt i,size;
12d62866d3SStefano Zampini     PetscInt *dirdofs_idxs;
13d62866d3SStefano Zampini 
14d62866d3SStefano Zampini     size = 0;
15d62866d3SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
16d62866d3SStefano Zampini       if (graph->count[i] && graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) size++;
17d62866d3SStefano Zampini     }
18d62866d3SStefano Zampini 
199566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(size,&dirdofs_idxs));
20d62866d3SStefano Zampini     size = 0;
21d62866d3SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
22d62866d3SStefano Zampini       if (graph->count[i] && graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) dirdofs_idxs[size++] = i;
23d62866d3SStefano Zampini     }
249566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF,size,dirdofs_idxs,PETSC_OWN_POINTER,&graph->dirdofsB));
259566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)graph->dirdofsB));
26d62866d3SStefano Zampini   }
27d62866d3SStefano Zampini   *dirdofs = graph->dirdofsB;
28d62866d3SStefano Zampini   PetscFunctionReturn(0);
29d62866d3SStefano Zampini }
30d62866d3SStefano Zampini 
311b968477SStefano Zampini PetscErrorCode PCBDDCGraphGetDirichletDofs(PCBDDCGraph graph, IS* dirdofs)
321b968477SStefano Zampini {
331b968477SStefano Zampini   PetscFunctionBegin;
34a07ea27aSStefano Zampini   if (graph->dirdofs) {
359566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)graph->dirdofs));
36a07ea27aSStefano Zampini   } else if (graph->has_dirichlet) {
37a07ea27aSStefano Zampini     PetscInt i,size;
38a07ea27aSStefano Zampini     PetscInt *dirdofs_idxs;
39a07ea27aSStefano Zampini 
401b968477SStefano Zampini     size = 0;
411b968477SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
421b968477SStefano Zampini       if (graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) size++;
431b968477SStefano Zampini     }
441b968477SStefano Zampini 
459566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(size,&dirdofs_idxs));
461b968477SStefano Zampini     size = 0;
471b968477SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
481b968477SStefano Zampini       if (graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) dirdofs_idxs[size++] = i;
491b968477SStefano Zampini     }
509566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)graph->l2gmap),size,dirdofs_idxs,PETSC_OWN_POINTER,&graph->dirdofs));
519566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)graph->dirdofs));
521b968477SStefano Zampini   }
53a07ea27aSStefano Zampini   *dirdofs = graph->dirdofs;
541b968477SStefano Zampini   PetscFunctionReturn(0);
551b968477SStefano Zampini }
561b968477SStefano Zampini 
57e49050b4SStefano Zampini PetscErrorCode PCBDDCGraphASCIIView(PCBDDCGraph graph, PetscInt verbosity_level, PetscViewer viewer)
58674ae819SStefano Zampini {
592b510759SStefano Zampini   PetscInt       i,j,tabs;
6093fb5fd3SStefano Zampini   PetscInt*      queue_in_global_numbering;
61674ae819SStefano Zampini 
62674ae819SStefano Zampini   PetscFunctionBegin;
639566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushSynchronized(viewer));
649566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIGetTab(viewer,&tabs));
659566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"--------------------------------------------------\n"));
669566063dSJacob Faibussowitsch   PetscCall(PetscViewerFlush(viewer));
679566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"Local BDDC graph for subdomain %04d\n",PetscGlobalRank));
6863a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"Number of vertices %" PetscInt_FMT "\n",graph->nvtxs));
6963a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"Number of local subdomains %" PetscInt_FMT "\n",graph->n_local_subs ? graph->n_local_subs : 1));
7063a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"Custom minimal size %" PetscInt_FMT "\n",graph->custom_minimal_size));
71aefa1729SStefano Zampini   if (graph->maxcount != PETSC_MAX_INT) {
7263a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"Max count %" PetscInt_FMT "\n",graph->maxcount));
73aefa1729SStefano Zampini   }
7463a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"Topological two dim? %s (set %s)\n",PetscBools[graph->twodim],PetscBools[graph->twodimset]));
75d4890cceSStefano Zampini   if (verbosity_level > 2) {
76674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
7763a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"%" PetscInt_FMT ":\n",i));
7863a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"   which_dof: %" PetscInt_FMT "\n",graph->which_dof[i]));
7963a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"   special_dof: %" PetscInt_FMT "\n",graph->special_dof[i]));
8063a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"   neighbours: %" PetscInt_FMT "\n",graph->count[i]));
819566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
82674ae819SStefano Zampini       if (graph->count[i]) {
839566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"     set of neighbours:"));
84674ae819SStefano Zampini         for (j=0;j<graph->count[i];j++) {
8563a3b9bcSJacob Faibussowitsch           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer," %" PetscInt_FMT,graph->neighbours_set[i][j]));
86674ae819SStefano Zampini         }
879566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"\n"));
88674ae819SStefano Zampini       }
899566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISetTab(viewer,tabs));
909566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
9151b0f6cfSStefano Zampini       if (graph->mirrors) {
9263a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"   mirrors: %" PetscInt_FMT "\n",graph->mirrors[i]));
9351b0f6cfSStefano Zampini         if (graph->mirrors[i]) {
949566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
959566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"     set of mirrors:"));
9651b0f6cfSStefano Zampini           for (j=0;j<graph->mirrors[i];j++) {
9763a3b9bcSJacob Faibussowitsch             PetscCall(PetscViewerASCIISynchronizedPrintf(viewer," %" PetscInt_FMT,graph->mirrors_set[i][j]));
9851b0f6cfSStefano Zampini           }
999566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"\n"));
1009566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIISetTab(viewer,tabs));
1019566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
10251b0f6cfSStefano Zampini         }
10351b0f6cfSStefano Zampini       }
104d4890cceSStefano Zampini       if (verbosity_level > 3) {
1051633d1f0SStefano Zampini         if (graph->xadj) {
1069566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"   local adj list:"));
1079566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
108674ae819SStefano Zampini           for (j=graph->xadj[i];j<graph->xadj[i+1];j++) {
10963a3b9bcSJacob Faibussowitsch             PetscCall(PetscViewerASCIISynchronizedPrintf(viewer," %" PetscInt_FMT,graph->adjncy[j]));
110674ae819SStefano Zampini           }
1119566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"\n"));
1129566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIISetTab(viewer,tabs));
1139566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
114ec1c413dSStefano Zampini         } else {
1159566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"   no adj info\n"));
116674ae819SStefano Zampini         }
117e49050b4SStefano Zampini       }
118531609faSStefano Zampini       if (graph->n_local_subs) {
11963a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"   local sub id: %" PetscInt_FMT "\n",graph->local_subs[i]));
120531609faSStefano Zampini       }
12163a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"   interface subset id: %" PetscInt_FMT "\n",graph->subset[i]));
122674ae819SStefano Zampini       if (graph->subset[i] && graph->subset_ncc) {
12363a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"   ncc for subset: %" PetscInt_FMT "\n",graph->subset_ncc[graph->subset[i]-1]));
124674ae819SStefano Zampini       }
125674ae819SStefano Zampini     }
126e49050b4SStefano Zampini   }
12763a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"Total number of connected components %" PetscInt_FMT "\n",graph->ncc));
1289566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(graph->cptr[graph->ncc],&queue_in_global_numbering));
1299566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_in_global_numbering));
130674ae819SStefano Zampini   for (i=0;i<graph->ncc;i++) {
1311ce3d396SStefano Zampini     PetscInt node_num=graph->queue[graph->cptr[i]];
1321ce3d396SStefano Zampini     PetscBool printcc = PETSC_FALSE;
13363a3b9bcSJacob 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]));
1349566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE));
135674ae819SStefano Zampini     for (j=0;j<graph->count[node_num];j++) {
13663a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer," %" PetscInt_FMT,graph->neighbours_set[node_num][j]));
137674ae819SStefano Zampini     }
138d4890cceSStefano Zampini     if (verbosity_level > 1) {
1399566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"):"));
140c33722ecSStefano 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)) {
1411ce3d396SStefano Zampini         printcc = PETSC_TRUE;
1421ce3d396SStefano Zampini       }
1431ce3d396SStefano Zampini       if (printcc) {
144674ae819SStefano Zampini         for (j=graph->cptr[i];j<graph->cptr[i+1];j++) {
14563a3b9bcSJacob Faibussowitsch           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer," %" PetscInt_FMT " (%" PetscInt_FMT ")",graph->queue[j],queue_in_global_numbering[j]));
146674ae819SStefano Zampini         }
1471ce3d396SStefano Zampini       }
148d4890cceSStefano Zampini     } else {
1499566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,")"));
150d4890cceSStefano Zampini     }
1519566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"\n"));
1529566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISetTab(viewer,tabs));
1539566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE));
154674ae819SStefano Zampini   }
1559566063dSJacob Faibussowitsch   PetscCall(PetscFree(queue_in_global_numbering));
1569566063dSJacob Faibussowitsch   PetscCall(PetscViewerFlush(viewer));
157674ae819SStefano Zampini   PetscFunctionReturn(0);
158674ae819SStefano Zampini }
159674ae819SStefano Zampini 
160c8272957SStefano Zampini PetscErrorCode PCBDDCGraphRestoreCandidatesIS(PCBDDCGraph graph, PetscInt *n_faces, IS *FacesIS[], PetscInt *n_edges, IS *EdgesIS[], IS *VerticesIS)
161c8272957SStefano Zampini {
162c8272957SStefano Zampini   PetscInt       i;
163c8272957SStefano Zampini 
164c8272957SStefano Zampini   PetscFunctionBegin;
165c8272957SStefano Zampini   if (n_faces) {
166c8272957SStefano Zampini     if (FacesIS) {
167c8272957SStefano Zampini       for (i=0;i<*n_faces;i++) {
1689566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&((*FacesIS)[i])));
169c8272957SStefano Zampini       }
1709566063dSJacob Faibussowitsch       PetscCall(PetscFree(*FacesIS));
171c8272957SStefano Zampini     }
172c8272957SStefano Zampini     *n_faces = 0;
173c8272957SStefano Zampini   }
174c8272957SStefano Zampini   if (n_edges) {
175c8272957SStefano Zampini     if (EdgesIS) {
176c8272957SStefano Zampini       for (i=0;i<*n_edges;i++) {
1779566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&((*EdgesIS)[i])));
178c8272957SStefano Zampini       }
1799566063dSJacob Faibussowitsch       PetscCall(PetscFree(*EdgesIS));
180c8272957SStefano Zampini     }
181c8272957SStefano Zampini     *n_edges = 0;
182c8272957SStefano Zampini   }
183c8272957SStefano Zampini   if (VerticesIS) {
1849566063dSJacob Faibussowitsch     PetscCall(ISDestroy(VerticesIS));
185c8272957SStefano Zampini   }
186c8272957SStefano Zampini   PetscFunctionReturn(0);
187c8272957SStefano Zampini }
188c8272957SStefano Zampini 
189a873d5d0SStefano Zampini PetscErrorCode PCBDDCGraphGetCandidatesIS(PCBDDCGraph graph, PetscInt *n_faces, IS *FacesIS[], PetscInt *n_edges, IS *EdgesIS[], IS *VerticesIS)
190674ae819SStefano Zampini {
191674ae819SStefano Zampini   IS             *ISForFaces,*ISForEdges,ISForVertices;
192be12c134Sstefano_zampini   PetscInt       i,nfc,nec,nvc,*idx,*mark;
193674ae819SStefano Zampini 
194674ae819SStefano Zampini   PetscFunctionBegin;
1959566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(graph->ncc,&mark));
196674ae819SStefano Zampini   /* loop on ccs to evalute number of faces, edges and vertices */
197674ae819SStefano Zampini   nfc = 0;
198674ae819SStefano Zampini   nec = 0;
199674ae819SStefano Zampini   nvc = 0;
200674ae819SStefano Zampini   for (i=0;i<graph->ncc;i++) {
2011e1f2d93SStefano Zampini     PetscInt repdof = graph->queue[graph->cptr[i]];
202be12c134Sstefano_zampini     if (graph->cptr[i+1]-graph->cptr[i] > graph->custom_minimal_size && graph->count[repdof] < graph->maxcount) {
203be12c134Sstefano_zampini       if (!graph->twodim && graph->count[repdof] == 1 && graph->special_dof[repdof] != PCBDDCGRAPH_NEUMANN_MARK) {
204674ae819SStefano Zampini         nfc++;
205be12c134Sstefano_zampini         mark[i] = 2;
206be12c134Sstefano_zampini       } else {
207674ae819SStefano Zampini         nec++;
208be12c134Sstefano_zampini         mark[i] = 1;
209674ae819SStefano Zampini       }
210674ae819SStefano Zampini     } else {
211674ae819SStefano Zampini       nvc += graph->cptr[i+1]-graph->cptr[i];
212674ae819SStefano Zampini     }
213674ae819SStefano Zampini   }
214adfc4fb2SStefano Zampini 
215674ae819SStefano Zampini   /* allocate IS arrays for faces, edges. Vertices need a single index set. */
216cf5a6209SStefano Zampini   if (FacesIS) {
2179566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nfc,&ISForFaces));
218cf5a6209SStefano Zampini   }
219cf5a6209SStefano Zampini   if (EdgesIS) {
2209566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nec,&ISForEdges));
221cf5a6209SStefano Zampini   }
222cf5a6209SStefano Zampini   if (VerticesIS) {
2239566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nvc,&idx));
224cf5a6209SStefano Zampini   }
225cf5a6209SStefano Zampini 
226674ae819SStefano Zampini   /* loop on ccs to compute index sets for faces and edges */
227acc113dbSStefano Zampini   if (!graph->queue_sorted) {
228acc113dbSStefano Zampini     PetscInt *queue_global;
229acc113dbSStefano Zampini 
2309566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(graph->cptr[graph->ncc],&queue_global));
2319566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_global));
232acc113dbSStefano Zampini     for (i=0;i<graph->ncc;i++) {
2339566063dSJacob Faibussowitsch       PetscCall(PetscSortIntWithArray(graph->cptr[i+1]-graph->cptr[i],&queue_global[graph->cptr[i]],&graph->queue[graph->cptr[i]]));
234acc113dbSStefano Zampini     }
2359566063dSJacob Faibussowitsch     PetscCall(PetscFree(queue_global));
236acc113dbSStefano Zampini     graph->queue_sorted = PETSC_TRUE;
237acc113dbSStefano Zampini   }
238674ae819SStefano Zampini   nfc = 0;
239674ae819SStefano Zampini   nec = 0;
240674ae819SStefano Zampini   for (i=0;i<graph->ncc;i++) {
241be12c134Sstefano_zampini     if (mark[i] == 2) {
242cf5a6209SStefano Zampini       if (FacesIS) {
2439566063dSJacob Faibussowitsch         PetscCall(ISCreateGeneral(PETSC_COMM_SELF,graph->cptr[i+1]-graph->cptr[i],&graph->queue[graph->cptr[i]],PETSC_USE_POINTER,&ISForFaces[nfc]));
244cf5a6209SStefano Zampini       }
245674ae819SStefano Zampini       nfc++;
246be12c134Sstefano_zampini     } else if (mark[i] == 1) {
247cf5a6209SStefano Zampini       if (EdgesIS) {
2489566063dSJacob Faibussowitsch         PetscCall(ISCreateGeneral(PETSC_COMM_SELF,graph->cptr[i+1]-graph->cptr[i],&graph->queue[graph->cptr[i]],PETSC_USE_POINTER,&ISForEdges[nec]));
249cf5a6209SStefano Zampini       }
250674ae819SStefano Zampini       nec++;
251674ae819SStefano Zampini     }
252674ae819SStefano Zampini   }
253be12c134Sstefano_zampini 
254674ae819SStefano Zampini   /* index set for vertices */
255cf5a6209SStefano Zampini   if (VerticesIS) {
256b8ffe317SStefano Zampini     nvc = 0;
257674ae819SStefano Zampini     for (i=0;i<graph->ncc;i++) {
258be12c134Sstefano_zampini       if (!mark[i]) {
259adfc4fb2SStefano Zampini         PetscInt j;
260adfc4fb2SStefano Zampini 
261674ae819SStefano Zampini         for (j=graph->cptr[i];j<graph->cptr[i+1];j++) {
262674ae819SStefano Zampini           idx[nvc]=graph->queue[j];
263674ae819SStefano Zampini           nvc++;
264674ae819SStefano Zampini         }
265674ae819SStefano Zampini       }
266674ae819SStefano Zampini     }
267674ae819SStefano Zampini     /* sort vertex set (by local ordering) */
2689566063dSJacob Faibussowitsch     PetscCall(PetscSortInt(nvc,idx));
2699566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF,nvc,idx,PETSC_OWN_POINTER,&ISForVertices));
270674ae819SStefano Zampini   }
2719566063dSJacob Faibussowitsch   PetscCall(PetscFree(mark));
272be12c134Sstefano_zampini 
273674ae819SStefano Zampini   /* get back info */
274a873d5d0SStefano Zampini   if (n_faces)       *n_faces = nfc;
275c8272957SStefano Zampini   if (FacesIS)       *FacesIS = ISForFaces;
276a873d5d0SStefano Zampini   if (n_edges)       *n_edges = nec;
277c8272957SStefano Zampini   if (EdgesIS)       *EdgesIS = ISForEdges;
278c8272957SStefano Zampini   if (VerticesIS) *VerticesIS = ISForVertices;
279674ae819SStefano Zampini   PetscFunctionReturn(0);
280674ae819SStefano Zampini }
281674ae819SStefano Zampini 
282674ae819SStefano Zampini PetscErrorCode PCBDDCGraphComputeConnectedComponents(PCBDDCGraph graph)
283674ae819SStefano Zampini {
2844a060362SStefano Zampini   PetscBool      adapt_interface_reduced;
285674ae819SStefano Zampini   MPI_Comm       interface_comm;
2864a060362SStefano Zampini   PetscMPIInt    size;
2878e4af1ccSStefano Zampini   PetscInt       i;
2881c7a958bSStefano Zampini   PetscBT        cornerp;
289674ae819SStefano Zampini 
290674ae819SStefano Zampini   PetscFunctionBegin;
291674ae819SStefano Zampini   /* compute connected components locally */
2929566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)(graph->l2gmap),&interface_comm));
2939566063dSJacob Faibussowitsch   PetscCall(PCBDDCGraphComputeConnectedComponentsLocal(graph));
2941c7a958bSStefano Zampini 
2951c7a958bSStefano Zampini   cornerp = NULL;
2961c7a958bSStefano Zampini   if (graph->active_coords) { /* face based corner selection */
2974f819b78SStefano Zampini     PetscBT   excluded;
2981c7a958bSStefano Zampini     PetscReal *wdist;
2991c7a958bSStefano Zampini     PetscInt  n_neigh,*neigh,*n_shared,**shared;
3001c7a958bSStefano Zampini     PetscInt  maxc, ns;
3011c7a958bSStefano Zampini 
3029566063dSJacob Faibussowitsch     PetscCall(PetscBTCreate(graph->nvtxs,&cornerp));
3039566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetInfo(graph->l2gmap,&n_neigh,&neigh,&n_shared,&shared));
3041c7a958bSStefano Zampini     for (ns = 1, maxc = 0; ns < n_neigh; ns++) maxc = PetscMax(maxc,n_shared[ns]);
3059566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(maxc*graph->cdim,&wdist));
3069566063dSJacob Faibussowitsch     PetscCall(PetscBTCreate(maxc,&excluded));
3071c7a958bSStefano Zampini 
3081c7a958bSStefano Zampini     for (ns = 1; ns < n_neigh; ns++) { /* first proc is self */
3091c7a958bSStefano Zampini       PetscReal *anchor,mdist;
3104f819b78SStefano Zampini       PetscInt  fst,j,k,d,cdim = graph->cdim,n = n_shared[ns];
31151ab8ad6SStefano Zampini       PetscInt  point1,point2,point3,point4;
3121c7a958bSStefano Zampini 
3131c7a958bSStefano Zampini       /* import coordinates on shared interface */
3149566063dSJacob Faibussowitsch       PetscCall(PetscBTMemzero(n,excluded));
3154f819b78SStefano Zampini       for (j=0,fst=-1,k=0;j<n;j++) {
3164f819b78SStefano Zampini         PetscBool skip = PETSC_FALSE;
3174f819b78SStefano Zampini         for (d=0;d<cdim;d++) {
3184f819b78SStefano Zampini           PetscReal c = graph->coords[shared[ns][j]*cdim+d];
3194f819b78SStefano Zampini           skip = (PetscBool)(skip || c == PETSC_MAX_REAL);
3204f819b78SStefano Zampini           wdist[k++] = c;
3214f819b78SStefano Zampini         }
3224f819b78SStefano Zampini         if (skip) {
3239566063dSJacob Faibussowitsch           PetscCall(PetscBTSet(excluded,j));
3244f819b78SStefano Zampini         } else if (fst == -1) fst = j;
3254f819b78SStefano Zampini       }
3264f819b78SStefano Zampini       if (fst == -1) continue;
3271c7a958bSStefano Zampini 
32851ab8ad6SStefano Zampini       /* the dofs are sorted by global numbering, so each rank starts from the same id
32951ab8ad6SStefano Zampini          and it will detect the same corners from the given set */
3301c7a958bSStefano Zampini 
3311c7a958bSStefano Zampini       /* find the farthest point from the starting one */
33251ab8ad6SStefano Zampini       anchor = wdist + fst*cdim;
3331c7a958bSStefano Zampini       mdist  = -1.0;
3344f819b78SStefano Zampini       point1 = fst;
3354f819b78SStefano Zampini       for (j=fst;j<n;j++) {
3361c7a958bSStefano Zampini         PetscReal dist = 0.0;
3371c7a958bSStefano Zampini 
3384f819b78SStefano Zampini         if (PetscUnlikely(PetscBTLookup(excluded,j))) continue;
3391c7a958bSStefano Zampini         for (d=0;d<cdim;d++) dist += (wdist[j*cdim+d]-anchor[d])*(wdist[j*cdim+d]-anchor[d]);
3401c7a958bSStefano Zampini         if (dist > mdist) { mdist = dist; point1 = j; }
3411c7a958bSStefano Zampini       }
3421c7a958bSStefano Zampini 
3431c7a958bSStefano Zampini       /* find the farthest point from point1 */
3441c7a958bSStefano Zampini       anchor = wdist + point1*cdim;
3451c7a958bSStefano Zampini       mdist  = -1.0;
3464f819b78SStefano Zampini       point2 = point1;
3474f819b78SStefano Zampini       for (j=fst;j<n;j++) {
3481c7a958bSStefano Zampini         PetscReal dist = 0.0;
3491c7a958bSStefano Zampini 
3504f819b78SStefano Zampini         if (PetscUnlikely(PetscBTLookup(excluded,j))) continue;
3511c7a958bSStefano Zampini         for (d=0;d<cdim;d++) dist += (wdist[j*cdim+d]-anchor[d])*(wdist[j*cdim+d]-anchor[d]);
3521c7a958bSStefano Zampini         if (dist > mdist) { mdist = dist; point2 = j; }
3531c7a958bSStefano Zampini       }
3541c7a958bSStefano Zampini 
3551c7a958bSStefano Zampini       /* find the third point maximizing the triangle area */
3561c7a958bSStefano Zampini       point3 = point2;
3571c7a958bSStefano Zampini       if (cdim > 2) {
3581c7a958bSStefano Zampini         PetscReal a = 0.0;
3591c7a958bSStefano Zampini 
3601c7a958bSStefano Zampini         for (d=0;d<cdim;d++) a += (wdist[point1*cdim+d]-wdist[point2*cdim+d])*(wdist[point1*cdim+d]-wdist[point2*cdim+d]);
361fff73eb4SStefano Zampini         a = PetscSqrtReal(a);
3621c7a958bSStefano Zampini         mdist = -1.0;
3634f819b78SStefano Zampini         for (j=fst;j<n;j++) {
364fff73eb4SStefano Zampini           PetscReal area,b = 0.0, c = 0.0,s;
3651c7a958bSStefano Zampini 
3664f819b78SStefano Zampini           if (PetscUnlikely(PetscBTLookup(excluded,j))) continue;
3671c7a958bSStefano Zampini           for (d=0;d<cdim;d++) {
3681c7a958bSStefano Zampini             b += (wdist[point1*cdim+d]-wdist[j*cdim+d])*(wdist[point1*cdim+d]-wdist[j*cdim+d]);
3691c7a958bSStefano Zampini             c += (wdist[point2*cdim+d]-wdist[j*cdim+d])*(wdist[point2*cdim+d]-wdist[j*cdim+d]);
3701c7a958bSStefano Zampini           }
371fff73eb4SStefano Zampini           b = PetscSqrtReal(b);
372fff73eb4SStefano Zampini           c = PetscSqrtReal(c);
373fff73eb4SStefano Zampini           s = 0.5*(a+b+c);
3744f819b78SStefano Zampini 
3754f819b78SStefano Zampini           /* Heron's formula, area squared */
3764f819b78SStefano Zampini           area = s*(s-a)*(s-b)*(s-c);
3771c7a958bSStefano Zampini           if (area > mdist) { mdist = area; point3 = j; }
3781c7a958bSStefano Zampini         }
3791c7a958bSStefano Zampini       }
3801c7a958bSStefano Zampini 
38151ab8ad6SStefano Zampini       /* find the farthest point from point3 different from point1 and point2 */
38251ab8ad6SStefano Zampini       anchor = wdist + point3*cdim;
38351ab8ad6SStefano Zampini       mdist  = -1.0;
38451ab8ad6SStefano Zampini       point4 = point3;
38551ab8ad6SStefano Zampini       for (j=fst;j<n;j++) {
38651ab8ad6SStefano Zampini         PetscReal dist = 0.0;
38751ab8ad6SStefano Zampini 
38851ab8ad6SStefano Zampini         if (PetscUnlikely(PetscBTLookup(excluded,j)) || j == point1 || j == point2 || j == point3) continue;
38951ab8ad6SStefano Zampini         for (d=0;d<cdim;d++) dist += (wdist[j*cdim+d]-anchor[d])*(wdist[j*cdim+d]-anchor[d]);
39051ab8ad6SStefano Zampini         if (dist > mdist) { mdist = dist; point4 = j; }
39151ab8ad6SStefano Zampini       }
39251ab8ad6SStefano Zampini 
3939566063dSJacob Faibussowitsch       PetscCall(PetscBTSet(cornerp,shared[ns][point1]));
3949566063dSJacob Faibussowitsch       PetscCall(PetscBTSet(cornerp,shared[ns][point2]));
3959566063dSJacob Faibussowitsch       PetscCall(PetscBTSet(cornerp,shared[ns][point3]));
39651ab8ad6SStefano Zampini       PetscCall(PetscBTSet(cornerp,shared[ns][point4]));
3974f819b78SStefano Zampini 
3981c7a958bSStefano Zampini       /* all dofs having the same coordinates will be primal */
3994f819b78SStefano Zampini       for (j=fst;j<n;j++) {
40051ab8ad6SStefano Zampini         PetscBool same[] = {PETSC_TRUE,PETSC_TRUE,PETSC_TRUE,PETSC_TRUE};
4011c7a958bSStefano Zampini 
4024f819b78SStefano Zampini         if (PetscUnlikely(PetscBTLookup(excluded,j))) continue;
4031c7a958bSStefano Zampini         for (d=0;d<cdim;d++) {
4041c7a958bSStefano Zampini           same[0] = (PetscBool)(same[0] && (PetscAbsReal(wdist[j*cdim + d]-wdist[point1*cdim+d]) < PETSC_SMALL));
4051c7a958bSStefano Zampini           same[1] = (PetscBool)(same[1] && (PetscAbsReal(wdist[j*cdim + d]-wdist[point2*cdim+d]) < PETSC_SMALL));
4061c7a958bSStefano Zampini           same[2] = (PetscBool)(same[2] && (PetscAbsReal(wdist[j*cdim + d]-wdist[point3*cdim+d]) < PETSC_SMALL));
40751ab8ad6SStefano Zampini           same[3] = (PetscBool)(same[3] && (PetscAbsReal(wdist[j*cdim + d]-wdist[point4*cdim+d]) < PETSC_SMALL));
4081c7a958bSStefano Zampini         }
40951ab8ad6SStefano Zampini         if (same[0] || same[1] || same[2] || same[3]) {
4109566063dSJacob Faibussowitsch           PetscCall(PetscBTSet(cornerp,shared[ns][j]));
4111c7a958bSStefano Zampini         }
4121c7a958bSStefano Zampini       }
4131c7a958bSStefano Zampini     }
4149566063dSJacob Faibussowitsch     PetscCall(PetscBTDestroy(&excluded));
4159566063dSJacob Faibussowitsch     PetscCall(PetscFree(wdist));
4169566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingRestoreInfo(graph->l2gmap,&n_neigh,&neigh,&n_shared,&shared));
4171c7a958bSStefano Zampini   }
4181c7a958bSStefano Zampini 
419674ae819SStefano Zampini   /* check consistency of connected components among neighbouring subdomains -> it adapt them in case it is needed */
4209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(interface_comm,&size));
4214a060362SStefano Zampini   adapt_interface_reduced = PETSC_FALSE;
4224a060362SStefano Zampini   if (size > 1) {
4234a060362SStefano Zampini     PetscInt i;
4241c7a958bSStefano Zampini     PetscBool adapt_interface = cornerp ? PETSC_TRUE : PETSC_FALSE;
4251c7a958bSStefano Zampini     for (i=0;i<graph->n_subsets && !adapt_interface;i++) {
426674ae819SStefano Zampini       /* We are not sure that on a given subset of the local interface,
427674ae819SStefano Zampini          with two connected components, the latters be the same among sharing subdomains */
4281c7a958bSStefano Zampini       if (graph->subset_ncc[i] > 1) adapt_interface = PETSC_TRUE;
429674ae819SStefano Zampini     }
4301c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce(&adapt_interface,&adapt_interface_reduced,1,MPIU_BOOL,MPI_LOR,interface_comm));
4314a060362SStefano Zampini   }
432674ae819SStefano Zampini 
433674ae819SStefano Zampini   if (graph->n_subsets && adapt_interface_reduced) {
434ec1c413dSStefano Zampini     PetscBT     subset_cc_adapt;
435ec1c413dSStefano Zampini     MPI_Request *send_requests,*recv_requests;
436ec1c413dSStefano Zampini     PetscInt    *send_buffer,*recv_buffer;
437ec1c413dSStefano Zampini     PetscInt    sum_requests,start_of_recv,start_of_send;
438ec1c413dSStefano Zampini     PetscInt    *cum_recv_counts;
439ec1c413dSStefano Zampini     PetscInt    *labels;
440b3cdcd63SStefano Zampini     PetscInt    ncc,cum_queue,mss,mns,j,k,s;
4418875e3e6SStefano Zampini     PetscInt    **refine_buffer=NULL,*private_labels = NULL;
44249e09c70SStefano Zampini     PetscBool   *subset_has_corn,*recv_buffer_bool,*send_buffer_bool;
4435b1b9aeaSStefano Zampini 
4449566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(graph->n_subsets,&subset_has_corn));
4451c7a958bSStefano Zampini     if (cornerp) {
4461c7a958bSStefano Zampini       for (i=0;i<graph->n_subsets;i++) {
4471c7a958bSStefano Zampini         for (j=0;j<graph->subset_size[i];j++) {
4481c7a958bSStefano Zampini           if (PetscBTLookup(cornerp,graph->subset_idxs[i][j])) {
4491c7a958bSStefano Zampini             subset_has_corn[i] = PETSC_TRUE;
4501c7a958bSStefano Zampini             break;
4511c7a958bSStefano Zampini           }
4521c7a958bSStefano Zampini         }
4531c7a958bSStefano Zampini       }
4541c7a958bSStefano Zampini     }
4559566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(graph->nvtxs,&labels));
4569566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(labels,graph->nvtxs));
4571c7a958bSStefano Zampini     for (i=0,k=0;i<graph->ncc;i++) {
4581c7a958bSStefano Zampini       PetscInt s = 1;
4591c7a958bSStefano Zampini       for (j=graph->cptr[i];j<graph->cptr[i+1];j++) {
4601c7a958bSStefano Zampini         if (cornerp && PetscBTLookup(cornerp,graph->queue[j])) {
4611c7a958bSStefano Zampini           labels[graph->queue[j]] = k+s;
4621c7a958bSStefano Zampini           s += 1;
4631c7a958bSStefano Zampini         } else {
4641c7a958bSStefano Zampini           labels[graph->queue[j]] = k;
4651c7a958bSStefano Zampini         }
4661c7a958bSStefano Zampini       }
4671c7a958bSStefano Zampini       k += s;
4681c7a958bSStefano Zampini     }
469b3cdcd63SStefano Zampini 
470674ae819SStefano Zampini     /* allocate some space */
4719566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(graph->n_subsets+1,&cum_recv_counts));
4729566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(cum_recv_counts,graph->n_subsets+1));
473b3cdcd63SStefano Zampini 
474674ae819SStefano Zampini     /* first count how many neighbours per connected component I will receive from */
475674ae819SStefano Zampini     cum_recv_counts[0] = 0;
476b3cdcd63SStefano 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]];
4779566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(graph->n_subsets,&send_buffer_bool));
4789566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(cum_recv_counts[graph->n_subsets],&recv_buffer_bool));
4799566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(cum_recv_counts[graph->n_subsets],&send_requests,cum_recv_counts[graph->n_subsets],&recv_requests));
480674ae819SStefano Zampini     for (i=0;i<cum_recv_counts[graph->n_subsets];i++) {
481674ae819SStefano Zampini       send_requests[i] = MPI_REQUEST_NULL;
482674ae819SStefano Zampini       recv_requests[i] = MPI_REQUEST_NULL;
483674ae819SStefano Zampini     }
484b3cdcd63SStefano Zampini 
485ec1c413dSStefano Zampini     /* exchange with my neighbours the number of my connected components on the subset of interface */
486674ae819SStefano Zampini     sum_requests = 0;
487674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
48849e09c70SStefano Zampini       send_buffer_bool[i] = (PetscBool)(graph->subset_ncc[i] > 1 || subset_has_corn[i]);
48949e09c70SStefano Zampini     }
49049e09c70SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
491ec1c413dSStefano Zampini       PetscMPIInt neigh,tag;
492b3cdcd63SStefano Zampini       PetscInt    count,*neighs;
493ec1c413dSStefano Zampini 
494b3cdcd63SStefano Zampini       count  = graph->count[graph->subset_idxs[i][0]];
495b3cdcd63SStefano Zampini       neighs = graph->neighbours_set[graph->subset_idxs[i][0]];
4969566063dSJacob Faibussowitsch       PetscCall(PetscMPIIntCast(2*graph->subset_ref_node[i],&tag));
497b3cdcd63SStefano Zampini       for (k=0;k<count;k++) {
4981c7a958bSStefano Zampini 
4999566063dSJacob Faibussowitsch         PetscCall(PetscMPIIntCast(neighs[k],&neigh));
5009566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Isend(send_buffer_bool + i,           1,MPIU_BOOL,neigh,tag,interface_comm,&send_requests[sum_requests]));
5019566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Irecv(recv_buffer_bool + sum_requests,1,MPIU_BOOL,neigh,tag,interface_comm,&recv_requests[sum_requests]));
502674ae819SStefano Zampini         sum_requests++;
503674ae819SStefano Zampini       }
504674ae819SStefano Zampini     }
5059566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(sum_requests,recv_requests,MPI_STATUSES_IGNORE));
5069566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(sum_requests,send_requests,MPI_STATUSES_IGNORE));
507b3cdcd63SStefano Zampini 
508b3cdcd63SStefano Zampini     /* determine the subsets I have to adapt (those having more than 1 cc) */
5099566063dSJacob Faibussowitsch     PetscCall(PetscBTCreate(graph->n_subsets,&subset_cc_adapt));
5109566063dSJacob Faibussowitsch     PetscCall(PetscBTMemzero(graph->n_subsets,subset_cc_adapt));
511674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
5121c7a958bSStefano Zampini       if (graph->subset_ncc[i] > 1 || subset_has_corn[i]) {
5139566063dSJacob Faibussowitsch         PetscCall(PetscBTSet(subset_cc_adapt,i));
514b3cdcd63SStefano Zampini         continue;
515b3cdcd63SStefano Zampini       }
516674ae819SStefano Zampini       for (j=cum_recv_counts[i];j<cum_recv_counts[i+1];j++) {
5171c7a958bSStefano Zampini          if (recv_buffer_bool[j]) {
5189566063dSJacob Faibussowitsch           PetscCall(PetscBTSet(subset_cc_adapt,i));
519674ae819SStefano Zampini           break;
520674ae819SStefano Zampini         }
521674ae819SStefano Zampini       }
522674ae819SStefano Zampini     }
5239566063dSJacob Faibussowitsch     PetscCall(PetscFree(send_buffer_bool));
5249566063dSJacob Faibussowitsch     PetscCall(PetscFree(recv_buffer_bool));
5259566063dSJacob Faibussowitsch     PetscCall(PetscFree(subset_has_corn));
526b3cdcd63SStefano Zampini 
527ec1c413dSStefano Zampini     /* determine send/recv buffers sizes */
528ec1c413dSStefano Zampini     j = 0;
529b3cdcd63SStefano Zampini     mss = 0;
530674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
531ec1c413dSStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i)) {
532b3cdcd63SStefano Zampini         j  += graph->subset_size[i];
533b3cdcd63SStefano Zampini         mss = PetscMax(graph->subset_size[i],mss);
534674ae819SStefano Zampini       }
535674ae819SStefano Zampini     }
536ec1c413dSStefano Zampini     k = 0;
537b3cdcd63SStefano Zampini     mns = 0;
538674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
539ec1c413dSStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i)) {
540b3cdcd63SStefano Zampini         k  += (cum_recv_counts[i+1]-cum_recv_counts[i])*graph->subset_size[i];
541b3cdcd63SStefano Zampini         mns = PetscMax(cum_recv_counts[i+1]-cum_recv_counts[i],mns);
542674ae819SStefano Zampini       }
543674ae819SStefano Zampini     }
5449566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(j,&send_buffer,k,&recv_buffer));
545ec1c413dSStefano Zampini 
546b3cdcd63SStefano Zampini     /* fill send buffer (order matters: subset_idxs ordered by global ordering) */
547ec1c413dSStefano Zampini     j = 0;
548b3cdcd63SStefano Zampini     for (i=0;i<graph->n_subsets;i++)
549b3cdcd63SStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i))
550b3cdcd63SStefano Zampini         for (k=0;k<graph->subset_size[i];k++)
551b3cdcd63SStefano Zampini           send_buffer[j++] = labels[graph->subset_idxs[i][k]];
552ec1c413dSStefano Zampini 
553674ae819SStefano Zampini     /* now exchange the data */
554674ae819SStefano Zampini     start_of_recv = 0;
555674ae819SStefano Zampini     start_of_send = 0;
556674ae819SStefano Zampini     sum_requests  = 0;
557674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
558ec1c413dSStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i)) {
559ec1c413dSStefano Zampini         PetscMPIInt neigh,tag;
560b3cdcd63SStefano Zampini         PetscInt    size_of_send = graph->subset_size[i];
561ec1c413dSStefano Zampini 
562b3cdcd63SStefano Zampini         j    = graph->subset_idxs[i][0];
5639566063dSJacob Faibussowitsch         PetscCall(PetscMPIIntCast(2*graph->subset_ref_node[i]+1,&tag));
564674ae819SStefano Zampini         for (k=0;k<graph->count[j];k++) {
5659566063dSJacob Faibussowitsch           PetscCall(PetscMPIIntCast(graph->neighbours_set[j][k],&neigh));
5669566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Isend(&send_buffer[start_of_send],size_of_send,MPIU_INT,neigh,tag,interface_comm,&send_requests[sum_requests]));
5679566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Irecv(&recv_buffer[start_of_recv],size_of_send,MPIU_INT,neigh,tag,interface_comm,&recv_requests[sum_requests]));
568ec1c413dSStefano Zampini           start_of_recv += size_of_send;
569674ae819SStefano Zampini           sum_requests++;
570674ae819SStefano Zampini         }
571674ae819SStefano Zampini         start_of_send += size_of_send;
572674ae819SStefano Zampini       }
573674ae819SStefano Zampini     }
5749566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(sum_requests,recv_requests,MPI_STATUSES_IGNORE));
575d52c4852SStefano Zampini 
576b3cdcd63SStefano Zampini     /* refine connected components */
577674ae819SStefano Zampini     start_of_recv = 0;
578b3cdcd63SStefano Zampini     /* allocate some temporary space */
579b3cdcd63SStefano Zampini     if (mss) {
5809566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(mss,&refine_buffer));
5819566063dSJacob Faibussowitsch       PetscCall(PetscMalloc2(mss*(mns+1),&refine_buffer[0],mss,&private_labels));
582b3cdcd63SStefano Zampini     }
583b3cdcd63SStefano Zampini     ncc = 0;
584b3cdcd63SStefano Zampini     cum_queue = 0;
585b3cdcd63SStefano Zampini     graph->cptr[0] = 0;
586674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
587ec1c413dSStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i)) {
588b3cdcd63SStefano Zampini         PetscInt subset_counter = 0;
589b3cdcd63SStefano Zampini         PetscInt sharingprocs = cum_recv_counts[i+1]-cum_recv_counts[i]+1; /* count myself */
590b3cdcd63SStefano Zampini         PetscInt buffer_size = graph->subset_size[i];
591ec1c413dSStefano Zampini 
592b3cdcd63SStefano Zampini         /* compute pointers */
593b3cdcd63SStefano Zampini         for (j=1;j<buffer_size;j++) refine_buffer[j] = refine_buffer[j-1] + sharingprocs;
594b3cdcd63SStefano Zampini         /* analyze contributions from subdomains that share the i-th subset
5951c7a958bSStefano Zampini            The structure of refine_buffer is suitable to find intersections of ccs among sharingprocs.
596b3cdcd63SStefano 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)
597b3cdcd63SStefano Zampini            sharing procs connected components:
598ec1c413dSStefano Zampini              neigh 0: [0 1 4], [2 3], labels [4,7]  (2 connected components)
599ec1c413dSStefano Zampini              neigh 1: [0 1], [2 3 4], labels [3 2]  (2 connected components)
600ec1c413dSStefano Zampini              neigh 2: [0 4], [1], [2 3], labels [1 5 6] (3 connected components)
601b3cdcd63SStefano Zampini            refine_buffer will be filled as:
602ec1c413dSStefano Zampini              [ 4, 3, 1;
603ec1c413dSStefano Zampini                4, 2, 1;
604ec1c413dSStefano Zampini                7, 2, 6;
605ec1c413dSStefano Zampini                4, 3, 5;
606ec1c413dSStefano Zampini                7, 2, 6; ];
607b3cdcd63SStefano Zampini            The connected components in local ordering are [0], [1], [2 3], [4] */
608b3cdcd63SStefano Zampini         /* fill temp_buffer */
609b3cdcd63SStefano Zampini         for (k=0;k<buffer_size;k++) refine_buffer[k][0] = labels[graph->subset_idxs[i][k]];
610b3cdcd63SStefano Zampini         for (j=0;j<sharingprocs-1;j++) {
611b3cdcd63SStefano Zampini           for (k=0;k<buffer_size;k++) refine_buffer[k][j+1] = recv_buffer[start_of_recv+k];
612b3cdcd63SStefano Zampini           start_of_recv += buffer_size;
613674ae819SStefano Zampini         }
6149566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(private_labels,buffer_size));
615b3cdcd63SStefano Zampini         for (j=0;j<buffer_size;j++) {
616b3cdcd63SStefano Zampini           if (!private_labels[j]) { /* found a new cc  */
617ec1c413dSStefano Zampini             PetscBool same_set;
618ec1c413dSStefano Zampini 
619b3cdcd63SStefano Zampini             graph->cptr[ncc] = cum_queue;
620b3cdcd63SStefano Zampini             ncc++;
621b3cdcd63SStefano Zampini             subset_counter++;
622b3cdcd63SStefano Zampini             private_labels[j] = subset_counter;
623b3cdcd63SStefano Zampini             graph->queue[cum_queue++] = graph->subset_idxs[i][j];
624b3cdcd63SStefano Zampini             for (k=j+1;k<buffer_size;k++) { /* check for other nodes in new cc */
625674ae819SStefano Zampini               same_set = PETSC_TRUE;
626b3cdcd63SStefano Zampini               for (s=0;s<sharingprocs;s++) {
627b3cdcd63SStefano Zampini                 if (refine_buffer[j][s] != refine_buffer[k][s]) {
628674ae819SStefano Zampini                   same_set = PETSC_FALSE;
629674ae819SStefano Zampini                   break;
630674ae819SStefano Zampini                 }
631674ae819SStefano Zampini               }
632674ae819SStefano Zampini               if (same_set) {
633b3cdcd63SStefano Zampini                 private_labels[k] = subset_counter;
634b3cdcd63SStefano Zampini                 graph->queue[cum_queue++] = graph->subset_idxs[i][k];
635674ae819SStefano Zampini               }
636674ae819SStefano Zampini             }
637674ae819SStefano Zampini           }
638674ae819SStefano Zampini         }
639b3cdcd63SStefano Zampini         graph->cptr[ncc]     = cum_queue;
640b3cdcd63SStefano Zampini         graph->subset_ncc[i] = subset_counter;
641b3cdcd63SStefano Zampini         graph->queue_sorted  = PETSC_FALSE;
642b3cdcd63SStefano Zampini       } else { /* this subset does not need to be adapted */
6439566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(graph->queue+cum_queue,graph->subset_idxs[i],graph->subset_size[i]));
644b3cdcd63SStefano Zampini         ncc++;
645b3cdcd63SStefano Zampini         cum_queue += graph->subset_size[i];
646b3cdcd63SStefano Zampini         graph->cptr[ncc] = cum_queue;
647674ae819SStefano Zampini       }
648674ae819SStefano Zampini     }
649b3cdcd63SStefano Zampini     graph->cptr[ncc] = cum_queue;
650b3cdcd63SStefano Zampini     graph->ncc       = ncc;
651b3cdcd63SStefano Zampini     if (mss) {
6529566063dSJacob Faibussowitsch       PetscCall(PetscFree2(refine_buffer[0],private_labels));
6539566063dSJacob Faibussowitsch       PetscCall(PetscFree(refine_buffer));
654b3cdcd63SStefano Zampini     }
6559566063dSJacob Faibussowitsch     PetscCall(PetscFree(labels));
6569566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Waitall(sum_requests,send_requests,MPI_STATUSES_IGNORE));
6579566063dSJacob Faibussowitsch     PetscCall(PetscFree2(send_requests,recv_requests));
6589566063dSJacob Faibussowitsch     PetscCall(PetscFree2(send_buffer,recv_buffer));
6599566063dSJacob Faibussowitsch     PetscCall(PetscFree(cum_recv_counts));
6609566063dSJacob Faibussowitsch     PetscCall(PetscBTDestroy(&subset_cc_adapt));
661674ae819SStefano Zampini   }
6629566063dSJacob Faibussowitsch   PetscCall(PetscBTDestroy(&cornerp));
6638e4af1ccSStefano Zampini 
664cc0a5675Sstefano_zampini   /* Determine if we are in 2D or 3D */
665cc0a5675Sstefano_zampini   if (!graph->twodimset) {
666cc0a5675Sstefano_zampini     PetscBool twodim = PETSC_TRUE;
667cc0a5675Sstefano_zampini     for (i=0;i<graph->ncc;i++) {
668cc0a5675Sstefano_zampini       PetscInt repdof = graph->queue[graph->cptr[i]];
669cc0a5675Sstefano_zampini       PetscInt ccsize = graph->cptr[i+1]-graph->cptr[i];
670cc0a5675Sstefano_zampini       if (graph->count[repdof] > 1 && ccsize > graph->custom_minimal_size) {
671cc0a5675Sstefano_zampini         twodim = PETSC_FALSE;
672cc0a5675Sstefano_zampini         break;
673cc0a5675Sstefano_zampini       }
674cc0a5675Sstefano_zampini     }
6751c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce(&twodim,&graph->twodim,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)graph->l2gmap)));
676cc0a5675Sstefano_zampini     graph->twodimset = PETSC_TRUE;
677cc0a5675Sstefano_zampini   }
678674ae819SStefano Zampini   PetscFunctionReturn(0);
679674ae819SStefano Zampini }
680674ae819SStefano Zampini 
6819fbee547SJacob Faibussowitsch static inline PetscErrorCode PCBDDCGraphComputeCC_Private(PCBDDCGraph graph,PetscInt pid,PetscInt* queue_tip,PetscInt n_prev,PetscInt* n_added)
682b3cdcd63SStefano Zampini {
683b3cdcd63SStefano Zampini   PetscInt       i,j,n;
684b3cdcd63SStefano Zampini   PetscInt       *xadj = graph->xadj,*adjncy = graph->adjncy;
685b3cdcd63SStefano Zampini   PetscBT        touched = graph->touched;
686d895b288Sstefano_zampini   PetscBool      havecsr = (PetscBool)(!!xadj);
687b3cdcd63SStefano Zampini   PetscBool      havesubs = (PetscBool)(!!graph->n_local_subs);
688b3cdcd63SStefano Zampini 
689b3cdcd63SStefano Zampini   PetscFunctionBegin;
690b3cdcd63SStefano Zampini   n = 0;
691b3cdcd63SStefano Zampini   if (havecsr && !havesubs) {
692b3cdcd63SStefano Zampini     for (i=-n_prev;i<0;i++) {
693b3cdcd63SStefano Zampini       PetscInt start_dof = queue_tip[i];
69454fffbccSStefano 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 */
69554fffbccSStefano Zampini       if (xadj[start_dof+1]-xadj[start_dof] == 1 && adjncy[xadj[start_dof]] < 0) {
69654fffbccSStefano Zampini         for (j=0;j<graph->subset_size[pid-1];j++) { /* pid \in [1,graph->n_subsets] */
69754fffbccSStefano Zampini           PetscInt dof = graph->subset_idxs[pid-1][j];
69854fffbccSStefano Zampini           if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid) {
6999566063dSJacob Faibussowitsch             PetscCall(PetscBTSet(touched,dof));
70054fffbccSStefano Zampini             queue_tip[n] = dof;
70154fffbccSStefano Zampini             n++;
70254fffbccSStefano Zampini           }
70354fffbccSStefano Zampini         }
70454fffbccSStefano Zampini       } else {
705b3cdcd63SStefano Zampini         for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
706b3cdcd63SStefano Zampini           PetscInt dof = adjncy[j];
707b3cdcd63SStefano Zampini           if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid) {
7089566063dSJacob Faibussowitsch             PetscCall(PetscBTSet(touched,dof));
709b3cdcd63SStefano Zampini             queue_tip[n] = dof;
710b3cdcd63SStefano Zampini             n++;
711b3cdcd63SStefano Zampini           }
712b3cdcd63SStefano Zampini         }
713b3cdcd63SStefano Zampini       }
71454fffbccSStefano Zampini     }
715b3cdcd63SStefano Zampini   } else if (havecsr && havesubs) {
716b3cdcd63SStefano Zampini     PetscInt sid = graph->local_subs[queue_tip[-n_prev]];
717b3cdcd63SStefano Zampini     for (i=-n_prev;i<0;i++) {
718b3cdcd63SStefano Zampini       PetscInt start_dof = queue_tip[i];
71954fffbccSStefano 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 */
72054fffbccSStefano Zampini       if (xadj[start_dof+1]-xadj[start_dof] == 1 && adjncy[xadj[start_dof]] < 0) {
72154fffbccSStefano Zampini         for (j=0;j<graph->subset_size[pid-1];j++) { /* pid \in [1,graph->n_subsets] */
72254fffbccSStefano Zampini           PetscInt dof = graph->subset_idxs[pid-1][j];
72354fffbccSStefano Zampini           if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid && graph->local_subs[dof] == sid) {
7249566063dSJacob Faibussowitsch             PetscCall(PetscBTSet(touched,dof));
72554fffbccSStefano Zampini             queue_tip[n] = dof;
72654fffbccSStefano Zampini             n++;
72754fffbccSStefano Zampini           }
72854fffbccSStefano Zampini         }
72954fffbccSStefano Zampini       } else {
730b3cdcd63SStefano Zampini         for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
731b3cdcd63SStefano Zampini           PetscInt dof = adjncy[j];
732b3cdcd63SStefano Zampini           if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid && graph->local_subs[dof] == sid) {
7339566063dSJacob Faibussowitsch             PetscCall(PetscBTSet(touched,dof));
734b3cdcd63SStefano Zampini             queue_tip[n] = dof;
735b3cdcd63SStefano Zampini             n++;
736b3cdcd63SStefano Zampini           }
737b3cdcd63SStefano Zampini         }
738b3cdcd63SStefano Zampini       }
73954fffbccSStefano Zampini     }
7401c7a958bSStefano Zampini   } else if (havesubs) { /* sub info only */
741b3cdcd63SStefano Zampini     PetscInt sid = graph->local_subs[queue_tip[-n_prev]];
742b3cdcd63SStefano Zampini     for (j=0;j<graph->subset_size[pid-1];j++) { /* pid \in [1,graph->n_subsets] */
743b3cdcd63SStefano Zampini       PetscInt dof = graph->subset_idxs[pid-1][j];
744b3cdcd63SStefano Zampini       if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid && graph->local_subs[dof] == sid) {
7459566063dSJacob Faibussowitsch         PetscCall(PetscBTSet(touched,dof));
746b3cdcd63SStefano Zampini         queue_tip[n] = dof;
747b3cdcd63SStefano Zampini         n++;
748b3cdcd63SStefano Zampini       }
749b3cdcd63SStefano Zampini     }
7501c7a958bSStefano Zampini   } else {
7511c7a958bSStefano Zampini     for (j=0;j<graph->subset_size[pid-1];j++) { /* pid \in [1,graph->n_subsets] */
7521c7a958bSStefano Zampini       PetscInt dof = graph->subset_idxs[pid-1][j];
7531c7a958bSStefano Zampini       if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid) {
7549566063dSJacob Faibussowitsch         PetscCall(PetscBTSet(touched,dof));
7551c7a958bSStefano Zampini         queue_tip[n] = dof;
7561c7a958bSStefano Zampini         n++;
7571c7a958bSStefano Zampini       }
7581c7a958bSStefano Zampini     }
759b3cdcd63SStefano Zampini   }
760b3cdcd63SStefano Zampini   *n_added = n;
761b3cdcd63SStefano Zampini   PetscFunctionReturn(0);
762b3cdcd63SStefano Zampini }
763b3cdcd63SStefano Zampini 
764674ae819SStefano Zampini PetscErrorCode PCBDDCGraphComputeConnectedComponentsLocal(PCBDDCGraph graph)
765674ae819SStefano Zampini {
766b3cdcd63SStefano Zampini   PetscInt       ncc,cum_queue,n;
767b3cdcd63SStefano Zampini   PetscMPIInt    commsize;
768674ae819SStefano Zampini 
769674ae819SStefano Zampini   PetscFunctionBegin;
77028b400f6SJacob Faibussowitsch   PetscCheck(graph->setupcalled,PetscObjectComm((PetscObject)graph->l2gmap),PETSC_ERR_ORDER,"PCBDDCGraphSetUp should be called first");
771b3cdcd63SStefano Zampini   /* quiet return if there isn't any local info */
7721633d1f0SStefano Zampini   if (!graph->xadj && !graph->n_local_subs) {
773674ae819SStefano Zampini     PetscFunctionReturn(0);
774674ae819SStefano Zampini   }
775674ae819SStefano Zampini 
776674ae819SStefano Zampini   /* reset any previous search of connected components */
7779566063dSJacob Faibussowitsch   PetscCall(PetscBTMemzero(graph->nvtxs,graph->touched));
7789566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)graph->l2gmap),&commsize));
7794b2aedd3SStefano Zampini   if (commsize > graph->commsizelimit) {
780b3cdcd63SStefano Zampini     PetscInt i;
781674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
7820cf82fd6SStefano Zampini       if (graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK || !graph->count[i]) {
7839566063dSJacob Faibussowitsch         PetscCall(PetscBTSet(graph->touched,i));
784674ae819SStefano Zampini       }
7854a060362SStefano Zampini     }
786b3cdcd63SStefano Zampini   }
787674ae819SStefano Zampini 
788674ae819SStefano Zampini   /* begin search for connected components */
789674ae819SStefano Zampini   cum_queue = 0;
790674ae819SStefano Zampini   ncc = 0;
791674ae819SStefano Zampini   for (n=0;n<graph->n_subsets;n++) {
792b3cdcd63SStefano Zampini     PetscInt pid = n+1;  /* partition labeled by 0 is discarded */
793b3cdcd63SStefano Zampini     PetscInt found = 0,prev = 0,first = 0,ncc_pid = 0;
794b3cdcd63SStefano Zampini     while (found != graph->subset_size[n]) {
795b3cdcd63SStefano Zampini       PetscInt added = 0;
796b3cdcd63SStefano Zampini       if (!prev) { /* search for new starting dof */
797b3cdcd63SStefano Zampini         while (PetscBTLookup(graph->touched,graph->subset_idxs[n][first])) first++;
7989566063dSJacob Faibussowitsch         PetscCall(PetscBTSet(graph->touched,graph->subset_idxs[n][first]));
799b3cdcd63SStefano Zampini         graph->queue[cum_queue] = graph->subset_idxs[n][first];
800674ae819SStefano Zampini         graph->cptr[ncc] = cum_queue;
801b3cdcd63SStefano Zampini         prev = 1;
802b3cdcd63SStefano Zampini         cum_queue++;
803b3cdcd63SStefano Zampini         found++;
804674ae819SStefano Zampini         ncc_pid++;
805b3cdcd63SStefano Zampini         ncc++;
806674ae819SStefano Zampini       }
8079566063dSJacob Faibussowitsch       PetscCall(PCBDDCGraphComputeCC_Private(graph,pid,graph->queue + cum_queue,prev,&added));
808b3cdcd63SStefano Zampini       if (!added) {
809674ae819SStefano Zampini         graph->subset_ncc[n] = ncc_pid;
810b3cdcd63SStefano Zampini         graph->cptr[ncc] = cum_queue;
811b3cdcd63SStefano Zampini       }
812b3cdcd63SStefano Zampini       prev = added;
813b3cdcd63SStefano Zampini       found += added;
814b3cdcd63SStefano Zampini       cum_queue += added;
815b3cdcd63SStefano Zampini       if (added && found == graph->subset_size[n]) {
816b3cdcd63SStefano Zampini         graph->subset_ncc[n] = ncc_pid;
817b3cdcd63SStefano Zampini         graph->cptr[ncc] = cum_queue;
818b3cdcd63SStefano Zampini       }
819b3cdcd63SStefano Zampini     }
820674ae819SStefano Zampini   }
821674ae819SStefano Zampini   graph->ncc = ncc;
822acc113dbSStefano Zampini   graph->queue_sorted = PETSC_FALSE;
823674ae819SStefano Zampini   PetscFunctionReturn(0);
824674ae819SStefano Zampini }
825674ae819SStefano Zampini 
826674ae819SStefano Zampini PetscErrorCode PCBDDCGraphSetUp(PCBDDCGraph graph, PetscInt custom_minimal_size, IS neumann_is, IS dirichlet_is, PetscInt n_ISForDofs, IS ISForDofs[], IS custom_primal_vertices)
827674ae819SStefano Zampini {
8287a0e7b2cSstefano_zampini   IS             subset,subset_n;
829674ae819SStefano Zampini   MPI_Comm       comm;
830674ae819SStefano Zampini   const PetscInt *is_indices;
831dc456d91SStefano Zampini   PetscInt       n_neigh,*neigh,*n_shared,**shared,*queue_global;
832674ae819SStefano Zampini   PetscInt       i,j,k,s,total_counts,nodes_touched,is_size;
833b3cdcd63SStefano Zampini   PetscMPIInt    commsize;
834cc0a5675Sstefano_zampini   PetscBool      same_set,mirrors_found;
835674ae819SStefano Zampini 
836674ae819SStefano Zampini   PetscFunctionBegin;
8377a0e7b2cSstefano_zampini   PetscValidLogicalCollectiveInt(graph->l2gmap,custom_minimal_size,2);
8387a0e7b2cSstefano_zampini   if (neumann_is) {
8397a0e7b2cSstefano_zampini     PetscValidHeaderSpecific(neumann_is,IS_CLASSID,3);
8407a0e7b2cSstefano_zampini     PetscCheckSameComm(graph->l2gmap,1,neumann_is,3);
8417a0e7b2cSstefano_zampini   }
842a07ea27aSStefano Zampini   graph->has_dirichlet = PETSC_FALSE;
843a07ea27aSStefano Zampini   if (dirichlet_is) {
8447a0e7b2cSstefano_zampini     PetscValidHeaderSpecific(dirichlet_is,IS_CLASSID,4);
845a07ea27aSStefano Zampini     PetscCheckSameComm(graph->l2gmap,1,dirichlet_is,4);
846a07ea27aSStefano Zampini     graph->has_dirichlet = PETSC_TRUE;
847a07ea27aSStefano Zampini   }
8487a0e7b2cSstefano_zampini   PetscValidLogicalCollectiveInt(graph->l2gmap,n_ISForDofs,5);
8497a0e7b2cSstefano_zampini   for (i=0;i<n_ISForDofs;i++) {
8507a0e7b2cSstefano_zampini     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,6);
8517a0e7b2cSstefano_zampini     PetscCheckSameComm(graph->l2gmap,1,ISForDofs[i],6);
8527a0e7b2cSstefano_zampini   }
8537a0e7b2cSstefano_zampini   if (custom_primal_vertices) {
854ab8c8b98SStefano Zampini     PetscValidHeaderSpecific(custom_primal_vertices,IS_CLASSID,7);
8557a0e7b2cSstefano_zampini     PetscCheckSameComm(graph->l2gmap,1,custom_primal_vertices,7);
8567a0e7b2cSstefano_zampini   }
8579566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)(graph->l2gmap),&comm));
8589566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm,&commsize));
859b3cdcd63SStefano Zampini 
860674ae819SStefano Zampini   /* custom_minimal_size */
86114f95afaSStefano Zampini   graph->custom_minimal_size = custom_minimal_size;
862674ae819SStefano Zampini   /* get info l2gmap and allocate work vectors  */
8639566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetInfo(graph->l2gmap,&n_neigh,&neigh,&n_shared,&shared));
8642635a1d4SStefano Zampini   /* check if we have any local periodic nodes (periodic BCs) */
8652635a1d4SStefano Zampini   mirrors_found = PETSC_FALSE;
8662635a1d4SStefano Zampini   if (graph->nvtxs && n_neigh) {
8672635a1d4SStefano Zampini     for (i=0; i<n_shared[0]; i++) graph->count[shared[0][i]] += 1;
8682635a1d4SStefano Zampini     for (i=0; i<n_shared[0]; i++) {
8692635a1d4SStefano Zampini       if (graph->count[shared[0][i]] > 1) {
8702635a1d4SStefano Zampini         mirrors_found = PETSC_TRUE;
8712635a1d4SStefano Zampini         break;
8722635a1d4SStefano Zampini       }
8732635a1d4SStefano Zampini     }
8742635a1d4SStefano Zampini   }
87551b0f6cfSStefano Zampini   /* compute local mirrors (if any) */
87651b0f6cfSStefano Zampini   if (mirrors_found) {
8777a0e7b2cSstefano_zampini     IS       to,from;
87851b0f6cfSStefano Zampini     PetscInt *local_indices,*global_indices;
8797a0e7b2cSstefano_zampini 
8809566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF,graph->nvtxs,0,1,&to));
8819566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingApplyIS(graph->l2gmap,to,&from));
88251b0f6cfSStefano Zampini     /* get arrays of local and global indices */
8839566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(graph->nvtxs,&local_indices));
8849566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(to,(const PetscInt**)&is_indices));
8859566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(local_indices,is_indices,graph->nvtxs));
8869566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(to,(const PetscInt**)&is_indices));
8879566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(graph->nvtxs,&global_indices));
8889566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(from,(const PetscInt**)&is_indices));
8899566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(global_indices,is_indices,graph->nvtxs));
8909566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(from,(const PetscInt**)&is_indices));
89151b0f6cfSStefano Zampini     /* allocate space for mirrors */
8929566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(graph->nvtxs,&graph->mirrors,graph->nvtxs,&graph->mirrors_set));
8939566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(graph->mirrors,graph->nvtxs));
8940a545947SLisandro Dalcin     graph->mirrors_set[0] = NULL;
89551b0f6cfSStefano Zampini 
89651b0f6cfSStefano Zampini     k=0;
89751b0f6cfSStefano Zampini     for (i=0;i<n_shared[0];i++) {
89851b0f6cfSStefano Zampini       j=shared[0][i];
89951b0f6cfSStefano Zampini       if (graph->count[j] > 1) {
90051b0f6cfSStefano Zampini         graph->mirrors[j]++;
90151b0f6cfSStefano Zampini         k++;
90251b0f6cfSStefano Zampini       }
90351b0f6cfSStefano Zampini     }
90451b0f6cfSStefano Zampini     /* allocate space for set of mirrors */
9059566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(k,&graph->mirrors_set[0]));
90651b0f6cfSStefano Zampini     for (i=1;i<graph->nvtxs;i++)
90751b0f6cfSStefano Zampini       graph->mirrors_set[i]=graph->mirrors_set[i-1]+graph->mirrors[i-1];
90851b0f6cfSStefano Zampini 
90951b0f6cfSStefano Zampini     /* fill arrays */
9109566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(graph->mirrors,graph->nvtxs));
91151b0f6cfSStefano Zampini     for (j=0;j<n_shared[0];j++) {
91251b0f6cfSStefano Zampini       i=shared[0][j];
91351b0f6cfSStefano Zampini       if (graph->count[i] > 1)
91451b0f6cfSStefano Zampini         graph->mirrors_set[i][graph->mirrors[i]++]=global_indices[i];
91551b0f6cfSStefano Zampini     }
9169566063dSJacob Faibussowitsch     PetscCall(PetscSortIntWithArray(graph->nvtxs,global_indices,local_indices));
91751b0f6cfSStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
91851b0f6cfSStefano Zampini       if (graph->mirrors[i] > 0) {
9199566063dSJacob Faibussowitsch         PetscCall(PetscFindInt(graph->mirrors_set[i][0],graph->nvtxs,global_indices,&k));
92051b0f6cfSStefano Zampini         j = global_indices[k];
92151b0f6cfSStefano Zampini         while (k > 0 && global_indices[k-1] == j) k--;
92251b0f6cfSStefano Zampini         for (j=0;j<graph->mirrors[i];j++) {
92351b0f6cfSStefano Zampini           graph->mirrors_set[i][j]=local_indices[k+j];
92451b0f6cfSStefano Zampini         }
9259566063dSJacob Faibussowitsch         PetscCall(PetscSortInt(graph->mirrors[i],graph->mirrors_set[i]));
92651b0f6cfSStefano Zampini       }
92751b0f6cfSStefano Zampini     }
9289566063dSJacob Faibussowitsch     PetscCall(PetscFree(local_indices));
9299566063dSJacob Faibussowitsch     PetscCall(PetscFree(global_indices));
9309566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&to));
9319566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&from));
9327a0e7b2cSstefano_zampini   }
9339566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(graph->count,graph->nvtxs));
93451b0f6cfSStefano Zampini 
935674ae819SStefano Zampini   /* Count total number of neigh per node */
936674ae819SStefano Zampini   k = 0;
937674ae819SStefano Zampini   for (i=1;i<n_neigh;i++) {
938674ae819SStefano Zampini     k += n_shared[i];
939674ae819SStefano Zampini     for (j=0;j<n_shared[i];j++) {
940674ae819SStefano Zampini       graph->count[shared[i][j]] += 1;
941674ae819SStefano Zampini     }
942674ae819SStefano Zampini   }
943674ae819SStefano Zampini   /* Allocate space for storing the set of neighbours for each node */
944674ae819SStefano Zampini   if (graph->nvtxs) {
9459566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(k,&graph->neighbours_set[0]));
946674ae819SStefano Zampini   }
947674ae819SStefano Zampini   for (i=1;i<graph->nvtxs;i++) { /* dont count myself */
948674ae819SStefano Zampini     graph->neighbours_set[i]=graph->neighbours_set[i-1]+graph->count[i-1];
949674ae819SStefano Zampini   }
950674ae819SStefano Zampini   /* Get information for sharing subdomains */
9519566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(graph->count,graph->nvtxs));
952674ae819SStefano Zampini   for (i=1;i<n_neigh;i++) { /* dont count myself */
953674ae819SStefano Zampini     s = n_shared[i];
954674ae819SStefano Zampini     for (j=0;j<s;j++) {
955674ae819SStefano Zampini       k = shared[i][j];
956674ae819SStefano Zampini       graph->neighbours_set[k][graph->count[k]] = neigh[i];
957674ae819SStefano Zampini       graph->count[k] += 1;
958674ae819SStefano Zampini     }
959674ae819SStefano Zampini   }
960674ae819SStefano Zampini   /* sort set of sharing subdomains */
961674ae819SStefano Zampini   for (i=0;i<graph->nvtxs;i++) {
9629566063dSJacob Faibussowitsch     PetscCall(PetscSortRemoveDupsInt(&graph->count[i],graph->neighbours_set[i]));
963674ae819SStefano Zampini   }
9647fb0e2dbSStefano Zampini   /* free memory allocated by ISLocalToGlobalMappingGetInfo */
9659566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingRestoreInfo(graph->l2gmap,&n_neigh,&neigh,&n_shared,&shared));
9667fb0e2dbSStefano Zampini 
96767c9da69SStefano Zampini   /*
96867c9da69SStefano Zampini      Get info for dofs splitting
9695777c63cSStefano Zampini      User can specify just a subset; an additional field is considered as a complementary field
97067c9da69SStefano Zampini   */
9710c85b387SStefano Zampini   for (i=0,k=0;i<n_ISForDofs;i++) {
9720c85b387SStefano Zampini     PetscInt bs;
9730c85b387SStefano Zampini 
9749566063dSJacob Faibussowitsch     PetscCall(ISGetBlockSize(ISForDofs[i],&bs));
9750c85b387SStefano Zampini     k   += bs;
9760c85b387SStefano Zampini   }
9770c85b387SStefano Zampini   for (i=0;i<graph->nvtxs;i++) graph->which_dof[i] = k; /* by default a dof belongs to the complement set */
9780c85b387SStefano Zampini   for (i=0,k=0;i<n_ISForDofs;i++) {
9790c85b387SStefano Zampini     PetscInt bs;
9800c85b387SStefano Zampini 
9819566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(ISForDofs[i],&is_size));
9829566063dSJacob Faibussowitsch     PetscCall(ISGetBlockSize(ISForDofs[i],&bs));
9839566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(ISForDofs[i],(const PetscInt**)&is_indices));
9840c85b387SStefano Zampini     for (j=0;j<is_size/bs;j++) {
9850c85b387SStefano Zampini       PetscInt b;
9860c85b387SStefano Zampini 
9870c85b387SStefano Zampini       for (b=0;b<bs;b++) {
9880c85b387SStefano Zampini         PetscInt jj = bs*j + b;
9890c85b387SStefano Zampini 
9900c85b387SStefano Zampini         if (is_indices[jj] > -1 && is_indices[jj] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
9910c85b387SStefano Zampini           graph->which_dof[is_indices[jj]] = k+b;
9920c85b387SStefano Zampini         }
993674ae819SStefano Zampini       }
99467c9da69SStefano Zampini     }
9959566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(ISForDofs[i],(const PetscInt**)&is_indices));
9960c85b387SStefano Zampini     k   += bs;
9975777c63cSStefano Zampini   }
9985777c63cSStefano Zampini 
999674ae819SStefano Zampini   /* Take into account Neumann nodes */
1000674ae819SStefano Zampini   if (neumann_is) {
10019566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(neumann_is,&is_size));
10029566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(neumann_is,(const PetscInt**)&is_indices));
1003674ae819SStefano Zampini     for (i=0;i<is_size;i++) {
1004d50ed60aSStefano Zampini       if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
10057a0e7b2cSstefano_zampini         graph->special_dof[is_indices[i]] = PCBDDCGRAPH_NEUMANN_MARK;
1006674ae819SStefano Zampini       }
10073c629590SStefano Zampini     }
10089566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(neumann_is,(const PetscInt**)&is_indices));
1009674ae819SStefano Zampini   }
10107a0e7b2cSstefano_zampini   /* Take into account Dirichlet nodes (they overwrite any neumann boundary mark previously set) */
1011674ae819SStefano Zampini   if (dirichlet_is) {
10129566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(dirichlet_is,&is_size));
10139566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(dirichlet_is,(const PetscInt**)&is_indices));
1014674ae819SStefano Zampini     for (i=0;i<is_size;i++) {
1015d50ed60aSStefano Zampini       if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
10167a0e7b2cSstefano_zampini         if (commsize > graph->commsizelimit) { /* dirichlet nodes treated as internal */
10179566063dSJacob Faibussowitsch           PetscCall(PetscBTSet(graph->touched,is_indices[i]));
10187a0e7b2cSstefano_zampini           graph->subset[is_indices[i]] = 0;
10197a0e7b2cSstefano_zampini         }
10207a0e7b2cSstefano_zampini         graph->special_dof[is_indices[i]] = PCBDDCGRAPH_DIRICHLET_MARK;
1021674ae819SStefano Zampini       }
10223c629590SStefano Zampini     }
10239566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(dirichlet_is,(const PetscInt**)&is_indices));
10247fb0e2dbSStefano Zampini   }
102508a5cf49SStefano Zampini   /* mark local periodic nodes (if any) and adapt CSR graph (if any) */
102651b0f6cfSStefano Zampini   if (graph->mirrors) {
102751b0f6cfSStefano Zampini     for (i=0;i<graph->nvtxs;i++)
102851b0f6cfSStefano Zampini       if (graph->mirrors[i])
10290cf82fd6SStefano Zampini         graph->special_dof[i] = PCBDDCGRAPH_LOCAL_PERIODIC_MARK;
103051b0f6cfSStefano Zampini 
10311633d1f0SStefano Zampini     if (graph->xadj) {
103208a5cf49SStefano Zampini       PetscInt *new_xadj,*new_adjncy;
103351b0f6cfSStefano Zampini       /* sort CSR graph */
103460d4fc61SSatish Balay       for (i=0;i<graph->nvtxs;i++) {
10359566063dSJacob Faibussowitsch         PetscCall(PetscSortInt(graph->xadj[i+1]-graph->xadj[i],&graph->adjncy[graph->xadj[i]]));
103660d4fc61SSatish Balay       }
103751b0f6cfSStefano Zampini       /* adapt local CSR graph in case of local periodicity */
103851b0f6cfSStefano Zampini       k = 0;
103951b0f6cfSStefano Zampini       for (i=0;i<graph->nvtxs;i++)
104051b0f6cfSStefano Zampini         for (j=graph->xadj[i];j<graph->xadj[i+1];j++)
104151b0f6cfSStefano Zampini           k += graph->mirrors[graph->adjncy[j]];
104251b0f6cfSStefano Zampini 
10439566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(graph->nvtxs+1,&new_xadj));
10449566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(k+graph->xadj[graph->nvtxs],&new_adjncy));
104551b0f6cfSStefano Zampini       new_xadj[0] = 0;
104651b0f6cfSStefano Zampini       for (i=0;i<graph->nvtxs;i++) {
104751b0f6cfSStefano Zampini         k = graph->xadj[i+1]-graph->xadj[i];
10489566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(&new_adjncy[new_xadj[i]],&graph->adjncy[graph->xadj[i]],k));
104951b0f6cfSStefano Zampini         new_xadj[i+1] = new_xadj[i]+k;
105051b0f6cfSStefano Zampini         for (j=graph->xadj[i];j<graph->xadj[i+1];j++) {
105151b0f6cfSStefano Zampini           k = graph->mirrors[graph->adjncy[j]];
10529566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(&new_adjncy[new_xadj[i+1]],graph->mirrors_set[graph->adjncy[j]],k));
105351b0f6cfSStefano Zampini           new_xadj[i+1] += k;
105451b0f6cfSStefano Zampini         }
105551b0f6cfSStefano Zampini         k = new_xadj[i+1]-new_xadj[i];
10569566063dSJacob Faibussowitsch         PetscCall(PetscSortRemoveDupsInt(&k,&new_adjncy[new_xadj[i]]));
105751b0f6cfSStefano Zampini         new_xadj[i+1] = new_xadj[i]+k;
105851b0f6cfSStefano Zampini       }
105951b0f6cfSStefano Zampini       /* set new CSR into graph */
10609566063dSJacob Faibussowitsch       PetscCall(PetscFree(graph->xadj));
10619566063dSJacob Faibussowitsch       PetscCall(PetscFree(graph->adjncy));
106251b0f6cfSStefano Zampini       graph->xadj = new_xadj;
106351b0f6cfSStefano Zampini       graph->adjncy = new_adjncy;
106451b0f6cfSStefano Zampini     }
106508a5cf49SStefano Zampini   }
106651b0f6cfSStefano Zampini 
106717eb1463SStefano Zampini   /* mark special nodes (if any) -> each will become a single node equivalence class */
1068674ae819SStefano Zampini   if (custom_primal_vertices) {
10699566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(custom_primal_vertices,&is_size));
10709566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(custom_primal_vertices,(const PetscInt**)&is_indices));
10717a0e7b2cSstefano_zampini     for (i=0,j=0;i<is_size;i++) {
10727a0e7b2cSstefano_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 */
10737a0e7b2cSstefano_zampini         graph->special_dof[is_indices[i]] = PCBDDCGRAPH_SPECIAL_MARK-j;
10749b70a373SStefano Zampini         j++;
10759b70a373SStefano Zampini       }
10769b70a373SStefano Zampini     }
10779566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(custom_primal_vertices,(const PetscInt**)&is_indices));
107817eb1463SStefano Zampini   }
10799b70a373SStefano Zampini 
10804b2aedd3SStefano Zampini   /* mark interior nodes (if commsize > graph->commsizelimit) as touched and belonging to partition number 0 */
10814b2aedd3SStefano Zampini   if (commsize > graph->commsizelimit) {
1082674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
1083674ae819SStefano Zampini       if (!graph->count[i]) {
10849566063dSJacob Faibussowitsch         PetscCall(PetscBTSet(graph->touched,i));
1085674ae819SStefano Zampini         graph->subset[i] = 0;
1086674ae819SStefano Zampini       }
1087674ae819SStefano Zampini     }
1088b3cdcd63SStefano Zampini   }
1089b3cdcd63SStefano Zampini 
1090674ae819SStefano Zampini   /* init graph structure and compute default subsets */
1091674ae819SStefano Zampini   nodes_touched = 0;
1092674ae819SStefano Zampini   for (i=0;i<graph->nvtxs;i++) {
1093df48933bSStefano Zampini     if (PetscBTLookup(graph->touched,i)) {
1094674ae819SStefano Zampini       nodes_touched++;
1095674ae819SStefano Zampini     }
1096674ae819SStefano Zampini   }
1097674ae819SStefano Zampini   i = 0;
1098674ae819SStefano Zampini   graph->ncc = 0;
1099674ae819SStefano Zampini   total_counts = 0;
11007babac9bSStefano Zampini 
11017babac9bSStefano Zampini   /* allocated space for queues */
11024b2aedd3SStefano Zampini   if (commsize == graph->commsizelimit) {
11039566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(graph->nvtxs+1,&graph->cptr,graph->nvtxs,&graph->queue));
11047babac9bSStefano Zampini   } else {
11057babac9bSStefano Zampini     PetscInt nused = graph->nvtxs - nodes_touched;
11069566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nused+1,&graph->cptr,nused,&graph->queue));
11077babac9bSStefano Zampini   }
11087babac9bSStefano Zampini 
1109674ae819SStefano Zampini   while (nodes_touched<graph->nvtxs) {
1110674ae819SStefano Zampini     /*  find first untouched node in local ordering */
1111b3cdcd63SStefano Zampini     while (PetscBTLookup(graph->touched,i)) i++;
11129566063dSJacob Faibussowitsch     PetscCall(PetscBTSet(graph->touched,i));
1113674ae819SStefano Zampini     graph->subset[i] = graph->ncc+1;
1114674ae819SStefano Zampini     graph->cptr[graph->ncc] = total_counts;
1115674ae819SStefano Zampini     graph->queue[total_counts] = i;
1116674ae819SStefano Zampini     total_counts++;
1117674ae819SStefano Zampini     nodes_touched++;
1118674ae819SStefano Zampini     /* now find all other nodes having the same set of sharing subdomains */
1119674ae819SStefano Zampini     for (j=i+1;j<graph->nvtxs;j++) {
112074e413f5SStefano Zampini       /* check for same number of sharing subdomains, dof number and same special mark */
1121df48933bSStefano 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]) {
1122674ae819SStefano Zampini         /* check for same set of sharing subdomains */
1123674ae819SStefano Zampini         same_set = PETSC_TRUE;
1124674ae819SStefano Zampini         for (k=0;k<graph->count[j];k++) {
1125674ae819SStefano Zampini           if (graph->neighbours_set[i][k] != graph->neighbours_set[j][k]) {
1126674ae819SStefano Zampini             same_set = PETSC_FALSE;
1127674ae819SStefano Zampini           }
1128674ae819SStefano Zampini         }
1129ab8c8b98SStefano Zampini         /* I have found a friend of mine */
1130674ae819SStefano Zampini         if (same_set) {
11319566063dSJacob Faibussowitsch           PetscCall(PetscBTSet(graph->touched,j));
1132b3cdcd63SStefano Zampini           graph->subset[j] = graph->ncc+1;
1133674ae819SStefano Zampini           nodes_touched++;
1134674ae819SStefano Zampini           graph->queue[total_counts] = j;
1135674ae819SStefano Zampini           total_counts++;
1136674ae819SStefano Zampini         }
1137674ae819SStefano Zampini       }
1138674ae819SStefano Zampini     }
1139674ae819SStefano Zampini     graph->ncc++;
1140674ae819SStefano Zampini   }
1141b3cdcd63SStefano 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 */
1142674ae819SStefano Zampini   graph->n_subsets = graph->ncc;
11439566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(graph->n_subsets,&graph->subset_ncc));
1144674ae819SStefano Zampini   for (i=0;i<graph->n_subsets;i++) {
1145674ae819SStefano Zampini     graph->subset_ncc[i] = 1;
1146674ae819SStefano Zampini   }
1147674ae819SStefano Zampini   /* final pointer */
1148674ae819SStefano Zampini   graph->cptr[graph->ncc] = total_counts;
1149ec1c413dSStefano Zampini 
1150b3cdcd63SStefano Zampini   /* For consistency reasons (among neighbours), I need to sort (by global ordering) each connected component */
1151ec1c413dSStefano Zampini   /* Get a reference node (min index in global ordering) for each subset for tagging messages */
11529566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(graph->ncc,&graph->subset_ref_node));
11539566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(graph->cptr[graph->ncc],&queue_global));
11549566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_global));
1155ec1c413dSStefano Zampini   for (j=0;j<graph->ncc;j++) {
11569566063dSJacob Faibussowitsch     PetscCall(PetscSortIntWithArray(graph->cptr[j+1]-graph->cptr[j],&queue_global[graph->cptr[j]],&graph->queue[graph->cptr[j]]));
1157dc456d91SStefano Zampini     graph->subset_ref_node[j] = graph->queue[graph->cptr[j]];
1158f0321c90SStefano Zampini   }
11599566063dSJacob Faibussowitsch   PetscCall(PetscFree(queue_global));
1160acc113dbSStefano Zampini   graph->queue_sorted = PETSC_TRUE;
116148bebe3cSStefano Zampini 
1162b3cdcd63SStefano Zampini   /* save information on subsets (needed when analyzing the connected components) */
11632f566a09SStefano Zampini   if (graph->ncc) {
11649566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(graph->ncc,&graph->subset_size,graph->ncc,&graph->subset_idxs));
11659566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(graph->cptr[graph->ncc],&graph->subset_idxs[0]));
11669566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(graph->subset_idxs[0],graph->cptr[graph->ncc]));
1167ec1c413dSStefano Zampini     for (j=1;j<graph->ncc;j++) {
1168b3cdcd63SStefano Zampini       graph->subset_size[j-1] = graph->cptr[j] - graph->cptr[j-1];
1169b3cdcd63SStefano Zampini       graph->subset_idxs[j] = graph->subset_idxs[j-1] + graph->subset_size[j-1];
1170ec1c413dSStefano Zampini     }
1171b3cdcd63SStefano Zampini     graph->subset_size[graph->ncc-1] = graph->cptr[graph->ncc] - graph->cptr[graph->ncc-1];
11729566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(graph->subset_idxs[0],graph->queue,graph->cptr[graph->ncc]));
1173ec1c413dSStefano Zampini   }
1174b3cdcd63SStefano Zampini 
1175f0321c90SStefano Zampini   /* renumber reference nodes */
11769566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)(graph->l2gmap)),graph->ncc,graph->subset_ref_node,PETSC_COPY_VALUES,&subset_n));
11779566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApplyIS(graph->l2gmap,subset_n,&subset));
11789566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&subset_n));
11799566063dSJacob Faibussowitsch   PetscCall(ISRenumber(subset,NULL,NULL,&subset_n));
11809566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&subset));
11819566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(subset_n,&k));
118263a3b9bcSJacob Faibussowitsch   PetscCheck(k == graph->ncc,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid size of new subset! %" PetscInt_FMT " != %" PetscInt_FMT,k,graph->ncc);
11839566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(subset_n,&is_indices));
11849566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(graph->subset_ref_node,is_indices,graph->ncc));
11859566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(subset_n,&is_indices));
11869566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&subset_n));
1187ec1c413dSStefano Zampini 
1188ec1c413dSStefano Zampini   /* free workspace */
1189b3cdcd63SStefano Zampini   graph->setupcalled = PETSC_TRUE;
1190674ae819SStefano Zampini   PetscFunctionReturn(0);
1191674ae819SStefano Zampini }
1192674ae819SStefano Zampini 
1193ab8c8b98SStefano Zampini PetscErrorCode PCBDDCGraphResetCoords(PCBDDCGraph graph)
1194ab8c8b98SStefano Zampini {
1195ab8c8b98SStefano Zampini   PetscFunctionBegin;
1196ab8c8b98SStefano Zampini   if (!graph) PetscFunctionReturn(0);
11979566063dSJacob Faibussowitsch   PetscCall(PetscFree(graph->coords));
1198ab8c8b98SStefano Zampini   graph->cdim  = 0;
1199ab8c8b98SStefano Zampini   graph->cnloc = 0;
1200ab8c8b98SStefano Zampini   graph->cloc  = PETSC_FALSE;
1201ab8c8b98SStefano Zampini   PetscFunctionReturn(0);
1202ab8c8b98SStefano Zampini }
1203ab8c8b98SStefano Zampini 
1204674ae819SStefano Zampini PetscErrorCode PCBDDCGraphResetCSR(PCBDDCGraph graph)
1205674ae819SStefano Zampini {
1206674ae819SStefano Zampini   PetscFunctionBegin;
120707d44595Sstefano_zampini   if (!graph) PetscFunctionReturn(0);
1208a1dbd327SStefano Zampini   if (graph->freecsr) {
12099566063dSJacob Faibussowitsch     PetscCall(PetscFree(graph->xadj));
12109566063dSJacob Faibussowitsch     PetscCall(PetscFree(graph->adjncy));
1211a1dbd327SStefano Zampini   } else {
1212a1dbd327SStefano Zampini     graph->xadj = NULL;
1213a1dbd327SStefano Zampini     graph->adjncy = NULL;
1214a1dbd327SStefano Zampini   }
1215c8272957SStefano Zampini   graph->freecsr = PETSC_FALSE;
1216575ad6abSStefano Zampini   graph->nvtxs_csr = 0;
1217674ae819SStefano Zampini   PetscFunctionReturn(0);
1218674ae819SStefano Zampini }
1219674ae819SStefano Zampini 
1220674ae819SStefano Zampini PetscErrorCode PCBDDCGraphReset(PCBDDCGraph graph)
1221674ae819SStefano Zampini {
1222674ae819SStefano Zampini   PetscFunctionBegin;
122307d44595Sstefano_zampini   if (!graph) PetscFunctionReturn(0);
12249566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingDestroy(&graph->l2gmap));
12259566063dSJacob Faibussowitsch   PetscCall(PetscFree(graph->subset_ncc));
12269566063dSJacob Faibussowitsch   PetscCall(PetscFree(graph->subset_ref_node));
1227674ae819SStefano Zampini   if (graph->nvtxs) {
12289566063dSJacob Faibussowitsch     PetscCall(PetscFree(graph->neighbours_set[0]));
1229674ae819SStefano Zampini   }
12309566063dSJacob Faibussowitsch   PetscCall(PetscBTDestroy(&graph->touched));
1231d0609cedSBarry Smith   PetscCall(PetscFree5(graph->count,graph->neighbours_set,graph->subset,graph->which_dof,graph->special_dof));
12329566063dSJacob Faibussowitsch   PetscCall(PetscFree2(graph->cptr,graph->queue));
123351b0f6cfSStefano Zampini   if (graph->mirrors) {
12349566063dSJacob Faibussowitsch     PetscCall(PetscFree(graph->mirrors_set[0]));
123551b0f6cfSStefano Zampini   }
12369566063dSJacob Faibussowitsch   PetscCall(PetscFree2(graph->mirrors,graph->mirrors_set));
1237b3cdcd63SStefano Zampini   if (graph->subset_idxs) {
12389566063dSJacob Faibussowitsch     PetscCall(PetscFree(graph->subset_idxs[0]));
1239ec1c413dSStefano Zampini   }
12409566063dSJacob Faibussowitsch   PetscCall(PetscFree2(graph->subset_size,graph->subset_idxs));
12419566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&graph->dirdofs));
12429566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&graph->dirdofsB));
12434f1b2e48SStefano Zampini   if (graph->n_local_subs) {
12449566063dSJacob Faibussowitsch     PetscCall(PetscFree(graph->local_subs));
12454f1b2e48SStefano Zampini   }
1246a07ea27aSStefano Zampini   graph->has_dirichlet       = PETSC_FALSE;
1247cc0a5675Sstefano_zampini   graph->twodimset           = PETSC_FALSE;
1248cc0a5675Sstefano_zampini   graph->twodim              = PETSC_FALSE;
1249674ae819SStefano Zampini   graph->nvtxs               = 0;
12507babac9bSStefano Zampini   graph->nvtxs_global        = 0;
1251674ae819SStefano Zampini   graph->n_subsets           = 0;
1252674ae819SStefano Zampini   graph->custom_minimal_size = 1;
12534f1b2e48SStefano Zampini   graph->n_local_subs        = 0;
1254be12c134Sstefano_zampini   graph->maxcount            = PETSC_MAX_INT;
1255fb597685SStefano Zampini   graph->setupcalled         = PETSC_FALSE;
1256674ae819SStefano Zampini   PetscFunctionReturn(0);
1257674ae819SStefano Zampini }
1258674ae819SStefano Zampini 
1259be12c134Sstefano_zampini PetscErrorCode PCBDDCGraphInit(PCBDDCGraph graph, ISLocalToGlobalMapping l2gmap, PetscInt N, PetscInt maxcount)
1260674ae819SStefano Zampini {
1261674ae819SStefano Zampini   PetscInt       n;
1262674ae819SStefano Zampini 
1263674ae819SStefano Zampini   PetscFunctionBegin;
1264674ae819SStefano Zampini   PetscValidPointer(graph,1);
1265674ae819SStefano Zampini   PetscValidHeaderSpecific(l2gmap,IS_LTOGM_CLASSID,2);
12667babac9bSStefano Zampini   PetscValidLogicalCollectiveInt(l2gmap,N,3);
1267be12c134Sstefano_zampini   PetscValidLogicalCollectiveInt(l2gmap,maxcount,4);
12688e61c736SStefano Zampini   /* raise an error if already allocated */
126928b400f6SJacob Faibussowitsch   PetscCheck(!graph->nvtxs_global,PetscObjectComm((PetscObject)l2gmap),PETSC_ERR_PLIB,"BDDCGraph already initialized");
1270674ae819SStefano Zampini   /* set number of vertices */
12719566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)l2gmap));
1272674ae819SStefano Zampini   graph->l2gmap = l2gmap;
12739566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetSize(l2gmap,&n));
1274674ae819SStefano Zampini   graph->nvtxs = n;
12757fb0e2dbSStefano Zampini   graph->nvtxs_global = N;
1276674ae819SStefano Zampini   /* allocate used space */
12779566063dSJacob Faibussowitsch   PetscCall(PetscBTCreate(graph->nvtxs,&graph->touched));
12789566063dSJacob 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));
1279674ae819SStefano Zampini   /* zeroes memory */
12809566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(graph->count,graph->nvtxs));
12819566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(graph->subset,graph->nvtxs));
128263602bcaSStefano Zampini   /* use -1 as a default value for which_dof array */
128363602bcaSStefano Zampini   for (n=0;n<graph->nvtxs;n++) graph->which_dof[n] = -1;
12849566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(graph->special_dof,graph->nvtxs));
1285674ae819SStefano Zampini   /* zeroes first pointer to neighbour set */
1286674ae819SStefano Zampini   if (graph->nvtxs) {
12870a545947SLisandro Dalcin     graph->neighbours_set[0] = NULL;
1288674ae819SStefano Zampini   }
1289674ae819SStefano Zampini   /* zeroes workspace for values of ncc */
12900a545947SLisandro Dalcin   graph->subset_ncc = NULL;
12910a545947SLisandro Dalcin   graph->subset_ref_node = NULL;
1292be12c134Sstefano_zampini   /* maxcount for cc */
1293be12c134Sstefano_zampini   graph->maxcount = maxcount;
1294674ae819SStefano Zampini   PetscFunctionReturn(0);
1295674ae819SStefano Zampini }
1296674ae819SStefano Zampini 
1297674ae819SStefano Zampini PetscErrorCode PCBDDCGraphDestroy(PCBDDCGraph* graph)
1298674ae819SStefano Zampini {
1299674ae819SStefano Zampini   PetscFunctionBegin;
13009566063dSJacob Faibussowitsch   PetscCall(PCBDDCGraphResetCSR(*graph));
13019566063dSJacob Faibussowitsch   PetscCall(PCBDDCGraphResetCoords(*graph));
13029566063dSJacob Faibussowitsch   PetscCall(PCBDDCGraphReset(*graph));
13039566063dSJacob Faibussowitsch   PetscCall(PetscFree(*graph));
1304674ae819SStefano Zampini   PetscFunctionReturn(0);
1305674ae819SStefano Zampini }
1306674ae819SStefano Zampini 
1307674ae819SStefano Zampini PetscErrorCode PCBDDCGraphCreate(PCBDDCGraph *graph)
1308674ae819SStefano Zampini {
1309674ae819SStefano Zampini   PCBDDCGraph    new_graph;
1310674ae819SStefano Zampini 
1311674ae819SStefano Zampini   PetscFunctionBegin;
13129566063dSJacob Faibussowitsch   PetscCall(PetscNew(&new_graph));
1313674ae819SStefano Zampini   new_graph->custom_minimal_size = 1;
13144b2aedd3SStefano Zampini   new_graph->commsizelimit = 1;
1315674ae819SStefano Zampini   *graph = new_graph;
1316674ae819SStefano Zampini   PetscFunctionReturn(0);
1317674ae819SStefano Zampini }
1318