xref: /petsc/src/mat/impls/aij/seq/fdaij.c (revision 531e53bd989298f6d4735b3cd53afa5a4ec56058)
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       bs,nis=iscoloring->n;
15   PetscBool      isBAIJ;
16 
17   PetscFunctionBegin;
18   /* set default brows and bcols for speedup inserting the dense matrix into sparse Jacobian */
19   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
20   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&isBAIJ);CHKERRQ(ierr);
21   if (isBAIJ) { /* brows and bcols will not be used */
22     c->brows = 0;
23     c->bcols = 0;
24   } else { /* seqaij matrix */
25     /* bcols is chosen s.t. dy-array takes 50% of memory space as mat */
26     Mat_SeqAIJ *spA = (Mat_SeqAIJ*)mat->data;
27     PetscReal  mem;
28     PetscInt   nz,brows,bcols,m=mat->rmap->n;
29 
30     bs    = 1; /* only bs=1 is supported for SeqAIJ matrix */
31 
32     nz    = spA->nz;
33     mem   = nz*(sizeof(PetscScalar) + sizeof(PetscInt)) + 3*m*sizeof(PetscInt);
34     bcols = (PetscInt)(0.5*mem /(m*sizeof(PetscScalar)));
35     brows = 1000/bcols;
36     if (bcols > nis) bcols = nis;
37     if (brows == 0 || brows > m) brows = m;
38     c->brows = brows;
39     c->bcols = bcols;
40   }
41 
42   c->M       = mat->rmap->N/bs;   /* set total rows, columns and local rows */
43   c->N       = mat->cmap->N/bs;
44   c->m       = mat->rmap->N/bs;
45   c->rstart  = 0;
46   c->ncolors = nis;
47   c->ctype   = IS_COLORING_GHOSTED;
48   PetscFunctionReturn(0);
49 }
50 
51 #undef __FUNCT__
52 #define __FUNCT__ "MatFDColoringSetUp_SeqXAIJ"
53 PetscErrorCode MatFDColoringSetUp_SeqXAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c)
54 {
55   PetscErrorCode ierr;
56   PetscInt       i,n,nrows,mbs=c->m,j,k,m,ncols,col,nis=iscoloring->n,*rowhit,bs,bs2,*spidx,nz;
57   const PetscInt *is,*row,*ci,*cj;
58   IS             *isa;
59   PetscBool      isBAIJ;
60   PetscScalar    *A_val,**valaddrhit;
61   MatEntry       *Jentry;
62   PetscInt       *color_start;
63   PetscInt       bcols=c->bcols;
64 
65   PetscFunctionBegin;
66   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
67 
68   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
69   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&isBAIJ);CHKERRQ(ierr);
70   if (isBAIJ) {
71     Mat_SeqBAIJ *spA = (Mat_SeqBAIJ*)mat->data;
72     A_val = spA->a;
73     nz    = spA->nz;
74   } else {
75     Mat_SeqAIJ  *spA = (Mat_SeqAIJ*)mat->data;
76     A_val = spA->a;
77     nz    = spA->nz;
78     bs    = 1; /* only bs=1 is supported for SeqAIJ matrix */
79   }
80 
81   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr);
82   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr);
83   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr);
84   ierr       = PetscLogObjectMemory((PetscObject)c,3*nis*sizeof(PetscInt));CHKERRQ(ierr);
85 
86   ierr       = PetscMalloc(nz*sizeof(MatEntry),&Jentry);CHKERRQ(ierr);
87   ierr       = PetscLogObjectMemory((PetscObject)c,nz*sizeof(MatEntry));CHKERRQ(ierr);
88   c->matentry = Jentry;
89 
90   if (isBAIJ) {
91     ierr = MatGetColumnIJ_SeqBAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
92   } else {
93     ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
94   }
95 
96   ierr = PetscMalloc3(c->m,PetscInt,&rowhit,c->m,PetscScalar*,&valaddrhit,nis+1,PetscInt,&color_start);CHKERRQ(ierr);
97   ierr = PetscMemzero(rowhit,c->m*sizeof(PetscInt));CHKERRQ(ierr);
98 
99   nz = 0;
100   for (i=0; i<nis; i++) { /* loop over colors */
101     color_start[i] = nz;
102 
103     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
104     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
105 
106     c->ncolumns[i] = n;
107     if (n) {
108       ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr);
109       ierr = PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));CHKERRQ(ierr);
110       ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr);
111     } else {
112       c->columns[i] = 0;
113     }
114 
115     /* fast, crude version requires O(N*N) work */
116     bs2   = bs*bs;
117     nrows = 0;
118     for (j=0; j<n; j++) {  /* loop over columns */
119       col    = is[j];
120       row    = cj + ci[col];
121       m      = ci[col+1] - ci[col];
122       nrows += m;
123       for (k=0; k<m; k++) {  /* loop over columns marking them in rowhit */
124         rowhit[*row]       = col + 1;
125         valaddrhit[*row++] = &A_val[bs2*spidx[ci[col] + k]];
126       }
127     }
128     c->nrows[i] = nrows; /* total num of rows for this color */
129 
130     for (j=0; j<mbs; j++) { /* loop over rows */
131       if (rowhit[j]) {
132         Jentry[nz].row     = j;              /* local row index */
133         Jentry[nz].col     = rowhit[j] - 1;  /* local column index */
134         Jentry[nz].valaddr = valaddrhit[j];  /* address of mat value for this entry */
135         nz++;
136         rowhit[j] = 0.0;                     /* zero rowhit for reuse */
137       }
138     }
139     ierr = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr);
140   }
141   color_start[nis] = nz;
142 
143   /*---------- reorder Jentry for faster MatFDColoringApply() ------------*/
144   if (!isBAIJ && bcols > 1) {
145     PetscInt nbcols=0,brows=c->brows,*row_start,*nrows_new,nz_new,row_end;
146     MatEntry *Jentry_new;
147 
148     m = mbs;
149     if (brows < 1 || brows > m) brows = m;
150 
151     ierr = PetscMalloc(nz*sizeof(MatEntry),&Jentry_new);CHKERRQ(ierr);
152     ierr = PetscMalloc(bcols*sizeof(PetscInt),&row_start);CHKERRQ(ierr);
153     ierr = PetscMalloc(nis*sizeof(PetscInt),&nrows_new);CHKERRQ(ierr);
154 
155     nz_new  = 0;
156     for (i=0; i<nis; i+=bcols) { /* loop over colors */
157       if (i + bcols > nis) bcols = nis - i;
158 
159       row_end = brows;
160       if (row_end > mbs) row_end = mbs;
161       for (j=0; j<bcols; j++) row_start[j] = 0;
162       while (row_end <= mbs) { /* loop over block rows */
163         for (j=0; j<bcols; j++) {       /* loop over block columns */
164           nrows = c->nrows[i+j];
165           nz    = color_start[i+j];
166           while (row_start[j] < nrows) {
167             if (Jentry[nz].row >= row_end) {
168               color_start[i+j] = nz;
169               break;
170             } else { /* copy Jentry[nz] to Jentry_new[nz_new] */
171               Jentry_new[nz_new].row     = Jentry[nz].row + j*mbs; /* index in dy-array */
172               Jentry_new[nz_new].col     = Jentry[nz].col;
173               Jentry_new[nz_new].valaddr = Jentry[nz].valaddr;
174               nz_new++; nz++; row_start[j]++;
175             }
176           }
177         }
178         if (row_end == mbs) break;
179         row_end += brows;
180         if (row_end > mbs) row_end = mbs;
181       }
182       nrows_new[nbcols++] = nz_new;
183     }
184     for (i=nbcols-1; i>0; i--) nrows_new[i] -= nrows_new[i-1];
185     ierr = PetscFree(c->nrows);CHKERRQ(ierr);
186     c->nrows = nrows_new;
187 
188     ierr = PetscFree(Jentry);CHKERRQ(ierr);
189     c->matentry = Jentry_new;
190     ierr = PetscMalloc(c->bcols*mat->rmap->n*sizeof(PetscScalar),&c->dy);CHKERRQ(ierr);
191     ierr = PetscFree(row_start);CHKERRQ(ierr);
192   }
193   /*---------------------------------------*/
194 
195   if (isBAIJ) {
196     ierr = MatRestoreColumnIJ_SeqBAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
197     ierr = PetscMalloc(bs*mat->rmap->n*sizeof(PetscScalar),&c->dy);CHKERRQ(ierr);
198   } else {
199     ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
200   }
201   ierr = PetscFree3(rowhit,valaddrhit,color_start);CHKERRQ(ierr);
202   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
203 
204   ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),mat->rmap->n,PETSC_DETERMINE,0,NULL,&c->vscale);CHKERRQ(ierr);
205   PetscFunctionReturn(0);
206 }
207