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__ 10*f86b9fbaSHong Zhang #define __FUNCT__ "MatFDColoringSetUp_SeqXAIJ" 11*f86b9fbaSHong Zhang PetscErrorCode MatFDColoringSetUp_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; 20*f86b9fbaSHong Zhang PetscInt *color_start,nz_new,row_end,*row_start,*nrows_new; 21*f86b9fbaSHong Zhang PetscInt brows,bcols; 2279c1e64dSHong Zhang 2379c1e64dSHong Zhang PetscFunctionBegin; 2479c1e64dSHong Zhang ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr); 2552011a10SHong Zhang 2679c1e64dSHong Zhang ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 27525d23c0SHong Zhang ierr = PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&isBAIJ);CHKERRQ(ierr); 28b312708bSHong Zhang if (isBAIJ) { 29b312708bSHong Zhang Mat_SeqBAIJ *spA = (Mat_SeqBAIJ*)mat->data; 30b312708bSHong Zhang A_val = spA->a; 31b312708bSHong Zhang nz = spA->nz; 32b312708bSHong Zhang } else { 33b312708bSHong Zhang Mat_SeqAIJ *spA = (Mat_SeqAIJ*)mat->data; 34b312708bSHong Zhang A_val = spA->a; 35b312708bSHong Zhang nz = spA->nz; 36b312708bSHong Zhang bs = 1; /* only bs=1 is supported for SeqAIJ matrix */ 3779c1e64dSHong Zhang } 38b312708bSHong Zhang 3979c1e64dSHong Zhang N = mat->cmap->N/bs; 4079c1e64dSHong Zhang c->M = mat->rmap->N/bs; /* set total rows, columns and local rows */ 4179c1e64dSHong Zhang c->N = mat->cmap->N/bs; 4279c1e64dSHong Zhang c->m = mat->rmap->N/bs; 4379c1e64dSHong Zhang c->rstart = 0; 4479c1e64dSHong Zhang 4579c1e64dSHong Zhang c->ncolors = nis; 4679c1e64dSHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr); 4779c1e64dSHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr); 4879c1e64dSHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr); 494b2e90caSHong Zhang ierr = PetscLogObjectMemory((PetscObject)c,3*nis*sizeof(PetscInt));CHKERRQ(ierr); 509a09c712SHong Zhang 51b312708bSHong Zhang ierr = PetscMalloc(nz*sizeof(MatEntry),&Jentry);CHKERRQ(ierr); 52b312708bSHong Zhang ierr = PetscLogObjectMemory((PetscObject)c,nz*sizeof(MatEntry));CHKERRQ(ierr); 539a09c712SHong Zhang c->matentry = Jentry; 5479c1e64dSHong Zhang 55e2c857f8SHong Zhang ierr = PetscMalloc2(nis+1,PetscInt,&color_start,nis+1,PetscInt,&row_start);CHKERRQ(ierr); 56e2c857f8SHong Zhang 57525d23c0SHong Zhang if (isBAIJ) { 58525d23c0SHong Zhang ierr = MatGetColumnIJ_SeqBAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr); 59525d23c0SHong Zhang } else { 606378f32dSHong Zhang ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr); 61525d23c0SHong Zhang } 6279c1e64dSHong Zhang 634b2e90caSHong Zhang ierr = PetscMalloc2(c->m,PetscInt,&rowhit,c->m,PetscScalar*,&valaddrhit);CHKERRQ(ierr); 6479c1e64dSHong Zhang ierr = PetscMemzero(rowhit,c->m*sizeof(PetscInt));CHKERRQ(ierr); 6579c1e64dSHong Zhang 666a509798SHong Zhang // estimate bcol 676a509798SHong Zhang PetscReal mem; 686a509798SHong Zhang mem = nz*(sizeof(PetscScalar) + sizeof(PetscInt)) + 3*c->m*sizeof(PetscInt); 696a509798SHong Zhang bcols = (PetscInt)(0.5*mem /(c->m*sizeof(PetscScalar))); 706a509798SHong Zhang if (bcols > nis) bcols = nis; 716a509798SHong Zhang brows = 1000/bcols; 726a509798SHong Zhang 7317a7dee1SHong Zhang nz = 0; 7479c1e64dSHong Zhang for (i=0; i<nis; i++) { /* loop over colors */ 75e2c857f8SHong Zhang color_start[i] = nz; 76e2c857f8SHong Zhang 7779c1e64dSHong Zhang ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr); 7879c1e64dSHong Zhang ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr); 7979c1e64dSHong Zhang 8079c1e64dSHong Zhang c->ncolumns[i] = n; 8179c1e64dSHong Zhang if (n) { 8279c1e64dSHong Zhang ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr); 834b2e90caSHong Zhang ierr = PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));CHKERRQ(ierr); 8479c1e64dSHong Zhang ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr); 8579c1e64dSHong Zhang } else { 8679c1e64dSHong Zhang c->columns[i] = 0; 8779c1e64dSHong Zhang } 8879c1e64dSHong Zhang 8979c1e64dSHong Zhang /* fast, crude version requires O(N*N) work */ 909a09c712SHong Zhang bs2 = bs*bs; 91476e0d0aSHong Zhang nrows = 0; 9279c1e64dSHong Zhang for (j=0; j<n; j++) { /* loop over columns */ 9379c1e64dSHong Zhang col = is[j]; 944b2e90caSHong Zhang row = cj + ci[col]; 9579c1e64dSHong Zhang m = ci[col+1] - ci[col]; 964b2e90caSHong Zhang nrows += m; 9779c1e64dSHong Zhang for (k=0; k<m; k++) { /* loop over columns marking them in rowhit */ 984b2e90caSHong Zhang rowhit[*row] = col + 1; 994b2e90caSHong Zhang valaddrhit[*row++] = &A_val[bs2*spidx[ci[col] + k]]; 10079c1e64dSHong Zhang } 10179c1e64dSHong Zhang } 102476e0d0aSHong Zhang c->nrows[i] = nrows; /* total num of rows for this color */ 10379c1e64dSHong Zhang 10479c1e64dSHong Zhang for (j=0; j<N; j++) { /* loop over rows */ 10579c1e64dSHong Zhang if (rowhit[j]) { 1069a09c712SHong Zhang Jentry[nz].row = j; /* local row index */ 1079a09c712SHong Zhang Jentry[nz].col = rowhit[j] - 1; /* local column index */ 1089a09c712SHong Zhang Jentry[nz].valaddr = valaddrhit[j]; /* address of mat value for this entry */ 1099a09c712SHong Zhang nz++; 11079c1e64dSHong Zhang rowhit[j] = 0.0; /* zero rowhit for reuse */ 11179c1e64dSHong Zhang } 11279c1e64dSHong Zhang } 11379c1e64dSHong Zhang ierr = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr); 11479c1e64dSHong Zhang } 115e2c857f8SHong Zhang color_start[nis] = nz; 116e2c857f8SHong Zhang 117e2c857f8SHong Zhang // ---------- reorder Jentry ------------ 118*f86b9fbaSHong Zhang if (!isBAIJ && bcols > 1) { 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; 135e2c857f8SHong Zhang 136e2c857f8SHong Zhang row_end = brows; 1376a509798SHong Zhang m = mat->rmap->n; 1386a509798SHong Zhang if (row_end > m) row_end = m; 1396a509798SHong Zhang while (row_end <= m) { /* loop over block rows */ 140e2c857f8SHong Zhang for (j=0; j<bcols; j++) { /* loop over block columns */ 141e2c857f8SHong Zhang nrows = c->nrows[i+j]; 142e2c857f8SHong Zhang for (nz=color_start[i+j]; nz<color_start[i+j+1]; nz++) { /* for each Jentry */ 143e2c857f8SHong Zhang if (row_start[i+j] >= nrows) break; 144e2c857f8SHong Zhang if (Jentry[nz].row >= row_end) { 145e2c857f8SHong Zhang color_start[i+j] = nz; 146e2c857f8SHong Zhang break; 147e2c857f8SHong Zhang } else { 1486a509798SHong Zhang Jentry_new[nz_new].row = Jentry[nz].row + j*m; /* index in dy-array */ 149d880da65SHong Zhang Jentry_new[nz_new].col = Jentry[nz].col; 150e2c857f8SHong Zhang Jentry_new[nz_new].valaddr = Jentry[nz].valaddr; 151e2c857f8SHong Zhang nz_new++; 152e2c857f8SHong Zhang row_start[i+j]++; 153e2c857f8SHong Zhang } 154e2c857f8SHong Zhang } 155e2c857f8SHong Zhang } 1566a509798SHong Zhang if (row_end == m) break; 157e2c857f8SHong Zhang row_end += brows; 1586a509798SHong Zhang if (row_end > m) row_end = m; 159e2c857f8SHong Zhang } 1606a509798SHong Zhang nrows_new[nbcols++] = nz_new; 161e2c857f8SHong Zhang } 1626a509798SHong Zhang for (i=nbcols-1; i>0; i--) nrows_new[i] -= nrows_new[i-1]; 1636a509798SHong Zhang ierr = PetscFree(c->nrows);CHKERRQ(ierr); 1646a509798SHong Zhang c->nrows = nrows_new; 165e2c857f8SHong Zhang 1666a509798SHong Zhang ierr = PetscFree(Jentry);CHKERRQ(ierr); 1676a509798SHong Zhang c->matentry = Jentry_new; 168e2c857f8SHong Zhang ierr = PetscMalloc(c->bcols*mat->rmap->n*sizeof(PetscScalar),&c->dy);CHKERRQ(ierr); 169e2c857f8SHong Zhang } 170e2c857f8SHong Zhang //--------------------------------------- 17185740eacSHong Zhang 172525d23c0SHong Zhang if (isBAIJ) { 173525d23c0SHong Zhang ierr = MatRestoreColumnIJ_SeqBAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr); 174e3225e9dSHong Zhang ierr = PetscMalloc(bs*mat->rmap->n*sizeof(PetscScalar),&c->dy);CHKERRQ(ierr); 175525d23c0SHong Zhang } else { 1766378f32dSHong Zhang ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr); 177525d23c0SHong Zhang } 1784b2e90caSHong Zhang ierr = PetscFree2(rowhit,valaddrhit);CHKERRQ(ierr); 17979c1e64dSHong Zhang ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr); 180e2c857f8SHong Zhang ierr = PetscFree2(color_start,row_start);CHKERRQ(ierr); 181476e0d0aSHong Zhang 1824737c552SHong Zhang c->ctype = IS_COLORING_GHOSTED; 18352011a10SHong Zhang ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),mat->rmap->n,PETSC_DETERMINE,0,NULL,&c->vscale);CHKERRQ(ierr); 18479c1e64dSHong Zhang PetscFunctionReturn(0); 18579c1e64dSHong Zhang } 186