Actual source code: fdmpiaij.c

  1: #define PETSCMAT_DLL

 3:  #include src/mat/impls/aij/mpi/mpiaij.h

  5: EXTERN PetscErrorCode CreateColmap_MPIAIJ_Private(Mat);

  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,*is,n,nrows,j,k,m,*rows = 0,*A_ci,*A_cj,ncols,col;
 15:   PetscInt              nis = iscoloring->n,nctot,*cols,*B_ci,*B_cj;
 16:   PetscInt              *rowhit,M = mat->rmap.n,cstart = mat->cmap.rstart,cend = mat->cmap.rend,colb;
 17:   PetscInt              *columnsforrow,l;
 18:   IS                    *isa;
 19:   PetscTruth             done,flg;
 20:   ISLocalToGlobalMapping map = mat->mapping;
 21:   PetscInt               *ltog = (map ? map->indices : 0) ,ctype=c->ctype;

 24:   if (!mat->assembled) {
 25:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be assembled first; MatAssemblyBegin/End();");
 26:   }
 27:   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:   ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);
 30:   c->M             = mat->rmap.N;  /* set the global rows and columns and local rows */
 31:   c->N             = mat->cmap.N;
 32:   c->m             = mat->rmap.n;
 33:   c->rstart        = mat->rmap.rstart;

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

 43:   /* Allow access to data structures of local part of matrix */
 44:   if (!aij->colmap) {
 45:     CreateColmap_MPIAIJ_Private(mat);
 46:   }
 47:   /*
 48:       Calls the _SeqAIJ() version of these routines to make sure it does not 
 49:      get the reduced (by inodes) version of I and J
 50:   */
 51:   MatGetColumnIJ_SeqAIJ(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);
 52:   MatGetColumnIJ_SeqAIJ(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);
 53: 
 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);
 60:     c->ncolumns[i] = n;
 61:     if (n) {
 62:       PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);
 63:       PetscLogObjectMemory(c,n*sizeof(PetscInt));
 64:       PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));
 65:     } else {
 66:       c->columns[i]  = 0;
 67:     }

 69:     if (ctype == IS_COLORING_GLOBAL){
 70:       /* Determine the total (parallel) number of columns of this color */
 71:       MPI_Comm_size(mat->comm,&size);
 72:       PetscMalloc(2*size*sizeof(PetscInt*),&ncolsonproc);
 73:       disp = ncolsonproc + size;

 75:       nn   = (PetscMPIInt)n;
 76:       MPI_Allgather(&nn,1,MPI_INT,ncolsonproc,1,MPI_INT,mat->comm);
 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(is,n,MPIU_INT,cols,ncolsonproc,disp,MPIU_INT,mat->comm);
 90:       PetscFree(ncolsonproc);
 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 {
 97:       SETERRQ(PETSC_ERR_SUP,"Not provided for this MatFDColoring type");
 98:     }

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

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

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

226:   /* Optimize by adding the vscale, and scaleforrow[][] fields */
227:   /*
228:        vscale will contain the "diagonal" on processor scalings followed by the off processor
229:   */
230:   if (ctype == IS_COLORING_GLOBAL) {
231:     VecCreateGhost(mat->comm,aij->A->rmap.n,PETSC_DETERMINE,aij->B->cmap.n,aij->garray,&c->vscale);
232:     PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);
233:     for (k=0; k<c->ncolors; k++) {
234:       PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);
235:       for (l=0; l<c->nrows[k]; l++) {
236:         col = c->columnsforrow[k][l];
237:         if (col >= cstart && col < cend) {
238:           /* column is in diagonal block of matrix */
239:           colb = col - cstart;
240:         } else {
241:           /* column  is in "off-processor" part */
242: #if defined (PETSC_USE_CTABLE)
243:           PetscTableFind(aij->colmap,col+1,&colb);
244:           colb --;
245: #else
246:           colb = aij->colmap[col] - 1;
247: #endif
248:           colb += cend - cstart;
249:         }
250:         c->vscaleforrow[k][l] = colb;
251:       }
252:     }
253:   } else if (ctype == IS_COLORING_GHOSTED) {
254:     /* Get gtol mapping */
255:     PetscInt N = mat->cmap.N, *gtol;
256:     PetscMalloc((N+1)*sizeof(PetscInt),&gtol);
257:     for (i=0; i<N; i++) gtol[i] = -1;
258:     for (i=0; i<map->n; i++) gtol[ltog[i]] = i;
259: 
260:     c->vscale = 0; /* will be created in MatFDColoringApply() */
261:     PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);
262:     for (k=0; k<c->ncolors; k++) {
263:       PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);
264:       for (l=0; l<c->nrows[k]; l++) {
265:         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_SeqAIJ(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);
276:   MatRestoreColumnIJ_SeqAIJ(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);
277:   return(0);
278: }