xref: /petsc/src/mat/impls/aij/seq/fdaij.c (revision e3225e9dc71240658ff8694e906a72dc65d8b579)
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__
10040ebd07SHong Zhang #define __FUNCT__ "MatFDColoringCreate_SeqXAIJ"
11040ebd07SHong Zhang PetscErrorCode MatFDColoringCreate_SeqXAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c)
1279c1e64dSHong Zhang {
1379c1e64dSHong Zhang   PetscErrorCode ierr;
14b312708bSHong Zhang   PetscInt       i,n,nrows,N,j,k,m,ncols,col,nis=iscoloring->n,*rowhit,bs,bs2,*spidx,nz;
154b2e90caSHong Zhang   const PetscInt *is,*row,*ci,*cj;
1679c1e64dSHong Zhang   IS             *isa;
17525d23c0SHong Zhang   PetscBool      isBAIJ;
18b312708bSHong Zhang   PetscScalar    *A_val,**valaddrhit;
199a09c712SHong Zhang   MatEntry       *Jentry;
2079c1e64dSHong Zhang 
2179c1e64dSHong Zhang   PetscFunctionBegin;
2279c1e64dSHong Zhang   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
2352011a10SHong Zhang 
24040ebd07SHong Zhang   /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */
2579c1e64dSHong Zhang   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
26525d23c0SHong Zhang   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&isBAIJ);CHKERRQ(ierr);
27b312708bSHong Zhang   if (isBAIJ) {
28b312708bSHong Zhang     Mat_SeqBAIJ *spA = (Mat_SeqBAIJ*)mat->data;
29b312708bSHong Zhang     A_val = spA->a;
30b312708bSHong Zhang     nz    = spA->nz;
31b312708bSHong Zhang   } else {
32b312708bSHong Zhang     Mat_SeqAIJ  *spA = (Mat_SeqAIJ*)mat->data;
33b312708bSHong Zhang     A_val = spA->a;
34b312708bSHong Zhang     nz    = spA->nz;
35b312708bSHong Zhang     bs    = 1; /* only bs=1 is supported for SeqAIJ matrix */
3679c1e64dSHong Zhang   }
37b312708bSHong Zhang 
3879c1e64dSHong Zhang   N         = mat->cmap->N/bs;
3979c1e64dSHong Zhang   c->M      = mat->rmap->N/bs;   /* set total rows, columns and local rows */
4079c1e64dSHong Zhang   c->N      = mat->cmap->N/bs;
4179c1e64dSHong Zhang   c->m      = mat->rmap->N/bs;
4279c1e64dSHong Zhang   c->rstart = 0;
4379c1e64dSHong Zhang 
4479c1e64dSHong Zhang   c->ncolors = nis;
4579c1e64dSHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr);
4679c1e64dSHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr);
4779c1e64dSHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr);
484b2e90caSHong Zhang   ierr       = PetscLogObjectMemory((PetscObject)c,3*nis*sizeof(PetscInt));CHKERRQ(ierr);
499a09c712SHong Zhang 
50b312708bSHong Zhang   ierr       = PetscMalloc(nz*sizeof(MatEntry),&Jentry);CHKERRQ(ierr);
51b312708bSHong Zhang   ierr       = PetscLogObjectMemory((PetscObject)c,nz*sizeof(MatEntry));CHKERRQ(ierr);
529a09c712SHong Zhang   c->matentry = Jentry;
5379c1e64dSHong Zhang 
54525d23c0SHong Zhang   if (isBAIJ) {
55525d23c0SHong Zhang     ierr = MatGetColumnIJ_SeqBAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
56525d23c0SHong Zhang   } else {
576378f32dSHong Zhang     ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
58525d23c0SHong Zhang   }
5979c1e64dSHong Zhang 
604b2e90caSHong Zhang   ierr = PetscMalloc2(c->m,PetscInt,&rowhit,c->m,PetscScalar*,&valaddrhit);CHKERRQ(ierr);
6179c1e64dSHong Zhang   ierr = PetscMemzero(rowhit,c->m*sizeof(PetscInt));CHKERRQ(ierr);
6279c1e64dSHong Zhang 
6317a7dee1SHong Zhang   nz = 0;
6479c1e64dSHong Zhang   for (i=0; i<nis; i++) { /* loop over colors */
6579c1e64dSHong Zhang     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
6679c1e64dSHong Zhang     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
6779c1e64dSHong Zhang 
6879c1e64dSHong Zhang     c->ncolumns[i] = n;
6979c1e64dSHong Zhang     if (n) {
7079c1e64dSHong Zhang       ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr);
714b2e90caSHong Zhang       ierr = PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));CHKERRQ(ierr);
7279c1e64dSHong Zhang       ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr);
7379c1e64dSHong Zhang     } else {
7479c1e64dSHong Zhang       c->columns[i] = 0;
7579c1e64dSHong Zhang     }
7679c1e64dSHong Zhang 
7779c1e64dSHong Zhang     /* fast, crude version requires O(N*N) work */
789a09c712SHong Zhang     bs2   = bs*bs;
79476e0d0aSHong Zhang     nrows = 0;
8079c1e64dSHong Zhang     for (j=0; j<n; j++) {  /* loop over columns */
8179c1e64dSHong Zhang       col    = is[j];
824b2e90caSHong Zhang       row    = cj + ci[col];
8379c1e64dSHong Zhang       m      = ci[col+1] - ci[col];
844b2e90caSHong Zhang       nrows += m;
8579c1e64dSHong Zhang       for (k=0; k<m; k++) {  /* loop over columns marking them in rowhit */
864b2e90caSHong Zhang         rowhit[*row]       = col + 1;
874b2e90caSHong Zhang         valaddrhit[*row++] = &A_val[bs2*spidx[ci[col] + k]];
8879c1e64dSHong Zhang       }
8979c1e64dSHong Zhang     }
90476e0d0aSHong Zhang     c->nrows[i] = nrows; /* total num of rows for this color */
9179c1e64dSHong Zhang 
9279c1e64dSHong Zhang     for (j=0; j<N; j++) { /* loop over rows */
9379c1e64dSHong Zhang       if (rowhit[j]) {
949a09c712SHong Zhang         Jentry[nz].row     = j;              /* local row index */
959a09c712SHong Zhang         Jentry[nz].col     = rowhit[j] - 1;  /* local column index */
969a09c712SHong Zhang         Jentry[nz].valaddr = valaddrhit[j];  /* address of mat value for this entry */
979a09c712SHong Zhang         nz++;
9879c1e64dSHong Zhang         rowhit[j] = 0.0;                     /* zero rowhit for reuse */
9979c1e64dSHong Zhang       }
10079c1e64dSHong Zhang     }
10179c1e64dSHong Zhang     ierr = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr);
10279c1e64dSHong Zhang   }
10385740eacSHong Zhang 
104525d23c0SHong Zhang   if (isBAIJ) {
105525d23c0SHong Zhang     ierr = MatRestoreColumnIJ_SeqBAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
106*e3225e9dSHong Zhang     ierr = PetscMalloc(bs*mat->rmap->n*sizeof(PetscScalar),&c->dy);CHKERRQ(ierr);
107525d23c0SHong Zhang   } else {
1086378f32dSHong Zhang     ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
109525d23c0SHong Zhang   }
1104b2e90caSHong Zhang   ierr = PetscFree2(rowhit,valaddrhit);CHKERRQ(ierr);
11179c1e64dSHong Zhang   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
112476e0d0aSHong Zhang 
1134737c552SHong Zhang   c->ctype = IS_COLORING_GHOSTED;
11452011a10SHong Zhang   ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),mat->rmap->n,PETSC_DETERMINE,0,NULL,&c->vscale);CHKERRQ(ierr);
11579c1e64dSHong Zhang   PetscFunctionReturn(0);
11679c1e64dSHong Zhang }
117