xref: /petsc/src/ksp/pc/impls/bddc/bddcgraph.c (revision a873d5d0a1a3c1e54cac90d59b4ecebe1ebedcd4)
1674ae819SStefano Zampini #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__
6674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphASCIIView"
7e49050b4SStefano Zampini PetscErrorCode PCBDDCGraphASCIIView(PCBDDCGraph graph, PetscInt verbosity_level, PetscViewer viewer)
8674ae819SStefano Zampini {
92b510759SStefano Zampini   PetscInt       i,j,tabs;
1093fb5fd3SStefano Zampini   PetscInt*      queue_in_global_numbering;
11674ae819SStefano Zampini   PetscErrorCode ierr;
12674ae819SStefano Zampini 
13674ae819SStefano Zampini   PetscFunctionBegin;
14674ae819SStefano Zampini   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
152b510759SStefano Zampini   ierr = PetscViewerASCIIGetTab(viewer,&tabs);CHKERRQ(ierr);
16674ae819SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"--------------------------------------------------\n");CHKERRQ(ierr);
17674ae819SStefano Zampini   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
18674ae819SStefano Zampini   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Local BDDC graph for subdomain %04d\n",PetscGlobalRank);CHKERRQ(ierr);
19674ae819SStefano Zampini   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Number of vertices %d\n",graph->nvtxs);CHKERRQ(ierr);
20e49050b4SStefano Zampini   if (verbosity_level > 1) {
21674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
22674ae819SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"%d:\n",i);CHKERRQ(ierr);
23674ae819SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   which_dof: %d\n",graph->which_dof[i]);CHKERRQ(ierr);
2474e413f5SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   special_dof: %d\n",graph->special_dof[i]);CHKERRQ(ierr);
25674ae819SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   neighbours: %d\n",graph->count[i]);CHKERRQ(ierr);
262b510759SStefano Zampini       ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
27674ae819SStefano Zampini       if (graph->count[i]) {
28674ae819SStefano Zampini         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"     set of neighbours:");CHKERRQ(ierr);
29674ae819SStefano Zampini         for (j=0;j<graph->count[i];j++) {
30674ae819SStefano Zampini           ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d",graph->neighbours_set[i][j]);CHKERRQ(ierr);
31674ae819SStefano Zampini         }
32674ae819SStefano Zampini         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
33674ae819SStefano Zampini       }
342b510759SStefano Zampini       ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr);
352b510759SStefano Zampini       ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
3651b0f6cfSStefano Zampini       if (graph->mirrors) {
3751b0f6cfSStefano Zampini         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   mirrors: %d\n",graph->mirrors[i]);CHKERRQ(ierr);
3851b0f6cfSStefano Zampini         if (graph->mirrors[i]) {
392b510759SStefano Zampini           ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
4051b0f6cfSStefano Zampini           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"     set of mirrors:");CHKERRQ(ierr);
4151b0f6cfSStefano Zampini           for (j=0;j<graph->mirrors[i];j++) {
4251b0f6cfSStefano Zampini             ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d",graph->mirrors_set[i][j]);CHKERRQ(ierr);
4351b0f6cfSStefano Zampini           }
4451b0f6cfSStefano Zampini           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
452b510759SStefano Zampini           ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr);
462b510759SStefano Zampini           ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
4751b0f6cfSStefano Zampini         }
4851b0f6cfSStefano Zampini       }
49e49050b4SStefano Zampini       if (verbosity_level > 2) {
50674ae819SStefano Zampini         if (graph->xadj && graph->adjncy) {
51674ae819SStefano Zampini           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   local adj list:");CHKERRQ(ierr);
522b510759SStefano Zampini           ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
53674ae819SStefano Zampini           for (j=graph->xadj[i];j<graph->xadj[i+1];j++) {
54674ae819SStefano Zampini             ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d",graph->adjncy[j]);CHKERRQ(ierr);
55674ae819SStefano Zampini           }
56674ae819SStefano Zampini           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
572b510759SStefano Zampini           ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr);
582b510759SStefano Zampini           ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
59674ae819SStefano Zampini         }
60e49050b4SStefano Zampini       }
61674ae819SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   interface subset id: %d\n",graph->subset[i]);CHKERRQ(ierr);
62674ae819SStefano Zampini       if (graph->subset[i] && graph->subset_ncc) {
63674ae819SStefano Zampini         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"   ncc for subset: %d\n",graph->subset_ncc[graph->subset[i]-1]);CHKERRQ(ierr);
64674ae819SStefano Zampini       }
65674ae819SStefano Zampini     }
66e49050b4SStefano Zampini   }
67674ae819SStefano Zampini   ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Total number of connected components %d\n",graph->ncc);CHKERRQ(ierr);
68785e854fSJed Brown   ierr = PetscMalloc1(graph->cptr[graph->ncc],&queue_in_global_numbering);CHKERRQ(ierr);
6993fb5fd3SStefano Zampini   ierr = ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_in_global_numbering);CHKERRQ(ierr);
70674ae819SStefano Zampini   for (i=0;i<graph->ncc;i++) {
711ce3d396SStefano Zampini     PetscInt node_num=graph->queue[graph->cptr[i]];
721ce3d396SStefano Zampini     PetscBool printcc = PETSC_FALSE;
73674ae819SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  %d (neighs:",i);CHKERRQ(ierr);
742b510759SStefano Zampini     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
75674ae819SStefano Zampini     for (j=0;j<graph->count[node_num];j++) {
76674ae819SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d",graph->neighbours_set[node_num][j]);CHKERRQ(ierr);
77674ae819SStefano Zampini     }
78674ae819SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"):");CHKERRQ(ierr);
791ce3d396SStefano Zampini     if (verbosity_level > 1) {
801ce3d396SStefano Zampini       printcc = PETSC_TRUE;
81e635b14bSStefano Zampini     } else if (graph->count[node_num] > 1 || (graph->count[node_num] == 1 && graph->special_dof[node_num] == PCBDDCGRAPH_NEUMANN_MARK)) {
821ce3d396SStefano Zampini       printcc = PETSC_TRUE;
831ce3d396SStefano Zampini     }
841ce3d396SStefano Zampini     if (printcc) {
85674ae819SStefano Zampini       for (j=graph->cptr[i];j<graph->cptr[i+1];j++) {
8693fb5fd3SStefano Zampini         ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d (%d)",graph->queue[j],queue_in_global_numbering[j]);CHKERRQ(ierr);
87674ae819SStefano Zampini       }
881ce3d396SStefano Zampini     }
89674ae819SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
902b510759SStefano Zampini     ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr);
912b510759SStefano Zampini     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
92674ae819SStefano Zampini   }
9393fb5fd3SStefano Zampini   ierr = PetscFree(queue_in_global_numbering);CHKERRQ(ierr);
94e49050b4SStefano Zampini   if (graph->custom_minimal_size > 1 && verbosity_level > 1) {
95674ae819SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Custom minimal size %d\n",graph->custom_minimal_size);CHKERRQ(ierr);
96674ae819SStefano Zampini   }
97674ae819SStefano Zampini   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
98674ae819SStefano Zampini   PetscFunctionReturn(0);
99674ae819SStefano Zampini }
100674ae819SStefano Zampini 
101674ae819SStefano Zampini #undef __FUNCT__
102674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphGetCandidatesIS"
103*a873d5d0SStefano Zampini PetscErrorCode PCBDDCGraphGetCandidatesIS(PCBDDCGraph graph, PetscInt *n_faces, IS *FacesIS[], PetscInt *n_edges, IS *EdgesIS[], IS *VerticesIS)
104674ae819SStefano Zampini {
105674ae819SStefano Zampini   IS             *ISForFaces,*ISForEdges,ISForVertices;
106674ae819SStefano Zampini   PetscInt       i,j,nfc,nec,nvc,*idx;
107674ae819SStefano Zampini   PetscBool      twodim_flag;
108674ae819SStefano Zampini   PetscErrorCode ierr;
109674ae819SStefano Zampini 
110674ae819SStefano Zampini   PetscFunctionBegin;
111674ae819SStefano Zampini   /* loop on ccs to evalute number of faces, edges and vertices */
112674ae819SStefano Zampini   nfc = 0;
113674ae819SStefano Zampini   nec = 0;
114674ae819SStefano Zampini   nvc = 0;
115674ae819SStefano Zampini   twodim_flag = PETSC_FALSE;
116674ae819SStefano Zampini   for (i=0;i<graph->ncc;i++) {
1171e1f2d93SStefano Zampini     PetscInt repdof = graph->queue[graph->cptr[i]];
118674ae819SStefano Zampini     if (graph->cptr[i+1]-graph->cptr[i] > graph->custom_minimal_size) {
119e635b14bSStefano Zampini       if (graph->count[repdof] == 1 && graph->special_dof[repdof] != PCBDDCGRAPH_NEUMANN_MARK) {
120674ae819SStefano Zampini         nfc++;
121674ae819SStefano Zampini       } else { /* note that nec will be zero in 2d */
122674ae819SStefano Zampini         nec++;
123674ae819SStefano Zampini       }
124674ae819SStefano Zampini     } else {
125674ae819SStefano Zampini       nvc += graph->cptr[i+1]-graph->cptr[i];
126674ae819SStefano Zampini     }
127674ae819SStefano Zampini   }
128aa9fcddfSStefano Zampini   j=0;
129aa9fcddfSStefano Zampini   ierr = MPI_Allreduce(&nec,&j,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)graph->l2gmap));CHKERRQ(ierr);
130aa9fcddfSStefano Zampini   if (!j) { /* we are in a 2D case -> no faces, only edges */
131674ae819SStefano Zampini     nec = nfc;
132674ae819SStefano Zampini     nfc = 0;
133674ae819SStefano Zampini     twodim_flag = PETSC_TRUE;
134674ae819SStefano Zampini   }
135674ae819SStefano Zampini   /* allocate IS arrays for faces, edges. Vertices need a single index set. */
136785e854fSJed Brown   ierr = PetscMalloc1(nfc,&ISForFaces);CHKERRQ(ierr);
137785e854fSJed Brown   ierr = PetscMalloc1(nec,&ISForEdges);CHKERRQ(ierr);
138785e854fSJed Brown   ierr = PetscMalloc1(nvc,&idx);CHKERRQ(ierr);
139674ae819SStefano Zampini   /* loop on ccs to compute index sets for faces and edges */
140674ae819SStefano Zampini   nfc = 0;
141674ae819SStefano Zampini   nec = 0;
142674ae819SStefano Zampini   for (i=0;i<graph->ncc;i++) {
1431e1f2d93SStefano Zampini     PetscInt repdof = graph->queue[graph->cptr[i]];
144674ae819SStefano Zampini     if (graph->cptr[i+1]-graph->cptr[i] > graph->custom_minimal_size) {
145e635b14bSStefano Zampini       if (graph->count[repdof] == 1 && graph->special_dof[repdof] != PCBDDCGRAPH_NEUMANN_MARK) {
146674ae819SStefano Zampini         if (twodim_flag) {
147674ae819SStefano 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);
148674ae819SStefano Zampini           nec++;
149674ae819SStefano Zampini         } else {
150674ae819SStefano 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);
151674ae819SStefano Zampini           nfc++;
152674ae819SStefano Zampini         }
153674ae819SStefano Zampini       } else {
154674ae819SStefano 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);
155674ae819SStefano Zampini         nec++;
156674ae819SStefano Zampini       }
157674ae819SStefano Zampini     }
158674ae819SStefano Zampini   }
159674ae819SStefano Zampini   /* index set for vertices */
160*a873d5d0SStefano Zampini   if (nvc) {
161b8ffe317SStefano Zampini     nvc = 0;
162674ae819SStefano Zampini     for (i=0;i<graph->ncc;i++) {
163674ae819SStefano Zampini       if (graph->cptr[i+1]-graph->cptr[i] <= graph->custom_minimal_size) {
164674ae819SStefano Zampini         for (j=graph->cptr[i];j<graph->cptr[i+1];j++) {
165674ae819SStefano Zampini           idx[nvc]=graph->queue[j];
166674ae819SStefano Zampini           nvc++;
167674ae819SStefano Zampini         }
168674ae819SStefano Zampini       }
169674ae819SStefano Zampini     }
170674ae819SStefano Zampini     /* sort vertex set (by local ordering) */
171674ae819SStefano Zampini     ierr = PetscSortInt(nvc,idx);CHKERRQ(ierr);
172674ae819SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,nvc,idx,PETSC_OWN_POINTER,&ISForVertices);CHKERRQ(ierr);
173674ae819SStefano Zampini   }
174674ae819SStefano Zampini   /* get back info */
175*a873d5d0SStefano Zampini   if (n_faces) *n_faces = nfc;
176d7796978SStefano Zampini   if (FacesIS) {
177674ae819SStefano Zampini     *FacesIS = ISForFaces;
178d7796978SStefano Zampini   } else {
179d7796978SStefano Zampini     for (i=0;i<nfc;i++) {
180d7796978SStefano Zampini       ierr = ISDestroy(&ISForFaces[i]);CHKERRQ(ierr);
181d7796978SStefano Zampini     }
182d7796978SStefano Zampini     ierr = PetscFree(ISForFaces);CHKERRQ(ierr);
183d7796978SStefano Zampini   }
184*a873d5d0SStefano Zampini   if (n_edges) *n_edges = nec;
185d7796978SStefano Zampini   if (EdgesIS) {
186674ae819SStefano Zampini     *EdgesIS = ISForEdges;
187d7796978SStefano Zampini   } else {
188d7796978SStefano Zampini     for (i=0;i<nec;i++) {
189d7796978SStefano Zampini       ierr = ISDestroy(&ISForEdges[i]);CHKERRQ(ierr);
190d7796978SStefano Zampini     }
191d7796978SStefano Zampini     ierr = PetscFree(ISForEdges);CHKERRQ(ierr);
192d7796978SStefano Zampini   }
193d7796978SStefano Zampini   if (VerticesIS) {
194674ae819SStefano Zampini     *VerticesIS = ISForVertices;
195d7796978SStefano Zampini   } else {
196d7796978SStefano Zampini     ierr = ISDestroy(&ISForVertices);CHKERRQ(ierr);
197d7796978SStefano Zampini   }
198674ae819SStefano Zampini   PetscFunctionReturn(0);
199674ae819SStefano Zampini }
200674ae819SStefano Zampini 
201674ae819SStefano Zampini #undef __FUNCT__
202674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphComputeConnectedComponents"
203674ae819SStefano Zampini PetscErrorCode PCBDDCGraphComputeConnectedComponents(PCBDDCGraph graph)
204674ae819SStefano Zampini {
2055b1b9aeaSStefano Zampini   PetscInt       i,adapt_interface,adapt_interface_reduced;
206674ae819SStefano Zampini   MPI_Comm       interface_comm;
207674ae819SStefano Zampini   PetscErrorCode ierr;
208674ae819SStefano Zampini 
209674ae819SStefano Zampini   PetscFunctionBegin;
210674ae819SStefano Zampini   /* compute connected components locally */
211674ae819SStefano Zampini   ierr = PetscObjectGetComm((PetscObject)(graph->l2gmap),&interface_comm);CHKERRQ(ierr);
212674ae819SStefano Zampini   ierr = PCBDDCGraphComputeConnectedComponentsLocal(graph);CHKERRQ(ierr);
213674ae819SStefano Zampini   /* check consistency of connected components among neighbouring subdomains -> it adapt them in case it is needed */
214674ae819SStefano Zampini   adapt_interface = 0;
215674ae819SStefano Zampini   adapt_interface_reduced = 0;
216674ae819SStefano Zampini   for (i=0;i<graph->n_subsets;i++) {
217674ae819SStefano Zampini     /* We are not sure that on a given subset of the local interface,
218674ae819SStefano Zampini        with two connected components, the latters be the same among sharing subdomains */
219674ae819SStefano Zampini     if (graph->subset_ncc[i] > 1) {
220674ae819SStefano Zampini       adapt_interface=1;
221674ae819SStefano Zampini       break;
222674ae819SStefano Zampini     }
223674ae819SStefano Zampini   }
224674ae819SStefano Zampini   ierr = MPI_Allreduce(&adapt_interface,&adapt_interface_reduced,1,MPIU_INT,MPI_LOR,interface_comm);CHKERRQ(ierr);
225674ae819SStefano Zampini 
226674ae819SStefano Zampini   if (graph->n_subsets && adapt_interface_reduced) {
2275b1b9aeaSStefano Zampini     MPI_Request *send_requests;
2285b1b9aeaSStefano Zampini     MPI_Request *recv_requests;
2295b1b9aeaSStefano Zampini     PetscInt    *aux_new_xadj,*new_xadj,*new_adjncy,**temp_buffer;
2305b1b9aeaSStefano Zampini     PetscInt    *old_xadj,*old_adjncy;
2315b1b9aeaSStefano Zampini     PetscInt    j,k,s,sum_requests,buffer_size,size_of_recv,temp_buffer_size;
2325b1b9aeaSStefano Zampini     PetscMPIInt rank,neigh,tag,mpi_buffer_size;
2335b1b9aeaSStefano Zampini     PetscInt    *cum_recv_counts,*subset_to_nodes_indices,*recv_buffer_subset,*nodes_to_temp_buffer_indices;
2345b1b9aeaSStefano Zampini     PetscInt    *send_buffer,*recv_buffer,*queue_in_global_numbering,*sizes_of_sends,*add_to_subset;
2355b1b9aeaSStefano Zampini     PetscInt    start_of_recv,start_of_send,size_of_send,global_subset_counter,ins_val;
2365b1b9aeaSStefano Zampini     PetscBool   *subset_cc_adapt,same_set;
2375b1b9aeaSStefano Zampini 
238674ae819SStefano Zampini     /* Retrict adjacency graph using information from previously computed connected components */
239785e854fSJed Brown     ierr = PetscMalloc1(graph->nvtxs,&aux_new_xadj);CHKERRQ(ierr);
240674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
241674ae819SStefano Zampini       aux_new_xadj[i]=1;
242674ae819SStefano Zampini     }
243674ae819SStefano Zampini     for (i=0;i<graph->ncc;i++) {
244674ae819SStefano Zampini       k = graph->cptr[i+1]-graph->cptr[i];
245674ae819SStefano Zampini       for (j=0;j<k;j++) {
246674ae819SStefano Zampini         aux_new_xadj[graph->queue[graph->cptr[i]+j]]=k;
247674ae819SStefano Zampini       }
248674ae819SStefano Zampini     }
249674ae819SStefano Zampini     j = 0;
250674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
251674ae819SStefano Zampini       j += aux_new_xadj[i];
252674ae819SStefano Zampini     }
253785e854fSJed Brown     ierr = PetscMalloc1((graph->nvtxs+1),&new_xadj);CHKERRQ(ierr);
254785e854fSJed Brown     ierr = PetscMalloc1(j,&new_adjncy);CHKERRQ(ierr);
255674ae819SStefano Zampini     new_xadj[0]=0;
256674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
257674ae819SStefano Zampini       new_xadj[i+1]=new_xadj[i]+aux_new_xadj[i];
258674ae819SStefano Zampini       if (aux_new_xadj[i]==1) {
259674ae819SStefano Zampini         new_adjncy[new_xadj[i]]=i;
260674ae819SStefano Zampini       }
261674ae819SStefano Zampini     }
262674ae819SStefano Zampini     ierr = PetscFree(aux_new_xadj);CHKERRQ(ierr);
263674ae819SStefano Zampini     for (i=0;i<graph->ncc;i++) {
264674ae819SStefano Zampini       k = graph->cptr[i+1]-graph->cptr[i];
265674ae819SStefano Zampini       for (j=0;j<k;j++) {
266674ae819SStefano Zampini         ierr = PetscMemcpy(&new_adjncy[new_xadj[graph->queue[graph->cptr[i]+j]]],&graph->queue[graph->cptr[i]],k*sizeof(PetscInt));CHKERRQ(ierr);
267674ae819SStefano Zampini       }
268674ae819SStefano Zampini     }
2695b1b9aeaSStefano Zampini     /* set temporarly new CSR into graph */
2705b1b9aeaSStefano Zampini     old_xadj = graph->xadj;
2715b1b9aeaSStefano Zampini     old_adjncy = graph->adjncy;
272674ae819SStefano Zampini     graph->xadj = new_xadj;
273674ae819SStefano Zampini     graph->adjncy = new_adjncy;
274674ae819SStefano Zampini     /* allocate some space */
275674ae819SStefano Zampini     ierr = MPI_Comm_rank(interface_comm,&rank);CHKERRQ(ierr);
276785e854fSJed Brown     ierr = PetscMalloc1((graph->n_subsets+1),&cum_recv_counts);CHKERRQ(ierr);
277674ae819SStefano Zampini     ierr = PetscMemzero(cum_recv_counts,(graph->n_subsets+1)*sizeof(*cum_recv_counts));CHKERRQ(ierr);
278785e854fSJed Brown     ierr = PetscMalloc1(graph->n_subsets,&subset_to_nodes_indices);CHKERRQ(ierr);
279674ae819SStefano Zampini     /* first count how many neighbours per connected component I will receive from */
280674ae819SStefano Zampini     cum_recv_counts[0]=0;
281674ae819SStefano Zampini     for (i=1;i<graph->n_subsets+1;i++) {
282674ae819SStefano Zampini       j = 0;
283674ae819SStefano Zampini       while (graph->subset[j] != i) {
284674ae819SStefano Zampini         j++;
285674ae819SStefano Zampini       }
286674ae819SStefano Zampini       subset_to_nodes_indices[i-1]=j;
287674ae819SStefano Zampini       /* We don't want sends/recvs_to/from_self -> here I don't count myself  */
288674ae819SStefano Zampini       cum_recv_counts[i]=cum_recv_counts[i-1]+graph->count[j];
289674ae819SStefano Zampini     }
290785e854fSJed Brown     ierr = PetscMalloc1(2*cum_recv_counts[graph->n_subsets],&recv_buffer_subset);CHKERRQ(ierr);
291785e854fSJed Brown     ierr = PetscMalloc1(cum_recv_counts[graph->n_subsets],&send_requests);CHKERRQ(ierr);
292785e854fSJed Brown     ierr = PetscMalloc1(cum_recv_counts[graph->n_subsets],&recv_requests);CHKERRQ(ierr);
293674ae819SStefano Zampini     for (i=0;i<cum_recv_counts[graph->n_subsets];i++) {
294674ae819SStefano Zampini       send_requests[i]=MPI_REQUEST_NULL;
295674ae819SStefano Zampini       recv_requests[i]=MPI_REQUEST_NULL;
296674ae819SStefano Zampini     }
297674ae819SStefano Zampini     /* exchange with my neighbours the number of my connected components on the shared interface */
298674ae819SStefano Zampini     sum_requests = 0;
299674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
300674ae819SStefano Zampini       j = subset_to_nodes_indices[i];
3013a5b03d1SStefano Zampini       ierr = PetscMPIIntCast(graph->subset_ref_node[i],&tag);CHKERRQ(ierr);
302674ae819SStefano Zampini       for (k=0;k<graph->count[j];k++) {
303674ae819SStefano Zampini         ierr = PetscMPIIntCast(graph->neighbours_set[j][k],&neigh);CHKERRQ(ierr);
304674ae819SStefano Zampini         ierr = MPI_Isend(&graph->subset_ncc[i],1,MPIU_INT,neigh,tag,interface_comm,&send_requests[sum_requests]);CHKERRQ(ierr);
305674ae819SStefano Zampini         ierr = MPI_Irecv(&recv_buffer_subset[sum_requests],1,MPIU_INT,neigh,tag,interface_comm,&recv_requests[sum_requests]);CHKERRQ(ierr);
306674ae819SStefano Zampini         sum_requests++;
307674ae819SStefano Zampini       }
308674ae819SStefano Zampini     }
309674ae819SStefano Zampini     ierr = MPI_Waitall(sum_requests,recv_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
310674ae819SStefano Zampini     ierr = MPI_Waitall(sum_requests,send_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
311674ae819SStefano Zampini     /* determine the connected component I need to adapt */
312785e854fSJed Brown     ierr = PetscMalloc1(graph->n_subsets,&subset_cc_adapt);CHKERRQ(ierr);
313674ae819SStefano Zampini     ierr = PetscMemzero(subset_cc_adapt,graph->n_subsets*sizeof(*subset_cc_adapt));CHKERRQ(ierr);
314674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
315674ae819SStefano Zampini       for (j=cum_recv_counts[i];j<cum_recv_counts[i+1];j++){
316674ae819SStefano Zampini         /* The first condition is natural (someone has a different number of ccs than me), the second one is just to be safe */
317674ae819SStefano Zampini         if (graph->subset_ncc[i] != recv_buffer_subset[j] || graph->subset_ncc[i] > 1) {
318674ae819SStefano Zampini           subset_cc_adapt[i] = PETSC_TRUE;
319674ae819SStefano Zampini           break;
320674ae819SStefano Zampini         }
321674ae819SStefano Zampini       }
322674ae819SStefano Zampini     }
323674ae819SStefano Zampini     buffer_size = 0;
324674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
325674ae819SStefano Zampini       if (subset_cc_adapt[i]) {
326674ae819SStefano Zampini         for (j=i;j<graph->ncc;j++) {
327674ae819SStefano Zampini           if (graph->subset[graph->queue[graph->cptr[j]]] == i+1) { /* WARNING -> subset values goes from 1 to graph->n_subsets included */
328674ae819SStefano Zampini             buffer_size += 1 + graph->cptr[j+1]-graph->cptr[j];
329674ae819SStefano Zampini           }
330674ae819SStefano Zampini         }
331674ae819SStefano Zampini       }
332674ae819SStefano Zampini     }
333785e854fSJed Brown     ierr = PetscMalloc1(buffer_size,&send_buffer);CHKERRQ(ierr);
334674ae819SStefano Zampini     /* now get from neighbours their ccs (in global numbering) and adapt them (in case it is needed) */
335785e854fSJed Brown     ierr = PetscMalloc1(graph->cptr[graph->ncc],&queue_in_global_numbering);CHKERRQ(ierr);
336674ae819SStefano Zampini     ierr = ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_in_global_numbering);CHKERRQ(ierr);
337674ae819SStefano Zampini     /* determine how much data to send (size of each queue plus the global indices) and communicate it to neighbours */
338785e854fSJed Brown     ierr = PetscMalloc1(graph->n_subsets,&sizes_of_sends);CHKERRQ(ierr);
339674ae819SStefano Zampini     ierr = PetscMemzero(sizes_of_sends,graph->n_subsets*sizeof(*sizes_of_sends));CHKERRQ(ierr);
340674ae819SStefano Zampini     sum_requests = 0;
341674ae819SStefano Zampini     start_of_send = 0;
342674ae819SStefano Zampini     start_of_recv = cum_recv_counts[graph->n_subsets];
343674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
344674ae819SStefano Zampini       if (subset_cc_adapt[i]) {
345674ae819SStefano Zampini         size_of_send = 0;
346674ae819SStefano Zampini         for (j=i;j<graph->ncc;j++) {
347674ae819SStefano Zampini           if (graph->subset[graph->queue[graph->cptr[j]]] == i+1) { /* WARNING -> subset values goes from 1 to graph->n_subsets included */
348674ae819SStefano Zampini             send_buffer[start_of_send+size_of_send]=graph->cptr[j+1]-graph->cptr[j];
349674ae819SStefano Zampini             size_of_send += 1;
350674ae819SStefano Zampini             ierr = PetscMemcpy(&send_buffer[start_of_send+size_of_send],
351674ae819SStefano Zampini                                &queue_in_global_numbering[graph->cptr[j]],
352674ae819SStefano Zampini                                (graph->cptr[j+1]-graph->cptr[j])*sizeof(*send_buffer));CHKERRQ(ierr);
353674ae819SStefano Zampini             size_of_send = size_of_send+graph->cptr[j+1]-graph->cptr[j];
354674ae819SStefano Zampini           }
355674ae819SStefano Zampini         }
356674ae819SStefano Zampini         j = subset_to_nodes_indices[i];
357674ae819SStefano Zampini         sizes_of_sends[i] = size_of_send;
3583a5b03d1SStefano Zampini         ierr = PetscMPIIntCast(graph->subset_ref_node[i]+1,&tag);CHKERRQ(ierr);
359674ae819SStefano Zampini         for (k=0;k<graph->count[j];k++) {
360674ae819SStefano Zampini           ierr = PetscMPIIntCast(graph->neighbours_set[j][k],&neigh);CHKERRQ(ierr);
361674ae819SStefano Zampini           ierr = MPI_Isend(&sizes_of_sends[i],1,MPIU_INT,neigh,tag,interface_comm,&send_requests[sum_requests]);CHKERRQ(ierr);
362674ae819SStefano Zampini           ierr = MPI_Irecv(&recv_buffer_subset[sum_requests+start_of_recv],1,MPIU_INT,neigh,tag,interface_comm,&recv_requests[sum_requests]);CHKERRQ(ierr);
363674ae819SStefano Zampini           sum_requests++;
364674ae819SStefano Zampini         }
365674ae819SStefano Zampini         start_of_send += size_of_send;
366674ae819SStefano Zampini       }
367674ae819SStefano Zampini     }
368674ae819SStefano Zampini     ierr = MPI_Waitall(sum_requests,send_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
369674ae819SStefano Zampini     ierr = MPI_Waitall(sum_requests,recv_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
370674ae819SStefano Zampini     buffer_size = 0;
371674ae819SStefano Zampini     for (k=0;k<sum_requests;k++) {
372674ae819SStefano Zampini       buffer_size += recv_buffer_subset[start_of_recv+k];
373674ae819SStefano Zampini     }
374785e854fSJed Brown     ierr = PetscMalloc1(buffer_size,&recv_buffer);CHKERRQ(ierr);
375674ae819SStefano Zampini     /* now exchange the data */
376674ae819SStefano Zampini     start_of_recv = 0;
377674ae819SStefano Zampini     start_of_send = 0;
378674ae819SStefano Zampini     sum_requests = 0;
379674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
380674ae819SStefano Zampini       if (subset_cc_adapt[i]) {
381674ae819SStefano Zampini         size_of_send = sizes_of_sends[i];
382674ae819SStefano Zampini         j = subset_to_nodes_indices[i];
3833a5b03d1SStefano Zampini         ierr = PetscMPIIntCast(graph->subset_ref_node[i]+2,&tag);CHKERRQ(ierr);
384674ae819SStefano Zampini         for (k=0;k<graph->count[j];k++) {
385674ae819SStefano Zampini           ierr = PetscMPIIntCast(graph->neighbours_set[j][k],&neigh);CHKERRQ(ierr);
386674ae819SStefano Zampini           ierr = PetscMPIIntCast(size_of_send,&mpi_buffer_size);CHKERRQ(ierr);
387674ae819SStefano Zampini           ierr = MPI_Isend(&send_buffer[start_of_send],mpi_buffer_size,MPIU_INT,neigh,tag,interface_comm,&send_requests[sum_requests]);CHKERRQ(ierr);
388674ae819SStefano Zampini           size_of_recv = recv_buffer_subset[cum_recv_counts[graph->n_subsets]+sum_requests];
389674ae819SStefano Zampini           ierr = PetscMPIIntCast(size_of_recv,&mpi_buffer_size);CHKERRQ(ierr);
390674ae819SStefano Zampini           ierr = MPI_Irecv(&recv_buffer[start_of_recv],mpi_buffer_size,MPIU_INT,neigh,tag,interface_comm,&recv_requests[sum_requests]);CHKERRQ(ierr);
391674ae819SStefano Zampini           start_of_recv += size_of_recv;
392674ae819SStefano Zampini           sum_requests++;
393674ae819SStefano Zampini         }
394674ae819SStefano Zampini         start_of_send += size_of_send;
395674ae819SStefano Zampini       }
396674ae819SStefano Zampini     }
397674ae819SStefano Zampini     ierr = MPI_Waitall(sum_requests,recv_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
398674ae819SStefano Zampini     ierr = MPI_Waitall(sum_requests,send_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
399674ae819SStefano Zampini     for (j=0;j<buffer_size;) {
400674ae819SStefano Zampini        ierr = ISGlobalToLocalMappingApply(graph->l2gmap,IS_GTOLM_MASK,recv_buffer[j],&recv_buffer[j+1],&recv_buffer[j],&recv_buffer[j+1]);CHKERRQ(ierr);
40151b0f6cfSStefano Zampini        /* we need to adapt the output of GlobalToLocal mapping if there are mirrored nodes */
40251b0f6cfSStefano Zampini        if (graph->mirrors) {
40351b0f6cfSStefano Zampini          PetscBool mirrored_found=PETSC_FALSE;
40451b0f6cfSStefano Zampini          for (k=0;k<recv_buffer[j];k++) {
40551b0f6cfSStefano Zampini            if (graph->mirrors[recv_buffer[j+k+1]]) {
40651b0f6cfSStefano Zampini              mirrored_found=PETSC_TRUE;
40751b0f6cfSStefano Zampini              recv_buffer[j+k+1]=graph->mirrors_set[recv_buffer[j+k+1]][0];
40851b0f6cfSStefano Zampini            }
40951b0f6cfSStefano Zampini          }
41051b0f6cfSStefano Zampini          if (mirrored_found) {
41151b0f6cfSStefano Zampini            ierr = PetscSortInt(recv_buffer[j],&recv_buffer[j+1]);CHKERRQ(ierr);
41251b0f6cfSStefano Zampini            k=0;
41351b0f6cfSStefano Zampini            while (k<recv_buffer[j]) {
41451b0f6cfSStefano Zampini              for (s=1;s<graph->mirrors[recv_buffer[j+1+k]];s++) {
41551b0f6cfSStefano Zampini                recv_buffer[j+1+k+s] = graph->mirrors_set[recv_buffer[j+1+k]][s];
41651b0f6cfSStefano Zampini              }
41751b0f6cfSStefano Zampini              k+=graph->mirrors[recv_buffer[j+1+k]]+s;
41851b0f6cfSStefano Zampini            }
41951b0f6cfSStefano Zampini          }
42051b0f6cfSStefano Zampini        }
421674ae819SStefano Zampini        k = recv_buffer[j]+1;
422674ae819SStefano Zampini        j += k;
423674ae819SStefano Zampini     }
424674ae819SStefano Zampini     sum_requests = cum_recv_counts[graph->n_subsets];
425674ae819SStefano Zampini     start_of_recv = 0;
426785e854fSJed Brown     ierr = PetscMalloc1(graph->nvtxs,&nodes_to_temp_buffer_indices);CHKERRQ(ierr);
427674ae819SStefano Zampini     global_subset_counter = 0;
428674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
429674ae819SStefano Zampini       if (subset_cc_adapt[i]) {
430674ae819SStefano Zampini         temp_buffer_size = 0;
431674ae819SStefano Zampini         /* find nodes on the shared interface we need to adapt */
432674ae819SStefano Zampini         for (j=0;j<graph->nvtxs;j++) {
433674ae819SStefano Zampini           if (graph->subset[j]==i+1) {
434674ae819SStefano Zampini             nodes_to_temp_buffer_indices[j] = temp_buffer_size;
435674ae819SStefano Zampini             temp_buffer_size++;
436674ae819SStefano Zampini           } else {
437674ae819SStefano Zampini             nodes_to_temp_buffer_indices[j] = -1;
438674ae819SStefano Zampini           }
439674ae819SStefano Zampini         }
440674ae819SStefano Zampini         /* allocate some temporary space */
441785e854fSJed Brown         ierr = PetscMalloc1(temp_buffer_size,&temp_buffer);CHKERRQ(ierr);
442785e854fSJed Brown         ierr = PetscMalloc1(temp_buffer_size*(cum_recv_counts[i+1]-cum_recv_counts[i]),&temp_buffer[0]);CHKERRQ(ierr);
443674ae819SStefano Zampini         ierr = PetscMemzero(temp_buffer[0],temp_buffer_size*(cum_recv_counts[i+1]-cum_recv_counts[i])*sizeof(PetscInt));CHKERRQ(ierr);
444674ae819SStefano Zampini         for (j=1;j<temp_buffer_size;j++) {
445674ae819SStefano Zampini           temp_buffer[j] = temp_buffer[j-1]+cum_recv_counts[i+1]-cum_recv_counts[i];
446674ae819SStefano Zampini         }
447674ae819SStefano Zampini         /* analyze contributions from neighbouring subdomains for i-th conn comp
448674ae819SStefano Zampini            temp buffer structure:
449674ae819SStefano Zampini            supposing part of the interface has dimension 5 (for example with global dofs 0,1,2,3,4)
450674ae819SStefano Zampini            3 neighs procs with structured connected components:
451674ae819SStefano Zampini              neigh 0: [0 1 4], [2 3];  (2 connected components)
452674ae819SStefano Zampini              neigh 1: [0 1], [2 3 4];  (2 connected components)
453674ae819SStefano Zampini              neigh 2: [0 4], [1], [2 3]; (3 connected components)
454674ae819SStefano Zampini            tempbuffer (row-oriented) will be filled as:
455674ae819SStefano Zampini              [ 0, 0, 0;
456674ae819SStefano Zampini                0, 0, 1;
457674ae819SStefano Zampini                1, 1, 2;
458674ae819SStefano Zampini                1, 1, 2;
459674ae819SStefano Zampini                0, 1, 0; ];
460674ae819SStefano Zampini            This way we can simply find intersections of ccs among neighs.
461674ae819SStefano Zampini            For the example above, the graph->subset array will be modified to reproduce the following 4 connected components [0], [1], [2 3], [4];
462674ae819SStefano Zampini                                                                                                                                    */
463674ae819SStefano Zampini         for (j=0;j<cum_recv_counts[i+1]-cum_recv_counts[i];j++) {
464674ae819SStefano Zampini           ins_val = 0;
465674ae819SStefano Zampini           size_of_recv = recv_buffer_subset[sum_requests];  /* total size of recv from neighs */
466674ae819SStefano Zampini           for (buffer_size=0;buffer_size<size_of_recv;) {  /* loop until all data from neighs has been taken into account */
467674ae819SStefano Zampini             for (k=1;k<recv_buffer[buffer_size+start_of_recv]+1;k++) { /* filling properly temp_buffer using data from a single recv */
468674ae819SStefano Zampini               temp_buffer[nodes_to_temp_buffer_indices[recv_buffer[start_of_recv+buffer_size+k]]][j] = ins_val;
469674ae819SStefano Zampini             }
470674ae819SStefano Zampini             buffer_size += k;
471674ae819SStefano Zampini             ins_val++;
472674ae819SStefano Zampini           }
473674ae819SStefano Zampini           start_of_recv += size_of_recv;
474674ae819SStefano Zampini           sum_requests++;
475674ae819SStefano Zampini         }
476785e854fSJed Brown         ierr = PetscMalloc1(temp_buffer_size,&add_to_subset);CHKERRQ(ierr);
477674ae819SStefano Zampini         ierr = PetscMemzero(add_to_subset,temp_buffer_size*sizeof(*add_to_subset));CHKERRQ(ierr);
478674ae819SStefano Zampini         for (j=0;j<temp_buffer_size;j++) {
479674ae819SStefano Zampini           if (!add_to_subset[j]) { /* found a new cc  */
480674ae819SStefano Zampini             global_subset_counter++;
481674ae819SStefano Zampini             add_to_subset[j] = global_subset_counter;
482674ae819SStefano Zampini             for (k=j+1;k<temp_buffer_size;k++) { /* check for other nodes in new cc */
483674ae819SStefano Zampini               same_set = PETSC_TRUE;
484674ae819SStefano Zampini               for (s=0;s<cum_recv_counts[i+1]-cum_recv_counts[i];s++) {
485674ae819SStefano Zampini                 if (temp_buffer[j][s]!=temp_buffer[k][s]) {
486674ae819SStefano Zampini                   same_set = PETSC_FALSE;
487674ae819SStefano Zampini                   break;
488674ae819SStefano Zampini                 }
489674ae819SStefano Zampini               }
490674ae819SStefano Zampini               if (same_set) {
491674ae819SStefano Zampini                 add_to_subset[k] = global_subset_counter;
492674ae819SStefano Zampini               }
493674ae819SStefano Zampini             }
494674ae819SStefano Zampini           }
495674ae819SStefano Zampini         }
496674ae819SStefano Zampini         /* insert new data in subset array */
497674ae819SStefano Zampini         temp_buffer_size = 0;
498674ae819SStefano Zampini         for (j=0;j<graph->nvtxs;j++) {
499674ae819SStefano Zampini           if (graph->subset[j]==i+1) {
500674ae819SStefano Zampini             graph->subset[j] = graph->n_subsets+add_to_subset[temp_buffer_size];
501674ae819SStefano Zampini             temp_buffer_size++;
502674ae819SStefano Zampini           }
503674ae819SStefano Zampini         }
504674ae819SStefano Zampini         ierr = PetscFree(temp_buffer[0]);CHKERRQ(ierr);
505674ae819SStefano Zampini         ierr = PetscFree(temp_buffer);CHKERRQ(ierr);
506674ae819SStefano Zampini         ierr = PetscFree(add_to_subset);CHKERRQ(ierr);
507674ae819SStefano Zampini       }
508674ae819SStefano Zampini     }
509674ae819SStefano Zampini     ierr = PetscFree(nodes_to_temp_buffer_indices);CHKERRQ(ierr);
510674ae819SStefano Zampini     ierr = PetscFree(sizes_of_sends);CHKERRQ(ierr);
511674ae819SStefano Zampini     ierr = PetscFree(send_requests);CHKERRQ(ierr);
512674ae819SStefano Zampini     ierr = PetscFree(recv_requests);CHKERRQ(ierr);
513674ae819SStefano Zampini     ierr = PetscFree(recv_buffer);CHKERRQ(ierr);
514674ae819SStefano Zampini     ierr = PetscFree(recv_buffer_subset);CHKERRQ(ierr);
515674ae819SStefano Zampini     ierr = PetscFree(send_buffer);CHKERRQ(ierr);
516674ae819SStefano Zampini     ierr = PetscFree(cum_recv_counts);CHKERRQ(ierr);
517674ae819SStefano Zampini     ierr = PetscFree(subset_to_nodes_indices);CHKERRQ(ierr);
518674ae819SStefano Zampini     ierr = PetscFree(subset_cc_adapt);CHKERRQ(ierr);
519674ae819SStefano Zampini     /* We are ready to find for connected components consistent among neighbouring subdomains */
520674ae819SStefano Zampini     if (global_subset_counter) {
521df48933bSStefano Zampini       ierr = PetscBTMemzero(graph->nvtxs,graph->touched);CHKERRQ(ierr);
522674ae819SStefano Zampini       global_subset_counter = 0;
523674ae819SStefano Zampini       for (i=0;i<graph->nvtxs;i++) {
524df48933bSStefano Zampini         if (graph->subset[i] && !PetscBTLookup(graph->touched,i)) {
525674ae819SStefano Zampini           global_subset_counter++;
526674ae819SStefano Zampini           for (j=i+1;j<graph->nvtxs;j++) {
527df48933bSStefano Zampini             if (!PetscBTLookup(graph->touched,j) && graph->subset[j]==graph->subset[i]) {
528674ae819SStefano Zampini               graph->subset[j] = global_subset_counter;
529df48933bSStefano Zampini               ierr = PetscBTSet(graph->touched,j);CHKERRQ(ierr);
530674ae819SStefano Zampini             }
531674ae819SStefano Zampini           }
532674ae819SStefano Zampini           graph->subset[i] = global_subset_counter;
533df48933bSStefano Zampini           ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
534674ae819SStefano Zampini         }
535674ae819SStefano Zampini       }
536674ae819SStefano Zampini       /* refine connected components locally */
537674ae819SStefano Zampini       ierr = PCBDDCGraphComputeConnectedComponentsLocal(graph);CHKERRQ(ierr);
538674ae819SStefano Zampini     }
5395b1b9aeaSStefano Zampini     /* restore original CSR graph of dofs */
5405b1b9aeaSStefano Zampini     ierr = PetscFree(new_xadj);CHKERRQ(ierr);
5415b1b9aeaSStefano Zampini     ierr = PetscFree(new_adjncy);CHKERRQ(ierr);
5425b1b9aeaSStefano Zampini     graph->xadj = old_xadj;
5435b1b9aeaSStefano Zampini     graph->adjncy = old_adjncy;
544674ae819SStefano Zampini     ierr = PetscFree(queue_in_global_numbering);CHKERRQ(ierr);
545674ae819SStefano Zampini   }
546674ae819SStefano Zampini   PetscFunctionReturn(0);
547674ae819SStefano Zampini }
548674ae819SStefano Zampini 
549674ae819SStefano Zampini /* The following code has been adapted from function IsConnectedSubdomain contained
550674ae819SStefano Zampini    in source file contig.c of METIS library (version 5.0.1)
551674ae819SStefano Zampini    It finds connected components for each subset  */
552674ae819SStefano Zampini #undef __FUNCT__
553674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphComputeConnectedComponentsLocal"
554674ae819SStefano Zampini PetscErrorCode PCBDDCGraphComputeConnectedComponentsLocal(PCBDDCGraph graph)
555674ae819SStefano Zampini {
556674ae819SStefano Zampini   PetscInt       i,j,k,first,last,nleft,ncc,pid,cum_queue,n,ncc_pid;
557674ae819SStefano Zampini   PetscInt       *queue_global;
558674ae819SStefano Zampini   PetscErrorCode ierr;
559674ae819SStefano Zampini 
560674ae819SStefano Zampini   PetscFunctionBegin;
561674ae819SStefano Zampini   /* quiet return if no csr info is available */
562674ae819SStefano Zampini   if (!graph->xadj || !graph->adjncy) {
563674ae819SStefano Zampini     PetscFunctionReturn(0);
564674ae819SStefano Zampini   }
565674ae819SStefano Zampini 
566674ae819SStefano Zampini   /* reset any previous search of connected components */
567df48933bSStefano Zampini   ierr = PetscBTMemzero(graph->nvtxs,graph->touched);CHKERRQ(ierr);
568674ae819SStefano Zampini   graph->n_subsets = 0;
569674ae819SStefano Zampini   for (i=0;i<graph->nvtxs;i++) {
5700cf82fd6SStefano Zampini     if (graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK || !graph->count[i]) {
571df48933bSStefano Zampini       ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
572674ae819SStefano Zampini       graph->subset[i] = 0;
573674ae819SStefano Zampini     }
574674ae819SStefano Zampini     graph->n_subsets = PetscMax(graph->n_subsets,graph->subset[i]);
575674ae819SStefano Zampini   }
576674ae819SStefano Zampini   ierr = PetscFree(graph->subset_ncc);CHKERRQ(ierr);
577785e854fSJed Brown   ierr = PetscMalloc1(graph->n_subsets,&graph->subset_ncc);CHKERRQ(ierr);
578674ae819SStefano Zampini   ierr = PetscMemzero(graph->subset_ncc,graph->n_subsets*sizeof(*graph->subset_ncc));CHKERRQ(ierr);
579674ae819SStefano Zampini   ierr = PetscMemzero(graph->cptr,(graph->nvtxs+1)*sizeof(*graph->cptr));CHKERRQ(ierr);
580674ae819SStefano Zampini   ierr = PetscMemzero(graph->queue,graph->nvtxs*sizeof(*graph->queue));CHKERRQ(ierr);
581674ae819SStefano Zampini 
582674ae819SStefano Zampini   /* begin search for connected components */
583674ae819SStefano Zampini   cum_queue = 0;
584674ae819SStefano Zampini   ncc = 0;
585674ae819SStefano Zampini   for (n=0;n<graph->n_subsets;n++) {
586674ae819SStefano Zampini     pid = n+1;  /* partition labeled by 0 is discarded */
587674ae819SStefano Zampini     nleft = 0;
588674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
589674ae819SStefano Zampini       if (graph->subset[i] == pid) {
590674ae819SStefano Zampini         nleft++;
591674ae819SStefano Zampini       }
592674ae819SStefano Zampini     }
593674ae819SStefano Zampini     for (i=0; i<graph->nvtxs; i++) {
594674ae819SStefano Zampini       if (graph->subset[i] == pid) {
595674ae819SStefano Zampini         break;
596674ae819SStefano Zampini       }
597674ae819SStefano Zampini     }
598df48933bSStefano Zampini     ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
599674ae819SStefano Zampini     graph->queue[cum_queue] = i;
600674ae819SStefano Zampini     first = 0;
601674ae819SStefano Zampini     last = 1;
602674ae819SStefano Zampini     graph->cptr[ncc] = cum_queue;
603674ae819SStefano Zampini     ncc_pid = 0;
604674ae819SStefano Zampini     while (first != nleft) {
605674ae819SStefano Zampini       if (first == last) {
606674ae819SStefano Zampini         graph->cptr[++ncc] = first+cum_queue;
607674ae819SStefano Zampini         ncc_pid++;
608df48933bSStefano Zampini         for (i=0; i<graph->nvtxs; i++) { /* TODO-> use a while! */
609df48933bSStefano Zampini           if (graph->subset[i] == pid && !PetscBTLookup(graph->touched,i)) {
610674ae819SStefano Zampini             break;
611674ae819SStefano Zampini           }
612674ae819SStefano Zampini         }
613674ae819SStefano Zampini         graph->queue[cum_queue+last] = i;
614674ae819SStefano Zampini         last++;
615df48933bSStefano Zampini         ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
616674ae819SStefano Zampini       }
617674ae819SStefano Zampini       i = graph->queue[cum_queue+first];
618674ae819SStefano Zampini       first++;
619674ae819SStefano Zampini       for (j=graph->xadj[i];j<graph->xadj[i+1];j++) {
620674ae819SStefano Zampini         k = graph->adjncy[j];
621df48933bSStefano Zampini         if (graph->subset[k] == pid && !PetscBTLookup(graph->touched,k)) {
622674ae819SStefano Zampini           graph->queue[cum_queue+last] = k;
623674ae819SStefano Zampini           last++;
624df48933bSStefano Zampini           ierr = PetscBTSet(graph->touched,k);CHKERRQ(ierr);
625674ae819SStefano Zampini         }
626674ae819SStefano Zampini       }
627674ae819SStefano Zampini     }
628674ae819SStefano Zampini     graph->cptr[++ncc] = first+cum_queue;
629674ae819SStefano Zampini     ncc_pid++;
630674ae819SStefano Zampini     cum_queue = graph->cptr[ncc];
631674ae819SStefano Zampini     graph->subset_ncc[n] = ncc_pid;
632674ae819SStefano Zampini   }
633674ae819SStefano Zampini   graph->ncc = ncc;
634674ae819SStefano Zampini   /* For consistency among neighbours, I need to sort (by global ordering) each connected component */
635785e854fSJed Brown   ierr = PetscMalloc1(graph->cptr[graph->ncc],&queue_global);CHKERRQ(ierr);
636674ae819SStefano Zampini   ierr = ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_global);CHKERRQ(ierr);
637674ae819SStefano Zampini   for (i=0;i<graph->ncc;i++) {
638674ae819SStefano Zampini     ierr = PetscSortIntWithArray(graph->cptr[i+1]-graph->cptr[i],&queue_global[graph->cptr[i]],&graph->queue[graph->cptr[i]]);CHKERRQ(ierr);
639674ae819SStefano Zampini   }
640674ae819SStefano Zampini   ierr = PetscFree(queue_global);CHKERRQ(ierr);
641674ae819SStefano Zampini   PetscFunctionReturn(0);
642674ae819SStefano Zampini }
643674ae819SStefano Zampini 
644674ae819SStefano Zampini #undef __FUNCT__
645674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphSetUp"
646674ae819SStefano Zampini PetscErrorCode PCBDDCGraphSetUp(PCBDDCGraph graph, PetscInt custom_minimal_size, IS neumann_is, IS dirichlet_is, PetscInt n_ISForDofs, IS ISForDofs[], IS custom_primal_vertices)
647674ae819SStefano Zampini {
648674ae819SStefano Zampini   VecScatter     scatter_ctx;
649674ae819SStefano Zampini   Vec            local_vec,local_vec2,global_vec;
650674ae819SStefano Zampini   IS             to,from;
651674ae819SStefano Zampini   MPI_Comm       comm;
652674ae819SStefano Zampini   PetscScalar    *array,*array2;
653674ae819SStefano Zampini   const PetscInt *is_indices;
6543a5b03d1SStefano Zampini   PetscInt       n_neigh,*neigh,*n_shared,**shared,*queue_global;
655674ae819SStefano Zampini   PetscInt       i,j,k,s,total_counts,nodes_touched,is_size;
656674ae819SStefano Zampini   PetscErrorCode ierr;
65751b0f6cfSStefano Zampini   PetscBool      same_set,mirrors_found;
658674ae819SStefano Zampini 
659674ae819SStefano Zampini   PetscFunctionBegin;
660674ae819SStefano Zampini   ierr = PetscObjectGetComm((PetscObject)(graph->l2gmap),&comm);CHKERRQ(ierr);
661674ae819SStefano Zampini   /* custom_minimal_size */
662674ae819SStefano Zampini   graph->custom_minimal_size = PetscMax(graph->custom_minimal_size,custom_minimal_size);
663674ae819SStefano Zampini   /* get info l2gmap and allocate work vectors  */
664674ae819SStefano Zampini   ierr = ISLocalToGlobalMappingGetInfo(graph->l2gmap,&n_neigh,&neigh,&n_shared,&shared);CHKERRQ(ierr);
665674ae819SStefano Zampini   ierr = ISLocalToGlobalMappingGetIndices(graph->l2gmap,&is_indices);CHKERRQ(ierr);
666674ae819SStefano Zampini   j = 0;
667674ae819SStefano Zampini   for (i=0;i<graph->nvtxs;i++) {
668674ae819SStefano Zampini     j = PetscMax(j,is_indices[i]);
669674ae819SStefano Zampini   }
670674ae819SStefano Zampini   ierr = MPI_Allreduce(&j,&i,1,MPIU_INT,MPI_MAX,comm);CHKERRQ(ierr);
671674ae819SStefano Zampini   i++;
672674ae819SStefano Zampini   ierr = VecCreate(PETSC_COMM_SELF,&local_vec);CHKERRQ(ierr);
673674ae819SStefano Zampini   ierr = VecSetSizes(local_vec,PETSC_DECIDE,graph->nvtxs);CHKERRQ(ierr);
674674ae819SStefano Zampini   ierr = VecSetType(local_vec,VECSTANDARD);CHKERRQ(ierr);
675674ae819SStefano Zampini   ierr = VecDuplicate(local_vec,&local_vec2);CHKERRQ(ierr);
676674ae819SStefano Zampini   ierr = VecCreate(comm,&global_vec);CHKERRQ(ierr);
677674ae819SStefano Zampini   ierr = VecSetSizes(global_vec,PETSC_DECIDE,i);CHKERRQ(ierr);
678674ae819SStefano Zampini   ierr = VecSetType(global_vec,VECSTANDARD);CHKERRQ(ierr);
679674ae819SStefano Zampini   ierr = ISCreateStride(PETSC_COMM_SELF,graph->nvtxs,0,1,&to);CHKERRQ(ierr);
680674ae819SStefano Zampini   ierr = ISLocalToGlobalMappingApplyIS(graph->l2gmap,to,&from);CHKERRQ(ierr);
681674ae819SStefano Zampini   ierr = VecScatterCreate(global_vec,from,local_vec,to,&scatter_ctx);CHKERRQ(ierr);
68251b0f6cfSStefano Zampini 
68351b0f6cfSStefano Zampini   /* get local periodic nodes */
68451b0f6cfSStefano Zampini   mirrors_found = PETSC_FALSE;
6859881197aSStefano Zampini   if (graph->nvtxs && n_neigh) {
68649eeff8cSStefano Zampini     for (i=0; i<n_shared[0]; i++) graph->count[shared[0][i]] += 1;
68751b0f6cfSStefano Zampini     for (i=0; i<n_shared[0]; i++) {
68851b0f6cfSStefano Zampini       if (graph->count[shared[0][i]] > 1) {
68951b0f6cfSStefano Zampini         mirrors_found = PETSC_TRUE;
69051b0f6cfSStefano Zampini         break;
69151b0f6cfSStefano Zampini       }
69251b0f6cfSStefano Zampini     }
69349eeff8cSStefano Zampini   }
69451b0f6cfSStefano Zampini   /* compute local mirrors (if any) */
69551b0f6cfSStefano Zampini   if (mirrors_found) {
69651b0f6cfSStefano Zampini     PetscInt *local_indices,*global_indices;
69751b0f6cfSStefano Zampini     /* get arrays of local and global indices */
698785e854fSJed Brown     ierr = PetscMalloc1(graph->nvtxs,&local_indices);CHKERRQ(ierr);
69951b0f6cfSStefano Zampini     ierr = ISGetIndices(to,(const PetscInt**)&is_indices);CHKERRQ(ierr);
70051b0f6cfSStefano Zampini     ierr = PetscMemcpy(local_indices,is_indices,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
70151b0f6cfSStefano Zampini     ierr = ISRestoreIndices(to,(const PetscInt**)&is_indices);CHKERRQ(ierr);
702785e854fSJed Brown     ierr = PetscMalloc1(graph->nvtxs,&global_indices);CHKERRQ(ierr);
70351b0f6cfSStefano Zampini     ierr = ISGetIndices(from,(const PetscInt**)&is_indices);CHKERRQ(ierr);
70451b0f6cfSStefano Zampini     ierr = PetscMemcpy(global_indices,is_indices,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
70551b0f6cfSStefano Zampini     ierr = ISRestoreIndices(from,(const PetscInt**)&is_indices);CHKERRQ(ierr);
70651b0f6cfSStefano Zampini     /* allocate space for mirrors */
7078e5aaad5SJed Brown     ierr = PetscMalloc2(graph->nvtxs,&graph->mirrors,graph->nvtxs,&graph->mirrors_set);CHKERRQ(ierr);
70851b0f6cfSStefano Zampini     ierr = PetscMemzero(graph->mirrors,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
70951b0f6cfSStefano Zampini     graph->mirrors_set[0] = 0;
71051b0f6cfSStefano Zampini 
71151b0f6cfSStefano Zampini     k=0;
71251b0f6cfSStefano Zampini     for (i=0;i<n_shared[0];i++) {
71351b0f6cfSStefano Zampini       j=shared[0][i];
71451b0f6cfSStefano Zampini       if (graph->count[j] > 1) {
71551b0f6cfSStefano Zampini         graph->mirrors[j]++;
71651b0f6cfSStefano Zampini         k++;
71751b0f6cfSStefano Zampini       }
71851b0f6cfSStefano Zampini     }
71951b0f6cfSStefano Zampini     /* allocate space for set of mirrors */
720785e854fSJed Brown     ierr = PetscMalloc1(k,&graph->mirrors_set[0]);CHKERRQ(ierr);
72151b0f6cfSStefano Zampini     for (i=1;i<graph->nvtxs;i++)
72251b0f6cfSStefano Zampini       graph->mirrors_set[i]=graph->mirrors_set[i-1]+graph->mirrors[i-1];
72351b0f6cfSStefano Zampini 
72451b0f6cfSStefano Zampini     /* fill arrays */
72551b0f6cfSStefano Zampini     ierr = PetscMemzero(graph->mirrors,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
72651b0f6cfSStefano Zampini     for (j=0;j<n_shared[0];j++) {
72751b0f6cfSStefano Zampini       i=shared[0][j];
72851b0f6cfSStefano Zampini       if (graph->count[i] > 1)
72951b0f6cfSStefano Zampini         graph->mirrors_set[i][graph->mirrors[i]++]=global_indices[i];
73051b0f6cfSStefano Zampini     }
73151b0f6cfSStefano Zampini     ierr = PetscSortIntWithArray(graph->nvtxs,global_indices,local_indices);CHKERRQ(ierr);
73251b0f6cfSStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
73351b0f6cfSStefano Zampini       if (graph->mirrors[i] > 0) {
73451b0f6cfSStefano Zampini         ierr = PetscFindInt(graph->mirrors_set[i][0],graph->nvtxs,global_indices,&k);CHKERRQ(ierr);
73551b0f6cfSStefano Zampini         j = global_indices[k];
73651b0f6cfSStefano Zampini         while ( k > 0 && global_indices[k-1] == j) k--;
73751b0f6cfSStefano Zampini         for (j=0;j<graph->mirrors[i];j++) {
73851b0f6cfSStefano Zampini           graph->mirrors_set[i][j]=local_indices[k+j];
73951b0f6cfSStefano Zampini         }
74051b0f6cfSStefano Zampini         ierr = PetscSortInt(graph->mirrors[i],graph->mirrors_set[i]);CHKERRQ(ierr);
74151b0f6cfSStefano Zampini       }
74251b0f6cfSStefano Zampini     }
74351b0f6cfSStefano Zampini     ierr = PetscFree(local_indices);CHKERRQ(ierr);
74451b0f6cfSStefano Zampini     ierr = PetscFree(global_indices);CHKERRQ(ierr);
74551b0f6cfSStefano Zampini   }
74651b0f6cfSStefano Zampini   ierr = PetscMemzero(graph->count,graph->nvtxs*sizeof(*graph->count));CHKERRQ(ierr);
747674ae819SStefano Zampini   ierr = ISDestroy(&to);CHKERRQ(ierr);
748674ae819SStefano Zampini   ierr = ISDestroy(&from);CHKERRQ(ierr);
74951b0f6cfSStefano Zampini 
750674ae819SStefano Zampini   /* Count total number of neigh per node */
751674ae819SStefano Zampini   k=0;
752674ae819SStefano Zampini   for (i=1;i<n_neigh;i++) {
753674ae819SStefano Zampini     k += n_shared[i];
754674ae819SStefano Zampini     for (j=0;j<n_shared[i];j++) {
755674ae819SStefano Zampini       graph->count[shared[i][j]] += 1;
756674ae819SStefano Zampini     }
757674ae819SStefano Zampini   }
758674ae819SStefano Zampini   /* Allocate space for storing the set of neighbours for each node */
759674ae819SStefano Zampini   if (graph->nvtxs) {
760785e854fSJed Brown     ierr = PetscMalloc1(k,&graph->neighbours_set[0]);CHKERRQ(ierr);
761674ae819SStefano Zampini   }
762674ae819SStefano Zampini   for (i=1;i<graph->nvtxs;i++) { /* dont count myself */
763674ae819SStefano Zampini     graph->neighbours_set[i]=graph->neighbours_set[i-1]+graph->count[i-1];
764674ae819SStefano Zampini   }
765674ae819SStefano Zampini   /* Get information for sharing subdomains */
766674ae819SStefano Zampini   ierr = PetscMemzero(graph->count,graph->nvtxs*sizeof(*graph->count));CHKERRQ(ierr);
767674ae819SStefano Zampini   for (i=1;i<n_neigh;i++) { /* dont count myself */
768674ae819SStefano Zampini     s = n_shared[i];
769674ae819SStefano Zampini     for (j=0;j<s;j++) {
770674ae819SStefano Zampini       k = shared[i][j];
771674ae819SStefano Zampini       graph->neighbours_set[k][graph->count[k]] = neigh[i];
772674ae819SStefano Zampini       graph->count[k] += 1;
773674ae819SStefano Zampini     }
774674ae819SStefano Zampini   }
775674ae819SStefano Zampini   /* sort set of sharing subdomains */
776674ae819SStefano Zampini   for (i=0;i<graph->nvtxs;i++) {
77793fb5fd3SStefano Zampini     ierr = PetscSortRemoveDupsInt(&graph->count[i],graph->neighbours_set[i]);CHKERRQ(ierr);
778674ae819SStefano Zampini   }
77967c9da69SStefano Zampini   /*
78067c9da69SStefano Zampini      Get info for dofs splitting
78167c9da69SStefano Zampini      User can specify only a subset; an additional field is considered as a complementary set
78267c9da69SStefano Zampini   */
78367c9da69SStefano Zampini   for (i=0;i<graph->nvtxs;i++) {
78467c9da69SStefano Zampini     graph->which_dof[i] = n_ISForDofs; /* by default a dof belongs to the complement set */
78567c9da69SStefano Zampini   }
786674ae819SStefano Zampini   for (i=0;i<n_ISForDofs;i++) {
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;
792674ae819SStefano Zampini       }
79367c9da69SStefano Zampini     }
794674ae819SStefano Zampini     ierr = ISRestoreIndices(ISForDofs[i],(const PetscInt**)&is_indices);CHKERRQ(ierr);
795674ae819SStefano Zampini   }
796674ae819SStefano Zampini   /* Take into account Neumann nodes */
797674ae819SStefano Zampini   ierr = VecSet(local_vec,0.0);CHKERRQ(ierr);
798674ae819SStefano Zampini   ierr = VecSet(local_vec2,0.0);CHKERRQ(ierr);
799674ae819SStefano Zampini   if (neumann_is) {
800674ae819SStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
80182ba6b80SStefano Zampini     ierr = ISGetLocalSize(neumann_is,&is_size);CHKERRQ(ierr);
802674ae819SStefano Zampini     ierr = ISGetIndices(neumann_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
803674ae819SStefano Zampini     for (i=0;i<is_size;i++) {
804d50ed60aSStefano Zampini       if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
805674ae819SStefano Zampini         array[is_indices[i]] = 1.0;
806674ae819SStefano Zampini       }
8073c629590SStefano Zampini     }
808674ae819SStefano Zampini     ierr = ISRestoreIndices(neumann_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
809674ae819SStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
810674ae819SStefano Zampini   }
811674ae819SStefano Zampini   /* Neumann nodes: impose consistency among neighbours */
812674ae819SStefano Zampini   ierr = VecSet(global_vec,0.0);CHKERRQ(ierr);
813674ae819SStefano Zampini   ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
814674ae819SStefano Zampini   ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
815674ae819SStefano Zampini   ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
816674ae819SStefano Zampini   ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
817674ae819SStefano Zampini   ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
818674ae819SStefano Zampini   for (i=0;i<graph->nvtxs;i++) {
8193c629590SStefano Zampini     if (PetscRealPart(array[i]) > 0.1) {
8200cf82fd6SStefano Zampini       graph->special_dof[i] = PCBDDCGRAPH_NEUMANN_MARK;
821674ae819SStefano Zampini     }
822674ae819SStefano Zampini   }
823674ae819SStefano Zampini   ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
824674ae819SStefano Zampini   /* Take into account Dirichlet nodes */
825674ae819SStefano Zampini   ierr = VecSet(local_vec2,0.0);CHKERRQ(ierr);
826674ae819SStefano Zampini   if (dirichlet_is) {
827674ae819SStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
828674ae819SStefano Zampini     ierr = VecGetArray(local_vec2,&array2);CHKERRQ(ierr);
82982ba6b80SStefano Zampini     ierr = ISGetLocalSize(dirichlet_is,&is_size);CHKERRQ(ierr);
830674ae819SStefano Zampini     ierr = ISGetIndices(dirichlet_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
831674ae819SStefano Zampini     for (i=0;i<is_size;i++){
832d50ed60aSStefano Zampini       if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
833674ae819SStefano Zampini         k = is_indices[i];
834df48933bSStefano Zampini         if (graph->count[k] && !PetscBTLookup(graph->touched,k)) {
8353c629590SStefano Zampini           if (PetscRealPart(array[k]) > 0.1) {
8363c629590SStefano 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);
837674ae819SStefano Zampini           }
838674ae819SStefano Zampini           array2[k] = 1.0;
839674ae819SStefano Zampini         }
840674ae819SStefano Zampini       }
8413c629590SStefano Zampini     }
842674ae819SStefano Zampini     ierr = ISRestoreIndices(dirichlet_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
843674ae819SStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
844674ae819SStefano Zampini     ierr = VecRestoreArray(local_vec2,&array2);CHKERRQ(ierr);
845674ae819SStefano Zampini   }
846674ae819SStefano Zampini   /* Dirichlet nodes: impose consistency among neighbours */
847674ae819SStefano Zampini   ierr = VecSet(global_vec,0.0);CHKERRQ(ierr);
848674ae819SStefano Zampini   ierr = VecScatterBegin(scatter_ctx,local_vec2,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
849674ae819SStefano Zampini   ierr = VecScatterEnd(scatter_ctx,local_vec2,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
850674ae819SStefano Zampini   ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
851674ae819SStefano Zampini   ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
852674ae819SStefano Zampini   ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
853674ae819SStefano Zampini   for (i=0;i<graph->nvtxs;i++) {
8543c629590SStefano Zampini     if (PetscRealPart(array[i]) > 0.1) {
855df48933bSStefano Zampini       ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
856674ae819SStefano Zampini       graph->subset[i] = 0; /* dirichlet nodes treated as internal -> is it ok? */
8570cf82fd6SStefano Zampini       graph->special_dof[i] = PCBDDCGRAPH_DIRICHLET_MARK;
858674ae819SStefano Zampini     }
859674ae819SStefano Zampini   }
860674ae819SStefano Zampini   ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
86151b0f6cfSStefano Zampini 
86208a5cf49SStefano Zampini   /* mark local periodic nodes (if any) and adapt CSR graph (if any) */
86351b0f6cfSStefano Zampini   if (graph->mirrors) {
86451b0f6cfSStefano Zampini     for (i=0;i<graph->nvtxs;i++)
86551b0f6cfSStefano Zampini       if (graph->mirrors[i])
8660cf82fd6SStefano Zampini         graph->special_dof[i] = PCBDDCGRAPH_LOCAL_PERIODIC_MARK;
86751b0f6cfSStefano Zampini 
86808a5cf49SStefano Zampini     if (graph->xadj && graph->adjncy) {
86908a5cf49SStefano Zampini       PetscInt *new_xadj,*new_adjncy;
87051b0f6cfSStefano Zampini       /* sort CSR graph */
87151b0f6cfSStefano Zampini       for (i=0;i<graph->nvtxs;i++)
87251b0f6cfSStefano Zampini         ierr = PetscSortInt(graph->xadj[i+1]-graph->xadj[i],&graph->adjncy[graph->xadj[i]]);CHKERRQ(ierr);
87351b0f6cfSStefano Zampini 
87451b0f6cfSStefano Zampini       /* adapt local CSR graph in case of local periodicity */
87551b0f6cfSStefano Zampini       k=0;
87651b0f6cfSStefano Zampini       for (i=0;i<graph->nvtxs;i++)
87751b0f6cfSStefano Zampini         for (j=graph->xadj[i];j<graph->xadj[i+1];j++)
87851b0f6cfSStefano Zampini           k += graph->mirrors[graph->adjncy[j]];
87951b0f6cfSStefano Zampini 
880785e854fSJed Brown       ierr = PetscMalloc1((graph->nvtxs+1),&new_xadj);CHKERRQ(ierr);
881785e854fSJed Brown       ierr = PetscMalloc1((k+graph->xadj[graph->nvtxs]),&new_adjncy);CHKERRQ(ierr);
88251b0f6cfSStefano Zampini       new_xadj[0]=0;
88351b0f6cfSStefano Zampini       for (i=0;i<graph->nvtxs;i++) {
88451b0f6cfSStefano Zampini         k = graph->xadj[i+1]-graph->xadj[i];
88551b0f6cfSStefano Zampini         ierr = PetscMemcpy(&new_adjncy[new_xadj[i]],&graph->adjncy[graph->xadj[i]],k*sizeof(PetscInt));CHKERRQ(ierr);
88651b0f6cfSStefano Zampini         new_xadj[i+1]=new_xadj[i]+k;
88751b0f6cfSStefano Zampini         for (j=graph->xadj[i];j<graph->xadj[i+1];j++) {
88851b0f6cfSStefano Zampini           k = graph->mirrors[graph->adjncy[j]];
88951b0f6cfSStefano Zampini           ierr = PetscMemcpy(&new_adjncy[new_xadj[i+1]],graph->mirrors_set[graph->adjncy[j]],k*sizeof(PetscInt));CHKERRQ(ierr);
89051b0f6cfSStefano Zampini           new_xadj[i+1]+=k;
89151b0f6cfSStefano Zampini         }
89251b0f6cfSStefano Zampini         k = new_xadj[i+1]-new_xadj[i];
89351b0f6cfSStefano Zampini         ierr = PetscSortRemoveDupsInt(&k,&new_adjncy[new_xadj[i]]);CHKERRQ(ierr);
89451b0f6cfSStefano Zampini         new_xadj[i+1]=new_xadj[i]+k;
89551b0f6cfSStefano Zampini       }
89651b0f6cfSStefano Zampini       /* set new CSR into graph */
89751b0f6cfSStefano Zampini       ierr = PetscFree(graph->xadj);CHKERRQ(ierr);
89851b0f6cfSStefano Zampini       ierr = PetscFree(graph->adjncy);CHKERRQ(ierr);
89951b0f6cfSStefano Zampini       graph->xadj = new_xadj;
90051b0f6cfSStefano Zampini       graph->adjncy = new_adjncy;
90151b0f6cfSStefano Zampini     }
90208a5cf49SStefano Zampini   }
90351b0f6cfSStefano Zampini 
904674ae819SStefano Zampini   /* mark special nodes -> each will become a single node equivalence class */
9059b70a373SStefano Zampini   ierr = VecSet(local_vec,0.0);CHKERRQ(ierr);
906674ae819SStefano Zampini   if (custom_primal_vertices) {
9079b70a373SStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
90863602bcaSStefano Zampini     ierr = ISGetLocalSize(custom_primal_vertices,&is_size);CHKERRQ(ierr);
909674ae819SStefano Zampini     ierr = ISGetIndices(custom_primal_vertices,(const PetscInt**)&is_indices);CHKERRQ(ierr);
910674ae819SStefano Zampini     for (i=0;i<is_size;i++){
9119b70a373SStefano Zampini       if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
9129b70a373SStefano Zampini         array[is_indices[i]] = 1.0;
9139b70a373SStefano Zampini       }
914674ae819SStefano Zampini     }
915674ae819SStefano Zampini     ierr = ISRestoreIndices(custom_primal_vertices,(const PetscInt**)&is_indices);CHKERRQ(ierr);
9169b70a373SStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
917674ae819SStefano Zampini   }
9189b70a373SStefano Zampini   /* special nodes: impose consistency among neighbours */
9199b70a373SStefano Zampini   ierr = VecSet(global_vec,0.0);CHKERRQ(ierr);
9209b70a373SStefano Zampini   ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9219b70a373SStefano Zampini   ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9229b70a373SStefano Zampini   ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9239b70a373SStefano Zampini   ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9249b70a373SStefano Zampini   ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
9259b70a373SStefano Zampini   j = 0;
9269b70a373SStefano Zampini   for (i=0;i<graph->nvtxs;i++) {
9279b70a373SStefano Zampini     if (PetscRealPart(array[i]) > 0.1 && graph->special_dof[i] != PCBDDCGRAPH_DIRICHLET_MARK) {
9289b70a373SStefano Zampini       graph->special_dof[i] = PCBDDCGRAPH_SPECIAL_MARK-j;
9299b70a373SStefano Zampini       j++;
9309b70a373SStefano Zampini     }
9319b70a373SStefano Zampini   }
9329b70a373SStefano Zampini   ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
9339b70a373SStefano Zampini 
934674ae819SStefano Zampini   /* mark interior nodes as touched and belonging to partition number 0 */
935674ae819SStefano Zampini   for (i=0;i<graph->nvtxs;i++) {
936674ae819SStefano Zampini     if (!graph->count[i]) {
937df48933bSStefano Zampini       ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
938674ae819SStefano Zampini       graph->subset[i] = 0;
939674ae819SStefano Zampini     }
940674ae819SStefano Zampini   }
941674ae819SStefano Zampini   /* init graph structure and compute default subsets */
942674ae819SStefano Zampini   nodes_touched=0;
943674ae819SStefano Zampini   for (i=0;i<graph->nvtxs;i++) {
944df48933bSStefano Zampini     if (PetscBTLookup(graph->touched,i)) {
945674ae819SStefano Zampini       nodes_touched++;
946674ae819SStefano Zampini     }
947674ae819SStefano Zampini   }
948674ae819SStefano Zampini   i = 0;
949674ae819SStefano Zampini   graph->ncc = 0;
950674ae819SStefano Zampini   total_counts = 0;
951674ae819SStefano Zampini   while (nodes_touched<graph->nvtxs) {
952674ae819SStefano Zampini     /*  find first untouched node in local ordering */
953df48933bSStefano Zampini     while (PetscBTLookup(graph->touched,i)) {
954674ae819SStefano Zampini       i++;
955674ae819SStefano Zampini     }
956df48933bSStefano Zampini     ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
957674ae819SStefano Zampini     graph->subset[i] = graph->ncc+1;
958674ae819SStefano Zampini     graph->cptr[graph->ncc] = total_counts;
959674ae819SStefano Zampini     graph->queue[total_counts] = i;
960674ae819SStefano Zampini     total_counts++;
961674ae819SStefano Zampini     nodes_touched++;
962674ae819SStefano Zampini     /* now find all other nodes having the same set of sharing subdomains */
963674ae819SStefano Zampini     for (j=i+1;j<graph->nvtxs;j++) {
96474e413f5SStefano Zampini       /* check for same number of sharing subdomains, dof number and same special mark */
965df48933bSStefano 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]) {
966674ae819SStefano Zampini         /* check for same set of sharing subdomains */
967674ae819SStefano Zampini         same_set=PETSC_TRUE;
968674ae819SStefano Zampini         for (k=0;k<graph->count[j];k++){
969674ae819SStefano Zampini           if (graph->neighbours_set[i][k]!=graph->neighbours_set[j][k]) {
970674ae819SStefano Zampini             same_set=PETSC_FALSE;
971674ae819SStefano Zampini           }
972674ae819SStefano Zampini         }
973674ae819SStefano Zampini         /* I found a friend of mine */
974674ae819SStefano Zampini         if (same_set) {
975674ae819SStefano Zampini           graph->subset[j]=graph->ncc+1;
976df48933bSStefano Zampini           ierr = PetscBTSet(graph->touched,j);CHKERRQ(ierr);
977674ae819SStefano Zampini           nodes_touched++;
978674ae819SStefano Zampini           graph->queue[total_counts] = j;
979674ae819SStefano Zampini           total_counts++;
980674ae819SStefano Zampini         }
981674ae819SStefano Zampini       }
982674ae819SStefano Zampini     }
983674ae819SStefano Zampini     graph->ncc++;
984674ae819SStefano Zampini   }
985674ae819SStefano Zampini   /* set default number of subsets (at this point no info on csr graph has been taken into account, so n_subsets = ncc */
986674ae819SStefano Zampini   graph->n_subsets = graph->ncc;
987785e854fSJed Brown   ierr = PetscMalloc1(graph->n_subsets,&graph->subset_ncc);CHKERRQ(ierr);
988674ae819SStefano Zampini   for (i=0;i<graph->n_subsets;i++) {
989674ae819SStefano Zampini     graph->subset_ncc[i] = 1;
990674ae819SStefano Zampini   }
991674ae819SStefano Zampini   /* final pointer */
992674ae819SStefano Zampini   graph->cptr[graph->ncc] = total_counts;
993674ae819SStefano Zampini   /* free memory allocated by ISLocalToGlobalMappingGetInfo */
994674ae819SStefano Zampini   ierr = ISLocalToGlobalMappingRestoreInfo(graph->l2gmap,&n_neigh,&neigh,&n_shared,&shared);CHKERRQ(ierr);
9953a5b03d1SStefano Zampini   /* get a reference node (min index in global ordering) for each subset */
9963a5b03d1SStefano Zampini   ierr = PetscMalloc1(graph->ncc,&graph->subset_ref_node);CHKERRQ(ierr);
9973a5b03d1SStefano Zampini   ierr = PetscMalloc1(graph->cptr[graph->ncc],&queue_global);CHKERRQ(ierr);
9983a5b03d1SStefano Zampini   ierr = ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_global);CHKERRQ(ierr);
9993a5b03d1SStefano Zampini   for (i=0;i<graph->ncc;i++) {
10003a5b03d1SStefano Zampini     PetscInt minval = queue_global[graph->cptr[i]];
10013a5b03d1SStefano Zampini     for (j=graph->cptr[i]+1;j<graph->cptr[i+1];j++) {
10023a5b03d1SStefano Zampini       minval = PetscMin(minval,queue_global[j]);
10033a5b03d1SStefano Zampini     }
10043a5b03d1SStefano Zampini     graph->subset_ref_node[i] = minval;
10053a5b03d1SStefano Zampini   }
10063a5b03d1SStefano Zampini   ierr = PetscFree(queue_global);CHKERRQ(ierr);
1007674ae819SStefano Zampini   /* free objects */
1008674ae819SStefano Zampini   ierr = VecDestroy(&local_vec);CHKERRQ(ierr);
1009674ae819SStefano Zampini   ierr = VecDestroy(&local_vec2);CHKERRQ(ierr);
1010674ae819SStefano Zampini   ierr = VecDestroy(&global_vec);CHKERRQ(ierr);
1011674ae819SStefano Zampini   ierr = VecScatterDestroy(&scatter_ctx);CHKERRQ(ierr);
1012674ae819SStefano Zampini   PetscFunctionReturn(0);
1013674ae819SStefano Zampini }
1014674ae819SStefano Zampini 
1015674ae819SStefano Zampini #undef __FUNCT__
1016674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphResetCSR"
1017674ae819SStefano Zampini PetscErrorCode PCBDDCGraphResetCSR(PCBDDCGraph graph)
1018674ae819SStefano Zampini {
1019674ae819SStefano Zampini   PetscErrorCode ierr;
1020674ae819SStefano Zampini 
1021674ae819SStefano Zampini   PetscFunctionBegin;
1022674ae819SStefano Zampini   ierr = PetscFree(graph->xadj);CHKERRQ(ierr);
1023674ae819SStefano Zampini   ierr = PetscFree(graph->adjncy);CHKERRQ(ierr);
1024575ad6abSStefano Zampini   graph->nvtxs_csr = 0;
1025674ae819SStefano Zampini   PetscFunctionReturn(0);
1026674ae819SStefano Zampini }
1027674ae819SStefano Zampini 
1028674ae819SStefano Zampini #undef __FUNCT__
1029674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphReset"
1030674ae819SStefano Zampini PetscErrorCode PCBDDCGraphReset(PCBDDCGraph graph)
1031674ae819SStefano Zampini {
1032674ae819SStefano Zampini   PetscErrorCode ierr;
1033674ae819SStefano Zampini 
1034674ae819SStefano Zampini   PetscFunctionBegin;
1035674ae819SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&graph->l2gmap);CHKERRQ(ierr);
1036674ae819SStefano Zampini   ierr = PetscFree(graph->subset_ncc);CHKERRQ(ierr);
10373a5b03d1SStefano Zampini   ierr = PetscFree(graph->subset_ref_node);CHKERRQ(ierr);
1038674ae819SStefano Zampini   if (graph->nvtxs) {
1039674ae819SStefano Zampini     ierr = PetscFree(graph->neighbours_set[0]);CHKERRQ(ierr);
1040674ae819SStefano Zampini   }
1041df48933bSStefano Zampini   ierr = PetscBTDestroy(&graph->touched);CHKERRQ(ierr);
1042df48933bSStefano Zampini   ierr = PetscFree7(graph->count,
1043674ae819SStefano Zampini                     graph->neighbours_set,
1044674ae819SStefano Zampini                     graph->subset,
1045674ae819SStefano Zampini                     graph->which_dof,
1046674ae819SStefano Zampini                     graph->cptr,
1047df48933bSStefano Zampini                     graph->queue,
1048df48933bSStefano Zampini                     graph->special_dof);CHKERRQ(ierr);
104951b0f6cfSStefano Zampini   if (graph->mirrors) {
105051b0f6cfSStefano Zampini     ierr = PetscFree(graph->mirrors_set[0]);CHKERRQ(ierr);
105151b0f6cfSStefano Zampini   }
105251b0f6cfSStefano Zampini   ierr = PetscFree2(graph->mirrors,graph->mirrors_set);CHKERRQ(ierr);
1053674ae819SStefano Zampini   graph->nvtxs = 0;
1054674ae819SStefano Zampini   graph->n_subsets = 0;
1055674ae819SStefano Zampini   graph->custom_minimal_size = 1;
1056674ae819SStefano Zampini   PetscFunctionReturn(0);
1057674ae819SStefano Zampini }
1058674ae819SStefano Zampini 
1059674ae819SStefano Zampini #undef __FUNCT__
1060674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphInit"
1061674ae819SStefano Zampini PetscErrorCode PCBDDCGraphInit(PCBDDCGraph graph, ISLocalToGlobalMapping l2gmap)
1062674ae819SStefano Zampini {
1063674ae819SStefano Zampini   PetscInt       n;
1064674ae819SStefano Zampini   PetscErrorCode ierr;
1065674ae819SStefano Zampini 
1066674ae819SStefano Zampini   PetscFunctionBegin;
1067674ae819SStefano Zampini   PetscValidPointer(graph,1);
1068674ae819SStefano Zampini   PetscValidHeaderSpecific(l2gmap,IS_LTOGM_CLASSID,2);
10698e61c736SStefano Zampini   /* raise an error if already allocated */
1070674ae819SStefano Zampini   if (graph->nvtxs) {
10718e61c736SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)l2gmap),PETSC_ERR_PLIB,"BDDCGraph already initialized");
1072674ae819SStefano Zampini   }
1073674ae819SStefano Zampini   /* set number of vertices */
1074674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)l2gmap);CHKERRQ(ierr);
1075674ae819SStefano Zampini   graph->l2gmap = l2gmap;
1076674ae819SStefano Zampini   ierr = ISLocalToGlobalMappingGetSize(l2gmap,&n);CHKERRQ(ierr);
1077674ae819SStefano Zampini   graph->nvtxs = n;
1078674ae819SStefano Zampini   /* allocate used space */
1079df48933bSStefano Zampini   ierr = PetscBTCreate(graph->nvtxs,&graph->touched);CHKERRQ(ierr);
10808e5aaad5SJed Brown   ierr = PetscMalloc7(graph->nvtxs,&graph->count,
10818e5aaad5SJed Brown                       graph->nvtxs,&graph->neighbours_set,
10828e5aaad5SJed Brown                       graph->nvtxs,&graph->subset,
10838e5aaad5SJed Brown                       graph->nvtxs,&graph->which_dof,
10848e5aaad5SJed Brown                       graph->nvtxs+1,&graph->cptr,
10858e5aaad5SJed Brown                       graph->nvtxs,&graph->queue,
10868e5aaad5SJed Brown                       graph->nvtxs,&graph->special_dof);CHKERRQ(ierr);
1087674ae819SStefano Zampini   /* zeroes memory */
1088674ae819SStefano Zampini   ierr = PetscMemzero(graph->count,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
1089674ae819SStefano Zampini   ierr = PetscMemzero(graph->subset,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
109063602bcaSStefano Zampini   /* use -1 as a default value for which_dof array */
109163602bcaSStefano Zampini   for (n=0;n<graph->nvtxs;n++) graph->which_dof[n] = -1;
1092674ae819SStefano Zampini   ierr = PetscMemzero(graph->cptr,(graph->nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr);
1093674ae819SStefano Zampini   ierr = PetscMemzero(graph->queue,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
109474e413f5SStefano Zampini   ierr = PetscMemzero(graph->special_dof,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
1095674ae819SStefano Zampini   /* zeroes first pointer to neighbour set */
1096674ae819SStefano Zampini   if (graph->nvtxs) {
1097674ae819SStefano Zampini     graph->neighbours_set[0] = 0;
1098674ae819SStefano Zampini   }
1099674ae819SStefano Zampini   /* zeroes workspace for values of ncc */
1100674ae819SStefano Zampini   graph->subset_ncc = 0;
11013a5b03d1SStefano Zampini   graph->subset_ref_node = 0;
1102674ae819SStefano Zampini   PetscFunctionReturn(0);
1103674ae819SStefano Zampini }
1104674ae819SStefano Zampini 
1105674ae819SStefano Zampini #undef __FUNCT__
1106674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphDestroy"
1107674ae819SStefano Zampini PetscErrorCode PCBDDCGraphDestroy(PCBDDCGraph* graph)
1108674ae819SStefano Zampini {
1109674ae819SStefano Zampini   PetscErrorCode ierr;
1110674ae819SStefano Zampini 
1111674ae819SStefano Zampini   PetscFunctionBegin;
1112674ae819SStefano Zampini   ierr = PCBDDCGraphReset(*graph);CHKERRQ(ierr);
1113674ae819SStefano Zampini   ierr = PetscFree(*graph);CHKERRQ(ierr);
1114674ae819SStefano Zampini   PetscFunctionReturn(0);
1115674ae819SStefano Zampini }
1116674ae819SStefano Zampini 
1117674ae819SStefano Zampini #undef __FUNCT__
1118674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphCreate"
1119674ae819SStefano Zampini PetscErrorCode PCBDDCGraphCreate(PCBDDCGraph *graph)
1120674ae819SStefano Zampini {
1121674ae819SStefano Zampini   PCBDDCGraph    new_graph;
1122674ae819SStefano Zampini   PetscErrorCode ierr;
1123674ae819SStefano Zampini 
1124674ae819SStefano Zampini   PetscFunctionBegin;
1125674ae819SStefano Zampini   ierr = PetscMalloc(sizeof(*new_graph),&new_graph);CHKERRQ(ierr);
1126674ae819SStefano Zampini   /* local to global mapping of dofs */
1127674ae819SStefano Zampini   new_graph->l2gmap = 0;
1128674ae819SStefano Zampini   /* vertex size */
1129674ae819SStefano Zampini   new_graph->nvtxs = 0;
1130674ae819SStefano Zampini   new_graph->n_subsets = 0;
1131674ae819SStefano Zampini   new_graph->custom_minimal_size = 1;
1132674ae819SStefano Zampini   /* zeroes ponters */
113351b0f6cfSStefano Zampini   new_graph->mirrors = 0;
113451b0f6cfSStefano Zampini   new_graph->mirrors_set = 0;
1135674ae819SStefano Zampini   new_graph->neighbours_set = 0;
1136674ae819SStefano Zampini   new_graph->subset = 0;
1137674ae819SStefano Zampini   new_graph->which_dof = 0;
113874e413f5SStefano Zampini   new_graph->special_dof = 0;
1139674ae819SStefano Zampini   new_graph->cptr = 0;
1140674ae819SStefano Zampini   new_graph->queue = 0;
1141674ae819SStefano Zampini   new_graph->count = 0;
1142674ae819SStefano Zampini   new_graph->subset_ncc = 0;
11433a5b03d1SStefano Zampini   new_graph->subset_ref_node = 0;
1144674ae819SStefano Zampini   new_graph->touched = 0;
1145674ae819SStefano Zampini   /* zeroes pointers to csr graph of local nodes connectivity (optional data) */
1146575ad6abSStefano Zampini   new_graph->nvtxs_csr = 0;
1147674ae819SStefano Zampini   new_graph->xadj = 0;
1148674ae819SStefano Zampini   new_graph->adjncy = 0;
1149674ae819SStefano Zampini   *graph = new_graph;
1150674ae819SStefano Zampini   PetscFunctionReturn(0);
1151674ae819SStefano Zampini }
1152