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; 19*e2c857f8SHong Zhang MatEntry *Jentry,*Jentry_new; 20*e2c857f8SHong Zhang PetscInt *color_start,brows=0,bcols=1,nz_new,row_end,*row_start; 2179c1e64dSHong Zhang 2279c1e64dSHong Zhang PetscFunctionBegin; 2379c1e64dSHong Zhang ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr); 2452011a10SHong Zhang 25040ebd07SHong Zhang /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */ 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); 49*e2c857f8SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt),&c->nrows_new);CHKERRQ(ierr); 504b2e90caSHong Zhang ierr = PetscLogObjectMemory((PetscObject)c,3*nis*sizeof(PetscInt));CHKERRQ(ierr); 519a09c712SHong Zhang 52b312708bSHong Zhang ierr = PetscMalloc(nz*sizeof(MatEntry),&Jentry);CHKERRQ(ierr); 53b312708bSHong Zhang ierr = PetscLogObjectMemory((PetscObject)c,nz*sizeof(MatEntry));CHKERRQ(ierr); 549a09c712SHong Zhang c->matentry = Jentry; 5579c1e64dSHong Zhang 56*e2c857f8SHong Zhang ierr = PetscMalloc2(nis+1,PetscInt,&color_start,nis+1,PetscInt,&row_start);CHKERRQ(ierr); 57*e2c857f8SHong Zhang ierr = PetscMalloc(nz*sizeof(MatEntry),&Jentry_new);CHKERRQ(ierr); 58*e2c857f8SHong Zhang c->matentry_new = Jentry_new; 59*e2c857f8SHong Zhang for (i=0; i<nis; i++) row_start[i] = 0; 60*e2c857f8SHong Zhang c->brows = brows; 61*e2c857f8SHong Zhang c->bcols = bcols; 62*e2c857f8SHong Zhang 63525d23c0SHong Zhang if (isBAIJ) { 64525d23c0SHong Zhang ierr = MatGetColumnIJ_SeqBAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr); 65525d23c0SHong Zhang } else { 666378f32dSHong Zhang ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr); 67525d23c0SHong Zhang } 6879c1e64dSHong Zhang 694b2e90caSHong Zhang ierr = PetscMalloc2(c->m,PetscInt,&rowhit,c->m,PetscScalar*,&valaddrhit);CHKERRQ(ierr); 7079c1e64dSHong Zhang ierr = PetscMemzero(rowhit,c->m*sizeof(PetscInt));CHKERRQ(ierr); 7179c1e64dSHong Zhang 7217a7dee1SHong Zhang nz = 0; 7379c1e64dSHong Zhang for (i=0; i<nis; i++) { /* loop over colors */ 74*e2c857f8SHong Zhang color_start[i] = nz; 75*e2c857f8SHong 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 } 114*e2c857f8SHong Zhang color_start[nis] = nz; 115*e2c857f8SHong Zhang 116*e2c857f8SHong Zhang // ---------- reorder Jentry ------------ 117*e2c857f8SHong Zhang ierr = PetscOptionsInt("-brows","The number of block rows","",brows,&brows,NULL);CHKERRQ(ierr); 118*e2c857f8SHong Zhang if (!isBAIJ && brows) { 119*e2c857f8SHong Zhang PetscInt nbcols = 0; 120*e2c857f8SHong Zhang ierr = PetscOptionsInt("-bcols","The number of block columns","",bcols,&bcols,NULL);CHKERRQ(ierr); 121*e2c857f8SHong Zhang if (bcols > nis) bcols = nis; 122*e2c857f8SHong Zhang c->brows = brows; 123*e2c857f8SHong Zhang c->bcols = bcols; 124*e2c857f8SHong Zhang 125*e2c857f8SHong Zhang nz_new = 0; 126*e2c857f8SHong Zhang for (i=0; i<nis; i+=bcols) { /* loop over colors */ 127*e2c857f8SHong Zhang if (i + bcols > nis) bcols = nis - i; 128*e2c857f8SHong Zhang printf(" color %d, bcols %d\n",i,bcols); 129*e2c857f8SHong Zhang 130*e2c857f8SHong Zhang row_end = brows; 131*e2c857f8SHong Zhang if (row_end > mat->rmap->n) row_end = mat->rmap->n; 132*e2c857f8SHong Zhang while (row_end <= mat->rmap->n) { /* loop over block rows */ 133*e2c857f8SHong Zhang for (j=0; j<bcols; j++) { /* loop over block columns */ 134*e2c857f8SHong Zhang nrows = c->nrows[i+j]; 135*e2c857f8SHong Zhang for (nz=color_start[i+j]; nz<color_start[i+j+1]; nz++) { /* for each Jentry */ 136*e2c857f8SHong Zhang if (row_start[i+j] >= nrows) break; 137*e2c857f8SHong Zhang if (Jentry[nz].row >= row_end) { 138*e2c857f8SHong Zhang color_start[i+j] = nz; 139*e2c857f8SHong Zhang break; 140*e2c857f8SHong Zhang } else { 141*e2c857f8SHong Zhang Jentry_new[nz_new].row = Jentry[nz].row; 142*e2c857f8SHong Zhang Jentry_new[nz_new].col = j; // j-th column in bcols 143*e2c857f8SHong Zhang Jentry_new[nz_new].valaddr = Jentry[nz].valaddr; 144*e2c857f8SHong Zhang nz_new++; 145*e2c857f8SHong Zhang row_start[i+j]++; 146*e2c857f8SHong Zhang } 147*e2c857f8SHong Zhang } 148*e2c857f8SHong Zhang } 149*e2c857f8SHong Zhang if (row_end == mat->rmap->n) break; 150*e2c857f8SHong Zhang row_end += brows; 151*e2c857f8SHong Zhang if (row_end > mat->rmap->n) row_end = mat->rmap->n; 152*e2c857f8SHong Zhang } 153*e2c857f8SHong Zhang c->nrows_new[nbcols++] = nz_new; 154*e2c857f8SHong Zhang } 155*e2c857f8SHong Zhang for (i=nbcols-1; i>0; i--) c->nrows_new[i] -= c->nrows_new[i-1]; 156*e2c857f8SHong Zhang 157*e2c857f8SHong Zhang ierr = PetscMalloc(c->bcols*mat->rmap->n*sizeof(PetscScalar),&c->dy);CHKERRQ(ierr); 158*e2c857f8SHong Zhang } 159*e2c857f8SHong Zhang //--------------------------------------- 16085740eacSHong Zhang 161525d23c0SHong Zhang if (isBAIJ) { 162525d23c0SHong Zhang ierr = MatRestoreColumnIJ_SeqBAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr); 163e3225e9dSHong Zhang ierr = PetscMalloc(bs*mat->rmap->n*sizeof(PetscScalar),&c->dy);CHKERRQ(ierr); 164525d23c0SHong Zhang } else { 1656378f32dSHong Zhang ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr); 166525d23c0SHong Zhang } 1674b2e90caSHong Zhang ierr = PetscFree2(rowhit,valaddrhit);CHKERRQ(ierr); 16879c1e64dSHong Zhang ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr); 169*e2c857f8SHong Zhang ierr = PetscFree2(color_start,row_start);CHKERRQ(ierr); 170476e0d0aSHong Zhang 1714737c552SHong Zhang c->ctype = IS_COLORING_GHOSTED; 17252011a10SHong Zhang ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),mat->rmap->n,PETSC_DETERMINE,0,NULL,&c->vscale);CHKERRQ(ierr); 17379c1e64dSHong Zhang PetscFunctionReturn(0); 17479c1e64dSHong Zhang } 175