xref: /petsc/src/mat/impls/aij/mpi/fdmpiaij.c (revision 3bb1ff401821b9e2ae019d3e61bc8ab4bd4e59d5)
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;
14afcb2eb5SJed Brown   const PetscInt         *is,*A_ci,*A_cj,*B_ci,*B_cj,*rows = 0,*ltog;
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;
21afcb2eb5SJed Brown   PetscInt               ctype=c->ctype;
22a64fbb32SBarry Smith 
233a40ed3dSBarry Smith   PetscFunctionBegin;
24ce94432eSBarry Smith   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Matrix must be assembled first; MatAssemblyBegin/End();");
25ce94432eSBarry Smith   if (ctype == IS_COLORING_GHOSTED && !map) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_INCOMP,"When using ghosted differencing matrix must have local to global mapping provided with MatSetLocalToGlobalMapping");
26522c5e43SBarry Smith 
27afcb2eb5SJed Brown   if (map) {ierr = ISLocalToGlobalMappingGetIndices(map,&ltog);CHKERRQ(ierr);}
28afcb2eb5SJed Brown   else     ltog = NULL;
29b9617806SBarry Smith   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
303acb8795SBarry Smith 
31f6d58c54SBarry Smith   M         = mat->rmap->n;
32f6d58c54SBarry Smith   cstart    = mat->cmap->rstart;
33f6d58c54SBarry Smith   cend      = mat->cmap->rend;
34f6d58c54SBarry Smith   c->M      = mat->rmap->N;         /* set the global rows and columns and local rows */
35f6d58c54SBarry Smith   c->N      = mat->cmap->N;
36f6d58c54SBarry Smith   c->m      = mat->rmap->n;
37f6d58c54SBarry Smith   c->rstart = mat->rmap->rstart;
38005c665bSBarry Smith 
39a64fbb32SBarry Smith   c->ncolors = nis;
40b1d57f15SBarry Smith   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr);
41b1d57f15SBarry Smith   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr);
42b1d57f15SBarry Smith   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr);
43b1d57f15SBarry Smith   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->rows);CHKERRQ(ierr);
44b1d57f15SBarry Smith   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->columnsforrow);CHKERRQ(ierr);
45*3bb1ff40SBarry Smith   ierr       = PetscLogObjectMemory((PetscObject)c,5*nis*sizeof(PetscInt));CHKERRQ(ierr);
466eaac0f3SBarry Smith 
476eaac0f3SBarry Smith   /* Allow access to data structures of local part of matrix */
486eaac0f3SBarry Smith   if (!aij->colmap) {
49ab9863d7SBarry Smith     ierr = MatCreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr);
506eaac0f3SBarry Smith   }
513acb8795SBarry Smith   ierr = MatGetColumnIJ(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);CHKERRQ(ierr);
523acb8795SBarry Smith   ierr = MatGetColumnIJ(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);CHKERRQ(ierr);
536eaac0f3SBarry Smith 
54b1d57f15SBarry Smith   ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr);
55b1d57f15SBarry Smith   ierr = PetscMalloc((M+1)*sizeof(PetscInt),&columnsforrow);CHKERRQ(ierr);
566eaac0f3SBarry Smith 
57a64fbb32SBarry Smith   for (i=0; i<nis; i++) {
58b9b97703SBarry Smith     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
59a64fbb32SBarry Smith     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
602205254eSKarl Rupp 
61a64fbb32SBarry Smith     c->ncolumns[i] = n;
62a64fbb32SBarry Smith     if (n) {
63b1d57f15SBarry Smith       ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr);
64*3bb1ff40SBarry Smith       ierr = PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));CHKERRQ(ierr);
65b1d57f15SBarry Smith       ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr);
66a64fbb32SBarry Smith     } else {
67a64fbb32SBarry Smith       c->columns[i] = 0;
68a64fbb32SBarry Smith     }
69a64fbb32SBarry Smith 
708ee2e534SBarry Smith     if (ctype == IS_COLORING_GLOBAL) {
716eaac0f3SBarry Smith       /* Determine the total (parallel) number of columns of this color */
72ce94432eSBarry Smith       ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);CHKERRQ(ierr);
73687f1162SBarry Smith       ierr = PetscMalloc2(size,PetscMPIInt,&ncolsonproc,size,PetscMPIInt,&disp);CHKERRQ(ierr);
74b8f8c88eSHong Zhang 
754dc2109aSBarry Smith       ierr  = PetscMPIIntCast(n,&nn);CHKERRQ(ierr);
76ce94432eSBarry Smith       ierr  = MPI_Allgather(&nn,1,MPI_INT,ncolsonproc,1,MPI_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
772205254eSKarl Rupp       nctot = 0; for (j=0; j<size; j++) nctot += ncolsonproc[j];
783a7fca6bSBarry Smith       if (!nctot) {
79ae15b995SBarry Smith         ierr = PetscInfo(mat,"Coloring of matrix has some unneeded colors with no corresponding rows\n");CHKERRQ(ierr);
803a7fca6bSBarry Smith       }
816eaac0f3SBarry Smith 
826eaac0f3SBarry Smith       disp[0] = 0;
836eaac0f3SBarry Smith       for (j=1; j<size; j++) {
846eaac0f3SBarry Smith         disp[j] = disp[j-1] + ncolsonproc[j-1];
856eaac0f3SBarry Smith       }
866eaac0f3SBarry Smith 
876eaac0f3SBarry Smith       /* Get complete list of columns for color on each processor */
88b1d57f15SBarry Smith       ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr);
89ce94432eSBarry Smith       ierr = MPI_Allgatherv((void*)is,n,MPIU_INT,cols,ncolsonproc,disp,MPIU_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
901d79065fSBarry Smith       ierr = PetscFree2(ncolsonproc,disp);CHKERRQ(ierr);
91b8f8c88eSHong Zhang     } else if (ctype == IS_COLORING_GHOSTED) {
92b8f8c88eSHong Zhang       /* Determine local number of columns of this color on this process, including ghost points */
93b8f8c88eSHong Zhang       nctot = n;
94b8f8c88eSHong Zhang       ierr  = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr);
95b8f8c88eSHong Zhang       ierr  = PetscMemcpy(cols,is,n*sizeof(PetscInt));CHKERRQ(ierr);
96f23aa3ddSBarry Smith     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not provided for this MatFDColoring type");
976eaac0f3SBarry Smith 
986eaac0f3SBarry Smith     /*
996eaac0f3SBarry Smith        Mark all rows affect by these columns
1006eaac0f3SBarry Smith     */
101b8f8c88eSHong Zhang     /* Temporary option to allow for debugging/testing */
10290d69ab7SBarry Smith     flg  = PETSC_FALSE;
1030298fd71SBarry Smith     ierr = PetscOptionsGetBool(NULL,"-matfdcoloring_slow",&flg,NULL);CHKERRQ(ierr);
104f158e583SBarry Smith     if (!flg) { /*-----------------------------------------------------------------------------*/
105f158e583SBarry Smith       /* crude, fast version */
106b1d57f15SBarry Smith       ierr = PetscMemzero(rowhit,M*sizeof(PetscInt));CHKERRQ(ierr);
107a64fbb32SBarry Smith       /* loop over columns*/
1086eaac0f3SBarry Smith       for (j=0; j<nctot; j++) {
109b8f8c88eSHong Zhang         if (ctype == IS_COLORING_GHOSTED) {
110b8f8c88eSHong Zhang           col = ltog[cols[j]];
111b8f8c88eSHong Zhang         } else {
1126eaac0f3SBarry Smith           col = cols[j];
113b8f8c88eSHong Zhang         }
1146eaac0f3SBarry Smith         if (col >= cstart && col < cend) {
1156eaac0f3SBarry Smith           /* column is in diagonal block of matrix */
1166eaac0f3SBarry Smith           rows = A_cj + A_ci[col-cstart];
1176eaac0f3SBarry Smith           m    = A_ci[col-cstart+1] - A_ci[col-cstart];
1186eaac0f3SBarry Smith         } else {
119aa482453SBarry Smith #if defined(PETSC_USE_CTABLE)
120cb9801acSJed Brown           ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr);
121fa46199cSSatish Balay           colb--;
122b3d2dc96SSatish Balay #else
1236eaac0f3SBarry Smith           colb = aij->colmap[col] - 1;
124b3d2dc96SSatish Balay #endif
1256eaac0f3SBarry Smith           if (colb == -1) {
1266eaac0f3SBarry Smith             m = 0;
1276eaac0f3SBarry Smith           } else {
1286eaac0f3SBarry Smith             rows = B_cj + B_ci[colb];
1296eaac0f3SBarry Smith             m    = B_ci[colb+1] - B_ci[colb];
1306eaac0f3SBarry Smith           }
1316eaac0f3SBarry Smith         }
132a64fbb32SBarry Smith         /* loop over columns marking them in rowhit */
133a64fbb32SBarry Smith         for (k=0; k<m; k++) {
134a64fbb32SBarry Smith           rowhit[*rows++] = col + 1;
135a64fbb32SBarry Smith         }
136a64fbb32SBarry Smith       }
1376eaac0f3SBarry Smith 
138a64fbb32SBarry Smith       /* count the number of hits */
139a64fbb32SBarry Smith       nrows = 0;
1406eaac0f3SBarry Smith       for (j=0; j<M; j++) {
141a64fbb32SBarry Smith         if (rowhit[j]) nrows++;
142a64fbb32SBarry Smith       }
143a64fbb32SBarry Smith       c->nrows[i] = nrows;
144b1d57f15SBarry Smith       ierr        = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr);
145b1d57f15SBarry Smith       ierr        = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr);
146*3bb1ff40SBarry Smith       ierr        = PetscLogObjectMemory((PetscObject)c,2*(nrows+1)*sizeof(PetscInt));CHKERRQ(ierr);
147a64fbb32SBarry Smith       nrows       = 0;
1486eaac0f3SBarry Smith       for (j=0; j<M; j++) {
149a64fbb32SBarry Smith         if (rowhit[j]) {
150a64fbb32SBarry Smith           c->rows[i][nrows]          = j;
151a64fbb32SBarry Smith           c->columnsforrow[i][nrows] = rowhit[j] - 1;
152a64fbb32SBarry Smith           nrows++;
153a64fbb32SBarry Smith         }
154a64fbb32SBarry Smith       }
155a64fbb32SBarry Smith     } else { /*-------------------------------------------------------------------------------*/
156f158e583SBarry Smith       /* slow version, using rowhit as a linked list */
157b1d57f15SBarry Smith       PetscInt currentcol,fm,mfm;
1586eaac0f3SBarry Smith       rowhit[M] = M;
159a64fbb32SBarry Smith       nrows     = 0;
160a64fbb32SBarry Smith       /* loop over columns*/
1616eaac0f3SBarry Smith       for (j=0; j<nctot; j++) {
162b8f8c88eSHong Zhang         if (ctype == IS_COLORING_GHOSTED) {
163b8f8c88eSHong Zhang           col = ltog[cols[j]];
164b8f8c88eSHong Zhang         } else {
1656eaac0f3SBarry Smith           col = cols[j];
166b8f8c88eSHong Zhang         }
1676eaac0f3SBarry Smith         if (col >= cstart && col < cend) {
1686eaac0f3SBarry Smith           /* column is in diagonal block of matrix */
1696eaac0f3SBarry Smith           rows = A_cj + A_ci[col-cstart];
1706eaac0f3SBarry Smith           m    = A_ci[col-cstart+1] - A_ci[col-cstart];
1716eaac0f3SBarry Smith         } else {
172aa482453SBarry Smith #if defined(PETSC_USE_CTABLE)
1730f5bd95cSBarry Smith           ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr);
174fa46199cSSatish Balay           colb--;
175b3d2dc96SSatish Balay #else
1766eaac0f3SBarry Smith           colb = aij->colmap[col] - 1;
177b3d2dc96SSatish Balay #endif
1786eaac0f3SBarry Smith           if (colb == -1) {
1796eaac0f3SBarry Smith             m = 0;
1806eaac0f3SBarry Smith           } else {
1816eaac0f3SBarry Smith             rows = B_cj + B_ci[colb];
1826eaac0f3SBarry Smith             m    = B_ci[colb+1] - B_ci[colb];
1836eaac0f3SBarry Smith           }
1846eaac0f3SBarry Smith         }
185b8f8c88eSHong Zhang 
186a64fbb32SBarry Smith         /* loop over columns marking them in rowhit */
1876eaac0f3SBarry Smith         fm = M;    /* fm points to first entry in linked list */
188a64fbb32SBarry Smith         for (k=0; k<m; k++) {
189a64fbb32SBarry Smith           currentcol = *rows++;
190a64fbb32SBarry Smith           /* is it already in the list? */
191a64fbb32SBarry Smith           do {
192a64fbb32SBarry Smith             mfm = fm;
193a64fbb32SBarry Smith             fm  = rowhit[fm];
194a64fbb32SBarry Smith           } while (fm < currentcol);
195a64fbb32SBarry Smith           /* not in list so add it */
196a64fbb32SBarry Smith           if (fm != currentcol) {
197a64fbb32SBarry Smith             nrows++;
198a64fbb32SBarry Smith             columnsforrow[currentcol] = col;
199a64fbb32SBarry Smith             /* next three lines insert new entry into linked list */
200a64fbb32SBarry Smith             rowhit[mfm]        = currentcol;
201a64fbb32SBarry Smith             rowhit[currentcol] = fm;
202a64fbb32SBarry Smith             fm                 = currentcol;
203a64fbb32SBarry Smith             /* fm points to present position in list since we know the columns are sorted */
204f23aa3ddSBarry Smith           } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid coloring of matrix detected");
205a64fbb32SBarry Smith         }
206a64fbb32SBarry Smith       }
207a64fbb32SBarry Smith       c->nrows[i] = nrows;
2082205254eSKarl Rupp 
209b1d57f15SBarry Smith       ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr);
210b1d57f15SBarry Smith       ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr);
211*3bb1ff40SBarry Smith       ierr = PetscLogObjectMemory((PetscObject)c,(nrows+1)*sizeof(PetscInt));CHKERRQ(ierr);
212a64fbb32SBarry Smith       /* now store the linked list of rows into c->rows[i] */
213a64fbb32SBarry Smith       nrows = 0;
2146eaac0f3SBarry Smith       fm    = rowhit[M];
215a64fbb32SBarry Smith       do {
216a64fbb32SBarry Smith         c->rows[i][nrows]            = fm;
217a64fbb32SBarry Smith         c->columnsforrow[i][nrows++] = columnsforrow[fm];
218a64fbb32SBarry Smith         fm                           = rowhit[fm];
2196eaac0f3SBarry Smith       } while (fm < M);
2206eaac0f3SBarry Smith     } /* ---------------------------------------------------------------------------------------*/
221606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
2226eaac0f3SBarry Smith   }
22330b34957SBarry Smith 
22430b34957SBarry Smith   /* Optimize by adding the vscale, and scaleforrow[][] fields */
22530b34957SBarry Smith   /*
22630b34957SBarry Smith        vscale will contain the "diagonal" on processor scalings followed by the off processor
22730b34957SBarry Smith   */
2288ee2e534SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
229ce94432eSBarry Smith     ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),aij->A->rmap->n,PETSC_DETERMINE,aij->B->cmap->n,aij->garray,&c->vscale);CHKERRQ(ierr);
230b1d57f15SBarry Smith     ierr = PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);CHKERRQ(ierr);
23130b34957SBarry Smith     for (k=0; k<c->ncolors; k++) {
232b1d57f15SBarry Smith       ierr = PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);CHKERRQ(ierr);
23330b34957SBarry Smith       for (l=0; l<c->nrows[k]; l++) {
23430b34957SBarry Smith         col = c->columnsforrow[k][l];
23530b34957SBarry Smith         if (col >= cstart && col < cend) {
23630b34957SBarry Smith           /* column is in diagonal block of matrix */
23730b34957SBarry Smith           colb = col - cstart;
23830b34957SBarry Smith         } else {
23930b34957SBarry Smith           /* column  is in "off-processor" part */
24030b34957SBarry Smith #if defined(PETSC_USE_CTABLE)
24130b34957SBarry Smith           ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr);
24230b34957SBarry Smith           colb--;
24330b34957SBarry Smith #else
24430b34957SBarry Smith           colb = aij->colmap[col] - 1;
24530b34957SBarry Smith #endif
24630b34957SBarry Smith           colb += cend - cstart;
24730b34957SBarry Smith         }
24830b34957SBarry Smith         c->vscaleforrow[k][l] = colb;
24930b34957SBarry Smith       }
25030b34957SBarry Smith     }
251b8f8c88eSHong Zhang   } else if (ctype == IS_COLORING_GHOSTED) {
252b8f8c88eSHong Zhang     /* Get gtol mapping */
253afcb2eb5SJed Brown     PetscInt N = mat->cmap->N,nlocal,*gtol;
254b8f8c88eSHong Zhang     ierr = PetscMalloc((N+1)*sizeof(PetscInt),&gtol);CHKERRQ(ierr);
255b8f8c88eSHong Zhang     for (i=0; i<N; i++) gtol[i] = -1;
256afcb2eb5SJed Brown     ierr = ISLocalToGlobalMappingGetSize(map,&nlocal);CHKERRQ(ierr);
257afcb2eb5SJed Brown     for (i=0; i<nlocal; i++) gtol[ltog[i]] = i;
258b8f8c88eSHong Zhang 
259b8f8c88eSHong Zhang     c->vscale = 0; /* will be created in MatFDColoringApply() */
260b8f8c88eSHong Zhang     ierr      = PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);CHKERRQ(ierr);
261b8f8c88eSHong Zhang     for (k=0; k<c->ncolors; k++) {
262b8f8c88eSHong Zhang       ierr = PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);CHKERRQ(ierr);
263b8f8c88eSHong Zhang       for (l=0; l<c->nrows[k]; l++) {
264b8f8c88eSHong Zhang         col = c->columnsforrow[k][l];      /* global column index */
2652205254eSKarl Rupp 
266b8f8c88eSHong Zhang         c->vscaleforrow[k][l] = gtol[col]; /* local column index */
267b8f8c88eSHong Zhang       }
268b8f8c88eSHong Zhang     }
269b8f8c88eSHong Zhang     ierr = PetscFree(gtol);CHKERRQ(ierr);
270b8f8c88eSHong Zhang   }
271b9617806SBarry Smith   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
27230b34957SBarry Smith 
273606d414cSSatish Balay   ierr = PetscFree(rowhit);CHKERRQ(ierr);
274606d414cSSatish Balay   ierr = PetscFree(columnsforrow);CHKERRQ(ierr);
2753acb8795SBarry Smith   ierr = MatRestoreColumnIJ(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);CHKERRQ(ierr);
2763acb8795SBarry Smith   ierr = MatRestoreColumnIJ(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);CHKERRQ(ierr);
277afcb2eb5SJed Brown   if (map) {ierr = ISLocalToGlobalMappingRestoreIndices(map,&ltog);CHKERRQ(ierr);}
2783a40ed3dSBarry Smith   PetscFunctionReturn(0);
279a64fbb32SBarry Smith }
280a64fbb32SBarry Smith 
281b9617806SBarry Smith 
282b9617806SBarry Smith 
283b9617806SBarry Smith 
284b9617806SBarry Smith 
285b9617806SBarry Smith 
286