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