Actual source code: fdmpiaij.c

petsc-3.4.2 2013-07-02
  2: #include <../src/mat/impls/aij/mpi/mpiaij.h>

  4: extern PetscErrorCode MatCreateColmap_MPIAIJ_Private(Mat);

  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,ncols,col;
 14:   const PetscInt         *is,*A_ci,*A_cj,*B_ci,*B_cj,*rows = 0,*ltog;
 15:   PetscInt               nis = iscoloring->n,nctot,*cols;
 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               ctype=c->ctype;

 24:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Matrix must be assembled first; MatAssemblyBegin/End();");
 25:   if (ctype == IS_COLORING_GHOSTED && !map) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_INCOMP,"When using ghosted differencing matrix must have local to global mapping provided with MatSetLocalToGlobalMapping");

 27:   if (map) {ISLocalToGlobalMappingGetIndices(map,&ltog);}
 28:   else     ltog = NULL;
 29:   ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);

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

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

 47:   /* Allow access to data structures of local part of matrix */
 48:   if (!aij->colmap) {
 49:     MatCreateColmap_MPIAIJ_Private(mat);
 50:   }
 51:   MatGetColumnIJ(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);
 52:   MatGetColumnIJ(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);

 54:   PetscMalloc((M+1)*sizeof(PetscInt),&rowhit);
 55:   PetscMalloc((M+1)*sizeof(PetscInt),&columnsforrow);

 57:   for (i=0; i<nis; i++) {
 58:     ISGetLocalSize(isa[i],&n);
 59:     ISGetIndices(isa[i],&is);

 61:     c->ncolumns[i] = n;
 62:     if (n) {
 63:       PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);
 64:       PetscLogObjectMemory(c,n*sizeof(PetscInt));
 65:       PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));
 66:     } else {
 67:       c->columns[i] = 0;
 68:     }

 70:     if (ctype == IS_COLORING_GLOBAL) {
 71:       /* Determine the total (parallel) number of columns of this color */
 72:       MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);
 73:       PetscMalloc2(size,PetscMPIInt,&ncolsonproc,size,PetscMPIInt,&disp);

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

 82:       disp[0] = 0;
 83:       for (j=1; j<size; j++) {
 84:         disp[j] = disp[j-1] + ncolsonproc[j-1];
 85:       }

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

 98:     /*
 99:        Mark all rows affect by these columns
100:     */
101:     /* Temporary option to allow for debugging/testing */
102:     flg  = PETSC_FALSE;
103:     PetscOptionsGetBool(NULL,"-matfdcoloring_slow",&flg,NULL);
104:     if (!flg) { /*-----------------------------------------------------------------------------*/
105:       /* crude, fast version */
106:       PetscMemzero(rowhit,M*sizeof(PetscInt));
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:           PetscTableFind(aij->colmap,col+1,&colb);
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:       }

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:       PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);
145:       PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);
146:       PetscLogObjectMemory(c,2*(nrows+1)*sizeof(PetscInt));
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:           PetscTableFind(aij->colmap,col+1,&colb);
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:         }

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 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid coloring of matrix detected");
205:         }
206:       }
207:       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(PetscObjectComm((PetscObject)mat),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,nlocal,*gtol;
254:     PetscMalloc((N+1)*sizeof(PetscInt),&gtol);
255:     for (i=0; i<N; i++) gtol[i] = -1;
256:     ISLocalToGlobalMappingGetSize(map,&nlocal);
257:     for (i=0; i<nlocal; i++) gtol[ltog[i]] = i;

259:     c->vscale = 0; /* will be created in MatFDColoringApply() */
260:     PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);
261:     for (k=0; k<c->ncolors; k++) {
262:       PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);
263:       for (l=0; l<c->nrows[k]; l++) {
264:         col = c->columnsforrow[k][l];      /* global column index */

266:         c->vscaleforrow[k][l] = gtol[col]; /* local column index */
267:       }
268:     }
269:     PetscFree(gtol);
270:   }
271:   ISColoringRestoreIS(iscoloring,&isa);

273:   PetscFree(rowhit);
274:   PetscFree(columnsforrow);
275:   MatRestoreColumnIJ(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);
276:   MatRestoreColumnIJ(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);
277:   if (map) {ISLocalToGlobalMappingRestoreIndices(map,&ltog);}
278:   return(0);
279: }