xref: /petsc/src/ksp/pc/impls/bddc/bddcgraph.c (revision a07ea27a0d0aa5839ea22e119d3de928ef75918f)
1 #include <petsc/private/petscimpl.h>
2 #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
3 #include <../src/ksp/pc/impls/bddc/bddcstructs.h>
4 
5 #undef __FUNCT__
6 #define __FUNCT__ "PCBDDCGraphGetDirichletDofs"
7 PetscErrorCode PCBDDCGraphGetDirichletDofs(PCBDDCGraph graph, IS* dirdofs)
8 {
9   PetscErrorCode ierr;
10 
11   PetscFunctionBegin;
12   if (graph->dirdofs) {
13     ierr = PetscObjectReference((PetscObject)graph->dirdofs);CHKERRQ(ierr);
14   } else if (graph->has_dirichlet) {
15     PetscInt i,size;
16     PetscInt *dirdofs_idxs;
17 
18     size = 0;
19     for (i=0;i<graph->nvtxs;i++) {
20       if (graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) size++;
21     }
22 
23     ierr = PetscMalloc1(size,&dirdofs_idxs);CHKERRQ(ierr);
24     size = 0;
25     for (i=0;i<graph->nvtxs;i++) {
26       if (graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) dirdofs_idxs[size++] = i;
27     }
28     ierr = ISCreateGeneral(PetscObjectComm((PetscObject)graph->l2gmap),size,dirdofs_idxs,PETSC_OWN_POINTER,&graph->dirdofs);CHKERRQ(ierr);
29     ierr = PetscObjectReference((PetscObject)graph->dirdofs);CHKERRQ(ierr);
30   }
31   *dirdofs = graph->dirdofs;
32   PetscFunctionReturn(0);
33 }
34 
35 #undef __FUNCT__
36 #define __FUNCT__ "PCBDDCGraphASCIIView"
37 PetscErrorCode PCBDDCGraphASCIIView(PCBDDCGraph graph, PetscInt verbosity_level, PetscViewer viewer)
38 {
39   PetscInt       i,j,tabs;
40   PetscInt*      queue_in_global_numbering;
41   PetscErrorCode ierr;
42 
43   PetscFunctionBegin;
44   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
45   ierr = PetscViewerASCIIGetTab(viewer,&tabs);CHKERRQ(ierr);
46   ierr = PetscViewerASCIIPrintf(viewer,"--------------------------------------------------\n");CHKERRQ(ierr);
47   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
48   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Local BDDC graph for subdomain %04d\n",PetscGlobalRank);CHKERRQ(ierr);
49   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Number of vertices %d\n",graph->nvtxs);CHKERRQ(ierr);
50   if (verbosity_level > 1) {
51     for (i=0;i<graph->nvtxs;i++) {
52       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"%d:\n",i);CHKERRQ(ierr);
53       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   which_dof: %d\n",graph->which_dof[i]);CHKERRQ(ierr);
54       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   special_dof: %d\n",graph->special_dof[i]);CHKERRQ(ierr);
55       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   neighbours: %d\n",graph->count[i]);CHKERRQ(ierr);
56       ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
57       if (graph->count[i]) {
58         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"     set of neighbours:");CHKERRQ(ierr);
59         for (j=0;j<graph->count[i];j++) {
60           ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d",graph->neighbours_set[i][j]);CHKERRQ(ierr);
61         }
62         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
63       }
64       ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr);
65       ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
66       if (graph->mirrors) {
67         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   mirrors: %d\n",graph->mirrors[i]);CHKERRQ(ierr);
68         if (graph->mirrors[i]) {
69           ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
70           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"     set of mirrors:");CHKERRQ(ierr);
71           for (j=0;j<graph->mirrors[i];j++) {
72             ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d",graph->mirrors_set[i][j]);CHKERRQ(ierr);
73           }
74           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
75           ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr);
76           ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
77         }
78       }
79       if (verbosity_level > 2) {
80         if (graph->xadj && graph->adjncy) {
81           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   local adj list:");CHKERRQ(ierr);
82           ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
83           for (j=graph->xadj[i];j<graph->xadj[i+1];j++) {
84             ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d",graph->adjncy[j]);CHKERRQ(ierr);
85           }
86           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
87           ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr);
88           ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
89         } else {
90           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   no adj info\n");CHKERRQ(ierr);
91         }
92       }
93       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   interface subset id: %d\n",graph->subset[i]);CHKERRQ(ierr);
94       if (graph->subset[i] && graph->subset_ncc) {
95         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   ncc for subset: %d\n",graph->subset_ncc[graph->subset[i]-1]);CHKERRQ(ierr);
96       }
97     }
98   }
99   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Total number of connected components %d\n",graph->ncc);CHKERRQ(ierr);
100   ierr = PetscMalloc1(graph->cptr[graph->ncc],&queue_in_global_numbering);CHKERRQ(ierr);
101   ierr = ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_in_global_numbering);CHKERRQ(ierr);
102   for (i=0;i<graph->ncc;i++) {
103     PetscInt node_num=graph->queue[graph->cptr[i]];
104     PetscBool printcc = PETSC_FALSE;
105     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  %d (neighs:",i);CHKERRQ(ierr);
106     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
107     for (j=0;j<graph->count[node_num];j++) {
108       ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d",graph->neighbours_set[node_num][j]);CHKERRQ(ierr);
109     }
110     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"):");CHKERRQ(ierr);
111     if (verbosity_level > 1) {
112       printcc = PETSC_TRUE;
113     } else if (graph->count[node_num] > 1 || (graph->count[node_num] == 1 && graph->special_dof[node_num] == PCBDDCGRAPH_NEUMANN_MARK)) {
114       printcc = PETSC_TRUE;
115     }
116     if (printcc) {
117       for (j=graph->cptr[i];j<graph->cptr[i+1];j++) {
118         ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d (%d)",graph->queue[j],queue_in_global_numbering[j]);CHKERRQ(ierr);
119       }
120     }
121     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
122     ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr);
123     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
124   }
125   ierr = PetscFree(queue_in_global_numbering);CHKERRQ(ierr);
126   if (graph->custom_minimal_size > 1 && verbosity_level > 1) {
127     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Custom minimal size %d\n",graph->custom_minimal_size);CHKERRQ(ierr);
128   }
129   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
130   PetscFunctionReturn(0);
131 }
132 
133 #undef __FUNCT__
134 #define __FUNCT__ "PCBDDCGraphGetCandidatesIS"
135 PetscErrorCode PCBDDCGraphGetCandidatesIS(PCBDDCGraph graph, PetscInt *n_faces, IS *FacesIS[], PetscInt *n_edges, IS *EdgesIS[], IS *VerticesIS)
136 {
137   IS             *ISForFaces,*ISForEdges,ISForVertices;
138   PetscInt       i,nfc,nec,nvc,*idx;
139   PetscErrorCode ierr;
140 
141   PetscFunctionBegin;
142   /* loop on ccs to evalute number of faces, edges and vertices */
143   nfc = 0;
144   nec = 0;
145   nvc = 0;
146   for (i=0;i<graph->ncc;i++) {
147     PetscInt repdof = graph->queue[graph->cptr[i]];
148     if (graph->cptr[i+1]-graph->cptr[i] > graph->custom_minimal_size) {
149       if (graph->count[repdof] == 1 && graph->special_dof[repdof] != PCBDDCGRAPH_NEUMANN_MARK) {
150         nfc++;
151       } else { /* note that nec will be zero in 2d */
152         nec++;
153       }
154     } else {
155       nvc += graph->cptr[i+1]-graph->cptr[i];
156     }
157   }
158 
159   /* check if we are in 2D or 3D */
160   if (graph->twodim) { /* we are in a 2D case -> edges are shared by 2 subregions and faces don't exist */
161     nec = nfc;
162     nfc = 0;
163   }
164 
165   /* allocate IS arrays for faces, edges. Vertices need a single index set. */
166   if (FacesIS) {
167     ierr = PetscMalloc1(nfc,&ISForFaces);CHKERRQ(ierr);
168   }
169   if (EdgesIS) {
170     ierr = PetscMalloc1(nec,&ISForEdges);CHKERRQ(ierr);
171   }
172   if (VerticesIS) {
173     ierr = PetscMalloc1(nvc,&idx);CHKERRQ(ierr);
174   }
175 
176   /* loop on ccs to compute index sets for faces and edges */
177   if (!graph->queue_sorted) {
178     PetscInt *queue_global;
179 
180     ierr = PetscMalloc1(graph->cptr[graph->ncc],&queue_global);CHKERRQ(ierr);
181     ierr = ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_global);CHKERRQ(ierr);
182     for (i=0;i<graph->ncc;i++) {
183       ierr = PetscSortIntWithArray(graph->cptr[i+1]-graph->cptr[i],&queue_global[graph->cptr[i]],&graph->queue[graph->cptr[i]]);CHKERRQ(ierr);
184     }
185     ierr = PetscFree(queue_global);CHKERRQ(ierr);
186     graph->queue_sorted = PETSC_TRUE;
187   }
188   nfc = 0;
189   nec = 0;
190   for (i=0;i<graph->ncc;i++) {
191     PetscInt repdof = graph->queue[graph->cptr[i]];
192     if (graph->cptr[i+1]-graph->cptr[i] > graph->custom_minimal_size) {
193       if (graph->count[repdof] == 1 && graph->special_dof[repdof] != PCBDDCGRAPH_NEUMANN_MARK) {
194         if (graph->twodim) {
195           if (EdgesIS) {
196             ierr = ISCreateGeneral(PETSC_COMM_SELF,graph->cptr[i+1]-graph->cptr[i],&graph->queue[graph->cptr[i]],PETSC_COPY_VALUES,&ISForEdges[nec]);CHKERRQ(ierr);
197           }
198           nec++;
199         } else {
200           if (FacesIS) {
201             ierr = ISCreateGeneral(PETSC_COMM_SELF,graph->cptr[i+1]-graph->cptr[i],&graph->queue[graph->cptr[i]],PETSC_COPY_VALUES,&ISForFaces[nfc]);CHKERRQ(ierr);
202           }
203           nfc++;
204         }
205       } else {
206         if (EdgesIS) {
207           ierr = ISCreateGeneral(PETSC_COMM_SELF,graph->cptr[i+1]-graph->cptr[i],&graph->queue[graph->cptr[i]],PETSC_COPY_VALUES,&ISForEdges[nec]);CHKERRQ(ierr);
208         }
209         nec++;
210       }
211     }
212   }
213   /* index set for vertices */
214   if (VerticesIS) {
215     nvc = 0;
216     for (i=0;i<graph->ncc;i++) {
217       if (graph->cptr[i+1]-graph->cptr[i] <= graph->custom_minimal_size) {
218         PetscInt j;
219 
220         for (j=graph->cptr[i];j<graph->cptr[i+1];j++) {
221           idx[nvc]=graph->queue[j];
222           nvc++;
223         }
224       }
225     }
226     /* sort vertex set (by local ordering) */
227     ierr = PetscSortInt(nvc,idx);CHKERRQ(ierr);
228     ierr = ISCreateGeneral(PETSC_COMM_SELF,nvc,idx,PETSC_OWN_POINTER,&ISForVertices);CHKERRQ(ierr);
229   }
230   /* get back info */
231   if (n_faces) *n_faces = nfc;
232   if (FacesIS) {
233     *FacesIS = ISForFaces;
234   }
235   if (n_edges) *n_edges = nec;
236   if (EdgesIS) {
237     *EdgesIS = ISForEdges;
238   }
239   if (VerticesIS) {
240     *VerticesIS = ISForVertices;
241   }
242   PetscFunctionReturn(0);
243 }
244 
245 #undef __FUNCT__
246 #define __FUNCT__ "PCBDDCGraphComputeConnectedComponents"
247 PetscErrorCode PCBDDCGraphComputeConnectedComponents(PCBDDCGraph graph)
248 {
249   PetscBool      adapt_interface_reduced;
250   MPI_Comm       interface_comm;
251   PetscMPIInt    size;
252   PetscInt       i;
253   PetscBool      twodim;
254   PetscErrorCode ierr;
255 
256   PetscFunctionBegin;
257   /* compute connected components locally */
258   ierr = PetscObjectGetComm((PetscObject)(graph->l2gmap),&interface_comm);CHKERRQ(ierr);
259   ierr = PCBDDCGraphComputeConnectedComponentsLocal(graph);CHKERRQ(ierr);
260   /* check consistency of connected components among neighbouring subdomains -> it adapt them in case it is needed */
261   ierr = MPI_Comm_size(interface_comm,&size);CHKERRQ(ierr);
262   adapt_interface_reduced = PETSC_FALSE;
263   if (size > 1) {
264     PetscInt i;
265     PetscBool adapt_interface = PETSC_FALSE;
266     for (i=0;i<graph->n_subsets;i++) {
267       /* We are not sure that on a given subset of the local interface,
268          with two connected components, the latters be the same among sharing subdomains */
269       if (graph->subset_ncc[i] > 1) {
270         adapt_interface = PETSC_TRUE;
271         break;
272       }
273     }
274     ierr = MPI_Allreduce(&adapt_interface,&adapt_interface_reduced,1,MPIU_BOOL,MPI_LOR,interface_comm);CHKERRQ(ierr);
275   }
276 
277   if (graph->n_subsets && adapt_interface_reduced) {
278     /* old csr stuff */
279     PetscInt    *aux_new_xadj,*new_xadj,*new_adjncy;
280     PetscInt    *old_xadj,*old_adjncy;
281     PetscBT     subset_cc_adapt;
282     /* MPI */
283     MPI_Request *send_requests,*recv_requests;
284     PetscInt    *send_buffer,*recv_buffer;
285     PetscInt    sum_requests,start_of_recv,start_of_send;
286     PetscInt    *cum_recv_counts;
287     /* others */
288     PetscInt    *labels;
289     PetscInt    global_subset_counter;
290     PetscInt    j,k,s;
291 
292     /* Retrict adjacency graph using information from previously computed connected components */
293     ierr = PetscMalloc1(graph->nvtxs,&aux_new_xadj);CHKERRQ(ierr);
294     for (i=0;i<graph->nvtxs;i++) {
295       aux_new_xadj[i]=1;
296     }
297     for (i=0;i<graph->ncc;i++) {
298       k = graph->cptr[i+1]-graph->cptr[i];
299       for (j=0;j<k;j++) {
300         aux_new_xadj[graph->queue[graph->cptr[i]+j]]=k;
301       }
302     }
303     j = 0;
304     for (i=0;i<graph->nvtxs;i++) {
305       j += aux_new_xadj[i];
306     }
307     ierr = PetscMalloc1(graph->nvtxs+1,&new_xadj);CHKERRQ(ierr);
308     ierr = PetscMalloc1(j,&new_adjncy);CHKERRQ(ierr);
309     new_xadj[0]=0;
310     for (i=0;i<graph->nvtxs;i++) {
311       new_xadj[i+1]=new_xadj[i]+aux_new_xadj[i];
312       if (aux_new_xadj[i]==1) {
313         new_adjncy[new_xadj[i]]=i;
314       }
315     }
316     ierr = PetscFree(aux_new_xadj);CHKERRQ(ierr);
317     ierr = PetscMalloc1(graph->nvtxs,&labels);CHKERRQ(ierr);
318     ierr = PetscMemzero(labels,graph->nvtxs*sizeof(*labels));CHKERRQ(ierr);
319     for (i=0;i<graph->ncc;i++) {
320       k = graph->cptr[i+1]-graph->cptr[i];
321       for (j=0;j<k;j++) {
322         ierr = PetscMemcpy(&new_adjncy[new_xadj[graph->queue[graph->cptr[i]+j]]],&graph->queue[graph->cptr[i]],k*sizeof(PetscInt));CHKERRQ(ierr);
323         labels[graph->queue[graph->cptr[i]+j]] = i;
324       }
325     }
326     /* set temporarly new CSR into graph */
327     old_xadj = graph->xadj;
328     old_adjncy = graph->adjncy;
329     graph->xadj = new_xadj;
330     graph->adjncy = new_adjncy;
331     /* allocate some space */
332     ierr = PetscMalloc1(graph->n_subsets+1,&cum_recv_counts);CHKERRQ(ierr);
333     ierr = PetscMemzero(cum_recv_counts,(graph->n_subsets+1)*sizeof(*cum_recv_counts));CHKERRQ(ierr);
334     /* first count how many neighbours per connected component I will receive from */
335     cum_recv_counts[0]=0;
336     for (i=0;i<graph->n_subsets;i++) {
337       cum_recv_counts[i+1]=cum_recv_counts[i]+graph->count[graph->subsets[i][0]];
338     }
339     ierr = PetscMalloc1(cum_recv_counts[graph->n_subsets],&recv_buffer);CHKERRQ(ierr);
340     ierr = PetscMalloc2(cum_recv_counts[graph->n_subsets],&send_requests,cum_recv_counts[graph->n_subsets],&recv_requests);CHKERRQ(ierr);
341     for (i=0;i<cum_recv_counts[graph->n_subsets];i++) {
342       send_requests[i]=MPI_REQUEST_NULL;
343       recv_requests[i]=MPI_REQUEST_NULL;
344     }
345     /* exchange with my neighbours the number of my connected components on the subset of interface */
346     sum_requests = 0;
347     for (i=0;i<graph->n_subsets;i++) {
348       PetscMPIInt neigh,tag;
349 
350       j = graph->subsets[i][0];
351       ierr = PetscMPIIntCast(2*graph->subset_ref_node[i],&tag);CHKERRQ(ierr);
352       for (k=0;k<graph->count[j];k++) {
353         ierr = PetscMPIIntCast(graph->neighbours_set[j][k],&neigh);CHKERRQ(ierr);
354         ierr = MPI_Isend(&graph->subset_ncc[i],1,MPIU_INT,neigh,tag,interface_comm,&send_requests[sum_requests]);CHKERRQ(ierr);
355         ierr = MPI_Irecv(&recv_buffer[sum_requests],1,MPIU_INT,neigh,tag,interface_comm,&recv_requests[sum_requests]);CHKERRQ(ierr);
356         sum_requests++;
357       }
358     }
359     ierr = MPI_Waitall(sum_requests,recv_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
360     ierr = MPI_Waitall(sum_requests,send_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
361     /* determine the subsets I have to adapt */
362     ierr = PetscBTCreate(graph->n_subsets,&subset_cc_adapt);CHKERRQ(ierr);
363     ierr = PetscBTMemzero(graph->n_subsets,subset_cc_adapt);CHKERRQ(ierr);
364     for (i=0;i<graph->n_subsets;i++) {
365       for (j=cum_recv_counts[i];j<cum_recv_counts[i+1];j++){
366         /* adapt if more than one cc is present */
367          if (graph->subset_ncc[i] > 1 || recv_buffer[j] > 1) {
368           ierr = PetscBTSet(subset_cc_adapt,i);
369           break;
370         }
371       }
372     }
373     ierr = PetscFree(recv_buffer);CHKERRQ(ierr);
374     /* determine send/recv buffers sizes */
375     j = 0;
376     for (i=0;i<graph->n_subsets;i++) {
377       if (PetscBTLookup(subset_cc_adapt,i)) {
378         j += graph->subsets_size[i];
379       }
380     }
381     k = 0;
382     for (i=0;i<graph->n_subsets;i++) {
383       if (PetscBTLookup(subset_cc_adapt,i)) {
384         k += (cum_recv_counts[i+1]-cum_recv_counts[i])*graph->subsets_size[i];
385       }
386     }
387     ierr = PetscMalloc2(j,&send_buffer,k,&recv_buffer);CHKERRQ(ierr);
388 
389     /* fill send buffer */
390     j = 0;
391     for (i=0;i<graph->n_subsets;i++) {
392       if (PetscBTLookup(subset_cc_adapt,i)) {
393         for (k=0;k<graph->subsets_size[i];k++) {
394           send_buffer[j++] = labels[graph->subsets[i][k]];
395         }
396       }
397     }
398     ierr = PetscFree(labels);CHKERRQ(ierr);
399 
400     /* now exchange the data */
401     start_of_recv = 0;
402     start_of_send = 0;
403     sum_requests = 0;
404     for (i=0;i<graph->n_subsets;i++) {
405       if (PetscBTLookup(subset_cc_adapt,i)) {
406         PetscMPIInt neigh,tag;
407         PetscInt    size_of_send = graph->subsets_size[i];
408 
409         j = graph->subsets[i][0];
410         ierr = PetscMPIIntCast(2*graph->subset_ref_node[i]+1,&tag);CHKERRQ(ierr);
411         for (k=0;k<graph->count[j];k++) {
412           ierr = PetscMPIIntCast(graph->neighbours_set[j][k],&neigh);CHKERRQ(ierr);
413           ierr = MPI_Isend(&send_buffer[start_of_send],size_of_send,MPIU_INT,neigh,tag,interface_comm,&send_requests[sum_requests]);CHKERRQ(ierr);
414           ierr = MPI_Irecv(&recv_buffer[start_of_recv],size_of_send,MPIU_INT,neigh,tag,interface_comm,&recv_requests[sum_requests]);CHKERRQ(ierr);
415           start_of_recv += size_of_send;
416           sum_requests++;
417         }
418         start_of_send += size_of_send;
419       }
420     }
421     ierr = MPI_Waitall(sum_requests,recv_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
422     ierr = MPI_Waitall(sum_requests,send_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
423     ierr = PetscFree2(send_requests,recv_requests);CHKERRQ(ierr);
424 
425     start_of_recv = 0;
426     global_subset_counter = 0;
427     for (i=0;i<graph->n_subsets;i++) {
428       if (PetscBTLookup(subset_cc_adapt,i)) {
429         PetscInt **temp_buffer,temp_buffer_size,*private_label;
430 
431         /* allocate some temporary space */
432         temp_buffer_size = graph->subsets_size[i];
433         ierr = PetscMalloc1(temp_buffer_size,&temp_buffer);CHKERRQ(ierr);
434         ierr = PetscMalloc1(temp_buffer_size*(cum_recv_counts[i+1]-cum_recv_counts[i]),&temp_buffer[0]);CHKERRQ(ierr);
435         ierr = PetscMemzero(temp_buffer[0],temp_buffer_size*(cum_recv_counts[i+1]-cum_recv_counts[i])*sizeof(PetscInt));CHKERRQ(ierr);
436         for (j=1;j<temp_buffer_size;j++) {
437           temp_buffer[j] = temp_buffer[j-1]+cum_recv_counts[i+1]-cum_recv_counts[i];
438         }
439         /* analyze contributions from neighbouring subdomains for i-th subset
440            temp buffer structure:
441            supposing part of the interface has dimension 5
442            e.g with global dofs 0,1,2,3,4, locally ordered on the current process as 0,4,3,1,2
443            3 neighs procs having connected components:
444              neigh 0: [0 1 4], [2 3], labels [4,7]  (2 connected components)
445              neigh 1: [0 1], [2 3 4], labels [3 2]  (2 connected components)
446              neigh 2: [0 4], [1], [2 3], labels [1 5 6] (3 connected components)
447            tempbuffer (row-oriented) will be filled as:
448              [ 4, 3, 1;
449                4, 2, 1;
450                7, 2, 6;
451                4, 3, 5;
452                7, 2, 6; ];
453            This way we can simply find intersections of ccs among neighs.
454            The graph->subset array will need to be modified. The output for the example is [0], [1], [2 3], [4];
455                                                                                                                                    */
456         for (j=0;j<cum_recv_counts[i+1]-cum_recv_counts[i];j++) {
457           for (k=0;k<temp_buffer_size;k++) {
458             temp_buffer[k][j] = recv_buffer[start_of_recv+k];
459           }
460           start_of_recv += temp_buffer_size;
461         }
462         ierr = PetscMalloc1(temp_buffer_size,&private_label);CHKERRQ(ierr);
463         ierr = PetscMemzero(private_label,temp_buffer_size*sizeof(*private_label));CHKERRQ(ierr);
464         for (j=0;j<temp_buffer_size;j++) {
465           if (!private_label[j]) { /* found a new cc  */
466             PetscBool same_set;
467 
468             global_subset_counter++;
469             private_label[j] = global_subset_counter;
470             for (k=j+1;k<temp_buffer_size;k++) { /* check for other nodes in new cc */
471               same_set = PETSC_TRUE;
472               for (s=0;s<cum_recv_counts[i+1]-cum_recv_counts[i];s++) {
473                 if (temp_buffer[j][s] != temp_buffer[k][s]) {
474                   same_set = PETSC_FALSE;
475                   break;
476                 }
477               }
478               if (same_set) {
479                 private_label[k] = global_subset_counter;
480               }
481             }
482           }
483         }
484         /* insert private labels in graph structure */
485         for (j=0;j<graph->subsets_size[i];j++) {
486           graph->subset[graph->subsets[i][j]] = graph->n_subsets+private_label[j];
487         }
488         ierr = PetscFree(private_label);CHKERRQ(ierr);
489         ierr = PetscFree(temp_buffer[0]);CHKERRQ(ierr);
490         ierr = PetscFree(temp_buffer);CHKERRQ(ierr);
491       }
492     }
493     ierr = PetscFree2(send_buffer,recv_buffer);CHKERRQ(ierr);
494     ierr = PetscFree(cum_recv_counts);CHKERRQ(ierr);
495     ierr = PetscBTDestroy(&subset_cc_adapt);CHKERRQ(ierr);
496     /* We are ready to find for connected components which are consistent among neighbouring subdomains */
497     if (global_subset_counter) {
498       ierr = PetscBTMemzero(graph->nvtxs,graph->touched);CHKERRQ(ierr);
499       global_subset_counter = 0;
500       for (i=0;i<graph->nvtxs;i++) {
501         if (graph->subset[i] && !PetscBTLookup(graph->touched,i)) {
502           global_subset_counter++;
503           for (j=i+1;j<graph->nvtxs;j++) {
504             if (!PetscBTLookup(graph->touched,j) && graph->subset[j]==graph->subset[i]) {
505               graph->subset[j] = global_subset_counter;
506               ierr = PetscBTSet(graph->touched,j);CHKERRQ(ierr);
507             }
508           }
509           graph->subset[i] = global_subset_counter;
510           ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
511         }
512       }
513       /* refine connected components locally */
514       ierr = PCBDDCGraphComputeConnectedComponentsLocal(graph);CHKERRQ(ierr);
515     }
516     /* restore original CSR graph of dofs */
517     ierr = PetscFree(new_xadj);CHKERRQ(ierr);
518     ierr = PetscFree(new_adjncy);CHKERRQ(ierr);
519     graph->xadj = old_xadj;
520     graph->adjncy = old_adjncy;
521   }
522 
523   /* Determine if we are in 2D or 3D */
524   twodim  = PETSC_TRUE;
525   for (i=0;i<graph->ncc;i++) {
526     PetscInt repdof = graph->queue[graph->cptr[i]];
527     if (graph->cptr[i+1]-graph->cptr[i] > graph->custom_minimal_size) {
528       if (graph->count[repdof] > 1 || graph->special_dof[repdof] == PCBDDCGRAPH_NEUMANN_MARK) {
529         twodim = PETSC_FALSE;
530         break;
531       }
532     }
533   }
534   ierr = MPI_Allreduce(&twodim,&graph->twodim,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)graph->l2gmap));CHKERRQ(ierr);
535   PetscFunctionReturn(0);
536 }
537 
538 /* The following code has been adapted from function IsConnectedSubdomain contained
539    in source file contig.c of METIS library (version 5.0.1)
540    It finds connected components for each subset  */
541 #undef __FUNCT__
542 #define __FUNCT__ "PCBDDCGraphComputeConnectedComponentsLocal"
543 PetscErrorCode PCBDDCGraphComputeConnectedComponentsLocal(PCBDDCGraph graph)
544 {
545   PetscInt       i,j,k,first,last,nleft,ncc,pid,cum_queue,n,ncc_pid;
546   PetscMPIInt    size;
547   PetscErrorCode ierr;
548 
549   PetscFunctionBegin;
550   /* quiet return if no csr info is available */
551   if (!graph->xadj || !graph->adjncy) {
552     PetscFunctionReturn(0);
553   }
554 
555   /* reset any previous search of connected components */
556   ierr = PetscBTMemzero(graph->nvtxs,graph->touched);CHKERRQ(ierr);
557   graph->n_subsets = 0;
558   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)graph->l2gmap),&size);CHKERRQ(ierr);
559   if (size == 1) {
560     if (graph->nvtxs) {
561       graph->n_subsets = 1;
562       for (i=0;i<graph->nvtxs;i++) {
563         graph->subset[i] = 1;
564       }
565     }
566   } else {
567     for (i=0;i<graph->nvtxs;i++) {
568       if (graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK || !graph->count[i]) {
569         ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
570         graph->subset[i] = 0;
571       }
572       graph->n_subsets = PetscMax(graph->n_subsets,graph->subset[i]);
573     }
574   }
575   ierr = PetscFree(graph->subset_ncc);CHKERRQ(ierr);
576   ierr = PetscMalloc1(graph->n_subsets,&graph->subset_ncc);CHKERRQ(ierr);
577 
578   /* begin search for connected components */
579   cum_queue = 0;
580   ncc = 0;
581   for (n=0;n<graph->n_subsets;n++) {
582     pid = n+1;  /* partition labeled by 0 is discarded */
583     nleft = 0;
584     for (i=0;i<graph->nvtxs;i++) {
585       if (graph->subset[i] == pid) {
586         nleft++;
587       }
588     }
589     for (i=0; i<graph->nvtxs; i++) {
590       if (graph->subset[i] == pid) {
591         break;
592       }
593     }
594     ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
595     graph->queue[cum_queue] = i;
596     first = 0;
597     last = 1;
598     graph->cptr[ncc] = cum_queue;
599     ncc_pid = 0;
600     while (first != nleft) {
601       if (first == last) {
602         graph->cptr[++ncc] = first+cum_queue;
603         ncc_pid++;
604         for (i=0; i<graph->nvtxs; i++) { /* TODO-> use a while! */
605           if (graph->subset[i] == pid && !PetscBTLookup(graph->touched,i)) {
606             break;
607           }
608         }
609         graph->queue[cum_queue+last] = i;
610         last++;
611         ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
612       }
613       i = graph->queue[cum_queue+first];
614       first++;
615       for (j=graph->xadj[i];j<graph->xadj[i+1];j++) {
616         k = graph->adjncy[j];
617         if (graph->subset[k] == pid && !PetscBTLookup(graph->touched,k)) {
618           graph->queue[cum_queue+last] = k;
619           last++;
620           ierr = PetscBTSet(graph->touched,k);CHKERRQ(ierr);
621         }
622       }
623     }
624     graph->cptr[++ncc] = first+cum_queue;
625     ncc_pid++;
626     cum_queue = graph->cptr[ncc];
627     graph->subset_ncc[n] = ncc_pid;
628   }
629   graph->ncc = ncc;
630   graph->queue_sorted = PETSC_FALSE;
631   PetscFunctionReturn(0);
632 }
633 
634 #undef __FUNCT__
635 #define __FUNCT__ "PCBDDCGraphSetUp"
636 PetscErrorCode PCBDDCGraphSetUp(PCBDDCGraph graph, PetscInt custom_minimal_size, IS neumann_is, IS dirichlet_is, PetscInt n_ISForDofs, IS ISForDofs[], IS custom_primal_vertices)
637 {
638   VecScatter     scatter_ctx;
639   Vec            local_vec,local_vec2,global_vec;
640   IS             to,from;
641   MPI_Comm       comm;
642   PetscScalar    *array,*array2;
643   const PetscInt *is_indices;
644   PetscInt       n_neigh,*neigh,*n_shared,**shared,*queue_global,*subset_ref_node_global;
645   PetscInt       i,j,k,s,total_counts,nodes_touched,is_size;
646   PetscMPIInt    size;
647   PetscBool      same_set,mirrors_found;
648   PetscErrorCode ierr;
649 
650   PetscFunctionBegin;
651   graph->has_dirichlet = PETSC_FALSE;
652   if (dirichlet_is) {
653     PetscCheckSameComm(graph->l2gmap,1,dirichlet_is,4);
654     graph->has_dirichlet = PETSC_TRUE;
655   }
656   ierr = PetscObjectGetComm((PetscObject)(graph->l2gmap),&comm);CHKERRQ(ierr);
657   /* custom_minimal_size */
658   graph->custom_minimal_size = PetscMax(graph->custom_minimal_size,custom_minimal_size);
659   /* get info l2gmap and allocate work vectors  */
660   ierr = ISLocalToGlobalMappingGetInfo(graph->l2gmap,&n_neigh,&neigh,&n_shared,&shared);CHKERRQ(ierr);
661   /* check if we have any local periodic nodes (periodic BCs) */
662   mirrors_found = PETSC_FALSE;
663   if (graph->nvtxs && n_neigh) {
664     for (i=0; i<n_shared[0]; i++) graph->count[shared[0][i]] += 1;
665     for (i=0; i<n_shared[0]; i++) {
666       if (graph->count[shared[0][i]] > 1) {
667         mirrors_found = PETSC_TRUE;
668         break;
669       }
670     }
671   }
672   /* create some workspace objects */
673   local_vec = NULL;
674   local_vec2 = NULL;
675   global_vec = NULL;
676   to = NULL;
677   from = NULL;
678   scatter_ctx = NULL;
679   if (n_ISForDofs || dirichlet_is || neumann_is || custom_primal_vertices) {
680     ierr = VecCreate(PETSC_COMM_SELF,&local_vec);CHKERRQ(ierr);
681     ierr = VecSetSizes(local_vec,PETSC_DECIDE,graph->nvtxs);CHKERRQ(ierr);
682     ierr = VecSetType(local_vec,VECSTANDARD);CHKERRQ(ierr);
683     ierr = VecDuplicate(local_vec,&local_vec2);CHKERRQ(ierr);
684     ierr = VecCreate(comm,&global_vec);CHKERRQ(ierr);
685     ierr = VecSetSizes(global_vec,PETSC_DECIDE,graph->nvtxs_global);CHKERRQ(ierr);
686     ierr = VecSetType(global_vec,VECSTANDARD);CHKERRQ(ierr);
687     ierr = ISCreateStride(PETSC_COMM_SELF,graph->nvtxs,0,1,&to);CHKERRQ(ierr);
688     ierr = ISLocalToGlobalMappingApplyIS(graph->l2gmap,to,&from);CHKERRQ(ierr);
689     ierr = VecScatterCreate(global_vec,from,local_vec,to,&scatter_ctx);CHKERRQ(ierr);
690   } else if (mirrors_found) {
691     ierr = ISCreateStride(PETSC_COMM_SELF,graph->nvtxs,0,1,&to);CHKERRQ(ierr);
692     ierr = ISLocalToGlobalMappingApplyIS(graph->l2gmap,to,&from);CHKERRQ(ierr);
693   }
694   /* compute local mirrors (if any) */
695   if (mirrors_found) {
696     PetscInt *local_indices,*global_indices;
697     /* get arrays of local and global indices */
698     ierr = PetscMalloc1(graph->nvtxs,&local_indices);CHKERRQ(ierr);
699     ierr = ISGetIndices(to,(const PetscInt**)&is_indices);CHKERRQ(ierr);
700     ierr = PetscMemcpy(local_indices,is_indices,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
701     ierr = ISRestoreIndices(to,(const PetscInt**)&is_indices);CHKERRQ(ierr);
702     ierr = PetscMalloc1(graph->nvtxs,&global_indices);CHKERRQ(ierr);
703     ierr = ISGetIndices(from,(const PetscInt**)&is_indices);CHKERRQ(ierr);
704     ierr = PetscMemcpy(global_indices,is_indices,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
705     ierr = ISRestoreIndices(from,(const PetscInt**)&is_indices);CHKERRQ(ierr);
706     /* allocate space for mirrors */
707     ierr = PetscMalloc2(graph->nvtxs,&graph->mirrors,graph->nvtxs,&graph->mirrors_set);CHKERRQ(ierr);
708     ierr = PetscMemzero(graph->mirrors,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
709     graph->mirrors_set[0] = 0;
710 
711     k=0;
712     for (i=0;i<n_shared[0];i++) {
713       j=shared[0][i];
714       if (graph->count[j] > 1) {
715         graph->mirrors[j]++;
716         k++;
717       }
718     }
719     /* allocate space for set of mirrors */
720     ierr = PetscMalloc1(k,&graph->mirrors_set[0]);CHKERRQ(ierr);
721     for (i=1;i<graph->nvtxs;i++)
722       graph->mirrors_set[i]=graph->mirrors_set[i-1]+graph->mirrors[i-1];
723 
724     /* fill arrays */
725     ierr = PetscMemzero(graph->mirrors,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
726     for (j=0;j<n_shared[0];j++) {
727       i=shared[0][j];
728       if (graph->count[i] > 1)
729         graph->mirrors_set[i][graph->mirrors[i]++]=global_indices[i];
730     }
731     ierr = PetscSortIntWithArray(graph->nvtxs,global_indices,local_indices);CHKERRQ(ierr);
732     for (i=0;i<graph->nvtxs;i++) {
733       if (graph->mirrors[i] > 0) {
734         ierr = PetscFindInt(graph->mirrors_set[i][0],graph->nvtxs,global_indices,&k);CHKERRQ(ierr);
735         j = global_indices[k];
736         while ( k > 0 && global_indices[k-1] == j) k--;
737         for (j=0;j<graph->mirrors[i];j++) {
738           graph->mirrors_set[i][j]=local_indices[k+j];
739         }
740         ierr = PetscSortInt(graph->mirrors[i],graph->mirrors_set[i]);CHKERRQ(ierr);
741       }
742     }
743     ierr = PetscFree(local_indices);CHKERRQ(ierr);
744     ierr = PetscFree(global_indices);CHKERRQ(ierr);
745   }
746   ierr = PetscMemzero(graph->count,graph->nvtxs*sizeof(*graph->count));CHKERRQ(ierr);
747   ierr = ISDestroy(&to);CHKERRQ(ierr);
748   ierr = ISDestroy(&from);CHKERRQ(ierr);
749 
750   /* Count total number of neigh per node */
751   k=0;
752   for (i=1;i<n_neigh;i++) {
753     k += n_shared[i];
754     for (j=0;j<n_shared[i];j++) {
755       graph->count[shared[i][j]] += 1;
756     }
757   }
758   /* Allocate space for storing the set of neighbours for each node */
759   if (graph->nvtxs) {
760     ierr = PetscMalloc1(k,&graph->neighbours_set[0]);CHKERRQ(ierr);
761   }
762   for (i=1;i<graph->nvtxs;i++) { /* dont count myself */
763     graph->neighbours_set[i]=graph->neighbours_set[i-1]+graph->count[i-1];
764   }
765   /* Get information for sharing subdomains */
766   ierr = PetscMemzero(graph->count,graph->nvtxs*sizeof(*graph->count));CHKERRQ(ierr);
767   for (i=1;i<n_neigh;i++) { /* dont count myself */
768     s = n_shared[i];
769     for (j=0;j<s;j++) {
770       k = shared[i][j];
771       graph->neighbours_set[k][graph->count[k]] = neigh[i];
772       graph->count[k] += 1;
773     }
774   }
775   /* sort set of sharing subdomains */
776   for (i=0;i<graph->nvtxs;i++) {
777     ierr = PetscSortRemoveDupsInt(&graph->count[i],graph->neighbours_set[i]);CHKERRQ(ierr);
778   }
779   /* free memory allocated by ISLocalToGlobalMappingGetInfo */
780   ierr = ISLocalToGlobalMappingRestoreInfo(graph->l2gmap,&n_neigh,&neigh,&n_shared,&shared);CHKERRQ(ierr);
781 
782   /*
783      Get info for dofs splitting
784      User can specify just a subset; an additional field is considered as a complementary field
785   */
786   for (i=0;i<graph->nvtxs;i++) {
787     graph->which_dof[i] = n_ISForDofs; /* by default a dof belongs to the complement set */
788   }
789   if (n_ISForDofs) {
790     ierr = VecSet(local_vec,-1.0);CHKERRQ(ierr);
791   }
792   for (i=0;i<n_ISForDofs;i++) {
793     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
794     ierr = ISGetLocalSize(ISForDofs[i],&is_size);CHKERRQ(ierr);
795     ierr = ISGetIndices(ISForDofs[i],(const PetscInt**)&is_indices);CHKERRQ(ierr);
796     for (j=0;j<is_size;j++) {
797       if (is_indices[j] > -1 && is_indices[j] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
798         graph->which_dof[is_indices[j]] = i;
799         array[is_indices[j]] = 1.*i;
800       }
801     }
802     ierr = ISRestoreIndices(ISForDofs[i],(const PetscInt**)&is_indices);CHKERRQ(ierr);
803     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
804   }
805   /* Check consistency among neighbours */
806   if (n_ISForDofs) {
807     ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
808     ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
809     ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
810     ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
811     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
812     ierr = VecGetArray(local_vec2,&array2);CHKERRQ(ierr);
813     for (i=0;i<graph->nvtxs;i++){
814       PetscInt field1,field2;
815 
816       field1 = (PetscInt)PetscRealPart(array[i]);
817       field2 = (PetscInt)PetscRealPart(array2[i]);
818       if (field1 != field2) {
819         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);
820       }
821     }
822     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
823     ierr = VecRestoreArray(local_vec2,&array2);CHKERRQ(ierr);
824   }
825   if (neumann_is || dirichlet_is) {
826     /* Take into account Neumann nodes */
827     ierr = VecSet(local_vec,0.0);CHKERRQ(ierr);
828     ierr = VecSet(local_vec2,0.0);CHKERRQ(ierr);
829     if (neumann_is) {
830       ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
831       ierr = ISGetLocalSize(neumann_is,&is_size);CHKERRQ(ierr);
832       ierr = ISGetIndices(neumann_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
833       for (i=0;i<is_size;i++) {
834         if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
835           array[is_indices[i]] = 1.0;
836         }
837       }
838       ierr = ISRestoreIndices(neumann_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
839       ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
840     }
841     /* Neumann nodes: impose consistency among neighbours */
842     ierr = VecSet(global_vec,0.0);CHKERRQ(ierr);
843     ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
844     ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
845     ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
846     ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
847     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
848     for (i=0;i<graph->nvtxs;i++) {
849       if (PetscRealPart(array[i]) > 0.1) {
850         graph->special_dof[i] = PCBDDCGRAPH_NEUMANN_MARK;
851       }
852     }
853     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
854     /* Take into account Dirichlet nodes */
855     ierr = VecSet(local_vec2,0.0);CHKERRQ(ierr);
856     if (dirichlet_is) {
857       ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
858       ierr = VecGetArray(local_vec2,&array2);CHKERRQ(ierr);
859       ierr = ISGetLocalSize(dirichlet_is,&is_size);CHKERRQ(ierr);
860       ierr = ISGetIndices(dirichlet_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
861       for (i=0;i<is_size;i++){
862         if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
863           k = is_indices[i];
864           if (graph->count[k] && !PetscBTLookup(graph->touched,k)) {
865             if (PetscRealPart(array[k]) > 0.1) {
866               SETERRQ1(comm,PETSC_ERR_USER,"BDDC cannot have boundary nodes which are marked Neumann and Dirichlet at the same time! Local node %d is wrong\n",k);
867             }
868             array2[k] = 1.0;
869           }
870         }
871       }
872       ierr = ISRestoreIndices(dirichlet_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
873       ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
874       ierr = VecRestoreArray(local_vec2,&array2);CHKERRQ(ierr);
875     }
876     /* Dirichlet nodes: impose consistency among neighbours */
877     ierr = VecSet(global_vec,0.0);CHKERRQ(ierr);
878     ierr = VecScatterBegin(scatter_ctx,local_vec2,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
879     ierr = VecScatterEnd(scatter_ctx,local_vec2,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
880     ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
881     ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
882     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
883     for (i=0;i<graph->nvtxs;i++) {
884       if (PetscRealPart(array[i]) > 0.1) {
885         ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
886         graph->subset[i] = 0; /* dirichlet nodes treated as internal -> is it ok? */
887         graph->special_dof[i] = PCBDDCGRAPH_DIRICHLET_MARK;
888       }
889     }
890     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
891   }
892   /* mark local periodic nodes (if any) and adapt CSR graph (if any) */
893   if (graph->mirrors) {
894     for (i=0;i<graph->nvtxs;i++)
895       if (graph->mirrors[i])
896         graph->special_dof[i] = PCBDDCGRAPH_LOCAL_PERIODIC_MARK;
897 
898     if (graph->xadj && graph->adjncy) {
899       PetscInt *new_xadj,*new_adjncy;
900       /* sort CSR graph */
901       for (i=0;i<graph->nvtxs;i++)
902         ierr = PetscSortInt(graph->xadj[i+1]-graph->xadj[i],&graph->adjncy[graph->xadj[i]]);CHKERRQ(ierr);
903 
904       /* adapt local CSR graph in case of local periodicity */
905       k=0;
906       for (i=0;i<graph->nvtxs;i++)
907         for (j=graph->xadj[i];j<graph->xadj[i+1];j++)
908           k += graph->mirrors[graph->adjncy[j]];
909 
910       ierr = PetscMalloc1(graph->nvtxs+1,&new_xadj);CHKERRQ(ierr);
911       ierr = PetscMalloc1(k+graph->xadj[graph->nvtxs],&new_adjncy);CHKERRQ(ierr);
912       new_xadj[0]=0;
913       for (i=0;i<graph->nvtxs;i++) {
914         k = graph->xadj[i+1]-graph->xadj[i];
915         ierr = PetscMemcpy(&new_adjncy[new_xadj[i]],&graph->adjncy[graph->xadj[i]],k*sizeof(PetscInt));CHKERRQ(ierr);
916         new_xadj[i+1]=new_xadj[i]+k;
917         for (j=graph->xadj[i];j<graph->xadj[i+1];j++) {
918           k = graph->mirrors[graph->adjncy[j]];
919           ierr = PetscMemcpy(&new_adjncy[new_xadj[i+1]],graph->mirrors_set[graph->adjncy[j]],k*sizeof(PetscInt));CHKERRQ(ierr);
920           new_xadj[i+1]+=k;
921         }
922         k = new_xadj[i+1]-new_xadj[i];
923         ierr = PetscSortRemoveDupsInt(&k,&new_adjncy[new_xadj[i]]);CHKERRQ(ierr);
924         new_xadj[i+1]=new_xadj[i]+k;
925       }
926       /* set new CSR into graph */
927       ierr = PetscFree(graph->xadj);CHKERRQ(ierr);
928       ierr = PetscFree(graph->adjncy);CHKERRQ(ierr);
929       graph->xadj = new_xadj;
930       graph->adjncy = new_adjncy;
931     }
932   }
933 
934   /* mark special nodes (if any) -> each will become a single node equivalence class */
935   if (custom_primal_vertices) {
936     ierr = VecSet(local_vec,0.0);CHKERRQ(ierr);
937     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
938     ierr = ISGetLocalSize(custom_primal_vertices,&is_size);CHKERRQ(ierr);
939     ierr = ISGetIndices(custom_primal_vertices,(const PetscInt**)&is_indices);CHKERRQ(ierr);
940     for (i=0;i<is_size;i++){
941       if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
942         array[is_indices[i]] = 1.0;
943       }
944     }
945     ierr = ISRestoreIndices(custom_primal_vertices,(const PetscInt**)&is_indices);CHKERRQ(ierr);
946     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
947     /* special nodes: impose consistency among neighbours */
948     ierr = VecSet(global_vec,0.0);CHKERRQ(ierr);
949     ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
950     ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
951     ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
952     ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
953     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
954     j = 0;
955     for (i=0;i<graph->nvtxs;i++) {
956       if (PetscRealPart(array[i]) > 0.1 && graph->special_dof[i] != PCBDDCGRAPH_DIRICHLET_MARK) {
957         graph->special_dof[i] = PCBDDCGRAPH_SPECIAL_MARK-j;
958         j++;
959       }
960     }
961     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
962   }
963 
964   /* mark interior nodes as touched and belonging to partition number 0 */
965   for (i=0;i<graph->nvtxs;i++) {
966     if (!graph->count[i]) {
967       ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
968       graph->subset[i] = 0;
969     }
970   }
971   /* init graph structure and compute default subsets */
972   nodes_touched=0;
973   for (i=0;i<graph->nvtxs;i++) {
974     if (PetscBTLookup(graph->touched,i)) {
975       nodes_touched++;
976     }
977   }
978   i = 0;
979   graph->ncc = 0;
980   total_counts = 0;
981 
982   /* allocated space for queues */
983   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
984   if (size == 1) {
985     ierr = PetscMalloc2(graph->nvtxs+1,&graph->cptr,graph->nvtxs,&graph->queue);CHKERRQ(ierr);
986   } else {
987     PetscInt nused = graph->nvtxs - nodes_touched;
988     ierr = PetscMalloc2(nused+1,&graph->cptr,nused,&graph->queue);CHKERRQ(ierr);
989   }
990 
991   while (nodes_touched<graph->nvtxs) {
992     /*  find first untouched node in local ordering */
993     while (PetscBTLookup(graph->touched,i)) {
994       i++;
995     }
996     ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
997     graph->subset[i] = graph->ncc+1;
998     graph->cptr[graph->ncc] = total_counts;
999     graph->queue[total_counts] = i;
1000     total_counts++;
1001     nodes_touched++;
1002     /* now find all other nodes having the same set of sharing subdomains */
1003     for (j=i+1;j<graph->nvtxs;j++) {
1004       /* check for same number of sharing subdomains, dof number and same special mark */
1005       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]) {
1006         /* check for same set of sharing subdomains */
1007         same_set=PETSC_TRUE;
1008         for (k=0;k<graph->count[j];k++){
1009           if (graph->neighbours_set[i][k]!=graph->neighbours_set[j][k]) {
1010             same_set=PETSC_FALSE;
1011           }
1012         }
1013         /* I found a friend of mine */
1014         if (same_set) {
1015           graph->subset[j]=graph->ncc+1;
1016           ierr = PetscBTSet(graph->touched,j);CHKERRQ(ierr);
1017           nodes_touched++;
1018           graph->queue[total_counts] = j;
1019           total_counts++;
1020         }
1021       }
1022     }
1023     graph->ncc++;
1024   }
1025   /* set default number of subsets (at this point no info on csr graph has been taken into account, so n_subsets = ncc */
1026   graph->n_subsets = graph->ncc;
1027   ierr = PetscMalloc1(graph->n_subsets,&graph->subset_ncc);CHKERRQ(ierr);
1028   for (i=0;i<graph->n_subsets;i++) {
1029     graph->subset_ncc[i] = 1;
1030   }
1031   /* final pointer */
1032   graph->cptr[graph->ncc] = total_counts;
1033 
1034   /* save information on subsets (needed if we have to adapt the connected components later) */
1035   /* For consistency reasons (among neighbours), I need to sort (by global ordering) each subset */
1036   /* Get a reference node (min index in global ordering) for each subset for tagging messages */
1037   ierr = PetscMalloc2(graph->ncc,&subset_ref_node_global,graph->cptr[graph->ncc],&queue_global);CHKERRQ(ierr);
1038   ierr = ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_global);CHKERRQ(ierr);
1039   for (j=0;j<graph->ncc;j++) {
1040     ierr = PetscSortIntWithArray(graph->cptr[j+1]-graph->cptr[j],&queue_global[graph->cptr[j]],&graph->queue[graph->cptr[j]]);CHKERRQ(ierr);
1041     subset_ref_node_global[j] = graph->queue[graph->cptr[j]];
1042   }
1043   graph->queue_sorted = PETSC_TRUE;
1044   if (graph->ncc) {
1045     ierr = PetscMalloc2(graph->ncc,&graph->subsets_size,graph->ncc,&graph->subsets);CHKERRQ(ierr);
1046     ierr = PetscMalloc1(graph->cptr[graph->ncc],&graph->subsets[0]);CHKERRQ(ierr);
1047     ierr = PetscMemzero(graph->subsets[0],graph->cptr[graph->ncc]*sizeof(PetscInt));CHKERRQ(ierr);
1048     for (j=1;j<graph->ncc;j++) {
1049       graph->subsets_size[j-1] = graph->cptr[j] - graph->cptr[j-1];
1050       graph->subsets[j] = graph->subsets[j-1] + graph->subsets_size[j-1];
1051     }
1052     graph->subsets_size[graph->ncc-1] = graph->cptr[graph->ncc] - graph->cptr[graph->ncc-1];
1053     for (j=0;j<graph->ncc;j++) {
1054       ierr = PetscMemcpy(graph->subsets[j],&graph->queue[graph->cptr[j]],graph->subsets_size[j]*sizeof(PetscInt));CHKERRQ(ierr);
1055     }
1056   }
1057   /* renumber reference nodes */
1058   ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)(graph->l2gmap)),graph->l2gmap,graph->ncc,subset_ref_node_global,NULL,&k,&graph->subset_ref_node);CHKERRQ(ierr);
1059 
1060   /* free workspace */
1061   ierr = PetscFree2(subset_ref_node_global,queue_global);CHKERRQ(ierr);
1062   ierr = VecDestroy(&local_vec);CHKERRQ(ierr);
1063   ierr = VecDestroy(&local_vec2);CHKERRQ(ierr);
1064   ierr = VecDestroy(&global_vec);CHKERRQ(ierr);
1065   ierr = VecScatterDestroy(&scatter_ctx);CHKERRQ(ierr);
1066   PetscFunctionReturn(0);
1067 }
1068 
1069 #undef __FUNCT__
1070 #define __FUNCT__ "PCBDDCGraphResetCSR"
1071 PetscErrorCode PCBDDCGraphResetCSR(PCBDDCGraph graph)
1072 {
1073   PetscErrorCode ierr;
1074 
1075   PetscFunctionBegin;
1076   ierr = PetscFree(graph->xadj);CHKERRQ(ierr);
1077   ierr = PetscFree(graph->adjncy);CHKERRQ(ierr);
1078   graph->nvtxs_csr = 0;
1079   PetscFunctionReturn(0);
1080 }
1081 
1082 #undef __FUNCT__
1083 #define __FUNCT__ "PCBDDCGraphReset"
1084 PetscErrorCode PCBDDCGraphReset(PCBDDCGraph graph)
1085 {
1086   PetscErrorCode ierr;
1087 
1088   PetscFunctionBegin;
1089   ierr = ISLocalToGlobalMappingDestroy(&graph->l2gmap);CHKERRQ(ierr);
1090   ierr = PetscFree(graph->subset_ncc);CHKERRQ(ierr);
1091   ierr = PetscFree(graph->subset_ref_node);CHKERRQ(ierr);
1092   if (graph->nvtxs) {
1093     ierr = PetscFree(graph->neighbours_set[0]);CHKERRQ(ierr);
1094   }
1095   ierr = PetscBTDestroy(&graph->touched);CHKERRQ(ierr);
1096   ierr = PetscFree5(graph->count,
1097                     graph->neighbours_set,
1098                     graph->subset,
1099                     graph->which_dof,
1100                     graph->special_dof);CHKERRQ(ierr);
1101   ierr = PetscFree2(graph->cptr,graph->queue);CHKERRQ(ierr);
1102   if (graph->mirrors) {
1103     ierr = PetscFree(graph->mirrors_set[0]);CHKERRQ(ierr);
1104   }
1105   ierr = PetscFree2(graph->mirrors,graph->mirrors_set);CHKERRQ(ierr);
1106   if (graph->subsets) {
1107     ierr = PetscFree(graph->subsets[0]);CHKERRQ(ierr);
1108   }
1109   ierr = PetscFree2(graph->subsets_size,graph->subsets);CHKERRQ(ierr);
1110   ierr = ISDestroy(&graph->dirdofs);CHKERRQ(ierr);
1111   graph->has_dirichlet = PETSC_FALSE;
1112   graph->nvtxs = 0;
1113   graph->nvtxs_global = 0;
1114   graph->n_subsets = 0;
1115   graph->custom_minimal_size = 1;
1116   PetscFunctionReturn(0);
1117 }
1118 
1119 #undef __FUNCT__
1120 #define __FUNCT__ "PCBDDCGraphInit"
1121 PetscErrorCode PCBDDCGraphInit(PCBDDCGraph graph, ISLocalToGlobalMapping l2gmap, PetscInt N)
1122 {
1123   PetscInt       n;
1124   PetscErrorCode ierr;
1125 
1126   PetscFunctionBegin;
1127   PetscValidPointer(graph,1);
1128   PetscValidHeaderSpecific(l2gmap,IS_LTOGM_CLASSID,2);
1129   PetscValidLogicalCollectiveInt(l2gmap,N,3);
1130   /* raise an error if already allocated */
1131   if (graph->nvtxs_global) {
1132     SETERRQ(PetscObjectComm((PetscObject)l2gmap),PETSC_ERR_PLIB,"BDDCGraph already initialized");
1133   }
1134   /* set number of vertices */
1135   ierr = PetscObjectReference((PetscObject)l2gmap);CHKERRQ(ierr);
1136   graph->l2gmap = l2gmap;
1137   ierr = ISLocalToGlobalMappingGetSize(l2gmap,&n);CHKERRQ(ierr);
1138   graph->nvtxs = n;
1139   graph->nvtxs_global = N;
1140   /* allocate used space */
1141   ierr = PetscBTCreate(graph->nvtxs,&graph->touched);CHKERRQ(ierr);
1142   ierr = PetscMalloc5(graph->nvtxs,&graph->count,
1143                       graph->nvtxs,&graph->neighbours_set,
1144                       graph->nvtxs,&graph->subset,
1145                       graph->nvtxs,&graph->which_dof,
1146                       graph->nvtxs,&graph->special_dof);CHKERRQ(ierr);
1147   /* zeroes memory */
1148   ierr = PetscMemzero(graph->count,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
1149   ierr = PetscMemzero(graph->subset,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
1150   /* use -1 as a default value for which_dof array */
1151   for (n=0;n<graph->nvtxs;n++) graph->which_dof[n] = -1;
1152   ierr = PetscMemzero(graph->special_dof,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
1153   /* zeroes first pointer to neighbour set */
1154   if (graph->nvtxs) {
1155     graph->neighbours_set[0] = 0;
1156   }
1157   /* zeroes workspace for values of ncc */
1158   graph->subset_ncc = 0;
1159   graph->subset_ref_node = 0;
1160   PetscFunctionReturn(0);
1161 }
1162 
1163 #undef __FUNCT__
1164 #define __FUNCT__ "PCBDDCGraphDestroy"
1165 PetscErrorCode PCBDDCGraphDestroy(PCBDDCGraph* graph)
1166 {
1167   PetscErrorCode ierr;
1168 
1169   PetscFunctionBegin;
1170   ierr = PCBDDCGraphReset(*graph);CHKERRQ(ierr);
1171   ierr = PetscFree(*graph);CHKERRQ(ierr);
1172   PetscFunctionReturn(0);
1173 }
1174 
1175 #undef __FUNCT__
1176 #define __FUNCT__ "PCBDDCGraphCreate"
1177 PetscErrorCode PCBDDCGraphCreate(PCBDDCGraph *graph)
1178 {
1179   PCBDDCGraph    new_graph;
1180   PetscErrorCode ierr;
1181 
1182   PetscFunctionBegin;
1183   ierr = PetscNew(&new_graph);CHKERRQ(ierr);
1184   new_graph->custom_minimal_size = 1;
1185   *graph = new_graph;
1186   PetscFunctionReturn(0);
1187 }
1188