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