xref: /petsc/src/mat/impls/aij/seq/fdaij.c (revision 531e53bd989298f6d4735b3cd53afa5a4ec56058)
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__
1093dfae19SHong Zhang #define __FUNCT__ "MatFDColoringCreate_SeqXAIJ"
1193dfae19SHong Zhang PetscErrorCode MatFDColoringCreate_SeqXAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c)
1293dfae19SHong Zhang {
1393dfae19SHong Zhang   PetscErrorCode ierr;
14*531e53bdSHong Zhang   PetscInt       bs,nis=iscoloring->n;
1593dfae19SHong Zhang   PetscBool      isBAIJ;
1693dfae19SHong Zhang 
1793dfae19SHong Zhang   PetscFunctionBegin;
18*531e53bdSHong Zhang   /* set default brows and bcols for speedup inserting the dense matrix into sparse Jacobian */
1993dfae19SHong Zhang   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
2093dfae19SHong Zhang   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&isBAIJ);CHKERRQ(ierr);
21*531e53bdSHong Zhang   if (isBAIJ) { /* brows and bcols will not be used */
22*531e53bdSHong Zhang     c->brows = 0;
23*531e53bdSHong Zhang     c->bcols = 0;
24*531e53bdSHong Zhang   } else { /* seqaij matrix */
25*531e53bdSHong Zhang     /* bcols is chosen s.t. dy-array takes 50% of memory space as mat */
2693dfae19SHong Zhang     Mat_SeqAIJ *spA = (Mat_SeqAIJ*)mat->data;
27*531e53bdSHong Zhang     PetscReal  mem;
28*531e53bdSHong Zhang     PetscInt   nz,brows,bcols,m=mat->rmap->n;
29*531e53bdSHong Zhang 
3093dfae19SHong Zhang     bs    = 1; /* only bs=1 is supported for SeqAIJ matrix */
31*531e53bdSHong Zhang 
32*531e53bdSHong Zhang     nz    = spA->nz;
33*531e53bdSHong Zhang     mem   = nz*(sizeof(PetscScalar) + sizeof(PetscInt)) + 3*m*sizeof(PetscInt);
34*531e53bdSHong Zhang     bcols = (PetscInt)(0.5*mem /(m*sizeof(PetscScalar)));
35*531e53bdSHong Zhang     brows = 1000/bcols;
36*531e53bdSHong Zhang     if (bcols > nis) bcols = nis;
37*531e53bdSHong Zhang     if (brows == 0 || brows > m) brows = m;
38*531e53bdSHong Zhang     c->brows = brows;
39*531e53bdSHong Zhang     c->bcols = bcols;
4093dfae19SHong Zhang   }
41*531e53bdSHong Zhang 
4293dfae19SHong Zhang   c->M       = mat->rmap->N/bs;   /* set total rows, columns and local rows */
4393dfae19SHong Zhang   c->N       = mat->cmap->N/bs;
4493dfae19SHong Zhang   c->m       = mat->rmap->N/bs;
4593dfae19SHong Zhang   c->rstart  = 0;
4693dfae19SHong Zhang   c->ncolors = nis;
4793dfae19SHong Zhang   c->ctype   = IS_COLORING_GHOSTED;
4893dfae19SHong Zhang   PetscFunctionReturn(0);
4993dfae19SHong Zhang }
5093dfae19SHong Zhang 
5193dfae19SHong Zhang #undef __FUNCT__
52f86b9fbaSHong Zhang #define __FUNCT__ "MatFDColoringSetUp_SeqXAIJ"
53f86b9fbaSHong Zhang PetscErrorCode MatFDColoringSetUp_SeqXAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c)
5479c1e64dSHong Zhang {
5579c1e64dSHong Zhang   PetscErrorCode ierr;
56c8a9c622SHong Zhang   PetscInt       i,n,nrows,mbs=c->m,j,k,m,ncols,col,nis=iscoloring->n,*rowhit,bs,bs2,*spidx,nz;
574b2e90caSHong Zhang   const PetscInt *is,*row,*ci,*cj;
5879c1e64dSHong Zhang   IS             *isa;
59525d23c0SHong Zhang   PetscBool      isBAIJ;
60b312708bSHong Zhang   PetscScalar    *A_val,**valaddrhit;
61*531e53bdSHong Zhang   MatEntry       *Jentry;
62*531e53bdSHong Zhang   PetscInt       *color_start;
6393dfae19SHong Zhang   PetscInt       bcols=c->bcols;
6479c1e64dSHong Zhang 
6579c1e64dSHong Zhang   PetscFunctionBegin;
6679c1e64dSHong Zhang   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
6752011a10SHong Zhang 
6879c1e64dSHong Zhang   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
69525d23c0SHong Zhang   ierr = PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&isBAIJ);CHKERRQ(ierr);
70b312708bSHong Zhang   if (isBAIJ) {
71b312708bSHong Zhang     Mat_SeqBAIJ *spA = (Mat_SeqBAIJ*)mat->data;
72b312708bSHong Zhang     A_val = spA->a;
73b312708bSHong Zhang     nz    = spA->nz;
74b312708bSHong Zhang   } else {
75b312708bSHong Zhang     Mat_SeqAIJ  *spA = (Mat_SeqAIJ*)mat->data;
76b312708bSHong Zhang     A_val = spA->a;
77b312708bSHong Zhang     nz    = spA->nz;
78b312708bSHong Zhang     bs    = 1; /* only bs=1 is supported for SeqAIJ matrix */
7979c1e64dSHong Zhang   }
80b312708bSHong Zhang 
8179c1e64dSHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr);
8279c1e64dSHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr);
8379c1e64dSHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr);
844b2e90caSHong Zhang   ierr       = PetscLogObjectMemory((PetscObject)c,3*nis*sizeof(PetscInt));CHKERRQ(ierr);
859a09c712SHong Zhang 
86b312708bSHong Zhang   ierr       = PetscMalloc(nz*sizeof(MatEntry),&Jentry);CHKERRQ(ierr);
87b312708bSHong Zhang   ierr       = PetscLogObjectMemory((PetscObject)c,nz*sizeof(MatEntry));CHKERRQ(ierr);
889a09c712SHong Zhang   c->matentry = Jentry;
8979c1e64dSHong Zhang 
90525d23c0SHong Zhang   if (isBAIJ) {
91525d23c0SHong Zhang     ierr = MatGetColumnIJ_SeqBAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
92525d23c0SHong Zhang   } else {
936378f32dSHong Zhang     ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
94525d23c0SHong Zhang   }
9579c1e64dSHong Zhang 
9693dfae19SHong Zhang   ierr = PetscMalloc3(c->m,PetscInt,&rowhit,c->m,PetscScalar*,&valaddrhit,nis+1,PetscInt,&color_start);CHKERRQ(ierr);
9779c1e64dSHong Zhang   ierr = PetscMemzero(rowhit,c->m*sizeof(PetscInt));CHKERRQ(ierr);
9879c1e64dSHong Zhang 
9917a7dee1SHong Zhang   nz = 0;
10079c1e64dSHong Zhang   for (i=0; i<nis; i++) { /* loop over colors */
101e2c857f8SHong Zhang     color_start[i] = nz;
102e2c857f8SHong Zhang 
10379c1e64dSHong Zhang     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
10479c1e64dSHong Zhang     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
10579c1e64dSHong Zhang 
10679c1e64dSHong Zhang     c->ncolumns[i] = n;
10779c1e64dSHong Zhang     if (n) {
10879c1e64dSHong Zhang       ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr);
1094b2e90caSHong Zhang       ierr = PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));CHKERRQ(ierr);
11079c1e64dSHong Zhang       ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr);
11179c1e64dSHong Zhang     } else {
11279c1e64dSHong Zhang       c->columns[i] = 0;
11379c1e64dSHong Zhang     }
11479c1e64dSHong Zhang 
11579c1e64dSHong Zhang     /* fast, crude version requires O(N*N) work */
1169a09c712SHong Zhang     bs2   = bs*bs;
117476e0d0aSHong Zhang     nrows = 0;
11879c1e64dSHong Zhang     for (j=0; j<n; j++) {  /* loop over columns */
11979c1e64dSHong Zhang       col    = is[j];
1204b2e90caSHong Zhang       row    = cj + ci[col];
12179c1e64dSHong Zhang       m      = ci[col+1] - ci[col];
1224b2e90caSHong Zhang       nrows += m;
12379c1e64dSHong Zhang       for (k=0; k<m; k++) {  /* loop over columns marking them in rowhit */
1244b2e90caSHong Zhang         rowhit[*row]       = col + 1;
1254b2e90caSHong Zhang         valaddrhit[*row++] = &A_val[bs2*spidx[ci[col] + k]];
12679c1e64dSHong Zhang       }
12779c1e64dSHong Zhang     }
128476e0d0aSHong Zhang     c->nrows[i] = nrows; /* total num of rows for this color */
12979c1e64dSHong Zhang 
130c8a9c622SHong Zhang     for (j=0; j<mbs; j++) { /* loop over rows */
13179c1e64dSHong Zhang       if (rowhit[j]) {
1329a09c712SHong Zhang         Jentry[nz].row     = j;              /* local row index */
1339a09c712SHong Zhang         Jentry[nz].col     = rowhit[j] - 1;  /* local column index */
1349a09c712SHong Zhang         Jentry[nz].valaddr = valaddrhit[j];  /* address of mat value for this entry */
1359a09c712SHong Zhang         nz++;
13679c1e64dSHong Zhang         rowhit[j] = 0.0;                     /* zero rowhit for reuse */
13779c1e64dSHong Zhang       }
13879c1e64dSHong Zhang     }
13979c1e64dSHong Zhang     ierr = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr);
14079c1e64dSHong Zhang   }
141e2c857f8SHong Zhang   color_start[nis] = nz;
142e2c857f8SHong Zhang 
143c8a9c622SHong Zhang   /*---------- reorder Jentry for faster MatFDColoringApply() ------------*/
144f86b9fbaSHong Zhang   if (!isBAIJ && bcols > 1) {
145*531e53bdSHong Zhang     PetscInt nbcols=0,brows=c->brows,*row_start,*nrows_new,nz_new,row_end;
146*531e53bdSHong Zhang     MatEntry *Jentry_new;
14793dfae19SHong Zhang 
148c8a9c622SHong Zhang     m = mbs;
149*531e53bdSHong Zhang     if (brows < 1 || brows > m) brows = m;
1506a509798SHong Zhang 
1516a509798SHong Zhang     ierr = PetscMalloc(nz*sizeof(MatEntry),&Jentry_new);CHKERRQ(ierr);
15293dfae19SHong Zhang     ierr = PetscMalloc(bcols*sizeof(PetscInt),&row_start);CHKERRQ(ierr);
1536a509798SHong Zhang     ierr = PetscMalloc(nis*sizeof(PetscInt),&nrows_new);CHKERRQ(ierr);
1546a509798SHong Zhang 
155e2c857f8SHong Zhang     nz_new  = 0;
156e2c857f8SHong Zhang     for (i=0; i<nis; i+=bcols) { /* loop over colors */
157e2c857f8SHong Zhang       if (i + bcols > nis) bcols = nis - i;
158e2c857f8SHong Zhang 
159e2c857f8SHong Zhang       row_end = brows;
160c8a9c622SHong Zhang       if (row_end > mbs) row_end = mbs;
16193dfae19SHong Zhang       for (j=0; j<bcols; j++) row_start[j] = 0;
162c8a9c622SHong Zhang       while (row_end <= mbs) { /* loop over block rows */
163e2c857f8SHong Zhang         for (j=0; j<bcols; j++) {       /* loop over block columns */
164e2c857f8SHong Zhang           nrows = c->nrows[i+j];
165c8a9c622SHong Zhang           nz    = color_start[i+j];
166c8a9c622SHong Zhang           while (row_start[j] < nrows) {
167e2c857f8SHong Zhang             if (Jentry[nz].row >= row_end) {
168e2c857f8SHong Zhang               color_start[i+j] = nz;
169e2c857f8SHong Zhang               break;
170c8a9c622SHong Zhang             } else { /* copy Jentry[nz] to Jentry_new[nz_new] */
171c8a9c622SHong Zhang               Jentry_new[nz_new].row     = Jentry[nz].row + j*mbs; /* index in dy-array */
172d880da65SHong Zhang               Jentry_new[nz_new].col     = Jentry[nz].col;
173e2c857f8SHong Zhang               Jentry_new[nz_new].valaddr = Jentry[nz].valaddr;
174c8a9c622SHong Zhang               nz_new++; nz++; row_start[j]++;
175e2c857f8SHong Zhang             }
176e2c857f8SHong Zhang           }
177e2c857f8SHong Zhang         }
178c8a9c622SHong Zhang         if (row_end == mbs) break;
179e2c857f8SHong Zhang         row_end += brows;
180c8a9c622SHong Zhang         if (row_end > mbs) row_end = mbs;
181e2c857f8SHong Zhang       }
1826a509798SHong Zhang       nrows_new[nbcols++] = nz_new;
183e2c857f8SHong Zhang     }
1846a509798SHong Zhang     for (i=nbcols-1; i>0; i--) nrows_new[i] -= nrows_new[i-1];
1856a509798SHong Zhang     ierr = PetscFree(c->nrows);CHKERRQ(ierr);
1866a509798SHong Zhang     c->nrows = nrows_new;
187e2c857f8SHong Zhang 
1886a509798SHong Zhang     ierr = PetscFree(Jentry);CHKERRQ(ierr);
1896a509798SHong Zhang     c->matentry = Jentry_new;
190e2c857f8SHong Zhang     ierr = PetscMalloc(c->bcols*mat->rmap->n*sizeof(PetscScalar),&c->dy);CHKERRQ(ierr);
19193dfae19SHong Zhang     ierr = PetscFree(row_start);CHKERRQ(ierr);
192e2c857f8SHong Zhang   }
193c8a9c622SHong Zhang   /*---------------------------------------*/
19485740eacSHong Zhang 
195525d23c0SHong Zhang   if (isBAIJ) {
196525d23c0SHong Zhang     ierr = MatRestoreColumnIJ_SeqBAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
197e3225e9dSHong Zhang     ierr = PetscMalloc(bs*mat->rmap->n*sizeof(PetscScalar),&c->dy);CHKERRQ(ierr);
198525d23c0SHong Zhang   } else {
1996378f32dSHong Zhang     ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
200525d23c0SHong Zhang   }
20193dfae19SHong Zhang   ierr = PetscFree3(rowhit,valaddrhit,color_start);CHKERRQ(ierr);
20279c1e64dSHong Zhang   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
203476e0d0aSHong Zhang 
20452011a10SHong Zhang   ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),mat->rmap->n,PETSC_DETERMINE,0,NULL,&c->vscale);CHKERRQ(ierr);
20579c1e64dSHong Zhang   PetscFunctionReturn(0);
20679c1e64dSHong Zhang }
207