xref: /petsc/src/mat/impls/aij/seq/fdaij.c (revision 59507afe4a127d2e4766b8e63aa06e2242b5ef04)
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,*Jentry_new;
20   PetscInt       *color_start,brows=0,bcols=1,nz_new,row_end,*row_start;
21 
22   PetscFunctionBegin;
23   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
24 
25   /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */
26   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
27   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&isBAIJ);CHKERRQ(ierr);
28   if (isBAIJ) {
29     Mat_SeqBAIJ *spA = (Mat_SeqBAIJ*)mat->data;
30     A_val = spA->a;
31     nz    = spA->nz;
32   } else {
33     Mat_SeqAIJ  *spA = (Mat_SeqAIJ*)mat->data;
34     A_val = spA->a;
35     nz    = spA->nz;
36     bs    = 1; /* only bs=1 is supported for SeqAIJ matrix */
37   }
38 
39   N         = mat->cmap->N/bs;
40   c->M      = mat->rmap->N/bs;   /* set total rows, columns and local rows */
41   c->N      = mat->cmap->N/bs;
42   c->m      = mat->rmap->N/bs;
43   c->rstart = 0;
44 
45   c->ncolors = nis;
46   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr);
47   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr);
48   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr);
49   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->nrows_new);CHKERRQ(ierr);
50   ierr       = PetscLogObjectMemory((PetscObject)c,3*nis*sizeof(PetscInt));CHKERRQ(ierr);
51 
52   ierr       = PetscMalloc(nz*sizeof(MatEntry),&Jentry);CHKERRQ(ierr);
53   ierr       = PetscLogObjectMemory((PetscObject)c,nz*sizeof(MatEntry));CHKERRQ(ierr);
54   c->matentry = Jentry;
55 
56   ierr       = PetscMalloc2(nis+1,PetscInt,&color_start,nis+1,PetscInt,&row_start);CHKERRQ(ierr);
57   ierr       = PetscMalloc(nz*sizeof(MatEntry),&Jentry_new);CHKERRQ(ierr);
58   c->matentry_new = Jentry_new;
59   for (i=0; i<nis; i++) row_start[i] = 0;
60   c->brows = brows;
61   c->bcols = bcols;
62 
63   if (isBAIJ) {
64     ierr = MatGetColumnIJ_SeqBAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
65   } else {
66     ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
67   }
68 
69   ierr = PetscMalloc2(c->m,PetscInt,&rowhit,c->m,PetscScalar*,&valaddrhit);CHKERRQ(ierr);
70   ierr = PetscMemzero(rowhit,c->m*sizeof(PetscInt));CHKERRQ(ierr);
71 
72   nz = 0;
73   for (i=0; i<nis; i++) { /* loop over colors */
74     color_start[i] = nz;
75 
76     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
77     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
78 
79     c->ncolumns[i] = n;
80     if (n) {
81       ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr);
82       ierr = PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));CHKERRQ(ierr);
83       ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr);
84     } else {
85       c->columns[i] = 0;
86     }
87 
88     /* fast, crude version requires O(N*N) work */
89     bs2   = bs*bs;
90     nrows = 0;
91     for (j=0; j<n; j++) {  /* loop over columns */
92       col    = is[j];
93       row    = cj + ci[col];
94       m      = ci[col+1] - ci[col];
95       nrows += m;
96       for (k=0; k<m; k++) {  /* loop over columns marking them in rowhit */
97         rowhit[*row]       = col + 1;
98         valaddrhit[*row++] = &A_val[bs2*spidx[ci[col] + k]];
99       }
100     }
101     c->nrows[i] = nrows; /* total num of rows for this color */
102 
103     for (j=0; j<N; j++) { /* loop over rows */
104       if (rowhit[j]) {
105         Jentry[nz].row     = j;              /* local row index */
106         Jentry[nz].col     = rowhit[j] - 1;  /* local column index */
107         Jentry[nz].valaddr = valaddrhit[j];  /* address of mat value for this entry */
108         nz++;
109         rowhit[j] = 0.0;                     /* zero rowhit for reuse */
110       }
111     }
112     ierr = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr);
113   }
114   color_start[nis] = nz;
115 
116   // ---------- reorder Jentry ------------
117   ierr = PetscOptionsInt("-brows","The number of block rows","",brows,&brows,NULL);CHKERRQ(ierr);
118   if (!isBAIJ && brows) {
119     PetscInt nbcols = 0;
120     ierr = PetscOptionsInt("-bcols","The number of block columns","",bcols,&bcols,NULL);CHKERRQ(ierr);
121     if (bcols > nis) bcols = nis;
122     c->brows = brows;
123     c->bcols = bcols;
124 
125     nz_new  = 0;
126     for (i=0; i<nis; i+=bcols) { /* loop over colors */
127       if (i + bcols > nis) bcols = nis - i;
128       printf(" color %d, bcols %d\n",i,bcols);
129 
130       row_end = brows;
131       if (row_end > mat->rmap->n) row_end = mat->rmap->n;
132       while (row_end <= mat->rmap->n) { /* loop over block rows */
133         for (j=0; j<bcols; j++) {       /* loop over block columns */
134           nrows = c->nrows[i+j];
135           for (nz=color_start[i+j]; nz<color_start[i+j+1]; nz++) { /* for each Jentry */
136             if (row_start[i+j] >= nrows) break;
137             if (Jentry[nz].row >= row_end) {
138               color_start[i+j] = nz;
139               break;
140             } else {
141               Jentry_new[nz_new].row     = Jentry[nz].row;
142               Jentry_new[nz_new].col     = j; // j-th column in bcols
143               Jentry_new[nz_new].valaddr = Jentry[nz].valaddr;
144               nz_new++;
145               row_start[i+j]++;
146             }
147           }
148         }
149         if (row_end == mat->rmap->n) break;
150         row_end += brows;
151         if (row_end > mat->rmap->n) row_end = mat->rmap->n;
152       }
153       c->nrows_new[nbcols++] = nz_new;
154     }
155     for (i=nbcols-1; i>0; i--) c->nrows_new[i] -= c->nrows_new[i-1];
156 
157     ierr = PetscMalloc(c->bcols*mat->rmap->n*sizeof(PetscScalar),&c->dy);CHKERRQ(ierr);
158   }
159   //---------------------------------------
160 
161   if (isBAIJ) {
162     ierr = MatRestoreColumnIJ_SeqBAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
163     ierr = PetscMalloc(bs*mat->rmap->n*sizeof(PetscScalar),&c->dy);CHKERRQ(ierr);
164   } else {
165     ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
166   }
167   ierr = PetscFree2(rowhit,valaddrhit);CHKERRQ(ierr);
168   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
169   ierr = PetscFree2(color_start,row_start);CHKERRQ(ierr);
170 
171   c->ctype = IS_COLORING_GHOSTED;
172   ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),mat->rmap->n,PETSC_DETERMINE,0,NULL,&c->vscale);CHKERRQ(ierr);
173   PetscFunctionReturn(0);
174 }
175