1674ae819SStefano Zampini #include <petsc-private/petscimpl.h> 2ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddcprivate.h> 3ab5c6b0cSJed Brown #include <../src/ksp/pc/impls/bddc/bddcstructs.h> 4674ae819SStefano Zampini 5674ae819SStefano Zampini #undef __FUNCT__ 6674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphASCIIView" 7e49050b4SStefano Zampini PetscErrorCode PCBDDCGraphASCIIView(PCBDDCGraph graph, PetscInt verbosity_level, PetscViewer viewer) 8674ae819SStefano Zampini { 92b510759SStefano Zampini PetscInt i,j,tabs; 1093fb5fd3SStefano Zampini PetscInt* queue_in_global_numbering; 11674ae819SStefano Zampini PetscErrorCode ierr; 12674ae819SStefano Zampini 13674ae819SStefano Zampini PetscFunctionBegin; 14674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 152b510759SStefano Zampini ierr = PetscViewerASCIIGetTab(viewer,&tabs);CHKERRQ(ierr); 16674ae819SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"--------------------------------------------------\n");CHKERRQ(ierr); 17674ae819SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 18674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Local BDDC graph for subdomain %04d\n",PetscGlobalRank);CHKERRQ(ierr); 19674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Number of vertices %d\n",graph->nvtxs);CHKERRQ(ierr); 20e49050b4SStefano Zampini if (verbosity_level > 1) { 21674ae819SStefano Zampini for (i=0;i<graph->nvtxs;i++) { 22674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"%d:\n",i);CHKERRQ(ierr); 23674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer," which_dof: %d\n",graph->which_dof[i]);CHKERRQ(ierr); 2474e413f5SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer," special_dof: %d\n",graph->special_dof[i]);CHKERRQ(ierr); 25674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer," neighbours: %d\n",graph->count[i]);CHKERRQ(ierr); 262b510759SStefano Zampini ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 27674ae819SStefano Zampini if (graph->count[i]) { 28674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer," set of neighbours:");CHKERRQ(ierr); 29674ae819SStefano Zampini for (j=0;j<graph->count[i];j++) { 30674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d",graph->neighbours_set[i][j]);CHKERRQ(ierr); 31674ae819SStefano Zampini } 32674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr); 33674ae819SStefano Zampini } 342b510759SStefano Zampini ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr); 352b510759SStefano Zampini ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 3651b0f6cfSStefano Zampini if (graph->mirrors) { 3751b0f6cfSStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer," mirrors: %d\n",graph->mirrors[i]);CHKERRQ(ierr); 3851b0f6cfSStefano Zampini if (graph->mirrors[i]) { 392b510759SStefano Zampini ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 4051b0f6cfSStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer," set of mirrors:");CHKERRQ(ierr); 4151b0f6cfSStefano Zampini for (j=0;j<graph->mirrors[i];j++) { 4251b0f6cfSStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d",graph->mirrors_set[i][j]);CHKERRQ(ierr); 4351b0f6cfSStefano Zampini } 4451b0f6cfSStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr); 452b510759SStefano Zampini ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr); 462b510759SStefano Zampini ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 4751b0f6cfSStefano Zampini } 4851b0f6cfSStefano Zampini } 49e49050b4SStefano Zampini if (verbosity_level > 2) { 50674ae819SStefano Zampini if (graph->xadj && graph->adjncy) { 51674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer," local adj list:");CHKERRQ(ierr); 522b510759SStefano Zampini ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 53674ae819SStefano Zampini for (j=graph->xadj[i];j<graph->xadj[i+1];j++) { 54674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d",graph->adjncy[j]);CHKERRQ(ierr); 55674ae819SStefano Zampini } 56674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr); 572b510759SStefano Zampini ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr); 582b510759SStefano Zampini ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 59674ae819SStefano Zampini } 60e49050b4SStefano Zampini } 61674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer," interface subset id: %d\n",graph->subset[i]);CHKERRQ(ierr); 62674ae819SStefano Zampini if (graph->subset[i] && graph->subset_ncc) { 63674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer," ncc for subset: %d\n",graph->subset_ncc[graph->subset[i]-1]);CHKERRQ(ierr); 64674ae819SStefano Zampini } 65674ae819SStefano Zampini } 66e49050b4SStefano Zampini } 67674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Total number of connected components %d\n",graph->ncc);CHKERRQ(ierr); 68785e854fSJed Brown ierr = PetscMalloc1(graph->cptr[graph->ncc],&queue_in_global_numbering);CHKERRQ(ierr); 6993fb5fd3SStefano Zampini ierr = ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_in_global_numbering);CHKERRQ(ierr); 70674ae819SStefano Zampini for (i=0;i<graph->ncc;i++) { 711ce3d396SStefano Zampini PetscInt node_num=graph->queue[graph->cptr[i]]; 721ce3d396SStefano Zampini PetscBool printcc = PETSC_FALSE; 73674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d (neighs:",i);CHKERRQ(ierr); 742b510759SStefano Zampini ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 75674ae819SStefano Zampini for (j=0;j<graph->count[node_num];j++) { 76674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d",graph->neighbours_set[node_num][j]);CHKERRQ(ierr); 77674ae819SStefano Zampini } 78674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"):");CHKERRQ(ierr); 791ce3d396SStefano Zampini if (verbosity_level > 1) { 801ce3d396SStefano Zampini printcc = PETSC_TRUE; 81e635b14bSStefano Zampini } else if (graph->count[node_num] > 1 || (graph->count[node_num] == 1 && graph->special_dof[node_num] == PCBDDCGRAPH_NEUMANN_MARK)) { 821ce3d396SStefano Zampini printcc = PETSC_TRUE; 831ce3d396SStefano Zampini } 841ce3d396SStefano Zampini if (printcc) { 85674ae819SStefano Zampini for (j=graph->cptr[i];j<graph->cptr[i+1];j++) { 8693fb5fd3SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d (%d)",graph->queue[j],queue_in_global_numbering[j]);CHKERRQ(ierr); 87674ae819SStefano Zampini } 881ce3d396SStefano Zampini } 89674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr); 902b510759SStefano Zampini ierr = PetscViewerASCIISetTab(viewer,tabs);CHKERRQ(ierr); 912b510759SStefano Zampini ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 92674ae819SStefano Zampini } 9393fb5fd3SStefano Zampini ierr = PetscFree(queue_in_global_numbering);CHKERRQ(ierr); 94e49050b4SStefano Zampini if (graph->custom_minimal_size > 1 && verbosity_level > 1) { 95674ae819SStefano Zampini ierr = PetscViewerASCIISynchronizedPrintf(viewer,"Custom minimal size %d\n",graph->custom_minimal_size);CHKERRQ(ierr); 96674ae819SStefano Zampini } 97674ae819SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 98674ae819SStefano Zampini PetscFunctionReturn(0); 99674ae819SStefano Zampini } 100674ae819SStefano Zampini 101674ae819SStefano Zampini #undef __FUNCT__ 102674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphGetCandidatesIS" 103674ae819SStefano 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) 104674ae819SStefano Zampini { 105674ae819SStefano Zampini IS *ISForFaces,*ISForEdges,ISForVertices; 106674ae819SStefano Zampini PetscInt i,j,nfc,nec,nvc,*idx; 107674ae819SStefano Zampini PetscBool twodim_flag; 108674ae819SStefano Zampini PetscErrorCode ierr; 109674ae819SStefano Zampini 110674ae819SStefano Zampini PetscFunctionBegin; 111674ae819SStefano Zampini /* loop on ccs to evalute number of faces, edges and vertices */ 112674ae819SStefano Zampini nfc = 0; 113674ae819SStefano Zampini nec = 0; 114674ae819SStefano Zampini nvc = 0; 115674ae819SStefano Zampini twodim_flag = PETSC_FALSE; 116674ae819SStefano Zampini for (i=0;i<graph->ncc;i++) { 1171e1f2d93SStefano Zampini PetscInt repdof = graph->queue[graph->cptr[i]]; 118674ae819SStefano Zampini if (graph->cptr[i+1]-graph->cptr[i] > graph->custom_minimal_size) { 119e635b14bSStefano Zampini if (graph->count[repdof] == 1 && graph->special_dof[repdof] != PCBDDCGRAPH_NEUMANN_MARK) { 120674ae819SStefano Zampini nfc++; 121674ae819SStefano Zampini } else { /* note that nec will be zero in 2d */ 122674ae819SStefano Zampini nec++; 123674ae819SStefano Zampini } 124674ae819SStefano Zampini } else { 125674ae819SStefano Zampini nvc += graph->cptr[i+1]-graph->cptr[i]; 126674ae819SStefano Zampini } 127674ae819SStefano Zampini } 128aa9fcddfSStefano Zampini j=0; 129aa9fcddfSStefano Zampini ierr = MPI_Allreduce(&nec,&j,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)graph->l2gmap));CHKERRQ(ierr); 130aa9fcddfSStefano Zampini if (!j) { /* we are in a 2D case -> no faces, only edges */ 131674ae819SStefano Zampini nec = nfc; 132674ae819SStefano Zampini nfc = 0; 133674ae819SStefano Zampini twodim_flag = PETSC_TRUE; 134674ae819SStefano Zampini } 135674ae819SStefano Zampini /* allocate IS arrays for faces, edges. Vertices need a single index set. */ 136674ae819SStefano Zampini ISForFaces = 0; 137674ae819SStefano Zampini ISForEdges = 0; 138674ae819SStefano Zampini ISForVertices = 0; 1391ee8691aSStefano Zampini if (use_faces && nfc) { 140785e854fSJed Brown ierr = PetscMalloc1(nfc,&ISForFaces);CHKERRQ(ierr); 141674ae819SStefano Zampini } 142b8ffe317SStefano Zampini if (use_edges && nec) { 143785e854fSJed Brown ierr = PetscMalloc1(nec,&ISForEdges);CHKERRQ(ierr); 144674ae819SStefano Zampini } 145b8ffe317SStefano Zampini if (use_vertices && nvc) { 146785e854fSJed Brown ierr = PetscMalloc1(nvc,&idx);CHKERRQ(ierr); 147674ae819SStefano Zampini } 148674ae819SStefano Zampini /* loop on ccs to compute index sets for faces and edges */ 149674ae819SStefano Zampini nfc = 0; 150674ae819SStefano Zampini nec = 0; 151674ae819SStefano Zampini for (i=0;i<graph->ncc;i++) { 1521e1f2d93SStefano Zampini PetscInt repdof = graph->queue[graph->cptr[i]]; 153674ae819SStefano Zampini if (graph->cptr[i+1]-graph->cptr[i] > graph->custom_minimal_size) { 154e635b14bSStefano Zampini if (graph->count[repdof] == 1 && graph->special_dof[repdof] != PCBDDCGRAPH_NEUMANN_MARK) { 155674ae819SStefano Zampini if (twodim_flag) { 156674ae819SStefano Zampini if (use_edges) { 157674ae819SStefano 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); 158674ae819SStefano Zampini nec++; 159674ae819SStefano Zampini } 160674ae819SStefano Zampini } else { 161674ae819SStefano Zampini if (use_faces) { 162674ae819SStefano 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); 163674ae819SStefano Zampini nfc++; 164674ae819SStefano Zampini } 165674ae819SStefano Zampini } 166674ae819SStefano Zampini } else { 167674ae819SStefano Zampini if (use_edges) { 168674ae819SStefano 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); 169674ae819SStefano Zampini nec++; 170674ae819SStefano Zampini } 171674ae819SStefano Zampini } 172674ae819SStefano Zampini } 173674ae819SStefano Zampini } 174674ae819SStefano Zampini /* index set for vertices */ 175b8ffe317SStefano Zampini if (use_vertices && nvc) { 176b8ffe317SStefano Zampini nvc = 0; 177674ae819SStefano Zampini for (i=0;i<graph->ncc;i++) { 178674ae819SStefano Zampini if (graph->cptr[i+1]-graph->cptr[i] <= graph->custom_minimal_size) { 179674ae819SStefano Zampini for (j=graph->cptr[i];j<graph->cptr[i+1];j++) { 180674ae819SStefano Zampini idx[nvc]=graph->queue[j]; 181674ae819SStefano Zampini nvc++; 182674ae819SStefano Zampini } 183674ae819SStefano Zampini } 184674ae819SStefano Zampini } 185674ae819SStefano Zampini /* sort vertex set (by local ordering) */ 186674ae819SStefano Zampini ierr = PetscSortInt(nvc,idx);CHKERRQ(ierr); 187674ae819SStefano Zampini ierr = ISCreateGeneral(PETSC_COMM_SELF,nvc,idx,PETSC_OWN_POINTER,&ISForVertices);CHKERRQ(ierr); 188674ae819SStefano Zampini } 189674ae819SStefano Zampini /* get back info */ 190674ae819SStefano Zampini *n_faces = nfc; 191*d7796978SStefano Zampini if (FacesIS) { 192674ae819SStefano Zampini *FacesIS = ISForFaces; 193*d7796978SStefano Zampini } else { 194*d7796978SStefano Zampini for (i=0;i<nfc;i++) { 195*d7796978SStefano Zampini ierr = ISDestroy(&ISForFaces[i]);CHKERRQ(ierr); 196*d7796978SStefano Zampini } 197*d7796978SStefano Zampini ierr = PetscFree(ISForFaces);CHKERRQ(ierr); 198*d7796978SStefano Zampini } 199674ae819SStefano Zampini *n_edges = nec; 200*d7796978SStefano Zampini if (EdgesIS) { 201674ae819SStefano Zampini *EdgesIS = ISForEdges; 202*d7796978SStefano Zampini } else { 203*d7796978SStefano Zampini for (i=0;i<nec;i++) { 204*d7796978SStefano Zampini ierr = ISDestroy(&ISForEdges[i]);CHKERRQ(ierr); 205*d7796978SStefano Zampini } 206*d7796978SStefano Zampini ierr = PetscFree(ISForEdges);CHKERRQ(ierr); 207*d7796978SStefano Zampini } 208*d7796978SStefano Zampini if (VerticesIS) { 209674ae819SStefano Zampini *VerticesIS = ISForVertices; 210*d7796978SStefano Zampini } else { 211*d7796978SStefano Zampini ierr = ISDestroy(&ISForVertices);CHKERRQ(ierr); 212*d7796978SStefano Zampini } 213674ae819SStefano Zampini PetscFunctionReturn(0); 214674ae819SStefano Zampini } 215674ae819SStefano Zampini 216674ae819SStefano Zampini #undef __FUNCT__ 217674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphComputeConnectedComponents" 218674ae819SStefano Zampini PetscErrorCode PCBDDCGraphComputeConnectedComponents(PCBDDCGraph graph) 219674ae819SStefano Zampini { 2205b1b9aeaSStefano Zampini PetscInt i,adapt_interface,adapt_interface_reduced; 221674ae819SStefano Zampini MPI_Comm interface_comm; 222674ae819SStefano Zampini PetscErrorCode ierr; 223674ae819SStefano Zampini 224674ae819SStefano Zampini PetscFunctionBegin; 225674ae819SStefano Zampini /* compute connected components locally */ 226674ae819SStefano Zampini ierr = PetscObjectGetComm((PetscObject)(graph->l2gmap),&interface_comm);CHKERRQ(ierr); 227674ae819SStefano Zampini ierr = PCBDDCGraphComputeConnectedComponentsLocal(graph);CHKERRQ(ierr); 228674ae819SStefano Zampini /* check consistency of connected components among neighbouring subdomains -> it adapt them in case it is needed */ 229674ae819SStefano Zampini adapt_interface = 0; 230674ae819SStefano Zampini adapt_interface_reduced = 0; 231674ae819SStefano Zampini for (i=0;i<graph->n_subsets;i++) { 232674ae819SStefano Zampini /* We are not sure that on a given subset of the local interface, 233674ae819SStefano Zampini with two connected components, the latters be the same among sharing subdomains */ 234674ae819SStefano Zampini if (graph->subset_ncc[i] > 1) { 235674ae819SStefano Zampini adapt_interface=1; 236674ae819SStefano Zampini break; 237674ae819SStefano Zampini } 238674ae819SStefano Zampini } 239674ae819SStefano Zampini ierr = MPI_Allreduce(&adapt_interface,&adapt_interface_reduced,1,MPIU_INT,MPI_LOR,interface_comm);CHKERRQ(ierr); 240674ae819SStefano Zampini 241674ae819SStefano Zampini if (graph->n_subsets && adapt_interface_reduced) { 2425b1b9aeaSStefano Zampini MPI_Request *send_requests; 2435b1b9aeaSStefano Zampini MPI_Request *recv_requests; 2445b1b9aeaSStefano Zampini PetscInt *aux_new_xadj,*new_xadj,*new_adjncy,**temp_buffer; 2455b1b9aeaSStefano Zampini PetscInt *old_xadj,*old_adjncy; 2465b1b9aeaSStefano Zampini PetscInt j,k,s,sum_requests,buffer_size,size_of_recv,temp_buffer_size; 2475b1b9aeaSStefano Zampini PetscMPIInt rank,neigh,tag,mpi_buffer_size; 2485b1b9aeaSStefano Zampini PetscInt *cum_recv_counts,*subset_to_nodes_indices,*recv_buffer_subset,*nodes_to_temp_buffer_indices; 2495b1b9aeaSStefano Zampini PetscInt *send_buffer,*recv_buffer,*queue_in_global_numbering,*sizes_of_sends,*add_to_subset; 2505b1b9aeaSStefano Zampini PetscInt start_of_recv,start_of_send,size_of_send,global_subset_counter,ins_val; 2515b1b9aeaSStefano Zampini PetscBool *subset_cc_adapt,same_set; 2525b1b9aeaSStefano Zampini 253674ae819SStefano Zampini /* Retrict adjacency graph using information from previously computed connected components */ 254785e854fSJed Brown ierr = PetscMalloc1(graph->nvtxs,&aux_new_xadj);CHKERRQ(ierr); 255674ae819SStefano Zampini for (i=0;i<graph->nvtxs;i++) { 256674ae819SStefano Zampini aux_new_xadj[i]=1; 257674ae819SStefano Zampini } 258674ae819SStefano Zampini for (i=0;i<graph->ncc;i++) { 259674ae819SStefano Zampini k = graph->cptr[i+1]-graph->cptr[i]; 260674ae819SStefano Zampini for (j=0;j<k;j++) { 261674ae819SStefano Zampini aux_new_xadj[graph->queue[graph->cptr[i]+j]]=k; 262674ae819SStefano Zampini } 263674ae819SStefano Zampini } 264674ae819SStefano Zampini j = 0; 265674ae819SStefano Zampini for (i=0;i<graph->nvtxs;i++) { 266674ae819SStefano Zampini j += aux_new_xadj[i]; 267674ae819SStefano Zampini } 268785e854fSJed Brown ierr = PetscMalloc1((graph->nvtxs+1),&new_xadj);CHKERRQ(ierr); 269785e854fSJed Brown ierr = PetscMalloc1(j,&new_adjncy);CHKERRQ(ierr); 270674ae819SStefano Zampini new_xadj[0]=0; 271674ae819SStefano Zampini for (i=0;i<graph->nvtxs;i++) { 272674ae819SStefano Zampini new_xadj[i+1]=new_xadj[i]+aux_new_xadj[i]; 273674ae819SStefano Zampini if (aux_new_xadj[i]==1) { 274674ae819SStefano Zampini new_adjncy[new_xadj[i]]=i; 275674ae819SStefano Zampini } 276674ae819SStefano Zampini } 277674ae819SStefano Zampini ierr = PetscFree(aux_new_xadj);CHKERRQ(ierr); 278674ae819SStefano Zampini for (i=0;i<graph->ncc;i++) { 279674ae819SStefano Zampini k = graph->cptr[i+1]-graph->cptr[i]; 280674ae819SStefano Zampini for (j=0;j<k;j++) { 281674ae819SStefano Zampini ierr = PetscMemcpy(&new_adjncy[new_xadj[graph->queue[graph->cptr[i]+j]]],&graph->queue[graph->cptr[i]],k*sizeof(PetscInt));CHKERRQ(ierr); 282674ae819SStefano Zampini } 283674ae819SStefano Zampini } 2845b1b9aeaSStefano Zampini /* set temporarly new CSR into graph */ 2855b1b9aeaSStefano Zampini old_xadj = graph->xadj; 2865b1b9aeaSStefano Zampini old_adjncy = graph->adjncy; 287674ae819SStefano Zampini graph->xadj = new_xadj; 288674ae819SStefano Zampini graph->adjncy = new_adjncy; 289674ae819SStefano Zampini /* allocate some space */ 290674ae819SStefano Zampini ierr = MPI_Comm_rank(interface_comm,&rank);CHKERRQ(ierr); 291785e854fSJed Brown ierr = PetscMalloc1((graph->n_subsets+1),&cum_recv_counts);CHKERRQ(ierr); 292674ae819SStefano Zampini ierr = PetscMemzero(cum_recv_counts,(graph->n_subsets+1)*sizeof(*cum_recv_counts));CHKERRQ(ierr); 293785e854fSJed Brown ierr = PetscMalloc1(graph->n_subsets,&subset_to_nodes_indices);CHKERRQ(ierr); 294674ae819SStefano Zampini /* first count how many neighbours per connected component I will receive from */ 295674ae819SStefano Zampini cum_recv_counts[0]=0; 296674ae819SStefano Zampini for (i=1;i<graph->n_subsets+1;i++) { 297674ae819SStefano Zampini j = 0; 298674ae819SStefano Zampini while (graph->subset[j] != i) { 299674ae819SStefano Zampini j++; 300674ae819SStefano Zampini } 301674ae819SStefano Zampini subset_to_nodes_indices[i-1]=j; 302674ae819SStefano Zampini /* We don't want sends/recvs_to/from_self -> here I don't count myself */ 303674ae819SStefano Zampini cum_recv_counts[i]=cum_recv_counts[i-1]+graph->count[j]; 304674ae819SStefano Zampini } 305785e854fSJed Brown ierr = PetscMalloc1(2*cum_recv_counts[graph->n_subsets],&recv_buffer_subset);CHKERRQ(ierr); 306785e854fSJed Brown ierr = PetscMalloc1(cum_recv_counts[graph->n_subsets],&send_requests);CHKERRQ(ierr); 307785e854fSJed Brown ierr = PetscMalloc1(cum_recv_counts[graph->n_subsets],&recv_requests);CHKERRQ(ierr); 308674ae819SStefano Zampini for (i=0;i<cum_recv_counts[graph->n_subsets];i++) { 309674ae819SStefano Zampini send_requests[i]=MPI_REQUEST_NULL; 310674ae819SStefano Zampini recv_requests[i]=MPI_REQUEST_NULL; 311674ae819SStefano Zampini } 312674ae819SStefano Zampini /* exchange with my neighbours the number of my connected components on the shared interface */ 313674ae819SStefano Zampini sum_requests = 0; 314674ae819SStefano Zampini for (i=0;i<graph->n_subsets;i++) { 315674ae819SStefano Zampini j = subset_to_nodes_indices[i]; 3163a5b03d1SStefano Zampini ierr = PetscMPIIntCast(graph->subset_ref_node[i],&tag);CHKERRQ(ierr); 317674ae819SStefano Zampini for (k=0;k<graph->count[j];k++) { 318674ae819SStefano Zampini ierr = PetscMPIIntCast(graph->neighbours_set[j][k],&neigh);CHKERRQ(ierr); 319674ae819SStefano Zampini ierr = MPI_Isend(&graph->subset_ncc[i],1,MPIU_INT,neigh,tag,interface_comm,&send_requests[sum_requests]);CHKERRQ(ierr); 320674ae819SStefano Zampini ierr = MPI_Irecv(&recv_buffer_subset[sum_requests],1,MPIU_INT,neigh,tag,interface_comm,&recv_requests[sum_requests]);CHKERRQ(ierr); 321674ae819SStefano Zampini sum_requests++; 322674ae819SStefano Zampini } 323674ae819SStefano Zampini } 324674ae819SStefano Zampini ierr = MPI_Waitall(sum_requests,recv_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 325674ae819SStefano Zampini ierr = MPI_Waitall(sum_requests,send_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 326674ae819SStefano Zampini /* determine the connected component I need to adapt */ 327785e854fSJed Brown ierr = PetscMalloc1(graph->n_subsets,&subset_cc_adapt);CHKERRQ(ierr); 328674ae819SStefano Zampini ierr = PetscMemzero(subset_cc_adapt,graph->n_subsets*sizeof(*subset_cc_adapt));CHKERRQ(ierr); 329674ae819SStefano Zampini for (i=0;i<graph->n_subsets;i++) { 330674ae819SStefano Zampini for (j=cum_recv_counts[i];j<cum_recv_counts[i+1];j++){ 331674ae819SStefano Zampini /* The first condition is natural (someone has a different number of ccs than me), the second one is just to be safe */ 332674ae819SStefano Zampini if (graph->subset_ncc[i] != recv_buffer_subset[j] || graph->subset_ncc[i] > 1) { 333674ae819SStefano Zampini subset_cc_adapt[i] = PETSC_TRUE; 334674ae819SStefano Zampini break; 335674ae819SStefano Zampini } 336674ae819SStefano Zampini } 337674ae819SStefano Zampini } 338674ae819SStefano Zampini buffer_size = 0; 339674ae819SStefano Zampini for (i=0;i<graph->n_subsets;i++) { 340674ae819SStefano Zampini if (subset_cc_adapt[i]) { 341674ae819SStefano Zampini for (j=i;j<graph->ncc;j++) { 342674ae819SStefano Zampini if (graph->subset[graph->queue[graph->cptr[j]]] == i+1) { /* WARNING -> subset values goes from 1 to graph->n_subsets included */ 343674ae819SStefano Zampini buffer_size += 1 + graph->cptr[j+1]-graph->cptr[j]; 344674ae819SStefano Zampini } 345674ae819SStefano Zampini } 346674ae819SStefano Zampini } 347674ae819SStefano Zampini } 348785e854fSJed Brown ierr = PetscMalloc1(buffer_size,&send_buffer);CHKERRQ(ierr); 349674ae819SStefano Zampini /* now get from neighbours their ccs (in global numbering) and adapt them (in case it is needed) */ 350785e854fSJed Brown ierr = PetscMalloc1(graph->cptr[graph->ncc],&queue_in_global_numbering);CHKERRQ(ierr); 351674ae819SStefano Zampini ierr = ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_in_global_numbering);CHKERRQ(ierr); 352674ae819SStefano Zampini /* determine how much data to send (size of each queue plus the global indices) and communicate it to neighbours */ 353785e854fSJed Brown ierr = PetscMalloc1(graph->n_subsets,&sizes_of_sends);CHKERRQ(ierr); 354674ae819SStefano Zampini ierr = PetscMemzero(sizes_of_sends,graph->n_subsets*sizeof(*sizes_of_sends));CHKERRQ(ierr); 355674ae819SStefano Zampini sum_requests = 0; 356674ae819SStefano Zampini start_of_send = 0; 357674ae819SStefano Zampini start_of_recv = cum_recv_counts[graph->n_subsets]; 358674ae819SStefano Zampini for (i=0;i<graph->n_subsets;i++) { 359674ae819SStefano Zampini if (subset_cc_adapt[i]) { 360674ae819SStefano Zampini size_of_send = 0; 361674ae819SStefano Zampini for (j=i;j<graph->ncc;j++) { 362674ae819SStefano Zampini if (graph->subset[graph->queue[graph->cptr[j]]] == i+1) { /* WARNING -> subset values goes from 1 to graph->n_subsets included */ 363674ae819SStefano Zampini send_buffer[start_of_send+size_of_send]=graph->cptr[j+1]-graph->cptr[j]; 364674ae819SStefano Zampini size_of_send += 1; 365674ae819SStefano Zampini ierr = PetscMemcpy(&send_buffer[start_of_send+size_of_send], 366674ae819SStefano Zampini &queue_in_global_numbering[graph->cptr[j]], 367674ae819SStefano Zampini (graph->cptr[j+1]-graph->cptr[j])*sizeof(*send_buffer));CHKERRQ(ierr); 368674ae819SStefano Zampini size_of_send = size_of_send+graph->cptr[j+1]-graph->cptr[j]; 369674ae819SStefano Zampini } 370674ae819SStefano Zampini } 371674ae819SStefano Zampini j = subset_to_nodes_indices[i]; 372674ae819SStefano Zampini sizes_of_sends[i] = size_of_send; 3733a5b03d1SStefano Zampini ierr = PetscMPIIntCast(graph->subset_ref_node[i]+1,&tag);CHKERRQ(ierr); 374674ae819SStefano Zampini for (k=0;k<graph->count[j];k++) { 375674ae819SStefano Zampini ierr = PetscMPIIntCast(graph->neighbours_set[j][k],&neigh);CHKERRQ(ierr); 376674ae819SStefano Zampini ierr = MPI_Isend(&sizes_of_sends[i],1,MPIU_INT,neigh,tag,interface_comm,&send_requests[sum_requests]);CHKERRQ(ierr); 377674ae819SStefano 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); 378674ae819SStefano Zampini sum_requests++; 379674ae819SStefano Zampini } 380674ae819SStefano Zampini start_of_send += size_of_send; 381674ae819SStefano Zampini } 382674ae819SStefano Zampini } 383674ae819SStefano Zampini ierr = MPI_Waitall(sum_requests,send_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 384674ae819SStefano Zampini ierr = MPI_Waitall(sum_requests,recv_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 385674ae819SStefano Zampini buffer_size = 0; 386674ae819SStefano Zampini for (k=0;k<sum_requests;k++) { 387674ae819SStefano Zampini buffer_size += recv_buffer_subset[start_of_recv+k]; 388674ae819SStefano Zampini } 389785e854fSJed Brown ierr = PetscMalloc1(buffer_size,&recv_buffer);CHKERRQ(ierr); 390674ae819SStefano Zampini /* now exchange the data */ 391674ae819SStefano Zampini start_of_recv = 0; 392674ae819SStefano Zampini start_of_send = 0; 393674ae819SStefano Zampini sum_requests = 0; 394674ae819SStefano Zampini for (i=0;i<graph->n_subsets;i++) { 395674ae819SStefano Zampini if (subset_cc_adapt[i]) { 396674ae819SStefano Zampini size_of_send = sizes_of_sends[i]; 397674ae819SStefano Zampini j = subset_to_nodes_indices[i]; 3983a5b03d1SStefano Zampini ierr = PetscMPIIntCast(graph->subset_ref_node[i]+2,&tag);CHKERRQ(ierr); 399674ae819SStefano Zampini for (k=0;k<graph->count[j];k++) { 400674ae819SStefano Zampini ierr = PetscMPIIntCast(graph->neighbours_set[j][k],&neigh);CHKERRQ(ierr); 401674ae819SStefano Zampini ierr = PetscMPIIntCast(size_of_send,&mpi_buffer_size);CHKERRQ(ierr); 402674ae819SStefano Zampini ierr = MPI_Isend(&send_buffer[start_of_send],mpi_buffer_size,MPIU_INT,neigh,tag,interface_comm,&send_requests[sum_requests]);CHKERRQ(ierr); 403674ae819SStefano Zampini size_of_recv = recv_buffer_subset[cum_recv_counts[graph->n_subsets]+sum_requests]; 404674ae819SStefano Zampini ierr = PetscMPIIntCast(size_of_recv,&mpi_buffer_size);CHKERRQ(ierr); 405674ae819SStefano Zampini ierr = MPI_Irecv(&recv_buffer[start_of_recv],mpi_buffer_size,MPIU_INT,neigh,tag,interface_comm,&recv_requests[sum_requests]);CHKERRQ(ierr); 406674ae819SStefano Zampini start_of_recv += size_of_recv; 407674ae819SStefano Zampini sum_requests++; 408674ae819SStefano Zampini } 409674ae819SStefano Zampini start_of_send += size_of_send; 410674ae819SStefano Zampini } 411674ae819SStefano Zampini } 412674ae819SStefano Zampini ierr = MPI_Waitall(sum_requests,recv_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 413674ae819SStefano Zampini ierr = MPI_Waitall(sum_requests,send_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr); 414674ae819SStefano Zampini for (j=0;j<buffer_size;) { 415674ae819SStefano Zampini ierr = ISGlobalToLocalMappingApply(graph->l2gmap,IS_GTOLM_MASK,recv_buffer[j],&recv_buffer[j+1],&recv_buffer[j],&recv_buffer[j+1]);CHKERRQ(ierr); 41651b0f6cfSStefano Zampini /* we need to adapt the output of GlobalToLocal mapping if there are mirrored nodes */ 41751b0f6cfSStefano Zampini if (graph->mirrors) { 41851b0f6cfSStefano Zampini PetscBool mirrored_found=PETSC_FALSE; 41951b0f6cfSStefano Zampini for (k=0;k<recv_buffer[j];k++) { 42051b0f6cfSStefano Zampini if (graph->mirrors[recv_buffer[j+k+1]]) { 42151b0f6cfSStefano Zampini mirrored_found=PETSC_TRUE; 42251b0f6cfSStefano Zampini recv_buffer[j+k+1]=graph->mirrors_set[recv_buffer[j+k+1]][0]; 42351b0f6cfSStefano Zampini } 42451b0f6cfSStefano Zampini } 42551b0f6cfSStefano Zampini if (mirrored_found) { 42651b0f6cfSStefano Zampini ierr = PetscSortInt(recv_buffer[j],&recv_buffer[j+1]);CHKERRQ(ierr); 42751b0f6cfSStefano Zampini k=0; 42851b0f6cfSStefano Zampini while (k<recv_buffer[j]) { 42951b0f6cfSStefano Zampini for (s=1;s<graph->mirrors[recv_buffer[j+1+k]];s++) { 43051b0f6cfSStefano Zampini recv_buffer[j+1+k+s] = graph->mirrors_set[recv_buffer[j+1+k]][s]; 43151b0f6cfSStefano Zampini } 43251b0f6cfSStefano Zampini k+=graph->mirrors[recv_buffer[j+1+k]]+s; 43351b0f6cfSStefano Zampini } 43451b0f6cfSStefano Zampini } 43551b0f6cfSStefano Zampini } 436674ae819SStefano Zampini k = recv_buffer[j]+1; 437674ae819SStefano Zampini j += k; 438674ae819SStefano Zampini } 439674ae819SStefano Zampini sum_requests = cum_recv_counts[graph->n_subsets]; 440674ae819SStefano Zampini start_of_recv = 0; 441785e854fSJed Brown ierr = PetscMalloc1(graph->nvtxs,&nodes_to_temp_buffer_indices);CHKERRQ(ierr); 442674ae819SStefano Zampini global_subset_counter = 0; 443674ae819SStefano Zampini for (i=0;i<graph->n_subsets;i++) { 444674ae819SStefano Zampini if (subset_cc_adapt[i]) { 445674ae819SStefano Zampini temp_buffer_size = 0; 446674ae819SStefano Zampini /* find nodes on the shared interface we need to adapt */ 447674ae819SStefano Zampini for (j=0;j<graph->nvtxs;j++) { 448674ae819SStefano Zampini if (graph->subset[j]==i+1) { 449674ae819SStefano Zampini nodes_to_temp_buffer_indices[j] = temp_buffer_size; 450674ae819SStefano Zampini temp_buffer_size++; 451674ae819SStefano Zampini } else { 452674ae819SStefano Zampini nodes_to_temp_buffer_indices[j] = -1; 453674ae819SStefano Zampini } 454674ae819SStefano Zampini } 455674ae819SStefano Zampini /* allocate some temporary space */ 456785e854fSJed Brown ierr = PetscMalloc1(temp_buffer_size,&temp_buffer);CHKERRQ(ierr); 457785e854fSJed Brown ierr = PetscMalloc1(temp_buffer_size*(cum_recv_counts[i+1]-cum_recv_counts[i]),&temp_buffer[0]);CHKERRQ(ierr); 458674ae819SStefano Zampini ierr = PetscMemzero(temp_buffer[0],temp_buffer_size*(cum_recv_counts[i+1]-cum_recv_counts[i])*sizeof(PetscInt));CHKERRQ(ierr); 459674ae819SStefano Zampini for (j=1;j<temp_buffer_size;j++) { 460674ae819SStefano Zampini temp_buffer[j] = temp_buffer[j-1]+cum_recv_counts[i+1]-cum_recv_counts[i]; 461674ae819SStefano Zampini } 462674ae819SStefano Zampini /* analyze contributions from neighbouring subdomains for i-th conn comp 463674ae819SStefano Zampini temp buffer structure: 464674ae819SStefano Zampini supposing part of the interface has dimension 5 (for example with global dofs 0,1,2,3,4) 465674ae819SStefano Zampini 3 neighs procs with structured connected components: 466674ae819SStefano Zampini neigh 0: [0 1 4], [2 3]; (2 connected components) 467674ae819SStefano Zampini neigh 1: [0 1], [2 3 4]; (2 connected components) 468674ae819SStefano Zampini neigh 2: [0 4], [1], [2 3]; (3 connected components) 469674ae819SStefano Zampini tempbuffer (row-oriented) will be filled as: 470674ae819SStefano Zampini [ 0, 0, 0; 471674ae819SStefano Zampini 0, 0, 1; 472674ae819SStefano Zampini 1, 1, 2; 473674ae819SStefano Zampini 1, 1, 2; 474674ae819SStefano Zampini 0, 1, 0; ]; 475674ae819SStefano Zampini This way we can simply find intersections of ccs among neighs. 476674ae819SStefano Zampini For the example above, the graph->subset array will be modified to reproduce the following 4 connected components [0], [1], [2 3], [4]; 477674ae819SStefano Zampini */ 478674ae819SStefano Zampini for (j=0;j<cum_recv_counts[i+1]-cum_recv_counts[i];j++) { 479674ae819SStefano Zampini ins_val = 0; 480674ae819SStefano Zampini size_of_recv = recv_buffer_subset[sum_requests]; /* total size of recv from neighs */ 481674ae819SStefano Zampini for (buffer_size=0;buffer_size<size_of_recv;) { /* loop until all data from neighs has been taken into account */ 482674ae819SStefano Zampini for (k=1;k<recv_buffer[buffer_size+start_of_recv]+1;k++) { /* filling properly temp_buffer using data from a single recv */ 483674ae819SStefano Zampini temp_buffer[nodes_to_temp_buffer_indices[recv_buffer[start_of_recv+buffer_size+k]]][j] = ins_val; 484674ae819SStefano Zampini } 485674ae819SStefano Zampini buffer_size += k; 486674ae819SStefano Zampini ins_val++; 487674ae819SStefano Zampini } 488674ae819SStefano Zampini start_of_recv += size_of_recv; 489674ae819SStefano Zampini sum_requests++; 490674ae819SStefano Zampini } 491785e854fSJed Brown ierr = PetscMalloc1(temp_buffer_size,&add_to_subset);CHKERRQ(ierr); 492674ae819SStefano Zampini ierr = PetscMemzero(add_to_subset,temp_buffer_size*sizeof(*add_to_subset));CHKERRQ(ierr); 493674ae819SStefano Zampini for (j=0;j<temp_buffer_size;j++) { 494674ae819SStefano Zampini if (!add_to_subset[j]) { /* found a new cc */ 495674ae819SStefano Zampini global_subset_counter++; 496674ae819SStefano Zampini add_to_subset[j] = global_subset_counter; 497674ae819SStefano Zampini for (k=j+1;k<temp_buffer_size;k++) { /* check for other nodes in new cc */ 498674ae819SStefano Zampini same_set = PETSC_TRUE; 499674ae819SStefano Zampini for (s=0;s<cum_recv_counts[i+1]-cum_recv_counts[i];s++) { 500674ae819SStefano Zampini if (temp_buffer[j][s]!=temp_buffer[k][s]) { 501674ae819SStefano Zampini same_set = PETSC_FALSE; 502674ae819SStefano Zampini break; 503674ae819SStefano Zampini } 504674ae819SStefano Zampini } 505674ae819SStefano Zampini if (same_set) { 506674ae819SStefano Zampini add_to_subset[k] = global_subset_counter; 507674ae819SStefano Zampini } 508674ae819SStefano Zampini } 509674ae819SStefano Zampini } 510674ae819SStefano Zampini } 511674ae819SStefano Zampini /* insert new data in subset array */ 512674ae819SStefano Zampini temp_buffer_size = 0; 513674ae819SStefano Zampini for (j=0;j<graph->nvtxs;j++) { 514674ae819SStefano Zampini if (graph->subset[j]==i+1) { 515674ae819SStefano Zampini graph->subset[j] = graph->n_subsets+add_to_subset[temp_buffer_size]; 516674ae819SStefano Zampini temp_buffer_size++; 517674ae819SStefano Zampini } 518674ae819SStefano Zampini } 519674ae819SStefano Zampini ierr = PetscFree(temp_buffer[0]);CHKERRQ(ierr); 520674ae819SStefano Zampini ierr = PetscFree(temp_buffer);CHKERRQ(ierr); 521674ae819SStefano Zampini ierr = PetscFree(add_to_subset);CHKERRQ(ierr); 522674ae819SStefano Zampini } 523674ae819SStefano Zampini } 524674ae819SStefano Zampini ierr = PetscFree(nodes_to_temp_buffer_indices);CHKERRQ(ierr); 525674ae819SStefano Zampini ierr = PetscFree(sizes_of_sends);CHKERRQ(ierr); 526674ae819SStefano Zampini ierr = PetscFree(send_requests);CHKERRQ(ierr); 527674ae819SStefano Zampini ierr = PetscFree(recv_requests);CHKERRQ(ierr); 528674ae819SStefano Zampini ierr = PetscFree(recv_buffer);CHKERRQ(ierr); 529674ae819SStefano Zampini ierr = PetscFree(recv_buffer_subset);CHKERRQ(ierr); 530674ae819SStefano Zampini ierr = PetscFree(send_buffer);CHKERRQ(ierr); 531674ae819SStefano Zampini ierr = PetscFree(cum_recv_counts);CHKERRQ(ierr); 532674ae819SStefano Zampini ierr = PetscFree(subset_to_nodes_indices);CHKERRQ(ierr); 533674ae819SStefano Zampini ierr = PetscFree(subset_cc_adapt);CHKERRQ(ierr); 534674ae819SStefano Zampini /* We are ready to find for connected components consistent among neighbouring subdomains */ 535674ae819SStefano Zampini if (global_subset_counter) { 536df48933bSStefano Zampini ierr = PetscBTMemzero(graph->nvtxs,graph->touched);CHKERRQ(ierr); 537674ae819SStefano Zampini global_subset_counter = 0; 538674ae819SStefano Zampini for (i=0;i<graph->nvtxs;i++) { 539df48933bSStefano Zampini if (graph->subset[i] && !PetscBTLookup(graph->touched,i)) { 540674ae819SStefano Zampini global_subset_counter++; 541674ae819SStefano Zampini for (j=i+1;j<graph->nvtxs;j++) { 542df48933bSStefano Zampini if (!PetscBTLookup(graph->touched,j) && graph->subset[j]==graph->subset[i]) { 543674ae819SStefano Zampini graph->subset[j] = global_subset_counter; 544df48933bSStefano Zampini ierr = PetscBTSet(graph->touched,j);CHKERRQ(ierr); 545674ae819SStefano Zampini } 546674ae819SStefano Zampini } 547674ae819SStefano Zampini graph->subset[i] = global_subset_counter; 548df48933bSStefano Zampini ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr); 549674ae819SStefano Zampini } 550674ae819SStefano Zampini } 551674ae819SStefano Zampini /* refine connected components locally */ 552674ae819SStefano Zampini ierr = PCBDDCGraphComputeConnectedComponentsLocal(graph);CHKERRQ(ierr); 553674ae819SStefano Zampini } 5545b1b9aeaSStefano Zampini /* restore original CSR graph of dofs */ 5555b1b9aeaSStefano Zampini ierr = PetscFree(new_xadj);CHKERRQ(ierr); 5565b1b9aeaSStefano Zampini ierr = PetscFree(new_adjncy);CHKERRQ(ierr); 5575b1b9aeaSStefano Zampini graph->xadj = old_xadj; 5585b1b9aeaSStefano Zampini graph->adjncy = old_adjncy; 559674ae819SStefano Zampini ierr = PetscFree(queue_in_global_numbering);CHKERRQ(ierr); 560674ae819SStefano Zampini } 561674ae819SStefano Zampini PetscFunctionReturn(0); 562674ae819SStefano Zampini } 563674ae819SStefano Zampini 564674ae819SStefano Zampini /* The following code has been adapted from function IsConnectedSubdomain contained 565674ae819SStefano Zampini in source file contig.c of METIS library (version 5.0.1) 566674ae819SStefano Zampini It finds connected components for each subset */ 567674ae819SStefano Zampini #undef __FUNCT__ 568674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphComputeConnectedComponentsLocal" 569674ae819SStefano Zampini PetscErrorCode PCBDDCGraphComputeConnectedComponentsLocal(PCBDDCGraph graph) 570674ae819SStefano Zampini { 571674ae819SStefano Zampini PetscInt i,j,k,first,last,nleft,ncc,pid,cum_queue,n,ncc_pid; 572674ae819SStefano Zampini PetscInt *queue_global; 573674ae819SStefano Zampini PetscErrorCode ierr; 574674ae819SStefano Zampini 575674ae819SStefano Zampini PetscFunctionBegin; 576674ae819SStefano Zampini /* quiet return if no csr info is available */ 577674ae819SStefano Zampini if (!graph->xadj || !graph->adjncy) { 578674ae819SStefano Zampini PetscFunctionReturn(0); 579674ae819SStefano Zampini } 580674ae819SStefano Zampini 581674ae819SStefano Zampini /* reset any previous search of connected components */ 582df48933bSStefano Zampini ierr = PetscBTMemzero(graph->nvtxs,graph->touched);CHKERRQ(ierr); 583674ae819SStefano Zampini graph->n_subsets = 0; 584674ae819SStefano Zampini for (i=0;i<graph->nvtxs;i++) { 5850cf82fd6SStefano Zampini if (graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK || !graph->count[i]) { 586df48933bSStefano Zampini ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr); 587674ae819SStefano Zampini graph->subset[i] = 0; 588674ae819SStefano Zampini } 589674ae819SStefano Zampini graph->n_subsets = PetscMax(graph->n_subsets,graph->subset[i]); 590674ae819SStefano Zampini } 591674ae819SStefano Zampini ierr = PetscFree(graph->subset_ncc);CHKERRQ(ierr); 592785e854fSJed Brown ierr = PetscMalloc1(graph->n_subsets,&graph->subset_ncc);CHKERRQ(ierr); 593674ae819SStefano Zampini ierr = PetscMemzero(graph->subset_ncc,graph->n_subsets*sizeof(*graph->subset_ncc));CHKERRQ(ierr); 594674ae819SStefano Zampini ierr = PetscMemzero(graph->cptr,(graph->nvtxs+1)*sizeof(*graph->cptr));CHKERRQ(ierr); 595674ae819SStefano Zampini ierr = PetscMemzero(graph->queue,graph->nvtxs*sizeof(*graph->queue));CHKERRQ(ierr); 596674ae819SStefano Zampini 597674ae819SStefano Zampini /* begin search for connected components */ 598674ae819SStefano Zampini cum_queue = 0; 599674ae819SStefano Zampini ncc = 0; 600674ae819SStefano Zampini for (n=0;n<graph->n_subsets;n++) { 601674ae819SStefano Zampini pid = n+1; /* partition labeled by 0 is discarded */ 602674ae819SStefano Zampini nleft = 0; 603674ae819SStefano Zampini for (i=0;i<graph->nvtxs;i++) { 604674ae819SStefano Zampini if (graph->subset[i] == pid) { 605674ae819SStefano Zampini nleft++; 606674ae819SStefano Zampini } 607674ae819SStefano Zampini } 608674ae819SStefano Zampini for (i=0; i<graph->nvtxs; i++) { 609674ae819SStefano Zampini if (graph->subset[i] == pid) { 610674ae819SStefano Zampini break; 611674ae819SStefano Zampini } 612674ae819SStefano Zampini } 613df48933bSStefano Zampini ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr); 614674ae819SStefano Zampini graph->queue[cum_queue] = i; 615674ae819SStefano Zampini first = 0; 616674ae819SStefano Zampini last = 1; 617674ae819SStefano Zampini graph->cptr[ncc] = cum_queue; 618674ae819SStefano Zampini ncc_pid = 0; 619674ae819SStefano Zampini while (first != nleft) { 620674ae819SStefano Zampini if (first == last) { 621674ae819SStefano Zampini graph->cptr[++ncc] = first+cum_queue; 622674ae819SStefano Zampini ncc_pid++; 623df48933bSStefano Zampini for (i=0; i<graph->nvtxs; i++) { /* TODO-> use a while! */ 624df48933bSStefano Zampini if (graph->subset[i] == pid && !PetscBTLookup(graph->touched,i)) { 625674ae819SStefano Zampini break; 626674ae819SStefano Zampini } 627674ae819SStefano Zampini } 628674ae819SStefano Zampini graph->queue[cum_queue+last] = i; 629674ae819SStefano Zampini last++; 630df48933bSStefano Zampini ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr); 631674ae819SStefano Zampini } 632674ae819SStefano Zampini i = graph->queue[cum_queue+first]; 633674ae819SStefano Zampini first++; 634674ae819SStefano Zampini for (j=graph->xadj[i];j<graph->xadj[i+1];j++) { 635674ae819SStefano Zampini k = graph->adjncy[j]; 636df48933bSStefano Zampini if (graph->subset[k] == pid && !PetscBTLookup(graph->touched,k)) { 637674ae819SStefano Zampini graph->queue[cum_queue+last] = k; 638674ae819SStefano Zampini last++; 639df48933bSStefano Zampini ierr = PetscBTSet(graph->touched,k);CHKERRQ(ierr); 640674ae819SStefano Zampini } 641674ae819SStefano Zampini } 642674ae819SStefano Zampini } 643674ae819SStefano Zampini graph->cptr[++ncc] = first+cum_queue; 644674ae819SStefano Zampini ncc_pid++; 645674ae819SStefano Zampini cum_queue = graph->cptr[ncc]; 646674ae819SStefano Zampini graph->subset_ncc[n] = ncc_pid; 647674ae819SStefano Zampini } 648674ae819SStefano Zampini graph->ncc = ncc; 649674ae819SStefano Zampini /* For consistency among neighbours, I need to sort (by global ordering) each connected component */ 650785e854fSJed Brown ierr = PetscMalloc1(graph->cptr[graph->ncc],&queue_global);CHKERRQ(ierr); 651674ae819SStefano Zampini ierr = ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_global);CHKERRQ(ierr); 652674ae819SStefano Zampini for (i=0;i<graph->ncc;i++) { 653674ae819SStefano Zampini ierr = PetscSortIntWithArray(graph->cptr[i+1]-graph->cptr[i],&queue_global[graph->cptr[i]],&graph->queue[graph->cptr[i]]);CHKERRQ(ierr); 654674ae819SStefano Zampini } 655674ae819SStefano Zampini ierr = PetscFree(queue_global);CHKERRQ(ierr); 656674ae819SStefano Zampini PetscFunctionReturn(0); 657674ae819SStefano Zampini } 658674ae819SStefano Zampini 659674ae819SStefano Zampini #undef __FUNCT__ 660674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphSetUp" 661674ae819SStefano Zampini PetscErrorCode PCBDDCGraphSetUp(PCBDDCGraph graph, PetscInt custom_minimal_size, IS neumann_is, IS dirichlet_is, PetscInt n_ISForDofs, IS ISForDofs[], IS custom_primal_vertices) 662674ae819SStefano Zampini { 663674ae819SStefano Zampini VecScatter scatter_ctx; 664674ae819SStefano Zampini Vec local_vec,local_vec2,global_vec; 665674ae819SStefano Zampini IS to,from; 666674ae819SStefano Zampini MPI_Comm comm; 667674ae819SStefano Zampini PetscScalar *array,*array2; 668674ae819SStefano Zampini const PetscInt *is_indices; 6693a5b03d1SStefano Zampini PetscInt n_neigh,*neigh,*n_shared,**shared,*queue_global; 670674ae819SStefano Zampini PetscInt i,j,k,s,total_counts,nodes_touched,is_size; 671674ae819SStefano Zampini PetscErrorCode ierr; 67251b0f6cfSStefano Zampini PetscBool same_set,mirrors_found; 673674ae819SStefano Zampini 674674ae819SStefano Zampini PetscFunctionBegin; 675674ae819SStefano Zampini ierr = PetscObjectGetComm((PetscObject)(graph->l2gmap),&comm);CHKERRQ(ierr); 676674ae819SStefano Zampini /* custom_minimal_size */ 677674ae819SStefano Zampini graph->custom_minimal_size = PetscMax(graph->custom_minimal_size,custom_minimal_size); 678674ae819SStefano Zampini /* get info l2gmap and allocate work vectors */ 679674ae819SStefano Zampini ierr = ISLocalToGlobalMappingGetInfo(graph->l2gmap,&n_neigh,&neigh,&n_shared,&shared);CHKERRQ(ierr); 680674ae819SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(graph->l2gmap,&is_indices);CHKERRQ(ierr); 681674ae819SStefano Zampini j = 0; 682674ae819SStefano Zampini for (i=0;i<graph->nvtxs;i++) { 683674ae819SStefano Zampini j = PetscMax(j,is_indices[i]); 684674ae819SStefano Zampini } 685674ae819SStefano Zampini ierr = MPI_Allreduce(&j,&i,1,MPIU_INT,MPI_MAX,comm);CHKERRQ(ierr); 686674ae819SStefano Zampini i++; 687674ae819SStefano Zampini ierr = VecCreate(PETSC_COMM_SELF,&local_vec);CHKERRQ(ierr); 688674ae819SStefano Zampini ierr = VecSetSizes(local_vec,PETSC_DECIDE,graph->nvtxs);CHKERRQ(ierr); 689674ae819SStefano Zampini ierr = VecSetType(local_vec,VECSTANDARD);CHKERRQ(ierr); 690674ae819SStefano Zampini ierr = VecDuplicate(local_vec,&local_vec2);CHKERRQ(ierr); 691674ae819SStefano Zampini ierr = VecCreate(comm,&global_vec);CHKERRQ(ierr); 692674ae819SStefano Zampini ierr = VecSetSizes(global_vec,PETSC_DECIDE,i);CHKERRQ(ierr); 693674ae819SStefano Zampini ierr = VecSetType(global_vec,VECSTANDARD);CHKERRQ(ierr); 694674ae819SStefano Zampini ierr = ISCreateStride(PETSC_COMM_SELF,graph->nvtxs,0,1,&to);CHKERRQ(ierr); 695674ae819SStefano Zampini ierr = ISLocalToGlobalMappingApplyIS(graph->l2gmap,to,&from);CHKERRQ(ierr); 696674ae819SStefano Zampini ierr = VecScatterCreate(global_vec,from,local_vec,to,&scatter_ctx);CHKERRQ(ierr); 69751b0f6cfSStefano Zampini 69851b0f6cfSStefano Zampini /* get local periodic nodes */ 69951b0f6cfSStefano Zampini mirrors_found = PETSC_FALSE; 7009881197aSStefano Zampini if (graph->nvtxs && n_neigh) { 70149eeff8cSStefano Zampini for (i=0; i<n_shared[0]; i++) graph->count[shared[0][i]] += 1; 70251b0f6cfSStefano Zampini for (i=0; i<n_shared[0]; i++) { 70351b0f6cfSStefano Zampini if (graph->count[shared[0][i]] > 1) { 70451b0f6cfSStefano Zampini mirrors_found = PETSC_TRUE; 70551b0f6cfSStefano Zampini break; 70651b0f6cfSStefano Zampini } 70751b0f6cfSStefano Zampini } 70849eeff8cSStefano Zampini } 70951b0f6cfSStefano Zampini /* compute local mirrors (if any) */ 71051b0f6cfSStefano Zampini if (mirrors_found) { 71151b0f6cfSStefano Zampini PetscInt *local_indices,*global_indices; 71251b0f6cfSStefano Zampini /* get arrays of local and global indices */ 713785e854fSJed Brown ierr = PetscMalloc1(graph->nvtxs,&local_indices);CHKERRQ(ierr); 71451b0f6cfSStefano Zampini ierr = ISGetIndices(to,(const PetscInt**)&is_indices);CHKERRQ(ierr); 71551b0f6cfSStefano Zampini ierr = PetscMemcpy(local_indices,is_indices,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr); 71651b0f6cfSStefano Zampini ierr = ISRestoreIndices(to,(const PetscInt**)&is_indices);CHKERRQ(ierr); 717785e854fSJed Brown ierr = PetscMalloc1(graph->nvtxs,&global_indices);CHKERRQ(ierr); 71851b0f6cfSStefano Zampini ierr = ISGetIndices(from,(const PetscInt**)&is_indices);CHKERRQ(ierr); 71951b0f6cfSStefano Zampini ierr = PetscMemcpy(global_indices,is_indices,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr); 72051b0f6cfSStefano Zampini ierr = ISRestoreIndices(from,(const PetscInt**)&is_indices);CHKERRQ(ierr); 72151b0f6cfSStefano Zampini /* allocate space for mirrors */ 7228e5aaad5SJed Brown ierr = PetscMalloc2(graph->nvtxs,&graph->mirrors,graph->nvtxs,&graph->mirrors_set);CHKERRQ(ierr); 72351b0f6cfSStefano Zampini ierr = PetscMemzero(graph->mirrors,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr); 72451b0f6cfSStefano Zampini graph->mirrors_set[0] = 0; 72551b0f6cfSStefano Zampini 72651b0f6cfSStefano Zampini k=0; 72751b0f6cfSStefano Zampini for (i=0;i<n_shared[0];i++) { 72851b0f6cfSStefano Zampini j=shared[0][i]; 72951b0f6cfSStefano Zampini if (graph->count[j] > 1) { 73051b0f6cfSStefano Zampini graph->mirrors[j]++; 73151b0f6cfSStefano Zampini k++; 73251b0f6cfSStefano Zampini } 73351b0f6cfSStefano Zampini } 73451b0f6cfSStefano Zampini /* allocate space for set of mirrors */ 735785e854fSJed Brown ierr = PetscMalloc1(k,&graph->mirrors_set[0]);CHKERRQ(ierr); 73651b0f6cfSStefano Zampini for (i=1;i<graph->nvtxs;i++) 73751b0f6cfSStefano Zampini graph->mirrors_set[i]=graph->mirrors_set[i-1]+graph->mirrors[i-1]; 73851b0f6cfSStefano Zampini 73951b0f6cfSStefano Zampini /* fill arrays */ 74051b0f6cfSStefano Zampini ierr = PetscMemzero(graph->mirrors,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr); 74151b0f6cfSStefano Zampini for (j=0;j<n_shared[0];j++) { 74251b0f6cfSStefano Zampini i=shared[0][j]; 74351b0f6cfSStefano Zampini if (graph->count[i] > 1) 74451b0f6cfSStefano Zampini graph->mirrors_set[i][graph->mirrors[i]++]=global_indices[i]; 74551b0f6cfSStefano Zampini } 74651b0f6cfSStefano Zampini ierr = PetscSortIntWithArray(graph->nvtxs,global_indices,local_indices);CHKERRQ(ierr); 74751b0f6cfSStefano Zampini for (i=0;i<graph->nvtxs;i++) { 74851b0f6cfSStefano Zampini if (graph->mirrors[i] > 0) { 74951b0f6cfSStefano Zampini ierr = PetscFindInt(graph->mirrors_set[i][0],graph->nvtxs,global_indices,&k);CHKERRQ(ierr); 75051b0f6cfSStefano Zampini j = global_indices[k]; 75151b0f6cfSStefano Zampini while ( k > 0 && global_indices[k-1] == j) k--; 75251b0f6cfSStefano Zampini for (j=0;j<graph->mirrors[i];j++) { 75351b0f6cfSStefano Zampini graph->mirrors_set[i][j]=local_indices[k+j]; 75451b0f6cfSStefano Zampini } 75551b0f6cfSStefano Zampini ierr = PetscSortInt(graph->mirrors[i],graph->mirrors_set[i]);CHKERRQ(ierr); 75651b0f6cfSStefano Zampini } 75751b0f6cfSStefano Zampini } 75851b0f6cfSStefano Zampini ierr = PetscFree(local_indices);CHKERRQ(ierr); 75951b0f6cfSStefano Zampini ierr = PetscFree(global_indices);CHKERRQ(ierr); 76051b0f6cfSStefano Zampini } 76151b0f6cfSStefano Zampini ierr = PetscMemzero(graph->count,graph->nvtxs*sizeof(*graph->count));CHKERRQ(ierr); 762674ae819SStefano Zampini ierr = ISDestroy(&to);CHKERRQ(ierr); 763674ae819SStefano Zampini ierr = ISDestroy(&from);CHKERRQ(ierr); 76451b0f6cfSStefano Zampini 765674ae819SStefano Zampini /* Count total number of neigh per node */ 766674ae819SStefano Zampini k=0; 767674ae819SStefano Zampini for (i=1;i<n_neigh;i++) { 768674ae819SStefano Zampini k += n_shared[i]; 769674ae819SStefano Zampini for (j=0;j<n_shared[i];j++) { 770674ae819SStefano Zampini graph->count[shared[i][j]] += 1; 771674ae819SStefano Zampini } 772674ae819SStefano Zampini } 773674ae819SStefano Zampini /* Allocate space for storing the set of neighbours for each node */ 774674ae819SStefano Zampini if (graph->nvtxs) { 775785e854fSJed Brown ierr = PetscMalloc1(k,&graph->neighbours_set[0]);CHKERRQ(ierr); 776674ae819SStefano Zampini } 777674ae819SStefano Zampini for (i=1;i<graph->nvtxs;i++) { /* dont count myself */ 778674ae819SStefano Zampini graph->neighbours_set[i]=graph->neighbours_set[i-1]+graph->count[i-1]; 779674ae819SStefano Zampini } 780674ae819SStefano Zampini /* Get information for sharing subdomains */ 781674ae819SStefano Zampini ierr = PetscMemzero(graph->count,graph->nvtxs*sizeof(*graph->count));CHKERRQ(ierr); 782674ae819SStefano Zampini for (i=1;i<n_neigh;i++) { /* dont count myself */ 783674ae819SStefano Zampini s = n_shared[i]; 784674ae819SStefano Zampini for (j=0;j<s;j++) { 785674ae819SStefano Zampini k = shared[i][j]; 786674ae819SStefano Zampini graph->neighbours_set[k][graph->count[k]] = neigh[i]; 787674ae819SStefano Zampini graph->count[k] += 1; 788674ae819SStefano Zampini } 789674ae819SStefano Zampini } 790674ae819SStefano Zampini /* sort set of sharing subdomains */ 791674ae819SStefano Zampini for (i=0;i<graph->nvtxs;i++) { 79293fb5fd3SStefano Zampini ierr = PetscSortRemoveDupsInt(&graph->count[i],graph->neighbours_set[i]);CHKERRQ(ierr); 793674ae819SStefano Zampini } 79467c9da69SStefano Zampini /* 79567c9da69SStefano Zampini Get info for dofs splitting 79667c9da69SStefano Zampini User can specify only a subset; an additional field is considered as a complementary set 79767c9da69SStefano Zampini */ 79867c9da69SStefano Zampini for (i=0;i<graph->nvtxs;i++) { 79967c9da69SStefano Zampini graph->which_dof[i] = n_ISForDofs; /* by default a dof belongs to the complement set */ 80067c9da69SStefano Zampini } 801674ae819SStefano Zampini for (i=0;i<n_ISForDofs;i++) { 80263602bcaSStefano Zampini ierr = ISGetLocalSize(ISForDofs[i],&is_size);CHKERRQ(ierr); 803674ae819SStefano Zampini ierr = ISGetIndices(ISForDofs[i],(const PetscInt**)&is_indices);CHKERRQ(ierr); 804674ae819SStefano Zampini for (j=0;j<is_size;j++) { 805d50ed60aSStefano Zampini if (is_indices[j] > -1 && is_indices[j] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */ 806674ae819SStefano Zampini graph->which_dof[is_indices[j]] = i; 807674ae819SStefano Zampini } 80867c9da69SStefano Zampini } 809674ae819SStefano Zampini ierr = ISRestoreIndices(ISForDofs[i],(const PetscInt**)&is_indices);CHKERRQ(ierr); 810674ae819SStefano Zampini } 811674ae819SStefano Zampini /* Take into account Neumann nodes */ 812674ae819SStefano Zampini ierr = VecSet(local_vec,0.0);CHKERRQ(ierr); 813674ae819SStefano Zampini ierr = VecSet(local_vec2,0.0);CHKERRQ(ierr); 814674ae819SStefano Zampini if (neumann_is) { 815674ae819SStefano Zampini ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr); 81682ba6b80SStefano Zampini ierr = ISGetLocalSize(neumann_is,&is_size);CHKERRQ(ierr); 817674ae819SStefano Zampini ierr = ISGetIndices(neumann_is,(const PetscInt**)&is_indices);CHKERRQ(ierr); 818674ae819SStefano Zampini for (i=0;i<is_size;i++) { 819d50ed60aSStefano Zampini if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */ 820674ae819SStefano Zampini array[is_indices[i]] = 1.0; 821674ae819SStefano Zampini } 8223c629590SStefano Zampini } 823674ae819SStefano Zampini ierr = ISRestoreIndices(neumann_is,(const PetscInt**)&is_indices);CHKERRQ(ierr); 824674ae819SStefano Zampini ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr); 825674ae819SStefano Zampini } 826674ae819SStefano Zampini /* Neumann nodes: impose consistency among neighbours */ 827674ae819SStefano Zampini ierr = VecSet(global_vec,0.0);CHKERRQ(ierr); 828674ae819SStefano Zampini ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 829674ae819SStefano Zampini ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 830674ae819SStefano Zampini ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 831674ae819SStefano Zampini ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 832674ae819SStefano Zampini ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr); 833674ae819SStefano Zampini for (i=0;i<graph->nvtxs;i++) { 8343c629590SStefano Zampini if (PetscRealPart(array[i]) > 0.1) { 8350cf82fd6SStefano Zampini graph->special_dof[i] = PCBDDCGRAPH_NEUMANN_MARK; 836674ae819SStefano Zampini } 837674ae819SStefano Zampini } 838674ae819SStefano Zampini ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr); 839674ae819SStefano Zampini /* Take into account Dirichlet nodes */ 840674ae819SStefano Zampini ierr = VecSet(local_vec2,0.0);CHKERRQ(ierr); 841674ae819SStefano Zampini if (dirichlet_is) { 842674ae819SStefano Zampini ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr); 843674ae819SStefano Zampini ierr = VecGetArray(local_vec2,&array2);CHKERRQ(ierr); 84482ba6b80SStefano Zampini ierr = ISGetLocalSize(dirichlet_is,&is_size);CHKERRQ(ierr); 845674ae819SStefano Zampini ierr = ISGetIndices(dirichlet_is,(const PetscInt**)&is_indices);CHKERRQ(ierr); 846674ae819SStefano Zampini for (i=0;i<is_size;i++){ 847d50ed60aSStefano Zampini if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */ 848674ae819SStefano Zampini k = is_indices[i]; 849df48933bSStefano Zampini if (graph->count[k] && !PetscBTLookup(graph->touched,k)) { 8503c629590SStefano Zampini if (PetscRealPart(array[k]) > 0.1) { 8513c629590SStefano 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); 852674ae819SStefano Zampini } 853674ae819SStefano Zampini array2[k] = 1.0; 854674ae819SStefano Zampini } 855674ae819SStefano Zampini } 8563c629590SStefano Zampini } 857674ae819SStefano Zampini ierr = ISRestoreIndices(dirichlet_is,(const PetscInt**)&is_indices);CHKERRQ(ierr); 858674ae819SStefano Zampini ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr); 859674ae819SStefano Zampini ierr = VecRestoreArray(local_vec2,&array2);CHKERRQ(ierr); 860674ae819SStefano Zampini } 861674ae819SStefano Zampini /* Dirichlet nodes: impose consistency among neighbours */ 862674ae819SStefano Zampini ierr = VecSet(global_vec,0.0);CHKERRQ(ierr); 863674ae819SStefano Zampini ierr = VecScatterBegin(scatter_ctx,local_vec2,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 864674ae819SStefano Zampini ierr = VecScatterEnd(scatter_ctx,local_vec2,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 865674ae819SStefano Zampini ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 866674ae819SStefano Zampini ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 867674ae819SStefano Zampini ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr); 868674ae819SStefano Zampini for (i=0;i<graph->nvtxs;i++) { 8693c629590SStefano Zampini if (PetscRealPart(array[i]) > 0.1) { 870df48933bSStefano Zampini ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr); 871674ae819SStefano Zampini graph->subset[i] = 0; /* dirichlet nodes treated as internal -> is it ok? */ 8720cf82fd6SStefano Zampini graph->special_dof[i] = PCBDDCGRAPH_DIRICHLET_MARK; 873674ae819SStefano Zampini } 874674ae819SStefano Zampini } 875674ae819SStefano Zampini ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr); 87651b0f6cfSStefano Zampini 87708a5cf49SStefano Zampini /* mark local periodic nodes (if any) and adapt CSR graph (if any) */ 87851b0f6cfSStefano Zampini if (graph->mirrors) { 87951b0f6cfSStefano Zampini for (i=0;i<graph->nvtxs;i++) 88051b0f6cfSStefano Zampini if (graph->mirrors[i]) 8810cf82fd6SStefano Zampini graph->special_dof[i] = PCBDDCGRAPH_LOCAL_PERIODIC_MARK; 88251b0f6cfSStefano Zampini 88308a5cf49SStefano Zampini if (graph->xadj && graph->adjncy) { 88408a5cf49SStefano Zampini PetscInt *new_xadj,*new_adjncy; 88551b0f6cfSStefano Zampini /* sort CSR graph */ 88651b0f6cfSStefano Zampini for (i=0;i<graph->nvtxs;i++) 88751b0f6cfSStefano Zampini ierr = PetscSortInt(graph->xadj[i+1]-graph->xadj[i],&graph->adjncy[graph->xadj[i]]);CHKERRQ(ierr); 88851b0f6cfSStefano Zampini 88951b0f6cfSStefano Zampini /* adapt local CSR graph in case of local periodicity */ 89051b0f6cfSStefano Zampini k=0; 89151b0f6cfSStefano Zampini for (i=0;i<graph->nvtxs;i++) 89251b0f6cfSStefano Zampini for (j=graph->xadj[i];j<graph->xadj[i+1];j++) 89351b0f6cfSStefano Zampini k += graph->mirrors[graph->adjncy[j]]; 89451b0f6cfSStefano Zampini 895785e854fSJed Brown ierr = PetscMalloc1((graph->nvtxs+1),&new_xadj);CHKERRQ(ierr); 896785e854fSJed Brown ierr = PetscMalloc1((k+graph->xadj[graph->nvtxs]),&new_adjncy);CHKERRQ(ierr); 89751b0f6cfSStefano Zampini new_xadj[0]=0; 89851b0f6cfSStefano Zampini for (i=0;i<graph->nvtxs;i++) { 89951b0f6cfSStefano Zampini k = graph->xadj[i+1]-graph->xadj[i]; 90051b0f6cfSStefano Zampini ierr = PetscMemcpy(&new_adjncy[new_xadj[i]],&graph->adjncy[graph->xadj[i]],k*sizeof(PetscInt));CHKERRQ(ierr); 90151b0f6cfSStefano Zampini new_xadj[i+1]=new_xadj[i]+k; 90251b0f6cfSStefano Zampini for (j=graph->xadj[i];j<graph->xadj[i+1];j++) { 90351b0f6cfSStefano Zampini k = graph->mirrors[graph->adjncy[j]]; 90451b0f6cfSStefano Zampini ierr = PetscMemcpy(&new_adjncy[new_xadj[i+1]],graph->mirrors_set[graph->adjncy[j]],k*sizeof(PetscInt));CHKERRQ(ierr); 90551b0f6cfSStefano Zampini new_xadj[i+1]+=k; 90651b0f6cfSStefano Zampini } 90751b0f6cfSStefano Zampini k = new_xadj[i+1]-new_xadj[i]; 90851b0f6cfSStefano Zampini ierr = PetscSortRemoveDupsInt(&k,&new_adjncy[new_xadj[i]]);CHKERRQ(ierr); 90951b0f6cfSStefano Zampini new_xadj[i+1]=new_xadj[i]+k; 91051b0f6cfSStefano Zampini } 91151b0f6cfSStefano Zampini /* set new CSR into graph */ 91251b0f6cfSStefano Zampini ierr = PetscFree(graph->xadj);CHKERRQ(ierr); 91351b0f6cfSStefano Zampini ierr = PetscFree(graph->adjncy);CHKERRQ(ierr); 91451b0f6cfSStefano Zampini graph->xadj = new_xadj; 91551b0f6cfSStefano Zampini graph->adjncy = new_adjncy; 91651b0f6cfSStefano Zampini } 91708a5cf49SStefano Zampini } 91851b0f6cfSStefano Zampini 919674ae819SStefano Zampini /* mark special nodes -> each will become a single node equivalence class */ 9209b70a373SStefano Zampini ierr = VecSet(local_vec,0.0);CHKERRQ(ierr); 921674ae819SStefano Zampini if (custom_primal_vertices) { 9229b70a373SStefano Zampini ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr); 92363602bcaSStefano Zampini ierr = ISGetLocalSize(custom_primal_vertices,&is_size);CHKERRQ(ierr); 924674ae819SStefano Zampini ierr = ISGetIndices(custom_primal_vertices,(const PetscInt**)&is_indices);CHKERRQ(ierr); 925674ae819SStefano Zampini for (i=0;i<is_size;i++){ 9269b70a373SStefano Zampini if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */ 9279b70a373SStefano Zampini array[is_indices[i]] = 1.0; 9289b70a373SStefano Zampini } 929674ae819SStefano Zampini } 930674ae819SStefano Zampini ierr = ISRestoreIndices(custom_primal_vertices,(const PetscInt**)&is_indices);CHKERRQ(ierr); 9319b70a373SStefano Zampini ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr); 932674ae819SStefano Zampini } 9339b70a373SStefano Zampini /* special nodes: impose consistency among neighbours */ 9349b70a373SStefano Zampini ierr = VecSet(global_vec,0.0);CHKERRQ(ierr); 9359b70a373SStefano Zampini ierr = VecScatterBegin(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 9369b70a373SStefano Zampini ierr = VecScatterEnd(scatter_ctx,local_vec,global_vec,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 9379b70a373SStefano Zampini ierr = VecScatterBegin(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9389b70a373SStefano Zampini ierr = VecScatterEnd(scatter_ctx,global_vec,local_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9399b70a373SStefano Zampini ierr = VecGetArray(local_vec,&array);CHKERRQ(ierr); 9409b70a373SStefano Zampini j = 0; 9419b70a373SStefano Zampini for (i=0;i<graph->nvtxs;i++) { 9429b70a373SStefano Zampini if (PetscRealPart(array[i]) > 0.1 && graph->special_dof[i] != PCBDDCGRAPH_DIRICHLET_MARK) { 9439b70a373SStefano Zampini graph->special_dof[i] = PCBDDCGRAPH_SPECIAL_MARK-j; 9449b70a373SStefano Zampini j++; 9459b70a373SStefano Zampini } 9469b70a373SStefano Zampini } 9479b70a373SStefano Zampini ierr = VecRestoreArray(local_vec,&array);CHKERRQ(ierr); 9489b70a373SStefano Zampini 949674ae819SStefano Zampini /* mark interior nodes as touched and belonging to partition number 0 */ 950674ae819SStefano Zampini for (i=0;i<graph->nvtxs;i++) { 951674ae819SStefano Zampini if (!graph->count[i]) { 952df48933bSStefano Zampini ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr); 953674ae819SStefano Zampini graph->subset[i] = 0; 954674ae819SStefano Zampini } 955674ae819SStefano Zampini } 956674ae819SStefano Zampini /* init graph structure and compute default subsets */ 957674ae819SStefano Zampini nodes_touched=0; 958674ae819SStefano Zampini for (i=0;i<graph->nvtxs;i++) { 959df48933bSStefano Zampini if (PetscBTLookup(graph->touched,i)) { 960674ae819SStefano Zampini nodes_touched++; 961674ae819SStefano Zampini } 962674ae819SStefano Zampini } 963674ae819SStefano Zampini i = 0; 964674ae819SStefano Zampini graph->ncc = 0; 965674ae819SStefano Zampini total_counts = 0; 966674ae819SStefano Zampini while (nodes_touched<graph->nvtxs) { 967674ae819SStefano Zampini /* find first untouched node in local ordering */ 968df48933bSStefano Zampini while (PetscBTLookup(graph->touched,i)) { 969674ae819SStefano Zampini i++; 970674ae819SStefano Zampini } 971df48933bSStefano Zampini ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr); 972674ae819SStefano Zampini graph->subset[i] = graph->ncc+1; 973674ae819SStefano Zampini graph->cptr[graph->ncc] = total_counts; 974674ae819SStefano Zampini graph->queue[total_counts] = i; 975674ae819SStefano Zampini total_counts++; 976674ae819SStefano Zampini nodes_touched++; 977674ae819SStefano Zampini /* now find all other nodes having the same set of sharing subdomains */ 978674ae819SStefano Zampini for (j=i+1;j<graph->nvtxs;j++) { 97974e413f5SStefano Zampini /* check for same number of sharing subdomains, dof number and same special mark */ 980df48933bSStefano 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]) { 981674ae819SStefano Zampini /* check for same set of sharing subdomains */ 982674ae819SStefano Zampini same_set=PETSC_TRUE; 983674ae819SStefano Zampini for (k=0;k<graph->count[j];k++){ 984674ae819SStefano Zampini if (graph->neighbours_set[i][k]!=graph->neighbours_set[j][k]) { 985674ae819SStefano Zampini same_set=PETSC_FALSE; 986674ae819SStefano Zampini } 987674ae819SStefano Zampini } 988674ae819SStefano Zampini /* I found a friend of mine */ 989674ae819SStefano Zampini if (same_set) { 990674ae819SStefano Zampini graph->subset[j]=graph->ncc+1; 991df48933bSStefano Zampini ierr = PetscBTSet(graph->touched,j);CHKERRQ(ierr); 992674ae819SStefano Zampini nodes_touched++; 993674ae819SStefano Zampini graph->queue[total_counts] = j; 994674ae819SStefano Zampini total_counts++; 995674ae819SStefano Zampini } 996674ae819SStefano Zampini } 997674ae819SStefano Zampini } 998674ae819SStefano Zampini graph->ncc++; 999674ae819SStefano Zampini } 1000674ae819SStefano Zampini /* set default number of subsets (at this point no info on csr graph has been taken into account, so n_subsets = ncc */ 1001674ae819SStefano Zampini graph->n_subsets = graph->ncc; 1002785e854fSJed Brown ierr = PetscMalloc1(graph->n_subsets,&graph->subset_ncc);CHKERRQ(ierr); 1003674ae819SStefano Zampini for (i=0;i<graph->n_subsets;i++) { 1004674ae819SStefano Zampini graph->subset_ncc[i] = 1; 1005674ae819SStefano Zampini } 1006674ae819SStefano Zampini /* final pointer */ 1007674ae819SStefano Zampini graph->cptr[graph->ncc] = total_counts; 1008674ae819SStefano Zampini /* free memory allocated by ISLocalToGlobalMappingGetInfo */ 1009674ae819SStefano Zampini ierr = ISLocalToGlobalMappingRestoreInfo(graph->l2gmap,&n_neigh,&neigh,&n_shared,&shared);CHKERRQ(ierr); 10103a5b03d1SStefano Zampini /* get a reference node (min index in global ordering) for each subset */ 10113a5b03d1SStefano Zampini ierr = PetscMalloc1(graph->ncc,&graph->subset_ref_node);CHKERRQ(ierr); 10123a5b03d1SStefano Zampini ierr = PetscMalloc1(graph->cptr[graph->ncc],&queue_global);CHKERRQ(ierr); 10133a5b03d1SStefano Zampini ierr = ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_global);CHKERRQ(ierr); 10143a5b03d1SStefano Zampini for (i=0;i<graph->ncc;i++) { 10153a5b03d1SStefano Zampini PetscInt minval = queue_global[graph->cptr[i]]; 10163a5b03d1SStefano Zampini for (j=graph->cptr[i]+1;j<graph->cptr[i+1];j++) { 10173a5b03d1SStefano Zampini minval = PetscMin(minval,queue_global[j]); 10183a5b03d1SStefano Zampini } 10193a5b03d1SStefano Zampini graph->subset_ref_node[i] = minval; 10203a5b03d1SStefano Zampini } 10213a5b03d1SStefano Zampini ierr = PetscFree(queue_global);CHKERRQ(ierr); 1022674ae819SStefano Zampini /* free objects */ 1023674ae819SStefano Zampini ierr = VecDestroy(&local_vec);CHKERRQ(ierr); 1024674ae819SStefano Zampini ierr = VecDestroy(&local_vec2);CHKERRQ(ierr); 1025674ae819SStefano Zampini ierr = VecDestroy(&global_vec);CHKERRQ(ierr); 1026674ae819SStefano Zampini ierr = VecScatterDestroy(&scatter_ctx);CHKERRQ(ierr); 1027674ae819SStefano Zampini PetscFunctionReturn(0); 1028674ae819SStefano Zampini } 1029674ae819SStefano Zampini 1030674ae819SStefano Zampini #undef __FUNCT__ 1031674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphResetCSR" 1032674ae819SStefano Zampini PetscErrorCode PCBDDCGraphResetCSR(PCBDDCGraph graph) 1033674ae819SStefano Zampini { 1034674ae819SStefano Zampini PetscErrorCode ierr; 1035674ae819SStefano Zampini 1036674ae819SStefano Zampini PetscFunctionBegin; 1037674ae819SStefano Zampini ierr = PetscFree(graph->xadj);CHKERRQ(ierr); 1038674ae819SStefano Zampini ierr = PetscFree(graph->adjncy);CHKERRQ(ierr); 1039575ad6abSStefano Zampini graph->nvtxs_csr = 0; 1040674ae819SStefano Zampini PetscFunctionReturn(0); 1041674ae819SStefano Zampini } 1042674ae819SStefano Zampini 1043674ae819SStefano Zampini #undef __FUNCT__ 1044674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphReset" 1045674ae819SStefano Zampini PetscErrorCode PCBDDCGraphReset(PCBDDCGraph graph) 1046674ae819SStefano Zampini { 1047674ae819SStefano Zampini PetscErrorCode ierr; 1048674ae819SStefano Zampini 1049674ae819SStefano Zampini PetscFunctionBegin; 1050674ae819SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&graph->l2gmap);CHKERRQ(ierr); 1051674ae819SStefano Zampini ierr = PetscFree(graph->subset_ncc);CHKERRQ(ierr); 10523a5b03d1SStefano Zampini ierr = PetscFree(graph->subset_ref_node);CHKERRQ(ierr); 1053674ae819SStefano Zampini if (graph->nvtxs) { 1054674ae819SStefano Zampini ierr = PetscFree(graph->neighbours_set[0]);CHKERRQ(ierr); 1055674ae819SStefano Zampini } 1056df48933bSStefano Zampini ierr = PetscBTDestroy(&graph->touched);CHKERRQ(ierr); 1057df48933bSStefano Zampini ierr = PetscFree7(graph->count, 1058674ae819SStefano Zampini graph->neighbours_set, 1059674ae819SStefano Zampini graph->subset, 1060674ae819SStefano Zampini graph->which_dof, 1061674ae819SStefano Zampini graph->cptr, 1062df48933bSStefano Zampini graph->queue, 1063df48933bSStefano Zampini graph->special_dof);CHKERRQ(ierr); 106451b0f6cfSStefano Zampini if (graph->mirrors) { 106551b0f6cfSStefano Zampini ierr = PetscFree(graph->mirrors_set[0]);CHKERRQ(ierr); 106651b0f6cfSStefano Zampini } 106751b0f6cfSStefano Zampini ierr = PetscFree2(graph->mirrors,graph->mirrors_set);CHKERRQ(ierr); 1068674ae819SStefano Zampini graph->nvtxs = 0; 1069674ae819SStefano Zampini graph->n_subsets = 0; 1070674ae819SStefano Zampini graph->custom_minimal_size = 1; 1071674ae819SStefano Zampini PetscFunctionReturn(0); 1072674ae819SStefano Zampini } 1073674ae819SStefano Zampini 1074674ae819SStefano Zampini #undef __FUNCT__ 1075674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphInit" 1076674ae819SStefano Zampini PetscErrorCode PCBDDCGraphInit(PCBDDCGraph graph, ISLocalToGlobalMapping l2gmap) 1077674ae819SStefano Zampini { 1078674ae819SStefano Zampini PetscInt n; 1079674ae819SStefano Zampini PetscErrorCode ierr; 1080674ae819SStefano Zampini 1081674ae819SStefano Zampini PetscFunctionBegin; 1082674ae819SStefano Zampini PetscValidPointer(graph,1); 1083674ae819SStefano Zampini PetscValidHeaderSpecific(l2gmap,IS_LTOGM_CLASSID,2); 10848e61c736SStefano Zampini /* raise an error if already allocated */ 1085674ae819SStefano Zampini if (graph->nvtxs) { 10868e61c736SStefano Zampini SETERRQ(PetscObjectComm((PetscObject)l2gmap),PETSC_ERR_PLIB,"BDDCGraph already initialized"); 1087674ae819SStefano Zampini } 1088674ae819SStefano Zampini /* set number of vertices */ 1089674ae819SStefano Zampini ierr = PetscObjectReference((PetscObject)l2gmap);CHKERRQ(ierr); 1090674ae819SStefano Zampini graph->l2gmap = l2gmap; 1091674ae819SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(l2gmap,&n);CHKERRQ(ierr); 1092674ae819SStefano Zampini graph->nvtxs = n; 1093674ae819SStefano Zampini /* allocate used space */ 1094df48933bSStefano Zampini ierr = PetscBTCreate(graph->nvtxs,&graph->touched);CHKERRQ(ierr); 10958e5aaad5SJed Brown ierr = PetscMalloc7(graph->nvtxs,&graph->count, 10968e5aaad5SJed Brown graph->nvtxs,&graph->neighbours_set, 10978e5aaad5SJed Brown graph->nvtxs,&graph->subset, 10988e5aaad5SJed Brown graph->nvtxs,&graph->which_dof, 10998e5aaad5SJed Brown graph->nvtxs+1,&graph->cptr, 11008e5aaad5SJed Brown graph->nvtxs,&graph->queue, 11018e5aaad5SJed Brown graph->nvtxs,&graph->special_dof);CHKERRQ(ierr); 1102674ae819SStefano Zampini /* zeroes memory */ 1103674ae819SStefano Zampini ierr = PetscMemzero(graph->count,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr); 1104674ae819SStefano Zampini ierr = PetscMemzero(graph->subset,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr); 110563602bcaSStefano Zampini /* use -1 as a default value for which_dof array */ 110663602bcaSStefano Zampini for (n=0;n<graph->nvtxs;n++) graph->which_dof[n] = -1; 1107674ae819SStefano Zampini ierr = PetscMemzero(graph->cptr,(graph->nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr); 1108674ae819SStefano Zampini ierr = PetscMemzero(graph->queue,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr); 110974e413f5SStefano Zampini ierr = PetscMemzero(graph->special_dof,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr); 1110674ae819SStefano Zampini /* zeroes first pointer to neighbour set */ 1111674ae819SStefano Zampini if (graph->nvtxs) { 1112674ae819SStefano Zampini graph->neighbours_set[0] = 0; 1113674ae819SStefano Zampini } 1114674ae819SStefano Zampini /* zeroes workspace for values of ncc */ 1115674ae819SStefano Zampini graph->subset_ncc = 0; 11163a5b03d1SStefano Zampini graph->subset_ref_node = 0; 1117674ae819SStefano Zampini PetscFunctionReturn(0); 1118674ae819SStefano Zampini } 1119674ae819SStefano Zampini 1120674ae819SStefano Zampini #undef __FUNCT__ 1121674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphDestroy" 1122674ae819SStefano Zampini PetscErrorCode PCBDDCGraphDestroy(PCBDDCGraph* graph) 1123674ae819SStefano Zampini { 1124674ae819SStefano Zampini PetscErrorCode ierr; 1125674ae819SStefano Zampini 1126674ae819SStefano Zampini PetscFunctionBegin; 1127674ae819SStefano Zampini ierr = PCBDDCGraphReset(*graph);CHKERRQ(ierr); 1128674ae819SStefano Zampini ierr = PetscFree(*graph);CHKERRQ(ierr); 1129674ae819SStefano Zampini PetscFunctionReturn(0); 1130674ae819SStefano Zampini } 1131674ae819SStefano Zampini 1132674ae819SStefano Zampini #undef __FUNCT__ 1133674ae819SStefano Zampini #define __FUNCT__ "PCBDDCGraphCreate" 1134674ae819SStefano Zampini PetscErrorCode PCBDDCGraphCreate(PCBDDCGraph *graph) 1135674ae819SStefano Zampini { 1136674ae819SStefano Zampini PCBDDCGraph new_graph; 1137674ae819SStefano Zampini PetscErrorCode ierr; 1138674ae819SStefano Zampini 1139674ae819SStefano Zampini PetscFunctionBegin; 1140674ae819SStefano Zampini ierr = PetscMalloc(sizeof(*new_graph),&new_graph);CHKERRQ(ierr); 1141674ae819SStefano Zampini /* local to global mapping of dofs */ 1142674ae819SStefano Zampini new_graph->l2gmap = 0; 1143674ae819SStefano Zampini /* vertex size */ 1144674ae819SStefano Zampini new_graph->nvtxs = 0; 1145674ae819SStefano Zampini new_graph->n_subsets = 0; 1146674ae819SStefano Zampini new_graph->custom_minimal_size = 1; 1147674ae819SStefano Zampini /* zeroes ponters */ 114851b0f6cfSStefano Zampini new_graph->mirrors = 0; 114951b0f6cfSStefano Zampini new_graph->mirrors_set = 0; 1150674ae819SStefano Zampini new_graph->neighbours_set = 0; 1151674ae819SStefano Zampini new_graph->subset = 0; 1152674ae819SStefano Zampini new_graph->which_dof = 0; 115374e413f5SStefano Zampini new_graph->special_dof = 0; 1154674ae819SStefano Zampini new_graph->cptr = 0; 1155674ae819SStefano Zampini new_graph->queue = 0; 1156674ae819SStefano Zampini new_graph->count = 0; 1157674ae819SStefano Zampini new_graph->subset_ncc = 0; 11583a5b03d1SStefano Zampini new_graph->subset_ref_node = 0; 1159674ae819SStefano Zampini new_graph->touched = 0; 1160674ae819SStefano Zampini /* zeroes pointers to csr graph of local nodes connectivity (optional data) */ 1161575ad6abSStefano Zampini new_graph->nvtxs_csr = 0; 1162674ae819SStefano Zampini new_graph->xadj = 0; 1163674ae819SStefano Zampini new_graph->adjncy = 0; 1164674ae819SStefano Zampini *graph = new_graph; 1165674ae819SStefano Zampini PetscFunctionReturn(0); 1166674ae819SStefano Zampini } 1167