Actual source code: fdaij.c

  1: /*$Id: fdaij.c,v 1.40 2001/06/21 21:16:21 bsmith Exp $*/

 3:  #include src/mat/impls/aij/seq/aij.h
 4:  #include src/vec/vecimpl.h

  6: EXTERN int MatGetColumnIJ_SeqAIJ(Mat,int,PetscTruth,int*,int**,int**,PetscTruth*);
  7: EXTERN int MatRestoreColumnIJ_SeqAIJ(Mat,int,PetscTruth,int*,int**,int**,PetscTruth*);

  9: int MatFDColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c)
 10: {
 11:   int        i,*is,n,nrows,N = mat->N,j,k,m,*rows,ierr,*ci,*cj,ncols,col;
 12:   int        nis = iscoloring->n,*rowhit,*columnsforrow,l;
 13:   IS         *isa;
 14:   PetscTruth done,flg;

 17:   if (!mat->assembled) {
 18:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be assembled by calls to MatAssemblyBegin/End();");
 19:   }

 21:   ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);
 22:   c->M       = mat->M;  /* set total rows, columns and local rows */
 23:   c->N       = mat->N;
 24:   c->m       = mat->M;
 25:   c->rstart  = 0;

 27:   c->ncolors = nis;
 28:   ierr       = PetscMalloc(nis*sizeof(int),&c->ncolumns);
 29:   ierr       = PetscMalloc(nis*sizeof(int*),&c->columns);
 30:   ierr       = PetscMalloc(nis*sizeof(int),&c->nrows);
 31:   ierr       = PetscMalloc(nis*sizeof(int*),&c->rows);
 32:   ierr       = PetscMalloc(nis*sizeof(int*),&c->columnsforrow);

 34:   /*
 35:       Calls the _SeqAIJ() version of these routines to make sure it does not 
 36:      get the reduced (by inodes) version of I and J
 37:   */
 38:   MatGetColumnIJ_SeqAIJ(mat,0,PETSC_FALSE,&ncols,&ci,&cj,&done);

 40:   /*
 41:      Temporary option to allow for debugging/testing
 42:   */
 43:   PetscOptionsHasName(PETSC_NULL,"-matfdcoloring_slow",&flg);

 45:   PetscMalloc((N+1)*sizeof(int),&rowhit);
 46:   PetscMalloc((N+1)*sizeof(int),&columnsforrow);

 48:   for (i=0; i<nis; i++) {
 49:     ISGetLocalSize(isa[i],&n);
 50:     ISGetIndices(isa[i],&is);
 51:     c->ncolumns[i] = n;
 52:     if (n) {
 53:       PetscMalloc(n*sizeof(int),&c->columns[i]);
 54:       PetscMemcpy(c->columns[i],is,n*sizeof(int));
 55:     } else {
 56:       c->columns[i]  = 0;
 57:     }

 59:     if (!flg) { /* ------------------------------------------------------------------------------*/
 60:       /* fast, crude version requires O(N*N) work */
 61:       PetscMemzero(rowhit,N*sizeof(int));
 62:       /* loop over columns*/
 63:       for (j=0; j<n; j++) {
 64:         col  = is[j];
 65:         rows = cj + ci[col];
 66:         m    = ci[col+1] - ci[col];
 67:         /* loop over columns marking them in rowhit */
 68:         for (k=0; k<m; k++) {
 69:           rowhit[*rows++] = col + 1;
 70:         }
 71:       }
 72:       /* count the number of hits */
 73:       nrows = 0;
 74:       for (j=0; j<N; j++) {
 75:         if (rowhit[j]) nrows++;
 76:       }
 77:       c->nrows[i] = nrows;
 78:       ierr        = PetscMalloc((nrows+1)*sizeof(int),&c->rows[i]);
 79:       ierr        = PetscMalloc((nrows+1)*sizeof(int),&c->columnsforrow[i]);
 80:       nrows       = 0;
 81:       for (j=0; j<N; j++) {
 82:         if (rowhit[j]) {
 83:           c->rows[i][nrows]          = j;
 84:           c->columnsforrow[i][nrows] = rowhit[j] - 1;
 85:           nrows++;
 86:         }
 87:       }
 88:     } else {  /*-------------------------------------------------------------------------------*/
 89:       /* slow version, using rowhit as a linked list */
 90:       int currentcol,fm,mfm;
 91:       rowhit[N] = N;
 92:       nrows     = 0;
 93:       /* loop over columns */
 94:       for (j=0; j<n; j++) {
 95:         col   = is[j];
 96:         rows  = cj + ci[col];
 97:         m     = ci[col+1] - ci[col];
 98:         /* loop over columns marking them in rowhit */
 99:         fm    = N; /* fm points to first entry in linked list */
100:         for (k=0; k<m; k++) {
101:           currentcol = *rows++;
102:           /* is it already in the list? */
103:           do {
104:             mfm  = fm;
105:             fm   = rowhit[fm];
106:           } while (fm < currentcol);
107:           /* not in list so add it */
108:           if (fm != currentcol) {
109:             nrows++;
110:             columnsforrow[currentcol] = col;
111:             /* next three lines insert new entry into linked list */
112:             rowhit[mfm]               = currentcol;
113:             rowhit[currentcol]        = fm;
114:             fm                        = currentcol;
115:             /* fm points to present position in list since we know the columns are sorted */
116:           } else {
117:             SETERRQ(PETSC_ERR_PLIB,"Detected invalid coloring");
118:           }

120:         }
121:       }
122:       c->nrows[i] = nrows;
123:       ierr        = PetscMalloc((nrows+1)*sizeof(int),&c->rows[i]);
124:       ierr        = PetscMalloc((nrows+1)*sizeof(int),&c->columnsforrow[i]);
125:       /* now store the linked list of rows into c->rows[i] */
126:       nrows       = 0;
127:       fm          = rowhit[N];
128:       do {
129:         c->rows[i][nrows]            = fm;
130:         c->columnsforrow[i][nrows++] = columnsforrow[fm];
131:         fm                           = rowhit[fm];
132:       } while (fm < N);
133:     } /* ---------------------------------------------------------------------------------------*/
134:     ISRestoreIndices(isa[i],&is);
135:   }
136:   MatRestoreColumnIJ_SeqAIJ(mat,0,PETSC_FALSE,&ncols,&ci,&cj,&done);

138:   PetscFree(rowhit);
139:   PetscFree(columnsforrow);

141:   /* Optimize by adding the vscale, and scaleforrow[][] fields */
142:   /*
143:        see the version for MPIAIJ
144:   */
145:   VecCreateGhost(mat->comm,mat->m,PETSC_DETERMINE,0,PETSC_NULL,&c->vscale);CHKERRQ(ierr)
146:   PetscMalloc(c->ncolors*sizeof(int*),&c->vscaleforrow);
147:   for (k=0; k<c->ncolors; k++) {
148:     PetscMalloc((c->nrows[k]+1)*sizeof(int),&c->vscaleforrow[k]);
149:     for (l=0; l<c->nrows[k]; l++) {
150:       col = c->columnsforrow[k][l];
151:       c->vscaleforrow[k][l] = col;
152:     }
153:   }
154:   ISColoringRestoreIS(iscoloring,&isa);
155:   return(0);
156: }

158: /*
159:      Makes a longer coloring[] array and calls the usual code with that
160: */
161: int MatColoringPatch_SeqAIJ_Inode(Mat mat,int nin,int ncolors,int *coloring,ISColoring *iscoloring)
162: {
163:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)mat->data;
164:   int        n = mat->n,ierr,m = a->inode.node_count,j,*ns = a->inode.size,row;
165:   int        *colorused,i,*newcolor;

168:   PetscMalloc((n+1)*sizeof(int),&newcolor);
169:   /* loop over inodes, marking a color for each column*/
170:   row = 0;
171:   for (i=0; i<m; i++){
172:     for (j=0; j<ns[i]; j++) {
173:       newcolor[row++] = coloring[i] + j*ncolors;
174:     }
175:   }

177:   /* eliminate unneeded colors */
178:   PetscMalloc(5*ncolors*sizeof(int),&colorused);
179:   PetscMemzero(colorused,5*ncolors*sizeof(int));
180:   for (i=0; i<n; i++) {
181:     colorused[newcolor[i]] = 1;
182:   }

184:   for (i=1; i<5*ncolors; i++) {
185:     colorused[i] += colorused[i-1];
186:   }
187:   ncolors = colorused[5*ncolors-1];
188:   for (i=0; i<n; i++) {
189:     newcolor[i] = colorused[newcolor[i]];
190:   }
191:   PetscFree(colorused);
192:   PetscFree(coloring);
193:   ISColoringCreate(mat->comm,n,newcolor,iscoloring);

195:   return(0);
196: }