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__ 1093dfae19SHong Zhang #define __FUNCT__ "MatFDColoringCreate_SeqXAIJ" 1193dfae19SHong Zhang PetscErrorCode MatFDColoringCreate_SeqXAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c) 1293dfae19SHong Zhang { 1393dfae19SHong Zhang PetscErrorCode ierr; 14a8971b87SHong Zhang PetscInt bs,nis=iscoloring->n,m=mat->rmap->n; 1593dfae19SHong Zhang PetscBool isBAIJ; 1693dfae19SHong Zhang 1793dfae19SHong Zhang PetscFunctionBegin; 18531e53bdSHong Zhang /* set default brows and bcols for speedup inserting the dense matrix into sparse Jacobian */ 1993dfae19SHong Zhang ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 2093dfae19SHong Zhang ierr = PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&isBAIJ);CHKERRQ(ierr); 21a8971b87SHong Zhang if (isBAIJ) { 22a8971b87SHong Zhang c->brows = m; 23a8971b87SHong Zhang c->bcols = 1; 24531e53bdSHong Zhang } else { /* seqaij matrix */ 25531e53bdSHong Zhang /* bcols is chosen s.t. dy-array takes 50% of memory space as mat */ 2693dfae19SHong Zhang Mat_SeqAIJ *spA = (Mat_SeqAIJ*)mat->data; 27531e53bdSHong Zhang PetscReal mem; 28a8971b87SHong Zhang PetscInt nz,brows,bcols; 29531e53bdSHong Zhang 3093dfae19SHong Zhang bs = 1; /* only bs=1 is supported for SeqAIJ matrix */ 31531e53bdSHong Zhang 32531e53bdSHong Zhang nz = spA->nz; 33531e53bdSHong Zhang mem = nz*(sizeof(PetscScalar) + sizeof(PetscInt)) + 3*m*sizeof(PetscInt); 34531e53bdSHong Zhang bcols = (PetscInt)(0.5*mem /(m*sizeof(PetscScalar))); 35531e53bdSHong Zhang brows = 1000/bcols; 36531e53bdSHong Zhang if (bcols > nis) bcols = nis; 37531e53bdSHong Zhang if (brows == 0 || brows > m) brows = m; 38531e53bdSHong Zhang c->brows = brows; 39531e53bdSHong Zhang c->bcols = bcols; 4093dfae19SHong Zhang } 41531e53bdSHong Zhang 4293dfae19SHong Zhang c->M = mat->rmap->N/bs; /* set total rows, columns and local rows */ 4393dfae19SHong Zhang c->N = mat->cmap->N/bs; 4493dfae19SHong Zhang c->m = mat->rmap->N/bs; 4593dfae19SHong Zhang c->rstart = 0; 4693dfae19SHong Zhang c->ncolors = nis; 4793dfae19SHong Zhang c->ctype = IS_COLORING_GHOSTED; 4893dfae19SHong Zhang PetscFunctionReturn(0); 4993dfae19SHong Zhang } 5093dfae19SHong Zhang 51*b3e1f37bSHong Zhang /* 52*b3e1f37bSHong Zhang Reorder Jentry such that blocked brows*bols of entries from dense matrix are inserted into sparse Jacobian 53*b3e1f37bSHong Zhang Input Parameters: 54*b3e1f37bSHong Zhang + mat - the matrix containing the nonzero structure of the Jacobian 55*b3e1f37bSHong Zhang . color - the coloring context 56*b3e1f37bSHong Zhang - nz - number of local non-zeros in mat 57*b3e1f37bSHong Zhang */ 5893dfae19SHong Zhang #undef __FUNCT__ 59a8971b87SHong Zhang #define __FUNCT__ "MatFDColoringSetUpBlocked_AIJ_Private" 60a8971b87SHong Zhang PetscErrorCode MatFDColoringSetUpBlocked_AIJ_Private(Mat mat,MatFDColoring c,PetscInt nz) 6179c1e64dSHong Zhang { 6279c1e64dSHong Zhang PetscErrorCode ierr; 63*b3e1f37bSHong Zhang PetscInt i,j,nrows,nbcols,brows=c->brows,bcols=c->bcols,mbs=c->m,nis=c->ncolors; 64*b3e1f37bSHong Zhang PetscInt *color_start,*row_start,*nrows_new,nz_new,row_end; 65a8971b87SHong Zhang MatEntry *Jentry_new,*Jentry=c->matentry; 6679c1e64dSHong Zhang 6779c1e64dSHong Zhang PetscFunctionBegin; 68a8971b87SHong Zhang if (brows < 1 || brows > mbs) brows = mbs; 69*b3e1f37bSHong Zhang ierr = PetscMalloc2(nis+1,PetscInt,&color_start,bcols,PetscInt,&row_start);CHKERRQ(ierr); 70a8971b87SHong Zhang color_start[0] = 0; 71a8971b87SHong Zhang for (i=0; i<nis; i++) color_start[i+1] = c->nrows[i] + color_start[i]; 726a509798SHong Zhang 736a509798SHong Zhang ierr = PetscMalloc(nz*sizeof(MatEntry),&Jentry_new);CHKERRQ(ierr); 746a509798SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt),&nrows_new);CHKERRQ(ierr); 75*b3e1f37bSHong Zhang ierr = PetscMalloc(bcols*mat->rmap->n*sizeof(PetscScalar),&c->dy);CHKERRQ(ierr); 76*b3e1f37bSHong Zhang ierr = PetscLogObjectMemory((PetscObject)c,bcols*mat->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 776a509798SHong Zhang 78e2c857f8SHong Zhang nz_new = 0; 79*b3e1f37bSHong Zhang nbcols = 0; 80e2c857f8SHong Zhang for (i=0; i<nis; i+=bcols) { /* loop over colors */ 81e2c857f8SHong Zhang if (i + bcols > nis) bcols = nis - i; 82e2c857f8SHong Zhang 83e2c857f8SHong Zhang row_end = brows; 84c8a9c622SHong Zhang if (row_end > mbs) row_end = mbs; 8593dfae19SHong Zhang for (j=0; j<bcols; j++) row_start[j] = 0; 86c8a9c622SHong Zhang while (row_end <= mbs) { /* loop over block rows */ 87e2c857f8SHong Zhang for (j=0; j<bcols; j++) { /* loop over block columns */ 88e2c857f8SHong Zhang nrows = c->nrows[i+j]; 89c8a9c622SHong Zhang nz = color_start[i+j]; 90c8a9c622SHong Zhang while (row_start[j] < nrows) { 91e2c857f8SHong Zhang if (Jentry[nz].row >= row_end) { 92e2c857f8SHong Zhang color_start[i+j] = nz; 93e2c857f8SHong Zhang break; 94c8a9c622SHong Zhang } else { /* copy Jentry[nz] to Jentry_new[nz_new] */ 95c8a9c622SHong Zhang Jentry_new[nz_new].row = Jentry[nz].row + j*mbs; /* index in dy-array */ 96d880da65SHong Zhang Jentry_new[nz_new].col = Jentry[nz].col; 97e2c857f8SHong Zhang Jentry_new[nz_new].valaddr = Jentry[nz].valaddr; 98c8a9c622SHong Zhang nz_new++; nz++; row_start[j]++; 99e2c857f8SHong Zhang } 100e2c857f8SHong Zhang } 101e2c857f8SHong Zhang } 102c8a9c622SHong Zhang if (row_end == mbs) break; 103e2c857f8SHong Zhang row_end += brows; 104c8a9c622SHong Zhang if (row_end > mbs) row_end = mbs; 105e2c857f8SHong Zhang } 1066a509798SHong Zhang nrows_new[nbcols++] = nz_new; 107e2c857f8SHong Zhang } 108*b3e1f37bSHong Zhang ierr = PetscFree2(color_start,row_start);CHKERRQ(ierr); 109*b3e1f37bSHong Zhang 1106a509798SHong Zhang for (i=nbcols-1; i>0; i--) nrows_new[i] -= nrows_new[i-1]; 1116a509798SHong Zhang ierr = PetscFree(c->nrows);CHKERRQ(ierr); 1126a509798SHong Zhang ierr = PetscFree(Jentry);CHKERRQ(ierr); 113*b3e1f37bSHong Zhang c->nrows = nrows_new; 1146a509798SHong Zhang c->matentry = Jentry_new; 115a8971b87SHong Zhang PetscFunctionReturn(0); 116e2c857f8SHong Zhang } 117a8971b87SHong Zhang 118a8971b87SHong Zhang #undef __FUNCT__ 119a8971b87SHong Zhang #define __FUNCT__ "MatFDColoringSetUp_SeqXAIJ" 120a8971b87SHong Zhang PetscErrorCode MatFDColoringSetUp_SeqXAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c) 121a8971b87SHong Zhang { 122a8971b87SHong Zhang PetscErrorCode ierr; 123a8971b87SHong Zhang PetscInt i,n,nrows,mbs=c->m,j,k,m,ncols,col,nis=iscoloring->n,*rowhit,bs,bs2,*spidx,nz; 124a8971b87SHong Zhang const PetscInt *is,*row,*ci,*cj; 125a8971b87SHong Zhang IS *isa; 126a8971b87SHong Zhang PetscBool isBAIJ; 127a8971b87SHong Zhang PetscScalar *A_val,**valaddrhit; 128a8971b87SHong Zhang MatEntry *Jentry; 129a8971b87SHong Zhang 130a8971b87SHong Zhang PetscFunctionBegin; 131a8971b87SHong Zhang ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr); 132a8971b87SHong Zhang 133a8971b87SHong Zhang ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 134a8971b87SHong Zhang ierr = PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&isBAIJ);CHKERRQ(ierr); 135a8971b87SHong Zhang if (isBAIJ) { 136a8971b87SHong Zhang Mat_SeqBAIJ *spA = (Mat_SeqBAIJ*)mat->data; 137a8971b87SHong Zhang A_val = spA->a; 138a8971b87SHong Zhang nz = spA->nz; 139a8971b87SHong Zhang } else { 140a8971b87SHong Zhang Mat_SeqAIJ *spA = (Mat_SeqAIJ*)mat->data; 141a8971b87SHong Zhang A_val = spA->a; 142a8971b87SHong Zhang nz = spA->nz; 143a8971b87SHong Zhang bs = 1; /* only bs=1 is supported for SeqAIJ matrix */ 144a8971b87SHong Zhang } 145a8971b87SHong Zhang 146a8971b87SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr); 147a8971b87SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr); 148a8971b87SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr); 149a8971b87SHong Zhang ierr = PetscLogObjectMemory((PetscObject)c,3*nis*sizeof(PetscInt));CHKERRQ(ierr); 150a8971b87SHong Zhang 151a8971b87SHong Zhang ierr = PetscMalloc(nz*sizeof(MatEntry),&Jentry);CHKERRQ(ierr); 152a8971b87SHong Zhang ierr = PetscLogObjectMemory((PetscObject)c,nz*sizeof(MatEntry));CHKERRQ(ierr); 153a8971b87SHong Zhang c->matentry = Jentry; 154a8971b87SHong Zhang 155a8971b87SHong Zhang if (isBAIJ) { 156a8971b87SHong Zhang ierr = MatGetColumnIJ_SeqBAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr); 157a8971b87SHong Zhang } else { 158a8971b87SHong Zhang ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr); 159a8971b87SHong Zhang } 160a8971b87SHong Zhang 161a8971b87SHong Zhang ierr = PetscMalloc2(c->m,PetscInt,&rowhit,c->m,PetscScalar*,&valaddrhit);CHKERRQ(ierr); 162a8971b87SHong Zhang ierr = PetscMemzero(rowhit,c->m*sizeof(PetscInt));CHKERRQ(ierr); 163a8971b87SHong Zhang 164a8971b87SHong Zhang nz = 0; 165a8971b87SHong Zhang for (i=0; i<nis; i++) { /* loop over colors */ 166a8971b87SHong Zhang ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr); 167a8971b87SHong Zhang ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr); 168a8971b87SHong Zhang 169a8971b87SHong Zhang c->ncolumns[i] = n; 170a8971b87SHong Zhang if (n) { 171a8971b87SHong Zhang ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr); 172a8971b87SHong Zhang ierr = PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));CHKERRQ(ierr); 173a8971b87SHong Zhang ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr); 174a8971b87SHong Zhang } else { 175a8971b87SHong Zhang c->columns[i] = 0; 176a8971b87SHong Zhang } 177a8971b87SHong Zhang 178a8971b87SHong Zhang /* fast, crude version requires O(N*N) work */ 179a8971b87SHong Zhang bs2 = bs*bs; 180a8971b87SHong Zhang nrows = 0; 181a8971b87SHong Zhang for (j=0; j<n; j++) { /* loop over columns */ 182a8971b87SHong Zhang col = is[j]; 183a8971b87SHong Zhang row = cj + ci[col]; 184a8971b87SHong Zhang m = ci[col+1] - ci[col]; 185a8971b87SHong Zhang nrows += m; 186a8971b87SHong Zhang for (k=0; k<m; k++) { /* loop over columns marking them in rowhit */ 187a8971b87SHong Zhang rowhit[*row] = col + 1; 188a8971b87SHong Zhang valaddrhit[*row++] = &A_val[bs2*spidx[ci[col] + k]]; 189a8971b87SHong Zhang } 190a8971b87SHong Zhang } 191a8971b87SHong Zhang c->nrows[i] = nrows; /* total num of rows for this color */ 192a8971b87SHong Zhang 193a8971b87SHong Zhang for (j=0; j<mbs; j++) { /* loop over rows */ 194a8971b87SHong Zhang if (rowhit[j]) { 195a8971b87SHong Zhang Jentry[nz].row = j; /* local row index */ 196a8971b87SHong Zhang Jentry[nz].col = rowhit[j] - 1; /* local column index */ 197a8971b87SHong Zhang Jentry[nz].valaddr = valaddrhit[j]; /* address of mat value for this entry */ 198a8971b87SHong Zhang nz++; 199a8971b87SHong Zhang rowhit[j] = 0.0; /* zero rowhit for reuse */ 200a8971b87SHong Zhang } 201a8971b87SHong Zhang } 202a8971b87SHong Zhang ierr = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr); 203a8971b87SHong Zhang } 204a8971b87SHong Zhang 205a8971b87SHong Zhang if (c->bcols > 1) { /* reorder Jentry for faster MatFDColoringApply() */ 206a8971b87SHong Zhang ierr = MatFDColoringSetUpBlocked_AIJ_Private(mat,c,nz);CHKERRQ(ierr); 207a8971b87SHong Zhang } 20885740eacSHong Zhang 209525d23c0SHong Zhang if (isBAIJ) { 210525d23c0SHong Zhang ierr = MatRestoreColumnIJ_SeqBAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr); 211e3225e9dSHong Zhang ierr = PetscMalloc(bs*mat->rmap->n*sizeof(PetscScalar),&c->dy);CHKERRQ(ierr); 212*b3e1f37bSHong Zhang ierr = PetscLogObjectMemory((PetscObject)c,bs*mat->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 213525d23c0SHong Zhang } else { 2146378f32dSHong Zhang ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr); 215525d23c0SHong Zhang } 216a8971b87SHong Zhang ierr = PetscFree2(rowhit,valaddrhit);CHKERRQ(ierr); 21779c1e64dSHong Zhang ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr); 218476e0d0aSHong Zhang 21952011a10SHong Zhang ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),mat->rmap->n,PETSC_DETERMINE,0,NULL,&c->vscale);CHKERRQ(ierr); 220*b3e1f37bSHong Zhang #if defined(PETSC_USE_INFO) 221*b3e1f37bSHong Zhang ierr = PetscInfo3(c,"ncolors %D, brows %D and bcols %D are used.\n",c->ncolors,c->brows,c->bcols);CHKERRQ(ierr); 222*b3e1f37bSHong Zhang #endif 22379c1e64dSHong Zhang PetscFunctionReturn(0); 22479c1e64dSHong Zhang } 225