xref: /petsc/src/mat/impls/aij/mpi/fdmpiaij.c (revision 69bb7ac94e38f481a514b28cb10ff7ece56113f9) !
1 #define PETSCMAT_DLL
2 
3 #include "src/mat/impls/aij/mpi/mpiaij.h"
4 
5 EXTERN PetscErrorCode CreateColmap_MPIAIJ_Private(Mat);
6 EXTERN PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscInt*,PetscInt*[],PetscInt*[],PetscTruth*);
7 EXTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscInt*,PetscInt*[],PetscInt*[],PetscTruth*);
8 
9 #undef __FUNCT__
10 #define __FUNCT__ "MatFDColoringCreate_MPIAIJ"
11 PetscErrorCode MatFDColoringCreate_MPIAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c)
12 {
13   Mat_MPIAIJ     *aij = (Mat_MPIAIJ*)mat->data;
14   PetscErrorCode ierr;
15   PetscMPIInt    size,*ncolsonproc,*disp,nn;
16   PetscInt       i,*is,n,nrows,j,k,m,*rows = 0,*A_ci,*A_cj,ncols,col;
17   PetscInt       nis = iscoloring->n,nctot,*cols,*B_ci,*B_cj;
18   PetscInt       *rowhit,M = mat->rmap.n,cstart = mat->cmap.rstart,cend = mat->cmap.rend,colb;
19   PetscInt       *columnsforrow,l;
20   IS             *isa;
21   PetscTruth     done,flg;
22   ISLocalToGlobalMapping map = mat->mapping;
23   PetscInt       *ltog = map->indices,ctype=c->ctype;
24 
25   PetscFunctionBegin;
26   if (!mat->assembled) {
27     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be assembled first; MatAssemblyBegin/End();");
28   }
29 
30   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
31   c->M             = mat->rmap.N;  /* set the global rows and columns and local rows */
32   c->N             = mat->cmap.N;
33   c->m             = mat->rmap.n;
34   c->rstart        = mat->rmap.rstart;
35 
36   c->ncolors       = nis;
37   ierr             = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr);
38   ierr             = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr);
39   ierr             = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr);
40   ierr             = PetscMalloc(nis*sizeof(PetscInt*),&c->rows);CHKERRQ(ierr);
41   ierr             = PetscMalloc(nis*sizeof(PetscInt*),&c->columnsforrow);CHKERRQ(ierr);
42   ierr = PetscLogObjectMemory(c,5*nis*sizeof(PetscInt));CHKERRQ(ierr);
43 
44   /* Allow access to data structures of local part of matrix */
45   if (!aij->colmap) {
46     ierr = CreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr);
47   }
48   /*
49       Calls the _SeqAIJ() version of these routines to make sure it does not
50      get the reduced (by inodes) version of I and J
51   */
52   ierr = MatGetColumnIJ_SeqAIJ(aij->A,0,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);CHKERRQ(ierr);
53   ierr = MatGetColumnIJ_SeqAIJ(aij->B,0,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);CHKERRQ(ierr);
54 
55   ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr);
56   ierr = PetscMalloc((M+1)*sizeof(PetscInt),&columnsforrow);CHKERRQ(ierr);
57 
58   for (i=0; i<nis; i++) {
59     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
60     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
61     c->ncolumns[i] = n;
62     if (n) {
63       ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr);
64       ierr = PetscLogObjectMemory(c,n*sizeof(PetscInt));CHKERRQ(ierr);
65       ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr);
66     } else {
67       c->columns[i]  = 0;
68     }
69 
70     if (ctype == IS_COLORING_LOCAL){
71       /* Determine the total (parallel) number of columns of this color */
72       ierr = MPI_Comm_size(mat->comm,&size);CHKERRQ(ierr);
73       ierr = PetscMalloc(2*size*sizeof(PetscInt*),&ncolsonproc);CHKERRQ(ierr);
74       disp = ncolsonproc + size;
75 
76       nn   = (PetscMPIInt)n;
77       ierr = MPI_Allgather(&nn,1,MPI_INT,ncolsonproc,1,MPI_INT,mat->comm);CHKERRQ(ierr);
78       nctot = 0; for (j=0; j<size; j++) {nctot += ncolsonproc[j];}
79       if (!nctot) {
80         ierr = PetscInfo(mat,"Coloring of matrix has some unneeded colors with no corresponding rows\n");CHKERRQ(ierr);
81       }
82 
83       disp[0] = 0;
84       for (j=1; j<size; j++) {
85         disp[j] = disp[j-1] + ncolsonproc[j-1];
86       }
87 
88       /* Get complete list of columns for color on each processor */
89       ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr);
90       ierr = MPI_Allgatherv(is,n,MPIU_INT,cols,ncolsonproc,disp,MPIU_INT,mat->comm);CHKERRQ(ierr);
91       ierr = PetscFree(ncolsonproc);CHKERRQ(ierr);
92     } else if (ctype == IS_COLORING_GHOSTED){
93       /* Determine local number of columns of this color on this process, including ghost points */
94       nctot = n;
95       ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr);
96       ierr = PetscMemcpy(cols,is,n*sizeof(PetscInt));CHKERRQ(ierr);
97     } else {
98       SETERRQ(PETSC_ERR_SUP,"Not provided for this MatFDColoring type");
99     }
100 
101     /*
102        Mark all rows affect by these columns
103     */
104     /* Temporary option to allow for debugging/testing */
105     ierr = PetscOptionsHasName(PETSC_NULL,"-matfdcoloring_slow",&flg);CHKERRQ(ierr);
106     if (!flg) {/*-----------------------------------------------------------------------------*/
107       /* crude, fast version */
108       ierr = PetscMemzero(rowhit,M*sizeof(PetscInt));CHKERRQ(ierr);
109       /* loop over columns*/
110       for (j=0; j<nctot; j++) {
111         if (ctype == IS_COLORING_GHOSTED) {
112           col = ltog[cols[j]];
113         } else {
114           col  = cols[j];
115         }
116         if (col >= cstart && col < cend) {
117           /* column is in diagonal block of matrix */
118           rows = A_cj + A_ci[col-cstart];
119           m    = A_ci[col-cstart+1] - A_ci[col-cstart];
120         } else {
121 #if defined (PETSC_USE_CTABLE)
122           ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr)
123 	  colb --;
124 #else
125           colb = aij->colmap[col] - 1;
126 #endif
127           if (colb == -1) {
128             m = 0;
129           } else {
130             rows = B_cj + B_ci[colb];
131             m    = B_ci[colb+1] - B_ci[colb];
132           }
133         }
134         /* loop over columns marking them in rowhit */
135         for (k=0; k<m; k++) {
136           rowhit[*rows++] = col + 1;
137         }
138       }
139 
140       /* count the number of hits */
141       nrows = 0;
142       for (j=0; j<M; j++) {
143         if (rowhit[j]) nrows++;
144       }
145       c->nrows[i]         = nrows;
146       ierr                = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr);
147       ierr                = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr);
148       ierr = PetscLogObjectMemory(c,2*(nrows+1)*sizeof(PetscInt));CHKERRQ(ierr);
149       nrows = 0;
150       for (j=0; j<M; j++) {
151         if (rowhit[j]) {
152           c->rows[i][nrows]           = j;
153           c->columnsforrow[i][nrows] = rowhit[j] - 1;
154           nrows++;
155         }
156       }
157     } else {/*-------------------------------------------------------------------------------*/
158       /* slow version, using rowhit as a linked list */
159       PetscInt currentcol,fm,mfm;
160       rowhit[M] = M;
161       nrows     = 0;
162       /* loop over columns*/
163       for (j=0; j<nctot; j++) {
164         if (ctype == IS_COLORING_GHOSTED) {
165           col = ltog[cols[j]];
166         } else {
167           col  = cols[j];
168         }
169         if (col >= cstart && col < cend) {
170           /* column is in diagonal block of matrix */
171           rows = A_cj + A_ci[col-cstart];
172           m    = A_ci[col-cstart+1] - A_ci[col-cstart];
173         } else {
174 #if defined (PETSC_USE_CTABLE)
175           ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr);
176           colb --;
177 #else
178           colb = aij->colmap[col] - 1;
179 #endif
180           if (colb == -1) {
181             m = 0;
182           } else {
183             rows = B_cj + B_ci[colb];
184             m    = B_ci[colb+1] - B_ci[colb];
185           }
186         }
187 
188         /* loop over columns marking them in rowhit */
189         fm    = M; /* fm points to first entry in linked list */
190         for (k=0; k<m; k++) {
191           currentcol = *rows++;
192 	  /* is it already in the list? */
193           do {
194             mfm  = fm;
195             fm   = rowhit[fm];
196           } while (fm < currentcol);
197           /* not in list so add it */
198           if (fm != currentcol) {
199             nrows++;
200             columnsforrow[currentcol] = col;
201             /* next three lines insert new entry into linked list */
202             rowhit[mfm]               = currentcol;
203             rowhit[currentcol]        = fm;
204             fm                        = currentcol;
205             /* fm points to present position in list since we know the columns are sorted */
206           } else {
207             SETERRQ(PETSC_ERR_PLIB,"Invalid coloring of matrix detected");
208           }
209         }
210       }
211       c->nrows[i]         = nrows;
212       ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr);
213       ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr);
214       ierr = PetscLogObjectMemory(c,(nrows+1)*sizeof(PetscInt));CHKERRQ(ierr);
215       /* now store the linked list of rows into c->rows[i] */
216       nrows = 0;
217       fm    = rowhit[M];
218       do {
219         c->rows[i][nrows]            = fm;
220         c->columnsforrow[i][nrows++] = columnsforrow[fm];
221         fm                           = rowhit[fm];
222       } while (fm < M);
223     } /* ---------------------------------------------------------------------------------------*/
224     ierr = PetscFree(cols);CHKERRQ(ierr);
225   }
226 
227   /* Optimize by adding the vscale, and scaleforrow[][] fields */
228   /*
229        vscale will contain the "diagonal" on processor scalings followed by the off processor
230   */
231   if (ctype == IS_COLORING_LOCAL) {
232     ierr = VecCreateGhost(mat->comm,aij->A->rmap.n,PETSC_DETERMINE,aij->B->cmap.n,aij->garray,&c->vscale);CHKERRQ(ierr);
233     ierr = PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);CHKERRQ(ierr);
234     for (k=0; k<c->ncolors; k++) {
235       ierr = PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);CHKERRQ(ierr);
236       for (l=0; l<c->nrows[k]; l++) {
237         col = c->columnsforrow[k][l];
238         if (col >= cstart && col < cend) {
239           /* column is in diagonal block of matrix */
240           colb = col - cstart;
241         } else {
242           /* column  is in "off-processor" part */
243 #if defined (PETSC_USE_CTABLE)
244           ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr);
245           colb --;
246 #else
247           colb = aij->colmap[col] - 1;
248 #endif
249           colb += cend - cstart;
250         }
251         c->vscaleforrow[k][l] = colb;
252       }
253     }
254   } else if (ctype == IS_COLORING_GHOSTED) {
255     /* Get gtol mapping */
256     PetscInt N = mat->cmap.N, *gtol;
257     ierr = PetscMalloc((N+1)*sizeof(PetscInt),&gtol);CHKERRQ(ierr);
258     for (i=0; i<N; i++) gtol[i] = -1;
259     for (i=0; i<map->n; i++) gtol[ltog[i]] = i;
260 
261     c->vscale = 0; /* will be created in MatFDColoringApply() */
262     ierr = PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);CHKERRQ(ierr);
263     for (k=0; k<c->ncolors; k++) {
264       ierr = PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);CHKERRQ(ierr);
265       for (l=0; l<c->nrows[k]; l++) {
266         col = c->columnsforrow[k][l];      /* global column index */
267         c->vscaleforrow[k][l] = gtol[col]; /* local column index */
268       }
269     }
270     ierr = PetscFree(gtol);CHKERRQ(ierr);
271   }
272   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
273 
274   ierr = PetscFree(rowhit);CHKERRQ(ierr);
275   ierr = PetscFree(columnsforrow);CHKERRQ(ierr);
276   ierr = MatRestoreColumnIJ_SeqAIJ(aij->A,0,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);CHKERRQ(ierr);
277   ierr = MatRestoreColumnIJ_SeqAIJ(aij->B,0,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);CHKERRQ(ierr);
278   PetscFunctionReturn(0);
279 }
280 
281 
282 
283 
284 
285 
286