xref: /petsc/src/ksp/pc/impls/bddc/bddcgraph.c (revision 7a0e7b2c98bf737ab9b4b75e6b535ea57e2ca273)
1af0996ceSBarry Smith #include <petsc/private/petscimpl.h>
2ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
3ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddcstructs.h>
4674ae819SStefano Zampini 
5674ae819SStefano Zampini #undef __FUNCT__
6d62866d3SStefano Zampini #define __FUNCT__ "PCBDDCGraphGetDirichletDofsB"
7d62866d3SStefano Zampini PetscErrorCode PCBDDCGraphGetDirichletDofsB(PCBDDCGraph graph, IS* dirdofs)
8d62866d3SStefano Zampini {
9d62866d3SStefano Zampini   PetscErrorCode ierr;
10d62866d3SStefano Zampini 
11d62866d3SStefano Zampini   PetscFunctionBegin;
12d62866d3SStefano Zampini   if (graph->dirdofsB) {
13d62866d3SStefano Zampini     ierr = PetscObjectReference((PetscObject)graph->dirdofsB);CHKERRQ(ierr);
14d62866d3SStefano Zampini   } else if (graph->has_dirichlet) {
15d62866d3SStefano Zampini     PetscInt i,size;
16d62866d3SStefano Zampini     PetscInt *dirdofs_idxs;
17d62866d3SStefano Zampini 
18d62866d3SStefano Zampini     size = 0;
19d62866d3SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
20d62866d3SStefano Zampini       if (graph->count[i] && graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) size++;
21d62866d3SStefano Zampini     }
22d62866d3SStefano Zampini 
23d62866d3SStefano Zampini     ierr = PetscMalloc1(size,&dirdofs_idxs);CHKERRQ(ierr);
24d62866d3SStefano Zampini     size = 0;
25d62866d3SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
26d62866d3SStefano Zampini       if (graph->count[i] && graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) dirdofs_idxs[size++] = i;
27d62866d3SStefano Zampini     }
28d62866d3SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,size,dirdofs_idxs,PETSC_OWN_POINTER,&graph->dirdofsB);CHKERRQ(ierr);
29d62866d3SStefano Zampini     ierr = PetscObjectReference((PetscObject)graph->dirdofsB);CHKERRQ(ierr);
30d62866d3SStefano Zampini   }
31d62866d3SStefano Zampini   *dirdofs = graph->dirdofsB;
32d62866d3SStefano Zampini   PetscFunctionReturn(0);
33d62866d3SStefano Zampini }
34d62866d3SStefano Zampini 
35d62866d3SStefano Zampini #undef __FUNCT__
361b968477SStefano Zampini #define __FUNCT__ "PCBDDCGraphGetDirichletDofs"
371b968477SStefano Zampini PetscErrorCode PCBDDCGraphGetDirichletDofs(PCBDDCGraph graph, IS* dirdofs)
381b968477SStefano Zampini {
391b968477SStefano Zampini   PetscErrorCode ierr;
401b968477SStefano Zampini 
411b968477SStefano Zampini   PetscFunctionBegin;
42a07ea27aSStefano Zampini   if (graph->dirdofs) {
43a07ea27aSStefano Zampini     ierr = PetscObjectReference((PetscObject)graph->dirdofs);CHKERRQ(ierr);
44a07ea27aSStefano Zampini   } else if (graph->has_dirichlet) {
45a07ea27aSStefano Zampini     PetscInt i,size;
46a07ea27aSStefano Zampini     PetscInt *dirdofs_idxs;
47a07ea27aSStefano Zampini 
481b968477SStefano Zampini     size = 0;
491b968477SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
501b968477SStefano Zampini       if (graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) size++;
511b968477SStefano Zampini     }
521b968477SStefano Zampini 
531b968477SStefano Zampini     ierr = PetscMalloc1(size,&dirdofs_idxs);CHKERRQ(ierr);
541b968477SStefano Zampini     size = 0;
551b968477SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
561b968477SStefano Zampini       if (graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) dirdofs_idxs[size++] = i;
571b968477SStefano Zampini     }
58a07ea27aSStefano Zampini     ierr = ISCreateGeneral(PetscObjectComm((PetscObject)graph->l2gmap),size,dirdofs_idxs,PETSC_OWN_POINTER,&graph->dirdofs);CHKERRQ(ierr);
59a07ea27aSStefano Zampini     ierr = PetscObjectReference((PetscObject)graph->dirdofs);CHKERRQ(ierr);
601b968477SStefano Zampini   }
61a07ea27aSStefano Zampini   *dirdofs = graph->dirdofs;
621b968477SStefano Zampini   PetscFunctionReturn(0);
631b968477SStefano Zampini }
641b968477SStefano Zampini 
651b968477SStefano Zampini #undef __FUNCT__
66674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphASCIIView"
67e49050b4SStefano Zampini PetscErrorCode PCBDDCGraphASCIIView(PCBDDCGraph graph, PetscInt verbosity_level, PetscViewer viewer)
68674ae819SStefano Zampini {
692b510759SStefano Zampini   PetscInt       i,j,tabs;
7093fb5fd3SStefano Zampini   PetscInt*      queue_in_global_numbering;
71674ae819SStefano Zampini   PetscErrorCode ierr;
72674ae819SStefano Zampini 
73674ae819SStefano Zampini   PetscFunctionBegin;
741575c14dSBarry Smith   ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
752b510759SStefano Zampini   ierr = PetscViewerASCIIGetTab(viewer,&tabs);CHKERRQ(ierr);
76674ae819SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"--------------------------------------------------\n");CHKERRQ(ierr);
77674ae819SStefano Zampini   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
78674ae819SStefano Zampini   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Local BDDC graph for subdomain %04d\n",PetscGlobalRank);CHKERRQ(ierr);
79674ae819SStefano Zampini   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Number of vertices %d\n",graph->nvtxs);CHKERRQ(ierr);
8014f95afaSStefano Zampini   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Custom minimal size %d\n",graph->custom_minimal_size);CHKERRQ(ierr);
81be12c134Sstefano_zampini   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Max count %d\n",graph->maxcount);CHKERRQ(ierr);
8248bebe3cSStefano Zampini   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Topological two dim? %d\n",graph->twodim);CHKERRQ(ierr);
83d4890cceSStefano Zampini   if (verbosity_level > 2) {
84674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
85674ae819SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"%d:\n",i);CHKERRQ(ierr);
86674ae819SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   which_dof: %d\n",graph->which_dof[i]);CHKERRQ(ierr);
8774e413f5SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   special_dof: %d\n",graph->special_dof[i]);CHKERRQ(ierr);
88674ae819SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   neighbours: %d\n",graph->count[i]);CHKERRQ(ierr);
892b510759SStefano Zampini       ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
90674ae819SStefano Zampini       if (graph->count[i]) {
91674ae819SStefano Zampini         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"     set of neighbours:");CHKERRQ(ierr);
92674ae819SStefano Zampini         for (j=0;j<graph->count[i];j++) {
93674ae819SStefano Zampini           ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d",graph->neighbours_set[i][j]);CHKERRQ(ierr);
94674ae819SStefano Zampini         }
95674ae819SStefano Zampini         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
96674ae819SStefano Zampini       }
972b510759SStefano Zampini       ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr);
982b510759SStefano Zampini       ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
9951b0f6cfSStefano Zampini       if (graph->mirrors) {
10051b0f6cfSStefano Zampini         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   mirrors: %d\n",graph->mirrors[i]);CHKERRQ(ierr);
10151b0f6cfSStefano Zampini         if (graph->mirrors[i]) {
1022b510759SStefano Zampini           ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
10351b0f6cfSStefano Zampini           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"     set of mirrors:");CHKERRQ(ierr);
10451b0f6cfSStefano Zampini           for (j=0;j<graph->mirrors[i];j++) {
10551b0f6cfSStefano Zampini             ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d",graph->mirrors_set[i][j]);CHKERRQ(ierr);
10651b0f6cfSStefano Zampini           }
10751b0f6cfSStefano Zampini           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
1082b510759SStefano Zampini           ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr);
1092b510759SStefano Zampini           ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
11051b0f6cfSStefano Zampini         }
11151b0f6cfSStefano Zampini       }
112d4890cceSStefano Zampini       if (verbosity_level > 3) {
1131633d1f0SStefano Zampini         if (graph->xadj) {
114674ae819SStefano Zampini           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   local adj list:");CHKERRQ(ierr);
1152b510759SStefano Zampini           ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
116674ae819SStefano Zampini           for (j=graph->xadj[i];j<graph->xadj[i+1];j++) {
117674ae819SStefano Zampini             ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d",graph->adjncy[j]);CHKERRQ(ierr);
118674ae819SStefano Zampini           }
119674ae819SStefano Zampini           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
1202b510759SStefano Zampini           ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr);
1212b510759SStefano Zampini           ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
122ec1c413dSStefano Zampini         } else {
123ec1c413dSStefano Zampini           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   no adj info\n");CHKERRQ(ierr);
124674ae819SStefano Zampini         }
125e49050b4SStefano Zampini       }
126531609faSStefano Zampini       if (graph->n_local_subs) {
127531609faSStefano Zampini         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   local sub id: %d\n",graph->local_subs[i]);CHKERRQ(ierr);
128531609faSStefano Zampini       }
129674ae819SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   interface subset id: %d\n",graph->subset[i]);CHKERRQ(ierr);
130674ae819SStefano Zampini       if (graph->subset[i] && graph->subset_ncc) {
131674ae819SStefano Zampini         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   ncc for subset: %d\n",graph->subset_ncc[graph->subset[i]-1]);CHKERRQ(ierr);
132674ae819SStefano Zampini       }
133674ae819SStefano Zampini     }
134e49050b4SStefano Zampini   }
135674ae819SStefano Zampini   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Total number of connected components %d\n",graph->ncc);CHKERRQ(ierr);
136785e854fSJed Brown   ierr = PetscMalloc1(graph->cptr[graph->ncc],&queue_in_global_numbering);CHKERRQ(ierr);
13793fb5fd3SStefano Zampini   ierr = ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_in_global_numbering);CHKERRQ(ierr);
138674ae819SStefano Zampini   for (i=0;i<graph->ncc;i++) {
1391ce3d396SStefano Zampini     PetscInt node_num=graph->queue[graph->cptr[i]];
1401ce3d396SStefano Zampini     PetscBool printcc = PETSC_FALSE;
14108b6f892Sstefano_zampini     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  cc %d (size %d, fid %d, neighs:",i,graph->cptr[i+1]-graph->cptr[i],graph->which_dof[node_num]);CHKERRQ(ierr);
1422b510759SStefano Zampini     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
143674ae819SStefano Zampini     for (j=0;j<graph->count[node_num];j++) {
144674ae819SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d",graph->neighbours_set[node_num][j]);CHKERRQ(ierr);
145674ae819SStefano Zampini     }
146d4890cceSStefano Zampini     if (verbosity_level > 1) {
147674ae819SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"):");CHKERRQ(ierr);
148d4890cceSStefano Zampini       if (graph->twodim || graph->count[node_num] > 1 || (graph->count[node_num] == 1 && graph->special_dof[node_num] == PCBDDCGRAPH_NEUMANN_MARK)) {
1491ce3d396SStefano Zampini         printcc = PETSC_TRUE;
1501ce3d396SStefano Zampini       }
1511ce3d396SStefano Zampini       if (printcc) {
152674ae819SStefano Zampini         for (j=graph->cptr[i];j<graph->cptr[i+1];j++) {
15393fb5fd3SStefano Zampini           ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d (%d)",graph->queue[j],queue_in_global_numbering[j]);CHKERRQ(ierr);
154674ae819SStefano Zampini         }
1551ce3d396SStefano Zampini       }
156d4890cceSStefano Zampini     } else {
157d4890cceSStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,")");CHKERRQ(ierr);
158d4890cceSStefano Zampini     }
159674ae819SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
1602b510759SStefano Zampini     ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr);
1612b510759SStefano Zampini     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
162674ae819SStefano Zampini   }
16393fb5fd3SStefano Zampini   ierr = PetscFree(queue_in_global_numbering);CHKERRQ(ierr);
164674ae819SStefano Zampini   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
165674ae819SStefano Zampini   PetscFunctionReturn(0);
166674ae819SStefano Zampini }
167674ae819SStefano Zampini 
168674ae819SStefano Zampini #undef __FUNCT__
169c8272957SStefano Zampini #define __FUNCT__ "PCBDDCGraphRestoreCandidatesIS"
170c8272957SStefano Zampini PetscErrorCode PCBDDCGraphRestoreCandidatesIS(PCBDDCGraph graph, PetscInt *n_faces, IS *FacesIS[], PetscInt *n_edges, IS *EdgesIS[], IS *VerticesIS)
171c8272957SStefano Zampini {
172c8272957SStefano Zampini   PetscInt       i;
173c8272957SStefano Zampini   PetscErrorCode ierr;
174c8272957SStefano Zampini 
175c8272957SStefano Zampini   PetscFunctionBegin;
176c8272957SStefano Zampini   if (n_faces) {
177c8272957SStefano Zampini     if (FacesIS) {
178c8272957SStefano Zampini       for (i=0;i<*n_faces;i++) {
179c8272957SStefano Zampini         ierr = ISDestroy(&((*FacesIS)[i]));CHKERRQ(ierr);
180c8272957SStefano Zampini       }
181c8272957SStefano Zampini       ierr = PetscFree(*FacesIS);CHKERRQ(ierr);
182c8272957SStefano Zampini     }
183c8272957SStefano Zampini     *n_faces = 0;
184c8272957SStefano Zampini   }
185c8272957SStefano Zampini   if (n_edges) {
186c8272957SStefano Zampini     if (EdgesIS) {
187c8272957SStefano Zampini       for (i=0;i<*n_edges;i++) {
188c8272957SStefano Zampini         ierr = ISDestroy(&((*EdgesIS)[i]));CHKERRQ(ierr);
189c8272957SStefano Zampini       }
190c8272957SStefano Zampini       ierr = PetscFree(*EdgesIS);CHKERRQ(ierr);
191c8272957SStefano Zampini     }
192c8272957SStefano Zampini     *n_edges = 0;
193c8272957SStefano Zampini   }
194c8272957SStefano Zampini   if (VerticesIS) {
195c8272957SStefano Zampini     ierr = ISDestroy(VerticesIS);CHKERRQ(ierr);
196c8272957SStefano Zampini   }
197c8272957SStefano Zampini   PetscFunctionReturn(0);
198c8272957SStefano Zampini }
199c8272957SStefano Zampini 
200c8272957SStefano Zampini #undef __FUNCT__
201674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphGetCandidatesIS"
202a873d5d0SStefano Zampini PetscErrorCode PCBDDCGraphGetCandidatesIS(PCBDDCGraph graph, PetscInt *n_faces, IS *FacesIS[], PetscInt *n_edges, IS *EdgesIS[], IS *VerticesIS)
203674ae819SStefano Zampini {
204674ae819SStefano Zampini   IS             *ISForFaces,*ISForEdges,ISForVertices;
205be12c134Sstefano_zampini   PetscInt       i,nfc,nec,nvc,*idx,*mark;
206674ae819SStefano Zampini   PetscErrorCode ierr;
207674ae819SStefano Zampini 
208674ae819SStefano Zampini   PetscFunctionBegin;
209be12c134Sstefano_zampini   ierr = PetscCalloc1(graph->ncc,&mark);CHKERRQ(ierr);
210674ae819SStefano Zampini   /* loop on ccs to evalute number of faces, edges and vertices */
211674ae819SStefano Zampini   nfc = 0;
212674ae819SStefano Zampini   nec = 0;
213674ae819SStefano Zampini   nvc = 0;
214674ae819SStefano Zampini   for (i=0;i<graph->ncc;i++) {
2151e1f2d93SStefano Zampini     PetscInt repdof = graph->queue[graph->cptr[i]];
216be12c134Sstefano_zampini     if (graph->cptr[i+1]-graph->cptr[i] > graph->custom_minimal_size && graph->count[repdof] < graph->maxcount) {
217be12c134Sstefano_zampini       if (!graph->twodim && graph->count[repdof] == 1 && graph->special_dof[repdof] != PCBDDCGRAPH_NEUMANN_MARK) {
218674ae819SStefano Zampini         nfc++;
219be12c134Sstefano_zampini         mark[i] = 2;
220be12c134Sstefano_zampini       } else {
221674ae819SStefano Zampini         nec++;
222be12c134Sstefano_zampini         mark[i] = 1;
223674ae819SStefano Zampini       }
224674ae819SStefano Zampini     } else {
225674ae819SStefano Zampini       nvc += graph->cptr[i+1]-graph->cptr[i];
226674ae819SStefano Zampini     }
227674ae819SStefano Zampini   }
228adfc4fb2SStefano Zampini 
229674ae819SStefano Zampini   /* allocate IS arrays for faces, edges. Vertices need a single index set. */
230cf5a6209SStefano Zampini   if (FacesIS) {
231785e854fSJed Brown     ierr = PetscMalloc1(nfc,&ISForFaces);CHKERRQ(ierr);
232cf5a6209SStefano Zampini   }
233cf5a6209SStefano Zampini   if (EdgesIS) {
234785e854fSJed Brown     ierr = PetscMalloc1(nec,&ISForEdges);CHKERRQ(ierr);
235cf5a6209SStefano Zampini   }
236cf5a6209SStefano Zampini   if (VerticesIS) {
237785e854fSJed Brown     ierr = PetscMalloc1(nvc,&idx);CHKERRQ(ierr);
238cf5a6209SStefano Zampini   }
239cf5a6209SStefano Zampini 
240674ae819SStefano Zampini   /* loop on ccs to compute index sets for faces and edges */
241acc113dbSStefano Zampini   if (!graph->queue_sorted) {
242acc113dbSStefano Zampini     PetscInt *queue_global;
243acc113dbSStefano Zampini 
244acc113dbSStefano Zampini     ierr = PetscMalloc1(graph->cptr[graph->ncc],&queue_global);CHKERRQ(ierr);
245acc113dbSStefano Zampini     ierr = ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_global);CHKERRQ(ierr);
246acc113dbSStefano Zampini     for (i=0;i<graph->ncc;i++) {
247acc113dbSStefano Zampini       ierr = PetscSortIntWithArray(graph->cptr[i+1]-graph->cptr[i],&queue_global[graph->cptr[i]],&graph->queue[graph->cptr[i]]);CHKERRQ(ierr);
248acc113dbSStefano Zampini     }
249acc113dbSStefano Zampini     ierr = PetscFree(queue_global);CHKERRQ(ierr);
250acc113dbSStefano Zampini     graph->queue_sorted = PETSC_TRUE;
251acc113dbSStefano Zampini   }
252674ae819SStefano Zampini   nfc = 0;
253674ae819SStefano Zampini   nec = 0;
254674ae819SStefano Zampini   for (i=0;i<graph->ncc;i++) {
255be12c134Sstefano_zampini     if (mark[i] == 2) {
256cf5a6209SStefano Zampini       if (FacesIS) {
257c8272957SStefano Zampini         ierr = ISCreateGeneral(PETSC_COMM_SELF,graph->cptr[i+1]-graph->cptr[i],&graph->queue[graph->cptr[i]],PETSC_USE_POINTER,&ISForFaces[nfc]);CHKERRQ(ierr);
258cf5a6209SStefano Zampini       }
259674ae819SStefano Zampini       nfc++;
260be12c134Sstefano_zampini     } else if (mark[i] == 1) {
261cf5a6209SStefano Zampini       if (EdgesIS) {
262c8272957SStefano Zampini         ierr = ISCreateGeneral(PETSC_COMM_SELF,graph->cptr[i+1]-graph->cptr[i],&graph->queue[graph->cptr[i]],PETSC_USE_POINTER,&ISForEdges[nec]);CHKERRQ(ierr);
263cf5a6209SStefano Zampini       }
264674ae819SStefano Zampini       nec++;
265674ae819SStefano Zampini     }
266674ae819SStefano Zampini   }
267be12c134Sstefano_zampini 
268674ae819SStefano Zampini   /* index set for vertices */
269cf5a6209SStefano Zampini   if (VerticesIS) {
270b8ffe317SStefano Zampini     nvc = 0;
271674ae819SStefano Zampini     for (i=0;i<graph->ncc;i++) {
272be12c134Sstefano_zampini       if (!mark[i]) {
273adfc4fb2SStefano Zampini         PetscInt j;
274adfc4fb2SStefano Zampini 
275674ae819SStefano Zampini         for (j=graph->cptr[i];j<graph->cptr[i+1];j++) {
276674ae819SStefano Zampini           idx[nvc]=graph->queue[j];
277674ae819SStefano Zampini           nvc++;
278674ae819SStefano Zampini         }
279674ae819SStefano Zampini       }
280674ae819SStefano Zampini     }
281674ae819SStefano Zampini     /* sort vertex set (by local ordering) */
282674ae819SStefano Zampini     ierr = PetscSortInt(nvc,idx);CHKERRQ(ierr);
283674ae819SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,nvc,idx,PETSC_OWN_POINTER,&ISForVertices);CHKERRQ(ierr);
284674ae819SStefano Zampini   }
285be12c134Sstefano_zampini   ierr = PetscFree(mark);CHKERRQ(ierr);
286be12c134Sstefano_zampini 
287674ae819SStefano Zampini   /* get back info */
288a873d5d0SStefano Zampini   if (n_faces)       *n_faces = nfc;
289c8272957SStefano Zampini   if (FacesIS)       *FacesIS = ISForFaces;
290a873d5d0SStefano Zampini   if (n_edges)       *n_edges = nec;
291c8272957SStefano Zampini   if (EdgesIS)       *EdgesIS = ISForEdges;
292c8272957SStefano Zampini   if (VerticesIS) *VerticesIS = ISForVertices;
293674ae819SStefano Zampini   PetscFunctionReturn(0);
294674ae819SStefano Zampini }
295674ae819SStefano Zampini 
296674ae819SStefano Zampini #undef __FUNCT__
297674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphComputeConnectedComponents"
298674ae819SStefano Zampini PetscErrorCode PCBDDCGraphComputeConnectedComponents(PCBDDCGraph graph)
299674ae819SStefano Zampini {
3004a060362SStefano Zampini   PetscBool      adapt_interface_reduced;
301674ae819SStefano Zampini   MPI_Comm       interface_comm;
3024a060362SStefano Zampini   PetscMPIInt    size;
3038e4af1ccSStefano Zampini   PetscInt       i;
304674ae819SStefano Zampini   PetscErrorCode ierr;
305674ae819SStefano Zampini 
306674ae819SStefano Zampini   PetscFunctionBegin;
307674ae819SStefano Zampini   /* compute connected components locally */
308674ae819SStefano Zampini   ierr = PetscObjectGetComm((PetscObject)(graph->l2gmap),&interface_comm);CHKERRQ(ierr);
309674ae819SStefano Zampini   ierr = PCBDDCGraphComputeConnectedComponentsLocal(graph);CHKERRQ(ierr);
310674ae819SStefano Zampini   /* check consistency of connected components among neighbouring subdomains -> it adapt them in case it is needed */
3114a060362SStefano Zampini   ierr = MPI_Comm_size(interface_comm,&size);CHKERRQ(ierr);
3124a060362SStefano Zampini   adapt_interface_reduced = PETSC_FALSE;
3134a060362SStefano Zampini   if (size > 1) {
3144a060362SStefano Zampini     PetscInt i;
3154a060362SStefano Zampini     PetscBool adapt_interface = PETSC_FALSE;
316674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
317674ae819SStefano Zampini       /* We are not sure that on a given subset of the local interface,
318674ae819SStefano Zampini          with two connected components, the latters be the same among sharing subdomains */
319674ae819SStefano Zampini       if (graph->subset_ncc[i] > 1) {
3204a060362SStefano Zampini         adapt_interface = PETSC_TRUE;
321674ae819SStefano Zampini         break;
322674ae819SStefano Zampini       }
323674ae819SStefano Zampini     }
324b2566f29SBarry Smith     ierr = MPIU_Allreduce(&adapt_interface,&adapt_interface_reduced,1,MPIU_BOOL,MPI_LOR,interface_comm);CHKERRQ(ierr);
3254a060362SStefano Zampini   }
326674ae819SStefano Zampini 
327674ae819SStefano Zampini   if (graph->n_subsets && adapt_interface_reduced) {
328ec1c413dSStefano Zampini     PetscBT     subset_cc_adapt;
329ec1c413dSStefano Zampini     MPI_Request *send_requests,*recv_requests;
330ec1c413dSStefano Zampini     PetscInt    *send_buffer,*recv_buffer;
331ec1c413dSStefano Zampini     PetscInt    sum_requests,start_of_recv,start_of_send;
332ec1c413dSStefano Zampini     PetscInt    *cum_recv_counts;
333ec1c413dSStefano Zampini     PetscInt    *labels;
334b3cdcd63SStefano Zampini     PetscInt    ncc,cum_queue,mss,mns,j,k,s;
3358875e3e6SStefano Zampini     PetscInt    **refine_buffer=NULL,*private_labels = NULL;
3365b1b9aeaSStefano Zampini 
337ec1c413dSStefano Zampini     ierr = PetscMalloc1(graph->nvtxs,&labels);CHKERRQ(ierr);
338ec1c413dSStefano Zampini     ierr = PetscMemzero(labels,graph->nvtxs*sizeof(*labels));CHKERRQ(ierr);
339b3cdcd63SStefano Zampini     for (i=0;i<graph->ncc;i++)
340b3cdcd63SStefano Zampini       for (j=graph->cptr[i];j<graph->cptr[i+1];j++)
341b3cdcd63SStefano Zampini         labels[graph->queue[j]] = i;
342b3cdcd63SStefano Zampini 
343674ae819SStefano Zampini     /* allocate some space */
344854ce69bSBarry Smith     ierr = PetscMalloc1(graph->n_subsets+1,&cum_recv_counts);CHKERRQ(ierr);
345674ae819SStefano Zampini     ierr = PetscMemzero(cum_recv_counts,(graph->n_subsets+1)*sizeof(*cum_recv_counts));CHKERRQ(ierr);
346b3cdcd63SStefano Zampini 
347674ae819SStefano Zampini     /* first count how many neighbours per connected component I will receive from */
348674ae819SStefano Zampini     cum_recv_counts[0] = 0;
349b3cdcd63SStefano 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]];
350ec1c413dSStefano Zampini     ierr = PetscMalloc1(cum_recv_counts[graph->n_subsets],&recv_buffer);CHKERRQ(ierr);
351ec1c413dSStefano Zampini     ierr = PetscMalloc2(cum_recv_counts[graph->n_subsets],&send_requests,cum_recv_counts[graph->n_subsets],&recv_requests);CHKERRQ(ierr);
352674ae819SStefano Zampini     for (i=0;i<cum_recv_counts[graph->n_subsets];i++) {
353674ae819SStefano Zampini       send_requests[i] = MPI_REQUEST_NULL;
354674ae819SStefano Zampini       recv_requests[i] = MPI_REQUEST_NULL;
355674ae819SStefano Zampini     }
356b3cdcd63SStefano Zampini 
357ec1c413dSStefano Zampini     /* exchange with my neighbours the number of my connected components on the subset of interface */
358674ae819SStefano Zampini     sum_requests = 0;
359674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
360ec1c413dSStefano Zampini       PetscMPIInt neigh,tag;
361b3cdcd63SStefano Zampini       PetscInt    count,*neighs;
362ec1c413dSStefano Zampini 
363b3cdcd63SStefano Zampini       count = graph->count[graph->subset_idxs[i][0]];
364b3cdcd63SStefano Zampini       neighs = graph->neighbours_set[graph->subset_idxs[i][0]];
365ec1c413dSStefano Zampini       ierr = PetscMPIIntCast(2*graph->subset_ref_node[i],&tag);CHKERRQ(ierr);
366b3cdcd63SStefano Zampini       for (k=0;k<count;k++) {
367b3cdcd63SStefano Zampini         ierr = PetscMPIIntCast(neighs[k],&neigh);CHKERRQ(ierr);
368674ae819SStefano Zampini         ierr = MPI_Isend(&graph->subset_ncc[i],1,MPIU_INT,neigh,tag,interface_comm,&send_requests[sum_requests]);CHKERRQ(ierr);
369ec1c413dSStefano Zampini         ierr = MPI_Irecv(&recv_buffer[sum_requests],1,MPIU_INT,neigh,tag,interface_comm,&recv_requests[sum_requests]);CHKERRQ(ierr);
370674ae819SStefano Zampini         sum_requests++;
371674ae819SStefano Zampini       }
372674ae819SStefano Zampini     }
373674ae819SStefano Zampini     ierr = MPI_Waitall(sum_requests,recv_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
374674ae819SStefano Zampini     ierr = MPI_Waitall(sum_requests,send_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
375b3cdcd63SStefano Zampini 
376b3cdcd63SStefano Zampini     /* determine the subsets I have to adapt (those having more than 1 cc) */
377ec1c413dSStefano Zampini     ierr = PetscBTCreate(graph->n_subsets,&subset_cc_adapt);CHKERRQ(ierr);
378ec1c413dSStefano Zampini     ierr = PetscBTMemzero(graph->n_subsets,subset_cc_adapt);CHKERRQ(ierr);
379674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
380b3cdcd63SStefano Zampini       if (graph->subset_ncc[i] > 1) {
381b3cdcd63SStefano Zampini         ierr = PetscBTSet(subset_cc_adapt,i);CHKERRQ(ierr);
382b3cdcd63SStefano Zampini         continue;
383b3cdcd63SStefano Zampini       }
384674ae819SStefano Zampini       for (j=cum_recv_counts[i];j<cum_recv_counts[i+1];j++){
385b3cdcd63SStefano Zampini          if (recv_buffer[j] > 1) {
386f913dca9SStefano Zampini           ierr = PetscBTSet(subset_cc_adapt,i);CHKERRQ(ierr);
387674ae819SStefano Zampini           break;
388674ae819SStefano Zampini         }
389674ae819SStefano Zampini       }
390674ae819SStefano Zampini     }
391ec1c413dSStefano Zampini     ierr = PetscFree(recv_buffer);CHKERRQ(ierr);
392b3cdcd63SStefano Zampini 
393ec1c413dSStefano Zampini     /* determine send/recv buffers sizes */
394ec1c413dSStefano Zampini     j = 0;
395b3cdcd63SStefano Zampini     mss = 0;
396674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
397ec1c413dSStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i)) {
398b3cdcd63SStefano Zampini         j += graph->subset_size[i];
399b3cdcd63SStefano Zampini         mss = PetscMax(graph->subset_size[i],mss);
400674ae819SStefano Zampini       }
401674ae819SStefano Zampini     }
402ec1c413dSStefano Zampini     k = 0;
403b3cdcd63SStefano Zampini     mns = 0;
404674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
405ec1c413dSStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i)) {
406b3cdcd63SStefano Zampini         k += (cum_recv_counts[i+1]-cum_recv_counts[i])*graph->subset_size[i];
407b3cdcd63SStefano Zampini         mns = PetscMax(cum_recv_counts[i+1]-cum_recv_counts[i],mns);
408674ae819SStefano Zampini       }
409674ae819SStefano Zampini     }
410ec1c413dSStefano Zampini     ierr = PetscMalloc2(j,&send_buffer,k,&recv_buffer);CHKERRQ(ierr);
411ec1c413dSStefano Zampini 
412b3cdcd63SStefano Zampini     /* fill send buffer (order matters: subset_idxs ordered by global ordering) */
413ec1c413dSStefano Zampini     j = 0;
414b3cdcd63SStefano Zampini     for (i=0;i<graph->n_subsets;i++)
415b3cdcd63SStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i))
416b3cdcd63SStefano Zampini         for (k=0;k<graph->subset_size[i];k++)
417b3cdcd63SStefano Zampini           send_buffer[j++] = labels[graph->subset_idxs[i][k]];
418ec1c413dSStefano Zampini 
419674ae819SStefano Zampini     /* now exchange the data */
420674ae819SStefano Zampini     start_of_recv = 0;
421674ae819SStefano Zampini     start_of_send = 0;
422674ae819SStefano Zampini     sum_requests = 0;
423674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
424ec1c413dSStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i)) {
425ec1c413dSStefano Zampini         PetscMPIInt neigh,tag;
426b3cdcd63SStefano Zampini         PetscInt    size_of_send = graph->subset_size[i];
427ec1c413dSStefano Zampini 
428b3cdcd63SStefano Zampini         j = graph->subset_idxs[i][0];
429ec1c413dSStefano Zampini         ierr = PetscMPIIntCast(2*graph->subset_ref_node[i]+1,&tag);CHKERRQ(ierr);
430674ae819SStefano Zampini         for (k=0;k<graph->count[j];k++) {
431674ae819SStefano Zampini           ierr = PetscMPIIntCast(graph->neighbours_set[j][k],&neigh);CHKERRQ(ierr);
432ec1c413dSStefano Zampini           ierr = MPI_Isend(&send_buffer[start_of_send],size_of_send,MPIU_INT,neigh,tag,interface_comm,&send_requests[sum_requests]);CHKERRQ(ierr);
433ec1c413dSStefano Zampini           ierr = MPI_Irecv(&recv_buffer[start_of_recv],size_of_send,MPIU_INT,neigh,tag,interface_comm,&recv_requests[sum_requests]);CHKERRQ(ierr);
434ec1c413dSStefano Zampini           start_of_recv += size_of_send;
435674ae819SStefano Zampini           sum_requests++;
436674ae819SStefano Zampini         }
437674ae819SStefano Zampini         start_of_send += size_of_send;
438674ae819SStefano Zampini       }
439674ae819SStefano Zampini     }
440674ae819SStefano Zampini     ierr = MPI_Waitall(sum_requests,recv_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
441d52c4852SStefano Zampini 
442b3cdcd63SStefano Zampini     /* refine connected components */
443674ae819SStefano Zampini     start_of_recv = 0;
444b3cdcd63SStefano Zampini     /* allocate some temporary space */
445b3cdcd63SStefano Zampini     if (mss) {
446b3cdcd63SStefano Zampini       ierr = PetscMalloc1(mss,&refine_buffer);CHKERRQ(ierr);
447b3cdcd63SStefano Zampini       ierr = PetscMalloc2(mss*(mns+1),&refine_buffer[0],mss,&private_labels);CHKERRQ(ierr);
448b3cdcd63SStefano Zampini     }
449b3cdcd63SStefano Zampini     ncc = 0;
450b3cdcd63SStefano Zampini     cum_queue = 0;
451b3cdcd63SStefano Zampini     graph->cptr[0] = 0;
452674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
453ec1c413dSStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i)) {
454b3cdcd63SStefano Zampini         PetscInt subset_counter = 0;
455b3cdcd63SStefano Zampini         PetscInt sharingprocs = cum_recv_counts[i+1]-cum_recv_counts[i]+1; /* count myself */
456b3cdcd63SStefano Zampini         PetscInt buffer_size = graph->subset_size[i];
457ec1c413dSStefano Zampini 
458b3cdcd63SStefano Zampini         /* compute pointers */
459b3cdcd63SStefano Zampini         for (j=1;j<buffer_size;j++) refine_buffer[j] = refine_buffer[j-1] + sharingprocs;
460b3cdcd63SStefano Zampini         /* analyze contributions from subdomains that share the i-th subset
461b3cdcd63SStefano Zampini            The stricture of refine_buffer is suitable to find intersections of ccs among sharingprocs.
462b3cdcd63SStefano 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)
463b3cdcd63SStefano Zampini            sharing procs connected components:
464ec1c413dSStefano Zampini              neigh 0: [0 1 4], [2 3], labels [4,7]  (2 connected components)
465ec1c413dSStefano Zampini              neigh 1: [0 1], [2 3 4], labels [3 2]  (2 connected components)
466ec1c413dSStefano Zampini              neigh 2: [0 4], [1], [2 3], labels [1 5 6] (3 connected components)
467b3cdcd63SStefano Zampini            refine_buffer will be filled as:
468ec1c413dSStefano Zampini              [ 4, 3, 1;
469ec1c413dSStefano Zampini                4, 2, 1;
470ec1c413dSStefano Zampini                7, 2, 6;
471ec1c413dSStefano Zampini                4, 3, 5;
472ec1c413dSStefano Zampini                7, 2, 6; ];
473b3cdcd63SStefano Zampini            The connected components in local ordering are [0], [1], [2 3], [4] */
474b3cdcd63SStefano Zampini         /* fill temp_buffer */
475b3cdcd63SStefano Zampini         for (k=0;k<buffer_size;k++) refine_buffer[k][0] = labels[graph->subset_idxs[i][k]];
476b3cdcd63SStefano Zampini         for (j=0;j<sharingprocs-1;j++) {
477b3cdcd63SStefano Zampini           for (k=0;k<buffer_size;k++) refine_buffer[k][j+1] = recv_buffer[start_of_recv+k];
478b3cdcd63SStefano Zampini           start_of_recv += buffer_size;
479674ae819SStefano Zampini         }
480b3cdcd63SStefano Zampini         ierr = PetscMemzero(private_labels,buffer_size*sizeof(PetscInt));CHKERRQ(ierr);
481b3cdcd63SStefano Zampini         for (j=0;j<buffer_size;j++) {
482b3cdcd63SStefano Zampini           if (!private_labels[j]) { /* found a new cc  */
483ec1c413dSStefano Zampini             PetscBool same_set;
484ec1c413dSStefano Zampini 
485b3cdcd63SStefano Zampini             graph->cptr[ncc] = cum_queue;
486b3cdcd63SStefano Zampini             ncc++;
487b3cdcd63SStefano Zampini             subset_counter++;
488b3cdcd63SStefano Zampini             private_labels[j] = subset_counter;
489b3cdcd63SStefano Zampini             graph->queue[cum_queue++] = graph->subset_idxs[i][j];
490b3cdcd63SStefano Zampini             for (k=j+1;k<buffer_size;k++) { /* check for other nodes in new cc */
491674ae819SStefano Zampini               same_set = PETSC_TRUE;
492b3cdcd63SStefano Zampini               for (s=0;s<sharingprocs;s++) {
493b3cdcd63SStefano Zampini                 if (refine_buffer[j][s] != refine_buffer[k][s]) {
494674ae819SStefano Zampini                   same_set = PETSC_FALSE;
495674ae819SStefano Zampini                   break;
496674ae819SStefano Zampini                 }
497674ae819SStefano Zampini               }
498674ae819SStefano Zampini               if (same_set) {
499b3cdcd63SStefano Zampini                 private_labels[k] = subset_counter;
500b3cdcd63SStefano Zampini                 graph->queue[cum_queue++] = graph->subset_idxs[i][k];
501674ae819SStefano Zampini               }
502674ae819SStefano Zampini             }
503674ae819SStefano Zampini           }
504674ae819SStefano Zampini         }
505b3cdcd63SStefano Zampini         graph->cptr[ncc] = cum_queue;
506b3cdcd63SStefano Zampini         graph->subset_ncc[i] = subset_counter;
507b3cdcd63SStefano Zampini         graph->queue_sorted = PETSC_FALSE;
508b3cdcd63SStefano Zampini       } else { /* this subset does not need to be adapted */
509b3cdcd63SStefano Zampini         ierr = PetscMemcpy(graph->queue+cum_queue,graph->subset_idxs[i],graph->subset_size[i]*sizeof(PetscInt));CHKERRQ(ierr);
510b3cdcd63SStefano Zampini         ncc++;
511b3cdcd63SStefano Zampini         cum_queue += graph->subset_size[i];
512b3cdcd63SStefano Zampini         graph->cptr[ncc] = cum_queue;
513674ae819SStefano Zampini       }
514674ae819SStefano Zampini     }
515b3cdcd63SStefano Zampini     graph->cptr[ncc] = cum_queue;
516b3cdcd63SStefano Zampini     graph->ncc = ncc;
517b3cdcd63SStefano Zampini     if (mss) {
518b3cdcd63SStefano Zampini       ierr = PetscFree2(refine_buffer[0],private_labels);CHKERRQ(ierr);
519b3cdcd63SStefano Zampini       ierr = PetscFree(refine_buffer);CHKERRQ(ierr);
520b3cdcd63SStefano Zampini     }
521b3cdcd63SStefano Zampini     ierr = PetscFree(labels);CHKERRQ(ierr);
522b3cdcd63SStefano Zampini     ierr = MPI_Waitall(sum_requests,send_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
523b3cdcd63SStefano Zampini     ierr = PetscFree2(send_requests,recv_requests);CHKERRQ(ierr);
52445b2a5aaSStefano Zampini     ierr = PetscFree2(send_buffer,recv_buffer);CHKERRQ(ierr);
525674ae819SStefano Zampini     ierr = PetscFree(cum_recv_counts);CHKERRQ(ierr);
526ec1c413dSStefano Zampini     ierr = PetscBTDestroy(&subset_cc_adapt);CHKERRQ(ierr);
527674ae819SStefano Zampini   }
5288e4af1ccSStefano Zampini 
529674ae819SStefano Zampini   PetscFunctionReturn(0);
530674ae819SStefano Zampini }
531674ae819SStefano Zampini 
532b3cdcd63SStefano Zampini 
533b3cdcd63SStefano Zampini #undef __FUNCT__
534b3cdcd63SStefano Zampini #define __FUNCT__ "PCBDDCGraphComputeCC_Private"
535b3cdcd63SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCGraphComputeCC_Private(PCBDDCGraph graph,PetscInt pid,PetscInt* queue_tip,PetscInt n_prev,PetscInt* n_added)
536b3cdcd63SStefano Zampini {
537b3cdcd63SStefano Zampini   PetscInt       i,j,n;
538b3cdcd63SStefano Zampini   PetscInt       *xadj = graph->xadj,*adjncy = graph->adjncy;
539b3cdcd63SStefano Zampini   PetscBT        touched = graph->touched;
540d895b288Sstefano_zampini   PetscBool      havecsr = (PetscBool)(!!xadj);
541b3cdcd63SStefano Zampini   PetscBool      havesubs = (PetscBool)(!!graph->n_local_subs);
542b3cdcd63SStefano Zampini   PetscErrorCode ierr;
543b3cdcd63SStefano Zampini 
544b3cdcd63SStefano Zampini   PetscFunctionBegin;
545b3cdcd63SStefano Zampini   n = 0;
546b3cdcd63SStefano Zampini   if (havecsr && !havesubs) {
547b3cdcd63SStefano Zampini     for (i=-n_prev;i<0;i++) {
548b3cdcd63SStefano Zampini       PetscInt start_dof = queue_tip[i];
54954fffbccSStefano 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 */
55054fffbccSStefano Zampini       if (xadj[start_dof+1]-xadj[start_dof] == 1 && adjncy[xadj[start_dof]] < 0) {
55154fffbccSStefano Zampini         for (j=0;j<graph->subset_size[pid-1];j++) { /* pid \in [1,graph->n_subsets] */
55254fffbccSStefano Zampini           PetscInt dof = graph->subset_idxs[pid-1][j];
55354fffbccSStefano Zampini           if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid) {
55454fffbccSStefano Zampini             ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
55554fffbccSStefano Zampini             queue_tip[n] = dof;
55654fffbccSStefano Zampini             n++;
55754fffbccSStefano Zampini           }
55854fffbccSStefano Zampini         }
55954fffbccSStefano Zampini       } else {
560b3cdcd63SStefano Zampini         for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
561b3cdcd63SStefano Zampini           PetscInt dof = adjncy[j];
562b3cdcd63SStefano Zampini           if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid) {
563b3cdcd63SStefano Zampini             ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
564b3cdcd63SStefano Zampini             queue_tip[n] = dof;
565b3cdcd63SStefano Zampini             n++;
566b3cdcd63SStefano Zampini           }
567b3cdcd63SStefano Zampini         }
568b3cdcd63SStefano Zampini       }
56954fffbccSStefano Zampini     }
570b3cdcd63SStefano Zampini   } else if (havecsr && havesubs) {
571b3cdcd63SStefano Zampini     PetscInt sid = graph->local_subs[queue_tip[-n_prev]];
572b3cdcd63SStefano Zampini     for (i=-n_prev;i<0;i++) {
573b3cdcd63SStefano Zampini       PetscInt start_dof = queue_tip[i];
57454fffbccSStefano 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 */
57554fffbccSStefano Zampini       if (xadj[start_dof+1]-xadj[start_dof] == 1 && adjncy[xadj[start_dof]] < 0) {
57654fffbccSStefano Zampini         for (j=0;j<graph->subset_size[pid-1];j++) { /* pid \in [1,graph->n_subsets] */
57754fffbccSStefano Zampini           PetscInt dof = graph->subset_idxs[pid-1][j];
57854fffbccSStefano Zampini           if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid && graph->local_subs[dof] == sid) {
57954fffbccSStefano Zampini             ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
58054fffbccSStefano Zampini             queue_tip[n] = dof;
58154fffbccSStefano Zampini             n++;
58254fffbccSStefano Zampini           }
58354fffbccSStefano Zampini         }
58454fffbccSStefano Zampini       } else {
585b3cdcd63SStefano Zampini         for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
586b3cdcd63SStefano Zampini           PetscInt dof = adjncy[j];
587b3cdcd63SStefano Zampini           if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid && graph->local_subs[dof] == sid) {
588b3cdcd63SStefano Zampini             ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
589b3cdcd63SStefano Zampini             queue_tip[n] = dof;
590b3cdcd63SStefano Zampini             n++;
591b3cdcd63SStefano Zampini           }
592b3cdcd63SStefano Zampini         }
593b3cdcd63SStefano Zampini       }
59454fffbccSStefano Zampini     }
595b3cdcd63SStefano Zampini   } else { /* sub info only */
596b3cdcd63SStefano Zampini     PetscInt sid = graph->local_subs[queue_tip[-n_prev]];
597b3cdcd63SStefano Zampini     for (j=0;j<graph->subset_size[pid-1];j++) { /* pid \in [1,graph->n_subsets] */
598b3cdcd63SStefano Zampini       PetscInt dof = graph->subset_idxs[pid-1][j];
599b3cdcd63SStefano Zampini       if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid && graph->local_subs[dof] == sid) {
600b3cdcd63SStefano Zampini         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
601b3cdcd63SStefano Zampini         queue_tip[n] = dof;
602b3cdcd63SStefano Zampini         n++;
603b3cdcd63SStefano Zampini       }
604b3cdcd63SStefano Zampini     }
605b3cdcd63SStefano Zampini   }
606b3cdcd63SStefano Zampini   *n_added = n;
607b3cdcd63SStefano Zampini   PetscFunctionReturn(0);
608b3cdcd63SStefano Zampini }
609b3cdcd63SStefano Zampini 
610674ae819SStefano Zampini #undef __FUNCT__
611674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphComputeConnectedComponentsLocal"
612674ae819SStefano Zampini PetscErrorCode PCBDDCGraphComputeConnectedComponentsLocal(PCBDDCGraph graph)
613674ae819SStefano Zampini {
614b3cdcd63SStefano Zampini   PetscInt       ncc,cum_queue,n;
615b3cdcd63SStefano Zampini   PetscMPIInt    commsize;
616674ae819SStefano Zampini   PetscErrorCode ierr;
617674ae819SStefano Zampini 
618674ae819SStefano Zampini   PetscFunctionBegin;
619b3cdcd63SStefano Zampini   if (!graph->setupcalled) SETERRQ(PetscObjectComm((PetscObject)graph->l2gmap),PETSC_ERR_ORDER,"PCBDDCGraphSetUp should be called first");
620b3cdcd63SStefano Zampini   /* quiet return if there isn't any local info */
6211633d1f0SStefano Zampini   if (!graph->xadj && !graph->n_local_subs) {
622674ae819SStefano Zampini     PetscFunctionReturn(0);
623674ae819SStefano Zampini   }
624674ae819SStefano Zampini 
625674ae819SStefano Zampini   /* reset any previous search of connected components */
626df48933bSStefano Zampini   ierr = PetscBTMemzero(graph->nvtxs,graph->touched);CHKERRQ(ierr);
627b3cdcd63SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)graph->l2gmap),&commsize);CHKERRQ(ierr);
6284b2aedd3SStefano Zampini   if (commsize > graph->commsizelimit) {
629b3cdcd63SStefano Zampini     PetscInt i;
630674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
6310cf82fd6SStefano Zampini       if (graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK || !graph->count[i]) {
632df48933bSStefano Zampini         ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
633674ae819SStefano Zampini       }
6344a060362SStefano Zampini     }
635b3cdcd63SStefano Zampini   }
636674ae819SStefano Zampini 
637674ae819SStefano Zampini   /* begin search for connected components */
638674ae819SStefano Zampini   cum_queue = 0;
639674ae819SStefano Zampini   ncc = 0;
640674ae819SStefano Zampini   for (n=0;n<graph->n_subsets;n++) {
641b3cdcd63SStefano Zampini     PetscInt pid = n+1;  /* partition labeled by 0 is discarded */
642b3cdcd63SStefano Zampini     PetscInt found = 0,prev = 0,first = 0,ncc_pid = 0;
643b3cdcd63SStefano Zampini     while (found != graph->subset_size[n]) {
644b3cdcd63SStefano Zampini       PetscInt added = 0;
645b3cdcd63SStefano Zampini       if (!prev) { /* search for new starting dof */
646b3cdcd63SStefano Zampini         while (PetscBTLookup(graph->touched,graph->subset_idxs[n][first])) first++;
647b3cdcd63SStefano Zampini         ierr = PetscBTSet(graph->touched,graph->subset_idxs[n][first]);CHKERRQ(ierr);
648b3cdcd63SStefano Zampini         graph->queue[cum_queue] = graph->subset_idxs[n][first];
649674ae819SStefano Zampini         graph->cptr[ncc] = cum_queue;
650b3cdcd63SStefano Zampini         prev = 1;
651b3cdcd63SStefano Zampini         cum_queue++;
652b3cdcd63SStefano Zampini         found++;
653674ae819SStefano Zampini         ncc_pid++;
654b3cdcd63SStefano Zampini         ncc++;
655674ae819SStefano Zampini       }
656b3cdcd63SStefano Zampini       ierr = PCBDDCGraphComputeCC_Private(graph,pid,graph->queue + cum_queue,prev,&added);CHKERRQ(ierr);
657b3cdcd63SStefano Zampini       if (!added) {
658674ae819SStefano Zampini         graph->subset_ncc[n] = ncc_pid;
659b3cdcd63SStefano Zampini         graph->cptr[ncc] = cum_queue;
660b3cdcd63SStefano Zampini       }
661b3cdcd63SStefano Zampini       prev = added;
662b3cdcd63SStefano Zampini       found += added;
663b3cdcd63SStefano Zampini       cum_queue += added;
664b3cdcd63SStefano Zampini       if (added && found == graph->subset_size[n]) {
665b3cdcd63SStefano Zampini         graph->subset_ncc[n] = ncc_pid;
666b3cdcd63SStefano Zampini         graph->cptr[ncc] = cum_queue;
667b3cdcd63SStefano Zampini       }
668b3cdcd63SStefano Zampini     }
669674ae819SStefano Zampini   }
670674ae819SStefano Zampini   graph->ncc = ncc;
671acc113dbSStefano Zampini   graph->queue_sorted = PETSC_FALSE;
672674ae819SStefano Zampini   PetscFunctionReturn(0);
673674ae819SStefano Zampini }
674674ae819SStefano Zampini 
675674ae819SStefano Zampini #undef __FUNCT__
676674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphSetUp"
677674ae819SStefano Zampini PetscErrorCode PCBDDCGraphSetUp(PCBDDCGraph graph, PetscInt custom_minimal_size, IS neumann_is, IS dirichlet_is, PetscInt n_ISForDofs, IS ISForDofs[], IS custom_primal_vertices)
678674ae819SStefano Zampini {
679*7a0e7b2cSstefano_zampini   IS             subset,subset_n;
680674ae819SStefano Zampini   MPI_Comm       comm;
681674ae819SStefano Zampini   const PetscInt *is_indices;
682dc456d91SStefano Zampini   PetscInt       n_neigh,*neigh,*n_shared,**shared,*queue_global;
683674ae819SStefano Zampini   PetscInt       i,j,k,s,total_counts,nodes_touched,is_size;
684b3cdcd63SStefano Zampini   PetscMPIInt    commsize;
6853ebc787cSStefano Zampini   PetscBool      same_set,mirrors_found,twodim;
6867babac9bSStefano Zampini   PetscErrorCode ierr;
687674ae819SStefano Zampini 
688674ae819SStefano Zampini   PetscFunctionBegin;
689*7a0e7b2cSstefano_zampini   PetscValidLogicalCollectiveInt(graph->l2gmap,custom_minimal_size,2);
690*7a0e7b2cSstefano_zampini   if (neumann_is) {
691*7a0e7b2cSstefano_zampini     PetscValidHeaderSpecific(neumann_is,IS_CLASSID,3);
692*7a0e7b2cSstefano_zampini     PetscCheckSameComm(graph->l2gmap,1,neumann_is,3);
693*7a0e7b2cSstefano_zampini   }
694a07ea27aSStefano Zampini   graph->has_dirichlet = PETSC_FALSE;
695a07ea27aSStefano Zampini   if (dirichlet_is) {
696*7a0e7b2cSstefano_zampini     PetscValidHeaderSpecific(dirichlet_is,IS_CLASSID,4);
697a07ea27aSStefano Zampini     PetscCheckSameComm(graph->l2gmap,1,dirichlet_is,4);
698a07ea27aSStefano Zampini     graph->has_dirichlet = PETSC_TRUE;
699a07ea27aSStefano Zampini   }
700*7a0e7b2cSstefano_zampini   PetscValidLogicalCollectiveInt(graph->l2gmap,n_ISForDofs,5);
701*7a0e7b2cSstefano_zampini   for (i=0;i<n_ISForDofs;i++) {
702*7a0e7b2cSstefano_zampini     PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,6);
703*7a0e7b2cSstefano_zampini     PetscCheckSameComm(graph->l2gmap,1,ISForDofs[i],6);
704*7a0e7b2cSstefano_zampini   }
705*7a0e7b2cSstefano_zampini   if (custom_primal_vertices) {
706*7a0e7b2cSstefano_zampini     PetscValidHeaderSpecific(custom_primal_vertices,IS_CLASSID,6);
707*7a0e7b2cSstefano_zampini     PetscCheckSameComm(graph->l2gmap,1,custom_primal_vertices,7);
708*7a0e7b2cSstefano_zampini   }
709674ae819SStefano Zampini   ierr = PetscObjectGetComm((PetscObject)(graph->l2gmap),&comm);CHKERRQ(ierr);
710b3cdcd63SStefano Zampini   ierr = MPI_Comm_size(comm,&commsize);CHKERRQ(ierr);
711b3cdcd63SStefano Zampini 
712674ae819SStefano Zampini   /* custom_minimal_size */
71314f95afaSStefano Zampini   graph->custom_minimal_size = custom_minimal_size;
714674ae819SStefano Zampini   /* get info l2gmap and allocate work vectors  */
715674ae819SStefano Zampini   ierr = ISLocalToGlobalMappingGetInfo(graph->l2gmap,&n_neigh,&neigh,&n_shared,&shared);CHKERRQ(ierr);
7162635a1d4SStefano Zampini   /* check if we have any local periodic nodes (periodic BCs) */
7172635a1d4SStefano Zampini   mirrors_found = PETSC_FALSE;
7182635a1d4SStefano Zampini   if (graph->nvtxs && n_neigh) {
7192635a1d4SStefano Zampini     for (i=0; i<n_shared[0]; i++) graph->count[shared[0][i]] += 1;
7202635a1d4SStefano Zampini     for (i=0; i<n_shared[0]; i++) {
7212635a1d4SStefano Zampini       if (graph->count[shared[0][i]] > 1) {
7222635a1d4SStefano Zampini         mirrors_found = PETSC_TRUE;
7232635a1d4SStefano Zampini         break;
7242635a1d4SStefano Zampini       }
7252635a1d4SStefano Zampini     }
7262635a1d4SStefano Zampini   }
72751b0f6cfSStefano Zampini   /* compute local mirrors (if any) */
72851b0f6cfSStefano Zampini   if (mirrors_found) {
729*7a0e7b2cSstefano_zampini     IS       to,from;
73051b0f6cfSStefano Zampini     PetscInt *local_indices,*global_indices;
731*7a0e7b2cSstefano_zampini 
732*7a0e7b2cSstefano_zampini     ierr = ISCreateStride(PETSC_COMM_SELF,graph->nvtxs,0,1,&to);CHKERRQ(ierr);
733*7a0e7b2cSstefano_zampini     ierr = ISLocalToGlobalMappingApplyIS(graph->l2gmap,to,&from);CHKERRQ(ierr);
73451b0f6cfSStefano Zampini     /* get arrays of local and global indices */
735785e854fSJed Brown     ierr = PetscMalloc1(graph->nvtxs,&local_indices);CHKERRQ(ierr);
73651b0f6cfSStefano Zampini     ierr = ISGetIndices(to,(const PetscInt**)&is_indices);CHKERRQ(ierr);
73751b0f6cfSStefano Zampini     ierr = PetscMemcpy(local_indices,is_indices,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
73851b0f6cfSStefano Zampini     ierr = ISRestoreIndices(to,(const PetscInt**)&is_indices);CHKERRQ(ierr);
739785e854fSJed Brown     ierr = PetscMalloc1(graph->nvtxs,&global_indices);CHKERRQ(ierr);
74051b0f6cfSStefano Zampini     ierr = ISGetIndices(from,(const PetscInt**)&is_indices);CHKERRQ(ierr);
74151b0f6cfSStefano Zampini     ierr = PetscMemcpy(global_indices,is_indices,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
74251b0f6cfSStefano Zampini     ierr = ISRestoreIndices(from,(const PetscInt**)&is_indices);CHKERRQ(ierr);
74351b0f6cfSStefano Zampini     /* allocate space for mirrors */
7448e5aaad5SJed Brown     ierr = PetscMalloc2(graph->nvtxs,&graph->mirrors,graph->nvtxs,&graph->mirrors_set);CHKERRQ(ierr);
74551b0f6cfSStefano Zampini     ierr = PetscMemzero(graph->mirrors,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
74651b0f6cfSStefano Zampini     graph->mirrors_set[0] = 0;
74751b0f6cfSStefano Zampini 
74851b0f6cfSStefano Zampini     k=0;
74951b0f6cfSStefano Zampini     for (i=0;i<n_shared[0];i++) {
75051b0f6cfSStefano Zampini       j=shared[0][i];
75151b0f6cfSStefano Zampini       if (graph->count[j] > 1) {
75251b0f6cfSStefano Zampini         graph->mirrors[j]++;
75351b0f6cfSStefano Zampini         k++;
75451b0f6cfSStefano Zampini       }
75551b0f6cfSStefano Zampini     }
75651b0f6cfSStefano Zampini     /* allocate space for set of mirrors */
757785e854fSJed Brown     ierr = PetscMalloc1(k,&graph->mirrors_set[0]);CHKERRQ(ierr);
75851b0f6cfSStefano Zampini     for (i=1;i<graph->nvtxs;i++)
75951b0f6cfSStefano Zampini       graph->mirrors_set[i]=graph->mirrors_set[i-1]+graph->mirrors[i-1];
76051b0f6cfSStefano Zampini 
76151b0f6cfSStefano Zampini     /* fill arrays */
76251b0f6cfSStefano Zampini     ierr = PetscMemzero(graph->mirrors,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
76351b0f6cfSStefano Zampini     for (j=0;j<n_shared[0];j++) {
76451b0f6cfSStefano Zampini       i=shared[0][j];
76551b0f6cfSStefano Zampini       if (graph->count[i] > 1)
76651b0f6cfSStefano Zampini         graph->mirrors_set[i][graph->mirrors[i]++]=global_indices[i];
76751b0f6cfSStefano Zampini     }
76851b0f6cfSStefano Zampini     ierr = PetscSortIntWithArray(graph->nvtxs,global_indices,local_indices);CHKERRQ(ierr);
76951b0f6cfSStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
77051b0f6cfSStefano Zampini       if (graph->mirrors[i] > 0) {
77151b0f6cfSStefano Zampini         ierr = PetscFindInt(graph->mirrors_set[i][0],graph->nvtxs,global_indices,&k);CHKERRQ(ierr);
77251b0f6cfSStefano Zampini         j = global_indices[k];
77351b0f6cfSStefano Zampini         while ( k > 0 && global_indices[k-1] == j) k--;
77451b0f6cfSStefano Zampini         for (j=0;j<graph->mirrors[i];j++) {
77551b0f6cfSStefano Zampini           graph->mirrors_set[i][j]=local_indices[k+j];
77651b0f6cfSStefano Zampini         }
77751b0f6cfSStefano Zampini         ierr = PetscSortInt(graph->mirrors[i],graph->mirrors_set[i]);CHKERRQ(ierr);
77851b0f6cfSStefano Zampini       }
77951b0f6cfSStefano Zampini     }
78051b0f6cfSStefano Zampini     ierr = PetscFree(local_indices);CHKERRQ(ierr);
78151b0f6cfSStefano Zampini     ierr = PetscFree(global_indices);CHKERRQ(ierr);
782674ae819SStefano Zampini     ierr = ISDestroy(&to);CHKERRQ(ierr);
783674ae819SStefano Zampini     ierr = ISDestroy(&from);CHKERRQ(ierr);
784*7a0e7b2cSstefano_zampini   }
785*7a0e7b2cSstefano_zampini   ierr = PetscMemzero(graph->count,graph->nvtxs*sizeof(*graph->count));CHKERRQ(ierr);
78651b0f6cfSStefano Zampini 
787674ae819SStefano Zampini   /* Count total number of neigh per node */
788674ae819SStefano Zampini   k = 0;
789674ae819SStefano Zampini   for (i=1;i<n_neigh;i++) {
790674ae819SStefano Zampini     k += n_shared[i];
791674ae819SStefano Zampini     for (j=0;j<n_shared[i];j++) {
792674ae819SStefano Zampini       graph->count[shared[i][j]] += 1;
793674ae819SStefano Zampini     }
794674ae819SStefano Zampini   }
795674ae819SStefano Zampini   /* Allocate space for storing the set of neighbours for each node */
796674ae819SStefano Zampini   if (graph->nvtxs) {
797785e854fSJed Brown     ierr = PetscMalloc1(k,&graph->neighbours_set[0]);CHKERRQ(ierr);
798674ae819SStefano Zampini   }
799674ae819SStefano Zampini   for (i=1;i<graph->nvtxs;i++) { /* dont count myself */
800674ae819SStefano Zampini     graph->neighbours_set[i]=graph->neighbours_set[i-1]+graph->count[i-1];
801674ae819SStefano Zampini   }
802674ae819SStefano Zampini   /* Get information for sharing subdomains */
803674ae819SStefano Zampini   ierr = PetscMemzero(graph->count,graph->nvtxs*sizeof(*graph->count));CHKERRQ(ierr);
804674ae819SStefano Zampini   for (i=1;i<n_neigh;i++) { /* dont count myself */
805674ae819SStefano Zampini     s = n_shared[i];
806674ae819SStefano Zampini     for (j=0;j<s;j++) {
807674ae819SStefano Zampini       k = shared[i][j];
808674ae819SStefano Zampini       graph->neighbours_set[k][graph->count[k]] = neigh[i];
809674ae819SStefano Zampini       graph->count[k] += 1;
810674ae819SStefano Zampini     }
811674ae819SStefano Zampini   }
812674ae819SStefano Zampini   /* sort set of sharing subdomains */
813674ae819SStefano Zampini   for (i=0;i<graph->nvtxs;i++) {
81493fb5fd3SStefano Zampini     ierr = PetscSortRemoveDupsInt(&graph->count[i],graph->neighbours_set[i]);CHKERRQ(ierr);
815674ae819SStefano Zampini   }
8167fb0e2dbSStefano Zampini   /* free memory allocated by ISLocalToGlobalMappingGetInfo */
8177fb0e2dbSStefano Zampini   ierr = ISLocalToGlobalMappingRestoreInfo(graph->l2gmap,&n_neigh,&neigh,&n_shared,&shared);CHKERRQ(ierr);
8187fb0e2dbSStefano Zampini 
81967c9da69SStefano Zampini   /*
82067c9da69SStefano Zampini      Get info for dofs splitting
8215777c63cSStefano Zampini      User can specify just a subset; an additional field is considered as a complementary field
82267c9da69SStefano Zampini   */
823b3cdcd63SStefano Zampini   for (i=0;i<graph->nvtxs;i++) graph->which_dof[i] = n_ISForDofs; /* by default a dof belongs to the complement set */
824674ae819SStefano Zampini   for (i=0;i<n_ISForDofs;i++) {
82563602bcaSStefano Zampini     ierr = ISGetLocalSize(ISForDofs[i],&is_size);CHKERRQ(ierr);
826674ae819SStefano Zampini     ierr = ISGetIndices(ISForDofs[i],(const PetscInt**)&is_indices);CHKERRQ(ierr);
827674ae819SStefano Zampini     for (j=0;j<is_size;j++) {
828d50ed60aSStefano Zampini       if (is_indices[j] > -1 && is_indices[j] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
829674ae819SStefano Zampini         graph->which_dof[is_indices[j]] = i;
830674ae819SStefano Zampini       }
83167c9da69SStefano Zampini     }
832674ae819SStefano Zampini     ierr = ISRestoreIndices(ISForDofs[i],(const PetscInt**)&is_indices);CHKERRQ(ierr);
8335777c63cSStefano Zampini   }
8345777c63cSStefano Zampini 
835674ae819SStefano Zampini   /* Take into account Neumann nodes */
836674ae819SStefano Zampini   if (neumann_is) {
83782ba6b80SStefano Zampini     ierr = ISGetLocalSize(neumann_is,&is_size);CHKERRQ(ierr);
838674ae819SStefano Zampini     ierr = ISGetIndices(neumann_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
839674ae819SStefano Zampini     for (i=0;i<is_size;i++) {
840d50ed60aSStefano Zampini       if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
841*7a0e7b2cSstefano_zampini         graph->special_dof[is_indices[i]] = PCBDDCGRAPH_NEUMANN_MARK;
842674ae819SStefano Zampini       }
8433c629590SStefano Zampini     }
844674ae819SStefano Zampini     ierr = ISRestoreIndices(neumann_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
845674ae819SStefano Zampini   }
846*7a0e7b2cSstefano_zampini   /* Take into account Dirichlet nodes (they overwrite any neumann boundary mark previously set) */
847674ae819SStefano Zampini   if (dirichlet_is) {
84882ba6b80SStefano Zampini     ierr = ISGetLocalSize(dirichlet_is,&is_size);CHKERRQ(ierr);
849674ae819SStefano Zampini     ierr = ISGetIndices(dirichlet_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
850674ae819SStefano Zampini     for (i=0;i<is_size;i++){
851d50ed60aSStefano Zampini       if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
852*7a0e7b2cSstefano_zampini         if (commsize > graph->commsizelimit) { /* dirichlet nodes treated as internal */
853*7a0e7b2cSstefano_zampini           ierr = PetscBTSet(graph->touched,is_indices[i]);CHKERRQ(ierr);
854*7a0e7b2cSstefano_zampini           graph->subset[is_indices[i]] = 0;
855*7a0e7b2cSstefano_zampini         }
856*7a0e7b2cSstefano_zampini         graph->special_dof[is_indices[i]] = PCBDDCGRAPH_DIRICHLET_MARK;
857674ae819SStefano Zampini       }
8583c629590SStefano Zampini     }
859674ae819SStefano Zampini     ierr = ISRestoreIndices(dirichlet_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
8607fb0e2dbSStefano Zampini   }
86108a5cf49SStefano Zampini   /* mark local periodic nodes (if any) and adapt CSR graph (if any) */
86251b0f6cfSStefano Zampini   if (graph->mirrors) {
86351b0f6cfSStefano Zampini     for (i=0;i<graph->nvtxs;i++)
86451b0f6cfSStefano Zampini       if (graph->mirrors[i])
8650cf82fd6SStefano Zampini         graph->special_dof[i] = PCBDDCGRAPH_LOCAL_PERIODIC_MARK;
86651b0f6cfSStefano Zampini 
8671633d1f0SStefano Zampini     if (graph->xadj) {
86808a5cf49SStefano Zampini       PetscInt *new_xadj,*new_adjncy;
86951b0f6cfSStefano Zampini       /* sort CSR graph */
87051b0f6cfSStefano Zampini       for (i=0;i<graph->nvtxs;i++)
87151b0f6cfSStefano Zampini         ierr = PetscSortInt(graph->xadj[i+1]-graph->xadj[i],&graph->adjncy[graph->xadj[i]]);CHKERRQ(ierr);
87251b0f6cfSStefano Zampini 
87351b0f6cfSStefano Zampini       /* adapt local CSR graph in case of local periodicity */
87451b0f6cfSStefano Zampini       k = 0;
87551b0f6cfSStefano Zampini       for (i=0;i<graph->nvtxs;i++)
87651b0f6cfSStefano Zampini         for (j=graph->xadj[i];j<graph->xadj[i+1];j++)
87751b0f6cfSStefano Zampini           k += graph->mirrors[graph->adjncy[j]];
87851b0f6cfSStefano Zampini 
879854ce69bSBarry Smith       ierr = PetscMalloc1(graph->nvtxs+1,&new_xadj);CHKERRQ(ierr);
880854ce69bSBarry Smith       ierr = PetscMalloc1(k+graph->xadj[graph->nvtxs],&new_adjncy);CHKERRQ(ierr);
88151b0f6cfSStefano Zampini       new_xadj[0] = 0;
88251b0f6cfSStefano Zampini       for (i=0;i<graph->nvtxs;i++) {
88351b0f6cfSStefano Zampini         k = graph->xadj[i+1]-graph->xadj[i];
88451b0f6cfSStefano Zampini         ierr = PetscMemcpy(&new_adjncy[new_xadj[i]],&graph->adjncy[graph->xadj[i]],k*sizeof(PetscInt));CHKERRQ(ierr);
88551b0f6cfSStefano Zampini         new_xadj[i+1] = new_xadj[i]+k;
88651b0f6cfSStefano Zampini         for (j=graph->xadj[i];j<graph->xadj[i+1];j++) {
88751b0f6cfSStefano Zampini           k = graph->mirrors[graph->adjncy[j]];
88851b0f6cfSStefano Zampini           ierr = PetscMemcpy(&new_adjncy[new_xadj[i+1]],graph->mirrors_set[graph->adjncy[j]],k*sizeof(PetscInt));CHKERRQ(ierr);
88951b0f6cfSStefano Zampini           new_xadj[i+1] += k;
89051b0f6cfSStefano Zampini         }
89151b0f6cfSStefano Zampini         k = new_xadj[i+1]-new_xadj[i];
89251b0f6cfSStefano Zampini         ierr = PetscSortRemoveDupsInt(&k,&new_adjncy[new_xadj[i]]);CHKERRQ(ierr);
89351b0f6cfSStefano Zampini         new_xadj[i+1] = new_xadj[i]+k;
89451b0f6cfSStefano Zampini       }
89551b0f6cfSStefano Zampini       /* set new CSR into graph */
89651b0f6cfSStefano Zampini       ierr = PetscFree(graph->xadj);CHKERRQ(ierr);
89751b0f6cfSStefano Zampini       ierr = PetscFree(graph->adjncy);CHKERRQ(ierr);
89851b0f6cfSStefano Zampini       graph->xadj = new_xadj;
89951b0f6cfSStefano Zampini       graph->adjncy = new_adjncy;
90051b0f6cfSStefano Zampini     }
90108a5cf49SStefano Zampini   }
90251b0f6cfSStefano Zampini 
90317eb1463SStefano Zampini   /* mark special nodes (if any) -> each will become a single node equivalence class */
904674ae819SStefano Zampini   if (custom_primal_vertices) {
90563602bcaSStefano Zampini     ierr = ISGetLocalSize(custom_primal_vertices,&is_size);CHKERRQ(ierr);
906674ae819SStefano Zampini     ierr = ISGetIndices(custom_primal_vertices,(const PetscInt**)&is_indices);CHKERRQ(ierr);
907*7a0e7b2cSstefano_zampini     for (i=0,j=0;i<is_size;i++){
908*7a0e7b2cSstefano_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 */
909*7a0e7b2cSstefano_zampini         graph->special_dof[is_indices[i]] = PCBDDCGRAPH_SPECIAL_MARK-j;
9109b70a373SStefano Zampini         j++;
9119b70a373SStefano Zampini       }
9129b70a373SStefano Zampini     }
913*7a0e7b2cSstefano_zampini     ierr = ISRestoreIndices(custom_primal_vertices,(const PetscInt**)&is_indices);CHKERRQ(ierr);
91417eb1463SStefano Zampini   }
9159b70a373SStefano Zampini 
9164b2aedd3SStefano Zampini   /* mark interior nodes (if commsize > graph->commsizelimit) as touched and belonging to partition number 0 */
9174b2aedd3SStefano Zampini   if (commsize > graph->commsizelimit) {
918674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
919674ae819SStefano Zampini       if (!graph->count[i]) {
920df48933bSStefano Zampini         ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
921674ae819SStefano Zampini         graph->subset[i] = 0;
922674ae819SStefano Zampini       }
923674ae819SStefano Zampini     }
924b3cdcd63SStefano Zampini   }
925b3cdcd63SStefano Zampini 
926674ae819SStefano Zampini   /* init graph structure and compute default subsets */
927674ae819SStefano Zampini   nodes_touched = 0;
928674ae819SStefano Zampini   for (i=0;i<graph->nvtxs;i++) {
929df48933bSStefano Zampini     if (PetscBTLookup(graph->touched,i)) {
930674ae819SStefano Zampini       nodes_touched++;
931674ae819SStefano Zampini     }
932674ae819SStefano Zampini   }
933674ae819SStefano Zampini   i = 0;
934674ae819SStefano Zampini   graph->ncc = 0;
935674ae819SStefano Zampini   total_counts = 0;
9367babac9bSStefano Zampini 
9377babac9bSStefano Zampini   /* allocated space for queues */
9384b2aedd3SStefano Zampini   if (commsize == graph->commsizelimit) {
9397babac9bSStefano Zampini     ierr = PetscMalloc2(graph->nvtxs+1,&graph->cptr,graph->nvtxs,&graph->queue);CHKERRQ(ierr);
9407babac9bSStefano Zampini   } else {
9417babac9bSStefano Zampini     PetscInt nused = graph->nvtxs - nodes_touched;
9427babac9bSStefano Zampini     ierr = PetscMalloc2(nused+1,&graph->cptr,nused,&graph->queue);CHKERRQ(ierr);
9437babac9bSStefano Zampini   }
9447babac9bSStefano Zampini 
945674ae819SStefano Zampini   while (nodes_touched<graph->nvtxs) {
946674ae819SStefano Zampini     /*  find first untouched node in local ordering */
947b3cdcd63SStefano Zampini     while (PetscBTLookup(graph->touched,i)) i++;
948df48933bSStefano Zampini     ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
949674ae819SStefano Zampini     graph->subset[i] = graph->ncc+1;
950674ae819SStefano Zampini     graph->cptr[graph->ncc] = total_counts;
951674ae819SStefano Zampini     graph->queue[total_counts] = i;
952674ae819SStefano Zampini     total_counts++;
953674ae819SStefano Zampini     nodes_touched++;
954674ae819SStefano Zampini     /* now find all other nodes having the same set of sharing subdomains */
955674ae819SStefano Zampini     for (j=i+1;j<graph->nvtxs;j++) {
95674e413f5SStefano Zampini       /* check for same number of sharing subdomains, dof number and same special mark */
957df48933bSStefano 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]) {
958674ae819SStefano Zampini         /* check for same set of sharing subdomains */
959674ae819SStefano Zampini         same_set = PETSC_TRUE;
960674ae819SStefano Zampini         for (k=0;k<graph->count[j];k++){
961674ae819SStefano Zampini           if (graph->neighbours_set[i][k] != graph->neighbours_set[j][k]) {
962674ae819SStefano Zampini             same_set = PETSC_FALSE;
963674ae819SStefano Zampini           }
964674ae819SStefano Zampini         }
965674ae819SStefano Zampini         /* I found a friend of mine */
966674ae819SStefano Zampini         if (same_set) {
967df48933bSStefano Zampini           ierr = PetscBTSet(graph->touched,j);CHKERRQ(ierr);
968b3cdcd63SStefano Zampini           graph->subset[j] = graph->ncc+1;
969674ae819SStefano Zampini           nodes_touched++;
970674ae819SStefano Zampini           graph->queue[total_counts] = j;
971674ae819SStefano Zampini           total_counts++;
972674ae819SStefano Zampini         }
973674ae819SStefano Zampini       }
974674ae819SStefano Zampini     }
975674ae819SStefano Zampini     graph->ncc++;
976674ae819SStefano Zampini   }
977b3cdcd63SStefano 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 */
978674ae819SStefano Zampini   graph->n_subsets = graph->ncc;
979785e854fSJed Brown   ierr = PetscMalloc1(graph->n_subsets,&graph->subset_ncc);CHKERRQ(ierr);
980674ae819SStefano Zampini   for (i=0;i<graph->n_subsets;i++) {
981674ae819SStefano Zampini     graph->subset_ncc[i] = 1;
982674ae819SStefano Zampini   }
983674ae819SStefano Zampini   /* final pointer */
984674ae819SStefano Zampini   graph->cptr[graph->ncc] = total_counts;
985ec1c413dSStefano Zampini 
986b3cdcd63SStefano Zampini   /* For consistency reasons (among neighbours), I need to sort (by global ordering) each connected component */
987ec1c413dSStefano Zampini   /* Get a reference node (min index in global ordering) for each subset for tagging messages */
988dc456d91SStefano Zampini   ierr = PetscMalloc1(graph->ncc,&graph->subset_ref_node);CHKERRQ(ierr);
989dc456d91SStefano Zampini   ierr = PetscMalloc1(graph->cptr[graph->ncc],&queue_global);CHKERRQ(ierr);
9903a5b03d1SStefano Zampini   ierr = ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_global);CHKERRQ(ierr);
991ec1c413dSStefano Zampini   for (j=0;j<graph->ncc;j++) {
992ec1c413dSStefano Zampini     ierr = PetscSortIntWithArray(graph->cptr[j+1]-graph->cptr[j],&queue_global[graph->cptr[j]],&graph->queue[graph->cptr[j]]);CHKERRQ(ierr);
993dc456d91SStefano Zampini     graph->subset_ref_node[j] = graph->queue[graph->cptr[j]];
994f0321c90SStefano Zampini   }
995dc456d91SStefano Zampini   ierr = PetscFree(queue_global);CHKERRQ(ierr);
996acc113dbSStefano Zampini   graph->queue_sorted = PETSC_TRUE;
99748bebe3cSStefano Zampini 
998b3cdcd63SStefano Zampini   /* save information on subsets (needed when analyzing the connected components) */
9992f566a09SStefano Zampini   if (graph->ncc) {
1000b3cdcd63SStefano Zampini     ierr = PetscMalloc2(graph->ncc,&graph->subset_size,graph->ncc,&graph->subset_idxs);CHKERRQ(ierr);
1001b3cdcd63SStefano Zampini     ierr = PetscMalloc1(graph->cptr[graph->ncc],&graph->subset_idxs[0]);CHKERRQ(ierr);
1002b3cdcd63SStefano Zampini     ierr = PetscMemzero(graph->subset_idxs[0],graph->cptr[graph->ncc]*sizeof(PetscInt));CHKERRQ(ierr);
1003ec1c413dSStefano Zampini     for (j=1;j<graph->ncc;j++) {
1004b3cdcd63SStefano Zampini       graph->subset_size[j-1] = graph->cptr[j] - graph->cptr[j-1];
1005b3cdcd63SStefano Zampini       graph->subset_idxs[j] = graph->subset_idxs[j-1] + graph->subset_size[j-1];
1006ec1c413dSStefano Zampini     }
1007b3cdcd63SStefano Zampini     graph->subset_size[graph->ncc-1] = graph->cptr[graph->ncc] - graph->cptr[graph->ncc-1];
1008b3cdcd63SStefano Zampini     ierr = PetscMemcpy(graph->subset_idxs[0],graph->queue,graph->cptr[graph->ncc]*sizeof(PetscInt));CHKERRQ(ierr);
1009ec1c413dSStefano Zampini   }
1010b3cdcd63SStefano Zampini 
1011f0321c90SStefano Zampini   /* renumber reference nodes */
1012dc456d91SStefano Zampini   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)(graph->l2gmap)),graph->ncc,graph->subset_ref_node,PETSC_COPY_VALUES,&subset_n);CHKERRQ(ierr);
1013dc456d91SStefano Zampini   ierr = ISLocalToGlobalMappingApplyIS(graph->l2gmap,subset_n,&subset);CHKERRQ(ierr);
1014dc456d91SStefano Zampini   ierr = ISDestroy(&subset_n);CHKERRQ(ierr);
10156583bcc1SStefano Zampini   ierr = ISRenumber(subset,NULL,NULL,&subset_n);CHKERRQ(ierr);
1016dc456d91SStefano Zampini   ierr = ISDestroy(&subset);CHKERRQ(ierr);
1017dc456d91SStefano Zampini   ierr = ISGetLocalSize(subset_n,&k);CHKERRQ(ierr);
10186c4ed002SBarry Smith   if (k != graph->ncc) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid size of new subset! %D != %D",k,graph->ncc);
1019dc456d91SStefano Zampini   ierr = ISGetIndices(subset_n,&is_indices);CHKERRQ(ierr);
1020dc456d91SStefano Zampini   ierr = PetscMemcpy(graph->subset_ref_node,is_indices,graph->ncc*sizeof(PetscInt));CHKERRQ(ierr);
1021dc456d91SStefano Zampini   ierr = ISRestoreIndices(subset_n,&is_indices);CHKERRQ(ierr);
1022dc456d91SStefano Zampini   ierr = ISDestroy(&subset_n);CHKERRQ(ierr);
1023ec1c413dSStefano Zampini 
10243ebc787cSStefano Zampini   /* Determine if we are in 2D or 3D */
10253ebc787cSStefano Zampini   twodim = PETSC_TRUE;
10263ebc787cSStefano Zampini   for (i=0;i<graph->ncc;i++) {
10273ebc787cSStefano Zampini     PetscInt repdof = graph->queue[graph->cptr[i]];
102848bebe3cSStefano Zampini     if (graph->count[repdof] > 1) {
10293ebc787cSStefano Zampini       twodim = PETSC_FALSE;
10303ebc787cSStefano Zampini       break;
10313ebc787cSStefano Zampini     }
10323ebc787cSStefano Zampini   }
10333ebc787cSStefano Zampini   ierr = MPIU_Allreduce(&twodim,&graph->twodim,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)graph->l2gmap));CHKERRQ(ierr);
10343ebc787cSStefano Zampini 
1035ec1c413dSStefano Zampini   /* free workspace */
1036b3cdcd63SStefano Zampini   graph->setupcalled = PETSC_TRUE;
1037674ae819SStefano Zampini   PetscFunctionReturn(0);
1038674ae819SStefano Zampini }
1039674ae819SStefano Zampini 
1040674ae819SStefano Zampini #undef __FUNCT__
1041674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphResetCSR"
1042674ae819SStefano Zampini PetscErrorCode PCBDDCGraphResetCSR(PCBDDCGraph graph)
1043674ae819SStefano Zampini {
1044674ae819SStefano Zampini   PetscErrorCode ierr;
1045674ae819SStefano Zampini 
1046674ae819SStefano Zampini   PetscFunctionBegin;
1047a1dbd327SStefano Zampini   if (graph->freecsr) {
1048674ae819SStefano Zampini     ierr = PetscFree(graph->xadj);CHKERRQ(ierr);
1049674ae819SStefano Zampini     ierr = PetscFree(graph->adjncy);CHKERRQ(ierr);
1050a1dbd327SStefano Zampini   } else {
1051a1dbd327SStefano Zampini     graph->xadj = NULL;
1052a1dbd327SStefano Zampini     graph->adjncy = NULL;
1053a1dbd327SStefano Zampini   }
1054c8272957SStefano Zampini   graph->freecsr = PETSC_FALSE;
1055575ad6abSStefano Zampini   graph->nvtxs_csr = 0;
1056674ae819SStefano Zampini   PetscFunctionReturn(0);
1057674ae819SStefano Zampini }
1058674ae819SStefano Zampini 
1059674ae819SStefano Zampini #undef __FUNCT__
1060674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphReset"
1061674ae819SStefano Zampini PetscErrorCode PCBDDCGraphReset(PCBDDCGraph graph)
1062674ae819SStefano Zampini {
1063674ae819SStefano Zampini   PetscErrorCode ierr;
1064674ae819SStefano Zampini 
1065674ae819SStefano Zampini   PetscFunctionBegin;
1066674ae819SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&graph->l2gmap);CHKERRQ(ierr);
1067674ae819SStefano Zampini   ierr = PetscFree(graph->subset_ncc);CHKERRQ(ierr);
10683a5b03d1SStefano Zampini   ierr = PetscFree(graph->subset_ref_node);CHKERRQ(ierr);
1069674ae819SStefano Zampini   if (graph->nvtxs) {
1070674ae819SStefano Zampini     ierr = PetscFree(graph->neighbours_set[0]);CHKERRQ(ierr);
1071674ae819SStefano Zampini   }
1072df48933bSStefano Zampini   ierr = PetscBTDestroy(&graph->touched);CHKERRQ(ierr);
10737babac9bSStefano Zampini   ierr = PetscFree5(graph->count,
1074674ae819SStefano Zampini                     graph->neighbours_set,
1075674ae819SStefano Zampini                     graph->subset,
1076674ae819SStefano Zampini                     graph->which_dof,
1077df48933bSStefano Zampini                     graph->special_dof);CHKERRQ(ierr);
10787babac9bSStefano Zampini   ierr = PetscFree2(graph->cptr,graph->queue);CHKERRQ(ierr);
107951b0f6cfSStefano Zampini   if (graph->mirrors) {
108051b0f6cfSStefano Zampini     ierr = PetscFree(graph->mirrors_set[0]);CHKERRQ(ierr);
108151b0f6cfSStefano Zampini   }
108251b0f6cfSStefano Zampini   ierr = PetscFree2(graph->mirrors,graph->mirrors_set);CHKERRQ(ierr);
1083b3cdcd63SStefano Zampini   if (graph->subset_idxs) {
1084b3cdcd63SStefano Zampini     ierr = PetscFree(graph->subset_idxs[0]);CHKERRQ(ierr);
1085ec1c413dSStefano Zampini   }
1086b3cdcd63SStefano Zampini   ierr = PetscFree2(graph->subset_size,graph->subset_idxs);CHKERRQ(ierr);
1087a07ea27aSStefano Zampini   ierr = ISDestroy(&graph->dirdofs);CHKERRQ(ierr);
1088d62866d3SStefano Zampini   ierr = ISDestroy(&graph->dirdofsB);CHKERRQ(ierr);
10894f1b2e48SStefano Zampini   if (graph->n_local_subs) {
10904f1b2e48SStefano Zampini     ierr = PetscFree(graph->local_subs);CHKERRQ(ierr);
10914f1b2e48SStefano Zampini   }
1092a07ea27aSStefano Zampini   graph->has_dirichlet       = PETSC_FALSE;
1093674ae819SStefano Zampini   graph->nvtxs               = 0;
10947babac9bSStefano Zampini   graph->nvtxs_global        = 0;
1095674ae819SStefano Zampini   graph->n_subsets           = 0;
1096674ae819SStefano Zampini   graph->custom_minimal_size = 1;
10974f1b2e48SStefano Zampini   graph->n_local_subs        = 0;
1098be12c134Sstefano_zampini   graph->maxcount            = PETSC_MAX_INT;
1099fb597685SStefano Zampini   graph->setupcalled         = PETSC_FALSE;
1100674ae819SStefano Zampini   PetscFunctionReturn(0);
1101674ae819SStefano Zampini }
1102674ae819SStefano Zampini 
1103674ae819SStefano Zampini #undef __FUNCT__
1104674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphInit"
1105be12c134Sstefano_zampini PetscErrorCode PCBDDCGraphInit(PCBDDCGraph graph, ISLocalToGlobalMapping l2gmap, PetscInt N, PetscInt maxcount)
1106674ae819SStefano Zampini {
1107674ae819SStefano Zampini   PetscInt       n;
1108674ae819SStefano Zampini   PetscErrorCode ierr;
1109674ae819SStefano Zampini 
1110674ae819SStefano Zampini   PetscFunctionBegin;
1111674ae819SStefano Zampini   PetscValidPointer(graph,1);
1112674ae819SStefano Zampini   PetscValidHeaderSpecific(l2gmap,IS_LTOGM_CLASSID,2);
11137babac9bSStefano Zampini   PetscValidLogicalCollectiveInt(l2gmap,N,3);
1114be12c134Sstefano_zampini   PetscValidLogicalCollectiveInt(l2gmap,maxcount,4);
11158e61c736SStefano Zampini   /* raise an error if already allocated */
11166c4ed002SBarry Smith   if (graph->nvtxs_global) SETERRQ(PetscObjectComm((PetscObject)l2gmap),PETSC_ERR_PLIB,"BDDCGraph already initialized");
1117674ae819SStefano Zampini   /* set number of vertices */
1118674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)l2gmap);CHKERRQ(ierr);
1119674ae819SStefano Zampini   graph->l2gmap = l2gmap;
1120674ae819SStefano Zampini   ierr = ISLocalToGlobalMappingGetSize(l2gmap,&n);CHKERRQ(ierr);
1121674ae819SStefano Zampini   graph->nvtxs = n;
11227fb0e2dbSStefano Zampini   graph->nvtxs_global = N;
1123674ae819SStefano Zampini   /* allocate used space */
1124df48933bSStefano Zampini   ierr = PetscBTCreate(graph->nvtxs,&graph->touched);CHKERRQ(ierr);
11257babac9bSStefano Zampini   ierr = PetscMalloc5(graph->nvtxs,&graph->count,
11268e5aaad5SJed Brown                       graph->nvtxs,&graph->neighbours_set,
11278e5aaad5SJed Brown                       graph->nvtxs,&graph->subset,
11288e5aaad5SJed Brown                       graph->nvtxs,&graph->which_dof,
11298e5aaad5SJed Brown                       graph->nvtxs,&graph->special_dof);CHKERRQ(ierr);
1130674ae819SStefano Zampini   /* zeroes memory */
1131674ae819SStefano Zampini   ierr = PetscMemzero(graph->count,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
1132674ae819SStefano Zampini   ierr = PetscMemzero(graph->subset,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
113363602bcaSStefano Zampini   /* use -1 as a default value for which_dof array */
113463602bcaSStefano Zampini   for (n=0;n<graph->nvtxs;n++) graph->which_dof[n] = -1;
113574e413f5SStefano Zampini   ierr = PetscMemzero(graph->special_dof,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
1136674ae819SStefano Zampini   /* zeroes first pointer to neighbour set */
1137674ae819SStefano Zampini   if (graph->nvtxs) {
1138674ae819SStefano Zampini     graph->neighbours_set[0] = 0;
1139674ae819SStefano Zampini   }
1140674ae819SStefano Zampini   /* zeroes workspace for values of ncc */
1141674ae819SStefano Zampini   graph->subset_ncc = 0;
11423a5b03d1SStefano Zampini   graph->subset_ref_node = 0;
1143be12c134Sstefano_zampini   /* maxcount for cc */
1144be12c134Sstefano_zampini   graph->maxcount = maxcount;
1145674ae819SStefano Zampini   PetscFunctionReturn(0);
1146674ae819SStefano Zampini }
1147674ae819SStefano Zampini 
1148674ae819SStefano Zampini #undef __FUNCT__
1149674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphDestroy"
1150674ae819SStefano Zampini PetscErrorCode PCBDDCGraphDestroy(PCBDDCGraph* graph)
1151674ae819SStefano Zampini {
1152674ae819SStefano Zampini   PetscErrorCode ierr;
1153674ae819SStefano Zampini 
1154674ae819SStefano Zampini   PetscFunctionBegin;
1155f5fcb376SStefano Zampini   ierr = PCBDDCGraphResetCSR(*graph);CHKERRQ(ierr);
1156674ae819SStefano Zampini   ierr = PCBDDCGraphReset(*graph);CHKERRQ(ierr);
1157674ae819SStefano Zampini   ierr = PetscFree(*graph);CHKERRQ(ierr);
1158674ae819SStefano Zampini   PetscFunctionReturn(0);
1159674ae819SStefano Zampini }
1160674ae819SStefano Zampini 
1161674ae819SStefano Zampini #undef __FUNCT__
1162674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphCreate"
1163674ae819SStefano Zampini PetscErrorCode PCBDDCGraphCreate(PCBDDCGraph *graph)
1164674ae819SStefano Zampini {
1165674ae819SStefano Zampini   PCBDDCGraph    new_graph;
1166674ae819SStefano Zampini   PetscErrorCode ierr;
1167674ae819SStefano Zampini 
1168674ae819SStefano Zampini   PetscFunctionBegin;
1169854ce69bSBarry Smith   ierr = PetscNew(&new_graph);CHKERRQ(ierr);
1170674ae819SStefano Zampini   new_graph->custom_minimal_size = 1;
11714b2aedd3SStefano Zampini   new_graph->commsizelimit = 1;
1172674ae819SStefano Zampini   *graph = new_graph;
1173674ae819SStefano Zampini   PetscFunctionReturn(0);
1174674ae819SStefano Zampini }
1175