xref: /petsc/src/mat/impls/aij/mpi/fdmpiaij.c (revision 8f7157ef067da8bcc67a53842181bca9de0f5d87)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2a64fbb32SBarry Smith 
36eaac0f3SBarry Smith #include "src/mat/impls/aij/mpi/mpiaij.h"
4a64fbb32SBarry Smith 
5dfbe8321SBarry Smith EXTERN PetscErrorCode CreateColmap_MPIAIJ_Private(Mat);
6a64fbb32SBarry Smith 
74a2ae208SSatish Balay #undef __FUNCT__
84a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringCreate_MPIAIJ"
9dfbe8321SBarry Smith PetscErrorCode MatFDColoringCreate_MPIAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c)
10a64fbb32SBarry Smith {
116eaac0f3SBarry Smith   Mat_MPIAIJ            *aij = (Mat_MPIAIJ*)mat->data;
126849ba73SBarry Smith   PetscErrorCode        ierr;
13b1d57f15SBarry Smith   PetscMPIInt           size,*ncolsonproc,*disp,nn;
14b1d57f15SBarry Smith   PetscInt              i,*is,n,nrows,j,k,m,*rows = 0,*A_ci,*A_cj,ncols,col;
15b1d57f15SBarry Smith   PetscInt              nis = iscoloring->n,nctot,*cols,*B_ci,*B_cj;
16899cda47SBarry Smith   PetscInt              *rowhit,M = mat->rmap.n,cstart = mat->cmap.rstart,cend = mat->cmap.rend,colb;
17b1d57f15SBarry Smith   PetscInt              *columnsforrow,l;
18b9617806SBarry Smith   IS                    *isa;
19f1af5d2fSBarry Smith   PetscTruth             done,flg;
20b8f8c88eSHong Zhang   ISLocalToGlobalMapping map = mat->mapping;
217cc688d7SBarry Smith   PetscInt               *ltog = (map ? map->indices : 0) ,ctype=c->ctype;
22a64fbb32SBarry Smith 
233a40ed3dSBarry Smith   PetscFunctionBegin;
24522c5e43SBarry Smith   if (!mat->assembled) {
2529bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be assembled first; MatAssemblyBegin/End();");
26522c5e43SBarry Smith   }
277cc688d7SBarry Smith   if (ctype == IS_COLORING_GHOSTED && !map) SETERRQ(PETSC_ERR_ARG_INCOMP,"When using ghosted differencing matrix must have local to global mapping provided with MatSetLocalToGlobalMapping");
28522c5e43SBarry Smith 
29b9617806SBarry Smith   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
30899cda47SBarry Smith   c->M             = mat->rmap.N;  /* set the global rows and columns and local rows */
31899cda47SBarry Smith   c->N             = mat->cmap.N;
32899cda47SBarry Smith   c->m             = mat->rmap.n;
33899cda47SBarry Smith   c->rstart        = mat->rmap.rstart;
34005c665bSBarry Smith 
35a64fbb32SBarry Smith   c->ncolors       = nis;
36b1d57f15SBarry Smith   ierr             = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr);
37b1d57f15SBarry Smith   ierr             = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr);
38b1d57f15SBarry Smith   ierr             = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr);
39b1d57f15SBarry Smith   ierr             = PetscMalloc(nis*sizeof(PetscInt*),&c->rows);CHKERRQ(ierr);
40b1d57f15SBarry Smith   ierr             = PetscMalloc(nis*sizeof(PetscInt*),&c->columnsforrow);CHKERRQ(ierr);
4152e6d16bSBarry Smith   ierr = PetscLogObjectMemory(c,5*nis*sizeof(PetscInt));CHKERRQ(ierr);
426eaac0f3SBarry Smith 
436eaac0f3SBarry Smith   /* Allow access to data structures of local part of matrix */
446eaac0f3SBarry Smith   if (!aij->colmap) {
456eaac0f3SBarry Smith     ierr = CreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr);
466eaac0f3SBarry Smith   }
4743a90d84SBarry Smith   /*
4843a90d84SBarry Smith       Calls the _SeqAIJ() version of these routines to make sure it does not
4943a90d84SBarry Smith      get the reduced (by inodes) version of I and J
5043a90d84SBarry Smith   */
51*8f7157efSSatish Balay   ierr = MatGetColumnIJ_SeqAIJ(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);CHKERRQ(ierr);
52*8f7157efSSatish Balay   ierr = MatGetColumnIJ_SeqAIJ(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);
60a64fbb32SBarry Smith     c->ncolumns[i] = n;
61a64fbb32SBarry Smith     if (n) {
62b1d57f15SBarry Smith       ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr);
6352e6d16bSBarry Smith       ierr = PetscLogObjectMemory(c,n*sizeof(PetscInt));CHKERRQ(ierr);
64b1d57f15SBarry Smith       ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr);
65a64fbb32SBarry Smith     } else {
66a64fbb32SBarry Smith       c->columns[i]  = 0;
67a64fbb32SBarry Smith     }
68a64fbb32SBarry Smith 
698ee2e534SBarry Smith     if (ctype == IS_COLORING_GLOBAL){
706eaac0f3SBarry Smith       /* Determine the total (parallel) number of columns of this color */
71b8f8c88eSHong Zhang       ierr = MPI_Comm_size(mat->comm,&size);CHKERRQ(ierr);
72b8f8c88eSHong Zhang       ierr = PetscMalloc(2*size*sizeof(PetscInt*),&ncolsonproc);CHKERRQ(ierr);
73b8f8c88eSHong Zhang       disp = ncolsonproc + size;
74b8f8c88eSHong Zhang 
75b1d57f15SBarry Smith       nn   = (PetscMPIInt)n;
76b1d57f15SBarry Smith       ierr = MPI_Allgather(&nn,1,MPI_INT,ncolsonproc,1,MPI_INT,mat->comm);CHKERRQ(ierr);
776eaac0f3SBarry Smith       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);
89b1d57f15SBarry Smith       ierr = MPI_Allgatherv(is,n,MPIU_INT,cols,ncolsonproc,disp,MPIU_INT,mat->comm);CHKERRQ(ierr);
90b8f8c88eSHong Zhang       ierr = PetscFree(ncolsonproc);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);
96b8f8c88eSHong Zhang     } else {
97b8f8c88eSHong Zhang       SETERRQ(PETSC_ERR_SUP,"Not provided for this MatFDColoring type");
98b8f8c88eSHong Zhang     }
996eaac0f3SBarry Smith 
1006eaac0f3SBarry Smith     /*
1016eaac0f3SBarry Smith        Mark all rows affect by these columns
1026eaac0f3SBarry Smith     */
103b8f8c88eSHong Zhang     /* Temporary option to allow for debugging/testing */
104b8f8c88eSHong Zhang     ierr = PetscOptionsHasName(PETSC_NULL,"-matfdcoloring_slow",&flg);CHKERRQ(ierr);
105f158e583SBarry Smith     if (!flg) {/*-----------------------------------------------------------------------------*/
106f158e583SBarry Smith       /* crude, fast version */
107b1d57f15SBarry Smith       ierr = PetscMemzero(rowhit,M*sizeof(PetscInt));CHKERRQ(ierr);
108a64fbb32SBarry Smith       /* loop over columns*/
1096eaac0f3SBarry Smith       for (j=0; j<nctot; j++) {
110b8f8c88eSHong Zhang         if (ctype == IS_COLORING_GHOSTED) {
111b8f8c88eSHong Zhang           col = ltog[cols[j]];
112b8f8c88eSHong Zhang         } else {
1136eaac0f3SBarry Smith           col  = cols[j];
114b8f8c88eSHong Zhang         }
1156eaac0f3SBarry Smith         if (col >= cstart && col < cend) {
1166eaac0f3SBarry Smith           /* column is in diagonal block of matrix */
1176eaac0f3SBarry Smith           rows = A_cj + A_ci[col-cstart];
1186eaac0f3SBarry Smith           m    = A_ci[col-cstart+1] - A_ci[col-cstart];
1196eaac0f3SBarry Smith         } else {
120aa482453SBarry Smith #if defined (PETSC_USE_CTABLE)
1210f5bd95cSBarry Smith           ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr)
122fa46199cSSatish Balay 	  colb --;
123b3d2dc96SSatish Balay #else
1246eaac0f3SBarry Smith           colb = aij->colmap[col] - 1;
125b3d2dc96SSatish Balay #endif
1266eaac0f3SBarry Smith           if (colb == -1) {
1276eaac0f3SBarry Smith             m = 0;
1286eaac0f3SBarry Smith           } else {
1296eaac0f3SBarry Smith             rows = B_cj + B_ci[colb];
1306eaac0f3SBarry Smith             m    = B_ci[colb+1] - B_ci[colb];
1316eaac0f3SBarry Smith           }
1326eaac0f3SBarry Smith         }
133a64fbb32SBarry Smith         /* loop over columns marking them in rowhit */
134a64fbb32SBarry Smith         for (k=0; k<m; k++) {
135a64fbb32SBarry Smith           rowhit[*rows++] = col + 1;
136a64fbb32SBarry Smith         }
137a64fbb32SBarry Smith       }
1386eaac0f3SBarry Smith 
139a64fbb32SBarry Smith       /* count the number of hits */
140a64fbb32SBarry Smith       nrows = 0;
1416eaac0f3SBarry Smith       for (j=0; j<M; j++) {
142a64fbb32SBarry Smith         if (rowhit[j]) nrows++;
143a64fbb32SBarry Smith       }
144a64fbb32SBarry Smith       c->nrows[i]         = nrows;
145b1d57f15SBarry Smith       ierr                = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr);
146b1d57f15SBarry Smith       ierr                = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr);
14752e6d16bSBarry Smith       ierr = PetscLogObjectMemory(c,2*(nrows+1)*sizeof(PetscInt));CHKERRQ(ierr);
148a64fbb32SBarry Smith       nrows = 0;
1496eaac0f3SBarry Smith       for (j=0; j<M; j++) {
150a64fbb32SBarry Smith         if (rowhit[j]) {
151a64fbb32SBarry Smith           c->rows[i][nrows]           = j;
152a64fbb32SBarry Smith           c->columnsforrow[i][nrows] = rowhit[j] - 1;
153a64fbb32SBarry Smith           nrows++;
154a64fbb32SBarry Smith         }
155a64fbb32SBarry Smith       }
156a64fbb32SBarry Smith     } else {/*-------------------------------------------------------------------------------*/
157f158e583SBarry Smith       /* slow version, using rowhit as a linked list */
158b1d57f15SBarry Smith       PetscInt currentcol,fm,mfm;
1596eaac0f3SBarry Smith       rowhit[M] = M;
160a64fbb32SBarry Smith       nrows     = 0;
161a64fbb32SBarry Smith       /* loop over columns*/
1626eaac0f3SBarry Smith       for (j=0; j<nctot; j++) {
163b8f8c88eSHong Zhang         if (ctype == IS_COLORING_GHOSTED) {
164b8f8c88eSHong Zhang           col = ltog[cols[j]];
165b8f8c88eSHong Zhang         } else {
1666eaac0f3SBarry Smith           col  = cols[j];
167b8f8c88eSHong Zhang         }
1686eaac0f3SBarry Smith         if (col >= cstart && col < cend) {
1696eaac0f3SBarry Smith           /* column is in diagonal block of matrix */
1706eaac0f3SBarry Smith           rows = A_cj + A_ci[col-cstart];
1716eaac0f3SBarry Smith           m    = A_ci[col-cstart+1] - A_ci[col-cstart];
1726eaac0f3SBarry Smith         } else {
173aa482453SBarry Smith #if defined (PETSC_USE_CTABLE)
1740f5bd95cSBarry Smith           ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr);
175fa46199cSSatish Balay           colb --;
176b3d2dc96SSatish Balay #else
1776eaac0f3SBarry Smith           colb = aij->colmap[col] - 1;
178b3d2dc96SSatish Balay #endif
1796eaac0f3SBarry Smith           if (colb == -1) {
1806eaac0f3SBarry Smith             m = 0;
1816eaac0f3SBarry Smith           } else {
1826eaac0f3SBarry Smith             rows = B_cj + B_ci[colb];
1836eaac0f3SBarry Smith             m    = B_ci[colb+1] - B_ci[colb];
1846eaac0f3SBarry Smith           }
1856eaac0f3SBarry Smith         }
186b8f8c88eSHong Zhang 
187a64fbb32SBarry Smith         /* loop over columns marking them in rowhit */
1886eaac0f3SBarry Smith         fm    = M; /* fm points to first entry in linked list */
189a64fbb32SBarry Smith         for (k=0; k<m; k++) {
190a64fbb32SBarry Smith           currentcol = *rows++;
191a64fbb32SBarry Smith 	  /* is it already in the list? */
192a64fbb32SBarry Smith           do {
193a64fbb32SBarry Smith             mfm  = fm;
194a64fbb32SBarry Smith             fm   = rowhit[fm];
195a64fbb32SBarry Smith           } while (fm < currentcol);
196a64fbb32SBarry Smith           /* not in list so add it */
197a64fbb32SBarry Smith           if (fm != currentcol) {
198a64fbb32SBarry Smith             nrows++;
199a64fbb32SBarry Smith             columnsforrow[currentcol] = col;
200a64fbb32SBarry Smith             /* next three lines insert new entry into linked list */
201a64fbb32SBarry Smith             rowhit[mfm]               = currentcol;
202a64fbb32SBarry Smith             rowhit[currentcol]        = fm;
203a64fbb32SBarry Smith             fm                        = currentcol;
204a64fbb32SBarry Smith             /* fm points to present position in list since we know the columns are sorted */
205a64fbb32SBarry Smith           } else {
20629bbc08cSBarry Smith             SETERRQ(PETSC_ERR_PLIB,"Invalid coloring of matrix detected");
207a64fbb32SBarry Smith           }
208a64fbb32SBarry Smith         }
209a64fbb32SBarry Smith       }
210a64fbb32SBarry Smith       c->nrows[i]         = nrows;
211b1d57f15SBarry Smith       ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr);
212b1d57f15SBarry Smith       ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr);
21352e6d16bSBarry Smith       ierr = PetscLogObjectMemory(c,(nrows+1)*sizeof(PetscInt));CHKERRQ(ierr);
214a64fbb32SBarry Smith       /* now store the linked list of rows into c->rows[i] */
215a64fbb32SBarry Smith       nrows = 0;
2166eaac0f3SBarry Smith       fm    = rowhit[M];
217a64fbb32SBarry Smith       do {
218a64fbb32SBarry Smith         c->rows[i][nrows]            = fm;
219a64fbb32SBarry Smith         c->columnsforrow[i][nrows++] = columnsforrow[fm];
220a64fbb32SBarry Smith         fm                           = rowhit[fm];
2216eaac0f3SBarry Smith       } while (fm < M);
2226eaac0f3SBarry Smith     } /* ---------------------------------------------------------------------------------------*/
223606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
2246eaac0f3SBarry Smith   }
22530b34957SBarry Smith 
22630b34957SBarry Smith   /* Optimize by adding the vscale, and scaleforrow[][] fields */
22730b34957SBarry Smith   /*
22830b34957SBarry Smith        vscale will contain the "diagonal" on processor scalings followed by the off processor
22930b34957SBarry Smith   */
2308ee2e534SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
231b8f8c88eSHong Zhang     ierr = VecCreateGhost(mat->comm,aij->A->rmap.n,PETSC_DETERMINE,aij->B->cmap.n,aij->garray,&c->vscale);CHKERRQ(ierr);
232b1d57f15SBarry Smith     ierr = PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);CHKERRQ(ierr);
23330b34957SBarry Smith     for (k=0; k<c->ncolors; k++) {
234b1d57f15SBarry Smith       ierr = PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);CHKERRQ(ierr);
23530b34957SBarry Smith       for (l=0; l<c->nrows[k]; l++) {
23630b34957SBarry Smith         col = c->columnsforrow[k][l];
23730b34957SBarry Smith         if (col >= cstart && col < cend) {
23830b34957SBarry Smith           /* column is in diagonal block of matrix */
23930b34957SBarry Smith           colb = col - cstart;
24030b34957SBarry Smith         } else {
24130b34957SBarry Smith           /* column  is in "off-processor" part */
24230b34957SBarry Smith #if defined (PETSC_USE_CTABLE)
24330b34957SBarry Smith           ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr);
24430b34957SBarry Smith           colb --;
24530b34957SBarry Smith #else
24630b34957SBarry Smith           colb = aij->colmap[col] - 1;
24730b34957SBarry Smith #endif
24830b34957SBarry Smith           colb += cend - cstart;
24930b34957SBarry Smith         }
25030b34957SBarry Smith         c->vscaleforrow[k][l] = colb;
25130b34957SBarry Smith       }
25230b34957SBarry Smith     }
253b8f8c88eSHong Zhang   } else if (ctype == IS_COLORING_GHOSTED) {
254b8f8c88eSHong Zhang     /* Get gtol mapping */
255b8f8c88eSHong Zhang     PetscInt N = mat->cmap.N, *gtol;
256b8f8c88eSHong Zhang     ierr = PetscMalloc((N+1)*sizeof(PetscInt),&gtol);CHKERRQ(ierr);
257b8f8c88eSHong Zhang     for (i=0; i<N; i++) gtol[i] = -1;
258b8f8c88eSHong Zhang     for (i=0; i<map->n; i++) gtol[ltog[i]] = i;
259b8f8c88eSHong Zhang 
260b8f8c88eSHong Zhang     c->vscale = 0; /* will be created in MatFDColoringApply() */
261b8f8c88eSHong Zhang     ierr = PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);CHKERRQ(ierr);
262b8f8c88eSHong Zhang     for (k=0; k<c->ncolors; k++) {
263b8f8c88eSHong Zhang       ierr = PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);CHKERRQ(ierr);
264b8f8c88eSHong Zhang       for (l=0; l<c->nrows[k]; l++) {
265b8f8c88eSHong Zhang         col = c->columnsforrow[k][l];      /* global column index */
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);
275*8f7157efSSatish Balay   ierr = MatRestoreColumnIJ_SeqAIJ(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);CHKERRQ(ierr);
276*8f7157efSSatish Balay   ierr = MatRestoreColumnIJ_SeqAIJ(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);CHKERRQ(ierr);
2773a40ed3dSBarry Smith   PetscFunctionReturn(0);
278a64fbb32SBarry Smith }
279a64fbb32SBarry Smith 
280b9617806SBarry Smith 
281b9617806SBarry Smith 
282b9617806SBarry Smith 
283b9617806SBarry Smith 
284b9617806SBarry Smith 
285