1 2 #include <../src/mat/impls/aij/seq/aij.h> 3 #include <../src/mat/impls/baij/seq/baij.h> 4 5 /* 6 This routine is shared by SeqAIJ and SeqBAIJ matrices, 7 since it operators only on the nonzero structure of the elements or blocks. 8 */ 9 #undef __FUNCT__ 10 #define __FUNCT__ "MatFDColoringCreate_SeqXAIJ" 11 PetscErrorCode MatFDColoringCreate_SeqXAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c) 12 { 13 PetscErrorCode ierr; 14 PetscInt i,n,nrows,N,j,k,m,ncols,col,nis=iscoloring->n,*rowhit,bs,bs2,*spidx,nz; 15 const PetscInt *is,*row,*ci,*cj; 16 IS *isa; 17 PetscBool isBAIJ; 18 PetscScalar *A_val,**valaddrhit; 19 MatEntry *Jentry,*Jentry_new; 20 PetscInt *color_start,brows,bcols,nz_new,row_end,*row_start,*nrows_new; 21 22 PetscFunctionBegin; 23 ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr); 24 25 ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 26 ierr = PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&isBAIJ);CHKERRQ(ierr); 27 if (isBAIJ) { 28 Mat_SeqBAIJ *spA = (Mat_SeqBAIJ*)mat->data; 29 A_val = spA->a; 30 nz = spA->nz; 31 } else { 32 Mat_SeqAIJ *spA = (Mat_SeqAIJ*)mat->data; 33 A_val = spA->a; 34 nz = spA->nz; 35 bs = 1; /* only bs=1 is supported for SeqAIJ matrix */ 36 } 37 38 N = mat->cmap->N/bs; 39 c->M = mat->rmap->N/bs; /* set total rows, columns and local rows */ 40 c->N = mat->cmap->N/bs; 41 c->m = mat->rmap->N/bs; 42 c->rstart = 0; 43 44 c->ncolors = nis; 45 ierr = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr); 46 ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr); 47 ierr = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr); 48 ierr = PetscLogObjectMemory((PetscObject)c,3*nis*sizeof(PetscInt));CHKERRQ(ierr); 49 50 ierr = PetscMalloc(nz*sizeof(MatEntry),&Jentry);CHKERRQ(ierr); 51 ierr = PetscLogObjectMemory((PetscObject)c,nz*sizeof(MatEntry));CHKERRQ(ierr); 52 c->matentry = Jentry; 53 54 ierr = PetscMalloc2(nis+1,PetscInt,&color_start,nis+1,PetscInt,&row_start);CHKERRQ(ierr); 55 56 if (isBAIJ) { 57 ierr = MatGetColumnIJ_SeqBAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr); 58 } else { 59 ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr); 60 } 61 62 ierr = PetscMalloc2(c->m,PetscInt,&rowhit,c->m,PetscScalar*,&valaddrhit);CHKERRQ(ierr); 63 ierr = PetscMemzero(rowhit,c->m*sizeof(PetscInt));CHKERRQ(ierr); 64 65 // estimate bcol 66 PetscReal mem; 67 mem = nz*(sizeof(PetscScalar) + sizeof(PetscInt)) + 3*c->m*sizeof(PetscInt); 68 bcols = (PetscInt)(0.5*mem /(c->m*sizeof(PetscScalar))); 69 if (bcols > nis) bcols = nis; 70 brows = 1000/bcols; 71 72 nz = 0; 73 for (i=0; i<nis; i++) { /* loop over colors */ 74 color_start[i] = nz; 75 76 ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr); 77 ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr); 78 79 c->ncolumns[i] = n; 80 if (n) { 81 ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr); 82 ierr = PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));CHKERRQ(ierr); 83 ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr); 84 } else { 85 c->columns[i] = 0; 86 } 87 88 /* fast, crude version requires O(N*N) work */ 89 bs2 = bs*bs; 90 nrows = 0; 91 for (j=0; j<n; j++) { /* loop over columns */ 92 col = is[j]; 93 row = cj + ci[col]; 94 m = ci[col+1] - ci[col]; 95 nrows += m; 96 for (k=0; k<m; k++) { /* loop over columns marking them in rowhit */ 97 rowhit[*row] = col + 1; 98 valaddrhit[*row++] = &A_val[bs2*spidx[ci[col] + k]]; 99 } 100 } 101 c->nrows[i] = nrows; /* total num of rows for this color */ 102 103 for (j=0; j<N; j++) { /* loop over rows */ 104 if (rowhit[j]) { 105 Jentry[nz].row = j; /* local row index */ 106 Jentry[nz].col = rowhit[j] - 1; /* local column index */ 107 Jentry[nz].valaddr = valaddrhit[j]; /* address of mat value for this entry */ 108 nz++; 109 rowhit[j] = 0.0; /* zero rowhit for reuse */ 110 } 111 } 112 ierr = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr); 113 } 114 color_start[nis] = nz; 115 116 // ---------- reorder Jentry ------------ 117 ierr = PetscOptionsInt("-brows","The number of block rows","",brows,&brows,NULL);CHKERRQ(ierr); 118 if (!isBAIJ && brows) { 119 PetscInt nbcols = 0; 120 121 ierr = PetscMalloc(nz*sizeof(MatEntry),&Jentry_new);CHKERRQ(ierr); 122 for (i=0; i<nis; i++) row_start[i] = 0; 123 124 125 ierr = PetscOptionsInt("-bcols","The number of block columns","",bcols,&bcols,NULL);CHKERRQ(ierr); 126 if (bcols > nis) bcols = nis; 127 c->brows = brows; 128 c->bcols = bcols; 129 130 ierr = PetscMalloc(nis*sizeof(PetscInt),&nrows_new);CHKERRQ(ierr); 131 132 nz_new = 0; 133 for (i=0; i<nis; i+=bcols) { /* loop over colors */ 134 if (i + bcols > nis) bcols = nis - i; 135 //printf(" color %d, bcols %d\n",i,bcols); 136 137 row_end = brows; 138 m = mat->rmap->n; 139 if (row_end > m) row_end = m; 140 while (row_end <= m) { /* loop over block rows */ 141 for (j=0; j<bcols; j++) { /* loop over block columns */ 142 nrows = c->nrows[i+j]; 143 for (nz=color_start[i+j]; nz<color_start[i+j+1]; nz++) { /* for each Jentry */ 144 if (row_start[i+j] >= nrows) break; 145 if (Jentry[nz].row >= row_end) { 146 color_start[i+j] = nz; 147 break; 148 } else { 149 Jentry_new[nz_new].row = Jentry[nz].row + j*m; /* index in dy-array */ 150 Jentry_new[nz_new].col = Jentry[nz].col; 151 Jentry_new[nz_new].valaddr = Jentry[nz].valaddr; 152 nz_new++; 153 row_start[i+j]++; 154 } 155 } 156 } 157 if (row_end == m) break; 158 row_end += brows; 159 if (row_end > m) row_end = m; 160 } 161 nrows_new[nbcols++] = nz_new; 162 } 163 for (i=nbcols-1; i>0; i--) nrows_new[i] -= nrows_new[i-1]; 164 ierr = PetscFree(c->nrows);CHKERRQ(ierr); 165 c->nrows = nrows_new; 166 167 ierr = PetscFree(Jentry);CHKERRQ(ierr); 168 c->matentry = Jentry_new; 169 ierr = PetscMalloc(c->bcols*mat->rmap->n*sizeof(PetscScalar),&c->dy);CHKERRQ(ierr); 170 } 171 //--------------------------------------- 172 173 if (isBAIJ) { 174 ierr = MatRestoreColumnIJ_SeqBAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr); 175 ierr = PetscMalloc(bs*mat->rmap->n*sizeof(PetscScalar),&c->dy);CHKERRQ(ierr); 176 } else { 177 ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr); 178 } 179 ierr = PetscFree2(rowhit,valaddrhit);CHKERRQ(ierr); 180 ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr); 181 ierr = PetscFree2(color_start,row_start);CHKERRQ(ierr); 182 183 c->ctype = IS_COLORING_GHOSTED; 184 ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),mat->rmap->n,PETSC_DETERMINE,0,NULL,&c->vscale);CHKERRQ(ierr); 185 PetscFunctionReturn(0); 186 } 187