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