Actual source code: fdmpiaij.c

  2: #include <../src/mat/impls/aij/mpi/mpiaij.h>


  8: PetscErrorCode MatFDColoringCreate_MPIAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c)
  9: {
 10:   Mat_MPIAIJ            *aij = (Mat_MPIAIJ*)mat->data;
 11:   PetscErrorCode        ierr;
 12:   PetscMPIInt           size,*ncolsonproc,*disp,nn;
 13:   PetscInt              i,n,nrows,j,k,m,*rows = 0,*A_ci,*A_cj,ncols,col;
 14:   const PetscInt        *is;
 15:   PetscInt              nis = iscoloring->n,nctot,*cols,*B_ci,*B_cj;
 16:   PetscInt              *rowhit,M,cstart,cend,colb;
 17:   PetscInt              *columnsforrow,l;
 18:   IS                    *isa;
 19:   PetscBool              done,flg;
 20:   ISLocalToGlobalMapping map = mat->cmap->mapping;
 21:   PetscInt               *ltog = (map ? map->indices : (PetscInt*) PETSC_NULL) ,ctype=c->ctype;

 24:   if (!mat->assembled) SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be assembled first; MatAssemblyBegin/End();");
 25:   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:   ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);

 29:   M                = mat->rmap->n;
 30:   cstart           = mat->cmap->rstart;
 31:   cend             = mat->cmap->rend;
 32:   c->M             = mat->rmap->N;  /* set the global rows and columns and local rows */
 33:   c->N             = mat->cmap->N;
 34:   c->m             = mat->rmap->n;
 35:   c->rstart        = mat->rmap->rstart;

 37:   c->ncolors       = nis;
 38:   PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);
 39:   PetscMalloc(nis*sizeof(PetscInt*),&c->columns);
 40:   PetscMalloc(nis*sizeof(PetscInt),&c->nrows);
 41:   PetscMalloc(nis*sizeof(PetscInt*),&c->rows);
 42:   PetscMalloc(nis*sizeof(PetscInt*),&c->columnsforrow);
 43:   PetscLogObjectMemory(c,5*nis*sizeof(PetscInt));

 45:   /* Allow access to data structures of local part of matrix */
 46:   if (!aij->colmap) {
 47:     CreateColmap_MPIAIJ_Private(mat);
 48:   }
 49:   MatGetColumnIJ(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);
 50:   MatGetColumnIJ(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);
 51: 
 52:   PetscMalloc((M+1)*sizeof(PetscInt),&rowhit);
 53:   PetscMalloc((M+1)*sizeof(PetscInt),&columnsforrow);

 55:   for (i=0; i<nis; i++) {
 56:     ISGetLocalSize(isa[i],&n);
 57:     ISGetIndices(isa[i],&is);
 58:     c->ncolumns[i] = n;
 59:     if (n) {
 60:       PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);
 61:       PetscLogObjectMemory(c,n*sizeof(PetscInt));
 62:       PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));
 63:     } else {
 64:       c->columns[i]  = 0;
 65:     }

 67:     if (ctype == IS_COLORING_GLOBAL){
 68:       /* Determine the total (parallel) number of columns of this color */
 69:       MPI_Comm_size(((PetscObject)mat)->comm,&size);
 70:       PetscMalloc2(size,PetscMPIInt,&ncolsonproc,size,PetscMPIInt,&disp);

 72:       nn   = PetscMPIIntCast(n);
 73:       MPI_Allgather(&nn,1,MPI_INT,ncolsonproc,1,MPI_INT,((PetscObject)mat)->comm);
 74:       nctot = 0; for (j=0; j<size; j++) {nctot += ncolsonproc[j];}
 75:       if (!nctot) {
 76:         PetscInfo(mat,"Coloring of matrix has some unneeded colors with no corresponding rows\n");
 77:       }

 79:       disp[0] = 0;
 80:       for (j=1; j<size; j++) {
 81:         disp[j] = disp[j-1] + ncolsonproc[j-1];
 82:       }

 84:       /* Get complete list of columns for color on each processor */
 85:       PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);
 86:       MPI_Allgatherv((void*)is,n,MPIU_INT,cols,ncolsonproc,disp,MPIU_INT,((PetscObject)mat)->comm);
 87:       PetscFree2(ncolsonproc,disp);
 88:     } else if (ctype == IS_COLORING_GHOSTED){
 89:       /* Determine local number of columns of this color on this process, including ghost points */
 90:       nctot = n;
 91:       PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);
 92:       PetscMemcpy(cols,is,n*sizeof(PetscInt));
 93:     } else {
 94:       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not provided for this MatFDColoring type");
 95:     }

 97:     /*
 98:        Mark all rows affect by these columns
 99:     */
100:     /* Temporary option to allow for debugging/testing */
101:     flg  = PETSC_FALSE;
102:     PetscOptionsGetBool(PETSC_NULL,"-matfdcoloring_slow",&flg,PETSC_NULL);
103:     if (!flg) {/*-----------------------------------------------------------------------------*/
104:       /* crude, fast version */
105:       PetscMemzero(rowhit,M*sizeof(PetscInt));
106:       /* loop over columns*/
107:       for (j=0; j<nctot; j++) {
108:         if (ctype == IS_COLORING_GHOSTED) {
109:           col = ltog[cols[j]];
110:         } else {
111:           col  = cols[j];
112:         }
113:         if (col >= cstart && col < cend) {
114:           /* column is in diagonal block of matrix */
115:           rows = A_cj + A_ci[col-cstart];
116:           m    = A_ci[col-cstart+1] - A_ci[col-cstart];
117:         } else {
118: #if defined (PETSC_USE_CTABLE)
119:           PetscTableFind(aij->colmap,col+1,&colb);
120:           colb --;
121: #else
122:           colb = aij->colmap[col] - 1;
123: #endif
124:           if (colb == -1) {
125:             m = 0;
126:           } else {
127:             rows = B_cj + B_ci[colb];
128:             m    = B_ci[colb+1] - B_ci[colb];
129:           }
130:         }
131:         /* loop over columns marking them in rowhit */
132:         for (k=0; k<m; k++) {
133:           rowhit[*rows++] = col + 1;
134:         }
135:       }

137:       /* count the number of hits */
138:       nrows = 0;
139:       for (j=0; j<M; j++) {
140:         if (rowhit[j]) nrows++;
141:       }
142:       c->nrows[i]         = nrows;
143:       PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);
144:       PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);
145:       PetscLogObjectMemory(c,2*(nrows+1)*sizeof(PetscInt));
146:       nrows = 0;
147:       for (j=0; j<M; j++) {
148:         if (rowhit[j]) {
149:           c->rows[i][nrows]           = j;
150:           c->columnsforrow[i][nrows] = rowhit[j] - 1;
151:           nrows++;
152:         }
153:       }
154:     } else {/*-------------------------------------------------------------------------------*/
155:       /* slow version, using rowhit as a linked list */
156:       PetscInt currentcol,fm,mfm;
157:       rowhit[M] = M;
158:       nrows     = 0;
159:       /* loop over columns*/
160:       for (j=0; j<nctot; j++) {
161:         if (ctype == IS_COLORING_GHOSTED) {
162:           col = ltog[cols[j]];
163:         } else {
164:           col  = cols[j];
165:         }
166:         if (col >= cstart && col < cend) {
167:           /* column is in diagonal block of matrix */
168:           rows = A_cj + A_ci[col-cstart];
169:           m    = A_ci[col-cstart+1] - A_ci[col-cstart];
170:         } else {
171: #if defined (PETSC_USE_CTABLE)
172:           PetscTableFind(aij->colmap,col+1,&colb);
173:           colb --;
174: #else
175:           colb = aij->colmap[col] - 1;
176: #endif
177:           if (colb == -1) {
178:             m = 0;
179:           } else {
180:             rows = B_cj + B_ci[colb];
181:             m    = B_ci[colb+1] - B_ci[colb];
182:           }
183:         }

185:         /* loop over columns marking them in rowhit */
186:         fm    = M; /* fm points to first entry in linked list */
187:         for (k=0; k<m; k++) {
188:           currentcol = *rows++;
189:           /* is it already in the list? */
190:           do {
191:             mfm  = fm;
192:             fm   = rowhit[fm];
193:           } while (fm < currentcol);
194:           /* not in list so add it */
195:           if (fm != currentcol) {
196:             nrows++;
197:             columnsforrow[currentcol] = col;
198:             /* next three lines insert new entry into linked list */
199:             rowhit[mfm]               = currentcol;
200:             rowhit[currentcol]        = fm;
201:             fm                        = currentcol;
202:             /* fm points to present position in list since we know the columns are sorted */
203:           } else {
204:             SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid coloring of matrix detected");
205:           }
206:         }
207:       }
208:       c->nrows[i]         = nrows;
209:       PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);
210:       PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);
211:       PetscLogObjectMemory(c,(nrows+1)*sizeof(PetscInt));
212:       /* now store the linked list of rows into c->rows[i] */
213:       nrows = 0;
214:       fm    = rowhit[M];
215:       do {
216:         c->rows[i][nrows]            = fm;
217:         c->columnsforrow[i][nrows++] = columnsforrow[fm];
218:         fm                           = rowhit[fm];
219:       } while (fm < M);
220:     } /* ---------------------------------------------------------------------------------------*/
221:     PetscFree(cols);
222:   }

224:   /* Optimize by adding the vscale, and scaleforrow[][] fields */
225:   /*
226:        vscale will contain the "diagonal" on processor scalings followed by the off processor
227:   */
228:   if (ctype == IS_COLORING_GLOBAL) {
229:     VecCreateGhost(((PetscObject)mat)->comm,aij->A->rmap->n,PETSC_DETERMINE,aij->B->cmap->n,aij->garray,&c->vscale);
230:     PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);
231:     for (k=0; k<c->ncolors; k++) {
232:       PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);
233:       for (l=0; l<c->nrows[k]; l++) {
234:         col = c->columnsforrow[k][l];
235:         if (col >= cstart && col < cend) {
236:           /* column is in diagonal block of matrix */
237:           colb = col - cstart;
238:         } else {
239:           /* column  is in "off-processor" part */
240: #if defined (PETSC_USE_CTABLE)
241:           PetscTableFind(aij->colmap,col+1,&colb);
242:           colb --;
243: #else
244:           colb = aij->colmap[col] - 1;
245: #endif
246:           colb += cend - cstart;
247:         }
248:         c->vscaleforrow[k][l] = colb;
249:       }
250:     }
251:   } else if (ctype == IS_COLORING_GHOSTED) {
252:     /* Get gtol mapping */
253:     PetscInt N = mat->cmap->N, *gtol;
254:     PetscMalloc((N+1)*sizeof(PetscInt),&gtol);
255:     for (i=0; i<N; i++) gtol[i] = -1;
256:     for (i=0; i<map->n; i++) gtol[ltog[i]] = i;
257: 
258:     c->vscale = 0; /* will be created in MatFDColoringApply() */
259:     PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);
260:     for (k=0; k<c->ncolors; k++) {
261:       PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);
262:       for (l=0; l<c->nrows[k]; l++) {
263:         col = c->columnsforrow[k][l];      /* global column index */
264:         c->vscaleforrow[k][l] = gtol[col]; /* local column index */
265:       }
266:     }
267:     PetscFree(gtol);
268:   }
269:   ISColoringRestoreIS(iscoloring,&isa);

271:   PetscFree(rowhit);
272:   PetscFree(columnsforrow);
273:   MatRestoreColumnIJ(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);
274:   MatRestoreColumnIJ(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);
275:   return(0);
276: }