xref: /petsc/src/mat/impls/aij/mpi/fdmpiaij.c (revision 2205254efee3a00a594e5e2a3a70f74dcb40bc03)
1a64fbb32SBarry Smith 
2c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h>
3a64fbb32SBarry Smith 
4ab9863d7SBarry Smith extern PetscErrorCode MatCreateColmap_MPIAIJ_Private(Mat);
5a64fbb32SBarry Smith 
64a2ae208SSatish Balay #undef __FUNCT__
74a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringCreate_MPIAIJ"
8dfbe8321SBarry Smith PetscErrorCode MatFDColoringCreate_MPIAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c)
9a64fbb32SBarry Smith {
106eaac0f3SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*)mat->data;
116849ba73SBarry Smith   PetscErrorCode         ierr;
12b1d57f15SBarry Smith   PetscMPIInt            size,*ncolsonproc,*disp,nn;
131a83f524SJed Brown   PetscInt               i,n,nrows,j,k,m,ncols,col;
141a83f524SJed Brown   const PetscInt         *is,*A_ci,*A_cj,*B_ci,*B_cj,*rows = 0;
151a83f524SJed Brown   PetscInt               nis = iscoloring->n,nctot,*cols;
16f6d58c54SBarry Smith   PetscInt               *rowhit,M,cstart,cend,colb;
17b1d57f15SBarry Smith   PetscInt               *columnsforrow,l;
18b9617806SBarry Smith   IS                     *isa;
19ace3abfcSBarry Smith   PetscBool              done,flg;
20992144d0SBarry Smith   ISLocalToGlobalMapping map   = mat->cmap->mapping;
21e60d7cb6SSatish Balay   PetscInt               *ltog = (map ? map->indices : (PetscInt*) PETSC_NULL),ctype=c->ctype;
22a64fbb32SBarry Smith 
233a40ed3dSBarry Smith   PetscFunctionBegin;
24e7e72b3dSBarry Smith   if (!mat->assembled) SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be assembled first; MatAssemblyBegin/End();");
25e7e72b3dSBarry Smith   if (ctype == IS_COLORING_GHOSTED && !map) SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_ARG_INCOMP,"When using ghosted differencing matrix must have local to global mapping provided with MatSetLocalToGlobalMapping");
26522c5e43SBarry Smith 
27b9617806SBarry Smith   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
283acb8795SBarry Smith 
29f6d58c54SBarry Smith   M         = mat->rmap->n;
30f6d58c54SBarry Smith   cstart    = mat->cmap->rstart;
31f6d58c54SBarry Smith   cend      = mat->cmap->rend;
32f6d58c54SBarry Smith   c->M      = mat->rmap->N;         /* set the global rows and columns and local rows */
33f6d58c54SBarry Smith   c->N      = mat->cmap->N;
34f6d58c54SBarry Smith   c->m      = mat->rmap->n;
35f6d58c54SBarry Smith   c->rstart = mat->rmap->rstart;
36005c665bSBarry Smith 
37a64fbb32SBarry Smith   c->ncolors = nis;
38b1d57f15SBarry Smith   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr);
39b1d57f15SBarry Smith   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr);
40b1d57f15SBarry Smith   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr);
41b1d57f15SBarry Smith   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->rows);CHKERRQ(ierr);
42b1d57f15SBarry Smith   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->columnsforrow);CHKERRQ(ierr);
4352e6d16bSBarry Smith   ierr       = PetscLogObjectMemory(c,5*nis*sizeof(PetscInt));CHKERRQ(ierr);
446eaac0f3SBarry Smith 
456eaac0f3SBarry Smith   /* Allow access to data structures of local part of matrix */
466eaac0f3SBarry Smith   if (!aij->colmap) {
47ab9863d7SBarry Smith     ierr = MatCreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr);
486eaac0f3SBarry Smith   }
493acb8795SBarry Smith   ierr = MatGetColumnIJ(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);CHKERRQ(ierr);
503acb8795SBarry Smith   ierr = MatGetColumnIJ(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);CHKERRQ(ierr);
516eaac0f3SBarry Smith 
52b1d57f15SBarry Smith   ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr);
53b1d57f15SBarry Smith   ierr = PetscMalloc((M+1)*sizeof(PetscInt),&columnsforrow);CHKERRQ(ierr);
546eaac0f3SBarry Smith 
55a64fbb32SBarry Smith   for (i=0; i<nis; i++) {
56b9b97703SBarry Smith     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
57a64fbb32SBarry Smith     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
58*2205254eSKarl Rupp 
59a64fbb32SBarry Smith     c->ncolumns[i] = n;
60a64fbb32SBarry Smith     if (n) {
61b1d57f15SBarry Smith       ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr);
6252e6d16bSBarry Smith       ierr = PetscLogObjectMemory(c,n*sizeof(PetscInt));CHKERRQ(ierr);
63b1d57f15SBarry Smith       ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr);
64a64fbb32SBarry Smith     } else {
65a64fbb32SBarry Smith       c->columns[i] = 0;
66a64fbb32SBarry Smith     }
67a64fbb32SBarry Smith 
688ee2e534SBarry Smith     if (ctype == IS_COLORING_GLOBAL) {
696eaac0f3SBarry Smith       /* Determine the total (parallel) number of columns of this color */
707adad957SLisandro Dalcin       ierr = MPI_Comm_size(((PetscObject)mat)->comm,&size);CHKERRQ(ierr);
71687f1162SBarry Smith       ierr = PetscMalloc2(size,PetscMPIInt,&ncolsonproc,size,PetscMPIInt,&disp);CHKERRQ(ierr);
72b8f8c88eSHong Zhang 
734dc2109aSBarry Smith       ierr = PetscMPIIntCast(n,&nn);CHKERRQ(ierr);
747adad957SLisandro Dalcin       ierr = MPI_Allgather(&nn,1,MPI_INT,ncolsonproc,1,MPI_INT,((PetscObject)mat)->comm);CHKERRQ(ierr);
75*2205254eSKarl Rupp       nctot = 0; for (j=0; j<size; j++) nctot += ncolsonproc[j];
763a7fca6bSBarry Smith       if (!nctot) {
77ae15b995SBarry Smith         ierr = PetscInfo(mat,"Coloring of matrix has some unneeded colors with no corresponding rows\n");CHKERRQ(ierr);
783a7fca6bSBarry Smith       }
796eaac0f3SBarry Smith 
806eaac0f3SBarry Smith       disp[0] = 0;
816eaac0f3SBarry Smith       for (j=1; j<size; j++) {
826eaac0f3SBarry Smith         disp[j] = disp[j-1] + ncolsonproc[j-1];
836eaac0f3SBarry Smith       }
846eaac0f3SBarry Smith 
856eaac0f3SBarry Smith       /* Get complete list of columns for color on each processor */
86b1d57f15SBarry Smith       ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr);
875d0c19d7SBarry Smith       ierr = MPI_Allgatherv((void*)is,n,MPIU_INT,cols,ncolsonproc,disp,MPIU_INT,((PetscObject)mat)->comm);CHKERRQ(ierr);
881d79065fSBarry Smith       ierr = PetscFree2(ncolsonproc,disp);CHKERRQ(ierr);
89b8f8c88eSHong Zhang     } else if (ctype == IS_COLORING_GHOSTED) {
90b8f8c88eSHong Zhang       /* Determine local number of columns of this color on this process, including ghost points */
91b8f8c88eSHong Zhang       nctot = n;
92b8f8c88eSHong Zhang       ierr  = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr);
93b8f8c88eSHong Zhang       ierr  = PetscMemcpy(cols,is,n*sizeof(PetscInt));CHKERRQ(ierr);
94f23aa3ddSBarry Smith     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not provided for this MatFDColoring type");
956eaac0f3SBarry Smith 
966eaac0f3SBarry Smith     /*
976eaac0f3SBarry Smith        Mark all rows affect by these columns
986eaac0f3SBarry Smith     */
99b8f8c88eSHong Zhang     /* Temporary option to allow for debugging/testing */
10090d69ab7SBarry Smith     flg  = PETSC_FALSE;
101acfcf0e5SJed Brown     ierr = PetscOptionsGetBool(PETSC_NULL,"-matfdcoloring_slow",&flg,PETSC_NULL);CHKERRQ(ierr);
102f158e583SBarry Smith     if (!flg) { /*-----------------------------------------------------------------------------*/
103f158e583SBarry Smith       /* crude, fast version */
104b1d57f15SBarry Smith       ierr = PetscMemzero(rowhit,M*sizeof(PetscInt));CHKERRQ(ierr);
105a64fbb32SBarry Smith       /* loop over columns*/
1066eaac0f3SBarry Smith       for (j=0; j<nctot; j++) {
107b8f8c88eSHong Zhang         if (ctype == IS_COLORING_GHOSTED) {
108b8f8c88eSHong Zhang           col = ltog[cols[j]];
109b8f8c88eSHong Zhang         } else {
1106eaac0f3SBarry Smith           col  = cols[j];
111b8f8c88eSHong Zhang         }
1126eaac0f3SBarry Smith         if (col >= cstart && col < cend) {
1136eaac0f3SBarry Smith           /* column is in diagonal block of matrix */
1146eaac0f3SBarry Smith           rows = A_cj + A_ci[col-cstart];
1156eaac0f3SBarry Smith           m    = A_ci[col-cstart+1] - A_ci[col-cstart];
1166eaac0f3SBarry Smith         } else {
117aa482453SBarry Smith #if defined(PETSC_USE_CTABLE)
118cb9801acSJed Brown           ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr);
119fa46199cSSatish Balay           colb--;
120b3d2dc96SSatish Balay #else
1216eaac0f3SBarry Smith           colb = aij->colmap[col] - 1;
122b3d2dc96SSatish Balay #endif
1236eaac0f3SBarry Smith           if (colb == -1) {
1246eaac0f3SBarry Smith             m = 0;
1256eaac0f3SBarry Smith           } else {
1266eaac0f3SBarry Smith             rows = B_cj + B_ci[colb];
1276eaac0f3SBarry Smith             m    = B_ci[colb+1] - B_ci[colb];
1286eaac0f3SBarry Smith           }
1296eaac0f3SBarry Smith         }
130a64fbb32SBarry Smith         /* loop over columns marking them in rowhit */
131a64fbb32SBarry Smith         for (k=0; k<m; k++) {
132a64fbb32SBarry Smith           rowhit[*rows++] = col + 1;
133a64fbb32SBarry Smith         }
134a64fbb32SBarry Smith       }
1356eaac0f3SBarry Smith 
136a64fbb32SBarry Smith       /* count the number of hits */
137a64fbb32SBarry Smith       nrows = 0;
1386eaac0f3SBarry Smith       for (j=0; j<M; j++) {
139a64fbb32SBarry Smith         if (rowhit[j]) nrows++;
140a64fbb32SBarry Smith       }
141a64fbb32SBarry Smith       c->nrows[i] = nrows;
142b1d57f15SBarry Smith       ierr        = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr);
143b1d57f15SBarry Smith       ierr        = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr);
14452e6d16bSBarry Smith       ierr        = PetscLogObjectMemory(c,2*(nrows+1)*sizeof(PetscInt));CHKERRQ(ierr);
145a64fbb32SBarry Smith       nrows       = 0;
1466eaac0f3SBarry Smith       for (j=0; j<M; j++) {
147a64fbb32SBarry Smith         if (rowhit[j]) {
148a64fbb32SBarry Smith           c->rows[i][nrows]          = j;
149a64fbb32SBarry Smith           c->columnsforrow[i][nrows] = rowhit[j] - 1;
150a64fbb32SBarry Smith           nrows++;
151a64fbb32SBarry Smith         }
152a64fbb32SBarry Smith       }
153a64fbb32SBarry Smith     } else { /*-------------------------------------------------------------------------------*/
154f158e583SBarry Smith       /* slow version, using rowhit as a linked list */
155b1d57f15SBarry Smith       PetscInt currentcol,fm,mfm;
1566eaac0f3SBarry Smith       rowhit[M] = M;
157a64fbb32SBarry Smith       nrows     = 0;
158a64fbb32SBarry Smith       /* loop over columns*/
1596eaac0f3SBarry Smith       for (j=0; j<nctot; j++) {
160b8f8c88eSHong Zhang         if (ctype == IS_COLORING_GHOSTED) {
161b8f8c88eSHong Zhang           col = ltog[cols[j]];
162b8f8c88eSHong Zhang         } else {
1636eaac0f3SBarry Smith           col = cols[j];
164b8f8c88eSHong Zhang         }
1656eaac0f3SBarry Smith         if (col >= cstart && col < cend) {
1666eaac0f3SBarry Smith           /* column is in diagonal block of matrix */
1676eaac0f3SBarry Smith           rows = A_cj + A_ci[col-cstart];
1686eaac0f3SBarry Smith           m    = A_ci[col-cstart+1] - A_ci[col-cstart];
1696eaac0f3SBarry Smith         } else {
170aa482453SBarry Smith #if defined(PETSC_USE_CTABLE)
1710f5bd95cSBarry Smith           ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr);
172fa46199cSSatish Balay           colb--;
173b3d2dc96SSatish Balay #else
1746eaac0f3SBarry Smith           colb = aij->colmap[col] - 1;
175b3d2dc96SSatish Balay #endif
1766eaac0f3SBarry Smith           if (colb == -1) {
1776eaac0f3SBarry Smith             m = 0;
1786eaac0f3SBarry Smith           } else {
1796eaac0f3SBarry Smith             rows = B_cj + B_ci[colb];
1806eaac0f3SBarry Smith             m    = B_ci[colb+1] - B_ci[colb];
1816eaac0f3SBarry Smith           }
1826eaac0f3SBarry Smith         }
183b8f8c88eSHong Zhang 
184a64fbb32SBarry Smith         /* loop over columns marking them in rowhit */
1856eaac0f3SBarry Smith         fm = M;    /* fm points to first entry in linked list */
186a64fbb32SBarry Smith         for (k=0; k<m; k++) {
187a64fbb32SBarry Smith           currentcol = *rows++;
188a64fbb32SBarry Smith           /* is it already in the list? */
189a64fbb32SBarry Smith           do {
190a64fbb32SBarry Smith             mfm = fm;
191a64fbb32SBarry Smith             fm  = rowhit[fm];
192a64fbb32SBarry Smith           } while (fm < currentcol);
193a64fbb32SBarry Smith           /* not in list so add it */
194a64fbb32SBarry Smith           if (fm != currentcol) {
195a64fbb32SBarry Smith             nrows++;
196a64fbb32SBarry Smith             columnsforrow[currentcol] = col;
197a64fbb32SBarry Smith             /* next three lines insert new entry into linked list */
198a64fbb32SBarry Smith             rowhit[mfm]        = currentcol;
199a64fbb32SBarry Smith             rowhit[currentcol] = fm;
200a64fbb32SBarry Smith             fm                 = currentcol;
201a64fbb32SBarry Smith             /* fm points to present position in list since we know the columns are sorted */
202f23aa3ddSBarry Smith           } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid coloring of matrix detected");
203a64fbb32SBarry Smith         }
204a64fbb32SBarry Smith       }
205a64fbb32SBarry Smith       c->nrows[i] = nrows;
206*2205254eSKarl Rupp 
207b1d57f15SBarry Smith       ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr);
208b1d57f15SBarry Smith       ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr);
20952e6d16bSBarry Smith       ierr = PetscLogObjectMemory(c,(nrows+1)*sizeof(PetscInt));CHKERRQ(ierr);
210a64fbb32SBarry Smith       /* now store the linked list of rows into c->rows[i] */
211a64fbb32SBarry Smith       nrows = 0;
2126eaac0f3SBarry Smith       fm    = rowhit[M];
213a64fbb32SBarry Smith       do {
214a64fbb32SBarry Smith         c->rows[i][nrows]            = fm;
215a64fbb32SBarry Smith         c->columnsforrow[i][nrows++] = columnsforrow[fm];
216a64fbb32SBarry Smith         fm                           = rowhit[fm];
2176eaac0f3SBarry Smith       } while (fm < M);
2186eaac0f3SBarry Smith     } /* ---------------------------------------------------------------------------------------*/
219606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
2206eaac0f3SBarry Smith   }
22130b34957SBarry Smith 
22230b34957SBarry Smith   /* Optimize by adding the vscale, and scaleforrow[][] fields */
22330b34957SBarry Smith   /*
22430b34957SBarry Smith        vscale will contain the "diagonal" on processor scalings followed by the off processor
22530b34957SBarry Smith   */
2268ee2e534SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
227d0f46423SBarry Smith     ierr = VecCreateGhost(((PetscObject)mat)->comm,aij->A->rmap->n,PETSC_DETERMINE,aij->B->cmap->n,aij->garray,&c->vscale);CHKERRQ(ierr);
228b1d57f15SBarry Smith     ierr = PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);CHKERRQ(ierr);
22930b34957SBarry Smith     for (k=0; k<c->ncolors; k++) {
230b1d57f15SBarry Smith       ierr = PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);CHKERRQ(ierr);
23130b34957SBarry Smith       for (l=0; l<c->nrows[k]; l++) {
23230b34957SBarry Smith         col = c->columnsforrow[k][l];
23330b34957SBarry Smith         if (col >= cstart && col < cend) {
23430b34957SBarry Smith           /* column is in diagonal block of matrix */
23530b34957SBarry Smith           colb = col - cstart;
23630b34957SBarry Smith         } else {
23730b34957SBarry Smith           /* column  is in "off-processor" part */
23830b34957SBarry Smith #if defined(PETSC_USE_CTABLE)
23930b34957SBarry Smith           ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr);
24030b34957SBarry Smith           colb--;
24130b34957SBarry Smith #else
24230b34957SBarry Smith           colb = aij->colmap[col] - 1;
24330b34957SBarry Smith #endif
24430b34957SBarry Smith           colb += cend - cstart;
24530b34957SBarry Smith         }
24630b34957SBarry Smith         c->vscaleforrow[k][l] = colb;
24730b34957SBarry Smith       }
24830b34957SBarry Smith     }
249b8f8c88eSHong Zhang   } else if (ctype == IS_COLORING_GHOSTED) {
250b8f8c88eSHong Zhang     /* Get gtol mapping */
251d0f46423SBarry Smith     PetscInt N = mat->cmap->N, *gtol;
252b8f8c88eSHong Zhang     ierr = PetscMalloc((N+1)*sizeof(PetscInt),&gtol);CHKERRQ(ierr);
253b8f8c88eSHong Zhang     for (i=0; i<N; i++) gtol[i] = -1;
254b8f8c88eSHong Zhang     for (i=0; i<map->n; i++) gtol[ltog[i]] = i;
255b8f8c88eSHong Zhang 
256b8f8c88eSHong Zhang     c->vscale = 0; /* will be created in MatFDColoringApply() */
257b8f8c88eSHong Zhang     ierr      = PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);CHKERRQ(ierr);
258b8f8c88eSHong Zhang     for (k=0; k<c->ncolors; k++) {
259b8f8c88eSHong Zhang       ierr = PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);CHKERRQ(ierr);
260b8f8c88eSHong Zhang       for (l=0; l<c->nrows[k]; l++) {
261b8f8c88eSHong Zhang         col = c->columnsforrow[k][l];      /* global column index */
262*2205254eSKarl Rupp 
263b8f8c88eSHong Zhang         c->vscaleforrow[k][l] = gtol[col]; /* local column index */
264b8f8c88eSHong Zhang       }
265b8f8c88eSHong Zhang     }
266b8f8c88eSHong Zhang     ierr = PetscFree(gtol);CHKERRQ(ierr);
267b8f8c88eSHong Zhang   }
268b9617806SBarry Smith   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
26930b34957SBarry Smith 
270606d414cSSatish Balay   ierr = PetscFree(rowhit);CHKERRQ(ierr);
271606d414cSSatish Balay   ierr = PetscFree(columnsforrow);CHKERRQ(ierr);
2723acb8795SBarry Smith   ierr = MatRestoreColumnIJ(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);CHKERRQ(ierr);
2733acb8795SBarry Smith   ierr = MatRestoreColumnIJ(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);CHKERRQ(ierr);
2743a40ed3dSBarry Smith   PetscFunctionReturn(0);
275a64fbb32SBarry Smith }
276a64fbb32SBarry Smith 
277b9617806SBarry Smith 
278b9617806SBarry Smith 
279b9617806SBarry Smith 
280b9617806SBarry Smith 
281b9617806SBarry Smith 
282