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: /* only works if matrix has a full set of diagonal entries */
16: PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalSet_SeqAIJ(Mat Y,Vec D,InsertMode is)
17: {
19: Mat_SeqAIJ *aij = (Mat_SeqAIJ*) Y->data;
20: PetscInt i,*diag, m = Y->rmap.n;
21: PetscScalar *v,*aa = aij->a;
24: MatMarkDiagonal_SeqAIJ(Y);
25: diag = aij->diag;
26: VecGetArray(D,&v);
27: if (is == INSERT_VALUES) {
28: for (i=0; i<m; i++) {
29: aa[diag[i]] = v[i];
30: }
31: } else {
32: for (i=0; i<m; i++) {
33: aa[diag[i]] += v[i];
34: }
35: }
36: VecRestoreArray(D,&v);
37: return(0);
38: }
42: PetscErrorCode MatGetRowIJ_SeqAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *m,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
43: {
44: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
46: PetscInt i,ishift;
47:
49: *m = A->rmap.n;
50: if (!ia) return(0);
51: ishift = 0;
52: if (symmetric && !A->structurally_symmetric) {
53: MatToSymmetricIJ_SeqAIJ(A->rmap.n,a->i,a->j,ishift,oshift,ia,ja);
54: } else if (oshift == 1) {
55: PetscInt nz = a->i[A->rmap.n];
56: /* malloc space and add 1 to i and j indices */
57: PetscMalloc((A->rmap.n+1)*sizeof(PetscInt),ia);
58: PetscMalloc((nz+1)*sizeof(PetscInt),ja);
59: for (i=0; i<nz; i++) (*ja)[i] = a->j[i] + 1;
60: for (i=0; i<A->rmap.n+1; i++) (*ia)[i] = a->i[i] + 1;
61: } else {
62: *ia = a->i; *ja = a->j;
63: }
64: return(0);
65: }
69: PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *n,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
70: {
72:
74: if (!ia) return(0);
75: if ((symmetric && !A->structurally_symmetric) || oshift == 1) {
76: PetscFree(*ia);
77: PetscFree(*ja);
78: }
79: return(0);
80: }
84: PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
85: {
86: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
88: PetscInt i,*collengths,*cia,*cja,n = A->cmap.n,m = A->rmap.n;
89: PetscInt nz = a->i[m],row,*jj,mr,col;
90:
92: *nn = n;
93: if (!ia) return(0);
94: if (symmetric) {
95: MatToSymmetricIJ_SeqAIJ(A->rmap.n,a->i,a->j,0,oshift,ia,ja);
96: } else {
97: PetscMalloc((n+1)*sizeof(PetscInt),&collengths);
98: PetscMemzero(collengths,n*sizeof(PetscInt));
99: PetscMalloc((n+1)*sizeof(PetscInt),&cia);
100: PetscMalloc((nz+1)*sizeof(PetscInt),&cja);
101: jj = a->j;
102: for (i=0; i<nz; i++) {
103: collengths[jj[i]]++;
104: }
105: cia[0] = oshift;
106: for (i=0; i<n; i++) {
107: cia[i+1] = cia[i] + collengths[i];
108: }
109: PetscMemzero(collengths,n*sizeof(PetscInt));
110: jj = a->j;
111: for (row=0; row<m; row++) {
112: mr = a->i[row+1] - a->i[row];
113: for (i=0; i<mr; i++) {
114: col = *jj++;
115: cja[cia[col] + collengths[col]++ - oshift] = row + oshift;
116: }
117: }
118: PetscFree(collengths);
119: *ia = cia; *ja = cja;
120: }
121: return(0);
122: }
126: PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *n,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
127: {
131: if (!ia) return(0);
133: PetscFree(*ia);
134: PetscFree(*ja);
135:
136: return(0);
137: }
141: PetscErrorCode MatSetValuesRow_SeqAIJ(Mat A,PetscInt row,const PetscScalar v[])
142: {
143: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
144: PetscInt *ai = a->i;
148: PetscMemcpy(a->a+ai[row],v,(ai[row+1]-ai[row])*sizeof(PetscScalar));
149: return(0);
150: }
152: #define CHUNKSIZE 15
156: PetscErrorCode MatSetValues_SeqAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
157: {
158: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
159: PetscInt *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N;
160: PetscInt *imax = a->imax,*ai = a->i,*ailen = a->ilen;
162: PetscInt *aj = a->j,nonew = a->nonew,lastcol = -1;
163: PetscScalar *ap,value,*aa = a->a;
164: PetscTruth ignorezeroentries = ((a->ignorezeroentries && is == ADD_VALUES) ? PETSC_TRUE:PETSC_FALSE);
165: PetscTruth roworiented = a->roworiented;
168: for (k=0; k<m; k++) { /* loop over added rows */
169: row = im[k];
170: if (row < 0) continue;
171: #if defined(PETSC_USE_DEBUG)
172: if (row >= A->rmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap.n-1);
173: #endif
174: rp = aj + ai[row]; ap = aa + ai[row];
175: rmax = imax[row]; nrow = ailen[row];
176: low = 0;
177: high = nrow;
178: for (l=0; l<n; l++) { /* loop over added columns */
179: if (in[l] < 0) continue;
180: #if defined(PETSC_USE_DEBUG)
181: if (in[l] >= A->cmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->cmap.n-1);
182: #endif
183: col = in[l];
184: if (roworiented) {
185: value = v[l + k*n];
186: } else {
187: value = v[k + l*m];
188: }
189: if (value == 0.0 && ignorezeroentries && (is == ADD_VALUES)) continue;
191: if (col <= lastcol) low = 0; else high = nrow;
192: lastcol = col;
193: while (high-low > 5) {
194: t = (low+high)/2;
195: if (rp[t] > col) high = t;
196: else low = t;
197: }
198: for (i=low; i<high; i++) {
199: if (rp[i] > col) break;
200: if (rp[i] == col) {
201: if (is == ADD_VALUES) ap[i] += value;
202: else ap[i] = value;
203: goto noinsert;
204: }
205: }
206: if (value == 0.0 && ignorezeroentries) goto noinsert;
207: if (nonew == 1) goto noinsert;
208: if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero at (%D,%D) in the matrix",row,col);
209: MatSeqXAIJReallocateAIJ(a,1,nrow,row,col,rmax,aa,ai,aj,A->rmap.n,rp,ap,imax,nonew);
210: N = nrow++ - 1; a->nz++; high++;
211: /* shift up all the later entries in this row */
212: for (ii=N; ii>=i; ii--) {
213: rp[ii+1] = rp[ii];
214: ap[ii+1] = ap[ii];
215: }
216: rp[i] = col;
217: ap[i] = value;
218: noinsert:;
219: low = i + 1;
220: }
221: ailen[row] = nrow;
222: }
223: A->same_nonzero = PETSC_FALSE;
224: return(0);
225: }
230: PetscErrorCode MatGetValues_SeqAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],PetscScalar v[])
231: {
232: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
233: PetscInt *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
234: PetscInt *ai = a->i,*ailen = a->ilen;
235: PetscScalar *ap,*aa = a->a,zero = 0.0;
238: for (k=0; k<m; k++) { /* loop over rows */
239: row = im[k];
240: if (row < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",row);
241: if (row >= A->rmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap.n-1);
242: rp = aj + ai[row]; ap = aa + ai[row];
243: nrow = ailen[row];
244: for (l=0; l<n; l++) { /* loop over columns */
245: if (in[l] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %D",in[l]);
246: if (in[l] >= A->cmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->cmap.n-1);
247: col = in[l] ;
248: high = nrow; low = 0; /* assume unsorted */
249: while (high-low > 5) {
250: t = (low+high)/2;
251: if (rp[t] > col) high = t;
252: else low = t;
253: }
254: for (i=low; i<high; i++) {
255: if (rp[i] > col) break;
256: if (rp[i] == col) {
257: *v++ = ap[i];
258: goto finished;
259: }
260: }
261: *v++ = zero;
262: finished:;
263: }
264: }
265: return(0);
266: }
271: PetscErrorCode MatView_SeqAIJ_Binary(Mat A,PetscViewer viewer)
272: {
273: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
275: PetscInt i,*col_lens;
276: int fd;
279: PetscViewerBinaryGetDescriptor(viewer,&fd);
280: PetscMalloc((4+A->rmap.n)*sizeof(PetscInt),&col_lens);
281: col_lens[0] = MAT_FILE_COOKIE;
282: col_lens[1] = A->rmap.n;
283: col_lens[2] = A->cmap.n;
284: col_lens[3] = a->nz;
286: /* store lengths of each row and write (including header) to file */
287: for (i=0; i<A->rmap.n; i++) {
288: col_lens[4+i] = a->i[i+1] - a->i[i];
289: }
290: PetscBinaryWrite(fd,col_lens,4+A->rmap.n,PETSC_INT,PETSC_TRUE);
291: PetscFree(col_lens);
293: /* store column indices (zero start index) */
294: PetscBinaryWrite(fd,a->j,a->nz,PETSC_INT,PETSC_FALSE);
296: /* store nonzero values */
297: PetscBinaryWrite(fd,a->a,a->nz,PETSC_SCALAR,PETSC_FALSE);
298: return(0);
299: }
301: EXTERN PetscErrorCode MatSeqAIJFactorInfo_Matlab(Mat,PetscViewer);
305: PetscErrorCode MatView_SeqAIJ_ASCII(Mat A,PetscViewer viewer)
306: {
307: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
308: PetscErrorCode ierr;
309: PetscInt i,j,m = A->rmap.n,shift=0;
310: const char *name;
311: PetscViewerFormat format;
314: PetscObjectGetName((PetscObject)A,&name);
315: PetscViewerGetFormat(viewer,&format);
316: if (format == PETSC_VIEWER_ASCII_MATLAB) {
317: PetscInt nofinalvalue = 0;
318: if ((a->i[m] == a->i[m-1]) || (a->j[a->nz-1] != A->cmap.n-!shift)) {
319: nofinalvalue = 1;
320: }
321: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
322: PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",m,A->cmap.n);
323: PetscViewerASCIIPrintf(viewer,"%% Nonzeros = %D \n",a->nz);
324: PetscViewerASCIIPrintf(viewer,"zzz = zeros(%D,3);\n",a->nz+nofinalvalue);
325: PetscViewerASCIIPrintf(viewer,"zzz = [\n");
327: for (i=0; i<m; i++) {
328: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
329: #if defined(PETSC_USE_COMPLEX)
330: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e + %18.16ei \n",i+1,a->j[j]+!shift,PetscRealPart(a->a[j]),PetscImaginaryPart(a->a[j]));
331: #else
332: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",i+1,a->j[j]+!shift,a->a[j]);
333: #endif
334: }
335: }
336: if (nofinalvalue) {
337: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",m,A->cmap.n,0.0);
338: }
339: PetscViewerASCIIPrintf(viewer,"];\n %s = spconvert(zzz);\n",name);
340: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
341: } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) {
342: return(0);
343: } else if (format == PETSC_VIEWER_ASCII_COMMON) {
344: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
345: for (i=0; i<m; i++) {
346: PetscViewerASCIIPrintf(viewer,"row %D:",i);
347: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
348: #if defined(PETSC_USE_COMPLEX)
349: if (PetscImaginaryPart(a->a[j]) > 0.0 && PetscRealPart(a->a[j]) != 0.0) {
350: PetscViewerASCIIPrintf(viewer," (%D, %G + %G i)",a->j[j]+shift,PetscRealPart(a->a[j]),PetscImaginaryPart(a->a[j]));
351: } else if (PetscImaginaryPart(a->a[j]) < 0.0 && PetscRealPart(a->a[j]) != 0.0) {
352: PetscViewerASCIIPrintf(viewer," (%D, %G - %G i)",a->j[j]+shift,PetscRealPart(a->a[j]),-PetscImaginaryPart(a->a[j]));
353: } else if (PetscRealPart(a->a[j]) != 0.0) {
354: PetscViewerASCIIPrintf(viewer," (%D, %G) ",a->j[j]+shift,PetscRealPart(a->a[j]));
355: }
356: #else
357: if (a->a[j] != 0.0) {PetscViewerASCIIPrintf(viewer," (%D, %G) ",a->j[j]+shift,a->a[j]);}
358: #endif
359: }
360: PetscViewerASCIIPrintf(viewer,"\n");
361: }
362: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
363: } else if (format == PETSC_VIEWER_ASCII_SYMMODU) {
364: PetscInt nzd=0,fshift=1,*sptr;
365: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
366: PetscMalloc((m+1)*sizeof(PetscInt),&sptr);
367: for (i=0; i<m; i++) {
368: sptr[i] = nzd+1;
369: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
370: if (a->j[j] >= i) {
371: #if defined(PETSC_USE_COMPLEX)
372: if (PetscImaginaryPart(a->a[j]) != 0.0 || PetscRealPart(a->a[j]) != 0.0) nzd++;
373: #else
374: if (a->a[j] != 0.0) nzd++;
375: #endif
376: }
377: }
378: }
379: sptr[m] = nzd+1;
380: PetscViewerASCIIPrintf(viewer," %D %D\n\n",m,nzd);
381: for (i=0; i<m+1; i+=6) {
382: 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]);}
383: 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]);}
384: else if (i+2<m) {PetscViewerASCIIPrintf(viewer," %D %D %D %D\n",sptr[i],sptr[i+1],sptr[i+2],sptr[i+3]);}
385: else if (i+1<m) {PetscViewerASCIIPrintf(viewer," %D %D %D\n",sptr[i],sptr[i+1],sptr[i+2]);}
386: else if (i<m) {PetscViewerASCIIPrintf(viewer," %D %D\n",sptr[i],sptr[i+1]);}
387: else {PetscViewerASCIIPrintf(viewer," %D\n",sptr[i]);}
388: }
389: PetscViewerASCIIPrintf(viewer,"\n");
390: PetscFree(sptr);
391: for (i=0; i<m; i++) {
392: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
393: if (a->j[j] >= i) {PetscViewerASCIIPrintf(viewer," %D ",a->j[j]+fshift);}
394: }
395: PetscViewerASCIIPrintf(viewer,"\n");
396: }
397: PetscViewerASCIIPrintf(viewer,"\n");
398: for (i=0; i<m; i++) {
399: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
400: if (a->j[j] >= i) {
401: #if defined(PETSC_USE_COMPLEX)
402: if (PetscImaginaryPart(a->a[j]) != 0.0 || PetscRealPart(a->a[j]) != 0.0) {
403: PetscViewerASCIIPrintf(viewer," %18.16e %18.16e ",PetscRealPart(a->a[j]),PetscImaginaryPart(a->a[j]));
404: }
405: #else
406: if (a->a[j] != 0.0) {PetscViewerASCIIPrintf(viewer," %18.16e ",a->a[j]);}
407: #endif
408: }
409: }
410: PetscViewerASCIIPrintf(viewer,"\n");
411: }
412: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
413: } else if (format == PETSC_VIEWER_ASCII_DENSE) {
414: PetscInt cnt = 0,jcnt;
415: PetscScalar value;
417: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
418: for (i=0; i<m; i++) {
419: jcnt = 0;
420: for (j=0; j<A->cmap.n; j++) {
421: if (jcnt < a->i[i+1]-a->i[i] && j == a->j[cnt]) {
422: value = a->a[cnt++];
423: jcnt++;
424: } else {
425: value = 0.0;
426: }
427: #if defined(PETSC_USE_COMPLEX)
428: PetscViewerASCIIPrintf(viewer," %7.5e+%7.5e i ",PetscRealPart(value),PetscImaginaryPart(value));
429: #else
430: PetscViewerASCIIPrintf(viewer," %7.5e ",value);
431: #endif
432: }
433: PetscViewerASCIIPrintf(viewer,"\n");
434: }
435: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
436: } else {
437: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
438: for (i=0; i<m; i++) {
439: PetscViewerASCIIPrintf(viewer,"row %D:",i);
440: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
441: #if defined(PETSC_USE_COMPLEX)
442: if (PetscImaginaryPart(a->a[j]) > 0.0) {
443: PetscViewerASCIIPrintf(viewer," (%D, %G + %G i)",a->j[j]+shift,PetscRealPart(a->a[j]),PetscImaginaryPart(a->a[j]));
444: } else if (PetscImaginaryPart(a->a[j]) < 0.0) {
445: PetscViewerASCIIPrintf(viewer," (%D, %G - %G i)",a->j[j]+shift,PetscRealPart(a->a[j]),-PetscImaginaryPart(a->a[j]));
446: } else {
447: PetscViewerASCIIPrintf(viewer," (%D, %G) ",a->j[j]+shift,PetscRealPart(a->a[j]));
448: }
449: #else
450: PetscViewerASCIIPrintf(viewer," (%D, %G) ",a->j[j]+shift,a->a[j]);
451: #endif
452: }
453: PetscViewerASCIIPrintf(viewer,"\n");
454: }
455: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
456: }
457: PetscViewerFlush(viewer);
458: return(0);
459: }
463: PetscErrorCode MatView_SeqAIJ_Draw_Zoom(PetscDraw draw,void *Aa)
464: {
465: Mat A = (Mat) Aa;
466: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
467: PetscErrorCode ierr;
468: PetscInt i,j,m = A->rmap.n,color;
469: PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r,maxv = 0.0;
470: PetscViewer viewer;
471: PetscViewerFormat format;
474: PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);
475: PetscViewerGetFormat(viewer,&format);
477: PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);
478: /* loop over matrix elements drawing boxes */
480: if (format != PETSC_VIEWER_DRAW_CONTOUR) {
481: /* Blue for negative, Cyan for zero and Red for positive */
482: color = PETSC_DRAW_BLUE;
483: for (i=0; i<m; i++) {
484: y_l = m - i - 1.0; y_r = y_l + 1.0;
485: for (j=a->i[i]; j<a->i[i+1]; j++) {
486: x_l = a->j[j] ; x_r = x_l + 1.0;
487: #if defined(PETSC_USE_COMPLEX)
488: if (PetscRealPart(a->a[j]) >= 0.) continue;
489: #else
490: if (a->a[j] >= 0.) continue;
491: #endif
492: PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
493: }
494: }
495: color = PETSC_DRAW_CYAN;
496: for (i=0; i<m; i++) {
497: y_l = m - i - 1.0; y_r = y_l + 1.0;
498: for (j=a->i[i]; j<a->i[i+1]; j++) {
499: x_l = a->j[j]; x_r = x_l + 1.0;
500: if (a->a[j] != 0.) continue;
501: PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
502: }
503: }
504: color = PETSC_DRAW_RED;
505: for (i=0; i<m; i++) {
506: y_l = m - i - 1.0; y_r = y_l + 1.0;
507: for (j=a->i[i]; j<a->i[i+1]; j++) {
508: x_l = a->j[j]; x_r = x_l + 1.0;
509: #if defined(PETSC_USE_COMPLEX)
510: if (PetscRealPart(a->a[j]) <= 0.) continue;
511: #else
512: if (a->a[j] <= 0.) continue;
513: #endif
514: PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
515: }
516: }
517: } else {
518: /* use contour shading to indicate magnitude of values */
519: /* first determine max of all nonzero values */
520: PetscInt nz = a->nz,count;
521: PetscDraw popup;
522: PetscReal scale;
524: for (i=0; i<nz; i++) {
525: if (PetscAbsScalar(a->a[i]) > maxv) maxv = PetscAbsScalar(a->a[i]);
526: }
527: scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv;
528: PetscDrawGetPopup(draw,&popup);
529: if (popup) {PetscDrawScalePopup(popup,0.0,maxv);}
530: count = 0;
531: for (i=0; i<m; i++) {
532: y_l = m - i - 1.0; y_r = y_l + 1.0;
533: for (j=a->i[i]; j<a->i[i+1]; j++) {
534: x_l = a->j[j]; x_r = x_l + 1.0;
535: color = PETSC_DRAW_BASIC_COLORS + (PetscInt)(scale*PetscAbsScalar(a->a[count]));
536: PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
537: count++;
538: }
539: }
540: }
541: return(0);
542: }
546: PetscErrorCode MatView_SeqAIJ_Draw(Mat A,PetscViewer viewer)
547: {
549: PetscDraw draw;
550: PetscReal xr,yr,xl,yl,h,w;
551: PetscTruth isnull;
554: PetscViewerDrawGetDraw(viewer,0,&draw);
555: PetscDrawIsNull(draw,&isnull);
556: if (isnull) return(0);
558: PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);
559: xr = A->cmap.n; yr = A->rmap.n; h = yr/10.0; w = xr/10.0;
560: xr += w; yr += h; xl = -w; yl = -h;
561: PetscDrawSetCoordinates(draw,xl,yl,xr,yr);
562: PetscDrawZoom(draw,MatView_SeqAIJ_Draw_Zoom,A);
563: PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);
564: return(0);
565: }
569: PetscErrorCode MatView_SeqAIJ(Mat A,PetscViewer viewer)
570: {
572: PetscTruth issocket,iascii,isbinary,isdraw;
575: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);
576: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
577: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
578: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
579: if (iascii) {
580: MatView_SeqAIJ_ASCII(A,viewer);
581: #if defined(PETSC_USE_SOCKET_VIEWER)
582: } else if (issocket) {
583: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
584: PetscViewerSocketPutSparse_Private(viewer,A->rmap.n,A->cmap.n,a->nz,a->a,a->i,a->j);
585: #endif
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: /* diagonals may have moved, so kill the diagonal pointers */
638: if (fshift && a->diag) {
639: PetscFree(a->diag);
640: PetscLogObjectMemory(A,-(m+1)*sizeof(PetscInt));
641: a->diag = 0;
642: }
643: PetscInfo4(A,"Matrix size: %D X %D; storage space: %D unneeded,%D used\n",m,A->cmap.n,fshift,a->nz);
644: PetscInfo1(A,"Number of mallocs during MatSetValues() is %D\n",a->reallocs);
645: PetscInfo1(A,"Maximum nonzeros in any row is %D\n",rmax);
646: a->reallocs = 0;
647: A->info.nz_unneeded = (double)fshift;
648: a->rmax = rmax;
650: /* check for zero rows. If found a large number of zero rows, use CompressedRow functions */
651: Mat_CheckCompressedRow(A,&a->compressedrow,a->i,m,ratio);
652: A->same_nonzero = PETSC_TRUE;
654: MatAssemblyEnd_Inode(A,mode);
655: return(0);
656: }
660: PetscErrorCode MatRealPart_SeqAIJ(Mat A)
661: {
662: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
663: PetscInt i,nz = a->nz;
664: PetscScalar *aa = a->a;
667: for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]);
668: return(0);
669: }
673: PetscErrorCode MatImaginaryPart_SeqAIJ(Mat A)
674: {
675: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
676: PetscInt i,nz = a->nz;
677: PetscScalar *aa = a->a;
680: for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]);
681: return(0);
682: }
686: PetscErrorCode MatZeroEntries_SeqAIJ(Mat A)
687: {
688: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
692: PetscMemzero(a->a,(a->i[A->rmap.n])*sizeof(PetscScalar));
693: return(0);
694: }
698: PetscErrorCode MatDestroy_SeqAIJ(Mat A)
699: {
700: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
704: #if defined(PETSC_USE_LOG)
705: PetscLogObjectState((PetscObject)A,"Rows=%D, Cols=%D, NZ=%D",A->rmap.n,A->cmap.n,a->nz);
706: #endif
707: if (a->freedata){
708: MatSeqXAIJFreeAIJ(a->singlemalloc,&a->a,&a->j,&a->i);
709: }
710: if (a->row) {
711: ISDestroy(a->row);
712: }
713: if (a->col) {
714: ISDestroy(a->col);
715: }
716: PetscFree(a->diag);
717: PetscFree2(a->imax,a->ilen);
718: PetscFree(a->idiag);
719: PetscFree(a->solve_work);
720: if (a->icol) {ISDestroy(a->icol);}
721: PetscFree(a->saved_values);
722: if (a->coloring) {ISColoringDestroy(a->coloring);}
723: PetscFree(a->xtoy);
724: if (a->compressedrow.use){PetscFree(a->compressedrow.i);}
726: MatDestroy_Inode(A);
728: PetscFree(a);
730: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSeqAIJSetColumnIndices_C","",PETSC_NULL);
731: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatStoreValues_C","",PETSC_NULL);
732: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatRetrieveValues_C","",PETSC_NULL);
733: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatConvert_seqaij_seqsbaij_C","",PETSC_NULL);
734: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatConvert_seqaij_seqbaij_C","",PETSC_NULL);
735: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatIsTranspose_C","",PETSC_NULL);
736: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSeqAIJSetPreallocation_C","",PETSC_NULL);
737: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSeqAIJSetPreallocationCSR_C","",PETSC_NULL);
738: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatReorderForNonzeroDiagonal_C","",PETSC_NULL);
739: return(0);
740: }
744: PetscErrorCode MatCompress_SeqAIJ(Mat A)
745: {
747: return(0);
748: }
752: PetscErrorCode MatSetOption_SeqAIJ(Mat A,MatOption op)
753: {
754: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
758: switch (op) {
759: case MAT_ROW_ORIENTED:
760: a->roworiented = PETSC_TRUE;
761: break;
762: case MAT_KEEP_ZEROED_ROWS:
763: a->keepzeroedrows = PETSC_TRUE;
764: break;
765: case MAT_COLUMN_ORIENTED:
766: a->roworiented = PETSC_FALSE;
767: break;
768: case MAT_COLUMNS_SORTED:
769: a->sorted = PETSC_TRUE;
770: break;
771: case MAT_COLUMNS_UNSORTED:
772: a->sorted = PETSC_FALSE;
773: break;
774: case MAT_NO_NEW_NONZERO_LOCATIONS:
775: a->nonew = 1;
776: break;
777: case MAT_NEW_NONZERO_LOCATION_ERR:
778: a->nonew = -1;
779: break;
780: case MAT_NEW_NONZERO_ALLOCATION_ERR:
781: a->nonew = -2;
782: break;
783: case MAT_YES_NEW_NONZERO_LOCATIONS:
784: a->nonew = 0;
785: break;
786: case MAT_IGNORE_ZERO_ENTRIES:
787: a->ignorezeroentries = PETSC_TRUE;
788: break;
789: case MAT_USE_COMPRESSEDROW:
790: a->compressedrow.use = PETSC_TRUE;
791: break;
792: case MAT_DO_NOT_USE_COMPRESSEDROW:
793: a->compressedrow.use = PETSC_FALSE;
794: break;
795: case MAT_ROWS_SORTED:
796: case MAT_ROWS_UNSORTED:
797: case MAT_YES_NEW_DIAGONALS:
798: case MAT_IGNORE_OFF_PROC_ENTRIES:
799: case MAT_USE_HASH_TABLE:
800: PetscInfo(A,"Option ignored\n");
801: break;
802: case MAT_NO_NEW_DIAGONALS:
803: SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
804: default:
805: break;
806: }
807: MatSetOption_Inode(A,op);
808: return(0);
809: }
813: PetscErrorCode MatGetDiagonal_SeqAIJ(Mat A,Vec v)
814: {
815: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
817: PetscInt i,j,n;
818: PetscScalar *x,zero = 0.0;
821: VecSet(v,zero);
822: VecGetArray(v,&x);
823: VecGetLocalSize(v,&n);
824: if (n != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
825: for (i=0; i<A->rmap.n; i++) {
826: for (j=a->i[i]; j<a->i[i+1]; j++) {
827: if (a->j[j] == i) {
828: x[i] = a->a[j];
829: break;
830: }
831: }
832: }
833: VecRestoreArray(v,&x);
834: return(0);
835: }
839: PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
840: {
841: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
842: PetscScalar *x,*y;
843: PetscErrorCode ierr;
844: PetscInt m = A->rmap.n;
845: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
846: PetscScalar *v,alpha;
847: PetscInt n,i,*idx,*ii,*ridx=PETSC_NULL;
848: Mat_CompressedRow cprow = a->compressedrow;
849: PetscTruth usecprow = cprow.use;
850: #endif
853: if (zz != yy) {VecCopy(zz,yy);}
854: VecGetArray(xx,&x);
855: VecGetArray(yy,&y);
857: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
858: fortranmulttransposeaddaij_(&m,x,a->i,a->j,a->a,y);
859: #else
860: if (usecprow){
861: m = cprow.nrows;
862: ii = cprow.i;
863: ridx = cprow.rindex;
864: } else {
865: ii = a->i;
866: }
867: for (i=0; i<m; i++) {
868: idx = a->j + ii[i] ;
869: v = a->a + ii[i] ;
870: n = ii[i+1] - ii[i];
871: if (usecprow){
872: alpha = x[ridx[i]];
873: } else {
874: alpha = x[i];
875: }
876: while (n-->0) {y[*idx++] += alpha * *v++;}
877: }
878: #endif
879: PetscLogFlops(2*a->nz);
880: VecRestoreArray(xx,&x);
881: VecRestoreArray(yy,&y);
882: return(0);
883: }
887: PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec xx,Vec yy)
888: {
889: PetscScalar zero = 0.0;
893: VecSet(yy,zero);
894: MatMultTransposeAdd_SeqAIJ(A,xx,yy,yy);
895: return(0);
896: }
901: PetscErrorCode MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
902: {
903: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
904: PetscScalar *x,*y,*aa;
906: PetscInt m=A->rmap.n,*aj,*ii;
907: PetscInt n,i,j,*ridx=PETSC_NULL;
908: PetscScalar sum;
909: PetscTruth usecprow=a->compressedrow.use;
910: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
911: PetscInt jrow;
912: #endif
914: #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
915: #pragma disjoint(*x,*y,*aa)
916: #endif
919: VecGetArray(xx,&x);
920: VecGetArray(yy,&y);
921: aj = a->j;
922: aa = a->a;
923: ii = a->i;
924: if (usecprow){ /* use compressed row format */
925: m = a->compressedrow.nrows;
926: ii = a->compressedrow.i;
927: ridx = a->compressedrow.rindex;
928: for (i=0; i<m; i++){
929: n = ii[i+1] - ii[i];
930: aj = a->j + ii[i];
931: aa = a->a + ii[i];
932: sum = 0.0;
933: for (j=0; j<n; j++) sum += (*aa++)*x[*aj++];
934: y[*ridx++] = sum;
935: }
936: } else { /* do not use compressed row format */
937: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
938: fortranmultaij_(&m,x,ii,aj,aa,y);
939: #else
940: for (i=0; i<m; i++) {
941: jrow = ii[i];
942: n = ii[i+1] - jrow;
943: sum = 0.0;
944: for (j=0; j<n; j++) {
945: sum += aa[jrow]*x[aj[jrow]]; jrow++;
946: }
947: y[i] = sum;
948: }
949: #endif
950: }
951: PetscLogFlops(2*a->nz - m);
952: VecRestoreArray(xx,&x);
953: VecRestoreArray(yy,&y);
954: return(0);
955: }
959: PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
960: {
961: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
962: PetscScalar *x,*y,*z,*aa;
964: PetscInt m = A->rmap.n,*aj,*ii;
965: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
966: PetscInt n,i,jrow,j,*ridx=PETSC_NULL;
967: PetscScalar sum;
968: PetscTruth usecprow=a->compressedrow.use;
969: #endif
972: VecGetArray(xx,&x);
973: VecGetArray(yy,&y);
974: if (zz != yy) {
975: VecGetArray(zz,&z);
976: } else {
977: z = y;
978: }
980: aj = a->j;
981: aa = a->a;
982: ii = a->i;
983: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
984: fortranmultaddaij_(&m,x,ii,aj,aa,y,z);
985: #else
986: if (usecprow){ /* use compressed row format */
987: if (zz != yy){
988: PetscMemcpy(z,y,m*sizeof(PetscScalar));
989: }
990: m = a->compressedrow.nrows;
991: ii = a->compressedrow.i;
992: ridx = a->compressedrow.rindex;
993: for (i=0; i<m; i++){
994: n = ii[i+1] - ii[i];
995: aj = a->j + ii[i];
996: aa = a->a + ii[i];
997: sum = y[*ridx];
998: for (j=0; j<n; j++) sum += (*aa++)*x[*aj++];
999: z[*ridx++] = sum;
1000: }
1001: } else { /* do not use compressed row format */
1002: for (i=0; i<m; i++) {
1003: jrow = ii[i];
1004: n = ii[i+1] - jrow;
1005: sum = y[i];
1006: for (j=0; j<n; j++) {
1007: sum += aa[jrow]*x[aj[jrow]]; jrow++;
1008: }
1009: z[i] = sum;
1010: }
1011: }
1012: #endif
1013: PetscLogFlops(2*a->nz);
1014: VecRestoreArray(xx,&x);
1015: VecRestoreArray(yy,&y);
1016: if (zz != yy) {
1017: VecRestoreArray(zz,&z);
1018: }
1019: return(0);
1020: }
1022: /*
1023: Adds diagonal pointers to sparse matrix structure.
1024: */
1027: PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat A)
1028: {
1029: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1031: PetscInt i,j,*diag,m = A->rmap.n;
1034: if (a->diag) return(0);
1036: PetscMalloc((m+1)*sizeof(PetscInt),&diag);
1037: PetscLogObjectMemory(A,(m+1)*sizeof(PetscInt));
1038: for (i=0; i<A->rmap.n; i++) {
1039: diag[i] = a->i[i+1];
1040: for (j=a->i[i]; j<a->i[i+1]; j++) {
1041: if (a->j[j] == i) {
1042: diag[i] = j;
1043: break;
1044: }
1045: }
1046: }
1047: a->diag = diag;
1048: return(0);
1049: }
1051: /*
1052: Checks for missing diagonals
1053: */
1056: PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat A)
1057: {
1058: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1060: PetscInt *diag,*jj = a->j,i;
1063: if (A->rmap.n > 0 && !jj) SETERRQ(PETSC_ERR_PLIB,"Matrix has no entries therefor is missing diagonal");
1064: MatMarkDiagonal_SeqAIJ(A);
1065: diag = a->diag;
1066: for (i=0; i<A->rmap.n; i++) {
1067: if (jj[diag[i]] != i) {
1068: SETERRQ1(PETSC_ERR_PLIB,"Matrix is missing diagonal number %D",i);
1069: }
1070: }
1071: return(0);
1072: }
1076: PetscErrorCode MatRelax_SeqAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
1077: {
1078: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1079: PetscScalar *x,d,*xs,sum,*t,scale,*idiag=0,*mdiag;
1080: const PetscScalar *v = a->a, *b, *bs,*xb, *ts;
1081: PetscErrorCode ierr;
1082: PetscInt n = A->cmap.n,m = A->rmap.n,i;
1083: const PetscInt *idx,*diag;
1086: its = its*lits;
1087: if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
1089: if (!a->diag) {MatMarkDiagonal_SeqAIJ(A);}
1090: diag = a->diag;
1091: if (!a->idiag) {
1092: PetscMalloc(3*m*sizeof(PetscScalar),&a->idiag);
1093: a->ssor = a->idiag + m;
1094: mdiag = a->ssor + m;
1096: v = a->a;
1098: /* this is wrong when fshift omega changes each iteration */
1099: if (omega == 1.0 && !fshift) {
1100: for (i=0; i<m; i++) {
1101: mdiag[i] = v[diag[i]];
1102: a->idiag[i] = 1.0/v[diag[i]];
1103: }
1104: PetscLogFlops(m);
1105: } else {
1106: for (i=0; i<m; i++) {
1107: mdiag[i] = v[diag[i]];
1108: a->idiag[i] = omega/(fshift + v[diag[i]]);
1109: }
1110: PetscLogFlops(2*m);
1111: }
1112: }
1113: t = a->ssor;
1114: idiag = a->idiag;
1115: mdiag = a->idiag + 2*m;
1117: VecGetArray(xx,&x);
1118: if (xx != bb) {
1119: VecGetArray(bb,(PetscScalar**)&b);
1120: } else {
1121: b = x;
1122: }
1124: /* We count flops by assuming the upper triangular and lower triangular parts have the same number of nonzeros */
1125: xs = x;
1126: if (flag == SOR_APPLY_UPPER) {
1127: /* apply (U + D/omega) to the vector */
1128: bs = b;
1129: for (i=0; i<m; i++) {
1130: d = fshift + a->a[diag[i]];
1131: n = a->i[i+1] - diag[i] - 1;
1132: idx = a->j + diag[i] + 1;
1133: v = a->a + diag[i] + 1;
1134: sum = b[i]*d/omega;
1135: SPARSEDENSEDOT(sum,bs,v,idx,n);
1136: x[i] = sum;
1137: }
1138: VecRestoreArray(xx,&x);
1139: if (bb != xx) {VecRestoreArray(bb,(PetscScalar**)&b);}
1140: PetscLogFlops(a->nz);
1141: return(0);
1142: }
1145: /* Let A = L + U + D; where L is lower trianglar,
1146: U is upper triangular, E is diagonal; This routine applies
1148: (L + E)^{-1} A (U + E)^{-1}
1150: to a vector efficiently using Eisenstat's trick. This is for
1151: the case of SSOR preconditioner, so E is D/omega where omega
1152: is the relaxation factor.
1153: */
1155: if (flag == SOR_APPLY_LOWER) {
1156: SETERRQ(PETSC_ERR_SUP,"SOR_APPLY_LOWER is not implemented");
1157: } else if (flag & SOR_EISENSTAT) {
1158: /* Let A = L + U + D; where L is lower trianglar,
1159: U is upper triangular, E is diagonal; This routine applies
1161: (L + E)^{-1} A (U + E)^{-1}
1163: to a vector efficiently using Eisenstat's trick. This is for
1164: the case of SSOR preconditioner, so E is D/omega where omega
1165: is the relaxation factor.
1166: */
1167: scale = (2.0/omega) - 1.0;
1169: /* x = (E + U)^{-1} b */
1170: for (i=m-1; i>=0; i--) {
1171: n = a->i[i+1] - diag[i] - 1;
1172: idx = a->j + diag[i] + 1;
1173: v = a->a + diag[i] + 1;
1174: sum = b[i];
1175: SPARSEDENSEMDOT(sum,xs,v,idx,n);
1176: x[i] = sum*idiag[i];
1177: }
1179: /* t = b - (2*E - D)x */
1180: v = a->a;
1181: for (i=0; i<m; i++) { t[i] = b[i] - scale*(v[*diag++])*x[i]; }
1183: /* t = (E + L)^{-1}t */
1184: ts = t;
1185: diag = a->diag;
1186: for (i=0; i<m; i++) {
1187: n = diag[i] - a->i[i];
1188: idx = a->j + a->i[i];
1189: v = a->a + a->i[i];
1190: sum = t[i];
1191: SPARSEDENSEMDOT(sum,ts,v,idx,n);
1192: t[i] = sum*idiag[i];
1193: /* x = x + t */
1194: x[i] += t[i];
1195: }
1197: PetscLogFlops(6*m-1 + 2*a->nz);
1198: VecRestoreArray(xx,&x);
1199: if (bb != xx) {VecRestoreArray(bb,(PetscScalar**)&b);}
1200: return(0);
1201: }
1202: if (flag & SOR_ZERO_INITIAL_GUESS) {
1203: if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
1204: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1205: fortranrelaxaijforwardzero_(&m,&omega,x,a->i,a->j,(PetscInt*)diag,idiag,a->a,(void*)b);
1206: #else
1207: for (i=0; i<m; i++) {
1208: n = diag[i] - a->i[i];
1209: idx = a->j + a->i[i];
1210: v = a->a + a->i[i];
1211: sum = b[i];
1212: SPARSEDENSEMDOT(sum,xs,v,idx,n);
1213: x[i] = sum*idiag[i];
1214: }
1215: #endif
1216: xb = x;
1217: PetscLogFlops(a->nz);
1218: } else xb = b;
1219: if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
1220: (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
1221: for (i=0; i<m; i++) {
1222: x[i] *= mdiag[i];
1223: }
1224: PetscLogFlops(m);
1225: }
1226: if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
1227: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1228: fortranrelaxaijbackwardzero_(&m,&omega,x,a->i,a->j,(PetscInt*)diag,idiag,a->a,(void*)xb);
1229: #else
1230: for (i=m-1; i>=0; i--) {
1231: n = a->i[i+1] - diag[i] - 1;
1232: idx = a->j + diag[i] + 1;
1233: v = a->a + diag[i] + 1;
1234: sum = xb[i];
1235: SPARSEDENSEMDOT(sum,xs,v,idx,n);
1236: x[i] = sum*idiag[i];
1237: }
1238: #endif
1239: PetscLogFlops(a->nz);
1240: }
1241: its--;
1242: }
1243: while (its--) {
1244: if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
1245: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1246: fortranrelaxaijforward_(&m,&omega,x,a->i,a->j,(PetscInt*)diag,a->a,(void*)b);
1247: #else
1248: for (i=0; i<m; i++) {
1249: d = fshift + a->a[diag[i]];
1250: n = a->i[i+1] - a->i[i];
1251: idx = a->j + a->i[i];
1252: v = a->a + a->i[i];
1253: sum = b[i];
1254: SPARSEDENSEMDOT(sum,xs,v,idx,n);
1255: x[i] = (1. - omega)*x[i] + (sum + mdiag[i]*x[i])*idiag[i];
1256: }
1257: #endif
1258: PetscLogFlops(a->nz);
1259: }
1260: if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
1261: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1262: fortranrelaxaijbackward_(&m,&omega,x,a->i,a->j,(PetscInt*)diag,a->a,(void*)b);
1263: #else
1264: for (i=m-1; i>=0; i--) {
1265: d = fshift + a->a[diag[i]];
1266: n = a->i[i+1] - a->i[i];
1267: idx = a->j + a->i[i];
1268: v = a->a + a->i[i];
1269: sum = b[i];
1270: SPARSEDENSEMDOT(sum,xs,v,idx,n);
1271: x[i] = (1. - omega)*x[i] + (sum + mdiag[i]*x[i])*idiag[i];
1272: }
1273: #endif
1274: PetscLogFlops(a->nz);
1275: }
1276: }
1277: VecRestoreArray(xx,&x);
1278: if (bb != xx) {VecRestoreArray(bb,(PetscScalar**)&b);}
1279: return(0);
1280: }
1284: PetscErrorCode MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,MatInfo *info)
1285: {
1286: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1289: info->rows_global = (double)A->rmap.n;
1290: info->columns_global = (double)A->cmap.n;
1291: info->rows_local = (double)A->rmap.n;
1292: info->columns_local = (double)A->cmap.n;
1293: info->block_size = 1.0;
1294: info->nz_allocated = (double)a->maxnz;
1295: info->nz_used = (double)a->nz;
1296: info->nz_unneeded = (double)(a->maxnz - a->nz);
1297: info->assemblies = (double)A->num_ass;
1298: info->mallocs = (double)a->reallocs;
1299: info->memory = A->mem;
1300: if (A->factor) {
1301: info->fill_ratio_given = A->info.fill_ratio_given;
1302: info->fill_ratio_needed = A->info.fill_ratio_needed;
1303: info->factor_mallocs = A->info.factor_mallocs;
1304: } else {
1305: info->fill_ratio_given = 0;
1306: info->fill_ratio_needed = 0;
1307: info->factor_mallocs = 0;
1308: }
1309: return(0);
1310: }
1314: PetscErrorCode MatZeroRows_SeqAIJ(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag)
1315: {
1316: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1317: PetscInt i,m = A->rmap.n - 1;
1321: if (a->keepzeroedrows) {
1322: for (i=0; i<N; i++) {
1323: if (rows[i] < 0 || rows[i] > m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range", rows[i]);
1324: PetscMemzero(&a->a[a->i[rows[i]]],a->ilen[rows[i]]*sizeof(PetscScalar));
1325: }
1326: if (diag != 0.0) {
1327: MatMissingDiagonal_SeqAIJ(A);
1328: MatMarkDiagonal_SeqAIJ(A);
1329: for (i=0; i<N; i++) {
1330: a->a[a->diag[rows[i]]] = diag;
1331: }
1332: }
1333: A->same_nonzero = PETSC_TRUE;
1334: } else {
1335: if (diag != 0.0) {
1336: for (i=0; i<N; i++) {
1337: if (rows[i] < 0 || rows[i] > m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range", rows[i]);
1338: if (a->ilen[rows[i]] > 0) {
1339: a->ilen[rows[i]] = 1;
1340: a->a[a->i[rows[i]]] = diag;
1341: a->j[a->i[rows[i]]] = rows[i];
1342: } else { /* in case row was completely empty */
1343: MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],&diag,INSERT_VALUES);
1344: }
1345: }
1346: } else {
1347: for (i=0; i<N; i++) {
1348: if (rows[i] < 0 || rows[i] > m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range", rows[i]);
1349: a->ilen[rows[i]] = 0;
1350: }
1351: }
1352: A->same_nonzero = PETSC_FALSE;
1353: }
1354: MatAssemblyEnd_SeqAIJ(A,MAT_FINAL_ASSEMBLY);
1355: return(0);
1356: }
1360: PetscErrorCode MatGetRow_SeqAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
1361: {
1362: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1363: PetscInt *itmp;
1366: if (row < 0 || row >= A->rmap.n) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Row %D out of range",row);
1368: *nz = a->i[row+1] - a->i[row];
1369: if (v) *v = a->a + a->i[row];
1370: if (idx) {
1371: itmp = a->j + a->i[row];
1372: if (*nz) {
1373: *idx = itmp;
1374: }
1375: else *idx = 0;
1376: }
1377: return(0);
1378: }
1380: /* remove this function? */
1383: PetscErrorCode MatRestoreRow_SeqAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
1384: {
1386: return(0);
1387: }
1391: PetscErrorCode MatNorm_SeqAIJ(Mat A,NormType type,PetscReal *nrm)
1392: {
1393: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1394: PetscScalar *v = a->a;
1395: PetscReal sum = 0.0;
1397: PetscInt i,j;
1400: if (type == NORM_FROBENIUS) {
1401: for (i=0; i<a->nz; i++) {
1402: #if defined(PETSC_USE_COMPLEX)
1403: sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1404: #else
1405: sum += (*v)*(*v); v++;
1406: #endif
1407: }
1408: *nrm = sqrt(sum);
1409: } else if (type == NORM_1) {
1410: PetscReal *tmp;
1411: PetscInt *jj = a->j;
1412: PetscMalloc((A->cmap.n+1)*sizeof(PetscReal),&tmp);
1413: PetscMemzero(tmp,A->cmap.n*sizeof(PetscReal));
1414: *nrm = 0.0;
1415: for (j=0; j<a->nz; j++) {
1416: tmp[*jj++] += PetscAbsScalar(*v); v++;
1417: }
1418: for (j=0; j<A->cmap.n; j++) {
1419: if (tmp[j] > *nrm) *nrm = tmp[j];
1420: }
1421: PetscFree(tmp);
1422: } else if (type == NORM_INFINITY) {
1423: *nrm = 0.0;
1424: for (j=0; j<A->rmap.n; j++) {
1425: v = a->a + a->i[j];
1426: sum = 0.0;
1427: for (i=0; i<a->i[j+1]-a->i[j]; i++) {
1428: sum += PetscAbsScalar(*v); v++;
1429: }
1430: if (sum > *nrm) *nrm = sum;
1431: }
1432: } else {
1433: SETERRQ(PETSC_ERR_SUP,"No support for two norm");
1434: }
1435: return(0);
1436: }
1440: PetscErrorCode MatTranspose_SeqAIJ(Mat A,Mat *B)
1441: {
1442: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1443: Mat C;
1445: PetscInt i,*aj = a->j,*ai = a->i,m = A->rmap.n,len,*col;
1446: PetscScalar *array = a->a;
1449: if (!B && m != A->cmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Square matrix only for in-place");
1450: PetscMalloc((1+A->cmap.n)*sizeof(PetscInt),&col);
1451: PetscMemzero(col,(1+A->cmap.n)*sizeof(PetscInt));
1452:
1453: for (i=0; i<ai[m]; i++) col[aj[i]] += 1;
1454: MatCreate(A->comm,&C);
1455: MatSetSizes(C,A->cmap.n,m,A->cmap.n,m);
1456: MatSetType(C,A->type_name);
1457: MatSeqAIJSetPreallocation_SeqAIJ(C,0,col);
1458: PetscFree(col);
1459: for (i=0; i<m; i++) {
1460: len = ai[i+1]-ai[i];
1461: MatSetValues_SeqAIJ(C,len,aj,1,&i,array,INSERT_VALUES);
1462: array += len;
1463: aj += len;
1464: }
1466: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1467: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
1469: if (B) {
1470: *B = C;
1471: } else {
1472: MatHeaderCopy(A,C);
1473: }
1474: return(0);
1475: }
1480: PetscErrorCode PETSCMAT_DLLEXPORT MatIsTranspose_SeqAIJ(Mat A,Mat B,PetscReal tol,PetscTruth *f)
1481: {
1482: Mat_SeqAIJ *aij = (Mat_SeqAIJ *) A->data,*bij = (Mat_SeqAIJ*) A->data;
1483: PetscInt *adx,*bdx,*aii,*bii,*aptr,*bptr; PetscScalar *va,*vb;
1485: PetscInt ma,na,mb,nb, i;
1488: bij = (Mat_SeqAIJ *) B->data;
1489:
1490: MatGetSize(A,&ma,&na);
1491: MatGetSize(B,&mb,&nb);
1492: if (ma!=nb || na!=mb){
1493: *f = PETSC_FALSE;
1494: return(0);
1495: }
1496: aii = aij->i; bii = bij->i;
1497: adx = aij->j; bdx = bij->j;
1498: va = aij->a; vb = bij->a;
1499: PetscMalloc(ma*sizeof(PetscInt),&aptr);
1500: PetscMalloc(mb*sizeof(PetscInt),&bptr);
1501: for (i=0; i<ma; i++) aptr[i] = aii[i];
1502: for (i=0; i<mb; i++) bptr[i] = bii[i];
1504: *f = PETSC_TRUE;
1505: for (i=0; i<ma; i++) {
1506: while (aptr[i]<aii[i+1]) {
1507: PetscInt idc,idr;
1508: PetscScalar vc,vr;
1509: /* column/row index/value */
1510: idc = adx[aptr[i]];
1511: idr = bdx[bptr[idc]];
1512: vc = va[aptr[i]];
1513: vr = vb[bptr[idc]];
1514: if (i!=idr || PetscAbsScalar(vc-vr) > tol) {
1515: *f = PETSC_FALSE;
1516: goto done;
1517: } else {
1518: aptr[i]++;
1519: if (B || i!=idc) bptr[idc]++;
1520: }
1521: }
1522: }
1523: done:
1524: PetscFree(aptr);
1525: if (B) {
1526: PetscFree(bptr);
1527: }
1528: return(0);
1529: }
1534: PetscErrorCode MatIsSymmetric_SeqAIJ(Mat A,PetscReal tol,PetscTruth *f)
1535: {
1538: MatIsTranspose_SeqAIJ(A,A,tol,f);
1539: return(0);
1540: }
1544: PetscErrorCode MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr)
1545: {
1546: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1547: PetscScalar *l,*r,x,*v;
1549: PetscInt i,j,m = A->rmap.n,n = A->cmap.n,M,nz = a->nz,*jj;
1552: if (ll) {
1553: /* The local size is used so that VecMPI can be passed to this routine
1554: by MatDiagonalScale_MPIAIJ */
1555: VecGetLocalSize(ll,&m);
1556: if (m != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length");
1557: VecGetArray(ll,&l);
1558: v = a->a;
1559: for (i=0; i<m; i++) {
1560: x = l[i];
1561: M = a->i[i+1] - a->i[i];
1562: for (j=0; j<M; j++) { (*v++) *= x;}
1563: }
1564: VecRestoreArray(ll,&l);
1565: PetscLogFlops(nz);
1566: }
1567: if (rr) {
1568: VecGetLocalSize(rr,&n);
1569: if (n != A->cmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vector wrong length");
1570: VecGetArray(rr,&r);
1571: v = a->a; jj = a->j;
1572: for (i=0; i<nz; i++) {
1573: (*v++) *= r[*jj++];
1574: }
1575: VecRestoreArray(rr,&r);
1576: PetscLogFlops(nz);
1577: }
1578: return(0);
1579: }
1583: PetscErrorCode MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,PetscInt csize,MatReuse scall,Mat *B)
1584: {
1585: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*c;
1587: PetscInt *smap,i,k,kstart,kend,oldcols = A->cmap.n,*lens;
1588: PetscInt row,mat_i,*mat_j,tcol,first,step,*mat_ilen,sum,lensi;
1589: PetscInt *irow,*icol,nrows,ncols;
1590: PetscInt *starts,*j_new,*i_new,*aj = a->j,*ai = a->i,ii,*ailen = a->ilen;
1591: PetscScalar *a_new,*mat_a;
1592: Mat C;
1593: PetscTruth stride;
1596: ISSorted(isrow,(PetscTruth*)&i);
1597: if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"ISrow is not sorted");
1598: ISSorted(iscol,(PetscTruth*)&i);
1599: if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IScol is not sorted");
1601: ISGetIndices(isrow,&irow);
1602: ISGetLocalSize(isrow,&nrows);
1603: ISGetLocalSize(iscol,&ncols);
1605: ISStrideGetInfo(iscol,&first,&step);
1606: ISStride(iscol,&stride);
1607: if (stride && step == 1) {
1608: /* special case of contiguous rows */
1609: PetscMalloc((2*nrows+1)*sizeof(PetscInt),&lens);
1610: starts = lens + nrows;
1611: /* loop over new rows determining lens and starting points */
1612: for (i=0; i<nrows; i++) {
1613: kstart = ai[irow[i]];
1614: kend = kstart + ailen[irow[i]];
1615: for (k=kstart; k<kend; k++) {
1616: if (aj[k] >= first) {
1617: starts[i] = k;
1618: break;
1619: }
1620: }
1621: sum = 0;
1622: while (k < kend) {
1623: if (aj[k++] >= first+ncols) break;
1624: sum++;
1625: }
1626: lens[i] = sum;
1627: }
1628: /* create submatrix */
1629: if (scall == MAT_REUSE_MATRIX) {
1630: PetscInt n_cols,n_rows;
1631: MatGetSize(*B,&n_rows,&n_cols);
1632: if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size");
1633: MatZeroEntries(*B);
1634: C = *B;
1635: } else {
1636: MatCreate(A->comm,&C);
1637: MatSetSizes(C,nrows,ncols,PETSC_DETERMINE,PETSC_DETERMINE);
1638: MatSetType(C,A->type_name);
1639: MatSeqAIJSetPreallocation_SeqAIJ(C,0,lens);
1640: }
1641: c = (Mat_SeqAIJ*)C->data;
1643: /* loop over rows inserting into submatrix */
1644: a_new = c->a;
1645: j_new = c->j;
1646: i_new = c->i;
1648: for (i=0; i<nrows; i++) {
1649: ii = starts[i];
1650: lensi = lens[i];
1651: for (k=0; k<lensi; k++) {
1652: *j_new++ = aj[ii+k] - first;
1653: }
1654: PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(PetscScalar));
1655: a_new += lensi;
1656: i_new[i+1] = i_new[i] + lensi;
1657: c->ilen[i] = lensi;
1658: }
1659: PetscFree(lens);
1660: } else {
1661: ISGetIndices(iscol,&icol);
1662: PetscMalloc((1+oldcols)*sizeof(PetscInt),&smap);
1663:
1664: PetscMalloc((1+nrows)*sizeof(PetscInt),&lens);
1665: PetscMemzero(smap,oldcols*sizeof(PetscInt));
1666: for (i=0; i<ncols; i++) smap[icol[i]] = i+1;
1667: /* determine lens of each row */
1668: for (i=0; i<nrows; i++) {
1669: kstart = ai[irow[i]];
1670: kend = kstart + a->ilen[irow[i]];
1671: lens[i] = 0;
1672: for (k=kstart; k<kend; k++) {
1673: if (smap[aj[k]]) {
1674: lens[i]++;
1675: }
1676: }
1677: }
1678: /* Create and fill new matrix */
1679: if (scall == MAT_REUSE_MATRIX) {
1680: PetscTruth equal;
1682: c = (Mat_SeqAIJ *)((*B)->data);
1683: if ((*B)->rmap.n != nrows || (*B)->cmap.n != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size");
1684: PetscMemcmp(c->ilen,lens,(*B)->rmap.n*sizeof(PetscInt),&equal);
1685: if (!equal) {
1686: SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros");
1687: }
1688: PetscMemzero(c->ilen,(*B)->rmap.n*sizeof(PetscInt));
1689: C = *B;
1690: } else {
1691: MatCreate(A->comm,&C);
1692: MatSetSizes(C,nrows,ncols,PETSC_DETERMINE,PETSC_DETERMINE);
1693: MatSetType(C,A->type_name);
1694: MatSeqAIJSetPreallocation_SeqAIJ(C,0,lens);
1695: }
1696: c = (Mat_SeqAIJ *)(C->data);
1697: for (i=0; i<nrows; i++) {
1698: row = irow[i];
1699: kstart = ai[row];
1700: kend = kstart + a->ilen[row];
1701: mat_i = c->i[i];
1702: mat_j = c->j + mat_i;
1703: mat_a = c->a + mat_i;
1704: mat_ilen = c->ilen + i;
1705: for (k=kstart; k<kend; k++) {
1706: if ((tcol=smap[a->j[k]])) {
1707: *mat_j++ = tcol - 1;
1708: *mat_a++ = a->a[k];
1709: (*mat_ilen)++;
1711: }
1712: }
1713: }
1714: /* Free work space */
1715: ISRestoreIndices(iscol,&icol);
1716: PetscFree(smap);
1717: PetscFree(lens);
1718: }
1719: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1720: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
1722: ISRestoreIndices(isrow,&irow);
1723: *B = C;
1724: return(0);
1725: }
1727: /*
1728: */
1731: PetscErrorCode MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,MatFactorInfo *info)
1732: {
1733: Mat_SeqAIJ *a = (Mat_SeqAIJ*)inA->data;
1735: Mat outA;
1736: PetscTruth row_identity,col_identity;
1739: if (info->levels != 0) SETERRQ(PETSC_ERR_SUP,"Only levels=0 supported for in-place ilu");
1740: ISIdentity(row,&row_identity);
1741: ISIdentity(col,&col_identity);
1742: if (!row_identity || !col_identity) {
1743: SETERRQ(PETSC_ERR_ARG_WRONG,"Row and column permutations must be identity for in-place ILU");
1744: }
1746: outA = inA;
1747: inA->factor = FACTOR_LU;
1748: a->row = row;
1749: a->col = col;
1750: PetscObjectReference((PetscObject)row);
1751: PetscObjectReference((PetscObject)col);
1753: /* Create the inverse permutation so that it can be used in MatLUFactorNumeric() */
1754: if (a->icol) {ISDestroy(a->icol);} /* need to remove old one */
1755: ISInvertPermutation(col,PETSC_DECIDE,&a->icol);
1756: PetscLogObjectParent(inA,a->icol);
1758: if (!a->solve_work) { /* this matrix may have been factored before */
1759: PetscMalloc((inA->rmap.n+1)*sizeof(PetscScalar),&a->solve_work);
1760: }
1762: if (!a->diag) {
1763: MatMarkDiagonal_SeqAIJ(inA);
1764: }
1765: MatLUFactorNumeric_SeqAIJ(inA,info,&outA);
1766: return(0);
1767: }
1769: #include petscblaslapack.h
1772: PetscErrorCode MatScale_SeqAIJ(Mat inA,PetscScalar alpha)
1773: {
1774: Mat_SeqAIJ *a = (Mat_SeqAIJ*)inA->data;
1775: PetscBLASInt bnz = (PetscBLASInt)a->nz,one = 1;
1776: PetscScalar oalpha = alpha;
1781: BLASscal_(&bnz,&oalpha,a->a,&one);
1782: PetscLogFlops(a->nz);
1783: return(0);
1784: }
1788: PetscErrorCode MatGetSubMatrices_SeqAIJ(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1789: {
1791: PetscInt i;
1794: if (scall == MAT_INITIAL_MATRIX) {
1795: PetscMalloc((n+1)*sizeof(Mat),B);
1796: }
1798: for (i=0; i<n; i++) {
1799: MatGetSubMatrix_SeqAIJ(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);
1800: }
1801: return(0);
1802: }
1806: PetscErrorCode MatIncreaseOverlap_SeqAIJ(Mat A,PetscInt is_max,IS is[],PetscInt ov)
1807: {
1808: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1810: PetscInt row,i,j,k,l,m,n,*idx,*nidx,isz,val;
1811: PetscInt start,end,*ai,*aj;
1812: PetscBT table;
1815: m = A->rmap.n;
1816: ai = a->i;
1817: aj = a->j;
1819: if (ov < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"illegal negative overlap value used");
1821: PetscMalloc((m+1)*sizeof(PetscInt),&nidx);
1822: PetscBTCreate(m,table);
1824: for (i=0; i<is_max; i++) {
1825: /* Initialize the two local arrays */
1826: isz = 0;
1827: PetscBTMemzero(m,table);
1828:
1829: /* Extract the indices, assume there can be duplicate entries */
1830: ISGetIndices(is[i],&idx);
1831: ISGetLocalSize(is[i],&n);
1832:
1833: /* Enter these into the temp arrays. I.e., mark table[row], enter row into new index */
1834: for (j=0; j<n ; ++j){
1835: if(!PetscBTLookupSet(table,idx[j])) { nidx[isz++] = idx[j];}
1836: }
1837: ISRestoreIndices(is[i],&idx);
1838: ISDestroy(is[i]);
1839:
1840: k = 0;
1841: for (j=0; j<ov; j++){ /* for each overlap */
1842: n = isz;
1843: for (; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */
1844: row = nidx[k];
1845: start = ai[row];
1846: end = ai[row+1];
1847: for (l = start; l<end ; l++){
1848: val = aj[l] ;
1849: if (!PetscBTLookupSet(table,val)) {nidx[isz++] = val;}
1850: }
1851: }
1852: }
1853: ISCreateGeneral(PETSC_COMM_SELF,isz,nidx,(is+i));
1854: }
1855: PetscBTDestroy(table);
1856: PetscFree(nidx);
1857: return(0);
1858: }
1860: /* -------------------------------------------------------------- */
1863: PetscErrorCode MatPermute_SeqAIJ(Mat A,IS rowp,IS colp,Mat *B)
1864: {
1865: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1867: PetscInt i,nz,m = A->rmap.n,n = A->cmap.n,*col;
1868: PetscInt *row,*cnew,j,*lens;
1869: IS icolp,irowp;
1870: PetscInt *cwork;
1871: PetscScalar *vwork;
1874: ISInvertPermutation(rowp,PETSC_DECIDE,&irowp);
1875: ISGetIndices(irowp,&row);
1876: ISInvertPermutation(colp,PETSC_DECIDE,&icolp);
1877: ISGetIndices(icolp,&col);
1878:
1879: /* determine lengths of permuted rows */
1880: PetscMalloc((m+1)*sizeof(PetscInt),&lens);
1881: for (i=0; i<m; i++) {
1882: lens[row[i]] = a->i[i+1] - a->i[i];
1883: }
1884: MatCreate(A->comm,B);
1885: MatSetSizes(*B,m,n,m,n);
1886: MatSetType(*B,A->type_name);
1887: MatSeqAIJSetPreallocation_SeqAIJ(*B,0,lens);
1888: PetscFree(lens);
1890: PetscMalloc(n*sizeof(PetscInt),&cnew);
1891: for (i=0; i<m; i++) {
1892: MatGetRow_SeqAIJ(A,i,&nz,&cwork,&vwork);
1893: for (j=0; j<nz; j++) { cnew[j] = col[cwork[j]];}
1894: MatSetValues_SeqAIJ(*B,1,&row[i],nz,cnew,vwork,INSERT_VALUES);
1895: MatRestoreRow_SeqAIJ(A,i,&nz,&cwork,&vwork);
1896: }
1897: PetscFree(cnew);
1898: (*B)->assembled = PETSC_FALSE;
1899: MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);
1900: MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);
1901: ISRestoreIndices(irowp,&row);
1902: ISRestoreIndices(icolp,&col);
1903: ISDestroy(irowp);
1904: ISDestroy(icolp);
1905: return(0);
1906: }
1910: PetscErrorCode MatPrintHelp_SeqAIJ(Mat A)
1911: {
1912: static PetscTruth called = PETSC_FALSE;
1913: MPI_Comm comm = A->comm;
1914: PetscErrorCode ierr;
1917: MatPrintHelp_Inode(A);
1918: if (called) {return(0);} else called = PETSC_TRUE;
1919: (*PetscHelpPrintf)(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n");
1920: (*PetscHelpPrintf)(comm," -mat_lu_pivotthreshold <threshold>: Set pivoting threshold\n");
1921: (*PetscHelpPrintf)(comm," -mat_aij_oneindex: internal indices begin at 1 instead of the default 0.\n");
1922: (*PetscHelpPrintf)(comm," -mat_no_compressedrow: Do not use compressedrow\n");
1923: return(0);
1924: }
1928: PetscErrorCode MatCopy_SeqAIJ(Mat A,Mat B,MatStructure str)
1929: {
1933: /* If the two matrices have the same copy implementation, use fast copy. */
1934: if (str == SAME_NONZERO_PATTERN && (A->ops->copy == B->ops->copy)) {
1935: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1936: Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data;
1938: if (a->i[A->rmap.n] != b->i[B->rmap.n]) {
1939: SETERRQ(PETSC_ERR_ARG_INCOMP,"Number of nonzeros in two matrices are different");
1940: }
1941: PetscMemcpy(b->a,a->a,(a->i[A->rmap.n])*sizeof(PetscScalar));
1942: } else {
1943: MatCopy_Basic(A,B,str);
1944: }
1945: return(0);
1946: }
1950: PetscErrorCode MatSetUpPreallocation_SeqAIJ(Mat A)
1951: {
1955: MatSeqAIJSetPreallocation_SeqAIJ(A,PETSC_DEFAULT,0);
1956: return(0);
1957: }
1961: PetscErrorCode MatGetArray_SeqAIJ(Mat A,PetscScalar *array[])
1962: {
1963: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1965: *array = a->a;
1966: return(0);
1967: }
1971: PetscErrorCode MatRestoreArray_SeqAIJ(Mat A,PetscScalar *array[])
1972: {
1974: return(0);
1975: }
1979: PetscErrorCode MatFDColoringApply_SeqAIJ(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
1980: {
1981: PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void *))coloring->f;
1983: PetscInt k,N,start,end,l,row,col,srow,**vscaleforrow,m1,m2;
1984: PetscScalar dx,*y,*xx,*w3_array;
1985: PetscScalar *vscale_array;
1986: PetscReal epsilon = coloring->error_rel,umin = coloring->umin;
1987: Vec w1,w2,w3;
1988: void *fctx = coloring->fctx;
1989: PetscTruth flg;
1992: if (!coloring->w1) {
1993: VecDuplicate(x1,&coloring->w1);
1994: PetscLogObjectParent(coloring,coloring->w1);
1995: VecDuplicate(x1,&coloring->w2);
1996: PetscLogObjectParent(coloring,coloring->w2);
1997: VecDuplicate(x1,&coloring->w3);
1998: PetscLogObjectParent(coloring,coloring->w3);
1999: }
2000: w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
2002: MatSetUnfactored(J);
2003: PetscOptionsHasName(coloring->prefix,"-mat_fd_coloring_dont_rezero",&flg);
2004: if (flg) {
2005: PetscInfo(coloring,"Not calling MatZeroEntries()\n");
2006: } else {
2007: PetscTruth assembled;
2008: MatAssembled(J,&assembled);
2009: if (assembled) {
2010: MatZeroEntries(J);
2011: }
2012: }
2014: VecGetOwnershipRange(x1,&start,&end);
2015: VecGetSize(x1,&N);
2017: /*
2018: This is a horrible, horrible, hack. See DMMGComputeJacobian_Multigrid() it inproperly sets
2019: coloring->F for the coarser grids from the finest
2020: */
2021: if (coloring->F) {
2022: VecGetLocalSize(coloring->F,&m1);
2023: VecGetLocalSize(w1,&m2);
2024: if (m1 != m2) {
2025: coloring->F = 0;
2026: }
2027: }
2029: if (coloring->F) {
2030: w1 = coloring->F;
2031: coloring->F = 0;
2032: } else {
2033: PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);
2034: (*f)(sctx,x1,w1,fctx);
2035: PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);
2036: }
2038: /*
2039: Compute all the scale factors and share with other processors
2040: */
2041: VecGetArray(x1,&xx);xx = xx - start;
2042: VecGetArray(coloring->vscale,&vscale_array);vscale_array = vscale_array - start;
2043: for (k=0; k<coloring->ncolors; k++) {
2044: /*
2045: Loop over each column associated with color adding the
2046: perturbation to the vector w3.
2047: */
2048: for (l=0; l<coloring->ncolumns[k]; l++) {
2049: col = coloring->columns[k][l]; /* column of the matrix we are probing for */
2050: dx = xx[col];
2051: if (dx == 0.0) dx = 1.0;
2052: #if !defined(PETSC_USE_COMPLEX)
2053: if (dx < umin && dx >= 0.0) dx = umin;
2054: else if (dx < 0.0 && dx > -umin) dx = -umin;
2055: #else
2056: if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin;
2057: else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
2058: #endif
2059: dx *= epsilon;
2060: vscale_array[col] = 1.0/dx;
2061: }
2062: }
2063: vscale_array = vscale_array + start;VecRestoreArray(coloring->vscale,&vscale_array);
2064: VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);
2065: VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);
2067: /* VecView(coloring->vscale,PETSC_VIEWER_STDOUT_WORLD);
2068: VecView(x1,PETSC_VIEWER_STDOUT_WORLD);*/
2070: if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
2071: else vscaleforrow = coloring->columnsforrow;
2073: VecGetArray(coloring->vscale,&vscale_array);
2074: /*
2075: Loop over each color
2076: */
2077: for (k=0; k<coloring->ncolors; k++) {
2078: coloring->currentcolor = k;
2079: VecCopy(x1,w3);
2080: VecGetArray(w3,&w3_array);w3_array = w3_array - start;
2081: /*
2082: Loop over each column associated with color adding the
2083: perturbation to the vector w3.
2084: */
2085: for (l=0; l<coloring->ncolumns[k]; l++) {
2086: col = coloring->columns[k][l]; /* column of the matrix we are probing for */
2087: dx = xx[col];
2088: if (dx == 0.0) dx = 1.0;
2089: #if !defined(PETSC_USE_COMPLEX)
2090: if (dx < umin && dx >= 0.0) dx = umin;
2091: else if (dx < 0.0 && dx > -umin) dx = -umin;
2092: #else
2093: if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin;
2094: else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
2095: #endif
2096: dx *= epsilon;
2097: if (!PetscAbsScalar(dx)) SETERRQ(PETSC_ERR_PLIB,"Computed 0 differencing parameter");
2098: w3_array[col] += dx;
2099: }
2100: w3_array = w3_array + start; VecRestoreArray(w3,&w3_array);
2102: /*
2103: Evaluate function at x1 + dx (here dx is a vector of perturbations)
2104: */
2106: PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);
2107: (*f)(sctx,w3,w2,fctx);
2108: PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);
2109: VecAXPY(w2,-1.0,w1);
2111: /*
2112: Loop over rows of vector, putting results into Jacobian matrix
2113: */
2114: VecGetArray(w2,&y);
2115: for (l=0; l<coloring->nrows[k]; l++) {
2116: row = coloring->rows[k][l];
2117: col = coloring->columnsforrow[k][l];
2118: y[row] *= vscale_array[vscaleforrow[k][l]];
2119: srow = row + start;
2120: MatSetValues_SeqAIJ(J,1,&srow,1,&col,y+row,INSERT_VALUES);
2121: }
2122: VecRestoreArray(w2,&y);
2123: }
2124: coloring->currentcolor = k;
2125: VecRestoreArray(coloring->vscale,&vscale_array);
2126: xx = xx + start; VecRestoreArray(x1,&xx);
2127: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
2128: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
2129: return(0);
2130: }
2132: #include petscblaslapack.h
2135: PetscErrorCode MatAXPY_SeqAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str)
2136: {
2138: PetscInt i;
2139: Mat_SeqAIJ *x = (Mat_SeqAIJ *)X->data,*y = (Mat_SeqAIJ *)Y->data;
2140: PetscBLASInt one=1,bnz = (PetscBLASInt)x->nz;
2143: if (str == SAME_NONZERO_PATTERN) {
2144: PetscScalar alpha = a;
2145: BLASaxpy_(&bnz,&alpha,x->a,&one,y->a,&one);
2146: } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
2147: if (y->xtoy && y->XtoY != X) {
2148: PetscFree(y->xtoy);
2149: MatDestroy(y->XtoY);
2150: }
2151: if (!y->xtoy) { /* get xtoy */
2152: MatAXPYGetxtoy_Private(X->rmap.n,x->i,x->j,PETSC_NULL, y->i,y->j,PETSC_NULL, &y->xtoy);
2153: y->XtoY = X;
2154: }
2155: for (i=0; i<x->nz; i++) y->a[y->xtoy[i]] += a*(x->a[i]);
2156: PetscInfo3(0,"ratio of nnz(X)/nnz(Y): %d/%d = %G\n",x->nz,y->nz,(PetscReal)(x->nz)/y->nz);
2157: } else {
2158: MatAXPY_Basic(Y,a,X,str);
2159: }
2160: return(0);
2161: }
2165: PetscErrorCode MatSetBlockSize_SeqAIJ(Mat A,PetscInt bs)
2166: {
2168: return(0);
2169: }
2173: PetscErrorCode PETSCMAT_DLLEXPORT MatConjugate_SeqAIJ(Mat mat)
2174: {
2175: #if defined(PETSC_USE_COMPLEX)
2176: Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
2177: PetscInt i,nz;
2178: PetscScalar *a;
2181: nz = aij->nz;
2182: a = aij->a;
2183: for (i=0; i<nz; i++) {
2184: a[i] = PetscConj(a[i]);
2185: }
2186: #else
2188: #endif
2189: return(0);
2190: }
2194: PetscErrorCode MatGetRowMax_SeqAIJ(Mat A,Vec v)
2195: {
2196: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
2198: PetscInt i,j,m = A->rmap.n,*ai,*aj,ncols,n;
2199: PetscReal atmp;
2200: PetscScalar *x,zero = 0.0;
2201: MatScalar *aa;
2204: if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2205: aa = a->a;
2206: ai = a->i;
2207: aj = a->j;
2209: VecSet(v,zero);
2210: VecGetArray(v,&x);
2211: VecGetLocalSize(v,&n);
2212: if (n != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
2213: for (i=0; i<m; i++) {
2214: ncols = ai[1] - ai[0]; ai++;
2215: for (j=0; j<ncols; j++){
2216: atmp = PetscAbsScalar(*aa); aa++;
2217: if (PetscAbsScalar(x[i]) < atmp) x[i] = atmp;
2218: aj++;
2219: }
2220: }
2221: VecRestoreArray(v,&x);
2222: return(0);
2223: }
2225: /* -------------------------------------------------------------------*/
2226: static struct _MatOps MatOps_Values = {MatSetValues_SeqAIJ,
2227: MatGetRow_SeqAIJ,
2228: MatRestoreRow_SeqAIJ,
2229: MatMult_SeqAIJ,
2230: /* 4*/ MatMultAdd_SeqAIJ,
2231: MatMultTranspose_SeqAIJ,
2232: MatMultTransposeAdd_SeqAIJ,
2233: MatSolve_SeqAIJ,
2234: MatSolveAdd_SeqAIJ,
2235: MatSolveTranspose_SeqAIJ,
2236: /*10*/ MatSolveTransposeAdd_SeqAIJ,
2237: MatLUFactor_SeqAIJ,
2238: 0,
2239: MatRelax_SeqAIJ,
2240: MatTranspose_SeqAIJ,
2241: /*15*/ MatGetInfo_SeqAIJ,
2242: MatEqual_SeqAIJ,
2243: MatGetDiagonal_SeqAIJ,
2244: MatDiagonalScale_SeqAIJ,
2245: MatNorm_SeqAIJ,
2246: /*20*/ 0,
2247: MatAssemblyEnd_SeqAIJ,
2248: MatCompress_SeqAIJ,
2249: MatSetOption_SeqAIJ,
2250: MatZeroEntries_SeqAIJ,
2251: /*25*/ MatZeroRows_SeqAIJ,
2252: MatLUFactorSymbolic_SeqAIJ,
2253: MatLUFactorNumeric_SeqAIJ,
2254: MatCholeskyFactorSymbolic_SeqAIJ,
2255: MatCholeskyFactorNumeric_SeqAIJ,
2256: /*30*/ MatSetUpPreallocation_SeqAIJ,
2257: MatILUFactorSymbolic_SeqAIJ,
2258: MatICCFactorSymbolic_SeqAIJ,
2259: MatGetArray_SeqAIJ,
2260: MatRestoreArray_SeqAIJ,
2261: /*35*/ MatDuplicate_SeqAIJ,
2262: 0,
2263: 0,
2264: MatILUFactor_SeqAIJ,
2265: 0,
2266: /*40*/ MatAXPY_SeqAIJ,
2267: MatGetSubMatrices_SeqAIJ,
2268: MatIncreaseOverlap_SeqAIJ,
2269: MatGetValues_SeqAIJ,
2270: MatCopy_SeqAIJ,
2271: /*45*/ MatPrintHelp_SeqAIJ,
2272: MatScale_SeqAIJ,
2273: 0,
2274: MatDiagonalSet_SeqAIJ,
2275: MatILUDTFactor_SeqAIJ,
2276: /*50*/ MatSetBlockSize_SeqAIJ,
2277: MatGetRowIJ_SeqAIJ,
2278: MatRestoreRowIJ_SeqAIJ,
2279: MatGetColumnIJ_SeqAIJ,
2280: MatRestoreColumnIJ_SeqAIJ,
2281: /*55*/ MatFDColoringCreate_SeqAIJ,
2282: 0,
2283: 0,
2284: MatPermute_SeqAIJ,
2285: 0,
2286: /*60*/ 0,
2287: MatDestroy_SeqAIJ,
2288: MatView_SeqAIJ,
2289: 0,
2290: 0,
2291: /*65*/ 0,
2292: 0,
2293: 0,
2294: 0,
2295: 0,
2296: /*70*/ MatGetRowMax_SeqAIJ,
2297: 0,
2298: MatSetColoring_SeqAIJ,
2299: #if defined(PETSC_HAVE_ADIC)
2300: MatSetValuesAdic_SeqAIJ,
2301: #else
2302: 0,
2303: #endif
2304: MatSetValuesAdifor_SeqAIJ,
2305: /*75*/ MatFDColoringApply_SeqAIJ,
2306: 0,
2307: 0,
2308: 0,
2309: 0,
2310: /*80*/ 0,
2311: 0,
2312: 0,
2313: 0,
2314: MatLoad_SeqAIJ,
2315: /*85*/ MatIsSymmetric_SeqAIJ,
2316: 0,
2317: 0,
2318: 0,
2319: 0,
2320: /*90*/ MatMatMult_SeqAIJ_SeqAIJ,
2321: MatMatMultSymbolic_SeqAIJ_SeqAIJ,
2322: MatMatMultNumeric_SeqAIJ_SeqAIJ,
2323: MatPtAP_Basic,
2324: MatPtAPSymbolic_SeqAIJ,
2325: /*95*/ MatPtAPNumeric_SeqAIJ,
2326: MatMatMultTranspose_SeqAIJ_SeqAIJ,
2327: MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ,
2328: MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ,
2329: MatPtAPSymbolic_SeqAIJ_SeqAIJ,
2330: /*100*/MatPtAPNumeric_SeqAIJ_SeqAIJ,
2331: 0,
2332: 0,
2333: MatConjugate_SeqAIJ,
2334: 0,
2335: /*105*/MatSetValuesRow_SeqAIJ,
2336: MatRealPart_SeqAIJ,
2337: MatImaginaryPart_SeqAIJ,
2338: 0,
2339: 0,
2340: /*110*/MatMatSolve_SeqAIJ
2341: };
2346: PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetColumnIndices_SeqAIJ(Mat mat,PetscInt *indices)
2347: {
2348: Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
2349: PetscInt i,nz,n;
2353: nz = aij->maxnz;
2354: n = mat->cmap.n;
2355: for (i=0; i<nz; i++) {
2356: aij->j[i] = indices[i];
2357: }
2358: aij->nz = nz;
2359: for (i=0; i<n; i++) {
2360: aij->ilen[i] = aij->imax[i];
2361: }
2363: return(0);
2364: }
2369: /*@
2370: MatSeqAIJSetColumnIndices - Set the column indices for all the rows
2371: in the matrix.
2373: Input Parameters:
2374: + mat - the SeqAIJ matrix
2375: - indices - the column indices
2377: Level: advanced
2379: Notes:
2380: This can be called if you have precomputed the nonzero structure of the
2381: matrix and want to provide it to the matrix object to improve the performance
2382: of the MatSetValues() operation.
2384: You MUST have set the correct numbers of nonzeros per row in the call to
2385: MatCreateSeqAIJ(), and the columns indices MUST be sorted.
2387: MUST be called before any calls to MatSetValues();
2389: The indices should start with zero, not one.
2391: @*/
2392: PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetColumnIndices(Mat mat,PetscInt *indices)
2393: {
2394: PetscErrorCode ierr,(*f)(Mat,PetscInt *);
2399: PetscObjectQueryFunction((PetscObject)mat,"MatSeqAIJSetColumnIndices_C",(void (**)(void))&f);
2400: if (f) {
2401: (*f)(mat,indices);
2402: } else {
2403: SETERRQ(PETSC_ERR_SUP,"Wrong type of matrix to set column indices");
2404: }
2405: return(0);
2406: }
2408: /* ----------------------------------------------------------------------------------------*/
2413: PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues_SeqAIJ(Mat mat)
2414: {
2415: Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
2417: size_t nz = aij->i[mat->rmap.n];
2420: if (aij->nonew != 1) {
2421: SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
2422: }
2424: /* allocate space for values if not already there */
2425: if (!aij->saved_values) {
2426: PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);
2427: }
2429: /* copy values over */
2430: PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));
2431: return(0);
2432: }
2437: /*@
2438: MatStoreValues - Stashes a copy of the matrix values; this allows, for
2439: example, reuse of the linear part of a Jacobian, while recomputing the
2440: nonlinear portion.
2442: Collect on Mat
2444: Input Parameters:
2445: . mat - the matrix (currently only AIJ matrices support this option)
2447: Level: advanced
2449: Common Usage, with SNESSolve():
2450: $ Create Jacobian matrix
2451: $ Set linear terms into matrix
2452: $ Apply boundary conditions to matrix, at this time matrix must have
2453: $ final nonzero structure (i.e. setting the nonlinear terms and applying
2454: $ boundary conditions again will not change the nonzero structure
2455: $ MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS);
2456: $ MatStoreValues(mat);
2457: $ Call SNESSetJacobian() with matrix
2458: $ In your Jacobian routine
2459: $ MatRetrieveValues(mat);
2460: $ Set nonlinear terms in matrix
2461:
2462: Common Usage without SNESSolve(), i.e. when you handle nonlinear solve yourself:
2463: $ // build linear portion of Jacobian
2464: $ MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS);
2465: $ MatStoreValues(mat);
2466: $ loop over nonlinear iterations
2467: $ MatRetrieveValues(mat);
2468: $ // call MatSetValues(mat,...) to set nonliner portion of Jacobian
2469: $ // call MatAssemblyBegin/End() on matrix
2470: $ Solve linear system with Jacobian
2471: $ endloop
2473: Notes:
2474: Matrix must already be assemblied before calling this routine
2475: Must set the matrix option MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS); before
2476: calling this routine.
2478: When this is called multiple times it overwrites the previous set of stored values
2479: and does not allocated additional space.
2481: .seealso: MatRetrieveValues()
2483: @*/
2484: PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues(Mat mat)
2485: {
2486: PetscErrorCode ierr,(*f)(Mat);
2490: if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2491: if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2493: PetscObjectQueryFunction((PetscObject)mat,"MatStoreValues_C",(void (**)(void))&f);
2494: if (f) {
2495: (*f)(mat);
2496: } else {
2497: SETERRQ(PETSC_ERR_SUP,"Wrong type of matrix to store values");
2498: }
2499: return(0);
2500: }
2505: PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues_SeqAIJ(Mat mat)
2506: {
2507: Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
2509: PetscInt nz = aij->i[mat->rmap.n];
2512: if (aij->nonew != 1) {
2513: SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
2514: }
2515: if (!aij->saved_values) {
2516: SETERRQ(PETSC_ERR_ORDER,"Must call MatStoreValues(A);first");
2517: }
2518: /* copy values over */
2519: PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));
2520: return(0);
2521: }
2526: /*@
2527: MatRetrieveValues - Retrieves the copy of the matrix values; this allows, for
2528: example, reuse of the linear part of a Jacobian, while recomputing the
2529: nonlinear portion.
2531: Collect on Mat
2533: Input Parameters:
2534: . mat - the matrix (currently on AIJ matrices support this option)
2536: Level: advanced
2538: .seealso: MatStoreValues()
2540: @*/
2541: PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues(Mat mat)
2542: {
2543: PetscErrorCode ierr,(*f)(Mat);
2547: if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2548: if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2550: PetscObjectQueryFunction((PetscObject)mat,"MatRetrieveValues_C",(void (**)(void))&f);
2551: if (f) {
2552: (*f)(mat);
2553: } else {
2554: SETERRQ(PETSC_ERR_SUP,"Wrong type of matrix to retrieve values");
2555: }
2556: return(0);
2557: }
2560: /* --------------------------------------------------------------------------------*/
2563: /*@C
2564: MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format
2565: (the default parallel PETSc format). For good matrix assembly performance
2566: the user should preallocate the matrix storage by setting the parameter nz
2567: (or the array nnz). By setting these parameters accurately, performance
2568: during matrix assembly can be increased by more than a factor of 50.
2570: Collective on MPI_Comm
2572: Input Parameters:
2573: + comm - MPI communicator, set to PETSC_COMM_SELF
2574: . m - number of rows
2575: . n - number of columns
2576: . nz - number of nonzeros per row (same for all rows)
2577: - nnz - array containing the number of nonzeros in the various rows
2578: (possibly different for each row) or PETSC_NULL
2580: Output Parameter:
2581: . A - the matrix
2583: Notes:
2584: If nnz is given then nz is ignored
2586: The AIJ format (also called the Yale sparse matrix format or
2587: compressed row storage), is fully compatible with standard Fortran 77
2588: storage. That is, the stored row and column indices can begin at
2589: either one (as in Fortran) or zero. See the users' manual for details.
2591: Specify the preallocated storage with either nz or nnz (not both).
2592: Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
2593: allocation. For large problems you MUST preallocate memory or you
2594: will get TERRIBLE performance, see the users' manual chapter on matrices.
2596: By default, this format uses inodes (identical nodes) when possible, to
2597: improve numerical efficiency of matrix-vector products and solves. We
2598: search for consecutive rows with the same nonzero structure, thereby
2599: reusing matrix information to achieve increased efficiency.
2601: Options Database Keys:
2602: + -mat_no_inode - Do not use inodes
2603: . -mat_inode_limit <limit> - Sets inode limit (max limit=5)
2604: - -mat_aij_oneindex - Internally use indexing starting at 1
2605: rather than 0. Note that when calling MatSetValues(),
2606: the user still MUST index entries starting at 0!
2608: Level: intermediate
2610: .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues(), MatSeqAIJSetColumnIndices(), MatCreateSeqAIJWithArrays()
2612: @*/
2613: PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJ(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A)
2614: {
2618: MatCreate(comm,A);
2619: MatSetSizes(*A,m,n,m,n);
2620: MatSetType(*A,MATSEQAIJ);
2621: MatSeqAIJSetPreallocation_SeqAIJ(*A,nz,(PetscInt*)nnz);
2622: return(0);
2623: }
2627: /*@C
2628: MatSeqAIJSetPreallocation - For good matrix assembly performance
2629: the user should preallocate the matrix storage by setting the parameter nz
2630: (or the array nnz). By setting these parameters accurately, performance
2631: during matrix assembly can be increased by more than a factor of 50.
2633: Collective on MPI_Comm
2635: Input Parameters:
2636: + B - The matrix
2637: . nz - number of nonzeros per row (same for all rows)
2638: - nnz - array containing the number of nonzeros in the various rows
2639: (possibly different for each row) or PETSC_NULL
2641: Notes:
2642: If nnz is given then nz is ignored
2644: The AIJ format (also called the Yale sparse matrix format or
2645: compressed row storage), is fully compatible with standard Fortran 77
2646: storage. That is, the stored row and column indices can begin at
2647: either one (as in Fortran) or zero. See the users' manual for details.
2649: Specify the preallocated storage with either nz or nnz (not both).
2650: Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
2651: allocation. For large problems you MUST preallocate memory or you
2652: will get TERRIBLE performance, see the users' manual chapter on matrices.
2654: Developers: Use nz of MAT_SKIP_ALLOCATION to not allocate any space for the matrix
2655: entries or columns indices
2657: By default, this format uses inodes (identical nodes) when possible, to
2658: improve numerical efficiency of matrix-vector products and solves. We
2659: search for consecutive rows with the same nonzero structure, thereby
2660: reusing matrix information to achieve increased efficiency.
2662: Options Database Keys:
2663: + -mat_no_inode - Do not use inodes
2664: . -mat_inode_limit <limit> - Sets inode limit (max limit=5)
2665: - -mat_aij_oneindex - Internally use indexing starting at 1
2666: rather than 0. Note that when calling MatSetValues(),
2667: the user still MUST index entries starting at 0!
2669: Level: intermediate
2671: .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues(), MatSeqAIJSetColumnIndices(), MatCreateSeqAIJWithArrays()
2673: @*/
2674: PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocation(Mat B,PetscInt nz,const PetscInt nnz[])
2675: {
2676: PetscErrorCode ierr,(*f)(Mat,PetscInt,const PetscInt[]);
2679: PetscObjectQueryFunction((PetscObject)B,"MatSeqAIJSetPreallocation_C",(void (**)(void))&f);
2680: if (f) {
2681: (*f)(B,nz,nnz);
2682: }
2683: return(0);
2684: }
2689: PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocation_SeqAIJ(Mat B,PetscInt nz,PetscInt *nnz)
2690: {
2691: Mat_SeqAIJ *b;
2692: PetscTruth skipallocation = PETSC_FALSE;
2694: PetscInt i;
2697:
2698: if (nz == MAT_SKIP_ALLOCATION) {
2699: skipallocation = PETSC_TRUE;
2700: nz = 0;
2701: }
2703: B->rmap.bs = B->cmap.bs = 1;
2704: PetscMapInitialize(B->comm,&B->rmap);
2705: PetscMapInitialize(B->comm,&B->cmap);
2707: if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
2708: if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz);
2709: if (nnz) {
2710: for (i=0; i<B->rmap.n; i++) {
2711: if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]);
2712: if (nnz[i] > B->cmap.n) SETERRQ3(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be greater than row length: local row %d value %d rowlength %d",i,nnz[i],B->cmap.n);
2713: }
2714: }
2716: B->preallocated = PETSC_TRUE;
2717: b = (Mat_SeqAIJ*)B->data;
2719: if (!skipallocation) {
2720: PetscMalloc2(B->rmap.n,PetscInt,&b->imax,B->rmap.n,PetscInt,&b->ilen);
2721: if (!nnz) {
2722: if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 10;
2723: else if (nz <= 0) nz = 1;
2724: for (i=0; i<B->rmap.n; i++) b->imax[i] = nz;
2725: nz = nz*B->rmap.n;
2726: } else {
2727: nz = 0;
2728: for (i=0; i<B->rmap.n; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
2729: }
2731: /* b->ilen will count nonzeros in each row so far. */
2732: for (i=0; i<B->rmap.n; i++) { b->ilen[i] = 0;}
2734: /* allocate the matrix space */
2735: PetscMalloc3(nz,PetscScalar,&b->a,nz,PetscInt,&b->j,B->rmap.n+1,PetscInt,&b->i);
2736: b->i[0] = 0;
2737: for (i=1; i<B->rmap.n+1; i++) {
2738: b->i[i] = b->i[i-1] + b->imax[i-1];
2739: }
2740: b->singlemalloc = PETSC_TRUE;
2741: b->freedata = PETSC_TRUE;
2742: } else {
2743: b->freedata = PETSC_FALSE;
2744: }
2746: b->nz = 0;
2747: b->maxnz = nz;
2748: B->info.nz_unneeded = (double)b->maxnz;
2749: return(0);
2750: }
2753: #undef __FUNCT__
2755: /*@C
2756: MatSeqAIJSetPreallocationCSR - Allocates memory for a sparse sequential matrix in AIJ format.
2758: Input Parameters:
2759: + B - the matrix
2760: . i - the indices into j for the start of each row (starts with zero)
2761: . j - the column indices for each row (starts with zero) these must be sorted for each row
2762: - v - optional values in the matrix
2764: Contributed by: Lisandro Dalchin
2766: Level: developer
2768: .keywords: matrix, aij, compressed row, sparse, sequential
2770: .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatSeqAIJSetPreallocation(), MatCreateSeqAIJ(), SeqAIJ
2771: @*/
2772: PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat B,const PetscInt i[],const PetscInt j[],const PetscScalar v[])
2773: {
2774: PetscErrorCode (*f)(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
2779: PetscObjectQueryFunction((PetscObject)B,"MatSeqAIJSetPreallocationCSR_C",(void (**)(void))&f);
2780: if (f) {
2781: (*f)(B,i,j,v);
2782: }
2783: return(0);
2784: }
2787: #undef __FUNCT__
2789: PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocationCSR_SeqAIJ(Mat B,const PetscInt I[],const PetscInt J[],const PetscScalar v[])
2790: {
2791: PetscInt i;
2792: PetscInt m,n;
2793: PetscInt nz;
2794: PetscInt *nnz, nz_max = 0;
2795: PetscScalar *values;
2799: MatGetSize(B, &m, &n);
2801: if (I[0]) {
2802: SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "I[0] must be 0 it is %D", I[0]);
2803: }
2804: PetscMalloc((m+1) * sizeof(PetscInt), &nnz);
2805: for(i = 0; i < m; i++) {
2806: nz = I[i+1]- I[i];
2807: nz_max = PetscMax(nz_max, nz);
2808: if (nz < 0) {
2809: SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE, "Local row %D has a negative number of columns %D", i, nnz);
2810: }
2811: nnz[i] = nz;
2812: }
2813: MatSeqAIJSetPreallocation(B, 0, nnz);
2814: PetscFree(nnz);
2816: if (v) {
2817: values = (PetscScalar*) v;
2818: } else {
2819: PetscMalloc((nz_max+1)*sizeof(PetscScalar), &values);
2820: PetscMemzero(values, nz_max*sizeof(PetscScalar));
2821: }
2823: MatSetOption(B,MAT_COLUMNS_SORTED);
2825: for(i = 0; i < m; i++) {
2826: nz = I[i+1] - I[i];
2827: MatSetValues_SeqAIJ(B, 1, &i, nz, J+I[i], values + (v ? I[i] : 0), INSERT_VALUES);
2828: }
2830: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
2831: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
2832: MatSetOption(B,MAT_COLUMNS_UNSORTED);
2834: if (!v) {
2835: PetscFree(values);
2836: }
2837: return(0);
2838: }
2841: /*MC
2842: MATSEQAIJ - MATSEQAIJ = "seqaij" - A matrix type to be used for sequential sparse matrices,
2843: based on compressed sparse row format.
2845: Options Database Keys:
2846: . -mat_type seqaij - sets the matrix type to "seqaij" during a call to MatSetFromOptions()
2848: Level: beginner
2850: .seealso: MatCreateSeqAIJ(), MatSetFromOptions(), MatSetType(), MatCreate(), MatType
2851: M*/
2856: PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqAIJ(Mat B)
2857: {
2858: Mat_SeqAIJ *b;
2860: PetscMPIInt size;
2863: MPI_Comm_size(B->comm,&size);
2864: if (size > 1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Comm must be of size 1");
2866: PetscNew(Mat_SeqAIJ,&b);
2867: B->data = (void*)b;
2868: PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
2869: B->factor = 0;
2870: B->mapping = 0;
2871: b->row = 0;
2872: b->col = 0;
2873: b->icol = 0;
2874: b->reallocs = 0;
2875: b->sorted = PETSC_FALSE;
2876: b->ignorezeroentries = PETSC_FALSE;
2877: b->roworiented = PETSC_TRUE;
2878: b->nonew = 0;
2879: b->diag = 0;
2880: b->solve_work = 0;
2881: B->spptr = 0;
2882: b->saved_values = 0;
2883: b->idiag = 0;
2884: b->ssor = 0;
2885: b->keepzeroedrows = PETSC_FALSE;
2886: b->xtoy = 0;
2887: b->XtoY = 0;
2888: b->compressedrow.use = PETSC_FALSE;
2889: b->compressedrow.nrows = B->rmap.n;
2890: b->compressedrow.i = PETSC_NULL;
2891: b->compressedrow.rindex = PETSC_NULL;
2892: b->compressedrow.checked = PETSC_FALSE;
2893: B->same_nonzero = PETSC_FALSE;
2895: PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);
2897: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqAIJSetColumnIndices_C",
2898: "MatSeqAIJSetColumnIndices_SeqAIJ",
2899: MatSeqAIJSetColumnIndices_SeqAIJ);
2900: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
2901: "MatStoreValues_SeqAIJ",
2902: MatStoreValues_SeqAIJ);
2903: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
2904: "MatRetrieveValues_SeqAIJ",
2905: MatRetrieveValues_SeqAIJ);
2906: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_seqsbaij_C",
2907: "MatConvert_SeqAIJ_SeqSBAIJ",
2908: MatConvert_SeqAIJ_SeqSBAIJ);
2909: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_seqbaij_C",
2910: "MatConvert_SeqAIJ_SeqBAIJ",
2911: MatConvert_SeqAIJ_SeqBAIJ);
2912: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatIsTranspose_C",
2913: "MatIsTranspose_SeqAIJ",
2914: MatIsTranspose_SeqAIJ);
2915: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqAIJSetPreallocation_C",
2916: "MatSeqAIJSetPreallocation_SeqAIJ",
2917: MatSeqAIJSetPreallocation_SeqAIJ);
2918: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqAIJSetPreallocationCSR_C",
2919: "MatSeqAIJSetPreallocationCSR_SeqAIJ",
2920: MatSeqAIJSetPreallocationCSR_SeqAIJ);
2921: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatReorderForNonzeroDiagonal_C",
2922: "MatReorderForNonzeroDiagonal_SeqAIJ",
2923: MatReorderForNonzeroDiagonal_SeqAIJ);
2924: MatCreate_Inode(B);
2925: return(0);
2926: }
2931: PetscErrorCode MatDuplicate_SeqAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
2932: {
2933: Mat C;
2934: Mat_SeqAIJ *c,*a = (Mat_SeqAIJ*)A->data;
2936: PetscInt i,m = A->rmap.n;
2939: *B = 0;
2940: MatCreate(A->comm,&C);
2941: MatSetSizes(C,A->rmap.n,A->cmap.n,A->rmap.n,A->cmap.n);
2942: MatSetType(C,A->type_name);
2943: PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));
2944:
2945: PetscMapCopy(A->comm,&A->rmap,&C->rmap);
2946: PetscMapCopy(A->comm,&A->cmap,&C->cmap);
2948: c = (Mat_SeqAIJ*)C->data;
2950: C->factor = A->factor;
2952: c->row = 0;
2953: c->col = 0;
2954: c->icol = 0;
2955: c->reallocs = 0;
2957: C->assembled = PETSC_TRUE;
2958:
2959: PetscMalloc2(m,PetscInt,&c->imax,m,PetscInt,&c->ilen);
2960: for (i=0; i<m; i++) {
2961: c->imax[i] = a->imax[i];
2962: c->ilen[i] = a->ilen[i];
2963: }
2965: /* allocate the matrix space */
2966: PetscMalloc3(a->i[m],PetscScalar,&c->a,a->i[m],PetscInt,&c->j,m+1,PetscInt,&c->i);
2967: c->singlemalloc = PETSC_TRUE;
2968: PetscMemcpy(c->i,a->i,(m+1)*sizeof(PetscInt));
2969: if (m > 0) {
2970: PetscMemcpy(c->j,a->j,(a->i[m])*sizeof(PetscInt));
2971: if (cpvalues == MAT_COPY_VALUES) {
2972: PetscMemcpy(c->a,a->a,(a->i[m])*sizeof(PetscScalar));
2973: } else {
2974: PetscMemzero(c->a,(a->i[m])*sizeof(PetscScalar));
2975: }
2976: }
2978: c->sorted = a->sorted;
2979: c->ignorezeroentries = a->ignorezeroentries;
2980: c->roworiented = a->roworiented;
2981: c->nonew = a->nonew;
2982: if (a->diag) {
2983: PetscMalloc((m+1)*sizeof(PetscInt),&c->diag);
2984: PetscLogObjectMemory(C,(m+1)*sizeof(PetscInt));
2985: for (i=0; i<m; i++) {
2986: c->diag[i] = a->diag[i];
2987: }
2988: } else c->diag = 0;
2989: c->solve_work = 0;
2990: c->saved_values = 0;
2991: c->idiag = 0;
2992: c->ssor = 0;
2993: c->keepzeroedrows = a->keepzeroedrows;
2994: c->freedata = PETSC_TRUE;
2995: c->xtoy = 0;
2996: c->XtoY = 0;
2998: c->nz = a->nz;
2999: c->maxnz = a->maxnz;
3000: C->preallocated = PETSC_TRUE;
3002: c->compressedrow.use = a->compressedrow.use;
3003: c->compressedrow.nrows = a->compressedrow.nrows;
3004: c->compressedrow.checked = a->compressedrow.checked;
3005: if ( a->compressedrow.checked && a->compressedrow.use){
3006: i = a->compressedrow.nrows;
3007: PetscMalloc((2*i+1)*sizeof(PetscInt),&c->compressedrow.i);
3008: c->compressedrow.rindex = c->compressedrow.i + i + 1;
3009: PetscMemcpy(c->compressedrow.i,a->compressedrow.i,(i+1)*sizeof(PetscInt));
3010: PetscMemcpy(c->compressedrow.rindex,a->compressedrow.rindex,i*sizeof(PetscInt));
3011: } else {
3012: c->compressedrow.use = PETSC_FALSE;
3013: c->compressedrow.i = PETSC_NULL;
3014: c->compressedrow.rindex = PETSC_NULL;
3015: }
3016: C->same_nonzero = A->same_nonzero;
3017: MatDuplicate_Inode(A,cpvalues,&C);
3019: *B = C;
3020: PetscFListDuplicate(A->qlist,&C->qlist);
3021: return(0);
3022: }
3026: PetscErrorCode MatLoad_SeqAIJ(PetscViewer viewer, MatType type,Mat *A)
3027: {
3028: Mat_SeqAIJ *a;
3029: Mat B;
3031: PetscInt i,sum,nz,header[4],*rowlengths = 0,M,N;
3032: int fd;
3033: PetscMPIInt size;
3034: MPI_Comm comm;
3035:
3037: PetscObjectGetComm((PetscObject)viewer,&comm);
3038: MPI_Comm_size(comm,&size);
3039: if (size > 1) SETERRQ(PETSC_ERR_ARG_SIZ,"view must have one processor");
3040: PetscViewerBinaryGetDescriptor(viewer,&fd);
3041: PetscBinaryRead(fd,header,4,PETSC_INT);
3042: if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object in file");
3043: M = header[1]; N = header[2]; nz = header[3];
3045: if (nz < 0) {
3046: SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format on disk,cannot load as SeqAIJ");
3047: }
3049: /* read in row lengths */
3050: PetscMalloc(M*sizeof(PetscInt),&rowlengths);
3051: PetscBinaryRead(fd,rowlengths,M,PETSC_INT);
3053: /* check if sum of rowlengths is same as nz */
3054: for (i=0,sum=0; i< M; i++) sum +=rowlengths[i];
3055: if (sum != nz) SETERRQ2(PETSC_ERR_FILE_READ,"Inconsistant matrix data in file. no-nonzeros = %d, sum-row-lengths = %d\n",nz,sum);
3057: /* create our matrix */
3058: MatCreate(comm,&B);
3059: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,M,N);
3060: MatSetType(B,type);
3061: MatSeqAIJSetPreallocation_SeqAIJ(B,0,rowlengths);
3062: a = (Mat_SeqAIJ*)B->data;
3064: /* read in column indices and adjust for Fortran indexing*/
3065: PetscBinaryRead(fd,a->j,nz,PETSC_INT);
3067: /* read in nonzero values */
3068: PetscBinaryRead(fd,a->a,nz,PETSC_SCALAR);
3070: /* set matrix "i" values */
3071: a->i[0] = 0;
3072: for (i=1; i<= M; i++) {
3073: a->i[i] = a->i[i-1] + rowlengths[i-1];
3074: a->ilen[i-1] = rowlengths[i-1];
3075: }
3076: PetscFree(rowlengths);
3078: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
3079: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
3080: *A = B;
3081: return(0);
3082: }
3086: PetscErrorCode MatEqual_SeqAIJ(Mat A,Mat B,PetscTruth* flg)
3087: {
3088: Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data,*b = (Mat_SeqAIJ *)B->data;
3092: /* If the matrix dimensions are not equal,or no of nonzeros */
3093: if ((A->rmap.n != B->rmap.n) || (A->cmap.n != B->cmap.n) ||(a->nz != b->nz)) {
3094: *flg = PETSC_FALSE;
3095: return(0);
3096: }
3097:
3098: /* if the a->i are the same */
3099: PetscMemcmp(a->i,b->i,(A->rmap.n+1)*sizeof(PetscInt),flg);
3100: if (!*flg) return(0);
3101:
3102: /* if a->j are the same */
3103: PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),flg);
3104: if (!*flg) return(0);
3105:
3106: /* if a->a are the same */
3107: PetscMemcmp(a->a,b->a,(a->nz)*sizeof(PetscScalar),flg);
3109: return(0);
3110:
3111: }
3115: /*@
3116: MatCreateSeqAIJWithArrays - Creates an sequential AIJ matrix using matrix elements (in CSR format)
3117: provided by the user.
3119: Coolective on MPI_Comm
3121: Input Parameters:
3122: + comm - must be an MPI communicator of size 1
3123: . m - number of rows
3124: . n - number of columns
3125: . i - row indices
3126: . j - column indices
3127: - a - matrix values
3129: Output Parameter:
3130: . mat - the matrix
3132: Level: intermediate
3134: Notes:
3135: The i, j, and a arrays are not copied by this routine, the user must free these arrays
3136: once the matrix is destroyed
3138: You cannot set new nonzero locations into this matrix, that will generate an error.
3140: The i and j indices are 0 based
3142: .seealso: MatCreate(), MatCreateMPIAIJ(), MatCreateSeqAIJ()
3144: @*/
3145: PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJWithArrays(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt* i,PetscInt*j,PetscScalar *a,Mat *mat)
3146: {
3148: PetscInt ii;
3149: Mat_SeqAIJ *aij;
3152: if (i[0]) {
3153: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"i (row indices) must start with 0");
3154: }
3155: MatCreate(comm,mat);
3156: MatSetSizes(*mat,m,n,m,n);
3157: MatSetType(*mat,MATSEQAIJ);
3158: MatSeqAIJSetPreallocation_SeqAIJ(*mat,MAT_SKIP_ALLOCATION,0);
3159: aij = (Mat_SeqAIJ*)(*mat)->data;
3160: PetscMalloc2(m,PetscInt,&aij->imax,m,PetscInt,&aij->ilen);
3162: aij->i = i;
3163: aij->j = j;
3164: aij->a = a;
3165: aij->singlemalloc = PETSC_FALSE;
3166: aij->nonew = -1; /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/
3167: aij->freedata = PETSC_FALSE;
3169: for (ii=0; ii<m; ii++) {
3170: aij->ilen[ii] = aij->imax[ii] = i[ii+1] - i[ii];
3171: #if defined(PETSC_USE_DEBUG)
3172: if (i[ii+1] - i[ii] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Negative row length in i (row indices) row = %d length = %d",ii,i[ii+1] - i[ii]);
3173: #endif
3174: }
3175: #if defined(PETSC_USE_DEBUG)
3176: for (ii=0; ii<aij->i[m]; ii++) {
3177: if (j[ii] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Negative column index at location = %d index = %d",ii,j[ii]);
3178: if (j[ii] > n - 1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column index to large at location = %d index = %d",ii,j[ii]);
3179: }
3180: #endif
3182: MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
3183: MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
3184: return(0);
3185: }
3189: PetscErrorCode MatSetColoring_SeqAIJ(Mat A,ISColoring coloring)
3190: {
3192: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
3195: if (coloring->ctype == IS_COLORING_LOCAL) {
3196: ISColoringReference(coloring);
3197: a->coloring = coloring;
3198: } else if (coloring->ctype == IS_COLORING_GHOSTED) {
3199: PetscInt i,*larray;
3200: ISColoring ocoloring;
3201: ISColoringValue *colors;
3203: /* set coloring for diagonal portion */
3204: PetscMalloc((A->cmap.n+1)*sizeof(PetscInt),&larray);
3205: for (i=0; i<A->cmap.n; i++) {
3206: larray[i] = i;
3207: }
3208: ISGlobalToLocalMappingApply(A->mapping,IS_GTOLM_MASK,A->cmap.n,larray,PETSC_NULL,larray);
3209: PetscMalloc((A->cmap.n+1)*sizeof(ISColoringValue),&colors);
3210: for (i=0; i<A->cmap.n; i++) {
3211: colors[i] = coloring->colors[larray[i]];
3212: }
3213: PetscFree(larray);
3214: ISColoringCreate(PETSC_COMM_SELF,coloring->n,A->cmap.n,colors,&ocoloring);
3215: a->coloring = ocoloring;
3216: }
3217: return(0);
3218: }
3220: #if defined(PETSC_HAVE_ADIC)
3222: #include "adic/ad_utils.h"
3227: PetscErrorCode MatSetValuesAdic_SeqAIJ(Mat A,void *advalues)
3228: {
3229: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
3230: PetscInt m = A->rmap.n,*ii = a->i,*jj = a->j,nz,i,j,nlen;
3231: PetscScalar *v = a->a,*values = ((PetscScalar*)advalues)+1;
3232: ISColoringValue *color;
3235: if (!a->coloring) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Coloring not set for matrix");
3236: nlen = PetscADGetDerivTypeSize()/sizeof(PetscScalar);
3237: color = a->coloring->colors;
3238: /* loop over rows */
3239: for (i=0; i<m; i++) {
3240: nz = ii[i+1] - ii[i];
3241: /* loop over columns putting computed value into matrix */
3242: for (j=0; j<nz; j++) {
3243: *v++ = values[color[*jj++]];
3244: }
3245: values += nlen; /* jump to next row of derivatives */
3246: }
3247: return(0);
3248: }
3249: #endif
3253: PetscErrorCode MatSetValuesAdifor_SeqAIJ(Mat A,PetscInt nl,void *advalues)
3254: {
3255: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
3256: PetscInt m = A->rmap.n,*ii = a->i,*jj = a->j,nz,i,j;
3257: PetscScalar *v = a->a,*values = (PetscScalar *)advalues;
3258: ISColoringValue *color;
3261: if (!a->coloring) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Coloring not set for matrix");
3262: color = a->coloring->colors;
3263: /* loop over rows */
3264: for (i=0; i<m; i++) {
3265: nz = ii[i+1] - ii[i];
3266: /* loop over columns putting computed value into matrix */
3267: for (j=0; j<nz; j++) {
3268: *v++ = values[color[*jj++]];
3269: }
3270: values += nl; /* jump to next row of derivatives */
3271: }
3272: return(0);
3273: }
3275: /*
3276: Special version for direct calls from Fortran
3277: */
3278: #if defined(PETSC_HAVE_FORTRAN_CAPS)
3279: #define matsetvaluesseqaij_ MATSETVALUESSEQAIJ
3280: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
3281: #define matsetvaluesseqaij_ matsetvaluesseqaij
3282: #endif
3284: /* Change these macros so can be used in void function */
3285: #undef CHKERRQ
3286: #define CHKERRQ(ierr) CHKERRABORT(A->comm,ierr)
3287: #undef SETERRQ2
3288: #define SETERRQ2(ierr,b,c,d) CHKERRABORT(A->comm,ierr)
3293: void PETSCMAT_DLLEXPORT matsetvaluesseqaij_(Mat *AA,PetscInt *mm,const PetscInt im[],PetscInt *nn,const PetscInt in[],const PetscScalar v[],InsertMode *isis)
3294: {
3295: Mat A = *AA;
3296: PetscInt m = *mm, n = *nn;
3297: InsertMode is = *isis;
3298: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
3299: PetscInt *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N;
3300: PetscInt *imax,*ai,*ailen;
3302: PetscInt *aj,nonew = a->nonew,lastcol = -1;
3303: PetscScalar *ap,value,*aa;
3304: PetscTruth ignorezeroentries = ((a->ignorezeroentries && is == ADD_VALUES) ? PETSC_TRUE:PETSC_FALSE);
3305: PetscTruth roworiented = a->roworiented;
3308: MatPreallocated(A);
3309: imax = a->imax;
3310: ai = a->i;
3311: ailen = a->ilen;
3312: aj = a->j;
3313: aa = a->a;
3315: for (k=0; k<m; k++) { /* loop over added rows */
3316: row = im[k];
3317: if (row < 0) continue;
3318: #if defined(PETSC_USE_DEBUG)
3319: if (row >= A->rmap.n) SETERRABORT(A->comm,PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
3320: #endif
3321: rp = aj + ai[row]; ap = aa + ai[row];
3322: rmax = imax[row]; nrow = ailen[row];
3323: low = 0;
3324: high = nrow;
3325: for (l=0; l<n; l++) { /* loop over added columns */
3326: if (in[l] < 0) continue;
3327: #if defined(PETSC_USE_DEBUG)
3328: if (in[l] >= A->cmap.n) SETERRABORT(A->comm,PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
3329: #endif
3330: col = in[l];
3331: if (roworiented) {
3332: value = v[l + k*n];
3333: } else {
3334: value = v[k + l*m];
3335: }
3336: if (value == 0.0 && ignorezeroentries && (is == ADD_VALUES)) continue;
3338: if (col <= lastcol) low = 0; else high = nrow;
3339: lastcol = col;
3340: while (high-low > 5) {
3341: t = (low+high)/2;
3342: if (rp[t] > col) high = t;
3343: else low = t;
3344: }
3345: for (i=low; i<high; i++) {
3346: if (rp[i] > col) break;
3347: if (rp[i] == col) {
3348: if (is == ADD_VALUES) ap[i] += value;
3349: else ap[i] = value;
3350: goto noinsert;
3351: }
3352: }
3353: if (value == 0.0 && ignorezeroentries) goto noinsert;
3354: if (nonew == 1) goto noinsert;
3355: if (nonew == -1) SETERRABORT(A->comm,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
3356: MatSeqXAIJReallocateAIJ(a,1,nrow,row,col,rmax,aa,ai,aj,A->rmap.n,rp,ap,imax,nonew);
3357: N = nrow++ - 1; a->nz++; high++;
3358: /* shift up all the later entries in this row */
3359: for (ii=N; ii>=i; ii--) {
3360: rp[ii+1] = rp[ii];
3361: ap[ii+1] = ap[ii];
3362: }
3363: rp[i] = col;
3364: ap[i] = value;
3365: noinsert:;
3366: low = i + 1;
3367: }
3368: ailen[row] = nrow;
3369: }
3370: A->same_nonzero = PETSC_FALSE;
3371: PetscFunctionReturnVoid();
3372: }