Actual source code: aij.c
1: #define PETSCMAT_DLL
3: /*
4: Defines the basic matrix operations for the AIJ (compressed row)
5: matrix storage format.
6: */
8: #include src/mat/impls/aij/seq/aij.h
9: #include src/inline/spops.h
10: #include src/inline/dot.h
11: #include petscbt.h
15: PetscErrorCode MatDiagonalSet_SeqAIJ(Mat Y,Vec D,InsertMode is)
16: {
18: Mat_SeqAIJ *aij = (Mat_SeqAIJ*) Y->data;
19: PetscInt i,*diag, m = Y->rmap.n;
20: PetscScalar *v,*aa = aij->a;
21: PetscTruth missing;
24: if (Y->assembled) {
25: MatMissingDiagonal_SeqAIJ(Y,&missing,PETSC_NULL);
26: if (!missing) {
27: diag = aij->diag;
28: VecGetArray(D,&v);
29: if (is == INSERT_VALUES) {
30: for (i=0; i<m; i++) {
31: aa[diag[i]] = v[i];
32: }
33: } else {
34: for (i=0; i<m; i++) {
35: aa[diag[i]] += v[i];
36: }
37: }
38: VecRestoreArray(D,&v);
39: return(0);
40: }
41: }
42: MatDiagonalSet_Default(Y,D,is);
43: return(0);
44: }
48: PetscErrorCode MatGetRowIJ_SeqAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscTruth blockcompressed,PetscInt *m,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
49: {
50: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
52: PetscInt i,ishift;
53:
55: *m = A->rmap.n;
56: if (!ia) return(0);
57: ishift = 0;
58: if (symmetric && !A->structurally_symmetric) {
59: MatToSymmetricIJ_SeqAIJ(A->rmap.n,a->i,a->j,ishift,oshift,ia,ja);
60: } else if (oshift == 1) {
61: PetscInt nz = a->i[A->rmap.n];
62: /* malloc space and add 1 to i and j indices */
63: PetscMalloc((A->rmap.n+1)*sizeof(PetscInt),ia);
64: PetscMalloc((nz+1)*sizeof(PetscInt),ja);
65: for (i=0; i<nz; i++) (*ja)[i] = a->j[i] + 1;
66: for (i=0; i<A->rmap.n+1; i++) (*ia)[i] = a->i[i] + 1;
67: } else {
68: *ia = a->i; *ja = a->j;
69: }
70: return(0);
71: }
75: PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscTruth blockcompressed,PetscInt *n,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
76: {
78:
80: if (!ia) return(0);
81: if ((symmetric && !A->structurally_symmetric) || oshift == 1) {
82: PetscFree(*ia);
83: PetscFree(*ja);
84: }
85: return(0);
86: }
90: PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscTruth blockcompressed,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
91: {
92: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
94: PetscInt i,*collengths,*cia,*cja,n = A->cmap.n,m = A->rmap.n;
95: PetscInt nz = a->i[m],row,*jj,mr,col;
96:
98: *nn = n;
99: if (!ia) return(0);
100: if (symmetric) {
101: MatToSymmetricIJ_SeqAIJ(A->rmap.n,a->i,a->j,0,oshift,ia,ja);
102: } else {
103: PetscMalloc((n+1)*sizeof(PetscInt),&collengths);
104: PetscMemzero(collengths,n*sizeof(PetscInt));
105: PetscMalloc((n+1)*sizeof(PetscInt),&cia);
106: PetscMalloc((nz+1)*sizeof(PetscInt),&cja);
107: jj = a->j;
108: for (i=0; i<nz; i++) {
109: collengths[jj[i]]++;
110: }
111: cia[0] = oshift;
112: for (i=0; i<n; i++) {
113: cia[i+1] = cia[i] + collengths[i];
114: }
115: PetscMemzero(collengths,n*sizeof(PetscInt));
116: jj = a->j;
117: for (row=0; row<m; row++) {
118: mr = a->i[row+1] - a->i[row];
119: for (i=0; i<mr; i++) {
120: col = *jj++;
121: cja[cia[col] + collengths[col]++ - oshift] = row + oshift;
122: }
123: }
124: PetscFree(collengths);
125: *ia = cia; *ja = cja;
126: }
127: return(0);
128: }
132: PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscTruth blockcompressed,PetscInt *n,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
133: {
137: if (!ia) return(0);
139: PetscFree(*ia);
140: PetscFree(*ja);
141:
142: return(0);
143: }
147: PetscErrorCode MatSetValuesRow_SeqAIJ(Mat A,PetscInt row,const PetscScalar v[])
148: {
149: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
150: PetscInt *ai = a->i;
154: PetscMemcpy(a->a+ai[row],v,(ai[row+1]-ai[row])*sizeof(PetscScalar));
155: return(0);
156: }
158: #define CHUNKSIZE 15
162: PetscErrorCode MatSetValues_SeqAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
163: {
164: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
165: PetscInt *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N;
166: PetscInt *imax = a->imax,*ai = a->i,*ailen = a->ilen;
168: PetscInt *aj = a->j,nonew = a->nonew,lastcol = -1;
169: PetscScalar *ap,value,*aa = a->a;
170: PetscTruth ignorezeroentries = ((a->ignorezeroentries && is == ADD_VALUES) ? PETSC_TRUE:PETSC_FALSE);
171: PetscTruth roworiented = a->roworiented;
174: for (k=0; k<m; k++) { /* loop over added rows */
175: row = im[k];
176: if (row < 0) continue;
177: #if defined(PETSC_USE_DEBUG)
178: if (row >= A->rmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap.n-1);
179: #endif
180: rp = aj + ai[row]; ap = aa + ai[row];
181: rmax = imax[row]; nrow = ailen[row];
182: low = 0;
183: high = nrow;
184: for (l=0; l<n; l++) { /* loop over added columns */
185: if (in[l] < 0) continue;
186: #if defined(PETSC_USE_DEBUG)
187: if (in[l] >= A->cmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->cmap.n-1);
188: #endif
189: col = in[l];
190: if (roworiented) {
191: value = v[l + k*n];
192: } else {
193: value = v[k + l*m];
194: }
195: if (value == 0.0 && ignorezeroentries && (is == ADD_VALUES)) continue;
197: if (col <= lastcol) low = 0; else high = nrow;
198: lastcol = col;
199: while (high-low > 5) {
200: t = (low+high)/2;
201: if (rp[t] > col) high = t;
202: else low = t;
203: }
204: for (i=low; i<high; i++) {
205: if (rp[i] > col) break;
206: if (rp[i] == col) {
207: if (is == ADD_VALUES) ap[i] += value;
208: else ap[i] = value;
209: goto noinsert;
210: }
211: }
212: if (value == 0.0 && ignorezeroentries) goto noinsert;
213: if (nonew == 1) goto noinsert;
214: if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero at (%D,%D) in the matrix",row,col);
215: MatSeqXAIJReallocateAIJ(A,A->rmap.n,1,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar);
216: N = nrow++ - 1; a->nz++; high++;
217: /* shift up all the later entries in this row */
218: for (ii=N; ii>=i; ii--) {
219: rp[ii+1] = rp[ii];
220: ap[ii+1] = ap[ii];
221: }
222: rp[i] = col;
223: ap[i] = value;
224: noinsert:;
225: low = i + 1;
226: }
227: ailen[row] = nrow;
228: }
229: A->same_nonzero = PETSC_FALSE;
230: return(0);
231: }
236: PetscErrorCode MatGetValues_SeqAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],PetscScalar v[])
237: {
238: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
239: PetscInt *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
240: PetscInt *ai = a->i,*ailen = a->ilen;
241: PetscScalar *ap,*aa = a->a,zero = 0.0;
244: for (k=0; k<m; k++) { /* loop over rows */
245: row = im[k];
246: if (row < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",row);
247: if (row >= A->rmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap.n-1);
248: rp = aj + ai[row]; ap = aa + ai[row];
249: nrow = ailen[row];
250: for (l=0; l<n; l++) { /* loop over columns */
251: if (in[l] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %D",in[l]);
252: if (in[l] >= A->cmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->cmap.n-1);
253: col = in[l] ;
254: high = nrow; low = 0; /* assume unsorted */
255: while (high-low > 5) {
256: t = (low+high)/2;
257: if (rp[t] > col) high = t;
258: else low = t;
259: }
260: for (i=low; i<high; i++) {
261: if (rp[i] > col) break;
262: if (rp[i] == col) {
263: *v++ = ap[i];
264: goto finished;
265: }
266: }
267: *v++ = zero;
268: finished:;
269: }
270: }
271: return(0);
272: }
277: PetscErrorCode MatView_SeqAIJ_Binary(Mat A,PetscViewer viewer)
278: {
279: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
281: PetscInt i,*col_lens;
282: int fd;
285: PetscViewerBinaryGetDescriptor(viewer,&fd);
286: PetscMalloc((4+A->rmap.n)*sizeof(PetscInt),&col_lens);
287: col_lens[0] = MAT_FILE_COOKIE;
288: col_lens[1] = A->rmap.n;
289: col_lens[2] = A->cmap.n;
290: col_lens[3] = a->nz;
292: /* store lengths of each row and write (including header) to file */
293: for (i=0; i<A->rmap.n; i++) {
294: col_lens[4+i] = a->i[i+1] - a->i[i];
295: }
296: PetscBinaryWrite(fd,col_lens,4+A->rmap.n,PETSC_INT,PETSC_TRUE);
297: PetscFree(col_lens);
299: /* store column indices (zero start index) */
300: PetscBinaryWrite(fd,a->j,a->nz,PETSC_INT,PETSC_FALSE);
302: /* store nonzero values */
303: PetscBinaryWrite(fd,a->a,a->nz,PETSC_SCALAR,PETSC_FALSE);
304: return(0);
305: }
307: EXTERN PetscErrorCode MatSeqAIJFactorInfo_Matlab(Mat,PetscViewer);
311: PetscErrorCode MatView_SeqAIJ_ASCII(Mat A,PetscViewer viewer)
312: {
313: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
314: PetscErrorCode ierr;
315: PetscInt i,j,m = A->rmap.n,shift=0;
316: const char *name;
317: PetscViewerFormat format;
320: PetscObjectGetName((PetscObject)A,&name);
321: PetscViewerGetFormat(viewer,&format);
322: if (format == PETSC_VIEWER_ASCII_MATLAB) {
323: PetscInt nofinalvalue = 0;
324: if ((a->i[m] == a->i[m-1]) || (a->j[a->nz-1] != A->cmap.n-!shift)) {
325: nofinalvalue = 1;
326: }
327: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
328: PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",m,A->cmap.n);
329: PetscViewerASCIIPrintf(viewer,"%% Nonzeros = %D \n",a->nz);
330: PetscViewerASCIIPrintf(viewer,"zzz = zeros(%D,3);\n",a->nz+nofinalvalue);
331: PetscViewerASCIIPrintf(viewer,"zzz = [\n");
333: for (i=0; i<m; i++) {
334: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
335: #if defined(PETSC_USE_COMPLEX)
336: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e + %18.16ei \n",i+1,a->j[j]+!shift,PetscRealPart(a->a[j]),PetscImaginaryPart(a->a[j]));
337: #else
338: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",i+1,a->j[j]+!shift,a->a[j]);
339: #endif
340: }
341: }
342: if (nofinalvalue) {
343: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",m,A->cmap.n,0.0);
344: }
345: PetscViewerASCIIPrintf(viewer,"];\n %s = spconvert(zzz);\n",name);
346: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
347: } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) {
348: return(0);
349: } else if (format == PETSC_VIEWER_ASCII_COMMON) {
350: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
351: for (i=0; i<m; i++) {
352: PetscViewerASCIIPrintf(viewer,"row %D:",i);
353: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
354: #if defined(PETSC_USE_COMPLEX)
355: if (PetscImaginaryPart(a->a[j]) > 0.0 && PetscRealPart(a->a[j]) != 0.0) {
356: PetscViewerASCIIPrintf(viewer," (%D, %G + %G i)",a->j[j]+shift,PetscRealPart(a->a[j]),PetscImaginaryPart(a->a[j]));
357: } else if (PetscImaginaryPart(a->a[j]) < 0.0 && PetscRealPart(a->a[j]) != 0.0) {
358: PetscViewerASCIIPrintf(viewer," (%D, %G - %G i)",a->j[j]+shift,PetscRealPart(a->a[j]),-PetscImaginaryPart(a->a[j]));
359: } else if (PetscRealPart(a->a[j]) != 0.0) {
360: PetscViewerASCIIPrintf(viewer," (%D, %G) ",a->j[j]+shift,PetscRealPart(a->a[j]));
361: }
362: #else
363: if (a->a[j] != 0.0) {PetscViewerASCIIPrintf(viewer," (%D, %G) ",a->j[j]+shift,a->a[j]);}
364: #endif
365: }
366: PetscViewerASCIIPrintf(viewer,"\n");
367: }
368: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
369: } else if (format == PETSC_VIEWER_ASCII_SYMMODU) {
370: PetscInt nzd=0,fshift=1,*sptr;
371: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
372: PetscMalloc((m+1)*sizeof(PetscInt),&sptr);
373: for (i=0; i<m; i++) {
374: sptr[i] = nzd+1;
375: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
376: if (a->j[j] >= i) {
377: #if defined(PETSC_USE_COMPLEX)
378: if (PetscImaginaryPart(a->a[j]) != 0.0 || PetscRealPart(a->a[j]) != 0.0) nzd++;
379: #else
380: if (a->a[j] != 0.0) nzd++;
381: #endif
382: }
383: }
384: }
385: sptr[m] = nzd+1;
386: PetscViewerASCIIPrintf(viewer," %D %D\n\n",m,nzd);
387: for (i=0; i<m+1; i+=6) {
388: if (i+4<m) {PetscViewerASCIIPrintf(viewer," %D %D %D %D %D %D\n",sptr[i],sptr[i+1],sptr[i+2],sptr[i+3],sptr[i+4],sptr[i+5]);}
389: else if (i+3<m) {PetscViewerASCIIPrintf(viewer," %D %D %D %D %D\n",sptr[i],sptr[i+1],sptr[i+2],sptr[i+3],sptr[i+4]);}
390: else if (i+2<m) {PetscViewerASCIIPrintf(viewer," %D %D %D %D\n",sptr[i],sptr[i+1],sptr[i+2],sptr[i+3]);}
391: else if (i+1<m) {PetscViewerASCIIPrintf(viewer," %D %D %D\n",sptr[i],sptr[i+1],sptr[i+2]);}
392: else if (i<m) {PetscViewerASCIIPrintf(viewer," %D %D\n",sptr[i],sptr[i+1]);}
393: else {PetscViewerASCIIPrintf(viewer," %D\n",sptr[i]);}
394: }
395: PetscViewerASCIIPrintf(viewer,"\n");
396: PetscFree(sptr);
397: for (i=0; i<m; i++) {
398: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
399: if (a->j[j] >= i) {PetscViewerASCIIPrintf(viewer," %D ",a->j[j]+fshift);}
400: }
401: PetscViewerASCIIPrintf(viewer,"\n");
402: }
403: PetscViewerASCIIPrintf(viewer,"\n");
404: for (i=0; i<m; i++) {
405: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
406: if (a->j[j] >= i) {
407: #if defined(PETSC_USE_COMPLEX)
408: if (PetscImaginaryPart(a->a[j]) != 0.0 || PetscRealPart(a->a[j]) != 0.0) {
409: PetscViewerASCIIPrintf(viewer," %18.16e %18.16e ",PetscRealPart(a->a[j]),PetscImaginaryPart(a->a[j]));
410: }
411: #else
412: if (a->a[j] != 0.0) {PetscViewerASCIIPrintf(viewer," %18.16e ",a->a[j]);}
413: #endif
414: }
415: }
416: PetscViewerASCIIPrintf(viewer,"\n");
417: }
418: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
419: } else if (format == PETSC_VIEWER_ASCII_DENSE) {
420: PetscInt cnt = 0,jcnt;
421: PetscScalar value;
423: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
424: for (i=0; i<m; i++) {
425: jcnt = 0;
426: for (j=0; j<A->cmap.n; j++) {
427: if (jcnt < a->i[i+1]-a->i[i] && j == a->j[cnt]) {
428: value = a->a[cnt++];
429: jcnt++;
430: } else {
431: value = 0.0;
432: }
433: #if defined(PETSC_USE_COMPLEX)
434: PetscViewerASCIIPrintf(viewer," %7.5e+%7.5e i ",PetscRealPart(value),PetscImaginaryPart(value));
435: #else
436: PetscViewerASCIIPrintf(viewer," %7.5e ",value);
437: #endif
438: }
439: PetscViewerASCIIPrintf(viewer,"\n");
440: }
441: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
442: } else {
443: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
444: for (i=0; i<m; i++) {
445: PetscViewerASCIIPrintf(viewer,"row %D:",i);
446: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
447: #if defined(PETSC_USE_COMPLEX)
448: if (PetscImaginaryPart(a->a[j]) > 0.0) {
449: PetscViewerASCIIPrintf(viewer," (%D, %G + %G i)",a->j[j]+shift,PetscRealPart(a->a[j]),PetscImaginaryPart(a->a[j]));
450: } else if (PetscImaginaryPart(a->a[j]) < 0.0) {
451: PetscViewerASCIIPrintf(viewer," (%D, %G - %G i)",a->j[j]+shift,PetscRealPart(a->a[j]),-PetscImaginaryPart(a->a[j]));
452: } else {
453: PetscViewerASCIIPrintf(viewer," (%D, %G) ",a->j[j]+shift,PetscRealPart(a->a[j]));
454: }
455: #else
456: PetscViewerASCIIPrintf(viewer," (%D, %G) ",a->j[j]+shift,a->a[j]);
457: #endif
458: }
459: PetscViewerASCIIPrintf(viewer,"\n");
460: }
461: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
462: }
463: PetscViewerFlush(viewer);
464: return(0);
465: }
469: PetscErrorCode MatView_SeqAIJ_Draw_Zoom(PetscDraw draw,void *Aa)
470: {
471: Mat A = (Mat) Aa;
472: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
473: PetscErrorCode ierr;
474: PetscInt i,j,m = A->rmap.n,color;
475: PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r,maxv = 0.0;
476: PetscViewer viewer;
477: PetscViewerFormat format;
480: PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);
481: PetscViewerGetFormat(viewer,&format);
483: PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);
484: /* loop over matrix elements drawing boxes */
486: if (format != PETSC_VIEWER_DRAW_CONTOUR) {
487: /* Blue for negative, Cyan for zero and Red for positive */
488: color = PETSC_DRAW_BLUE;
489: for (i=0; i<m; i++) {
490: y_l = m - i - 1.0; y_r = y_l + 1.0;
491: for (j=a->i[i]; j<a->i[i+1]; j++) {
492: x_l = a->j[j] ; x_r = x_l + 1.0;
493: #if defined(PETSC_USE_COMPLEX)
494: if (PetscRealPart(a->a[j]) >= 0.) continue;
495: #else
496: if (a->a[j] >= 0.) continue;
497: #endif
498: PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
499: }
500: }
501: color = PETSC_DRAW_CYAN;
502: for (i=0; i<m; i++) {
503: y_l = m - i - 1.0; y_r = y_l + 1.0;
504: for (j=a->i[i]; j<a->i[i+1]; j++) {
505: x_l = a->j[j]; x_r = x_l + 1.0;
506: if (a->a[j] != 0.) continue;
507: PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
508: }
509: }
510: color = PETSC_DRAW_RED;
511: for (i=0; i<m; i++) {
512: y_l = m - i - 1.0; y_r = y_l + 1.0;
513: for (j=a->i[i]; j<a->i[i+1]; j++) {
514: x_l = a->j[j]; x_r = x_l + 1.0;
515: #if defined(PETSC_USE_COMPLEX)
516: if (PetscRealPart(a->a[j]) <= 0.) continue;
517: #else
518: if (a->a[j] <= 0.) continue;
519: #endif
520: PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
521: }
522: }
523: } else {
524: /* use contour shading to indicate magnitude of values */
525: /* first determine max of all nonzero values */
526: PetscInt nz = a->nz,count;
527: PetscDraw popup;
528: PetscReal scale;
530: for (i=0; i<nz; i++) {
531: if (PetscAbsScalar(a->a[i]) > maxv) maxv = PetscAbsScalar(a->a[i]);
532: }
533: scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv;
534: PetscDrawGetPopup(draw,&popup);
535: if (popup) {PetscDrawScalePopup(popup,0.0,maxv);}
536: count = 0;
537: for (i=0; i<m; i++) {
538: y_l = m - i - 1.0; y_r = y_l + 1.0;
539: for (j=a->i[i]; j<a->i[i+1]; j++) {
540: x_l = a->j[j]; x_r = x_l + 1.0;
541: color = PETSC_DRAW_BASIC_COLORS + (PetscInt)(scale*PetscAbsScalar(a->a[count]));
542: PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
543: count++;
544: }
545: }
546: }
547: return(0);
548: }
552: PetscErrorCode MatView_SeqAIJ_Draw(Mat A,PetscViewer viewer)
553: {
555: PetscDraw draw;
556: PetscReal xr,yr,xl,yl,h,w;
557: PetscTruth isnull;
560: PetscViewerDrawGetDraw(viewer,0,&draw);
561: PetscDrawIsNull(draw,&isnull);
562: if (isnull) return(0);
564: PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);
565: xr = A->cmap.n; yr = A->rmap.n; h = yr/10.0; w = xr/10.0;
566: xr += w; yr += h; xl = -w; yl = -h;
567: PetscDrawSetCoordinates(draw,xl,yl,xr,yr);
568: PetscDrawZoom(draw,MatView_SeqAIJ_Draw_Zoom,A);
569: PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);
570: return(0);
571: }
575: PetscErrorCode MatView_SeqAIJ(Mat A,PetscViewer viewer)
576: {
578: PetscTruth iascii,isbinary,isdraw;
581: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
582: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
583: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
584: if (iascii) {
585: MatView_SeqAIJ_ASCII(A,viewer);
586: } else if (isbinary) {
587: MatView_SeqAIJ_Binary(A,viewer);
588: } else if (isdraw) {
589: MatView_SeqAIJ_Draw(A,viewer);
590: } else {
591: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by SeqAIJ matrices",((PetscObject)viewer)->type_name);
592: }
593: MatView_Inode(A,viewer);
594: return(0);
595: }
599: PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode)
600: {
601: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
603: PetscInt fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
604: PetscInt m = A->rmap.n,*ip,N,*ailen = a->ilen,rmax = 0;
605: PetscScalar *aa = a->a,*ap;
606: PetscReal ratio=0.6;
609: if (mode == MAT_FLUSH_ASSEMBLY) return(0);
611: if (m) rmax = ailen[0]; /* determine row with most nonzeros */
612: for (i=1; i<m; i++) {
613: /* move each row back by the amount of empty slots (fshift) before it*/
614: fshift += imax[i-1] - ailen[i-1];
615: rmax = PetscMax(rmax,ailen[i]);
616: if (fshift) {
617: ip = aj + ai[i] ;
618: ap = aa + ai[i] ;
619: N = ailen[i];
620: for (j=0; j<N; j++) {
621: ip[j-fshift] = ip[j];
622: ap[j-fshift] = ap[j];
623: }
624: }
625: ai[i] = ai[i-1] + ailen[i-1];
626: }
627: if (m) {
628: fshift += imax[m-1] - ailen[m-1];
629: ai[m] = ai[m-1] + ailen[m-1];
630: }
631: /* reset ilen and imax for each row */
632: for (i=0; i<m; i++) {
633: ailen[i] = imax[i] = ai[i+1] - ai[i];
634: }
635: a->nz = ai[m];
637: MatMarkDiagonal_SeqAIJ(A);
638: PetscInfo4(A,"Matrix size: %D X %D; storage space: %D unneeded,%D used\n",m,A->cmap.n,fshift,a->nz);
639: PetscInfo1(A,"Number of mallocs during MatSetValues() is %D\n",a->reallocs);
640: PetscInfo1(A,"Maximum nonzeros in any row is %D\n",rmax);
642: a->reallocs = 0;
643: A->info.nz_unneeded = (double)fshift;
644: a->rmax = rmax;
646: /* check for zero rows. If found a large number of zero rows, use CompressedRow functions */
647: Mat_CheckCompressedRow(A,&a->compressedrow,a->i,m,ratio);
648: A->same_nonzero = PETSC_TRUE;
650: MatAssemblyEnd_Inode(A,mode);
651: return(0);
652: }
656: PetscErrorCode MatRealPart_SeqAIJ(Mat A)
657: {
658: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
659: PetscInt i,nz = a->nz;
660: PetscScalar *aa = a->a;
663: for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]);
664: return(0);
665: }
669: PetscErrorCode MatImaginaryPart_SeqAIJ(Mat A)
670: {
671: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
672: PetscInt i,nz = a->nz;
673: PetscScalar *aa = a->a;
676: for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]);
677: return(0);
678: }
682: PetscErrorCode MatZeroEntries_SeqAIJ(Mat A)
683: {
684: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
688: PetscMemzero(a->a,(a->i[A->rmap.n])*sizeof(PetscScalar));
689: return(0);
690: }
694: PetscErrorCode MatDestroy_SeqAIJ(Mat A)
695: {
696: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
700: #if defined(PETSC_USE_LOG)
701: PetscLogObjectState((PetscObject)A,"Rows=%D, Cols=%D, NZ=%D",A->rmap.n,A->cmap.n,a->nz);
702: #endif
703: MatSeqXAIJFreeAIJ(A,&a->a,&a->j,&a->i);
704: if (a->row) {
705: ISDestroy(a->row);
706: }
707: if (a->col) {
708: ISDestroy(a->col);
709: }
710: PetscFree(a->diag);
711: PetscFree2(a->imax,a->ilen);
712: PetscFree(a->idiag);
713: PetscFree(a->solve_work);
714: if (a->icol) {ISDestroy(a->icol);}
715: PetscFree(a->saved_values);
716: if (a->coloring) {ISColoringDestroy(a->coloring);}
717: PetscFree(a->xtoy);
718: if (a->XtoY) {MatDestroy(a->XtoY);}
719: if (a->compressedrow.use){PetscFree(a->compressedrow.i);}
721: MatDestroy_Inode(A);
723: PetscFree(a);
725: PetscObjectChangeTypeName((PetscObject)A,0);
726: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSeqAIJSetColumnIndices_C","",PETSC_NULL);
727: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatStoreValues_C","",PETSC_NULL);
728: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatRetrieveValues_C","",PETSC_NULL);
729: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatConvert_seqaij_seqsbaij_C","",PETSC_NULL);
730: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatConvert_seqaij_seqbaij_C","",PETSC_NULL);
731: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatConvert_seqaij_seqcsrperm_C","",PETSC_NULL);
732: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatIsTranspose_C","",PETSC_NULL);
733: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSeqAIJSetPreallocation_C","",PETSC_NULL);
734: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSeqAIJSetPreallocationCSR_C","",PETSC_NULL);
735: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatReorderForNonzeroDiagonal_C","",PETSC_NULL);
736: return(0);
737: }
741: PetscErrorCode MatCompress_SeqAIJ(Mat A)
742: {
744: return(0);
745: }
749: PetscErrorCode MatSetOption_SeqAIJ(Mat A,MatOption op)
750: {
751: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
755: switch (op) {
756: case MAT_ROW_ORIENTED:
757: a->roworiented = PETSC_TRUE;
758: break;
759: case MAT_KEEP_ZEROED_ROWS:
760: a->keepzeroedrows = PETSC_TRUE;
761: break;
762: case MAT_COLUMN_ORIENTED:
763: a->roworiented = PETSC_FALSE;
764: break;
765: case MAT_COLUMNS_SORTED:
766: a->sorted = PETSC_TRUE;
767: break;
768: case MAT_COLUMNS_UNSORTED:
769: a->sorted = PETSC_FALSE;
770: break;
771: case MAT_NO_NEW_NONZERO_LOCATIONS:
772: a->nonew = 1;
773: break;
774: case MAT_NEW_NONZERO_LOCATION_ERR:
775: a->nonew = -1;
776: break;
777: case MAT_NEW_NONZERO_ALLOCATION_ERR:
778: a->nonew = -2;
779: break;
780: case MAT_YES_NEW_NONZERO_LOCATIONS:
781: a->nonew = 0;
782: break;
783: case MAT_IGNORE_ZERO_ENTRIES:
784: a->ignorezeroentries = PETSC_TRUE;
785: break;
786: case MAT_USE_COMPRESSEDROW:
787: a->compressedrow.use = PETSC_TRUE;
788: break;
789: case MAT_DO_NOT_USE_COMPRESSEDROW:
790: a->compressedrow.use = PETSC_FALSE;
791: break;
792: case MAT_ROWS_SORTED:
793: case MAT_ROWS_UNSORTED:
794: case MAT_YES_NEW_DIAGONALS:
795: case MAT_IGNORE_OFF_PROC_ENTRIES:
796: case MAT_USE_HASH_TABLE:
797: PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
798: break;
799: case MAT_NO_NEW_DIAGONALS:
800: SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
801: default:
802: break;
803: }
804: MatSetOption_Inode(A,op);
805: return(0);
806: }
810: PetscErrorCode MatGetDiagonal_SeqAIJ(Mat A,Vec v)
811: {
812: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
814: PetscInt i,j,n;
815: PetscScalar *x,zero = 0.0;
818: VecSet(v,zero);
819: VecGetArray(v,&x);
820: VecGetLocalSize(v,&n);
821: if (n != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
822: for (i=0; i<A->rmap.n; i++) {
823: for (j=a->i[i]; j<a->i[i+1]; j++) {
824: if (a->j[j] == i) {
825: x[i] = a->a[j];
826: break;
827: }
828: }
829: }
830: VecRestoreArray(v,&x);
831: return(0);
832: }
836: PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
837: {
838: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
839: PetscScalar *x,*y;
840: PetscErrorCode ierr;
841: PetscInt m = A->rmap.n;
842: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
843: PetscScalar *v,alpha;
844: PetscInt n,i,*idx,*ii,*ridx=PETSC_NULL;
845: Mat_CompressedRow cprow = a->compressedrow;
846: PetscTruth usecprow = cprow.use;
847: #endif
850: if (zz != yy) {VecCopy(zz,yy);}
851: VecGetArray(xx,&x);
852: VecGetArray(yy,&y);
854: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
855: fortranmulttransposeaddaij_(&m,x,a->i,a->j,a->a,y);
856: #else
857: if (usecprow){
858: m = cprow.nrows;
859: ii = cprow.i;
860: ridx = cprow.rindex;
861: } else {
862: ii = a->i;
863: }
864: for (i=0; i<m; i++) {
865: idx = a->j + ii[i] ;
866: v = a->a + ii[i] ;
867: n = ii[i+1] - ii[i];
868: if (usecprow){
869: alpha = x[ridx[i]];
870: } else {
871: alpha = x[i];
872: }
873: while (n-->0) {y[*idx++] += alpha * *v++;}
874: }
875: #endif
876: PetscLogFlops(2*a->nz);
877: VecRestoreArray(xx,&x);
878: VecRestoreArray(yy,&y);
879: return(0);
880: }
884: PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec xx,Vec yy)
885: {
886: PetscScalar zero = 0.0;
890: VecSet(yy,zero);
891: MatMultTransposeAdd_SeqAIJ(A,xx,yy,yy);
892: return(0);
893: }
898: PetscErrorCode MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
899: {
900: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
901: PetscScalar *x,*y,*aa;
903: PetscInt m=A->rmap.n,*aj,*ii;
904: PetscInt n,i,j,nonzerorow=0,*ridx=PETSC_NULL;
905: PetscScalar sum;
906: PetscTruth usecprow=a->compressedrow.use;
907: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
908: PetscInt jrow;
909: #endif
911: #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
912: #pragma disjoint(*x,*y,*aa)
913: #endif
916: VecGetArray(xx,&x);
917: VecGetArray(yy,&y);
918: aj = a->j;
919: aa = a->a;
920: ii = a->i;
921: if (usecprow){ /* use compressed row format */
922: m = a->compressedrow.nrows;
923: ii = a->compressedrow.i;
924: ridx = a->compressedrow.rindex;
925: for (i=0; i<m; i++){
926: n = ii[i+1] - ii[i];
927: aj = a->j + ii[i];
928: aa = a->a + ii[i];
929: sum = 0.0;
930: nonzerorow += (n>0);
931: for (j=0; j<n; j++) sum += (*aa++)*x[*aj++];
932: y[*ridx++] = sum;
933: }
934: } else { /* do not use compressed row format */
935: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
936: fortranmultaij_(&m,x,ii,aj,aa,y);
937: #else
938: for (i=0; i<m; i++) {
939: jrow = ii[i];
940: n = ii[i+1] - jrow;
941: sum = 0.0;
942: nonzerorow += (n>0);
943: for (j=0; j<n; j++) {
944: sum += aa[jrow]*x[aj[jrow]]; jrow++;
945: }
946: y[i] = sum;
947: }
948: #endif
949: }
950: PetscLogFlops(2*a->nz - nonzerorow);
951: VecRestoreArray(xx,&x);
952: VecRestoreArray(yy,&y);
953: return(0);
954: }
958: PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
959: {
960: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
961: PetscScalar *x,*y,*z,*aa;
963: PetscInt m = A->rmap.n,*aj,*ii;
964: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
965: PetscInt n,i,jrow,j,*ridx=PETSC_NULL;
966: PetscScalar sum;
967: PetscTruth usecprow=a->compressedrow.use;
968: #endif
971: VecGetArray(xx,&x);
972: VecGetArray(yy,&y);
973: if (zz != yy) {
974: VecGetArray(zz,&z);
975: } else {
976: z = y;
977: }
979: aj = a->j;
980: aa = a->a;
981: ii = a->i;
982: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
983: fortranmultaddaij_(&m,x,ii,aj,aa,y,z);
984: #else
985: if (usecprow){ /* use compressed row format */
986: if (zz != yy){
987: PetscMemcpy(z,y,m*sizeof(PetscScalar));
988: }
989: m = a->compressedrow.nrows;
990: ii = a->compressedrow.i;
991: ridx = a->compressedrow.rindex;
992: for (i=0; i<m; i++){
993: n = ii[i+1] - ii[i];
994: aj = a->j + ii[i];
995: aa = a->a + ii[i];
996: sum = y[*ridx];
997: for (j=0; j<n; j++) sum += (*aa++)*x[*aj++];
998: z[*ridx++] = sum;
999: }
1000: } else { /* do not use compressed row format */
1001: for (i=0; i<m; i++) {
1002: jrow = ii[i];
1003: n = ii[i+1] - jrow;
1004: sum = y[i];
1005: for (j=0; j<n; j++) {
1006: sum += aa[jrow]*x[aj[jrow]]; jrow++;
1007: }
1008: z[i] = sum;
1009: }
1010: }
1011: #endif
1012: PetscLogFlops(2*a->nz);
1013: VecRestoreArray(xx,&x);
1014: VecRestoreArray(yy,&y);
1015: if (zz != yy) {
1016: VecRestoreArray(zz,&z);
1017: }
1018: return(0);
1019: }
1021: /*
1022: Adds diagonal pointers to sparse matrix structure.
1023: */
1026: PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat A)
1027: {
1028: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1030: PetscInt i,j,m = A->rmap.n;
1033: if (!a->diag) {
1034: PetscMalloc(m*sizeof(PetscInt),&a->diag);
1035: PetscLogObjectMemory(A, m*sizeof(PetscInt));
1036: }
1037: for (i=0; i<A->rmap.n; i++) {
1038: a->diag[i] = a->i[i+1];
1039: for (j=a->i[i]; j<a->i[i+1]; j++) {
1040: if (a->j[j] == i) {
1041: a->diag[i] = j;
1042: break;
1043: }
1044: }
1045: }
1046: return(0);
1047: }
1049: /*
1050: Checks for missing diagonals
1051: */
1054: PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat A,PetscTruth *missing,PetscInt *d)
1055: {
1056: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1057: PetscInt *diag,*jj = a->j,i;
1060: *missing = PETSC_FALSE;
1061: if (A->rmap.n > 0 && !jj) {
1062: *missing = PETSC_TRUE;
1063: if (d) *d = 0;
1064: PetscInfo(A,"Matrix has no entries therefor is missing diagonal");
1065: } else {
1066: diag = a->diag;
1067: for (i=0; i<A->rmap.n; i++) {
1068: if (jj[diag[i]] != i) {
1069: *missing = PETSC_TRUE;
1070: if (d) *d = i;
1071: PetscInfo1(A,"Matrix is missing diagonal number %D",i);
1072: }
1073: }
1074: }
1075: return(0);
1076: }
1080: PetscErrorCode MatRelax_SeqAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
1081: {
1082: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1083: PetscScalar *x,d,*xs,sum,*t,scale,*idiag=0,*mdiag;
1084: const PetscScalar *v = a->a, *b, *bs,*xb, *ts;
1085: PetscErrorCode ierr;
1086: PetscInt n = A->cmap.n,m = A->rmap.n,i;
1087: const PetscInt *idx,*diag;
1090: its = its*lits;
1091: if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
1093: diag = a->diag;
1094: if (!a->idiag) {
1095: PetscMalloc(3*m*sizeof(PetscScalar),&a->idiag);
1096: PetscLogObjectMemory(A, 3*m*sizeof(PetscScalar));
1097: a->ssor = a->idiag + m;
1098: mdiag = a->ssor + m;
1099: v = a->a;
1101: /* this is wrong when fshift omega changes each iteration */
1102: if (omega == 1.0 && !fshift) {
1103: for (i=0; i<m; i++) {
1104: mdiag[i] = v[diag[i]];
1105: if (!PetscAbsScalar(mdiag[i])) SETERRQ1(PETSC_ERR_ARG_INCOMP,"Zero diagonal on row %D",i);
1106: a->idiag[i] = 1.0/v[diag[i]];
1107: }
1108: PetscLogFlops(m);
1109: } else {
1110: for (i=0; i<m; i++) {
1111: mdiag[i] = v[diag[i]];
1112: a->idiag[i] = omega/(fshift + v[diag[i]]);
1113: }
1114: PetscLogFlops(2*m);
1115: }
1116: }
1117: t = a->ssor;
1118: idiag = a->idiag;
1119: mdiag = a->idiag + 2*m;
1121: VecGetArray(xx,&x);
1122: if (xx != bb) {
1123: VecGetArray(bb,(PetscScalar**)&b);
1124: } else {
1125: b = x;
1126: }
1128: /* We count flops by assuming the upper triangular and lower triangular parts have the same number of nonzeros */
1129: xs = x;
1130: if (flag == SOR_APPLY_UPPER) {
1131: /* apply (U + D/omega) to the vector */
1132: bs = b;
1133: for (i=0; i<m; i++) {
1134: d = fshift + a->a[diag[i]];
1135: n = a->i[i+1] - diag[i] - 1;
1136: idx = a->j + diag[i] + 1;
1137: v = a->a + diag[i] + 1;
1138: sum = b[i]*d/omega;
1139: SPARSEDENSEDOT(sum,bs,v,idx,n);
1140: x[i] = sum;
1141: }
1142: VecRestoreArray(xx,&x);
1143: if (bb != xx) {VecRestoreArray(bb,(PetscScalar**)&b);}
1144: PetscLogFlops(a->nz);
1145: return(0);
1146: }
1149: /* Let A = L + U + D; where L is lower trianglar,
1150: U is upper triangular, E is diagonal; This routine applies
1152: (L + E)^{-1} A (U + E)^{-1}
1154: to a vector efficiently using Eisenstat's trick. This is for
1155: the case of SSOR preconditioner, so E is D/omega where omega
1156: is the relaxation factor.
1157: */
1159: if (flag == SOR_APPLY_LOWER) {
1160: SETERRQ(PETSC_ERR_SUP,"SOR_APPLY_LOWER is not implemented");
1161: } else if (flag & SOR_EISENSTAT) {
1162: /* Let A = L + U + D; where L is lower trianglar,
1163: U is upper triangular, E is diagonal; This routine applies
1165: (L + E)^{-1} A (U + E)^{-1}
1167: to a vector efficiently using Eisenstat's trick. This is for
1168: the case of SSOR preconditioner, so E is D/omega where omega
1169: is the relaxation factor.
1170: */
1171: scale = (2.0/omega) - 1.0;
1173: /* x = (E + U)^{-1} b */
1174: for (i=m-1; i>=0; i--) {
1175: n = a->i[i+1] - diag[i] - 1;
1176: idx = a->j + diag[i] + 1;
1177: v = a->a + diag[i] + 1;
1178: sum = b[i];
1179: SPARSEDENSEMDOT(sum,xs,v,idx,n);
1180: x[i] = sum*idiag[i];
1181: }
1183: /* t = b - (2*E - D)x */
1184: v = a->a;
1185: for (i=0; i<m; i++) { t[i] = b[i] - scale*(v[*diag++])*x[i]; }
1187: /* t = (E + L)^{-1}t */
1188: ts = t;
1189: diag = a->diag;
1190: for (i=0; i<m; i++) {
1191: n = diag[i] - a->i[i];
1192: idx = a->j + a->i[i];
1193: v = a->a + a->i[i];
1194: sum = t[i];
1195: SPARSEDENSEMDOT(sum,ts,v,idx,n);
1196: t[i] = sum*idiag[i];
1197: /* x = x + t */
1198: x[i] += t[i];
1199: }
1201: PetscLogFlops(6*m-1 + 2*a->nz);
1202: VecRestoreArray(xx,&x);
1203: if (bb != xx) {VecRestoreArray(bb,(PetscScalar**)&b);}
1204: return(0);
1205: }
1206: if (flag & SOR_ZERO_INITIAL_GUESS) {
1207: if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
1208: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1209: fortranrelaxaijforwardzero_(&m,&omega,x,a->i,a->j,(PetscInt*)diag,idiag,a->a,(void*)b);
1210: #else
1211: for (i=0; i<m; i++) {
1212: n = diag[i] - a->i[i];
1213: idx = a->j + a->i[i];
1214: v = a->a + a->i[i];
1215: sum = b[i];
1216: SPARSEDENSEMDOT(sum,xs,v,idx,n);
1217: x[i] = sum*idiag[i];
1218: }
1219: #endif
1220: xb = x;
1221: PetscLogFlops(a->nz);
1222: } else xb = b;
1223: if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
1224: (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
1225: for (i=0; i<m; i++) {
1226: x[i] *= mdiag[i];
1227: }
1228: PetscLogFlops(m);
1229: }
1230: if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
1231: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1232: fortranrelaxaijbackwardzero_(&m,&omega,x,a->i,a->j,(PetscInt*)diag,idiag,a->a,(void*)xb);
1233: #else
1234: for (i=m-1; i>=0; i--) {
1235: n = a->i[i+1] - diag[i] - 1;
1236: idx = a->j + diag[i] + 1;
1237: v = a->a + diag[i] + 1;
1238: sum = xb[i];
1239: SPARSEDENSEMDOT(sum,xs,v,idx,n);
1240: x[i] = sum*idiag[i];
1241: }
1242: #endif
1243: PetscLogFlops(a->nz);
1244: }
1245: its--;
1246: }
1247: while (its--) {
1248: if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
1249: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1250: fortranrelaxaijforward_(&m,&omega,x,a->i,a->j,(PetscInt*)diag,a->a,(void*)b);
1251: #else
1252: for (i=0; i<m; i++) {
1253: n = a->i[i+1] - a->i[i];
1254: idx = a->j + a->i[i];
1255: v = a->a + a->i[i];
1256: sum = b[i];
1257: SPARSEDENSEMDOT(sum,xs,v,idx,n);
1258: x[i] = (1. - omega)*x[i] + (sum + mdiag[i]*x[i])*idiag[i];
1259: }
1260: #endif
1261: PetscLogFlops(a->nz);
1262: }
1263: if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
1264: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1265: fortranrelaxaijbackward_(&m,&omega,x,a->i,a->j,(PetscInt*)diag,a->a,(void*)b);
1266: #else
1267: for (i=m-1; i>=0; i--) {
1268: n = a->i[i+1] - a->i[i];
1269: idx = a->j + a->i[i];
1270: v = a->a + a->i[i];
1271: sum = b[i];
1272: SPARSEDENSEMDOT(sum,xs,v,idx,n);
1273: x[i] = (1. - omega)*x[i] + (sum + mdiag[i]*x[i])*idiag[i];
1274: }
1275: #endif
1276: PetscLogFlops(a->nz);
1277: }
1278: }
1279: VecRestoreArray(xx,&x);
1280: if (bb != xx) {VecRestoreArray(bb,(PetscScalar**)&b);}
1281: return(0);
1282: }
1286: PetscErrorCode MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,MatInfo *info)
1287: {
1288: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1291: info->rows_global = (double)A->rmap.n;
1292: info->columns_global = (double)A->cmap.n;
1293: info->rows_local = (double)A->rmap.n;
1294: info->columns_local = (double)A->cmap.n;
1295: info->block_size = 1.0;
1296: info->nz_allocated = (double)a->maxnz;
1297: info->nz_used = (double)a->nz;
1298: info->nz_unneeded = (double)(a->maxnz - a->nz);
1299: info->assemblies = (double)A->num_ass;
1300: info->mallocs = (double)a->reallocs;
1301: info->memory = A->mem;
1302: if (A->factor) {
1303: info->fill_ratio_given = A->info.fill_ratio_given;
1304: info->fill_ratio_needed = A->info.fill_ratio_needed;
1305: info->factor_mallocs = A->info.factor_mallocs;
1306: } else {
1307: info->fill_ratio_given = 0;
1308: info->fill_ratio_needed = 0;
1309: info->factor_mallocs = 0;
1310: }
1311: return(0);
1312: }
1316: PetscErrorCode MatZeroRows_SeqAIJ(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag)
1317: {
1318: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1319: PetscInt i,m = A->rmap.n - 1,d;
1321: PetscTruth missing;
1324: if (a->keepzeroedrows) {
1325: for (i=0; i<N; i++) {
1326: if (rows[i] < 0 || rows[i] > m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range", rows[i]);
1327: PetscMemzero(&a->a[a->i[rows[i]]],a->ilen[rows[i]]*sizeof(PetscScalar));
1328: }
1329: if (diag != 0.0) {
1330: MatMissingDiagonal_SeqAIJ(A,&missing,&d);
1331: if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry in row %D",d);
1332: for (i=0; i<N; i++) {
1333: a->a[a->diag[rows[i]]] = diag;
1334: }
1335: }
1336: A->same_nonzero = PETSC_TRUE;
1337: } else {
1338: if (diag != 0.0) {
1339: for (i=0; i<N; i++) {
1340: if (rows[i] < 0 || rows[i] > m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range", rows[i]);
1341: if (a->ilen[rows[i]] > 0) {
1342: a->ilen[rows[i]] = 1;
1343: a->a[a->i[rows[i]]] = diag;
1344: a->j[a->i[rows[i]]] = rows[i];
1345: } else { /* in case row was completely empty */
1346: MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],&diag,INSERT_VALUES);
1347: }
1348: }
1349: } else {
1350: for (i=0; i<N; i++) {
1351: if (rows[i] < 0 || rows[i] > m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range", rows[i]);
1352: a->ilen[rows[i]] = 0;
1353: }
1354: }
1355: A->same_nonzero = PETSC_FALSE;
1356: }
1357: MatAssemblyEnd_SeqAIJ(A,M