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 bs,nis=iscoloring->n,m=mat->rmap->n; 15 PetscBool isBAIJ; 16 17 PetscFunctionBegin; 18 /* set default brows and bcols for speedup inserting the dense matrix into sparse Jacobian */ 19 ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 20 ierr = PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&isBAIJ);CHKERRQ(ierr); 21 if (isBAIJ) { 22 c->brows = m; 23 c->bcols = 1; 24 } else { /* seqaij matrix */ 25 /* bcols is chosen s.t. dy-array takes 50% of memory space as mat */ 26 Mat_SeqAIJ *spA = (Mat_SeqAIJ*)mat->data; 27 PetscReal mem; 28 PetscInt nz,brows,bcols; 29 30 bs = 1; /* only bs=1 is supported for SeqAIJ matrix */ 31 32 nz = spA->nz; 33 mem = nz*(sizeof(PetscScalar) + sizeof(PetscInt)) + 3*m*sizeof(PetscInt); 34 bcols = (PetscInt)(0.5*mem /(m*sizeof(PetscScalar))); 35 brows = 1000/bcols; 36 if (bcols > nis) bcols = nis; 37 if (brows == 0 || brows > m) brows = m; 38 c->brows = brows; 39 c->bcols = bcols; 40 } 41 42 c->M = mat->rmap->N/bs; /* set total rows, columns and local rows */ 43 c->N = mat->cmap->N/bs; 44 c->m = mat->rmap->N/bs; 45 c->rstart = 0; 46 c->ncolors = nis; 47 c->ctype = IS_COLORING_GHOSTED; 48 PetscFunctionReturn(0); 49 } 50 51 /* 52 Reorder Jentry such that blocked brows*bols of entries from dense matrix are inserted into sparse Jacobian 53 Input Parameters: 54 + mat - the matrix containing the nonzero structure of the Jacobian 55 . color - the coloring context 56 - nz - number of local non-zeros in mat 57 */ 58 #undef __FUNCT__ 59 #define __FUNCT__ "MatFDColoringSetUpBlocked_AIJ_Private" 60 PetscErrorCode MatFDColoringSetUpBlocked_AIJ_Private(Mat mat,MatFDColoring c,PetscInt nz) 61 { 62 PetscErrorCode ierr; 63 PetscInt i,j,nrows,nbcols,brows=c->brows,bcols=c->bcols,mbs=c->m,nis=c->ncolors; 64 PetscInt *color_start,*row_start,*nrows_new,nz_new,row_end; 65 MatEntry *Jentry_new,*Jentry=c->matentry; 66 67 PetscFunctionBegin; 68 if (brows < 1 || brows > mbs) brows = mbs; 69 ierr = PetscMalloc2(nis+1,PetscInt,&color_start,bcols,PetscInt,&row_start);CHKERRQ(ierr); 70 color_start[0] = 0; 71 for (i=0; i<nis; i++) color_start[i+1] = c->nrows[i] + color_start[i]; 72 73 ierr = PetscMalloc(nz*sizeof(MatEntry),&Jentry_new);CHKERRQ(ierr); 74 ierr = PetscMalloc(nis*sizeof(PetscInt),&nrows_new);CHKERRQ(ierr); 75 ierr = PetscMalloc(bcols*mat->rmap->n*sizeof(PetscScalar),&c->dy);CHKERRQ(ierr); 76 ierr = PetscLogObjectMemory((PetscObject)c,bcols*mat->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 77 78 nz_new = 0; 79 nbcols = 0; 80 for (i=0; i<nis; i+=bcols) { /* loop over colors */ 81 if (i + bcols > nis) bcols = nis - i; 82 83 row_end = brows; 84 if (row_end > mbs) row_end = mbs; 85 for (j=0; j<bcols; j++) row_start[j] = 0; 86 while (row_end <= mbs) { /* loop over block rows */ 87 for (j=0; j<bcols; j++) { /* loop over block columns */ 88 nrows = c->nrows[i+j]; 89 nz = color_start[i+j]; 90 while (row_start[j] < nrows) { 91 if (Jentry[nz].row >= row_end) { 92 color_start[i+j] = nz; 93 break; 94 } else { /* copy Jentry[nz] to Jentry_new[nz_new] */ 95 Jentry_new[nz_new].row = Jentry[nz].row + j*mbs; /* index in dy-array */ 96 Jentry_new[nz_new].col = Jentry[nz].col; 97 Jentry_new[nz_new].valaddr = Jentry[nz].valaddr; 98 nz_new++; nz++; row_start[j]++; 99 } 100 } 101 } 102 if (row_end == mbs) break; 103 row_end += brows; 104 if (row_end > mbs) row_end = mbs; 105 } 106 nrows_new[nbcols++] = nz_new; 107 } 108 ierr = PetscFree2(color_start,row_start);CHKERRQ(ierr); 109 110 for (i=nbcols-1; i>0; i--) nrows_new[i] -= nrows_new[i-1]; 111 ierr = PetscFree(c->nrows);CHKERRQ(ierr); 112 ierr = PetscFree(Jentry);CHKERRQ(ierr); 113 c->nrows = nrows_new; 114 c->matentry = Jentry_new; 115 PetscFunctionReturn(0); 116 } 117 118 #undef __FUNCT__ 119 #define __FUNCT__ "MatFDColoringSetUp_SeqXAIJ" 120 PetscErrorCode MatFDColoringSetUp_SeqXAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c) 121 { 122 PetscErrorCode ierr; 123 PetscInt i,n,nrows,mbs=c->m,j,k,m,ncols,col,nis=iscoloring->n,*rowhit,bs,bs2,*spidx,nz; 124 const PetscInt *is,*row,*ci,*cj; 125 IS *isa; 126 PetscBool isBAIJ; 127 PetscScalar *A_val,**valaddrhit; 128 MatEntry *Jentry; 129 130 PetscFunctionBegin; 131 ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr); 132 133 ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 134 ierr = PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&isBAIJ);CHKERRQ(ierr); 135 if (isBAIJ) { 136 Mat_SeqBAIJ *spA = (Mat_SeqBAIJ*)mat->data; 137 A_val = spA->a; 138 nz = spA->nz; 139 } else { 140 Mat_SeqAIJ *spA = (Mat_SeqAIJ*)mat->data; 141 A_val = spA->a; 142 nz = spA->nz; 143 bs = 1; /* only bs=1 is supported for SeqAIJ matrix */ 144 } 145 146 ierr = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr); 147 ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr); 148 ierr = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr); 149 ierr = PetscLogObjectMemory((PetscObject)c,3*nis*sizeof(PetscInt));CHKERRQ(ierr); 150 151 ierr = PetscMalloc(nz*sizeof(MatEntry),&Jentry);CHKERRQ(ierr); 152 ierr = PetscLogObjectMemory((PetscObject)c,nz*sizeof(MatEntry));CHKERRQ(ierr); 153 c->matentry = Jentry; 154 155 if (isBAIJ) { 156 ierr = MatGetColumnIJ_SeqBAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr); 157 } else { 158 ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr); 159 } 160 161 ierr = PetscMalloc2(c->m,PetscInt,&rowhit,c->m,PetscScalar*,&valaddrhit);CHKERRQ(ierr); 162 ierr = PetscMemzero(rowhit,c->m*sizeof(PetscInt));CHKERRQ(ierr); 163 164 nz = 0; 165 for (i=0; i<nis; i++) { /* loop over colors */ 166 ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr); 167 ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr); 168 169 c->ncolumns[i] = n; 170 if (n) { 171 ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr); 172 ierr = PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));CHKERRQ(ierr); 173 ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr); 174 } else { 175 c->columns[i] = 0; 176 } 177 178 /* fast, crude version requires O(N*N) work */ 179 bs2 = bs*bs; 180 nrows = 0; 181 for (j=0; j<n; j++) { /* loop over columns */ 182 col = is[j]; 183 row = cj + ci[col]; 184 m = ci[col+1] - ci[col]; 185 nrows += m; 186 for (k=0; k<m; k++) { /* loop over columns marking them in rowhit */ 187 rowhit[*row] = col + 1; 188 valaddrhit[*row++] = &A_val[bs2*spidx[ci[col] + k]]; 189 } 190 } 191 c->nrows[i] = nrows; /* total num of rows for this color */ 192 193 for (j=0; j<mbs; j++) { /* loop over rows */ 194 if (rowhit[j]) { 195 Jentry[nz].row = j; /* local row index */ 196 Jentry[nz].col = rowhit[j] - 1; /* local column index */ 197 Jentry[nz].valaddr = valaddrhit[j]; /* address of mat value for this entry */ 198 nz++; 199 rowhit[j] = 0.0; /* zero rowhit for reuse */ 200 } 201 } 202 ierr = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr); 203 } 204 205 if (c->bcols > 1) { /* reorder Jentry for faster MatFDColoringApply() */ 206 ierr = MatFDColoringSetUpBlocked_AIJ_Private(mat,c,nz);CHKERRQ(ierr); 207 } 208 209 if (isBAIJ) { 210 ierr = MatRestoreColumnIJ_SeqBAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr); 211 ierr = PetscMalloc(bs*mat->rmap->n*sizeof(PetscScalar),&c->dy);CHKERRQ(ierr); 212 ierr = PetscLogObjectMemory((PetscObject)c,bs*mat->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 213 } else { 214 ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr); 215 } 216 ierr = PetscFree2(rowhit,valaddrhit);CHKERRQ(ierr); 217 ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr); 218 219 ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),mat->rmap->n,PETSC_DETERMINE,0,NULL,&c->vscale);CHKERRQ(ierr); 220 #if defined(PETSC_USE_INFO) 221 ierr = PetscInfo3(c,"ncolors %D, brows %D and bcols %D are used.\n",c->ncolors,c->brows,c->bcols);CHKERRQ(ierr); 222 #endif 223 PetscFunctionReturn(0); 224 } 225