xref: /petsc/src/mat/utils/overlapsplit.c (revision bac5b06fd7acf2dbdffb3ef782ccfe7a95cc790d)
1*bac5b06fSFande Kong /*
2*bac5b06fSFande Kong  * overlapsplit.c: increase the overlap of a 'big' subdomain across several processor cores
3*bac5b06fSFande Kong  *
4*bac5b06fSFande Kong  * Author: Fande Kong <fdkong.jd@gmail.com>
5*bac5b06fSFande Kong  */
6*bac5b06fSFande 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;
24*bac5b06fSFande 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);
41*bac5b06fSFande Kong #if 0
42*bac5b06fSFande Kong   ierr = MPI_Barrier(gcomm);CHKERRQ(ierr);
43*bac5b06fSFande Kong   ierr = PetscPrintf(PETSC_COMM_SELF,"before mat->ops->increaseoverlap\n");CHKERRQ(ierr);
44ca2fc57aSFande Kong #endif
452452736bSFande Kong   /*increase overlap on each individual subdomain*/
462452736bSFande Kong   ierr = (*mat->ops->increaseoverlap)(mat,1,is,ov);CHKERRQ(ierr);
47*bac5b06fSFande Kong #if 0
48ca2fc57aSFande Kong   ierr = PetscPrintf(gcomm,"after mat->ops->increaseoverlap \n");CHKERRQ(ierr);
49a69400e0SFande Kong   ierr = ISView(*is,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
50*bac5b06fSFande Kong   ierr = MPI_Barrier(gcomm);CHKERRQ(ierr);
51ca2fc57aSFande Kong #endif
522452736bSFande Kong   /*compare communicators */
532452736bSFande Kong   ierr = MPI_Comm_compare(gcomm,scomm,&issamecomm);CHKERRQ(ierr);
542452736bSFande Kong   /* if the sub-communicator is the same as the global communicator,
552452736bSFande Kong    * user does not want to use a sub-communicator
562452736bSFande Kong    * */
57a69400e0SFande Kong   if(issamecomm == MPI_IDENT || issamecomm == MPI_CONGRUENT) PetscFunctionReturn(0);
582452736bSFande Kong   /* if the sub-communicator is petsc_comm_self,
592452736bSFande Kong    * user also does not care the sub-communicator
602452736bSFande Kong    * */
612452736bSFande Kong   ierr = MPI_Comm_compare(scomm,PETSC_COMM_SELF,&issamecomm);CHKERRQ(ierr);
62a69400e0SFande Kong   if(issamecomm == MPI_IDENT || issamecomm == MPI_CONGRUENT){PetscFunctionReturn(0);}
63ca2fc57aSFande Kong   /*local rank, size in a sub-communicator  */
642452736bSFande Kong   ierr = MPI_Comm_rank(scomm,&srank);CHKERRQ(ierr);
652452736bSFande Kong   ierr = MPI_Comm_size(scomm,&ssize);CHKERRQ(ierr);
66a69400e0SFande Kong   ierr = MPI_Comm_rank(gcomm,&grank);CHKERRQ(ierr);
67ca2fc57aSFande Kong   /*create a new IS based on sub-communicator
68ca2fc57aSFande Kong    * since the old IS is often based on petsc_comm_self
692452736bSFande Kong    * */
702452736bSFande Kong   ierr = ISGetLocalSize(*is,&nindx);CHKERRQ(ierr);
712452736bSFande Kong   ierr = PetscCalloc1(nindx,&indices_sc);CHKERRQ(ierr);
722452736bSFande Kong   ierr = ISGetIndices(*is,&indices);CHKERRQ(ierr);
732452736bSFande Kong   ierr = PetscMemcpy(indices_sc,indices,sizeof(PetscInt)*nindx);CHKERRQ(ierr);
742452736bSFande Kong   ierr = ISRestoreIndices(*is,&indices);CHKERRQ(ierr);
752452736bSFande Kong   /*we do not need any more*/
762452736bSFande Kong   ierr = ISDestroy(is);CHKERRQ(ierr);
772452736bSFande Kong   /*create a index set based on the sub communicator  */
782452736bSFande Kong   ierr = ISCreateGeneral(scomm,nindx,indices_sc,PETSC_OWN_POINTER,&is_sc);CHKERRQ(ierr);
792452736bSFande Kong   /*gather all indices within  the sub communicator*/
802452736bSFande Kong   ierr = ISAllGather(is_sc,&allis_sc);CHKERRQ(ierr);
81*bac5b06fSFande Kong #if 0
82a69400e0SFande Kong   ierr = ISView(allis_sc,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
83a69400e0SFande Kong   ierr = MPI_Barrier(gcomm);CHKERRQ(ierr);
84a69400e0SFande Kong #endif
852452736bSFande Kong   ierr = ISDestroy(&is_sc);CHKERRQ(ierr);
862452736bSFande Kong   /* gather local sizes */
872452736bSFande Kong   ierr = PetscMalloc1(ssize,&localsizes_sc);CHKERRQ(ierr);
88ca2fc57aSFande Kong   /*get individual local sizes for all index sets*/
892452736bSFande Kong   ierr = MPI_Gather(&nindx,1,MPIU_INT,localsizes_sc,1,MPIU_INT,0,scomm);CHKERRQ(ierr);
90*bac5b06fSFande Kong #if 0
91a69400e0SFande Kong   if(!srank){
92a69400e0SFande Kong 	for(i=0; i<ssize; i++){
93a69400e0SFande Kong 	  ierr = PetscPrintf(PETSC_COMM_SELF," localsize[%d]: %d \n",i,localsizes_sc[i]);CHKERRQ(ierr);
94a69400e0SFande Kong 	}
95a69400e0SFande Kong   }
96a69400e0SFande Kong   ierr = MPI_Barrier(gcomm);CHKERRQ(ierr);
97a69400e0SFande Kong #endif
98ca2fc57aSFande Kong   /*only root does these computations */
992452736bSFande Kong   if(!srank){
1002452736bSFande Kong    /*get local size for the big index set*/
1012452736bSFande Kong    ierr = ISGetLocalSize(allis_sc,&localsize);CHKERRQ(ierr);
1022452736bSFande Kong    ierr = PetscCalloc2(localsize,&indices_ov,localsize,&sources_sc);CHKERRQ(ierr);
1032452736bSFande Kong    ierr = PetscCalloc2(localsize,&indices_ov_rd,localsize,&sources_sc_rd);CHKERRQ(ierr);
1042452736bSFande Kong    ierr = ISGetIndices(allis_sc,&indices);CHKERRQ(ierr);
1052452736bSFande Kong    ierr = PetscMemcpy(indices_ov,indices,sizeof(PetscInt)*localsize);CHKERRQ(ierr);
1062452736bSFande Kong    ierr = ISRestoreIndices(allis_sc,&indices);CHKERRQ(ierr);
107ca2fc57aSFande Kong    /*we do not need it any more */
1082452736bSFande Kong    ierr = ISDestroy(&allis_sc);CHKERRQ(ierr);
1092452736bSFande Kong    /*assign corresponding sources */
1102452736bSFande Kong    localsize_tmp = 0;
1112452736bSFande Kong    for(k=0; k<ssize; k++){
1122452736bSFande Kong      for(i=0; i<localsizes_sc[k]; i++){
1132452736bSFande Kong        sources_sc[localsize_tmp++] = k;
1142452736bSFande Kong      }
1152452736bSFande Kong    }
1162452736bSFande Kong    /*record where indices come from */
1172452736bSFande Kong    ierr = PetscSortIntWithArray(localsize,indices_ov,sources_sc);CHKERRQ(ierr);
118a69400e0SFande Kong #if 0
119a69400e0SFande Kong    ierr = PetscIntView(localsize,indices_ov,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
120a69400e0SFande Kong    ierr = PetscIntView(localsize,sources_sc,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
121a69400e0SFande Kong #endif
122ca2fc57aSFande Kong    /*count local sizes for reduced indices */
1232452736bSFande Kong    ierr = PetscMemzero(localsizes_sc,sizeof(PetscInt)*ssize);CHKERRQ(ierr);
124ca2fc57aSFande Kong    /*initialize the first entity*/
1252452736bSFande Kong    if(localsize){
1262452736bSFande Kong 	 indices_ov_rd[0] = indices_ov[0];
1272452736bSFande Kong 	 sources_sc_rd[0] = sources_sc[0];
1282452736bSFande Kong 	 localsizes_sc[sources_sc[0]]++;
1292452736bSFande Kong    }
130ca2fc57aSFande Kong    localsize_tmp = 1;
1312452736bSFande Kong    /*remove duplicate integers */
1322452736bSFande Kong    for(i=1; i<localsize; i++){
1332452736bSFande Kong 	 if(indices_ov[i] != indices_ov[i-1]){
1342452736bSFande Kong 	   indices_ov_rd[localsize_tmp]   = indices_ov[i];
1352452736bSFande Kong 	   sources_sc_rd[localsize_tmp++] = sources_sc[i];
1362452736bSFande Kong 	   localsizes_sc[sources_sc[i]]++;
1372452736bSFande Kong 	 }
1382452736bSFande Kong    }
1392452736bSFande Kong    ierr = PetscFree2(indices_ov,sources_sc);CHKERRQ(ierr);
1402452736bSFande Kong    ierr = PetscCalloc1(ssize+1,&localoffsets);CHKERRQ(ierr);
1412452736bSFande Kong    for(k=0; k<ssize; k++){
142a69400e0SFande Kong 	 localoffsets[k+1] = localoffsets[k] + localsizes_sc[k];
1432452736bSFande Kong    }
144ca2fc57aSFande Kong    /*construct a star forest to send data back */
1452452736bSFande Kong    nleaves = localoffsets[ssize];
1462452736bSFande Kong    ierr = PetscMemzero(localoffsets,(ssize+1)*sizeof(PetscInt));CHKERRQ(ierr);
1472452736bSFande Kong    nroots  = localsizes_sc[srank];
1482452736bSFande Kong    ierr = PetscCalloc1(nleaves,&remote);CHKERRQ(ierr);
1492452736bSFande Kong    for(i=0; i<nleaves; i++){
150ca2fc57aSFande Kong 	 remote[i].rank  = sources_sc_rd[i];
151ca2fc57aSFande Kong 	 remote[i].index = localoffsets[sources_sc_rd[i]]++;
1522452736bSFande Kong    }
153ca2fc57aSFande Kong    ierr = PetscFree(localoffsets);CHKERRQ(ierr);
154a69400e0SFande Kong #if 0
155a69400e0SFande Kong    if(grank==2){
156a69400e0SFande Kong 	 ierr = PetscIntView(localsize_tmp,indices_ov_rd,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
157a69400e0SFande Kong    }
158a69400e0SFande Kong #endif
1592452736bSFande Kong   }else{
1602452736bSFande Kong    ierr = ISDestroy(&allis_sc);CHKERRQ(ierr);
161*bac5b06fSFande Kong    /*Allocate a 'zero' pointer */
162*bac5b06fSFande Kong    ierr = PetscCalloc1(0,&remote);CHKERRQ(ierr);
1632452736bSFande Kong    nleaves = 0;
1642452736bSFande Kong    indices_ov_rd = 0;
1652452736bSFande Kong    sources_sc_rd = 0;
1662452736bSFande Kong   }
1672452736bSFande Kong   /*scatter sizes to everybody */
1682452736bSFande Kong   ierr = MPI_Scatter(localsizes_sc,1, MPIU_INT,&nroots,1, MPIU_INT,0,scomm);CHKERRQ(ierr);
169ca2fc57aSFande Kong   /*free memory */
170ca2fc57aSFande Kong   ierr = PetscFree(localsizes_sc);CHKERRQ(ierr);
1712452736bSFande Kong   ierr = PetscCalloc1(nroots,&indices_recv);CHKERRQ(ierr);
172a69400e0SFande Kong   /*ierr = MPI_Comm_dup(scomm,&dcomm);CHKERRQ(ierr);*/
1732452736bSFande Kong   /*set data back to every body */
1742452736bSFande Kong   ierr = PetscSFCreate(scomm,&sf);CHKERRQ(ierr);
1752452736bSFande Kong   ierr = PetscSFSetType(sf,PETSCSFBASIC);CHKERRQ(ierr);
1762452736bSFande Kong   ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr);
1772452736bSFande Kong   ierr = PetscSFSetGraph(sf,nroots,nleaves,PETSC_NULL,PETSC_OWN_POINTER,remote,PETSC_OWN_POINTER);CHKERRQ(ierr);
178a69400e0SFande Kong #if 0
179a69400e0SFande Kong   ierr = PetscSFView(sf,PETSC_NULL);CHKERRQ(ierr);
180a69400e0SFande Kong   ierr = MPI_Barrier(gcomm);CHKERRQ(ierr);
181a69400e0SFande Kong #endif
182a69400e0SFande Kong   ierr = PetscSFReduceBegin(sf,MPIU_INT,indices_ov_rd,indices_recv,MPIU_REPLACE);CHKERRQ(ierr);
183a69400e0SFande Kong   ierr = PetscSFReduceEnd(sf,MPIU_INT,indices_ov_rd,indices_recv,MPIU_REPLACE);CHKERRQ(ierr);
1842452736bSFande Kong   ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
1852452736bSFande Kong   /* free memory */
1862452736bSFande Kong   ierr = PetscFree2(indices_ov_rd,sources_sc_rd);CHKERRQ(ierr);
1872452736bSFande Kong   /*create a index set*/
1882452736bSFande Kong   ierr = ISCreateGeneral(scomm,nroots,indices_recv,PETSC_OWN_POINTER,&is_sc);CHKERRQ(ierr);
189*bac5b06fSFande Kong #if 0
190bdeb5f4eSFande Kong   ierr = ISView(is_sc,PETSC_NULL);CHKERRQ(ierr);
191a69400e0SFande Kong   ierr = MPI_Barrier(gcomm);CHKERRQ(ierr);
192a69400e0SFande Kong #endif
1932452736bSFande Kong   /*construct a parallel submatrix */
1942452736bSFande Kong   ierr = PetscCalloc1(1,&smat);CHKERRQ(ierr);
195a69400e0SFande Kong #if 0
196a69400e0SFande Kong   ierr = ISView(allis_sc,PETSC_NULL);CHKERRQ(ierr);
197a69400e0SFande Kong   ierr = MPI_Barrier(gcomm);CHKERRQ(ierr);
198a69400e0SFande Kong   //ierr = ISView(is_sc,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
199a69400e0SFande Kong   MPI_Comm   comm1, comm2;
200a69400e0SFande Kong   ierr = PetscObjectGetComm((PetscObject)is_sc,&comm1);CHKERRQ(ierr);
201a69400e0SFande Kong   ierr = PetscObjectGetComm((PetscObject)allis_sc,&comm2);CHKERRQ(ierr);
202a69400e0SFande Kong   /*ierr = PetscCommDuplicate(comm1,&comm2,NULL);CHKERRQ(ierr);*/
203a69400e0SFande Kong   /*ierr = MPI_Comm_dup(comm1,&comm2);CHKERRQ(ierr);*/
204a69400e0SFande Kong   ierr = MPI_Comm_compare(comm2,comm1,&issamecomm);CHKERRQ(ierr);
205a69400e0SFande Kong   if(issamecomm == MPI_IDENT){
206a69400e0SFande Kong     ierr=PetscPrintf(gcomm,"the same communicator \n");CHKERRQ(ierr);
207a69400e0SFande Kong   }else{
208a69400e0SFande Kong   	ierr=PetscPrintf(gcomm,"different communicator \n");CHKERRQ(ierr);
209a69400e0SFande Kong   }
210a69400e0SFande Kong #endif
211*bac5b06fSFande Kong #if 0
212a69400e0SFande Kong   ierr = ISView(is_sc,PETSC_NULL);CHKERRQ(ierr);
213a69400e0SFande Kong #endif
214*bac5b06fSFande Kong   /*Get square sub matrices */
215*bac5b06fSFande Kong   ierr = MatGetSubMatricesMPI(mat,1,&is_sc,&is_sc,MAT_INITIAL_MATRIX,&smat);CHKERRQ(ierr);
2162452736bSFande Kong   /* we do not need them any more */
2172452736bSFande Kong   ierr = ISDestroy(&allis_sc);CHKERRQ(ierr);
218*bac5b06fSFande Kong #if 0
219a69400e0SFande Kong   ierr = MatView(smat[0],PETSC_NULL);CHKERRQ(ierr);
220*bac5b06fSFande Kong   ierr = MPI_Barrier(gcomm);CHKERRQ(ierr);
221a69400e0SFande Kong #endif
2222452736bSFande Kong   /*create a partitioner to repartition the sub-matrix*/
2232452736bSFande Kong   ierr = MatPartitioningCreate(scomm,&part);CHKERRQ(ierr);
2242452736bSFande Kong   ierr = MatPartitioningSetAdjacency(part,smat[0]);CHKERRQ(ierr);
2252452736bSFande Kong #if PETSC_HAVE_PARMETIS
2262452736bSFande Kong   /* if there exists a ParMETIS installation, we try to use ParMETIS
2272452736bSFande Kong    * because a repartition routine possibly work better
2282452736bSFande Kong    * */
2292452736bSFande Kong   ierr = MatPartitioningSetType(part,MATPARTITIONINGPARMETIS);CHKERRQ(ierr);
2302452736bSFande Kong   /*try to use reparition function, instead of partition function */
2312452736bSFande Kong   ierr = MatPartitioningParmetisSetRepartition(part);CHKERRQ(ierr);
2322452736bSFande Kong #else
2332452736bSFande Kong   /*we at least provide a default partitioner to rebalance the computation  */
2342452736bSFande Kong   ierr = MatPartitioningSetType(part,MATPARTITIONINGAVERAGE);CHKERRQ(ierr);
2352452736bSFande Kong #endif
2362452736bSFande Kong   /*user can pick up any partitioner by using an option*/
2372452736bSFande Kong   ierr = MatPartitioningSetFromOptions(part);CHKERRQ(ierr);
2382452736bSFande Kong   /* apply partition */
2392452736bSFande Kong   ierr = MatPartitioningApply(part,&partitioning);CHKERRQ(ierr);
2402452736bSFande Kong   ierr = MatPartitioningDestroy(&part);CHKERRQ(ierr);
241ca2fc57aSFande Kong   ierr = MatDestroy(&(smat[0]));CHKERRQ(ierr);
2422452736bSFande Kong   ierr = PetscFree(smat);CHKERRQ(ierr);
243*bac5b06fSFande Kong #if 0
244bdeb5f4eSFande Kong   ierr = ISView(partitioning,PETSC_NULL);CHKERRQ(ierr);
245*bac5b06fSFande Kong   ierr = MPI_Barrier(gcomm);CHKERRQ(ierr);
246bdeb5f4eSFande Kong #endif
2472452736bSFande Kong   /* get local rows including  overlap */
248bdeb5f4eSFande Kong   ierr = ISBuildTwoSided(partitioning,is_sc,is);CHKERRQ(ierr);
249*bac5b06fSFande Kong #if 0
250*bac5b06fSFande Kong   ierr = ISView(*is,PETSC_NULL);CHKERRQ(ierr);
251*bac5b06fSFande Kong   ierr = MPI_Barrier(gcomm);CHKERRQ(ierr);
252*bac5b06fSFande Kong #endif
253bdeb5f4eSFande Kong   /* destroy */
254bdeb5f4eSFande Kong   ierr = ISDestroy(&is_sc);CHKERRQ(ierr);
2552452736bSFande Kong   PetscFunctionReturn(0);
2562452736bSFande Kong }
2572452736bSFande Kong 
2582452736bSFande Kong 
259