xref: /petsc/src/mat/impls/aij/seq/fdaij.c (revision e3225e9dc71240658ff8694e906a72dc65d8b579)
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       i,n,nrows,N,j,k,m,ncols,col,nis=iscoloring->n,*rowhit,bs,bs2,*spidx,nz;
15   const PetscInt *is,*row,*ci,*cj;
16   IS             *isa;
17   PetscBool      isBAIJ;
18   PetscScalar    *A_val,**valaddrhit;
19   MatEntry       *Jentry;
20 
21   PetscFunctionBegin;
22   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
23 
24   /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */
25   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
26   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&isBAIJ);CHKERRQ(ierr);
27   if (isBAIJ) {
28     Mat_SeqBAIJ *spA = (Mat_SeqBAIJ*)mat->data;
29     A_val = spA->a;
30     nz    = spA->nz;
31   } else {
32     Mat_SeqAIJ  *spA = (Mat_SeqAIJ*)mat->data;
33     A_val = spA->a;
34     nz    = spA->nz;
35     bs    = 1; /* only bs=1 is supported for SeqAIJ matrix */
36   }
37 
38   N         = mat->cmap->N/bs;
39   c->M      = mat->rmap->N/bs;   /* set total rows, columns and local rows */
40   c->N      = mat->cmap->N/bs;
41   c->m      = mat->rmap->N/bs;
42   c->rstart = 0;
43 
44   c->ncolors = nis;
45   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr);
46   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr);
47   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr);
48   ierr       = PetscLogObjectMemory((PetscObject)c,3*nis*sizeof(PetscInt));CHKERRQ(ierr);
49 
50   ierr       = PetscMalloc(nz*sizeof(MatEntry),&Jentry);CHKERRQ(ierr);
51   ierr       = PetscLogObjectMemory((PetscObject)c,nz*sizeof(MatEntry));CHKERRQ(ierr);
52   c->matentry = Jentry;
53 
54   if (isBAIJ) {
55     ierr = MatGetColumnIJ_SeqBAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
56   } else {
57     ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
58   }
59 
60   ierr = PetscMalloc2(c->m,PetscInt,&rowhit,c->m,PetscScalar*,&valaddrhit);CHKERRQ(ierr);
61   ierr = PetscMemzero(rowhit,c->m*sizeof(PetscInt));CHKERRQ(ierr);
62 
63   nz = 0;
64   for (i=0; i<nis; i++) { /* loop over colors */
65     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
66     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
67 
68     c->ncolumns[i] = n;
69     if (n) {
70       ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr);
71       ierr = PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));CHKERRQ(ierr);
72       ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr);
73     } else {
74       c->columns[i] = 0;
75     }
76 
77     /* fast, crude version requires O(N*N) work */
78     bs2   = bs*bs;
79     nrows = 0;
80     for (j=0; j<n; j++) {  /* loop over columns */
81       col    = is[j];
82       row    = cj + ci[col];
83       m      = ci[col+1] - ci[col];
84       nrows += m;
85       for (k=0; k<m; k++) {  /* loop over columns marking them in rowhit */
86         rowhit[*row]       = col + 1;
87         valaddrhit[*row++] = &A_val[bs2*spidx[ci[col] + k]];
88       }
89     }
90     c->nrows[i] = nrows; /* total num of rows for this color */
91 
92     for (j=0; j<N; j++) { /* loop over rows */
93       if (rowhit[j]) {
94         Jentry[nz].row     = j;              /* local row index */
95         Jentry[nz].col     = rowhit[j] - 1;  /* local column index */
96         Jentry[nz].valaddr = valaddrhit[j];  /* address of mat value for this entry */
97         nz++;
98         rowhit[j] = 0.0;                     /* zero rowhit for reuse */
99       }
100     }
101     ierr = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr);
102   }
103 
104   if (isBAIJ) {
105     ierr = MatRestoreColumnIJ_SeqBAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
106     ierr = PetscMalloc(bs*mat->rmap->n*sizeof(PetscScalar),&c->dy);CHKERRQ(ierr);
107   } else {
108     ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
109   }
110   ierr = PetscFree2(rowhit,valaddrhit);CHKERRQ(ierr);
111   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
112 
113   c->ctype = IS_COLORING_GHOSTED;
114   ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),mat->rmap->n,PETSC_DETERMINE,0,NULL,&c->vscale);CHKERRQ(ierr);
115   PetscFunctionReturn(0);
116 }
117