xref: /petsc/src/ksp/pc/impls/bddc/bddcgraph.c (revision fb59768521fcfea99af49e491227530eb032fa67)
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);
81d4890cceSStefano 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       }
110d4890cceSStefano 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     }
144d4890cceSStefano Zampini     if (verbosity_level > 1) {
145674ae819SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"):");CHKERRQ(ierr);
146d4890cceSStefano 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       }
154d4890cceSStefano Zampini     } else {
155d4890cceSStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,")");CHKERRQ(ierr);
156d4890cceSStefano 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__
167c8272957SStefano Zampini #define __FUNCT__ "PCBDDCGraphRestoreCandidatesIS"
168c8272957SStefano Zampini PetscErrorCode PCBDDCGraphRestoreCandidatesIS(PCBDDCGraph graph, PetscInt *n_faces, IS *FacesIS[], PetscInt *n_edges, IS *EdgesIS[], IS *VerticesIS)
169c8272957SStefano Zampini {
170c8272957SStefano Zampini   PetscInt       i;
171c8272957SStefano Zampini   PetscErrorCode ierr;
172c8272957SStefano Zampini 
173c8272957SStefano Zampini   PetscFunctionBegin;
174c8272957SStefano Zampini   if (n_faces) {
175c8272957SStefano Zampini     if (FacesIS) {
176c8272957SStefano Zampini       for (i=0;i<*n_faces;i++) {
177c8272957SStefano Zampini         ierr = ISDestroy(&((*FacesIS)[i]));CHKERRQ(ierr);
178c8272957SStefano Zampini       }
179c8272957SStefano Zampini       ierr = PetscFree(*FacesIS);CHKERRQ(ierr);
180c8272957SStefano Zampini     }
181c8272957SStefano Zampini     *n_faces = 0;
182c8272957SStefano Zampini   }
183c8272957SStefano Zampini   if (n_edges) {
184c8272957SStefano Zampini     if (EdgesIS) {
185c8272957SStefano Zampini       for (i=0;i<*n_edges;i++) {
186c8272957SStefano Zampini         ierr = ISDestroy(&((*EdgesIS)[i]));CHKERRQ(ierr);
187c8272957SStefano Zampini       }
188c8272957SStefano Zampini       ierr = PetscFree(*EdgesIS);CHKERRQ(ierr);
189c8272957SStefano Zampini     }
190c8272957SStefano Zampini     *n_edges = 0;
191c8272957SStefano Zampini   }
192c8272957SStefano Zampini   if (VerticesIS) {
193c8272957SStefano Zampini     ierr = ISDestroy(VerticesIS);CHKERRQ(ierr);
194c8272957SStefano Zampini   }
195c8272957SStefano Zampini   PetscFunctionReturn(0);
196c8272957SStefano Zampini }
197c8272957SStefano Zampini 
198c8272957SStefano Zampini #undef __FUNCT__
199674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphGetCandidatesIS"
200a873d5d0SStefano Zampini PetscErrorCode PCBDDCGraphGetCandidatesIS(PCBDDCGraph graph, PetscInt *n_faces, IS *FacesIS[], PetscInt *n_edges, IS *EdgesIS[], IS *VerticesIS)
201674ae819SStefano Zampini {
202674ae819SStefano Zampini   IS             *ISForFaces,*ISForEdges,ISForVertices;
203adfc4fb2SStefano Zampini   PetscInt       i,nfc,nec,nvc,*idx;
204674ae819SStefano Zampini   PetscErrorCode ierr;
205674ae819SStefano Zampini 
206674ae819SStefano Zampini   PetscFunctionBegin;
207674ae819SStefano Zampini   /* loop on ccs to evalute number of faces, edges and vertices */
208674ae819SStefano Zampini   nfc = 0;
209674ae819SStefano Zampini   nec = 0;
210674ae819SStefano Zampini   nvc = 0;
211674ae819SStefano Zampini   for (i=0;i<graph->ncc;i++) {
2121e1f2d93SStefano Zampini     PetscInt repdof = graph->queue[graph->cptr[i]];
213674ae819SStefano Zampini     if (graph->cptr[i+1]-graph->cptr[i] > graph->custom_minimal_size) {
214e635b14bSStefano Zampini       if (graph->count[repdof] == 1 && graph->special_dof[repdof] != PCBDDCGRAPH_NEUMANN_MARK) {
215674ae819SStefano Zampini         nfc++;
216674ae819SStefano Zampini       } else { /* note that nec will be zero in 2d */
217674ae819SStefano Zampini         nec++;
218674ae819SStefano Zampini       }
219674ae819SStefano Zampini     } else {
220674ae819SStefano Zampini       nvc += graph->cptr[i+1]-graph->cptr[i];
221674ae819SStefano Zampini     }
222674ae819SStefano Zampini   }
223adfc4fb2SStefano Zampini 
2248e4af1ccSStefano Zampini   /* check if we are in 2D or 3D */
2258e4af1ccSStefano Zampini   if (graph->twodim) { /* we are in a 2D case -> edges are shared by 2 subregions and faces don't exist */
226674ae819SStefano Zampini     nec = nfc;
227674ae819SStefano Zampini     nfc = 0;
228674ae819SStefano Zampini   }
229adfc4fb2SStefano Zampini 
230674ae819SStefano Zampini   /* allocate IS arrays for faces, edges. Vertices need a single index set. */
231cf5a6209SStefano Zampini   if (FacesIS) {
232785e854fSJed Brown     ierr = PetscMalloc1(nfc,&ISForFaces);CHKERRQ(ierr);
233cf5a6209SStefano Zampini   }
234cf5a6209SStefano Zampini   if (EdgesIS) {
235785e854fSJed Brown     ierr = PetscMalloc1(nec,&ISForEdges);CHKERRQ(ierr);
236cf5a6209SStefano Zampini   }
237cf5a6209SStefano Zampini   if (VerticesIS) {
238785e854fSJed Brown     ierr = PetscMalloc1(nvc,&idx);CHKERRQ(ierr);
239cf5a6209SStefano Zampini   }
240cf5a6209SStefano Zampini 
241674ae819SStefano Zampini   /* loop on ccs to compute index sets for faces and edges */
242acc113dbSStefano Zampini   if (!graph->queue_sorted) {
243acc113dbSStefano Zampini     PetscInt *queue_global;
244acc113dbSStefano Zampini 
245acc113dbSStefano Zampini     ierr = PetscMalloc1(graph->cptr[graph->ncc],&queue_global);CHKERRQ(ierr);
246acc113dbSStefano Zampini     ierr = ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_global);CHKERRQ(ierr);
247acc113dbSStefano Zampini     for (i=0;i<graph->ncc;i++) {
248acc113dbSStefano Zampini       ierr = PetscSortIntWithArray(graph->cptr[i+1]-graph->cptr[i],&queue_global[graph->cptr[i]],&graph->queue[graph->cptr[i]]);CHKERRQ(ierr);
249acc113dbSStefano Zampini     }
250acc113dbSStefano Zampini     ierr = PetscFree(queue_global);CHKERRQ(ierr);
251acc113dbSStefano Zampini     graph->queue_sorted = PETSC_TRUE;
252acc113dbSStefano Zampini   }
253674ae819SStefano Zampini   nfc = 0;
254674ae819SStefano Zampini   nec = 0;
255674ae819SStefano Zampini   for (i=0;i<graph->ncc;i++) {
2561e1f2d93SStefano Zampini     PetscInt repdof = graph->queue[graph->cptr[i]];
257674ae819SStefano Zampini     if (graph->cptr[i+1]-graph->cptr[i] > graph->custom_minimal_size) {
258e635b14bSStefano Zampini       if (graph->count[repdof] == 1 && graph->special_dof[repdof] != PCBDDCGRAPH_NEUMANN_MARK) {
2598e4af1ccSStefano Zampini         if (graph->twodim) {
260cf5a6209SStefano Zampini           if (EdgesIS) {
261c8272957SStefano Zampini             ierr = ISCreateGeneral(PETSC_COMM_SELF,graph->cptr[i+1]-graph->cptr[i],&graph->queue[graph->cptr[i]],PETSC_USE_POINTER,&ISForEdges[nec]);CHKERRQ(ierr);
262cf5a6209SStefano Zampini           }
263674ae819SStefano Zampini           nec++;
264674ae819SStefano Zampini         } else {
265cf5a6209SStefano Zampini           if (FacesIS) {
266c8272957SStefano Zampini             ierr = ISCreateGeneral(PETSC_COMM_SELF,graph->cptr[i+1]-graph->cptr[i],&graph->queue[graph->cptr[i]],PETSC_USE_POINTER,&ISForFaces[nfc]);CHKERRQ(ierr);
267cf5a6209SStefano Zampini           }
268674ae819SStefano Zampini           nfc++;
269674ae819SStefano Zampini         }
270674ae819SStefano Zampini       } else {
271cf5a6209SStefano Zampini         if (EdgesIS) {
272c8272957SStefano Zampini           ierr = ISCreateGeneral(PETSC_COMM_SELF,graph->cptr[i+1]-graph->cptr[i],&graph->queue[graph->cptr[i]],PETSC_USE_POINTER,&ISForEdges[nec]);CHKERRQ(ierr);
273cf5a6209SStefano Zampini         }
274674ae819SStefano Zampini         nec++;
275674ae819SStefano Zampini       }
276674ae819SStefano Zampini     }
277674ae819SStefano Zampini   }
278674ae819SStefano Zampini   /* index set for vertices */
279cf5a6209SStefano Zampini   if (VerticesIS) {
280b8ffe317SStefano Zampini     nvc = 0;
281674ae819SStefano Zampini     for (i=0;i<graph->ncc;i++) {
282674ae819SStefano Zampini       if (graph->cptr[i+1]-graph->cptr[i] <= graph->custom_minimal_size) {
283adfc4fb2SStefano Zampini         PetscInt j;
284adfc4fb2SStefano Zampini 
285674ae819SStefano Zampini         for (j=graph->cptr[i];j<graph->cptr[i+1];j++) {
286674ae819SStefano Zampini           idx[nvc]=graph->queue[j];
287674ae819SStefano Zampini           nvc++;
288674ae819SStefano Zampini         }
289674ae819SStefano Zampini       }
290674ae819SStefano Zampini     }
291674ae819SStefano Zampini     /* sort vertex set (by local ordering) */
292674ae819SStefano Zampini     ierr = PetscSortInt(nvc,idx);CHKERRQ(ierr);
293674ae819SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,nvc,idx,PETSC_OWN_POINTER,&ISForVertices);CHKERRQ(ierr);
294674ae819SStefano Zampini   }
295674ae819SStefano Zampini   /* get back info */
296a873d5d0SStefano Zampini   if (n_faces) *n_faces = nfc;
297c8272957SStefano Zampini   if (FacesIS) *FacesIS = ISForFaces;
298a873d5d0SStefano Zampini   if (n_edges) *n_edges = nec;
299c8272957SStefano Zampini   if (EdgesIS) *EdgesIS = ISForEdges;
300c8272957SStefano Zampini   if (VerticesIS) *VerticesIS = ISForVertices;
301674ae819SStefano Zampini   PetscFunctionReturn(0);
302674ae819SStefano Zampini }
303674ae819SStefano Zampini 
304674ae819SStefano Zampini #undef __FUNCT__
305674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphComputeConnectedComponents"
306674ae819SStefano Zampini PetscErrorCode PCBDDCGraphComputeConnectedComponents(PCBDDCGraph graph)
307674ae819SStefano Zampini {
3084a060362SStefano Zampini   PetscBool      adapt_interface_reduced;
309674ae819SStefano Zampini   MPI_Comm       interface_comm;
3104a060362SStefano Zampini   PetscMPIInt    size;
3118e4af1ccSStefano Zampini   PetscInt       i;
3128e4af1ccSStefano Zampini   PetscBool      twodim;
313674ae819SStefano Zampini   PetscErrorCode ierr;
314674ae819SStefano Zampini 
315674ae819SStefano Zampini   PetscFunctionBegin;
316674ae819SStefano Zampini   /* compute connected components locally */
317674ae819SStefano Zampini   ierr = PetscObjectGetComm((PetscObject)(graph->l2gmap),&interface_comm);CHKERRQ(ierr);
318674ae819SStefano Zampini   ierr = PCBDDCGraphComputeConnectedComponentsLocal(graph);CHKERRQ(ierr);
319674ae819SStefano Zampini   /* check consistency of connected components among neighbouring subdomains -> it adapt them in case it is needed */
3204a060362SStefano Zampini   ierr = MPI_Comm_size(interface_comm,&size);CHKERRQ(ierr);
3214a060362SStefano Zampini   adapt_interface_reduced = PETSC_FALSE;
3224a060362SStefano Zampini   if (size > 1) {
3234a060362SStefano Zampini     PetscInt i;
3244a060362SStefano Zampini     PetscBool adapt_interface = PETSC_FALSE;
325674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
326674ae819SStefano Zampini       /* We are not sure that on a given subset of the local interface,
327674ae819SStefano Zampini          with two connected components, the latters be the same among sharing subdomains */
328674ae819SStefano Zampini       if (graph->subset_ncc[i] > 1) {
3294a060362SStefano Zampini         adapt_interface = PETSC_TRUE;
330674ae819SStefano Zampini         break;
331674ae819SStefano Zampini       }
332674ae819SStefano Zampini     }
333b2566f29SBarry Smith     ierr = MPIU_Allreduce(&adapt_interface,&adapt_interface_reduced,1,MPIU_BOOL,MPI_LOR,interface_comm);CHKERRQ(ierr);
3344a060362SStefano Zampini   }
335674ae819SStefano Zampini 
336674ae819SStefano Zampini   if (graph->n_subsets && adapt_interface_reduced) {
337ec1c413dSStefano Zampini     PetscBT     subset_cc_adapt;
338ec1c413dSStefano Zampini     MPI_Request *send_requests,*recv_requests;
339ec1c413dSStefano Zampini     PetscInt    *send_buffer,*recv_buffer;
340ec1c413dSStefano Zampini     PetscInt    sum_requests,start_of_recv,start_of_send;
341ec1c413dSStefano Zampini     PetscInt    *cum_recv_counts;
342ec1c413dSStefano Zampini     PetscInt    *labels;
343b3cdcd63SStefano Zampini     PetscInt    ncc,cum_queue,mss,mns,j,k,s;
3448875e3e6SStefano Zampini     PetscInt    **refine_buffer=NULL,*private_labels = NULL;
3455b1b9aeaSStefano Zampini 
346ec1c413dSStefano Zampini     ierr = PetscMalloc1(graph->nvtxs,&labels);CHKERRQ(ierr);
347ec1c413dSStefano Zampini     ierr = PetscMemzero(labels,graph->nvtxs*sizeof(*labels));CHKERRQ(ierr);
348b3cdcd63SStefano Zampini     for (i=0;i<graph->ncc;i++)
349b3cdcd63SStefano Zampini       for (j=graph->cptr[i];j<graph->cptr[i+1];j++)
350b3cdcd63SStefano Zampini         labels[graph->queue[j]] = i;
351b3cdcd63SStefano Zampini 
352674ae819SStefano Zampini     /* allocate some space */
353854ce69bSBarry Smith     ierr = PetscMalloc1(graph->n_subsets+1,&cum_recv_counts);CHKERRQ(ierr);
354674ae819SStefano Zampini     ierr = PetscMemzero(cum_recv_counts,(graph->n_subsets+1)*sizeof(*cum_recv_counts));CHKERRQ(ierr);
355b3cdcd63SStefano Zampini 
356674ae819SStefano Zampini     /* first count how many neighbours per connected component I will receive from */
357674ae819SStefano Zampini     cum_recv_counts[0] = 0;
358b3cdcd63SStefano 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]];
359ec1c413dSStefano Zampini     ierr = PetscMalloc1(cum_recv_counts[graph->n_subsets],&recv_buffer);CHKERRQ(ierr);
360ec1c413dSStefano Zampini     ierr = PetscMalloc2(cum_recv_counts[graph->n_subsets],&send_requests,cum_recv_counts[graph->n_subsets],&recv_requests);CHKERRQ(ierr);
361674ae819SStefano Zampini     for (i=0;i<cum_recv_counts[graph->n_subsets];i++) {
362674ae819SStefano Zampini       send_requests[i] = MPI_REQUEST_NULL;
363674ae819SStefano Zampini       recv_requests[i] = MPI_REQUEST_NULL;
364674ae819SStefano Zampini     }
365b3cdcd63SStefano Zampini 
366ec1c413dSStefano Zampini     /* exchange with my neighbours the number of my connected components on the subset of interface */
367674ae819SStefano Zampini     sum_requests = 0;
368674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
369ec1c413dSStefano Zampini       PetscMPIInt neigh,tag;
370b3cdcd63SStefano Zampini       PetscInt    count,*neighs;
371ec1c413dSStefano Zampini 
372b3cdcd63SStefano Zampini       count = graph->count[graph->subset_idxs[i][0]];
373b3cdcd63SStefano Zampini       neighs = graph->neighbours_set[graph->subset_idxs[i][0]];
374ec1c413dSStefano Zampini       ierr = PetscMPIIntCast(2*graph->subset_ref_node[i],&tag);CHKERRQ(ierr);
375b3cdcd63SStefano Zampini       for (k=0;k<count;k++) {
376b3cdcd63SStefano Zampini         ierr = PetscMPIIntCast(neighs[k],&neigh);CHKERRQ(ierr);
377674ae819SStefano Zampini         ierr = MPI_Isend(&graph->subset_ncc[i],1,MPIU_INT,neigh,tag,interface_comm,&send_requests[sum_requests]);CHKERRQ(ierr);
378ec1c413dSStefano Zampini         ierr = MPI_Irecv(&recv_buffer[sum_requests],1,MPIU_INT,neigh,tag,interface_comm,&recv_requests[sum_requests]);CHKERRQ(ierr);
379674ae819SStefano Zampini         sum_requests++;
380674ae819SStefano Zampini       }
381674ae819SStefano Zampini     }
382674ae819SStefano Zampini     ierr = MPI_Waitall(sum_requests,recv_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
383674ae819SStefano Zampini     ierr = MPI_Waitall(sum_requests,send_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
384b3cdcd63SStefano Zampini 
385b3cdcd63SStefano Zampini     /* determine the subsets I have to adapt (those having more than 1 cc) */
386ec1c413dSStefano Zampini     ierr = PetscBTCreate(graph->n_subsets,&subset_cc_adapt);CHKERRQ(ierr);
387ec1c413dSStefano Zampini     ierr = PetscBTMemzero(graph->n_subsets,subset_cc_adapt);CHKERRQ(ierr);
388674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
389b3cdcd63SStefano Zampini       if (graph->subset_ncc[i] > 1) {
390b3cdcd63SStefano Zampini         ierr = PetscBTSet(subset_cc_adapt,i);CHKERRQ(ierr);
391b3cdcd63SStefano Zampini         continue;
392b3cdcd63SStefano Zampini       }
393674ae819SStefano Zampini       for (j=cum_recv_counts[i];j<cum_recv_counts[i+1];j++){
394b3cdcd63SStefano Zampini          if (recv_buffer[j] > 1) {
395ec1c413dSStefano Zampini           ierr = PetscBTSet(subset_cc_adapt,i);
396674ae819SStefano Zampini           break;
397674ae819SStefano Zampini         }
398674ae819SStefano Zampini       }
399674ae819SStefano Zampini     }
400ec1c413dSStefano Zampini     ierr = PetscFree(recv_buffer);CHKERRQ(ierr);
401b3cdcd63SStefano Zampini 
402ec1c413dSStefano Zampini     /* determine send/recv buffers sizes */
403ec1c413dSStefano Zampini     j = 0;
404b3cdcd63SStefano Zampini     mss = 0;
405674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
406ec1c413dSStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i)) {
407b3cdcd63SStefano Zampini         j += graph->subset_size[i];
408b3cdcd63SStefano Zampini         mss = PetscMax(graph->subset_size[i],mss);
409674ae819SStefano Zampini       }
410674ae819SStefano Zampini     }
411ec1c413dSStefano Zampini     k = 0;
412b3cdcd63SStefano Zampini     mns = 0;
413674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
414ec1c413dSStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i)) {
415b3cdcd63SStefano Zampini         k += (cum_recv_counts[i+1]-cum_recv_counts[i])*graph->subset_size[i];
416b3cdcd63SStefano Zampini         mns = PetscMax(cum_recv_counts[i+1]-cum_recv_counts[i],mns);
417674ae819SStefano Zampini       }
418674ae819SStefano Zampini     }
419ec1c413dSStefano Zampini     ierr = PetscMalloc2(j,&send_buffer,k,&recv_buffer);CHKERRQ(ierr);
420ec1c413dSStefano Zampini 
421b3cdcd63SStefano Zampini     /* fill send buffer (order matters: subset_idxs ordered by global ordering) */
422ec1c413dSStefano Zampini     j = 0;
423b3cdcd63SStefano Zampini     for (i=0;i<graph->n_subsets;i++)
424b3cdcd63SStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i))
425b3cdcd63SStefano Zampini         for (k=0;k<graph->subset_size[i];k++)
426b3cdcd63SStefano Zampini           send_buffer[j++] = labels[graph->subset_idxs[i][k]];
427ec1c413dSStefano Zampini 
428674ae819SStefano Zampini     /* now exchange the data */
429674ae819SStefano Zampini     start_of_recv = 0;
430674ae819SStefano Zampini     start_of_send = 0;
431674ae819SStefano Zampini     sum_requests = 0;
432674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
433ec1c413dSStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i)) {
434ec1c413dSStefano Zampini         PetscMPIInt neigh,tag;
435b3cdcd63SStefano Zampini         PetscInt    size_of_send = graph->subset_size[i];
436ec1c413dSStefano Zampini 
437b3cdcd63SStefano Zampini         j = graph->subset_idxs[i][0];
438ec1c413dSStefano Zampini         ierr = PetscMPIIntCast(2*graph->subset_ref_node[i]+1,&tag);CHKERRQ(ierr);
439674ae819SStefano Zampini         for (k=0;k<graph->count[j];k++) {
440674ae819SStefano Zampini           ierr = PetscMPIIntCast(graph->neighbours_set[j][k],&neigh);CHKERRQ(ierr);
441ec1c413dSStefano Zampini           ierr = MPI_Isend(&send_buffer[start_of_send],size_of_send,MPIU_INT,neigh,tag,interface_comm,&send_requests[sum_requests]);CHKERRQ(ierr);
442ec1c413dSStefano Zampini           ierr = MPI_Irecv(&recv_buffer[start_of_recv],size_of_send,MPIU_INT,neigh,tag,interface_comm,&recv_requests[sum_requests]);CHKERRQ(ierr);
443ec1c413dSStefano Zampini           start_of_recv += size_of_send;
444674ae819SStefano Zampini           sum_requests++;
445674ae819SStefano Zampini         }
446674ae819SStefano Zampini         start_of_send += size_of_send;
447674ae819SStefano Zampini       }
448674ae819SStefano Zampini     }
449674ae819SStefano Zampini     ierr = MPI_Waitall(sum_requests,recv_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
450d52c4852SStefano Zampini 
451b3cdcd63SStefano Zampini     /* refine connected components */
452674ae819SStefano Zampini     start_of_recv = 0;
453b3cdcd63SStefano Zampini     /* allocate some temporary space */
454b3cdcd63SStefano Zampini     if (mss) {
455b3cdcd63SStefano Zampini       ierr = PetscMalloc1(mss,&refine_buffer);CHKERRQ(ierr);
456b3cdcd63SStefano Zampini       ierr = PetscMalloc2(mss*(mns+1),&refine_buffer[0],mss,&private_labels);CHKERRQ(ierr);
457b3cdcd63SStefano Zampini     }
458b3cdcd63SStefano Zampini     ncc = 0;
459b3cdcd63SStefano Zampini     cum_queue = 0;
460b3cdcd63SStefano Zampini     graph->cptr[0] = 0;
461674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
462ec1c413dSStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i)) {
463b3cdcd63SStefano Zampini         PetscInt subset_counter = 0;
464b3cdcd63SStefano Zampini         PetscInt sharingprocs = cum_recv_counts[i+1]-cum_recv_counts[i]+1; /* count myself */
465b3cdcd63SStefano Zampini         PetscInt buffer_size = graph->subset_size[i];
466ec1c413dSStefano Zampini 
467b3cdcd63SStefano Zampini         /* compute pointers */
468b3cdcd63SStefano Zampini         for (j=1;j<buffer_size;j++) refine_buffer[j] = refine_buffer[j-1] + sharingprocs;
469b3cdcd63SStefano Zampini         /* analyze contributions from subdomains that share the i-th subset
470b3cdcd63SStefano Zampini            The stricture of refine_buffer is suitable to find intersections of ccs among sharingprocs.
471b3cdcd63SStefano 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)
472b3cdcd63SStefano Zampini            sharing procs connected components:
473ec1c413dSStefano Zampini              neigh 0: [0 1 4], [2 3], labels [4,7]  (2 connected components)
474ec1c413dSStefano Zampini              neigh 1: [0 1], [2 3 4], labels [3 2]  (2 connected components)
475ec1c413dSStefano Zampini              neigh 2: [0 4], [1], [2 3], labels [1 5 6] (3 connected components)
476b3cdcd63SStefano Zampini            refine_buffer will be filled as:
477ec1c413dSStefano Zampini              [ 4, 3, 1;
478ec1c413dSStefano Zampini                4, 2, 1;
479ec1c413dSStefano Zampini                7, 2, 6;
480ec1c413dSStefano Zampini                4, 3, 5;
481ec1c413dSStefano Zampini                7, 2, 6; ];
482b3cdcd63SStefano Zampini            The connected components in local ordering are [0], [1], [2 3], [4] */
483b3cdcd63SStefano Zampini         /* fill temp_buffer */
484b3cdcd63SStefano Zampini         for (k=0;k<buffer_size;k++) refine_buffer[k][0] = labels[graph->subset_idxs[i][k]];
485b3cdcd63SStefano Zampini         for (j=0;j<sharingprocs-1;j++) {
486b3cdcd63SStefano Zampini           for (k=0;k<buffer_size;k++) refine_buffer[k][j+1] = recv_buffer[start_of_recv+k];
487b3cdcd63SStefano Zampini           start_of_recv += buffer_size;
488674ae819SStefano Zampini         }
489b3cdcd63SStefano Zampini         ierr = PetscMemzero(private_labels,buffer_size*sizeof(PetscInt));CHKERRQ(ierr);
490b3cdcd63SStefano Zampini         for (j=0;j<buffer_size;j++) {
491b3cdcd63SStefano Zampini           if (!private_labels[j]) { /* found a new cc  */
492ec1c413dSStefano Zampini             PetscBool same_set;
493ec1c413dSStefano Zampini 
494b3cdcd63SStefano Zampini             graph->cptr[ncc] = cum_queue;
495b3cdcd63SStefano Zampini             ncc++;
496b3cdcd63SStefano Zampini             subset_counter++;
497b3cdcd63SStefano Zampini             private_labels[j] = subset_counter;
498b3cdcd63SStefano Zampini             graph->queue[cum_queue++] = graph->subset_idxs[i][j];
499b3cdcd63SStefano Zampini             for (k=j+1;k<buffer_size;k++) { /* check for other nodes in new cc */
500674ae819SStefano Zampini               same_set = PETSC_TRUE;
501b3cdcd63SStefano Zampini               for (s=0;s<sharingprocs;s++) {
502b3cdcd63SStefano Zampini                 if (refine_buffer[j][s] != refine_buffer[k][s]) {
503674ae819SStefano Zampini                   same_set = PETSC_FALSE;
504674ae819SStefano Zampini                   break;
505674ae819SStefano Zampini                 }
506674ae819SStefano Zampini               }
507674ae819SStefano Zampini               if (same_set) {
508b3cdcd63SStefano Zampini                 private_labels[k] = subset_counter;
509b3cdcd63SStefano Zampini                 graph->queue[cum_queue++] = graph->subset_idxs[i][k];
510674ae819SStefano Zampini               }
511674ae819SStefano Zampini             }
512674ae819SStefano Zampini           }
513674ae819SStefano Zampini         }
514b3cdcd63SStefano Zampini         graph->cptr[ncc] = cum_queue;
515b3cdcd63SStefano Zampini         graph->subset_ncc[i] = subset_counter;
516b3cdcd63SStefano Zampini         graph->queue_sorted = PETSC_FALSE;
517b3cdcd63SStefano Zampini       } else { /* this subset does not need to be adapted */
518b3cdcd63SStefano Zampini         ierr = PetscMemcpy(graph->queue+cum_queue,graph->subset_idxs[i],graph->subset_size[i]*sizeof(PetscInt));CHKERRQ(ierr);
519b3cdcd63SStefano Zampini         ncc++;
520b3cdcd63SStefano Zampini         cum_queue += graph->subset_size[i];
521b3cdcd63SStefano Zampini         graph->cptr[ncc] = cum_queue;
522674ae819SStefano Zampini       }
523674ae819SStefano Zampini     }
524b3cdcd63SStefano Zampini     graph->cptr[ncc] = cum_queue;
525b3cdcd63SStefano Zampini     graph->ncc = ncc;
526b3cdcd63SStefano Zampini     if (mss) {
527b3cdcd63SStefano Zampini       ierr = PetscFree2(refine_buffer[0],private_labels);CHKERRQ(ierr);
528b3cdcd63SStefano Zampini       ierr = PetscFree(refine_buffer);CHKERRQ(ierr);
529b3cdcd63SStefano Zampini     }
530b3cdcd63SStefano Zampini     ierr = PetscFree(labels);CHKERRQ(ierr);
531b3cdcd63SStefano Zampini     ierr = MPI_Waitall(sum_requests,send_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
532b3cdcd63SStefano Zampini     ierr = PetscFree2(send_requests,recv_requests);CHKERRQ(ierr);
53345b2a5aaSStefano Zampini     ierr = PetscFree2(send_buffer,recv_buffer);CHKERRQ(ierr);
534674ae819SStefano Zampini     ierr = PetscFree(cum_recv_counts);CHKERRQ(ierr);
535ec1c413dSStefano Zampini     ierr = PetscBTDestroy(&subset_cc_adapt);CHKERRQ(ierr);
536674ae819SStefano Zampini   }
5378e4af1ccSStefano Zampini 
5388e4af1ccSStefano Zampini   /* Determine if we are in 2D or 3D */
5398e4af1ccSStefano Zampini   twodim  = PETSC_TRUE;
5408e4af1ccSStefano Zampini   for (i=0;i<graph->ncc;i++) {
5418e4af1ccSStefano Zampini     PetscInt repdof = graph->queue[graph->cptr[i]];
5428e4af1ccSStefano Zampini     if (graph->cptr[i+1]-graph->cptr[i] > graph->custom_minimal_size) {
5438e4af1ccSStefano Zampini       if (graph->count[repdof] > 1 || graph->special_dof[repdof] == PCBDDCGRAPH_NEUMANN_MARK) {
5448e4af1ccSStefano Zampini         twodim = PETSC_FALSE;
5458e4af1ccSStefano Zampini         break;
5468e4af1ccSStefano Zampini       }
5478e4af1ccSStefano Zampini     }
5488e4af1ccSStefano Zampini   }
549b2566f29SBarry Smith   ierr = MPIU_Allreduce(&twodim,&graph->twodim,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)graph->l2gmap));CHKERRQ(ierr);
550674ae819SStefano Zampini   PetscFunctionReturn(0);
551674ae819SStefano Zampini }
552674ae819SStefano Zampini 
553b3cdcd63SStefano Zampini 
554b3cdcd63SStefano Zampini #undef __FUNCT__
555b3cdcd63SStefano Zampini #define __FUNCT__ "PCBDDCGraphComputeCC_Private"
556b3cdcd63SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCGraphComputeCC_Private(PCBDDCGraph graph,PetscInt pid,PetscInt* queue_tip,PetscInt n_prev,PetscInt* n_added)
557b3cdcd63SStefano Zampini {
558b3cdcd63SStefano Zampini   PetscInt       i,j,n;
559b3cdcd63SStefano Zampini   PetscInt       *xadj = graph->xadj,*adjncy = graph->adjncy;
560b3cdcd63SStefano Zampini   PetscBT        touched = graph->touched;
561b3cdcd63SStefano Zampini   PetscBool      havecsr = (PetscBool)(xadj && adjncy);
562b3cdcd63SStefano Zampini   PetscBool      havesubs = (PetscBool)(!!graph->n_local_subs);
563b3cdcd63SStefano Zampini   PetscErrorCode ierr;
564b3cdcd63SStefano Zampini 
565b3cdcd63SStefano Zampini   PetscFunctionBegin;
566b3cdcd63SStefano Zampini   n = 0;
567b3cdcd63SStefano Zampini   if (havecsr && !havesubs) {
568b3cdcd63SStefano Zampini     for (i=-n_prev;i<0;i++) {
569b3cdcd63SStefano Zampini       PetscInt start_dof = queue_tip[i];
570b3cdcd63SStefano Zampini       for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
571b3cdcd63SStefano Zampini         PetscInt dof = adjncy[j];
572b3cdcd63SStefano Zampini         if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid) {
573b3cdcd63SStefano Zampini           ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
574b3cdcd63SStefano Zampini           queue_tip[n] = dof;
575b3cdcd63SStefano Zampini           n++;
576b3cdcd63SStefano Zampini         }
577b3cdcd63SStefano Zampini       }
578b3cdcd63SStefano Zampini     }
579b3cdcd63SStefano Zampini   } else if (havecsr && havesubs) {
580b3cdcd63SStefano Zampini     PetscInt sid = graph->local_subs[queue_tip[-n_prev]];
581b3cdcd63SStefano Zampini     for (i=-n_prev;i<0;i++) {
582b3cdcd63SStefano Zampini       PetscInt start_dof = queue_tip[i];
583b3cdcd63SStefano Zampini       for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
584b3cdcd63SStefano Zampini         PetscInt dof = adjncy[j];
585b3cdcd63SStefano Zampini         if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid && graph->local_subs[dof] == sid) {
586b3cdcd63SStefano Zampini           ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
587b3cdcd63SStefano Zampini           queue_tip[n] = dof;
588b3cdcd63SStefano Zampini           n++;
589b3cdcd63SStefano Zampini         }
590b3cdcd63SStefano Zampini       }
591b3cdcd63SStefano Zampini     }
592b3cdcd63SStefano Zampini   } else { /* sub info only */
593b3cdcd63SStefano Zampini     PetscInt sid = graph->local_subs[queue_tip[-n_prev]];
594b3cdcd63SStefano Zampini     for (j=0;j<graph->subset_size[pid-1];j++) { /* pid \in [1,graph->n_subsets] */
595b3cdcd63SStefano Zampini       PetscInt dof = graph->subset_idxs[pid-1][j];
596b3cdcd63SStefano Zampini       if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid && graph->local_subs[dof] == sid) {
597b3cdcd63SStefano Zampini         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
598b3cdcd63SStefano Zampini         queue_tip[n] = dof;
599b3cdcd63SStefano Zampini         n++;
600b3cdcd63SStefano Zampini       }
601b3cdcd63SStefano Zampini     }
602b3cdcd63SStefano Zampini   }
603b3cdcd63SStefano Zampini   *n_added = n;
604b3cdcd63SStefano Zampini   PetscFunctionReturn(0);
605b3cdcd63SStefano Zampini }
606b3cdcd63SStefano Zampini 
607674ae819SStefano Zampini #undef __FUNCT__
608674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphComputeConnectedComponentsLocal"
609674ae819SStefano Zampini PetscErrorCode PCBDDCGraphComputeConnectedComponentsLocal(PCBDDCGraph graph)
610674ae819SStefano Zampini {
611b3cdcd63SStefano Zampini   PetscInt       ncc,cum_queue,n;
612b3cdcd63SStefano Zampini   PetscMPIInt    commsize;
613674ae819SStefano Zampini   PetscErrorCode ierr;
614674ae819SStefano Zampini 
615674ae819SStefano Zampini   PetscFunctionBegin;
616b3cdcd63SStefano Zampini   if (!graph->setupcalled) SETERRQ(PetscObjectComm((PetscObject)graph->l2gmap),PETSC_ERR_ORDER,"PCBDDCGraphSetUp should be called first");
617b3cdcd63SStefano Zampini   /* quiet return if there isn't any local info */
6184f1b2e48SStefano Zampini   if ((!graph->xadj || !graph->adjncy) && !graph->n_local_subs) {
619674ae819SStefano Zampini     PetscFunctionReturn(0);
620674ae819SStefano Zampini   }
621674ae819SStefano Zampini 
622674ae819SStefano Zampini   /* reset any previous search of connected components */
623df48933bSStefano Zampini   ierr = PetscBTMemzero(graph->nvtxs,graph->touched);CHKERRQ(ierr);
624b3cdcd63SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)graph->l2gmap),&commsize);CHKERRQ(ierr);
6254b2aedd3SStefano Zampini   if (commsize > graph->commsizelimit) {
626b3cdcd63SStefano Zampini     PetscInt i;
627674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
6280cf82fd6SStefano Zampini       if (graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK || !graph->count[i]) {
629df48933bSStefano Zampini         ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
630674ae819SStefano Zampini       }
6314a060362SStefano Zampini     }
632b3cdcd63SStefano Zampini   }
633674ae819SStefano Zampini 
634674ae819SStefano Zampini   /* begin search for connected components */
635674ae819SStefano Zampini   cum_queue = 0;
636674ae819SStefano Zampini   ncc = 0;
637674ae819SStefano Zampini   for (n=0;n<graph->n_subsets;n++) {
638b3cdcd63SStefano Zampini     PetscInt pid = n+1;  /* partition labeled by 0 is discarded */
639b3cdcd63SStefano Zampini     PetscInt found = 0,prev = 0,first = 0,ncc_pid = 0;
640b3cdcd63SStefano Zampini     while (found != graph->subset_size[n]) {
641b3cdcd63SStefano Zampini       PetscInt added = 0;
642b3cdcd63SStefano Zampini       if (!prev) { /* search for new starting dof */
643b3cdcd63SStefano Zampini         while (PetscBTLookup(graph->touched,graph->subset_idxs[n][first])) first++;
644b3cdcd63SStefano Zampini         ierr = PetscBTSet(graph->touched,graph->subset_idxs[n][first]);CHKERRQ(ierr);
645b3cdcd63SStefano Zampini         graph->queue[cum_queue] = graph->subset_idxs[n][first];
646674ae819SStefano Zampini         graph->cptr[ncc] = cum_queue;
647b3cdcd63SStefano Zampini         prev = 1;
648b3cdcd63SStefano Zampini         cum_queue++;
649b3cdcd63SStefano Zampini         found++;
650674ae819SStefano Zampini         ncc_pid++;
651b3cdcd63SStefano Zampini         ncc++;
652674ae819SStefano Zampini       }
653b3cdcd63SStefano Zampini       ierr = PCBDDCGraphComputeCC_Private(graph,pid,graph->queue + cum_queue,prev,&added);CHKERRQ(ierr);
654b3cdcd63SStefano Zampini       if (!added) {
655674ae819SStefano Zampini         graph->subset_ncc[n] = ncc_pid;
656b3cdcd63SStefano Zampini         graph->cptr[ncc] = cum_queue;
657b3cdcd63SStefano Zampini       }
658b3cdcd63SStefano Zampini       prev = added;
659b3cdcd63SStefano Zampini       found += added;
660b3cdcd63SStefano Zampini       cum_queue += added;
661b3cdcd63SStefano Zampini       if (added && found == graph->subset_size[n]) {
662b3cdcd63SStefano Zampini         graph->subset_ncc[n] = ncc_pid;
663b3cdcd63SStefano Zampini         graph->cptr[ncc] = cum_queue;
664b3cdcd63SStefano Zampini       }
665b3cdcd63SStefano Zampini     }
666674ae819SStefano Zampini   }
667674ae819SStefano Zampini   graph->ncc = ncc;
668acc113dbSStefano Zampini   graph->queue_sorted = PETSC_FALSE;
669674ae819SStefano Zampini   PetscFunctionReturn(0);
670674ae819SStefano Zampini }
671674ae819SStefano Zampini 
672674ae819SStefano Zampini #undef __FUNCT__
673674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphSetUp"
674674ae819SStefano Zampini PetscErrorCode PCBDDCGraphSetUp(PCBDDCGraph graph, PetscInt custom_minimal_size, IS neumann_is, IS dirichlet_is, PetscInt n_ISForDofs, IS ISForDofs[], IS custom_primal_vertices)
675674ae819SStefano Zampini {
676674ae819SStefano Zampini   VecScatter     scatter_ctx;
677674ae819SStefano Zampini   Vec            local_vec,local_vec2,global_vec;
678dc456d91SStefano Zampini   IS             to,from,subset,subset_n;
679674ae819SStefano Zampini   MPI_Comm       comm;
680674ae819SStefano Zampini   PetscScalar    *array,*array2;
681674ae819SStefano Zampini   const PetscInt *is_indices;
682dc456d91SStefano Zampini   PetscInt       n_neigh,*neigh,*n_shared,**shared,*queue_global;
683674ae819SStefano Zampini   PetscInt       i,j,k,s,total_counts,nodes_touched,is_size;
684b3cdcd63SStefano Zampini   PetscMPIInt    commsize;
68551b0f6cfSStefano Zampini   PetscBool      same_set,mirrors_found;
6867babac9bSStefano Zampini   PetscErrorCode ierr;
687674ae819SStefano Zampini 
688674ae819SStefano Zampini   PetscFunctionBegin;
689a07ea27aSStefano Zampini   graph->has_dirichlet = PETSC_FALSE;
690a07ea27aSStefano Zampini   if (dirichlet_is) {
691a07ea27aSStefano Zampini     PetscCheckSameComm(graph->l2gmap,1,dirichlet_is,4);
692a07ea27aSStefano Zampini     graph->has_dirichlet = PETSC_TRUE;
693a07ea27aSStefano Zampini   }
694674ae819SStefano Zampini   ierr = PetscObjectGetComm((PetscObject)(graph->l2gmap),&comm);CHKERRQ(ierr);
695b3cdcd63SStefano Zampini   ierr = MPI_Comm_size(comm,&commsize);CHKERRQ(ierr);
696b3cdcd63SStefano Zampini 
697674ae819SStefano Zampini   /* custom_minimal_size */
69814f95afaSStefano Zampini   graph->custom_minimal_size = custom_minimal_size;
699674ae819SStefano Zampini   /* get info l2gmap and allocate work vectors  */
700674ae819SStefano Zampini   ierr = ISLocalToGlobalMappingGetInfo(graph->l2gmap,&n_neigh,&neigh,&n_shared,&shared);CHKERRQ(ierr);
7012635a1d4SStefano Zampini   /* check if we have any local periodic nodes (periodic BCs) */
7022635a1d4SStefano Zampini   mirrors_found = PETSC_FALSE;
7032635a1d4SStefano Zampini   if (graph->nvtxs && n_neigh) {
7042635a1d4SStefano Zampini     for (i=0; i<n_shared[0]; i++) graph->count[shared[0][i]] += 1;
7052635a1d4SStefano Zampini     for (i=0; i<n_shared[0]; i++) {
7062635a1d4SStefano Zampini       if (graph->count[shared[0][i]] > 1) {
7072635a1d4SStefano Zampini         mirrors_found = PETSC_TRUE;
7082635a1d4SStefano Zampini         break;
7092635a1d4SStefano Zampini       }
7102635a1d4SStefano Zampini     }
7112635a1d4SStefano Zampini   }
7122635a1d4SStefano Zampini   /* create some workspace objects */
7132635a1d4SStefano Zampini   local_vec = NULL;
7142635a1d4SStefano Zampini   local_vec2 = NULL;
7152635a1d4SStefano Zampini   global_vec = NULL;
7162635a1d4SStefano Zampini   to = NULL;
7172635a1d4SStefano Zampini   from = NULL;
7182635a1d4SStefano Zampini   scatter_ctx = NULL;
7192635a1d4SStefano Zampini   if (n_ISForDofs || dirichlet_is || neumann_is || custom_primal_vertices) {
720674ae819SStefano Zampini     ierr = VecCreate(PETSC_COMM_SELF,&local_vec);CHKERRQ(ierr);
721674ae819SStefano Zampini     ierr = VecSetSizes(local_vec,PETSC_DECIDE,graph->nvtxs);CHKERRQ(ierr);
722674ae819SStefano Zampini     ierr = VecSetType(local_vec,VECSTANDARD);CHKERRQ(ierr);
723674ae819SStefano Zampini     ierr = VecDuplicate(local_vec,&local_vec2);CHKERRQ(ierr);
724674ae819SStefano Zampini     ierr = VecCreate(comm,&global_vec);CHKERRQ(ierr);
7257fb0e2dbSStefano Zampini     ierr = VecSetSizes(global_vec,PETSC_DECIDE,graph->nvtxs_global);CHKERRQ(ierr);
726674ae819SStefano Zampini     ierr = VecSetType(global_vec,VECSTANDARD);CHKERRQ(ierr);
727674ae819SStefano Zampini     ierr = ISCreateStride(PETSC_COMM_SELF,graph->nvtxs,0,1,&to);CHKERRQ(ierr);
728674ae819SStefano Zampini     ierr = ISLocalToGlobalMappingApplyIS(graph->l2gmap,to,&from);CHKERRQ(ierr);
729674ae819SStefano Zampini     ierr = VecScatterCreate(global_vec,from,local_vec,to,&scatter_ctx);CHKERRQ(ierr);
7302635a1d4SStefano Zampini   } else if (mirrors_found) {
7312635a1d4SStefano Zampini     ierr = ISCreateStride(PETSC_COMM_SELF,graph->nvtxs,0,1,&to);CHKERRQ(ierr);
7322635a1d4SStefano Zampini     ierr = ISLocalToGlobalMappingApplyIS(graph->l2gmap,to,&from);CHKERRQ(ierr);
73349eeff8cSStefano Zampini   }
73451b0f6cfSStefano Zampini   /* compute local mirrors (if any) */
73551b0f6cfSStefano Zampini   if (mirrors_found) {
73651b0f6cfSStefano Zampini     PetscInt *local_indices,*global_indices;
73751b0f6cfSStefano Zampini     /* get arrays of local and global indices */
738785e854fSJed Brown     ierr = PetscMalloc1(graph->nvtxs,&local_indices);CHKERRQ(ierr);
73951b0f6cfSStefano Zampini     ierr = ISGetIndices(to,(const PetscInt**)&is_indices);CHKERRQ(ierr);
74051b0f6cfSStefano Zampini     ierr = PetscMemcpy(local_indices,is_indices,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
74151b0f6cfSStefano Zampini     ierr = ISRestoreIndices(to,(const PetscInt**)&is_indices);CHKERRQ(ierr);
742785e854fSJed Brown     ierr = PetscMalloc1(graph->nvtxs,&global_indices);CHKERRQ(ierr);
74351b0f6cfSStefano Zampini     ierr = ISGetIndices(from,(const PetscInt**)&is_indices);CHKERRQ(ierr);
74451b0f6cfSStefano Zampini     ierr = PetscMemcpy(global_indices,is_indices,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
74551b0f6cfSStefano Zampini     ierr = ISRestoreIndices(from,(const PetscInt**)&is_indices);CHKERRQ(ierr);
74651b0f6cfSStefano Zampini     /* allocate space for mirrors */
7478e5aaad5SJed Brown     ierr = PetscMalloc2(graph->nvtxs,&graph->mirrors,graph->nvtxs,&graph->mirrors_set);CHKERRQ(ierr);
74851b0f6cfSStefano Zampini     ierr = PetscMemzero(graph->mirrors,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
74951b0f6cfSStefano Zampini     graph->mirrors_set[0] = 0;
75051b0f6cfSStefano Zampini 
75151b0f6cfSStefano Zampini     k=0;
75251b0f6cfSStefano Zampini     for (i=0;i<n_shared[0];i++) {
75351b0f6cfSStefano Zampini       j=shared[0][i];
75451b0f6cfSStefano Zampini       if (graph->count[j] > 1) {
75551b0f6cfSStefano Zampini         graph->mirrors[j]++;
75651b0f6cfSStefano Zampini         k++;
75751b0f6cfSStefano Zampini       }
75851b0f6cfSStefano Zampini     }
75951b0f6cfSStefano Zampini     /* allocate space for set of mirrors */
760785e854fSJed Brown     ierr = PetscMalloc1(k,&graph->mirrors_set[0]);CHKERRQ(ierr);
76151b0f6cfSStefano Zampini     for (i=1;i<graph->nvtxs;i++)
76251b0f6cfSStefano Zampini       graph->mirrors_set[i]=graph->mirrors_set[i-1]+graph->mirrors[i-1];
76351b0f6cfSStefano Zampini 
76451b0f6cfSStefano Zampini     /* fill arrays */
76551b0f6cfSStefano Zampini     ierr = PetscMemzero(graph->mirrors,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
76651b0f6cfSStefano Zampini     for (j=0;j<n_shared[0];j++) {
76751b0f6cfSStefano Zampini       i=shared[0][j];
76851b0f6cfSStefano Zampini       if (graph->count[i] > 1)
76951b0f6cfSStefano Zampini         graph->mirrors_set[i][graph->mirrors[i]++]=global_indices[i];
77051b0f6cfSStefano Zampini     }
77151b0f6cfSStefano Zampini     ierr = PetscSortIntWithArray(graph->nvtxs,global_indices,local_indices);CHKERRQ(ierr);
77251b0f6cfSStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
77351b0f6cfSStefano Zampini       if (graph->mirrors[i] > 0) {
77451b0f6cfSStefano Zampini         ierr = PetscFindInt(graph->mirrors_set[i][0],graph->nvtxs,global_indices,&k);CHKERRQ(ierr);
77551b0f6cfSStefano Zampini         j = global_indices[k];
77651b0f6cfSStefano Zampini         while ( k > 0 && global_indices[k-1] == j) k--;
77751b0f6cfSStefano Zampini         for (j=0;j<graph->mirrors[i];j++) {
77851b0f6cfSStefano Zampini           graph->mirrors_set[i][j]=local_indices[k+j];
77951b0f6cfSStefano Zampini         }
78051b0f6cfSStefano Zampini         ierr = PetscSortInt(graph->mirrors[i],graph->mirrors_set[i]);CHKERRQ(ierr);
78151b0f6cfSStefano Zampini       }
78251b0f6cfSStefano Zampini     }
78351b0f6cfSStefano Zampini     ierr = PetscFree(local_indices);CHKERRQ(ierr);
78451b0f6cfSStefano Zampini     ierr = PetscFree(global_indices);CHKERRQ(ierr);
78551b0f6cfSStefano Zampini   }
78651b0f6cfSStefano Zampini   ierr = PetscMemzero(graph->count,graph->nvtxs*sizeof(*graph->count));CHKERRQ(ierr);
787674ae819SStefano Zampini   ierr = ISDestroy(&to);CHKERRQ(ierr);
788674ae819SStefano Zampini   ierr = ISDestroy(&from);CHKERRQ(ierr);
78951b0f6cfSStefano Zampini 
790674ae819SStefano Zampini   /* Count total number of neigh per node */
791674ae819SStefano Zampini   k = 0;
792674ae819SStefano Zampini   for (i=1;i<n_neigh;i++) {
793674ae819SStefano Zampini     k += n_shared[i];
794674ae819SStefano Zampini     for (j=0;j<n_shared[i];j++) {
795674ae819SStefano Zampini       graph->count[shared[i][j]] += 1;
796674ae819SStefano Zampini     }
797674ae819SStefano Zampini   }
798674ae819SStefano Zampini   /* Allocate space for storing the set of neighbours for each node */
799674ae819SStefano Zampini   if (graph->nvtxs) {
800785e854fSJed Brown     ierr = PetscMalloc1(k,&graph->neighbours_set[0]);CHKERRQ(ierr);
801674ae819SStefano Zampini   }
802674ae819SStefano Zampini   for (i=1;i<graph->nvtxs;i++) { /* dont count myself */
803674ae819SStefano Zampini     graph->neighbours_set[i]=graph->neighbours_set[i-1]+graph->count[i-1];
804674ae819SStefano Zampini   }
805674ae819SStefano Zampini   /* Get information for sharing subdomains */
806674ae819SStefano Zampini   ierr = PetscMemzero(graph->count,graph->nvtxs*sizeof(*graph->count));CHKERRQ(ierr);
807674ae819SStefano Zampini   for (i=1;i<n_neigh;i++) { /* dont count myself */
808674ae819SStefano Zampini     s = n_shared[i];
809674ae819SStefano Zampini     for (j=0;j<s;j++) {
810674ae819SStefano Zampini       k = shared[i][j];
811674ae819SStefano Zampini       graph->neighbours_set[k][graph->count[k]] = neigh[i];
812674ae819SStefano Zampini       graph->count[k] += 1;
813674ae819SStefano Zampini     }
814674ae819SStefano Zampini   }
815674ae819SStefano Zampini   /* sort set of sharing subdomains */
816674ae819SStefano Zampini   for (i=0;i<graph->nvtxs;i++) {
81793fb5fd3SStefano Zampini     ierr = PetscSortRemoveDupsInt(&graph->count[i],graph->neighbours_set[i]);CHKERRQ(ierr);
818674ae819SStefano Zampini   }
8197fb0e2dbSStefano Zampini   /* free memory allocated by ISLocalToGlobalMappingGetInfo */
8207fb0e2dbSStefano Zampini   ierr = ISLocalToGlobalMappingRestoreInfo(graph->l2gmap,&n_neigh,&neigh,&n_shared,&shared);CHKERRQ(ierr);
8217fb0e2dbSStefano Zampini 
82267c9da69SStefano Zampini   /*
82367c9da69SStefano Zampini      Get info for dofs splitting
8245777c63cSStefano Zampini      User can specify just a subset; an additional field is considered as a complementary field
82567c9da69SStefano Zampini   */
826b3cdcd63SStefano Zampini   for (i=0;i<graph->nvtxs;i++) graph->which_dof[i] = n_ISForDofs; /* by default a dof belongs to the complement set */
8275777c63cSStefano Zampini   if (n_ISForDofs) {
8285777c63cSStefano Zampini     ierr = VecSet(local_vec,-1.0);CHKERRQ(ierr);
8295777c63cSStefano Zampini   }
830674ae819SStefano Zampini   for (i=0;i<n_ISForDofs;i++) {
8315777c63cSStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
83263602bcaSStefano Zampini     ierr = ISGetLocalSize(ISForDofs[i],&is_size);CHKERRQ(ierr);
833674ae819SStefano Zampini     ierr = ISGetIndices(ISForDofs[i],(const PetscInt**)&is_indices);CHKERRQ(ierr);
834674ae819SStefano Zampini     for (j=0;j<is_size;j++) {
835d50ed60aSStefano Zampini       if (is_indices[j] > -1 && is_indices[j] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
836674ae819SStefano Zampini         graph->which_dof[is_indices[j]] = i;
8375777c63cSStefano Zampini         array[is_indices[j]] = 1.*i;
838674ae819SStefano Zampini       }
83967c9da69SStefano Zampini     }
840674ae819SStefano Zampini     ierr = ISRestoreIndices(ISForDofs[i],(const PetscInt**)&is_indices);CHKERRQ(ierr);
8415777c63cSStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
8425777c63cSStefano Zampini   }
8435777c63cSStefano Zampini   /* Check consistency among neighbours */
8445777c63cSStefano Zampini   if (n_ISForDofs) {
8455777c63cSStefano Zampini     ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8465777c63cSStefano Zampini     ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8475777c63cSStefano Zampini     ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8485777c63cSStefano Zampini     ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8495777c63cSStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
8505777c63cSStefano Zampini     ierr = VecGetArray(local_vec2,&array2);CHKERRQ(ierr);
8515777c63cSStefano Zampini     for (i=0;i<graph->nvtxs;i++){
8525777c63cSStefano Zampini       PetscInt field1,field2;
8535777c63cSStefano Zampini 
8545777c63cSStefano Zampini       field1 = (PetscInt)PetscRealPart(array[i]);
8555777c63cSStefano Zampini       field2 = (PetscInt)PetscRealPart(array2[i]);
8566c4ed002SBarry 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);
8575777c63cSStefano Zampini     }
8585777c63cSStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
8595777c63cSStefano Zampini     ierr = VecRestoreArray(local_vec2,&array2);CHKERRQ(ierr);
860674ae819SStefano Zampini   }
861674ae819SStefano Zampini   /* Take into account Neumann nodes */
862674ae819SStefano Zampini   if (neumann_is) {
8637e5be6c5SStefano Zampini     ierr = VecSet(local_vec,0.0);CHKERRQ(ierr);
864674ae819SStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
86582ba6b80SStefano Zampini     ierr = ISGetLocalSize(neumann_is,&is_size);CHKERRQ(ierr);
866674ae819SStefano Zampini     ierr = ISGetIndices(neumann_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
867674ae819SStefano Zampini     for (i=0;i<is_size;i++) {
868d50ed60aSStefano Zampini       if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
869674ae819SStefano Zampini         array[is_indices[i]] = 1.0;
870674ae819SStefano Zampini       }
8713c629590SStefano Zampini     }
872674ae819SStefano Zampini     ierr = ISRestoreIndices(neumann_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
873674ae819SStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
874674ae819SStefano Zampini     /* Neumann nodes: impose consistency among neighbours */
875674ae819SStefano Zampini     ierr = VecSet(global_vec,0.0);CHKERRQ(ierr);
876674ae819SStefano Zampini     ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
877674ae819SStefano Zampini     ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
878674ae819SStefano Zampini     ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
879674ae819SStefano Zampini     ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
880674ae819SStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
881674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
8823c629590SStefano Zampini       if (PetscRealPart(array[i]) > 0.1) {
8830cf82fd6SStefano Zampini         graph->special_dof[i] = PCBDDCGRAPH_NEUMANN_MARK;
884674ae819SStefano Zampini       }
885674ae819SStefano Zampini     }
886674ae819SStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
8877e5be6c5SStefano Zampini   }
8887e5be6c5SStefano Zampini   /* Take into account Dirichlet nodes (they overwrite any Neumann boundary node previously set) */
889674ae819SStefano Zampini   if (dirichlet_is) {
8907e5be6c5SStefano Zampini     ierr = VecSet(local_vec,0.0);CHKERRQ(ierr);
891674ae819SStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
89282ba6b80SStefano Zampini     ierr = ISGetLocalSize(dirichlet_is,&is_size);CHKERRQ(ierr);
893674ae819SStefano Zampini     ierr = ISGetIndices(dirichlet_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
894674ae819SStefano Zampini     for (i=0;i<is_size;i++){
895d50ed60aSStefano Zampini       if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
8967e5be6c5SStefano Zampini         array[is_indices[i]] = 1.0;
897674ae819SStefano Zampini       }
8983c629590SStefano Zampini     }
899674ae819SStefano Zampini     ierr = ISRestoreIndices(dirichlet_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
900674ae819SStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
901674ae819SStefano Zampini     /* Dirichlet nodes: impose consistency among neighbours */
902674ae819SStefano Zampini     ierr = VecSet(global_vec,0.0);CHKERRQ(ierr);
9037e5be6c5SStefano Zampini     ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9047e5be6c5SStefano Zampini     ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
905674ae819SStefano Zampini     ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
906674ae819SStefano Zampini     ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
907674ae819SStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
908674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
9093c629590SStefano Zampini       if (PetscRealPart(array[i]) > 0.1) {
9104b2aedd3SStefano Zampini         if (commsize > graph->commsizelimit) { /* dirichlet nodes treated as internal */
911df48933bSStefano Zampini           ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
912b3cdcd63SStefano Zampini           graph->subset[i] = 0;
913b3cdcd63SStefano Zampini         }
9140cf82fd6SStefano Zampini         graph->special_dof[i] = PCBDDCGRAPH_DIRICHLET_MARK;
915674ae819SStefano Zampini       }
916674ae819SStefano Zampini     }
917674ae819SStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
9187fb0e2dbSStefano Zampini   }
91908a5cf49SStefano Zampini   /* mark local periodic nodes (if any) and adapt CSR graph (if any) */
92051b0f6cfSStefano Zampini   if (graph->mirrors) {
92151b0f6cfSStefano Zampini     for (i=0;i<graph->nvtxs;i++)
92251b0f6cfSStefano Zampini       if (graph->mirrors[i])
9230cf82fd6SStefano Zampini         graph->special_dof[i] = PCBDDCGRAPH_LOCAL_PERIODIC_MARK;
92451b0f6cfSStefano Zampini 
92508a5cf49SStefano Zampini     if (graph->xadj && graph->adjncy) {
92608a5cf49SStefano Zampini       PetscInt *new_xadj,*new_adjncy;
92751b0f6cfSStefano Zampini       /* sort CSR graph */
92851b0f6cfSStefano Zampini       for (i=0;i<graph->nvtxs;i++)
92951b0f6cfSStefano Zampini         ierr = PetscSortInt(graph->xadj[i+1]-graph->xadj[i],&graph->adjncy[graph->xadj[i]]);CHKERRQ(ierr);
93051b0f6cfSStefano Zampini 
93151b0f6cfSStefano Zampini       /* adapt local CSR graph in case of local periodicity */
93251b0f6cfSStefano Zampini       k = 0;
93351b0f6cfSStefano Zampini       for (i=0;i<graph->nvtxs;i++)
93451b0f6cfSStefano Zampini         for (j=graph->xadj[i];j<graph->xadj[i+1];j++)
93551b0f6cfSStefano Zampini           k += graph->mirrors[graph->adjncy[j]];
93651b0f6cfSStefano Zampini 
937854ce69bSBarry Smith       ierr = PetscMalloc1(graph->nvtxs+1,&new_xadj);CHKERRQ(ierr);
938854ce69bSBarry Smith       ierr = PetscMalloc1(k+graph->xadj[graph->nvtxs],&new_adjncy);CHKERRQ(ierr);
93951b0f6cfSStefano Zampini       new_xadj[0] = 0;
94051b0f6cfSStefano Zampini       for (i=0;i<graph->nvtxs;i++) {
94151b0f6cfSStefano Zampini         k = graph->xadj[i+1]-graph->xadj[i];
94251b0f6cfSStefano Zampini         ierr = PetscMemcpy(&new_adjncy[new_xadj[i]],&graph->adjncy[graph->xadj[i]],k*sizeof(PetscInt));CHKERRQ(ierr);
94351b0f6cfSStefano Zampini         new_xadj[i+1] = new_xadj[i]+k;
94451b0f6cfSStefano Zampini         for (j=graph->xadj[i];j<graph->xadj[i+1];j++) {
94551b0f6cfSStefano Zampini           k = graph->mirrors[graph->adjncy[j]];
94651b0f6cfSStefano Zampini           ierr = PetscMemcpy(&new_adjncy[new_xadj[i+1]],graph->mirrors_set[graph->adjncy[j]],k*sizeof(PetscInt));CHKERRQ(ierr);
94751b0f6cfSStefano Zampini           new_xadj[i+1] += k;
94851b0f6cfSStefano Zampini         }
94951b0f6cfSStefano Zampini         k = new_xadj[i+1]-new_xadj[i];
95051b0f6cfSStefano Zampini         ierr = PetscSortRemoveDupsInt(&k,&new_adjncy[new_xadj[i]]);CHKERRQ(ierr);
95151b0f6cfSStefano Zampini         new_xadj[i+1] = new_xadj[i]+k;
95251b0f6cfSStefano Zampini       }
95351b0f6cfSStefano Zampini       /* set new CSR into graph */
95451b0f6cfSStefano Zampini       ierr = PetscFree(graph->xadj);CHKERRQ(ierr);
95551b0f6cfSStefano Zampini       ierr = PetscFree(graph->adjncy);CHKERRQ(ierr);
95651b0f6cfSStefano Zampini       graph->xadj = new_xadj;
95751b0f6cfSStefano Zampini       graph->adjncy = new_adjncy;
95851b0f6cfSStefano Zampini     }
95908a5cf49SStefano Zampini   }
96051b0f6cfSStefano Zampini 
96117eb1463SStefano Zampini   /* mark special nodes (if any) -> each will become a single node equivalence class */
962674ae819SStefano Zampini   if (custom_primal_vertices) {
96317eb1463SStefano Zampini     ierr = VecSet(local_vec,0.0);CHKERRQ(ierr);
9649b70a373SStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
96563602bcaSStefano Zampini     ierr = ISGetLocalSize(custom_primal_vertices,&is_size);CHKERRQ(ierr);
966674ae819SStefano Zampini     ierr = ISGetIndices(custom_primal_vertices,(const PetscInt**)&is_indices);CHKERRQ(ierr);
967674ae819SStefano Zampini     for (i=0;i<is_size;i++){
9689b70a373SStefano Zampini       if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
9699b70a373SStefano Zampini         array[is_indices[i]] = 1.0;
9709b70a373SStefano Zampini       }
971674ae819SStefano Zampini     }
972674ae819SStefano Zampini     ierr = ISRestoreIndices(custom_primal_vertices,(const PetscInt**)&is_indices);CHKERRQ(ierr);
9739b70a373SStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
9749b70a373SStefano Zampini     /* special nodes: impose consistency among neighbours */
9759b70a373SStefano Zampini     ierr = VecSet(global_vec,0.0);CHKERRQ(ierr);
9769b70a373SStefano Zampini     ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9779b70a373SStefano Zampini     ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9789b70a373SStefano Zampini     ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9799b70a373SStefano Zampini     ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9809b70a373SStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
9819b70a373SStefano Zampini     j = 0;
9829b70a373SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
9839b70a373SStefano Zampini       if (PetscRealPart(array[i]) > 0.1 && graph->special_dof[i] != PCBDDCGRAPH_DIRICHLET_MARK) {
9849b70a373SStefano Zampini         graph->special_dof[i] = PCBDDCGRAPH_SPECIAL_MARK-j;
9859b70a373SStefano Zampini         j++;
9869b70a373SStefano Zampini       }
9879b70a373SStefano Zampini     }
9889b70a373SStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
98917eb1463SStefano Zampini   }
9909b70a373SStefano Zampini 
9914b2aedd3SStefano Zampini   /* mark interior nodes (if commsize > graph->commsizelimit) as touched and belonging to partition number 0 */
9924b2aedd3SStefano Zampini   if (commsize > graph->commsizelimit) {
993674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
994674ae819SStefano Zampini       if (!graph->count[i]) {
995df48933bSStefano Zampini         ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
996674ae819SStefano Zampini         graph->subset[i] = 0;
997674ae819SStefano Zampini       }
998674ae819SStefano Zampini     }
999b3cdcd63SStefano Zampini   }
1000b3cdcd63SStefano Zampini 
1001674ae819SStefano Zampini   /* init graph structure and compute default subsets */
1002674ae819SStefano Zampini   nodes_touched = 0;
1003674ae819SStefano Zampini   for (i=0;i<graph->nvtxs;i++) {
1004df48933bSStefano Zampini     if (PetscBTLookup(graph->touched,i)) {
1005674ae819SStefano Zampini       nodes_touched++;
1006674ae819SStefano Zampini     }
1007674ae819SStefano Zampini   }
1008674ae819SStefano Zampini   i = 0;
1009674ae819SStefano Zampini   graph->ncc = 0;
1010674ae819SStefano Zampini   total_counts = 0;
10117babac9bSStefano Zampini 
10127babac9bSStefano Zampini   /* allocated space for queues */
10134b2aedd3SStefano Zampini   if (commsize == graph->commsizelimit) {
10147babac9bSStefano Zampini     ierr = PetscMalloc2(graph->nvtxs+1,&graph->cptr,graph->nvtxs,&graph->queue);CHKERRQ(ierr);
10157babac9bSStefano Zampini   } else {
10167babac9bSStefano Zampini     PetscInt nused = graph->nvtxs - nodes_touched;
10177babac9bSStefano Zampini     ierr = PetscMalloc2(nused+1,&graph->cptr,nused,&graph->queue);CHKERRQ(ierr);
10187babac9bSStefano Zampini   }
10197babac9bSStefano Zampini 
1020674ae819SStefano Zampini   while (nodes_touched<graph->nvtxs) {
1021674ae819SStefano Zampini     /*  find first untouched node in local ordering */
1022b3cdcd63SStefano Zampini     while (PetscBTLookup(graph->touched,i)) i++;
1023df48933bSStefano Zampini     ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
1024674ae819SStefano Zampini     graph->subset[i] = graph->ncc+1;
1025674ae819SStefano Zampini     graph->cptr[graph->ncc] = total_counts;
1026674ae819SStefano Zampini     graph->queue[total_counts] = i;
1027674ae819SStefano Zampini     total_counts++;
1028674ae819SStefano Zampini     nodes_touched++;
1029674ae819SStefano Zampini     /* now find all other nodes having the same set of sharing subdomains */
1030674ae819SStefano Zampini     for (j=i+1;j<graph->nvtxs;j++) {
103174e413f5SStefano Zampini       /* check for same number of sharing subdomains, dof number and same special mark */
1032df48933bSStefano 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]) {
1033674ae819SStefano Zampini         /* check for same set of sharing subdomains */
1034674ae819SStefano Zampini         same_set = PETSC_TRUE;
1035674ae819SStefano Zampini         for (k=0;k<graph->count[j];k++){
1036674ae819SStefano Zampini           if (graph->neighbours_set[i][k] != graph->neighbours_set[j][k]) {
1037674ae819SStefano Zampini             same_set = PETSC_FALSE;
1038674ae819SStefano Zampini           }
1039674ae819SStefano Zampini         }
1040674ae819SStefano Zampini         /* I found a friend of mine */
1041674ae819SStefano Zampini         if (same_set) {
1042df48933bSStefano Zampini           ierr = PetscBTSet(graph->touched,j);CHKERRQ(ierr);
1043b3cdcd63SStefano Zampini           graph->subset[j] = graph->ncc+1;
1044674ae819SStefano Zampini           nodes_touched++;
1045674ae819SStefano Zampini           graph->queue[total_counts] = j;
1046674ae819SStefano Zampini           total_counts++;
1047674ae819SStefano Zampini         }
1048674ae819SStefano Zampini       }
1049674ae819SStefano Zampini     }
1050674ae819SStefano Zampini     graph->ncc++;
1051674ae819SStefano Zampini   }
1052b3cdcd63SStefano 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 */
1053674ae819SStefano Zampini   graph->n_subsets = graph->ncc;
1054785e854fSJed Brown   ierr = PetscMalloc1(graph->n_subsets,&graph->subset_ncc);CHKERRQ(ierr);
1055674ae819SStefano Zampini   for (i=0;i<graph->n_subsets;i++) {
1056674ae819SStefano Zampini     graph->subset_ncc[i] = 1;
1057674ae819SStefano Zampini   }
1058674ae819SStefano Zampini   /* final pointer */
1059674ae819SStefano Zampini   graph->cptr[graph->ncc] = total_counts;
1060ec1c413dSStefano Zampini 
1061b3cdcd63SStefano Zampini   /* For consistency reasons (among neighbours), I need to sort (by global ordering) each connected component */
1062ec1c413dSStefano Zampini   /* Get a reference node (min index in global ordering) for each subset for tagging messages */
1063dc456d91SStefano Zampini   ierr = PetscMalloc1(graph->ncc,&graph->subset_ref_node);CHKERRQ(ierr);
1064dc456d91SStefano Zampini   ierr = PetscMalloc1(graph->cptr[graph->ncc],&queue_global);CHKERRQ(ierr);
10653a5b03d1SStefano Zampini   ierr = ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_global);CHKERRQ(ierr);
1066ec1c413dSStefano Zampini   for (j=0;j<graph->ncc;j++) {
1067ec1c413dSStefano Zampini     ierr = PetscSortIntWithArray(graph->cptr[j+1]-graph->cptr[j],&queue_global[graph->cptr[j]],&graph->queue[graph->cptr[j]]);CHKERRQ(ierr);
1068dc456d91SStefano Zampini     graph->subset_ref_node[j] = graph->queue[graph->cptr[j]];
1069f0321c90SStefano Zampini   }
1070dc456d91SStefano Zampini   ierr = PetscFree(queue_global);CHKERRQ(ierr);
1071acc113dbSStefano Zampini   graph->queue_sorted = PETSC_TRUE;
1072b3cdcd63SStefano Zampini   /* save information on subsets (needed when analyzing the connected components) */
10732f566a09SStefano Zampini   if (graph->ncc) {
1074b3cdcd63SStefano Zampini     ierr = PetscMalloc2(graph->ncc,&graph->subset_size,graph->ncc,&graph->subset_idxs);CHKERRQ(ierr);
1075b3cdcd63SStefano Zampini     ierr = PetscMalloc1(graph->cptr[graph->ncc],&graph->subset_idxs[0]);CHKERRQ(ierr);
1076b3cdcd63SStefano Zampini     ierr = PetscMemzero(graph->subset_idxs[0],graph->cptr[graph->ncc]*sizeof(PetscInt));CHKERRQ(ierr);
1077ec1c413dSStefano Zampini     for (j=1;j<graph->ncc;j++) {
1078b3cdcd63SStefano Zampini       graph->subset_size[j-1] = graph->cptr[j] - graph->cptr[j-1];
1079b3cdcd63SStefano Zampini       graph->subset_idxs[j] = graph->subset_idxs[j-1] + graph->subset_size[j-1];
1080ec1c413dSStefano Zampini     }
1081b3cdcd63SStefano Zampini     graph->subset_size[graph->ncc-1] = graph->cptr[graph->ncc] - graph->cptr[graph->ncc-1];
1082b3cdcd63SStefano Zampini     ierr = PetscMemcpy(graph->subset_idxs[0],graph->queue,graph->cptr[graph->ncc]*sizeof(PetscInt));CHKERRQ(ierr);
1083ec1c413dSStefano Zampini   }
1084b3cdcd63SStefano Zampini 
1085f0321c90SStefano Zampini   /* renumber reference nodes */
1086dc456d91SStefano Zampini   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)(graph->l2gmap)),graph->ncc,graph->subset_ref_node,PETSC_COPY_VALUES,&subset_n);CHKERRQ(ierr);
1087dc456d91SStefano Zampini   ierr = ISLocalToGlobalMappingApplyIS(graph->l2gmap,subset_n,&subset);CHKERRQ(ierr);
1088dc456d91SStefano Zampini   ierr = ISDestroy(&subset_n);CHKERRQ(ierr);
10896583bcc1SStefano Zampini   ierr = ISRenumber(subset,NULL,NULL,&subset_n);CHKERRQ(ierr);
1090dc456d91SStefano Zampini   ierr = ISDestroy(&subset);CHKERRQ(ierr);
1091dc456d91SStefano Zampini   ierr = ISGetLocalSize(subset_n,&k);CHKERRQ(ierr);
10926c4ed002SBarry Smith   if (k != graph->ncc) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid size of new subset! %D != %D",k,graph->ncc);
1093dc456d91SStefano Zampini   ierr = ISGetIndices(subset_n,&is_indices);CHKERRQ(ierr);
1094dc456d91SStefano Zampini   ierr = PetscMemcpy(graph->subset_ref_node,is_indices,graph->ncc*sizeof(PetscInt));CHKERRQ(ierr);
1095dc456d91SStefano Zampini   ierr = ISRestoreIndices(subset_n,&is_indices);CHKERRQ(ierr);
1096dc456d91SStefano Zampini   ierr = ISDestroy(&subset_n);CHKERRQ(ierr);
1097ec1c413dSStefano Zampini 
1098ec1c413dSStefano Zampini   /* free workspace */
1099674ae819SStefano Zampini   ierr = VecDestroy(&local_vec);CHKERRQ(ierr);
1100674ae819SStefano Zampini   ierr = VecDestroy(&local_vec2);CHKERRQ(ierr);
1101674ae819SStefano Zampini   ierr = VecDestroy(&global_vec);CHKERRQ(ierr);
1102674ae819SStefano Zampini   ierr = VecScatterDestroy(&scatter_ctx);CHKERRQ(ierr);
1103b3cdcd63SStefano Zampini   graph->setupcalled = PETSC_TRUE;
1104674ae819SStefano Zampini   PetscFunctionReturn(0);
1105674ae819SStefano Zampini }
1106674ae819SStefano Zampini 
1107674ae819SStefano Zampini #undef __FUNCT__
1108674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphResetCSR"
1109674ae819SStefano Zampini PetscErrorCode PCBDDCGraphResetCSR(PCBDDCGraph graph)
1110674ae819SStefano Zampini {
1111674ae819SStefano Zampini   PetscErrorCode ierr;
1112674ae819SStefano Zampini 
1113674ae819SStefano Zampini   PetscFunctionBegin;
1114a1dbd327SStefano Zampini   if (graph->freecsr) {
1115674ae819SStefano Zampini     ierr = PetscFree(graph->xadj);CHKERRQ(ierr);
1116674ae819SStefano Zampini     ierr = PetscFree(graph->adjncy);CHKERRQ(ierr);
1117a1dbd327SStefano Zampini   } else {
1118a1dbd327SStefano Zampini     graph->xadj = NULL;
1119a1dbd327SStefano Zampini     graph->adjncy = NULL;
1120a1dbd327SStefano Zampini   }
1121c8272957SStefano Zampini   graph->freecsr = PETSC_FALSE;
1122575ad6abSStefano Zampini   graph->nvtxs_csr = 0;
1123674ae819SStefano Zampini   PetscFunctionReturn(0);
1124674ae819SStefano Zampini }
1125674ae819SStefano Zampini 
1126674ae819SStefano Zampini #undef __FUNCT__
1127674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphReset"
1128674ae819SStefano Zampini PetscErrorCode PCBDDCGraphReset(PCBDDCGraph graph)
1129674ae819SStefano Zampini {
1130674ae819SStefano Zampini   PetscErrorCode ierr;
1131674ae819SStefano Zampini 
1132674ae819SStefano Zampini   PetscFunctionBegin;
1133*fb597685SStefano Zampini   ierr = PCBDDCGraphResetCSR(graph);CHKERRQ(ierr);
1134674ae819SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&graph->l2gmap);CHKERRQ(ierr);
1135674ae819SStefano Zampini   ierr = PetscFree(graph->subset_ncc);CHKERRQ(ierr);
11363a5b03d1SStefano Zampini   ierr = PetscFree(graph->subset_ref_node);CHKERRQ(ierr);
1137674ae819SStefano Zampini   if (graph->nvtxs) {
1138674ae819SStefano Zampini     ierr = PetscFree(graph->neighbours_set[0]);CHKERRQ(ierr);
1139674ae819SStefano Zampini   }
1140df48933bSStefano Zampini   ierr = PetscBTDestroy(&graph->touched);CHKERRQ(ierr);
11417babac9bSStefano Zampini   ierr = PetscFree5(graph->count,
1142674ae819SStefano Zampini                     graph->neighbours_set,
1143674ae819SStefano Zampini                     graph->subset,
1144674ae819SStefano Zampini                     graph->which_dof,
1145df48933bSStefano Zampini                     graph->special_dof);CHKERRQ(ierr);
11467babac9bSStefano Zampini   ierr = PetscFree2(graph->cptr,graph->queue);CHKERRQ(ierr);
114751b0f6cfSStefano Zampini   if (graph->mirrors) {
114851b0f6cfSStefano Zampini     ierr = PetscFree(graph->mirrors_set[0]);CHKERRQ(ierr);
114951b0f6cfSStefano Zampini   }
115051b0f6cfSStefano Zampini   ierr = PetscFree2(graph->mirrors,graph->mirrors_set);CHKERRQ(ierr);
1151b3cdcd63SStefano Zampini   if (graph->subset_idxs) {
1152b3cdcd63SStefano Zampini     ierr = PetscFree(graph->subset_idxs[0]);CHKERRQ(ierr);
1153ec1c413dSStefano Zampini   }
1154b3cdcd63SStefano Zampini   ierr = PetscFree2(graph->subset_size,graph->subset_idxs);CHKERRQ(ierr);
1155a07ea27aSStefano Zampini   ierr = ISDestroy(&graph->dirdofs);CHKERRQ(ierr);
1156d62866d3SStefano Zampini   ierr = ISDestroy(&graph->dirdofsB);CHKERRQ(ierr);
11574f1b2e48SStefano Zampini   if (graph->n_local_subs) {
11584f1b2e48SStefano Zampini     ierr = PetscFree(graph->local_subs);CHKERRQ(ierr);
11594f1b2e48SStefano Zampini   }
1160a07ea27aSStefano Zampini   graph->has_dirichlet = PETSC_FALSE;
1161674ae819SStefano Zampini   graph->nvtxs = 0;
11627babac9bSStefano Zampini   graph->nvtxs_global = 0;
1163674ae819SStefano Zampini   graph->n_subsets = 0;
1164674ae819SStefano Zampini   graph->custom_minimal_size = 1;
11654f1b2e48SStefano Zampini   graph->n_local_subs = 0;
1166*fb597685SStefano Zampini   graph->setupcalled = PETSC_FALSE;
1167674ae819SStefano Zampini   PetscFunctionReturn(0);
1168674ae819SStefano Zampini }
1169674ae819SStefano Zampini 
1170674ae819SStefano Zampini #undef __FUNCT__
1171674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphInit"
11727fb0e2dbSStefano Zampini PetscErrorCode PCBDDCGraphInit(PCBDDCGraph graph, ISLocalToGlobalMapping l2gmap, PetscInt N)
1173674ae819SStefano Zampini {
1174674ae819SStefano Zampini   PetscInt       n;
1175674ae819SStefano Zampini   PetscErrorCode ierr;
1176674ae819SStefano Zampini 
1177674ae819SStefano Zampini   PetscFunctionBegin;
1178674ae819SStefano Zampini   PetscValidPointer(graph,1);
1179674ae819SStefano Zampini   PetscValidHeaderSpecific(l2gmap,IS_LTOGM_CLASSID,2);
11807babac9bSStefano Zampini   PetscValidLogicalCollectiveInt(l2gmap,N,3);
11818e61c736SStefano Zampini   /* raise an error if already allocated */
11826c4ed002SBarry Smith   if (graph->nvtxs_global) SETERRQ(PetscObjectComm((PetscObject)l2gmap),PETSC_ERR_PLIB,"BDDCGraph already initialized");
1183674ae819SStefano Zampini   /* set number of vertices */
1184674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)l2gmap);CHKERRQ(ierr);
1185674ae819SStefano Zampini   graph->l2gmap = l2gmap;
1186674ae819SStefano Zampini   ierr = ISLocalToGlobalMappingGetSize(l2gmap,&n);CHKERRQ(ierr);
1187674ae819SStefano Zampini   graph->nvtxs = n;
11887fb0e2dbSStefano Zampini   graph->nvtxs_global = N;
1189674ae819SStefano Zampini   /* allocate used space */
1190df48933bSStefano Zampini   ierr = PetscBTCreate(graph->nvtxs,&graph->touched);CHKERRQ(ierr);
11917babac9bSStefano Zampini   ierr = PetscMalloc5(graph->nvtxs,&graph->count,
11928e5aaad5SJed Brown                       graph->nvtxs,&graph->neighbours_set,
11938e5aaad5SJed Brown                       graph->nvtxs,&graph->subset,
11948e5aaad5SJed Brown                       graph->nvtxs,&graph->which_dof,
11958e5aaad5SJed Brown                       graph->nvtxs,&graph->special_dof);CHKERRQ(ierr);
1196674ae819SStefano Zampini   /* zeroes memory */
1197674ae819SStefano Zampini   ierr = PetscMemzero(graph->count,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
1198674ae819SStefano Zampini   ierr = PetscMemzero(graph->subset,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
119963602bcaSStefano Zampini   /* use -1 as a default value for which_dof array */
120063602bcaSStefano Zampini   for (n=0;n<graph->nvtxs;n++) graph->which_dof[n] = -1;
120174e413f5SStefano Zampini   ierr = PetscMemzero(graph->special_dof,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
1202674ae819SStefano Zampini   /* zeroes first pointer to neighbour set */
1203674ae819SStefano Zampini   if (graph->nvtxs) {
1204674ae819SStefano Zampini     graph->neighbours_set[0] = 0;
1205674ae819SStefano Zampini   }
1206674ae819SStefano Zampini   /* zeroes workspace for values of ncc */
1207674ae819SStefano Zampini   graph->subset_ncc = 0;
12083a5b03d1SStefano Zampini   graph->subset_ref_node = 0;
1209a1dbd327SStefano Zampini   /* default flag for csr */
1210c8272957SStefano Zampini   graph->freecsr = PETSC_FALSE;
1211674ae819SStefano Zampini   PetscFunctionReturn(0);
1212674ae819SStefano Zampini }
1213674ae819SStefano Zampini 
1214674ae819SStefano Zampini #undef __FUNCT__
1215674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphDestroy"
1216674ae819SStefano Zampini PetscErrorCode PCBDDCGraphDestroy(PCBDDCGraph* graph)
1217674ae819SStefano Zampini {
1218674ae819SStefano Zampini   PetscErrorCode ierr;
1219674ae819SStefano Zampini 
1220674ae819SStefano Zampini   PetscFunctionBegin;
1221674ae819SStefano Zampini   ierr = PCBDDCGraphReset(*graph);CHKERRQ(ierr);
1222674ae819SStefano Zampini   ierr = PetscFree(*graph);CHKERRQ(ierr);
1223674ae819SStefano Zampini   PetscFunctionReturn(0);
1224674ae819SStefano Zampini }
1225674ae819SStefano Zampini 
1226674ae819SStefano Zampini #undef __FUNCT__
1227674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphCreate"
1228674ae819SStefano Zampini PetscErrorCode PCBDDCGraphCreate(PCBDDCGraph *graph)
1229674ae819SStefano Zampini {
1230674ae819SStefano Zampini   PCBDDCGraph    new_graph;
1231674ae819SStefano Zampini   PetscErrorCode ierr;
1232674ae819SStefano Zampini 
1233674ae819SStefano Zampini   PetscFunctionBegin;
1234854ce69bSBarry Smith   ierr = PetscNew(&new_graph);CHKERRQ(ierr);
1235674ae819SStefano Zampini   new_graph->custom_minimal_size = 1;
12364b2aedd3SStefano Zampini   new_graph->commsizelimit = 1;
1237674ae819SStefano Zampini   *graph = new_graph;
1238674ae819SStefano Zampini   PetscFunctionReturn(0);
1239674ae819SStefano Zampini }
1240