1 /* 2 * overlapsplit.c: increase the overlap of a 'big' subdomain across several processor cores 3 * 4 * Author: Fande Kong <fdkong.jd@gmail.com> 5 */ 6 7 #include <petscsf.h> 8 #include <petsc/private/matimpl.h> 9 10 11 12 /* 13 * Increase overlap for the sub-matrix across sub communicator 14 * sub-matrix could be a graph or numerical matrix 15 * */ 16 PetscErrorCode MatIncreaseOverlapSplit_Single(Mat mat,IS *is,PetscInt ov) 17 { 18 PetscInt i,nindx,*indices_sc,*indices_ov,localsize,*localsizes_sc,localsize_tmp; 19 PetscInt *indices_ov_rd,nroots,nleaves,*localoffsets,*indices_recv,*sources_sc,*sources_sc_rd; 20 const PetscInt *indices; 21 PetscMPIInt srank,ssize,issamecomm,k,grank; 22 IS is_sc,allis_sc,partitioning; 23 MPI_Comm gcomm,dcomm,scomm; 24 PetscSF sf; 25 PetscSFNode *remote; 26 Mat *smat; 27 MatPartitioning part; 28 PetscErrorCode ierr; 29 30 PetscFunctionBegin; 31 /* get a sub communicator before call individual MatIncreaseOverlap 32 * since the sub communicator may be changed. 33 * */ 34 ierr = PetscObjectGetComm((PetscObject)(*is),&dcomm);CHKERRQ(ierr); 35 /* make a copy before the original one is deleted */ 36 ierr = PetscCommDuplicate(dcomm,&scomm,NULL);CHKERRQ(ierr); 37 /* get a global communicator, where mat should be a global matrix */ 38 ierr = PetscObjectGetComm((PetscObject)mat,&gcomm);CHKERRQ(ierr); 39 ierr = (*mat->ops->increaseoverlap)(mat,1,is,ov);CHKERRQ(ierr); 40 ierr = MPI_Comm_compare(gcomm,scomm,&issamecomm);CHKERRQ(ierr); 41 /* if the sub-communicator is the same as the global communicator, 42 * user does not want to use a sub-communicator 43 * */ 44 if(issamecomm == MPI_IDENT || issamecomm == MPI_CONGRUENT){ 45 ierr = PetscCommDestroy(&scomm);CHKERRQ(ierr); 46 PetscFunctionReturn(0); 47 } 48 /* if the sub-communicator is petsc_comm_self, 49 * user also does not care the sub-communicator 50 * */ 51 ierr = MPI_Comm_compare(scomm,PETSC_COMM_SELF,&issamecomm);CHKERRQ(ierr); 52 if(issamecomm == MPI_IDENT || issamecomm == MPI_CONGRUENT){ 53 ierr = PetscCommDestroy(&scomm);CHKERRQ(ierr); 54 PetscFunctionReturn(0); 55 } 56 ierr = MPI_Comm_rank(scomm,&srank);CHKERRQ(ierr); 57 ierr = MPI_Comm_size(scomm,&ssize);CHKERRQ(ierr); 58 ierr = MPI_Comm_rank(gcomm,&grank);CHKERRQ(ierr); 59 /* create a new IS based on sub-communicator 60 * since the old IS is often based on petsc_comm_self 61 * */ 62 ierr = ISGetLocalSize(*is,&nindx);CHKERRQ(ierr); 63 ierr = PetscCalloc1(nindx,&indices_sc);CHKERRQ(ierr); 64 ierr = ISGetIndices(*is,&indices);CHKERRQ(ierr); 65 ierr = PetscMemcpy(indices_sc,indices,sizeof(PetscInt)*nindx);CHKERRQ(ierr); 66 ierr = ISRestoreIndices(*is,&indices);CHKERRQ(ierr); 67 /* we do not need any more */ 68 ierr = ISDestroy(is);CHKERRQ(ierr); 69 /* create a index set based on the sub communicator */ 70 ierr = ISCreateGeneral(scomm,nindx,indices_sc,PETSC_OWN_POINTER,&is_sc);CHKERRQ(ierr); 71 /* gather all indices within the sub communicator */ 72 ierr = ISAllGather(is_sc,&allis_sc);CHKERRQ(ierr); 73 ierr = ISDestroy(&is_sc);CHKERRQ(ierr); 74 /* gather local sizes */ 75 ierr = PetscMalloc1(ssize,&localsizes_sc);CHKERRQ(ierr); 76 /* get individual local sizes for all index sets */ 77 ierr = MPI_Gather(&nindx,1,MPIU_INT,localsizes_sc,1,MPIU_INT,0,scomm);CHKERRQ(ierr); 78 /* only root does these computations */ 79 if(!srank){ 80 /* get local size for the big index set */ 81 ierr = ISGetLocalSize(allis_sc,&localsize);CHKERRQ(ierr); 82 ierr = PetscCalloc2(localsize,&indices_ov,localsize,&sources_sc);CHKERRQ(ierr); 83 ierr = PetscCalloc2(localsize,&indices_ov_rd,localsize,&sources_sc_rd);CHKERRQ(ierr); 84 ierr = ISGetIndices(allis_sc,&indices);CHKERRQ(ierr); 85 ierr = PetscMemcpy(indices_ov,indices,sizeof(PetscInt)*localsize);CHKERRQ(ierr); 86 ierr = ISRestoreIndices(allis_sc,&indices);CHKERRQ(ierr); 87 ierr = ISDestroy(&allis_sc);CHKERRQ(ierr); 88 /* assign corresponding sources */ 89 localsize_tmp = 0; 90 for(k=0; k<ssize; k++){ 91 for(i=0; i<localsizes_sc[k]; i++){ 92 sources_sc[localsize_tmp++] = k; 93 } 94 } 95 /* record where indices come from */ 96 ierr = PetscSortIntWithArray(localsize,indices_ov,sources_sc);CHKERRQ(ierr); 97 /* count local sizes for reduced indices */ 98 ierr = PetscMemzero(localsizes_sc,sizeof(PetscInt)*ssize);CHKERRQ(ierr); 99 /* initialize the first entity */ 100 if(localsize){ 101 indices_ov_rd[0] = indices_ov[0]; 102 sources_sc_rd[0] = sources_sc[0]; 103 localsizes_sc[sources_sc[0]]++; 104 } 105 localsize_tmp = 1; 106 /* remove duplicate integers */ 107 for(i=1; i<localsize; i++){ 108 if(indices_ov[i] != indices_ov[i-1]){ 109 indices_ov_rd[localsize_tmp] = indices_ov[i]; 110 sources_sc_rd[localsize_tmp++] = sources_sc[i]; 111 localsizes_sc[sources_sc[i]]++; 112 } 113 } 114 ierr = PetscFree2(indices_ov,sources_sc);CHKERRQ(ierr); 115 ierr = PetscCalloc1(ssize+1,&localoffsets);CHKERRQ(ierr); 116 for(k=0; k<ssize; k++){ 117 localoffsets[k+1] = localoffsets[k] + localsizes_sc[k]; 118 } 119 nleaves = localoffsets[ssize]; 120 ierr = PetscMemzero(localoffsets,(ssize+1)*sizeof(PetscInt));CHKERRQ(ierr); 121 nroots = localsizes_sc[srank]; 122 ierr = PetscCalloc1(nleaves,&remote);CHKERRQ(ierr); 123 for(i=0; i<nleaves; i++){ 124 remote[i].rank = sources_sc_rd[i]; 125 remote[i].index = localoffsets[sources_sc_rd[i]]++; 126 } 127 ierr = PetscFree(localoffsets);CHKERRQ(ierr); 128 }else{ 129 ierr = ISDestroy(&allis_sc);CHKERRQ(ierr); 130 /* Allocate a 'zero' pointer to avoid using uninitialized variable */ 131 ierr = PetscCalloc1(0,&remote);CHKERRQ(ierr); 132 nleaves = 0; 133 indices_ov_rd = 0; 134 sources_sc_rd = 0; 135 } 136 /* scatter sizes to everybody */ 137 ierr = MPI_Scatter(localsizes_sc,1, MPIU_INT,&nroots,1, MPIU_INT,0,scomm);CHKERRQ(ierr); 138 ierr = PetscFree(localsizes_sc);CHKERRQ(ierr); 139 ierr = PetscCalloc1(nroots,&indices_recv);CHKERRQ(ierr); 140 /* set data back to every body */ 141 ierr = PetscSFCreate(scomm,&sf);CHKERRQ(ierr); 142 ierr = PetscSFSetType(sf,PETSCSFBASIC);CHKERRQ(ierr); 143 ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr); 144 ierr = PetscSFSetGraph(sf,nroots,nleaves,NULL,PETSC_OWN_POINTER,remote,PETSC_OWN_POINTER);CHKERRQ(ierr); 145 ierr = PetscSFReduceBegin(sf,MPIU_INT,indices_ov_rd,indices_recv,MPIU_REPLACE);CHKERRQ(ierr); 146 ierr = PetscSFReduceEnd(sf,MPIU_INT,indices_ov_rd,indices_recv,MPIU_REPLACE);CHKERRQ(ierr); 147 ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 148 ierr = PetscFree2(indices_ov_rd,sources_sc_rd);CHKERRQ(ierr); 149 ierr = ISCreateGeneral(scomm,nroots,indices_recv,PETSC_OWN_POINTER,&is_sc);CHKERRQ(ierr); 150 ierr = MatCreateSubMatricesMPI(mat,1,&is_sc,&is_sc,MAT_INITIAL_MATRIX,&smat);CHKERRQ(ierr); 151 ierr = ISDestroy(&allis_sc);CHKERRQ(ierr); 152 /* create a partitioner to repartition the sub-matrix */ 153 ierr = MatPartitioningCreate(scomm,&part);CHKERRQ(ierr); 154 ierr = MatPartitioningSetAdjacency(part,smat[0]);CHKERRQ(ierr); 155 #if defined(PETSC_HAVE_PARMETIS) 156 /* if there exists a ParMETIS installation, we try to use ParMETIS 157 * because a repartition routine possibly work better 158 * */ 159 ierr = MatPartitioningSetType(part,MATPARTITIONINGPARMETIS);CHKERRQ(ierr); 160 /* try to use reparition function, instead of partition function */ 161 ierr = MatPartitioningParmetisSetRepartition(part);CHKERRQ(ierr); 162 #else 163 /* we at least provide a default partitioner to rebalance the computation */ 164 ierr = MatPartitioningSetType(part,MATPARTITIONINGAVERAGE);CHKERRQ(ierr); 165 #endif 166 /* user can pick up any partitioner by using an option */ 167 ierr = MatPartitioningSetFromOptions(part);CHKERRQ(ierr); 168 ierr = MatPartitioningApply(part,&partitioning);CHKERRQ(ierr); 169 ierr = MatPartitioningDestroy(&part);CHKERRQ(ierr); 170 ierr = MatDestroy(&(smat[0]));CHKERRQ(ierr); 171 ierr = PetscFree(smat);CHKERRQ(ierr); 172 /* get local rows including overlap */ 173 ierr = ISBuildTwoSided(partitioning,is_sc,is);CHKERRQ(ierr); 174 ierr = ISDestroy(&is_sc);CHKERRQ(ierr); 175 ierr = ISDestroy(&partitioning);CHKERRQ(ierr); 176 ierr = PetscCommDestroy(&scomm);CHKERRQ(ierr); 177 PetscFunctionReturn(0); 178 } 179 180 181