xref: /petsc/src/ksp/pc/impls/bddc/bddcgraph.c (revision be12c13469a3bc9ab17b9f62f1726d40dec0aeb3)
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);
81*be12c134Sstefano_zampini   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Max count %d\n",graph->maxcount);CHKERRQ(ierr);
82d4890cceSStefano Zampini   if (verbosity_level > 2) {
83674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
84674ae819SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"%d:\n",i);CHKERRQ(ierr);
85674ae819SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   which_dof: %d\n",graph->which_dof[i]);CHKERRQ(ierr);
8674e413f5SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   special_dof: %d\n",graph->special_dof[i]);CHKERRQ(ierr);
87674ae819SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   neighbours: %d\n",graph->count[i]);CHKERRQ(ierr);
882b510759SStefano Zampini       ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
89674ae819SStefano Zampini       if (graph->count[i]) {
90674ae819SStefano Zampini         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"     set of neighbours:");CHKERRQ(ierr);
91674ae819SStefano Zampini         for (j=0;j<graph->count[i];j++) {
92674ae819SStefano Zampini           ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d",graph->neighbours_set[i][j]);CHKERRQ(ierr);
93674ae819SStefano Zampini         }
94674ae819SStefano Zampini         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
95674ae819SStefano Zampini       }
962b510759SStefano Zampini       ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr);
972b510759SStefano Zampini       ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
9851b0f6cfSStefano Zampini       if (graph->mirrors) {
9951b0f6cfSStefano Zampini         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   mirrors: %d\n",graph->mirrors[i]);CHKERRQ(ierr);
10051b0f6cfSStefano Zampini         if (graph->mirrors[i]) {
1012b510759SStefano Zampini           ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
10251b0f6cfSStefano Zampini           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"     set of mirrors:");CHKERRQ(ierr);
10351b0f6cfSStefano Zampini           for (j=0;j<graph->mirrors[i];j++) {
10451b0f6cfSStefano Zampini             ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d",graph->mirrors_set[i][j]);CHKERRQ(ierr);
10551b0f6cfSStefano Zampini           }
10651b0f6cfSStefano Zampini           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
1072b510759SStefano Zampini           ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr);
1082b510759SStefano Zampini           ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
10951b0f6cfSStefano Zampini         }
11051b0f6cfSStefano Zampini       }
111d4890cceSStefano Zampini       if (verbosity_level > 3) {
112674ae819SStefano Zampini         if (graph->xadj && graph->adjncy) {
113674ae819SStefano Zampini           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   local adj list:");CHKERRQ(ierr);
1142b510759SStefano Zampini           ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
115674ae819SStefano Zampini           for (j=graph->xadj[i];j<graph->xadj[i+1];j++) {
116674ae819SStefano Zampini             ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d",graph->adjncy[j]);CHKERRQ(ierr);
117674ae819SStefano Zampini           }
118674ae819SStefano Zampini           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
1192b510759SStefano Zampini           ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr);
1202b510759SStefano Zampini           ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
121ec1c413dSStefano Zampini         } else {
122ec1c413dSStefano Zampini           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   no adj info\n");CHKERRQ(ierr);
123674ae819SStefano Zampini         }
124e49050b4SStefano Zampini       }
125531609faSStefano Zampini       if (graph->n_local_subs) {
126531609faSStefano Zampini         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   local sub id: %d\n",graph->local_subs[i]);CHKERRQ(ierr);
127531609faSStefano Zampini       }
128674ae819SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   interface subset id: %d\n",graph->subset[i]);CHKERRQ(ierr);
129674ae819SStefano Zampini       if (graph->subset[i] && graph->subset_ncc) {
130674ae819SStefano Zampini         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   ncc for subset: %d\n",graph->subset_ncc[graph->subset[i]-1]);CHKERRQ(ierr);
131674ae819SStefano Zampini       }
132674ae819SStefano Zampini     }
133e49050b4SStefano Zampini   }
134674ae819SStefano Zampini   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Total number of connected components %d\n",graph->ncc);CHKERRQ(ierr);
135785e854fSJed Brown   ierr = PetscMalloc1(graph->cptr[graph->ncc],&queue_in_global_numbering);CHKERRQ(ierr);
13693fb5fd3SStefano Zampini   ierr = ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_in_global_numbering);CHKERRQ(ierr);
137674ae819SStefano Zampini   for (i=0;i<graph->ncc;i++) {
1381ce3d396SStefano Zampini     PetscInt node_num=graph->queue[graph->cptr[i]];
1391ce3d396SStefano Zampini     PetscBool printcc = PETSC_FALSE;
140674ae819SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  %d (neighs:",i);CHKERRQ(ierr);
1412b510759SStefano Zampini     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
142674ae819SStefano Zampini     for (j=0;j<graph->count[node_num];j++) {
143674ae819SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d",graph->neighbours_set[node_num][j]);CHKERRQ(ierr);
144674ae819SStefano Zampini     }
145d4890cceSStefano Zampini     if (verbosity_level > 1) {
146674ae819SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"):");CHKERRQ(ierr);
147d4890cceSStefano Zampini       if (graph->twodim || graph->count[node_num] > 1 || (graph->count[node_num] == 1 && graph->special_dof[node_num] == PCBDDCGRAPH_NEUMANN_MARK)) {
1481ce3d396SStefano Zampini         printcc = PETSC_TRUE;
1491ce3d396SStefano Zampini       }
1501ce3d396SStefano Zampini       if (printcc) {
151674ae819SStefano Zampini         for (j=graph->cptr[i];j<graph->cptr[i+1];j++) {
15293fb5fd3SStefano Zampini           ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d (%d)",graph->queue[j],queue_in_global_numbering[j]);CHKERRQ(ierr);
153674ae819SStefano Zampini         }
1541ce3d396SStefano Zampini       }
155d4890cceSStefano Zampini     } else {
156d4890cceSStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,")");CHKERRQ(ierr);
157d4890cceSStefano Zampini     }
158674ae819SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
1592b510759SStefano Zampini     ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr);
1602b510759SStefano Zampini     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
161674ae819SStefano Zampini   }
16293fb5fd3SStefano Zampini   ierr = PetscFree(queue_in_global_numbering);CHKERRQ(ierr);
163674ae819SStefano Zampini   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
164674ae819SStefano Zampini   PetscFunctionReturn(0);
165674ae819SStefano Zampini }
166674ae819SStefano Zampini 
167674ae819SStefano Zampini #undef __FUNCT__
168c8272957SStefano Zampini #define __FUNCT__ "PCBDDCGraphRestoreCandidatesIS"
169c8272957SStefano Zampini PetscErrorCode PCBDDCGraphRestoreCandidatesIS(PCBDDCGraph graph, PetscInt *n_faces, IS *FacesIS[], PetscInt *n_edges, IS *EdgesIS[], IS *VerticesIS)
170c8272957SStefano Zampini {
171c8272957SStefano Zampini   PetscInt       i;
172c8272957SStefano Zampini   PetscErrorCode ierr;
173c8272957SStefano Zampini 
174c8272957SStefano Zampini   PetscFunctionBegin;
175c8272957SStefano Zampini   if (n_faces) {
176c8272957SStefano Zampini     if (FacesIS) {
177c8272957SStefano Zampini       for (i=0;i<*n_faces;i++) {
178c8272957SStefano Zampini         ierr = ISDestroy(&((*FacesIS)[i]));CHKERRQ(ierr);
179c8272957SStefano Zampini       }
180c8272957SStefano Zampini       ierr = PetscFree(*FacesIS);CHKERRQ(ierr);
181c8272957SStefano Zampini     }
182c8272957SStefano Zampini     *n_faces = 0;
183c8272957SStefano Zampini   }
184c8272957SStefano Zampini   if (n_edges) {
185c8272957SStefano Zampini     if (EdgesIS) {
186c8272957SStefano Zampini       for (i=0;i<*n_edges;i++) {
187c8272957SStefano Zampini         ierr = ISDestroy(&((*EdgesIS)[i]));CHKERRQ(ierr);
188c8272957SStefano Zampini       }
189c8272957SStefano Zampini       ierr = PetscFree(*EdgesIS);CHKERRQ(ierr);
190c8272957SStefano Zampini     }
191c8272957SStefano Zampini     *n_edges = 0;
192c8272957SStefano Zampini   }
193c8272957SStefano Zampini   if (VerticesIS) {
194c8272957SStefano Zampini     ierr = ISDestroy(VerticesIS);CHKERRQ(ierr);
195c8272957SStefano Zampini   }
196c8272957SStefano Zampini   PetscFunctionReturn(0);
197c8272957SStefano Zampini }
198c8272957SStefano Zampini 
199c8272957SStefano Zampini #undef __FUNCT__
200674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphGetCandidatesIS"
201a873d5d0SStefano Zampini PetscErrorCode PCBDDCGraphGetCandidatesIS(PCBDDCGraph graph, PetscInt *n_faces, IS *FacesIS[], PetscInt *n_edges, IS *EdgesIS[], IS *VerticesIS)
202674ae819SStefano Zampini {
203674ae819SStefano Zampini   IS             *ISForFaces,*ISForEdges,ISForVertices;
204*be12c134Sstefano_zampini   PetscInt       i,nfc,nec,nvc,*idx,*mark;
205674ae819SStefano Zampini   PetscErrorCode ierr;
206674ae819SStefano Zampini 
207674ae819SStefano Zampini   PetscFunctionBegin;
208*be12c134Sstefano_zampini   ierr = PetscCalloc1(graph->ncc,&mark);CHKERRQ(ierr);
209674ae819SStefano Zampini   /* loop on ccs to evalute number of faces, edges and vertices */
210674ae819SStefano Zampini   nfc = 0;
211674ae819SStefano Zampini   nec = 0;
212674ae819SStefano Zampini   nvc = 0;
213674ae819SStefano Zampini   for (i=0;i<graph->ncc;i++) {
2141e1f2d93SStefano Zampini     PetscInt repdof = graph->queue[graph->cptr[i]];
215*be12c134Sstefano_zampini     if (graph->cptr[i+1]-graph->cptr[i] > graph->custom_minimal_size && graph->count[repdof] < graph->maxcount) {
216*be12c134Sstefano_zampini       if (!graph->twodim && graph->count[repdof] == 1 && graph->special_dof[repdof] != PCBDDCGRAPH_NEUMANN_MARK) {
217674ae819SStefano Zampini         nfc++;
218*be12c134Sstefano_zampini         mark[i] = 2;
219*be12c134Sstefano_zampini       } else {
220674ae819SStefano Zampini         nec++;
221*be12c134Sstefano_zampini         mark[i] = 1;
222674ae819SStefano Zampini       }
223674ae819SStefano Zampini     } else {
224674ae819SStefano Zampini       nvc += graph->cptr[i+1]-graph->cptr[i];
225674ae819SStefano Zampini     }
226674ae819SStefano Zampini   }
227adfc4fb2SStefano Zampini 
228674ae819SStefano Zampini   /* allocate IS arrays for faces, edges. Vertices need a single index set. */
229cf5a6209SStefano Zampini   if (FacesIS) {
230785e854fSJed Brown     ierr = PetscMalloc1(nfc,&ISForFaces);CHKERRQ(ierr);
231cf5a6209SStefano Zampini   }
232cf5a6209SStefano Zampini   if (EdgesIS) {
233785e854fSJed Brown     ierr = PetscMalloc1(nec,&ISForEdges);CHKERRQ(ierr);
234cf5a6209SStefano Zampini   }
235cf5a6209SStefano Zampini   if (VerticesIS) {
236785e854fSJed Brown     ierr = PetscMalloc1(nvc,&idx);CHKERRQ(ierr);
237cf5a6209SStefano Zampini   }
238cf5a6209SStefano Zampini 
239674ae819SStefano Zampini   /* loop on ccs to compute index sets for faces and edges */
240acc113dbSStefano Zampini   if (!graph->queue_sorted) {
241acc113dbSStefano Zampini     PetscInt *queue_global;
242acc113dbSStefano Zampini 
243acc113dbSStefano Zampini     ierr = PetscMalloc1(graph->cptr[graph->ncc],&queue_global);CHKERRQ(ierr);
244acc113dbSStefano Zampini     ierr = ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_global);CHKERRQ(ierr);
245acc113dbSStefano Zampini     for (i=0;i<graph->ncc;i++) {
246acc113dbSStefano Zampini       ierr = PetscSortIntWithArray(graph->cptr[i+1]-graph->cptr[i],&queue_global[graph->cptr[i]],&graph->queue[graph->cptr[i]]);CHKERRQ(ierr);
247acc113dbSStefano Zampini     }
248acc113dbSStefano Zampini     ierr = PetscFree(queue_global);CHKERRQ(ierr);
249acc113dbSStefano Zampini     graph->queue_sorted = PETSC_TRUE;
250acc113dbSStefano Zampini   }
251674ae819SStefano Zampini   nfc = 0;
252674ae819SStefano Zampini   nec = 0;
253674ae819SStefano Zampini   for (i=0;i<graph->ncc;i++) {
2541e1f2d93SStefano Zampini     PetscInt repdof = graph->queue[graph->cptr[i]];
255*be12c134Sstefano_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++;
260*be12c134Sstefano_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   }
267*be12c134Sstefano_zampini 
268674ae819SStefano Zampini   /* index set for vertices */
269cf5a6209SStefano Zampini   if (VerticesIS) {
270b8ffe317SStefano Zampini     nvc = 0;
271674ae819SStefano Zampini     for (i=0;i<graph->ncc;i++) {
272*be12c134Sstefano_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   }
285*be12c134Sstefano_zampini   ierr = PetscFree(mark);CHKERRQ(ierr);
286*be12c134Sstefano_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;
3048e4af1ccSStefano Zampini   PetscBool      twodim;
305674ae819SStefano Zampini   PetscErrorCode ierr;
306674ae819SStefano Zampini 
307674ae819SStefano Zampini   PetscFunctionBegin;
308674ae819SStefano Zampini   /* compute connected components locally */
309674ae819SStefano Zampini   ierr = PetscObjectGetComm((PetscObject)(graph->l2gmap),&interface_comm);CHKERRQ(ierr);
310674ae819SStefano Zampini   ierr = PCBDDCGraphComputeConnectedComponentsLocal(graph);CHKERRQ(ierr);
311674ae819SStefano Zampini   /* check consistency of connected components among neighbouring subdomains -> it adapt them in case it is needed */
3124a060362SStefano Zampini   ierr = MPI_Comm_size(interface_comm,&size);CHKERRQ(ierr);
3134a060362SStefano Zampini   adapt_interface_reduced = PETSC_FALSE;
3144a060362SStefano Zampini   if (size > 1) {
3154a060362SStefano Zampini     PetscInt i;
3164a060362SStefano Zampini     PetscBool adapt_interface = PETSC_FALSE;
317674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
318674ae819SStefano Zampini       /* We are not sure that on a given subset of the local interface,
319674ae819SStefano Zampini          with two connected components, the latters be the same among sharing subdomains */
320674ae819SStefano Zampini       if (graph->subset_ncc[i] > 1) {
3214a060362SStefano Zampini         adapt_interface = PETSC_TRUE;
322674ae819SStefano Zampini         break;
323674ae819SStefano Zampini       }
324674ae819SStefano Zampini     }
325b2566f29SBarry Smith     ierr = MPIU_Allreduce(&adapt_interface,&adapt_interface_reduced,1,MPIU_BOOL,MPI_LOR,interface_comm);CHKERRQ(ierr);
3264a060362SStefano Zampini   }
327674ae819SStefano Zampini 
328674ae819SStefano Zampini   if (graph->n_subsets && adapt_interface_reduced) {
329ec1c413dSStefano Zampini     PetscBT     subset_cc_adapt;
330ec1c413dSStefano Zampini     MPI_Request *send_requests,*recv_requests;
331ec1c413dSStefano Zampini     PetscInt    *send_buffer,*recv_buffer;
332ec1c413dSStefano Zampini     PetscInt    sum_requests,start_of_recv,start_of_send;
333ec1c413dSStefano Zampini     PetscInt    *cum_recv_counts;
334ec1c413dSStefano Zampini     PetscInt    *labels;
335b3cdcd63SStefano Zampini     PetscInt    ncc,cum_queue,mss,mns,j,k,s;
3368875e3e6SStefano Zampini     PetscInt    **refine_buffer=NULL,*private_labels = NULL;
3375b1b9aeaSStefano Zampini 
338ec1c413dSStefano Zampini     ierr = PetscMalloc1(graph->nvtxs,&labels);CHKERRQ(ierr);
339ec1c413dSStefano Zampini     ierr = PetscMemzero(labels,graph->nvtxs*sizeof(*labels));CHKERRQ(ierr);
340b3cdcd63SStefano Zampini     for (i=0;i<graph->ncc;i++)
341b3cdcd63SStefano Zampini       for (j=graph->cptr[i];j<graph->cptr[i+1];j++)
342b3cdcd63SStefano Zampini         labels[graph->queue[j]] = i;
343b3cdcd63SStefano Zampini 
344674ae819SStefano Zampini     /* allocate some space */
345854ce69bSBarry Smith     ierr = PetscMalloc1(graph->n_subsets+1,&cum_recv_counts);CHKERRQ(ierr);
346674ae819SStefano Zampini     ierr = PetscMemzero(cum_recv_counts,(graph->n_subsets+1)*sizeof(*cum_recv_counts));CHKERRQ(ierr);
347b3cdcd63SStefano Zampini 
348674ae819SStefano Zampini     /* first count how many neighbours per connected component I will receive from */
349674ae819SStefano Zampini     cum_recv_counts[0] = 0;
350b3cdcd63SStefano 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]];
351ec1c413dSStefano Zampini     ierr = PetscMalloc1(cum_recv_counts[graph->n_subsets],&recv_buffer);CHKERRQ(ierr);
352ec1c413dSStefano Zampini     ierr = PetscMalloc2(cum_recv_counts[graph->n_subsets],&send_requests,cum_recv_counts[graph->n_subsets],&recv_requests);CHKERRQ(ierr);
353674ae819SStefano Zampini     for (i=0;i<cum_recv_counts[graph->n_subsets];i++) {
354674ae819SStefano Zampini       send_requests[i] = MPI_REQUEST_NULL;
355674ae819SStefano Zampini       recv_requests[i] = MPI_REQUEST_NULL;
356674ae819SStefano Zampini     }
357b3cdcd63SStefano Zampini 
358ec1c413dSStefano Zampini     /* exchange with my neighbours the number of my connected components on the subset of interface */
359674ae819SStefano Zampini     sum_requests = 0;
360674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
361ec1c413dSStefano Zampini       PetscMPIInt neigh,tag;
362b3cdcd63SStefano Zampini       PetscInt    count,*neighs;
363ec1c413dSStefano Zampini 
364b3cdcd63SStefano Zampini       count = graph->count[graph->subset_idxs[i][0]];
365b3cdcd63SStefano Zampini       neighs = graph->neighbours_set[graph->subset_idxs[i][0]];
366ec1c413dSStefano Zampini       ierr = PetscMPIIntCast(2*graph->subset_ref_node[i],&tag);CHKERRQ(ierr);
367b3cdcd63SStefano Zampini       for (k=0;k<count;k++) {
368b3cdcd63SStefano Zampini         ierr = PetscMPIIntCast(neighs[k],&neigh);CHKERRQ(ierr);
369674ae819SStefano Zampini         ierr = MPI_Isend(&graph->subset_ncc[i],1,MPIU_INT,neigh,tag,interface_comm,&send_requests[sum_requests]);CHKERRQ(ierr);
370ec1c413dSStefano Zampini         ierr = MPI_Irecv(&recv_buffer[sum_requests],1,MPIU_INT,neigh,tag,interface_comm,&recv_requests[sum_requests]);CHKERRQ(ierr);
371674ae819SStefano Zampini         sum_requests++;
372674ae819SStefano Zampini       }
373674ae819SStefano Zampini     }
374674ae819SStefano Zampini     ierr = MPI_Waitall(sum_requests,recv_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
375674ae819SStefano Zampini     ierr = MPI_Waitall(sum_requests,send_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
376b3cdcd63SStefano Zampini 
377b3cdcd63SStefano Zampini     /* determine the subsets I have to adapt (those having more than 1 cc) */
378ec1c413dSStefano Zampini     ierr = PetscBTCreate(graph->n_subsets,&subset_cc_adapt);CHKERRQ(ierr);
379ec1c413dSStefano Zampini     ierr = PetscBTMemzero(graph->n_subsets,subset_cc_adapt);CHKERRQ(ierr);
380674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
381b3cdcd63SStefano Zampini       if (graph->subset_ncc[i] > 1) {
382b3cdcd63SStefano Zampini         ierr = PetscBTSet(subset_cc_adapt,i);CHKERRQ(ierr);
383b3cdcd63SStefano Zampini         continue;
384b3cdcd63SStefano Zampini       }
385674ae819SStefano Zampini       for (j=cum_recv_counts[i];j<cum_recv_counts[i+1];j++){
386b3cdcd63SStefano Zampini          if (recv_buffer[j] > 1) {
387ec1c413dSStefano Zampini           ierr = PetscBTSet(subset_cc_adapt,i);
388674ae819SStefano Zampini           break;
389674ae819SStefano Zampini         }
390674ae819SStefano Zampini       }
391674ae819SStefano Zampini     }
392ec1c413dSStefano Zampini     ierr = PetscFree(recv_buffer);CHKERRQ(ierr);
393b3cdcd63SStefano Zampini 
394ec1c413dSStefano Zampini     /* determine send/recv buffers sizes */
395ec1c413dSStefano Zampini     j = 0;
396b3cdcd63SStefano Zampini     mss = 0;
397674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
398ec1c413dSStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i)) {
399b3cdcd63SStefano Zampini         j += graph->subset_size[i];
400b3cdcd63SStefano Zampini         mss = PetscMax(graph->subset_size[i],mss);
401674ae819SStefano Zampini       }
402674ae819SStefano Zampini     }
403ec1c413dSStefano Zampini     k = 0;
404b3cdcd63SStefano Zampini     mns = 0;
405674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
406ec1c413dSStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i)) {
407b3cdcd63SStefano Zampini         k += (cum_recv_counts[i+1]-cum_recv_counts[i])*graph->subset_size[i];
408b3cdcd63SStefano Zampini         mns = PetscMax(cum_recv_counts[i+1]-cum_recv_counts[i],mns);
409674ae819SStefano Zampini       }
410674ae819SStefano Zampini     }
411ec1c413dSStefano Zampini     ierr = PetscMalloc2(j,&send_buffer,k,&recv_buffer);CHKERRQ(ierr);
412ec1c413dSStefano Zampini 
413b3cdcd63SStefano Zampini     /* fill send buffer (order matters: subset_idxs ordered by global ordering) */
414ec1c413dSStefano Zampini     j = 0;
415b3cdcd63SStefano Zampini     for (i=0;i<graph->n_subsets;i++)
416b3cdcd63SStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i))
417b3cdcd63SStefano Zampini         for (k=0;k<graph->subset_size[i];k++)
418b3cdcd63SStefano Zampini           send_buffer[j++] = labels[graph->subset_idxs[i][k]];
419ec1c413dSStefano Zampini 
420674ae819SStefano Zampini     /* now exchange the data */
421674ae819SStefano Zampini     start_of_recv = 0;
422674ae819SStefano Zampini     start_of_send = 0;
423674ae819SStefano Zampini     sum_requests = 0;
424674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
425ec1c413dSStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i)) {
426ec1c413dSStefano Zampini         PetscMPIInt neigh,tag;
427b3cdcd63SStefano Zampini         PetscInt    size_of_send = graph->subset_size[i];
428ec1c413dSStefano Zampini 
429b3cdcd63SStefano Zampini         j = graph->subset_idxs[i][0];
430ec1c413dSStefano Zampini         ierr = PetscMPIIntCast(2*graph->subset_ref_node[i]+1,&tag);CHKERRQ(ierr);
431674ae819SStefano Zampini         for (k=0;k<graph->count[j];k++) {
432674ae819SStefano Zampini           ierr = PetscMPIIntCast(graph->neighbours_set[j][k],&neigh);CHKERRQ(ierr);
433ec1c413dSStefano Zampini           ierr = MPI_Isend(&send_buffer[start_of_send],size_of_send,MPIU_INT,neigh,tag,interface_comm,&send_requests[sum_requests]);CHKERRQ(ierr);
434ec1c413dSStefano Zampini           ierr = MPI_Irecv(&recv_buffer[start_of_recv],size_of_send,MPIU_INT,neigh,tag,interface_comm,&recv_requests[sum_requests]);CHKERRQ(ierr);
435ec1c413dSStefano Zampini           start_of_recv += size_of_send;
436674ae819SStefano Zampini           sum_requests++;
437674ae819SStefano Zampini         }
438674ae819SStefano Zampini         start_of_send += size_of_send;
439674ae819SStefano Zampini       }
440674ae819SStefano Zampini     }
441674ae819SStefano Zampini     ierr = MPI_Waitall(sum_requests,recv_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
442d52c4852SStefano Zampini 
443b3cdcd63SStefano Zampini     /* refine connected components */
444674ae819SStefano Zampini     start_of_recv = 0;
445b3cdcd63SStefano Zampini     /* allocate some temporary space */
446b3cdcd63SStefano Zampini     if (mss) {
447b3cdcd63SStefano Zampini       ierr = PetscMalloc1(mss,&refine_buffer);CHKERRQ(ierr);
448b3cdcd63SStefano Zampini       ierr = PetscMalloc2(mss*(mns+1),&refine_buffer[0],mss,&private_labels);CHKERRQ(ierr);
449b3cdcd63SStefano Zampini     }
450b3cdcd63SStefano Zampini     ncc = 0;
451b3cdcd63SStefano Zampini     cum_queue = 0;
452b3cdcd63SStefano Zampini     graph->cptr[0] = 0;
453674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
454ec1c413dSStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i)) {
455b3cdcd63SStefano Zampini         PetscInt subset_counter = 0;
456b3cdcd63SStefano Zampini         PetscInt sharingprocs = cum_recv_counts[i+1]-cum_recv_counts[i]+1; /* count myself */
457b3cdcd63SStefano Zampini         PetscInt buffer_size = graph->subset_size[i];
458ec1c413dSStefano Zampini 
459b3cdcd63SStefano Zampini         /* compute pointers */
460b3cdcd63SStefano Zampini         for (j=1;j<buffer_size;j++) refine_buffer[j] = refine_buffer[j-1] + sharingprocs;
461b3cdcd63SStefano Zampini         /* analyze contributions from subdomains that share the i-th subset
462b3cdcd63SStefano Zampini            The stricture of refine_buffer is suitable to find intersections of ccs among sharingprocs.
463b3cdcd63SStefano 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)
464b3cdcd63SStefano Zampini            sharing procs connected components:
465ec1c413dSStefano Zampini              neigh 0: [0 1 4], [2 3], labels [4,7]  (2 connected components)
466ec1c413dSStefano Zampini              neigh 1: [0 1], [2 3 4], labels [3 2]  (2 connected components)
467ec1c413dSStefano Zampini              neigh 2: [0 4], [1], [2 3], labels [1 5 6] (3 connected components)
468b3cdcd63SStefano Zampini            refine_buffer will be filled as:
469ec1c413dSStefano Zampini              [ 4, 3, 1;
470ec1c413dSStefano Zampini                4, 2, 1;
471ec1c413dSStefano Zampini                7, 2, 6;
472ec1c413dSStefano Zampini                4, 3, 5;
473ec1c413dSStefano Zampini                7, 2, 6; ];
474b3cdcd63SStefano Zampini            The connected components in local ordering are [0], [1], [2 3], [4] */
475b3cdcd63SStefano Zampini         /* fill temp_buffer */
476b3cdcd63SStefano Zampini         for (k=0;k<buffer_size;k++) refine_buffer[k][0] = labels[graph->subset_idxs[i][k]];
477b3cdcd63SStefano Zampini         for (j=0;j<sharingprocs-1;j++) {
478b3cdcd63SStefano Zampini           for (k=0;k<buffer_size;k++) refine_buffer[k][j+1] = recv_buffer[start_of_recv+k];
479b3cdcd63SStefano Zampini           start_of_recv += buffer_size;
480674ae819SStefano Zampini         }
481b3cdcd63SStefano Zampini         ierr = PetscMemzero(private_labels,buffer_size*sizeof(PetscInt));CHKERRQ(ierr);
482b3cdcd63SStefano Zampini         for (j=0;j<buffer_size;j++) {
483b3cdcd63SStefano Zampini           if (!private_labels[j]) { /* found a new cc  */
484ec1c413dSStefano Zampini             PetscBool same_set;
485ec1c413dSStefano Zampini 
486b3cdcd63SStefano Zampini             graph->cptr[ncc] = cum_queue;
487b3cdcd63SStefano Zampini             ncc++;
488b3cdcd63SStefano Zampini             subset_counter++;
489b3cdcd63SStefano Zampini             private_labels[j] = subset_counter;
490b3cdcd63SStefano Zampini             graph->queue[cum_queue++] = graph->subset_idxs[i][j];
491b3cdcd63SStefano Zampini             for (k=j+1;k<buffer_size;k++) { /* check for other nodes in new cc */
492674ae819SStefano Zampini               same_set = PETSC_TRUE;
493b3cdcd63SStefano Zampini               for (s=0;s<sharingprocs;s++) {
494b3cdcd63SStefano Zampini                 if (refine_buffer[j][s] != refine_buffer[k][s]) {
495674ae819SStefano Zampini                   same_set = PETSC_FALSE;
496674ae819SStefano Zampini                   break;
497674ae819SStefano Zampini                 }
498674ae819SStefano Zampini               }
499674ae819SStefano Zampini               if (same_set) {
500b3cdcd63SStefano Zampini                 private_labels[k] = subset_counter;
501b3cdcd63SStefano Zampini                 graph->queue[cum_queue++] = graph->subset_idxs[i][k];
502674ae819SStefano Zampini               }
503674ae819SStefano Zampini             }
504674ae819SStefano Zampini           }
505674ae819SStefano Zampini         }
506b3cdcd63SStefano Zampini         graph->cptr[ncc] = cum_queue;
507b3cdcd63SStefano Zampini         graph->subset_ncc[i] = subset_counter;
508b3cdcd63SStefano Zampini         graph->queue_sorted = PETSC_FALSE;
509b3cdcd63SStefano Zampini       } else { /* this subset does not need to be adapted */
510b3cdcd63SStefano Zampini         ierr = PetscMemcpy(graph->queue+cum_queue,graph->subset_idxs[i],graph->subset_size[i]*sizeof(PetscInt));CHKERRQ(ierr);
511b3cdcd63SStefano Zampini         ncc++;
512b3cdcd63SStefano Zampini         cum_queue += graph->subset_size[i];
513b3cdcd63SStefano Zampini         graph->cptr[ncc] = cum_queue;
514674ae819SStefano Zampini       }
515674ae819SStefano Zampini     }
516b3cdcd63SStefano Zampini     graph->cptr[ncc] = cum_queue;
517b3cdcd63SStefano Zampini     graph->ncc = ncc;
518b3cdcd63SStefano Zampini     if (mss) {
519b3cdcd63SStefano Zampini       ierr = PetscFree2(refine_buffer[0],private_labels);CHKERRQ(ierr);
520b3cdcd63SStefano Zampini       ierr = PetscFree(refine_buffer);CHKERRQ(ierr);
521b3cdcd63SStefano Zampini     }
522b3cdcd63SStefano Zampini     ierr = PetscFree(labels);CHKERRQ(ierr);
523b3cdcd63SStefano Zampini     ierr = MPI_Waitall(sum_requests,send_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
524b3cdcd63SStefano Zampini     ierr = PetscFree2(send_requests,recv_requests);CHKERRQ(ierr);
52545b2a5aaSStefano Zampini     ierr = PetscFree2(send_buffer,recv_buffer);CHKERRQ(ierr);
526674ae819SStefano Zampini     ierr = PetscFree(cum_recv_counts);CHKERRQ(ierr);
527ec1c413dSStefano Zampini     ierr = PetscBTDestroy(&subset_cc_adapt);CHKERRQ(ierr);
528674ae819SStefano Zampini   }
5298e4af1ccSStefano Zampini 
5308e4af1ccSStefano Zampini   /* Determine if we are in 2D or 3D */
5318e4af1ccSStefano Zampini   twodim  = PETSC_TRUE;
5328e4af1ccSStefano Zampini   for (i=0;i<graph->ncc;i++) {
5338e4af1ccSStefano Zampini     PetscInt repdof = graph->queue[graph->cptr[i]];
5348e4af1ccSStefano Zampini     if (graph->cptr[i+1]-graph->cptr[i] > graph->custom_minimal_size) {
5358e4af1ccSStefano Zampini       if (graph->count[repdof] > 1 || graph->special_dof[repdof] == PCBDDCGRAPH_NEUMANN_MARK) {
5368e4af1ccSStefano Zampini         twodim = PETSC_FALSE;
5378e4af1ccSStefano Zampini         break;
5388e4af1ccSStefano Zampini       }
5398e4af1ccSStefano Zampini     }
5408e4af1ccSStefano Zampini   }
541b2566f29SBarry Smith   ierr = MPIU_Allreduce(&twodim,&graph->twodim,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)graph->l2gmap));CHKERRQ(ierr);
542674ae819SStefano Zampini   PetscFunctionReturn(0);
543674ae819SStefano Zampini }
544674ae819SStefano Zampini 
545b3cdcd63SStefano Zampini 
546b3cdcd63SStefano Zampini #undef __FUNCT__
547b3cdcd63SStefano Zampini #define __FUNCT__ "PCBDDCGraphComputeCC_Private"
548b3cdcd63SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCGraphComputeCC_Private(PCBDDCGraph graph,PetscInt pid,PetscInt* queue_tip,PetscInt n_prev,PetscInt* n_added)
549b3cdcd63SStefano Zampini {
550b3cdcd63SStefano Zampini   PetscInt       i,j,n;
551b3cdcd63SStefano Zampini   PetscInt       *xadj = graph->xadj,*adjncy = graph->adjncy;
552b3cdcd63SStefano Zampini   PetscBT        touched = graph->touched;
553b3cdcd63SStefano Zampini   PetscBool      havecsr = (PetscBool)(xadj && adjncy);
554b3cdcd63SStefano Zampini   PetscBool      havesubs = (PetscBool)(!!graph->n_local_subs);
555b3cdcd63SStefano Zampini   PetscErrorCode ierr;
556b3cdcd63SStefano Zampini 
557b3cdcd63SStefano Zampini   PetscFunctionBegin;
558b3cdcd63SStefano Zampini   n = 0;
559b3cdcd63SStefano Zampini   if (havecsr && !havesubs) {
560b3cdcd63SStefano Zampini     for (i=-n_prev;i<0;i++) {
561b3cdcd63SStefano Zampini       PetscInt start_dof = queue_tip[i];
562b3cdcd63SStefano Zampini       for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
563b3cdcd63SStefano Zampini         PetscInt dof = adjncy[j];
564b3cdcd63SStefano Zampini         if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid) {
565b3cdcd63SStefano Zampini           ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
566b3cdcd63SStefano Zampini           queue_tip[n] = dof;
567b3cdcd63SStefano Zampini           n++;
568b3cdcd63SStefano Zampini         }
569b3cdcd63SStefano Zampini       }
570b3cdcd63SStefano Zampini     }
571b3cdcd63SStefano Zampini   } else if (havecsr && havesubs) {
572b3cdcd63SStefano Zampini     PetscInt sid = graph->local_subs[queue_tip[-n_prev]];
573b3cdcd63SStefano Zampini     for (i=-n_prev;i<0;i++) {
574b3cdcd63SStefano Zampini       PetscInt start_dof = queue_tip[i];
575b3cdcd63SStefano Zampini       for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
576b3cdcd63SStefano Zampini         PetscInt dof = adjncy[j];
577b3cdcd63SStefano Zampini         if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid && graph->local_subs[dof] == sid) {
578b3cdcd63SStefano Zampini           ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
579b3cdcd63SStefano Zampini           queue_tip[n] = dof;
580b3cdcd63SStefano Zampini           n++;
581b3cdcd63SStefano Zampini         }
582b3cdcd63SStefano Zampini       }
583b3cdcd63SStefano Zampini     }
584b3cdcd63SStefano Zampini   } else { /* sub info only */
585b3cdcd63SStefano Zampini     PetscInt sid = graph->local_subs[queue_tip[-n_prev]];
586b3cdcd63SStefano Zampini     for (j=0;j<graph->subset_size[pid-1];j++) { /* pid \in [1,graph->n_subsets] */
587b3cdcd63SStefano Zampini       PetscInt dof = graph->subset_idxs[pid-1][j];
588b3cdcd63SStefano Zampini       if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid && graph->local_subs[dof] == sid) {
589b3cdcd63SStefano Zampini         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
590b3cdcd63SStefano Zampini         queue_tip[n] = dof;
591b3cdcd63SStefano Zampini         n++;
592b3cdcd63SStefano Zampini       }
593b3cdcd63SStefano Zampini     }
594b3cdcd63SStefano Zampini   }
595b3cdcd63SStefano Zampini   *n_added = n;
596b3cdcd63SStefano Zampini   PetscFunctionReturn(0);
597b3cdcd63SStefano Zampini }
598b3cdcd63SStefano Zampini 
599674ae819SStefano Zampini #undef __FUNCT__
600674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphComputeConnectedComponentsLocal"
601674ae819SStefano Zampini PetscErrorCode PCBDDCGraphComputeConnectedComponentsLocal(PCBDDCGraph graph)
602674ae819SStefano Zampini {
603b3cdcd63SStefano Zampini   PetscInt       ncc,cum_queue,n;
604b3cdcd63SStefano Zampini   PetscMPIInt    commsize;
605674ae819SStefano Zampini   PetscErrorCode ierr;
606674ae819SStefano Zampini 
607674ae819SStefano Zampini   PetscFunctionBegin;
608b3cdcd63SStefano Zampini   if (!graph->setupcalled) SETERRQ(PetscObjectComm((PetscObject)graph->l2gmap),PETSC_ERR_ORDER,"PCBDDCGraphSetUp should be called first");
609b3cdcd63SStefano Zampini   /* quiet return if there isn't any local info */
6104f1b2e48SStefano Zampini   if ((!graph->xadj || !graph->adjncy) && !graph->n_local_subs) {
611674ae819SStefano Zampini     PetscFunctionReturn(0);
612674ae819SStefano Zampini   }
613674ae819SStefano Zampini 
614674ae819SStefano Zampini   /* reset any previous search of connected components */
615df48933bSStefano Zampini   ierr = PetscBTMemzero(graph->nvtxs,graph->touched);CHKERRQ(ierr);
616b3cdcd63SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)graph->l2gmap),&commsize);CHKERRQ(ierr);
6174b2aedd3SStefano Zampini   if (commsize > graph->commsizelimit) {
618b3cdcd63SStefano Zampini     PetscInt i;
619674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
6200cf82fd6SStefano Zampini       if (graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK || !graph->count[i]) {
621df48933bSStefano Zampini         ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
622674ae819SStefano Zampini       }
6234a060362SStefano Zampini     }
624b3cdcd63SStefano Zampini   }
625674ae819SStefano Zampini 
626674ae819SStefano Zampini   /* begin search for connected components */
627674ae819SStefano Zampini   cum_queue = 0;
628674ae819SStefano Zampini   ncc = 0;
629674ae819SStefano Zampini   for (n=0;n<graph->n_subsets;n++) {
630b3cdcd63SStefano Zampini     PetscInt pid = n+1;  /* partition labeled by 0 is discarded */
631b3cdcd63SStefano Zampini     PetscInt found = 0,prev = 0,first = 0,ncc_pid = 0;
632b3cdcd63SStefano Zampini     while (found != graph->subset_size[n]) {
633b3cdcd63SStefano Zampini       PetscInt added = 0;
634b3cdcd63SStefano Zampini       if (!prev) { /* search for new starting dof */
635b3cdcd63SStefano Zampini         while (PetscBTLookup(graph->touched,graph->subset_idxs[n][first])) first++;
636b3cdcd63SStefano Zampini         ierr = PetscBTSet(graph->touched,graph->subset_idxs[n][first]);CHKERRQ(ierr);
637b3cdcd63SStefano Zampini         graph->queue[cum_queue] = graph->subset_idxs[n][first];
638674ae819SStefano Zampini         graph->cptr[ncc] = cum_queue;
639b3cdcd63SStefano Zampini         prev = 1;
640b3cdcd63SStefano Zampini         cum_queue++;
641b3cdcd63SStefano Zampini         found++;
642674ae819SStefano Zampini         ncc_pid++;
643b3cdcd63SStefano Zampini         ncc++;
644674ae819SStefano Zampini       }
645b3cdcd63SStefano Zampini       ierr = PCBDDCGraphComputeCC_Private(graph,pid,graph->queue + cum_queue,prev,&added);CHKERRQ(ierr);
646b3cdcd63SStefano Zampini       if (!added) {
647674ae819SStefano Zampini         graph->subset_ncc[n] = ncc_pid;
648b3cdcd63SStefano Zampini         graph->cptr[ncc] = cum_queue;
649b3cdcd63SStefano Zampini       }
650b3cdcd63SStefano Zampini       prev = added;
651b3cdcd63SStefano Zampini       found += added;
652b3cdcd63SStefano Zampini       cum_queue += added;
653b3cdcd63SStefano Zampini       if (added && found == graph->subset_size[n]) {
654b3cdcd63SStefano Zampini         graph->subset_ncc[n] = ncc_pid;
655b3cdcd63SStefano Zampini         graph->cptr[ncc] = cum_queue;
656b3cdcd63SStefano Zampini       }
657b3cdcd63SStefano Zampini     }
658674ae819SStefano Zampini   }
659674ae819SStefano Zampini   graph->ncc = ncc;
660acc113dbSStefano Zampini   graph->queue_sorted = PETSC_FALSE;
661674ae819SStefano Zampini   PetscFunctionReturn(0);
662674ae819SStefano Zampini }
663674ae819SStefano Zampini 
664674ae819SStefano Zampini #undef __FUNCT__
665674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphSetUp"
666674ae819SStefano Zampini PetscErrorCode PCBDDCGraphSetUp(PCBDDCGraph graph, PetscInt custom_minimal_size, IS neumann_is, IS dirichlet_is, PetscInt n_ISForDofs, IS ISForDofs[], IS custom_primal_vertices)
667674ae819SStefano Zampini {
668674ae819SStefano Zampini   VecScatter     scatter_ctx;
669674ae819SStefano Zampini   Vec            local_vec,local_vec2,global_vec;
670dc456d91SStefano Zampini   IS             to,from,subset,subset_n;
671674ae819SStefano Zampini   MPI_Comm       comm;
672674ae819SStefano Zampini   PetscScalar    *array,*array2;
673674ae819SStefano Zampini   const PetscInt *is_indices;
674dc456d91SStefano Zampini   PetscInt       n_neigh,*neigh,*n_shared,**shared,*queue_global;
675674ae819SStefano Zampini   PetscInt       i,j,k,s,total_counts,nodes_touched,is_size;
676b3cdcd63SStefano Zampini   PetscMPIInt    commsize;
67751b0f6cfSStefano Zampini   PetscBool      same_set,mirrors_found;
6787babac9bSStefano Zampini   PetscErrorCode ierr;
679674ae819SStefano Zampini 
680674ae819SStefano Zampini   PetscFunctionBegin;
681a07ea27aSStefano Zampini   graph->has_dirichlet = PETSC_FALSE;
682a07ea27aSStefano Zampini   if (dirichlet_is) {
683a07ea27aSStefano Zampini     PetscCheckSameComm(graph->l2gmap,1,dirichlet_is,4);
684a07ea27aSStefano Zampini     graph->has_dirichlet = PETSC_TRUE;
685a07ea27aSStefano Zampini   }
686674ae819SStefano Zampini   ierr = PetscObjectGetComm((PetscObject)(graph->l2gmap),&comm);CHKERRQ(ierr);
687b3cdcd63SStefano Zampini   ierr = MPI_Comm_size(comm,&commsize);CHKERRQ(ierr);
688b3cdcd63SStefano Zampini 
689674ae819SStefano Zampini   /* custom_minimal_size */
69014f95afaSStefano Zampini   graph->custom_minimal_size = custom_minimal_size;
691674ae819SStefano Zampini   /* get info l2gmap and allocate work vectors  */
692674ae819SStefano Zampini   ierr = ISLocalToGlobalMappingGetInfo(graph->l2gmap,&n_neigh,&neigh,&n_shared,&shared);CHKERRQ(ierr);
6932635a1d4SStefano Zampini   /* check if we have any local periodic nodes (periodic BCs) */
6942635a1d4SStefano Zampini   mirrors_found = PETSC_FALSE;
6952635a1d4SStefano Zampini   if (graph->nvtxs && n_neigh) {
6962635a1d4SStefano Zampini     for (i=0; i<n_shared[0]; i++) graph->count[shared[0][i]] += 1;
6972635a1d4SStefano Zampini     for (i=0; i<n_shared[0]; i++) {
6982635a1d4SStefano Zampini       if (graph->count[shared[0][i]] > 1) {
6992635a1d4SStefano Zampini         mirrors_found = PETSC_TRUE;
7002635a1d4SStefano Zampini         break;
7012635a1d4SStefano Zampini       }
7022635a1d4SStefano Zampini     }
7032635a1d4SStefano Zampini   }
7042635a1d4SStefano Zampini   /* create some workspace objects */
7052635a1d4SStefano Zampini   local_vec = NULL;
7062635a1d4SStefano Zampini   local_vec2 = NULL;
7072635a1d4SStefano Zampini   global_vec = NULL;
7082635a1d4SStefano Zampini   to = NULL;
7092635a1d4SStefano Zampini   from = NULL;
7102635a1d4SStefano Zampini   scatter_ctx = NULL;
7112635a1d4SStefano Zampini   if (n_ISForDofs || dirichlet_is || neumann_is || custom_primal_vertices) {
712674ae819SStefano Zampini     ierr = VecCreate(PETSC_COMM_SELF,&local_vec);CHKERRQ(ierr);
713674ae819SStefano Zampini     ierr = VecSetSizes(local_vec,PETSC_DECIDE,graph->nvtxs);CHKERRQ(ierr);
714674ae819SStefano Zampini     ierr = VecSetType(local_vec,VECSTANDARD);CHKERRQ(ierr);
715674ae819SStefano Zampini     ierr = VecDuplicate(local_vec,&local_vec2);CHKERRQ(ierr);
716674ae819SStefano Zampini     ierr = VecCreate(comm,&global_vec);CHKERRQ(ierr);
7177fb0e2dbSStefano Zampini     ierr = VecSetSizes(global_vec,PETSC_DECIDE,graph->nvtxs_global);CHKERRQ(ierr);
718674ae819SStefano Zampini     ierr = VecSetType(global_vec,VECSTANDARD);CHKERRQ(ierr);
719674ae819SStefano Zampini     ierr = ISCreateStride(PETSC_COMM_SELF,graph->nvtxs,0,1,&to);CHKERRQ(ierr);
720674ae819SStefano Zampini     ierr = ISLocalToGlobalMappingApplyIS(graph->l2gmap,to,&from);CHKERRQ(ierr);
721674ae819SStefano Zampini     ierr = VecScatterCreate(global_vec,from,local_vec,to,&scatter_ctx);CHKERRQ(ierr);
7222635a1d4SStefano Zampini   } else if (mirrors_found) {
7232635a1d4SStefano Zampini     ierr = ISCreateStride(PETSC_COMM_SELF,graph->nvtxs,0,1,&to);CHKERRQ(ierr);
7242635a1d4SStefano Zampini     ierr = ISLocalToGlobalMappingApplyIS(graph->l2gmap,to,&from);CHKERRQ(ierr);
72549eeff8cSStefano Zampini   }
72651b0f6cfSStefano Zampini   /* compute local mirrors (if any) */
72751b0f6cfSStefano Zampini   if (mirrors_found) {
72851b0f6cfSStefano Zampini     PetscInt *local_indices,*global_indices;
72951b0f6cfSStefano Zampini     /* get arrays of local and global indices */
730785e854fSJed Brown     ierr = PetscMalloc1(graph->nvtxs,&local_indices);CHKERRQ(ierr);
73151b0f6cfSStefano Zampini     ierr = ISGetIndices(to,(const PetscInt**)&is_indices);CHKERRQ(ierr);
73251b0f6cfSStefano Zampini     ierr = PetscMemcpy(local_indices,is_indices,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
73351b0f6cfSStefano Zampini     ierr = ISRestoreIndices(to,(const PetscInt**)&is_indices);CHKERRQ(ierr);
734785e854fSJed Brown     ierr = PetscMalloc1(graph->nvtxs,&global_indices);CHKERRQ(ierr);
73551b0f6cfSStefano Zampini     ierr = ISGetIndices(from,(const PetscInt**)&is_indices);CHKERRQ(ierr);
73651b0f6cfSStefano Zampini     ierr = PetscMemcpy(global_indices,is_indices,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
73751b0f6cfSStefano Zampini     ierr = ISRestoreIndices(from,(const PetscInt**)&is_indices);CHKERRQ(ierr);
73851b0f6cfSStefano Zampini     /* allocate space for mirrors */
7398e5aaad5SJed Brown     ierr = PetscMalloc2(graph->nvtxs,&graph->mirrors,graph->nvtxs,&graph->mirrors_set);CHKERRQ(ierr);
74051b0f6cfSStefano Zampini     ierr = PetscMemzero(graph->mirrors,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
74151b0f6cfSStefano Zampini     graph->mirrors_set[0] = 0;
74251b0f6cfSStefano Zampini 
74351b0f6cfSStefano Zampini     k=0;
74451b0f6cfSStefano Zampini     for (i=0;i<n_shared[0];i++) {
74551b0f6cfSStefano Zampini       j=shared[0][i];
74651b0f6cfSStefano Zampini       if (graph->count[j] > 1) {
74751b0f6cfSStefano Zampini         graph->mirrors[j]++;
74851b0f6cfSStefano Zampini         k++;
74951b0f6cfSStefano Zampini       }
75051b0f6cfSStefano Zampini     }
75151b0f6cfSStefano Zampini     /* allocate space for set of mirrors */
752785e854fSJed Brown     ierr = PetscMalloc1(k,&graph->mirrors_set[0]);CHKERRQ(ierr);
75351b0f6cfSStefano Zampini     for (i=1;i<graph->nvtxs;i++)
75451b0f6cfSStefano Zampini       graph->mirrors_set[i]=graph->mirrors_set[i-1]+graph->mirrors[i-1];
75551b0f6cfSStefano Zampini 
75651b0f6cfSStefano Zampini     /* fill arrays */
75751b0f6cfSStefano Zampini     ierr = PetscMemzero(graph->mirrors,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
75851b0f6cfSStefano Zampini     for (j=0;j<n_shared[0];j++) {
75951b0f6cfSStefano Zampini       i=shared[0][j];
76051b0f6cfSStefano Zampini       if (graph->count[i] > 1)
76151b0f6cfSStefano Zampini         graph->mirrors_set[i][graph->mirrors[i]++]=global_indices[i];
76251b0f6cfSStefano Zampini     }
76351b0f6cfSStefano Zampini     ierr = PetscSortIntWithArray(graph->nvtxs,global_indices,local_indices);CHKERRQ(ierr);
76451b0f6cfSStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
76551b0f6cfSStefano Zampini       if (graph->mirrors[i] > 0) {
76651b0f6cfSStefano Zampini         ierr = PetscFindInt(graph->mirrors_set[i][0],graph->nvtxs,global_indices,&k);CHKERRQ(ierr);
76751b0f6cfSStefano Zampini         j = global_indices[k];
76851b0f6cfSStefano Zampini         while ( k > 0 && global_indices[k-1] == j) k--;
76951b0f6cfSStefano Zampini         for (j=0;j<graph->mirrors[i];j++) {
77051b0f6cfSStefano Zampini           graph->mirrors_set[i][j]=local_indices[k+j];
77151b0f6cfSStefano Zampini         }
77251b0f6cfSStefano Zampini         ierr = PetscSortInt(graph->mirrors[i],graph->mirrors_set[i]);CHKERRQ(ierr);
77351b0f6cfSStefano Zampini       }
77451b0f6cfSStefano Zampini     }
77551b0f6cfSStefano Zampini     ierr = PetscFree(local_indices);CHKERRQ(ierr);
77651b0f6cfSStefano Zampini     ierr = PetscFree(global_indices);CHKERRQ(ierr);
77751b0f6cfSStefano Zampini   }
77851b0f6cfSStefano Zampini   ierr = PetscMemzero(graph->count,graph->nvtxs*sizeof(*graph->count));CHKERRQ(ierr);
779674ae819SStefano Zampini   ierr = ISDestroy(&to);CHKERRQ(ierr);
780674ae819SStefano Zampini   ierr = ISDestroy(&from);CHKERRQ(ierr);
78151b0f6cfSStefano Zampini 
782674ae819SStefano Zampini   /* Count total number of neigh per node */
783674ae819SStefano Zampini   k = 0;
784674ae819SStefano Zampini   for (i=1;i<n_neigh;i++) {
785674ae819SStefano Zampini     k += n_shared[i];
786674ae819SStefano Zampini     for (j=0;j<n_shared[i];j++) {
787674ae819SStefano Zampini       graph->count[shared[i][j]] += 1;
788674ae819SStefano Zampini     }
789674ae819SStefano Zampini   }
790674ae819SStefano Zampini   /* Allocate space for storing the set of neighbours for each node */
791674ae819SStefano Zampini   if (graph->nvtxs) {
792785e854fSJed Brown     ierr = PetscMalloc1(k,&graph->neighbours_set[0]);CHKERRQ(ierr);
793674ae819SStefano Zampini   }
794674ae819SStefano Zampini   for (i=1;i<graph->nvtxs;i++) { /* dont count myself */
795674ae819SStefano Zampini     graph->neighbours_set[i]=graph->neighbours_set[i-1]+graph->count[i-1];
796674ae819SStefano Zampini   }
797674ae819SStefano Zampini   /* Get information for sharing subdomains */
798674ae819SStefano Zampini   ierr = PetscMemzero(graph->count,graph->nvtxs*sizeof(*graph->count));CHKERRQ(ierr);
799674ae819SStefano Zampini   for (i=1;i<n_neigh;i++) { /* dont count myself */
800674ae819SStefano Zampini     s = n_shared[i];
801674ae819SStefano Zampini     for (j=0;j<s;j++) {
802674ae819SStefano Zampini       k = shared[i][j];
803674ae819SStefano Zampini       graph->neighbours_set[k][graph->count[k]] = neigh[i];
804674ae819SStefano Zampini       graph->count[k] += 1;
805674ae819SStefano Zampini     }
806674ae819SStefano Zampini   }
807674ae819SStefano Zampini   /* sort set of sharing subdomains */
808674ae819SStefano Zampini   for (i=0;i<graph->nvtxs;i++) {
80993fb5fd3SStefano Zampini     ierr = PetscSortRemoveDupsInt(&graph->count[i],graph->neighbours_set[i]);CHKERRQ(ierr);
810674ae819SStefano Zampini   }
8117fb0e2dbSStefano Zampini   /* free memory allocated by ISLocalToGlobalMappingGetInfo */
8127fb0e2dbSStefano Zampini   ierr = ISLocalToGlobalMappingRestoreInfo(graph->l2gmap,&n_neigh,&neigh,&n_shared,&shared);CHKERRQ(ierr);
8137fb0e2dbSStefano Zampini 
81467c9da69SStefano Zampini   /*
81567c9da69SStefano Zampini      Get info for dofs splitting
8165777c63cSStefano Zampini      User can specify just a subset; an additional field is considered as a complementary field
81767c9da69SStefano Zampini   */
818b3cdcd63SStefano Zampini   for (i=0;i<graph->nvtxs;i++) graph->which_dof[i] = n_ISForDofs; /* by default a dof belongs to the complement set */
8195777c63cSStefano Zampini   if (n_ISForDofs) {
8205777c63cSStefano Zampini     ierr = VecSet(local_vec,-1.0);CHKERRQ(ierr);
8215777c63cSStefano Zampini   }
822674ae819SStefano Zampini   for (i=0;i<n_ISForDofs;i++) {
8235777c63cSStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
82463602bcaSStefano Zampini     ierr = ISGetLocalSize(ISForDofs[i],&is_size);CHKERRQ(ierr);
825674ae819SStefano Zampini     ierr = ISGetIndices(ISForDofs[i],(const PetscInt**)&is_indices);CHKERRQ(ierr);
826674ae819SStefano Zampini     for (j=0;j<is_size;j++) {
827d50ed60aSStefano Zampini       if (is_indices[j] > -1 && is_indices[j] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
828674ae819SStefano Zampini         graph->which_dof[is_indices[j]] = i;
8295777c63cSStefano Zampini         array[is_indices[j]] = 1.*i;
830674ae819SStefano Zampini       }
83167c9da69SStefano Zampini     }
832674ae819SStefano Zampini     ierr = ISRestoreIndices(ISForDofs[i],(const PetscInt**)&is_indices);CHKERRQ(ierr);
8335777c63cSStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
8345777c63cSStefano Zampini   }
8355777c63cSStefano Zampini   /* Check consistency among neighbours */
8365777c63cSStefano Zampini   if (n_ISForDofs) {
8375777c63cSStefano Zampini     ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8385777c63cSStefano Zampini     ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8395777c63cSStefano Zampini     ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8405777c63cSStefano Zampini     ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8415777c63cSStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
8425777c63cSStefano Zampini     ierr = VecGetArray(local_vec2,&array2);CHKERRQ(ierr);
8435777c63cSStefano Zampini     for (i=0;i<graph->nvtxs;i++){
8445777c63cSStefano Zampini       PetscInt field1,field2;
8455777c63cSStefano Zampini 
8465777c63cSStefano Zampini       field1 = (PetscInt)PetscRealPart(array[i]);
8475777c63cSStefano Zampini       field2 = (PetscInt)PetscRealPart(array2[i]);
8486c4ed002SBarry Smith       if (field1 != field2) SETERRQ3(comm,PETSC_ERR_USER,"Local node %D have been assigned two different field ids %D and %D at the same time\n",i,field1,field2);
8495777c63cSStefano Zampini     }
8505777c63cSStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
8515777c63cSStefano Zampini     ierr = VecRestoreArray(local_vec2,&array2);CHKERRQ(ierr);
852674ae819SStefano Zampini   }
853674ae819SStefano Zampini   /* Take into account Neumann nodes */
854674ae819SStefano Zampini   if (neumann_is) {
8557e5be6c5SStefano Zampini     ierr = VecSet(local_vec,0.0);CHKERRQ(ierr);
856674ae819SStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
85782ba6b80SStefano Zampini     ierr = ISGetLocalSize(neumann_is,&is_size);CHKERRQ(ierr);
858674ae819SStefano Zampini     ierr = ISGetIndices(neumann_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
859674ae819SStefano Zampini     for (i=0;i<is_size;i++) {
860d50ed60aSStefano Zampini       if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
861674ae819SStefano Zampini         array[is_indices[i]] = 1.0;
862674ae819SStefano Zampini       }
8633c629590SStefano Zampini     }
864674ae819SStefano Zampini     ierr = ISRestoreIndices(neumann_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
865674ae819SStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
866674ae819SStefano Zampini     /* Neumann nodes: impose consistency among neighbours */
867674ae819SStefano Zampini     ierr = VecSet(global_vec,0.0);CHKERRQ(ierr);
868674ae819SStefano Zampini     ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
869674ae819SStefano Zampini     ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
870674ae819SStefano Zampini     ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
871674ae819SStefano Zampini     ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
872674ae819SStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
873674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
8743c629590SStefano Zampini       if (PetscRealPart(array[i]) > 0.1) {
8750cf82fd6SStefano Zampini         graph->special_dof[i] = PCBDDCGRAPH_NEUMANN_MARK;
876674ae819SStefano Zampini       }
877674ae819SStefano Zampini     }
878674ae819SStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
8797e5be6c5SStefano Zampini   }
8807e5be6c5SStefano Zampini   /* Take into account Dirichlet nodes (they overwrite any Neumann boundary node previously set) */
881674ae819SStefano Zampini   if (dirichlet_is) {
8827e5be6c5SStefano Zampini     ierr = VecSet(local_vec,0.0);CHKERRQ(ierr);
883674ae819SStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
88482ba6b80SStefano Zampini     ierr = ISGetLocalSize(dirichlet_is,&is_size);CHKERRQ(ierr);
885674ae819SStefano Zampini     ierr = ISGetIndices(dirichlet_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
886674ae819SStefano Zampini     for (i=0;i<is_size;i++){
887d50ed60aSStefano Zampini       if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
8887e5be6c5SStefano Zampini         array[is_indices[i]] = 1.0;
889674ae819SStefano Zampini       }
8903c629590SStefano Zampini     }
891674ae819SStefano Zampini     ierr = ISRestoreIndices(dirichlet_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
892674ae819SStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
893674ae819SStefano Zampini     /* Dirichlet nodes: impose consistency among neighbours */
894674ae819SStefano Zampini     ierr = VecSet(global_vec,0.0);CHKERRQ(ierr);
8957e5be6c5SStefano Zampini     ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8967e5be6c5SStefano Zampini     ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
897674ae819SStefano Zampini     ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
898674ae819SStefano Zampini     ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
899674ae819SStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
900674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
9013c629590SStefano Zampini       if (PetscRealPart(array[i]) > 0.1) {
9024b2aedd3SStefano Zampini         if (commsize > graph->commsizelimit) { /* dirichlet nodes treated as internal */
903df48933bSStefano Zampini           ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
904b3cdcd63SStefano Zampini           graph->subset[i] = 0;
905b3cdcd63SStefano Zampini         }
9060cf82fd6SStefano Zampini         graph->special_dof[i] = PCBDDCGRAPH_DIRICHLET_MARK;
907674ae819SStefano Zampini       }
908674ae819SStefano Zampini     }
909674ae819SStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
9107fb0e2dbSStefano Zampini   }
91108a5cf49SStefano Zampini   /* mark local periodic nodes (if any) and adapt CSR graph (if any) */
91251b0f6cfSStefano Zampini   if (graph->mirrors) {
91351b0f6cfSStefano Zampini     for (i=0;i<graph->nvtxs;i++)
91451b0f6cfSStefano Zampini       if (graph->mirrors[i])
9150cf82fd6SStefano Zampini         graph->special_dof[i] = PCBDDCGRAPH_LOCAL_PERIODIC_MARK;
91651b0f6cfSStefano Zampini 
91708a5cf49SStefano Zampini     if (graph->xadj && graph->adjncy) {
91808a5cf49SStefano Zampini       PetscInt *new_xadj,*new_adjncy;
91951b0f6cfSStefano Zampini       /* sort CSR graph */
92051b0f6cfSStefano Zampini       for (i=0;i<graph->nvtxs;i++)
92151b0f6cfSStefano Zampini         ierr = PetscSortInt(graph->xadj[i+1]-graph->xadj[i],&graph->adjncy[graph->xadj[i]]);CHKERRQ(ierr);
92251b0f6cfSStefano Zampini 
92351b0f6cfSStefano Zampini       /* adapt local CSR graph in case of local periodicity */
92451b0f6cfSStefano Zampini       k = 0;
92551b0f6cfSStefano Zampini       for (i=0;i<graph->nvtxs;i++)
92651b0f6cfSStefano Zampini         for (j=graph->xadj[i];j<graph->xadj[i+1];j++)
92751b0f6cfSStefano Zampini           k += graph->mirrors[graph->adjncy[j]];
92851b0f6cfSStefano Zampini 
929854ce69bSBarry Smith       ierr = PetscMalloc1(graph->nvtxs+1,&new_xadj);CHKERRQ(ierr);
930854ce69bSBarry Smith       ierr = PetscMalloc1(k+graph->xadj[graph->nvtxs],&new_adjncy);CHKERRQ(ierr);
93151b0f6cfSStefano Zampini       new_xadj[0] = 0;
93251b0f6cfSStefano Zampini       for (i=0;i<graph->nvtxs;i++) {
93351b0f6cfSStefano Zampini         k = graph->xadj[i+1]-graph->xadj[i];
93451b0f6cfSStefano Zampini         ierr = PetscMemcpy(&new_adjncy[new_xadj[i]],&graph->adjncy[graph->xadj[i]],k*sizeof(PetscInt));CHKERRQ(ierr);
93551b0f6cfSStefano Zampini         new_xadj[i+1] = new_xadj[i]+k;
93651b0f6cfSStefano Zampini         for (j=graph->xadj[i];j<graph->xadj[i+1];j++) {
93751b0f6cfSStefano Zampini           k = graph->mirrors[graph->adjncy[j]];
93851b0f6cfSStefano Zampini           ierr = PetscMemcpy(&new_adjncy[new_xadj[i+1]],graph->mirrors_set[graph->adjncy[j]],k*sizeof(PetscInt));CHKERRQ(ierr);
93951b0f6cfSStefano Zampini           new_xadj[i+1] += k;
94051b0f6cfSStefano Zampini         }
94151b0f6cfSStefano Zampini         k = new_xadj[i+1]-new_xadj[i];
94251b0f6cfSStefano Zampini         ierr = PetscSortRemoveDupsInt(&k,&new_adjncy[new_xadj[i]]);CHKERRQ(ierr);
94351b0f6cfSStefano Zampini         new_xadj[i+1] = new_xadj[i]+k;
94451b0f6cfSStefano Zampini       }
94551b0f6cfSStefano Zampini       /* set new CSR into graph */
94651b0f6cfSStefano Zampini       ierr = PetscFree(graph->xadj);CHKERRQ(ierr);
94751b0f6cfSStefano Zampini       ierr = PetscFree(graph->adjncy);CHKERRQ(ierr);
94851b0f6cfSStefano Zampini       graph->xadj = new_xadj;
94951b0f6cfSStefano Zampini       graph->adjncy = new_adjncy;
95051b0f6cfSStefano Zampini     }
95108a5cf49SStefano Zampini   }
95251b0f6cfSStefano Zampini 
95317eb1463SStefano Zampini   /* mark special nodes (if any) -> each will become a single node equivalence class */
954674ae819SStefano Zampini   if (custom_primal_vertices) {
95517eb1463SStefano Zampini     ierr = VecSet(local_vec,0.0);CHKERRQ(ierr);
9569b70a373SStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
95763602bcaSStefano Zampini     ierr = ISGetLocalSize(custom_primal_vertices,&is_size);CHKERRQ(ierr);
958674ae819SStefano Zampini     ierr = ISGetIndices(custom_primal_vertices,(const PetscInt**)&is_indices);CHKERRQ(ierr);
959674ae819SStefano Zampini     for (i=0;i<is_size;i++){
9609b70a373SStefano Zampini       if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
9619b70a373SStefano Zampini         array[is_indices[i]] = 1.0;
9629b70a373SStefano Zampini       }
963674ae819SStefano Zampini     }
964674ae819SStefano Zampini     ierr = ISRestoreIndices(custom_primal_vertices,(const PetscInt**)&is_indices);CHKERRQ(ierr);
9659b70a373SStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
9669b70a373SStefano Zampini     /* special nodes: impose consistency among neighbours */
9679b70a373SStefano Zampini     ierr = VecSet(global_vec,0.0);CHKERRQ(ierr);
9689b70a373SStefano Zampini     ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9699b70a373SStefano Zampini     ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9709b70a373SStefano Zampini     ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9719b70a373SStefano Zampini     ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9729b70a373SStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
9739b70a373SStefano Zampini     j = 0;
9749b70a373SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
9759b70a373SStefano Zampini       if (PetscRealPart(array[i]) > 0.1 && graph->special_dof[i] != PCBDDCGRAPH_DIRICHLET_MARK) {
9769b70a373SStefano Zampini         graph->special_dof[i] = PCBDDCGRAPH_SPECIAL_MARK-j;
9779b70a373SStefano Zampini         j++;
9789b70a373SStefano Zampini       }
9799b70a373SStefano Zampini     }
9809b70a373SStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
98117eb1463SStefano Zampini   }
9829b70a373SStefano Zampini 
9834b2aedd3SStefano Zampini   /* mark interior nodes (if commsize > graph->commsizelimit) as touched and belonging to partition number 0 */
9844b2aedd3SStefano Zampini   if (commsize > graph->commsizelimit) {
985674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
986674ae819SStefano Zampini       if (!graph->count[i]) {
987df48933bSStefano Zampini         ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
988674ae819SStefano Zampini         graph->subset[i] = 0;
989674ae819SStefano Zampini       }
990674ae819SStefano Zampini     }
991b3cdcd63SStefano Zampini   }
992b3cdcd63SStefano Zampini 
993674ae819SStefano Zampini   /* init graph structure and compute default subsets */
994674ae819SStefano Zampini   nodes_touched = 0;
995674ae819SStefano Zampini   for (i=0;i<graph->nvtxs;i++) {
996df48933bSStefano Zampini     if (PetscBTLookup(graph->touched,i)) {
997674ae819SStefano Zampini       nodes_touched++;
998674ae819SStefano Zampini     }
999674ae819SStefano Zampini   }
1000674ae819SStefano Zampini   i = 0;
1001674ae819SStefano Zampini   graph->ncc = 0;
1002674ae819SStefano Zampini   total_counts = 0;
10037babac9bSStefano Zampini 
10047babac9bSStefano Zampini   /* allocated space for queues */
10054b2aedd3SStefano Zampini   if (commsize == graph->commsizelimit) {
10067babac9bSStefano Zampini     ierr = PetscMalloc2(graph->nvtxs+1,&graph->cptr,graph->nvtxs,&graph->queue);CHKERRQ(ierr);
10077babac9bSStefano Zampini   } else {
10087babac9bSStefano Zampini     PetscInt nused = graph->nvtxs - nodes_touched;
10097babac9bSStefano Zampini     ierr = PetscMalloc2(nused+1,&graph->cptr,nused,&graph->queue);CHKERRQ(ierr);
10107babac9bSStefano Zampini   }
10117babac9bSStefano Zampini 
1012674ae819SStefano Zampini   while (nodes_touched<graph->nvtxs) {
1013674ae819SStefano Zampini     /*  find first untouched node in local ordering */
1014b3cdcd63SStefano Zampini     while (PetscBTLookup(graph->touched,i)) i++;
1015df48933bSStefano Zampini     ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
1016674ae819SStefano Zampini     graph->subset[i] = graph->ncc+1;
1017674ae819SStefano Zampini     graph->cptr[graph->ncc] = total_counts;
1018674ae819SStefano Zampini     graph->queue[total_counts] = i;
1019674ae819SStefano Zampini     total_counts++;
1020674ae819SStefano Zampini     nodes_touched++;
1021674ae819SStefano Zampini     /* now find all other nodes having the same set of sharing subdomains */
1022674ae819SStefano Zampini     for (j=i+1;j<graph->nvtxs;j++) {
102374e413f5SStefano Zampini       /* check for same number of sharing subdomains, dof number and same special mark */
1024df48933bSStefano 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]) {
1025674ae819SStefano Zampini         /* check for same set of sharing subdomains */
1026674ae819SStefano Zampini         same_set = PETSC_TRUE;
1027674ae819SStefano Zampini         for (k=0;k<graph->count[j];k++){
1028674ae819SStefano Zampini           if (graph->neighbours_set[i][k] != graph->neighbours_set[j][k]) {
1029674ae819SStefano Zampini             same_set = PETSC_FALSE;
1030674ae819SStefano Zampini           }
1031674ae819SStefano Zampini         }
1032674ae819SStefano Zampini         /* I found a friend of mine */
1033674ae819SStefano Zampini         if (same_set) {
1034df48933bSStefano Zampini           ierr = PetscBTSet(graph->touched,j);CHKERRQ(ierr);
1035b3cdcd63SStefano Zampini           graph->subset[j] = graph->ncc+1;
1036674ae819SStefano Zampini           nodes_touched++;
1037674ae819SStefano Zampini           graph->queue[total_counts] = j;
1038674ae819SStefano Zampini           total_counts++;
1039674ae819SStefano Zampini         }
1040674ae819SStefano Zampini       }
1041674ae819SStefano Zampini     }
1042674ae819SStefano Zampini     graph->ncc++;
1043674ae819SStefano Zampini   }
1044b3cdcd63SStefano 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 */
1045674ae819SStefano Zampini   graph->n_subsets = graph->ncc;
1046785e854fSJed Brown   ierr = PetscMalloc1(graph->n_subsets,&graph->subset_ncc);CHKERRQ(ierr);
1047674ae819SStefano Zampini   for (i=0;i<graph->n_subsets;i++) {
1048674ae819SStefano Zampini     graph->subset_ncc[i] = 1;
1049674ae819SStefano Zampini   }
1050674ae819SStefano Zampini   /* final pointer */
1051674ae819SStefano Zampini   graph->cptr[graph->ncc] = total_counts;
1052ec1c413dSStefano Zampini 
1053b3cdcd63SStefano Zampini   /* For consistency reasons (among neighbours), I need to sort (by global ordering) each connected component */
1054ec1c413dSStefano Zampini   /* Get a reference node (min index in global ordering) for each subset for tagging messages */
1055dc456d91SStefano Zampini   ierr = PetscMalloc1(graph->ncc,&graph->subset_ref_node);CHKERRQ(ierr);
1056dc456d91SStefano Zampini   ierr = PetscMalloc1(graph->cptr[graph->ncc],&queue_global);CHKERRQ(ierr);
10573a5b03d1SStefano Zampini   ierr = ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_global);CHKERRQ(ierr);
1058ec1c413dSStefano Zampini   for (j=0;j<graph->ncc;j++) {
1059ec1c413dSStefano Zampini     ierr = PetscSortIntWithArray(graph->cptr[j+1]-graph->cptr[j],&queue_global[graph->cptr[j]],&graph->queue[graph->cptr[j]]);CHKERRQ(ierr);
1060dc456d91SStefano Zampini     graph->subset_ref_node[j] = graph->queue[graph->cptr[j]];
1061f0321c90SStefano Zampini   }
1062dc456d91SStefano Zampini   ierr = PetscFree(queue_global);CHKERRQ(ierr);
1063acc113dbSStefano Zampini   graph->queue_sorted = PETSC_TRUE;
1064b3cdcd63SStefano Zampini   /* save information on subsets (needed when analyzing the connected components) */
10652f566a09SStefano Zampini   if (graph->ncc) {
1066b3cdcd63SStefano Zampini     ierr = PetscMalloc2(graph->ncc,&graph->subset_size,graph->ncc,&graph->subset_idxs);CHKERRQ(ierr);
1067b3cdcd63SStefano Zampini     ierr = PetscMalloc1(graph->cptr[graph->ncc],&graph->subset_idxs[0]);CHKERRQ(ierr);
1068b3cdcd63SStefano Zampini     ierr = PetscMemzero(graph->subset_idxs[0],graph->cptr[graph->ncc]*sizeof(PetscInt));CHKERRQ(ierr);
1069ec1c413dSStefano Zampini     for (j=1;j<graph->ncc;j++) {
1070b3cdcd63SStefano Zampini       graph->subset_size[j-1] = graph->cptr[j] - graph->cptr[j-1];
1071b3cdcd63SStefano Zampini       graph->subset_idxs[j] = graph->subset_idxs[j-1] + graph->subset_size[j-1];
1072ec1c413dSStefano Zampini     }
1073b3cdcd63SStefano Zampini     graph->subset_size[graph->ncc-1] = graph->cptr[graph->ncc] - graph->cptr[graph->ncc-1];
1074b3cdcd63SStefano Zampini     ierr = PetscMemcpy(graph->subset_idxs[0],graph->queue,graph->cptr[graph->ncc]*sizeof(PetscInt));CHKERRQ(ierr);
1075ec1c413dSStefano Zampini   }
1076b3cdcd63SStefano Zampini 
1077f0321c90SStefano Zampini   /* renumber reference nodes */
1078dc456d91SStefano Zampini   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)(graph->l2gmap)),graph->ncc,graph->subset_ref_node,PETSC_COPY_VALUES,&subset_n);CHKERRQ(ierr);
1079dc456d91SStefano Zampini   ierr = ISLocalToGlobalMappingApplyIS(graph->l2gmap,subset_n,&subset);CHKERRQ(ierr);
1080dc456d91SStefano Zampini   ierr = ISDestroy(&subset_n);CHKERRQ(ierr);
10816583bcc1SStefano Zampini   ierr = ISRenumber(subset,NULL,NULL,&subset_n);CHKERRQ(ierr);
1082dc456d91SStefano Zampini   ierr = ISDestroy(&subset);CHKERRQ(ierr);
1083dc456d91SStefano Zampini   ierr = ISGetLocalSize(subset_n,&k);CHKERRQ(ierr);
10846c4ed002SBarry Smith   if (k != graph->ncc) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid size of new subset! %D != %D",k,graph->ncc);
1085dc456d91SStefano Zampini   ierr = ISGetIndices(subset_n,&is_indices);CHKERRQ(ierr);
1086dc456d91SStefano Zampini   ierr = PetscMemcpy(graph->subset_ref_node,is_indices,graph->ncc*sizeof(PetscInt));CHKERRQ(ierr);
1087dc456d91SStefano Zampini   ierr = ISRestoreIndices(subset_n,&is_indices);CHKERRQ(ierr);
1088dc456d91SStefano Zampini   ierr = ISDestroy(&subset_n);CHKERRQ(ierr);
1089ec1c413dSStefano Zampini 
1090ec1c413dSStefano Zampini   /* free workspace */
1091674ae819SStefano Zampini   ierr = VecDestroy(&local_vec);CHKERRQ(ierr);
1092674ae819SStefano Zampini   ierr = VecDestroy(&local_vec2);CHKERRQ(ierr);
1093674ae819SStefano Zampini   ierr = VecDestroy(&global_vec);CHKERRQ(ierr);
1094674ae819SStefano Zampini   ierr = VecScatterDestroy(&scatter_ctx);CHKERRQ(ierr);
1095b3cdcd63SStefano Zampini   graph->setupcalled = PETSC_TRUE;
1096674ae819SStefano Zampini   PetscFunctionReturn(0);
1097674ae819SStefano Zampini }
1098674ae819SStefano Zampini 
1099674ae819SStefano Zampini #undef __FUNCT__
1100674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphResetCSR"
1101674ae819SStefano Zampini PetscErrorCode PCBDDCGraphResetCSR(PCBDDCGraph graph)
1102674ae819SStefano Zampini {
1103674ae819SStefano Zampini   PetscErrorCode ierr;
1104674ae819SStefano Zampini 
1105674ae819SStefano Zampini   PetscFunctionBegin;
1106a1dbd327SStefano Zampini   if (graph->freecsr) {
1107674ae819SStefano Zampini     ierr = PetscFree(graph->xadj);CHKERRQ(ierr);
1108674ae819SStefano Zampini     ierr = PetscFree(graph->adjncy);CHKERRQ(ierr);
1109a1dbd327SStefano Zampini   } else {
1110a1dbd327SStefano Zampini     graph->xadj = NULL;
1111a1dbd327SStefano Zampini     graph->adjncy = NULL;
1112a1dbd327SStefano Zampini   }
1113c8272957SStefano Zampini   graph->freecsr = PETSC_FALSE;
1114575ad6abSStefano Zampini   graph->nvtxs_csr = 0;
1115674ae819SStefano Zampini   PetscFunctionReturn(0);
1116674ae819SStefano Zampini }
1117674ae819SStefano Zampini 
1118674ae819SStefano Zampini #undef __FUNCT__
1119674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphReset"
1120674ae819SStefano Zampini PetscErrorCode PCBDDCGraphReset(PCBDDCGraph graph)
1121674ae819SStefano Zampini {
1122674ae819SStefano Zampini   PetscErrorCode ierr;
1123674ae819SStefano Zampini 
1124674ae819SStefano Zampini   PetscFunctionBegin;
1125674ae819SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&graph->l2gmap);CHKERRQ(ierr);
1126674ae819SStefano Zampini   ierr = PetscFree(graph->subset_ncc);CHKERRQ(ierr);
11273a5b03d1SStefano Zampini   ierr = PetscFree(graph->subset_ref_node);CHKERRQ(ierr);
1128674ae819SStefano Zampini   if (graph->nvtxs) {
1129674ae819SStefano Zampini     ierr = PetscFree(graph->neighbours_set[0]);CHKERRQ(ierr);
1130674ae819SStefano Zampini   }
1131df48933bSStefano Zampini   ierr = PetscBTDestroy(&graph->touched);CHKERRQ(ierr);
11327babac9bSStefano Zampini   ierr = PetscFree5(graph->count,
1133674ae819SStefano Zampini                     graph->neighbours_set,
1134674ae819SStefano Zampini                     graph->subset,
1135674ae819SStefano Zampini                     graph->which_dof,
1136df48933bSStefano Zampini                     graph->special_dof);CHKERRQ(ierr);
11377babac9bSStefano Zampini   ierr = PetscFree2(graph->cptr,graph->queue);CHKERRQ(ierr);
113851b0f6cfSStefano Zampini   if (graph->mirrors) {
113951b0f6cfSStefano Zampini     ierr = PetscFree(graph->mirrors_set[0]);CHKERRQ(ierr);
114051b0f6cfSStefano Zampini   }
114151b0f6cfSStefano Zampini   ierr = PetscFree2(graph->mirrors,graph->mirrors_set);CHKERRQ(ierr);
1142b3cdcd63SStefano Zampini   if (graph->subset_idxs) {
1143b3cdcd63SStefano Zampini     ierr = PetscFree(graph->subset_idxs[0]);CHKERRQ(ierr);
1144ec1c413dSStefano Zampini   }
1145b3cdcd63SStefano Zampini   ierr = PetscFree2(graph->subset_size,graph->subset_idxs);CHKERRQ(ierr);
1146a07ea27aSStefano Zampini   ierr = ISDestroy(&graph->dirdofs);CHKERRQ(ierr);
1147d62866d3SStefano Zampini   ierr = ISDestroy(&graph->dirdofsB);CHKERRQ(ierr);
11484f1b2e48SStefano Zampini   if (graph->n_local_subs) {
11494f1b2e48SStefano Zampini     ierr = PetscFree(graph->local_subs);CHKERRQ(ierr);
11504f1b2e48SStefano Zampini   }
1151a07ea27aSStefano Zampini   graph->has_dirichlet       = PETSC_FALSE;
1152674ae819SStefano Zampini   graph->nvtxs               = 0;
11537babac9bSStefano Zampini   graph->nvtxs_global        = 0;
1154674ae819SStefano Zampini   graph->n_subsets           = 0;
1155674ae819SStefano Zampini   graph->custom_minimal_size = 1;
11564f1b2e48SStefano Zampini   graph->n_local_subs        = 0;
1157*be12c134Sstefano_zampini   graph->maxcount            = PETSC_MAX_INT;
1158fb597685SStefano Zampini   graph->setupcalled         = PETSC_FALSE;
1159674ae819SStefano Zampini   PetscFunctionReturn(0);
1160674ae819SStefano Zampini }
1161674ae819SStefano Zampini 
1162674ae819SStefano Zampini #undef __FUNCT__
1163674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphInit"
1164*be12c134Sstefano_zampini PetscErrorCode PCBDDCGraphInit(PCBDDCGraph graph, ISLocalToGlobalMapping l2gmap, PetscInt N, PetscInt maxcount)
1165674ae819SStefano Zampini {
1166674ae819SStefano Zampini   PetscInt       n;
1167674ae819SStefano Zampini   PetscErrorCode ierr;
1168674ae819SStefano Zampini 
1169674ae819SStefano Zampini   PetscFunctionBegin;
1170674ae819SStefano Zampini   PetscValidPointer(graph,1);
1171674ae819SStefano Zampini   PetscValidHeaderSpecific(l2gmap,IS_LTOGM_CLASSID,2);
11727babac9bSStefano Zampini   PetscValidLogicalCollectiveInt(l2gmap,N,3);
1173*be12c134Sstefano_zampini   PetscValidLogicalCollectiveInt(l2gmap,maxcount,4);
11748e61c736SStefano Zampini   /* raise an error if already allocated */
11756c4ed002SBarry Smith   if (graph->nvtxs_global) SETERRQ(PetscObjectComm((PetscObject)l2gmap),PETSC_ERR_PLIB,"BDDCGraph already initialized");
1176674ae819SStefano Zampini   /* set number of vertices */
1177674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)l2gmap);CHKERRQ(ierr);
1178674ae819SStefano Zampini   graph->l2gmap = l2gmap;
1179674ae819SStefano Zampini   ierr = ISLocalToGlobalMappingGetSize(l2gmap,&n);CHKERRQ(ierr);
1180674ae819SStefano Zampini   graph->nvtxs = n;
11817fb0e2dbSStefano Zampini   graph->nvtxs_global = N;
1182674ae819SStefano Zampini   /* allocate used space */
1183df48933bSStefano Zampini   ierr = PetscBTCreate(graph->nvtxs,&graph->touched);CHKERRQ(ierr);
11847babac9bSStefano Zampini   ierr = PetscMalloc5(graph->nvtxs,&graph->count,
11858e5aaad5SJed Brown                       graph->nvtxs,&graph->neighbours_set,
11868e5aaad5SJed Brown                       graph->nvtxs,&graph->subset,
11878e5aaad5SJed Brown                       graph->nvtxs,&graph->which_dof,
11888e5aaad5SJed Brown                       graph->nvtxs,&graph->special_dof);CHKERRQ(ierr);
1189674ae819SStefano Zampini   /* zeroes memory */
1190674ae819SStefano Zampini   ierr = PetscMemzero(graph->count,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
1191674ae819SStefano Zampini   ierr = PetscMemzero(graph->subset,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
119263602bcaSStefano Zampini   /* use -1 as a default value for which_dof array */
119363602bcaSStefano Zampini   for (n=0;n<graph->nvtxs;n++) graph->which_dof[n] = -1;
119474e413f5SStefano Zampini   ierr = PetscMemzero(graph->special_dof,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
1195674ae819SStefano Zampini   /* zeroes first pointer to neighbour set */
1196674ae819SStefano Zampini   if (graph->nvtxs) {
1197674ae819SStefano Zampini     graph->neighbours_set[0] = 0;
1198674ae819SStefano Zampini   }
1199674ae819SStefano Zampini   /* zeroes workspace for values of ncc */
1200674ae819SStefano Zampini   graph->subset_ncc = 0;
12013a5b03d1SStefano Zampini   graph->subset_ref_node = 0;
1202a1dbd327SStefano Zampini   /* default flag for csr */
1203c8272957SStefano Zampini   graph->freecsr = PETSC_FALSE;
1204*be12c134Sstefano_zampini   /* maxcount for cc */
1205*be12c134Sstefano_zampini   graph->maxcount = maxcount;
1206674ae819SStefano Zampini   PetscFunctionReturn(0);
1207674ae819SStefano Zampini }
1208674ae819SStefano Zampini 
1209674ae819SStefano Zampini #undef __FUNCT__
1210674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphDestroy"
1211674ae819SStefano Zampini PetscErrorCode PCBDDCGraphDestroy(PCBDDCGraph* graph)
1212674ae819SStefano Zampini {
1213674ae819SStefano Zampini   PetscErrorCode ierr;
1214674ae819SStefano Zampini 
1215674ae819SStefano Zampini   PetscFunctionBegin;
1216674ae819SStefano Zampini   ierr = PCBDDCGraphReset(*graph);CHKERRQ(ierr);
1217674ae819SStefano Zampini   ierr = PetscFree(*graph);CHKERRQ(ierr);
1218674ae819SStefano Zampini   PetscFunctionReturn(0);
1219674ae819SStefano Zampini }
1220674ae819SStefano Zampini 
1221674ae819SStefano Zampini #undef __FUNCT__
1222674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphCreate"
1223674ae819SStefano Zampini PetscErrorCode PCBDDCGraphCreate(PCBDDCGraph *graph)
1224674ae819SStefano Zampini {
1225674ae819SStefano Zampini   PCBDDCGraph    new_graph;
1226674ae819SStefano Zampini   PetscErrorCode ierr;
1227674ae819SStefano Zampini 
1228674ae819SStefano Zampini   PetscFunctionBegin;
1229854ce69bSBarry Smith   ierr = PetscNew(&new_graph);CHKERRQ(ierr);
1230674ae819SStefano Zampini   new_graph->custom_minimal_size = 1;
12314b2aedd3SStefano Zampini   new_graph->commsizelimit = 1;
1232674ae819SStefano Zampini   *graph = new_graph;
1233674ae819SStefano Zampini   PetscFunctionReturn(0);
1234674ae819SStefano Zampini }
1235