Actual source code: fdaij.c

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

  6: /*
  7:     This routine is shared by AIJ and BAIJ matrices, since it operators only on the nonzero structure of the elements or blocks.
  8:     This is why it has the ugly code with the MatGetBlockSize()
  9: */
 10: PetscErrorCode MatFDColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c)
 11: {
 13:   PetscInt       i,n,nrows,N,j,k,m,ncols,col;
 14:   const PetscInt *is,*rows,*ci,*cj;
 15:   PetscInt       nis = iscoloring->n,*rowhit,*columnsforrow,l,bs = 1;
 16:   IS             *isa;
 17:   PetscBool      done,flg = PETSC_FALSE;
 18:   PetscBool      flg1,flg2;

 21:   if (!mat->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be assembled by calls to MatAssemblyBegin/End();");

 23:   ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);
 24:   /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */
 25:   PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);
 26:   PetscObjectTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);
 27:   if (flg1 || flg2) {
 28:     MatGetBlockSize(mat,&bs);
 29:   }

 31:   N         = mat->cmap->N/bs;
 32:   c->M      = mat->rmap->N/bs;   /* set total rows, columns and local rows */
 33:   c->N      = mat->cmap->N/bs;
 34:   c->m      = mat->rmap->N/bs;
 35:   c->rstart = 0;

 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);

 44:   MatGetColumnIJ(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&done);
 45:   if (!done) SETERRQ1(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"MatGetColumnIJ() not supported for matrix type %s",((PetscObject)mat)->type_name);

 47:   /*
 48:      Temporary option to allow for debugging/testing
 49:   */
 50:   PetscOptionsGetBool(NULL,"-matfdcoloring_slow",&flg,NULL);

 52:   PetscMalloc((N+1)*sizeof(PetscInt),&rowhit);
 53:   PetscMalloc((N+1)*sizeof(PetscInt),&columnsforrow);

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

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

 67:     if (!flg) { /* ------------------------------------------------------------------------------*/
 68:       /* fast, crude version requires O(N*N) work */
 69:       PetscMemzero(rowhit,N*sizeof(PetscInt));
 70:       /* loop over columns*/
 71:       for (j=0; j<n; j++) {
 72:         col  = is[j];
 73:         rows = cj + ci[col];
 74:         m    = ci[col+1] - ci[col];
 75:         /* loop over columns marking them in rowhit */
 76:         for (k=0; k<m; k++) {
 77:           rowhit[*rows++] = col + 1;
 78:         }
 79:       }
 80:       /* count the number of hits */
 81:       nrows = 0;
 82:       for (j=0; j<N; j++) {
 83:         if (rowhit[j]) nrows++;
 84:       }
 85:       c->nrows[i] = nrows;
 86:       PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);
 87:       PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);
 88:       nrows       = 0;
 89:       for (j=0; j<N; j++) {
 90:         if (rowhit[j]) {
 91:           c->rows[i][nrows]          = j;
 92:           c->columnsforrow[i][nrows] = rowhit[j] - 1;
 93:           nrows++;
 94:         }
 95:       }
 96:     } else {  /*-------------------------------------------------------------------------------*/
 97:       /* slow version, using rowhit as a linked list */
 98:       PetscInt currentcol,fm,mfm;
 99:       rowhit[N] = N;
100:       nrows     = 0;
101:       /* loop over columns */
102:       for (j=0; j<n; j++) {
103:         col  = is[j];
104:         rows = cj + ci[col];
105:         m    = ci[col+1] - ci[col];
106:         /* loop over columns marking them in rowhit */
107:         fm = N;    /* fm points to first entry in linked list */
108:         for (k=0; k<m; k++) {
109:           currentcol = *rows++;
110:           /* is it already in the list? */
111:           do {
112:             mfm = fm;
113:             fm  = rowhit[fm];
114:           } while (fm < currentcol);
115:           /* not in list so add it */
116:           if (fm != currentcol) {
117:             nrows++;
118:             columnsforrow[currentcol] = col;
119:             /* next three lines insert new entry into linked list */
120:             rowhit[mfm]        = currentcol;
121:             rowhit[currentcol] = fm;
122:             fm                 = currentcol;
123:             /* fm points to present position in list since we know the columns are sorted */
124:           } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Detected invalid coloring");
125:         }
126:       }
127:       c->nrows[i] = nrows;
128:       PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);
129:       PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);
130:       /* now store the linked list of rows into c->rows[i] */
131:       nrows = 0;
132:       fm    = rowhit[N];
133:       do {
134:         c->rows[i][nrows]            = fm;
135:         c->columnsforrow[i][nrows++] = columnsforrow[fm];
136:         fm                           = rowhit[fm];
137:       } while (fm < N);
138:     } /* ---------------------------------------------------------------------------------------*/
139:     ISRestoreIndices(isa[i],&is);
140:   }
141:   MatRestoreColumnIJ(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&done);

143:   PetscFree(rowhit);
144:   PetscFree(columnsforrow);

146:   /* Optimize by adding the vscale, and scaleforrow[][] fields */
147:   /*
148:        see the version for MPIAIJ
149:   */
150:   VecCreateGhost(PetscObjectComm((PetscObject)mat),mat->rmap->n,PETSC_DETERMINE,0,NULL,&c->vscale);
151:   PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);
152:   for (k=0; k<c->ncolors; k++) {
153:     PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);
154:     for (l=0; l<c->nrows[k]; l++) {
155:       col = c->columnsforrow[k][l];

157:       c->vscaleforrow[k][l] = col;
158:     }
159:   }
160:   ISColoringRestoreIS(iscoloring,&isa);
161:   return(0);
162: }