xref: /petsc/src/ksp/pc/impls/bddc/bddcgraph.c (revision b3cdcd63525199af32e7b457f094c7aca4893208)
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);
81e49050b4SStefano Zampini   if (verbosity_level > 1) {
82674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
83674ae819SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"%d:\n",i);CHKERRQ(ierr);
84674ae819SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   which_dof: %d\n",graph->which_dof[i]);CHKERRQ(ierr);
8574e413f5SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   special_dof: %d\n",graph->special_dof[i]);CHKERRQ(ierr);
86674ae819SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   neighbours: %d\n",graph->count[i]);CHKERRQ(ierr);
872b510759SStefano Zampini       ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
88674ae819SStefano Zampini       if (graph->count[i]) {
89674ae819SStefano Zampini         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"     set of neighbours:");CHKERRQ(ierr);
90674ae819SStefano Zampini         for (j=0;j<graph->count[i];j++) {
91674ae819SStefano Zampini           ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d",graph->neighbours_set[i][j]);CHKERRQ(ierr);
92674ae819SStefano Zampini         }
93674ae819SStefano Zampini         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
94674ae819SStefano Zampini       }
952b510759SStefano Zampini       ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr);
962b510759SStefano Zampini       ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
9751b0f6cfSStefano Zampini       if (graph->mirrors) {
9851b0f6cfSStefano Zampini         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   mirrors: %d\n",graph->mirrors[i]);CHKERRQ(ierr);
9951b0f6cfSStefano Zampini         if (graph->mirrors[i]) {
1002b510759SStefano Zampini           ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
10151b0f6cfSStefano Zampini           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"     set of mirrors:");CHKERRQ(ierr);
10251b0f6cfSStefano Zampini           for (j=0;j<graph->mirrors[i];j++) {
10351b0f6cfSStefano Zampini             ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d",graph->mirrors_set[i][j]);CHKERRQ(ierr);
10451b0f6cfSStefano Zampini           }
10551b0f6cfSStefano Zampini           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
1062b510759SStefano Zampini           ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr);
1072b510759SStefano Zampini           ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
10851b0f6cfSStefano Zampini         }
10951b0f6cfSStefano Zampini       }
110e49050b4SStefano Zampini       if (verbosity_level > 2) {
111674ae819SStefano Zampini         if (graph->xadj && graph->adjncy) {
112674ae819SStefano Zampini           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   local adj list:");CHKERRQ(ierr);
1132b510759SStefano Zampini           ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
114674ae819SStefano Zampini           for (j=graph->xadj[i];j<graph->xadj[i+1];j++) {
115674ae819SStefano Zampini             ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d",graph->adjncy[j]);CHKERRQ(ierr);
116674ae819SStefano Zampini           }
117674ae819SStefano Zampini           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
1182b510759SStefano Zampini           ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr);
1192b510759SStefano Zampini           ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
120ec1c413dSStefano Zampini         } else {
121ec1c413dSStefano Zampini           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   no adj info\n");CHKERRQ(ierr);
122674ae819SStefano Zampini         }
123e49050b4SStefano Zampini       }
124531609faSStefano Zampini       if (graph->n_local_subs) {
125531609faSStefano Zampini         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   local sub id: %d\n",graph->local_subs[i]);CHKERRQ(ierr);
126531609faSStefano Zampini       }
127674ae819SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   interface subset id: %d\n",graph->subset[i]);CHKERRQ(ierr);
128674ae819SStefano Zampini       if (graph->subset[i] && graph->subset_ncc) {
129674ae819SStefano Zampini         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   ncc for subset: %d\n",graph->subset_ncc[graph->subset[i]-1]);CHKERRQ(ierr);
130674ae819SStefano Zampini       }
131674ae819SStefano Zampini     }
132e49050b4SStefano Zampini   }
133674ae819SStefano Zampini   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Total number of connected components %d\n",graph->ncc);CHKERRQ(ierr);
134785e854fSJed Brown   ierr = PetscMalloc1(graph->cptr[graph->ncc],&queue_in_global_numbering);CHKERRQ(ierr);
13593fb5fd3SStefano Zampini   ierr = ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_in_global_numbering);CHKERRQ(ierr);
136674ae819SStefano Zampini   for (i=0;i<graph->ncc;i++) {
1371ce3d396SStefano Zampini     PetscInt node_num=graph->queue[graph->cptr[i]];
1381ce3d396SStefano Zampini     PetscBool printcc = PETSC_FALSE;
139674ae819SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  %d (neighs:",i);CHKERRQ(ierr);
1402b510759SStefano Zampini     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
141674ae819SStefano Zampini     for (j=0;j<graph->count[node_num];j++) {
142674ae819SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d",graph->neighbours_set[node_num][j]);CHKERRQ(ierr);
143674ae819SStefano Zampini     }
144674ae819SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"):");CHKERRQ(ierr);
145e2291ee8SStefano Zampini     if (verbosity_level > 1 || graph->twodim) {
1461ce3d396SStefano Zampini       printcc = PETSC_TRUE;
147e635b14bSStefano Zampini     } else if (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     }
155674ae819SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
1562b510759SStefano Zampini     ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr);
1572b510759SStefano Zampini     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
158674ae819SStefano Zampini   }
15993fb5fd3SStefano Zampini   ierr = PetscFree(queue_in_global_numbering);CHKERRQ(ierr);
160674ae819SStefano Zampini   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
161674ae819SStefano Zampini   PetscFunctionReturn(0);
162674ae819SStefano Zampini }
163674ae819SStefano Zampini 
164674ae819SStefano Zampini #undef __FUNCT__
165674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphGetCandidatesIS"
166a873d5d0SStefano Zampini PetscErrorCode PCBDDCGraphGetCandidatesIS(PCBDDCGraph graph, PetscInt *n_faces, IS *FacesIS[], PetscInt *n_edges, IS *EdgesIS[], IS *VerticesIS)
167674ae819SStefano Zampini {
168674ae819SStefano Zampini   IS             *ISForFaces,*ISForEdges,ISForVertices;
169adfc4fb2SStefano Zampini   PetscInt       i,nfc,nec,nvc,*idx;
170674ae819SStefano Zampini   PetscErrorCode ierr;
171674ae819SStefano Zampini 
172674ae819SStefano Zampini   PetscFunctionBegin;
173674ae819SStefano Zampini   /* loop on ccs to evalute number of faces, edges and vertices */
174674ae819SStefano Zampini   nfc = 0;
175674ae819SStefano Zampini   nec = 0;
176674ae819SStefano Zampini   nvc = 0;
177674ae819SStefano Zampini   for (i=0;i<graph->ncc;i++) {
1781e1f2d93SStefano Zampini     PetscInt repdof = graph->queue[graph->cptr[i]];
179674ae819SStefano Zampini     if (graph->cptr[i+1]-graph->cptr[i] > graph->custom_minimal_size) {
180e635b14bSStefano Zampini       if (graph->count[repdof] == 1 && graph->special_dof[repdof] != PCBDDCGRAPH_NEUMANN_MARK) {
181674ae819SStefano Zampini         nfc++;
182674ae819SStefano Zampini       } else { /* note that nec will be zero in 2d */
183674ae819SStefano Zampini         nec++;
184674ae819SStefano Zampini       }
185674ae819SStefano Zampini     } else {
186674ae819SStefano Zampini       nvc += graph->cptr[i+1]-graph->cptr[i];
187674ae819SStefano Zampini     }
188674ae819SStefano Zampini   }
189adfc4fb2SStefano Zampini 
1908e4af1ccSStefano Zampini   /* check if we are in 2D or 3D */
1918e4af1ccSStefano Zampini   if (graph->twodim) { /* we are in a 2D case -> edges are shared by 2 subregions and faces don't exist */
192674ae819SStefano Zampini     nec = nfc;
193674ae819SStefano Zampini     nfc = 0;
194674ae819SStefano Zampini   }
195adfc4fb2SStefano Zampini 
196674ae819SStefano Zampini   /* allocate IS arrays for faces, edges. Vertices need a single index set. */
197cf5a6209SStefano Zampini   if (FacesIS) {
198785e854fSJed Brown     ierr = PetscMalloc1(nfc,&ISForFaces);CHKERRQ(ierr);
199cf5a6209SStefano Zampini   }
200cf5a6209SStefano Zampini   if (EdgesIS) {
201785e854fSJed Brown     ierr = PetscMalloc1(nec,&ISForEdges);CHKERRQ(ierr);
202cf5a6209SStefano Zampini   }
203cf5a6209SStefano Zampini   if (VerticesIS) {
204785e854fSJed Brown     ierr = PetscMalloc1(nvc,&idx);CHKERRQ(ierr);
205cf5a6209SStefano Zampini   }
206cf5a6209SStefano Zampini 
207674ae819SStefano Zampini   /* loop on ccs to compute index sets for faces and edges */
208acc113dbSStefano Zampini   if (!graph->queue_sorted) {
209acc113dbSStefano Zampini     PetscInt *queue_global;
210acc113dbSStefano Zampini 
211acc113dbSStefano Zampini     ierr = PetscMalloc1(graph->cptr[graph->ncc],&queue_global);CHKERRQ(ierr);
212acc113dbSStefano Zampini     ierr = ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_global);CHKERRQ(ierr);
213acc113dbSStefano Zampini     for (i=0;i<graph->ncc;i++) {
214acc113dbSStefano Zampini       ierr = PetscSortIntWithArray(graph->cptr[i+1]-graph->cptr[i],&queue_global[graph->cptr[i]],&graph->queue[graph->cptr[i]]);CHKERRQ(ierr);
215acc113dbSStefano Zampini     }
216acc113dbSStefano Zampini     ierr = PetscFree(queue_global);CHKERRQ(ierr);
217acc113dbSStefano Zampini     graph->queue_sorted = PETSC_TRUE;
218acc113dbSStefano Zampini   }
219674ae819SStefano Zampini   nfc = 0;
220674ae819SStefano Zampini   nec = 0;
221674ae819SStefano Zampini   for (i=0;i<graph->ncc;i++) {
2221e1f2d93SStefano Zampini     PetscInt repdof = graph->queue[graph->cptr[i]];
223674ae819SStefano Zampini     if (graph->cptr[i+1]-graph->cptr[i] > graph->custom_minimal_size) {
224e635b14bSStefano Zampini       if (graph->count[repdof] == 1 && graph->special_dof[repdof] != PCBDDCGRAPH_NEUMANN_MARK) {
2258e4af1ccSStefano Zampini         if (graph->twodim) {
226cf5a6209SStefano Zampini           if (EdgesIS) {
227674ae819SStefano Zampini             ierr = ISCreateGeneral(PETSC_COMM_SELF,graph->cptr[i+1]-graph->cptr[i],&graph->queue[graph->cptr[i]],PETSC_COPY_VALUES,&ISForEdges[nec]);CHKERRQ(ierr);
228cf5a6209SStefano Zampini           }
229674ae819SStefano Zampini           nec++;
230674ae819SStefano Zampini         } else {
231cf5a6209SStefano Zampini           if (FacesIS) {
232674ae819SStefano Zampini             ierr = ISCreateGeneral(PETSC_COMM_SELF,graph->cptr[i+1]-graph->cptr[i],&graph->queue[graph->cptr[i]],PETSC_COPY_VALUES,&ISForFaces[nfc]);CHKERRQ(ierr);
233cf5a6209SStefano Zampini           }
234674ae819SStefano Zampini           nfc++;
235674ae819SStefano Zampini         }
236674ae819SStefano Zampini       } else {
237cf5a6209SStefano Zampini         if (EdgesIS) {
238674ae819SStefano Zampini           ierr = ISCreateGeneral(PETSC_COMM_SELF,graph->cptr[i+1]-graph->cptr[i],&graph->queue[graph->cptr[i]],PETSC_COPY_VALUES,&ISForEdges[nec]);CHKERRQ(ierr);
239cf5a6209SStefano Zampini         }
240674ae819SStefano Zampini         nec++;
241674ae819SStefano Zampini       }
242674ae819SStefano Zampini     }
243674ae819SStefano Zampini   }
244674ae819SStefano Zampini   /* index set for vertices */
245cf5a6209SStefano Zampini   if (VerticesIS) {
246b8ffe317SStefano Zampini     nvc = 0;
247674ae819SStefano Zampini     for (i=0;i<graph->ncc;i++) {
248674ae819SStefano Zampini       if (graph->cptr[i+1]-graph->cptr[i] <= graph->custom_minimal_size) {
249adfc4fb2SStefano Zampini         PetscInt j;
250adfc4fb2SStefano Zampini 
251674ae819SStefano Zampini         for (j=graph->cptr[i];j<graph->cptr[i+1];j++) {
252674ae819SStefano Zampini           idx[nvc]=graph->queue[j];
253674ae819SStefano Zampini           nvc++;
254674ae819SStefano Zampini         }
255674ae819SStefano Zampini       }
256674ae819SStefano Zampini     }
257674ae819SStefano Zampini     /* sort vertex set (by local ordering) */
258674ae819SStefano Zampini     ierr = PetscSortInt(nvc,idx);CHKERRQ(ierr);
259674ae819SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,nvc,idx,PETSC_OWN_POINTER,&ISForVertices);CHKERRQ(ierr);
260674ae819SStefano Zampini   }
261674ae819SStefano Zampini   /* get back info */
262a873d5d0SStefano Zampini   if (n_faces) *n_faces = nfc;
263d7796978SStefano Zampini   if (FacesIS) {
264674ae819SStefano Zampini     *FacesIS = ISForFaces;
265d7796978SStefano Zampini   }
266a873d5d0SStefano Zampini   if (n_edges) *n_edges = nec;
267d7796978SStefano Zampini   if (EdgesIS) {
268674ae819SStefano Zampini     *EdgesIS = ISForEdges;
269d7796978SStefano Zampini   }
270d7796978SStefano Zampini   if (VerticesIS) {
271674ae819SStefano Zampini     *VerticesIS = ISForVertices;
272d7796978SStefano Zampini   }
273674ae819SStefano Zampini   PetscFunctionReturn(0);
274674ae819SStefano Zampini }
275674ae819SStefano Zampini 
276674ae819SStefano Zampini #undef __FUNCT__
277674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphComputeConnectedComponents"
278674ae819SStefano Zampini PetscErrorCode PCBDDCGraphComputeConnectedComponents(PCBDDCGraph graph)
279674ae819SStefano Zampini {
2804a060362SStefano Zampini   PetscBool      adapt_interface_reduced;
281674ae819SStefano Zampini   MPI_Comm       interface_comm;
2824a060362SStefano Zampini   PetscMPIInt    size;
2838e4af1ccSStefano Zampini   PetscInt       i;
2848e4af1ccSStefano Zampini   PetscBool      twodim;
285674ae819SStefano Zampini   PetscErrorCode ierr;
286674ae819SStefano Zampini 
287674ae819SStefano Zampini   PetscFunctionBegin;
288674ae819SStefano Zampini   /* compute connected components locally */
289674ae819SStefano Zampini   ierr = PetscObjectGetComm((PetscObject)(graph->l2gmap),&interface_comm);CHKERRQ(ierr);
290674ae819SStefano Zampini   ierr = PCBDDCGraphComputeConnectedComponentsLocal(graph);CHKERRQ(ierr);
291674ae819SStefano Zampini   /* check consistency of connected components among neighbouring subdomains -> it adapt them in case it is needed */
2924a060362SStefano Zampini   ierr = MPI_Comm_size(interface_comm,&size);CHKERRQ(ierr);
2934a060362SStefano Zampini   adapt_interface_reduced = PETSC_FALSE;
2944a060362SStefano Zampini   if (size > 1) {
2954a060362SStefano Zampini     PetscInt i;
2964a060362SStefano Zampini     PetscBool adapt_interface = PETSC_FALSE;
297674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
298674ae819SStefano Zampini       /* We are not sure that on a given subset of the local interface,
299674ae819SStefano Zampini          with two connected components, the latters be the same among sharing subdomains */
300674ae819SStefano Zampini       if (graph->subset_ncc[i] > 1) {
3014a060362SStefano Zampini         adapt_interface = PETSC_TRUE;
302674ae819SStefano Zampini         break;
303674ae819SStefano Zampini       }
304674ae819SStefano Zampini     }
305b2566f29SBarry Smith     ierr = MPIU_Allreduce(&adapt_interface,&adapt_interface_reduced,1,MPIU_BOOL,MPI_LOR,interface_comm);CHKERRQ(ierr);
3064a060362SStefano Zampini   }
307674ae819SStefano Zampini 
308674ae819SStefano Zampini   if (graph->n_subsets && adapt_interface_reduced) {
309ec1c413dSStefano Zampini     PetscBT     subset_cc_adapt;
310ec1c413dSStefano Zampini     MPI_Request *send_requests,*recv_requests;
311ec1c413dSStefano Zampini     PetscInt    *send_buffer,*recv_buffer;
312ec1c413dSStefano Zampini     PetscInt    sum_requests,start_of_recv,start_of_send;
313ec1c413dSStefano Zampini     PetscInt    *cum_recv_counts;
314ec1c413dSStefano Zampini     PetscInt    *labels;
315*b3cdcd63SStefano Zampini     PetscInt    ncc,cum_queue,mss,mns,j,k,s;
316*b3cdcd63SStefano Zampini     PetscInt    **refine_buffer,*private_labels;
3175b1b9aeaSStefano Zampini 
318ec1c413dSStefano Zampini     ierr = PetscMalloc1(graph->nvtxs,&labels);CHKERRQ(ierr);
319ec1c413dSStefano Zampini     ierr = PetscMemzero(labels,graph->nvtxs*sizeof(*labels));CHKERRQ(ierr);
320*b3cdcd63SStefano Zampini     for (i=0;i<graph->ncc;i++)
321*b3cdcd63SStefano Zampini       for (j=graph->cptr[i];j<graph->cptr[i+1];j++)
322*b3cdcd63SStefano Zampini         labels[graph->queue[j]] = i;
323*b3cdcd63SStefano Zampini 
324674ae819SStefano Zampini     /* allocate some space */
325854ce69bSBarry Smith     ierr = PetscMalloc1(graph->n_subsets+1,&cum_recv_counts);CHKERRQ(ierr);
326674ae819SStefano Zampini     ierr = PetscMemzero(cum_recv_counts,(graph->n_subsets+1)*sizeof(*cum_recv_counts));CHKERRQ(ierr);
327*b3cdcd63SStefano Zampini 
328674ae819SStefano Zampini     /* first count how many neighbours per connected component I will receive from */
329674ae819SStefano Zampini     cum_recv_counts[0] = 0;
330*b3cdcd63SStefano 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]];
331ec1c413dSStefano Zampini     ierr = PetscMalloc1(cum_recv_counts[graph->n_subsets],&recv_buffer);CHKERRQ(ierr);
332ec1c413dSStefano Zampini     ierr = PetscMalloc2(cum_recv_counts[graph->n_subsets],&send_requests,cum_recv_counts[graph->n_subsets],&recv_requests);CHKERRQ(ierr);
333674ae819SStefano Zampini     for (i=0;i<cum_recv_counts[graph->n_subsets];i++) {
334674ae819SStefano Zampini       send_requests[i] = MPI_REQUEST_NULL;
335674ae819SStefano Zampini       recv_requests[i] = MPI_REQUEST_NULL;
336674ae819SStefano Zampini     }
337*b3cdcd63SStefano Zampini 
338ec1c413dSStefano Zampini     /* exchange with my neighbours the number of my connected components on the subset of interface */
339674ae819SStefano Zampini     sum_requests = 0;
340674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
341ec1c413dSStefano Zampini       PetscMPIInt neigh,tag;
342*b3cdcd63SStefano Zampini       PetscInt    count,*neighs;
343ec1c413dSStefano Zampini 
344*b3cdcd63SStefano Zampini       count = graph->count[graph->subset_idxs[i][0]];
345*b3cdcd63SStefano Zampini       neighs = graph->neighbours_set[graph->subset_idxs[i][0]];
346ec1c413dSStefano Zampini       ierr = PetscMPIIntCast(2*graph->subset_ref_node[i],&tag);CHKERRQ(ierr);
347*b3cdcd63SStefano Zampini       for (k=0;k<count;k++) {
348*b3cdcd63SStefano Zampini         ierr = PetscMPIIntCast(neighs[k],&neigh);CHKERRQ(ierr);
349674ae819SStefano Zampini         ierr = MPI_Isend(&graph->subset_ncc[i],1,MPIU_INT,neigh,tag,interface_comm,&send_requests[sum_requests]);CHKERRQ(ierr);
350ec1c413dSStefano Zampini         ierr = MPI_Irecv(&recv_buffer[sum_requests],1,MPIU_INT,neigh,tag,interface_comm,&recv_requests[sum_requests]);CHKERRQ(ierr);
351674ae819SStefano Zampini         sum_requests++;
352674ae819SStefano Zampini       }
353674ae819SStefano Zampini     }
354674ae819SStefano Zampini     ierr = MPI_Waitall(sum_requests,recv_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
355674ae819SStefano Zampini     ierr = MPI_Waitall(sum_requests,send_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
356*b3cdcd63SStefano Zampini 
357*b3cdcd63SStefano Zampini     /* determine the subsets I have to adapt (those having more than 1 cc) */
358ec1c413dSStefano Zampini     ierr = PetscBTCreate(graph->n_subsets,&subset_cc_adapt);CHKERRQ(ierr);
359ec1c413dSStefano Zampini     ierr = PetscBTMemzero(graph->n_subsets,subset_cc_adapt);CHKERRQ(ierr);
360674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
361*b3cdcd63SStefano Zampini       if (graph->subset_ncc[i] > 1) {
362*b3cdcd63SStefano Zampini         ierr = PetscBTSet(subset_cc_adapt,i);CHKERRQ(ierr);
363*b3cdcd63SStefano Zampini         continue;
364*b3cdcd63SStefano Zampini       }
365674ae819SStefano Zampini       for (j=cum_recv_counts[i];j<cum_recv_counts[i+1];j++){
366*b3cdcd63SStefano Zampini          if (recv_buffer[j] > 1) {
367ec1c413dSStefano Zampini           ierr = PetscBTSet(subset_cc_adapt,i);
368674ae819SStefano Zampini           break;
369674ae819SStefano Zampini         }
370674ae819SStefano Zampini       }
371674ae819SStefano Zampini     }
372ec1c413dSStefano Zampini     ierr = PetscFree(recv_buffer);CHKERRQ(ierr);
373*b3cdcd63SStefano Zampini 
374ec1c413dSStefano Zampini     /* determine send/recv buffers sizes */
375ec1c413dSStefano Zampini     j = 0;
376*b3cdcd63SStefano Zampini     mss = 0;
377674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
378ec1c413dSStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i)) {
379*b3cdcd63SStefano Zampini         j += graph->subset_size[i];
380*b3cdcd63SStefano Zampini         mss = PetscMax(graph->subset_size[i],mss);
381674ae819SStefano Zampini       }
382674ae819SStefano Zampini     }
383ec1c413dSStefano Zampini     k = 0;
384*b3cdcd63SStefano Zampini     mns = 0;
385674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
386ec1c413dSStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i)) {
387*b3cdcd63SStefano Zampini         k += (cum_recv_counts[i+1]-cum_recv_counts[i])*graph->subset_size[i];
388*b3cdcd63SStefano Zampini         mns = PetscMax(cum_recv_counts[i+1]-cum_recv_counts[i],mns);
389674ae819SStefano Zampini       }
390674ae819SStefano Zampini     }
391ec1c413dSStefano Zampini     ierr = PetscMalloc2(j,&send_buffer,k,&recv_buffer);CHKERRQ(ierr);
392ec1c413dSStefano Zampini 
393*b3cdcd63SStefano Zampini     /* fill send buffer (order matters: subset_idxs ordered by global ordering) */
394ec1c413dSStefano Zampini     j = 0;
395*b3cdcd63SStefano Zampini     for (i=0;i<graph->n_subsets;i++)
396*b3cdcd63SStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i))
397*b3cdcd63SStefano Zampini         for (k=0;k<graph->subset_size[i];k++)
398*b3cdcd63SStefano Zampini           send_buffer[j++] = labels[graph->subset_idxs[i][k]];
399ec1c413dSStefano Zampini 
400674ae819SStefano Zampini     /* now exchange the data */
401674ae819SStefano Zampini     start_of_recv = 0;
402674ae819SStefano Zampini     start_of_send = 0;
403674ae819SStefano Zampini     sum_requests = 0;
404674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
405ec1c413dSStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i)) {
406ec1c413dSStefano Zampini         PetscMPIInt neigh,tag;
407*b3cdcd63SStefano Zampini         PetscInt    size_of_send = graph->subset_size[i];
408ec1c413dSStefano Zampini 
409*b3cdcd63SStefano Zampini         j = graph->subset_idxs[i][0];
410ec1c413dSStefano Zampini         ierr = PetscMPIIntCast(2*graph->subset_ref_node[i]+1,&tag);CHKERRQ(ierr);
411674ae819SStefano Zampini         for (k=0;k<graph->count[j];k++) {
412674ae819SStefano Zampini           ierr = PetscMPIIntCast(graph->neighbours_set[j][k],&neigh);CHKERRQ(ierr);
413ec1c413dSStefano Zampini           ierr = MPI_Isend(&send_buffer[start_of_send],size_of_send,MPIU_INT,neigh,tag,interface_comm,&send_requests[sum_requests]);CHKERRQ(ierr);
414ec1c413dSStefano Zampini           ierr = MPI_Irecv(&recv_buffer[start_of_recv],size_of_send,MPIU_INT,neigh,tag,interface_comm,&recv_requests[sum_requests]);CHKERRQ(ierr);
415ec1c413dSStefano Zampini           start_of_recv += size_of_send;
416674ae819SStefano Zampini           sum_requests++;
417674ae819SStefano Zampini         }
418674ae819SStefano Zampini         start_of_send += size_of_send;
419674ae819SStefano Zampini       }
420674ae819SStefano Zampini     }
421674ae819SStefano Zampini     ierr = MPI_Waitall(sum_requests,recv_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
422d52c4852SStefano Zampini 
423*b3cdcd63SStefano Zampini     /* refine connected components */
424674ae819SStefano Zampini     start_of_recv = 0;
425*b3cdcd63SStefano Zampini     /* allocate some temporary space */
426*b3cdcd63SStefano Zampini     if (mss) {
427*b3cdcd63SStefano Zampini       ierr = PetscMalloc1(mss,&refine_buffer);CHKERRQ(ierr);
428*b3cdcd63SStefano Zampini       ierr = PetscMalloc2(mss*(mns+1),&refine_buffer[0],mss,&private_labels);CHKERRQ(ierr);
429*b3cdcd63SStefano Zampini     }
430*b3cdcd63SStefano Zampini     ncc = 0;
431*b3cdcd63SStefano Zampini     cum_queue = 0;
432*b3cdcd63SStefano Zampini     graph->cptr[0] = 0;
433674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
434ec1c413dSStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i)) {
435*b3cdcd63SStefano Zampini         PetscInt subset_counter = 0;
436*b3cdcd63SStefano Zampini         PetscInt sharingprocs = cum_recv_counts[i+1]-cum_recv_counts[i]+1; /* count myself */
437*b3cdcd63SStefano Zampini         PetscInt buffer_size = graph->subset_size[i];
438ec1c413dSStefano Zampini 
439*b3cdcd63SStefano Zampini         /* compute pointers */
440*b3cdcd63SStefano Zampini         for (j=1;j<buffer_size;j++) refine_buffer[j] = refine_buffer[j-1] + sharingprocs;
441*b3cdcd63SStefano Zampini         /* analyze contributions from subdomains that share the i-th subset
442*b3cdcd63SStefano Zampini            The stricture of refine_buffer is suitable to find intersections of ccs among sharingprocs.
443*b3cdcd63SStefano 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)
444*b3cdcd63SStefano Zampini            sharing procs connected components:
445ec1c413dSStefano Zampini              neigh 0: [0 1 4], [2 3], labels [4,7]  (2 connected components)
446ec1c413dSStefano Zampini              neigh 1: [0 1], [2 3 4], labels [3 2]  (2 connected components)
447ec1c413dSStefano Zampini              neigh 2: [0 4], [1], [2 3], labels [1 5 6] (3 connected components)
448*b3cdcd63SStefano Zampini            refine_buffer will be filled as:
449ec1c413dSStefano Zampini              [ 4, 3, 1;
450ec1c413dSStefano Zampini                4, 2, 1;
451ec1c413dSStefano Zampini                7, 2, 6;
452ec1c413dSStefano Zampini                4, 3, 5;
453ec1c413dSStefano Zampini                7, 2, 6; ];
454*b3cdcd63SStefano Zampini            The connected components in local ordering are [0], [1], [2 3], [4] */
455*b3cdcd63SStefano Zampini         /* fill temp_buffer */
456*b3cdcd63SStefano Zampini         for (k=0;k<buffer_size;k++) refine_buffer[k][0] = labels[graph->subset_idxs[i][k]];
457*b3cdcd63SStefano Zampini         for (j=0;j<sharingprocs-1;j++) {
458*b3cdcd63SStefano Zampini           for (k=0;k<buffer_size;k++) refine_buffer[k][j+1] = recv_buffer[start_of_recv+k];
459*b3cdcd63SStefano Zampini           start_of_recv += buffer_size;
460674ae819SStefano Zampini         }
461*b3cdcd63SStefano Zampini         ierr = PetscMemzero(private_labels,buffer_size*sizeof(PetscInt));CHKERRQ(ierr);
462*b3cdcd63SStefano Zampini         for (j=0;j<buffer_size;j++) {
463*b3cdcd63SStefano Zampini           if (!private_labels[j]) { /* found a new cc  */
464ec1c413dSStefano Zampini             PetscBool same_set;
465ec1c413dSStefano Zampini 
466*b3cdcd63SStefano Zampini             graph->cptr[ncc] = cum_queue;
467*b3cdcd63SStefano Zampini             ncc++;
468*b3cdcd63SStefano Zampini             subset_counter++;
469*b3cdcd63SStefano Zampini             private_labels[j] = subset_counter;
470*b3cdcd63SStefano Zampini             graph->queue[cum_queue++] = graph->subset_idxs[i][j];
471*b3cdcd63SStefano Zampini             for (k=j+1;k<buffer_size;k++) { /* check for other nodes in new cc */
472674ae819SStefano Zampini               same_set = PETSC_TRUE;
473*b3cdcd63SStefano Zampini               for (s=0;s<sharingprocs;s++) {
474*b3cdcd63SStefano Zampini                 if (refine_buffer[j][s] != refine_buffer[k][s]) {
475674ae819SStefano Zampini                   same_set = PETSC_FALSE;
476674ae819SStefano Zampini                   break;
477674ae819SStefano Zampini                 }
478674ae819SStefano Zampini               }
479674ae819SStefano Zampini               if (same_set) {
480*b3cdcd63SStefano Zampini                 private_labels[k] = subset_counter;
481*b3cdcd63SStefano Zampini                 graph->queue[cum_queue++] = graph->subset_idxs[i][k];
482674ae819SStefano Zampini               }
483674ae819SStefano Zampini             }
484674ae819SStefano Zampini           }
485674ae819SStefano Zampini         }
486*b3cdcd63SStefano Zampini         graph->cptr[ncc] = cum_queue;
487*b3cdcd63SStefano Zampini         graph->subset_ncc[i] = subset_counter;
488*b3cdcd63SStefano Zampini         graph->queue_sorted = PETSC_FALSE;
489*b3cdcd63SStefano Zampini       } else { /* this subset does not need to be adapted */
490*b3cdcd63SStefano Zampini         ierr = PetscMemcpy(graph->queue+cum_queue,graph->subset_idxs[i],graph->subset_size[i]*sizeof(PetscInt));CHKERRQ(ierr);
491*b3cdcd63SStefano Zampini         ncc++;
492*b3cdcd63SStefano Zampini         cum_queue += graph->subset_size[i];
493*b3cdcd63SStefano Zampini         graph->cptr[ncc] = cum_queue;
494674ae819SStefano Zampini       }
495674ae819SStefano Zampini     }
496*b3cdcd63SStefano Zampini     graph->cptr[ncc] = cum_queue;
497*b3cdcd63SStefano Zampini     graph->ncc = ncc;
498*b3cdcd63SStefano Zampini     if (mss) {
499*b3cdcd63SStefano Zampini       ierr = PetscFree2(refine_buffer[0],private_labels);CHKERRQ(ierr);
500*b3cdcd63SStefano Zampini       ierr = PetscFree(refine_buffer);CHKERRQ(ierr);
501*b3cdcd63SStefano Zampini     }
502*b3cdcd63SStefano Zampini     ierr = PetscFree(labels);CHKERRQ(ierr);
503*b3cdcd63SStefano Zampini     ierr = MPI_Waitall(sum_requests,send_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
504*b3cdcd63SStefano Zampini     ierr = PetscFree2(send_requests,recv_requests);CHKERRQ(ierr);
50545b2a5aaSStefano Zampini     ierr = PetscFree2(send_buffer,recv_buffer);CHKERRQ(ierr);
506674ae819SStefano Zampini     ierr = PetscFree(cum_recv_counts);CHKERRQ(ierr);
507ec1c413dSStefano Zampini     ierr = PetscBTDestroy(&subset_cc_adapt);CHKERRQ(ierr);
508674ae819SStefano Zampini   }
5098e4af1ccSStefano Zampini 
5108e4af1ccSStefano Zampini   /* Determine if we are in 2D or 3D */
5118e4af1ccSStefano Zampini   twodim  = PETSC_TRUE;
5128e4af1ccSStefano Zampini   for (i=0;i<graph->ncc;i++) {
5138e4af1ccSStefano Zampini     PetscInt repdof = graph->queue[graph->cptr[i]];
5148e4af1ccSStefano Zampini     if (graph->cptr[i+1]-graph->cptr[i] > graph->custom_minimal_size) {
5158e4af1ccSStefano Zampini       if (graph->count[repdof] > 1 || graph->special_dof[repdof] == PCBDDCGRAPH_NEUMANN_MARK) {
5168e4af1ccSStefano Zampini         twodim = PETSC_FALSE;
5178e4af1ccSStefano Zampini         break;
5188e4af1ccSStefano Zampini       }
5198e4af1ccSStefano Zampini     }
5208e4af1ccSStefano Zampini   }
521b2566f29SBarry Smith   ierr = MPIU_Allreduce(&twodim,&graph->twodim,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)graph->l2gmap));CHKERRQ(ierr);
522674ae819SStefano Zampini   PetscFunctionReturn(0);
523674ae819SStefano Zampini }
524674ae819SStefano Zampini 
525*b3cdcd63SStefano Zampini 
526*b3cdcd63SStefano Zampini #undef __FUNCT__
527*b3cdcd63SStefano Zampini #define __FUNCT__ "PCBDDCGraphComputeCC_Private"
528*b3cdcd63SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCGraphComputeCC_Private(PCBDDCGraph graph,PetscInt pid,PetscInt* queue_tip,PetscInt n_prev,PetscInt* n_added)
529*b3cdcd63SStefano Zampini {
530*b3cdcd63SStefano Zampini   PetscInt       i,j,n;
531*b3cdcd63SStefano Zampini   PetscInt       *xadj = graph->xadj,*adjncy = graph->adjncy;
532*b3cdcd63SStefano Zampini   PetscBT        touched = graph->touched;
533*b3cdcd63SStefano Zampini   PetscBool      havecsr = (PetscBool)(xadj && adjncy);
534*b3cdcd63SStefano Zampini   PetscBool      havesubs = (PetscBool)(!!graph->n_local_subs);
535*b3cdcd63SStefano Zampini   PetscErrorCode ierr;
536*b3cdcd63SStefano Zampini 
537*b3cdcd63SStefano Zampini   PetscFunctionBegin;
538*b3cdcd63SStefano Zampini   n = 0;
539*b3cdcd63SStefano Zampini   if (havecsr && !havesubs) {
540*b3cdcd63SStefano Zampini     for (i=-n_prev;i<0;i++) {
541*b3cdcd63SStefano Zampini       PetscInt start_dof = queue_tip[i];
542*b3cdcd63SStefano Zampini       for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
543*b3cdcd63SStefano Zampini         PetscInt dof = adjncy[j];
544*b3cdcd63SStefano Zampini         if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid) {
545*b3cdcd63SStefano Zampini           ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
546*b3cdcd63SStefano Zampini           queue_tip[n] = dof;
547*b3cdcd63SStefano Zampini           n++;
548*b3cdcd63SStefano Zampini         }
549*b3cdcd63SStefano Zampini       }
550*b3cdcd63SStefano Zampini     }
551*b3cdcd63SStefano Zampini   } else if (havecsr && havesubs) {
552*b3cdcd63SStefano Zampini     PetscInt sid = graph->local_subs[queue_tip[-n_prev]];
553*b3cdcd63SStefano Zampini     for (i=-n_prev;i<0;i++) {
554*b3cdcd63SStefano Zampini       PetscInt start_dof = queue_tip[i];
555*b3cdcd63SStefano Zampini       for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
556*b3cdcd63SStefano Zampini         PetscInt dof = adjncy[j];
557*b3cdcd63SStefano Zampini         if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid && graph->local_subs[dof] == sid) {
558*b3cdcd63SStefano Zampini           ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
559*b3cdcd63SStefano Zampini           queue_tip[n] = dof;
560*b3cdcd63SStefano Zampini           n++;
561*b3cdcd63SStefano Zampini         }
562*b3cdcd63SStefano Zampini       }
563*b3cdcd63SStefano Zampini     }
564*b3cdcd63SStefano Zampini   } else { /* sub info only */
565*b3cdcd63SStefano Zampini     PetscInt sid = graph->local_subs[queue_tip[-n_prev]];
566*b3cdcd63SStefano Zampini     for (j=0;j<graph->subset_size[pid-1];j++) { /* pid \in [1,graph->n_subsets] */
567*b3cdcd63SStefano Zampini       PetscInt dof = graph->subset_idxs[pid-1][j];
568*b3cdcd63SStefano Zampini       if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid && graph->local_subs[dof] == sid) {
569*b3cdcd63SStefano Zampini         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
570*b3cdcd63SStefano Zampini         queue_tip[n] = dof;
571*b3cdcd63SStefano Zampini         n++;
572*b3cdcd63SStefano Zampini       }
573*b3cdcd63SStefano Zampini     }
574*b3cdcd63SStefano Zampini   }
575*b3cdcd63SStefano Zampini   *n_added = n;
576*b3cdcd63SStefano Zampini   PetscFunctionReturn(0);
577*b3cdcd63SStefano Zampini }
578*b3cdcd63SStefano Zampini 
579674ae819SStefano Zampini #undef __FUNCT__
580674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphComputeConnectedComponentsLocal"
581674ae819SStefano Zampini PetscErrorCode PCBDDCGraphComputeConnectedComponentsLocal(PCBDDCGraph graph)
582674ae819SStefano Zampini {
583*b3cdcd63SStefano Zampini   PetscInt       ncc,cum_queue,n;
584*b3cdcd63SStefano Zampini   PetscMPIInt    commsize;
585674ae819SStefano Zampini   PetscErrorCode ierr;
586674ae819SStefano Zampini 
587674ae819SStefano Zampini   PetscFunctionBegin;
588*b3cdcd63SStefano Zampini   if (!graph->setupcalled) SETERRQ(PetscObjectComm((PetscObject)graph->l2gmap),PETSC_ERR_ORDER,"PCBDDCGraphSetUp should be called first");
589*b3cdcd63SStefano Zampini   /* quiet return if there isn't any local info */
5904f1b2e48SStefano Zampini   if ((!graph->xadj || !graph->adjncy) && !graph->n_local_subs) {
591674ae819SStefano Zampini     PetscFunctionReturn(0);
592674ae819SStefano Zampini   }
593674ae819SStefano Zampini 
594674ae819SStefano Zampini   /* reset any previous search of connected components */
595df48933bSStefano Zampini   ierr = PetscBTMemzero(graph->nvtxs,graph->touched);CHKERRQ(ierr);
596*b3cdcd63SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)graph->l2gmap),&commsize);CHKERRQ(ierr);
597*b3cdcd63SStefano Zampini   if (commsize > 1) {
598*b3cdcd63SStefano Zampini     PetscInt i;
599674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
6000cf82fd6SStefano Zampini       if (graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK || !graph->count[i]) {
601df48933bSStefano Zampini         ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
602674ae819SStefano Zampini       }
6034a060362SStefano Zampini     }
604*b3cdcd63SStefano Zampini   }
605674ae819SStefano Zampini 
606674ae819SStefano Zampini   /* begin search for connected components */
607674ae819SStefano Zampini   cum_queue = 0;
608674ae819SStefano Zampini   ncc = 0;
609674ae819SStefano Zampini   for (n=0;n<graph->n_subsets;n++) {
610*b3cdcd63SStefano Zampini     PetscInt pid = n+1;  /* partition labeled by 0 is discarded */
611*b3cdcd63SStefano Zampini     PetscInt found = 0,prev = 0,first = 0,ncc_pid = 0;
612*b3cdcd63SStefano Zampini     while (found != graph->subset_size[n]) {
613*b3cdcd63SStefano Zampini       PetscInt added = 0;
614*b3cdcd63SStefano Zampini       if (!prev) { /* search for new starting dof */
615*b3cdcd63SStefano Zampini         while (PetscBTLookup(graph->touched,graph->subset_idxs[n][first])) first++;
616*b3cdcd63SStefano Zampini         ierr = PetscBTSet(graph->touched,graph->subset_idxs[n][first]);CHKERRQ(ierr);
617*b3cdcd63SStefano Zampini         graph->queue[cum_queue] = graph->subset_idxs[n][first];
618674ae819SStefano Zampini         graph->cptr[ncc] = cum_queue;
619*b3cdcd63SStefano Zampini         prev = 1;
620*b3cdcd63SStefano Zampini         cum_queue++;
621*b3cdcd63SStefano Zampini         found++;
622674ae819SStefano Zampini         ncc_pid++;
623*b3cdcd63SStefano Zampini         ncc++;
624674ae819SStefano Zampini       }
625*b3cdcd63SStefano Zampini       ierr = PCBDDCGraphComputeCC_Private(graph,pid,graph->queue + cum_queue,prev,&added);CHKERRQ(ierr);
626*b3cdcd63SStefano Zampini       if (!added) {
627674ae819SStefano Zampini         graph->subset_ncc[n] = ncc_pid;
628*b3cdcd63SStefano Zampini         graph->cptr[ncc] = cum_queue;
629*b3cdcd63SStefano Zampini       }
630*b3cdcd63SStefano Zampini       prev = added;
631*b3cdcd63SStefano Zampini       found += added;
632*b3cdcd63SStefano Zampini       cum_queue += added;
633*b3cdcd63SStefano Zampini       if (added && found == graph->subset_size[n]) {
634*b3cdcd63SStefano Zampini         graph->subset_ncc[n] = ncc_pid;
635*b3cdcd63SStefano Zampini         graph->cptr[ncc] = cum_queue;
636*b3cdcd63SStefano Zampini       }
637*b3cdcd63SStefano Zampini     }
638674ae819SStefano Zampini   }
639674ae819SStefano Zampini   graph->ncc = ncc;
640acc113dbSStefano Zampini   graph->queue_sorted = PETSC_FALSE;
641674ae819SStefano Zampini   PetscFunctionReturn(0);
642674ae819SStefano Zampini }
643674ae819SStefano Zampini 
644674ae819SStefano Zampini #undef __FUNCT__
645674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphSetUp"
646674ae819SStefano Zampini PetscErrorCode PCBDDCGraphSetUp(PCBDDCGraph graph, PetscInt custom_minimal_size, IS neumann_is, IS dirichlet_is, PetscInt n_ISForDofs, IS ISForDofs[], IS custom_primal_vertices)
647674ae819SStefano Zampini {
648674ae819SStefano Zampini   VecScatter     scatter_ctx;
649674ae819SStefano Zampini   Vec            local_vec,local_vec2,global_vec;
650dc456d91SStefano Zampini   IS             to,from,subset,subset_n;
651674ae819SStefano Zampini   MPI_Comm       comm;
652674ae819SStefano Zampini   PetscScalar    *array,*array2;
653674ae819SStefano Zampini   const PetscInt *is_indices;
654dc456d91SStefano Zampini   PetscInt       n_neigh,*neigh,*n_shared,**shared,*queue_global;
655674ae819SStefano Zampini   PetscInt       i,j,k,s,total_counts,nodes_touched,is_size;
656*b3cdcd63SStefano Zampini   PetscMPIInt    commsize;
65751b0f6cfSStefano Zampini   PetscBool      same_set,mirrors_found;
6587babac9bSStefano Zampini   PetscErrorCode ierr;
659674ae819SStefano Zampini 
660674ae819SStefano Zampini   PetscFunctionBegin;
661a07ea27aSStefano Zampini   graph->has_dirichlet = PETSC_FALSE;
662a07ea27aSStefano Zampini   if (dirichlet_is) {
663a07ea27aSStefano Zampini     PetscCheckSameComm(graph->l2gmap,1,dirichlet_is,4);
664a07ea27aSStefano Zampini     graph->has_dirichlet = PETSC_TRUE;
665a07ea27aSStefano Zampini   }
666674ae819SStefano Zampini   ierr = PetscObjectGetComm((PetscObject)(graph->l2gmap),&comm);CHKERRQ(ierr);
667*b3cdcd63SStefano Zampini   ierr = MPI_Comm_size(comm,&commsize);CHKERRQ(ierr);
668*b3cdcd63SStefano Zampini 
669674ae819SStefano Zampini   /* custom_minimal_size */
67014f95afaSStefano Zampini   graph->custom_minimal_size = custom_minimal_size;
671674ae819SStefano Zampini   /* get info l2gmap and allocate work vectors  */
672674ae819SStefano Zampini   ierr = ISLocalToGlobalMappingGetInfo(graph->l2gmap,&n_neigh,&neigh,&n_shared,&shared);CHKERRQ(ierr);
6732635a1d4SStefano Zampini   /* check if we have any local periodic nodes (periodic BCs) */
6742635a1d4SStefano Zampini   mirrors_found = PETSC_FALSE;
6752635a1d4SStefano Zampini   if (graph->nvtxs && n_neigh) {
6762635a1d4SStefano Zampini     for (i=0; i<n_shared[0]; i++) graph->count[shared[0][i]] += 1;
6772635a1d4SStefano Zampini     for (i=0; i<n_shared[0]; i++) {
6782635a1d4SStefano Zampini       if (graph->count[shared[0][i]] > 1) {
6792635a1d4SStefano Zampini         mirrors_found = PETSC_TRUE;
6802635a1d4SStefano Zampini         break;
6812635a1d4SStefano Zampini       }
6822635a1d4SStefano Zampini     }
6832635a1d4SStefano Zampini   }
6842635a1d4SStefano Zampini   /* create some workspace objects */
6852635a1d4SStefano Zampini   local_vec = NULL;
6862635a1d4SStefano Zampini   local_vec2 = NULL;
6872635a1d4SStefano Zampini   global_vec = NULL;
6882635a1d4SStefano Zampini   to = NULL;
6892635a1d4SStefano Zampini   from = NULL;
6902635a1d4SStefano Zampini   scatter_ctx = NULL;
6912635a1d4SStefano Zampini   if (n_ISForDofs || dirichlet_is || neumann_is || custom_primal_vertices) {
692674ae819SStefano Zampini     ierr = VecCreate(PETSC_COMM_SELF,&local_vec);CHKERRQ(ierr);
693674ae819SStefano Zampini     ierr = VecSetSizes(local_vec,PETSC_DECIDE,graph->nvtxs);CHKERRQ(ierr);
694674ae819SStefano Zampini     ierr = VecSetType(local_vec,VECSTANDARD);CHKERRQ(ierr);
695674ae819SStefano Zampini     ierr = VecDuplicate(local_vec,&local_vec2);CHKERRQ(ierr);
696674ae819SStefano Zampini     ierr = VecCreate(comm,&global_vec);CHKERRQ(ierr);
6977fb0e2dbSStefano Zampini     ierr = VecSetSizes(global_vec,PETSC_DECIDE,graph->nvtxs_global);CHKERRQ(ierr);
698674ae819SStefano Zampini     ierr = VecSetType(global_vec,VECSTANDARD);CHKERRQ(ierr);
699674ae819SStefano Zampini     ierr = ISCreateStride(PETSC_COMM_SELF,graph->nvtxs,0,1,&to);CHKERRQ(ierr);
700674ae819SStefano Zampini     ierr = ISLocalToGlobalMappingApplyIS(graph->l2gmap,to,&from);CHKERRQ(ierr);
701674ae819SStefano Zampini     ierr = VecScatterCreate(global_vec,from,local_vec,to,&scatter_ctx);CHKERRQ(ierr);
7022635a1d4SStefano Zampini   } else if (mirrors_found) {
7032635a1d4SStefano Zampini     ierr = ISCreateStride(PETSC_COMM_SELF,graph->nvtxs,0,1,&to);CHKERRQ(ierr);
7042635a1d4SStefano Zampini     ierr = ISLocalToGlobalMappingApplyIS(graph->l2gmap,to,&from);CHKERRQ(ierr);
70549eeff8cSStefano Zampini   }
70651b0f6cfSStefano Zampini   /* compute local mirrors (if any) */
70751b0f6cfSStefano Zampini   if (mirrors_found) {
70851b0f6cfSStefano Zampini     PetscInt *local_indices,*global_indices;
70951b0f6cfSStefano Zampini     /* get arrays of local and global indices */
710785e854fSJed Brown     ierr = PetscMalloc1(graph->nvtxs,&local_indices);CHKERRQ(ierr);
71151b0f6cfSStefano Zampini     ierr = ISGetIndices(to,(const PetscInt**)&is_indices);CHKERRQ(ierr);
71251b0f6cfSStefano Zampini     ierr = PetscMemcpy(local_indices,is_indices,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
71351b0f6cfSStefano Zampini     ierr = ISRestoreIndices(to,(const PetscInt**)&is_indices);CHKERRQ(ierr);
714785e854fSJed Brown     ierr = PetscMalloc1(graph->nvtxs,&global_indices);CHKERRQ(ierr);
71551b0f6cfSStefano Zampini     ierr = ISGetIndices(from,(const PetscInt**)&is_indices);CHKERRQ(ierr);
71651b0f6cfSStefano Zampini     ierr = PetscMemcpy(global_indices,is_indices,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
71751b0f6cfSStefano Zampini     ierr = ISRestoreIndices(from,(const PetscInt**)&is_indices);CHKERRQ(ierr);
71851b0f6cfSStefano Zampini     /* allocate space for mirrors */
7198e5aaad5SJed Brown     ierr = PetscMalloc2(graph->nvtxs,&graph->mirrors,graph->nvtxs,&graph->mirrors_set);CHKERRQ(ierr);
72051b0f6cfSStefano Zampini     ierr = PetscMemzero(graph->mirrors,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
72151b0f6cfSStefano Zampini     graph->mirrors_set[0] = 0;
72251b0f6cfSStefano Zampini 
72351b0f6cfSStefano Zampini     k=0;
72451b0f6cfSStefano Zampini     for (i=0;i<n_shared[0];i++) {
72551b0f6cfSStefano Zampini       j=shared[0][i];
72651b0f6cfSStefano Zampini       if (graph->count[j] > 1) {
72751b0f6cfSStefano Zampini         graph->mirrors[j]++;
72851b0f6cfSStefano Zampini         k++;
72951b0f6cfSStefano Zampini       }
73051b0f6cfSStefano Zampini     }
73151b0f6cfSStefano Zampini     /* allocate space for set of mirrors */
732785e854fSJed Brown     ierr = PetscMalloc1(k,&graph->mirrors_set[0]);CHKERRQ(ierr);
73351b0f6cfSStefano Zampini     for (i=1;i<graph->nvtxs;i++)
73451b0f6cfSStefano Zampini       graph->mirrors_set[i]=graph->mirrors_set[i-1]+graph->mirrors[i-1];
73551b0f6cfSStefano Zampini 
73651b0f6cfSStefano Zampini     /* fill arrays */
73751b0f6cfSStefano Zampini     ierr = PetscMemzero(graph->mirrors,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
73851b0f6cfSStefano Zampini     for (j=0;j<n_shared[0];j++) {
73951b0f6cfSStefano Zampini       i=shared[0][j];
74051b0f6cfSStefano Zampini       if (graph->count[i] > 1)
74151b0f6cfSStefano Zampini         graph->mirrors_set[i][graph->mirrors[i]++]=global_indices[i];
74251b0f6cfSStefano Zampini     }
74351b0f6cfSStefano Zampini     ierr = PetscSortIntWithArray(graph->nvtxs,global_indices,local_indices);CHKERRQ(ierr);
74451b0f6cfSStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
74551b0f6cfSStefano Zampini       if (graph->mirrors[i] > 0) {
74651b0f6cfSStefano Zampini         ierr = PetscFindInt(graph->mirrors_set[i][0],graph->nvtxs,global_indices,&k);CHKERRQ(ierr);
74751b0f6cfSStefano Zampini         j = global_indices[k];
74851b0f6cfSStefano Zampini         while ( k > 0 && global_indices[k-1] == j) k--;
74951b0f6cfSStefano Zampini         for (j=0;j<graph->mirrors[i];j++) {
75051b0f6cfSStefano Zampini           graph->mirrors_set[i][j]=local_indices[k+j];
75151b0f6cfSStefano Zampini         }
75251b0f6cfSStefano Zampini         ierr = PetscSortInt(graph->mirrors[i],graph->mirrors_set[i]);CHKERRQ(ierr);
75351b0f6cfSStefano Zampini       }
75451b0f6cfSStefano Zampini     }
75551b0f6cfSStefano Zampini     ierr = PetscFree(local_indices);CHKERRQ(ierr);
75651b0f6cfSStefano Zampini     ierr = PetscFree(global_indices);CHKERRQ(ierr);
75751b0f6cfSStefano Zampini   }
75851b0f6cfSStefano Zampini   ierr = PetscMemzero(graph->count,graph->nvtxs*sizeof(*graph->count));CHKERRQ(ierr);
759674ae819SStefano Zampini   ierr = ISDestroy(&to);CHKERRQ(ierr);
760674ae819SStefano Zampini   ierr = ISDestroy(&from);CHKERRQ(ierr);
76151b0f6cfSStefano Zampini 
762674ae819SStefano Zampini   /* Count total number of neigh per node */
763674ae819SStefano Zampini   k = 0;
764674ae819SStefano Zampini   for (i=1;i<n_neigh;i++) {
765674ae819SStefano Zampini     k += n_shared[i];
766674ae819SStefano Zampini     for (j=0;j<n_shared[i];j++) {
767674ae819SStefano Zampini       graph->count[shared[i][j]] += 1;
768674ae819SStefano Zampini     }
769674ae819SStefano Zampini   }
770674ae819SStefano Zampini   /* Allocate space for storing the set of neighbours for each node */
771674ae819SStefano Zampini   if (graph->nvtxs) {
772785e854fSJed Brown     ierr = PetscMalloc1(k,&graph->neighbours_set[0]);CHKERRQ(ierr);
773674ae819SStefano Zampini   }
774674ae819SStefano Zampini   for (i=1;i<graph->nvtxs;i++) { /* dont count myself */
775674ae819SStefano Zampini     graph->neighbours_set[i]=graph->neighbours_set[i-1]+graph->count[i-1];
776674ae819SStefano Zampini   }
777674ae819SStefano Zampini   /* Get information for sharing subdomains */
778674ae819SStefano Zampini   ierr = PetscMemzero(graph->count,graph->nvtxs*sizeof(*graph->count));CHKERRQ(ierr);
779674ae819SStefano Zampini   for (i=1;i<n_neigh;i++) { /* dont count myself */
780674ae819SStefano Zampini     s = n_shared[i];
781674ae819SStefano Zampini     for (j=0;j<s;j++) {
782674ae819SStefano Zampini       k = shared[i][j];
783674ae819SStefano Zampini       graph->neighbours_set[k][graph->count[k]] = neigh[i];
784674ae819SStefano Zampini       graph->count[k] += 1;
785674ae819SStefano Zampini     }
786674ae819SStefano Zampini   }
787674ae819SStefano Zampini   /* sort set of sharing subdomains */
788674ae819SStefano Zampini   for (i=0;i<graph->nvtxs;i++) {
78993fb5fd3SStefano Zampini     ierr = PetscSortRemoveDupsInt(&graph->count[i],graph->neighbours_set[i]);CHKERRQ(ierr);
790674ae819SStefano Zampini   }
7917fb0e2dbSStefano Zampini   /* free memory allocated by ISLocalToGlobalMappingGetInfo */
7927fb0e2dbSStefano Zampini   ierr = ISLocalToGlobalMappingRestoreInfo(graph->l2gmap,&n_neigh,&neigh,&n_shared,&shared);CHKERRQ(ierr);
7937fb0e2dbSStefano Zampini 
79467c9da69SStefano Zampini   /*
79567c9da69SStefano Zampini      Get info for dofs splitting
7965777c63cSStefano Zampini      User can specify just a subset; an additional field is considered as a complementary field
79767c9da69SStefano Zampini   */
798*b3cdcd63SStefano Zampini   for (i=0;i<graph->nvtxs;i++) graph->which_dof[i] = n_ISForDofs; /* by default a dof belongs to the complement set */
7995777c63cSStefano Zampini   if (n_ISForDofs) {
8005777c63cSStefano Zampini     ierr = VecSet(local_vec,-1.0);CHKERRQ(ierr);
8015777c63cSStefano Zampini   }
802674ae819SStefano Zampini   for (i=0;i<n_ISForDofs;i++) {
8035777c63cSStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
80463602bcaSStefano Zampini     ierr = ISGetLocalSize(ISForDofs[i],&is_size);CHKERRQ(ierr);
805674ae819SStefano Zampini     ierr = ISGetIndices(ISForDofs[i],(const PetscInt**)&is_indices);CHKERRQ(ierr);
806674ae819SStefano Zampini     for (j=0;j<is_size;j++) {
807d50ed60aSStefano Zampini       if (is_indices[j] > -1 && is_indices[j] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
808674ae819SStefano Zampini         graph->which_dof[is_indices[j]] = i;
8095777c63cSStefano Zampini         array[is_indices[j]] = 1.*i;
810674ae819SStefano Zampini       }
81167c9da69SStefano Zampini     }
812674ae819SStefano Zampini     ierr = ISRestoreIndices(ISForDofs[i],(const PetscInt**)&is_indices);CHKERRQ(ierr);
8135777c63cSStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
8145777c63cSStefano Zampini   }
8155777c63cSStefano Zampini   /* Check consistency among neighbours */
8165777c63cSStefano Zampini   if (n_ISForDofs) {
8175777c63cSStefano Zampini     ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8185777c63cSStefano Zampini     ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8195777c63cSStefano Zampini     ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8205777c63cSStefano Zampini     ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8215777c63cSStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
8225777c63cSStefano Zampini     ierr = VecGetArray(local_vec2,&array2);CHKERRQ(ierr);
8235777c63cSStefano Zampini     for (i=0;i<graph->nvtxs;i++){
8245777c63cSStefano Zampini       PetscInt field1,field2;
8255777c63cSStefano Zampini 
8265777c63cSStefano Zampini       field1 = (PetscInt)PetscRealPart(array[i]);
8275777c63cSStefano Zampini       field2 = (PetscInt)PetscRealPart(array2[i]);
8286c4ed002SBarry 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);
8295777c63cSStefano Zampini     }
8305777c63cSStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
8315777c63cSStefano Zampini     ierr = VecRestoreArray(local_vec2,&array2);CHKERRQ(ierr);
832674ae819SStefano Zampini   }
8337fb0e2dbSStefano Zampini   if (neumann_is || dirichlet_is) {
834674ae819SStefano Zampini     /* Take into account Neumann nodes */
835674ae819SStefano Zampini     ierr = VecSet(local_vec,0.0);CHKERRQ(ierr);
836674ae819SStefano Zampini     ierr = VecSet(local_vec2,0.0);CHKERRQ(ierr);
837674ae819SStefano Zampini     if (neumann_is) {
838674ae819SStefano Zampini       ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
83982ba6b80SStefano Zampini       ierr = ISGetLocalSize(neumann_is,&is_size);CHKERRQ(ierr);
840674ae819SStefano Zampini       ierr = ISGetIndices(neumann_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
841674ae819SStefano Zampini       for (i=0;i<is_size;i++) {
842d50ed60aSStefano Zampini         if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
843674ae819SStefano Zampini           array[is_indices[i]] = 1.0;
844674ae819SStefano Zampini         }
8453c629590SStefano Zampini       }
846674ae819SStefano Zampini       ierr = ISRestoreIndices(neumann_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
847674ae819SStefano Zampini       ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
848674ae819SStefano Zampini     }
849674ae819SStefano Zampini     /* Neumann nodes: impose consistency among neighbours */
850674ae819SStefano Zampini     ierr = VecSet(global_vec,0.0);CHKERRQ(ierr);
851674ae819SStefano Zampini     ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
852674ae819SStefano Zampini     ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
853674ae819SStefano Zampini     ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
854674ae819SStefano Zampini     ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
855674ae819SStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
856674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
8573c629590SStefano Zampini       if (PetscRealPart(array[i]) > 0.1) {
8580cf82fd6SStefano Zampini         graph->special_dof[i] = PCBDDCGRAPH_NEUMANN_MARK;
859674ae819SStefano Zampini       }
860674ae819SStefano Zampini     }
861674ae819SStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
862674ae819SStefano Zampini     /* Take into account Dirichlet nodes */
863674ae819SStefano Zampini     ierr = VecSet(local_vec2,0.0);CHKERRQ(ierr);
864674ae819SStefano Zampini     if (dirichlet_is) {
865674ae819SStefano Zampini       ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
866674ae819SStefano Zampini       ierr = VecGetArray(local_vec2,&array2);CHKERRQ(ierr);
86782ba6b80SStefano Zampini       ierr = ISGetLocalSize(dirichlet_is,&is_size);CHKERRQ(ierr);
868674ae819SStefano Zampini       ierr = ISGetIndices(dirichlet_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           k = is_indices[i];
872e477d8c3SStefano Zampini           if (!PetscBTLookup(graph->touched,k)) {
873af25d912SStefano Zampini             if (PetscRealPart(array[k]) > 0.1) 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);
874674ae819SStefano Zampini             array2[k] = 1.0;
875674ae819SStefano Zampini           }
876674ae819SStefano Zampini         }
8773c629590SStefano Zampini       }
878674ae819SStefano Zampini       ierr = ISRestoreIndices(dirichlet_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
879674ae819SStefano Zampini       ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
880674ae819SStefano Zampini       ierr = VecRestoreArray(local_vec2,&array2);CHKERRQ(ierr);
881674ae819SStefano Zampini     }
882674ae819SStefano Zampini     /* Dirichlet nodes: impose consistency among neighbours */
883674ae819SStefano Zampini     ierr = VecSet(global_vec,0.0);CHKERRQ(ierr);
884674ae819SStefano Zampini     ierr = VecScatterBegin(scatter_ctx,local_vec2,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
885674ae819SStefano Zampini     ierr = VecScatterEnd(scatter_ctx,local_vec2,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
886674ae819SStefano Zampini     ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
887674ae819SStefano Zampini     ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
888674ae819SStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
889674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
8903c629590SStefano Zampini       if (PetscRealPart(array[i]) > 0.1) {
891*b3cdcd63SStefano Zampini         if (commsize > 1) { /* dirichlet nodes treated as internal */
892df48933bSStefano Zampini           ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
893*b3cdcd63SStefano Zampini           graph->subset[i] = 0;
894*b3cdcd63SStefano Zampini         }
8950cf82fd6SStefano Zampini         graph->special_dof[i] = PCBDDCGRAPH_DIRICHLET_MARK;
896674ae819SStefano Zampini       }
897674ae819SStefano Zampini     }
898674ae819SStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
8997fb0e2dbSStefano Zampini   }
90008a5cf49SStefano Zampini   /* mark local periodic nodes (if any) and adapt CSR graph (if any) */
90151b0f6cfSStefano Zampini   if (graph->mirrors) {
90251b0f6cfSStefano Zampini     for (i=0;i<graph->nvtxs;i++)
90351b0f6cfSStefano Zampini       if (graph->mirrors[i])
9040cf82fd6SStefano Zampini         graph->special_dof[i] = PCBDDCGRAPH_LOCAL_PERIODIC_MARK;
90551b0f6cfSStefano Zampini 
90608a5cf49SStefano Zampini     if (graph->xadj && graph->adjncy) {
90708a5cf49SStefano Zampini       PetscInt *new_xadj,*new_adjncy;
90851b0f6cfSStefano Zampini       /* sort CSR graph */
90951b0f6cfSStefano Zampini       for (i=0;i<graph->nvtxs;i++)
91051b0f6cfSStefano Zampini         ierr = PetscSortInt(graph->xadj[i+1]-graph->xadj[i],&graph->adjncy[graph->xadj[i]]);CHKERRQ(ierr);
91151b0f6cfSStefano Zampini 
91251b0f6cfSStefano Zampini       /* adapt local CSR graph in case of local periodicity */
91351b0f6cfSStefano Zampini       k = 0;
91451b0f6cfSStefano Zampini       for (i=0;i<graph->nvtxs;i++)
91551b0f6cfSStefano Zampini         for (j=graph->xadj[i];j<graph->xadj[i+1];j++)
91651b0f6cfSStefano Zampini           k += graph->mirrors[graph->adjncy[j]];
91751b0f6cfSStefano Zampini 
918854ce69bSBarry Smith       ierr = PetscMalloc1(graph->nvtxs+1,&new_xadj);CHKERRQ(ierr);
919854ce69bSBarry Smith       ierr = PetscMalloc1(k+graph->xadj[graph->nvtxs],&new_adjncy);CHKERRQ(ierr);
92051b0f6cfSStefano Zampini       new_xadj[0] = 0;
92151b0f6cfSStefano Zampini       for (i=0;i<graph->nvtxs;i++) {
92251b0f6cfSStefano Zampini         k = graph->xadj[i+1]-graph->xadj[i];
92351b0f6cfSStefano Zampini         ierr = PetscMemcpy(&new_adjncy[new_xadj[i]],&graph->adjncy[graph->xadj[i]],k*sizeof(PetscInt));CHKERRQ(ierr);
92451b0f6cfSStefano Zampini         new_xadj[i+1] = new_xadj[i]+k;
92551b0f6cfSStefano Zampini         for (j=graph->xadj[i];j<graph->xadj[i+1];j++) {
92651b0f6cfSStefano Zampini           k = graph->mirrors[graph->adjncy[j]];
92751b0f6cfSStefano Zampini           ierr = PetscMemcpy(&new_adjncy[new_xadj[i+1]],graph->mirrors_set[graph->adjncy[j]],k*sizeof(PetscInt));CHKERRQ(ierr);
92851b0f6cfSStefano Zampini           new_xadj[i+1] += k;
92951b0f6cfSStefano Zampini         }
93051b0f6cfSStefano Zampini         k = new_xadj[i+1]-new_xadj[i];
93151b0f6cfSStefano Zampini         ierr = PetscSortRemoveDupsInt(&k,&new_adjncy[new_xadj[i]]);CHKERRQ(ierr);
93251b0f6cfSStefano Zampini         new_xadj[i+1] = new_xadj[i]+k;
93351b0f6cfSStefano Zampini       }
93451b0f6cfSStefano Zampini       /* set new CSR into graph */
93551b0f6cfSStefano Zampini       ierr = PetscFree(graph->xadj);CHKERRQ(ierr);
93651b0f6cfSStefano Zampini       ierr = PetscFree(graph->adjncy);CHKERRQ(ierr);
93751b0f6cfSStefano Zampini       graph->xadj = new_xadj;
93851b0f6cfSStefano Zampini       graph->adjncy = new_adjncy;
93951b0f6cfSStefano Zampini     }
94008a5cf49SStefano Zampini   }
94151b0f6cfSStefano Zampini 
94217eb1463SStefano Zampini   /* mark special nodes (if any) -> each will become a single node equivalence class */
943674ae819SStefano Zampini   if (custom_primal_vertices) {
94417eb1463SStefano Zampini     ierr = VecSet(local_vec,0.0);CHKERRQ(ierr);
9459b70a373SStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
94663602bcaSStefano Zampini     ierr = ISGetLocalSize(custom_primal_vertices,&is_size);CHKERRQ(ierr);
947674ae819SStefano Zampini     ierr = ISGetIndices(custom_primal_vertices,(const PetscInt**)&is_indices);CHKERRQ(ierr);
948674ae819SStefano Zampini     for (i=0;i<is_size;i++){
9499b70a373SStefano Zampini       if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
9509b70a373SStefano Zampini         array[is_indices[i]] = 1.0;
9519b70a373SStefano Zampini       }
952674ae819SStefano Zampini     }
953674ae819SStefano Zampini     ierr = ISRestoreIndices(custom_primal_vertices,(const PetscInt**)&is_indices);CHKERRQ(ierr);
9549b70a373SStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
9559b70a373SStefano Zampini     /* special nodes: impose consistency among neighbours */
9569b70a373SStefano Zampini     ierr = VecSet(global_vec,0.0);CHKERRQ(ierr);
9579b70a373SStefano Zampini     ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9589b70a373SStefano Zampini     ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9599b70a373SStefano Zampini     ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9609b70a373SStefano Zampini     ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9619b70a373SStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
9629b70a373SStefano Zampini     j = 0;
9639b70a373SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
9649b70a373SStefano Zampini       if (PetscRealPart(array[i]) > 0.1 && graph->special_dof[i] != PCBDDCGRAPH_DIRICHLET_MARK) {
9659b70a373SStefano Zampini         graph->special_dof[i] = PCBDDCGRAPH_SPECIAL_MARK-j;
9669b70a373SStefano Zampini         j++;
9679b70a373SStefano Zampini       }
9689b70a373SStefano Zampini     }
9699b70a373SStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
97017eb1463SStefano Zampini   }
9719b70a373SStefano Zampini 
972*b3cdcd63SStefano Zampini   /* mark interior nodes (if commsize > 1) as touched and belonging to partition number 0 */
973*b3cdcd63SStefano Zampini   if (commsize > 1) {
974674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
975674ae819SStefano Zampini       if (!graph->count[i]) {
976df48933bSStefano Zampini         ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
977674ae819SStefano Zampini         graph->subset[i] = 0;
978674ae819SStefano Zampini       }
979674ae819SStefano Zampini     }
980*b3cdcd63SStefano Zampini   }
981*b3cdcd63SStefano Zampini 
982674ae819SStefano Zampini   /* init graph structure and compute default subsets */
983674ae819SStefano Zampini   nodes_touched = 0;
984674ae819SStefano Zampini   for (i=0;i<graph->nvtxs;i++) {
985df48933bSStefano Zampini     if (PetscBTLookup(graph->touched,i)) {
986674ae819SStefano Zampini       nodes_touched++;
987674ae819SStefano Zampini     }
988674ae819SStefano Zampini   }
989674ae819SStefano Zampini   i = 0;
990674ae819SStefano Zampini   graph->ncc = 0;
991674ae819SStefano Zampini   total_counts = 0;
9927babac9bSStefano Zampini 
9937babac9bSStefano Zampini   /* allocated space for queues */
994*b3cdcd63SStefano Zampini   if (commsize == 1) {
9957babac9bSStefano Zampini     ierr = PetscMalloc2(graph->nvtxs+1,&graph->cptr,graph->nvtxs,&graph->queue);CHKERRQ(ierr);
9967babac9bSStefano Zampini   } else {
9977babac9bSStefano Zampini     PetscInt nused = graph->nvtxs - nodes_touched;
9987babac9bSStefano Zampini     ierr = PetscMalloc2(nused+1,&graph->cptr,nused,&graph->queue);CHKERRQ(ierr);
9997babac9bSStefano Zampini   }
10007babac9bSStefano Zampini 
1001674ae819SStefano Zampini   while (nodes_touched<graph->nvtxs) {
1002674ae819SStefano Zampini     /*  find first untouched node in local ordering */
1003*b3cdcd63SStefano Zampini     while (PetscBTLookup(graph->touched,i)) i++;
1004df48933bSStefano Zampini     ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
1005674ae819SStefano Zampini     graph->subset[i] = graph->ncc+1;
1006674ae819SStefano Zampini     graph->cptr[graph->ncc] = total_counts;
1007674ae819SStefano Zampini     graph->queue[total_counts] = i;
1008674ae819SStefano Zampini     total_counts++;
1009674ae819SStefano Zampini     nodes_touched++;
1010674ae819SStefano Zampini     /* now find all other nodes having the same set of sharing subdomains */
1011674ae819SStefano Zampini     for (j=i+1;j<graph->nvtxs;j++) {
101274e413f5SStefano Zampini       /* check for same number of sharing subdomains, dof number and same special mark */
1013df48933bSStefano 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]) {
1014674ae819SStefano Zampini         /* check for same set of sharing subdomains */
1015674ae819SStefano Zampini         same_set = PETSC_TRUE;
1016674ae819SStefano Zampini         for (k=0;k<graph->count[j];k++){
1017674ae819SStefano Zampini           if (graph->neighbours_set[i][k] != graph->neighbours_set[j][k]) {
1018674ae819SStefano Zampini             same_set = PETSC_FALSE;
1019674ae819SStefano Zampini           }
1020674ae819SStefano Zampini         }
1021674ae819SStefano Zampini         /* I found a friend of mine */
1022674ae819SStefano Zampini         if (same_set) {
1023df48933bSStefano Zampini           ierr = PetscBTSet(graph->touched,j);CHKERRQ(ierr);
1024*b3cdcd63SStefano Zampini           graph->subset[j] = graph->ncc+1;
1025674ae819SStefano Zampini           nodes_touched++;
1026674ae819SStefano Zampini           graph->queue[total_counts] = j;
1027674ae819SStefano Zampini           total_counts++;
1028674ae819SStefano Zampini         }
1029674ae819SStefano Zampini       }
1030674ae819SStefano Zampini     }
1031674ae819SStefano Zampini     graph->ncc++;
1032674ae819SStefano Zampini   }
1033*b3cdcd63SStefano 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 */
1034674ae819SStefano Zampini   graph->n_subsets = graph->ncc;
1035785e854fSJed Brown   ierr = PetscMalloc1(graph->n_subsets,&graph->subset_ncc);CHKERRQ(ierr);
1036674ae819SStefano Zampini   for (i=0;i<graph->n_subsets;i++) {
1037674ae819SStefano Zampini     graph->subset_ncc[i] = 1;
1038674ae819SStefano Zampini   }
1039674ae819SStefano Zampini   /* final pointer */
1040674ae819SStefano Zampini   graph->cptr[graph->ncc] = total_counts;
1041ec1c413dSStefano Zampini 
1042*b3cdcd63SStefano Zampini   /* For consistency reasons (among neighbours), I need to sort (by global ordering) each connected component */
1043ec1c413dSStefano Zampini   /* Get a reference node (min index in global ordering) for each subset for tagging messages */
1044dc456d91SStefano Zampini   ierr = PetscMalloc1(graph->ncc,&graph->subset_ref_node);CHKERRQ(ierr);
1045dc456d91SStefano Zampini   ierr = PetscMalloc1(graph->cptr[graph->ncc],&queue_global);CHKERRQ(ierr);
10463a5b03d1SStefano Zampini   ierr = ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_global);CHKERRQ(ierr);
1047ec1c413dSStefano Zampini   for (j=0;j<graph->ncc;j++) {
1048ec1c413dSStefano Zampini     ierr = PetscSortIntWithArray(graph->cptr[j+1]-graph->cptr[j],&queue_global[graph->cptr[j]],&graph->queue[graph->cptr[j]]);CHKERRQ(ierr);
1049dc456d91SStefano Zampini     graph->subset_ref_node[j] = graph->queue[graph->cptr[j]];
1050f0321c90SStefano Zampini   }
1051dc456d91SStefano Zampini   ierr = PetscFree(queue_global);CHKERRQ(ierr);
1052acc113dbSStefano Zampini   graph->queue_sorted = PETSC_TRUE;
1053*b3cdcd63SStefano Zampini   /* save information on subsets (needed when analyzing the connected components) */
10542f566a09SStefano Zampini   if (graph->ncc) {
1055*b3cdcd63SStefano Zampini     ierr = PetscMalloc2(graph->ncc,&graph->subset_size,graph->ncc,&graph->subset_idxs);CHKERRQ(ierr);
1056*b3cdcd63SStefano Zampini     ierr = PetscMalloc1(graph->cptr[graph->ncc],&graph->subset_idxs[0]);CHKERRQ(ierr);
1057*b3cdcd63SStefano Zampini     ierr = PetscMemzero(graph->subset_idxs[0],graph->cptr[graph->ncc]*sizeof(PetscInt));CHKERRQ(ierr);
1058ec1c413dSStefano Zampini     for (j=1;j<graph->ncc;j++) {
1059*b3cdcd63SStefano Zampini       graph->subset_size[j-1] = graph->cptr[j] - graph->cptr[j-1];
1060*b3cdcd63SStefano Zampini       graph->subset_idxs[j] = graph->subset_idxs[j-1] + graph->subset_size[j-1];
1061ec1c413dSStefano Zampini     }
1062*b3cdcd63SStefano Zampini     graph->subset_size[graph->ncc-1] = graph->cptr[graph->ncc] - graph->cptr[graph->ncc-1];
1063*b3cdcd63SStefano Zampini     ierr = PetscMemcpy(graph->subset_idxs[0],graph->queue,graph->cptr[graph->ncc]*sizeof(PetscInt));CHKERRQ(ierr);
1064ec1c413dSStefano Zampini   }
1065*b3cdcd63SStefano Zampini 
1066f0321c90SStefano Zampini   /* renumber reference nodes */
1067dc456d91SStefano Zampini   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)(graph->l2gmap)),graph->ncc,graph->subset_ref_node,PETSC_COPY_VALUES,&subset_n);CHKERRQ(ierr);
1068dc456d91SStefano Zampini   ierr = ISLocalToGlobalMappingApplyIS(graph->l2gmap,subset_n,&subset);CHKERRQ(ierr);
1069dc456d91SStefano Zampini   ierr = ISDestroy(&subset_n);CHKERRQ(ierr);
1070dc456d91SStefano Zampini   ierr = PCBDDCSubsetNumbering(subset,NULL,NULL,&subset_n);CHKERRQ(ierr);
1071dc456d91SStefano Zampini   ierr = ISDestroy(&subset);CHKERRQ(ierr);
1072dc456d91SStefano Zampini   ierr = ISGetLocalSize(subset_n,&k);CHKERRQ(ierr);
10736c4ed002SBarry Smith   if (k != graph->ncc) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid size of new subset! %D != %D",k,graph->ncc);
1074dc456d91SStefano Zampini   ierr = ISGetIndices(subset_n,&is_indices);CHKERRQ(ierr);
1075dc456d91SStefano Zampini   ierr = PetscMemcpy(graph->subset_ref_node,is_indices,graph->ncc*sizeof(PetscInt));CHKERRQ(ierr);
1076dc456d91SStefano Zampini   ierr = ISRestoreIndices(subset_n,&is_indices);CHKERRQ(ierr);
1077dc456d91SStefano Zampini   ierr = ISDestroy(&subset_n);CHKERRQ(ierr);
1078ec1c413dSStefano Zampini 
1079ec1c413dSStefano Zampini   /* free workspace */
1080674ae819SStefano Zampini   ierr = VecDestroy(&local_vec);CHKERRQ(ierr);
1081674ae819SStefano Zampini   ierr = VecDestroy(&local_vec2);CHKERRQ(ierr);
1082674ae819SStefano Zampini   ierr = VecDestroy(&global_vec);CHKERRQ(ierr);
1083674ae819SStefano Zampini   ierr = VecScatterDestroy(&scatter_ctx);CHKERRQ(ierr);
1084*b3cdcd63SStefano Zampini   graph->setupcalled = PETSC_TRUE;
1085674ae819SStefano Zampini   PetscFunctionReturn(0);
1086674ae819SStefano Zampini }
1087674ae819SStefano Zampini 
1088674ae819SStefano Zampini #undef __FUNCT__
1089674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphResetCSR"
1090674ae819SStefano Zampini PetscErrorCode PCBDDCGraphResetCSR(PCBDDCGraph graph)
1091674ae819SStefano Zampini {
1092674ae819SStefano Zampini   PetscErrorCode ierr;
1093674ae819SStefano Zampini 
1094674ae819SStefano Zampini   PetscFunctionBegin;
1095674ae819SStefano Zampini   ierr = PetscFree(graph->xadj);CHKERRQ(ierr);
1096674ae819SStefano Zampini   ierr = PetscFree(graph->adjncy);CHKERRQ(ierr);
1097575ad6abSStefano Zampini   graph->nvtxs_csr = 0;
1098674ae819SStefano Zampini   PetscFunctionReturn(0);
1099674ae819SStefano Zampini }
1100674ae819SStefano Zampini 
1101674ae819SStefano Zampini #undef __FUNCT__
1102674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphReset"
1103674ae819SStefano Zampini PetscErrorCode PCBDDCGraphReset(PCBDDCGraph graph)
1104674ae819SStefano Zampini {
1105674ae819SStefano Zampini   PetscErrorCode ierr;
1106674ae819SStefano Zampini 
1107674ae819SStefano Zampini   PetscFunctionBegin;
1108*b3cdcd63SStefano Zampini   graph->setupcalled = PETSC_FALSE;
1109674ae819SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&graph->l2gmap);CHKERRQ(ierr);
1110674ae819SStefano Zampini   ierr = PetscFree(graph->subset_ncc);CHKERRQ(ierr);
11113a5b03d1SStefano Zampini   ierr = PetscFree(graph->subset_ref_node);CHKERRQ(ierr);
1112674ae819SStefano Zampini   if (graph->nvtxs) {
1113674ae819SStefano Zampini     ierr = PetscFree(graph->neighbours_set[0]);CHKERRQ(ierr);
1114674ae819SStefano Zampini   }
1115df48933bSStefano Zampini   ierr = PetscBTDestroy(&graph->touched);CHKERRQ(ierr);
11167babac9bSStefano Zampini   ierr = PetscFree5(graph->count,
1117674ae819SStefano Zampini                     graph->neighbours_set,
1118674ae819SStefano Zampini                     graph->subset,
1119674ae819SStefano Zampini                     graph->which_dof,
1120df48933bSStefano Zampini                     graph->special_dof);CHKERRQ(ierr);
11217babac9bSStefano Zampini   ierr = PetscFree2(graph->cptr,graph->queue);CHKERRQ(ierr);
112251b0f6cfSStefano Zampini   if (graph->mirrors) {
112351b0f6cfSStefano Zampini     ierr = PetscFree(graph->mirrors_set[0]);CHKERRQ(ierr);
112451b0f6cfSStefano Zampini   }
112551b0f6cfSStefano Zampini   ierr = PetscFree2(graph->mirrors,graph->mirrors_set);CHKERRQ(ierr);
1126*b3cdcd63SStefano Zampini   if (graph->subset_idxs) {
1127*b3cdcd63SStefano Zampini     ierr = PetscFree(graph->subset_idxs[0]);CHKERRQ(ierr);
1128ec1c413dSStefano Zampini   }
1129*b3cdcd63SStefano Zampini   ierr = PetscFree2(graph->subset_size,graph->subset_idxs);CHKERRQ(ierr);
1130a07ea27aSStefano Zampini   ierr = ISDestroy(&graph->dirdofs);CHKERRQ(ierr);
1131d62866d3SStefano Zampini   ierr = ISDestroy(&graph->dirdofsB);CHKERRQ(ierr);
11324f1b2e48SStefano Zampini   if (graph->n_local_subs) {
11334f1b2e48SStefano Zampini     ierr = PetscFree(graph->local_subs);CHKERRQ(ierr);
11344f1b2e48SStefano Zampini   }
1135a07ea27aSStefano Zampini   graph->has_dirichlet = PETSC_FALSE;
1136674ae819SStefano Zampini   graph->nvtxs = 0;
11377babac9bSStefano Zampini   graph->nvtxs_global = 0;
1138674ae819SStefano Zampini   graph->n_subsets = 0;
1139674ae819SStefano Zampini   graph->custom_minimal_size = 1;
11404f1b2e48SStefano Zampini   graph->n_local_subs = 0;
1141674ae819SStefano Zampini   PetscFunctionReturn(0);
1142674ae819SStefano Zampini }
1143674ae819SStefano Zampini 
1144674ae819SStefano Zampini #undef __FUNCT__
1145674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphInit"
11467fb0e2dbSStefano Zampini PetscErrorCode PCBDDCGraphInit(PCBDDCGraph graph, ISLocalToGlobalMapping l2gmap, PetscInt N)
1147674ae819SStefano Zampini {
1148674ae819SStefano Zampini   PetscInt       n;
1149674ae819SStefano Zampini   PetscErrorCode ierr;
1150674ae819SStefano Zampini 
1151674ae819SStefano Zampini   PetscFunctionBegin;
1152674ae819SStefano Zampini   PetscValidPointer(graph,1);
1153674ae819SStefano Zampini   PetscValidHeaderSpecific(l2gmap,IS_LTOGM_CLASSID,2);
11547babac9bSStefano Zampini   PetscValidLogicalCollectiveInt(l2gmap,N,3);
11558e61c736SStefano Zampini   /* raise an error if already allocated */
11566c4ed002SBarry Smith   if (graph->nvtxs_global) SETERRQ(PetscObjectComm((PetscObject)l2gmap),PETSC_ERR_PLIB,"BDDCGraph already initialized");
1157674ae819SStefano Zampini   /* set number of vertices */
1158674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)l2gmap);CHKERRQ(ierr);
1159674ae819SStefano Zampini   graph->l2gmap = l2gmap;
1160674ae819SStefano Zampini   ierr = ISLocalToGlobalMappingGetSize(l2gmap,&n);CHKERRQ(ierr);
1161674ae819SStefano Zampini   graph->nvtxs = n;
11627fb0e2dbSStefano Zampini   graph->nvtxs_global = N;
1163674ae819SStefano Zampini   /* allocate used space */
1164df48933bSStefano Zampini   ierr = PetscBTCreate(graph->nvtxs,&graph->touched);CHKERRQ(ierr);
11657babac9bSStefano Zampini   ierr = PetscMalloc5(graph->nvtxs,&graph->count,
11668e5aaad5SJed Brown                       graph->nvtxs,&graph->neighbours_set,
11678e5aaad5SJed Brown                       graph->nvtxs,&graph->subset,
11688e5aaad5SJed Brown                       graph->nvtxs,&graph->which_dof,
11698e5aaad5SJed Brown                       graph->nvtxs,&graph->special_dof);CHKERRQ(ierr);
1170674ae819SStefano Zampini   /* zeroes memory */
1171674ae819SStefano Zampini   ierr = PetscMemzero(graph->count,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
1172674ae819SStefano Zampini   ierr = PetscMemzero(graph->subset,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
117363602bcaSStefano Zampini   /* use -1 as a default value for which_dof array */
117463602bcaSStefano Zampini   for (n=0;n<graph->nvtxs;n++) graph->which_dof[n] = -1;
117574e413f5SStefano Zampini   ierr = PetscMemzero(graph->special_dof,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
1176674ae819SStefano Zampini   /* zeroes first pointer to neighbour set */
1177674ae819SStefano Zampini   if (graph->nvtxs) {
1178674ae819SStefano Zampini     graph->neighbours_set[0] = 0;
1179674ae819SStefano Zampini   }
1180674ae819SStefano Zampini   /* zeroes workspace for values of ncc */
1181674ae819SStefano Zampini   graph->subset_ncc = 0;
11823a5b03d1SStefano Zampini   graph->subset_ref_node = 0;
1183674ae819SStefano Zampini   PetscFunctionReturn(0);
1184674ae819SStefano Zampini }
1185674ae819SStefano Zampini 
1186674ae819SStefano Zampini #undef __FUNCT__
1187674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphDestroy"
1188674ae819SStefano Zampini PetscErrorCode PCBDDCGraphDestroy(PCBDDCGraph* graph)
1189674ae819SStefano Zampini {
1190674ae819SStefano Zampini   PetscErrorCode ierr;
1191674ae819SStefano Zampini 
1192674ae819SStefano Zampini   PetscFunctionBegin;
1193674ae819SStefano Zampini   ierr = PCBDDCGraphReset(*graph);CHKERRQ(ierr);
1194674ae819SStefano Zampini   ierr = PetscFree(*graph);CHKERRQ(ierr);
1195674ae819SStefano Zampini   PetscFunctionReturn(0);
1196674ae819SStefano Zampini }
1197674ae819SStefano Zampini 
1198674ae819SStefano Zampini #undef __FUNCT__
1199674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphCreate"
1200674ae819SStefano Zampini PetscErrorCode PCBDDCGraphCreate(PCBDDCGraph *graph)
1201674ae819SStefano Zampini {
1202674ae819SStefano Zampini   PCBDDCGraph    new_graph;
1203674ae819SStefano Zampini   PetscErrorCode ierr;
1204674ae819SStefano Zampini 
1205674ae819SStefano Zampini   PetscFunctionBegin;
1206854ce69bSBarry Smith   ierr = PetscNew(&new_graph);CHKERRQ(ierr);
1207674ae819SStefano Zampini   new_graph->custom_minimal_size = 1;
1208674ae819SStefano Zampini   *graph = new_graph;
1209674ae819SStefano Zampini   PetscFunctionReturn(0);
1210674ae819SStefano Zampini }
1211