xref: /petsc/src/mat/impls/aij/seq/fdaij.c (revision 040ebd07df7f48cdffac40673e8597e52ed39091)
170c3da92SBarry Smith 
2c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h>
3525d23c0SHong Zhang #include <../src/mat/impls/baij/seq/baij.h>
4525d23c0SHong Zhang 
5*040ebd07SHong Zhang /*
6*040ebd07SHong Zhang     This routine is shared by SeqAIJ and SeqBAIJ matrices,
7*040ebd07SHong Zhang     since it operators only on the nonzero structure of the elements or blocks.
8*040ebd07SHong Zhang */
9b582cc96SHong Zhang #undef __FUNCT__
10*040ebd07SHong Zhang #define __FUNCT__ "MatFDColoringCreate_SeqXAIJ"
11*040ebd07SHong Zhang PetscErrorCode MatFDColoringCreate_SeqXAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c)
1279c1e64dSHong Zhang {
1379c1e64dSHong Zhang   PetscErrorCode ierr;
1479c1e64dSHong Zhang   PetscInt       i,n,nrows,N,j,k,m,ncols,col;
154b2e90caSHong Zhang   const PetscInt *is,*row,*ci,*cj;
164b2e90caSHong Zhang   PetscInt       nis=iscoloring->n,*rowhit,bs,bs2,*spidx,nz;
1779c1e64dSHong Zhang   IS             *isa;
18525d23c0SHong Zhang   PetscBool      isBAIJ;
194b2e90caSHong Zhang   Mat_SeqAIJ     *spA = (Mat_SeqAIJ*)mat->data;
204b2e90caSHong Zhang   PetscScalar    *A_val=spA->a;
219a09c712SHong Zhang   PetscScalar    **valaddrhit;
229a09c712SHong Zhang   MatEntry       *Jentry;
2379c1e64dSHong Zhang 
2479c1e64dSHong Zhang   PetscFunctionBegin;
2579c1e64dSHong Zhang   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
2652011a10SHong Zhang 
27*040ebd07SHong Zhang   /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */
2879c1e64dSHong Zhang   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
29525d23c0SHong Zhang   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&isBAIJ);CHKERRQ(ierr);
30525d23c0SHong Zhang   if (!isBAIJ) {
3152011a10SHong Zhang     bs = 1; /* only bs=1 is supported for non SEQBAIJ matrix */
3279c1e64dSHong Zhang   }
3379c1e64dSHong Zhang   N         = mat->cmap->N/bs;
3479c1e64dSHong Zhang   c->M      = mat->rmap->N/bs;   /* set total rows, columns and local rows */
3579c1e64dSHong Zhang   c->N      = mat->cmap->N/bs;
3679c1e64dSHong Zhang   c->m      = mat->rmap->N/bs;
3779c1e64dSHong Zhang   c->rstart = 0;
3879c1e64dSHong Zhang 
3979c1e64dSHong Zhang   c->ncolors = nis;
4079c1e64dSHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr);
4179c1e64dSHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr);
4279c1e64dSHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr);
434b2e90caSHong Zhang   ierr       = PetscLogObjectMemory((PetscObject)c,3*nis*sizeof(PetscInt));CHKERRQ(ierr);
449a09c712SHong Zhang 
454b2e90caSHong Zhang   ierr       = PetscMalloc(spA->nz*sizeof(MatEntry),&Jentry);CHKERRQ(ierr);
464b2e90caSHong Zhang   ierr       = PetscLogObjectMemory((PetscObject)c,spA->nz*sizeof(MatEntry));CHKERRQ(ierr);
479a09c712SHong Zhang   c->matentry = Jentry;
4879c1e64dSHong Zhang 
49525d23c0SHong Zhang   if (isBAIJ) {
50525d23c0SHong Zhang     ierr = MatGetColumnIJ_SeqBAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
51525d23c0SHong Zhang   } else {
526378f32dSHong Zhang     ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
53525d23c0SHong Zhang   }
5479c1e64dSHong Zhang 
554b2e90caSHong Zhang   ierr = PetscMalloc2(c->m,PetscInt,&rowhit,c->m,PetscScalar*,&valaddrhit);CHKERRQ(ierr);
5679c1e64dSHong Zhang   ierr = PetscMemzero(rowhit,c->m*sizeof(PetscInt));CHKERRQ(ierr);
5779c1e64dSHong Zhang 
5817a7dee1SHong Zhang   nz = 0;
5979c1e64dSHong Zhang   for (i=0; i<nis; i++) { /* loop over colors */
6079c1e64dSHong Zhang     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
6179c1e64dSHong Zhang     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
6279c1e64dSHong Zhang 
6379c1e64dSHong Zhang     c->ncolumns[i] = n;
6479c1e64dSHong Zhang     if (n) {
6579c1e64dSHong Zhang       ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr);
664b2e90caSHong Zhang       ierr = PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));CHKERRQ(ierr);
6779c1e64dSHong Zhang       ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr);
6879c1e64dSHong Zhang     } else {
6979c1e64dSHong Zhang       c->columns[i] = 0;
7079c1e64dSHong Zhang     }
7179c1e64dSHong Zhang 
7279c1e64dSHong Zhang     /* fast, crude version requires O(N*N) work */
739a09c712SHong Zhang     bs2   = bs*bs;
74476e0d0aSHong Zhang     nrows = 0;
7579c1e64dSHong Zhang     for (j=0; j<n; j++) {  /* loop over columns */
7679c1e64dSHong Zhang       col    = is[j];
774b2e90caSHong Zhang       row    = cj + ci[col];
7879c1e64dSHong Zhang       m      = ci[col+1] - ci[col];
794b2e90caSHong Zhang       nrows += m;
8079c1e64dSHong Zhang       for (k=0; k<m; k++) {  /* loop over columns marking them in rowhit */
814b2e90caSHong Zhang         rowhit[*row]       = col + 1;
824b2e90caSHong Zhang         valaddrhit[*row++] = &A_val[bs2*spidx[ci[col] + k]];
8379c1e64dSHong Zhang       }
8479c1e64dSHong Zhang     }
85476e0d0aSHong Zhang     c->nrows[i] = nrows; /* total num of rows for this color */
8679c1e64dSHong Zhang 
8779c1e64dSHong Zhang     for (j=0; j<N; j++) { /* loop over rows */
8879c1e64dSHong Zhang       if (rowhit[j]) {
899a09c712SHong Zhang         Jentry[nz].row     = j;              /* local row index */
909a09c712SHong Zhang         Jentry[nz].col     = rowhit[j] - 1;  /* local column index */
919a09c712SHong Zhang         Jentry[nz].valaddr = valaddrhit[j];  /* address of mat value for this entry */
929a09c712SHong Zhang         nz++;
9379c1e64dSHong Zhang         rowhit[j] = 0.0;                     /* zero rowhit for reuse */
9479c1e64dSHong Zhang       }
9579c1e64dSHong Zhang     }
9679c1e64dSHong Zhang     ierr = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr);
9779c1e64dSHong Zhang   }
984b2e90caSHong Zhang   if (nz != spA->nz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"nz %d != mat->nz %d\n",nz,spA->nz);
9985740eacSHong Zhang 
100525d23c0SHong Zhang   if (isBAIJ) {
101525d23c0SHong Zhang     ierr = MatRestoreColumnIJ_SeqBAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
102525d23c0SHong Zhang   } else {
1036378f32dSHong Zhang     ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
104525d23c0SHong Zhang   }
1054b2e90caSHong Zhang   ierr = PetscFree2(rowhit,valaddrhit);CHKERRQ(ierr);
10679c1e64dSHong Zhang   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
107476e0d0aSHong Zhang 
1084737c552SHong Zhang   c->ctype = IS_COLORING_GHOSTED;
109b582cc96SHong Zhang   if (isBAIJ) {
110b7799381SHong Zhang     ierr = PetscMalloc(bs*mat->cmap->N*sizeof(PetscScalar),&c->dy);CHKERRQ(ierr);
111b582cc96SHong Zhang   }
11252011a10SHong Zhang   ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),mat->rmap->n,PETSC_DETERMINE,0,NULL,&c->vscale);CHKERRQ(ierr);
11379c1e64dSHong Zhang   PetscFunctionReturn(0);
11479c1e64dSHong Zhang }
115