xref: /petsc/src/mat/impls/aij/seq/fdaij.c (revision 5bdb020c42840d115a46dbe91df9860332cbbbfd)
170c3da92SBarry Smith 
2c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h>
3525d23c0SHong Zhang #include <../src/mat/impls/baij/seq/baij.h>
4af0996ceSBarry Smith #include <petsc/private/isimpl.h>
5525d23c0SHong Zhang 
6040ebd07SHong Zhang /*
7040ebd07SHong Zhang     This routine is shared by SeqAIJ and SeqBAIJ matrices,
8040ebd07SHong Zhang     since it operators only on the nonzero structure of the elements or blocks.
9040ebd07SHong Zhang */
10b582cc96SHong Zhang #undef __FUNCT__
1193dfae19SHong Zhang #define __FUNCT__ "MatFDColoringCreate_SeqXAIJ"
1293dfae19SHong Zhang PetscErrorCode MatFDColoringCreate_SeqXAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c)
1393dfae19SHong Zhang {
1493dfae19SHong Zhang   PetscErrorCode ierr;
15a8971b87SHong Zhang   PetscInt       bs,nis=iscoloring->n,m=mat->rmap->n;
1693dfae19SHong Zhang   PetscBool      isBAIJ;
1793dfae19SHong Zhang 
1893dfae19SHong Zhang   PetscFunctionBegin;
19531e53bdSHong Zhang   /* set default brows and bcols for speedup inserting the dense matrix into sparse Jacobian */
2093dfae19SHong Zhang   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
2193dfae19SHong Zhang   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&isBAIJ);CHKERRQ(ierr);
22a8971b87SHong Zhang   if (isBAIJ) {
23a8971b87SHong Zhang     c->brows = m;
24a8971b87SHong Zhang     c->bcols = 1;
25531e53bdSHong Zhang   } else { /* seqaij matrix */
26531e53bdSHong Zhang     /* bcols is chosen s.t. dy-array takes 50% of memory space as mat */
2793dfae19SHong Zhang     Mat_SeqAIJ *spA = (Mat_SeqAIJ*)mat->data;
28531e53bdSHong Zhang     PetscReal  mem;
29a8971b87SHong Zhang     PetscInt   nz,brows,bcols;
30531e53bdSHong Zhang 
3193dfae19SHong Zhang     bs    = 1; /* only bs=1 is supported for SeqAIJ matrix */
32531e53bdSHong Zhang 
33531e53bdSHong Zhang     nz    = spA->nz;
34531e53bdSHong Zhang     mem   = nz*(sizeof(PetscScalar) + sizeof(PetscInt)) + 3*m*sizeof(PetscInt);
35531e53bdSHong Zhang     bcols = (PetscInt)(0.5*mem /(m*sizeof(PetscScalar)));
36531e53bdSHong Zhang     brows = 1000/bcols;
37531e53bdSHong Zhang     if (bcols > nis) bcols = nis;
38531e53bdSHong Zhang     if (brows == 0 || brows > m) brows = m;
39531e53bdSHong Zhang     c->brows = brows;
40531e53bdSHong Zhang     c->bcols = bcols;
4193dfae19SHong Zhang   }
42531e53bdSHong Zhang 
4393dfae19SHong Zhang   c->M       = mat->rmap->N/bs;   /* set total rows, columns and local rows */
4493dfae19SHong Zhang   c->N       = mat->cmap->N/bs;
4593dfae19SHong Zhang   c->m       = mat->rmap->N/bs;
4693dfae19SHong Zhang   c->rstart  = 0;
4793dfae19SHong Zhang   c->ncolors = nis;
48*5bdb020cSBarry Smith   c->ctype   = IS_COLORING_LOCAL;
4993dfae19SHong Zhang   PetscFunctionReturn(0);
5093dfae19SHong Zhang }
5193dfae19SHong Zhang 
52b3e1f37bSHong Zhang /*
530df34763SHong Zhang  Reorder Jentry such that blocked brows*bols of entries from dense matrix are inserted into Jacobian for improved cache performance
54b3e1f37bSHong Zhang    Input Parameters:
55b3e1f37bSHong Zhang +  mat - the matrix containing the nonzero structure of the Jacobian
56b3e1f37bSHong Zhang .  color - the coloring context
57b3e1f37bSHong Zhang -  nz - number of local non-zeros in mat
58b3e1f37bSHong Zhang */
5993dfae19SHong Zhang #undef __FUNCT__
60a8971b87SHong Zhang #define __FUNCT__ "MatFDColoringSetUpBlocked_AIJ_Private"
61a8971b87SHong Zhang PetscErrorCode MatFDColoringSetUpBlocked_AIJ_Private(Mat mat,MatFDColoring c,PetscInt nz)
6279c1e64dSHong Zhang {
6379c1e64dSHong Zhang   PetscErrorCode ierr;
64b3e1f37bSHong Zhang   PetscInt       i,j,nrows,nbcols,brows=c->brows,bcols=c->bcols,mbs=c->m,nis=c->ncolors;
65b3e1f37bSHong Zhang   PetscInt       *color_start,*row_start,*nrows_new,nz_new,row_end;
6679c1e64dSHong Zhang 
6779c1e64dSHong Zhang   PetscFunctionBegin;
68a8971b87SHong Zhang   if (brows < 1 || brows > mbs) brows = mbs;
69dcca6d9dSJed Brown   ierr = PetscMalloc2(bcols+1,&color_start,bcols,&row_start);CHKERRQ(ierr);
70*5bdb020cSBarry Smith   ierr = PetscCalloc1(nis,&nrows_new);CHKERRQ(ierr);
71785e854fSJed Brown   ierr = PetscMalloc1(bcols*mat->rmap->n,&c->dy);CHKERRQ(ierr);
72b3e1f37bSHong Zhang   ierr = PetscLogObjectMemory((PetscObject)c,bcols*mat->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
736a509798SHong Zhang 
74e2c857f8SHong Zhang   nz_new = 0;
75b3e1f37bSHong Zhang   nbcols = 0;
76054951adSHong Zhang   color_start[bcols] = 0;
770df34763SHong Zhang 
780df34763SHong Zhang   if (c->htype[0] == 'd') { /* ----  c->htype == 'ds', use MatEntry --------*/
790df34763SHong Zhang     MatEntry       *Jentry_new,*Jentry=c->matentry;
80785e854fSJed Brown     ierr = PetscMalloc1(nz,&Jentry_new);CHKERRQ(ierr);
81e2c857f8SHong Zhang     for (i=0; i<nis; i+=bcols) { /* loop over colors */
82054951adSHong Zhang       if (i + bcols > nis) {
83054951adSHong Zhang         color_start[nis - i] = color_start[bcols];
84054951adSHong Zhang         bcols                = nis - i;
85054951adSHong Zhang       }
86054951adSHong Zhang 
87054951adSHong Zhang       color_start[0] = color_start[bcols];
88054951adSHong Zhang       for (j=0; j<bcols; j++) {
89054951adSHong Zhang         color_start[j+1] = c->nrows[i+j] + color_start[j];
90054951adSHong Zhang         row_start[j]     = 0;
91054951adSHong Zhang       }
92e2c857f8SHong Zhang 
93e2c857f8SHong Zhang       row_end = brows;
94c8a9c622SHong Zhang       if (row_end > mbs) row_end = mbs;
95054951adSHong Zhang 
96c8a9c622SHong Zhang       while (row_end <= mbs) {   /* loop over block rows */
97e2c857f8SHong Zhang         for (j=0; j<bcols; j++) {       /* loop over block columns */
98e2c857f8SHong Zhang           nrows = c->nrows[i+j];
99054951adSHong Zhang           nz    = color_start[j];
100c8a9c622SHong Zhang           while (row_start[j] < nrows) {
101e2c857f8SHong Zhang             if (Jentry[nz].row >= row_end) {
102054951adSHong Zhang               color_start[j] = nz;
103e2c857f8SHong Zhang               break;
104c8a9c622SHong Zhang             } else { /* copy Jentry[nz] to Jentry_new[nz_new] */
105c8a9c622SHong Zhang               Jentry_new[nz_new].row     = Jentry[nz].row + j*mbs; /* index in dy-array */
106d880da65SHong Zhang               Jentry_new[nz_new].col     = Jentry[nz].col;
107e2c857f8SHong Zhang               Jentry_new[nz_new].valaddr = Jentry[nz].valaddr;
108c8a9c622SHong Zhang               nz_new++; nz++; row_start[j]++;
109e2c857f8SHong Zhang             }
110e2c857f8SHong Zhang           }
111e2c857f8SHong Zhang         }
112c8a9c622SHong Zhang         if (row_end == mbs) break;
113e2c857f8SHong Zhang         row_end += brows;
114c8a9c622SHong Zhang         if (row_end > mbs) row_end = mbs;
115e2c857f8SHong Zhang       }
1166a509798SHong Zhang       nrows_new[nbcols++] = nz_new;
117e2c857f8SHong Zhang     }
1180df34763SHong Zhang     ierr = PetscFree(Jentry);CHKERRQ(ierr);
1190df34763SHong Zhang     c->matentry = Jentry_new;
1200df34763SHong Zhang   } else { /* ---------  c->htype == 'wp', use MatEntry2 ------------------*/
1210df34763SHong Zhang     MatEntry2      *Jentry2_new,*Jentry2=c->matentry2;
122785e854fSJed Brown     ierr = PetscMalloc1(nz,&Jentry2_new);CHKERRQ(ierr);
1230df34763SHong Zhang     for (i=0; i<nis; i+=bcols) { /* loop over colors */
1240df34763SHong Zhang       if (i + bcols > nis) {
1250df34763SHong Zhang         color_start[nis - i] = color_start[bcols];
1260df34763SHong Zhang         bcols                = nis - i;
1270df34763SHong Zhang       }
1280df34763SHong Zhang 
1290df34763SHong Zhang       color_start[0] = color_start[bcols];
1300df34763SHong Zhang       for (j=0; j<bcols; j++) {
1310df34763SHong Zhang         color_start[j+1] = c->nrows[i+j] + color_start[j];
1320df34763SHong Zhang         row_start[j]     = 0;
1330df34763SHong Zhang       }
1340df34763SHong Zhang 
1350df34763SHong Zhang       row_end = brows;
1360df34763SHong Zhang       if (row_end > mbs) row_end = mbs;
1370df34763SHong Zhang 
1380df34763SHong Zhang       while (row_end <= mbs) {   /* loop over block rows */
1390df34763SHong Zhang         for (j=0; j<bcols; j++) {       /* loop over block columns */
1400df34763SHong Zhang           nrows = c->nrows[i+j];
1410df34763SHong Zhang           nz    = color_start[j];
1420df34763SHong Zhang           while (row_start[j] < nrows) {
1430df34763SHong Zhang             if (Jentry2[nz].row >= row_end) {
1440df34763SHong Zhang               color_start[j] = nz;
1450df34763SHong Zhang               break;
1460df34763SHong Zhang             } else { /* copy Jentry2[nz] to Jentry2_new[nz_new] */
1470df34763SHong Zhang               Jentry2_new[nz_new].row     = Jentry2[nz].row + j*mbs; /* index in dy-array */
1480df34763SHong Zhang               Jentry2_new[nz_new].valaddr = Jentry2[nz].valaddr;
1490df34763SHong Zhang               nz_new++; nz++; row_start[j]++;
1500df34763SHong Zhang             }
1510df34763SHong Zhang           }
1520df34763SHong Zhang         }
1530df34763SHong Zhang         if (row_end == mbs) break;
1540df34763SHong Zhang         row_end += brows;
1550df34763SHong Zhang         if (row_end > mbs) row_end = mbs;
1560df34763SHong Zhang       }
1570df34763SHong Zhang       nrows_new[nbcols++] = nz_new;
1580df34763SHong Zhang     }
1590df34763SHong Zhang     ierr = PetscFree(Jentry2);CHKERRQ(ierr);
1600df34763SHong Zhang     c->matentry2 = Jentry2_new;
1610df34763SHong Zhang   } /* ---------------------------------------------*/
1620df34763SHong Zhang 
163b3e1f37bSHong Zhang   ierr = PetscFree2(color_start,row_start);CHKERRQ(ierr);
164b3e1f37bSHong Zhang 
1656a509798SHong Zhang   for (i=nbcols-1; i>0; i--) nrows_new[i] -= nrows_new[i-1];
1666a509798SHong Zhang   ierr = PetscFree(c->nrows);CHKERRQ(ierr);
167b3e1f37bSHong Zhang   c->nrows = nrows_new;
168a8971b87SHong Zhang   PetscFunctionReturn(0);
169e2c857f8SHong Zhang }
170a8971b87SHong Zhang 
171a8971b87SHong Zhang #undef __FUNCT__
172a8971b87SHong Zhang #define __FUNCT__ "MatFDColoringSetUp_SeqXAIJ"
173a8971b87SHong Zhang PetscErrorCode MatFDColoringSetUp_SeqXAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c)
174a8971b87SHong Zhang {
175a8971b87SHong Zhang   PetscErrorCode ierr;
176a8971b87SHong Zhang   PetscInt       i,n,nrows,mbs=c->m,j,k,m,ncols,col,nis=iscoloring->n,*rowhit,bs,bs2,*spidx,nz;
177a8971b87SHong Zhang   const PetscInt *is,*row,*ci,*cj;
178a8971b87SHong Zhang   IS             *isa;
179a8971b87SHong Zhang   PetscBool      isBAIJ;
180a8971b87SHong Zhang   PetscScalar    *A_val,**valaddrhit;
181a8971b87SHong Zhang   MatEntry       *Jentry;
1820df34763SHong Zhang   MatEntry2      *Jentry2;
183a8971b87SHong Zhang 
184a8971b87SHong Zhang   PetscFunctionBegin;
185a8971b87SHong Zhang   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
186a8971b87SHong Zhang 
187a8971b87SHong Zhang   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
188a8971b87SHong Zhang   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&isBAIJ);CHKERRQ(ierr);
189a8971b87SHong Zhang   if (isBAIJ) {
190a8971b87SHong Zhang     Mat_SeqBAIJ *spA = (Mat_SeqBAIJ*)mat->data;
191a8971b87SHong Zhang     A_val = spA->a;
192a8971b87SHong Zhang     nz    = spA->nz;
193a8971b87SHong Zhang   } else {
194a8971b87SHong Zhang     Mat_SeqAIJ  *spA = (Mat_SeqAIJ*)mat->data;
195a8971b87SHong Zhang     A_val = spA->a;
196a8971b87SHong Zhang     nz    = spA->nz;
197a8971b87SHong Zhang     bs    = 1; /* only bs=1 is supported for SeqAIJ matrix */
198a8971b87SHong Zhang   }
199a8971b87SHong Zhang 
200785e854fSJed Brown   ierr       = PetscMalloc1(nis,&c->ncolumns);CHKERRQ(ierr);
201785e854fSJed Brown   ierr       = PetscMalloc1(nis,&c->columns);CHKERRQ(ierr);
202785e854fSJed Brown   ierr       = PetscMalloc1(nis,&c->nrows);CHKERRQ(ierr);
203a8971b87SHong Zhang   ierr       = PetscLogObjectMemory((PetscObject)c,3*nis*sizeof(PetscInt));CHKERRQ(ierr);
204a8971b87SHong Zhang 
2050df34763SHong Zhang   if (c->htype[0] == 'd') {
206785e854fSJed Brown     ierr       = PetscMalloc1(nz,&Jentry);CHKERRQ(ierr);
207a8971b87SHong Zhang     ierr       = PetscLogObjectMemory((PetscObject)c,nz*sizeof(MatEntry));CHKERRQ(ierr);
208a8971b87SHong Zhang     c->matentry = Jentry;
2090df34763SHong Zhang   } else if (c->htype[0] == 'w') {
210785e854fSJed Brown     ierr       = PetscMalloc1(nz,&Jentry2);CHKERRQ(ierr);
2110df34763SHong Zhang     ierr       = PetscLogObjectMemory((PetscObject)c,nz*sizeof(MatEntry2));CHKERRQ(ierr);
2120df34763SHong Zhang     c->matentry2 = Jentry2;
2130df34763SHong Zhang   } else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"htype is not supported");
214a8971b87SHong Zhang 
215a8971b87SHong Zhang   if (isBAIJ) {
216a8971b87SHong Zhang     ierr = MatGetColumnIJ_SeqBAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
217a8971b87SHong Zhang   } else {
218a8971b87SHong Zhang     ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
219a8971b87SHong Zhang   }
220a8971b87SHong Zhang 
221dcca6d9dSJed Brown   ierr = PetscMalloc2(c->m,&rowhit,c->m,&valaddrhit);CHKERRQ(ierr);
222a8971b87SHong Zhang   ierr = PetscMemzero(rowhit,c->m*sizeof(PetscInt));CHKERRQ(ierr);
223a8971b87SHong Zhang 
224a8971b87SHong Zhang   nz = 0;
225a8971b87SHong Zhang   for (i=0; i<nis; i++) { /* loop over colors */
226a8971b87SHong Zhang     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
227a8971b87SHong Zhang     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
228a8971b87SHong Zhang 
229a8971b87SHong Zhang     c->ncolumns[i] = n;
230a8971b87SHong Zhang     if (n) {
231785e854fSJed Brown       ierr = PetscMalloc1(n,&c->columns[i]);CHKERRQ(ierr);
232a8971b87SHong Zhang       ierr = PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));CHKERRQ(ierr);
233a8971b87SHong Zhang       ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr);
234a8971b87SHong Zhang     } else {
235a8971b87SHong Zhang       c->columns[i] = 0;
236a8971b87SHong Zhang     }
237a8971b87SHong Zhang 
238a8971b87SHong Zhang     /* fast, crude version requires O(N*N) work */
239a8971b87SHong Zhang     bs2   = bs*bs;
240a8971b87SHong Zhang     nrows = 0;
241a8971b87SHong Zhang     for (j=0; j<n; j++) {  /* loop over columns */
242a8971b87SHong Zhang       col    = is[j];
243a8971b87SHong Zhang       row    = cj + ci[col];
244a8971b87SHong Zhang       m      = ci[col+1] - ci[col];
245a8971b87SHong Zhang       nrows += m;
246a8971b87SHong Zhang       for (k=0; k<m; k++) {  /* loop over columns marking them in rowhit */
247a8971b87SHong Zhang         rowhit[*row]       = col + 1;
248a8971b87SHong Zhang         valaddrhit[*row++] = &A_val[bs2*spidx[ci[col] + k]];
249a8971b87SHong Zhang       }
250a8971b87SHong Zhang     }
251a8971b87SHong Zhang     c->nrows[i] = nrows; /* total num of rows for this color */
252a8971b87SHong Zhang 
2530df34763SHong Zhang     if (c->htype[0] == 'd') {
254a8971b87SHong Zhang       for (j=0; j<mbs; j++) { /* loop over rows */
255a8971b87SHong Zhang         if (rowhit[j]) {
256a8971b87SHong Zhang           Jentry[nz].row     = j;              /* local row index */
257a8971b87SHong Zhang           Jentry[nz].col     = rowhit[j] - 1;  /* local column index */
258a8971b87SHong Zhang           Jentry[nz].valaddr = valaddrhit[j];  /* address of mat value for this entry */
259a8971b87SHong Zhang           nz++;
260a8971b87SHong Zhang           rowhit[j] = 0.0;                     /* zero rowhit for reuse */
261a8971b87SHong Zhang         }
262a8971b87SHong Zhang       }
2630df34763SHong Zhang     }  else { /* c->htype == 'wp' */
2640df34763SHong Zhang       for (j=0; j<mbs; j++) { /* loop over rows */
2650df34763SHong Zhang         if (rowhit[j]) {
2660df34763SHong Zhang           Jentry2[nz].row     = j;              /* local row index */
2670df34763SHong Zhang           Jentry2[nz].valaddr = valaddrhit[j];  /* address of mat value for this entry */
2680df34763SHong Zhang           nz++;
2690df34763SHong Zhang           rowhit[j] = 0.0;                     /* zero rowhit for reuse */
2700df34763SHong Zhang         }
2710df34763SHong Zhang       }
2720df34763SHong Zhang     }
273a8971b87SHong Zhang     ierr = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr);
274a8971b87SHong Zhang   }
275a8971b87SHong Zhang 
276a8971b87SHong Zhang   if (c->bcols > 1) {  /* reorder Jentry for faster MatFDColoringApply() */
277a8971b87SHong Zhang     ierr = MatFDColoringSetUpBlocked_AIJ_Private(mat,c,nz);CHKERRQ(ierr);
278a8971b87SHong Zhang   }
27985740eacSHong Zhang 
280525d23c0SHong Zhang   if (isBAIJ) {
281525d23c0SHong Zhang     ierr = MatRestoreColumnIJ_SeqBAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
282785e854fSJed Brown     ierr = PetscMalloc1(bs*mat->rmap->n,&c->dy);CHKERRQ(ierr);
283b3e1f37bSHong Zhang     ierr = PetscLogObjectMemory((PetscObject)c,bs*mat->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr);
284525d23c0SHong Zhang   } else {
2856378f32dSHong Zhang     ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
286525d23c0SHong Zhang   }
287a8971b87SHong Zhang   ierr = PetscFree2(rowhit,valaddrhit);CHKERRQ(ierr);
28879c1e64dSHong Zhang   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
289476e0d0aSHong Zhang 
29052011a10SHong Zhang   ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),mat->rmap->n,PETSC_DETERMINE,0,NULL,&c->vscale);CHKERRQ(ierr);
291b3e1f37bSHong Zhang   ierr = PetscInfo3(c,"ncolors %D, brows %D and bcols %D are used.\n",c->ncolors,c->brows,c->bcols);CHKERRQ(ierr);
29279c1e64dSHong Zhang   PetscFunctionReturn(0);
29379c1e64dSHong Zhang }
294