xref: /petsc/src/mat/utils/overlapsplit.c (revision e12b4729b15ad5de2ed17324d40445da9b32e630)
1bac5b06fSFande Kong /*
2bac5b06fSFande Kong  * overlapsplit.c: increase the overlap of a 'big' subdomain across several processor cores
3bac5b06fSFande Kong  *
4bac5b06fSFande Kong  * Author: Fande Kong <fdkong.jd@gmail.com>
5bac5b06fSFande Kong  */
6bac5b06fSFande Kong 
72452736bSFande Kong #include <petscsf.h>
82452736bSFande Kong #include <petsc/private/matimpl.h>
92452736bSFande Kong 
102452736bSFande Kong 
11ca2fc57aSFande Kong #undef __FUNCT__
12ca2fc57aSFande Kong #define __FUNCT__ "MatIncreaseOverlapSplit_Single"
13ca2fc57aSFande Kong 
142452736bSFande Kong /*
152452736bSFande Kong  * Increase overlap for the sub-matrix across sub communicator
162452736bSFande Kong  * sub-matrix could be a graph or numerical matrix
172452736bSFande Kong  * */
182452736bSFande Kong PetscErrorCode  MatIncreaseOverlapSplit_Single(Mat mat,IS *is,PetscInt ov)
192452736bSFande Kong {
202452736bSFande Kong   PetscInt         i,nindx,*indices_sc,*indices_ov,localsize,*localsizes_sc,localsize_tmp;
212452736bSFande Kong   PetscInt         *indices_ov_rd,nroots,nleaves,*localoffsets,*indices_recv,*sources_sc,*sources_sc_rd;
222452736bSFande Kong   const PetscInt   *indices;
23a69400e0SFande Kong   PetscMPIInt      srank,ssize,issamecomm,k,grank;
24bac5b06fSFande Kong   IS               is_sc,allis_sc,partitioning;
25a69400e0SFande Kong   MPI_Comm         gcomm,dcomm,scomm;
262452736bSFande Kong   PetscSF          sf;
272452736bSFande Kong   PetscSFNode      *remote;
282452736bSFande Kong   Mat              *smat;
292452736bSFande Kong   MatPartitioning  part;
302452736bSFande Kong   PetscErrorCode   ierr;
312452736bSFande Kong 
322452736bSFande Kong   PetscFunctionBegin;
332452736bSFande Kong   /* get a sub communicator before call individual MatIncreaseOverlap
342452736bSFande Kong    * since the sub communicator may be changed.
352452736bSFande Kong    * */
36a69400e0SFande Kong   ierr = PetscObjectGetComm((PetscObject)(*is),&dcomm);CHKERRQ(ierr);
37a69400e0SFande Kong   /* make a copy before the original one is deleted */
38a69400e0SFande Kong   ierr = PetscCommDuplicate(dcomm,&scomm,NULL);CHKERRQ(ierr);
39ca2fc57aSFande Kong   /* get a global communicator, where mat should be a global matrix  */
40ca2fc57aSFande Kong   ierr = PetscObjectGetComm((PetscObject)mat,&gcomm);CHKERRQ(ierr);
412452736bSFande Kong   ierr = (*mat->ops->increaseoverlap)(mat,1,is,ov);CHKERRQ(ierr);
422452736bSFande Kong   ierr = MPI_Comm_compare(gcomm,scomm,&issamecomm);CHKERRQ(ierr);
432452736bSFande Kong   /* if the sub-communicator is the same as the global communicator,
442452736bSFande Kong    * user does not want to use a sub-communicator
452452736bSFande Kong    * */
46a69400e0SFande Kong   if(issamecomm == MPI_IDENT || issamecomm == MPI_CONGRUENT) PetscFunctionReturn(0);
472452736bSFande Kong   /* if the sub-communicator is petsc_comm_self,
482452736bSFande Kong    * user also does not care the sub-communicator
492452736bSFande Kong    * */
502452736bSFande Kong   ierr = MPI_Comm_compare(scomm,PETSC_COMM_SELF,&issamecomm);CHKERRQ(ierr);
51a69400e0SFande Kong   if(issamecomm == MPI_IDENT || issamecomm == MPI_CONGRUENT){PetscFunctionReturn(0);}
522452736bSFande Kong   ierr = MPI_Comm_rank(scomm,&srank);CHKERRQ(ierr);
532452736bSFande Kong   ierr = MPI_Comm_size(scomm,&ssize);CHKERRQ(ierr);
54a69400e0SFande Kong   ierr = MPI_Comm_rank(gcomm,&grank);CHKERRQ(ierr);
55ca2fc57aSFande Kong   /* create a new IS based on sub-communicator
56ca2fc57aSFande Kong    * since the old IS is often based on petsc_comm_self
572452736bSFande Kong    * */
582452736bSFande Kong   ierr = ISGetLocalSize(*is,&nindx);CHKERRQ(ierr);
592452736bSFande Kong   ierr = PetscCalloc1(nindx,&indices_sc);CHKERRQ(ierr);
602452736bSFande Kong   ierr = ISGetIndices(*is,&indices);CHKERRQ(ierr);
612452736bSFande Kong   ierr = PetscMemcpy(indices_sc,indices,sizeof(PetscInt)*nindx);CHKERRQ(ierr);
622452736bSFande Kong   ierr = ISRestoreIndices(*is,&indices);CHKERRQ(ierr);
632452736bSFande Kong   /* we do not need any more */
642452736bSFande Kong   ierr = ISDestroy(is);CHKERRQ(ierr);
652452736bSFande Kong   /* create a index set based on the sub communicator  */
662452736bSFande Kong   ierr = ISCreateGeneral(scomm,nindx,indices_sc,PETSC_OWN_POINTER,&is_sc);CHKERRQ(ierr);
672452736bSFande Kong   /* gather all indices within  the sub communicator */
682452736bSFande Kong   ierr = ISAllGather(is_sc,&allis_sc);CHKERRQ(ierr);
692452736bSFande Kong   ierr = ISDestroy(&is_sc);CHKERRQ(ierr);
702452736bSFande Kong   /* gather local sizes */
712452736bSFande Kong   ierr = PetscMalloc1(ssize,&localsizes_sc);CHKERRQ(ierr);
72ca2fc57aSFande Kong   /* get individual local sizes for all index sets */
732452736bSFande Kong   ierr = MPI_Gather(&nindx,1,MPIU_INT,localsizes_sc,1,MPIU_INT,0,scomm);CHKERRQ(ierr);
74ca2fc57aSFande Kong   /* only root does these computations */
752452736bSFande Kong   if(!srank){
762452736bSFande Kong    /* get local size for the big index set */
772452736bSFande Kong    ierr = ISGetLocalSize(allis_sc,&localsize);CHKERRQ(ierr);
782452736bSFande Kong    ierr = PetscCalloc2(localsize,&indices_ov,localsize,&sources_sc);CHKERRQ(ierr);
792452736bSFande Kong    ierr = PetscCalloc2(localsize,&indices_ov_rd,localsize,&sources_sc_rd);CHKERRQ(ierr);
802452736bSFande Kong    ierr = ISGetIndices(allis_sc,&indices);CHKERRQ(ierr);
812452736bSFande Kong    ierr = PetscMemcpy(indices_ov,indices,sizeof(PetscInt)*localsize);CHKERRQ(ierr);
822452736bSFande Kong    ierr = ISRestoreIndices(allis_sc,&indices);CHKERRQ(ierr);
832452736bSFande Kong    ierr = ISDestroy(&allis_sc);CHKERRQ(ierr);
842452736bSFande Kong    /* assign corresponding sources */
852452736bSFande Kong    localsize_tmp = 0;
862452736bSFande Kong    for(k=0; k<ssize; k++){
872452736bSFande Kong      for(i=0; i<localsizes_sc[k]; i++){
882452736bSFande Kong        sources_sc[localsize_tmp++] = k;
892452736bSFande Kong      }
902452736bSFande Kong    }
912452736bSFande Kong    /* record where indices come from */
922452736bSFande Kong    ierr = PetscSortIntWithArray(localsize,indices_ov,sources_sc);CHKERRQ(ierr);
93ca2fc57aSFande Kong    /* count local sizes for reduced indices */
942452736bSFande Kong    ierr = PetscMemzero(localsizes_sc,sizeof(PetscInt)*ssize);CHKERRQ(ierr);
95ca2fc57aSFande Kong    /* initialize the first entity */
962452736bSFande Kong    if(localsize){
972452736bSFande Kong 	 indices_ov_rd[0] = indices_ov[0];
982452736bSFande Kong 	 sources_sc_rd[0] = sources_sc[0];
992452736bSFande Kong 	 localsizes_sc[sources_sc[0]]++;
1002452736bSFande Kong    }
101ca2fc57aSFande Kong    localsize_tmp = 1;
1022452736bSFande Kong    /* remove duplicate integers */
1032452736bSFande Kong    for(i=1; i<localsize; i++){
1042452736bSFande Kong 	 if(indices_ov[i] != indices_ov[i-1]){
1052452736bSFande Kong 	   indices_ov_rd[localsize_tmp]   = indices_ov[i];
1062452736bSFande Kong 	   sources_sc_rd[localsize_tmp++] = sources_sc[i];
1072452736bSFande Kong 	   localsizes_sc[sources_sc[i]]++;
1082452736bSFande Kong 	 }
1092452736bSFande Kong    }
1102452736bSFande Kong    ierr = PetscFree2(indices_ov,sources_sc);CHKERRQ(ierr);
1112452736bSFande Kong    ierr = PetscCalloc1(ssize+1,&localoffsets);CHKERRQ(ierr);
1122452736bSFande Kong    for(k=0; k<ssize; k++){
113a69400e0SFande Kong 	 localoffsets[k+1] = localoffsets[k] + localsizes_sc[k];
1142452736bSFande Kong    }
1152452736bSFande Kong    nleaves = localoffsets[ssize];
1162452736bSFande Kong    ierr = PetscMemzero(localoffsets,(ssize+1)*sizeof(PetscInt));CHKERRQ(ierr);
1172452736bSFande Kong    nroots  = localsizes_sc[srank];
1182452736bSFande Kong    ierr = PetscCalloc1(nleaves,&remote);CHKERRQ(ierr);
1192452736bSFande Kong    for(i=0; i<nleaves; i++){
120ca2fc57aSFande Kong 	 remote[i].rank  = sources_sc_rd[i];
121ca2fc57aSFande Kong 	 remote[i].index = localoffsets[sources_sc_rd[i]]++;
1222452736bSFande Kong    }
123ca2fc57aSFande Kong    ierr = PetscFree(localoffsets);CHKERRQ(ierr);
1242452736bSFande Kong   }else{
1252452736bSFande Kong    ierr = ISDestroy(&allis_sc);CHKERRQ(ierr);
126*e12b4729SFande Kong    /* Allocate a 'zero' pointer to avoid using uninitialized variable  */
127bac5b06fSFande Kong    ierr = PetscCalloc1(0,&remote);CHKERRQ(ierr);
1282452736bSFande Kong    nleaves = 0;
1292452736bSFande Kong    indices_ov_rd = 0;
1302452736bSFande Kong    sources_sc_rd = 0;
1312452736bSFande Kong   }
1322452736bSFande Kong   /* scatter sizes to everybody */
1332452736bSFande Kong   ierr = MPI_Scatter(localsizes_sc,1, MPIU_INT,&nroots,1, MPIU_INT,0,scomm);CHKERRQ(ierr);
134ca2fc57aSFande Kong   ierr = PetscFree(localsizes_sc);CHKERRQ(ierr);
1352452736bSFande Kong   ierr = PetscCalloc1(nroots,&indices_recv);CHKERRQ(ierr);
1362452736bSFande Kong   /* set data back to every body */
1372452736bSFande Kong   ierr = PetscSFCreate(scomm,&sf);CHKERRQ(ierr);
1382452736bSFande Kong   ierr = PetscSFSetType(sf,PETSCSFBASIC);CHKERRQ(ierr);
1392452736bSFande Kong   ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr);
1402452736bSFande Kong   ierr = PetscSFSetGraph(sf,nroots,nleaves,PETSC_NULL,PETSC_OWN_POINTER,remote,PETSC_OWN_POINTER);CHKERRQ(ierr);
141a69400e0SFande Kong   ierr = PetscSFReduceBegin(sf,MPIU_INT,indices_ov_rd,indices_recv,MPIU_REPLACE);CHKERRQ(ierr);
142a69400e0SFande Kong   ierr = PetscSFReduceEnd(sf,MPIU_INT,indices_ov_rd,indices_recv,MPIU_REPLACE);CHKERRQ(ierr);
1432452736bSFande Kong   ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
1442452736bSFande Kong   ierr = PetscFree2(indices_ov_rd,sources_sc_rd);CHKERRQ(ierr);
1452452736bSFande Kong   ierr = ISCreateGeneral(scomm,nroots,indices_recv,PETSC_OWN_POINTER,&is_sc);CHKERRQ(ierr);
146bac5b06fSFande Kong   ierr = MatGetSubMatricesMPI(mat,1,&is_sc,&is_sc,MAT_INITIAL_MATRIX,&smat);CHKERRQ(ierr);
1472452736bSFande Kong   ierr = ISDestroy(&allis_sc);CHKERRQ(ierr);
1482452736bSFande Kong   /* create a partitioner to repartition the sub-matrix */
1492452736bSFande Kong   ierr = MatPartitioningCreate(scomm,&part);CHKERRQ(ierr);
1502452736bSFande Kong   ierr = MatPartitioningSetAdjacency(part,smat[0]);CHKERRQ(ierr);
1512452736bSFande Kong #if PETSC_HAVE_PARMETIS
1522452736bSFande Kong   /* if there exists a ParMETIS installation, we try to use ParMETIS
1532452736bSFande Kong    * because a repartition routine possibly work better
1542452736bSFande Kong    * */
1552452736bSFande Kong   ierr = MatPartitioningSetType(part,MATPARTITIONINGPARMETIS);CHKERRQ(ierr);
1562452736bSFande Kong   /* try to use reparition function, instead of partition function */
1572452736bSFande Kong   ierr = MatPartitioningParmetisSetRepartition(part);CHKERRQ(ierr);
1582452736bSFande Kong #else
1592452736bSFande Kong   /* we at least provide a default partitioner to rebalance the computation  */
1602452736bSFande Kong   ierr = MatPartitioningSetType(part,MATPARTITIONINGAVERAGE);CHKERRQ(ierr);
1612452736bSFande Kong #endif
1622452736bSFande Kong   /* user can pick up any partitioner by using an option */
1632452736bSFande Kong   ierr = MatPartitioningSetFromOptions(part);CHKERRQ(ierr);
1642452736bSFande Kong   ierr = MatPartitioningApply(part,&partitioning);CHKERRQ(ierr);
1652452736bSFande Kong   ierr = MatPartitioningDestroy(&part);CHKERRQ(ierr);
166ca2fc57aSFande Kong   ierr = MatDestroy(&(smat[0]));CHKERRQ(ierr);
1672452736bSFande Kong   ierr = PetscFree(smat);CHKERRQ(ierr);
1682452736bSFande Kong   /* get local rows including  overlap */
169bdeb5f4eSFande Kong   ierr = ISBuildTwoSided(partitioning,is_sc,is);CHKERRQ(ierr);
170bdeb5f4eSFande Kong   ierr = ISDestroy(&is_sc);CHKERRQ(ierr);
17108902ef5SFande Kong   ierr = ISDestroy(&partitioning);CHKERRQ(ierr);
172e3c39df8SFande Kong   ierr = PetscCommDestroy(&scomm);CHKERRQ(ierr);
1732452736bSFande Kong   PetscFunctionReturn(0);
1742452736bSFande Kong }
1752452736bSFande Kong 
1762452736bSFande Kong 
177