xref: /petsc/src/mat/impls/aij/seq/fdaij.c (revision d880da652b9b0febbb366ced8120a3620a856c60)
170c3da92SBarry Smith 
2c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h>
3525d23c0SHong Zhang #include <../src/mat/impls/baij/seq/baij.h>
4525d23c0SHong Zhang 
5040ebd07SHong Zhang /*
6040ebd07SHong Zhang     This routine is shared by SeqAIJ and SeqBAIJ matrices,
7040ebd07SHong Zhang     since it operators only on the nonzero structure of the elements or blocks.
8040ebd07SHong Zhang */
9b582cc96SHong Zhang #undef __FUNCT__
10040ebd07SHong Zhang #define __FUNCT__ "MatFDColoringCreate_SeqXAIJ"
11040ebd07SHong Zhang PetscErrorCode MatFDColoringCreate_SeqXAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c)
1279c1e64dSHong Zhang {
1379c1e64dSHong Zhang   PetscErrorCode ierr;
14b312708bSHong Zhang   PetscInt       i,n,nrows,N,j,k,m,ncols,col,nis=iscoloring->n,*rowhit,bs,bs2,*spidx,nz;
154b2e90caSHong Zhang   const PetscInt *is,*row,*ci,*cj;
1679c1e64dSHong Zhang   IS             *isa;
17525d23c0SHong Zhang   PetscBool      isBAIJ;
18b312708bSHong Zhang   PetscScalar    *A_val,**valaddrhit;
19e2c857f8SHong Zhang   MatEntry       *Jentry,*Jentry_new;
206a509798SHong Zhang   PetscInt       *color_start,brows,bcols,nz_new,row_end,*row_start,*nrows_new;
2179c1e64dSHong Zhang 
2279c1e64dSHong Zhang   PetscFunctionBegin;
2379c1e64dSHong Zhang   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
2452011a10SHong Zhang 
2579c1e64dSHong Zhang   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
26525d23c0SHong Zhang   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&isBAIJ);CHKERRQ(ierr);
27b312708bSHong Zhang   if (isBAIJ) {
28b312708bSHong Zhang     Mat_SeqBAIJ *spA = (Mat_SeqBAIJ*)mat->data;
29b312708bSHong Zhang     A_val = spA->a;
30b312708bSHong Zhang     nz    = spA->nz;
31b312708bSHong Zhang   } else {
32b312708bSHong Zhang     Mat_SeqAIJ  *spA = (Mat_SeqAIJ*)mat->data;
33b312708bSHong Zhang     A_val = spA->a;
34b312708bSHong Zhang     nz    = spA->nz;
35b312708bSHong Zhang     bs    = 1; /* only bs=1 is supported for SeqAIJ matrix */
3679c1e64dSHong Zhang   }
37b312708bSHong Zhang 
3879c1e64dSHong Zhang   N         = mat->cmap->N/bs;
3979c1e64dSHong Zhang   c->M      = mat->rmap->N/bs;   /* set total rows, columns and local rows */
4079c1e64dSHong Zhang   c->N      = mat->cmap->N/bs;
4179c1e64dSHong Zhang   c->m      = mat->rmap->N/bs;
4279c1e64dSHong Zhang   c->rstart = 0;
4379c1e64dSHong Zhang 
4479c1e64dSHong Zhang   c->ncolors = nis;
4579c1e64dSHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr);
4679c1e64dSHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr);
4779c1e64dSHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr);
484b2e90caSHong Zhang   ierr       = PetscLogObjectMemory((PetscObject)c,3*nis*sizeof(PetscInt));CHKERRQ(ierr);
499a09c712SHong Zhang 
50b312708bSHong Zhang   ierr       = PetscMalloc(nz*sizeof(MatEntry),&Jentry);CHKERRQ(ierr);
51b312708bSHong Zhang   ierr       = PetscLogObjectMemory((PetscObject)c,nz*sizeof(MatEntry));CHKERRQ(ierr);
529a09c712SHong Zhang   c->matentry = Jentry;
5379c1e64dSHong Zhang 
54e2c857f8SHong Zhang   ierr       = PetscMalloc2(nis+1,PetscInt,&color_start,nis+1,PetscInt,&row_start);CHKERRQ(ierr);
55e2c857f8SHong Zhang 
56525d23c0SHong Zhang   if (isBAIJ) {
57525d23c0SHong Zhang     ierr = MatGetColumnIJ_SeqBAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
58525d23c0SHong Zhang   } else {
596378f32dSHong Zhang     ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
60525d23c0SHong Zhang   }
6179c1e64dSHong Zhang 
624b2e90caSHong Zhang   ierr = PetscMalloc2(c->m,PetscInt,&rowhit,c->m,PetscScalar*,&valaddrhit);CHKERRQ(ierr);
6379c1e64dSHong Zhang   ierr = PetscMemzero(rowhit,c->m*sizeof(PetscInt));CHKERRQ(ierr);
6479c1e64dSHong Zhang 
656a509798SHong Zhang   // estimate bcol
666a509798SHong Zhang   PetscReal mem;
676a509798SHong Zhang   mem = nz*(sizeof(PetscScalar) + sizeof(PetscInt)) + 3*c->m*sizeof(PetscInt);
686a509798SHong Zhang   bcols = (PetscInt)(0.5*mem /(c->m*sizeof(PetscScalar)));
696a509798SHong Zhang   if (bcols > nis) bcols = nis;
706a509798SHong Zhang   brows = 1000/bcols;
716a509798SHong Zhang 
7217a7dee1SHong Zhang   nz = 0;
7379c1e64dSHong Zhang   for (i=0; i<nis; i++) { /* loop over colors */
74e2c857f8SHong Zhang     color_start[i] = nz;
75e2c857f8SHong Zhang 
7679c1e64dSHong Zhang     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
7779c1e64dSHong Zhang     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
7879c1e64dSHong Zhang 
7979c1e64dSHong Zhang     c->ncolumns[i] = n;
8079c1e64dSHong Zhang     if (n) {
8179c1e64dSHong Zhang       ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr);
824b2e90caSHong Zhang       ierr = PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));CHKERRQ(ierr);
8379c1e64dSHong Zhang       ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr);
8479c1e64dSHong Zhang     } else {
8579c1e64dSHong Zhang       c->columns[i] = 0;
8679c1e64dSHong Zhang     }
8779c1e64dSHong Zhang 
8879c1e64dSHong Zhang     /* fast, crude version requires O(N*N) work */
899a09c712SHong Zhang     bs2   = bs*bs;
90476e0d0aSHong Zhang     nrows = 0;
9179c1e64dSHong Zhang     for (j=0; j<n; j++) {  /* loop over columns */
9279c1e64dSHong Zhang       col    = is[j];
934b2e90caSHong Zhang       row    = cj + ci[col];
9479c1e64dSHong Zhang       m      = ci[col+1] - ci[col];
954b2e90caSHong Zhang       nrows += m;
9679c1e64dSHong Zhang       for (k=0; k<m; k++) {  /* loop over columns marking them in rowhit */
974b2e90caSHong Zhang         rowhit[*row]       = col + 1;
984b2e90caSHong Zhang         valaddrhit[*row++] = &A_val[bs2*spidx[ci[col] + k]];
9979c1e64dSHong Zhang       }
10079c1e64dSHong Zhang     }
101476e0d0aSHong Zhang     c->nrows[i] = nrows; /* total num of rows for this color */
10279c1e64dSHong Zhang 
10379c1e64dSHong Zhang     for (j=0; j<N; j++) { /* loop over rows */
10479c1e64dSHong Zhang       if (rowhit[j]) {
1059a09c712SHong Zhang         Jentry[nz].row     = j;              /* local row index */
1069a09c712SHong Zhang         Jentry[nz].col     = rowhit[j] - 1;  /* local column index */
1079a09c712SHong Zhang         Jentry[nz].valaddr = valaddrhit[j];  /* address of mat value for this entry */
1089a09c712SHong Zhang         nz++;
10979c1e64dSHong Zhang         rowhit[j] = 0.0;                     /* zero rowhit for reuse */
11079c1e64dSHong Zhang       }
11179c1e64dSHong Zhang     }
11279c1e64dSHong Zhang     ierr = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr);
11379c1e64dSHong Zhang   }
114e2c857f8SHong Zhang   color_start[nis] = nz;
115e2c857f8SHong Zhang 
116e2c857f8SHong Zhang   // ---------- reorder Jentry ------------
117e2c857f8SHong Zhang   ierr = PetscOptionsInt("-brows","The number of block rows","",brows,&brows,NULL);CHKERRQ(ierr);
118e2c857f8SHong Zhang   if (!isBAIJ && brows) {
119e2c857f8SHong Zhang     PetscInt nbcols = 0;
1206a509798SHong Zhang 
1216a509798SHong Zhang     ierr = PetscMalloc(nz*sizeof(MatEntry),&Jentry_new);CHKERRQ(ierr);
1226a509798SHong Zhang     for (i=0; i<nis; i++) row_start[i] = 0;
1236a509798SHong Zhang 
1246a509798SHong Zhang 
125e2c857f8SHong Zhang     ierr = PetscOptionsInt("-bcols","The number of block columns","",bcols,&bcols,NULL);CHKERRQ(ierr);
126e2c857f8SHong Zhang     if (bcols > nis) bcols = nis;
127e2c857f8SHong Zhang     c->brows = brows;
128e2c857f8SHong Zhang     c->bcols = bcols;
129e2c857f8SHong Zhang 
1306a509798SHong Zhang     ierr = PetscMalloc(nis*sizeof(PetscInt),&nrows_new);CHKERRQ(ierr);
1316a509798SHong Zhang 
132e2c857f8SHong Zhang     nz_new  = 0;
133e2c857f8SHong Zhang     for (i=0; i<nis; i+=bcols) { /* loop over colors */
134e2c857f8SHong Zhang       if (i + bcols > nis) bcols = nis - i;
1356a509798SHong Zhang       //printf(" color %d, bcols %d\n",i,bcols);
136e2c857f8SHong Zhang 
137e2c857f8SHong Zhang       row_end = brows;
1386a509798SHong Zhang       m       = mat->rmap->n;
1396a509798SHong Zhang       if (row_end > m) row_end = m;
1406a509798SHong Zhang       while (row_end <= m) { /* loop over block rows */
141e2c857f8SHong Zhang         for (j=0; j<bcols; j++) {       /* loop over block columns */
142e2c857f8SHong Zhang           nrows = c->nrows[i+j];
143e2c857f8SHong Zhang           for (nz=color_start[i+j]; nz<color_start[i+j+1]; nz++) { /* for each Jentry */
144e2c857f8SHong Zhang             if (row_start[i+j] >= nrows) break;
145e2c857f8SHong Zhang             if (Jentry[nz].row >= row_end) {
146e2c857f8SHong Zhang               color_start[i+j] = nz;
147e2c857f8SHong Zhang               break;
148e2c857f8SHong Zhang             } else {
1496a509798SHong Zhang               Jentry_new[nz_new].row     = Jentry[nz].row + j*m; /* index in dy-array */
150*d880da65SHong Zhang               Jentry_new[nz_new].col     = Jentry[nz].col;
151e2c857f8SHong Zhang               Jentry_new[nz_new].valaddr = Jentry[nz].valaddr;
152e2c857f8SHong Zhang               nz_new++;
153e2c857f8SHong Zhang               row_start[i+j]++;
154e2c857f8SHong Zhang             }
155e2c857f8SHong Zhang           }
156e2c857f8SHong Zhang         }
1576a509798SHong Zhang         if (row_end == m) break;
158e2c857f8SHong Zhang         row_end += brows;
1596a509798SHong Zhang         if (row_end > m) row_end = m;
160e2c857f8SHong Zhang       }
1616a509798SHong Zhang       nrows_new[nbcols++] = nz_new;
162e2c857f8SHong Zhang     }
1636a509798SHong Zhang     for (i=nbcols-1; i>0; i--) nrows_new[i] -= nrows_new[i-1];
1646a509798SHong Zhang     ierr = PetscFree(c->nrows);CHKERRQ(ierr);
1656a509798SHong Zhang     c->nrows = nrows_new;
166e2c857f8SHong Zhang 
1676a509798SHong Zhang     ierr = PetscFree(Jentry);CHKERRQ(ierr);
1686a509798SHong Zhang     c->matentry = Jentry_new;
169e2c857f8SHong Zhang     ierr = PetscMalloc(c->bcols*mat->rmap->n*sizeof(PetscScalar),&c->dy);CHKERRQ(ierr);
170e2c857f8SHong Zhang   }
171e2c857f8SHong Zhang   //---------------------------------------
17285740eacSHong Zhang 
173525d23c0SHong Zhang   if (isBAIJ) {
174525d23c0SHong Zhang     ierr = MatRestoreColumnIJ_SeqBAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
175e3225e9dSHong Zhang     ierr = PetscMalloc(bs*mat->rmap->n*sizeof(PetscScalar),&c->dy);CHKERRQ(ierr);
176525d23c0SHong Zhang   } else {
1776378f32dSHong Zhang     ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
178525d23c0SHong Zhang   }
1794b2e90caSHong Zhang   ierr = PetscFree2(rowhit,valaddrhit);CHKERRQ(ierr);
18079c1e64dSHong Zhang   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
181e2c857f8SHong Zhang   ierr = PetscFree2(color_start,row_start);CHKERRQ(ierr);
182476e0d0aSHong Zhang 
1834737c552SHong Zhang   c->ctype = IS_COLORING_GHOSTED;
18452011a10SHong Zhang   ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),mat->rmap->n,PETSC_DETERMINE,0,NULL,&c->vscale);CHKERRQ(ierr);
18579c1e64dSHong Zhang   PetscFunctionReturn(0);
18679c1e64dSHong Zhang }
187