Actual source code: fdaij.c
petsc-3.5.4 2015-05-23
2: #include <../src/mat/impls/aij/seq/aij.h>
3: #include <../src/mat/impls/baij/seq/baij.h>
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: */
11: PetscErrorCode MatFDColoringCreate_SeqXAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c)
12: {
14: PetscInt bs,nis=iscoloring->n,m=mat->rmap->n;
15: PetscBool isBAIJ;
18: /* set default brows and bcols for speedup inserting the dense matrix into sparse Jacobian */
19: MatGetBlockSize(mat,&bs);
20: PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&isBAIJ);
21: if (isBAIJ) {
22: c->brows = m;
23: c->bcols = 1;
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;
30: bs = 1; /* only bs=1 is supported for SeqAIJ matrix */
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: }
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: return(0);
49: }
51: /*
52: Reorder Jentry such that blocked brows*bols of entries from dense matrix are inserted into Jacobian for improved cache performance
53: Input Parameters:
54: + mat - the matrix containing the nonzero structure of the Jacobian
55: . color - the coloring context
56: - nz - number of local non-zeros in mat
57: */
60: PetscErrorCode MatFDColoringSetUpBlocked_AIJ_Private(Mat mat,MatFDColoring c,PetscInt nz)
61: {
63: PetscInt i,j,nrows,nbcols,brows=c->brows,bcols=c->bcols,mbs=c->m,nis=c->ncolors;
64: PetscInt *color_start,*row_start,*nrows_new,nz_new,row_end;
67: if (brows < 1 || brows > mbs) brows = mbs;
68: PetscMalloc2(bcols+1,&color_start,bcols,&row_start);
69: PetscMalloc1(nis,&nrows_new);
70: PetscMalloc1(bcols*mat->rmap->n,&c->dy);
71: PetscLogObjectMemory((PetscObject)c,bcols*mat->rmap->n*sizeof(PetscScalar));
73: nz_new = 0;
74: nbcols = 0;
75: color_start[bcols] = 0;
77: if (c->htype[0] == 'd') { /* ---- c->htype == 'ds', use MatEntry --------*/
78: MatEntry *Jentry_new,*Jentry=c->matentry;
79: PetscMalloc1(nz,&Jentry_new);
80: for (i=0; i<nis; i+=bcols) { /* loop over colors */
81: if (i + bcols > nis) {
82: color_start[nis - i] = color_start[bcols];
83: bcols = nis - i;
84: }
86: color_start[0] = color_start[bcols];
87: for (j=0; j<bcols; j++) {
88: color_start[j+1] = c->nrows[i+j] + color_start[j];
89: row_start[j] = 0;
90: }
92: row_end = brows;
93: if (row_end > mbs) row_end = mbs;
95: while (row_end <= mbs) { /* loop over block rows */
96: for (j=0; j<bcols; j++) { /* loop over block columns */
97: nrows = c->nrows[i+j];
98: nz = color_start[j];
99: while (row_start[j] < nrows) {
100: if (Jentry[nz].row >= row_end) {
101: color_start[j] = nz;
102: break;
103: } else { /* copy Jentry[nz] to Jentry_new[nz_new] */
104: Jentry_new[nz_new].row = Jentry[nz].row + j*mbs; /* index in dy-array */
105: Jentry_new[nz_new].col = Jentry[nz].col;
106: Jentry_new[nz_new].valaddr = Jentry[nz].valaddr;
107: nz_new++; nz++; row_start[j]++;
108: }
109: }
110: }
111: if (row_end == mbs) break;
112: row_end += brows;
113: if (row_end > mbs) row_end = mbs;
114: }
115: nrows_new[nbcols++] = nz_new;
116: }
117: PetscFree(Jentry);
118: c->matentry = Jentry_new;
119: } else { /* --------- c->htype == 'wp', use MatEntry2 ------------------*/
120: MatEntry2 *Jentry2_new,*Jentry2=c->matentry2;
121: PetscMalloc1(nz,&Jentry2_new);
122: for (i=0; i<nis; i+=bcols) { /* loop over colors */
123: if (i + bcols > nis) {
124: color_start[nis - i] = color_start[bcols];
125: bcols = nis - i;
126: }
128: color_start[0] = color_start[bcols];
129: for (j=0; j<bcols; j++) {
130: color_start[j+1] = c->nrows[i+j] + color_start[j];
131: row_start[j] = 0;
132: }
134: row_end = brows;
135: if (row_end > mbs) row_end = mbs;
137: while (row_end <= mbs) { /* loop over block rows */
138: for (j=0; j<bcols; j++) { /* loop over block columns */
139: nrows = c->nrows[i+j];
140: nz = color_start[j];
141: while (row_start[j] < nrows) {
142: if (Jentry2[nz].row >= row_end) {
143: color_start[j] = nz;
144: break;
145: } else { /* copy Jentry2[nz] to Jentry2_new[nz_new] */
146: Jentry2_new[nz_new].row = Jentry2[nz].row + j*mbs; /* index in dy-array */
147: Jentry2_new[nz_new].valaddr = Jentry2[nz].valaddr;
148: nz_new++; nz++; row_start[j]++;
149: }
150: }
151: }
152: if (row_end == mbs) break;
153: row_end += brows;
154: if (row_end > mbs) row_end = mbs;
155: }
156: nrows_new[nbcols++] = nz_new;
157: }
158: PetscFree(Jentry2);
159: c->matentry2 = Jentry2_new;
160: } /* ---------------------------------------------*/
162: PetscFree2(color_start,row_start);
164: for (i=nbcols-1; i>0; i--) nrows_new[i] -= nrows_new[i-1];
165: PetscFree(c->nrows);
166: c->nrows = nrows_new;
167: return(0);
168: }
172: PetscErrorCode MatFDColoringSetUp_SeqXAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c)
173: {
175: PetscInt i,n,nrows,mbs=c->m,j,k,m,ncols,col,nis=iscoloring->n,*rowhit,bs,bs2,*spidx,nz;
176: const PetscInt *is,*row,*ci,*cj;
177: IS *isa;
178: PetscBool isBAIJ;
179: PetscScalar *A_val,**valaddrhit;
180: MatEntry *Jentry;
181: MatEntry2 *Jentry2;
184: ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);
186: MatGetBlockSize(mat,&bs);
187: PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&isBAIJ);
188: if (isBAIJ) {
189: Mat_SeqBAIJ *spA = (Mat_SeqBAIJ*)mat->data;
190: A_val = spA->a;
191: nz = spA->nz;
192: } else {
193: Mat_SeqAIJ *spA = (Mat_SeqAIJ*)mat->data;
194: A_val = spA->a;
195: nz = spA->nz;
196: bs = 1; /* only bs=1 is supported for SeqAIJ matrix */
197: }
199: PetscMalloc1(nis,&c->ncolumns);
200: PetscMalloc1(nis,&c->columns);
201: PetscMalloc1(nis,&c->nrows);
202: PetscLogObjectMemory((PetscObject)c,3*nis*sizeof(PetscInt));
204: if (c->htype[0] == 'd') {
205: PetscMalloc1(nz,&Jentry);
206: PetscLogObjectMemory((PetscObject)c,nz*sizeof(MatEntry));
207: c->matentry = Jentry;
208: } else if (c->htype[0] == 'w') {
209: PetscMalloc1(nz,&Jentry2);
210: PetscLogObjectMemory((PetscObject)c,nz*sizeof(MatEntry2));
211: c->matentry2 = Jentry2;
212: } else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"htype is not supported");
214: if (isBAIJ) {
215: MatGetColumnIJ_SeqBAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);
216: } else {
217: MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);
218: }
220: PetscMalloc2(c->m,&rowhit,c->m,&valaddrhit);
221: PetscMemzero(rowhit,c->m*sizeof(PetscInt));
223: nz = 0;
224: for (i=0; i<nis; i++) { /* loop over colors */
225: ISGetLocalSize(isa[i],&n);
226: ISGetIndices(isa[i],&is);
228: c->ncolumns[i] = n;
229: if (n) {
230: PetscMalloc1(n,&c->columns[i]);
231: PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));
232: PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));
233: } else {
234: c->columns[i] = 0;
235: }
237: /* fast, crude version requires O(N*N) work */
238: bs2 = bs*bs;
239: nrows = 0;
240: for (j=0; j<n; j++) { /* loop over columns */
241: col = is[j];
242: row = cj + ci[col];
243: m = ci[col+1] - ci[col];
244: nrows += m;
245: for (k=0; k<m; k++) { /* loop over columns marking them in rowhit */
246: rowhit[*row] = col + 1;
247: valaddrhit[*row++] = &A_val[bs2*spidx[ci[col] + k]];
248: }
249: }
250: c->nrows[i] = nrows; /* total num of rows for this color */
252: if (c->htype[0] == 'd') {
253: for (j=0; j<mbs; j++) { /* loop over rows */
254: if (rowhit[j]) {
255: Jentry[nz].row = j; /* local row index */
256: Jentry[nz].col = rowhit[j] - 1; /* local column index */
257: Jentry[nz].valaddr = valaddrhit[j]; /* address of mat value for this entry */
258: nz++;
259: rowhit[j] = 0.0; /* zero rowhit for reuse */
260: }
261: }
262: } else { /* c->htype == 'wp' */
263: for (j=0; j<mbs; j++) { /* loop over rows */
264: if (rowhit[j]) {
265: Jentry2[nz].row = j; /* local row index */
266: Jentry2[nz].valaddr = valaddrhit[j]; /* address of mat value for this entry */
267: nz++;
268: rowhit[j] = 0.0; /* zero rowhit for reuse */
269: }
270: }
271: }
272: ISRestoreIndices(isa[i],&is);
273: }
275: if (c->bcols > 1) { /* reorder Jentry for faster MatFDColoringApply() */
276: MatFDColoringSetUpBlocked_AIJ_Private(mat,c,nz);
277: }
279: if (isBAIJ) {
280: MatRestoreColumnIJ_SeqBAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);
281: PetscMalloc1(bs*mat->rmap->n,&c->dy);
282: PetscLogObjectMemory((PetscObject)c,bs*mat->rmap->n*sizeof(PetscScalar));
283: } else {
284: MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);
285: }
286: PetscFree2(rowhit,valaddrhit);
287: ISColoringRestoreIS(iscoloring,&isa);
289: VecCreateGhost(PetscObjectComm((PetscObject)mat),mat->rmap->n,PETSC_DETERMINE,0,NULL,&c->vscale);
290: PetscInfo3(c,"ncolors %D, brows %D and bcols %D are used.\n",c->ncolors,c->brows,c->bcols);
291: return(0);
292: }