xref: /petsc/src/ksp/pc/impls/bddc/bddcgraph.c (revision 08b6f8921c6dba6c65fa1b9fa8116dc04b43cd3a)
1af0996ceSBarry Smith #include <petsc/private/petscimpl.h>
2ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
3ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddcstructs.h>
4674ae819SStefano Zampini 
5674ae819SStefano Zampini #undef __FUNCT__
6d62866d3SStefano Zampini #define __FUNCT__ "PCBDDCGraphGetDirichletDofsB"
7d62866d3SStefano Zampini PetscErrorCode PCBDDCGraphGetDirichletDofsB(PCBDDCGraph graph, IS* dirdofs)
8d62866d3SStefano Zampini {
9d62866d3SStefano Zampini   PetscErrorCode ierr;
10d62866d3SStefano Zampini 
11d62866d3SStefano Zampini   PetscFunctionBegin;
12d62866d3SStefano Zampini   if (graph->dirdofsB) {
13d62866d3SStefano Zampini     ierr = PetscObjectReference((PetscObject)graph->dirdofsB);CHKERRQ(ierr);
14d62866d3SStefano Zampini   } else if (graph->has_dirichlet) {
15d62866d3SStefano Zampini     PetscInt i,size;
16d62866d3SStefano Zampini     PetscInt *dirdofs_idxs;
17d62866d3SStefano Zampini 
18d62866d3SStefano Zampini     size = 0;
19d62866d3SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
20d62866d3SStefano Zampini       if (graph->count[i] && graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) size++;
21d62866d3SStefano Zampini     }
22d62866d3SStefano Zampini 
23d62866d3SStefano Zampini     ierr = PetscMalloc1(size,&dirdofs_idxs);CHKERRQ(ierr);
24d62866d3SStefano Zampini     size = 0;
25d62866d3SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
26d62866d3SStefano Zampini       if (graph->count[i] && graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) dirdofs_idxs[size++] = i;
27d62866d3SStefano Zampini     }
28d62866d3SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,size,dirdofs_idxs,PETSC_OWN_POINTER,&graph->dirdofsB);CHKERRQ(ierr);
29d62866d3SStefano Zampini     ierr = PetscObjectReference((PetscObject)graph->dirdofsB);CHKERRQ(ierr);
30d62866d3SStefano Zampini   }
31d62866d3SStefano Zampini   *dirdofs = graph->dirdofsB;
32d62866d3SStefano Zampini   PetscFunctionReturn(0);
33d62866d3SStefano Zampini }
34d62866d3SStefano Zampini 
35d62866d3SStefano Zampini #undef __FUNCT__
361b968477SStefano Zampini #define __FUNCT__ "PCBDDCGraphGetDirichletDofs"
371b968477SStefano Zampini PetscErrorCode PCBDDCGraphGetDirichletDofs(PCBDDCGraph graph, IS* dirdofs)
381b968477SStefano Zampini {
391b968477SStefano Zampini   PetscErrorCode ierr;
401b968477SStefano Zampini 
411b968477SStefano Zampini   PetscFunctionBegin;
42a07ea27aSStefano Zampini   if (graph->dirdofs) {
43a07ea27aSStefano Zampini     ierr = PetscObjectReference((PetscObject)graph->dirdofs);CHKERRQ(ierr);
44a07ea27aSStefano Zampini   } else if (graph->has_dirichlet) {
45a07ea27aSStefano Zampini     PetscInt i,size;
46a07ea27aSStefano Zampini     PetscInt *dirdofs_idxs;
47a07ea27aSStefano Zampini 
481b968477SStefano Zampini     size = 0;
491b968477SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
501b968477SStefano Zampini       if (graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) size++;
511b968477SStefano Zampini     }
521b968477SStefano Zampini 
531b968477SStefano Zampini     ierr = PetscMalloc1(size,&dirdofs_idxs);CHKERRQ(ierr);
541b968477SStefano Zampini     size = 0;
551b968477SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
561b968477SStefano Zampini       if (graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) dirdofs_idxs[size++] = i;
571b968477SStefano Zampini     }
58a07ea27aSStefano Zampini     ierr = ISCreateGeneral(PetscObjectComm((PetscObject)graph->l2gmap),size,dirdofs_idxs,PETSC_OWN_POINTER,&graph->dirdofs);CHKERRQ(ierr);
59a07ea27aSStefano Zampini     ierr = PetscObjectReference((PetscObject)graph->dirdofs);CHKERRQ(ierr);
601b968477SStefano Zampini   }
61a07ea27aSStefano Zampini   *dirdofs = graph->dirdofs;
621b968477SStefano Zampini   PetscFunctionReturn(0);
631b968477SStefano Zampini }
641b968477SStefano Zampini 
651b968477SStefano Zampini #undef __FUNCT__
66674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphASCIIView"
67e49050b4SStefano Zampini PetscErrorCode PCBDDCGraphASCIIView(PCBDDCGraph graph, PetscInt verbosity_level, PetscViewer viewer)
68674ae819SStefano Zampini {
692b510759SStefano Zampini   PetscInt       i,j,tabs;
7093fb5fd3SStefano Zampini   PetscInt*      queue_in_global_numbering;
71674ae819SStefano Zampini   PetscErrorCode ierr;
72674ae819SStefano Zampini 
73674ae819SStefano Zampini   PetscFunctionBegin;
741575c14dSBarry Smith   ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
752b510759SStefano Zampini   ierr = PetscViewerASCIIGetTab(viewer,&tabs);CHKERRQ(ierr);
76674ae819SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"--------------------------------------------------\n");CHKERRQ(ierr);
77674ae819SStefano Zampini   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
78674ae819SStefano Zampini   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Local BDDC graph for subdomain %04d\n",PetscGlobalRank);CHKERRQ(ierr);
79674ae819SStefano Zampini   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Number of vertices %d\n",graph->nvtxs);CHKERRQ(ierr);
8014f95afaSStefano Zampini   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Custom minimal size %d\n",graph->custom_minimal_size);CHKERRQ(ierr);
81be12c134Sstefano_zampini   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Max count %d\n",graph->maxcount);CHKERRQ(ierr);
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;
140*08b6f892Sstefano_zampini     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  cc %d (size %d, fid %d, neighs:",i,graph->cptr[i+1]-graph->cptr[i],graph->which_dof[node_num]);CHKERRQ(ierr);
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;
204be12c134Sstefano_zampini   PetscInt       i,nfc,nec,nvc,*idx,*mark;
205674ae819SStefano Zampini   PetscErrorCode ierr;
206674ae819SStefano Zampini 
207674ae819SStefano Zampini   PetscFunctionBegin;
208be12c134Sstefano_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]];
215be12c134Sstefano_zampini     if (graph->cptr[i+1]-graph->cptr[i] > graph->custom_minimal_size && graph->count[repdof] < graph->maxcount) {
216be12c134Sstefano_zampini       if (!graph->twodim && graph->count[repdof] == 1 && graph->special_dof[repdof] != PCBDDCGRAPH_NEUMANN_MARK) {
217674ae819SStefano Zampini         nfc++;
218be12c134Sstefano_zampini         mark[i] = 2;
219be12c134Sstefano_zampini       } else {
220674ae819SStefano Zampini         nec++;
221be12c134Sstefano_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++) {
254be12c134Sstefano_zampini     if (mark[i] == 2) {
255cf5a6209SStefano Zampini       if (FacesIS) {
256c8272957SStefano 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);
257cf5a6209SStefano Zampini       }
258674ae819SStefano Zampini       nfc++;
259be12c134Sstefano_zampini     } else if (mark[i] == 1) {
260cf5a6209SStefano Zampini       if (EdgesIS) {
261c8272957SStefano 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);
262cf5a6209SStefano Zampini       }
263674ae819SStefano Zampini       nec++;
264674ae819SStefano Zampini     }
265674ae819SStefano Zampini   }
266be12c134Sstefano_zampini 
267674ae819SStefano Zampini   /* index set for vertices */
268cf5a6209SStefano Zampini   if (VerticesIS) {
269b8ffe317SStefano Zampini     nvc = 0;
270674ae819SStefano Zampini     for (i=0;i<graph->ncc;i++) {
271be12c134Sstefano_zampini       if (!mark[i]) {
272adfc4fb2SStefano Zampini         PetscInt j;
273adfc4fb2SStefano Zampini 
274674ae819SStefano Zampini         for (j=graph->cptr[i];j<graph->cptr[i+1];j++) {
275674ae819SStefano Zampini           idx[nvc]=graph->queue[j];
276674ae819SStefano Zampini           nvc++;
277674ae819SStefano Zampini         }
278674ae819SStefano Zampini       }
279674ae819SStefano Zampini     }
280674ae819SStefano Zampini     /* sort vertex set (by local ordering) */
281674ae819SStefano Zampini     ierr = PetscSortInt(nvc,idx);CHKERRQ(ierr);
282674ae819SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,nvc,idx,PETSC_OWN_POINTER,&ISForVertices);CHKERRQ(ierr);
283674ae819SStefano Zampini   }
284be12c134Sstefano_zampini   ierr = PetscFree(mark);CHKERRQ(ierr);
285be12c134Sstefano_zampini 
286674ae819SStefano Zampini   /* get back info */
287a873d5d0SStefano Zampini   if (n_faces)       *n_faces = nfc;
288c8272957SStefano Zampini   if (FacesIS)       *FacesIS = ISForFaces;
289a873d5d0SStefano Zampini   if (n_edges)       *n_edges = nec;
290c8272957SStefano Zampini   if (EdgesIS)       *EdgesIS = ISForEdges;
291c8272957SStefano Zampini   if (VerticesIS) *VerticesIS = ISForVertices;
292674ae819SStefano Zampini   PetscFunctionReturn(0);
293674ae819SStefano Zampini }
294674ae819SStefano Zampini 
295674ae819SStefano Zampini #undef __FUNCT__
296674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphComputeConnectedComponents"
297674ae819SStefano Zampini PetscErrorCode PCBDDCGraphComputeConnectedComponents(PCBDDCGraph graph)
298674ae819SStefano Zampini {
2994a060362SStefano Zampini   PetscBool      adapt_interface_reduced;
300674ae819SStefano Zampini   MPI_Comm       interface_comm;
3014a060362SStefano Zampini   PetscMPIInt    size;
3028e4af1ccSStefano Zampini   PetscInt       i;
303674ae819SStefano Zampini   PetscErrorCode ierr;
304674ae819SStefano Zampini 
305674ae819SStefano Zampini   PetscFunctionBegin;
306674ae819SStefano Zampini   /* compute connected components locally */
307674ae819SStefano Zampini   ierr = PetscObjectGetComm((PetscObject)(graph->l2gmap),&interface_comm);CHKERRQ(ierr);
308674ae819SStefano Zampini   ierr = PCBDDCGraphComputeConnectedComponentsLocal(graph);CHKERRQ(ierr);
309674ae819SStefano Zampini   /* check consistency of connected components among neighbouring subdomains -> it adapt them in case it is needed */
3104a060362SStefano Zampini   ierr = MPI_Comm_size(interface_comm,&size);CHKERRQ(ierr);
3114a060362SStefano Zampini   adapt_interface_reduced = PETSC_FALSE;
3124a060362SStefano Zampini   if (size > 1) {
3134a060362SStefano Zampini     PetscInt i;
3144a060362SStefano Zampini     PetscBool adapt_interface = PETSC_FALSE;
315674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
316674ae819SStefano Zampini       /* We are not sure that on a given subset of the local interface,
317674ae819SStefano Zampini          with two connected components, the latters be the same among sharing subdomains */
318674ae819SStefano Zampini       if (graph->subset_ncc[i] > 1) {
3194a060362SStefano Zampini         adapt_interface = PETSC_TRUE;
320674ae819SStefano Zampini         break;
321674ae819SStefano Zampini       }
322674ae819SStefano Zampini     }
323b2566f29SBarry Smith     ierr = MPIU_Allreduce(&adapt_interface,&adapt_interface_reduced,1,MPIU_BOOL,MPI_LOR,interface_comm);CHKERRQ(ierr);
3244a060362SStefano Zampini   }
325674ae819SStefano Zampini 
326674ae819SStefano Zampini   if (graph->n_subsets && adapt_interface_reduced) {
327ec1c413dSStefano Zampini     PetscBT     subset_cc_adapt;
328ec1c413dSStefano Zampini     MPI_Request *send_requests,*recv_requests;
329ec1c413dSStefano Zampini     PetscInt    *send_buffer,*recv_buffer;
330ec1c413dSStefano Zampini     PetscInt    sum_requests,start_of_recv,start_of_send;
331ec1c413dSStefano Zampini     PetscInt    *cum_recv_counts;
332ec1c413dSStefano Zampini     PetscInt    *labels;
333b3cdcd63SStefano Zampini     PetscInt    ncc,cum_queue,mss,mns,j,k,s;
3348875e3e6SStefano Zampini     PetscInt    **refine_buffer=NULL,*private_labels = NULL;
3355b1b9aeaSStefano Zampini 
336ec1c413dSStefano Zampini     ierr = PetscMalloc1(graph->nvtxs,&labels);CHKERRQ(ierr);
337ec1c413dSStefano Zampini     ierr = PetscMemzero(labels,graph->nvtxs*sizeof(*labels));CHKERRQ(ierr);
338b3cdcd63SStefano Zampini     for (i=0;i<graph->ncc;i++)
339b3cdcd63SStefano Zampini       for (j=graph->cptr[i];j<graph->cptr[i+1];j++)
340b3cdcd63SStefano Zampini         labels[graph->queue[j]] = i;
341b3cdcd63SStefano Zampini 
342674ae819SStefano Zampini     /* allocate some space */
343854ce69bSBarry Smith     ierr = PetscMalloc1(graph->n_subsets+1,&cum_recv_counts);CHKERRQ(ierr);
344674ae819SStefano Zampini     ierr = PetscMemzero(cum_recv_counts,(graph->n_subsets+1)*sizeof(*cum_recv_counts));CHKERRQ(ierr);
345b3cdcd63SStefano Zampini 
346674ae819SStefano Zampini     /* first count how many neighbours per connected component I will receive from */
347674ae819SStefano Zampini     cum_recv_counts[0] = 0;
348b3cdcd63SStefano 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]];
349ec1c413dSStefano Zampini     ierr = PetscMalloc1(cum_recv_counts[graph->n_subsets],&recv_buffer);CHKERRQ(ierr);
350ec1c413dSStefano Zampini     ierr = PetscMalloc2(cum_recv_counts[graph->n_subsets],&send_requests,cum_recv_counts[graph->n_subsets],&recv_requests);CHKERRQ(ierr);
351674ae819SStefano Zampini     for (i=0;i<cum_recv_counts[graph->n_subsets];i++) {
352674ae819SStefano Zampini       send_requests[i] = MPI_REQUEST_NULL;
353674ae819SStefano Zampini       recv_requests[i] = MPI_REQUEST_NULL;
354674ae819SStefano Zampini     }
355b3cdcd63SStefano Zampini 
356ec1c413dSStefano Zampini     /* exchange with my neighbours the number of my connected components on the subset of interface */
357674ae819SStefano Zampini     sum_requests = 0;
358674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
359ec1c413dSStefano Zampini       PetscMPIInt neigh,tag;
360b3cdcd63SStefano Zampini       PetscInt    count,*neighs;
361ec1c413dSStefano Zampini 
362b3cdcd63SStefano Zampini       count = graph->count[graph->subset_idxs[i][0]];
363b3cdcd63SStefano Zampini       neighs = graph->neighbours_set[graph->subset_idxs[i][0]];
364ec1c413dSStefano Zampini       ierr = PetscMPIIntCast(2*graph->subset_ref_node[i],&tag);CHKERRQ(ierr);
365b3cdcd63SStefano Zampini       for (k=0;k<count;k++) {
366b3cdcd63SStefano Zampini         ierr = PetscMPIIntCast(neighs[k],&neigh);CHKERRQ(ierr);
367674ae819SStefano Zampini         ierr = MPI_Isend(&graph->subset_ncc[i],1,MPIU_INT,neigh,tag,interface_comm,&send_requests[sum_requests]);CHKERRQ(ierr);
368ec1c413dSStefano Zampini         ierr = MPI_Irecv(&recv_buffer[sum_requests],1,MPIU_INT,neigh,tag,interface_comm,&recv_requests[sum_requests]);CHKERRQ(ierr);
369674ae819SStefano Zampini         sum_requests++;
370674ae819SStefano Zampini       }
371674ae819SStefano Zampini     }
372674ae819SStefano Zampini     ierr = MPI_Waitall(sum_requests,recv_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
373674ae819SStefano Zampini     ierr = MPI_Waitall(sum_requests,send_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
374b3cdcd63SStefano Zampini 
375b3cdcd63SStefano Zampini     /* determine the subsets I have to adapt (those having more than 1 cc) */
376ec1c413dSStefano Zampini     ierr = PetscBTCreate(graph->n_subsets,&subset_cc_adapt);CHKERRQ(ierr);
377ec1c413dSStefano Zampini     ierr = PetscBTMemzero(graph->n_subsets,subset_cc_adapt);CHKERRQ(ierr);
378674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
379b3cdcd63SStefano Zampini       if (graph->subset_ncc[i] > 1) {
380b3cdcd63SStefano Zampini         ierr = PetscBTSet(subset_cc_adapt,i);CHKERRQ(ierr);
381b3cdcd63SStefano Zampini         continue;
382b3cdcd63SStefano Zampini       }
383674ae819SStefano Zampini       for (j=cum_recv_counts[i];j<cum_recv_counts[i+1];j++){
384b3cdcd63SStefano Zampini          if (recv_buffer[j] > 1) {
385f913dca9SStefano Zampini           ierr = PetscBTSet(subset_cc_adapt,i);CHKERRQ(ierr);
386674ae819SStefano Zampini           break;
387674ae819SStefano Zampini         }
388674ae819SStefano Zampini       }
389674ae819SStefano Zampini     }
390ec1c413dSStefano Zampini     ierr = PetscFree(recv_buffer);CHKERRQ(ierr);
391b3cdcd63SStefano Zampini 
392ec1c413dSStefano Zampini     /* determine send/recv buffers sizes */
393ec1c413dSStefano Zampini     j = 0;
394b3cdcd63SStefano Zampini     mss = 0;
395674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
396ec1c413dSStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i)) {
397b3cdcd63SStefano Zampini         j += graph->subset_size[i];
398b3cdcd63SStefano Zampini         mss = PetscMax(graph->subset_size[i],mss);
399674ae819SStefano Zampini       }
400674ae819SStefano Zampini     }
401ec1c413dSStefano Zampini     k = 0;
402b3cdcd63SStefano Zampini     mns = 0;
403674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
404ec1c413dSStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i)) {
405b3cdcd63SStefano Zampini         k += (cum_recv_counts[i+1]-cum_recv_counts[i])*graph->subset_size[i];
406b3cdcd63SStefano Zampini         mns = PetscMax(cum_recv_counts[i+1]-cum_recv_counts[i],mns);
407674ae819SStefano Zampini       }
408674ae819SStefano Zampini     }
409ec1c413dSStefano Zampini     ierr = PetscMalloc2(j,&send_buffer,k,&recv_buffer);CHKERRQ(ierr);
410ec1c413dSStefano Zampini 
411b3cdcd63SStefano Zampini     /* fill send buffer (order matters: subset_idxs ordered by global ordering) */
412ec1c413dSStefano Zampini     j = 0;
413b3cdcd63SStefano Zampini     for (i=0;i<graph->n_subsets;i++)
414b3cdcd63SStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i))
415b3cdcd63SStefano Zampini         for (k=0;k<graph->subset_size[i];k++)
416b3cdcd63SStefano Zampini           send_buffer[j++] = labels[graph->subset_idxs[i][k]];
417ec1c413dSStefano Zampini 
418674ae819SStefano Zampini     /* now exchange the data */
419674ae819SStefano Zampini     start_of_recv = 0;
420674ae819SStefano Zampini     start_of_send = 0;
421674ae819SStefano Zampini     sum_requests = 0;
422674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
423ec1c413dSStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i)) {
424ec1c413dSStefano Zampini         PetscMPIInt neigh,tag;
425b3cdcd63SStefano Zampini         PetscInt    size_of_send = graph->subset_size[i];
426ec1c413dSStefano Zampini 
427b3cdcd63SStefano Zampini         j = graph->subset_idxs[i][0];
428ec1c413dSStefano Zampini         ierr = PetscMPIIntCast(2*graph->subset_ref_node[i]+1,&tag);CHKERRQ(ierr);
429674ae819SStefano Zampini         for (k=0;k<graph->count[j];k++) {
430674ae819SStefano Zampini           ierr = PetscMPIIntCast(graph->neighbours_set[j][k],&neigh);CHKERRQ(ierr);
431ec1c413dSStefano Zampini           ierr = MPI_Isend(&send_buffer[start_of_send],size_of_send,MPIU_INT,neigh,tag,interface_comm,&send_requests[sum_requests]);CHKERRQ(ierr);
432ec1c413dSStefano Zampini           ierr = MPI_Irecv(&recv_buffer[start_of_recv],size_of_send,MPIU_INT,neigh,tag,interface_comm,&recv_requests[sum_requests]);CHKERRQ(ierr);
433ec1c413dSStefano Zampini           start_of_recv += size_of_send;
434674ae819SStefano Zampini           sum_requests++;
435674ae819SStefano Zampini         }
436674ae819SStefano Zampini         start_of_send += size_of_send;
437674ae819SStefano Zampini       }
438674ae819SStefano Zampini     }
439674ae819SStefano Zampini     ierr = MPI_Waitall(sum_requests,recv_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
440d52c4852SStefano Zampini 
441b3cdcd63SStefano Zampini     /* refine connected components */
442674ae819SStefano Zampini     start_of_recv = 0;
443b3cdcd63SStefano Zampini     /* allocate some temporary space */
444b3cdcd63SStefano Zampini     if (mss) {
445b3cdcd63SStefano Zampini       ierr = PetscMalloc1(mss,&refine_buffer);CHKERRQ(ierr);
446b3cdcd63SStefano Zampini       ierr = PetscMalloc2(mss*(mns+1),&refine_buffer[0],mss,&private_labels);CHKERRQ(ierr);
447b3cdcd63SStefano Zampini     }
448b3cdcd63SStefano Zampini     ncc = 0;
449b3cdcd63SStefano Zampini     cum_queue = 0;
450b3cdcd63SStefano Zampini     graph->cptr[0] = 0;
451674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
452ec1c413dSStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i)) {
453b3cdcd63SStefano Zampini         PetscInt subset_counter = 0;
454b3cdcd63SStefano Zampini         PetscInt sharingprocs = cum_recv_counts[i+1]-cum_recv_counts[i]+1; /* count myself */
455b3cdcd63SStefano Zampini         PetscInt buffer_size = graph->subset_size[i];
456ec1c413dSStefano Zampini 
457b3cdcd63SStefano Zampini         /* compute pointers */
458b3cdcd63SStefano Zampini         for (j=1;j<buffer_size;j++) refine_buffer[j] = refine_buffer[j-1] + sharingprocs;
459b3cdcd63SStefano Zampini         /* analyze contributions from subdomains that share the i-th subset
460b3cdcd63SStefano Zampini            The stricture of refine_buffer is suitable to find intersections of ccs among sharingprocs.
461b3cdcd63SStefano 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)
462b3cdcd63SStefano Zampini            sharing procs connected components:
463ec1c413dSStefano Zampini              neigh 0: [0 1 4], [2 3], labels [4,7]  (2 connected components)
464ec1c413dSStefano Zampini              neigh 1: [0 1], [2 3 4], labels [3 2]  (2 connected components)
465ec1c413dSStefano Zampini              neigh 2: [0 4], [1], [2 3], labels [1 5 6] (3 connected components)
466b3cdcd63SStefano Zampini            refine_buffer will be filled as:
467ec1c413dSStefano Zampini              [ 4, 3, 1;
468ec1c413dSStefano Zampini                4, 2, 1;
469ec1c413dSStefano Zampini                7, 2, 6;
470ec1c413dSStefano Zampini                4, 3, 5;
471ec1c413dSStefano Zampini                7, 2, 6; ];
472b3cdcd63SStefano Zampini            The connected components in local ordering are [0], [1], [2 3], [4] */
473b3cdcd63SStefano Zampini         /* fill temp_buffer */
474b3cdcd63SStefano Zampini         for (k=0;k<buffer_size;k++) refine_buffer[k][0] = labels[graph->subset_idxs[i][k]];
475b3cdcd63SStefano Zampini         for (j=0;j<sharingprocs-1;j++) {
476b3cdcd63SStefano Zampini           for (k=0;k<buffer_size;k++) refine_buffer[k][j+1] = recv_buffer[start_of_recv+k];
477b3cdcd63SStefano Zampini           start_of_recv += buffer_size;
478674ae819SStefano Zampini         }
479b3cdcd63SStefano Zampini         ierr = PetscMemzero(private_labels,buffer_size*sizeof(PetscInt));CHKERRQ(ierr);
480b3cdcd63SStefano Zampini         for (j=0;j<buffer_size;j++) {
481b3cdcd63SStefano Zampini           if (!private_labels[j]) { /* found a new cc  */
482ec1c413dSStefano Zampini             PetscBool same_set;
483ec1c413dSStefano Zampini 
484b3cdcd63SStefano Zampini             graph->cptr[ncc] = cum_queue;
485b3cdcd63SStefano Zampini             ncc++;
486b3cdcd63SStefano Zampini             subset_counter++;
487b3cdcd63SStefano Zampini             private_labels[j] = subset_counter;
488b3cdcd63SStefano Zampini             graph->queue[cum_queue++] = graph->subset_idxs[i][j];
489b3cdcd63SStefano Zampini             for (k=j+1;k<buffer_size;k++) { /* check for other nodes in new cc */
490674ae819SStefano Zampini               same_set = PETSC_TRUE;
491b3cdcd63SStefano Zampini               for (s=0;s<sharingprocs;s++) {
492b3cdcd63SStefano Zampini                 if (refine_buffer[j][s] != refine_buffer[k][s]) {
493674ae819SStefano Zampini                   same_set = PETSC_FALSE;
494674ae819SStefano Zampini                   break;
495674ae819SStefano Zampini                 }
496674ae819SStefano Zampini               }
497674ae819SStefano Zampini               if (same_set) {
498b3cdcd63SStefano Zampini                 private_labels[k] = subset_counter;
499b3cdcd63SStefano Zampini                 graph->queue[cum_queue++] = graph->subset_idxs[i][k];
500674ae819SStefano Zampini               }
501674ae819SStefano Zampini             }
502674ae819SStefano Zampini           }
503674ae819SStefano Zampini         }
504b3cdcd63SStefano Zampini         graph->cptr[ncc] = cum_queue;
505b3cdcd63SStefano Zampini         graph->subset_ncc[i] = subset_counter;
506b3cdcd63SStefano Zampini         graph->queue_sorted = PETSC_FALSE;
507b3cdcd63SStefano Zampini       } else { /* this subset does not need to be adapted */
508b3cdcd63SStefano Zampini         ierr = PetscMemcpy(graph->queue+cum_queue,graph->subset_idxs[i],graph->subset_size[i]*sizeof(PetscInt));CHKERRQ(ierr);
509b3cdcd63SStefano Zampini         ncc++;
510b3cdcd63SStefano Zampini         cum_queue += graph->subset_size[i];
511b3cdcd63SStefano Zampini         graph->cptr[ncc] = cum_queue;
512674ae819SStefano Zampini       }
513674ae819SStefano Zampini     }
514b3cdcd63SStefano Zampini     graph->cptr[ncc] = cum_queue;
515b3cdcd63SStefano Zampini     graph->ncc = ncc;
516b3cdcd63SStefano Zampini     if (mss) {
517b3cdcd63SStefano Zampini       ierr = PetscFree2(refine_buffer[0],private_labels);CHKERRQ(ierr);
518b3cdcd63SStefano Zampini       ierr = PetscFree(refine_buffer);CHKERRQ(ierr);
519b3cdcd63SStefano Zampini     }
520b3cdcd63SStefano Zampini     ierr = PetscFree(labels);CHKERRQ(ierr);
521b3cdcd63SStefano Zampini     ierr = MPI_Waitall(sum_requests,send_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
522b3cdcd63SStefano Zampini     ierr = PetscFree2(send_requests,recv_requests);CHKERRQ(ierr);
52345b2a5aaSStefano Zampini     ierr = PetscFree2(send_buffer,recv_buffer);CHKERRQ(ierr);
524674ae819SStefano Zampini     ierr = PetscFree(cum_recv_counts);CHKERRQ(ierr);
525ec1c413dSStefano Zampini     ierr = PetscBTDestroy(&subset_cc_adapt);CHKERRQ(ierr);
526674ae819SStefano Zampini   }
5278e4af1ccSStefano Zampini 
528674ae819SStefano Zampini   PetscFunctionReturn(0);
529674ae819SStefano Zampini }
530674ae819SStefano Zampini 
531b3cdcd63SStefano Zampini 
532b3cdcd63SStefano Zampini #undef __FUNCT__
533b3cdcd63SStefano Zampini #define __FUNCT__ "PCBDDCGraphComputeCC_Private"
534b3cdcd63SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCGraphComputeCC_Private(PCBDDCGraph graph,PetscInt pid,PetscInt* queue_tip,PetscInt n_prev,PetscInt* n_added)
535b3cdcd63SStefano Zampini {
536b3cdcd63SStefano Zampini   PetscInt       i,j,n;
537b3cdcd63SStefano Zampini   PetscInt       *xadj = graph->xadj,*adjncy = graph->adjncy;
538b3cdcd63SStefano Zampini   PetscBT        touched = graph->touched;
539b3cdcd63SStefano Zampini   PetscBool      havecsr = (PetscBool)(xadj && adjncy);
540b3cdcd63SStefano Zampini   PetscBool      havesubs = (PetscBool)(!!graph->n_local_subs);
541b3cdcd63SStefano Zampini   PetscErrorCode ierr;
542b3cdcd63SStefano Zampini 
543b3cdcd63SStefano Zampini   PetscFunctionBegin;
544b3cdcd63SStefano Zampini   n = 0;
545b3cdcd63SStefano Zampini   if (havecsr && !havesubs) {
546b3cdcd63SStefano Zampini     for (i=-n_prev;i<0;i++) {
547b3cdcd63SStefano Zampini       PetscInt start_dof = queue_tip[i];
54854fffbccSStefano Zampini       /* we assume that if a dof has a size 1 adjacency list and the corresponding entry is negative, it is connected to all dofs */
54954fffbccSStefano Zampini       if (xadj[start_dof+1]-xadj[start_dof] == 1 && adjncy[xadj[start_dof]] < 0) {
55054fffbccSStefano Zampini         for (j=0;j<graph->subset_size[pid-1];j++) { /* pid \in [1,graph->n_subsets] */
55154fffbccSStefano Zampini           PetscInt dof = graph->subset_idxs[pid-1][j];
55254fffbccSStefano Zampini           if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid) {
55354fffbccSStefano Zampini             ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
55454fffbccSStefano Zampini             queue_tip[n] = dof;
55554fffbccSStefano Zampini             n++;
55654fffbccSStefano Zampini           }
55754fffbccSStefano Zampini         }
55854fffbccSStefano Zampini       } else {
559b3cdcd63SStefano Zampini         for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
560b3cdcd63SStefano Zampini           PetscInt dof = adjncy[j];
561b3cdcd63SStefano Zampini           if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid) {
562b3cdcd63SStefano Zampini             ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
563b3cdcd63SStefano Zampini             queue_tip[n] = dof;
564b3cdcd63SStefano Zampini             n++;
565b3cdcd63SStefano Zampini           }
566b3cdcd63SStefano Zampini         }
567b3cdcd63SStefano Zampini       }
56854fffbccSStefano Zampini     }
569b3cdcd63SStefano Zampini   } else if (havecsr && havesubs) {
570b3cdcd63SStefano Zampini     PetscInt sid = graph->local_subs[queue_tip[-n_prev]];
571b3cdcd63SStefano Zampini     for (i=-n_prev;i<0;i++) {
572b3cdcd63SStefano Zampini       PetscInt start_dof = queue_tip[i];
57354fffbccSStefano Zampini       /* we assume that if a dof has a size 1 adjacency list and the corresponding entry is negative, it is connected to all dofs belonging to the local sub */
57454fffbccSStefano Zampini       if (xadj[start_dof+1]-xadj[start_dof] == 1 && adjncy[xadj[start_dof]] < 0) {
57554fffbccSStefano Zampini         for (j=0;j<graph->subset_size[pid-1];j++) { /* pid \in [1,graph->n_subsets] */
57654fffbccSStefano Zampini           PetscInt dof = graph->subset_idxs[pid-1][j];
57754fffbccSStefano Zampini           if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid && graph->local_subs[dof] == sid) {
57854fffbccSStefano Zampini             ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
57954fffbccSStefano Zampini             queue_tip[n] = dof;
58054fffbccSStefano Zampini             n++;
58154fffbccSStefano Zampini           }
58254fffbccSStefano Zampini         }
58354fffbccSStefano Zampini       } else {
584b3cdcd63SStefano Zampini         for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
585b3cdcd63SStefano Zampini           PetscInt dof = adjncy[j];
586b3cdcd63SStefano Zampini           if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid && graph->local_subs[dof] == sid) {
587b3cdcd63SStefano Zampini             ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
588b3cdcd63SStefano Zampini             queue_tip[n] = dof;
589b3cdcd63SStefano Zampini             n++;
590b3cdcd63SStefano Zampini           }
591b3cdcd63SStefano Zampini         }
592b3cdcd63SStefano Zampini       }
59354fffbccSStefano Zampini     }
594b3cdcd63SStefano Zampini   } else { /* sub info only */
595b3cdcd63SStefano Zampini     PetscInt sid = graph->local_subs[queue_tip[-n_prev]];
596b3cdcd63SStefano Zampini     for (j=0;j<graph->subset_size[pid-1];j++) { /* pid \in [1,graph->n_subsets] */
597b3cdcd63SStefano Zampini       PetscInt dof = graph->subset_idxs[pid-1][j];
598b3cdcd63SStefano Zampini       if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid && graph->local_subs[dof] == sid) {
599b3cdcd63SStefano Zampini         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
600b3cdcd63SStefano Zampini         queue_tip[n] = dof;
601b3cdcd63SStefano Zampini         n++;
602b3cdcd63SStefano Zampini       }
603b3cdcd63SStefano Zampini     }
604b3cdcd63SStefano Zampini   }
605b3cdcd63SStefano Zampini   *n_added = n;
606b3cdcd63SStefano Zampini   PetscFunctionReturn(0);
607b3cdcd63SStefano Zampini }
608b3cdcd63SStefano Zampini 
609674ae819SStefano Zampini #undef __FUNCT__
610674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphComputeConnectedComponentsLocal"
611674ae819SStefano Zampini PetscErrorCode PCBDDCGraphComputeConnectedComponentsLocal(PCBDDCGraph graph)
612674ae819SStefano Zampini {
613b3cdcd63SStefano Zampini   PetscInt       ncc,cum_queue,n;
614b3cdcd63SStefano Zampini   PetscMPIInt    commsize;
615674ae819SStefano Zampini   PetscErrorCode ierr;
616674ae819SStefano Zampini 
617674ae819SStefano Zampini   PetscFunctionBegin;
618b3cdcd63SStefano Zampini   if (!graph->setupcalled) SETERRQ(PetscObjectComm((PetscObject)graph->l2gmap),PETSC_ERR_ORDER,"PCBDDCGraphSetUp should be called first");
619b3cdcd63SStefano Zampini   /* quiet return if there isn't any local info */
6204f1b2e48SStefano Zampini   if ((!graph->xadj || !graph->adjncy) && !graph->n_local_subs) {
621674ae819SStefano Zampini     PetscFunctionReturn(0);
622674ae819SStefano Zampini   }
623674ae819SStefano Zampini 
624674ae819SStefano Zampini   /* reset any previous search of connected components */
625df48933bSStefano Zampini   ierr = PetscBTMemzero(graph->nvtxs,graph->touched);CHKERRQ(ierr);
626b3cdcd63SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)graph->l2gmap),&commsize);CHKERRQ(ierr);
6274b2aedd3SStefano Zampini   if (commsize > graph->commsizelimit) {
628b3cdcd63SStefano Zampini     PetscInt i;
629674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
6300cf82fd6SStefano Zampini       if (graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK || !graph->count[i]) {
631df48933bSStefano Zampini         ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
632674ae819SStefano Zampini       }
6334a060362SStefano Zampini     }
634b3cdcd63SStefano Zampini   }
635674ae819SStefano Zampini 
636674ae819SStefano Zampini   /* begin search for connected components */
637674ae819SStefano Zampini   cum_queue = 0;
638674ae819SStefano Zampini   ncc = 0;
639674ae819SStefano Zampini   for (n=0;n<graph->n_subsets;n++) {
640b3cdcd63SStefano Zampini     PetscInt pid = n+1;  /* partition labeled by 0 is discarded */
641b3cdcd63SStefano Zampini     PetscInt found = 0,prev = 0,first = 0,ncc_pid = 0;
642b3cdcd63SStefano Zampini     while (found != graph->subset_size[n]) {
643b3cdcd63SStefano Zampini       PetscInt added = 0;
644b3cdcd63SStefano Zampini       if (!prev) { /* search for new starting dof */
645b3cdcd63SStefano Zampini         while (PetscBTLookup(graph->touched,graph->subset_idxs[n][first])) first++;
646b3cdcd63SStefano Zampini         ierr = PetscBTSet(graph->touched,graph->subset_idxs[n][first]);CHKERRQ(ierr);
647b3cdcd63SStefano Zampini         graph->queue[cum_queue] = graph->subset_idxs[n][first];
648674ae819SStefano Zampini         graph->cptr[ncc] = cum_queue;
649b3cdcd63SStefano Zampini         prev = 1;
650b3cdcd63SStefano Zampini         cum_queue++;
651b3cdcd63SStefano Zampini         found++;
652674ae819SStefano Zampini         ncc_pid++;
653b3cdcd63SStefano Zampini         ncc++;
654674ae819SStefano Zampini       }
655b3cdcd63SStefano Zampini       ierr = PCBDDCGraphComputeCC_Private(graph,pid,graph->queue + cum_queue,prev,&added);CHKERRQ(ierr);
656b3cdcd63SStefano Zampini       if (!added) {
657674ae819SStefano Zampini         graph->subset_ncc[n] = ncc_pid;
658b3cdcd63SStefano Zampini         graph->cptr[ncc] = cum_queue;
659b3cdcd63SStefano Zampini       }
660b3cdcd63SStefano Zampini       prev = added;
661b3cdcd63SStefano Zampini       found += added;
662b3cdcd63SStefano Zampini       cum_queue += added;
663b3cdcd63SStefano Zampini       if (added && found == graph->subset_size[n]) {
664b3cdcd63SStefano Zampini         graph->subset_ncc[n] = ncc_pid;
665b3cdcd63SStefano Zampini         graph->cptr[ncc] = cum_queue;
666b3cdcd63SStefano Zampini       }
667b3cdcd63SStefano Zampini     }
668674ae819SStefano Zampini   }
669674ae819SStefano Zampini   graph->ncc = ncc;
670acc113dbSStefano Zampini   graph->queue_sorted = PETSC_FALSE;
671674ae819SStefano Zampini   PetscFunctionReturn(0);
672674ae819SStefano Zampini }
673674ae819SStefano Zampini 
674674ae819SStefano Zampini #undef __FUNCT__
675674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphSetUp"
676674ae819SStefano Zampini PetscErrorCode PCBDDCGraphSetUp(PCBDDCGraph graph, PetscInt custom_minimal_size, IS neumann_is, IS dirichlet_is, PetscInt n_ISForDofs, IS ISForDofs[], IS custom_primal_vertices)
677674ae819SStefano Zampini {
678674ae819SStefano Zampini   VecScatter     scatter_ctx;
679674ae819SStefano Zampini   Vec            local_vec,local_vec2,global_vec;
680dc456d91SStefano Zampini   IS             to,from,subset,subset_n;
681674ae819SStefano Zampini   MPI_Comm       comm;
682674ae819SStefano Zampini   PetscScalar    *array,*array2;
683674ae819SStefano Zampini   const PetscInt *is_indices;
684dc456d91SStefano Zampini   PetscInt       n_neigh,*neigh,*n_shared,**shared,*queue_global;
685674ae819SStefano Zampini   PetscInt       i,j,k,s,total_counts,nodes_touched,is_size;
686b3cdcd63SStefano Zampini   PetscMPIInt    commsize;
6873ebc787cSStefano Zampini   PetscBool      same_set,mirrors_found,twodim;
6887babac9bSStefano Zampini   PetscErrorCode ierr;
689674ae819SStefano Zampini 
690674ae819SStefano Zampini   PetscFunctionBegin;
691a07ea27aSStefano Zampini   graph->has_dirichlet = PETSC_FALSE;
692a07ea27aSStefano Zampini   if (dirichlet_is) {
693a07ea27aSStefano Zampini     PetscCheckSameComm(graph->l2gmap,1,dirichlet_is,4);
694a07ea27aSStefano Zampini     graph->has_dirichlet = PETSC_TRUE;
695a07ea27aSStefano Zampini   }
696674ae819SStefano Zampini   ierr = PetscObjectGetComm((PetscObject)(graph->l2gmap),&comm);CHKERRQ(ierr);
697b3cdcd63SStefano Zampini   ierr = MPI_Comm_size(comm,&commsize);CHKERRQ(ierr);
698b3cdcd63SStefano Zampini 
699674ae819SStefano Zampini   /* custom_minimal_size */
70014f95afaSStefano Zampini   graph->custom_minimal_size = custom_minimal_size;
701674ae819SStefano Zampini   /* get info l2gmap and allocate work vectors  */
702674ae819SStefano Zampini   ierr = ISLocalToGlobalMappingGetInfo(graph->l2gmap,&n_neigh,&neigh,&n_shared,&shared);CHKERRQ(ierr);
7032635a1d4SStefano Zampini   /* check if we have any local periodic nodes (periodic BCs) */
7042635a1d4SStefano Zampini   mirrors_found = PETSC_FALSE;
7052635a1d4SStefano Zampini   if (graph->nvtxs && n_neigh) {
7062635a1d4SStefano Zampini     for (i=0; i<n_shared[0]; i++) graph->count[shared[0][i]] += 1;
7072635a1d4SStefano Zampini     for (i=0; i<n_shared[0]; i++) {
7082635a1d4SStefano Zampini       if (graph->count[shared[0][i]] > 1) {
7092635a1d4SStefano Zampini         mirrors_found = PETSC_TRUE;
7102635a1d4SStefano Zampini         break;
7112635a1d4SStefano Zampini       }
7122635a1d4SStefano Zampini     }
7132635a1d4SStefano Zampini   }
7142635a1d4SStefano Zampini   /* create some workspace objects */
7152635a1d4SStefano Zampini   local_vec = NULL;
7162635a1d4SStefano Zampini   local_vec2 = NULL;
7172635a1d4SStefano Zampini   global_vec = NULL;
7182635a1d4SStefano Zampini   to = NULL;
7192635a1d4SStefano Zampini   from = NULL;
7202635a1d4SStefano Zampini   scatter_ctx = NULL;
7212635a1d4SStefano Zampini   if (n_ISForDofs || dirichlet_is || neumann_is || custom_primal_vertices) {
722674ae819SStefano Zampini     ierr = VecCreate(PETSC_COMM_SELF,&local_vec);CHKERRQ(ierr);
723674ae819SStefano Zampini     ierr = VecSetSizes(local_vec,PETSC_DECIDE,graph->nvtxs);CHKERRQ(ierr);
724674ae819SStefano Zampini     ierr = VecSetType(local_vec,VECSTANDARD);CHKERRQ(ierr);
725674ae819SStefano Zampini     ierr = VecDuplicate(local_vec,&local_vec2);CHKERRQ(ierr);
726674ae819SStefano Zampini     ierr = VecCreate(comm,&global_vec);CHKERRQ(ierr);
7277fb0e2dbSStefano Zampini     ierr = VecSetSizes(global_vec,PETSC_DECIDE,graph->nvtxs_global);CHKERRQ(ierr);
728674ae819SStefano Zampini     ierr = VecSetType(global_vec,VECSTANDARD);CHKERRQ(ierr);
729674ae819SStefano Zampini     ierr = ISCreateStride(PETSC_COMM_SELF,graph->nvtxs,0,1,&to);CHKERRQ(ierr);
730674ae819SStefano Zampini     ierr = ISLocalToGlobalMappingApplyIS(graph->l2gmap,to,&from);CHKERRQ(ierr);
731674ae819SStefano Zampini     ierr = VecScatterCreate(global_vec,from,local_vec,to,&scatter_ctx);CHKERRQ(ierr);
7322635a1d4SStefano Zampini   } else if (mirrors_found) {
7332635a1d4SStefano Zampini     ierr = ISCreateStride(PETSC_COMM_SELF,graph->nvtxs,0,1,&to);CHKERRQ(ierr);
7342635a1d4SStefano Zampini     ierr = ISLocalToGlobalMappingApplyIS(graph->l2gmap,to,&from);CHKERRQ(ierr);
73549eeff8cSStefano Zampini   }
73651b0f6cfSStefano Zampini   /* compute local mirrors (if any) */
73751b0f6cfSStefano Zampini   if (mirrors_found) {
73851b0f6cfSStefano Zampini     PetscInt *local_indices,*global_indices;
73951b0f6cfSStefano Zampini     /* get arrays of local and global indices */
740785e854fSJed Brown     ierr = PetscMalloc1(graph->nvtxs,&local_indices);CHKERRQ(ierr);
74151b0f6cfSStefano Zampini     ierr = ISGetIndices(to,(const PetscInt**)&is_indices);CHKERRQ(ierr);
74251b0f6cfSStefano Zampini     ierr = PetscMemcpy(local_indices,is_indices,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
74351b0f6cfSStefano Zampini     ierr = ISRestoreIndices(to,(const PetscInt**)&is_indices);CHKERRQ(ierr);
744785e854fSJed Brown     ierr = PetscMalloc1(graph->nvtxs,&global_indices);CHKERRQ(ierr);
74551b0f6cfSStefano Zampini     ierr = ISGetIndices(from,(const PetscInt**)&is_indices);CHKERRQ(ierr);
74651b0f6cfSStefano Zampini     ierr = PetscMemcpy(global_indices,is_indices,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
74751b0f6cfSStefano Zampini     ierr = ISRestoreIndices(from,(const PetscInt**)&is_indices);CHKERRQ(ierr);
74851b0f6cfSStefano Zampini     /* allocate space for mirrors */
7498e5aaad5SJed Brown     ierr = PetscMalloc2(graph->nvtxs,&graph->mirrors,graph->nvtxs,&graph->mirrors_set);CHKERRQ(ierr);
75051b0f6cfSStefano Zampini     ierr = PetscMemzero(graph->mirrors,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
75151b0f6cfSStefano Zampini     graph->mirrors_set[0] = 0;
75251b0f6cfSStefano Zampini 
75351b0f6cfSStefano Zampini     k=0;
75451b0f6cfSStefano Zampini     for (i=0;i<n_shared[0];i++) {
75551b0f6cfSStefano Zampini       j=shared[0][i];
75651b0f6cfSStefano Zampini       if (graph->count[j] > 1) {
75751b0f6cfSStefano Zampini         graph->mirrors[j]++;
75851b0f6cfSStefano Zampini         k++;
75951b0f6cfSStefano Zampini       }
76051b0f6cfSStefano Zampini     }
76151b0f6cfSStefano Zampini     /* allocate space for set of mirrors */
762785e854fSJed Brown     ierr = PetscMalloc1(k,&graph->mirrors_set[0]);CHKERRQ(ierr);
76351b0f6cfSStefano Zampini     for (i=1;i<graph->nvtxs;i++)
76451b0f6cfSStefano Zampini       graph->mirrors_set[i]=graph->mirrors_set[i-1]+graph->mirrors[i-1];
76551b0f6cfSStefano Zampini 
76651b0f6cfSStefano Zampini     /* fill arrays */
76751b0f6cfSStefano Zampini     ierr = PetscMemzero(graph->mirrors,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
76851b0f6cfSStefano Zampini     for (j=0;j<n_shared[0];j++) {
76951b0f6cfSStefano Zampini       i=shared[0][j];
77051b0f6cfSStefano Zampini       if (graph->count[i] > 1)
77151b0f6cfSStefano Zampini         graph->mirrors_set[i][graph->mirrors[i]++]=global_indices[i];
77251b0f6cfSStefano Zampini     }
77351b0f6cfSStefano Zampini     ierr = PetscSortIntWithArray(graph->nvtxs,global_indices,local_indices);CHKERRQ(ierr);
77451b0f6cfSStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
77551b0f6cfSStefano Zampini       if (graph->mirrors[i] > 0) {
77651b0f6cfSStefano Zampini         ierr = PetscFindInt(graph->mirrors_set[i][0],graph->nvtxs,global_indices,&k);CHKERRQ(ierr);
77751b0f6cfSStefano Zampini         j = global_indices[k];
77851b0f6cfSStefano Zampini         while ( k > 0 && global_indices[k-1] == j) k--;
77951b0f6cfSStefano Zampini         for (j=0;j<graph->mirrors[i];j++) {
78051b0f6cfSStefano Zampini           graph->mirrors_set[i][j]=local_indices[k+j];
78151b0f6cfSStefano Zampini         }
78251b0f6cfSStefano Zampini         ierr = PetscSortInt(graph->mirrors[i],graph->mirrors_set[i]);CHKERRQ(ierr);
78351b0f6cfSStefano Zampini       }
78451b0f6cfSStefano Zampini     }
78551b0f6cfSStefano Zampini     ierr = PetscFree(local_indices);CHKERRQ(ierr);
78651b0f6cfSStefano Zampini     ierr = PetscFree(global_indices);CHKERRQ(ierr);
78751b0f6cfSStefano Zampini   }
78851b0f6cfSStefano Zampini   ierr = PetscMemzero(graph->count,graph->nvtxs*sizeof(*graph->count));CHKERRQ(ierr);
789674ae819SStefano Zampini   ierr = ISDestroy(&to);CHKERRQ(ierr);
790674ae819SStefano Zampini   ierr = ISDestroy(&from);CHKERRQ(ierr);
79151b0f6cfSStefano Zampini 
792674ae819SStefano Zampini   /* Count total number of neigh per node */
793674ae819SStefano Zampini   k = 0;
794674ae819SStefano Zampini   for (i=1;i<n_neigh;i++) {
795674ae819SStefano Zampini     k += n_shared[i];
796674ae819SStefano Zampini     for (j=0;j<n_shared[i];j++) {
797674ae819SStefano Zampini       graph->count[shared[i][j]] += 1;
798674ae819SStefano Zampini     }
799674ae819SStefano Zampini   }
800674ae819SStefano Zampini   /* Allocate space for storing the set of neighbours for each node */
801674ae819SStefano Zampini   if (graph->nvtxs) {
802785e854fSJed Brown     ierr = PetscMalloc1(k,&graph->neighbours_set[0]);CHKERRQ(ierr);
803674ae819SStefano Zampini   }
804674ae819SStefano Zampini   for (i=1;i<graph->nvtxs;i++) { /* dont count myself */
805674ae819SStefano Zampini     graph->neighbours_set[i]=graph->neighbours_set[i-1]+graph->count[i-1];
806674ae819SStefano Zampini   }
807674ae819SStefano Zampini   /* Get information for sharing subdomains */
808674ae819SStefano Zampini   ierr = PetscMemzero(graph->count,graph->nvtxs*sizeof(*graph->count));CHKERRQ(ierr);
809674ae819SStefano Zampini   for (i=1;i<n_neigh;i++) { /* dont count myself */
810674ae819SStefano Zampini     s = n_shared[i];
811674ae819SStefano Zampini     for (j=0;j<s;j++) {
812674ae819SStefano Zampini       k = shared[i][j];
813674ae819SStefano Zampini       graph->neighbours_set[k][graph->count[k]] = neigh[i];
814674ae819SStefano Zampini       graph->count[k] += 1;
815674ae819SStefano Zampini     }
816674ae819SStefano Zampini   }
817674ae819SStefano Zampini   /* sort set of sharing subdomains */
818674ae819SStefano Zampini   for (i=0;i<graph->nvtxs;i++) {
81993fb5fd3SStefano Zampini     ierr = PetscSortRemoveDupsInt(&graph->count[i],graph->neighbours_set[i]);CHKERRQ(ierr);
820674ae819SStefano Zampini   }
8217fb0e2dbSStefano Zampini   /* free memory allocated by ISLocalToGlobalMappingGetInfo */
8227fb0e2dbSStefano Zampini   ierr = ISLocalToGlobalMappingRestoreInfo(graph->l2gmap,&n_neigh,&neigh,&n_shared,&shared);CHKERRQ(ierr);
8237fb0e2dbSStefano Zampini 
82467c9da69SStefano Zampini   /*
82567c9da69SStefano Zampini      Get info for dofs splitting
8265777c63cSStefano Zampini      User can specify just a subset; an additional field is considered as a complementary field
82767c9da69SStefano Zampini   */
828b3cdcd63SStefano Zampini   for (i=0;i<graph->nvtxs;i++) graph->which_dof[i] = n_ISForDofs; /* by default a dof belongs to the complement set */
8295777c63cSStefano Zampini   if (n_ISForDofs) {
8305777c63cSStefano Zampini     ierr = VecSet(local_vec,-1.0);CHKERRQ(ierr);
8315777c63cSStefano Zampini   }
832674ae819SStefano Zampini   for (i=0;i<n_ISForDofs;i++) {
8335777c63cSStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
83463602bcaSStefano Zampini     ierr = ISGetLocalSize(ISForDofs[i],&is_size);CHKERRQ(ierr);
835674ae819SStefano Zampini     ierr = ISGetIndices(ISForDofs[i],(const PetscInt**)&is_indices);CHKERRQ(ierr);
836674ae819SStefano Zampini     for (j=0;j<is_size;j++) {
837d50ed60aSStefano Zampini       if (is_indices[j] > -1 && is_indices[j] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
838674ae819SStefano Zampini         graph->which_dof[is_indices[j]] = i;
8395777c63cSStefano Zampini         array[is_indices[j]] = 1.*i;
840674ae819SStefano Zampini       }
84167c9da69SStefano Zampini     }
842674ae819SStefano Zampini     ierr = ISRestoreIndices(ISForDofs[i],(const PetscInt**)&is_indices);CHKERRQ(ierr);
8435777c63cSStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
8445777c63cSStefano Zampini   }
8455777c63cSStefano Zampini   /* Check consistency among neighbours */
8465777c63cSStefano Zampini   if (n_ISForDofs) {
8475777c63cSStefano Zampini     ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8485777c63cSStefano Zampini     ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8495777c63cSStefano Zampini     ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8505777c63cSStefano Zampini     ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8515777c63cSStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
8525777c63cSStefano Zampini     ierr = VecGetArray(local_vec2,&array2);CHKERRQ(ierr);
8535777c63cSStefano Zampini     for (i=0;i<graph->nvtxs;i++){
8545777c63cSStefano Zampini       PetscInt field1,field2;
8555777c63cSStefano Zampini 
8565777c63cSStefano Zampini       field1 = (PetscInt)PetscRealPart(array[i]);
8575777c63cSStefano Zampini       field2 = (PetscInt)PetscRealPart(array2[i]);
8586c4ed002SBarry 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);
8595777c63cSStefano Zampini     }
8605777c63cSStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
8615777c63cSStefano Zampini     ierr = VecRestoreArray(local_vec2,&array2);CHKERRQ(ierr);
862674ae819SStefano Zampini   }
863674ae819SStefano Zampini   /* Take into account Neumann nodes */
864674ae819SStefano Zampini   if (neumann_is) {
8657e5be6c5SStefano Zampini     ierr = VecSet(local_vec,0.0);CHKERRQ(ierr);
866674ae819SStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
86782ba6b80SStefano Zampini     ierr = ISGetLocalSize(neumann_is,&is_size);CHKERRQ(ierr);
868674ae819SStefano Zampini     ierr = ISGetIndices(neumann_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
869674ae819SStefano Zampini     for (i=0;i<is_size;i++) {
870d50ed60aSStefano Zampini       if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
871674ae819SStefano Zampini         array[is_indices[i]] = 1.0;
872674ae819SStefano Zampini       }
8733c629590SStefano Zampini     }
874674ae819SStefano Zampini     ierr = ISRestoreIndices(neumann_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
875674ae819SStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
876674ae819SStefano Zampini     /* Neumann nodes: impose consistency among neighbours */
877674ae819SStefano Zampini     ierr = VecSet(global_vec,0.0);CHKERRQ(ierr);
878674ae819SStefano Zampini     ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
879674ae819SStefano Zampini     ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
880674ae819SStefano Zampini     ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
881674ae819SStefano Zampini     ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
882674ae819SStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
883674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
8843c629590SStefano Zampini       if (PetscRealPart(array[i]) > 0.1) {
8850cf82fd6SStefano Zampini         graph->special_dof[i] = PCBDDCGRAPH_NEUMANN_MARK;
886674ae819SStefano Zampini       }
887674ae819SStefano Zampini     }
888674ae819SStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
8897e5be6c5SStefano Zampini   }
8907e5be6c5SStefano Zampini   /* Take into account Dirichlet nodes (they overwrite any Neumann boundary node previously set) */
891674ae819SStefano Zampini   if (dirichlet_is) {
8927e5be6c5SStefano Zampini     ierr = VecSet(local_vec,0.0);CHKERRQ(ierr);
893674ae819SStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
89482ba6b80SStefano Zampini     ierr = ISGetLocalSize(dirichlet_is,&is_size);CHKERRQ(ierr);
895674ae819SStefano Zampini     ierr = ISGetIndices(dirichlet_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
896674ae819SStefano Zampini     for (i=0;i<is_size;i++){
897d50ed60aSStefano Zampini       if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
8987e5be6c5SStefano Zampini         array[is_indices[i]] = 1.0;
899674ae819SStefano Zampini       }
9003c629590SStefano Zampini     }
901674ae819SStefano Zampini     ierr = ISRestoreIndices(dirichlet_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
902674ae819SStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
903674ae819SStefano Zampini     /* Dirichlet nodes: impose consistency among neighbours */
904674ae819SStefano Zampini     ierr = VecSet(global_vec,0.0);CHKERRQ(ierr);
9057e5be6c5SStefano Zampini     ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9067e5be6c5SStefano Zampini     ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
907674ae819SStefano Zampini     ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
908674ae819SStefano Zampini     ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
909674ae819SStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
910674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
9113c629590SStefano Zampini       if (PetscRealPart(array[i]) > 0.1) {
9124b2aedd3SStefano Zampini         if (commsize > graph->commsizelimit) { /* dirichlet nodes treated as internal */
913df48933bSStefano Zampini           ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
914b3cdcd63SStefano Zampini           graph->subset[i] = 0;
915b3cdcd63SStefano Zampini         }
9160cf82fd6SStefano Zampini         graph->special_dof[i] = PCBDDCGRAPH_DIRICHLET_MARK;
917674ae819SStefano Zampini       }
918674ae819SStefano Zampini     }
919674ae819SStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
9207fb0e2dbSStefano Zampini   }
92108a5cf49SStefano Zampini   /* mark local periodic nodes (if any) and adapt CSR graph (if any) */
92251b0f6cfSStefano Zampini   if (graph->mirrors) {
92351b0f6cfSStefano Zampini     for (i=0;i<graph->nvtxs;i++)
92451b0f6cfSStefano Zampini       if (graph->mirrors[i])
9250cf82fd6SStefano Zampini         graph->special_dof[i] = PCBDDCGRAPH_LOCAL_PERIODIC_MARK;
92651b0f6cfSStefano Zampini 
92708a5cf49SStefano Zampini     if (graph->xadj && graph->adjncy) {
92808a5cf49SStefano Zampini       PetscInt *new_xadj,*new_adjncy;
92951b0f6cfSStefano Zampini       /* sort CSR graph */
93051b0f6cfSStefano Zampini       for (i=0;i<graph->nvtxs;i++)
93151b0f6cfSStefano Zampini         ierr = PetscSortInt(graph->xadj[i+1]-graph->xadj[i],&graph->adjncy[graph->xadj[i]]);CHKERRQ(ierr);
93251b0f6cfSStefano Zampini 
93351b0f6cfSStefano Zampini       /* adapt local CSR graph in case of local periodicity */
93451b0f6cfSStefano Zampini       k = 0;
93551b0f6cfSStefano Zampini       for (i=0;i<graph->nvtxs;i++)
93651b0f6cfSStefano Zampini         for (j=graph->xadj[i];j<graph->xadj[i+1];j++)
93751b0f6cfSStefano Zampini           k += graph->mirrors[graph->adjncy[j]];
93851b0f6cfSStefano Zampini 
939854ce69bSBarry Smith       ierr = PetscMalloc1(graph->nvtxs+1,&new_xadj);CHKERRQ(ierr);
940854ce69bSBarry Smith       ierr = PetscMalloc1(k+graph->xadj[graph->nvtxs],&new_adjncy);CHKERRQ(ierr);
94151b0f6cfSStefano Zampini       new_xadj[0] = 0;
94251b0f6cfSStefano Zampini       for (i=0;i<graph->nvtxs;i++) {
94351b0f6cfSStefano Zampini         k = graph->xadj[i+1]-graph->xadj[i];
94451b0f6cfSStefano Zampini         ierr = PetscMemcpy(&new_adjncy[new_xadj[i]],&graph->adjncy[graph->xadj[i]],k*sizeof(PetscInt));CHKERRQ(ierr);
94551b0f6cfSStefano Zampini         new_xadj[i+1] = new_xadj[i]+k;
94651b0f6cfSStefano Zampini         for (j=graph->xadj[i];j<graph->xadj[i+1];j++) {
94751b0f6cfSStefano Zampini           k = graph->mirrors[graph->adjncy[j]];
94851b0f6cfSStefano Zampini           ierr = PetscMemcpy(&new_adjncy[new_xadj[i+1]],graph->mirrors_set[graph->adjncy[j]],k*sizeof(PetscInt));CHKERRQ(ierr);
94951b0f6cfSStefano Zampini           new_xadj[i+1] += k;
95051b0f6cfSStefano Zampini         }
95151b0f6cfSStefano Zampini         k = new_xadj[i+1]-new_xadj[i];
95251b0f6cfSStefano Zampini         ierr = PetscSortRemoveDupsInt(&k,&new_adjncy[new_xadj[i]]);CHKERRQ(ierr);
95351b0f6cfSStefano Zampini         new_xadj[i+1] = new_xadj[i]+k;
95451b0f6cfSStefano Zampini       }
95551b0f6cfSStefano Zampini       /* set new CSR into graph */
95651b0f6cfSStefano Zampini       ierr = PetscFree(graph->xadj);CHKERRQ(ierr);
95751b0f6cfSStefano Zampini       ierr = PetscFree(graph->adjncy);CHKERRQ(ierr);
95851b0f6cfSStefano Zampini       graph->xadj = new_xadj;
95951b0f6cfSStefano Zampini       graph->adjncy = new_adjncy;
96051b0f6cfSStefano Zampini     }
96108a5cf49SStefano Zampini   }
96251b0f6cfSStefano Zampini 
96317eb1463SStefano Zampini   /* mark special nodes (if any) -> each will become a single node equivalence class */
964674ae819SStefano Zampini   if (custom_primal_vertices) {
96517eb1463SStefano Zampini     ierr = VecSet(local_vec,0.0);CHKERRQ(ierr);
9669b70a373SStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
96763602bcaSStefano Zampini     ierr = ISGetLocalSize(custom_primal_vertices,&is_size);CHKERRQ(ierr);
968674ae819SStefano Zampini     ierr = ISGetIndices(custom_primal_vertices,(const PetscInt**)&is_indices);CHKERRQ(ierr);
969674ae819SStefano Zampini     for (i=0;i<is_size;i++){
9709b70a373SStefano Zampini       if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
9719b70a373SStefano Zampini         array[is_indices[i]] = 1.0;
9729b70a373SStefano Zampini       }
973674ae819SStefano Zampini     }
974674ae819SStefano Zampini     ierr = ISRestoreIndices(custom_primal_vertices,(const PetscInt**)&is_indices);CHKERRQ(ierr);
9759b70a373SStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
9769b70a373SStefano Zampini     /* special nodes: impose consistency among neighbours */
9779b70a373SStefano Zampini     ierr = VecSet(global_vec,0.0);CHKERRQ(ierr);
9789b70a373SStefano Zampini     ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9799b70a373SStefano Zampini     ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9809b70a373SStefano Zampini     ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9819b70a373SStefano Zampini     ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9829b70a373SStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
9839b70a373SStefano Zampini     j = 0;
9849b70a373SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
9859b70a373SStefano Zampini       if (PetscRealPart(array[i]) > 0.1 && graph->special_dof[i] != PCBDDCGRAPH_DIRICHLET_MARK) {
9869b70a373SStefano Zampini         graph->special_dof[i] = PCBDDCGRAPH_SPECIAL_MARK-j;
9879b70a373SStefano Zampini         j++;
9889b70a373SStefano Zampini       }
9899b70a373SStefano Zampini     }
9909b70a373SStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
99117eb1463SStefano Zampini   }
9929b70a373SStefano Zampini 
9934b2aedd3SStefano Zampini   /* mark interior nodes (if commsize > graph->commsizelimit) as touched and belonging to partition number 0 */
9944b2aedd3SStefano Zampini   if (commsize > graph->commsizelimit) {
995674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
996674ae819SStefano Zampini       if (!graph->count[i]) {
997df48933bSStefano Zampini         ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
998674ae819SStefano Zampini         graph->subset[i] = 0;
999674ae819SStefano Zampini       }
1000674ae819SStefano Zampini     }
1001b3cdcd63SStefano Zampini   }
1002b3cdcd63SStefano Zampini 
1003674ae819SStefano Zampini   /* init graph structure and compute default subsets */
1004674ae819SStefano Zampini   nodes_touched = 0;
1005674ae819SStefano Zampini   for (i=0;i<graph->nvtxs;i++) {
1006df48933bSStefano Zampini     if (PetscBTLookup(graph->touched,i)) {
1007674ae819SStefano Zampini       nodes_touched++;
1008674ae819SStefano Zampini     }
1009674ae819SStefano Zampini   }
1010674ae819SStefano Zampini   i = 0;
1011674ae819SStefano Zampini   graph->ncc = 0;
1012674ae819SStefano Zampini   total_counts = 0;
10137babac9bSStefano Zampini 
10147babac9bSStefano Zampini   /* allocated space for queues */
10154b2aedd3SStefano Zampini   if (commsize == graph->commsizelimit) {
10167babac9bSStefano Zampini     ierr = PetscMalloc2(graph->nvtxs+1,&graph->cptr,graph->nvtxs,&graph->queue);CHKERRQ(ierr);
10177babac9bSStefano Zampini   } else {
10187babac9bSStefano Zampini     PetscInt nused = graph->nvtxs - nodes_touched;
10197babac9bSStefano Zampini     ierr = PetscMalloc2(nused+1,&graph->cptr,nused,&graph->queue);CHKERRQ(ierr);
10207babac9bSStefano Zampini   }
10217babac9bSStefano Zampini 
1022674ae819SStefano Zampini   while (nodes_touched<graph->nvtxs) {
1023674ae819SStefano Zampini     /*  find first untouched node in local ordering */
1024b3cdcd63SStefano Zampini     while (PetscBTLookup(graph->touched,i)) i++;
1025df48933bSStefano Zampini     ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
1026674ae819SStefano Zampini     graph->subset[i] = graph->ncc+1;
1027674ae819SStefano Zampini     graph->cptr[graph->ncc] = total_counts;
1028674ae819SStefano Zampini     graph->queue[total_counts] = i;
1029674ae819SStefano Zampini     total_counts++;
1030674ae819SStefano Zampini     nodes_touched++;
1031674ae819SStefano Zampini     /* now find all other nodes having the same set of sharing subdomains */
1032674ae819SStefano Zampini     for (j=i+1;j<graph->nvtxs;j++) {
103374e413f5SStefano Zampini       /* check for same number of sharing subdomains, dof number and same special mark */
1034df48933bSStefano 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]) {
1035674ae819SStefano Zampini         /* check for same set of sharing subdomains */
1036674ae819SStefano Zampini         same_set = PETSC_TRUE;
1037674ae819SStefano Zampini         for (k=0;k<graph->count[j];k++){
1038674ae819SStefano Zampini           if (graph->neighbours_set[i][k] != graph->neighbours_set[j][k]) {
1039674ae819SStefano Zampini             same_set = PETSC_FALSE;
1040674ae819SStefano Zampini           }
1041674ae819SStefano Zampini         }
1042674ae819SStefano Zampini         /* I found a friend of mine */
1043674ae819SStefano Zampini         if (same_set) {
1044df48933bSStefano Zampini           ierr = PetscBTSet(graph->touched,j);CHKERRQ(ierr);
1045b3cdcd63SStefano Zampini           graph->subset[j] = graph->ncc+1;
1046674ae819SStefano Zampini           nodes_touched++;
1047674ae819SStefano Zampini           graph->queue[total_counts] = j;
1048674ae819SStefano Zampini           total_counts++;
1049674ae819SStefano Zampini         }
1050674ae819SStefano Zampini       }
1051674ae819SStefano Zampini     }
1052674ae819SStefano Zampini     graph->ncc++;
1053674ae819SStefano Zampini   }
1054b3cdcd63SStefano 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 */
1055674ae819SStefano Zampini   graph->n_subsets = graph->ncc;
1056785e854fSJed Brown   ierr = PetscMalloc1(graph->n_subsets,&graph->subset_ncc);CHKERRQ(ierr);
1057674ae819SStefano Zampini   for (i=0;i<graph->n_subsets;i++) {
1058674ae819SStefano Zampini     graph->subset_ncc[i] = 1;
1059674ae819SStefano Zampini   }
1060674ae819SStefano Zampini   /* final pointer */
1061674ae819SStefano Zampini   graph->cptr[graph->ncc] = total_counts;
1062ec1c413dSStefano Zampini 
1063b3cdcd63SStefano Zampini   /* For consistency reasons (among neighbours), I need to sort (by global ordering) each connected component */
1064ec1c413dSStefano Zampini   /* Get a reference node (min index in global ordering) for each subset for tagging messages */
1065dc456d91SStefano Zampini   ierr = PetscMalloc1(graph->ncc,&graph->subset_ref_node);CHKERRQ(ierr);
1066dc456d91SStefano Zampini   ierr = PetscMalloc1(graph->cptr[graph->ncc],&queue_global);CHKERRQ(ierr);
10673a5b03d1SStefano Zampini   ierr = ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_global);CHKERRQ(ierr);
1068ec1c413dSStefano Zampini   for (j=0;j<graph->ncc;j++) {
1069ec1c413dSStefano Zampini     ierr = PetscSortIntWithArray(graph->cptr[j+1]-graph->cptr[j],&queue_global[graph->cptr[j]],&graph->queue[graph->cptr[j]]);CHKERRQ(ierr);
1070dc456d91SStefano Zampini     graph->subset_ref_node[j] = graph->queue[graph->cptr[j]];
1071f0321c90SStefano Zampini   }
1072dc456d91SStefano Zampini   ierr = PetscFree(queue_global);CHKERRQ(ierr);
1073acc113dbSStefano Zampini   graph->queue_sorted = PETSC_TRUE;
1074b3cdcd63SStefano Zampini   /* save information on subsets (needed when analyzing the connected components) */
10752f566a09SStefano Zampini   if (graph->ncc) {
1076b3cdcd63SStefano Zampini     ierr = PetscMalloc2(graph->ncc,&graph->subset_size,graph->ncc,&graph->subset_idxs);CHKERRQ(ierr);
1077b3cdcd63SStefano Zampini     ierr = PetscMalloc1(graph->cptr[graph->ncc],&graph->subset_idxs[0]);CHKERRQ(ierr);
1078b3cdcd63SStefano Zampini     ierr = PetscMemzero(graph->subset_idxs[0],graph->cptr[graph->ncc]*sizeof(PetscInt));CHKERRQ(ierr);
1079ec1c413dSStefano Zampini     for (j=1;j<graph->ncc;j++) {
1080b3cdcd63SStefano Zampini       graph->subset_size[j-1] = graph->cptr[j] - graph->cptr[j-1];
1081b3cdcd63SStefano Zampini       graph->subset_idxs[j] = graph->subset_idxs[j-1] + graph->subset_size[j-1];
1082ec1c413dSStefano Zampini     }
1083b3cdcd63SStefano Zampini     graph->subset_size[graph->ncc-1] = graph->cptr[graph->ncc] - graph->cptr[graph->ncc-1];
1084b3cdcd63SStefano Zampini     ierr = PetscMemcpy(graph->subset_idxs[0],graph->queue,graph->cptr[graph->ncc]*sizeof(PetscInt));CHKERRQ(ierr);
1085ec1c413dSStefano Zampini   }
1086b3cdcd63SStefano Zampini 
1087f0321c90SStefano Zampini   /* renumber reference nodes */
1088dc456d91SStefano Zampini   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)(graph->l2gmap)),graph->ncc,graph->subset_ref_node,PETSC_COPY_VALUES,&subset_n);CHKERRQ(ierr);
1089dc456d91SStefano Zampini   ierr = ISLocalToGlobalMappingApplyIS(graph->l2gmap,subset_n,&subset);CHKERRQ(ierr);
1090dc456d91SStefano Zampini   ierr = ISDestroy(&subset_n);CHKERRQ(ierr);
10916583bcc1SStefano Zampini   ierr = ISRenumber(subset,NULL,NULL,&subset_n);CHKERRQ(ierr);
1092dc456d91SStefano Zampini   ierr = ISDestroy(&subset);CHKERRQ(ierr);
1093dc456d91SStefano Zampini   ierr = ISGetLocalSize(subset_n,&k);CHKERRQ(ierr);
10946c4ed002SBarry Smith   if (k != graph->ncc) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid size of new subset! %D != %D",k,graph->ncc);
1095dc456d91SStefano Zampini   ierr = ISGetIndices(subset_n,&is_indices);CHKERRQ(ierr);
1096dc456d91SStefano Zampini   ierr = PetscMemcpy(graph->subset_ref_node,is_indices,graph->ncc*sizeof(PetscInt));CHKERRQ(ierr);
1097dc456d91SStefano Zampini   ierr = ISRestoreIndices(subset_n,&is_indices);CHKERRQ(ierr);
1098dc456d91SStefano Zampini   ierr = ISDestroy(&subset_n);CHKERRQ(ierr);
1099ec1c413dSStefano Zampini 
11003ebc787cSStefano Zampini   /* Determine if we are in 2D or 3D */
11013ebc787cSStefano Zampini   twodim  = PETSC_TRUE;
11023ebc787cSStefano Zampini   for (i=0;i<graph->ncc;i++) {
11033ebc787cSStefano Zampini     PetscInt repdof = graph->queue[graph->cptr[i]];
11043ebc787cSStefano Zampini     if (graph->cptr[i+1]-graph->cptr[i] > graph->custom_minimal_size) {
11053ebc787cSStefano Zampini       if (graph->count[repdof] > 1 || graph->special_dof[repdof] == PCBDDCGRAPH_NEUMANN_MARK) {
11063ebc787cSStefano Zampini         twodim = PETSC_FALSE;
11073ebc787cSStefano Zampini         break;
11083ebc787cSStefano Zampini       }
11093ebc787cSStefano Zampini     }
11103ebc787cSStefano Zampini   }
11113ebc787cSStefano Zampini   ierr = MPIU_Allreduce(&twodim,&graph->twodim,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)graph->l2gmap));CHKERRQ(ierr);
11123ebc787cSStefano Zampini 
1113ec1c413dSStefano Zampini   /* free workspace */
1114674ae819SStefano Zampini   ierr = VecDestroy(&local_vec);CHKERRQ(ierr);
1115674ae819SStefano Zampini   ierr = VecDestroy(&local_vec2);CHKERRQ(ierr);
1116674ae819SStefano Zampini   ierr = VecDestroy(&global_vec);CHKERRQ(ierr);
1117674ae819SStefano Zampini   ierr = VecScatterDestroy(&scatter_ctx);CHKERRQ(ierr);
1118b3cdcd63SStefano Zampini   graph->setupcalled = PETSC_TRUE;
1119674ae819SStefano Zampini   PetscFunctionReturn(0);
1120674ae819SStefano Zampini }
1121674ae819SStefano Zampini 
1122674ae819SStefano Zampini #undef __FUNCT__
1123674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphResetCSR"
1124674ae819SStefano Zampini PetscErrorCode PCBDDCGraphResetCSR(PCBDDCGraph graph)
1125674ae819SStefano Zampini {
1126674ae819SStefano Zampini   PetscErrorCode ierr;
1127674ae819SStefano Zampini 
1128674ae819SStefano Zampini   PetscFunctionBegin;
1129a1dbd327SStefano Zampini   if (graph->freecsr) {
1130674ae819SStefano Zampini     ierr = PetscFree(graph->xadj);CHKERRQ(ierr);
1131674ae819SStefano Zampini     ierr = PetscFree(graph->adjncy);CHKERRQ(ierr);
1132a1dbd327SStefano Zampini   } else {
1133a1dbd327SStefano Zampini     graph->xadj = NULL;
1134a1dbd327SStefano Zampini     graph->adjncy = NULL;
1135a1dbd327SStefano Zampini   }
1136c8272957SStefano Zampini   graph->freecsr = PETSC_FALSE;
1137575ad6abSStefano Zampini   graph->nvtxs_csr = 0;
1138674ae819SStefano Zampini   PetscFunctionReturn(0);
1139674ae819SStefano Zampini }
1140674ae819SStefano Zampini 
1141674ae819SStefano Zampini #undef __FUNCT__
1142674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphReset"
1143674ae819SStefano Zampini PetscErrorCode PCBDDCGraphReset(PCBDDCGraph graph)
1144674ae819SStefano Zampini {
1145674ae819SStefano Zampini   PetscErrorCode ierr;
1146674ae819SStefano Zampini 
1147674ae819SStefano Zampini   PetscFunctionBegin;
1148674ae819SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&graph->l2gmap);CHKERRQ(ierr);
1149674ae819SStefano Zampini   ierr = PetscFree(graph->subset_ncc);CHKERRQ(ierr);
11503a5b03d1SStefano Zampini   ierr = PetscFree(graph->subset_ref_node);CHKERRQ(ierr);
1151674ae819SStefano Zampini   if (graph->nvtxs) {
1152674ae819SStefano Zampini     ierr = PetscFree(graph->neighbours_set[0]);CHKERRQ(ierr);
1153674ae819SStefano Zampini   }
1154df48933bSStefano Zampini   ierr = PetscBTDestroy(&graph->touched);CHKERRQ(ierr);
11557babac9bSStefano Zampini   ierr = PetscFree5(graph->count,
1156674ae819SStefano Zampini                     graph->neighbours_set,
1157674ae819SStefano Zampini                     graph->subset,
1158674ae819SStefano Zampini                     graph->which_dof,
1159df48933bSStefano Zampini                     graph->special_dof);CHKERRQ(ierr);
11607babac9bSStefano Zampini   ierr = PetscFree2(graph->cptr,graph->queue);CHKERRQ(ierr);
116151b0f6cfSStefano Zampini   if (graph->mirrors) {
116251b0f6cfSStefano Zampini     ierr = PetscFree(graph->mirrors_set[0]);CHKERRQ(ierr);
116351b0f6cfSStefano Zampini   }
116451b0f6cfSStefano Zampini   ierr = PetscFree2(graph->mirrors,graph->mirrors_set);CHKERRQ(ierr);
1165b3cdcd63SStefano Zampini   if (graph->subset_idxs) {
1166b3cdcd63SStefano Zampini     ierr = PetscFree(graph->subset_idxs[0]);CHKERRQ(ierr);
1167ec1c413dSStefano Zampini   }
1168b3cdcd63SStefano Zampini   ierr = PetscFree2(graph->subset_size,graph->subset_idxs);CHKERRQ(ierr);
1169a07ea27aSStefano Zampini   ierr = ISDestroy(&graph->dirdofs);CHKERRQ(ierr);
1170d62866d3SStefano Zampini   ierr = ISDestroy(&graph->dirdofsB);CHKERRQ(ierr);
11714f1b2e48SStefano Zampini   if (graph->n_local_subs) {
11724f1b2e48SStefano Zampini     ierr = PetscFree(graph->local_subs);CHKERRQ(ierr);
11734f1b2e48SStefano Zampini   }
1174a07ea27aSStefano Zampini   graph->has_dirichlet       = PETSC_FALSE;
1175674ae819SStefano Zampini   graph->nvtxs               = 0;
11767babac9bSStefano Zampini   graph->nvtxs_global        = 0;
1177674ae819SStefano Zampini   graph->n_subsets           = 0;
1178674ae819SStefano Zampini   graph->custom_minimal_size = 1;
11794f1b2e48SStefano Zampini   graph->n_local_subs        = 0;
1180be12c134Sstefano_zampini   graph->maxcount            = PETSC_MAX_INT;
1181fb597685SStefano Zampini   graph->setupcalled         = PETSC_FALSE;
1182674ae819SStefano Zampini   PetscFunctionReturn(0);
1183674ae819SStefano Zampini }
1184674ae819SStefano Zampini 
1185674ae819SStefano Zampini #undef __FUNCT__
1186674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphInit"
1187be12c134Sstefano_zampini PetscErrorCode PCBDDCGraphInit(PCBDDCGraph graph, ISLocalToGlobalMapping l2gmap, PetscInt N, PetscInt maxcount)
1188674ae819SStefano Zampini {
1189674ae819SStefano Zampini   PetscInt       n;
1190674ae819SStefano Zampini   PetscErrorCode ierr;
1191674ae819SStefano Zampini 
1192674ae819SStefano Zampini   PetscFunctionBegin;
1193674ae819SStefano Zampini   PetscValidPointer(graph,1);
1194674ae819SStefano Zampini   PetscValidHeaderSpecific(l2gmap,IS_LTOGM_CLASSID,2);
11957babac9bSStefano Zampini   PetscValidLogicalCollectiveInt(l2gmap,N,3);
1196be12c134Sstefano_zampini   PetscValidLogicalCollectiveInt(l2gmap,maxcount,4);
11978e61c736SStefano Zampini   /* raise an error if already allocated */
11986c4ed002SBarry Smith   if (graph->nvtxs_global) SETERRQ(PetscObjectComm((PetscObject)l2gmap),PETSC_ERR_PLIB,"BDDCGraph already initialized");
1199674ae819SStefano Zampini   /* set number of vertices */
1200674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)l2gmap);CHKERRQ(ierr);
1201674ae819SStefano Zampini   graph->l2gmap = l2gmap;
1202674ae819SStefano Zampini   ierr = ISLocalToGlobalMappingGetSize(l2gmap,&n);CHKERRQ(ierr);
1203674ae819SStefano Zampini   graph->nvtxs = n;
12047fb0e2dbSStefano Zampini   graph->nvtxs_global = N;
1205674ae819SStefano Zampini   /* allocate used space */
1206df48933bSStefano Zampini   ierr = PetscBTCreate(graph->nvtxs,&graph->touched);CHKERRQ(ierr);
12077babac9bSStefano Zampini   ierr = PetscMalloc5(graph->nvtxs,&graph->count,
12088e5aaad5SJed Brown                       graph->nvtxs,&graph->neighbours_set,
12098e5aaad5SJed Brown                       graph->nvtxs,&graph->subset,
12108e5aaad5SJed Brown                       graph->nvtxs,&graph->which_dof,
12118e5aaad5SJed Brown                       graph->nvtxs,&graph->special_dof);CHKERRQ(ierr);
1212674ae819SStefano Zampini   /* zeroes memory */
1213674ae819SStefano Zampini   ierr = PetscMemzero(graph->count,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
1214674ae819SStefano Zampini   ierr = PetscMemzero(graph->subset,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
121563602bcaSStefano Zampini   /* use -1 as a default value for which_dof array */
121663602bcaSStefano Zampini   for (n=0;n<graph->nvtxs;n++) graph->which_dof[n] = -1;
121774e413f5SStefano Zampini   ierr = PetscMemzero(graph->special_dof,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
1218674ae819SStefano Zampini   /* zeroes first pointer to neighbour set */
1219674ae819SStefano Zampini   if (graph->nvtxs) {
1220674ae819SStefano Zampini     graph->neighbours_set[0] = 0;
1221674ae819SStefano Zampini   }
1222674ae819SStefano Zampini   /* zeroes workspace for values of ncc */
1223674ae819SStefano Zampini   graph->subset_ncc = 0;
12243a5b03d1SStefano Zampini   graph->subset_ref_node = 0;
1225be12c134Sstefano_zampini   /* maxcount for cc */
1226be12c134Sstefano_zampini   graph->maxcount = maxcount;
1227674ae819SStefano Zampini   PetscFunctionReturn(0);
1228674ae819SStefano Zampini }
1229674ae819SStefano Zampini 
1230674ae819SStefano Zampini #undef __FUNCT__
1231674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphDestroy"
1232674ae819SStefano Zampini PetscErrorCode PCBDDCGraphDestroy(PCBDDCGraph* graph)
1233674ae819SStefano Zampini {
1234674ae819SStefano Zampini   PetscErrorCode ierr;
1235674ae819SStefano Zampini 
1236674ae819SStefano Zampini   PetscFunctionBegin;
1237f5fcb376SStefano Zampini   ierr = PCBDDCGraphResetCSR(*graph);CHKERRQ(ierr);
1238674ae819SStefano Zampini   ierr = PCBDDCGraphReset(*graph);CHKERRQ(ierr);
1239674ae819SStefano Zampini   ierr = PetscFree(*graph);CHKERRQ(ierr);
1240674ae819SStefano Zampini   PetscFunctionReturn(0);
1241674ae819SStefano Zampini }
1242674ae819SStefano Zampini 
1243674ae819SStefano Zampini #undef __FUNCT__
1244674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphCreate"
1245674ae819SStefano Zampini PetscErrorCode PCBDDCGraphCreate(PCBDDCGraph *graph)
1246674ae819SStefano Zampini {
1247674ae819SStefano Zampini   PCBDDCGraph    new_graph;
1248674ae819SStefano Zampini   PetscErrorCode ierr;
1249674ae819SStefano Zampini 
1250674ae819SStefano Zampini   PetscFunctionBegin;
1251854ce69bSBarry Smith   ierr = PetscNew(&new_graph);CHKERRQ(ierr);
1252674ae819SStefano Zampini   new_graph->custom_minimal_size = 1;
12534b2aedd3SStefano Zampini   new_graph->commsizelimit = 1;
1254674ae819SStefano Zampini   *graph = new_graph;
1255674ae819SStefano Zampini   PetscFunctionReturn(0);
1256674ae819SStefano Zampini }
1257