xref: /petsc/src/ksp/pc/impls/bddc/bddcgraph.c (revision 45b2a5aa91e2ff93a2de5b69a3a548b602a9ce05)
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__
61b968477SStefano Zampini #define __FUNCT__ "PCBDDCGraphGetDirichletDofs"
71b968477SStefano Zampini PetscErrorCode PCBDDCGraphGetDirichletDofs(PCBDDCGraph graph,IS* dirdofs)
81b968477SStefano Zampini {
91b968477SStefano Zampini   PetscInt       i,size,sizer;
101b968477SStefano Zampini   PetscErrorCode ierr;
111b968477SStefano Zampini 
121b968477SStefano Zampini   PetscFunctionBegin;
131b968477SStefano Zampini   size = 0;
141b968477SStefano Zampini   for (i=0;i<graph->nvtxs;i++) {
151b968477SStefano Zampini     if (graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) size++;
161b968477SStefano Zampini   }
171b968477SStefano Zampini   ierr = MPI_Allreduce(&size,&sizer,1,MPIU_INT,MPI_LOR,PetscObjectComm((PetscObject)graph->l2gmap));CHKERRQ(ierr);
181b968477SStefano Zampini   if (sizer) {
191b968477SStefano Zampini     PetscInt *dirdofs_idxs;
201b968477SStefano Zampini 
211b968477SStefano Zampini     ierr = PetscMalloc1(size,&dirdofs_idxs);CHKERRQ(ierr);
221b968477SStefano Zampini     size = 0;
231b968477SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
241b968477SStefano Zampini       if (graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) dirdofs_idxs[size++] = i;
251b968477SStefano Zampini     }
261b968477SStefano Zampini     ierr = ISCreateGeneral(PetscObjectComm((PetscObject)graph->l2gmap),size,dirdofs_idxs,PETSC_OWN_POINTER,dirdofs);CHKERRQ(ierr);
271b968477SStefano Zampini   } else {
281b968477SStefano Zampini     *dirdofs = NULL;
291b968477SStefano Zampini   }
301b968477SStefano Zampini   PetscFunctionReturn(0);
311b968477SStefano Zampini }
321b968477SStefano Zampini 
331b968477SStefano Zampini #undef __FUNCT__
34674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphASCIIView"
35e49050b4SStefano Zampini PetscErrorCode PCBDDCGraphASCIIView(PCBDDCGraph graph, PetscInt verbosity_level, PetscViewer viewer)
36674ae819SStefano Zampini {
372b510759SStefano Zampini   PetscInt       i,j,tabs;
3893fb5fd3SStefano Zampini   PetscInt*      queue_in_global_numbering;
39674ae819SStefano Zampini   PetscErrorCode ierr;
40674ae819SStefano Zampini 
41674ae819SStefano Zampini   PetscFunctionBegin;
42674ae819SStefano Zampini   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
432b510759SStefano Zampini   ierr = PetscViewerASCIIGetTab(viewer,&tabs);CHKERRQ(ierr);
44674ae819SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"--------------------------------------------------\n");CHKERRQ(ierr);
45674ae819SStefano Zampini   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
46674ae819SStefano Zampini   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Local BDDC graph for subdomain %04d\n",PetscGlobalRank);CHKERRQ(ierr);
47674ae819SStefano Zampini   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Number of vertices %d\n",graph->nvtxs);CHKERRQ(ierr);
48e49050b4SStefano Zampini   if (verbosity_level > 1) {
49674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
50674ae819SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"%d:\n",i);CHKERRQ(ierr);
51674ae819SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   which_dof: %d\n",graph->which_dof[i]);CHKERRQ(ierr);
5274e413f5SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   special_dof: %d\n",graph->special_dof[i]);CHKERRQ(ierr);
53674ae819SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   neighbours: %d\n",graph->count[i]);CHKERRQ(ierr);
542b510759SStefano Zampini       ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
55674ae819SStefano Zampini       if (graph->count[i]) {
56674ae819SStefano Zampini         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"     set of neighbours:");CHKERRQ(ierr);
57674ae819SStefano Zampini         for (j=0;j<graph->count[i];j++) {
58674ae819SStefano Zampini           ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d",graph->neighbours_set[i][j]);CHKERRQ(ierr);
59674ae819SStefano Zampini         }
60674ae819SStefano Zampini         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
61674ae819SStefano Zampini       }
622b510759SStefano Zampini       ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr);
632b510759SStefano Zampini       ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
6451b0f6cfSStefano Zampini       if (graph->mirrors) {
6551b0f6cfSStefano Zampini         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   mirrors: %d\n",graph->mirrors[i]);CHKERRQ(ierr);
6651b0f6cfSStefano Zampini         if (graph->mirrors[i]) {
672b510759SStefano Zampini           ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
6851b0f6cfSStefano Zampini           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"     set of mirrors:");CHKERRQ(ierr);
6951b0f6cfSStefano Zampini           for (j=0;j<graph->mirrors[i];j++) {
7051b0f6cfSStefano Zampini             ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d",graph->mirrors_set[i][j]);CHKERRQ(ierr);
7151b0f6cfSStefano Zampini           }
7251b0f6cfSStefano Zampini           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
732b510759SStefano Zampini           ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr);
742b510759SStefano Zampini           ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
7551b0f6cfSStefano Zampini         }
7651b0f6cfSStefano Zampini       }
77e49050b4SStefano Zampini       if (verbosity_level > 2) {
78674ae819SStefano Zampini         if (graph->xadj && graph->adjncy) {
79674ae819SStefano Zampini           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   local adj list:");CHKERRQ(ierr);
802b510759SStefano Zampini           ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
81674ae819SStefano Zampini           for (j=graph->xadj[i];j<graph->xadj[i+1];j++) {
82674ae819SStefano Zampini             ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d",graph->adjncy[j]);CHKERRQ(ierr);
83674ae819SStefano Zampini           }
84674ae819SStefano Zampini           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
852b510759SStefano Zampini           ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr);
862b510759SStefano Zampini           ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
87ec1c413dSStefano Zampini         } else {
88ec1c413dSStefano Zampini           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   no adj info\n");CHKERRQ(ierr);
89674ae819SStefano Zampini         }
90e49050b4SStefano Zampini       }
91674ae819SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   interface subset id: %d\n",graph->subset[i]);CHKERRQ(ierr);
92674ae819SStefano Zampini       if (graph->subset[i] && graph->subset_ncc) {
93674ae819SStefano Zampini         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   ncc for subset: %d\n",graph->subset_ncc[graph->subset[i]-1]);CHKERRQ(ierr);
94674ae819SStefano Zampini       }
95674ae819SStefano Zampini     }
96e49050b4SStefano Zampini   }
97674ae819SStefano Zampini   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Total number of connected components %d\n",graph->ncc);CHKERRQ(ierr);
98785e854fSJed Brown   ierr = PetscMalloc1(graph->cptr[graph->ncc],&queue_in_global_numbering);CHKERRQ(ierr);
9993fb5fd3SStefano Zampini   ierr = ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_in_global_numbering);CHKERRQ(ierr);
100674ae819SStefano Zampini   for (i=0;i<graph->ncc;i++) {
1011ce3d396SStefano Zampini     PetscInt node_num=graph->queue[graph->cptr[i]];
1021ce3d396SStefano Zampini     PetscBool printcc = PETSC_FALSE;
103674ae819SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  %d (neighs:",i);CHKERRQ(ierr);
1042b510759SStefano Zampini     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
105674ae819SStefano Zampini     for (j=0;j<graph->count[node_num];j++) {
106674ae819SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d",graph->neighbours_set[node_num][j]);CHKERRQ(ierr);
107674ae819SStefano Zampini     }
108674ae819SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"):");CHKERRQ(ierr);
1091ce3d396SStefano Zampini     if (verbosity_level > 1) {
1101ce3d396SStefano Zampini       printcc = PETSC_TRUE;
111e635b14bSStefano Zampini     } else if (graph->count[node_num] > 1 || (graph->count[node_num] == 1 && graph->special_dof[node_num] == PCBDDCGRAPH_NEUMANN_MARK)) {
1121ce3d396SStefano Zampini       printcc = PETSC_TRUE;
1131ce3d396SStefano Zampini     }
1141ce3d396SStefano Zampini     if (printcc) {
115674ae819SStefano Zampini       for (j=graph->cptr[i];j<graph->cptr[i+1];j++) {
11693fb5fd3SStefano Zampini         ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d (%d)",graph->queue[j],queue_in_global_numbering[j]);CHKERRQ(ierr);
117674ae819SStefano Zampini       }
1181ce3d396SStefano Zampini     }
119674ae819SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
1202b510759SStefano Zampini     ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr);
1212b510759SStefano Zampini     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
122674ae819SStefano Zampini   }
12393fb5fd3SStefano Zampini   ierr = PetscFree(queue_in_global_numbering);CHKERRQ(ierr);
124e49050b4SStefano Zampini   if (graph->custom_minimal_size > 1 && verbosity_level > 1) {
125674ae819SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Custom minimal size %d\n",graph->custom_minimal_size);CHKERRQ(ierr);
126674ae819SStefano Zampini   }
127674ae819SStefano Zampini   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
128674ae819SStefano Zampini   PetscFunctionReturn(0);
129674ae819SStefano Zampini }
130674ae819SStefano Zampini 
131674ae819SStefano Zampini #undef __FUNCT__
132674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphGetCandidatesIS"
133a873d5d0SStefano Zampini PetscErrorCode PCBDDCGraphGetCandidatesIS(PCBDDCGraph graph, PetscInt *n_faces, IS *FacesIS[], PetscInt *n_edges, IS *EdgesIS[], IS *VerticesIS)
134674ae819SStefano Zampini {
135674ae819SStefano Zampini   IS             *ISForFaces,*ISForEdges,ISForVertices;
136adfc4fb2SStefano Zampini   PetscInt       i,nfc,nec,nvc,*idx;
137674ae819SStefano Zampini   PetscErrorCode ierr;
138674ae819SStefano Zampini 
139674ae819SStefano Zampini   PetscFunctionBegin;
140674ae819SStefano Zampini   /* loop on ccs to evalute number of faces, edges and vertices */
141674ae819SStefano Zampini   nfc = 0;
142674ae819SStefano Zampini   nec = 0;
143674ae819SStefano Zampini   nvc = 0;
144674ae819SStefano Zampini   for (i=0;i<graph->ncc;i++) {
1451e1f2d93SStefano Zampini     PetscInt repdof = graph->queue[graph->cptr[i]];
146674ae819SStefano Zampini     if (graph->cptr[i+1]-graph->cptr[i] > graph->custom_minimal_size) {
147e635b14bSStefano Zampini       if (graph->count[repdof] == 1 && graph->special_dof[repdof] != PCBDDCGRAPH_NEUMANN_MARK) {
148674ae819SStefano Zampini         nfc++;
149674ae819SStefano Zampini       } else { /* note that nec will be zero in 2d */
150674ae819SStefano Zampini         nec++;
151674ae819SStefano Zampini       }
152674ae819SStefano Zampini     } else {
153674ae819SStefano Zampini       nvc += graph->cptr[i+1]-graph->cptr[i];
154674ae819SStefano Zampini     }
155674ae819SStefano Zampini   }
156adfc4fb2SStefano Zampini 
1578e4af1ccSStefano Zampini   /* check if we are in 2D or 3D */
1588e4af1ccSStefano Zampini   if (graph->twodim) { /* we are in a 2D case -> edges are shared by 2 subregions and faces don't exist */
159674ae819SStefano Zampini     nec = nfc;
160674ae819SStefano Zampini     nfc = 0;
161674ae819SStefano Zampini   }
162adfc4fb2SStefano Zampini 
163674ae819SStefano Zampini   /* allocate IS arrays for faces, edges. Vertices need a single index set. */
164cf5a6209SStefano Zampini   if (FacesIS) {
165785e854fSJed Brown     ierr = PetscMalloc1(nfc,&ISForFaces);CHKERRQ(ierr);
166cf5a6209SStefano Zampini   }
167cf5a6209SStefano Zampini   if (EdgesIS) {
168785e854fSJed Brown     ierr = PetscMalloc1(nec,&ISForEdges);CHKERRQ(ierr);
169cf5a6209SStefano Zampini   }
170cf5a6209SStefano Zampini   if (VerticesIS) {
171785e854fSJed Brown     ierr = PetscMalloc1(nvc,&idx);CHKERRQ(ierr);
172cf5a6209SStefano Zampini   }
173cf5a6209SStefano Zampini 
174674ae819SStefano Zampini   /* loop on ccs to compute index sets for faces and edges */
175acc113dbSStefano Zampini   if (!graph->queue_sorted) {
176acc113dbSStefano Zampini     PetscInt *queue_global;
177acc113dbSStefano Zampini 
178acc113dbSStefano Zampini     ierr = PetscMalloc1(graph->cptr[graph->ncc],&queue_global);CHKERRQ(ierr);
179acc113dbSStefano Zampini     ierr = ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_global);CHKERRQ(ierr);
180acc113dbSStefano Zampini     for (i=0;i<graph->ncc;i++) {
181acc113dbSStefano Zampini       ierr = PetscSortIntWithArray(graph->cptr[i+1]-graph->cptr[i],&queue_global[graph->cptr[i]],&graph->queue[graph->cptr[i]]);CHKERRQ(ierr);
182acc113dbSStefano Zampini     }
183acc113dbSStefano Zampini     ierr = PetscFree(queue_global);CHKERRQ(ierr);
184acc113dbSStefano Zampini     graph->queue_sorted = PETSC_TRUE;
185acc113dbSStefano Zampini   }
186674ae819SStefano Zampini   nfc = 0;
187674ae819SStefano Zampini   nec = 0;
188674ae819SStefano Zampini   for (i=0;i<graph->ncc;i++) {
1891e1f2d93SStefano Zampini     PetscInt repdof = graph->queue[graph->cptr[i]];
190674ae819SStefano Zampini     if (graph->cptr[i+1]-graph->cptr[i] > graph->custom_minimal_size) {
191e635b14bSStefano Zampini       if (graph->count[repdof] == 1 && graph->special_dof[repdof] != PCBDDCGRAPH_NEUMANN_MARK) {
1928e4af1ccSStefano Zampini         if (graph->twodim) {
193cf5a6209SStefano Zampini           if (EdgesIS) {
194674ae819SStefano 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);
195cf5a6209SStefano Zampini           }
196674ae819SStefano Zampini           nec++;
197674ae819SStefano Zampini         } else {
198cf5a6209SStefano Zampini           if (FacesIS) {
199674ae819SStefano 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);
200cf5a6209SStefano Zampini           }
201674ae819SStefano Zampini           nfc++;
202674ae819SStefano Zampini         }
203674ae819SStefano Zampini       } else {
204cf5a6209SStefano Zampini         if (EdgesIS) {
205674ae819SStefano 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);
206cf5a6209SStefano Zampini         }
207674ae819SStefano Zampini         nec++;
208674ae819SStefano Zampini       }
209674ae819SStefano Zampini     }
210674ae819SStefano Zampini   }
211674ae819SStefano Zampini   /* index set for vertices */
212cf5a6209SStefano Zampini   if (VerticesIS) {
213b8ffe317SStefano Zampini     nvc = 0;
214674ae819SStefano Zampini     for (i=0;i<graph->ncc;i++) {
215674ae819SStefano Zampini       if (graph->cptr[i+1]-graph->cptr[i] <= graph->custom_minimal_size) {
216adfc4fb2SStefano Zampini         PetscInt j;
217adfc4fb2SStefano Zampini 
218674ae819SStefano Zampini         for (j=graph->cptr[i];j<graph->cptr[i+1];j++) {
219674ae819SStefano Zampini           idx[nvc]=graph->queue[j];
220674ae819SStefano Zampini           nvc++;
221674ae819SStefano Zampini         }
222674ae819SStefano Zampini       }
223674ae819SStefano Zampini     }
224674ae819SStefano Zampini     /* sort vertex set (by local ordering) */
225674ae819SStefano Zampini     ierr = PetscSortInt(nvc,idx);CHKERRQ(ierr);
226674ae819SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,nvc,idx,PETSC_OWN_POINTER,&ISForVertices);CHKERRQ(ierr);
227674ae819SStefano Zampini   }
228674ae819SStefano Zampini   /* get back info */
229a873d5d0SStefano Zampini   if (n_faces) *n_faces = nfc;
230d7796978SStefano Zampini   if (FacesIS) {
231674ae819SStefano Zampini     *FacesIS = ISForFaces;
232d7796978SStefano Zampini   }
233a873d5d0SStefano Zampini   if (n_edges) *n_edges = nec;
234d7796978SStefano Zampini   if (EdgesIS) {
235674ae819SStefano Zampini     *EdgesIS = ISForEdges;
236d7796978SStefano Zampini   }
237d7796978SStefano Zampini   if (VerticesIS) {
238674ae819SStefano Zampini     *VerticesIS = ISForVertices;
239d7796978SStefano Zampini   }
240674ae819SStefano Zampini   PetscFunctionReturn(0);
241674ae819SStefano Zampini }
242674ae819SStefano Zampini 
243674ae819SStefano Zampini #undef __FUNCT__
244674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphComputeConnectedComponents"
245674ae819SStefano Zampini PetscErrorCode PCBDDCGraphComputeConnectedComponents(PCBDDCGraph graph)
246674ae819SStefano Zampini {
2474a060362SStefano Zampini   PetscBool      adapt_interface_reduced;
248674ae819SStefano Zampini   MPI_Comm       interface_comm;
2494a060362SStefano Zampini   PetscMPIInt    size;
2508e4af1ccSStefano Zampini   PetscInt       i;
2518e4af1ccSStefano Zampini   PetscBool      twodim;
252674ae819SStefano Zampini   PetscErrorCode ierr;
253674ae819SStefano Zampini 
254674ae819SStefano Zampini   PetscFunctionBegin;
255674ae819SStefano Zampini   /* compute connected components locally */
256674ae819SStefano Zampini   ierr = PetscObjectGetComm((PetscObject)(graph->l2gmap),&interface_comm);CHKERRQ(ierr);
257674ae819SStefano Zampini   ierr = PCBDDCGraphComputeConnectedComponentsLocal(graph);CHKERRQ(ierr);
258674ae819SStefano Zampini   /* check consistency of connected components among neighbouring subdomains -> it adapt them in case it is needed */
2594a060362SStefano Zampini   ierr = MPI_Comm_size(interface_comm,&size);CHKERRQ(ierr);
2604a060362SStefano Zampini   adapt_interface_reduced = PETSC_FALSE;
2614a060362SStefano Zampini   if (size > 1) {
2624a060362SStefano Zampini     PetscInt i;
2634a060362SStefano Zampini     PetscBool adapt_interface = PETSC_FALSE;
264674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
265674ae819SStefano Zampini       /* We are not sure that on a given subset of the local interface,
266674ae819SStefano Zampini          with two connected components, the latters be the same among sharing subdomains */
267674ae819SStefano Zampini       if (graph->subset_ncc[i] > 1) {
2684a060362SStefano Zampini         adapt_interface = PETSC_TRUE;
269674ae819SStefano Zampini         break;
270674ae819SStefano Zampini       }
271674ae819SStefano Zampini     }
2724a060362SStefano Zampini     ierr = MPI_Allreduce(&adapt_interface,&adapt_interface_reduced,1,MPIU_BOOL,MPI_LOR,interface_comm);CHKERRQ(ierr);
2734a060362SStefano Zampini   }
274674ae819SStefano Zampini 
275674ae819SStefano Zampini   if (graph->n_subsets && adapt_interface_reduced) {
276ec1c413dSStefano Zampini     /* old csr stuff */
277ec1c413dSStefano Zampini     PetscInt    *aux_new_xadj,*new_xadj,*new_adjncy;
2785b1b9aeaSStefano Zampini     PetscInt    *old_xadj,*old_adjncy;
279ec1c413dSStefano Zampini     PetscBT     subset_cc_adapt;
280ec1c413dSStefano Zampini     /* MPI */
281ec1c413dSStefano Zampini     MPI_Request *send_requests,*recv_requests;
282ec1c413dSStefano Zampini     PetscInt    *send_buffer,*recv_buffer;
283ec1c413dSStefano Zampini     PetscInt    sum_requests,start_of_recv,start_of_send;
284ec1c413dSStefano Zampini     PetscInt    *cum_recv_counts;
285ec1c413dSStefano Zampini     /* others */
286ec1c413dSStefano Zampini     PetscInt    *labels;
287ec1c413dSStefano Zampini     PetscInt    global_subset_counter;
288ec1c413dSStefano Zampini     PetscInt    j,k,s;
2895b1b9aeaSStefano Zampini 
290674ae819SStefano Zampini     /* Retrict adjacency graph using information from previously computed connected components */
291785e854fSJed Brown     ierr = PetscMalloc1(graph->nvtxs,&aux_new_xadj);CHKERRQ(ierr);
292674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
293674ae819SStefano Zampini       aux_new_xadj[i]=1;
294674ae819SStefano Zampini     }
295674ae819SStefano Zampini     for (i=0;i<graph->ncc;i++) {
296674ae819SStefano Zampini       k = graph->cptr[i+1]-graph->cptr[i];
297674ae819SStefano Zampini       for (j=0;j<k;j++) {
298674ae819SStefano Zampini         aux_new_xadj[graph->queue[graph->cptr[i]+j]]=k;
299674ae819SStefano Zampini       }
300674ae819SStefano Zampini     }
301674ae819SStefano Zampini     j = 0;
302674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
303674ae819SStefano Zampini       j += aux_new_xadj[i];
304674ae819SStefano Zampini     }
305854ce69bSBarry Smith     ierr = PetscMalloc1(graph->nvtxs+1,&new_xadj);CHKERRQ(ierr);
306785e854fSJed Brown     ierr = PetscMalloc1(j,&new_adjncy);CHKERRQ(ierr);
307674ae819SStefano Zampini     new_xadj[0]=0;
308674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
309674ae819SStefano Zampini       new_xadj[i+1]=new_xadj[i]+aux_new_xadj[i];
310674ae819SStefano Zampini       if (aux_new_xadj[i]==1) {
311674ae819SStefano Zampini         new_adjncy[new_xadj[i]]=i;
312674ae819SStefano Zampini       }
313674ae819SStefano Zampini     }
314674ae819SStefano Zampini     ierr = PetscFree(aux_new_xadj);CHKERRQ(ierr);
315ec1c413dSStefano Zampini     ierr = PetscMalloc1(graph->nvtxs,&labels);CHKERRQ(ierr);
316ec1c413dSStefano Zampini     ierr = PetscMemzero(labels,graph->nvtxs*sizeof(*labels));CHKERRQ(ierr);
317674ae819SStefano Zampini     for (i=0;i<graph->ncc;i++) {
318674ae819SStefano Zampini       k = graph->cptr[i+1]-graph->cptr[i];
319674ae819SStefano Zampini       for (j=0;j<k;j++) {
320674ae819SStefano Zampini         ierr = PetscMemcpy(&new_adjncy[new_xadj[graph->queue[graph->cptr[i]+j]]],&graph->queue[graph->cptr[i]],k*sizeof(PetscInt));CHKERRQ(ierr);
321ec1c413dSStefano Zampini         labels[graph->queue[graph->cptr[i]+j]] = i;
322674ae819SStefano Zampini       }
323674ae819SStefano Zampini     }
3245b1b9aeaSStefano Zampini     /* set temporarly new CSR into graph */
3255b1b9aeaSStefano Zampini     old_xadj = graph->xadj;
3265b1b9aeaSStefano Zampini     old_adjncy = graph->adjncy;
327674ae819SStefano Zampini     graph->xadj = new_xadj;
328674ae819SStefano Zampini     graph->adjncy = new_adjncy;
329674ae819SStefano Zampini     /* allocate some space */
330854ce69bSBarry Smith     ierr = PetscMalloc1(graph->n_subsets+1,&cum_recv_counts);CHKERRQ(ierr);
331674ae819SStefano Zampini     ierr = PetscMemzero(cum_recv_counts,(graph->n_subsets+1)*sizeof(*cum_recv_counts));CHKERRQ(ierr);
332674ae819SStefano Zampini     /* first count how many neighbours per connected component I will receive from */
333674ae819SStefano Zampini     cum_recv_counts[0]=0;
334ec1c413dSStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
335ec1c413dSStefano Zampini       cum_recv_counts[i+1]=cum_recv_counts[i]+graph->count[graph->subsets[i][0]];
336674ae819SStefano Zampini     }
337ec1c413dSStefano Zampini     ierr = PetscMalloc1(cum_recv_counts[graph->n_subsets],&recv_buffer);CHKERRQ(ierr);
338ec1c413dSStefano Zampini     ierr = PetscMalloc2(cum_recv_counts[graph->n_subsets],&send_requests,cum_recv_counts[graph->n_subsets],&recv_requests);CHKERRQ(ierr);
339674ae819SStefano Zampini     for (i=0;i<cum_recv_counts[graph->n_subsets];i++) {
340674ae819SStefano Zampini       send_requests[i]=MPI_REQUEST_NULL;
341674ae819SStefano Zampini       recv_requests[i]=MPI_REQUEST_NULL;
342674ae819SStefano Zampini     }
343ec1c413dSStefano Zampini     /* exchange with my neighbours the number of my connected components on the subset of interface */
344674ae819SStefano Zampini     sum_requests = 0;
345674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
346ec1c413dSStefano Zampini       PetscMPIInt neigh,tag;
347ec1c413dSStefano Zampini 
348ec1c413dSStefano Zampini       j = graph->subsets[i][0];
349ec1c413dSStefano Zampini       ierr = PetscMPIIntCast(2*graph->subset_ref_node[i],&tag);CHKERRQ(ierr);
350674ae819SStefano Zampini       for (k=0;k<graph->count[j];k++) {
351674ae819SStefano Zampini         ierr = PetscMPIIntCast(graph->neighbours_set[j][k],&neigh);CHKERRQ(ierr);
352674ae819SStefano Zampini         ierr = MPI_Isend(&graph->subset_ncc[i],1,MPIU_INT,neigh,tag,interface_comm,&send_requests[sum_requests]);CHKERRQ(ierr);
353ec1c413dSStefano Zampini         ierr = MPI_Irecv(&recv_buffer[sum_requests],1,MPIU_INT,neigh,tag,interface_comm,&recv_requests[sum_requests]);CHKERRQ(ierr);
354674ae819SStefano Zampini         sum_requests++;
355674ae819SStefano Zampini       }
356674ae819SStefano Zampini     }
357674ae819SStefano Zampini     ierr = MPI_Waitall(sum_requests,recv_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
358674ae819SStefano Zampini     ierr = MPI_Waitall(sum_requests,send_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
359ec1c413dSStefano Zampini     /* determine the subsets I have to adapt */
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++) {
363674ae819SStefano Zampini       for (j=cum_recv_counts[i];j<cum_recv_counts[i+1];j++){
364ec1c413dSStefano Zampini         /* adapt if more than one cc is present */
365ec1c413dSStefano Zampini          if (graph->subset_ncc[i] > 1 || recv_buffer[j] > 1) {
366ec1c413dSStefano Zampini           ierr = PetscBTSet(subset_cc_adapt,i);
367674ae819SStefano Zampini           break;
368674ae819SStefano Zampini         }
369674ae819SStefano Zampini       }
370674ae819SStefano Zampini     }
371ec1c413dSStefano Zampini     ierr = PetscFree(recv_buffer);CHKERRQ(ierr);
372ec1c413dSStefano Zampini     /* determine send/recv buffers sizes */
373ec1c413dSStefano Zampini     j = 0;
374674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
375ec1c413dSStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i)) {
376ec1c413dSStefano Zampini         j += graph->subsets_size[i];
377674ae819SStefano Zampini       }
378674ae819SStefano Zampini     }
379ec1c413dSStefano Zampini     k = 0;
380674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
381ec1c413dSStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i)) {
382ec1c413dSStefano Zampini         k += (cum_recv_counts[i+1]-cum_recv_counts[i])*graph->subsets_size[i];
383674ae819SStefano Zampini       }
384674ae819SStefano Zampini     }
385ec1c413dSStefano Zampini     ierr = PetscMalloc2(j,&send_buffer,k,&recv_buffer);CHKERRQ(ierr);
386ec1c413dSStefano Zampini 
387ec1c413dSStefano Zampini     /* fill send buffer */
388ec1c413dSStefano Zampini     j = 0;
389ec1c413dSStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
390ec1c413dSStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i)) {
391ec1c413dSStefano Zampini         for (k=0;k<graph->subsets_size[i];k++) {
392ec1c413dSStefano Zampini           send_buffer[j++] = labels[graph->subsets[i][k]];
393674ae819SStefano Zampini         }
394674ae819SStefano Zampini       }
395674ae819SStefano Zampini     }
396ec1c413dSStefano Zampini     ierr = PetscFree(labels);CHKERRQ(ierr);
397ec1c413dSStefano Zampini 
398674ae819SStefano Zampini     /* now exchange the data */
399674ae819SStefano Zampini     start_of_recv = 0;
400674ae819SStefano Zampini     start_of_send = 0;
401674ae819SStefano Zampini     sum_requests = 0;
402674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
403ec1c413dSStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i)) {
404ec1c413dSStefano Zampini         PetscMPIInt neigh,tag;
405ec1c413dSStefano Zampini         PetscInt    size_of_send = graph->subsets_size[i];
406ec1c413dSStefano Zampini 
407ec1c413dSStefano Zampini         j = graph->subsets[i][0];
408ec1c413dSStefano Zampini         ierr = PetscMPIIntCast(2*graph->subset_ref_node[i]+1,&tag);CHKERRQ(ierr);
409674ae819SStefano Zampini         for (k=0;k<graph->count[j];k++) {
410674ae819SStefano Zampini           ierr = PetscMPIIntCast(graph->neighbours_set[j][k],&neigh);CHKERRQ(ierr);
411ec1c413dSStefano Zampini           ierr = MPI_Isend(&send_buffer[start_of_send],size_of_send,MPIU_INT,neigh,tag,interface_comm,&send_requests[sum_requests]);CHKERRQ(ierr);
412ec1c413dSStefano Zampini           ierr = MPI_Irecv(&recv_buffer[start_of_recv],size_of_send,MPIU_INT,neigh,tag,interface_comm,&recv_requests[sum_requests]);CHKERRQ(ierr);
413ec1c413dSStefano Zampini           start_of_recv += size_of_send;
414674ae819SStefano Zampini           sum_requests++;
415674ae819SStefano Zampini         }
416674ae819SStefano Zampini         start_of_send += size_of_send;
417674ae819SStefano Zampini       }
418674ae819SStefano Zampini     }
419674ae819SStefano Zampini     ierr = MPI_Waitall(sum_requests,recv_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
420674ae819SStefano Zampini     ierr = MPI_Waitall(sum_requests,send_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
421ec1c413dSStefano Zampini     ierr = PetscFree2(send_requests,recv_requests);CHKERRQ(ierr);
422d52c4852SStefano Zampini 
423674ae819SStefano Zampini     start_of_recv = 0;
424674ae819SStefano Zampini     global_subset_counter = 0;
425674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
426ec1c413dSStefano Zampini       if (PetscBTLookup(subset_cc_adapt,i)) {
427ec1c413dSStefano Zampini         PetscInt **temp_buffer,temp_buffer_size,*private_label;
428ec1c413dSStefano Zampini 
429674ae819SStefano Zampini         /* allocate some temporary space */
430ec1c413dSStefano Zampini         temp_buffer_size = graph->subsets_size[i];
431785e854fSJed Brown         ierr = PetscMalloc1(temp_buffer_size,&temp_buffer);CHKERRQ(ierr);
432785e854fSJed Brown         ierr = PetscMalloc1(temp_buffer_size*(cum_recv_counts[i+1]-cum_recv_counts[i]),&temp_buffer[0]);CHKERRQ(ierr);
433674ae819SStefano Zampini         ierr = PetscMemzero(temp_buffer[0],temp_buffer_size*(cum_recv_counts[i+1]-cum_recv_counts[i])*sizeof(PetscInt));CHKERRQ(ierr);
434674ae819SStefano Zampini         for (j=1;j<temp_buffer_size;j++) {
435674ae819SStefano Zampini           temp_buffer[j] = temp_buffer[j-1]+cum_recv_counts[i+1]-cum_recv_counts[i];
436674ae819SStefano Zampini         }
437ec1c413dSStefano Zampini         /* analyze contributions from neighbouring subdomains for i-th subset
438674ae819SStefano Zampini            temp buffer structure:
439ec1c413dSStefano Zampini            supposing part of the interface has dimension 5
440ec1c413dSStefano Zampini            e.g with global dofs 0,1,2,3,4, locally ordered on the current process as 0,4,3,1,2
441ec1c413dSStefano Zampini            3 neighs procs having connected components:
442ec1c413dSStefano Zampini              neigh 0: [0 1 4], [2 3], labels [4,7]  (2 connected components)
443ec1c413dSStefano Zampini              neigh 1: [0 1], [2 3 4], labels [3 2]  (2 connected components)
444ec1c413dSStefano Zampini              neigh 2: [0 4], [1], [2 3], labels [1 5 6] (3 connected components)
445674ae819SStefano Zampini            tempbuffer (row-oriented) will be filled as:
446ec1c413dSStefano Zampini              [ 4, 3, 1;
447ec1c413dSStefano Zampini                4, 2, 1;
448ec1c413dSStefano Zampini                7, 2, 6;
449ec1c413dSStefano Zampini                4, 3, 5;
450ec1c413dSStefano Zampini                7, 2, 6; ];
451674ae819SStefano Zampini            This way we can simply find intersections of ccs among neighs.
452ec1c413dSStefano Zampini            The graph->subset array will need to be modified. The output for the example is [0], [1], [2 3], [4];
453674ae819SStefano Zampini                                                                                                                                    */
454674ae819SStefano Zampini         for (j=0;j<cum_recv_counts[i+1]-cum_recv_counts[i];j++) {
455ec1c413dSStefano Zampini           for (k=0;k<temp_buffer_size;k++) {
456ec1c413dSStefano Zampini             temp_buffer[k][j] = recv_buffer[start_of_recv+k];
457674ae819SStefano Zampini           }
458ec1c413dSStefano Zampini           start_of_recv += temp_buffer_size;
459674ae819SStefano Zampini         }
460ec1c413dSStefano Zampini         ierr = PetscMalloc1(temp_buffer_size,&private_label);CHKERRQ(ierr);
461ec1c413dSStefano Zampini         ierr = PetscMemzero(private_label,temp_buffer_size*sizeof(*private_label));CHKERRQ(ierr);
462674ae819SStefano Zampini         for (j=0;j<temp_buffer_size;j++) {
463ec1c413dSStefano Zampini           if (!private_label[j]) { /* found a new cc  */
464ec1c413dSStefano Zampini             PetscBool same_set;
465ec1c413dSStefano Zampini 
466674ae819SStefano Zampini             global_subset_counter++;
467ec1c413dSStefano Zampini             private_label[j] = global_subset_counter;
468674ae819SStefano Zampini             for (k=j+1;k<temp_buffer_size;k++) { /* check for other nodes in new cc */
469674ae819SStefano Zampini               same_set = PETSC_TRUE;
470674ae819SStefano Zampini               for (s=0;s<cum_recv_counts[i+1]-cum_recv_counts[i];s++) {
471674ae819SStefano Zampini                 if (temp_buffer[j][s] != temp_buffer[k][s]) {
472674ae819SStefano Zampini                   same_set = PETSC_FALSE;
473674ae819SStefano Zampini                   break;
474674ae819SStefano Zampini                 }
475674ae819SStefano Zampini               }
476674ae819SStefano Zampini               if (same_set) {
477ec1c413dSStefano Zampini                 private_label[k] = global_subset_counter;
478674ae819SStefano Zampini               }
479674ae819SStefano Zampini             }
480674ae819SStefano Zampini           }
481674ae819SStefano Zampini         }
482ec1c413dSStefano Zampini         /* insert private labels in graph structure */
483ec1c413dSStefano Zampini         for (j=0;j<graph->subsets_size[i];j++) {
484ec1c413dSStefano Zampini           graph->subset[graph->subsets[i][j]] = graph->n_subsets+private_label[j];
485674ae819SStefano Zampini         }
486ec1c413dSStefano Zampini         ierr = PetscFree(private_label);CHKERRQ(ierr);
487674ae819SStefano Zampini         ierr = PetscFree(temp_buffer[0]);CHKERRQ(ierr);
488674ae819SStefano Zampini         ierr = PetscFree(temp_buffer);CHKERRQ(ierr);
489674ae819SStefano Zampini       }
490674ae819SStefano Zampini     }
491*45b2a5aaSStefano Zampini     ierr = PetscFree2(send_buffer,recv_buffer);CHKERRQ(ierr);
492674ae819SStefano Zampini     ierr = PetscFree(cum_recv_counts);CHKERRQ(ierr);
493ec1c413dSStefano Zampini     ierr = PetscBTDestroy(&subset_cc_adapt);CHKERRQ(ierr);
494ec1c413dSStefano Zampini     /* We are ready to find for connected components which are consistent among neighbouring subdomains */
495674ae819SStefano Zampini     if (global_subset_counter) {
496df48933bSStefano Zampini       ierr = PetscBTMemzero(graph->nvtxs,graph->touched);CHKERRQ(ierr);
497674ae819SStefano Zampini       global_subset_counter = 0;
498674ae819SStefano Zampini       for (i=0;i<graph->nvtxs;i++) {
499df48933bSStefano Zampini         if (graph->subset[i] && !PetscBTLookup(graph->touched,i)) {
500674ae819SStefano Zampini           global_subset_counter++;
501674ae819SStefano Zampini           for (j=i+1;j<graph->nvtxs;j++) {
502df48933bSStefano Zampini             if (!PetscBTLookup(graph->touched,j) && graph->subset[j]==graph->subset[i]) {
503674ae819SStefano Zampini               graph->subset[j] = global_subset_counter;
504df48933bSStefano Zampini               ierr = PetscBTSet(graph->touched,j);CHKERRQ(ierr);
505674ae819SStefano Zampini             }
506674ae819SStefano Zampini           }
507674ae819SStefano Zampini           graph->subset[i] = global_subset_counter;
508df48933bSStefano Zampini           ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
509674ae819SStefano Zampini         }
510674ae819SStefano Zampini       }
511674ae819SStefano Zampini       /* refine connected components locally */
512674ae819SStefano Zampini       ierr = PCBDDCGraphComputeConnectedComponentsLocal(graph);CHKERRQ(ierr);
513674ae819SStefano Zampini     }
5145b1b9aeaSStefano Zampini     /* restore original CSR graph of dofs */
5155b1b9aeaSStefano Zampini     ierr = PetscFree(new_xadj);CHKERRQ(ierr);
5165b1b9aeaSStefano Zampini     ierr = PetscFree(new_adjncy);CHKERRQ(ierr);
5175b1b9aeaSStefano Zampini     graph->xadj = old_xadj;
5185b1b9aeaSStefano Zampini     graph->adjncy = old_adjncy;
519674ae819SStefano Zampini   }
5208e4af1ccSStefano Zampini 
5218e4af1ccSStefano Zampini   /* Determine if we are in 2D or 3D */
5228e4af1ccSStefano Zampini   twodim  = PETSC_TRUE;
5238e4af1ccSStefano Zampini   for (i=0;i<graph->ncc;i++) {
5248e4af1ccSStefano Zampini     PetscInt repdof = graph->queue[graph->cptr[i]];
5258e4af1ccSStefano Zampini     if (graph->cptr[i+1]-graph->cptr[i] > graph->custom_minimal_size) {
5268e4af1ccSStefano Zampini       if (graph->count[repdof] > 1 || graph->special_dof[repdof] == PCBDDCGRAPH_NEUMANN_MARK) {
5278e4af1ccSStefano Zampini         twodim = PETSC_FALSE;
5288e4af1ccSStefano Zampini         break;
5298e4af1ccSStefano Zampini       }
5308e4af1ccSStefano Zampini     }
5318e4af1ccSStefano Zampini   }
5328e4af1ccSStefano Zampini   ierr = MPI_Allreduce(&twodim,&graph->twodim,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)graph->l2gmap));CHKERRQ(ierr);
533674ae819SStefano Zampini   PetscFunctionReturn(0);
534674ae819SStefano Zampini }
535674ae819SStefano Zampini 
536674ae819SStefano Zampini /* The following code has been adapted from function IsConnectedSubdomain contained
537674ae819SStefano Zampini    in source file contig.c of METIS library (version 5.0.1)
538674ae819SStefano Zampini    It finds connected components for each subset  */
539674ae819SStefano Zampini #undef __FUNCT__
540674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphComputeConnectedComponentsLocal"
541674ae819SStefano Zampini PetscErrorCode PCBDDCGraphComputeConnectedComponentsLocal(PCBDDCGraph graph)
542674ae819SStefano Zampini {
543674ae819SStefano Zampini   PetscInt       i,j,k,first,last,nleft,ncc,pid,cum_queue,n,ncc_pid;
5444a060362SStefano Zampini   PetscMPIInt    size;
545674ae819SStefano Zampini   PetscErrorCode ierr;
546674ae819SStefano Zampini 
547674ae819SStefano Zampini   PetscFunctionBegin;
548674ae819SStefano Zampini   /* quiet return if no csr info is available */
549674ae819SStefano Zampini   if (!graph->xadj || !graph->adjncy) {
550674ae819SStefano Zampini     PetscFunctionReturn(0);
551674ae819SStefano Zampini   }
552674ae819SStefano Zampini 
553674ae819SStefano Zampini   /* reset any previous search of connected components */
554df48933bSStefano Zampini   ierr = PetscBTMemzero(graph->nvtxs,graph->touched);CHKERRQ(ierr);
555674ae819SStefano Zampini   graph->n_subsets = 0;
5567babac9bSStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)graph->l2gmap),&size);CHKERRQ(ierr);
5574a060362SStefano Zampini   if (size == 1) {
558dc36bf05SStefano Zampini     if (graph->nvtxs) {
5594a060362SStefano Zampini       graph->n_subsets = 1;
5604a060362SStefano Zampini       for (i=0;i<graph->nvtxs;i++) {
5614a060362SStefano Zampini         graph->subset[i] = 1;
5624a060362SStefano Zampini       }
563dc36bf05SStefano Zampini     }
5644a060362SStefano Zampini   } else {
565674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
5660cf82fd6SStefano Zampini       if (graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK || !graph->count[i]) {
567df48933bSStefano Zampini         ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
568674ae819SStefano Zampini         graph->subset[i] = 0;
569674ae819SStefano Zampini       }
570674ae819SStefano Zampini       graph->n_subsets = PetscMax(graph->n_subsets,graph->subset[i]);
571674ae819SStefano Zampini     }
5724a060362SStefano Zampini   }
573674ae819SStefano Zampini   ierr = PetscFree(graph->subset_ncc);CHKERRQ(ierr);
574785e854fSJed Brown   ierr = PetscMalloc1(graph->n_subsets,&graph->subset_ncc);CHKERRQ(ierr);
575674ae819SStefano Zampini 
576674ae819SStefano Zampini   /* begin search for connected components */
577674ae819SStefano Zampini   cum_queue = 0;
578674ae819SStefano Zampini   ncc = 0;
579674ae819SStefano Zampini   for (n=0;n<graph->n_subsets;n++) {
580674ae819SStefano Zampini     pid = n+1;  /* partition labeled by 0 is discarded */
581674ae819SStefano Zampini     nleft = 0;
582674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
583674ae819SStefano Zampini       if (graph->subset[i] == pid) {
584674ae819SStefano Zampini         nleft++;
585674ae819SStefano Zampini       }
586674ae819SStefano Zampini     }
587674ae819SStefano Zampini     for (i=0; i<graph->nvtxs; i++) {
588674ae819SStefano Zampini       if (graph->subset[i] == pid) {
589674ae819SStefano Zampini         break;
590674ae819SStefano Zampini       }
591674ae819SStefano Zampini     }
592df48933bSStefano Zampini     ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
593674ae819SStefano Zampini     graph->queue[cum_queue] = i;
594674ae819SStefano Zampini     first = 0;
595674ae819SStefano Zampini     last = 1;
596674ae819SStefano Zampini     graph->cptr[ncc] = cum_queue;
597674ae819SStefano Zampini     ncc_pid = 0;
598674ae819SStefano Zampini     while (first != nleft) {
599674ae819SStefano Zampini       if (first == last) {
600674ae819SStefano Zampini         graph->cptr[++ncc] = first+cum_queue;
601674ae819SStefano Zampini         ncc_pid++;
602df48933bSStefano Zampini         for (i=0; i<graph->nvtxs; i++) { /* TODO-> use a while! */
603df48933bSStefano Zampini           if (graph->subset[i] == pid && !PetscBTLookup(graph->touched,i)) {
604674ae819SStefano Zampini             break;
605674ae819SStefano Zampini           }
606674ae819SStefano Zampini         }
607674ae819SStefano Zampini         graph->queue[cum_queue+last] = i;
608674ae819SStefano Zampini         last++;
609df48933bSStefano Zampini         ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
610674ae819SStefano Zampini       }
611674ae819SStefano Zampini       i = graph->queue[cum_queue+first];
612674ae819SStefano Zampini       first++;
613674ae819SStefano Zampini       for (j=graph->xadj[i];j<graph->xadj[i+1];j++) {
614674ae819SStefano Zampini         k = graph->adjncy[j];
615df48933bSStefano Zampini         if (graph->subset[k] == pid && !PetscBTLookup(graph->touched,k)) {
616674ae819SStefano Zampini           graph->queue[cum_queue+last] = k;
617674ae819SStefano Zampini           last++;
618df48933bSStefano Zampini           ierr = PetscBTSet(graph->touched,k);CHKERRQ(ierr);
619674ae819SStefano Zampini         }
620674ae819SStefano Zampini       }
621674ae819SStefano Zampini     }
622674ae819SStefano Zampini     graph->cptr[++ncc] = first+cum_queue;
623674ae819SStefano Zampini     ncc_pid++;
624674ae819SStefano Zampini     cum_queue = graph->cptr[ncc];
625674ae819SStefano Zampini     graph->subset_ncc[n] = ncc_pid;
626674ae819SStefano Zampini   }
627674ae819SStefano Zampini   graph->ncc = ncc;
628acc113dbSStefano Zampini   graph->queue_sorted = PETSC_FALSE;
629674ae819SStefano Zampini   PetscFunctionReturn(0);
630674ae819SStefano Zampini }
631674ae819SStefano Zampini 
632674ae819SStefano Zampini #undef __FUNCT__
633674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphSetUp"
634674ae819SStefano Zampini PetscErrorCode PCBDDCGraphSetUp(PCBDDCGraph graph, PetscInt custom_minimal_size, IS neumann_is, IS dirichlet_is, PetscInt n_ISForDofs, IS ISForDofs[], IS custom_primal_vertices)
635674ae819SStefano Zampini {
636674ae819SStefano Zampini   VecScatter     scatter_ctx;
637674ae819SStefano Zampini   Vec            local_vec,local_vec2,global_vec;
638674ae819SStefano Zampini   IS             to,from;
639674ae819SStefano Zampini   MPI_Comm       comm;
640674ae819SStefano Zampini   PetscScalar    *array,*array2;
641674ae819SStefano Zampini   const PetscInt *is_indices;
642f0321c90SStefano Zampini   PetscInt       n_neigh,*neigh,*n_shared,**shared,*queue_global,*subset_ref_node_global;
643674ae819SStefano Zampini   PetscInt       i,j,k,s,total_counts,nodes_touched,is_size;
6447babac9bSStefano Zampini   PetscMPIInt    size;
64551b0f6cfSStefano Zampini   PetscBool      same_set,mirrors_found;
6467babac9bSStefano Zampini   PetscErrorCode ierr;
647674ae819SStefano Zampini 
648674ae819SStefano Zampini   PetscFunctionBegin;
649674ae819SStefano Zampini   ierr = PetscObjectGetComm((PetscObject)(graph->l2gmap),&comm);CHKERRQ(ierr);
650674ae819SStefano Zampini   /* custom_minimal_size */
651674ae819SStefano Zampini   graph->custom_minimal_size = PetscMax(graph->custom_minimal_size,custom_minimal_size);
652674ae819SStefano Zampini   /* get info l2gmap and allocate work vectors  */
653674ae819SStefano Zampini   ierr = ISLocalToGlobalMappingGetInfo(graph->l2gmap,&n_neigh,&neigh,&n_shared,&shared);CHKERRQ(ierr);
6542635a1d4SStefano Zampini   /* check if we have any local periodic nodes (periodic BCs) */
6552635a1d4SStefano Zampini   mirrors_found = PETSC_FALSE;
6562635a1d4SStefano Zampini   if (graph->nvtxs && n_neigh) {
6572635a1d4SStefano Zampini     for (i=0; i<n_shared[0]; i++) graph->count[shared[0][i]] += 1;
6582635a1d4SStefano Zampini     for (i=0; i<n_shared[0]; i++) {
6592635a1d4SStefano Zampini       if (graph->count[shared[0][i]] > 1) {
6602635a1d4SStefano Zampini         mirrors_found = PETSC_TRUE;
6612635a1d4SStefano Zampini         break;
6622635a1d4SStefano Zampini       }
6632635a1d4SStefano Zampini     }
6642635a1d4SStefano Zampini   }
6652635a1d4SStefano Zampini   /* create some workspace objects */
6662635a1d4SStefano Zampini   local_vec = NULL;
6672635a1d4SStefano Zampini   local_vec2 = NULL;
6682635a1d4SStefano Zampini   global_vec = NULL;
6692635a1d4SStefano Zampini   to = NULL;
6702635a1d4SStefano Zampini   from = NULL;
6712635a1d4SStefano Zampini   scatter_ctx = NULL;
6722635a1d4SStefano Zampini   if (n_ISForDofs || dirichlet_is || neumann_is || custom_primal_vertices) {
673674ae819SStefano Zampini     ierr = VecCreate(PETSC_COMM_SELF,&local_vec);CHKERRQ(ierr);
674674ae819SStefano Zampini     ierr = VecSetSizes(local_vec,PETSC_DECIDE,graph->nvtxs);CHKERRQ(ierr);
675674ae819SStefano Zampini     ierr = VecSetType(local_vec,VECSTANDARD);CHKERRQ(ierr);
676674ae819SStefano Zampini     ierr = VecDuplicate(local_vec,&local_vec2);CHKERRQ(ierr);
677674ae819SStefano Zampini     ierr = VecCreate(comm,&global_vec);CHKERRQ(ierr);
6787fb0e2dbSStefano Zampini     ierr = VecSetSizes(global_vec,PETSC_DECIDE,graph->nvtxs_global);CHKERRQ(ierr);
679674ae819SStefano Zampini     ierr = VecSetType(global_vec,VECSTANDARD);CHKERRQ(ierr);
680674ae819SStefano Zampini     ierr = ISCreateStride(PETSC_COMM_SELF,graph->nvtxs,0,1,&to);CHKERRQ(ierr);
681674ae819SStefano Zampini     ierr = ISLocalToGlobalMappingApplyIS(graph->l2gmap,to,&from);CHKERRQ(ierr);
682674ae819SStefano Zampini     ierr = VecScatterCreate(global_vec,from,local_vec,to,&scatter_ctx);CHKERRQ(ierr);
6832635a1d4SStefano Zampini   } else if (mirrors_found) {
6842635a1d4SStefano Zampini     ierr = ISCreateStride(PETSC_COMM_SELF,graph->nvtxs,0,1,&to);CHKERRQ(ierr);
6852635a1d4SStefano Zampini     ierr = ISLocalToGlobalMappingApplyIS(graph->l2gmap,to,&from);CHKERRQ(ierr);
68649eeff8cSStefano Zampini   }
68751b0f6cfSStefano Zampini   /* compute local mirrors (if any) */
68851b0f6cfSStefano Zampini   if (mirrors_found) {
68951b0f6cfSStefano Zampini     PetscInt *local_indices,*global_indices;
69051b0f6cfSStefano Zampini     /* get arrays of local and global indices */
691785e854fSJed Brown     ierr = PetscMalloc1(graph->nvtxs,&local_indices);CHKERRQ(ierr);
69251b0f6cfSStefano Zampini     ierr = ISGetIndices(to,(const PetscInt**)&is_indices);CHKERRQ(ierr);
69351b0f6cfSStefano Zampini     ierr = PetscMemcpy(local_indices,is_indices,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
69451b0f6cfSStefano Zampini     ierr = ISRestoreIndices(to,(const PetscInt**)&is_indices);CHKERRQ(ierr);
695785e854fSJed Brown     ierr = PetscMalloc1(graph->nvtxs,&global_indices);CHKERRQ(ierr);
69651b0f6cfSStefano Zampini     ierr = ISGetIndices(from,(const PetscInt**)&is_indices);CHKERRQ(ierr);
69751b0f6cfSStefano Zampini     ierr = PetscMemcpy(global_indices,is_indices,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
69851b0f6cfSStefano Zampini     ierr = ISRestoreIndices(from,(const PetscInt**)&is_indices);CHKERRQ(ierr);
69951b0f6cfSStefano Zampini     /* allocate space for mirrors */
7008e5aaad5SJed Brown     ierr = PetscMalloc2(graph->nvtxs,&graph->mirrors,graph->nvtxs,&graph->mirrors_set);CHKERRQ(ierr);
70151b0f6cfSStefano Zampini     ierr = PetscMemzero(graph->mirrors,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
70251b0f6cfSStefano Zampini     graph->mirrors_set[0] = 0;
70351b0f6cfSStefano Zampini 
70451b0f6cfSStefano Zampini     k=0;
70551b0f6cfSStefano Zampini     for (i=0;i<n_shared[0];i++) {
70651b0f6cfSStefano Zampini       j=shared[0][i];
70751b0f6cfSStefano Zampini       if (graph->count[j] > 1) {
70851b0f6cfSStefano Zampini         graph->mirrors[j]++;
70951b0f6cfSStefano Zampini         k++;
71051b0f6cfSStefano Zampini       }
71151b0f6cfSStefano Zampini     }
71251b0f6cfSStefano Zampini     /* allocate space for set of mirrors */
713785e854fSJed Brown     ierr = PetscMalloc1(k,&graph->mirrors_set[0]);CHKERRQ(ierr);
71451b0f6cfSStefano Zampini     for (i=1;i<graph->nvtxs;i++)
71551b0f6cfSStefano Zampini       graph->mirrors_set[i]=graph->mirrors_set[i-1]+graph->mirrors[i-1];
71651b0f6cfSStefano Zampini 
71751b0f6cfSStefano Zampini     /* fill arrays */
71851b0f6cfSStefano Zampini     ierr = PetscMemzero(graph->mirrors,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
71951b0f6cfSStefano Zampini     for (j=0;j<n_shared[0];j++) {
72051b0f6cfSStefano Zampini       i=shared[0][j];
72151b0f6cfSStefano Zampini       if (graph->count[i] > 1)
72251b0f6cfSStefano Zampini         graph->mirrors_set[i][graph->mirrors[i]++]=global_indices[i];
72351b0f6cfSStefano Zampini     }
72451b0f6cfSStefano Zampini     ierr = PetscSortIntWithArray(graph->nvtxs,global_indices,local_indices);CHKERRQ(ierr);
72551b0f6cfSStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
72651b0f6cfSStefano Zampini       if (graph->mirrors[i] > 0) {
72751b0f6cfSStefano Zampini         ierr = PetscFindInt(graph->mirrors_set[i][0],graph->nvtxs,global_indices,&k);CHKERRQ(ierr);
72851b0f6cfSStefano Zampini         j = global_indices[k];
72951b0f6cfSStefano Zampini         while ( k > 0 && global_indices[k-1] == j) k--;
73051b0f6cfSStefano Zampini         for (j=0;j<graph->mirrors[i];j++) {
73151b0f6cfSStefano Zampini           graph->mirrors_set[i][j]=local_indices[k+j];
73251b0f6cfSStefano Zampini         }
73351b0f6cfSStefano Zampini         ierr = PetscSortInt(graph->mirrors[i],graph->mirrors_set[i]);CHKERRQ(ierr);
73451b0f6cfSStefano Zampini       }
73551b0f6cfSStefano Zampini     }
73651b0f6cfSStefano Zampini     ierr = PetscFree(local_indices);CHKERRQ(ierr);
73751b0f6cfSStefano Zampini     ierr = PetscFree(global_indices);CHKERRQ(ierr);
73851b0f6cfSStefano Zampini   }
73951b0f6cfSStefano Zampini   ierr = PetscMemzero(graph->count,graph->nvtxs*sizeof(*graph->count));CHKERRQ(ierr);
740674ae819SStefano Zampini   ierr = ISDestroy(&to);CHKERRQ(ierr);
741674ae819SStefano Zampini   ierr = ISDestroy(&from);CHKERRQ(ierr);
74251b0f6cfSStefano Zampini 
743674ae819SStefano Zampini   /* Count total number of neigh per node */
744674ae819SStefano Zampini   k=0;
745674ae819SStefano Zampini   for (i=1;i<n_neigh;i++) {
746674ae819SStefano Zampini     k += n_shared[i];
747674ae819SStefano Zampini     for (j=0;j<n_shared[i];j++) {
748674ae819SStefano Zampini       graph->count[shared[i][j]] += 1;
749674ae819SStefano Zampini     }
750674ae819SStefano Zampini   }
751674ae819SStefano Zampini   /* Allocate space for storing the set of neighbours for each node */
752674ae819SStefano Zampini   if (graph->nvtxs) {
753785e854fSJed Brown     ierr = PetscMalloc1(k,&graph->neighbours_set[0]);CHKERRQ(ierr);
754674ae819SStefano Zampini   }
755674ae819SStefano Zampini   for (i=1;i<graph->nvtxs;i++) { /* dont count myself */
756674ae819SStefano Zampini     graph->neighbours_set[i]=graph->neighbours_set[i-1]+graph->count[i-1];
757674ae819SStefano Zampini   }
758674ae819SStefano Zampini   /* Get information for sharing subdomains */
759674ae819SStefano Zampini   ierr = PetscMemzero(graph->count,graph->nvtxs*sizeof(*graph->count));CHKERRQ(ierr);
760674ae819SStefano Zampini   for (i=1;i<n_neigh;i++) { /* dont count myself */
761674ae819SStefano Zampini     s = n_shared[i];
762674ae819SStefano Zampini     for (j=0;j<s;j++) {
763674ae819SStefano Zampini       k = shared[i][j];
764674ae819SStefano Zampini       graph->neighbours_set[k][graph->count[k]] = neigh[i];
765674ae819SStefano Zampini       graph->count[k] += 1;
766674ae819SStefano Zampini     }
767674ae819SStefano Zampini   }
768674ae819SStefano Zampini   /* sort set of sharing subdomains */
769674ae819SStefano Zampini   for (i=0;i<graph->nvtxs;i++) {
77093fb5fd3SStefano Zampini     ierr = PetscSortRemoveDupsInt(&graph->count[i],graph->neighbours_set[i]);CHKERRQ(ierr);
771674ae819SStefano Zampini   }
7727fb0e2dbSStefano Zampini   /* free memory allocated by ISLocalToGlobalMappingGetInfo */
7737fb0e2dbSStefano Zampini   ierr = ISLocalToGlobalMappingRestoreInfo(graph->l2gmap,&n_neigh,&neigh,&n_shared,&shared);CHKERRQ(ierr);
7747fb0e2dbSStefano Zampini 
77567c9da69SStefano Zampini   /*
77667c9da69SStefano Zampini      Get info for dofs splitting
7775777c63cSStefano Zampini      User can specify just a subset; an additional field is considered as a complementary field
77867c9da69SStefano Zampini   */
77967c9da69SStefano Zampini   for (i=0;i<graph->nvtxs;i++) {
78067c9da69SStefano Zampini     graph->which_dof[i] = n_ISForDofs; /* by default a dof belongs to the complement set */
78167c9da69SStefano Zampini   }
7825777c63cSStefano Zampini   if (n_ISForDofs) {
7835777c63cSStefano Zampini     ierr = VecSet(local_vec,-1.0);CHKERRQ(ierr);
7845777c63cSStefano Zampini   }
785674ae819SStefano Zampini   for (i=0;i<n_ISForDofs;i++) {
7865777c63cSStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
78763602bcaSStefano Zampini     ierr = ISGetLocalSize(ISForDofs[i],&is_size);CHKERRQ(ierr);
788674ae819SStefano Zampini     ierr = ISGetIndices(ISForDofs[i],(const PetscInt**)&is_indices);CHKERRQ(ierr);
789674ae819SStefano Zampini     for (j=0;j<is_size;j++) {
790d50ed60aSStefano Zampini       if (is_indices[j] > -1 && is_indices[j] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
791674ae819SStefano Zampini         graph->which_dof[is_indices[j]] = i;
7925777c63cSStefano Zampini         array[is_indices[j]] = 1.*i;
793674ae819SStefano Zampini       }
79467c9da69SStefano Zampini     }
795674ae819SStefano Zampini     ierr = ISRestoreIndices(ISForDofs[i],(const PetscInt**)&is_indices);CHKERRQ(ierr);
7965777c63cSStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
7975777c63cSStefano Zampini   }
7985777c63cSStefano Zampini   /* Check consistency among neighbours */
7995777c63cSStefano Zampini   if (n_ISForDofs) {
8005777c63cSStefano Zampini     ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8015777c63cSStefano Zampini     ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
8025777c63cSStefano Zampini     ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8035777c63cSStefano Zampini     ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8045777c63cSStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
8055777c63cSStefano Zampini     ierr = VecGetArray(local_vec2,&array2);CHKERRQ(ierr);
8065777c63cSStefano Zampini     for (i=0;i<graph->nvtxs;i++){
8075777c63cSStefano Zampini       PetscInt field1,field2;
8085777c63cSStefano Zampini 
8095777c63cSStefano Zampini       field1 = (PetscInt)PetscRealPart(array[i]);
8105777c63cSStefano Zampini       field2 = (PetscInt)PetscRealPart(array2[i]);
8115777c63cSStefano Zampini       if (field1 != field2) {
8125777c63cSStefano Zampini         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);
8135777c63cSStefano Zampini       }
8145777c63cSStefano Zampini     }
8155777c63cSStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
8165777c63cSStefano Zampini     ierr = VecRestoreArray(local_vec2,&array2);CHKERRQ(ierr);
817674ae819SStefano Zampini   }
8187fb0e2dbSStefano Zampini   if (neumann_is || dirichlet_is) {
819674ae819SStefano Zampini     /* Take into account Neumann nodes */
820674ae819SStefano Zampini     ierr = VecSet(local_vec,0.0);CHKERRQ(ierr);
821674ae819SStefano Zampini     ierr = VecSet(local_vec2,0.0);CHKERRQ(ierr);
822674ae819SStefano Zampini     if (neumann_is) {
823674ae819SStefano Zampini       ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
82482ba6b80SStefano Zampini       ierr = ISGetLocalSize(neumann_is,&is_size);CHKERRQ(ierr);
825674ae819SStefano Zampini       ierr = ISGetIndices(neumann_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
826674ae819SStefano Zampini       for (i=0;i<is_size;i++) {
827d50ed60aSStefano Zampini         if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
828674ae819SStefano Zampini           array[is_indices[i]] = 1.0;
829674ae819SStefano Zampini         }
8303c629590SStefano Zampini       }
831674ae819SStefano Zampini       ierr = ISRestoreIndices(neumann_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
832674ae819SStefano Zampini       ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
833674ae819SStefano Zampini     }
834674ae819SStefano Zampini     /* Neumann nodes: impose consistency among neighbours */
835674ae819SStefano Zampini     ierr = VecSet(global_vec,0.0);CHKERRQ(ierr);
836674ae819SStefano Zampini     ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
837674ae819SStefano Zampini     ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
838674ae819SStefano Zampini     ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
839674ae819SStefano Zampini     ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
840674ae819SStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
841674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
8423c629590SStefano Zampini       if (PetscRealPart(array[i]) > 0.1) {
8430cf82fd6SStefano Zampini         graph->special_dof[i] = PCBDDCGRAPH_NEUMANN_MARK;
844674ae819SStefano Zampini       }
845674ae819SStefano Zampini     }
846674ae819SStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
847674ae819SStefano Zampini     /* Take into account Dirichlet nodes */
848674ae819SStefano Zampini     ierr = VecSet(local_vec2,0.0);CHKERRQ(ierr);
849674ae819SStefano Zampini     if (dirichlet_is) {
850674ae819SStefano Zampini       ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
851674ae819SStefano Zampini       ierr = VecGetArray(local_vec2,&array2);CHKERRQ(ierr);
85282ba6b80SStefano Zampini       ierr = ISGetLocalSize(dirichlet_is,&is_size);CHKERRQ(ierr);
853674ae819SStefano Zampini       ierr = ISGetIndices(dirichlet_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
854674ae819SStefano Zampini       for (i=0;i<is_size;i++){
855d50ed60aSStefano Zampini         if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
856674ae819SStefano Zampini           k = is_indices[i];
857df48933bSStefano Zampini           if (graph->count[k] && !PetscBTLookup(graph->touched,k)) {
8583c629590SStefano Zampini             if (PetscRealPart(array[k]) > 0.1) {
8595777c63cSStefano Zampini               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);
860674ae819SStefano Zampini             }
861674ae819SStefano Zampini             array2[k] = 1.0;
862674ae819SStefano Zampini           }
863674ae819SStefano Zampini         }
8643c629590SStefano Zampini       }
865674ae819SStefano Zampini       ierr = ISRestoreIndices(dirichlet_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
866674ae819SStefano Zampini       ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
867674ae819SStefano Zampini       ierr = VecRestoreArray(local_vec2,&array2);CHKERRQ(ierr);
868674ae819SStefano Zampini     }
869674ae819SStefano Zampini     /* Dirichlet nodes: impose consistency among neighbours */
870674ae819SStefano Zampini     ierr = VecSet(global_vec,0.0);CHKERRQ(ierr);
871674ae819SStefano Zampini     ierr = VecScatterBegin(scatter_ctx,local_vec2,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
872674ae819SStefano Zampini     ierr = VecScatterEnd(scatter_ctx,local_vec2,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
873674ae819SStefano Zampini     ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
874674ae819SStefano Zampini     ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
875674ae819SStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
876674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
8773c629590SStefano Zampini       if (PetscRealPart(array[i]) > 0.1) {
878df48933bSStefano Zampini         ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
879674ae819SStefano Zampini         graph->subset[i] = 0; /* dirichlet nodes treated as internal -> is it ok? */
8800cf82fd6SStefano Zampini         graph->special_dof[i] = PCBDDCGRAPH_DIRICHLET_MARK;
881674ae819SStefano Zampini       }
882674ae819SStefano Zampini     }
883674ae819SStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
8847fb0e2dbSStefano Zampini   }
88508a5cf49SStefano Zampini   /* mark local periodic nodes (if any) and adapt CSR graph (if any) */
88651b0f6cfSStefano Zampini   if (graph->mirrors) {
88751b0f6cfSStefano Zampini     for (i=0;i<graph->nvtxs;i++)
88851b0f6cfSStefano Zampini       if (graph->mirrors[i])
8890cf82fd6SStefano Zampini         graph->special_dof[i] = PCBDDCGRAPH_LOCAL_PERIODIC_MARK;
89051b0f6cfSStefano Zampini 
89108a5cf49SStefano Zampini     if (graph->xadj && graph->adjncy) {
89208a5cf49SStefano Zampini       PetscInt *new_xadj,*new_adjncy;
89351b0f6cfSStefano Zampini       /* sort CSR graph */
89451b0f6cfSStefano Zampini       for (i=0;i<graph->nvtxs;i++)
89551b0f6cfSStefano Zampini         ierr = PetscSortInt(graph->xadj[i+1]-graph->xadj[i],&graph->adjncy[graph->xadj[i]]);CHKERRQ(ierr);
89651b0f6cfSStefano Zampini 
89751b0f6cfSStefano Zampini       /* adapt local CSR graph in case of local periodicity */
89851b0f6cfSStefano Zampini       k=0;
89951b0f6cfSStefano Zampini       for (i=0;i<graph->nvtxs;i++)
90051b0f6cfSStefano Zampini         for (j=graph->xadj[i];j<graph->xadj[i+1];j++)
90151b0f6cfSStefano Zampini           k += graph->mirrors[graph->adjncy[j]];
90251b0f6cfSStefano Zampini 
903854ce69bSBarry Smith       ierr = PetscMalloc1(graph->nvtxs+1,&new_xadj);CHKERRQ(ierr);
904854ce69bSBarry Smith       ierr = PetscMalloc1(k+graph->xadj[graph->nvtxs],&new_adjncy);CHKERRQ(ierr);
90551b0f6cfSStefano Zampini       new_xadj[0]=0;
90651b0f6cfSStefano Zampini       for (i=0;i<graph->nvtxs;i++) {
90751b0f6cfSStefano Zampini         k = graph->xadj[i+1]-graph->xadj[i];
90851b0f6cfSStefano Zampini         ierr = PetscMemcpy(&new_adjncy[new_xadj[i]],&graph->adjncy[graph->xadj[i]],k*sizeof(PetscInt));CHKERRQ(ierr);
90951b0f6cfSStefano Zampini         new_xadj[i+1]=new_xadj[i]+k;
91051b0f6cfSStefano Zampini         for (j=graph->xadj[i];j<graph->xadj[i+1];j++) {
91151b0f6cfSStefano Zampini           k = graph->mirrors[graph->adjncy[j]];
91251b0f6cfSStefano Zampini           ierr = PetscMemcpy(&new_adjncy[new_xadj[i+1]],graph->mirrors_set[graph->adjncy[j]],k*sizeof(PetscInt));CHKERRQ(ierr);
91351b0f6cfSStefano Zampini           new_xadj[i+1]+=k;
91451b0f6cfSStefano Zampini         }
91551b0f6cfSStefano Zampini         k = new_xadj[i+1]-new_xadj[i];
91651b0f6cfSStefano Zampini         ierr = PetscSortRemoveDupsInt(&k,&new_adjncy[new_xadj[i]]);CHKERRQ(ierr);
91751b0f6cfSStefano Zampini         new_xadj[i+1]=new_xadj[i]+k;
91851b0f6cfSStefano Zampini       }
91951b0f6cfSStefano Zampini       /* set new CSR into graph */
92051b0f6cfSStefano Zampini       ierr = PetscFree(graph->xadj);CHKERRQ(ierr);
92151b0f6cfSStefano Zampini       ierr = PetscFree(graph->adjncy);CHKERRQ(ierr);
92251b0f6cfSStefano Zampini       graph->xadj = new_xadj;
92351b0f6cfSStefano Zampini       graph->adjncy = new_adjncy;
92451b0f6cfSStefano Zampini     }
92508a5cf49SStefano Zampini   }
92651b0f6cfSStefano Zampini 
92717eb1463SStefano Zampini   /* mark special nodes (if any) -> each will become a single node equivalence class */
928674ae819SStefano Zampini   if (custom_primal_vertices) {
92917eb1463SStefano Zampini     ierr = VecSet(local_vec,0.0);CHKERRQ(ierr);
9309b70a373SStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
93163602bcaSStefano Zampini     ierr = ISGetLocalSize(custom_primal_vertices,&is_size);CHKERRQ(ierr);
932674ae819SStefano Zampini     ierr = ISGetIndices(custom_primal_vertices,(const PetscInt**)&is_indices);CHKERRQ(ierr);
933674ae819SStefano Zampini     for (i=0;i<is_size;i++){
9349b70a373SStefano Zampini       if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
9359b70a373SStefano Zampini         array[is_indices[i]] = 1.0;
9369b70a373SStefano Zampini       }
937674ae819SStefano Zampini     }
938674ae819SStefano Zampini     ierr = ISRestoreIndices(custom_primal_vertices,(const PetscInt**)&is_indices);CHKERRQ(ierr);
9399b70a373SStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
9409b70a373SStefano Zampini     /* special nodes: impose consistency among neighbours */
9419b70a373SStefano Zampini     ierr = VecSet(global_vec,0.0);CHKERRQ(ierr);
9429b70a373SStefano Zampini     ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9439b70a373SStefano Zampini     ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9449b70a373SStefano Zampini     ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9459b70a373SStefano Zampini     ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9469b70a373SStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
9479b70a373SStefano Zampini     j = 0;
9489b70a373SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
9499b70a373SStefano Zampini       if (PetscRealPart(array[i]) > 0.1 && graph->special_dof[i] != PCBDDCGRAPH_DIRICHLET_MARK) {
9509b70a373SStefano Zampini         graph->special_dof[i] = PCBDDCGRAPH_SPECIAL_MARK-j;
9519b70a373SStefano Zampini         j++;
9529b70a373SStefano Zampini       }
9539b70a373SStefano Zampini     }
9549b70a373SStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
95517eb1463SStefano Zampini   }
9569b70a373SStefano Zampini 
957674ae819SStefano Zampini   /* mark interior nodes as touched and belonging to partition number 0 */
958674ae819SStefano Zampini   for (i=0;i<graph->nvtxs;i++) {
959674ae819SStefano Zampini     if (!graph->count[i]) {
960df48933bSStefano Zampini       ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
961674ae819SStefano Zampini       graph->subset[i] = 0;
962674ae819SStefano Zampini     }
963674ae819SStefano Zampini   }
964674ae819SStefano Zampini   /* init graph structure and compute default subsets */
965674ae819SStefano Zampini   nodes_touched=0;
966674ae819SStefano Zampini   for (i=0;i<graph->nvtxs;i++) {
967df48933bSStefano Zampini     if (PetscBTLookup(graph->touched,i)) {
968674ae819SStefano Zampini       nodes_touched++;
969674ae819SStefano Zampini     }
970674ae819SStefano Zampini   }
971674ae819SStefano Zampini   i = 0;
972674ae819SStefano Zampini   graph->ncc = 0;
973674ae819SStefano Zampini   total_counts = 0;
9747babac9bSStefano Zampini 
9757babac9bSStefano Zampini   /* allocated space for queues */
9767babac9bSStefano Zampini   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
9777babac9bSStefano Zampini   if (size == 1) {
9787babac9bSStefano Zampini     ierr = PetscMalloc2(graph->nvtxs+1,&graph->cptr,graph->nvtxs,&graph->queue);CHKERRQ(ierr);
9797babac9bSStefano Zampini   } else {
9807babac9bSStefano Zampini     PetscInt nused = graph->nvtxs - nodes_touched;
9817babac9bSStefano Zampini     ierr = PetscMalloc2(nused+1,&graph->cptr,nused,&graph->queue);CHKERRQ(ierr);
9827babac9bSStefano Zampini   }
9837babac9bSStefano Zampini 
984674ae819SStefano Zampini   while (nodes_touched<graph->nvtxs) {
985674ae819SStefano Zampini     /*  find first untouched node in local ordering */
986df48933bSStefano Zampini     while (PetscBTLookup(graph->touched,i)) {
987674ae819SStefano Zampini       i++;
988674ae819SStefano Zampini     }
989df48933bSStefano Zampini     ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
990674ae819SStefano Zampini     graph->subset[i] = graph->ncc+1;
991674ae819SStefano Zampini     graph->cptr[graph->ncc] = total_counts;
992674ae819SStefano Zampini     graph->queue[total_counts] = i;
993674ae819SStefano Zampini     total_counts++;
994674ae819SStefano Zampini     nodes_touched++;
995674ae819SStefano Zampini     /* now find all other nodes having the same set of sharing subdomains */
996674ae819SStefano Zampini     for (j=i+1;j<graph->nvtxs;j++) {
99774e413f5SStefano Zampini       /* check for same number of sharing subdomains, dof number and same special mark */
998df48933bSStefano 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]) {
999674ae819SStefano Zampini         /* check for same set of sharing subdomains */
1000674ae819SStefano Zampini         same_set=PETSC_TRUE;
1001674ae819SStefano Zampini         for (k=0;k<graph->count[j];k++){
1002674ae819SStefano Zampini           if (graph->neighbours_set[i][k]!=graph->neighbours_set[j][k]) {
1003674ae819SStefano Zampini             same_set=PETSC_FALSE;
1004674ae819SStefano Zampini           }
1005674ae819SStefano Zampini         }
1006674ae819SStefano Zampini         /* I found a friend of mine */
1007674ae819SStefano Zampini         if (same_set) {
1008674ae819SStefano Zampini           graph->subset[j]=graph->ncc+1;
1009df48933bSStefano Zampini           ierr = PetscBTSet(graph->touched,j);CHKERRQ(ierr);
1010674ae819SStefano Zampini           nodes_touched++;
1011674ae819SStefano Zampini           graph->queue[total_counts] = j;
1012674ae819SStefano Zampini           total_counts++;
1013674ae819SStefano Zampini         }
1014674ae819SStefano Zampini       }
1015674ae819SStefano Zampini     }
1016674ae819SStefano Zampini     graph->ncc++;
1017674ae819SStefano Zampini   }
1018674ae819SStefano Zampini   /* set default number of subsets (at this point no info on csr graph has been taken into account, so n_subsets = ncc */
1019674ae819SStefano Zampini   graph->n_subsets = graph->ncc;
1020785e854fSJed Brown   ierr = PetscMalloc1(graph->n_subsets,&graph->subset_ncc);CHKERRQ(ierr);
1021674ae819SStefano Zampini   for (i=0;i<graph->n_subsets;i++) {
1022674ae819SStefano Zampini     graph->subset_ncc[i] = 1;
1023674ae819SStefano Zampini   }
1024674ae819SStefano Zampini   /* final pointer */
1025674ae819SStefano Zampini   graph->cptr[graph->ncc] = total_counts;
1026ec1c413dSStefano Zampini 
1027ec1c413dSStefano Zampini   /* save information on subsets (needed if we have to adapt the connected components later) */
1028ec1c413dSStefano Zampini   /* For consistency reasons (among neighbours), I need to sort (by global ordering) each subset */
1029ec1c413dSStefano Zampini   /* Get a reference node (min index in global ordering) for each subset for tagging messages */
1030ec1c413dSStefano Zampini   ierr = PetscMalloc2(graph->ncc,&subset_ref_node_global,graph->cptr[graph->ncc],&queue_global);CHKERRQ(ierr);
10313a5b03d1SStefano Zampini   ierr = ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_global);CHKERRQ(ierr);
1032ec1c413dSStefano Zampini   for (j=0;j<graph->ncc;j++) {
1033ec1c413dSStefano Zampini     ierr = PetscSortIntWithArray(graph->cptr[j+1]-graph->cptr[j],&queue_global[graph->cptr[j]],&graph->queue[graph->cptr[j]]);CHKERRQ(ierr);
1034ec1c413dSStefano Zampini     subset_ref_node_global[j] = graph->queue[graph->cptr[j]];
1035f0321c90SStefano Zampini   }
1036acc113dbSStefano Zampini   graph->queue_sorted = PETSC_TRUE;
10372f566a09SStefano Zampini   if (graph->ncc) {
1038ec1c413dSStefano Zampini     ierr = PetscMalloc2(graph->ncc,&graph->subsets_size,graph->ncc,&graph->subsets);CHKERRQ(ierr);
1039ec1c413dSStefano Zampini     ierr = PetscMalloc1(graph->cptr[graph->ncc],&graph->subsets[0]);CHKERRQ(ierr);
1040ec1c413dSStefano Zampini     ierr = PetscMemzero(graph->subsets[0],graph->cptr[graph->ncc]*sizeof(PetscInt));CHKERRQ(ierr);
1041ec1c413dSStefano Zampini     for (j=1;j<graph->ncc;j++) {
1042ec1c413dSStefano Zampini       graph->subsets_size[j-1] = graph->cptr[j] - graph->cptr[j-1];
1043ec1c413dSStefano Zampini       graph->subsets[j] = graph->subsets[j-1] + graph->subsets_size[j-1];
1044ec1c413dSStefano Zampini     }
1045ec1c413dSStefano Zampini     graph->subsets_size[graph->ncc-1] = graph->cptr[graph->ncc] - graph->cptr[graph->ncc-1];
1046ec1c413dSStefano Zampini     for (j=0;j<graph->ncc;j++) {
1047ec1c413dSStefano Zampini       ierr = PetscMemcpy(graph->subsets[j],&graph->queue[graph->cptr[j]],graph->subsets_size[j]*sizeof(PetscInt));CHKERRQ(ierr);
1048ec1c413dSStefano Zampini     }
10492f566a09SStefano Zampini   }
1050f0321c90SStefano Zampini   /* renumber reference nodes */
1051f0321c90SStefano Zampini   ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)(graph->l2gmap)),graph->l2gmap,graph->ncc,subset_ref_node_global,NULL,&k,&graph->subset_ref_node);CHKERRQ(ierr);
1052ec1c413dSStefano Zampini 
1053ec1c413dSStefano Zampini   /* free workspace */
1054ec1c413dSStefano Zampini   ierr = PetscFree2(subset_ref_node_global,queue_global);CHKERRQ(ierr);
1055674ae819SStefano Zampini   ierr = VecDestroy(&local_vec);CHKERRQ(ierr);
1056674ae819SStefano Zampini   ierr = VecDestroy(&local_vec2);CHKERRQ(ierr);
1057674ae819SStefano Zampini   ierr = VecDestroy(&global_vec);CHKERRQ(ierr);
1058674ae819SStefano Zampini   ierr = VecScatterDestroy(&scatter_ctx);CHKERRQ(ierr);
1059674ae819SStefano Zampini   PetscFunctionReturn(0);
1060674ae819SStefano Zampini }
1061674ae819SStefano Zampini 
1062674ae819SStefano Zampini #undef __FUNCT__
1063674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphResetCSR"
1064674ae819SStefano Zampini PetscErrorCode PCBDDCGraphResetCSR(PCBDDCGraph graph)
1065674ae819SStefano Zampini {
1066674ae819SStefano Zampini   PetscErrorCode ierr;
1067674ae819SStefano Zampini 
1068674ae819SStefano Zampini   PetscFunctionBegin;
1069674ae819SStefano Zampini   ierr = PetscFree(graph->xadj);CHKERRQ(ierr);
1070674ae819SStefano Zampini   ierr = PetscFree(graph->adjncy);CHKERRQ(ierr);
1071575ad6abSStefano Zampini   graph->nvtxs_csr = 0;
1072674ae819SStefano Zampini   PetscFunctionReturn(0);
1073674ae819SStefano Zampini }
1074674ae819SStefano Zampini 
1075674ae819SStefano Zampini #undef __FUNCT__
1076674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphReset"
1077674ae819SStefano Zampini PetscErrorCode PCBDDCGraphReset(PCBDDCGraph graph)
1078674ae819SStefano Zampini {
1079674ae819SStefano Zampini   PetscErrorCode ierr;
1080674ae819SStefano Zampini 
1081674ae819SStefano Zampini   PetscFunctionBegin;
1082674ae819SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&graph->l2gmap);CHKERRQ(ierr);
1083674ae819SStefano Zampini   ierr = PetscFree(graph->subset_ncc);CHKERRQ(ierr);
10843a5b03d1SStefano Zampini   ierr = PetscFree(graph->subset_ref_node);CHKERRQ(ierr);
1085674ae819SStefano Zampini   if (graph->nvtxs) {
1086674ae819SStefano Zampini     ierr = PetscFree(graph->neighbours_set[0]);CHKERRQ(ierr);
1087674ae819SStefano Zampini   }
1088df48933bSStefano Zampini   ierr = PetscBTDestroy(&graph->touched);CHKERRQ(ierr);
10897babac9bSStefano Zampini   ierr = PetscFree5(graph->count,
1090674ae819SStefano Zampini                     graph->neighbours_set,
1091674ae819SStefano Zampini                     graph->subset,
1092674ae819SStefano Zampini                     graph->which_dof,
1093df48933bSStefano Zampini                     graph->special_dof);CHKERRQ(ierr);
10947babac9bSStefano Zampini   ierr = PetscFree2(graph->cptr,graph->queue);CHKERRQ(ierr);
109551b0f6cfSStefano Zampini   if (graph->mirrors) {
109651b0f6cfSStefano Zampini     ierr = PetscFree(graph->mirrors_set[0]);CHKERRQ(ierr);
109751b0f6cfSStefano Zampini   }
109851b0f6cfSStefano Zampini   ierr = PetscFree2(graph->mirrors,graph->mirrors_set);CHKERRQ(ierr);
1099ec1c413dSStefano Zampini   if (graph->subsets) {
1100ec1c413dSStefano Zampini     ierr = PetscFree(graph->subsets[0]);CHKERRQ(ierr);
1101ec1c413dSStefano Zampini   }
1102ec1c413dSStefano Zampini   ierr = PetscFree2(graph->subsets_size,graph->subsets);CHKERRQ(ierr);
1103674ae819SStefano Zampini   graph->nvtxs = 0;
11047babac9bSStefano Zampini   graph->nvtxs_global = 0;
1105674ae819SStefano Zampini   graph->n_subsets = 0;
1106674ae819SStefano Zampini   graph->custom_minimal_size = 1;
1107674ae819SStefano Zampini   PetscFunctionReturn(0);
1108674ae819SStefano Zampini }
1109674ae819SStefano Zampini 
1110674ae819SStefano Zampini #undef __FUNCT__
1111674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphInit"
11127fb0e2dbSStefano Zampini PetscErrorCode PCBDDCGraphInit(PCBDDCGraph graph, ISLocalToGlobalMapping l2gmap, PetscInt N)
1113674ae819SStefano Zampini {
1114674ae819SStefano Zampini   PetscInt       n;
1115674ae819SStefano Zampini   PetscErrorCode ierr;
1116674ae819SStefano Zampini 
1117674ae819SStefano Zampini   PetscFunctionBegin;
1118674ae819SStefano Zampini   PetscValidPointer(graph,1);
1119674ae819SStefano Zampini   PetscValidHeaderSpecific(l2gmap,IS_LTOGM_CLASSID,2);
11207babac9bSStefano Zampini   PetscValidLogicalCollectiveInt(l2gmap,N,3);
11218e61c736SStefano Zampini   /* raise an error if already allocated */
11227babac9bSStefano Zampini   if (graph->nvtxs_global) {
11238e61c736SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)l2gmap),PETSC_ERR_PLIB,"BDDCGraph already initialized");
1124674ae819SStefano Zampini   }
1125674ae819SStefano Zampini   /* set number of vertices */
1126674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)l2gmap);CHKERRQ(ierr);
1127674ae819SStefano Zampini   graph->l2gmap = l2gmap;
1128674ae819SStefano Zampini   ierr = ISLocalToGlobalMappingGetSize(l2gmap,&n);CHKERRQ(ierr);
1129674ae819SStefano Zampini   graph->nvtxs = n;
11307fb0e2dbSStefano Zampini   graph->nvtxs_global = N;
1131674ae819SStefano Zampini   /* allocate used space */
1132df48933bSStefano Zampini   ierr = PetscBTCreate(graph->nvtxs,&graph->touched);CHKERRQ(ierr);
11337babac9bSStefano Zampini   ierr = PetscMalloc5(graph->nvtxs,&graph->count,
11348e5aaad5SJed Brown                       graph->nvtxs,&graph->neighbours_set,
11358e5aaad5SJed Brown                       graph->nvtxs,&graph->subset,
11368e5aaad5SJed Brown                       graph->nvtxs,&graph->which_dof,
11378e5aaad5SJed Brown                       graph->nvtxs,&graph->special_dof);CHKERRQ(ierr);
1138674ae819SStefano Zampini   /* zeroes memory */
1139674ae819SStefano Zampini   ierr = PetscMemzero(graph->count,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
1140674ae819SStefano Zampini   ierr = PetscMemzero(graph->subset,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
114163602bcaSStefano Zampini   /* use -1 as a default value for which_dof array */
114263602bcaSStefano Zampini   for (n=0;n<graph->nvtxs;n++) graph->which_dof[n] = -1;
114374e413f5SStefano Zampini   ierr = PetscMemzero(graph->special_dof,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
1144674ae819SStefano Zampini   /* zeroes first pointer to neighbour set */
1145674ae819SStefano Zampini   if (graph->nvtxs) {
1146674ae819SStefano Zampini     graph->neighbours_set[0] = 0;
1147674ae819SStefano Zampini   }
1148674ae819SStefano Zampini   /* zeroes workspace for values of ncc */
1149674ae819SStefano Zampini   graph->subset_ncc = 0;
11503a5b03d1SStefano Zampini   graph->subset_ref_node = 0;
1151674ae819SStefano Zampini   PetscFunctionReturn(0);
1152674ae819SStefano Zampini }
1153674ae819SStefano Zampini 
1154674ae819SStefano Zampini #undef __FUNCT__
1155674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphDestroy"
1156674ae819SStefano Zampini PetscErrorCode PCBDDCGraphDestroy(PCBDDCGraph* graph)
1157674ae819SStefano Zampini {
1158674ae819SStefano Zampini   PetscErrorCode ierr;
1159674ae819SStefano Zampini 
1160674ae819SStefano Zampini   PetscFunctionBegin;
1161674ae819SStefano Zampini   ierr = PCBDDCGraphReset(*graph);CHKERRQ(ierr);
1162674ae819SStefano Zampini   ierr = PetscFree(*graph);CHKERRQ(ierr);
1163674ae819SStefano Zampini   PetscFunctionReturn(0);
1164674ae819SStefano Zampini }
1165674ae819SStefano Zampini 
1166674ae819SStefano Zampini #undef __FUNCT__
1167674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphCreate"
1168674ae819SStefano Zampini PetscErrorCode PCBDDCGraphCreate(PCBDDCGraph *graph)
1169674ae819SStefano Zampini {
1170674ae819SStefano Zampini   PCBDDCGraph    new_graph;
1171674ae819SStefano Zampini   PetscErrorCode ierr;
1172674ae819SStefano Zampini 
1173674ae819SStefano Zampini   PetscFunctionBegin;
1174854ce69bSBarry Smith   ierr = PetscNew(&new_graph);CHKERRQ(ierr);
1175674ae819SStefano Zampini   new_graph->custom_minimal_size = 1;
1176674ae819SStefano Zampini   *graph = new_graph;
1177674ae819SStefano Zampini   PetscFunctionReturn(0);
1178674ae819SStefano Zampini }
1179