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; 14*a8971b87SHong 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); 21*a8971b87SHong Zhang if (isBAIJ) { 22*a8971b87SHong Zhang c->brows = m; 23*a8971b87SHong 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; 28*a8971b87SHong 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 5193dfae19SHong Zhang #undef __FUNCT__ 52*a8971b87SHong Zhang #define __FUNCT__ "MatFDColoringSetUpBlocked_AIJ_Private" 53*a8971b87SHong Zhang PetscErrorCode MatFDColoringSetUpBlocked_AIJ_Private(Mat mat,MatFDColoring c,PetscInt nz) 5479c1e64dSHong Zhang { 5579c1e64dSHong Zhang PetscErrorCode ierr; 56*a8971b87SHong Zhang PetscInt i,j,nrows,nbcols=0,brows=c->brows,bcols=c->bcols,mbs=c->m,*row_start,*nrows_new,nz_new,row_end,nis=c->ncolors; 57531e53bdSHong Zhang PetscInt *color_start; 58*a8971b87SHong Zhang MatEntry *Jentry_new,*Jentry=c->matentry; 5979c1e64dSHong Zhang 6079c1e64dSHong Zhang PetscFunctionBegin; 61*a8971b87SHong Zhang if (brows < 1 || brows > mbs) brows = mbs; 62*a8971b87SHong Zhang ierr = PetscMalloc((nis+1)*sizeof(PetscInt),&color_start);CHKERRQ(ierr); 63*a8971b87SHong Zhang color_start[0] = 0; 64*a8971b87SHong Zhang for (i=0; i<nis; i++) color_start[i+1] = c->nrows[i] + color_start[i]; 656a509798SHong Zhang 666a509798SHong Zhang ierr = PetscMalloc(nz*sizeof(MatEntry),&Jentry_new);CHKERRQ(ierr); 6793dfae19SHong Zhang ierr = PetscMalloc(bcols*sizeof(PetscInt),&row_start);CHKERRQ(ierr); 686a509798SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt),&nrows_new);CHKERRQ(ierr); 696a509798SHong Zhang 70e2c857f8SHong Zhang nz_new = 0; 71e2c857f8SHong Zhang for (i=0; i<nis; i+=bcols) { /* loop over colors */ 72e2c857f8SHong Zhang if (i + bcols > nis) bcols = nis - i; 73e2c857f8SHong Zhang 74e2c857f8SHong Zhang row_end = brows; 75c8a9c622SHong Zhang if (row_end > mbs) row_end = mbs; 7693dfae19SHong Zhang for (j=0; j<bcols; j++) row_start[j] = 0; 77c8a9c622SHong Zhang while (row_end <= mbs) { /* loop over block rows */ 78e2c857f8SHong Zhang for (j=0; j<bcols; j++) { /* loop over block columns */ 79e2c857f8SHong Zhang nrows = c->nrows[i+j]; 80c8a9c622SHong Zhang nz = color_start[i+j]; 81c8a9c622SHong Zhang while (row_start[j] < nrows) { 82e2c857f8SHong Zhang if (Jentry[nz].row >= row_end) { 83e2c857f8SHong Zhang color_start[i+j] = nz; 84e2c857f8SHong Zhang break; 85c8a9c622SHong Zhang } else { /* copy Jentry[nz] to Jentry_new[nz_new] */ 86c8a9c622SHong Zhang Jentry_new[nz_new].row = Jentry[nz].row + j*mbs; /* index in dy-array */ 87d880da65SHong Zhang Jentry_new[nz_new].col = Jentry[nz].col; 88e2c857f8SHong Zhang Jentry_new[nz_new].valaddr = Jentry[nz].valaddr; 89c8a9c622SHong Zhang nz_new++; nz++; row_start[j]++; 90e2c857f8SHong Zhang } 91e2c857f8SHong Zhang } 92e2c857f8SHong Zhang } 93c8a9c622SHong Zhang if (row_end == mbs) break; 94e2c857f8SHong Zhang row_end += brows; 95c8a9c622SHong Zhang if (row_end > mbs) row_end = mbs; 96e2c857f8SHong Zhang } 976a509798SHong Zhang nrows_new[nbcols++] = nz_new; 98e2c857f8SHong Zhang } 996a509798SHong Zhang for (i=nbcols-1; i>0; i--) nrows_new[i] -= nrows_new[i-1]; 1006a509798SHong Zhang ierr = PetscFree(c->nrows);CHKERRQ(ierr); 1016a509798SHong Zhang c->nrows = nrows_new; 102e2c857f8SHong Zhang 1036a509798SHong Zhang ierr = PetscFree(Jentry);CHKERRQ(ierr); 1046a509798SHong Zhang c->matentry = Jentry_new; 105e2c857f8SHong Zhang ierr = PetscMalloc(c->bcols*mat->rmap->n*sizeof(PetscScalar),&c->dy);CHKERRQ(ierr); 10693dfae19SHong Zhang ierr = PetscFree(row_start);CHKERRQ(ierr); 107*a8971b87SHong Zhang ierr = PetscFree(color_start);CHKERRQ(ierr); 108*a8971b87SHong Zhang PetscFunctionReturn(0); 109e2c857f8SHong Zhang } 110*a8971b87SHong Zhang 111*a8971b87SHong Zhang #undef __FUNCT__ 112*a8971b87SHong Zhang #define __FUNCT__ "MatFDColoringSetUp_SeqXAIJ" 113*a8971b87SHong Zhang PetscErrorCode MatFDColoringSetUp_SeqXAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c) 114*a8971b87SHong Zhang { 115*a8971b87SHong Zhang PetscErrorCode ierr; 116*a8971b87SHong Zhang PetscInt i,n,nrows,mbs=c->m,j,k,m,ncols,col,nis=iscoloring->n,*rowhit,bs,bs2,*spidx,nz; 117*a8971b87SHong Zhang const PetscInt *is,*row,*ci,*cj; 118*a8971b87SHong Zhang IS *isa; 119*a8971b87SHong Zhang PetscBool isBAIJ; 120*a8971b87SHong Zhang PetscScalar *A_val,**valaddrhit; 121*a8971b87SHong Zhang MatEntry *Jentry; 122*a8971b87SHong Zhang 123*a8971b87SHong Zhang PetscFunctionBegin; 124*a8971b87SHong Zhang ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr); 125*a8971b87SHong Zhang 126*a8971b87SHong Zhang ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 127*a8971b87SHong Zhang ierr = PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&isBAIJ);CHKERRQ(ierr); 128*a8971b87SHong Zhang if (isBAIJ) { 129*a8971b87SHong Zhang Mat_SeqBAIJ *spA = (Mat_SeqBAIJ*)mat->data; 130*a8971b87SHong Zhang A_val = spA->a; 131*a8971b87SHong Zhang nz = spA->nz; 132*a8971b87SHong Zhang } else { 133*a8971b87SHong Zhang Mat_SeqAIJ *spA = (Mat_SeqAIJ*)mat->data; 134*a8971b87SHong Zhang A_val = spA->a; 135*a8971b87SHong Zhang nz = spA->nz; 136*a8971b87SHong Zhang bs = 1; /* only bs=1 is supported for SeqAIJ matrix */ 137*a8971b87SHong Zhang } 138*a8971b87SHong Zhang 139*a8971b87SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr); 140*a8971b87SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr); 141*a8971b87SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr); 142*a8971b87SHong Zhang ierr = PetscLogObjectMemory((PetscObject)c,3*nis*sizeof(PetscInt));CHKERRQ(ierr); 143*a8971b87SHong Zhang 144*a8971b87SHong Zhang ierr = PetscMalloc(nz*sizeof(MatEntry),&Jentry);CHKERRQ(ierr); 145*a8971b87SHong Zhang ierr = PetscLogObjectMemory((PetscObject)c,nz*sizeof(MatEntry));CHKERRQ(ierr); 146*a8971b87SHong Zhang c->matentry = Jentry; 147*a8971b87SHong Zhang 148*a8971b87SHong Zhang if (isBAIJ) { 149*a8971b87SHong Zhang ierr = MatGetColumnIJ_SeqBAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr); 150*a8971b87SHong Zhang } else { 151*a8971b87SHong Zhang ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr); 152*a8971b87SHong Zhang } 153*a8971b87SHong Zhang 154*a8971b87SHong Zhang ierr = PetscMalloc2(c->m,PetscInt,&rowhit,c->m,PetscScalar*,&valaddrhit);CHKERRQ(ierr); 155*a8971b87SHong Zhang ierr = PetscMemzero(rowhit,c->m*sizeof(PetscInt));CHKERRQ(ierr); 156*a8971b87SHong Zhang 157*a8971b87SHong Zhang nz = 0; 158*a8971b87SHong Zhang for (i=0; i<nis; i++) { /* loop over colors */ 159*a8971b87SHong Zhang ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr); 160*a8971b87SHong Zhang ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr); 161*a8971b87SHong Zhang 162*a8971b87SHong Zhang c->ncolumns[i] = n; 163*a8971b87SHong Zhang if (n) { 164*a8971b87SHong Zhang ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr); 165*a8971b87SHong Zhang ierr = PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));CHKERRQ(ierr); 166*a8971b87SHong Zhang ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr); 167*a8971b87SHong Zhang } else { 168*a8971b87SHong Zhang c->columns[i] = 0; 169*a8971b87SHong Zhang } 170*a8971b87SHong Zhang 171*a8971b87SHong Zhang /* fast, crude version requires O(N*N) work */ 172*a8971b87SHong Zhang bs2 = bs*bs; 173*a8971b87SHong Zhang nrows = 0; 174*a8971b87SHong Zhang for (j=0; j<n; j++) { /* loop over columns */ 175*a8971b87SHong Zhang col = is[j]; 176*a8971b87SHong Zhang row = cj + ci[col]; 177*a8971b87SHong Zhang m = ci[col+1] - ci[col]; 178*a8971b87SHong Zhang nrows += m; 179*a8971b87SHong Zhang for (k=0; k<m; k++) { /* loop over columns marking them in rowhit */ 180*a8971b87SHong Zhang rowhit[*row] = col + 1; 181*a8971b87SHong Zhang valaddrhit[*row++] = &A_val[bs2*spidx[ci[col] + k]]; 182*a8971b87SHong Zhang } 183*a8971b87SHong Zhang } 184*a8971b87SHong Zhang c->nrows[i] = nrows; /* total num of rows for this color */ 185*a8971b87SHong Zhang 186*a8971b87SHong Zhang for (j=0; j<mbs; j++) { /* loop over rows */ 187*a8971b87SHong Zhang if (rowhit[j]) { 188*a8971b87SHong Zhang Jentry[nz].row = j; /* local row index */ 189*a8971b87SHong Zhang Jentry[nz].col = rowhit[j] - 1; /* local column index */ 190*a8971b87SHong Zhang Jentry[nz].valaddr = valaddrhit[j]; /* address of mat value for this entry */ 191*a8971b87SHong Zhang nz++; 192*a8971b87SHong Zhang rowhit[j] = 0.0; /* zero rowhit for reuse */ 193*a8971b87SHong Zhang } 194*a8971b87SHong Zhang } 195*a8971b87SHong Zhang ierr = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr); 196*a8971b87SHong Zhang } 197*a8971b87SHong Zhang 198*a8971b87SHong Zhang if (c->bcols > 1) { /* reorder Jentry for faster MatFDColoringApply() */ 199*a8971b87SHong Zhang ierr = MatFDColoringSetUpBlocked_AIJ_Private(mat,c,nz);CHKERRQ(ierr); 200*a8971b87SHong Zhang } 20185740eacSHong Zhang 202525d23c0SHong Zhang if (isBAIJ) { 203525d23c0SHong Zhang ierr = MatRestoreColumnIJ_SeqBAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr); 204e3225e9dSHong Zhang ierr = PetscMalloc(bs*mat->rmap->n*sizeof(PetscScalar),&c->dy);CHKERRQ(ierr); 205525d23c0SHong Zhang } else { 2066378f32dSHong Zhang ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr); 207525d23c0SHong Zhang } 208*a8971b87SHong Zhang ierr = PetscFree2(rowhit,valaddrhit);CHKERRQ(ierr); 20979c1e64dSHong Zhang ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr); 210476e0d0aSHong Zhang 21152011a10SHong Zhang ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),mat->rmap->n,PETSC_DETERMINE,0,NULL,&c->vscale);CHKERRQ(ierr); 21279c1e64dSHong Zhang PetscFunctionReturn(0); 21379c1e64dSHong Zhang } 214