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__ 61b968477SStefano Zampini #define __FUNCT__ "PCBDDCGraphGetDirichletDofs" 71b968477SStefano Zampini PetscErrorCode PCBDDCGraphGetDirichletDofs(PCBDDCGraph graph,IS* dirdofs) 81b968477SStefano Zampini { 91b968477SStefano Zampini PetscInt i,size,sizer; 101b968477SStefano Zampini PetscErrorCode ierr; 111b968477SStefano Zampini 121b968477SStefano Zampini PetscFunctionBegin; 131b968477SStefano Zampini size = 0; 141b968477SStefano Zampini for (i=0;i<graph->nvtxs;i++) { 151b968477SStefano Zampini if (graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) size++; 161b968477SStefano Zampini } 171b968477SStefano Zampini ierr = MPI_Allreduce(&size,&sizer,1,MPIU_INT,MPI_LOR,PetscObjectComm((PetscObject)graph->l2gmap));CHKERRQ(ierr); 181b968477SStefano Zampini if (sizer) { 191b968477SStefano Zampini PetscInt *dirdofs_idxs; 201b968477SStefano Zampini 211b968477SStefano Zampini ierr = PetscMalloc1(size,&dirdofs_idxs);CHKERRQ(ierr); 221b968477SStefano Zampini size = 0; 231b968477SStefano Zampini for (i=0;i<graph->nvtxs;i++) { 241b968477SStefano Zampini if (graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) dirdofs_idxs[size++] = i; 251b968477SStefano Zampini } 261b968477SStefano Zampini ierr = ISCreateGeneral(PetscObjectComm((PetscObject)graph->l2gmap),size,dirdofs_idxs,PETSC_OWN_POINTER,dirdofs);CHKERRQ(ierr); 271b968477SStefano Zampini } else { 281b968477SStefano Zampini *dirdofs = NULL; 291b968477SStefano Zampini } 301b968477SStefano Zampini PetscFunctionReturn(0); 311b968477SStefano Zampini } 321b968477SStefano Zampini 331b968477SStefano Zampini #undef __FUNCT__ 34674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphASCIIView" 35e49050b4SStefano Zampini PetscErrorCode PCBDDCGraphASCIIView(PCBDDCGraph graph, PetscInt verbosity_level, PetscViewer viewer) 36674ae819SStefano Zampini { 372b510759SStefano Zampini PetscInt i,j,tabs; 3893fb5fd3SStefano Zampini PetscInt* queue_in_global_numbering; 39674ae819SStefano Zampini PetscErrorCode ierr; 40674ae819SStefano Zampini 41674ae819SStefano Zampini PetscFunctionBegin; 42674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 432b510759SStefano Zampini ierr = PetscViewerASCIIGetTab(viewer,&tabs);CHKERRQ(ierr); 44674ae819SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"--------------------------------------------------\n");CHKERRQ(ierr); 45674ae819SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 46674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Local BDDC graph for subdomain %04d\n",PetscGlobalRank);CHKERRQ(ierr); 47674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Number of vertices %d\n",graph->nvtxs);CHKERRQ(ierr); 48e49050b4SStefano Zampini if (verbosity_level > 1) { 49674ae819SStefano Zampini for (i=0;i<graph->nvtxs;i++) { 50674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"%d:\n",i);CHKERRQ(ierr); 51674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer," which_dof: %d\n",graph->which_dof[i]);CHKERRQ(ierr); 5274e413f5SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer," special_dof: %d\n",graph->special_dof[i]);CHKERRQ(ierr); 53674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer," neighbours: %d\n",graph->count[i]);CHKERRQ(ierr); 542b510759SStefano Zampini ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 55674ae819SStefano Zampini if (graph->count[i]) { 56674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer," set of neighbours:");CHKERRQ(ierr); 57674ae819SStefano Zampini for (j=0;j<graph->count[i];j++) { 58674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d",graph->neighbours_set[i][j]);CHKERRQ(ierr); 59674ae819SStefano Zampini } 60674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr); 61674ae819SStefano Zampini } 622b510759SStefano Zampini ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr); 632b510759SStefano Zampini ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 6451b0f6cfSStefano Zampini if (graph->mirrors) { 6551b0f6cfSStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer," mirrors: %d\n",graph->mirrors[i]);CHKERRQ(ierr); 6651b0f6cfSStefano Zampini if (graph->mirrors[i]) { 672b510759SStefano Zampini ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 6851b0f6cfSStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer," set of mirrors:");CHKERRQ(ierr); 6951b0f6cfSStefano Zampini for (j=0;j<graph->mirrors[i];j++) { 7051b0f6cfSStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d",graph->mirrors_set[i][j]);CHKERRQ(ierr); 7151b0f6cfSStefano Zampini } 7251b0f6cfSStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr); 732b510759SStefano Zampini ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr); 742b510759SStefano Zampini ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 7551b0f6cfSStefano Zampini } 7651b0f6cfSStefano Zampini } 77e49050b4SStefano Zampini if (verbosity_level > 2) { 78674ae819SStefano Zampini if (graph->xadj && graph->adjncy) { 79674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer," local adj list:");CHKERRQ(ierr); 802b510759SStefano Zampini ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 81674ae819SStefano Zampini for (j=graph->xadj[i];j<graph->xadj[i+1];j++) { 82674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d",graph->adjncy[j]);CHKERRQ(ierr); 83674ae819SStefano Zampini } 84674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr); 852b510759SStefano Zampini ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr); 862b510759SStefano Zampini ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 87674ae819SStefano Zampini } 88e49050b4SStefano Zampini } 89674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer," interface subset id: %d\n",graph->subset[i]);CHKERRQ(ierr); 90674ae819SStefano Zampini if (graph->subset[i] && graph->subset_ncc) { 91674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer," ncc for subset: %d\n",graph->subset_ncc[graph->subset[i]-1]);CHKERRQ(ierr); 92674ae819SStefano Zampini } 93674ae819SStefano Zampini } 94e49050b4SStefano Zampini } 95674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Total number of connected components %d\n",graph->ncc);CHKERRQ(ierr); 96785e854fSJed Brown ierr = PetscMalloc1(graph->cptr[graph->ncc],&queue_in_global_numbering);CHKERRQ(ierr); 9793fb5fd3SStefano Zampini ierr = ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_in_global_numbering);CHKERRQ(ierr); 98674ae819SStefano Zampini for (i=0;i<graph->ncc;i++) { 991ce3d396SStefano Zampini PetscInt node_num=graph->queue[graph->cptr[i]]; 1001ce3d396SStefano Zampini PetscBool printcc = PETSC_FALSE; 101674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d (neighs:",i);CHKERRQ(ierr); 1022b510759SStefano Zampini ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 103674ae819SStefano Zampini for (j=0;j<graph->count[node_num];j++) { 104674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d",graph->neighbours_set[node_num][j]);CHKERRQ(ierr); 105674ae819SStefano Zampini } 106674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"):");CHKERRQ(ierr); 1071ce3d396SStefano Zampini if (verbosity_level > 1) { 1081ce3d396SStefano Zampini printcc = PETSC_TRUE; 109e635b14bSStefano Zampini } else if (graph->count[node_num] > 1 || (graph->count[node_num] == 1 && graph->special_dof[node_num] == PCBDDCGRAPH_NEUMANN_MARK)) { 1101ce3d396SStefano Zampini printcc = PETSC_TRUE; 1111ce3d396SStefano Zampini } 1121ce3d396SStefano Zampini if (printcc) { 113674ae819SStefano Zampini for (j=graph->cptr[i];j<graph->cptr[i+1];j++) { 11493fb5fd3SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d (%d)",graph->queue[j],queue_in_global_numbering[j]);CHKERRQ(ierr); 115674ae819SStefano Zampini } 1161ce3d396SStefano Zampini } 117674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr); 1182b510759SStefano Zampini ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr); 1192b510759SStefano Zampini ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 120674ae819SStefano Zampini } 12193fb5fd3SStefano Zampini ierr = PetscFree(queue_in_global_numbering);CHKERRQ(ierr); 122e49050b4SStefano Zampini if (graph->custom_minimal_size > 1 && verbosity_level > 1) { 123674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Custom minimal size %d\n",graph->custom_minimal_size);CHKERRQ(ierr); 124674ae819SStefano Zampini } 125674ae819SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 126674ae819SStefano Zampini PetscFunctionReturn(0); 127674ae819SStefano Zampini } 128674ae819SStefano Zampini 129674ae819SStefano Zampini #undef __FUNCT__ 130674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphGetCandidatesIS" 131a873d5d0SStefano Zampini PetscErrorCode PCBDDCGraphGetCandidatesIS(PCBDDCGraph graph, PetscInt *n_faces, IS *FacesIS[], PetscInt *n_edges, IS *EdgesIS[], IS *VerticesIS) 132674ae819SStefano Zampini { 133674ae819SStefano Zampini IS *ISForFaces,*ISForEdges,ISForVertices; 134674ae819SStefano Zampini PetscInt i,j,nfc,nec,nvc,*idx; 135674ae819SStefano Zampini PetscBool twodim_flag; 136674ae819SStefano Zampini PetscErrorCode ierr; 137674ae819SStefano Zampini 138674ae819SStefano Zampini PetscFunctionBegin; 139674ae819SStefano Zampini /* loop on ccs to evalute number of faces, edges and vertices */ 140674ae819SStefano Zampini nfc = 0; 141674ae819SStefano Zampini nec = 0; 142674ae819SStefano Zampini nvc = 0; 143674ae819SStefano Zampini twodim_flag = PETSC_FALSE; 144674ae819SStefano Zampini for (i=0;i<graph->ncc;i++) { 1451e1f2d93SStefano Zampini PetscInt repdof = graph->queue[graph->cptr[i]]; 146674ae819SStefano Zampini if (graph->cptr[i+1]-graph->cptr[i] > graph->custom_minimal_size) { 147e635b14bSStefano Zampini if (graph->count[repdof] == 1 && graph->special_dof[repdof] != PCBDDCGRAPH_NEUMANN_MARK) { 148674ae819SStefano Zampini nfc++; 149674ae819SStefano Zampini } else { /* note that nec will be zero in 2d */ 150674ae819SStefano Zampini nec++; 151674ae819SStefano Zampini } 152674ae819SStefano Zampini } else { 153674ae819SStefano Zampini nvc += graph->cptr[i+1]-graph->cptr[i]; 154674ae819SStefano Zampini } 155674ae819SStefano Zampini } 156aa9fcddfSStefano Zampini j=0; 157aa9fcddfSStefano Zampini ierr = MPI_Allreduce(&nec,&j,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)graph->l2gmap));CHKERRQ(ierr); 158aa9fcddfSStefano Zampini if (!j) { /* we are in a 2D case -> no faces, only edges */ 159674ae819SStefano Zampini nec = nfc; 160674ae819SStefano Zampini nfc = 0; 161674ae819SStefano Zampini twodim_flag = PETSC_TRUE; 162674ae819SStefano Zampini } 163674ae819SStefano Zampini /* allocate IS arrays for faces, edges. Vertices need a single index set. */ 164*cf5a6209SStefano Zampini if (FacesIS) { 165785e854fSJed Brown ierr = PetscMalloc1(nfc,&ISForFaces);CHKERRQ(ierr); 166*cf5a6209SStefano Zampini } 167*cf5a6209SStefano Zampini if (EdgesIS) { 168785e854fSJed Brown ierr = PetscMalloc1(nec,&ISForEdges);CHKERRQ(ierr); 169*cf5a6209SStefano Zampini } 170*cf5a6209SStefano Zampini if (VerticesIS) { 171785e854fSJed Brown ierr = PetscMalloc1(nvc,&idx);CHKERRQ(ierr); 172*cf5a6209SStefano Zampini } 173*cf5a6209SStefano Zampini 174674ae819SStefano Zampini /* loop on ccs to compute index sets for faces and edges */ 175674ae819SStefano Zampini nfc = 0; 176674ae819SStefano Zampini nec = 0; 177674ae819SStefano Zampini for (i=0;i<graph->ncc;i++) { 1781e1f2d93SStefano Zampini PetscInt repdof = graph->queue[graph->cptr[i]]; 179674ae819SStefano Zampini if (graph->cptr[i+1]-graph->cptr[i] > graph->custom_minimal_size) { 180e635b14bSStefano Zampini if (graph->count[repdof] == 1 && graph->special_dof[repdof] != PCBDDCGRAPH_NEUMANN_MARK) { 181674ae819SStefano Zampini if (twodim_flag) { 182*cf5a6209SStefano Zampini if (EdgesIS) { 183674ae819SStefano 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); 184*cf5a6209SStefano Zampini } 185674ae819SStefano Zampini nec++; 186674ae819SStefano Zampini } else { 187*cf5a6209SStefano Zampini if (FacesIS) { 188674ae819SStefano 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); 189*cf5a6209SStefano Zampini } 190674ae819SStefano Zampini nfc++; 191674ae819SStefano Zampini } 192674ae819SStefano Zampini } else { 193*cf5a6209SStefano Zampini if (EdgesIS) { 194674ae819SStefano Zampini ierr = ISCreateGeneral(PETSC_COMM_SELF,graph->cptr[i+1]-graph->cptr[i],&graph->queue[graph->cptr[i]],PETSC_COPY_VALUES,&ISForEdges[nec]);CHKERRQ(ierr); 195*cf5a6209SStefano Zampini } 196674ae819SStefano Zampini nec++; 197674ae819SStefano Zampini } 198674ae819SStefano Zampini } 199674ae819SStefano Zampini } 200674ae819SStefano Zampini /* index set for vertices */ 201*cf5a6209SStefano Zampini if (VerticesIS) { 202b8ffe317SStefano Zampini nvc = 0; 203674ae819SStefano Zampini for (i=0;i<graph->ncc;i++) { 204674ae819SStefano Zampini if (graph->cptr[i+1]-graph->cptr[i] <= graph->custom_minimal_size) { 205674ae819SStefano Zampini for (j=graph->cptr[i];j<graph->cptr[i+1];j++) { 206674ae819SStefano Zampini idx[nvc]=graph->queue[j]; 207674ae819SStefano Zampini nvc++; 208674ae819SStefano Zampini } 209674ae819SStefano Zampini } 210674ae819SStefano Zampini } 211674ae819SStefano Zampini /* sort vertex set (by local ordering) */ 212674ae819SStefano Zampini ierr = PetscSortInt(nvc,idx);CHKERRQ(ierr); 213674ae819SStefano Zampini ierr = ISCreateGeneral(PETSC_COMM_SELF,nvc,idx,PETSC_OWN_POINTER,&ISForVertices);CHKERRQ(ierr); 214674ae819SStefano Zampini } 215674ae819SStefano Zampini /* get back info */ 216a873d5d0SStefano Zampini if (n_faces) *n_faces = nfc; 217d7796978SStefano Zampini if (FacesIS) { 218674ae819SStefano Zampini *FacesIS = ISForFaces; 219d7796978SStefano Zampini } 220a873d5d0SStefano Zampini if (n_edges) *n_edges = nec; 221d7796978SStefano Zampini if (EdgesIS) { 222674ae819SStefano Zampini *EdgesIS = ISForEdges; 223d7796978SStefano Zampini } 224d7796978SStefano Zampini if (VerticesIS) { 225674ae819SStefano Zampini *VerticesIS = ISForVertices; 226d7796978SStefano Zampini } 227674ae819SStefano Zampini PetscFunctionReturn(0); 228674ae819SStefano Zampini } 229674ae819SStefano Zampini 230674ae819SStefano Zampini #undef __FUNCT__ 231674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphComputeConnectedComponents" 232674ae819SStefano Zampini PetscErrorCode PCBDDCGraphComputeConnectedComponents(PCBDDCGraph graph) 233674ae819SStefano Zampini { 2344a060362SStefano Zampini PetscBool adapt_interface_reduced; 235674ae819SStefano Zampini MPI_Comm interface_comm; 2364a060362SStefano Zampini PetscMPIInt size; 237674ae819SStefano Zampini PetscErrorCode ierr; 238674ae819SStefano Zampini 239674ae819SStefano Zampini PetscFunctionBegin; 240674ae819SStefano Zampini /* compute connected components locally */ 241674ae819SStefano Zampini ierr = PetscObjectGetComm((PetscObject)(graph->l2gmap),&interface_comm);CHKERRQ(ierr); 242674ae819SStefano Zampini ierr = PCBDDCGraphComputeConnectedComponentsLocal(graph);CHKERRQ(ierr); 243674ae819SStefano Zampini /* check consistency of connected components among neighbouring subdomains -> it adapt them in case it is needed */ 2444a060362SStefano Zampini ierr = MPI_Comm_size(interface_comm,&size);CHKERRQ(ierr); 2454a060362SStefano Zampini adapt_interface_reduced = PETSC_FALSE; 2464a060362SStefano Zampini if (size > 1) { 2474a060362SStefano Zampini PetscInt i; 2484a060362SStefano Zampini PetscBool adapt_interface = PETSC_FALSE; 249674ae819SStefano Zampini for (i=0;i<graph->n_subsets;i++) { 250674ae819SStefano Zampini /* We are not sure that on a given subset of the local interface, 251674ae819SStefano Zampini with two connected components, the latters be the same among sharing subdomains */ 252674ae819SStefano Zampini if (graph->subset_ncc[i] > 1) { 2534a060362SStefano Zampini adapt_interface = PETSC_TRUE; 254674ae819SStefano Zampini break; 255674ae819SStefano Zampini } 256674ae819SStefano Zampini } 2574a060362SStefano Zampini ierr = MPI_Allreduce(&adapt_interface,&adapt_interface_reduced,1,MPIU_BOOL,MPI_LOR,interface_comm);CHKERRQ(ierr); 2584a060362SStefano Zampini } 259674ae819SStefano Zampini 260674ae819SStefano Zampini if (graph->n_subsets && adapt_interface_reduced) { 2615b1b9aeaSStefano Zampini MPI_Request *send_requests; 2625b1b9aeaSStefano Zampini MPI_Request *recv_requests; 2635b1b9aeaSStefano Zampini PetscInt *aux_new_xadj,*new_xadj,*new_adjncy,**temp_buffer; 2645b1b9aeaSStefano Zampini PetscInt *old_xadj,*old_adjncy; 2654a060362SStefano Zampini PetscInt i,j,k,s,sum_requests,buffer_size,size_of_recv,temp_buffer_size; 2665b1b9aeaSStefano Zampini PetscMPIInt rank,neigh,tag,mpi_buffer_size; 2675b1b9aeaSStefano Zampini PetscInt *cum_recv_counts,*subset_to_nodes_indices,*recv_buffer_subset,*nodes_to_temp_buffer_indices; 2685b1b9aeaSStefano Zampini PetscInt *send_buffer,*recv_buffer,*queue_in_global_numbering,*sizes_of_sends,*add_to_subset; 2695b1b9aeaSStefano Zampini PetscInt start_of_recv,start_of_send,size_of_send,global_subset_counter,ins_val; 2705b1b9aeaSStefano Zampini PetscBool *subset_cc_adapt,same_set; 2715b1b9aeaSStefano Zampini 272674ae819SStefano Zampini /* Retrict adjacency graph using information from previously computed connected components */ 273785e854fSJed Brown ierr = PetscMalloc1(graph->nvtxs,&aux_new_xadj);CHKERRQ(ierr); 274674ae819SStefano Zampini for (i=0;i<graph->nvtxs;i++) { 275674ae819SStefano Zampini aux_new_xadj[i]=1; 276674ae819SStefano Zampini } 277674ae819SStefano Zampini for (i=0;i<graph->ncc;i++) { 278674ae819SStefano Zampini k = graph->cptr[i+1]-graph->cptr[i]; 279674ae819SStefano Zampini for (j=0;j<k;j++) { 280674ae819SStefano Zampini aux_new_xadj[graph->queue[graph->cptr[i]+j]]=k; 281674ae819SStefano Zampini } 282674ae819SStefano Zampini } 283674ae819SStefano Zampini j = 0; 284674ae819SStefano Zampini for (i=0;i<graph->nvtxs;i++) { 285674ae819SStefano Zampini j += aux_new_xadj[i]; 286674ae819SStefano Zampini } 287854ce69bSBarry Smith ierr = PetscMalloc1(graph->nvtxs+1,&new_xadj);CHKERRQ(ierr); 288785e854fSJed Brown ierr = PetscMalloc1(j,&new_adjncy);CHKERRQ(ierr); 289674ae819SStefano Zampini new_xadj[0]=0; 290674ae819SStefano Zampini for (i=0;i<graph->nvtxs;i++) { 291674ae819SStefano Zampini new_xadj[i+1]=new_xadj[i]+aux_new_xadj[i]; 292674ae819SStefano Zampini if (aux_new_xadj[i]==1) { 293674ae819SStefano Zampini new_adjncy[new_xadj[i]]=i; 294674ae819SStefano Zampini } 295674ae819SStefano Zampini } 296674ae819SStefano Zampini ierr = PetscFree(aux_new_xadj);CHKERRQ(ierr); 297674ae819SStefano Zampini for (i=0;i<graph->ncc;i++) { 298674ae819SStefano Zampini k = graph->cptr[i+1]-graph->cptr[i]; 299674ae819SStefano Zampini for (j=0;j<k;j++) { 300674ae819SStefano Zampini ierr = PetscMemcpy(&new_adjncy[new_xadj[graph->queue[graph->cptr[i]+j]]],&graph->queue[graph->cptr[i]],k*sizeof(PetscInt));CHKERRQ(ierr); 301674ae819SStefano Zampini } 302674ae819SStefano Zampini } 3035b1b9aeaSStefano Zampini /* set temporarly new CSR into graph */ 3045b1b9aeaSStefano Zampini old_xadj = graph->xadj; 3055b1b9aeaSStefano Zampini old_adjncy = graph->adjncy; 306674ae819SStefano Zampini graph->xadj = new_xadj; 307674ae819SStefano Zampini graph->adjncy = new_adjncy; 308674ae819SStefano Zampini /* allocate some space */ 309674ae819SStefano Zampini ierr = MPI_Comm_rank(interface_comm,&rank);CHKERRQ(ierr); 310854ce69bSBarry Smith ierr = PetscMalloc1(graph->n_subsets+1,&cum_recv_counts);CHKERRQ(ierr); 311674ae819SStefano Zampini ierr = PetscMemzero(cum_recv_counts,(graph->n_subsets+1)*sizeof(*cum_recv_counts));CHKERRQ(ierr); 312785e854fSJed Brown ierr = PetscMalloc1(graph->n_subsets,&subset_to_nodes_indices);CHKERRQ(ierr); 313674ae819SStefano Zampini /* first count how many neighbours per connected component I will receive from */ 314674ae819SStefano Zampini cum_recv_counts[0]=0; 315674ae819SStefano Zampini for (i=1;i<graph->n_subsets+1;i++) { 316674ae819SStefano Zampini j = 0; 317674ae819SStefano Zampini while (graph->subset[j] != i) { 318674ae819SStefano Zampini j++; 319674ae819SStefano Zampini } 320674ae819SStefano Zampini subset_to_nodes_indices[i-1]=j; 321674ae819SStefano Zampini /* We don't want sends/recvs_to/from_self -> here I don't count myself */ 322674ae819SStefano Zampini cum_recv_counts[i]=cum_recv_counts[i-1]+graph->count[j]; 323674ae819SStefano Zampini } 324785e854fSJed Brown ierr = PetscMalloc1(2*cum_recv_counts[graph->n_subsets],&recv_buffer_subset);CHKERRQ(ierr); 325785e854fSJed Brown ierr = PetscMalloc1(cum_recv_counts[graph->n_subsets],&send_requests);CHKERRQ(ierr); 326785e854fSJed Brown ierr = PetscMalloc1(cum_recv_counts[graph->n_subsets],&recv_requests);CHKERRQ(ierr); 327674ae819SStefano Zampini for (i=0;i<cum_recv_counts[graph->n_subsets];i++) { 328674ae819SStefano Zampini send_requests[i]=MPI_REQUEST_NULL; 329674ae819SStefano Zampini recv_requests[i]=MPI_REQUEST_NULL; 330674ae819SStefano Zampini } 331674ae819SStefano Zampini /* exchange with my neighbours the number of my connected components on the shared interface */ 332674ae819SStefano Zampini sum_requests = 0; 333674ae819SStefano Zampini for (i=0;i<graph->n_subsets;i++) { 334674ae819SStefano Zampini j = subset_to_nodes_indices[i]; 335f0321c90SStefano Zampini ierr = PetscMPIIntCast(3*graph->subset_ref_node[i],&tag);CHKERRQ(ierr); 336674ae819SStefano Zampini for (k=0;k<graph->count[j];k++) { 337674ae819SStefano Zampini ierr = PetscMPIIntCast(graph->neighbours_set[j][k],&neigh);CHKERRQ(ierr); 338674ae819SStefano Zampini ierr = MPI_Isend(&graph->subset_ncc[i],1,MPIU_INT,neigh,tag,interface_comm,&send_requests[sum_requests]);CHKERRQ(ierr); 339674ae819SStefano Zampini ierr = MPI_Irecv(&recv_buffer_subset[sum_requests],1,MPIU_INT,neigh,tag,interface_comm,&recv_requests[sum_requests]);CHKERRQ(ierr); 340674ae819SStefano Zampini sum_requests++; 341674ae819SStefano Zampini } 342674ae819SStefano Zampini } 343674ae819SStefano Zampini ierr = MPI_Waitall(sum_requests,recv_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 344674ae819SStefano Zampini ierr = MPI_Waitall(sum_requests,send_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 345674ae819SStefano Zampini /* determine the connected component I need to adapt */ 346785e854fSJed Brown ierr = PetscMalloc1(graph->n_subsets,&subset_cc_adapt);CHKERRQ(ierr); 347674ae819SStefano Zampini ierr = PetscMemzero(subset_cc_adapt,graph->n_subsets*sizeof(*subset_cc_adapt));CHKERRQ(ierr); 348674ae819SStefano Zampini for (i=0;i<graph->n_subsets;i++) { 349674ae819SStefano Zampini for (j=cum_recv_counts[i];j<cum_recv_counts[i+1];j++){ 350674ae819SStefano Zampini /* The first condition is natural (someone has a different number of ccs than me), the second one is just to be safe */ 351674ae819SStefano Zampini if (graph->subset_ncc[i] != recv_buffer_subset[j] || graph->subset_ncc[i] > 1) { 352674ae819SStefano Zampini subset_cc_adapt[i] = PETSC_TRUE; 353674ae819SStefano Zampini break; 354674ae819SStefano Zampini } 355674ae819SStefano Zampini } 356674ae819SStefano Zampini } 357674ae819SStefano Zampini buffer_size = 0; 358674ae819SStefano Zampini for (i=0;i<graph->n_subsets;i++) { 359674ae819SStefano Zampini if (subset_cc_adapt[i]) { 360674ae819SStefano Zampini for (j=i;j<graph->ncc;j++) { 361674ae819SStefano Zampini if (graph->subset[graph->queue[graph->cptr[j]]] == i+1) { /* WARNING -> subset values goes from 1 to graph->n_subsets included */ 362674ae819SStefano Zampini buffer_size += 1 + graph->cptr[j+1]-graph->cptr[j]; 363674ae819SStefano Zampini } 364674ae819SStefano Zampini } 365674ae819SStefano Zampini } 366674ae819SStefano Zampini } 367785e854fSJed Brown ierr = PetscMalloc1(buffer_size,&send_buffer);CHKERRQ(ierr); 368674ae819SStefano Zampini /* now get from neighbours their ccs (in global numbering) and adapt them (in case it is needed) */ 369785e854fSJed Brown ierr = PetscMalloc1(graph->cptr[graph->ncc],&queue_in_global_numbering);CHKERRQ(ierr); 370674ae819SStefano Zampini ierr = ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_in_global_numbering);CHKERRQ(ierr); 371674ae819SStefano Zampini /* determine how much data to send (size of each queue plus the global indices) and communicate it to neighbours */ 372785e854fSJed Brown ierr = PetscMalloc1(graph->n_subsets,&sizes_of_sends);CHKERRQ(ierr); 373674ae819SStefano Zampini ierr = PetscMemzero(sizes_of_sends,graph->n_subsets*sizeof(*sizes_of_sends));CHKERRQ(ierr); 374674ae819SStefano Zampini sum_requests = 0; 375674ae819SStefano Zampini start_of_send = 0; 376674ae819SStefano Zampini start_of_recv = cum_recv_counts[graph->n_subsets]; 377674ae819SStefano Zampini for (i=0;i<graph->n_subsets;i++) { 378674ae819SStefano Zampini if (subset_cc_adapt[i]) { 379674ae819SStefano Zampini size_of_send = 0; 380674ae819SStefano Zampini for (j=i;j<graph->ncc;j++) { 381674ae819SStefano Zampini if (graph->subset[graph->queue[graph->cptr[j]]] == i+1) { /* WARNING -> subset values goes from 1 to graph->n_subsets included */ 382674ae819SStefano Zampini send_buffer[start_of_send+size_of_send]=graph->cptr[j+1]-graph->cptr[j]; 383674ae819SStefano Zampini size_of_send += 1; 384674ae819SStefano Zampini ierr = PetscMemcpy(&send_buffer[start_of_send+size_of_send], 385674ae819SStefano Zampini &queue_in_global_numbering[graph->cptr[j]], 386674ae819SStefano Zampini (graph->cptr[j+1]-graph->cptr[j])*sizeof(*send_buffer));CHKERRQ(ierr); 387674ae819SStefano Zampini size_of_send = size_of_send+graph->cptr[j+1]-graph->cptr[j]; 388674ae819SStefano Zampini } 389674ae819SStefano Zampini } 390674ae819SStefano Zampini j = subset_to_nodes_indices[i]; 391674ae819SStefano Zampini sizes_of_sends[i] = size_of_send; 392f0321c90SStefano Zampini ierr = PetscMPIIntCast(3*graph->subset_ref_node[i]+1,&tag);CHKERRQ(ierr); 393674ae819SStefano Zampini for (k=0;k<graph->count[j];k++) { 394674ae819SStefano Zampini ierr = PetscMPIIntCast(graph->neighbours_set[j][k],&neigh);CHKERRQ(ierr); 395674ae819SStefano Zampini ierr = MPI_Isend(&sizes_of_sends[i],1,MPIU_INT,neigh,tag,interface_comm,&send_requests[sum_requests]);CHKERRQ(ierr); 396674ae819SStefano 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); 397674ae819SStefano Zampini sum_requests++; 398674ae819SStefano Zampini } 399674ae819SStefano Zampini start_of_send += size_of_send; 400674ae819SStefano Zampini } 401674ae819SStefano Zampini } 402674ae819SStefano Zampini ierr = MPI_Waitall(sum_requests,send_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 403674ae819SStefano Zampini ierr = MPI_Waitall(sum_requests,recv_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 404674ae819SStefano Zampini buffer_size = 0; 405674ae819SStefano Zampini for (k=0;k<sum_requests;k++) { 406674ae819SStefano Zampini buffer_size += recv_buffer_subset[start_of_recv+k]; 407674ae819SStefano Zampini } 408785e854fSJed Brown ierr = PetscMalloc1(buffer_size,&recv_buffer);CHKERRQ(ierr); 409674ae819SStefano Zampini /* now exchange the data */ 410674ae819SStefano Zampini start_of_recv = 0; 411674ae819SStefano Zampini start_of_send = 0; 412674ae819SStefano Zampini sum_requests = 0; 413674ae819SStefano Zampini for (i=0;i<graph->n_subsets;i++) { 414674ae819SStefano Zampini if (subset_cc_adapt[i]) { 415674ae819SStefano Zampini size_of_send = sizes_of_sends[i]; 416674ae819SStefano Zampini j = subset_to_nodes_indices[i]; 417f0321c90SStefano Zampini ierr = PetscMPIIntCast(3*graph->subset_ref_node[i]+2,&tag);CHKERRQ(ierr); 418674ae819SStefano Zampini for (k=0;k<graph->count[j];k++) { 419674ae819SStefano Zampini ierr = PetscMPIIntCast(graph->neighbours_set[j][k],&neigh);CHKERRQ(ierr); 420674ae819SStefano Zampini ierr = PetscMPIIntCast(size_of_send,&mpi_buffer_size);CHKERRQ(ierr); 421674ae819SStefano Zampini ierr = MPI_Isend(&send_buffer[start_of_send],mpi_buffer_size,MPIU_INT,neigh,tag,interface_comm,&send_requests[sum_requests]);CHKERRQ(ierr); 422674ae819SStefano Zampini size_of_recv = recv_buffer_subset[cum_recv_counts[graph->n_subsets]+sum_requests]; 423674ae819SStefano Zampini ierr = PetscMPIIntCast(size_of_recv,&mpi_buffer_size);CHKERRQ(ierr); 424674ae819SStefano Zampini ierr = MPI_Irecv(&recv_buffer[start_of_recv],mpi_buffer_size,MPIU_INT,neigh,tag,interface_comm,&recv_requests[sum_requests]);CHKERRQ(ierr); 425674ae819SStefano Zampini start_of_recv += size_of_recv; 426674ae819SStefano Zampini sum_requests++; 427674ae819SStefano Zampini } 428674ae819SStefano Zampini start_of_send += size_of_send; 429674ae819SStefano Zampini } 430674ae819SStefano Zampini } 431674ae819SStefano Zampini ierr = MPI_Waitall(sum_requests,recv_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 432674ae819SStefano Zampini ierr = MPI_Waitall(sum_requests,send_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 433674ae819SStefano Zampini for (j=0;j<buffer_size;) { 434d52c4852SStefano Zampini PetscInt rbs = recv_buffer[j]; 435d52c4852SStefano Zampini 436d52c4852SStefano Zampini ierr = ISGlobalToLocalMappingApply(graph->l2gmap,IS_GTOLM_DROP,rbs,&recv_buffer[j+1],&recv_buffer[j],&recv_buffer[j+1]);CHKERRQ(ierr); 437d52c4852SStefano Zampini if (rbs != recv_buffer[j]) { 438d52c4852SStefano Zampini SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Mapped received buffer size %d differs from input size %d\n",rbs,recv_buffer[j]); 439d52c4852SStefano Zampini } 44051b0f6cfSStefano Zampini /* we need to adapt the output of GlobalToLocal mapping if there are mirrored nodes */ 44151b0f6cfSStefano Zampini if (graph->mirrors) { 44251b0f6cfSStefano Zampini PetscBool mirrored_found=PETSC_FALSE; 44351b0f6cfSStefano Zampini for (k=0;k<recv_buffer[j];k++) { 44451b0f6cfSStefano Zampini if (graph->mirrors[recv_buffer[j+k+1]]) { 44551b0f6cfSStefano Zampini mirrored_found=PETSC_TRUE; 44651b0f6cfSStefano Zampini recv_buffer[j+k+1]=graph->mirrors_set[recv_buffer[j+k+1]][0]; 44751b0f6cfSStefano Zampini } 44851b0f6cfSStefano Zampini } 44951b0f6cfSStefano Zampini if (mirrored_found) { 45051b0f6cfSStefano Zampini ierr = PetscSortInt(recv_buffer[j],&recv_buffer[j+1]);CHKERRQ(ierr); 45151b0f6cfSStefano Zampini k=0; 45251b0f6cfSStefano Zampini while (k<recv_buffer[j]) { 453d52c4852SStefano Zampini PetscInt nmir = graph->mirrors[recv_buffer[j+1+k]]; 454d52c4852SStefano Zampini for (s=1;s<nmir;s++) { 45551b0f6cfSStefano Zampini recv_buffer[j+1+k+s] = graph->mirrors_set[recv_buffer[j+1+k]][s]; 45651b0f6cfSStefano Zampini } 457d52c4852SStefano Zampini k += nmir; 45851b0f6cfSStefano Zampini } 45951b0f6cfSStefano Zampini } 46051b0f6cfSStefano Zampini } 461674ae819SStefano Zampini k = recv_buffer[j]+1; 462674ae819SStefano Zampini j += k; 463674ae819SStefano Zampini } 464674ae819SStefano Zampini sum_requests = cum_recv_counts[graph->n_subsets]; 465674ae819SStefano Zampini start_of_recv = 0; 466785e854fSJed Brown ierr = PetscMalloc1(graph->nvtxs,&nodes_to_temp_buffer_indices);CHKERRQ(ierr); 467674ae819SStefano Zampini global_subset_counter = 0; 468674ae819SStefano Zampini for (i=0;i<graph->n_subsets;i++) { 469674ae819SStefano Zampini if (subset_cc_adapt[i]) { 470674ae819SStefano Zampini temp_buffer_size = 0; 471674ae819SStefano Zampini /* find nodes on the shared interface we need to adapt */ 472674ae819SStefano Zampini for (j=0;j<graph->nvtxs;j++) { 473674ae819SStefano Zampini if (graph->subset[j]==i+1) { 474674ae819SStefano Zampini nodes_to_temp_buffer_indices[j] = temp_buffer_size; 475674ae819SStefano Zampini temp_buffer_size++; 476674ae819SStefano Zampini } else { 477674ae819SStefano Zampini nodes_to_temp_buffer_indices[j] = -1; 478674ae819SStefano Zampini } 479674ae819SStefano Zampini } 480674ae819SStefano Zampini /* allocate some temporary space */ 481785e854fSJed Brown ierr = PetscMalloc1(temp_buffer_size,&temp_buffer);CHKERRQ(ierr); 482785e854fSJed Brown ierr = PetscMalloc1(temp_buffer_size*(cum_recv_counts[i+1]-cum_recv_counts[i]),&temp_buffer[0]);CHKERRQ(ierr); 483674ae819SStefano Zampini ierr = PetscMemzero(temp_buffer[0],temp_buffer_size*(cum_recv_counts[i+1]-cum_recv_counts[i])*sizeof(PetscInt));CHKERRQ(ierr); 484674ae819SStefano Zampini for (j=1;j<temp_buffer_size;j++) { 485674ae819SStefano Zampini temp_buffer[j] = temp_buffer[j-1]+cum_recv_counts[i+1]-cum_recv_counts[i]; 486674ae819SStefano Zampini } 487674ae819SStefano Zampini /* analyze contributions from neighbouring subdomains for i-th conn comp 488674ae819SStefano Zampini temp buffer structure: 489674ae819SStefano Zampini supposing part of the interface has dimension 5 (for example with global dofs 0,1,2,3,4) 490674ae819SStefano Zampini 3 neighs procs with structured connected components: 491674ae819SStefano Zampini neigh 0: [0 1 4], [2 3]; (2 connected components) 492674ae819SStefano Zampini neigh 1: [0 1], [2 3 4]; (2 connected components) 493674ae819SStefano Zampini neigh 2: [0 4], [1], [2 3]; (3 connected components) 494674ae819SStefano Zampini tempbuffer (row-oriented) will be filled as: 495674ae819SStefano Zampini [ 0, 0, 0; 496674ae819SStefano Zampini 0, 0, 1; 497674ae819SStefano Zampini 1, 1, 2; 498674ae819SStefano Zampini 1, 1, 2; 499674ae819SStefano Zampini 0, 1, 0; ]; 500674ae819SStefano Zampini This way we can simply find intersections of ccs among neighs. 501674ae819SStefano Zampini For the example above, the graph->subset array will be modified to reproduce the following 4 connected components [0], [1], [2 3], [4]; 502674ae819SStefano Zampini */ 503674ae819SStefano Zampini for (j=0;j<cum_recv_counts[i+1]-cum_recv_counts[i];j++) { 504674ae819SStefano Zampini ins_val = 0; 505674ae819SStefano Zampini size_of_recv = recv_buffer_subset[sum_requests]; /* total size of recv from neighs */ 506674ae819SStefano Zampini for (buffer_size=0;buffer_size<size_of_recv;) { /* loop until all data from neighs has been taken into account */ 507674ae819SStefano Zampini for (k=1;k<recv_buffer[buffer_size+start_of_recv]+1;k++) { /* filling properly temp_buffer using data from a single recv */ 508674ae819SStefano Zampini temp_buffer[nodes_to_temp_buffer_indices[recv_buffer[start_of_recv+buffer_size+k]]][j] = ins_val; 509674ae819SStefano Zampini } 510674ae819SStefano Zampini buffer_size += k; 511674ae819SStefano Zampini ins_val++; 512674ae819SStefano Zampini } 513674ae819SStefano Zampini start_of_recv += size_of_recv; 514674ae819SStefano Zampini sum_requests++; 515674ae819SStefano Zampini } 516785e854fSJed Brown ierr = PetscMalloc1(temp_buffer_size,&add_to_subset);CHKERRQ(ierr); 517674ae819SStefano Zampini ierr = PetscMemzero(add_to_subset,temp_buffer_size*sizeof(*add_to_subset));CHKERRQ(ierr); 518674ae819SStefano Zampini for (j=0;j<temp_buffer_size;j++) { 519674ae819SStefano Zampini if (!add_to_subset[j]) { /* found a new cc */ 520674ae819SStefano Zampini global_subset_counter++; 521674ae819SStefano Zampini add_to_subset[j] = global_subset_counter; 522674ae819SStefano Zampini for (k=j+1;k<temp_buffer_size;k++) { /* check for other nodes in new cc */ 523674ae819SStefano Zampini same_set = PETSC_TRUE; 524674ae819SStefano Zampini for (s=0;s<cum_recv_counts[i+1]-cum_recv_counts[i];s++) { 525674ae819SStefano Zampini if (temp_buffer[j][s]!=temp_buffer[k][s]) { 526674ae819SStefano Zampini same_set = PETSC_FALSE; 527674ae819SStefano Zampini break; 528674ae819SStefano Zampini } 529674ae819SStefano Zampini } 530674ae819SStefano Zampini if (same_set) { 531674ae819SStefano Zampini add_to_subset[k] = global_subset_counter; 532674ae819SStefano Zampini } 533674ae819SStefano Zampini } 534674ae819SStefano Zampini } 535674ae819SStefano Zampini } 536674ae819SStefano Zampini /* insert new data in subset array */ 537674ae819SStefano Zampini temp_buffer_size = 0; 538674ae819SStefano Zampini for (j=0;j<graph->nvtxs;j++) { 539674ae819SStefano Zampini if (graph->subset[j]==i+1) { 540674ae819SStefano Zampini graph->subset[j] = graph->n_subsets+add_to_subset[temp_buffer_size]; 541674ae819SStefano Zampini temp_buffer_size++; 542674ae819SStefano Zampini } 543674ae819SStefano Zampini } 544674ae819SStefano Zampini ierr = PetscFree(temp_buffer[0]);CHKERRQ(ierr); 545674ae819SStefano Zampini ierr = PetscFree(temp_buffer);CHKERRQ(ierr); 546674ae819SStefano Zampini ierr = PetscFree(add_to_subset);CHKERRQ(ierr); 547674ae819SStefano Zampini } 548674ae819SStefano Zampini } 549674ae819SStefano Zampini ierr = PetscFree(nodes_to_temp_buffer_indices);CHKERRQ(ierr); 550674ae819SStefano Zampini ierr = PetscFree(sizes_of_sends);CHKERRQ(ierr); 551674ae819SStefano Zampini ierr = PetscFree(send_requests);CHKERRQ(ierr); 552674ae819SStefano Zampini ierr = PetscFree(recv_requests);CHKERRQ(ierr); 553674ae819SStefano Zampini ierr = PetscFree(recv_buffer);CHKERRQ(ierr); 554674ae819SStefano Zampini ierr = PetscFree(recv_buffer_subset);CHKERRQ(ierr); 555674ae819SStefano Zampini ierr = PetscFree(send_buffer);CHKERRQ(ierr); 556674ae819SStefano Zampini ierr = PetscFree(cum_recv_counts);CHKERRQ(ierr); 557674ae819SStefano Zampini ierr = PetscFree(subset_to_nodes_indices);CHKERRQ(ierr); 558674ae819SStefano Zampini ierr = PetscFree(subset_cc_adapt);CHKERRQ(ierr); 559674ae819SStefano Zampini /* We are ready to find for connected components consistent among neighbouring subdomains */ 560674ae819SStefano Zampini if (global_subset_counter) { 561df48933bSStefano Zampini ierr = PetscBTMemzero(graph->nvtxs,graph->touched);CHKERRQ(ierr); 562674ae819SStefano Zampini global_subset_counter = 0; 563674ae819SStefano Zampini for (i=0;i<graph->nvtxs;i++) { 564df48933bSStefano Zampini if (graph->subset[i] && !PetscBTLookup(graph->touched,i)) { 565674ae819SStefano Zampini global_subset_counter++; 566674ae819SStefano Zampini for (j=i+1;j<graph->nvtxs;j++) { 567df48933bSStefano Zampini if (!PetscBTLookup(graph->touched,j) && graph->subset[j]==graph->subset[i]) { 568674ae819SStefano Zampini graph->subset[j] = global_subset_counter; 569df48933bSStefano Zampini ierr = PetscBTSet(graph->touched,j);CHKERRQ(ierr); 570674ae819SStefano Zampini } 571674ae819SStefano Zampini } 572674ae819SStefano Zampini graph->subset[i] = global_subset_counter; 573df48933bSStefano Zampini ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr); 574674ae819SStefano Zampini } 575674ae819SStefano Zampini } 576674ae819SStefano Zampini /* refine connected components locally */ 577674ae819SStefano Zampini ierr = PCBDDCGraphComputeConnectedComponentsLocal(graph);CHKERRQ(ierr); 578674ae819SStefano Zampini } 5795b1b9aeaSStefano Zampini /* restore original CSR graph of dofs */ 5805b1b9aeaSStefano Zampini ierr = PetscFree(new_xadj);CHKERRQ(ierr); 5815b1b9aeaSStefano Zampini ierr = PetscFree(new_adjncy);CHKERRQ(ierr); 5825b1b9aeaSStefano Zampini graph->xadj = old_xadj; 5835b1b9aeaSStefano Zampini graph->adjncy = old_adjncy; 584674ae819SStefano Zampini ierr = PetscFree(queue_in_global_numbering);CHKERRQ(ierr); 585674ae819SStefano Zampini } 586674ae819SStefano Zampini PetscFunctionReturn(0); 587674ae819SStefano Zampini } 588674ae819SStefano Zampini 589674ae819SStefano Zampini /* The following code has been adapted from function IsConnectedSubdomain contained 590674ae819SStefano Zampini in source file contig.c of METIS library (version 5.0.1) 591674ae819SStefano Zampini It finds connected components for each subset */ 592674ae819SStefano Zampini #undef __FUNCT__ 593674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphComputeConnectedComponentsLocal" 594674ae819SStefano Zampini PetscErrorCode PCBDDCGraphComputeConnectedComponentsLocal(PCBDDCGraph graph) 595674ae819SStefano Zampini { 596674ae819SStefano Zampini PetscInt i,j,k,first,last,nleft,ncc,pid,cum_queue,n,ncc_pid; 597674ae819SStefano Zampini PetscInt *queue_global; 5984a060362SStefano Zampini PetscMPIInt size; 599674ae819SStefano Zampini PetscErrorCode ierr; 600674ae819SStefano Zampini 601674ae819SStefano Zampini PetscFunctionBegin; 6024a060362SStefano Zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)graph->l2gmap),&size);CHKERRQ(ierr); 603674ae819SStefano Zampini /* quiet return if no csr info is available */ 604674ae819SStefano Zampini if (!graph->xadj || !graph->adjncy) { 605674ae819SStefano Zampini PetscFunctionReturn(0); 606674ae819SStefano Zampini } 607674ae819SStefano Zampini 608674ae819SStefano Zampini /* reset any previous search of connected components */ 609df48933bSStefano Zampini ierr = PetscBTMemzero(graph->nvtxs,graph->touched);CHKERRQ(ierr); 610674ae819SStefano Zampini graph->n_subsets = 0; 6114a060362SStefano Zampini if (size == 1) { 612dc36bf05SStefano Zampini if (graph->nvtxs) { 6134a060362SStefano Zampini graph->n_subsets = 1; 6144a060362SStefano Zampini for (i=0;i<graph->nvtxs;i++) { 6154a060362SStefano Zampini graph->subset[i] = 1; 6164a060362SStefano Zampini } 617dc36bf05SStefano Zampini } 6184a060362SStefano Zampini } else { 619674ae819SStefano Zampini for (i=0;i<graph->nvtxs;i++) { 6200cf82fd6SStefano Zampini if (graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK || !graph->count[i]) { 621df48933bSStefano Zampini ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr); 622674ae819SStefano Zampini graph->subset[i] = 0; 623674ae819SStefano Zampini } 624674ae819SStefano Zampini graph->n_subsets = PetscMax(graph->n_subsets,graph->subset[i]); 625674ae819SStefano Zampini } 6264a060362SStefano Zampini } 627674ae819SStefano Zampini ierr = PetscFree(graph->subset_ncc);CHKERRQ(ierr); 628785e854fSJed Brown ierr = PetscMalloc1(graph->n_subsets,&graph->subset_ncc);CHKERRQ(ierr); 629674ae819SStefano Zampini ierr = PetscMemzero(graph->subset_ncc,graph->n_subsets*sizeof(*graph->subset_ncc));CHKERRQ(ierr); 630674ae819SStefano Zampini ierr = PetscMemzero(graph->cptr,(graph->nvtxs+1)*sizeof(*graph->cptr));CHKERRQ(ierr); 631674ae819SStefano Zampini ierr = PetscMemzero(graph->queue,graph->nvtxs*sizeof(*graph->queue));CHKERRQ(ierr); 632674ae819SStefano Zampini 633674ae819SStefano Zampini /* begin search for connected components */ 634674ae819SStefano Zampini cum_queue = 0; 635674ae819SStefano Zampini ncc = 0; 636674ae819SStefano Zampini for (n=0;n<graph->n_subsets;n++) { 637674ae819SStefano Zampini pid = n+1; /* partition labeled by 0 is discarded */ 638674ae819SStefano Zampini nleft = 0; 639674ae819SStefano Zampini for (i=0;i<graph->nvtxs;i++) { 640674ae819SStefano Zampini if (graph->subset[i] == pid) { 641674ae819SStefano Zampini nleft++; 642674ae819SStefano Zampini } 643674ae819SStefano Zampini } 644674ae819SStefano Zampini for (i=0; i<graph->nvtxs; i++) { 645674ae819SStefano Zampini if (graph->subset[i] == pid) { 646674ae819SStefano Zampini break; 647674ae819SStefano Zampini } 648674ae819SStefano Zampini } 649df48933bSStefano Zampini ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr); 650674ae819SStefano Zampini graph->queue[cum_queue] = i; 651674ae819SStefano Zampini first = 0; 652674ae819SStefano Zampini last = 1; 653674ae819SStefano Zampini graph->cptr[ncc] = cum_queue; 654674ae819SStefano Zampini ncc_pid = 0; 655674ae819SStefano Zampini while (first != nleft) { 656674ae819SStefano Zampini if (first == last) { 657674ae819SStefano Zampini graph->cptr[++ncc] = first+cum_queue; 658674ae819SStefano Zampini ncc_pid++; 659df48933bSStefano Zampini for (i=0; i<graph->nvtxs; i++) { /* TODO-> use a while! */ 660df48933bSStefano Zampini if (graph->subset[i] == pid && !PetscBTLookup(graph->touched,i)) { 661674ae819SStefano Zampini break; 662674ae819SStefano Zampini } 663674ae819SStefano Zampini } 664674ae819SStefano Zampini graph->queue[cum_queue+last] = i; 665674ae819SStefano Zampini last++; 666df48933bSStefano Zampini ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr); 667674ae819SStefano Zampini } 668674ae819SStefano Zampini i = graph->queue[cum_queue+first]; 669674ae819SStefano Zampini first++; 670674ae819SStefano Zampini for (j=graph->xadj[i];j<graph->xadj[i+1];j++) { 671674ae819SStefano Zampini k = graph->adjncy[j]; 672df48933bSStefano Zampini if (graph->subset[k] == pid && !PetscBTLookup(graph->touched,k)) { 673674ae819SStefano Zampini graph->queue[cum_queue+last] = k; 674674ae819SStefano Zampini last++; 675df48933bSStefano Zampini ierr = PetscBTSet(graph->touched,k);CHKERRQ(ierr); 676674ae819SStefano Zampini } 677674ae819SStefano Zampini } 678674ae819SStefano Zampini } 679674ae819SStefano Zampini graph->cptr[++ncc] = first+cum_queue; 680674ae819SStefano Zampini ncc_pid++; 681674ae819SStefano Zampini cum_queue = graph->cptr[ncc]; 682674ae819SStefano Zampini graph->subset_ncc[n] = ncc_pid; 683674ae819SStefano Zampini } 684674ae819SStefano Zampini graph->ncc = ncc; 685674ae819SStefano Zampini /* For consistency among neighbours, I need to sort (by global ordering) each connected component */ 6864a060362SStefano Zampini if (size > 1) { 687785e854fSJed Brown ierr = PetscMalloc1(graph->cptr[graph->ncc],&queue_global);CHKERRQ(ierr); 688674ae819SStefano Zampini ierr = ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_global);CHKERRQ(ierr); 689674ae819SStefano Zampini for (i=0;i<graph->ncc;i++) { 690674ae819SStefano Zampini ierr = PetscSortIntWithArray(graph->cptr[i+1]-graph->cptr[i],&queue_global[graph->cptr[i]],&graph->queue[graph->cptr[i]]);CHKERRQ(ierr); 691674ae819SStefano Zampini } 692674ae819SStefano Zampini ierr = PetscFree(queue_global);CHKERRQ(ierr); 6934a060362SStefano Zampini } 694674ae819SStefano Zampini PetscFunctionReturn(0); 695674ae819SStefano Zampini } 696674ae819SStefano Zampini 697674ae819SStefano Zampini #undef __FUNCT__ 698674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphSetUp" 699674ae819SStefano Zampini PetscErrorCode PCBDDCGraphSetUp(PCBDDCGraph graph, PetscInt custom_minimal_size, IS neumann_is, IS dirichlet_is, PetscInt n_ISForDofs, IS ISForDofs[], IS custom_primal_vertices) 700674ae819SStefano Zampini { 701674ae819SStefano Zampini VecScatter scatter_ctx; 702674ae819SStefano Zampini Vec local_vec,local_vec2,global_vec; 703674ae819SStefano Zampini IS to,from; 704674ae819SStefano Zampini MPI_Comm comm; 705674ae819SStefano Zampini PetscScalar *array,*array2; 706674ae819SStefano Zampini const PetscInt *is_indices; 707f0321c90SStefano Zampini PetscInt n_neigh,*neigh,*n_shared,**shared,*queue_global,*subset_ref_node_global; 708674ae819SStefano Zampini PetscInt i,j,k,s,total_counts,nodes_touched,is_size; 709674ae819SStefano Zampini PetscErrorCode ierr; 71051b0f6cfSStefano Zampini PetscBool same_set,mirrors_found; 711674ae819SStefano Zampini 712674ae819SStefano Zampini PetscFunctionBegin; 713674ae819SStefano Zampini ierr = PetscObjectGetComm((PetscObject)(graph->l2gmap),&comm);CHKERRQ(ierr); 714674ae819SStefano Zampini /* custom_minimal_size */ 715674ae819SStefano Zampini graph->custom_minimal_size = PetscMax(graph->custom_minimal_size,custom_minimal_size); 716674ae819SStefano Zampini /* get info l2gmap and allocate work vectors */ 717674ae819SStefano Zampini ierr = ISLocalToGlobalMappingGetInfo(graph->l2gmap,&n_neigh,&neigh,&n_shared,&shared);CHKERRQ(ierr); 718674ae819SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(graph->l2gmap,&is_indices);CHKERRQ(ierr); 719674ae819SStefano Zampini j = 0; 720674ae819SStefano Zampini for (i=0;i<graph->nvtxs;i++) { 721674ae819SStefano Zampini j = PetscMax(j,is_indices[i]); 722674ae819SStefano Zampini } 72325a755b6SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(graph->l2gmap,&is_indices);CHKERRQ(ierr); 724674ae819SStefano Zampini ierr = MPI_Allreduce(&j,&i,1,MPIU_INT,MPI_MAX,comm);CHKERRQ(ierr); 725674ae819SStefano Zampini i++; 726674ae819SStefano Zampini ierr = VecCreate(PETSC_COMM_SELF,&local_vec);CHKERRQ(ierr); 727674ae819SStefano Zampini ierr = VecSetSizes(local_vec,PETSC_DECIDE,graph->nvtxs);CHKERRQ(ierr); 728674ae819SStefano Zampini ierr = VecSetType(local_vec,VECSTANDARD);CHKERRQ(ierr); 729674ae819SStefano Zampini ierr = VecDuplicate(local_vec,&local_vec2);CHKERRQ(ierr); 730674ae819SStefano Zampini ierr = VecCreate(comm,&global_vec);CHKERRQ(ierr); 731674ae819SStefano Zampini ierr = VecSetSizes(global_vec,PETSC_DECIDE,i);CHKERRQ(ierr); 732674ae819SStefano Zampini ierr = VecSetType(global_vec,VECSTANDARD);CHKERRQ(ierr); 733674ae819SStefano Zampini ierr = ISCreateStride(PETSC_COMM_SELF,graph->nvtxs,0,1,&to);CHKERRQ(ierr); 734674ae819SStefano Zampini ierr = ISLocalToGlobalMappingApplyIS(graph->l2gmap,to,&from);CHKERRQ(ierr); 735674ae819SStefano Zampini ierr = VecScatterCreate(global_vec,from,local_vec,to,&scatter_ctx);CHKERRQ(ierr); 73651b0f6cfSStefano Zampini 73751b0f6cfSStefano Zampini /* get local periodic nodes */ 73851b0f6cfSStefano Zampini mirrors_found = PETSC_FALSE; 7399881197aSStefano Zampini if (graph->nvtxs && n_neigh) { 74049eeff8cSStefano Zampini for (i=0; i<n_shared[0]; i++) graph->count[shared[0][i]] += 1; 74151b0f6cfSStefano Zampini for (i=0; i<n_shared[0]; i++) { 74251b0f6cfSStefano Zampini if (graph->count[shared[0][i]] > 1) { 74351b0f6cfSStefano Zampini mirrors_found = PETSC_TRUE; 74451b0f6cfSStefano Zampini break; 74551b0f6cfSStefano Zampini } 74651b0f6cfSStefano Zampini } 74749eeff8cSStefano Zampini } 74851b0f6cfSStefano Zampini /* compute local mirrors (if any) */ 74951b0f6cfSStefano Zampini if (mirrors_found) { 75051b0f6cfSStefano Zampini PetscInt *local_indices,*global_indices; 75151b0f6cfSStefano Zampini /* get arrays of local and global indices */ 752785e854fSJed Brown ierr = PetscMalloc1(graph->nvtxs,&local_indices);CHKERRQ(ierr); 75351b0f6cfSStefano Zampini ierr = ISGetIndices(to,(const PetscInt**)&is_indices);CHKERRQ(ierr); 75451b0f6cfSStefano Zampini ierr = PetscMemcpy(local_indices,is_indices,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr); 75551b0f6cfSStefano Zampini ierr = ISRestoreIndices(to,(const PetscInt**)&is_indices);CHKERRQ(ierr); 756785e854fSJed Brown ierr = PetscMalloc1(graph->nvtxs,&global_indices);CHKERRQ(ierr); 75751b0f6cfSStefano Zampini ierr = ISGetIndices(from,(const PetscInt**)&is_indices);CHKERRQ(ierr); 75851b0f6cfSStefano Zampini ierr = PetscMemcpy(global_indices,is_indices,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr); 75951b0f6cfSStefano Zampini ierr = ISRestoreIndices(from,(const PetscInt**)&is_indices);CHKERRQ(ierr); 76051b0f6cfSStefano Zampini /* allocate space for mirrors */ 7618e5aaad5SJed Brown ierr = PetscMalloc2(graph->nvtxs,&graph->mirrors,graph->nvtxs,&graph->mirrors_set);CHKERRQ(ierr); 76251b0f6cfSStefano Zampini ierr = PetscMemzero(graph->mirrors,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr); 76351b0f6cfSStefano Zampini graph->mirrors_set[0] = 0; 76451b0f6cfSStefano Zampini 76551b0f6cfSStefano Zampini k=0; 76651b0f6cfSStefano Zampini for (i=0;i<n_shared[0];i++) { 76751b0f6cfSStefano Zampini j=shared[0][i]; 76851b0f6cfSStefano Zampini if (graph->count[j] > 1) { 76951b0f6cfSStefano Zampini graph->mirrors[j]++; 77051b0f6cfSStefano Zampini k++; 77151b0f6cfSStefano Zampini } 77251b0f6cfSStefano Zampini } 77351b0f6cfSStefano Zampini /* allocate space for set of mirrors */ 774785e854fSJed Brown ierr = PetscMalloc1(k,&graph->mirrors_set[0]);CHKERRQ(ierr); 77551b0f6cfSStefano Zampini for (i=1;i<graph->nvtxs;i++) 77651b0f6cfSStefano Zampini graph->mirrors_set[i]=graph->mirrors_set[i-1]+graph->mirrors[i-1]; 77751b0f6cfSStefano Zampini 77851b0f6cfSStefano Zampini /* fill arrays */ 77951b0f6cfSStefano Zampini ierr = PetscMemzero(graph->mirrors,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr); 78051b0f6cfSStefano Zampini for (j=0;j<n_shared[0];j++) { 78151b0f6cfSStefano Zampini i=shared[0][j]; 78251b0f6cfSStefano Zampini if (graph->count[i] > 1) 78351b0f6cfSStefano Zampini graph->mirrors_set[i][graph->mirrors[i]++]=global_indices[i]; 78451b0f6cfSStefano Zampini } 78551b0f6cfSStefano Zampini ierr = PetscSortIntWithArray(graph->nvtxs,global_indices,local_indices);CHKERRQ(ierr); 78651b0f6cfSStefano Zampini for (i=0;i<graph->nvtxs;i++) { 78751b0f6cfSStefano Zampini if (graph->mirrors[i] > 0) { 78851b0f6cfSStefano Zampini ierr = PetscFindInt(graph->mirrors_set[i][0],graph->nvtxs,global_indices,&k);CHKERRQ(ierr); 78951b0f6cfSStefano Zampini j = global_indices[k]; 79051b0f6cfSStefano Zampini while ( k > 0 && global_indices[k-1] == j) k--; 79151b0f6cfSStefano Zampini for (j=0;j<graph->mirrors[i];j++) { 79251b0f6cfSStefano Zampini graph->mirrors_set[i][j]=local_indices[k+j]; 79351b0f6cfSStefano Zampini } 79451b0f6cfSStefano Zampini ierr = PetscSortInt(graph->mirrors[i],graph->mirrors_set[i]);CHKERRQ(ierr); 79551b0f6cfSStefano Zampini } 79651b0f6cfSStefano Zampini } 79751b0f6cfSStefano Zampini ierr = PetscFree(local_indices);CHKERRQ(ierr); 79851b0f6cfSStefano Zampini ierr = PetscFree(global_indices);CHKERRQ(ierr); 79951b0f6cfSStefano Zampini } 80051b0f6cfSStefano Zampini ierr = PetscMemzero(graph->count,graph->nvtxs*sizeof(*graph->count));CHKERRQ(ierr); 801674ae819SStefano Zampini ierr = ISDestroy(&to);CHKERRQ(ierr); 802674ae819SStefano Zampini ierr = ISDestroy(&from);CHKERRQ(ierr); 80351b0f6cfSStefano Zampini 804674ae819SStefano Zampini /* Count total number of neigh per node */ 805674ae819SStefano Zampini k=0; 806674ae819SStefano Zampini for (i=1;i<n_neigh;i++) { 807674ae819SStefano Zampini k += n_shared[i]; 808674ae819SStefano Zampini for (j=0;j<n_shared[i];j++) { 809674ae819SStefano Zampini graph->count[shared[i][j]] += 1; 810674ae819SStefano Zampini } 811674ae819SStefano Zampini } 812674ae819SStefano Zampini /* Allocate space for storing the set of neighbours for each node */ 813674ae819SStefano Zampini if (graph->nvtxs) { 814785e854fSJed Brown ierr = PetscMalloc1(k,&graph->neighbours_set[0]);CHKERRQ(ierr); 815674ae819SStefano Zampini } 816674ae819SStefano Zampini for (i=1;i<graph->nvtxs;i++) { /* dont count myself */ 817674ae819SStefano Zampini graph->neighbours_set[i]=graph->neighbours_set[i-1]+graph->count[i-1]; 818674ae819SStefano Zampini } 819674ae819SStefano Zampini /* Get information for sharing subdomains */ 820674ae819SStefano Zampini ierr = PetscMemzero(graph->count,graph->nvtxs*sizeof(*graph->count));CHKERRQ(ierr); 821674ae819SStefano Zampini for (i=1;i<n_neigh;i++) { /* dont count myself */ 822674ae819SStefano Zampini s = n_shared[i]; 823674ae819SStefano Zampini for (j=0;j<s;j++) { 824674ae819SStefano Zampini k = shared[i][j]; 825674ae819SStefano Zampini graph->neighbours_set[k][graph->count[k]] = neigh[i]; 826674ae819SStefano Zampini graph->count[k] += 1; 827674ae819SStefano Zampini } 828674ae819SStefano Zampini } 829674ae819SStefano Zampini /* sort set of sharing subdomains */ 830674ae819SStefano Zampini for (i=0;i<graph->nvtxs;i++) { 83193fb5fd3SStefano Zampini ierr = PetscSortRemoveDupsInt(&graph->count[i],graph->neighbours_set[i]);CHKERRQ(ierr); 832674ae819SStefano Zampini } 83367c9da69SStefano Zampini /* 83467c9da69SStefano Zampini Get info for dofs splitting 83567c9da69SStefano Zampini User can specify only a subset; an additional field is considered as a complementary set 83667c9da69SStefano Zampini */ 83767c9da69SStefano Zampini for (i=0;i<graph->nvtxs;i++) { 83867c9da69SStefano Zampini graph->which_dof[i] = n_ISForDofs; /* by default a dof belongs to the complement set */ 83967c9da69SStefano Zampini } 840674ae819SStefano Zampini for (i=0;i<n_ISForDofs;i++) { 84163602bcaSStefano Zampini ierr = ISGetLocalSize(ISForDofs[i],&is_size);CHKERRQ(ierr); 842674ae819SStefano Zampini ierr = ISGetIndices(ISForDofs[i],(const PetscInt**)&is_indices);CHKERRQ(ierr); 843674ae819SStefano Zampini for (j=0;j<is_size;j++) { 844d50ed60aSStefano Zampini if (is_indices[j] > -1 && is_indices[j] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */ 845674ae819SStefano Zampini graph->which_dof[is_indices[j]] = i; 846674ae819SStefano Zampini } 84767c9da69SStefano Zampini } 848674ae819SStefano Zampini ierr = ISRestoreIndices(ISForDofs[i],(const PetscInt**)&is_indices);CHKERRQ(ierr); 849674ae819SStefano Zampini } 850674ae819SStefano Zampini /* Take into account Neumann nodes */ 851674ae819SStefano Zampini ierr = VecSet(local_vec,0.0);CHKERRQ(ierr); 852674ae819SStefano Zampini ierr = VecSet(local_vec2,0.0);CHKERRQ(ierr); 853674ae819SStefano Zampini if (neumann_is) { 854674ae819SStefano Zampini ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr); 85582ba6b80SStefano Zampini ierr = ISGetLocalSize(neumann_is,&is_size);CHKERRQ(ierr); 856674ae819SStefano Zampini ierr = ISGetIndices(neumann_is,(const PetscInt**)&is_indices);CHKERRQ(ierr); 857674ae819SStefano Zampini for (i=0;i<is_size;i++) { 858d50ed60aSStefano Zampini if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */ 859674ae819SStefano Zampini array[is_indices[i]] = 1.0; 860674ae819SStefano Zampini } 8613c629590SStefano Zampini } 862674ae819SStefano Zampini ierr = ISRestoreIndices(neumann_is,(const PetscInt**)&is_indices);CHKERRQ(ierr); 863674ae819SStefano Zampini ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr); 864674ae819SStefano Zampini } 865674ae819SStefano Zampini /* Neumann nodes: impose consistency among neighbours */ 866674ae819SStefano Zampini ierr = VecSet(global_vec,0.0);CHKERRQ(ierr); 867674ae819SStefano Zampini ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 868674ae819SStefano Zampini ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 869674ae819SStefano Zampini ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 870674ae819SStefano Zampini ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 871674ae819SStefano Zampini ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr); 872674ae819SStefano Zampini for (i=0;i<graph->nvtxs;i++) { 8733c629590SStefano Zampini if (PetscRealPart(array[i]) > 0.1) { 8740cf82fd6SStefano Zampini graph->special_dof[i] = PCBDDCGRAPH_NEUMANN_MARK; 875674ae819SStefano Zampini } 876674ae819SStefano Zampini } 877674ae819SStefano Zampini ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr); 878674ae819SStefano Zampini /* Take into account Dirichlet nodes */ 879674ae819SStefano Zampini ierr = VecSet(local_vec2,0.0);CHKERRQ(ierr); 880674ae819SStefano Zampini if (dirichlet_is) { 881674ae819SStefano Zampini ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr); 882674ae819SStefano Zampini ierr = VecGetArray(local_vec2,&array2);CHKERRQ(ierr); 88382ba6b80SStefano Zampini ierr = ISGetLocalSize(dirichlet_is,&is_size);CHKERRQ(ierr); 884674ae819SStefano Zampini ierr = ISGetIndices(dirichlet_is,(const PetscInt**)&is_indices);CHKERRQ(ierr); 885674ae819SStefano Zampini for (i=0;i<is_size;i++){ 886d50ed60aSStefano Zampini if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */ 887674ae819SStefano Zampini k = is_indices[i]; 888df48933bSStefano Zampini if (graph->count[k] && !PetscBTLookup(graph->touched,k)) { 8893c629590SStefano Zampini if (PetscRealPart(array[k]) > 0.1) { 8903c629590SStefano 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); 891674ae819SStefano Zampini } 892674ae819SStefano Zampini array2[k] = 1.0; 893674ae819SStefano Zampini } 894674ae819SStefano Zampini } 8953c629590SStefano Zampini } 896674ae819SStefano Zampini ierr = ISRestoreIndices(dirichlet_is,(const PetscInt**)&is_indices);CHKERRQ(ierr); 897674ae819SStefano Zampini ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr); 898674ae819SStefano Zampini ierr = VecRestoreArray(local_vec2,&array2);CHKERRQ(ierr); 899674ae819SStefano Zampini } 900674ae819SStefano Zampini /* Dirichlet nodes: impose consistency among neighbours */ 901674ae819SStefano Zampini ierr = VecSet(global_vec,0.0);CHKERRQ(ierr); 902674ae819SStefano Zampini ierr = VecScatterBegin(scatter_ctx,local_vec2,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 903674ae819SStefano Zampini ierr = VecScatterEnd(scatter_ctx,local_vec2,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 904674ae819SStefano Zampini ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 905674ae819SStefano Zampini ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 906674ae819SStefano Zampini ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr); 907674ae819SStefano Zampini for (i=0;i<graph->nvtxs;i++) { 9083c629590SStefano Zampini if (PetscRealPart(array[i]) > 0.1) { 909df48933bSStefano Zampini ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr); 910674ae819SStefano Zampini graph->subset[i] = 0; /* dirichlet nodes treated as internal -> is it ok? */ 9110cf82fd6SStefano Zampini graph->special_dof[i] = PCBDDCGRAPH_DIRICHLET_MARK; 912674ae819SStefano Zampini } 913674ae819SStefano Zampini } 914674ae819SStefano Zampini ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr); 91551b0f6cfSStefano Zampini 91608a5cf49SStefano Zampini /* mark local periodic nodes (if any) and adapt CSR graph (if any) */ 91751b0f6cfSStefano Zampini if (graph->mirrors) { 91851b0f6cfSStefano Zampini for (i=0;i<graph->nvtxs;i++) 91951b0f6cfSStefano Zampini if (graph->mirrors[i]) 9200cf82fd6SStefano Zampini graph->special_dof[i] = PCBDDCGRAPH_LOCAL_PERIODIC_MARK; 92151b0f6cfSStefano Zampini 92208a5cf49SStefano Zampini if (graph->xadj && graph->adjncy) { 92308a5cf49SStefano Zampini PetscInt *new_xadj,*new_adjncy; 92451b0f6cfSStefano Zampini /* sort CSR graph */ 92551b0f6cfSStefano Zampini for (i=0;i<graph->nvtxs;i++) 92651b0f6cfSStefano Zampini ierr = PetscSortInt(graph->xadj[i+1]-graph->xadj[i],&graph->adjncy[graph->xadj[i]]);CHKERRQ(ierr); 92751b0f6cfSStefano Zampini 92851b0f6cfSStefano Zampini /* adapt local CSR graph in case of local periodicity */ 92951b0f6cfSStefano Zampini k=0; 93051b0f6cfSStefano Zampini for (i=0;i<graph->nvtxs;i++) 93151b0f6cfSStefano Zampini for (j=graph->xadj[i];j<graph->xadj[i+1];j++) 93251b0f6cfSStefano Zampini k += graph->mirrors[graph->adjncy[j]]; 93351b0f6cfSStefano Zampini 934854ce69bSBarry Smith ierr = PetscMalloc1(graph->nvtxs+1,&new_xadj);CHKERRQ(ierr); 935854ce69bSBarry Smith ierr = PetscMalloc1(k+graph->xadj[graph->nvtxs],&new_adjncy);CHKERRQ(ierr); 93651b0f6cfSStefano Zampini new_xadj[0]=0; 93751b0f6cfSStefano Zampini for (i=0;i<graph->nvtxs;i++) { 93851b0f6cfSStefano Zampini k = graph->xadj[i+1]-graph->xadj[i]; 93951b0f6cfSStefano Zampini ierr = PetscMemcpy(&new_adjncy[new_xadj[i]],&graph->adjncy[graph->xadj[i]],k*sizeof(PetscInt));CHKERRQ(ierr); 94051b0f6cfSStefano Zampini new_xadj[i+1]=new_xadj[i]+k; 94151b0f6cfSStefano Zampini for (j=graph->xadj[i];j<graph->xadj[i+1];j++) { 94251b0f6cfSStefano Zampini k = graph->mirrors[graph->adjncy[j]]; 94351b0f6cfSStefano Zampini ierr = PetscMemcpy(&new_adjncy[new_xadj[i+1]],graph->mirrors_set[graph->adjncy[j]],k*sizeof(PetscInt));CHKERRQ(ierr); 94451b0f6cfSStefano Zampini new_xadj[i+1]+=k; 94551b0f6cfSStefano Zampini } 94651b0f6cfSStefano Zampini k = new_xadj[i+1]-new_xadj[i]; 94751b0f6cfSStefano Zampini ierr = PetscSortRemoveDupsInt(&k,&new_adjncy[new_xadj[i]]);CHKERRQ(ierr); 94851b0f6cfSStefano Zampini new_xadj[i+1]=new_xadj[i]+k; 94951b0f6cfSStefano Zampini } 95051b0f6cfSStefano Zampini /* set new CSR into graph */ 95151b0f6cfSStefano Zampini ierr = PetscFree(graph->xadj);CHKERRQ(ierr); 95251b0f6cfSStefano Zampini ierr = PetscFree(graph->adjncy);CHKERRQ(ierr); 95351b0f6cfSStefano Zampini graph->xadj = new_xadj; 95451b0f6cfSStefano Zampini graph->adjncy = new_adjncy; 95551b0f6cfSStefano Zampini } 95608a5cf49SStefano Zampini } 95751b0f6cfSStefano Zampini 958674ae819SStefano Zampini /* mark special nodes -> each will become a single node equivalence class */ 9599b70a373SStefano Zampini ierr = VecSet(local_vec,0.0);CHKERRQ(ierr); 960674ae819SStefano Zampini if (custom_primal_vertices) { 9619b70a373SStefano Zampini ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr); 96263602bcaSStefano Zampini ierr = ISGetLocalSize(custom_primal_vertices,&is_size);CHKERRQ(ierr); 963674ae819SStefano Zampini ierr = ISGetIndices(custom_primal_vertices,(const PetscInt**)&is_indices);CHKERRQ(ierr); 964674ae819SStefano Zampini for (i=0;i<is_size;i++){ 9659b70a373SStefano Zampini if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */ 9669b70a373SStefano Zampini array[is_indices[i]] = 1.0; 9679b70a373SStefano Zampini } 968674ae819SStefano Zampini } 969674ae819SStefano Zampini ierr = ISRestoreIndices(custom_primal_vertices,(const PetscInt**)&is_indices);CHKERRQ(ierr); 9709b70a373SStefano Zampini ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr); 971674ae819SStefano Zampini } 9729b70a373SStefano Zampini /* special nodes: impose consistency among neighbours */ 9739b70a373SStefano Zampini ierr = VecSet(global_vec,0.0);CHKERRQ(ierr); 9749b70a373SStefano Zampini ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 9759b70a373SStefano Zampini ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 9769b70a373SStefano Zampini ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9779b70a373SStefano Zampini ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9789b70a373SStefano Zampini ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr); 9799b70a373SStefano Zampini j = 0; 9809b70a373SStefano Zampini for (i=0;i<graph->nvtxs;i++) { 9819b70a373SStefano Zampini if (PetscRealPart(array[i]) > 0.1 && graph->special_dof[i] != PCBDDCGRAPH_DIRICHLET_MARK) { 9829b70a373SStefano Zampini graph->special_dof[i] = PCBDDCGRAPH_SPECIAL_MARK-j; 9839b70a373SStefano Zampini j++; 9849b70a373SStefano Zampini } 9859b70a373SStefano Zampini } 9869b70a373SStefano Zampini ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr); 9879b70a373SStefano Zampini 988674ae819SStefano Zampini /* mark interior nodes as touched and belonging to partition number 0 */ 989674ae819SStefano Zampini for (i=0;i<graph->nvtxs;i++) { 990674ae819SStefano Zampini if (!graph->count[i]) { 991df48933bSStefano Zampini ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr); 992674ae819SStefano Zampini graph->subset[i] = 0; 993674ae819SStefano Zampini } 994674ae819SStefano Zampini } 995674ae819SStefano Zampini /* init graph structure and compute default subsets */ 996674ae819SStefano Zampini nodes_touched=0; 997674ae819SStefano Zampini for (i=0;i<graph->nvtxs;i++) { 998df48933bSStefano Zampini if (PetscBTLookup(graph->touched,i)) { 999674ae819SStefano Zampini nodes_touched++; 1000674ae819SStefano Zampini } 1001674ae819SStefano Zampini } 1002674ae819SStefano Zampini i = 0; 1003674ae819SStefano Zampini graph->ncc = 0; 1004674ae819SStefano Zampini total_counts = 0; 1005674ae819SStefano Zampini while (nodes_touched<graph->nvtxs) { 1006674ae819SStefano Zampini /* find first untouched node in local ordering */ 1007df48933bSStefano Zampini while (PetscBTLookup(graph->touched,i)) { 1008674ae819SStefano Zampini i++; 1009674ae819SStefano Zampini } 1010df48933bSStefano Zampini ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr); 1011674ae819SStefano Zampini graph->subset[i] = graph->ncc+1; 1012674ae819SStefano Zampini graph->cptr[graph->ncc] = total_counts; 1013674ae819SStefano Zampini graph->queue[total_counts] = i; 1014674ae819SStefano Zampini total_counts++; 1015674ae819SStefano Zampini nodes_touched++; 1016674ae819SStefano Zampini /* now find all other nodes having the same set of sharing subdomains */ 1017674ae819SStefano Zampini for (j=i+1;j<graph->nvtxs;j++) { 101874e413f5SStefano Zampini /* check for same number of sharing subdomains, dof number and same special mark */ 1019df48933bSStefano 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]) { 1020674ae819SStefano Zampini /* check for same set of sharing subdomains */ 1021674ae819SStefano Zampini same_set=PETSC_TRUE; 1022674ae819SStefano Zampini for (k=0;k<graph->count[j];k++){ 1023674ae819SStefano Zampini if (graph->neighbours_set[i][k]!=graph->neighbours_set[j][k]) { 1024674ae819SStefano Zampini same_set=PETSC_FALSE; 1025674ae819SStefano Zampini } 1026674ae819SStefano Zampini } 1027674ae819SStefano Zampini /* I found a friend of mine */ 1028674ae819SStefano Zampini if (same_set) { 1029674ae819SStefano Zampini graph->subset[j]=graph->ncc+1; 1030df48933bSStefano Zampini ierr = PetscBTSet(graph->touched,j);CHKERRQ(ierr); 1031674ae819SStefano Zampini nodes_touched++; 1032674ae819SStefano Zampini graph->queue[total_counts] = j; 1033674ae819SStefano Zampini total_counts++; 1034674ae819SStefano Zampini } 1035674ae819SStefano Zampini } 1036674ae819SStefano Zampini } 1037674ae819SStefano Zampini graph->ncc++; 1038674ae819SStefano Zampini } 1039674ae819SStefano Zampini /* set default number of subsets (at this point no info on csr graph has been taken into account, so n_subsets = ncc */ 1040674ae819SStefano Zampini graph->n_subsets = graph->ncc; 1041785e854fSJed Brown ierr = PetscMalloc1(graph->n_subsets,&graph->subset_ncc);CHKERRQ(ierr); 1042674ae819SStefano Zampini for (i=0;i<graph->n_subsets;i++) { 1043674ae819SStefano Zampini graph->subset_ncc[i] = 1; 1044674ae819SStefano Zampini } 1045674ae819SStefano Zampini /* final pointer */ 1046674ae819SStefano Zampini graph->cptr[graph->ncc] = total_counts; 1047674ae819SStefano Zampini /* free memory allocated by ISLocalToGlobalMappingGetInfo */ 1048674ae819SStefano Zampini ierr = ISLocalToGlobalMappingRestoreInfo(graph->l2gmap,&n_neigh,&neigh,&n_shared,&shared);CHKERRQ(ierr); 10493a5b03d1SStefano Zampini /* get a reference node (min index in global ordering) for each subset */ 1050f0321c90SStefano Zampini ierr = PetscMalloc1(graph->ncc,&subset_ref_node_global);CHKERRQ(ierr); 10513a5b03d1SStefano Zampini ierr = PetscMalloc1(graph->cptr[graph->ncc],&queue_global);CHKERRQ(ierr); 10523a5b03d1SStefano Zampini ierr = ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_global);CHKERRQ(ierr); 10533a5b03d1SStefano Zampini for (i=0;i<graph->ncc;i++) { 10543a5b03d1SStefano Zampini PetscInt minval = queue_global[graph->cptr[i]]; 1055f0321c90SStefano Zampini PetscInt minloc = graph->queue[graph->cptr[i]]; 10563a5b03d1SStefano Zampini for (j=graph->cptr[i]+1;j<graph->cptr[i+1];j++) { 1057f0321c90SStefano Zampini if (minval > queue_global[j]) { 1058f0321c90SStefano Zampini minval = queue_global[j]; 1059f0321c90SStefano Zampini minloc = graph->queue[j]; 10603a5b03d1SStefano Zampini } 10613a5b03d1SStefano Zampini } 1062f0321c90SStefano Zampini subset_ref_node_global[i] = minloc; 10634a060362SStefano Zampini /* For consistency among neighbours, I need to sort (by global ordering) each connected component */ 10644a060362SStefano Zampini ierr = PetscSortIntWithArray(graph->cptr[i+1]-graph->cptr[i],&queue_global[graph->cptr[i]],&graph->queue[graph->cptr[i]]);CHKERRQ(ierr); 1065f0321c90SStefano Zampini } 1066f0321c90SStefano Zampini /* renumber reference nodes */ 1067f0321c90SStefano Zampini ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)(graph->l2gmap)),graph->l2gmap,graph->ncc,subset_ref_node_global,NULL,&k,&graph->subset_ref_node);CHKERRQ(ierr); 1068f0321c90SStefano Zampini ierr = PetscFree(subset_ref_node_global);CHKERRQ(ierr); 10693a5b03d1SStefano Zampini ierr = PetscFree(queue_global);CHKERRQ(ierr); 1070674ae819SStefano Zampini /* free objects */ 1071674ae819SStefano Zampini ierr = VecDestroy(&local_vec);CHKERRQ(ierr); 1072674ae819SStefano Zampini ierr = VecDestroy(&local_vec2);CHKERRQ(ierr); 1073674ae819SStefano Zampini ierr = VecDestroy(&global_vec);CHKERRQ(ierr); 1074674ae819SStefano Zampini ierr = VecScatterDestroy(&scatter_ctx);CHKERRQ(ierr); 1075674ae819SStefano Zampini PetscFunctionReturn(0); 1076674ae819SStefano Zampini } 1077674ae819SStefano Zampini 1078674ae819SStefano Zampini #undef __FUNCT__ 1079674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphResetCSR" 1080674ae819SStefano Zampini PetscErrorCode PCBDDCGraphResetCSR(PCBDDCGraph graph) 1081674ae819SStefano Zampini { 1082674ae819SStefano Zampini PetscErrorCode ierr; 1083674ae819SStefano Zampini 1084674ae819SStefano Zampini PetscFunctionBegin; 1085674ae819SStefano Zampini ierr = PetscFree(graph->xadj);CHKERRQ(ierr); 1086674ae819SStefano Zampini ierr = PetscFree(graph->adjncy);CHKERRQ(ierr); 1087575ad6abSStefano Zampini graph->nvtxs_csr = 0; 1088674ae819SStefano Zampini PetscFunctionReturn(0); 1089674ae819SStefano Zampini } 1090674ae819SStefano Zampini 1091674ae819SStefano Zampini #undef __FUNCT__ 1092674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphReset" 1093674ae819SStefano Zampini PetscErrorCode PCBDDCGraphReset(PCBDDCGraph graph) 1094674ae819SStefano Zampini { 1095674ae819SStefano Zampini PetscErrorCode ierr; 1096674ae819SStefano Zampini 1097674ae819SStefano Zampini PetscFunctionBegin; 1098674ae819SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&graph->l2gmap);CHKERRQ(ierr); 1099674ae819SStefano Zampini ierr = PetscFree(graph->subset_ncc);CHKERRQ(ierr); 11003a5b03d1SStefano Zampini ierr = PetscFree(graph->subset_ref_node);CHKERRQ(ierr); 1101674ae819SStefano Zampini if (graph->nvtxs) { 1102674ae819SStefano Zampini ierr = PetscFree(graph->neighbours_set[0]);CHKERRQ(ierr); 1103674ae819SStefano Zampini } 1104df48933bSStefano Zampini ierr = PetscBTDestroy(&graph->touched);CHKERRQ(ierr); 1105df48933bSStefano Zampini ierr = PetscFree7(graph->count, 1106674ae819SStefano Zampini graph->neighbours_set, 1107674ae819SStefano Zampini graph->subset, 1108674ae819SStefano Zampini graph->which_dof, 1109674ae819SStefano Zampini graph->cptr, 1110df48933bSStefano Zampini graph->queue, 1111df48933bSStefano Zampini graph->special_dof);CHKERRQ(ierr); 111251b0f6cfSStefano Zampini if (graph->mirrors) { 111351b0f6cfSStefano Zampini ierr = PetscFree(graph->mirrors_set[0]);CHKERRQ(ierr); 111451b0f6cfSStefano Zampini } 111551b0f6cfSStefano Zampini ierr = PetscFree2(graph->mirrors,graph->mirrors_set);CHKERRQ(ierr); 1116674ae819SStefano Zampini graph->nvtxs = 0; 1117674ae819SStefano Zampini graph->n_subsets = 0; 1118674ae819SStefano Zampini graph->custom_minimal_size = 1; 1119674ae819SStefano Zampini PetscFunctionReturn(0); 1120674ae819SStefano Zampini } 1121674ae819SStefano Zampini 1122674ae819SStefano Zampini #undef __FUNCT__ 1123674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphInit" 1124674ae819SStefano Zampini PetscErrorCode PCBDDCGraphInit(PCBDDCGraph graph, ISLocalToGlobalMapping l2gmap) 1125674ae819SStefano Zampini { 1126674ae819SStefano Zampini PetscInt n; 1127674ae819SStefano Zampini PetscErrorCode ierr; 1128674ae819SStefano Zampini 1129674ae819SStefano Zampini PetscFunctionBegin; 1130674ae819SStefano Zampini PetscValidPointer(graph,1); 1131674ae819SStefano Zampini PetscValidHeaderSpecific(l2gmap,IS_LTOGM_CLASSID,2); 11328e61c736SStefano Zampini /* raise an error if already allocated */ 1133674ae819SStefano Zampini if (graph->nvtxs) { 11348e61c736SStefano Zampini SETERRQ(PetscObjectComm((PetscObject)l2gmap),PETSC_ERR_PLIB,"BDDCGraph already initialized"); 1135674ae819SStefano Zampini } 1136674ae819SStefano Zampini /* set number of vertices */ 1137674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)l2gmap);CHKERRQ(ierr); 1138674ae819SStefano Zampini graph->l2gmap = l2gmap; 1139674ae819SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(l2gmap,&n);CHKERRQ(ierr); 1140674ae819SStefano Zampini graph->nvtxs = n; 1141674ae819SStefano Zampini /* allocate used space */ 1142df48933bSStefano Zampini ierr = PetscBTCreate(graph->nvtxs,&graph->touched);CHKERRQ(ierr); 11438e5aaad5SJed Brown ierr = PetscMalloc7(graph->nvtxs,&graph->count, 11448e5aaad5SJed Brown graph->nvtxs,&graph->neighbours_set, 11458e5aaad5SJed Brown graph->nvtxs,&graph->subset, 11468e5aaad5SJed Brown graph->nvtxs,&graph->which_dof, 11478e5aaad5SJed Brown graph->nvtxs+1,&graph->cptr, 11488e5aaad5SJed Brown graph->nvtxs,&graph->queue, 11498e5aaad5SJed Brown graph->nvtxs,&graph->special_dof);CHKERRQ(ierr); 1150674ae819SStefano Zampini /* zeroes memory */ 1151674ae819SStefano Zampini ierr = PetscMemzero(graph->count,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr); 1152674ae819SStefano Zampini ierr = PetscMemzero(graph->subset,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr); 115363602bcaSStefano Zampini /* use -1 as a default value for which_dof array */ 115463602bcaSStefano Zampini for (n=0;n<graph->nvtxs;n++) graph->which_dof[n] = -1; 1155674ae819SStefano Zampini ierr = PetscMemzero(graph->cptr,(graph->nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr); 1156674ae819SStefano Zampini ierr = PetscMemzero(graph->queue,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr); 115774e413f5SStefano Zampini ierr = PetscMemzero(graph->special_dof,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr); 1158674ae819SStefano Zampini /* zeroes first pointer to neighbour set */ 1159674ae819SStefano Zampini if (graph->nvtxs) { 1160674ae819SStefano Zampini graph->neighbours_set[0] = 0; 1161674ae819SStefano Zampini } 1162674ae819SStefano Zampini /* zeroes workspace for values of ncc */ 1163674ae819SStefano Zampini graph->subset_ncc = 0; 11643a5b03d1SStefano Zampini graph->subset_ref_node = 0; 1165674ae819SStefano Zampini PetscFunctionReturn(0); 1166674ae819SStefano Zampini } 1167674ae819SStefano Zampini 1168674ae819SStefano Zampini #undef __FUNCT__ 1169674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphDestroy" 1170674ae819SStefano Zampini PetscErrorCode PCBDDCGraphDestroy(PCBDDCGraph* graph) 1171674ae819SStefano Zampini { 1172674ae819SStefano Zampini PetscErrorCode ierr; 1173674ae819SStefano Zampini 1174674ae819SStefano Zampini PetscFunctionBegin; 1175674ae819SStefano Zampini ierr = PCBDDCGraphReset(*graph);CHKERRQ(ierr); 1176674ae819SStefano Zampini ierr = PetscFree(*graph);CHKERRQ(ierr); 1177674ae819SStefano Zampini PetscFunctionReturn(0); 1178674ae819SStefano Zampini } 1179674ae819SStefano Zampini 1180674ae819SStefano Zampini #undef __FUNCT__ 1181674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphCreate" 1182674ae819SStefano Zampini PetscErrorCode PCBDDCGraphCreate(PCBDDCGraph *graph) 1183674ae819SStefano Zampini { 1184674ae819SStefano Zampini PCBDDCGraph new_graph; 1185674ae819SStefano Zampini PetscErrorCode ierr; 1186674ae819SStefano Zampini 1187674ae819SStefano Zampini PetscFunctionBegin; 1188854ce69bSBarry Smith ierr = PetscNew(&new_graph);CHKERRQ(ierr); 1189674ae819SStefano Zampini /* local to global mapping of dofs */ 1190674ae819SStefano Zampini new_graph->l2gmap = 0; 1191674ae819SStefano Zampini /* vertex size */ 1192674ae819SStefano Zampini new_graph->nvtxs = 0; 1193674ae819SStefano Zampini new_graph->n_subsets = 0; 1194674ae819SStefano Zampini new_graph->custom_minimal_size = 1; 1195674ae819SStefano Zampini /* zeroes ponters */ 119651b0f6cfSStefano Zampini new_graph->mirrors = 0; 119751b0f6cfSStefano Zampini new_graph->mirrors_set = 0; 1198674ae819SStefano Zampini new_graph->neighbours_set = 0; 1199674ae819SStefano Zampini new_graph->subset = 0; 1200674ae819SStefano Zampini new_graph->which_dof = 0; 120174e413f5SStefano Zampini new_graph->special_dof = 0; 1202674ae819SStefano Zampini new_graph->cptr = 0; 1203674ae819SStefano Zampini new_graph->queue = 0; 1204674ae819SStefano Zampini new_graph->count = 0; 1205674ae819SStefano Zampini new_graph->subset_ncc = 0; 12063a5b03d1SStefano Zampini new_graph->subset_ref_node = 0; 1207674ae819SStefano Zampini new_graph->touched = 0; 1208674ae819SStefano Zampini /* zeroes pointers to csr graph of local nodes connectivity (optional data) */ 1209575ad6abSStefano Zampini new_graph->nvtxs_csr = 0; 1210674ae819SStefano Zampini new_graph->xadj = 0; 1211674ae819SStefano Zampini new_graph->adjncy = 0; 1212674ae819SStefano Zampini *graph = new_graph; 1213674ae819SStefano Zampini PetscFunctionReturn(0); 1214674ae819SStefano Zampini } 1215