xref: /petsc/src/ksp/pc/impls/bddc/bddcgraph.c (revision d4890ccebc2b4fa3cd5558bb0f40abcf7b0bf636)
1af0996ceSBarry Smith #include <petsc/private/petscimpl.h>
2ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
3ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddcstructs.h>
4674ae819SStefano Zampini 
5674ae819SStefano Zampini #undef __FUNCT__
6d62866d3SStefano Zampini #define __FUNCT__ "PCBDDCGraphGetDirichletDofsB"
7d62866d3SStefano Zampini PetscErrorCode PCBDDCGraphGetDirichletDofsB(PCBDDCGraph graph, IS* dirdofs)
8d62866d3SStefano Zampini {
9d62866d3SStefano Zampini   PetscErrorCode ierr;
10d62866d3SStefano Zampini 
11d62866d3SStefano Zampini   PetscFunctionBegin;
12d62866d3SStefano Zampini   if (graph->dirdofsB) {
13d62866d3SStefano Zampini     ierr = PetscObjectReference((PetscObject)graph->dirdofsB);CHKERRQ(ierr);
14d62866d3SStefano Zampini   } else if (graph->has_dirichlet) {
15d62866d3SStefano Zampini     PetscInt i,size;
16d62866d3SStefano Zampini     PetscInt *dirdofs_idxs;
17d62866d3SStefano Zampini 
18d62866d3SStefano Zampini     size = 0;
19d62866d3SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
20d62866d3SStefano Zampini       if (graph->count[i] && graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) size++;
21d62866d3SStefano Zampini     }
22d62866d3SStefano Zampini 
23d62866d3SStefano Zampini     ierr = PetscMalloc1(size,&dirdofs_idxs);CHKERRQ(ierr);
24d62866d3SStefano Zampini     size = 0;
25d62866d3SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
26d62866d3SStefano Zampini       if (graph->count[i] && graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) dirdofs_idxs[size++] = i;
27d62866d3SStefano Zampini     }
28d62866d3SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,size,dirdofs_idxs,PETSC_OWN_POINTER,&graph->dirdofsB);CHKERRQ(ierr);
29d62866d3SStefano Zampini     ierr = PetscObjectReference((PetscObject)graph->dirdofsB);CHKERRQ(ierr);
30d62866d3SStefano Zampini   }
31d62866d3SStefano Zampini   *dirdofs = graph->dirdofsB;
32d62866d3SStefano Zampini   PetscFunctionReturn(0);
33d62866d3SStefano Zampini }
34d62866d3SStefano Zampini 
35d62866d3SStefano Zampini #undef __FUNCT__
361b968477SStefano Zampini #define __FUNCT__ "PCBDDCGraphGetDirichletDofs"
371b968477SStefano Zampini PetscErrorCode PCBDDCGraphGetDirichletDofs(PCBDDCGraph graph, IS* dirdofs)
381b968477SStefano Zampini {
391b968477SStefano Zampini   PetscErrorCode ierr;
401b968477SStefano Zampini 
411b968477SStefano Zampini   PetscFunctionBegin;
42a07ea27aSStefano Zampini   if (graph->dirdofs) {
43a07ea27aSStefano Zampini     ierr = PetscObjectReference((PetscObject)graph->dirdofs);CHKERRQ(ierr);
44a07ea27aSStefano Zampini   } else if (graph->has_dirichlet) {
45a07ea27aSStefano Zampini     PetscInt i,size;
46a07ea27aSStefano Zampini     PetscInt *dirdofs_idxs;
47a07ea27aSStefano Zampini 
481b968477SStefano Zampini     size = 0;
491b968477SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
501b968477SStefano Zampini       if (graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) size++;
511b968477SStefano Zampini     }
521b968477SStefano Zampini 
531b968477SStefano Zampini     ierr = PetscMalloc1(size,&dirdofs_idxs);CHKERRQ(ierr);
541b968477SStefano Zampini     size = 0;
551b968477SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
561b968477SStefano Zampini       if (graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) dirdofs_idxs[size++] = i;
571b968477SStefano Zampini     }
58a07ea27aSStefano Zampini     ierr = ISCreateGeneral(PetscObjectComm((PetscObject)graph->l2gmap),size,dirdofs_idxs,PETSC_OWN_POINTER,&graph->dirdofs);CHKERRQ(ierr);
59a07ea27aSStefano Zampini     ierr = PetscObjectReference((PetscObject)graph->dirdofs);CHKERRQ(ierr);
601b968477SStefano Zampini   }
61a07ea27aSStefano Zampini   *dirdofs = graph->dirdofs;
621b968477SStefano Zampini   PetscFunctionReturn(0);
631b968477SStefano Zampini }
641b968477SStefano Zampini 
651b968477SStefano Zampini #undef __FUNCT__
66674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphASCIIView"
67e49050b4SStefano Zampini PetscErrorCode PCBDDCGraphASCIIView(PCBDDCGraph graph, PetscInt verbosity_level, PetscViewer viewer)
68674ae819SStefano Zampini {
692b510759SStefano Zampini   PetscInt       i,j,tabs;
7093fb5fd3SStefano Zampini   PetscInt*      queue_in_global_numbering;
71674ae819SStefano Zampini   PetscErrorCode ierr;
72674ae819SStefano Zampini 
73674ae819SStefano Zampini   PetscFunctionBegin;
741575c14dSBarry Smith   ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
752b510759SStefano Zampini   ierr = PetscViewerASCIIGetTab(viewer,&tabs);CHKERRQ(ierr);
76674ae819SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"--------------------------------------------------\n");CHKERRQ(ierr);
77674ae819SStefano Zampini   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
78674ae819SStefano Zampini   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Local BDDC graph for subdomain %04d\n",PetscGlobalRank);CHKERRQ(ierr);
79674ae819SStefano Zampini   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Number of vertices %d\n",graph->nvtxs);CHKERRQ(ierr);
8014f95afaSStefano Zampini   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Custom minimal size %d\n",graph->custom_minimal_size);CHKERRQ(ierr);
81*d4890cceSStefano Zampini   if (verbosity_level > 2) {
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       }
110*d4890cceSStefano Zampini       if (verbosity_level > 3) {
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     }
144*d4890cceSStefano Zampini     if (verbosity_level > 1) {
145674ae819SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"):");CHKERRQ(ierr);
146*d4890cceSStefano Zampini       if (graph->twodim || graph->count[node_num] > 1 || (graph->count[node_num] == 1 && graph->special_dof[node_num] == PCBDDCGRAPH_NEUMANN_MARK)) {
1471ce3d396SStefano Zampini         printcc = PETSC_TRUE;
1481ce3d396SStefano Zampini       }
1491ce3d396SStefano Zampini       if (printcc) {
150674ae819SStefano Zampini         for (j=graph->cptr[i];j<graph->cptr[i+1];j++) {
15193fb5fd3SStefano Zampini           ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d (%d)",graph->queue[j],queue_in_global_numbering[j]);CHKERRQ(ierr);
152674ae819SStefano Zampini         }
1531ce3d396SStefano Zampini       }
154*d4890cceSStefano Zampini     } else {
155*d4890cceSStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,")");CHKERRQ(ierr);
156*d4890cceSStefano Zampini     }
157674ae819SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
1582b510759SStefano Zampini     ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr);
1592b510759SStefano Zampini     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
160674ae819SStefano Zampini   }
16193fb5fd3SStefano Zampini   ierr = PetscFree(queue_in_global_numbering);CHKERRQ(ierr);
162674ae819SStefano Zampini   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
163674ae819SStefano Zampini   PetscFunctionReturn(0);
164674ae819SStefano Zampini }
165674ae819SStefano Zampini 
166674ae819SStefano Zampini #undef __FUNCT__
167674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphGetCandidatesIS"
168a873d5d0SStefano Zampini PetscErrorCode PCBDDCGraphGetCandidatesIS(PCBDDCGraph graph, PetscInt *n_faces, IS *FacesIS[], PetscInt *n_edges, IS *EdgesIS[], IS *VerticesIS)
169674ae819SStefano Zampini {
170674ae819SStefano Zampini   IS             *ISForFaces,*ISForEdges,ISForVertices;
171adfc4fb2SStefano Zampini   PetscInt       i,nfc,nec,nvc,*idx;
172674ae819SStefano Zampini   PetscErrorCode ierr;
173674ae819SStefano Zampini 
174674ae819SStefano Zampini   PetscFunctionBegin;
175674ae819SStefano Zampini   /* loop on ccs to evalute number of faces, edges and vertices */
176674ae819SStefano Zampini   nfc = 0;
177674ae819SStefano Zampini   nec = 0;
178674ae819SStefano Zampini   nvc = 0;
179674ae819SStefano Zampini   for (i=0;i<graph->ncc;i++) {
1801e1f2d93SStefano Zampini     PetscInt repdof = graph->queue[graph->cptr[i]];
181674ae819SStefano Zampini     if (graph->cptr[i+1]-graph->cptr[i] > graph->custom_minimal_size) {
182e635b14bSStefano Zampini       if (graph->count[repdof] == 1 && graph->special_dof[repdof] != PCBDDCGRAPH_NEUMANN_MARK) {
183674ae819SStefano Zampini         nfc++;
184674ae819SStefano Zampini       } else { /* note that nec will be zero in 2d */
185674ae819SStefano Zampini         nec++;
186674ae819SStefano Zampini       }
187674ae819SStefano Zampini     } else {
188674ae819SStefano Zampini       nvc += graph->cptr[i+1]-graph->cptr[i];
189674ae819SStefano Zampini     }
190674ae819SStefano Zampini   }
191adfc4fb2SStefano Zampini 
1928e4af1ccSStefano Zampini   /* check if we are in 2D or 3D */
1938e4af1ccSStefano Zampini   if (graph->twodim) { /* we are in a 2D case -> edges are shared by 2 subregions and faces don't exist */
194674ae819SStefano Zampini     nec = nfc;
195674ae819SStefano Zampini     nfc = 0;
196674ae819SStefano Zampini   }
197adfc4fb2SStefano Zampini 
198674ae819SStefano Zampini   /* allocate IS arrays for faces, edges. Vertices need a single index set. */
199cf5a6209SStefano Zampini   if (FacesIS) {
200785e854fSJed Brown     ierr = PetscMalloc1(nfc,&ISForFaces);CHKERRQ(ierr);
201cf5a6209SStefano Zampini   }
202cf5a6209SStefano Zampini   if (EdgesIS) {
203785e854fSJed Brown     ierr = PetscMalloc1(nec,&ISForEdges);CHKERRQ(ierr);
204cf5a6209SStefano Zampini   }
205cf5a6209SStefano Zampini   if (VerticesIS) {
206785e854fSJed Brown     ierr = PetscMalloc1(nvc,&idx);CHKERRQ(ierr);
207cf5a6209SStefano Zampini   }
208cf5a6209SStefano Zampini 
209674ae819SStefano Zampini   /* loop on ccs to compute index sets for faces and edges */
210acc113dbSStefano Zampini   if (!graph->queue_sorted) {
211acc113dbSStefano Zampini     PetscInt *queue_global;
212acc113dbSStefano Zampini 
213acc113dbSStefano Zampini     ierr = PetscMalloc1(graph->cptr[graph->ncc],&queue_global);CHKERRQ(ierr);
214acc113dbSStefano Zampini     ierr = ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_global);CHKERRQ(ierr);
215acc113dbSStefano Zampini     for (i=0;i<graph->ncc;i++) {
216acc113dbSStefano Zampini       ierr = PetscSortIntWithArray(graph->cptr[i+1]-graph->cptr[i],&queue_global[graph->cptr[i]],&graph->queue[graph->cptr[i]]);CHKERRQ(ierr);
217acc113dbSStefano Zampini     }
218acc113dbSStefano Zampini     ierr = PetscFree(queue_global);CHKERRQ(ierr);
219acc113dbSStefano Zampini     graph->queue_sorted = PETSC_TRUE;
220acc113dbSStefano Zampini   }
221674ae819SStefano Zampini   nfc = 0;
222674ae819SStefano Zampini   nec = 0;
223674ae819SStefano Zampini   for (i=0;i<graph->ncc;i++) {
2241e1f2d93SStefano Zampini     PetscInt repdof = graph->queue[graph->cptr[i]];
225674ae819SStefano Zampini     if (graph->cptr[i+1]-graph->cptr[i] > graph->custom_minimal_size) {
226e635b14bSStefano Zampini       if (graph->count[repdof] == 1 && graph->special_dof[repdof] != PCBDDCGRAPH_NEUMANN_MARK) {
2278e4af1ccSStefano Zampini         if (graph->twodim) {
228cf5a6209SStefano Zampini           if (EdgesIS) {
229674ae819SStefano 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);
230cf5a6209SStefano Zampini           }
231674ae819SStefano Zampini           nec++;
232674ae819SStefano Zampini         } else {
233cf5a6209SStefano Zampini           if (FacesIS) {
234674ae819SStefano 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);
235cf5a6209SStefano Zampini           }
236674ae819SStefano Zampini           nfc++;
237674ae819SStefano Zampini         }
238674ae819SStefano Zampini       } else {
239cf5a6209SStefano Zampini         if (EdgesIS) {
240674ae819SStefano 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);
241cf5a6209SStefano Zampini         }
242674ae819SStefano Zampini         nec++;
243674ae819SStefano Zampini       }
244674ae819SStefano Zampini     }
245674ae819SStefano Zampini   }
246674ae819SStefano Zampini   /* index set for vertices */
247cf5a6209SStefano Zampini   if (VerticesIS) {
248b8ffe317SStefano Zampini     nvc = 0;
249674ae819SStefano Zampini     for (i=0;i<graph->ncc;i++) {
250674ae819SStefano Zampini       if (graph->cptr[i+1]-graph->cptr[i] <= graph->custom_minimal_size) {
251adfc4fb2SStefano Zampini         PetscInt j;
252adfc4fb2SStefano Zampini 
253674ae819SStefano Zampini         for (j=graph->cptr[i];j<graph->cptr[i+1];j++) {
254674ae819SStefano Zampini           idx[nvc]=graph->queue[j];
255674ae819SStefano Zampini           nvc++;
256674ae819SStefano Zampini         }
257674ae819SStefano Zampini       }
258674ae819SStefano Zampini     }
259674ae819SStefano Zampini     /* sort vertex set (by local ordering) */
260674ae819SStefano Zampini     ierr = PetscSortInt(nvc,idx);CHKERRQ(ierr);
261674ae819SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,nvc,idx,PETSC_OWN_POINTER,&ISForVertices);CHKERRQ(ierr);
262674ae819SStefano Zampini   }
263674ae819SStefano Zampini   /* get back info */
264a873d5d0SStefano Zampini   if (n_faces) *n_faces = nfc;
265d7796978SStefano Zampini   if (FacesIS) {
266674ae819SStefano Zampini     *FacesIS = ISForFaces;
267d7796978SStefano Zampini   }
268a873d5d0SStefano Zampini   if (n_edges) *n_edges = nec;
269d7796978SStefano Zampini   if (EdgesIS) {
270674ae819SStefano Zampini     *EdgesIS = ISForEdges;
271d7796978SStefano Zampini   }
272d7796978SStefano Zampini   if (VerticesIS) {
273674ae819SStefano Zampini     *VerticesIS = ISForVertices;
274d7796978SStefano Zampini   }
275674ae819SStefano Zampini   PetscFunctionReturn(0);
276674ae819SStefano Zampini }
277674ae819SStefano Zampini 
278674ae819SStefano Zampini #undef __FUNCT__
279674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphComputeConnectedComponents"
280674ae819SStefano Zampini PetscErrorCode PCBDDCGraphComputeConnectedComponents(PCBDDCGraph graph)
281674ae819SStefano Zampini {
2824a060362SStefano Zampini   PetscBool      adapt_interface_reduced;
283674ae819SStefano Zampini   MPI_Comm       interface_comm;
2844a060362SStefano Zampini   PetscMPIInt    size;
2858e4af1ccSStefano Zampini   PetscInt       i;
2868e4af1ccSStefano Zampini   PetscBool      twodim;
287674ae819SStefano Zampini   PetscErrorCode ierr;
288674ae819SStefano Zampini 
289674ae819SStefano Zampini   PetscFunctionBegin;
290674ae819SStefano Zampini   /* compute connected components locally */
291674ae819SStefano Zampini   ierr = PetscObjectGetComm((PetscObject)(graph->l2gmap),&interface_comm);CHKERRQ(ierr);
292674ae819SStefano Zampini   ierr = PCBDDCGraphComputeConnectedComponentsLocal(graph);CHKERRQ(ierr);
293674ae819SStefano Zampini   /* check consistency of connected components among neighbouring subdomains -> it adapt them in case it is needed */
2944a060362SStefano Zampini   ierr = MPI_Comm_size(interface_comm,&size);CHKERRQ(ierr);
2954a060362SStefano Zampini   adapt_interface_reduced = PETSC_FALSE;
2964a060362SStefano Zampini   if (size > 1) {
2974a060362SStefano Zampini     PetscInt i;
2984a060362SStefano Zampini     PetscBool adapt_interface = PETSC_FALSE;
299674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
300674ae819SStefano Zampini       /* We are not sure that on a given subset of the local interface,
301674ae819SStefano Zampini          with two connected components, the latters be the same among sharing subdomains */
302674ae819SStefano Zampini       if (graph->subset_ncc[i] > 1) {
3034a060362SStefano Zampini         adapt_interface = PETSC_TRUE;
304674ae819SStefano Zampini         break;
305674ae819SStefano Zampini       }
306674ae819SStefano Zampini     }
307b2566f29SBarry Smith     ierr = MPIU_Allreduce(&adapt_interface,&adapt_interface_reduced,1,MPIU_BOOL,MPI_LOR,interface_comm);CHKERRQ(ierr);
3084a060362SStefano Zampini   }
309674ae819SStefano Zampini 
310674ae819SStefano Zampini   if (graph->n_subsets && adapt_interface_reduced) {
311ec1c413dSStefano Zampini     PetscBT     subset_cc_adapt;
312ec1c413dSStefano Zampini     MPI_Request *send_requests,*recv_requests;
313ec1c413dSStefano Zampini     PetscInt    *send_buffer,*recv_buffer;
314ec1c413dSStefano Zampini     PetscInt    sum_requests,start_of_recv,start_of_send;
315ec1c413dSStefano Zampini     PetscInt    *cum_recv_counts;
316ec1c413dSStefano Zampini     PetscInt    *labels;
317b3cdcd63SStefano Zampini     PetscInt    ncc,cum_queue,mss,mns,j,k,s;
3188875e3e6SStefano Zampini     PetscInt    **refine_buffer=NULL,*private_labels = NULL;
3195b1b9aeaSStefano Zampini 
320ec1c413dSStefano Zampini     ierr = PetscMalloc1(graph->nvtxs,&labels);CHKERRQ(ierr);
321ec1c413dSStefano Zampini     ierr = PetscMemzero(labels,graph->nvtxs*sizeof(*labels));CHKERRQ(ierr);
322b3cdcd63SStefano Zampini     for (i=0;i<graph->ncc;i++)
323b3cdcd63SStefano Zampini       for (j=graph->cptr[i];j<graph->cptr[i+1];j++)
324b3cdcd63SStefano Zampini         labels[graph->queue[j]] = i;
325b3cdcd63SStefano Zampini 
326674ae819SStefano Zampini     /* allocate some space */
327854ce69bSBarry Smith     ierr = PetscMalloc1(graph->n_subsets+1,&cum_recv_counts);CHKERRQ(ierr);
328674ae819SStefano Zampini     ierr = PetscMemzero(cum_recv_counts,(graph->n_subsets+1)*sizeof(*cum_recv_counts));CHKERRQ(ierr);
329b3cdcd63SStefano Zampini 
330674ae819SStefano Zampini     /* first count how many neighbours per connected component I will receive from */
331674ae819SStefano Zampini     cum_recv_counts[0] = 0;
332b3cdcd63SStefano 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]];
333ec1c413dSStefano Zampini     ierr = PetscMalloc1(cum_recv_counts[graph->n_subsets],&recv_buffer);CHKERRQ(ierr);
334ec1c413dSStefano Zampini     ierr = PetscMalloc2(cum_recv_counts[graph->n_subsets],&send_requests,cum_recv_counts[graph->n_subsets],&recv_requests);CHKERRQ(ierr);
335674ae819SStefano Zampini     for (i=0;i<cum_recv_counts[graph->n_subsets];i++) {
336674ae819SStefano Zampini       send_requests[i] = MPI_REQUEST_NULL;
337674ae819SStefano Zampini       recv_requests[i] = MPI_REQUEST_NULL;
338674ae819SStefano Zampini     }
339b3cdcd63SStefano Zampini 
340ec1c413dSStefano Zampini     /* exchange with my neighbours the number of my connected components on the subset of interface */
341674ae819SStefano Zampini     sum_requests = 0;
342674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
343ec1c413dSStefano Zampini       PetscMPIInt neigh,tag;
344b3cdcd63SStefano Zampini       PetscInt    count,*neighs;
345ec1c413dSStefano Zampini 
346b3cdcd63SStefano Zampini       count = graph->count[graph->subset_idxs[i][0]];
347b3cdcd63SStefano Zampini       neighs = graph->neighbours_set[graph->subset_idxs[i][0]];
348ec1c413dSStefano Zampini       ierr = PetscMPIIntCast(2*graph->subset_ref_node[i],&tag);CHKERRQ(ierr);
349b3cdcd63SStefano Zampini       for (k=0;k<count;k++) {
350b3cdcd63SStefano Zampini         ierr = PetscMPIIntCast(neighs[k],&neigh);CHKERRQ(ierr);
351674ae819SStefano Zampini         ierr = MPI_Isend(&graph->subset_ncc[i],1,MPIU_INT,neigh,tag,interface_comm,&send_requests[sum_requests]);CHKERRQ(ierr);
352ec1c413dSStefano Zampini         ierr = MPI_Irecv(&recv_buffer[sum_requests],1,MPIU_INT,neigh,tag,interface_comm,&recv_requests[sum_requests]);CHKERRQ(ierr);
353674ae819SStefano Zampini         sum_requests++;
354674ae819SStefano Zampini       }
355674ae819SStefano Zampini     }
356674ae819SStefano Zampini     ierr = MPI_Waitall(sum_requests,recv_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
357674ae819SStefano Zampini     ierr = MPI_Waitall(sum_requests,send_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
358b3cdcd63SStefano Zampini 
359b3cdcd63SStefano Zampini     /* determine the subsets I have to adapt (those having more than 1 cc) */
360ec1c413dSStefano Zampini     ierr = PetscBTCreate(graph->n_subsets,&subset_cc_adapt);CHKERRQ(ierr);
361ec1c413dSStefano Zampini     ierr = PetscBTMemzero(graph->n_subsets,subset_cc_adapt);CHKERRQ(ierr);
362674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
363b3cdcd63SStefano Zampini       if (graph->subset_ncc[i] > 1) {
364b3cdcd63SStefano Zampini         ierr = PetscBTSet(subset_cc_adapt,i);CHKERRQ(ierr);
365b3cdcd63SStefano Zampini         continue;
366b3cdcd63SStefano Zampini       }
367674ae819SStefano Zampini       for (j=cum_recv_counts[i];j<cum_recv_counts[i+1];j++){
368b3cdcd63SStefano Zampini          if (recv_buffer[j] > 1) {
369ec1c413dSStefano Zampini           ierr = PetscBTSet(subset_cc_adapt,i);
370674ae819SStefano Zampini           break;
371674ae819SStefano Zampini         }
372674ae819SStefano Zampini       }
373674ae819SStefano Zampini     }
374ec1c413dSStefano Zampini     ierr = PetscFree(recv_buffer);CHKERRQ(ierr);
375b3cdcd63SStefano Zampini 
376ec1c413dSStefano Zampini     /* determine send/recv buffers sizes */
377ec1c413dSStefano Zampini     j = 0;
378b3cdcd63SStefano Zampini     mss = 0;
379674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
380ec1c413dSStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i)) {
381b3cdcd63SStefano Zampini         j += graph->subset_size[i];
382b3cdcd63SStefano Zampini         mss = PetscMax(graph->subset_size[i],mss);
383674ae819SStefano Zampini       }
384674ae819SStefano Zampini     }
385ec1c413dSStefano Zampini     k = 0;
386b3cdcd63SStefano Zampini     mns = 0;
387674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
388ec1c413dSStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i)) {
389b3cdcd63SStefano Zampini         k += (cum_recv_counts[i+1]-cum_recv_counts[i])*graph->subset_size[i];
390b3cdcd63SStefano Zampini         mns = PetscMax(cum_recv_counts[i+1]-cum_recv_counts[i],mns);
391674ae819SStefano Zampini       }
392674ae819SStefano Zampini     }
393ec1c413dSStefano Zampini     ierr = PetscMalloc2(j,&send_buffer,k,&recv_buffer);CHKERRQ(ierr);
394ec1c413dSStefano Zampini 
395b3cdcd63SStefano Zampini     /* fill send buffer (order matters: subset_idxs ordered by global ordering) */
396ec1c413dSStefano Zampini     j = 0;
397b3cdcd63SStefano Zampini     for (i=0;i<graph->n_subsets;i++)
398b3cdcd63SStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i))
399b3cdcd63SStefano Zampini         for (k=0;k<graph->subset_size[i];k++)
400b3cdcd63SStefano Zampini           send_buffer[j++] = labels[graph->subset_idxs[i][k]];
401ec1c413dSStefano Zampini 
402674ae819SStefano Zampini     /* now exchange the data */
403674ae819SStefano Zampini     start_of_recv = 0;
404674ae819SStefano Zampini     start_of_send = 0;
405674ae819SStefano Zampini     sum_requests = 0;
406674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
407ec1c413dSStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i)) {
408ec1c413dSStefano Zampini         PetscMPIInt neigh,tag;
409b3cdcd63SStefano Zampini         PetscInt    size_of_send = graph->subset_size[i];
410ec1c413dSStefano Zampini 
411b3cdcd63SStefano Zampini         j = graph->subset_idxs[i][0];
412ec1c413dSStefano Zampini         ierr = PetscMPIIntCast(2*graph->subset_ref_node[i]+1,&tag);CHKERRQ(ierr);
413674ae819SStefano Zampini         for (k=0;k<graph->count[j];k++) {
414674ae819SStefano Zampini           ierr = PetscMPIIntCast(graph->neighbours_set[j][k],&neigh);CHKERRQ(ierr);
415ec1c413dSStefano Zampini           ierr = MPI_Isend(&send_buffer[start_of_send],size_of_send,MPIU_INT,neigh,tag,interface_comm,&send_requests[sum_requests]);CHKERRQ(ierr);
416ec1c413dSStefano Zampini           ierr = MPI_Irecv(&recv_buffer[start_of_recv],size_of_send,MPIU_INT,neigh,tag,interface_comm,&recv_requests[sum_requests]);CHKERRQ(ierr);
417ec1c413dSStefano Zampini           start_of_recv += size_of_send;
418674ae819SStefano Zampini           sum_requests++;
419674ae819SStefano Zampini         }
420674ae819SStefano Zampini         start_of_send += size_of_send;
421674ae819SStefano Zampini       }
422674ae819SStefano Zampini     }
423674ae819SStefano Zampini     ierr = MPI_Waitall(sum_requests,recv_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
424d52c4852SStefano Zampini 
425b3cdcd63SStefano Zampini     /* refine connected components */
426674ae819SStefano Zampini     start_of_recv = 0;
427b3cdcd63SStefano Zampini     /* allocate some temporary space */
428b3cdcd63SStefano Zampini     if (mss) {
429b3cdcd63SStefano Zampini       ierr = PetscMalloc1(mss,&refine_buffer);CHKERRQ(ierr);
430b3cdcd63SStefano Zampini       ierr = PetscMalloc2(mss*(mns+1),&refine_buffer[0],mss,&private_labels);CHKERRQ(ierr);
431b3cdcd63SStefano Zampini     }
432b3cdcd63SStefano Zampini     ncc = 0;
433b3cdcd63SStefano Zampini     cum_queue = 0;
434b3cdcd63SStefano Zampini     graph->cptr[0] = 0;
435674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
436ec1c413dSStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i)) {
437b3cdcd63SStefano Zampini         PetscInt subset_counter = 0;
438b3cdcd63SStefano Zampini         PetscInt sharingprocs = cum_recv_counts[i+1]-cum_recv_counts[i]+1; /* count myself */
439b3cdcd63SStefano Zampini         PetscInt buffer_size = graph->subset_size[i];
440ec1c413dSStefano Zampini 
441b3cdcd63SStefano Zampini         /* compute pointers */
442b3cdcd63SStefano Zampini         for (j=1;j<buffer_size;j++) refine_buffer[j] = refine_buffer[j-1] + sharingprocs;
443b3cdcd63SStefano Zampini         /* analyze contributions from subdomains that share the i-th subset
444b3cdcd63SStefano Zampini            The stricture of refine_buffer is suitable to find intersections of ccs among sharingprocs.
445b3cdcd63SStefano 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)
446b3cdcd63SStefano Zampini            sharing procs connected components:
447ec1c413dSStefano Zampini              neigh 0: [0 1 4], [2 3], labels [4,7]  (2 connected components)
448ec1c413dSStefano Zampini              neigh 1: [0 1], [2 3 4], labels [3 2]  (2 connected components)
449ec1c413dSStefano Zampini              neigh 2: [0 4], [1], [2 3], labels [1 5 6] (3 connected components)
450b3cdcd63SStefano Zampini            refine_buffer will be filled as:
451ec1c413dSStefano Zampini              [ 4, 3, 1;
452ec1c413dSStefano Zampini                4, 2, 1;
453ec1c413dSStefano Zampini                7, 2, 6;
454ec1c413dSStefano Zampini                4, 3, 5;
455ec1c413dSStefano Zampini                7, 2, 6; ];
456b3cdcd63SStefano Zampini            The connected components in local ordering are [0], [1], [2 3], [4] */
457b3cdcd63SStefano Zampini         /* fill temp_buffer */
458b3cdcd63SStefano Zampini         for (k=0;k<buffer_size;k++) refine_buffer[k][0] = labels[graph->subset_idxs[i][k]];
459b3cdcd63SStefano Zampini         for (j=0;j<sharingprocs-1;j++) {
460b3cdcd63SStefano Zampini           for (k=0;k<buffer_size;k++) refine_buffer[k][j+1] = recv_buffer[start_of_recv+k];
461b3cdcd63SStefano Zampini           start_of_recv += buffer_size;
462674ae819SStefano Zampini         }
463b3cdcd63SStefano Zampini         ierr = PetscMemzero(private_labels,buffer_size*sizeof(PetscInt));CHKERRQ(ierr);
464b3cdcd63SStefano Zampini         for (j=0;j<buffer_size;j++) {
465b3cdcd63SStefano Zampini           if (!private_labels[j]) { /* found a new cc  */
466ec1c413dSStefano Zampini             PetscBool same_set;
467ec1c413dSStefano Zampini 
468b3cdcd63SStefano Zampini             graph->cptr[ncc] = cum_queue;
469b3cdcd63SStefano Zampini             ncc++;
470b3cdcd63SStefano Zampini             subset_counter++;
471b3cdcd63SStefano Zampini             private_labels[j] = subset_counter;
472b3cdcd63SStefano Zampini             graph->queue[cum_queue++] = graph->subset_idxs[i][j];
473b3cdcd63SStefano Zampini             for (k=j+1;k<buffer_size;k++) { /* check for other nodes in new cc */
474674ae819SStefano Zampini               same_set = PETSC_TRUE;
475b3cdcd63SStefano Zampini               for (s=0;s<sharingprocs;s++) {
476b3cdcd63SStefano Zampini                 if (refine_buffer[j][s] != refine_buffer[k][s]) {
477674ae819SStefano Zampini                   same_set = PETSC_FALSE;
478674ae819SStefano Zampini                   break;
479674ae819SStefano Zampini                 }
480674ae819SStefano Zampini               }
481674ae819SStefano Zampini               if (same_set) {
482b3cdcd63SStefano Zampini                 private_labels[k] = subset_counter;
483b3cdcd63SStefano Zampini                 graph->queue[cum_queue++] = graph->subset_idxs[i][k];
484674ae819SStefano Zampini               }
485674ae819SStefano Zampini             }
486674ae819SStefano Zampini           }
487674ae819SStefano Zampini         }
488b3cdcd63SStefano Zampini         graph->cptr[ncc] = cum_queue;
489b3cdcd63SStefano Zampini         graph->subset_ncc[i] = subset_counter;
490b3cdcd63SStefano Zampini         graph->queue_sorted = PETSC_FALSE;
491b3cdcd63SStefano Zampini       } else { /* this subset does not need to be adapted */
492b3cdcd63SStefano Zampini         ierr = PetscMemcpy(graph->queue+cum_queue,graph->subset_idxs[i],graph->subset_size[i]*sizeof(PetscInt));CHKERRQ(ierr);
493b3cdcd63SStefano Zampini         ncc++;
494b3cdcd63SStefano Zampini         cum_queue += graph->subset_size[i];
495b3cdcd63SStefano Zampini         graph->cptr[ncc] = cum_queue;
496674ae819SStefano Zampini       }
497674ae819SStefano Zampini     }
498b3cdcd63SStefano Zampini     graph->cptr[ncc] = cum_queue;
499b3cdcd63SStefano Zampini     graph->ncc = ncc;
500b3cdcd63SStefano Zampini     if (mss) {
501b3cdcd63SStefano Zampini       ierr = PetscFree2(refine_buffer[0],private_labels);CHKERRQ(ierr);
502b3cdcd63SStefano Zampini       ierr = PetscFree(refine_buffer);CHKERRQ(ierr);
503b3cdcd63SStefano Zampini     }
504b3cdcd63SStefano Zampini     ierr = PetscFree(labels);CHKERRQ(ierr);
505b3cdcd63SStefano Zampini     ierr = MPI_Waitall(sum_requests,send_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
506b3cdcd63SStefano Zampini     ierr = PetscFree2(send_requests,recv_requests);CHKERRQ(ierr);
50745b2a5aaSStefano Zampini     ierr = PetscFree2(send_buffer,recv_buffer);CHKERRQ(ierr);
508674ae819SStefano Zampini     ierr = PetscFree(cum_recv_counts);CHKERRQ(ierr);
509ec1c413dSStefano Zampini     ierr = PetscBTDestroy(&subset_cc_adapt);CHKERRQ(ierr);
510674ae819SStefano Zampini   }
5118e4af1ccSStefano Zampini 
5128e4af1ccSStefano Zampini   /* Determine if we are in 2D or 3D */
5138e4af1ccSStefano Zampini   twodim  = PETSC_TRUE;
5148e4af1ccSStefano Zampini   for (i=0;i<graph->ncc;i++) {
5158e4af1ccSStefano Zampini     PetscInt repdof = graph->queue[graph->cptr[i]];
5168e4af1ccSStefano Zampini     if (graph->cptr[i+1]-graph->cptr[i] > graph->custom_minimal_size) {
5178e4af1ccSStefano Zampini       if (graph->count[repdof] > 1 || graph->special_dof[repdof] == PCBDDCGRAPH_NEUMANN_MARK) {
5188e4af1ccSStefano Zampini         twodim = PETSC_FALSE;
5198e4af1ccSStefano Zampini         break;
5208e4af1ccSStefano Zampini       }
5218e4af1ccSStefano Zampini     }
5228e4af1ccSStefano Zampini   }
523b2566f29SBarry Smith   ierr = MPIU_Allreduce(&twodim,&graph->twodim,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)graph->l2gmap));CHKERRQ(ierr);
524674ae819SStefano Zampini   PetscFunctionReturn(0);
525674ae819SStefano Zampini }
526674ae819SStefano Zampini 
527b3cdcd63SStefano Zampini 
528b3cdcd63SStefano Zampini #undef __FUNCT__
529b3cdcd63SStefano Zampini #define __FUNCT__ "PCBDDCGraphComputeCC_Private"
530b3cdcd63SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCGraphComputeCC_Private(PCBDDCGraph graph,PetscInt pid,PetscInt* queue_tip,PetscInt n_prev,PetscInt* n_added)
531b3cdcd63SStefano Zampini {
532b3cdcd63SStefano Zampini   PetscInt       i,j,n;
533b3cdcd63SStefano Zampini   PetscInt       *xadj = graph->xadj,*adjncy = graph->adjncy;
534b3cdcd63SStefano Zampini   PetscBT        touched = graph->touched;
535b3cdcd63SStefano Zampini   PetscBool      havecsr = (PetscBool)(xadj && adjncy);
536b3cdcd63SStefano Zampini   PetscBool      havesubs = (PetscBool)(!!graph->n_local_subs);
537b3cdcd63SStefano Zampini   PetscErrorCode ierr;
538b3cdcd63SStefano Zampini 
539b3cdcd63SStefano Zampini   PetscFunctionBegin;
540b3cdcd63SStefano Zampini   n = 0;
541b3cdcd63SStefano Zampini   if (havecsr && !havesubs) {
542b3cdcd63SStefano Zampini     for (i=-n_prev;i<0;i++) {
543b3cdcd63SStefano Zampini       PetscInt start_dof = queue_tip[i];
544b3cdcd63SStefano Zampini       for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
545b3cdcd63SStefano Zampini         PetscInt dof = adjncy[j];
546b3cdcd63SStefano Zampini         if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid) {
547b3cdcd63SStefano Zampini           ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
548b3cdcd63SStefano Zampini           queue_tip[n] = dof;
549b3cdcd63SStefano Zampini           n++;
550b3cdcd63SStefano Zampini         }
551b3cdcd63SStefano Zampini       }
552b3cdcd63SStefano Zampini     }
553b3cdcd63SStefano Zampini   } else if (havecsr && havesubs) {
554b3cdcd63SStefano Zampini     PetscInt sid = graph->local_subs[queue_tip[-n_prev]];
555b3cdcd63SStefano Zampini     for (i=-n_prev;i<0;i++) {
556b3cdcd63SStefano Zampini       PetscInt start_dof = queue_tip[i];
557b3cdcd63SStefano Zampini       for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
558b3cdcd63SStefano Zampini         PetscInt dof = adjncy[j];
559b3cdcd63SStefano Zampini         if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid && graph->local_subs[dof] == sid) {
560b3cdcd63SStefano Zampini           ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
561b3cdcd63SStefano Zampini           queue_tip[n] = dof;
562b3cdcd63SStefano Zampini           n++;
563b3cdcd63SStefano Zampini         }
564b3cdcd63SStefano Zampini       }
565b3cdcd63SStefano Zampini     }
566b3cdcd63SStefano Zampini   } else { /* sub info only */
567b3cdcd63SStefano Zampini     PetscInt sid = graph->local_subs[queue_tip[-n_prev]];
568b3cdcd63SStefano Zampini     for (j=0;j<graph->subset_size[pid-1];j++) { /* pid \in [1,graph->n_subsets] */
569b3cdcd63SStefano Zampini       PetscInt dof = graph->subset_idxs[pid-1][j];
570b3cdcd63SStefano Zampini       if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid && graph->local_subs[dof] == sid) {
571b3cdcd63SStefano Zampini         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
572b3cdcd63SStefano Zampini         queue_tip[n] = dof;
573b3cdcd63SStefano Zampini         n++;
574b3cdcd63SStefano Zampini       }
575b3cdcd63SStefano Zampini     }
576b3cdcd63SStefano Zampini   }
577b3cdcd63SStefano Zampini   *n_added = n;
578b3cdcd63SStefano Zampini   PetscFunctionReturn(0);
579b3cdcd63SStefano Zampini }
580b3cdcd63SStefano Zampini 
581674ae819SStefano Zampini #undef __FUNCT__
582674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphComputeConnectedComponentsLocal"
583674ae819SStefano Zampini PetscErrorCode PCBDDCGraphComputeConnectedComponentsLocal(PCBDDCGraph graph)
584674ae819SStefano Zampini {
585b3cdcd63SStefano Zampini   PetscInt       ncc,cum_queue,n;
586b3cdcd63SStefano Zampini   PetscMPIInt    commsize;
587674ae819SStefano Zampini   PetscErrorCode ierr;
588674ae819SStefano Zampini 
589674ae819SStefano Zampini   PetscFunctionBegin;
590b3cdcd63SStefano Zampini   if (!graph->setupcalled) SETERRQ(PetscObjectComm((PetscObject)graph->l2gmap),PETSC_ERR_ORDER,"PCBDDCGraphSetUp should be called first");
591b3cdcd63SStefano Zampini   /* quiet return if there isn't any local info */
5924f1b2e48SStefano Zampini   if ((!graph->xadj || !graph->adjncy) && !graph->n_local_subs) {
593674ae819SStefano Zampini     PetscFunctionReturn(0);
594674ae819SStefano Zampini   }
595674ae819SStefano Zampini 
596674ae819SStefano Zampini   /* reset any previous search of connected components */
597df48933bSStefano Zampini   ierr = PetscBTMemzero(graph->nvtxs,graph->touched);CHKERRQ(ierr);
598b3cdcd63SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)graph->l2gmap),&commsize);CHKERRQ(ierr);
5994b2aedd3SStefano Zampini   if (commsize > graph->commsizelimit) {
600b3cdcd63SStefano Zampini     PetscInt i;
601674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
6020cf82fd6SStefano Zampini       if (graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK || !graph->count[i]) {
603df48933bSStefano Zampini         ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
604674ae819SStefano Zampini       }
6054a060362SStefano Zampini     }
606b3cdcd63SStefano Zampini   }
607674ae819SStefano Zampini 
608674ae819SStefano Zampini   /* begin search for connected components */
609674ae819SStefano Zampini   cum_queue = 0;
610674ae819SStefano Zampini   ncc = 0;
611674ae819SStefano Zampini   for (n=0;n<graph->n_subsets;n++) {
612b3cdcd63SStefano Zampini     PetscInt pid = n+1;  /* partition labeled by 0 is discarded */
613b3cdcd63SStefano Zampini     PetscInt found = 0,prev = 0,first = 0,ncc_pid = 0;
614b3cdcd63SStefano Zampini     while (found != graph->subset_size[n]) {
615b3cdcd63SStefano Zampini       PetscInt added = 0;
616b3cdcd63SStefano Zampini       if (!prev) { /* search for new starting dof */
617b3cdcd63SStefano Zampini         while (PetscBTLookup(graph->touched,graph->subset_idxs[n][first])) first++;
618b3cdcd63SStefano Zampini         ierr = PetscBTSet(graph->touched,graph->subset_idxs[n][first]);CHKERRQ(ierr);
619b3cdcd63SStefano Zampini         graph->queue[cum_queue] = graph->subset_idxs[n][first];
620674ae819SStefano Zampini         graph->cptr[ncc] = cum_queue;
621b3cdcd63SStefano Zampini         prev = 1;
622b3cdcd63SStefano Zampini         cum_queue++;
623b3cdcd63SStefano Zampini         found++;
624674ae819SStefano Zampini         ncc_pid++;
625b3cdcd63SStefano Zampini         ncc++;
626674ae819SStefano Zampini       }
627b3cdcd63SStefano Zampini       ierr = PCBDDCGraphComputeCC_Private(graph,pid,graph->queue + cum_queue,prev,&added);CHKERRQ(ierr);
628b3cdcd63SStefano Zampini       if (!added) {
629674ae819SStefano Zampini         graph->subset_ncc[n] = ncc_pid;
630b3cdcd63SStefano Zampini         graph->cptr[ncc] = cum_queue;
631b3cdcd63SStefano Zampini       }
632b3cdcd63SStefano Zampini       prev = added;
633b3cdcd63SStefano Zampini       found += added;
634b3cdcd63SStefano Zampini       cum_queue += added;
635b3cdcd63SStefano Zampini       if (added && found == graph->subset_size[n]) {
636b3cdcd63SStefano Zampini         graph->subset_ncc[n] = ncc_pid;
637b3cdcd63SStefano Zampini         graph->cptr[ncc] = cum_queue;
638b3cdcd63SStefano Zampini       }
639b3cdcd63SStefano Zampini     }
640674ae819SStefano Zampini   }
641674ae819SStefano Zampini   graph->ncc = ncc;
642acc113dbSStefano Zampini   graph->queue_sorted = PETSC_FALSE;
643674ae819SStefano Zampini   PetscFunctionReturn(0);
644674ae819SStefano Zampini }
645674ae819SStefano Zampini 
646674ae819SStefano Zampini #undef __FUNCT__
647674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphSetUp"
648674ae819SStefano Zampini PetscErrorCode PCBDDCGraphSetUp(PCBDDCGraph graph, PetscInt custom_minimal_size, IS neumann_is, IS dirichlet_is, PetscInt n_ISForDofs, IS ISForDofs[], IS custom_primal_vertices)
649674ae819SStefano Zampini {
650674ae819SStefano Zampini   VecScatter     scatter_ctx;
651674ae819SStefano Zampini   Vec            local_vec,local_vec2,global_vec;
652dc456d91SStefano Zampini   IS             to,from,subset,subset_n;
653674ae819SStefano Zampini   MPI_Comm       comm;
654674ae819SStefano Zampini   PetscScalar    *array,*array2;
655674ae819SStefano Zampini   const PetscInt *is_indices;
656dc456d91SStefano Zampini   PetscInt       n_neigh,*neigh,*n_shared,**shared,*queue_global;
657674ae819SStefano Zampini   PetscInt       i,j,k,s,total_counts,nodes_touched,is_size;
658b3cdcd63SStefano Zampini   PetscMPIInt    commsize;
65951b0f6cfSStefano Zampini   PetscBool      same_set,mirrors_found;
6607babac9bSStefano Zampini   PetscErrorCode ierr;
661674ae819SStefano Zampini 
662674ae819SStefano Zampini   PetscFunctionBegin;
663a07ea27aSStefano Zampini   graph->has_dirichlet = PETSC_FALSE;
664a07ea27aSStefano Zampini   if (dirichlet_is) {
665a07ea27aSStefano Zampini     PetscCheckSameComm(graph->l2gmap,1,dirichlet_is,4);
666a07ea27aSStefano Zampini     graph->has_dirichlet = PETSC_TRUE;
667a07ea27aSStefano Zampini   }
668674ae819SStefano Zampini   ierr = PetscObjectGetComm((PetscObject)(graph->l2gmap),&comm);CHKERRQ(ierr);
669b3cdcd63SStefano Zampini   ierr = MPI_Comm_size(comm,&commsize);CHKERRQ(ierr);
670b3cdcd63SStefano Zampini 
671674ae819SStefano Zampini   /* custom_minimal_size */
67214f95afaSStefano Zampini   graph->custom_minimal_size = custom_minimal_size;
673674ae819SStefano Zampini   /* get info l2gmap and allocate work vectors  */
674674ae819SStefano Zampini   ierr = ISLocalToGlobalMappingGetInfo(graph->l2gmap,&n_neigh,&neigh,&n_shared,&shared);CHKERRQ(ierr);
6752635a1d4SStefano Zampini   /* check if we have any local periodic nodes (periodic BCs) */
6762635a1d4SStefano Zampini   mirrors_found = PETSC_FALSE;
6772635a1d4SStefano Zampini   if (graph->nvtxs && n_neigh) {
6782635a1d4SStefano Zampini     for (i=0; i<n_shared[0]; i++) graph->count[shared[0][i]] += 1;
6792635a1d4SStefano Zampini     for (i=0; i<n_shared[0]; i++) {
6802635a1d4SStefano Zampini       if (graph->count[shared[0][i]] > 1) {
6812635a1d4SStefano Zampini         mirrors_found = PETSC_TRUE;
6822635a1d4SStefano Zampini         break;
6832635a1d4SStefano Zampini       }
6842635a1d4SStefano Zampini     }
6852635a1d4SStefano Zampini   }
6862635a1d4SStefano Zampini   /* create some workspace objects */
6872635a1d4SStefano Zampini   local_vec = NULL;
6882635a1d4SStefano Zampini   local_vec2 = NULL;
6892635a1d4SStefano Zampini   global_vec = NULL;
6902635a1d4SStefano Zampini   to = NULL;
6912635a1d4SStefano Zampini   from = NULL;
6922635a1d4SStefano Zampini   scatter_ctx = NULL;
6932635a1d4SStefano Zampini   if (n_ISForDofs || dirichlet_is || neumann_is || custom_primal_vertices) {
694674ae819SStefano Zampini     ierr = VecCreate(PETSC_COMM_SELF,&local_vec);CHKERRQ(ierr);
695674ae819SStefano Zampini     ierr = VecSetSizes(local_vec,PETSC_DECIDE,graph->nvtxs);CHKERRQ(ierr);
696674ae819SStefano Zampini     ierr = VecSetType(local_vec,VECSTANDARD);CHKERRQ(ierr);
697674ae819SStefano Zampini     ierr = VecDuplicate(local_vec,&local_vec2);CHKERRQ(ierr);
698674ae819SStefano Zampini     ierr = VecCreate(comm,&global_vec);CHKERRQ(ierr);
6997fb0e2dbSStefano Zampini     ierr = VecSetSizes(global_vec,PETSC_DECIDE,graph->nvtxs_global);CHKERRQ(ierr);
700674ae819SStefano Zampini     ierr = VecSetType(global_vec,VECSTANDARD);CHKERRQ(ierr);
701674ae819SStefano Zampini     ierr = ISCreateStride(PETSC_COMM_SELF,graph->nvtxs,0,1,&to);CHKERRQ(ierr);
702674ae819SStefano Zampini     ierr = ISLocalToGlobalMappingApplyIS(graph->l2gmap,to,&from);CHKERRQ(ierr);
703674ae819SStefano Zampini     ierr = VecScatterCreate(global_vec,from,local_vec,to,&scatter_ctx);CHKERRQ(ierr);
7042635a1d4SStefano Zampini   } else if (mirrors_found) {
7052635a1d4SStefano Zampini     ierr = ISCreateStride(PETSC_COMM_SELF,graph->nvtxs,0,1,&to);CHKERRQ(ierr);
7062635a1d4SStefano Zampini     ierr = ISLocalToGlobalMappingApplyIS(graph->l2gmap,to,&from);CHKERRQ(ierr);
70749eeff8cSStefano Zampini   }
70851b0f6cfSStefano Zampini   /* compute local mirrors (if any) */
70951b0f6cfSStefano Zampini   if (mirrors_found) {
71051b0f6cfSStefano Zampini     PetscInt *local_indices,*global_indices;
71151b0f6cfSStefano Zampini     /* get arrays of local and global indices */
712785e854fSJed Brown     ierr = PetscMalloc1(graph->nvtxs,&local_indices);CHKERRQ(ierr);
71351b0f6cfSStefano Zampini     ierr = ISGetIndices(to,(const PetscInt**)&is_indices);CHKERRQ(ierr);
71451b0f6cfSStefano Zampini     ierr = PetscMemcpy(local_indices,is_indices,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
71551b0f6cfSStefano Zampini     ierr = ISRestoreIndices(to,(const PetscInt**)&is_indices);CHKERRQ(ierr);
716785e854fSJed Brown     ierr = PetscMalloc1(graph->nvtxs,&global_indices);CHKERRQ(ierr);
71751b0f6cfSStefano Zampini     ierr = ISGetIndices(from,(const PetscInt**)&is_indices);CHKERRQ(ierr);
71851b0f6cfSStefano Zampini     ierr = PetscMemcpy(global_indices,is_indices,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
71951b0f6cfSStefano Zampini     ierr = ISRestoreIndices(from,(const PetscInt**)&is_indices);CHKERRQ(ierr);
72051b0f6cfSStefano Zampini     /* allocate space for mirrors */
7218e5aaad5SJed Brown     ierr = PetscMalloc2(graph->nvtxs,&graph->mirrors,graph->nvtxs,&graph->mirrors_set);CHKERRQ(ierr);
72251b0f6cfSStefano Zampini     ierr = PetscMemzero(graph->mirrors,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
72351b0f6cfSStefano Zampini     graph->mirrors_set[0] = 0;
72451b0f6cfSStefano Zampini 
72551b0f6cfSStefano Zampini     k=0;
72651b0f6cfSStefano Zampini     for (i=0;i<n_shared[0];i++) {
72751b0f6cfSStefano Zampini       j=shared[0][i];
72851b0f6cfSStefano Zampini       if (graph->count[j] > 1) {
72951b0f6cfSStefano Zampini         graph->mirrors[j]++;
73051b0f6cfSStefano Zampini         k++;
73151b0f6cfSStefano Zampini       }
73251b0f6cfSStefano Zampini     }
73351b0f6cfSStefano Zampini     /* allocate space for set of mirrors */
734785e854fSJed Brown     ierr = PetscMalloc1(k,&graph->mirrors_set[0]);CHKERRQ(ierr);
73551b0f6cfSStefano Zampini     for (i=1;i<graph->nvtxs;i++)
73651b0f6cfSStefano Zampini       graph->mirrors_set[i]=graph->mirrors_set[i-1]+graph->mirrors[i-1];
73751b0f6cfSStefano Zampini 
73851b0f6cfSStefano Zampini     /* fill arrays */
73951b0f6cfSStefano Zampini     ierr = PetscMemzero(graph->mirrors,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
74051b0f6cfSStefano Zampini     for (j=0;j<n_shared[0];j++) {
74151b0f6cfSStefano Zampini       i=shared[0][j];
74251b0f6cfSStefano Zampini       if (graph->count[i] > 1)
74351b0f6cfSStefano Zampini         graph->mirrors_set[i][graph->mirrors[i]++]=global_indices[i];
74451b0f6cfSStefano Zampini     }
74551b0f6cfSStefano Zampini     ierr = PetscSortIntWithArray(graph->nvtxs,global_indices,local_indices);CHKERRQ(ierr);
74651b0f6cfSStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
74751b0f6cfSStefano Zampini       if (graph->mirrors[i] > 0) {
74851b0f6cfSStefano Zampini         ierr = PetscFindInt(graph->mirrors_set[i][0],graph->nvtxs,global_indices,&k);CHKERRQ(ierr);
74951b0f6cfSStefano Zampini         j = global_indices[k];
75051b0f6cfSStefano Zampini         while ( k > 0 && global_indices[k-1] == j) k--;
75151b0f6cfSStefano Zampini         for (j=0;j<graph->mirrors[i];j++) {
75251b0f6cfSStefano Zampini           graph->mirrors_set[i][j]=local_indices[k+j];
75351b0f6cfSStefano Zampini         }
75451b0f6cfSStefano Zampini         ierr = PetscSortInt(graph->mirrors[i],graph->mirrors_set[i]);CHKERRQ(ierr);
75551b0f6cfSStefano Zampini       }
75651b0f6cfSStefano Zampini     }
75751b0f6cfSStefano Zampini     ierr = PetscFree(local_indices);CHKERRQ(ierr);
75851b0f6cfSStefano Zampini     ierr = PetscFree(global_indices);CHKERRQ(ierr);
75951b0f6cfSStefano Zampini   }
76051b0f6cfSStefano Zampini   ierr = PetscMemzero(graph->count,graph->nvtxs*sizeof(*graph->count));CHKERRQ(ierr);
761674ae819SStefano Zampini   ierr = ISDestroy(&to);CHKERRQ(ierr);
762674ae819SStefano Zampini   ierr = ISDestroy(&from);CHKERRQ(ierr);
76351b0f6cfSStefano Zampini 
764674ae819SStefano Zampini   /* Count total number of neigh per node */
765674ae819SStefano Zampini   k = 0;
766674ae819SStefano Zampini   for (i=1;i<n_neigh;i++) {
767674ae819SStefano Zampini     k += n_shared[i];
768674ae819SStefano Zampini     for (j=0;j<n_shared[i];j++) {
769674ae819SStefano Zampini       graph->count[shared[i][j]] += 1;
770674ae819SStefano Zampini     }
771674ae819SStefano Zampini   }
772674ae819SStefano Zampini   /* Allocate space for storing the set of neighbours for each node */
773674ae819SStefano Zampini   if (graph->nvtxs) {
774785e854fSJed Brown     ierr = PetscMalloc1(k,&graph->neighbours_set[0]);CHKERRQ(ierr);
775674ae819SStefano Zampini   }
776674ae819SStefano Zampini   for (i=1;i<graph->nvtxs;i++) { /* dont count myself */
777674ae819SStefano Zampini     graph->neighbours_set[i]=graph->neighbours_set[i-1]+graph->count[i-1];
778674ae819SStefano Zampini   }
779674ae819SStefano Zampini   /* Get information for sharing subdomains */
780674ae819SStefano Zampini   ierr = PetscMemzero(graph->count,graph->nvtxs*sizeof(*graph->count));CHKERRQ(ierr);
781674ae819SStefano Zampini   for (i=1;i<n_neigh;i++) { /* dont count myself */
782674ae819SStefano Zampini     s = n_shared[i];
783674ae819SStefano Zampini     for (j=0;j<s;j++) {
784674ae819SStefano Zampini       k = shared[i][j];
785674ae819SStefano Zampini       graph->neighbours_set[k][graph->count[k]] = neigh[i];
786674ae819SStefano Zampini       graph->count[k] += 1;
787674ae819SStefano Zampini     }
788674ae819SStefano Zampini   }
789674ae819SStefano Zampini   /* sort set of sharing subdomains */
790674ae819SStefano Zampini   for (i=0;i<graph->nvtxs;i++) {
79193fb5fd3SStefano Zampini     ierr = PetscSortRemoveDupsInt(&graph->count[i],graph->neighbours_set[i]);CHKERRQ(ierr);
792674ae819SStefano Zampini   }
7937fb0e2dbSStefano Zampini   /* free memory allocated by ISLocalToGlobalMappingGetInfo */
7947fb0e2dbSStefano Zampini   ierr = ISLocalToGlobalMappingRestoreInfo(graph->l2gmap,&n_neigh,&neigh,&n_shared,&shared);CHKERRQ(ierr);
7957fb0e2dbSStefano Zampini 
79667c9da69SStefano Zampini   /*
79767c9da69SStefano Zampini      Get info for dofs splitting
7985777c63cSStefano Zampini      User can specify just a subset; an additional field is considered as a complementary field
79967c9da69SStefano Zampini   */
800b3cdcd63SStefano Zampini   for (i=0;i<graph->nvtxs;i++) graph->which_dof[i] = n_ISForDofs; /* by default a dof belongs to the complement set */
8015777c63cSStefano Zampini   if (n_ISForDofs) {
8025777c63cSStefano Zampini     ierr = VecSet(local_vec,-1.0);CHKERRQ(ierr);
8035777c63cSStefano Zampini   }
804674ae819SStefano Zampini   for (i=0;i<n_ISForDofs;i++) {
8055777c63cSStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
80663602bcaSStefano Zampini     ierr = ISGetLocalSize(ISForDofs[i],&is_size);CHKERRQ(ierr);
807674ae819SStefano Zampini     ierr = ISGetIndices(ISForDofs[i],(const PetscInt**)&is_indices);CHKERRQ(ierr);
808674ae819SStefano Zampini     for (j=0;j<is_size;j++) {
809d50ed60aSStefano Zampini       if (is_indices[j] > -1 && is_indices[j] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
810674ae819SStefano Zampini         graph->which_dof[is_indices[j]] = i;
8115777c63cSStefano Zampini         array[is_indices[j]] = 1.*i;
812674ae819SStefano Zampini       }
81367c9da69SStefano Zampini     }
814674ae819SStefano Zampini     ierr = ISRestoreIndices(ISForDofs[i],(const PetscInt**)&is_indices);CHKERRQ(ierr);
8155777c63cSStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
8165777c63cSStefano Zampini   }
8175777c63cSStefano Zampini   /* Check consistency among neighbours */
8185777c63cSStefano Zampini   if (n_ISForDofs) {
8195777c63cSStefano Zampini     ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8205777c63cSStefano Zampini     ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8215777c63cSStefano Zampini     ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8225777c63cSStefano Zampini     ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8235777c63cSStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
8245777c63cSStefano Zampini     ierr = VecGetArray(local_vec2,&array2);CHKERRQ(ierr);
8255777c63cSStefano Zampini     for (i=0;i<graph->nvtxs;i++){
8265777c63cSStefano Zampini       PetscInt field1,field2;
8275777c63cSStefano Zampini 
8285777c63cSStefano Zampini       field1 = (PetscInt)PetscRealPart(array[i]);
8295777c63cSStefano Zampini       field2 = (PetscInt)PetscRealPart(array2[i]);
8306c4ed002SBarry 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);
8315777c63cSStefano Zampini     }
8325777c63cSStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
8335777c63cSStefano Zampini     ierr = VecRestoreArray(local_vec2,&array2);CHKERRQ(ierr);
834674ae819SStefano Zampini   }
8357fb0e2dbSStefano Zampini   if (neumann_is || dirichlet_is) {
836674ae819SStefano Zampini     /* Take into account Neumann nodes */
837674ae819SStefano Zampini     ierr = VecSet(local_vec,0.0);CHKERRQ(ierr);
838674ae819SStefano Zampini     ierr = VecSet(local_vec2,0.0);CHKERRQ(ierr);
839674ae819SStefano Zampini     if (neumann_is) {
840674ae819SStefano Zampini       ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
84182ba6b80SStefano Zampini       ierr = ISGetLocalSize(neumann_is,&is_size);CHKERRQ(ierr);
842674ae819SStefano Zampini       ierr = ISGetIndices(neumann_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
843674ae819SStefano Zampini       for (i=0;i<is_size;i++) {
844d50ed60aSStefano Zampini         if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
845674ae819SStefano Zampini           array[is_indices[i]] = 1.0;
846674ae819SStefano Zampini         }
8473c629590SStefano Zampini       }
848674ae819SStefano Zampini       ierr = ISRestoreIndices(neumann_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
849674ae819SStefano Zampini       ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
850674ae819SStefano Zampini     }
851674ae819SStefano Zampini     /* Neumann nodes: impose consistency among neighbours */
852674ae819SStefano Zampini     ierr = VecSet(global_vec,0.0);CHKERRQ(ierr);
853674ae819SStefano Zampini     ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
854674ae819SStefano Zampini     ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
855674ae819SStefano Zampini     ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
856674ae819SStefano Zampini     ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
857674ae819SStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
858674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
8593c629590SStefano Zampini       if (PetscRealPart(array[i]) > 0.1) {
8600cf82fd6SStefano Zampini         graph->special_dof[i] = PCBDDCGRAPH_NEUMANN_MARK;
861674ae819SStefano Zampini       }
862674ae819SStefano Zampini     }
863674ae819SStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
864674ae819SStefano Zampini     /* Take into account Dirichlet nodes */
865674ae819SStefano Zampini     ierr = VecSet(local_vec2,0.0);CHKERRQ(ierr);
866674ae819SStefano Zampini     if (dirichlet_is) {
867674ae819SStefano Zampini       ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
868674ae819SStefano Zampini       ierr = VecGetArray(local_vec2,&array2);CHKERRQ(ierr);
86982ba6b80SStefano Zampini       ierr = ISGetLocalSize(dirichlet_is,&is_size);CHKERRQ(ierr);
870674ae819SStefano Zampini       ierr = ISGetIndices(dirichlet_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
871674ae819SStefano Zampini       for (i=0;i<is_size;i++){
872d50ed60aSStefano Zampini         if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
873674ae819SStefano Zampini           k = is_indices[i];
874e477d8c3SStefano Zampini           if (!PetscBTLookup(graph->touched,k)) {
875af25d912SStefano 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);
876674ae819SStefano Zampini             array2[k] = 1.0;
877674ae819SStefano Zampini           }
878674ae819SStefano Zampini         }
8793c629590SStefano Zampini       }
880674ae819SStefano Zampini       ierr = ISRestoreIndices(dirichlet_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
881674ae819SStefano Zampini       ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
882674ae819SStefano Zampini       ierr = VecRestoreArray(local_vec2,&array2);CHKERRQ(ierr);
883674ae819SStefano Zampini     }
884674ae819SStefano Zampini     /* Dirichlet nodes: impose consistency among neighbours */
885674ae819SStefano Zampini     ierr = VecSet(global_vec,0.0);CHKERRQ(ierr);
886674ae819SStefano Zampini     ierr = VecScatterBegin(scatter_ctx,local_vec2,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
887674ae819SStefano Zampini     ierr = VecScatterEnd(scatter_ctx,local_vec2,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
888674ae819SStefano Zampini     ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
889674ae819SStefano Zampini     ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
890674ae819SStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
891674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
8923c629590SStefano Zampini       if (PetscRealPart(array[i]) > 0.1) {
8934b2aedd3SStefano Zampini         if (commsize > graph->commsizelimit) { /* dirichlet nodes treated as internal */
894df48933bSStefano Zampini           ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
895b3cdcd63SStefano Zampini           graph->subset[i] = 0;
896b3cdcd63SStefano Zampini         }
8970cf82fd6SStefano Zampini         graph->special_dof[i] = PCBDDCGRAPH_DIRICHLET_MARK;
898674ae819SStefano Zampini       }
899674ae819SStefano Zampini     }
900674ae819SStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
9017fb0e2dbSStefano Zampini   }
90208a5cf49SStefano Zampini   /* mark local periodic nodes (if any) and adapt CSR graph (if any) */
90351b0f6cfSStefano Zampini   if (graph->mirrors) {
90451b0f6cfSStefano Zampini     for (i=0;i<graph->nvtxs;i++)
90551b0f6cfSStefano Zampini       if (graph->mirrors[i])
9060cf82fd6SStefano Zampini         graph->special_dof[i] = PCBDDCGRAPH_LOCAL_PERIODIC_MARK;
90751b0f6cfSStefano Zampini 
90808a5cf49SStefano Zampini     if (graph->xadj && graph->adjncy) {
90908a5cf49SStefano Zampini       PetscInt *new_xadj,*new_adjncy;
91051b0f6cfSStefano Zampini       /* sort CSR graph */
91151b0f6cfSStefano Zampini       for (i=0;i<graph->nvtxs;i++)
91251b0f6cfSStefano Zampini         ierr = PetscSortInt(graph->xadj[i+1]-graph->xadj[i],&graph->adjncy[graph->xadj[i]]);CHKERRQ(ierr);
91351b0f6cfSStefano Zampini 
91451b0f6cfSStefano Zampini       /* adapt local CSR graph in case of local periodicity */
91551b0f6cfSStefano Zampini       k = 0;
91651b0f6cfSStefano Zampini       for (i=0;i<graph->nvtxs;i++)
91751b0f6cfSStefano Zampini         for (j=graph->xadj[i];j<graph->xadj[i+1];j++)
91851b0f6cfSStefano Zampini           k += graph->mirrors[graph->adjncy[j]];
91951b0f6cfSStefano Zampini 
920854ce69bSBarry Smith       ierr = PetscMalloc1(graph->nvtxs+1,&new_xadj);CHKERRQ(ierr);
921854ce69bSBarry Smith       ierr = PetscMalloc1(k+graph->xadj[graph->nvtxs],&new_adjncy);CHKERRQ(ierr);
92251b0f6cfSStefano Zampini       new_xadj[0] = 0;
92351b0f6cfSStefano Zampini       for (i=0;i<graph->nvtxs;i++) {
92451b0f6cfSStefano Zampini         k = graph->xadj[i+1]-graph->xadj[i];
92551b0f6cfSStefano Zampini         ierr = PetscMemcpy(&new_adjncy[new_xadj[i]],&graph->adjncy[graph->xadj[i]],k*sizeof(PetscInt));CHKERRQ(ierr);
92651b0f6cfSStefano Zampini         new_xadj[i+1] = new_xadj[i]+k;
92751b0f6cfSStefano Zampini         for (j=graph->xadj[i];j<graph->xadj[i+1];j++) {
92851b0f6cfSStefano Zampini           k = graph->mirrors[graph->adjncy[j]];
92951b0f6cfSStefano Zampini           ierr = PetscMemcpy(&new_adjncy[new_xadj[i+1]],graph->mirrors_set[graph->adjncy[j]],k*sizeof(PetscInt));CHKERRQ(ierr);
93051b0f6cfSStefano Zampini           new_xadj[i+1] += k;
93151b0f6cfSStefano Zampini         }
93251b0f6cfSStefano Zampini         k = new_xadj[i+1]-new_xadj[i];
93351b0f6cfSStefano Zampini         ierr = PetscSortRemoveDupsInt(&k,&new_adjncy[new_xadj[i]]);CHKERRQ(ierr);
93451b0f6cfSStefano Zampini         new_xadj[i+1] = new_xadj[i]+k;
93551b0f6cfSStefano Zampini       }
93651b0f6cfSStefano Zampini       /* set new CSR into graph */
93751b0f6cfSStefano Zampini       ierr = PetscFree(graph->xadj);CHKERRQ(ierr);
93851b0f6cfSStefano Zampini       ierr = PetscFree(graph->adjncy);CHKERRQ(ierr);
93951b0f6cfSStefano Zampini       graph->xadj = new_xadj;
94051b0f6cfSStefano Zampini       graph->adjncy = new_adjncy;
94151b0f6cfSStefano Zampini     }
94208a5cf49SStefano Zampini   }
94351b0f6cfSStefano Zampini 
94417eb1463SStefano Zampini   /* mark special nodes (if any) -> each will become a single node equivalence class */
945674ae819SStefano Zampini   if (custom_primal_vertices) {
94617eb1463SStefano Zampini     ierr = VecSet(local_vec,0.0);CHKERRQ(ierr);
9479b70a373SStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
94863602bcaSStefano Zampini     ierr = ISGetLocalSize(custom_primal_vertices,&is_size);CHKERRQ(ierr);
949674ae819SStefano Zampini     ierr = ISGetIndices(custom_primal_vertices,(const PetscInt**)&is_indices);CHKERRQ(ierr);
950674ae819SStefano Zampini     for (i=0;i<is_size;i++){
9519b70a373SStefano Zampini       if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
9529b70a373SStefano Zampini         array[is_indices[i]] = 1.0;
9539b70a373SStefano Zampini       }
954674ae819SStefano Zampini     }
955674ae819SStefano Zampini     ierr = ISRestoreIndices(custom_primal_vertices,(const PetscInt**)&is_indices);CHKERRQ(ierr);
9569b70a373SStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
9579b70a373SStefano Zampini     /* special nodes: impose consistency among neighbours */
9589b70a373SStefano Zampini     ierr = VecSet(global_vec,0.0);CHKERRQ(ierr);
9599b70a373SStefano Zampini     ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9609b70a373SStefano Zampini     ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9619b70a373SStefano Zampini     ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9629b70a373SStefano Zampini     ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9639b70a373SStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
9649b70a373SStefano Zampini     j = 0;
9659b70a373SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
9669b70a373SStefano Zampini       if (PetscRealPart(array[i]) > 0.1 && graph->special_dof[i] != PCBDDCGRAPH_DIRICHLET_MARK) {
9679b70a373SStefano Zampini         graph->special_dof[i] = PCBDDCGRAPH_SPECIAL_MARK-j;
9689b70a373SStefano Zampini         j++;
9699b70a373SStefano Zampini       }
9709b70a373SStefano Zampini     }
9719b70a373SStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
97217eb1463SStefano Zampini   }
9739b70a373SStefano Zampini 
9744b2aedd3SStefano Zampini   /* mark interior nodes (if commsize > graph->commsizelimit) as touched and belonging to partition number 0 */
9754b2aedd3SStefano Zampini   if (commsize > graph->commsizelimit) {
976674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
977674ae819SStefano Zampini       if (!graph->count[i]) {
978df48933bSStefano Zampini         ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
979674ae819SStefano Zampini         graph->subset[i] = 0;
980674ae819SStefano Zampini       }
981674ae819SStefano Zampini     }
982b3cdcd63SStefano Zampini   }
983b3cdcd63SStefano Zampini 
984674ae819SStefano Zampini   /* init graph structure and compute default subsets */
985674ae819SStefano Zampini   nodes_touched = 0;
986674ae819SStefano Zampini   for (i=0;i<graph->nvtxs;i++) {
987df48933bSStefano Zampini     if (PetscBTLookup(graph->touched,i)) {
988674ae819SStefano Zampini       nodes_touched++;
989674ae819SStefano Zampini     }
990674ae819SStefano Zampini   }
991674ae819SStefano Zampini   i = 0;
992674ae819SStefano Zampini   graph->ncc = 0;
993674ae819SStefano Zampini   total_counts = 0;
9947babac9bSStefano Zampini 
9957babac9bSStefano Zampini   /* allocated space for queues */
9964b2aedd3SStefano Zampini   if (commsize == graph->commsizelimit) {
9977babac9bSStefano Zampini     ierr = PetscMalloc2(graph->nvtxs+1,&graph->cptr,graph->nvtxs,&graph->queue);CHKERRQ(ierr);
9987babac9bSStefano Zampini   } else {
9997babac9bSStefano Zampini     PetscInt nused = graph->nvtxs - nodes_touched;
10007babac9bSStefano Zampini     ierr = PetscMalloc2(nused+1,&graph->cptr,nused,&graph->queue);CHKERRQ(ierr);
10017babac9bSStefano Zampini   }
10027babac9bSStefano Zampini 
1003674ae819SStefano Zampini   while (nodes_touched<graph->nvtxs) {
1004674ae819SStefano Zampini     /*  find first untouched node in local ordering */
1005b3cdcd63SStefano Zampini     while (PetscBTLookup(graph->touched,i)) i++;
1006df48933bSStefano Zampini     ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
1007674ae819SStefano Zampini     graph->subset[i] = graph->ncc+1;
1008674ae819SStefano Zampini     graph->cptr[graph->ncc] = total_counts;
1009674ae819SStefano Zampini     graph->queue[total_counts] = i;
1010674ae819SStefano Zampini     total_counts++;
1011674ae819SStefano Zampini     nodes_touched++;
1012674ae819SStefano Zampini     /* now find all other nodes having the same set of sharing subdomains */
1013674ae819SStefano Zampini     for (j=i+1;j<graph->nvtxs;j++) {
101474e413f5SStefano Zampini       /* check for same number of sharing subdomains, dof number and same special mark */
1015df48933bSStefano 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]) {
1016674ae819SStefano Zampini         /* check for same set of sharing subdomains */
1017674ae819SStefano Zampini         same_set = PETSC_TRUE;
1018674ae819SStefano Zampini         for (k=0;k<graph->count[j];k++){
1019674ae819SStefano Zampini           if (graph->neighbours_set[i][k] != graph->neighbours_set[j][k]) {
1020674ae819SStefano Zampini             same_set = PETSC_FALSE;
1021674ae819SStefano Zampini           }
1022674ae819SStefano Zampini         }
1023674ae819SStefano Zampini         /* I found a friend of mine */
1024674ae819SStefano Zampini         if (same_set) {
1025df48933bSStefano Zampini           ierr = PetscBTSet(graph->touched,j);CHKERRQ(ierr);
1026b3cdcd63SStefano Zampini           graph->subset[j] = graph->ncc+1;
1027674ae819SStefano Zampini           nodes_touched++;
1028674ae819SStefano Zampini           graph->queue[total_counts] = j;
1029674ae819SStefano Zampini           total_counts++;
1030674ae819SStefano Zampini         }
1031674ae819SStefano Zampini       }
1032674ae819SStefano Zampini     }
1033674ae819SStefano Zampini     graph->ncc++;
1034674ae819SStefano Zampini   }
1035b3cdcd63SStefano 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 */
1036674ae819SStefano Zampini   graph->n_subsets = graph->ncc;
1037785e854fSJed Brown   ierr = PetscMalloc1(graph->n_subsets,&graph->subset_ncc);CHKERRQ(ierr);
1038674ae819SStefano Zampini   for (i=0;i<graph->n_subsets;i++) {
1039674ae819SStefano Zampini     graph->subset_ncc[i] = 1;
1040674ae819SStefano Zampini   }
1041674ae819SStefano Zampini   /* final pointer */
1042674ae819SStefano Zampini   graph->cptr[graph->ncc] = total_counts;
1043ec1c413dSStefano Zampini 
1044b3cdcd63SStefano Zampini   /* For consistency reasons (among neighbours), I need to sort (by global ordering) each connected component */
1045ec1c413dSStefano Zampini   /* Get a reference node (min index in global ordering) for each subset for tagging messages */
1046dc456d91SStefano Zampini   ierr = PetscMalloc1(graph->ncc,&graph->subset_ref_node);CHKERRQ(ierr);
1047dc456d91SStefano Zampini   ierr = PetscMalloc1(graph->cptr[graph->ncc],&queue_global);CHKERRQ(ierr);
10483a5b03d1SStefano Zampini   ierr = ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_global);CHKERRQ(ierr);
1049ec1c413dSStefano Zampini   for (j=0;j<graph->ncc;j++) {
1050ec1c413dSStefano Zampini     ierr = PetscSortIntWithArray(graph->cptr[j+1]-graph->cptr[j],&queue_global[graph->cptr[j]],&graph->queue[graph->cptr[j]]);CHKERRQ(ierr);
1051dc456d91SStefano Zampini     graph->subset_ref_node[j] = graph->queue[graph->cptr[j]];
1052f0321c90SStefano Zampini   }
1053dc456d91SStefano Zampini   ierr = PetscFree(queue_global);CHKERRQ(ierr);
1054acc113dbSStefano Zampini   graph->queue_sorted = PETSC_TRUE;
1055b3cdcd63SStefano Zampini   /* save information on subsets (needed when analyzing the connected components) */
10562f566a09SStefano Zampini   if (graph->ncc) {
1057b3cdcd63SStefano Zampini     ierr = PetscMalloc2(graph->ncc,&graph->subset_size,graph->ncc,&graph->subset_idxs);CHKERRQ(ierr);
1058b3cdcd63SStefano Zampini     ierr = PetscMalloc1(graph->cptr[graph->ncc],&graph->subset_idxs[0]);CHKERRQ(ierr);
1059b3cdcd63SStefano Zampini     ierr = PetscMemzero(graph->subset_idxs[0],graph->cptr[graph->ncc]*sizeof(PetscInt));CHKERRQ(ierr);
1060ec1c413dSStefano Zampini     for (j=1;j<graph->ncc;j++) {
1061b3cdcd63SStefano Zampini       graph->subset_size[j-1] = graph->cptr[j] - graph->cptr[j-1];
1062b3cdcd63SStefano Zampini       graph->subset_idxs[j] = graph->subset_idxs[j-1] + graph->subset_size[j-1];
1063ec1c413dSStefano Zampini     }
1064b3cdcd63SStefano Zampini     graph->subset_size[graph->ncc-1] = graph->cptr[graph->ncc] - graph->cptr[graph->ncc-1];
1065b3cdcd63SStefano Zampini     ierr = PetscMemcpy(graph->subset_idxs[0],graph->queue,graph->cptr[graph->ncc]*sizeof(PetscInt));CHKERRQ(ierr);
1066ec1c413dSStefano Zampini   }
1067b3cdcd63SStefano Zampini 
1068f0321c90SStefano Zampini   /* renumber reference nodes */
1069dc456d91SStefano Zampini   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)(graph->l2gmap)),graph->ncc,graph->subset_ref_node,PETSC_COPY_VALUES,&subset_n);CHKERRQ(ierr);
1070dc456d91SStefano Zampini   ierr = ISLocalToGlobalMappingApplyIS(graph->l2gmap,subset_n,&subset);CHKERRQ(ierr);
1071dc456d91SStefano Zampini   ierr = ISDestroy(&subset_n);CHKERRQ(ierr);
10726583bcc1SStefano Zampini   ierr = ISRenumber(subset,NULL,NULL,&subset_n);CHKERRQ(ierr);
1073dc456d91SStefano Zampini   ierr = ISDestroy(&subset);CHKERRQ(ierr);
1074dc456d91SStefano Zampini   ierr = ISGetLocalSize(subset_n,&k);CHKERRQ(ierr);
10756c4ed002SBarry Smith   if (k != graph->ncc) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid size of new subset! %D != %D",k,graph->ncc);
1076dc456d91SStefano Zampini   ierr = ISGetIndices(subset_n,&is_indices);CHKERRQ(ierr);
1077dc456d91SStefano Zampini   ierr = PetscMemcpy(graph->subset_ref_node,is_indices,graph->ncc*sizeof(PetscInt));CHKERRQ(ierr);
1078dc456d91SStefano Zampini   ierr = ISRestoreIndices(subset_n,&is_indices);CHKERRQ(ierr);
1079dc456d91SStefano Zampini   ierr = ISDestroy(&subset_n);CHKERRQ(ierr);
1080ec1c413dSStefano Zampini 
1081ec1c413dSStefano Zampini   /* free workspace */
1082674ae819SStefano Zampini   ierr = VecDestroy(&local_vec);CHKERRQ(ierr);
1083674ae819SStefano Zampini   ierr = VecDestroy(&local_vec2);CHKERRQ(ierr);
1084674ae819SStefano Zampini   ierr = VecDestroy(&global_vec);CHKERRQ(ierr);
1085674ae819SStefano Zampini   ierr = VecScatterDestroy(&scatter_ctx);CHKERRQ(ierr);
1086b3cdcd63SStefano Zampini   graph->setupcalled = PETSC_TRUE;
1087674ae819SStefano Zampini   PetscFunctionReturn(0);
1088674ae819SStefano Zampini }
1089674ae819SStefano Zampini 
1090674ae819SStefano Zampini #undef __FUNCT__
1091674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphResetCSR"
1092674ae819SStefano Zampini PetscErrorCode PCBDDCGraphResetCSR(PCBDDCGraph graph)
1093674ae819SStefano Zampini {
1094674ae819SStefano Zampini   PetscErrorCode ierr;
1095674ae819SStefano Zampini 
1096674ae819SStefano Zampini   PetscFunctionBegin;
1097674ae819SStefano Zampini   ierr = PetscFree(graph->xadj);CHKERRQ(ierr);
1098674ae819SStefano Zampini   ierr = PetscFree(graph->adjncy);CHKERRQ(ierr);
1099575ad6abSStefano Zampini   graph->nvtxs_csr = 0;
1100674ae819SStefano Zampini   PetscFunctionReturn(0);
1101674ae819SStefano Zampini }
1102674ae819SStefano Zampini 
1103674ae819SStefano Zampini #undef __FUNCT__
1104674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphReset"
1105674ae819SStefano Zampini PetscErrorCode PCBDDCGraphReset(PCBDDCGraph graph)
1106674ae819SStefano Zampini {
1107674ae819SStefano Zampini   PetscErrorCode ierr;
1108674ae819SStefano Zampini 
1109674ae819SStefano Zampini   PetscFunctionBegin;
1110b3cdcd63SStefano Zampini   graph->setupcalled = PETSC_FALSE;
1111674ae819SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&graph->l2gmap);CHKERRQ(ierr);
1112674ae819SStefano Zampini   ierr = PetscFree(graph->subset_ncc);CHKERRQ(ierr);
11133a5b03d1SStefano Zampini   ierr = PetscFree(graph->subset_ref_node);CHKERRQ(ierr);
1114674ae819SStefano Zampini   if (graph->nvtxs) {
1115674ae819SStefano Zampini     ierr = PetscFree(graph->neighbours_set[0]);CHKERRQ(ierr);
1116674ae819SStefano Zampini   }
1117df48933bSStefano Zampini   ierr = PetscBTDestroy(&graph->touched);CHKERRQ(ierr);
11187babac9bSStefano Zampini   ierr = PetscFree5(graph->count,
1119674ae819SStefano Zampini                     graph->neighbours_set,
1120674ae819SStefano Zampini                     graph->subset,
1121674ae819SStefano Zampini                     graph->which_dof,
1122df48933bSStefano Zampini                     graph->special_dof);CHKERRQ(ierr);
11237babac9bSStefano Zampini   ierr = PetscFree2(graph->cptr,graph->queue);CHKERRQ(ierr);
112451b0f6cfSStefano Zampini   if (graph->mirrors) {
112551b0f6cfSStefano Zampini     ierr = PetscFree(graph->mirrors_set[0]);CHKERRQ(ierr);
112651b0f6cfSStefano Zampini   }
112751b0f6cfSStefano Zampini   ierr = PetscFree2(graph->mirrors,graph->mirrors_set);CHKERRQ(ierr);
1128b3cdcd63SStefano Zampini   if (graph->subset_idxs) {
1129b3cdcd63SStefano Zampini     ierr = PetscFree(graph->subset_idxs[0]);CHKERRQ(ierr);
1130ec1c413dSStefano Zampini   }
1131b3cdcd63SStefano Zampini   ierr = PetscFree2(graph->subset_size,graph->subset_idxs);CHKERRQ(ierr);
1132a07ea27aSStefano Zampini   ierr = ISDestroy(&graph->dirdofs);CHKERRQ(ierr);
1133d62866d3SStefano Zampini   ierr = ISDestroy(&graph->dirdofsB);CHKERRQ(ierr);
11344f1b2e48SStefano Zampini   if (graph->n_local_subs) {
11354f1b2e48SStefano Zampini     ierr = PetscFree(graph->local_subs);CHKERRQ(ierr);
11364f1b2e48SStefano Zampini   }
1137a07ea27aSStefano Zampini   graph->has_dirichlet = PETSC_FALSE;
1138674ae819SStefano Zampini   graph->nvtxs = 0;
11397babac9bSStefano Zampini   graph->nvtxs_global = 0;
1140674ae819SStefano Zampini   graph->n_subsets = 0;
1141674ae819SStefano Zampini   graph->custom_minimal_size = 1;
11424f1b2e48SStefano Zampini   graph->n_local_subs = 0;
1143674ae819SStefano Zampini   PetscFunctionReturn(0);
1144674ae819SStefano Zampini }
1145674ae819SStefano Zampini 
1146674ae819SStefano Zampini #undef __FUNCT__
1147674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphInit"
11487fb0e2dbSStefano Zampini PetscErrorCode PCBDDCGraphInit(PCBDDCGraph graph, ISLocalToGlobalMapping l2gmap, PetscInt N)
1149674ae819SStefano Zampini {
1150674ae819SStefano Zampini   PetscInt       n;
1151674ae819SStefano Zampini   PetscErrorCode ierr;
1152674ae819SStefano Zampini 
1153674ae819SStefano Zampini   PetscFunctionBegin;
1154674ae819SStefano Zampini   PetscValidPointer(graph,1);
1155674ae819SStefano Zampini   PetscValidHeaderSpecific(l2gmap,IS_LTOGM_CLASSID,2);
11567babac9bSStefano Zampini   PetscValidLogicalCollectiveInt(l2gmap,N,3);
11578e61c736SStefano Zampini   /* raise an error if already allocated */
11586c4ed002SBarry Smith   if (graph->nvtxs_global) SETERRQ(PetscObjectComm((PetscObject)l2gmap),PETSC_ERR_PLIB,"BDDCGraph already initialized");
1159674ae819SStefano Zampini   /* set number of vertices */
1160674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)l2gmap);CHKERRQ(ierr);
1161674ae819SStefano Zampini   graph->l2gmap = l2gmap;
1162674ae819SStefano Zampini   ierr = ISLocalToGlobalMappingGetSize(l2gmap,&n);CHKERRQ(ierr);
1163674ae819SStefano Zampini   graph->nvtxs = n;
11647fb0e2dbSStefano Zampini   graph->nvtxs_global = N;
1165674ae819SStefano Zampini   /* allocate used space */
1166df48933bSStefano Zampini   ierr = PetscBTCreate(graph->nvtxs,&graph->touched);CHKERRQ(ierr);
11677babac9bSStefano Zampini   ierr = PetscMalloc5(graph->nvtxs,&graph->count,
11688e5aaad5SJed Brown                       graph->nvtxs,&graph->neighbours_set,
11698e5aaad5SJed Brown                       graph->nvtxs,&graph->subset,
11708e5aaad5SJed Brown                       graph->nvtxs,&graph->which_dof,
11718e5aaad5SJed Brown                       graph->nvtxs,&graph->special_dof);CHKERRQ(ierr);
1172674ae819SStefano Zampini   /* zeroes memory */
1173674ae819SStefano Zampini   ierr = PetscMemzero(graph->count,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
1174674ae819SStefano Zampini   ierr = PetscMemzero(graph->subset,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
117563602bcaSStefano Zampini   /* use -1 as a default value for which_dof array */
117663602bcaSStefano Zampini   for (n=0;n<graph->nvtxs;n++) graph->which_dof[n] = -1;
117774e413f5SStefano Zampini   ierr = PetscMemzero(graph->special_dof,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
1178674ae819SStefano Zampini   /* zeroes first pointer to neighbour set */
1179674ae819SStefano Zampini   if (graph->nvtxs) {
1180674ae819SStefano Zampini     graph->neighbours_set[0] = 0;
1181674ae819SStefano Zampini   }
1182674ae819SStefano Zampini   /* zeroes workspace for values of ncc */
1183674ae819SStefano Zampini   graph->subset_ncc = 0;
11843a5b03d1SStefano Zampini   graph->subset_ref_node = 0;
1185674ae819SStefano Zampini   PetscFunctionReturn(0);
1186674ae819SStefano Zampini }
1187674ae819SStefano Zampini 
1188674ae819SStefano Zampini #undef __FUNCT__
1189674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphDestroy"
1190674ae819SStefano Zampini PetscErrorCode PCBDDCGraphDestroy(PCBDDCGraph* graph)
1191674ae819SStefano Zampini {
1192674ae819SStefano Zampini   PetscErrorCode ierr;
1193674ae819SStefano Zampini 
1194674ae819SStefano Zampini   PetscFunctionBegin;
1195674ae819SStefano Zampini   ierr = PCBDDCGraphReset(*graph);CHKERRQ(ierr);
1196674ae819SStefano Zampini   ierr = PetscFree(*graph);CHKERRQ(ierr);
1197674ae819SStefano Zampini   PetscFunctionReturn(0);
1198674ae819SStefano Zampini }
1199674ae819SStefano Zampini 
1200674ae819SStefano Zampini #undef __FUNCT__
1201674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphCreate"
1202674ae819SStefano Zampini PetscErrorCode PCBDDCGraphCreate(PCBDDCGraph *graph)
1203674ae819SStefano Zampini {
1204674ae819SStefano Zampini   PCBDDCGraph    new_graph;
1205674ae819SStefano Zampini   PetscErrorCode ierr;
1206674ae819SStefano Zampini 
1207674ae819SStefano Zampini   PetscFunctionBegin;
1208854ce69bSBarry Smith   ierr = PetscNew(&new_graph);CHKERRQ(ierr);
1209674ae819SStefano Zampini   new_graph->custom_minimal_size = 1;
12104b2aedd3SStefano Zampini   new_graph->commsizelimit = 1;
1211674ae819SStefano Zampini   *graph = new_graph;
1212674ae819SStefano Zampini   PetscFunctionReturn(0);
1213674ae819SStefano Zampini }
1214