xref: /petsc/src/ksp/pc/impls/bddc/bddcgraph.c (revision 67c9da694bc355307fef111169dbc0aca808e1d6)
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++) {
71674ae819SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  %d (neighs:",i);CHKERRQ(ierr);
722b510759SStefano Zampini     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
73674ae819SStefano Zampini     PetscInt node_num=graph->queue[graph->cptr[i]];
74674ae819SStefano Zampini     for (j=0;j<graph->count[node_num];j++) {
75674ae819SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d",graph->neighbours_set[node_num][j]);CHKERRQ(ierr);
76674ae819SStefano Zampini     }
77674ae819SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"):");CHKERRQ(ierr);
78674ae819SStefano Zampini     for (j=graph->cptr[i];j<graph->cptr[i+1];j++) {
7993fb5fd3SStefano Zampini       ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d (%d)",graph->queue[j],queue_in_global_numbering[j]);CHKERRQ(ierr);
80674ae819SStefano Zampini     }
81674ae819SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
822b510759SStefano Zampini     ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr);
832b510759SStefano Zampini     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
84674ae819SStefano Zampini   }
8593fb5fd3SStefano Zampini   ierr = PetscFree(queue_in_global_numbering);CHKERRQ(ierr);
86e49050b4SStefano Zampini   if (graph->custom_minimal_size > 1 && verbosity_level > 1) {
87674ae819SStefano Zampini     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Custom minimal size %d\n",graph->custom_minimal_size);CHKERRQ(ierr);
88674ae819SStefano Zampini   }
89674ae819SStefano Zampini   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
90674ae819SStefano Zampini   PetscFunctionReturn(0);
91674ae819SStefano Zampini }
92674ae819SStefano Zampini 
93674ae819SStefano Zampini #undef __FUNCT__
94674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphGetCandidatesIS"
95674ae819SStefano Zampini PetscErrorCode PCBDDCGraphGetCandidatesIS(PCBDDCGraph graph, PetscBool use_faces, PetscBool use_edges, PetscBool use_vertices, PetscInt *n_faces, IS *FacesIS[], PetscInt *n_edges, IS *EdgesIS[], IS *VerticesIS)
96674ae819SStefano Zampini {
97674ae819SStefano Zampini   IS             *ISForFaces,*ISForEdges,ISForVertices;
98674ae819SStefano Zampini   PetscInt       i,j,nfc,nec,nvc,*idx;
99674ae819SStefano Zampini   PetscBool      twodim_flag;
100674ae819SStefano Zampini   PetscErrorCode ierr;
101674ae819SStefano Zampini 
102674ae819SStefano Zampini   PetscFunctionBegin;
103674ae819SStefano Zampini   /* loop on ccs to evalute number of faces, edges and vertices */
104674ae819SStefano Zampini   nfc = 0;
105674ae819SStefano Zampini   nec = 0;
106674ae819SStefano Zampini   nvc = 0;
107674ae819SStefano Zampini   twodim_flag = PETSC_FALSE;
108674ae819SStefano Zampini   for (i=0;i<graph->ncc;i++) {
109674ae819SStefano Zampini     if (graph->cptr[i+1]-graph->cptr[i] > graph->custom_minimal_size) {
1100cf82fd6SStefano Zampini       if (graph->count[graph->queue[graph->cptr[i]]] == 1 && graph->special_dof[graph->queue[graph->cptr[i]]] != PCBDDCGRAPH_NEUMANN_MARK) {
111674ae819SStefano Zampini         nfc++;
112674ae819SStefano Zampini       } else { /* note that nec will be zero in 2d */
113674ae819SStefano Zampini         nec++;
114674ae819SStefano Zampini       }
115674ae819SStefano Zampini     } else {
116674ae819SStefano Zampini       nvc += graph->cptr[i+1]-graph->cptr[i];
117674ae819SStefano Zampini     }
118674ae819SStefano Zampini   }
119aa9fcddfSStefano Zampini   j=0;
120aa9fcddfSStefano Zampini   ierr = MPI_Allreduce(&nec,&j,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)graph->l2gmap));CHKERRQ(ierr);
121aa9fcddfSStefano Zampini   if (!j) { /* we are in a 2D case -> no faces, only edges */
122674ae819SStefano Zampini     nec = nfc;
123674ae819SStefano Zampini     nfc = 0;
124674ae819SStefano Zampini     twodim_flag = PETSC_TRUE;
125674ae819SStefano Zampini   }
126674ae819SStefano Zampini   /* allocate IS arrays for faces, edges. Vertices need a single index set. */
127674ae819SStefano Zampini   ISForFaces = 0;
128674ae819SStefano Zampini   ISForEdges = 0;
129674ae819SStefano Zampini   ISForVertices = 0;
1301ee8691aSStefano Zampini   if (use_faces && nfc) {
131785e854fSJed Brown     ierr = PetscMalloc1(nfc,&ISForFaces);CHKERRQ(ierr);
132674ae819SStefano Zampini   }
133b8ffe317SStefano Zampini   if (use_edges && nec) {
134785e854fSJed Brown     ierr = PetscMalloc1(nec,&ISForEdges);CHKERRQ(ierr);
135674ae819SStefano Zampini   }
136b8ffe317SStefano Zampini   if (use_vertices && nvc) {
137785e854fSJed Brown     ierr = PetscMalloc1(nvc,&idx);CHKERRQ(ierr);
138674ae819SStefano Zampini   }
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++) {
143674ae819SStefano Zampini     if (graph->cptr[i+1]-graph->cptr[i] > graph->custom_minimal_size) {
1440cf82fd6SStefano Zampini       if (graph->count[graph->queue[graph->cptr[i]]] == 1 && graph->special_dof[graph->queue[graph->cptr[i]]] != PCBDDCGRAPH_NEUMANN_MARK) {
145674ae819SStefano Zampini         if (twodim_flag) {
146674ae819SStefano Zampini           if (use_edges) {
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           }
150674ae819SStefano Zampini         } else {
151674ae819SStefano Zampini           if (use_faces) {
152674ae819SStefano 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);
153674ae819SStefano Zampini             nfc++;
154674ae819SStefano Zampini           }
155674ae819SStefano Zampini         }
156674ae819SStefano Zampini       } else {
157674ae819SStefano Zampini         if (use_edges) {
158674ae819SStefano 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);
159674ae819SStefano Zampini           nec++;
160674ae819SStefano Zampini         }
161674ae819SStefano Zampini       }
162674ae819SStefano Zampini     }
163674ae819SStefano Zampini   }
164674ae819SStefano Zampini   /* index set for vertices */
165b8ffe317SStefano Zampini   if (use_vertices && nvc) {
166b8ffe317SStefano Zampini     nvc = 0;
167674ae819SStefano Zampini     for (i=0;i<graph->ncc;i++) {
168674ae819SStefano Zampini       if (graph->cptr[i+1]-graph->cptr[i] <= graph->custom_minimal_size) {
169674ae819SStefano Zampini         for (j=graph->cptr[i];j<graph->cptr[i+1];j++) {
170674ae819SStefano Zampini           idx[nvc]=graph->queue[j];
171674ae819SStefano Zampini           nvc++;
172674ae819SStefano Zampini         }
173674ae819SStefano Zampini       }
174674ae819SStefano Zampini     }
175674ae819SStefano Zampini     /* sort vertex set (by local ordering) */
176674ae819SStefano Zampini     ierr = PetscSortInt(nvc,idx);CHKERRQ(ierr);
177674ae819SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,nvc,idx,PETSC_OWN_POINTER,&ISForVertices);CHKERRQ(ierr);
178674ae819SStefano Zampini   }
179674ae819SStefano Zampini   /* get back info */
180674ae819SStefano Zampini   *n_faces = nfc;
181674ae819SStefano Zampini   *FacesIS = ISForFaces;
182674ae819SStefano Zampini   *n_edges = nec;
183674ae819SStefano Zampini   *EdgesIS = ISForEdges;
184674ae819SStefano Zampini   *VerticesIS = ISForVertices;
185674ae819SStefano Zampini   PetscFunctionReturn(0);
186674ae819SStefano Zampini }
187674ae819SStefano Zampini 
188674ae819SStefano Zampini #undef __FUNCT__
189674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphComputeConnectedComponents"
190674ae819SStefano Zampini PetscErrorCode PCBDDCGraphComputeConnectedComponents(PCBDDCGraph graph)
191674ae819SStefano Zampini {
1925b1b9aeaSStefano Zampini   PetscInt       i,adapt_interface,adapt_interface_reduced;
193674ae819SStefano Zampini   MPI_Comm       interface_comm;
194674ae819SStefano Zampini   PetscErrorCode ierr;
195674ae819SStefano Zampini 
196674ae819SStefano Zampini   PetscFunctionBegin;
197674ae819SStefano Zampini   /* compute connected components locally */
198674ae819SStefano Zampini   ierr = PetscObjectGetComm((PetscObject)(graph->l2gmap),&interface_comm);CHKERRQ(ierr);
199674ae819SStefano Zampini   ierr = PCBDDCGraphComputeConnectedComponentsLocal(graph);CHKERRQ(ierr);
200674ae819SStefano Zampini   /* check consistency of connected components among neighbouring subdomains -> it adapt them in case it is needed */
201674ae819SStefano Zampini   adapt_interface = 0;
202674ae819SStefano Zampini   adapt_interface_reduced = 0;
203674ae819SStefano Zampini   for (i=0;i<graph->n_subsets;i++) {
204674ae819SStefano Zampini     /* We are not sure that on a given subset of the local interface,
205674ae819SStefano Zampini        with two connected components, the latters be the same among sharing subdomains */
206674ae819SStefano Zampini     if (graph->subset_ncc[i] > 1) {
207674ae819SStefano Zampini       adapt_interface=1;
208674ae819SStefano Zampini       break;
209674ae819SStefano Zampini     }
210674ae819SStefano Zampini   }
211674ae819SStefano Zampini   ierr = MPI_Allreduce(&adapt_interface,&adapt_interface_reduced,1,MPIU_INT,MPI_LOR,interface_comm);CHKERRQ(ierr);
212674ae819SStefano Zampini 
213674ae819SStefano Zampini   if (graph->n_subsets && adapt_interface_reduced) {
2145b1b9aeaSStefano Zampini     MPI_Request *send_requests;
2155b1b9aeaSStefano Zampini     MPI_Request *recv_requests;
2165b1b9aeaSStefano Zampini     PetscInt    *aux_new_xadj,*new_xadj,*new_adjncy,**temp_buffer;
2175b1b9aeaSStefano Zampini     PetscInt    *old_xadj,*old_adjncy;
2185b1b9aeaSStefano Zampini     PetscInt    j,k,s,sum_requests,buffer_size,size_of_recv,temp_buffer_size;
2195b1b9aeaSStefano Zampini     PetscMPIInt rank,neigh,tag,mpi_buffer_size;
2205b1b9aeaSStefano Zampini     PetscInt    *cum_recv_counts,*subset_to_nodes_indices,*recv_buffer_subset,*nodes_to_temp_buffer_indices;
2215b1b9aeaSStefano Zampini     PetscInt    *send_buffer,*recv_buffer,*queue_in_global_numbering,*sizes_of_sends,*add_to_subset;
2225b1b9aeaSStefano Zampini     PetscInt    start_of_recv,start_of_send,size_of_send,global_subset_counter,ins_val;
2235b1b9aeaSStefano Zampini     PetscBool   *subset_cc_adapt,same_set;
2245b1b9aeaSStefano Zampini 
225674ae819SStefano Zampini     /* Retrict adjacency graph using information from previously computed connected components */
226785e854fSJed Brown     ierr = PetscMalloc1(graph->nvtxs,&aux_new_xadj);CHKERRQ(ierr);
227674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
228674ae819SStefano Zampini       aux_new_xadj[i]=1;
229674ae819SStefano Zampini     }
230674ae819SStefano Zampini     for (i=0;i<graph->ncc;i++) {
231674ae819SStefano Zampini       k = graph->cptr[i+1]-graph->cptr[i];
232674ae819SStefano Zampini       for (j=0;j<k;j++) {
233674ae819SStefano Zampini         aux_new_xadj[graph->queue[graph->cptr[i]+j]]=k;
234674ae819SStefano Zampini       }
235674ae819SStefano Zampini     }
236674ae819SStefano Zampini     j = 0;
237674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
238674ae819SStefano Zampini       j += aux_new_xadj[i];
239674ae819SStefano Zampini     }
240785e854fSJed Brown     ierr = PetscMalloc1((graph->nvtxs+1),&new_xadj);CHKERRQ(ierr);
241785e854fSJed Brown     ierr = PetscMalloc1(j,&new_adjncy);CHKERRQ(ierr);
242674ae819SStefano Zampini     new_xadj[0]=0;
243674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
244674ae819SStefano Zampini       new_xadj[i+1]=new_xadj[i]+aux_new_xadj[i];
245674ae819SStefano Zampini       if (aux_new_xadj[i]==1) {
246674ae819SStefano Zampini         new_adjncy[new_xadj[i]]=i;
247674ae819SStefano Zampini       }
248674ae819SStefano Zampini     }
249674ae819SStefano Zampini     ierr = PetscFree(aux_new_xadj);CHKERRQ(ierr);
250674ae819SStefano Zampini     for (i=0;i<graph->ncc;i++) {
251674ae819SStefano Zampini       k = graph->cptr[i+1]-graph->cptr[i];
252674ae819SStefano Zampini       for (j=0;j<k;j++) {
253674ae819SStefano Zampini         ierr = PetscMemcpy(&new_adjncy[new_xadj[graph->queue[graph->cptr[i]+j]]],&graph->queue[graph->cptr[i]],k*sizeof(PetscInt));CHKERRQ(ierr);
254674ae819SStefano Zampini       }
255674ae819SStefano Zampini     }
2565b1b9aeaSStefano Zampini     /* set temporarly new CSR into graph */
2575b1b9aeaSStefano Zampini     old_xadj = graph->xadj;
2585b1b9aeaSStefano Zampini     old_adjncy = graph->adjncy;
259674ae819SStefano Zampini     graph->xadj = new_xadj;
260674ae819SStefano Zampini     graph->adjncy = new_adjncy;
261674ae819SStefano Zampini     /* allocate some space */
262674ae819SStefano Zampini     ierr = MPI_Comm_rank(interface_comm,&rank);CHKERRQ(ierr);
263785e854fSJed Brown     ierr = PetscMalloc1((graph->n_subsets+1),&cum_recv_counts);CHKERRQ(ierr);
264674ae819SStefano Zampini     ierr = PetscMemzero(cum_recv_counts,(graph->n_subsets+1)*sizeof(*cum_recv_counts));CHKERRQ(ierr);
265785e854fSJed Brown     ierr = PetscMalloc1(graph->n_subsets,&subset_to_nodes_indices);CHKERRQ(ierr);
266674ae819SStefano Zampini     /* first count how many neighbours per connected component I will receive from */
267674ae819SStefano Zampini     cum_recv_counts[0]=0;
268674ae819SStefano Zampini     for (i=1;i<graph->n_subsets+1;i++) {
269674ae819SStefano Zampini       j = 0;
270674ae819SStefano Zampini       while (graph->subset[j] != i) {
271674ae819SStefano Zampini         j++;
272674ae819SStefano Zampini       }
273674ae819SStefano Zampini       subset_to_nodes_indices[i-1]=j;
274674ae819SStefano Zampini       /* We don't want sends/recvs_to/from_self -> here I don't count myself  */
275674ae819SStefano Zampini       cum_recv_counts[i]=cum_recv_counts[i-1]+graph->count[j];
276674ae819SStefano Zampini     }
277785e854fSJed Brown     ierr = PetscMalloc1(2*cum_recv_counts[graph->n_subsets],&recv_buffer_subset);CHKERRQ(ierr);
278785e854fSJed Brown     ierr = PetscMalloc1(cum_recv_counts[graph->n_subsets],&send_requests);CHKERRQ(ierr);
279785e854fSJed Brown     ierr = PetscMalloc1(cum_recv_counts[graph->n_subsets],&recv_requests);CHKERRQ(ierr);
280674ae819SStefano Zampini     for (i=0;i<cum_recv_counts[graph->n_subsets];i++) {
281674ae819SStefano Zampini       send_requests[i]=MPI_REQUEST_NULL;
282674ae819SStefano Zampini       recv_requests[i]=MPI_REQUEST_NULL;
283674ae819SStefano Zampini     }
284674ae819SStefano Zampini     /* exchange with my neighbours the number of my connected components on the shared interface */
285674ae819SStefano Zampini     sum_requests = 0;
286674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
287674ae819SStefano Zampini       j = subset_to_nodes_indices[i];
2883a5b03d1SStefano Zampini       ierr = PetscMPIIntCast(graph->subset_ref_node[i],&tag);CHKERRQ(ierr);
289674ae819SStefano Zampini       for (k=0;k<graph->count[j];k++) {
290674ae819SStefano Zampini         ierr = PetscMPIIntCast(graph->neighbours_set[j][k],&neigh);CHKERRQ(ierr);
291674ae819SStefano Zampini         ierr = MPI_Isend(&graph->subset_ncc[i],1,MPIU_INT,neigh,tag,interface_comm,&send_requests[sum_requests]);CHKERRQ(ierr);
292674ae819SStefano Zampini         ierr = MPI_Irecv(&recv_buffer_subset[sum_requests],1,MPIU_INT,neigh,tag,interface_comm,&recv_requests[sum_requests]);CHKERRQ(ierr);
293674ae819SStefano Zampini         sum_requests++;
294674ae819SStefano Zampini       }
295674ae819SStefano Zampini     }
296674ae819SStefano Zampini     ierr = MPI_Waitall(sum_requests,recv_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
297674ae819SStefano Zampini     ierr = MPI_Waitall(sum_requests,send_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
298674ae819SStefano Zampini     /* determine the connected component I need to adapt */
299785e854fSJed Brown     ierr = PetscMalloc1(graph->n_subsets,&subset_cc_adapt);CHKERRQ(ierr);
300674ae819SStefano Zampini     ierr = PetscMemzero(subset_cc_adapt,graph->n_subsets*sizeof(*subset_cc_adapt));CHKERRQ(ierr);
301674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
302674ae819SStefano Zampini       for (j=cum_recv_counts[i];j<cum_recv_counts[i+1];j++){
303674ae819SStefano Zampini         /* The first condition is natural (someone has a different number of ccs than me), the second one is just to be safe */
304674ae819SStefano Zampini         if (graph->subset_ncc[i] != recv_buffer_subset[j] || graph->subset_ncc[i] > 1) {
305674ae819SStefano Zampini           subset_cc_adapt[i] = PETSC_TRUE;
306674ae819SStefano Zampini           break;
307674ae819SStefano Zampini         }
308674ae819SStefano Zampini       }
309674ae819SStefano Zampini     }
310674ae819SStefano Zampini     buffer_size = 0;
311674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
312674ae819SStefano Zampini       if (subset_cc_adapt[i]) {
313674ae819SStefano Zampini         for (j=i;j<graph->ncc;j++) {
314674ae819SStefano Zampini           if (graph->subset[graph->queue[graph->cptr[j]]] == i+1) { /* WARNING -> subset values goes from 1 to graph->n_subsets included */
315674ae819SStefano Zampini             buffer_size += 1 + graph->cptr[j+1]-graph->cptr[j];
316674ae819SStefano Zampini           }
317674ae819SStefano Zampini         }
318674ae819SStefano Zampini       }
319674ae819SStefano Zampini     }
320785e854fSJed Brown     ierr = PetscMalloc1(buffer_size,&send_buffer);CHKERRQ(ierr);
321674ae819SStefano Zampini     /* now get from neighbours their ccs (in global numbering) and adapt them (in case it is needed) */
322785e854fSJed Brown     ierr = PetscMalloc1(graph->cptr[graph->ncc],&queue_in_global_numbering);CHKERRQ(ierr);
323674ae819SStefano Zampini     ierr = ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_in_global_numbering);CHKERRQ(ierr);
324674ae819SStefano Zampini     /* determine how much data to send (size of each queue plus the global indices) and communicate it to neighbours */
325785e854fSJed Brown     ierr = PetscMalloc1(graph->n_subsets,&sizes_of_sends);CHKERRQ(ierr);
326674ae819SStefano Zampini     ierr = PetscMemzero(sizes_of_sends,graph->n_subsets*sizeof(*sizes_of_sends));CHKERRQ(ierr);
327674ae819SStefano Zampini     sum_requests = 0;
328674ae819SStefano Zampini     start_of_send = 0;
329674ae819SStefano Zampini     start_of_recv = cum_recv_counts[graph->n_subsets];
330674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
331674ae819SStefano Zampini       if (subset_cc_adapt[i]) {
332674ae819SStefano Zampini         size_of_send = 0;
333674ae819SStefano Zampini         for (j=i;j<graph->ncc;j++) {
334674ae819SStefano Zampini           if (graph->subset[graph->queue[graph->cptr[j]]] == i+1) { /* WARNING -> subset values goes from 1 to graph->n_subsets included */
335674ae819SStefano Zampini             send_buffer[start_of_send+size_of_send]=graph->cptr[j+1]-graph->cptr[j];
336674ae819SStefano Zampini             size_of_send += 1;
337674ae819SStefano Zampini             ierr = PetscMemcpy(&send_buffer[start_of_send+size_of_send],
338674ae819SStefano Zampini                                &queue_in_global_numbering[graph->cptr[j]],
339674ae819SStefano Zampini                                (graph->cptr[j+1]-graph->cptr[j])*sizeof(*send_buffer));CHKERRQ(ierr);
340674ae819SStefano Zampini             size_of_send = size_of_send+graph->cptr[j+1]-graph->cptr[j];
341674ae819SStefano Zampini           }
342674ae819SStefano Zampini         }
343674ae819SStefano Zampini         j = subset_to_nodes_indices[i];
344674ae819SStefano Zampini         sizes_of_sends[i] = size_of_send;
3453a5b03d1SStefano Zampini         ierr = PetscMPIIntCast(graph->subset_ref_node[i]+1,&tag);CHKERRQ(ierr);
346674ae819SStefano Zampini         for (k=0;k<graph->count[j];k++) {
347674ae819SStefano Zampini           ierr = PetscMPIIntCast(graph->neighbours_set[j][k],&neigh);CHKERRQ(ierr);
348674ae819SStefano Zampini           ierr = MPI_Isend(&sizes_of_sends[i],1,MPIU_INT,neigh,tag,interface_comm,&send_requests[sum_requests]);CHKERRQ(ierr);
349674ae819SStefano 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);
350674ae819SStefano Zampini           sum_requests++;
351674ae819SStefano Zampini         }
352674ae819SStefano Zampini         start_of_send += size_of_send;
353674ae819SStefano Zampini       }
354674ae819SStefano Zampini     }
355674ae819SStefano Zampini     ierr = MPI_Waitall(sum_requests,send_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
356674ae819SStefano Zampini     ierr = MPI_Waitall(sum_requests,recv_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
357674ae819SStefano Zampini     buffer_size = 0;
358674ae819SStefano Zampini     for (k=0;k<sum_requests;k++) {
359674ae819SStefano Zampini       buffer_size += recv_buffer_subset[start_of_recv+k];
360674ae819SStefano Zampini     }
361785e854fSJed Brown     ierr = PetscMalloc1(buffer_size,&recv_buffer);CHKERRQ(ierr);
362674ae819SStefano Zampini     /* now exchange the data */
363674ae819SStefano Zampini     start_of_recv = 0;
364674ae819SStefano Zampini     start_of_send = 0;
365674ae819SStefano Zampini     sum_requests = 0;
366674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
367674ae819SStefano Zampini       if (subset_cc_adapt[i]) {
368674ae819SStefano Zampini         size_of_send = sizes_of_sends[i];
369674ae819SStefano Zampini         j = subset_to_nodes_indices[i];
3703a5b03d1SStefano Zampini         ierr = PetscMPIIntCast(graph->subset_ref_node[i]+2,&tag);CHKERRQ(ierr);
371674ae819SStefano Zampini         for (k=0;k<graph->count[j];k++) {
372674ae819SStefano Zampini           ierr = PetscMPIIntCast(graph->neighbours_set[j][k],&neigh);CHKERRQ(ierr);
373674ae819SStefano Zampini           ierr = PetscMPIIntCast(size_of_send,&mpi_buffer_size);CHKERRQ(ierr);
374674ae819SStefano Zampini           ierr = MPI_Isend(&send_buffer[start_of_send],mpi_buffer_size,MPIU_INT,neigh,tag,interface_comm,&send_requests[sum_requests]);CHKERRQ(ierr);
375674ae819SStefano Zampini           size_of_recv = recv_buffer_subset[cum_recv_counts[graph->n_subsets]+sum_requests];
376674ae819SStefano Zampini           ierr = PetscMPIIntCast(size_of_recv,&mpi_buffer_size);CHKERRQ(ierr);
377674ae819SStefano Zampini           ierr = MPI_Irecv(&recv_buffer[start_of_recv],mpi_buffer_size,MPIU_INT,neigh,tag,interface_comm,&recv_requests[sum_requests]);CHKERRQ(ierr);
378674ae819SStefano Zampini           start_of_recv += size_of_recv;
379674ae819SStefano Zampini           sum_requests++;
380674ae819SStefano Zampini         }
381674ae819SStefano Zampini         start_of_send += size_of_send;
382674ae819SStefano Zampini       }
383674ae819SStefano Zampini     }
384674ae819SStefano Zampini     ierr = MPI_Waitall(sum_requests,recv_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
385674ae819SStefano Zampini     ierr = MPI_Waitall(sum_requests,send_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr);
386674ae819SStefano Zampini     for (j=0;j<buffer_size;) {
387674ae819SStefano Zampini        ierr = ISGlobalToLocalMappingApply(graph->l2gmap,IS_GTOLM_MASK,recv_buffer[j],&recv_buffer[j+1],&recv_buffer[j],&recv_buffer[j+1]);CHKERRQ(ierr);
38851b0f6cfSStefano Zampini        /* we need to adapt the output of GlobalToLocal mapping if there are mirrored nodes */
38951b0f6cfSStefano Zampini        if (graph->mirrors) {
39051b0f6cfSStefano Zampini          PetscBool mirrored_found=PETSC_FALSE;
39151b0f6cfSStefano Zampini          for (k=0;k<recv_buffer[j];k++) {
39251b0f6cfSStefano Zampini            if (graph->mirrors[recv_buffer[j+k+1]]) {
39351b0f6cfSStefano Zampini              mirrored_found=PETSC_TRUE;
39451b0f6cfSStefano Zampini              recv_buffer[j+k+1]=graph->mirrors_set[recv_buffer[j+k+1]][0];
39551b0f6cfSStefano Zampini            }
39651b0f6cfSStefano Zampini          }
39751b0f6cfSStefano Zampini          if (mirrored_found) {
39851b0f6cfSStefano Zampini            ierr = PetscSortInt(recv_buffer[j],&recv_buffer[j+1]);CHKERRQ(ierr);
39951b0f6cfSStefano Zampini            k=0;
40051b0f6cfSStefano Zampini            while (k<recv_buffer[j]) {
40151b0f6cfSStefano Zampini              for (s=1;s<graph->mirrors[recv_buffer[j+1+k]];s++) {
40251b0f6cfSStefano Zampini                recv_buffer[j+1+k+s] = graph->mirrors_set[recv_buffer[j+1+k]][s];
40351b0f6cfSStefano Zampini              }
40451b0f6cfSStefano Zampini              k+=graph->mirrors[recv_buffer[j+1+k]]+s;
40551b0f6cfSStefano Zampini            }
40651b0f6cfSStefano Zampini          }
40751b0f6cfSStefano Zampini        }
408674ae819SStefano Zampini        k = recv_buffer[j]+1;
409674ae819SStefano Zampini        j += k;
410674ae819SStefano Zampini     }
411674ae819SStefano Zampini     sum_requests = cum_recv_counts[graph->n_subsets];
412674ae819SStefano Zampini     start_of_recv = 0;
413785e854fSJed Brown     ierr = PetscMalloc1(graph->nvtxs,&nodes_to_temp_buffer_indices);CHKERRQ(ierr);
414674ae819SStefano Zampini     global_subset_counter = 0;
415674ae819SStefano Zampini     for (i=0;i<graph->n_subsets;i++) {
416674ae819SStefano Zampini       if (subset_cc_adapt[i]) {
417674ae819SStefano Zampini         temp_buffer_size = 0;
418674ae819SStefano Zampini         /* find nodes on the shared interface we need to adapt */
419674ae819SStefano Zampini         for (j=0;j<graph->nvtxs;j++) {
420674ae819SStefano Zampini           if (graph->subset[j]==i+1) {
421674ae819SStefano Zampini             nodes_to_temp_buffer_indices[j] = temp_buffer_size;
422674ae819SStefano Zampini             temp_buffer_size++;
423674ae819SStefano Zampini           } else {
424674ae819SStefano Zampini             nodes_to_temp_buffer_indices[j] = -1;
425674ae819SStefano Zampini           }
426674ae819SStefano Zampini         }
427674ae819SStefano Zampini         /* allocate some temporary space */
428785e854fSJed Brown         ierr = PetscMalloc1(temp_buffer_size,&temp_buffer);CHKERRQ(ierr);
429785e854fSJed Brown         ierr = PetscMalloc1(temp_buffer_size*(cum_recv_counts[i+1]-cum_recv_counts[i]),&temp_buffer[0]);CHKERRQ(ierr);
430674ae819SStefano Zampini         ierr = PetscMemzero(temp_buffer[0],temp_buffer_size*(cum_recv_counts[i+1]-cum_recv_counts[i])*sizeof(PetscInt));CHKERRQ(ierr);
431674ae819SStefano Zampini         for (j=1;j<temp_buffer_size;j++) {
432674ae819SStefano Zampini           temp_buffer[j] = temp_buffer[j-1]+cum_recv_counts[i+1]-cum_recv_counts[i];
433674ae819SStefano Zampini         }
434674ae819SStefano Zampini         /* analyze contributions from neighbouring subdomains for i-th conn comp
435674ae819SStefano Zampini            temp buffer structure:
436674ae819SStefano Zampini            supposing part of the interface has dimension 5 (for example with global dofs 0,1,2,3,4)
437674ae819SStefano Zampini            3 neighs procs with structured connected components:
438674ae819SStefano Zampini              neigh 0: [0 1 4], [2 3];  (2 connected components)
439674ae819SStefano Zampini              neigh 1: [0 1], [2 3 4];  (2 connected components)
440674ae819SStefano Zampini              neigh 2: [0 4], [1], [2 3]; (3 connected components)
441674ae819SStefano Zampini            tempbuffer (row-oriented) will be filled as:
442674ae819SStefano Zampini              [ 0, 0, 0;
443674ae819SStefano Zampini                0, 0, 1;
444674ae819SStefano Zampini                1, 1, 2;
445674ae819SStefano Zampini                1, 1, 2;
446674ae819SStefano Zampini                0, 1, 0; ];
447674ae819SStefano Zampini            This way we can simply find intersections of ccs among neighs.
448674ae819SStefano Zampini            For the example above, the graph->subset array will be modified to reproduce the following 4 connected components [0], [1], [2 3], [4];
449674ae819SStefano Zampini                                                                                                                                    */
450674ae819SStefano Zampini         for (j=0;j<cum_recv_counts[i+1]-cum_recv_counts[i];j++) {
451674ae819SStefano Zampini           ins_val = 0;
452674ae819SStefano Zampini           size_of_recv = recv_buffer_subset[sum_requests];  /* total size of recv from neighs */
453674ae819SStefano Zampini           for (buffer_size=0;buffer_size<size_of_recv;) {  /* loop until all data from neighs has been taken into account */
454674ae819SStefano Zampini             for (k=1;k<recv_buffer[buffer_size+start_of_recv]+1;k++) { /* filling properly temp_buffer using data from a single recv */
455674ae819SStefano Zampini               temp_buffer[nodes_to_temp_buffer_indices[recv_buffer[start_of_recv+buffer_size+k]]][j] = ins_val;
456674ae819SStefano Zampini             }
457674ae819SStefano Zampini             buffer_size += k;
458674ae819SStefano Zampini             ins_val++;
459674ae819SStefano Zampini           }
460674ae819SStefano Zampini           start_of_recv += size_of_recv;
461674ae819SStefano Zampini           sum_requests++;
462674ae819SStefano Zampini         }
463785e854fSJed Brown         ierr = PetscMalloc1(temp_buffer_size,&add_to_subset);CHKERRQ(ierr);
464674ae819SStefano Zampini         ierr = PetscMemzero(add_to_subset,temp_buffer_size*sizeof(*add_to_subset));CHKERRQ(ierr);
465674ae819SStefano Zampini         for (j=0;j<temp_buffer_size;j++) {
466674ae819SStefano Zampini           if (!add_to_subset[j]) { /* found a new cc  */
467674ae819SStefano Zampini             global_subset_counter++;
468674ae819SStefano Zampini             add_to_subset[j] = global_subset_counter;
469674ae819SStefano Zampini             for (k=j+1;k<temp_buffer_size;k++) { /* check for other nodes in new cc */
470674ae819SStefano Zampini               same_set = PETSC_TRUE;
471674ae819SStefano Zampini               for (s=0;s<cum_recv_counts[i+1]-cum_recv_counts[i];s++) {
472674ae819SStefano Zampini                 if (temp_buffer[j][s]!=temp_buffer[k][s]) {
473674ae819SStefano Zampini                   same_set = PETSC_FALSE;
474674ae819SStefano Zampini                   break;
475674ae819SStefano Zampini                 }
476674ae819SStefano Zampini               }
477674ae819SStefano Zampini               if (same_set) {
478674ae819SStefano Zampini                 add_to_subset[k] = global_subset_counter;
479674ae819SStefano Zampini               }
480674ae819SStefano Zampini             }
481674ae819SStefano Zampini           }
482674ae819SStefano Zampini         }
483674ae819SStefano Zampini         /* insert new data in subset array */
484674ae819SStefano Zampini         temp_buffer_size = 0;
485674ae819SStefano Zampini         for (j=0;j<graph->nvtxs;j++) {
486674ae819SStefano Zampini           if (graph->subset[j]==i+1) {
487674ae819SStefano Zampini             graph->subset[j] = graph->n_subsets+add_to_subset[temp_buffer_size];
488674ae819SStefano Zampini             temp_buffer_size++;
489674ae819SStefano Zampini           }
490674ae819SStefano Zampini         }
491674ae819SStefano Zampini         ierr = PetscFree(temp_buffer[0]);CHKERRQ(ierr);
492674ae819SStefano Zampini         ierr = PetscFree(temp_buffer);CHKERRQ(ierr);
493674ae819SStefano Zampini         ierr = PetscFree(add_to_subset);CHKERRQ(ierr);
494674ae819SStefano Zampini       }
495674ae819SStefano Zampini     }
496674ae819SStefano Zampini     ierr = PetscFree(nodes_to_temp_buffer_indices);CHKERRQ(ierr);
497674ae819SStefano Zampini     ierr = PetscFree(sizes_of_sends);CHKERRQ(ierr);
498674ae819SStefano Zampini     ierr = PetscFree(send_requests);CHKERRQ(ierr);
499674ae819SStefano Zampini     ierr = PetscFree(recv_requests);CHKERRQ(ierr);
500674ae819SStefano Zampini     ierr = PetscFree(recv_buffer);CHKERRQ(ierr);
501674ae819SStefano Zampini     ierr = PetscFree(recv_buffer_subset);CHKERRQ(ierr);
502674ae819SStefano Zampini     ierr = PetscFree(send_buffer);CHKERRQ(ierr);
503674ae819SStefano Zampini     ierr = PetscFree(cum_recv_counts);CHKERRQ(ierr);
504674ae819SStefano Zampini     ierr = PetscFree(subset_to_nodes_indices);CHKERRQ(ierr);
505674ae819SStefano Zampini     ierr = PetscFree(subset_cc_adapt);CHKERRQ(ierr);
506674ae819SStefano Zampini     /* We are ready to find for connected components consistent among neighbouring subdomains */
507674ae819SStefano Zampini     if (global_subset_counter) {
508df48933bSStefano Zampini       ierr = PetscBTMemzero(graph->nvtxs,graph->touched);CHKERRQ(ierr);
509674ae819SStefano Zampini       global_subset_counter = 0;
510674ae819SStefano Zampini       for (i=0;i<graph->nvtxs;i++) {
511df48933bSStefano Zampini         if (graph->subset[i] && !PetscBTLookup(graph->touched,i)) {
512674ae819SStefano Zampini           global_subset_counter++;
513674ae819SStefano Zampini           for (j=i+1;j<graph->nvtxs;j++) {
514df48933bSStefano Zampini             if (!PetscBTLookup(graph->touched,j) && graph->subset[j]==graph->subset[i]) {
515674ae819SStefano Zampini               graph->subset[j] = global_subset_counter;
516df48933bSStefano Zampini               ierr = PetscBTSet(graph->touched,j);CHKERRQ(ierr);
517674ae819SStefano Zampini             }
518674ae819SStefano Zampini           }
519674ae819SStefano Zampini           graph->subset[i] = global_subset_counter;
520df48933bSStefano Zampini           ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
521674ae819SStefano Zampini         }
522674ae819SStefano Zampini       }
523674ae819SStefano Zampini       /* refine connected components locally */
524674ae819SStefano Zampini       ierr = PCBDDCGraphComputeConnectedComponentsLocal(graph);CHKERRQ(ierr);
525674ae819SStefano Zampini     }
5265b1b9aeaSStefano Zampini     /* restore original CSR graph of dofs */
5275b1b9aeaSStefano Zampini     ierr = PetscFree(new_xadj);CHKERRQ(ierr);
5285b1b9aeaSStefano Zampini     ierr = PetscFree(new_adjncy);CHKERRQ(ierr);
5295b1b9aeaSStefano Zampini     graph->xadj = old_xadj;
5305b1b9aeaSStefano Zampini     graph->adjncy = old_adjncy;
531674ae819SStefano Zampini     ierr = PetscFree(queue_in_global_numbering);CHKERRQ(ierr);
532674ae819SStefano Zampini   }
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;
544674ae819SStefano Zampini   PetscInt       *queue_global;
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;
556674ae819SStefano Zampini   for (i=0;i<graph->nvtxs;i++) {
5570cf82fd6SStefano Zampini     if (graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK || !graph->count[i]) {
558df48933bSStefano Zampini       ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
559674ae819SStefano Zampini       graph->subset[i] = 0;
560674ae819SStefano Zampini     }
561674ae819SStefano Zampini     graph->n_subsets = PetscMax(graph->n_subsets,graph->subset[i]);
562674ae819SStefano Zampini   }
563674ae819SStefano Zampini   ierr = PetscFree(graph->subset_ncc);CHKERRQ(ierr);
564785e854fSJed Brown   ierr = PetscMalloc1(graph->n_subsets,&graph->subset_ncc);CHKERRQ(ierr);
565674ae819SStefano Zampini   ierr = PetscMemzero(graph->subset_ncc,graph->n_subsets*sizeof(*graph->subset_ncc));CHKERRQ(ierr);
566674ae819SStefano Zampini   ierr = PetscMemzero(graph->cptr,(graph->nvtxs+1)*sizeof(*graph->cptr));CHKERRQ(ierr);
567674ae819SStefano Zampini   ierr = PetscMemzero(graph->queue,graph->nvtxs*sizeof(*graph->queue));CHKERRQ(ierr);
568674ae819SStefano Zampini 
569674ae819SStefano Zampini   /* begin search for connected components */
570674ae819SStefano Zampini   cum_queue = 0;
571674ae819SStefano Zampini   ncc = 0;
572674ae819SStefano Zampini   for (n=0;n<graph->n_subsets;n++) {
573674ae819SStefano Zampini     pid = n+1;  /* partition labeled by 0 is discarded */
574674ae819SStefano Zampini     nleft = 0;
575674ae819SStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
576674ae819SStefano Zampini       if (graph->subset[i] == pid) {
577674ae819SStefano Zampini         nleft++;
578674ae819SStefano Zampini       }
579674ae819SStefano Zampini     }
580674ae819SStefano Zampini     for (i=0; i<graph->nvtxs; i++) {
581674ae819SStefano Zampini       if (graph->subset[i] == pid) {
582674ae819SStefano Zampini         break;
583674ae819SStefano Zampini       }
584674ae819SStefano Zampini     }
585df48933bSStefano Zampini     ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
586674ae819SStefano Zampini     graph->queue[cum_queue] = i;
587674ae819SStefano Zampini     first = 0;
588674ae819SStefano Zampini     last = 1;
589674ae819SStefano Zampini     graph->cptr[ncc] = cum_queue;
590674ae819SStefano Zampini     ncc_pid = 0;
591674ae819SStefano Zampini     while (first != nleft) {
592674ae819SStefano Zampini       if (first == last) {
593674ae819SStefano Zampini         graph->cptr[++ncc] = first+cum_queue;
594674ae819SStefano Zampini         ncc_pid++;
595df48933bSStefano Zampini         for (i=0; i<graph->nvtxs; i++) { /* TODO-> use a while! */
596df48933bSStefano Zampini           if (graph->subset[i] == pid && !PetscBTLookup(graph->touched,i)) {
597674ae819SStefano Zampini             break;
598674ae819SStefano Zampini           }
599674ae819SStefano Zampini         }
600674ae819SStefano Zampini         graph->queue[cum_queue+last] = i;
601674ae819SStefano Zampini         last++;
602df48933bSStefano Zampini         ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
603674ae819SStefano Zampini       }
604674ae819SStefano Zampini       i = graph->queue[cum_queue+first];
605674ae819SStefano Zampini       first++;
606674ae819SStefano Zampini       for (j=graph->xadj[i];j<graph->xadj[i+1];j++) {
607674ae819SStefano Zampini         k = graph->adjncy[j];
608df48933bSStefano Zampini         if (graph->subset[k] == pid && !PetscBTLookup(graph->touched,k)) {
609674ae819SStefano Zampini           graph->queue[cum_queue+last] = k;
610674ae819SStefano Zampini           last++;
611df48933bSStefano Zampini           ierr = PetscBTSet(graph->touched,k);CHKERRQ(ierr);
612674ae819SStefano Zampini         }
613674ae819SStefano Zampini       }
614674ae819SStefano Zampini     }
615674ae819SStefano Zampini     graph->cptr[++ncc] = first+cum_queue;
616674ae819SStefano Zampini     ncc_pid++;
617674ae819SStefano Zampini     cum_queue = graph->cptr[ncc];
618674ae819SStefano Zampini     graph->subset_ncc[n] = ncc_pid;
619674ae819SStefano Zampini   }
620674ae819SStefano Zampini   graph->ncc = ncc;
621674ae819SStefano Zampini   /* For consistency among neighbours, I need to sort (by global ordering) each connected component */
622785e854fSJed Brown   ierr = PetscMalloc1(graph->cptr[graph->ncc],&queue_global);CHKERRQ(ierr);
623674ae819SStefano Zampini   ierr = ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_global);CHKERRQ(ierr);
624674ae819SStefano Zampini   for (i=0;i<graph->ncc;i++) {
625674ae819SStefano Zampini     ierr = PetscSortIntWithArray(graph->cptr[i+1]-graph->cptr[i],&queue_global[graph->cptr[i]],&graph->queue[graph->cptr[i]]);CHKERRQ(ierr);
626674ae819SStefano Zampini   }
627674ae819SStefano Zampini   ierr = PetscFree(queue_global);CHKERRQ(ierr);
628674ae819SStefano Zampini   PetscFunctionReturn(0);
629674ae819SStefano Zampini }
630674ae819SStefano Zampini 
631674ae819SStefano Zampini #undef __FUNCT__
632674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphSetUp"
633674ae819SStefano Zampini PetscErrorCode PCBDDCGraphSetUp(PCBDDCGraph graph, PetscInt custom_minimal_size, IS neumann_is, IS dirichlet_is, PetscInt n_ISForDofs, IS ISForDofs[], IS custom_primal_vertices)
634674ae819SStefano Zampini {
635674ae819SStefano Zampini   VecScatter     scatter_ctx;
636674ae819SStefano Zampini   Vec            local_vec,local_vec2,global_vec;
637674ae819SStefano Zampini   IS             to,from;
638674ae819SStefano Zampini   MPI_Comm       comm;
639674ae819SStefano Zampini   PetscScalar    *array,*array2;
640674ae819SStefano Zampini   const PetscInt *is_indices;
6413a5b03d1SStefano Zampini   PetscInt       n_neigh,*neigh,*n_shared,**shared,*queue_global;
642674ae819SStefano Zampini   PetscInt       i,j,k,s,total_counts,nodes_touched,is_size;
643674ae819SStefano Zampini   PetscErrorCode ierr;
64451b0f6cfSStefano Zampini   PetscBool      same_set,mirrors_found;
645674ae819SStefano Zampini 
646674ae819SStefano Zampini   PetscFunctionBegin;
647674ae819SStefano Zampini   ierr = PetscObjectGetComm((PetscObject)(graph->l2gmap),&comm);CHKERRQ(ierr);
648674ae819SStefano Zampini   /* custom_minimal_size */
649674ae819SStefano Zampini   graph->custom_minimal_size = PetscMax(graph->custom_minimal_size,custom_minimal_size);
650674ae819SStefano Zampini   /* get info l2gmap and allocate work vectors  */
651674ae819SStefano Zampini   ierr = ISLocalToGlobalMappingGetInfo(graph->l2gmap,&n_neigh,&neigh,&n_shared,&shared);CHKERRQ(ierr);
652674ae819SStefano Zampini   ierr = ISLocalToGlobalMappingGetIndices(graph->l2gmap,&is_indices);CHKERRQ(ierr);
653674ae819SStefano Zampini   j = 0;
654674ae819SStefano Zampini   for (i=0;i<graph->nvtxs;i++) {
655674ae819SStefano Zampini     j = PetscMax(j,is_indices[i]);
656674ae819SStefano Zampini   }
657674ae819SStefano Zampini   ierr = MPI_Allreduce(&j,&i,1,MPIU_INT,MPI_MAX,comm);CHKERRQ(ierr);
658674ae819SStefano Zampini   i++;
659674ae819SStefano Zampini   ierr = VecCreate(PETSC_COMM_SELF,&local_vec);CHKERRQ(ierr);
660674ae819SStefano Zampini   ierr = VecSetSizes(local_vec,PETSC_DECIDE,graph->nvtxs);CHKERRQ(ierr);
661674ae819SStefano Zampini   ierr = VecSetType(local_vec,VECSTANDARD);CHKERRQ(ierr);
662674ae819SStefano Zampini   ierr = VecDuplicate(local_vec,&local_vec2);CHKERRQ(ierr);
663674ae819SStefano Zampini   ierr = VecCreate(comm,&global_vec);CHKERRQ(ierr);
664674ae819SStefano Zampini   ierr = VecSetSizes(global_vec,PETSC_DECIDE,i);CHKERRQ(ierr);
665674ae819SStefano Zampini   ierr = VecSetType(global_vec,VECSTANDARD);CHKERRQ(ierr);
666674ae819SStefano Zampini   ierr = ISCreateStride(PETSC_COMM_SELF,graph->nvtxs,0,1,&to);CHKERRQ(ierr);
667674ae819SStefano Zampini   ierr = ISLocalToGlobalMappingApplyIS(graph->l2gmap,to,&from);CHKERRQ(ierr);
668674ae819SStefano Zampini   ierr = VecScatterCreate(global_vec,from,local_vec,to,&scatter_ctx);CHKERRQ(ierr);
66951b0f6cfSStefano Zampini 
67051b0f6cfSStefano Zampini   /* get local periodic nodes */
67151b0f6cfSStefano Zampini   mirrors_found = PETSC_FALSE;
6729881197aSStefano Zampini   if (graph->nvtxs && n_neigh) {
67349eeff8cSStefano Zampini     for (i=0; i<n_shared[0]; i++) graph->count[shared[0][i]] += 1;
67451b0f6cfSStefano Zampini     for (i=0; i<n_shared[0]; i++) {
67551b0f6cfSStefano Zampini       if (graph->count[shared[0][i]] > 1) {
67651b0f6cfSStefano Zampini         mirrors_found = PETSC_TRUE;
67751b0f6cfSStefano Zampini         break;
67851b0f6cfSStefano Zampini       }
67951b0f6cfSStefano Zampini     }
68049eeff8cSStefano Zampini   }
68151b0f6cfSStefano Zampini   /* compute local mirrors (if any) */
68251b0f6cfSStefano Zampini   if (mirrors_found) {
68351b0f6cfSStefano Zampini     PetscInt *local_indices,*global_indices;
68451b0f6cfSStefano Zampini     /* get arrays of local and global indices */
685785e854fSJed Brown     ierr = PetscMalloc1(graph->nvtxs,&local_indices);CHKERRQ(ierr);
68651b0f6cfSStefano Zampini     ierr = ISGetIndices(to,(const PetscInt**)&is_indices);CHKERRQ(ierr);
68751b0f6cfSStefano Zampini     ierr = PetscMemcpy(local_indices,is_indices,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
68851b0f6cfSStefano Zampini     ierr = ISRestoreIndices(to,(const PetscInt**)&is_indices);CHKERRQ(ierr);
689785e854fSJed Brown     ierr = PetscMalloc1(graph->nvtxs,&global_indices);CHKERRQ(ierr);
69051b0f6cfSStefano Zampini     ierr = ISGetIndices(from,(const PetscInt**)&is_indices);CHKERRQ(ierr);
69151b0f6cfSStefano Zampini     ierr = PetscMemcpy(global_indices,is_indices,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
69251b0f6cfSStefano Zampini     ierr = ISRestoreIndices(from,(const PetscInt**)&is_indices);CHKERRQ(ierr);
69351b0f6cfSStefano Zampini     /* allocate space for mirrors */
6948e5aaad5SJed Brown     ierr = PetscMalloc2(graph->nvtxs,&graph->mirrors,graph->nvtxs,&graph->mirrors_set);CHKERRQ(ierr);
69551b0f6cfSStefano Zampini     ierr = PetscMemzero(graph->mirrors,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
69651b0f6cfSStefano Zampini     graph->mirrors_set[0] = 0;
69751b0f6cfSStefano Zampini 
69851b0f6cfSStefano Zampini     k=0;
69951b0f6cfSStefano Zampini     for (i=0;i<n_shared[0];i++) {
70051b0f6cfSStefano Zampini       j=shared[0][i];
70151b0f6cfSStefano Zampini       if (graph->count[j] > 1) {
70251b0f6cfSStefano Zampini         graph->mirrors[j]++;
70351b0f6cfSStefano Zampini         k++;
70451b0f6cfSStefano Zampini       }
70551b0f6cfSStefano Zampini     }
70651b0f6cfSStefano Zampini     /* allocate space for set of mirrors */
707785e854fSJed Brown     ierr = PetscMalloc1(k,&graph->mirrors_set[0]);CHKERRQ(ierr);
70851b0f6cfSStefano Zampini     for (i=1;i<graph->nvtxs;i++)
70951b0f6cfSStefano Zampini       graph->mirrors_set[i]=graph->mirrors_set[i-1]+graph->mirrors[i-1];
71051b0f6cfSStefano Zampini 
71151b0f6cfSStefano Zampini     /* fill arrays */
71251b0f6cfSStefano Zampini     ierr = PetscMemzero(graph->mirrors,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
71351b0f6cfSStefano Zampini     for (j=0;j<n_shared[0];j++) {
71451b0f6cfSStefano Zampini       i=shared[0][j];
71551b0f6cfSStefano Zampini       if (graph->count[i] > 1)
71651b0f6cfSStefano Zampini         graph->mirrors_set[i][graph->mirrors[i]++]=global_indices[i];
71751b0f6cfSStefano Zampini     }
71851b0f6cfSStefano Zampini     ierr = PetscSortIntWithArray(graph->nvtxs,global_indices,local_indices);CHKERRQ(ierr);
71951b0f6cfSStefano Zampini     for (i=0;i<graph->nvtxs;i++) {
72051b0f6cfSStefano Zampini       if (graph->mirrors[i] > 0) {
72151b0f6cfSStefano Zampini         ierr = PetscFindInt(graph->mirrors_set[i][0],graph->nvtxs,global_indices,&k);CHKERRQ(ierr);
72251b0f6cfSStefano Zampini         j = global_indices[k];
72351b0f6cfSStefano Zampini         while ( k > 0 && global_indices[k-1] == j) k--;
72451b0f6cfSStefano Zampini         for (j=0;j<graph->mirrors[i];j++) {
72551b0f6cfSStefano Zampini           graph->mirrors_set[i][j]=local_indices[k+j];
72651b0f6cfSStefano Zampini         }
72751b0f6cfSStefano Zampini         ierr = PetscSortInt(graph->mirrors[i],graph->mirrors_set[i]);CHKERRQ(ierr);
72851b0f6cfSStefano Zampini       }
72951b0f6cfSStefano Zampini     }
73051b0f6cfSStefano Zampini     ierr = PetscFree(local_indices);CHKERRQ(ierr);
73151b0f6cfSStefano Zampini     ierr = PetscFree(global_indices);CHKERRQ(ierr);
73251b0f6cfSStefano Zampini   }
73351b0f6cfSStefano Zampini   ierr = PetscMemzero(graph->count,graph->nvtxs*sizeof(*graph->count));CHKERRQ(ierr);
734674ae819SStefano Zampini   ierr = ISDestroy(&to);CHKERRQ(ierr);
735674ae819SStefano Zampini   ierr = ISDestroy(&from);CHKERRQ(ierr);
73651b0f6cfSStefano Zampini 
737674ae819SStefano Zampini   /* Count total number of neigh per node */
738674ae819SStefano Zampini   k=0;
739674ae819SStefano Zampini   for (i=1;i<n_neigh;i++) {
740674ae819SStefano Zampini     k += n_shared[i];
741674ae819SStefano Zampini     for (j=0;j<n_shared[i];j++) {
742674ae819SStefano Zampini       graph->count[shared[i][j]] += 1;
743674ae819SStefano Zampini     }
744674ae819SStefano Zampini   }
745674ae819SStefano Zampini   /* Allocate space for storing the set of neighbours for each node */
746674ae819SStefano Zampini   if (graph->nvtxs) {
747785e854fSJed Brown     ierr = PetscMalloc1(k,&graph->neighbours_set[0]);CHKERRQ(ierr);
748674ae819SStefano Zampini   }
749674ae819SStefano Zampini   for (i=1;i<graph->nvtxs;i++) { /* dont count myself */
750674ae819SStefano Zampini     graph->neighbours_set[i]=graph->neighbours_set[i-1]+graph->count[i-1];
751674ae819SStefano Zampini   }
752674ae819SStefano Zampini   /* Get information for sharing subdomains */
753674ae819SStefano Zampini   ierr = PetscMemzero(graph->count,graph->nvtxs*sizeof(*graph->count));CHKERRQ(ierr);
754674ae819SStefano Zampini   for (i=1;i<n_neigh;i++) { /* dont count myself */
755674ae819SStefano Zampini     s = n_shared[i];
756674ae819SStefano Zampini     for (j=0;j<s;j++) {
757674ae819SStefano Zampini       k = shared[i][j];
758674ae819SStefano Zampini       graph->neighbours_set[k][graph->count[k]] = neigh[i];
759674ae819SStefano Zampini       graph->count[k] += 1;
760674ae819SStefano Zampini     }
761674ae819SStefano Zampini   }
762674ae819SStefano Zampini   /* sort set of sharing subdomains */
763674ae819SStefano Zampini   for (i=0;i<graph->nvtxs;i++) {
76493fb5fd3SStefano Zampini     ierr = PetscSortRemoveDupsInt(&graph->count[i],graph->neighbours_set[i]);CHKERRQ(ierr);
765674ae819SStefano Zampini   }
766*67c9da69SStefano Zampini   /*
767*67c9da69SStefano Zampini      Get info for dofs splitting
768*67c9da69SStefano Zampini      User can specify only a subset; an additional field is considered as a complementary set
769*67c9da69SStefano Zampini   */
770*67c9da69SStefano Zampini   for (i=0;i<graph->nvtxs;i++) {
771*67c9da69SStefano Zampini     graph->which_dof[i] = n_ISForDofs; /* by default a dof belongs to the complement set */
772*67c9da69SStefano Zampini   }
773674ae819SStefano Zampini   for (i=0;i<n_ISForDofs;i++) {
77463602bcaSStefano Zampini     ierr = ISGetLocalSize(ISForDofs[i],&is_size);CHKERRQ(ierr);
775674ae819SStefano Zampini     ierr = ISGetIndices(ISForDofs[i],(const PetscInt**)&is_indices);CHKERRQ(ierr);
776674ae819SStefano Zampini     for (j=0;j<is_size;j++) {
777*67c9da69SStefano Zampini       if (is_indices[j] > 0 && is_indices[j] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */
778674ae819SStefano Zampini         graph->which_dof[is_indices[j]] = i;
779674ae819SStefano Zampini       }
780*67c9da69SStefano Zampini     }
781674ae819SStefano Zampini     ierr = ISRestoreIndices(ISForDofs[i],(const PetscInt**)&is_indices);CHKERRQ(ierr);
782674ae819SStefano Zampini   }
783674ae819SStefano Zampini   /* Take into account Neumann nodes */
784674ae819SStefano Zampini   ierr = VecSet(local_vec,0.0);CHKERRQ(ierr);
785674ae819SStefano Zampini   ierr = VecSet(local_vec2,0.0);CHKERRQ(ierr);
786674ae819SStefano Zampini   if (neumann_is) {
787674ae819SStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
78882ba6b80SStefano Zampini     ierr = ISGetLocalSize(neumann_is,&is_size);CHKERRQ(ierr);
789674ae819SStefano Zampini     ierr = ISGetIndices(neumann_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
790674ae819SStefano Zampini     for (i=0;i<is_size;i++) {
791674ae819SStefano Zampini       array[is_indices[i]] = 1.0;
792674ae819SStefano Zampini     }
793674ae819SStefano Zampini     ierr = ISRestoreIndices(neumann_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
794674ae819SStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
795674ae819SStefano Zampini   }
796674ae819SStefano Zampini   /* Neumann nodes: impose consistency among neighbours */
797674ae819SStefano Zampini   ierr = VecSet(global_vec,0.0);CHKERRQ(ierr);
798674ae819SStefano Zampini   ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
799674ae819SStefano Zampini   ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
800674ae819SStefano Zampini   ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
801674ae819SStefano Zampini   ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
802674ae819SStefano Zampini   ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
803674ae819SStefano Zampini   for (i=0;i<graph->nvtxs;i++) {
8045b08dc53SStefano Zampini     if (PetscRealPart(array[i]) > 0.0) {
8050cf82fd6SStefano Zampini       graph->special_dof[i] = PCBDDCGRAPH_NEUMANN_MARK;
806674ae819SStefano Zampini     }
807674ae819SStefano Zampini   }
808674ae819SStefano Zampini   ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
809674ae819SStefano Zampini   /* Take into account Dirichlet nodes */
810674ae819SStefano Zampini   ierr = VecSet(local_vec2,0.0);CHKERRQ(ierr);
811674ae819SStefano Zampini   if (dirichlet_is) {
812674ae819SStefano Zampini     ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
813674ae819SStefano Zampini     ierr = VecGetArray(local_vec2,&array2);CHKERRQ(ierr);
81482ba6b80SStefano Zampini     ierr = ISGetLocalSize(dirichlet_is,&is_size);CHKERRQ(ierr);
815674ae819SStefano Zampini     ierr = ISGetIndices(dirichlet_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
816674ae819SStefano Zampini     for (i=0;i<is_size;i++){
817674ae819SStefano Zampini       k = is_indices[i];
818df48933bSStefano Zampini       if (graph->count[k] && !PetscBTLookup(graph->touched,k)) {
8195b08dc53SStefano Zampini         if (PetscRealPart(array[k]) > 0.0) {
820674ae819SStefano Zampini           SETERRQ1(PETSC_COMM_SELF,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);
821674ae819SStefano Zampini         }
822674ae819SStefano Zampini         array2[k] = 1.0;
823674ae819SStefano Zampini       }
824674ae819SStefano Zampini     }
825674ae819SStefano Zampini     ierr = ISRestoreIndices(dirichlet_is,(const PetscInt**)&is_indices);CHKERRQ(ierr);
826674ae819SStefano Zampini     ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
827674ae819SStefano Zampini     ierr = VecRestoreArray(local_vec2,&array2);CHKERRQ(ierr);
828674ae819SStefano Zampini   }
829674ae819SStefano Zampini   /* Dirichlet nodes: impose consistency among neighbours */
830674ae819SStefano Zampini   ierr = VecSet(global_vec,0.0);CHKERRQ(ierr);
831674ae819SStefano Zampini   ierr = VecScatterBegin(scatter_ctx,local_vec2,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
832674ae819SStefano Zampini   ierr = VecScatterEnd(scatter_ctx,local_vec2,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
833674ae819SStefano Zampini   ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
834674ae819SStefano Zampini   ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
835674ae819SStefano Zampini   ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr);
836674ae819SStefano Zampini   for (i=0;i<graph->nvtxs;i++) {
8375b08dc53SStefano Zampini     if (PetscRealPart(array[i]) > 0.0) {
838df48933bSStefano Zampini       ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
839674ae819SStefano Zampini       graph->subset[i] = 0; /* dirichlet nodes treated as internal -> is it ok? */
8400cf82fd6SStefano Zampini       graph->special_dof[i] = PCBDDCGRAPH_DIRICHLET_MARK;
841674ae819SStefano Zampini     }
842674ae819SStefano Zampini   }
843674ae819SStefano Zampini   ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr);
84451b0f6cfSStefano Zampini 
84508a5cf49SStefano Zampini   /* mark local periodic nodes (if any) and adapt CSR graph (if any) */
84651b0f6cfSStefano Zampini   if (graph->mirrors) {
84751b0f6cfSStefano Zampini     for (i=0;i<graph->nvtxs;i++)
84851b0f6cfSStefano Zampini       if (graph->mirrors[i])
8490cf82fd6SStefano Zampini         graph->special_dof[i] = PCBDDCGRAPH_LOCAL_PERIODIC_MARK;
85051b0f6cfSStefano Zampini 
85108a5cf49SStefano Zampini     if (graph->xadj && graph->adjncy) {
85208a5cf49SStefano Zampini       PetscInt *new_xadj,*new_adjncy;
85351b0f6cfSStefano Zampini       /* sort CSR graph */
85451b0f6cfSStefano Zampini       for (i=0;i<graph->nvtxs;i++)
85551b0f6cfSStefano Zampini         ierr = PetscSortInt(graph->xadj[i+1]-graph->xadj[i],&graph->adjncy[graph->xadj[i]]);CHKERRQ(ierr);
85651b0f6cfSStefano Zampini 
85751b0f6cfSStefano Zampini       /* adapt local CSR graph in case of local periodicity */
85851b0f6cfSStefano Zampini       k=0;
85951b0f6cfSStefano Zampini       for (i=0;i<graph->nvtxs;i++)
86051b0f6cfSStefano Zampini         for (j=graph->xadj[i];j<graph->xadj[i+1];j++)
86151b0f6cfSStefano Zampini           k += graph->mirrors[graph->adjncy[j]];
86251b0f6cfSStefano Zampini 
863785e854fSJed Brown       ierr = PetscMalloc1((graph->nvtxs+1),&new_xadj);CHKERRQ(ierr);
864785e854fSJed Brown       ierr = PetscMalloc1((k+graph->xadj[graph->nvtxs]),&new_adjncy);CHKERRQ(ierr);
86551b0f6cfSStefano Zampini       new_xadj[0]=0;
86651b0f6cfSStefano Zampini       for (i=0;i<graph->nvtxs;i++) {
86751b0f6cfSStefano Zampini         k = graph->xadj[i+1]-graph->xadj[i];
86851b0f6cfSStefano Zampini         ierr = PetscMemcpy(&new_adjncy[new_xadj[i]],&graph->adjncy[graph->xadj[i]],k*sizeof(PetscInt));CHKERRQ(ierr);
86951b0f6cfSStefano Zampini         new_xadj[i+1]=new_xadj[i]+k;
87051b0f6cfSStefano Zampini         for (j=graph->xadj[i];j<graph->xadj[i+1];j++) {
87151b0f6cfSStefano Zampini           k = graph->mirrors[graph->adjncy[j]];
87251b0f6cfSStefano Zampini           ierr = PetscMemcpy(&new_adjncy[new_xadj[i+1]],graph->mirrors_set[graph->adjncy[j]],k*sizeof(PetscInt));CHKERRQ(ierr);
87351b0f6cfSStefano Zampini           new_xadj[i+1]+=k;
87451b0f6cfSStefano Zampini         }
87551b0f6cfSStefano Zampini         k = new_xadj[i+1]-new_xadj[i];
87651b0f6cfSStefano Zampini         ierr = PetscSortRemoveDupsInt(&k,&new_adjncy[new_xadj[i]]);CHKERRQ(ierr);
87751b0f6cfSStefano Zampini         new_xadj[i+1]=new_xadj[i]+k;
87851b0f6cfSStefano Zampini       }
87951b0f6cfSStefano Zampini       /* set new CSR into graph */
88051b0f6cfSStefano Zampini       ierr = PetscFree(graph->xadj);CHKERRQ(ierr);
88151b0f6cfSStefano Zampini       ierr = PetscFree(graph->adjncy);CHKERRQ(ierr);
88251b0f6cfSStefano Zampini       graph->xadj = new_xadj;
88351b0f6cfSStefano Zampini       graph->adjncy = new_adjncy;
88451b0f6cfSStefano Zampini     }
88508a5cf49SStefano Zampini   }
88651b0f6cfSStefano Zampini 
887674ae819SStefano Zampini   /* mark special nodes -> each will become a single node equivalence class */
888674ae819SStefano Zampini   if (custom_primal_vertices) {
88963602bcaSStefano Zampini     ierr = ISGetLocalSize(custom_primal_vertices,&is_size);CHKERRQ(ierr);
890674ae819SStefano Zampini     ierr = ISGetIndices(custom_primal_vertices,(const PetscInt**)&is_indices);CHKERRQ(ierr);
891674ae819SStefano Zampini     for (i=0;i<is_size;i++) {
8920cf82fd6SStefano Zampini       graph->special_dof[is_indices[i]] = PCBDDCGRAPH_SPECIAL_MARK-i;
893674ae819SStefano Zampini     }
894674ae819SStefano Zampini     ierr = ISRestoreIndices(custom_primal_vertices,(const PetscInt**)&is_indices);CHKERRQ(ierr);
895674ae819SStefano Zampini   }
896674ae819SStefano Zampini   /* mark interior nodes as touched and belonging to partition number 0 */
897674ae819SStefano Zampini   for (i=0;i<graph->nvtxs;i++) {
898674ae819SStefano Zampini     if (!graph->count[i]) {
899df48933bSStefano Zampini       ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
900674ae819SStefano Zampini       graph->subset[i] = 0;
901674ae819SStefano Zampini     }
902674ae819SStefano Zampini   }
903674ae819SStefano Zampini   /* init graph structure and compute default subsets */
904674ae819SStefano Zampini   nodes_touched=0;
905674ae819SStefano Zampini   for (i=0;i<graph->nvtxs;i++) {
906df48933bSStefano Zampini     if (PetscBTLookup(graph->touched,i)) {
907674ae819SStefano Zampini       nodes_touched++;
908674ae819SStefano Zampini     }
909674ae819SStefano Zampini   }
910674ae819SStefano Zampini   i = 0;
911674ae819SStefano Zampini   graph->ncc = 0;
912674ae819SStefano Zampini   total_counts = 0;
913674ae819SStefano Zampini   while (nodes_touched<graph->nvtxs) {
914674ae819SStefano Zampini     /*  find first untouched node in local ordering */
915df48933bSStefano Zampini     while (PetscBTLookup(graph->touched,i)) {
916674ae819SStefano Zampini       i++;
917674ae819SStefano Zampini     }
918df48933bSStefano Zampini     ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr);
919674ae819SStefano Zampini     graph->subset[i] = graph->ncc+1;
920674ae819SStefano Zampini     graph->cptr[graph->ncc] = total_counts;
921674ae819SStefano Zampini     graph->queue[total_counts] = i;
922674ae819SStefano Zampini     total_counts++;
923674ae819SStefano Zampini     nodes_touched++;
924674ae819SStefano Zampini     /* now find all other nodes having the same set of sharing subdomains */
925674ae819SStefano Zampini     for (j=i+1;j<graph->nvtxs;j++) {
92674e413f5SStefano Zampini       /* check for same number of sharing subdomains, dof number and same special mark */
927df48933bSStefano 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]) {
928674ae819SStefano Zampini         /* check for same set of sharing subdomains */
929674ae819SStefano Zampini         same_set=PETSC_TRUE;
930674ae819SStefano Zampini         for (k=0;k<graph->count[j];k++){
931674ae819SStefano Zampini           if (graph->neighbours_set[i][k]!=graph->neighbours_set[j][k]) {
932674ae819SStefano Zampini             same_set=PETSC_FALSE;
933674ae819SStefano Zampini           }
934674ae819SStefano Zampini         }
935674ae819SStefano Zampini         /* I found a friend of mine */
936674ae819SStefano Zampini         if (same_set) {
937674ae819SStefano Zampini           graph->subset[j]=graph->ncc+1;
938df48933bSStefano Zampini           ierr = PetscBTSet(graph->touched,j);CHKERRQ(ierr);
939674ae819SStefano Zampini           nodes_touched++;
940674ae819SStefano Zampini           graph->queue[total_counts] = j;
941674ae819SStefano Zampini           total_counts++;
942674ae819SStefano Zampini         }
943674ae819SStefano Zampini       }
944674ae819SStefano Zampini     }
945674ae819SStefano Zampini     graph->ncc++;
946674ae819SStefano Zampini   }
947674ae819SStefano Zampini   /* set default number of subsets (at this point no info on csr graph has been taken into account, so n_subsets = ncc */
948674ae819SStefano Zampini   graph->n_subsets = graph->ncc;
949785e854fSJed Brown   ierr = PetscMalloc1(graph->n_subsets,&graph->subset_ncc);CHKERRQ(ierr);
950674ae819SStefano Zampini   for (i=0;i<graph->n_subsets;i++) {
951674ae819SStefano Zampini     graph->subset_ncc[i] = 1;
952674ae819SStefano Zampini   }
953674ae819SStefano Zampini   /* final pointer */
954674ae819SStefano Zampini   graph->cptr[graph->ncc] = total_counts;
955674ae819SStefano Zampini   /* free memory allocated by ISLocalToGlobalMappingGetInfo */
956674ae819SStefano Zampini   ierr = ISLocalToGlobalMappingRestoreInfo(graph->l2gmap,&n_neigh,&neigh,&n_shared,&shared);CHKERRQ(ierr);
9573a5b03d1SStefano Zampini   /* get a reference node (min index in global ordering) for each subset */
9583a5b03d1SStefano Zampini   ierr = PetscMalloc1(graph->ncc,&graph->subset_ref_node);CHKERRQ(ierr);
9593a5b03d1SStefano Zampini   ierr = PetscMalloc1(graph->cptr[graph->ncc],&queue_global);CHKERRQ(ierr);
9603a5b03d1SStefano Zampini   ierr = ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_global);CHKERRQ(ierr);
9613a5b03d1SStefano Zampini   for (i=0;i<graph->ncc;i++) {
9623a5b03d1SStefano Zampini     PetscInt minval = queue_global[graph->cptr[i]];
9633a5b03d1SStefano Zampini     for (j=graph->cptr[i]+1;j<graph->cptr[i+1];j++) {
9643a5b03d1SStefano Zampini       minval = PetscMin(minval,queue_global[j]);
9653a5b03d1SStefano Zampini     }
9663a5b03d1SStefano Zampini     graph->subset_ref_node[i] = minval;
9673a5b03d1SStefano Zampini   }
9683a5b03d1SStefano Zampini   ierr = PetscFree(queue_global);CHKERRQ(ierr);
969674ae819SStefano Zampini   /* free objects */
970674ae819SStefano Zampini   ierr = VecDestroy(&local_vec);CHKERRQ(ierr);
971674ae819SStefano Zampini   ierr = VecDestroy(&local_vec2);CHKERRQ(ierr);
972674ae819SStefano Zampini   ierr = VecDestroy(&global_vec);CHKERRQ(ierr);
973674ae819SStefano Zampini   ierr = VecScatterDestroy(&scatter_ctx);CHKERRQ(ierr);
974674ae819SStefano Zampini   PetscFunctionReturn(0);
975674ae819SStefano Zampini }
976674ae819SStefano Zampini 
977674ae819SStefano Zampini #undef __FUNCT__
978674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphResetCSR"
979674ae819SStefano Zampini PetscErrorCode PCBDDCGraphResetCSR(PCBDDCGraph graph)
980674ae819SStefano Zampini {
981674ae819SStefano Zampini   PetscErrorCode ierr;
982674ae819SStefano Zampini 
983674ae819SStefano Zampini   PetscFunctionBegin;
984674ae819SStefano Zampini   ierr = PetscFree(graph->xadj);CHKERRQ(ierr);
985674ae819SStefano Zampini   ierr = PetscFree(graph->adjncy);CHKERRQ(ierr);
986575ad6abSStefano Zampini   graph->nvtxs_csr = 0;
987674ae819SStefano Zampini   PetscFunctionReturn(0);
988674ae819SStefano Zampini }
989674ae819SStefano Zampini 
990674ae819SStefano Zampini #undef __FUNCT__
991674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphReset"
992674ae819SStefano Zampini PetscErrorCode PCBDDCGraphReset(PCBDDCGraph graph)
993674ae819SStefano Zampini {
994674ae819SStefano Zampini   PetscErrorCode ierr;
995674ae819SStefano Zampini 
996674ae819SStefano Zampini   PetscFunctionBegin;
997674ae819SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&graph->l2gmap);CHKERRQ(ierr);
998674ae819SStefano Zampini   ierr = PetscFree(graph->subset_ncc);CHKERRQ(ierr);
9993a5b03d1SStefano Zampini   ierr = PetscFree(graph->subset_ref_node);CHKERRQ(ierr);
1000674ae819SStefano Zampini   if (graph->nvtxs) {
1001674ae819SStefano Zampini     ierr = PetscFree(graph->neighbours_set[0]);CHKERRQ(ierr);
1002674ae819SStefano Zampini   }
1003df48933bSStefano Zampini   ierr = PetscBTDestroy(&graph->touched);CHKERRQ(ierr);
1004df48933bSStefano Zampini   ierr = PetscFree7(graph->count,
1005674ae819SStefano Zampini                     graph->neighbours_set,
1006674ae819SStefano Zampini                     graph->subset,
1007674ae819SStefano Zampini                     graph->which_dof,
1008674ae819SStefano Zampini                     graph->cptr,
1009df48933bSStefano Zampini                     graph->queue,
1010df48933bSStefano Zampini                     graph->special_dof);CHKERRQ(ierr);
101151b0f6cfSStefano Zampini   if (graph->mirrors) {
101251b0f6cfSStefano Zampini     ierr = PetscFree(graph->mirrors_set[0]);CHKERRQ(ierr);
101351b0f6cfSStefano Zampini   }
101451b0f6cfSStefano Zampini   ierr = PetscFree2(graph->mirrors,graph->mirrors_set);CHKERRQ(ierr);
1015674ae819SStefano Zampini   graph->nvtxs = 0;
1016674ae819SStefano Zampini   graph->n_subsets = 0;
1017674ae819SStefano Zampini   graph->custom_minimal_size = 1;
1018674ae819SStefano Zampini   PetscFunctionReturn(0);
1019674ae819SStefano Zampini }
1020674ae819SStefano Zampini 
1021674ae819SStefano Zampini #undef __FUNCT__
1022674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphInit"
1023674ae819SStefano Zampini PetscErrorCode PCBDDCGraphInit(PCBDDCGraph graph, ISLocalToGlobalMapping l2gmap)
1024674ae819SStefano Zampini {
1025674ae819SStefano Zampini   PetscInt       n;
1026674ae819SStefano Zampini   PetscErrorCode ierr;
1027674ae819SStefano Zampini 
1028674ae819SStefano Zampini   PetscFunctionBegin;
1029674ae819SStefano Zampini   PetscValidPointer(graph,1);
1030674ae819SStefano Zampini   PetscValidHeaderSpecific(l2gmap,IS_LTOGM_CLASSID,2);
10318e61c736SStefano Zampini   /* raise an error if already allocated */
1032674ae819SStefano Zampini   if (graph->nvtxs) {
10338e61c736SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)l2gmap),PETSC_ERR_PLIB,"BDDCGraph already initialized");
1034674ae819SStefano Zampini   }
1035674ae819SStefano Zampini   /* set number of vertices */
1036674ae819SStefano Zampini   ierr = PetscObjectReference((PetscObject)l2gmap);CHKERRQ(ierr);
1037674ae819SStefano Zampini   graph->l2gmap = l2gmap;
1038674ae819SStefano Zampini   ierr = ISLocalToGlobalMappingGetSize(l2gmap,&n);CHKERRQ(ierr);
1039674ae819SStefano Zampini   graph->nvtxs = n;
1040674ae819SStefano Zampini   /* allocate used space */
1041df48933bSStefano Zampini   ierr = PetscBTCreate(graph->nvtxs,&graph->touched);CHKERRQ(ierr);
10428e5aaad5SJed Brown   ierr = PetscMalloc7(graph->nvtxs,&graph->count,
10438e5aaad5SJed Brown                       graph->nvtxs,&graph->neighbours_set,
10448e5aaad5SJed Brown                       graph->nvtxs,&graph->subset,
10458e5aaad5SJed Brown                       graph->nvtxs,&graph->which_dof,
10468e5aaad5SJed Brown                       graph->nvtxs+1,&graph->cptr,
10478e5aaad5SJed Brown                       graph->nvtxs,&graph->queue,
10488e5aaad5SJed Brown                       graph->nvtxs,&graph->special_dof);CHKERRQ(ierr);
1049674ae819SStefano Zampini   /* zeroes memory */
1050674ae819SStefano Zampini   ierr = PetscMemzero(graph->count,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
1051674ae819SStefano Zampini   ierr = PetscMemzero(graph->subset,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
105263602bcaSStefano Zampini   /* use -1 as a default value for which_dof array */
105363602bcaSStefano Zampini   for (n=0;n<graph->nvtxs;n++) graph->which_dof[n] = -1;
1054674ae819SStefano Zampini   ierr = PetscMemzero(graph->cptr,(graph->nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr);
1055674ae819SStefano Zampini   ierr = PetscMemzero(graph->queue,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
105674e413f5SStefano Zampini   ierr = PetscMemzero(graph->special_dof,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr);
1057674ae819SStefano Zampini   /* zeroes first pointer to neighbour set */
1058674ae819SStefano Zampini   if (graph->nvtxs) {
1059674ae819SStefano Zampini     graph->neighbours_set[0] = 0;
1060674ae819SStefano Zampini   }
1061674ae819SStefano Zampini   /* zeroes workspace for values of ncc */
1062674ae819SStefano Zampini   graph->subset_ncc = 0;
10633a5b03d1SStefano Zampini   graph->subset_ref_node = 0;
1064674ae819SStefano Zampini   PetscFunctionReturn(0);
1065674ae819SStefano Zampini }
1066674ae819SStefano Zampini 
1067674ae819SStefano Zampini #undef __FUNCT__
1068674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphDestroy"
1069674ae819SStefano Zampini PetscErrorCode PCBDDCGraphDestroy(PCBDDCGraph* graph)
1070674ae819SStefano Zampini {
1071674ae819SStefano Zampini   PetscErrorCode ierr;
1072674ae819SStefano Zampini 
1073674ae819SStefano Zampini   PetscFunctionBegin;
1074674ae819SStefano Zampini   ierr = PCBDDCGraphReset(*graph);CHKERRQ(ierr);
1075674ae819SStefano Zampini   ierr = PetscFree(*graph);CHKERRQ(ierr);
1076674ae819SStefano Zampini   PetscFunctionReturn(0);
1077674ae819SStefano Zampini }
1078674ae819SStefano Zampini 
1079674ae819SStefano Zampini #undef __FUNCT__
1080674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphCreate"
1081674ae819SStefano Zampini PetscErrorCode PCBDDCGraphCreate(PCBDDCGraph *graph)
1082674ae819SStefano Zampini {
1083674ae819SStefano Zampini   PCBDDCGraph    new_graph;
1084674ae819SStefano Zampini   PetscErrorCode ierr;
1085674ae819SStefano Zampini 
1086674ae819SStefano Zampini   PetscFunctionBegin;
1087674ae819SStefano Zampini   ierr = PetscMalloc(sizeof(*new_graph),&new_graph);CHKERRQ(ierr);
1088674ae819SStefano Zampini   /* local to global mapping of dofs */
1089674ae819SStefano Zampini   new_graph->l2gmap = 0;
1090674ae819SStefano Zampini   /* vertex size */
1091674ae819SStefano Zampini   new_graph->nvtxs = 0;
1092674ae819SStefano Zampini   new_graph->n_subsets = 0;
1093674ae819SStefano Zampini   new_graph->custom_minimal_size = 1;
1094674ae819SStefano Zampini   /* zeroes ponters */
109551b0f6cfSStefano Zampini   new_graph->mirrors = 0;
109651b0f6cfSStefano Zampini   new_graph->mirrors_set = 0;
1097674ae819SStefano Zampini   new_graph->neighbours_set = 0;
1098674ae819SStefano Zampini   new_graph->subset = 0;
1099674ae819SStefano Zampini   new_graph->which_dof = 0;
110074e413f5SStefano Zampini   new_graph->special_dof = 0;
1101674ae819SStefano Zampini   new_graph->cptr = 0;
1102674ae819SStefano Zampini   new_graph->queue = 0;
1103674ae819SStefano Zampini   new_graph->count = 0;
1104674ae819SStefano Zampini   new_graph->subset_ncc = 0;
11053a5b03d1SStefano Zampini   new_graph->subset_ref_node = 0;
1106674ae819SStefano Zampini   new_graph->touched = 0;
1107674ae819SStefano Zampini   /* zeroes pointers to csr graph of local nodes connectivity (optional data) */
1108575ad6abSStefano Zampini   new_graph->nvtxs_csr = 0;
1109674ae819SStefano Zampini   new_graph->xadj = 0;
1110674ae819SStefano Zampini   new_graph->adjncy = 0;
1111674ae819SStefano Zampini   *graph = new_graph;
1112674ae819SStefano Zampini   PetscFunctionReturn(0);
1113674ae819SStefano Zampini }
1114