Actual source code: aij.c
1: #define PETSCMAT_DLL
3: /*
4: Defines the basic matrix operations for the AIJ (compressed row)
5: matrix storage format.
6: */
8: #include src/mat/impls/aij/seq/aij.h
9: #include src/inline/spops.h
10: #include src/inline/dot.h
11: #include petscbt.h
15: PetscErrorCode MatDiagonalSet_SeqAIJ(Mat Y,Vec D,InsertMode is)
16: {
18: Mat_SeqAIJ *aij = (Mat_SeqAIJ*) Y->data;
19: PetscInt i,*diag, m = Y->rmap.n;
20: PetscScalar *v,*aa = aij->a;
21: PetscTruth missing;
24: if (Y->assembled) {
25: MatMissingDiagonal_SeqAIJ(Y,&missing,PETSC_NULL);
26: if (!missing) {
27: diag = aij->diag;
28: VecGetArray(D,&v);
29: if (is == INSERT_VALUES) {
30: for (i=0; i<m; i++) {
31: aa[diag[i]] = v[i];
32: }
33: } else {
34: for (i=0; i<m; i++) {
35: aa[diag[i]] += v[i];
36: }
37: }
38: VecRestoreArray(D,&v);
39: return(0);
40: }
41: }
42: MatDiagonalSet_Default(Y,D,is);
43: return(0);
44: }
48: PetscErrorCode MatGetRowIJ_SeqAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscTruth blockcompressed,PetscInt *m,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
49: {
50: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
52: PetscInt i,ishift;
53:
55: *m = A->rmap.n;
56: if (!ia) return(0);
57: ishift = 0;
58: if (symmetric && !A->structurally_symmetric) {
59: MatToSymmetricIJ_SeqAIJ(A->rmap.n,a->i,a->j,ishift,oshift,ia,ja);
60: } else if (oshift == 1) {
61: PetscInt nz = a->i[A->rmap.n];
62: /* malloc space and add 1 to i and j indices */
63: PetscMalloc((A->rmap.n+1)*sizeof(PetscInt),ia);
64: PetscMalloc((nz+1)*sizeof(PetscInt),ja);
65: for (i=0; i<nz; i++) (*ja)[i] = a->j[i] + 1;
66: for (i=0; i<A->rmap.n+1; i++) (*ia)[i] = a->i[i] + 1;
67: } else {
68: *ia = a->i; *ja = a->j;
69: }
70: return(0);
71: }
75: PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscTruth blockcompressed,PetscInt *n,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
76: {
78:
80: if (!ia) return(0);
81: if ((symmetric && !A->structurally_symmetric) || oshift == 1) {
82: PetscFree(*ia);
83: PetscFree(*ja);
84: }
85: return(0);
86: }
90: PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscTruth blockcompressed,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
91: {
92: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
94: PetscInt i,*collengths,*cia,*cja,n = A->cmap.n,m = A->rmap.n;
95: PetscInt nz = a->i[m],row,*jj,mr,col;
96:
98: *nn = n;
99: if (!ia) return(0);
100: if (symmetric) {
101: MatToSymmetricIJ_SeqAIJ(A->rmap.n,a->i,a->j,0,oshift,ia,ja);
102: } else {
103: PetscMalloc((n+1)*sizeof(PetscInt),&collengths);
104: PetscMemzero(collengths,n*sizeof(PetscInt));
105: PetscMalloc((n+1)*sizeof(PetscInt),&cia);
106: PetscMalloc((nz+1)*sizeof(PetscInt),&cja);
107: jj = a->j;
108: for (i=0; i<nz; i++) {
109: collengths[jj[i]]++;
110: }
111: cia[0] = oshift;
112: for (i=0; i<n; i++) {
113: cia[i+1] = cia[i] + collengths[i];
114: }
115: PetscMemzero(collengths,n*sizeof(PetscInt));
116: jj = a->j;
117: for (row=0; row<m; row++) {
118: mr = a->i[row+1] - a->i[row];
119: for (i=0; i<mr; i++) {
120: col = *jj++;
121: cja[cia[col] + collengths[col]++ - oshift] = row + oshift;
122: }
123: }
124: PetscFree(collengths);
125: *ia = cia; *ja = cja;
126: }
127: return(0);
128: }
132: PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscTruth blockcompressed,PetscInt *n,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
133: {
137: if (!ia) return(0);
139: PetscFree(*ia);
140: PetscFree(*ja);
141:
142: return(0);
143: }
147: PetscErrorCode MatSetValuesRow_SeqAIJ(Mat A,PetscInt row,const PetscScalar v[])
148: {
149: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
150: PetscInt *ai = a->i;
154: PetscMemcpy(a->a+ai[row],v,(ai[row+1]-ai[row])*sizeof(PetscScalar));
155: return(0);
156: }
158: #define CHUNKSIZE 15
162: PetscErrorCode MatSetValues_SeqAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
163: {
164: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
165: PetscInt *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N;
166: PetscInt *imax = a->imax,*ai = a->i,*ailen = a->ilen;
168: PetscInt *aj = a->j,nonew = a->nonew,lastcol = -1;
169: PetscScalar *ap,value,*aa = a->a;
170: PetscTruth ignorezeroentries = a->ignorezeroentries;
171: PetscTruth roworiented = a->roworiented;
174: for (k=0; k<m; k++) { /* loop over added rows */
175: row = im[k];
176: if (row < 0) continue;
177: #if defined(PETSC_USE_DEBUG)
178: if (row >= A->rmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap.n-1);
179: #endif
180: rp = aj + ai[row]; ap = aa + ai[row];
181: rmax = imax[row]; nrow = ailen[row];
182: low = 0;
183: high = nrow;
184: for (l=0; l<n; l++) { /* loop over added columns */
185: if (in[l] < 0) continue;
186: #if defined(PETSC_USE_DEBUG)
187: if (in[l] >= A->cmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->cmap.n-1);
188: #endif
189: col = in[l];
190: if (roworiented) {
191: value = v[l + k*n];
192: } else {
193: value = v[k + l*m];
194: }
195: if (value == 0.0 && ignorezeroentries && (is == ADD_VALUES)) continue;
197: if (col <= lastcol) low = 0; else high = nrow;
198: lastcol = col;
199: while (high-low > 5) {
200: t = (low+high)/2;
201: if (rp[t] > col) high = t;
202: else low = t;
203: }
204: for (i=low; i<high; i++) {
205: if (rp[i] > col) break;
206: if (rp[i] == col) {
207: if (is == ADD_VALUES) ap[i] += value;
208: else ap[i] = value;
209: goto noinsert;
210: }
211: }
212: if (value == 0.0 && ignorezeroentries) goto noinsert;
213: if (nonew == 1) goto noinsert;
214: if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero at (%D,%D) in the matrix",row,col);
215: MatSeqXAIJReallocateAIJ(A,A->rmap.n,1,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar);
216: N = nrow++ - 1; a->nz++; high++;
217: /* shift up all the later entries in this row */
218: for (ii=N; ii>=i; ii--) {
219: rp[ii+1] = rp[ii];
220: ap[ii+1] = ap[ii];
221: }
222: rp[i] = col;
223: ap[i] = value;
224: noinsert:;
225: low = i + 1;
226: }
227: ailen[row] = nrow;
228: }
229: A->same_nonzero = PETSC_FALSE;
230: return(0);
231: }
236: PetscErrorCode MatGetValues_SeqAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],PetscScalar v[])
237: {
238: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
239: PetscInt *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
240: PetscInt *ai = a->i,*ailen = a->ilen;
241: PetscScalar *ap,*aa = a->a,zero = 0.0;
244: for (k=0; k<m; k++) { /* loop over rows */
245: row = im[k];
246: if (row < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",row);
247: if (row >= A->rmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap.n-1);
248: rp = aj + ai[row]; ap = aa + ai[row];
249: nrow = ailen[row];
250: for (l=0; l<n; l++) { /* loop over columns */
251: if (in[l] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %D",in[l]);
252: if (in[l] >= A->cmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->cmap.n-1);
253: col = in[l] ;
254: high = nrow; low = 0; /* assume unsorted */
255: while (high-low > 5) {
256: t = (low+high)/2;
257: if (rp[t] > col) high = t;
258: else low = t;
259: }
260: for (i=low; i<high; i++) {
261: if (rp[i] > col) break;
262: if (rp[i] == col) {
263: *v++ = ap[i];
264: goto finished;
265: }
266: }
267: *v++ = zero;
268: finished:;
269: }
270: }
271: return(0);
272: }
277: PetscErrorCode MatView_SeqAIJ_Binary(Mat A,PetscViewer viewer)
278: {
279: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
281: PetscInt i,*col_lens;
282: int fd;
285: PetscViewerBinaryGetDescriptor(viewer,&fd);
286: PetscMalloc((4+A->rmap.n)*sizeof(PetscInt),&col_lens);
287: col_lens[0] = MAT_FILE_COOKIE;
288: col_lens[1] = A->rmap.n;
289: col_lens[2] = A->cmap.n;
290: col_lens[3] = a->nz;
292: /* store lengths of each row and write (including header) to file */
293: for (i=0; i<A->rmap.n; i++) {
294: col_lens[4+i] = a->i[i+1] - a->i[i];
295: }
296: PetscBinaryWrite(fd,col_lens,4+A->rmap.n,PETSC_INT,PETSC_TRUE);
297: PetscFree(col_lens);
299: /* store column indices (zero start index) */
300: PetscBinaryWrite(fd,a->j,a->nz,PETSC_INT,PETSC_FALSE);
302: /* store nonzero values */
303: PetscBinaryWrite(fd,a->a,a->nz,PETSC_SCALAR,PETSC_FALSE);
304: return(0);
305: }
307: EXTERN PetscErrorCode MatSeqAIJFactorInfo_Matlab(Mat,PetscViewer);
311: PetscErrorCode MatView_SeqAIJ_ASCII(Mat A,PetscViewer viewer)
312: {
313: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
314: PetscErrorCode ierr;
315: PetscInt i,j,m = A->rmap.n,shift=0;
316: const char *name;
317: PetscViewerFormat format;
320: PetscObjectGetName((PetscObject)A,&name);
321: PetscViewerGetFormat(viewer,&format);
322: if (format == PETSC_VIEWER_ASCII_MATLAB) {
323: PetscInt nofinalvalue = 0;
324: if ((a->i[m] == a->i[m-1]) || (a->j[a->nz-1] != A->cmap.n-!shift)) {
325: nofinalvalue = 1;
326: }
327: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
328: PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",m,A->cmap.n);
329: PetscViewerASCIIPrintf(viewer,"%% Nonzeros = %D \n",a->nz);
330: PetscViewerASCIIPrintf(viewer,"zzz = zeros(%D,3);\n",a->nz+nofinalvalue);
331: PetscViewerASCIIPrintf(viewer,"zzz = [\n");
333: for (i=0; i<m; i++) {
334: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
335: #if defined(PETSC_USE_COMPLEX)
336: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e + %18.16ei \n",i+1,a->j[j]+!shift,PetscRealPart(a->a[j]),PetscImaginaryPart(a->a[j]));
337: #else
338: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",i+1,a->j[j]+!shift,a->a[j]);
339: #endif
340: }
341: }
342: if (nofinalvalue) {
343: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",m,A->cmap.n,0.0);
344: }
345: PetscViewerASCIIPrintf(viewer,"];\n %s = spconvert(zzz);\n",name);
346: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
347: } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) {
348: return(0);
349: } else if (format == PETSC_VIEWER_ASCII_COMMON) {
350: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
351: for (i=0; i<m; i++) {
352: PetscViewerASCIIPrintf(viewer,"row %D:",i);
353: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
354: #if defined(PETSC_USE_COMPLEX)
355: if (PetscImaginaryPart(a->a[j]) > 0.0 && PetscRealPart(a->a[j]) != 0.0) {
356: PetscViewerASCIIPrintf(viewer," (%D, %G + %G i)",a->j[j]+shift,PetscRealPart(a->a[j]),PetscImaginaryPart(a->a[j]));
357: } else if (PetscImaginaryPart(a->a[j]) < 0.0 && PetscRealPart(a->a[j]) != 0.0) {
358: PetscViewerASCIIPrintf(viewer," (%D, %G - %G i)",a->j[j]+shift,PetscRealPart(a->a[j]),-PetscImaginaryPart(a->a[j]));
359: } else if (PetscRealPart(a->a[j]) != 0.0) {
360: PetscViewerASCIIPrintf(viewer," (%D, %G) ",a->j[j]+shift,PetscRealPart(a->a[j]));
361: }
362: #else
363: if (a->a[j] != 0.0) {PetscViewerASCIIPrintf(viewer," (%D, %G) ",a->j[j]+shift,a->a[j]);}
364: #endif
365: }
366: PetscViewerASCIIPrintf(viewer,"\n");
367: }
368: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
369: } else if (format == PETSC_VIEWER_ASCII_SYMMODU) {
370: PetscInt nzd=0,fshift=1,*sptr;
371: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
372: PetscMalloc((m+1)*sizeof(PetscInt),&sptr);
373: for (i=0; i<m; i++) {
374: sptr[i] = nzd+1;
375: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
376: if (a->j[j] >= i) {
377: #if defined(PETSC_USE_COMPLEX)
378: if (PetscImaginaryPart(a->a[j]) != 0.0 || PetscRealPart(a->a[j]) != 0.0) nzd++;
379: #else
380: if (a->a[j] != 0.0) nzd++;
381: #endif
382: }
383: }
384: }
385: sptr[m] = nzd+1;
386: PetscViewerASCIIPrintf(viewer," %D %D\n\n",m,nzd);
387: for (i=0; i<m+1; i+=6) {
388: if (i+4<m) {PetscViewerASCIIPrintf(viewer," %D %D %D %D %D %D\n",sptr[i],sptr[i+1],sptr[i+2],sptr[i+3],sptr[i+4],sptr[i+5]);}
389: else if (i+3<m) {PetscViewerASCIIPrintf(viewer," %D %D %D %D %D\n",sptr[i],sptr[i+1],sptr[i+2],sptr[i+3],sptr[i+4]);}
390: else if (i+2<m) {PetscViewerASCIIPrintf(viewer," %D %D %D %D\n",sptr[i],sptr[i+1],sptr[i+2],sptr[i+3]);}
391: else if (i+1<m) {PetscViewerASCIIPrintf(viewer," %D %D %D\n",sptr[i],sptr[i+1],sptr[i+2]);}
392: else if (i<m) {PetscViewerASCIIPrintf(viewer," %D %D\n",sptr[i],sptr[i+1]);}
393: else {PetscViewerASCIIPrintf(viewer," %D\n",sptr[i]);}
394: }
395: PetscViewerASCIIPrintf(viewer,"\n");
396: PetscFree(sptr);
397: for (i=0; i<m; i++) {
398: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
399: if (a->j[j] >= i) {PetscViewerASCIIPrintf(viewer," %D ",a->j[j]+fshift);}
400: }
401: PetscViewerASCIIPrintf(viewer,"\n");
402: }
403: PetscViewerASCIIPrintf(viewer,"\n");
404: for (i=0; i<m; i++) {
405: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
406: if (a->j[j] >= i) {
407: #if defined(PETSC_USE_COMPLEX)
408: if (PetscImaginaryPart(a->a[j]) != 0.0 || PetscRealPart(a->a[j]) != 0.0) {
409: PetscViewerASCIIPrintf(viewer," %18.16e %18.16e ",PetscRealPart(a->a[j]),PetscImaginaryPart(a->a[j]));
410: }
411: #else
412: if (a->a[j] != 0.0) {PetscViewerASCIIPrintf(viewer," %18.16e ",a->a[j]);}
413: #endif
414: }
415: }
416: PetscViewerASCIIPrintf(viewer,"\n");
417: }
418: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
419: } else if (format == PETSC_VIEWER_ASCII_DENSE) {
420: PetscInt cnt = 0,jcnt;
421: PetscScalar value;
423: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
424: for (i=0; i<m; i++) {
425: jcnt = 0;
426: for (j=0; j<A->cmap.n; j++) {
427: if (jcnt < a->i[i+1]-a->i[i] && j == a->j[cnt]) {
428: value = a->a[cnt++];
429: jcnt++;
430: } else {
431: value = 0.0;
432: }
433: #if defined(PETSC_USE_COMPLEX)
434: PetscViewerASCIIPrintf(viewer," %7.5e+%7.5e i ",PetscRealPart(value),PetscImaginaryPart(value));
435: #else
436: PetscViewerASCIIPrintf(viewer," %7.5e ",value);
437: #endif
438: }
439: PetscViewerASCIIPrintf(viewer,"\n");
440: }
441: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
442: } else {
443: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
444: for (i=0; i<m; i++) {
445: PetscViewerASCIIPrintf(viewer,"row %D:",i);
446: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
447: #if defined(PETSC_USE_COMPLEX)
448: if (PetscImaginaryPart(a->a[j]) > 0.0) {
449: PetscViewerASCIIPrintf(viewer," (%D, %G + %G i)",a->j[j]+shift,PetscRealPart(a->a[j]),PetscImaginaryPart(a->a[j]));
450: } else if (PetscImaginaryPart(a->a[j]) < 0.0) {
451: PetscViewerASCIIPrintf(viewer," (%D, %G - %G i)",a->j[j]+shift,PetscRealPart(a->a[j]),-PetscImaginaryPart(a->a[j]));
452: } else {
453: PetscViewerASCIIPrintf(viewer," (%D, %G) ",a->j[j]+shift,PetscRealPart(a->a[j]));
454: }
455: #else
456: PetscViewerASCIIPrintf(viewer," (%D, %G) ",a->j[j]+shift,a->a[j]);
457: #endif
458: }
459: PetscViewerASCIIPrintf(viewer,"\n");
460: }
461: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
462: }
463: PetscViewerFlush(viewer);
464: return(0);
465: }
469: PetscErrorCode MatView_SeqAIJ_Draw_Zoom(PetscDraw draw,void *Aa)
470: {
471: Mat A = (Mat) Aa;
472: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
473: PetscErrorCode ierr;
474: PetscInt i,j,m = A->rmap.n,color;
475: PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r,maxv = 0.0;
476: PetscViewer viewer;
477: PetscViewerFormat format;
480: PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);
481: PetscViewerGetFormat(viewer,&format);
483: PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);
484: /* loop over matrix elements drawing boxes */
486: if (format != PETSC_VIEWER_DRAW_CONTOUR) {
487: /* Blue for negative, Cyan for zero and Red for positive */
488: color = PETSC_DRAW_BLUE;
489: for (i=0; i<m; i++) {
490: y_l = m - i - 1.0; y_r = y_l + 1.0;
491: for (j=a->i[i]; j<a->i[i+1]; j++) {
492: x_l = a->j[j] ; x_r = x_l + 1.0;
493: #if defined(PETSC_USE_COMPLEX)
494: if (PetscRealPart(a->a[j]) >= 0.) continue;
495: #else
496: if (a->a[j] >= 0.) continue;
497: #endif
498: PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
499: }
500: }
501: color = PETSC_DRAW_CYAN;
502: for (i=0; i<m; i++) {
503: y_l = m - i - 1.0; y_r = y_l + 1.0;
504: for (j=a->i[i]; j<a->i[i+1]; j++) {
505: x_l = a->j[j]; x_r = x_l + 1.0;
506: if (a->a[j] != 0.) continue;
507: PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
508: }
509: }
510: color = PETSC_DRAW_RED;
511: for (i=0; i<m; i++) {
512: y_l = m - i - 1.0; y_r = y_l + 1.0;
513: for (j=a->i[i]; j<a->i[i+1]; j++) {
514: x_l = a->j[j]; x_r = x_l + 1.0;
515: #if defined(PETSC_USE_COMPLEX)
516: if (PetscRealPart(a->a[j]) <= 0.) continue;
517: #else
518: if (a->a[j] <= 0.) continue;
519: #endif
520: PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
521: }
522: }
523: } else {
524: /* use contour shading to indicate magnitude of values */
525: /* first determine max of all nonzero values */
526: PetscInt nz = a->nz,count;
527: PetscDraw popup;
528: PetscReal scale;
530: for (i=0; i<nz; i++) {
531: if (PetscAbsScalar(a->a[i]) > maxv) maxv = PetscAbsScalar(a->a[i]);
532: }
533: scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv;
534: PetscDrawGetPopup(draw,&popup);
535: if (popup) {PetscDrawScalePopup(popup,0.0,maxv);}
536: count = 0;
537: for (i=0; i<m; i++) {
538: y_l = m - i - 1.0; y_r = y_l + 1.0;
539: for (j=a->i[i]; j<a->i[i+1]; j++) {
540: x_l = a->j[j]; x_r = x_l + 1.0;
541: color = PETSC_DRAW_BASIC_COLORS + (PetscInt)(scale*PetscAbsScalar(a->a[count]));
542: PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
543: count++;
544: }
545: }
546: }
547: return(0);
548: }
552: PetscErrorCode MatView_SeqAIJ_Draw(Mat A,PetscViewer viewer)
553: {
555: PetscDraw draw;
556: PetscReal xr,yr,xl,yl,h,w;
557: PetscTruth isnull;
560: PetscViewerDrawGetDraw(viewer,0,&draw);
561: PetscDrawIsNull(draw,&isnull);
562: if (isnull) return(0);
564: PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);
565: xr = A->cmap.n; yr = A->rmap.n; h = yr/10.0; w = xr/10.0;
566: xr += w; yr += h; xl = -w; yl = -h;
567: PetscDrawSetCoordinates(draw,xl,yl,xr,yr);
568: PetscDrawZoom(draw,MatView_SeqAIJ_Draw_Zoom,A);
569: PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);
570: return(0);
571: }
575: PetscErrorCode MatView_SeqAIJ(Mat A,PetscViewer viewer)
576: {
578: PetscTruth iascii,isbinary,isdraw;
581: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
582: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
583: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
584: if (iascii) {
585: MatView_SeqAIJ_ASCII(A,viewer);
586: } else if (isbinary) {
587: MatView_SeqAIJ_Binary(A,viewer);
588: } else if (isdraw) {
589: MatView_SeqAIJ_Draw(A,viewer);
590: } else {
591: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by SeqAIJ matrices",((PetscObject)viewer)->type_name);
592: }
593: MatView_Inode(A,viewer);
594: return(0);
595: }
599: PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode)
600: {
601: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
603: PetscInt fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
604: PetscInt m = A->rmap.n,*ip,N,*ailen = a->ilen,rmax = 0;
605: PetscScalar *aa = a->a,*ap;
606: PetscReal ratio=0.6;
609: if (mode == MAT_FLUSH_ASSEMBLY) return(0);
611: if (m) rmax = ailen[0]; /* determine row with most nonzeros */
612: for (i=1; i<m; i++) {
613: /* move each row back by the amount of empty slots (fshift) before it*/
614: fshift += imax[i-1] - ailen[i-1];
615: rmax = PetscMax(rmax,ailen[i]);
616: if (fshift) {
617: ip = aj + ai[i] ;
618: ap = aa + ai[i] ;
619: N = ailen[i];
620: for (j=0; j<N; j++) {
621: ip[j-fshift] = ip[j];
622: ap[j-fshift] = ap[j];
623: }
624: }
625: ai[i] = ai[i-1] + ailen[i-1];
626: }
627: if (m) {
628: fshift += imax[m-1] - ailen[m-1];
629: ai[m] = ai[m-1] + ailen[m-1];
630: }
631: /* reset ilen and imax for each row */
632: for (i=0; i<m; i++) {
633: ailen[i] = imax[i] = ai[i+1] - ai[i];
634: }
635: a->nz = ai[m];
637: MatMarkDiagonal_SeqAIJ(A);
638: PetscInfo4(A,"Matrix size: %D X %D; storage space: %D unneeded,%D used\n",m,A->cmap.n,fshift,a->nz);
639: PetscInfo1(A,"Number of mallocs during MatSetValues() is %D\n",a->reallocs);
640: PetscInfo1(A,"Maximum nonzeros in any row is %D\n",rmax);
642: a->reallocs = 0;
643: A->info.nz_unneeded = (double)fshift;
644: a->rmax = rmax;
646: /* check for zero rows. If found a large number of zero rows, use CompressedRow functions */
647: Mat_CheckCompressedRow(A,&a->compressedrow,a->i,m,ratio);
648: A->same_nonzero = PETSC_TRUE;
650: MatAssemblyEnd_Inode(A,mode);
651: return(0);
652: }
656: PetscErrorCode MatRealPart_SeqAIJ(Mat A)
657: {
658: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
659: PetscInt i,nz = a->nz;
660: PetscScalar *aa = a->a;
663: for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]);
664: return(0);
665: }
669: PetscErrorCode MatImaginaryPart_SeqAIJ(Mat A)
670: {
671: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
672: PetscInt i,nz = a->nz;
673: PetscScalar *aa = a->a;
676: for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]);
677: return(0);
678: }
682: PetscErrorCode MatZeroEntries_SeqAIJ(Mat A)
683: {
684: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
688: PetscMemzero(a->a,(a->i[A->rmap.n])*sizeof(PetscScalar));
689: return(0);
690: }
694: PetscErrorCode MatDestroy_SeqAIJ(Mat A)
695: {
696: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
700: #if defined(PETSC_USE_LOG)
701: PetscLogObjectState((PetscObject)A,"Rows=%D, Cols=%D, NZ=%D",A->rmap.n,A->cmap.n,a->nz);
702: #endif
703: MatSeqXAIJFreeAIJ(A,&a->a,&a->j,&a->i);
704: if (a->row) {
705: ISDestroy(a->row);
706: }
707: if (a->col) {
708: ISDestroy(a->col);
709: }
710: PetscFree(a->diag);
711: PetscFree2(a->imax,a->ilen);
712: PetscFree(a->idiag);
713: PetscFree(a->solve_work);
714: if (a->icol) {ISDestroy(a->icol);}
715: PetscFree(a->saved_values);
716: if (a->coloring) {ISColoringDestroy(a->coloring);}
717: PetscFree(a->xtoy);
718: if (a->XtoY) {MatDestroy(a->XtoY);}
719: if (a->compressedrow.use){PetscFree(a->compressedrow.i);}
721: MatDestroy_Inode(A);
723: PetscFree(a);
725: PetscObjectChangeTypeName((PetscObject)A,0);
726: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSeqAIJSetColumnIndices_C","",PETSC_NULL);
727: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatStoreValues_C","",PETSC_NULL);
728: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatRetrieveValues_C","",PETSC_NULL);
729: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatConvert_seqaij_seqsbaij_C","",PETSC_NULL);
730: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatConvert_seqaij_seqbaij_C","",PETSC_NULL);
731: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatConvert_seqaij_seqcsrperm_C","",PETSC_NULL);
732: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatIsTranspose_C","",PETSC_NULL);
733: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSeqAIJSetPreallocation_C","",PETSC_NULL);
734: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSeqAIJSetPreallocationCSR_C","",PETSC_NULL);
735: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatReorderForNonzeroDiagonal_C","",PETSC_NULL);
736: return(0);
737: }
741: PetscErrorCode MatCompress_SeqAIJ(Mat A)
742: {
744: return(0);
745: }
749: PetscErrorCode MatSetOption_SeqAIJ(Mat A,MatOption op)
750: {
751: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
755: switch (op) {
756: case MAT_ROW_ORIENTED:
757: a->roworiented = PETSC_TRUE;
758: break;
759: case MAT_KEEP_ZEROED_ROWS:
760: a->keepzeroedrows = PETSC_TRUE;
761: break;
762: case MAT_COLUMN_ORIENTED:
763: a->roworiented = PETSC_FALSE;
764: break;
765: case MAT_COLUMNS_SORTED:
766: a->sorted = PETSC_TRUE;
767: break;
768: case MAT_COLUMNS_UNSORTED:
769: a->sorted = PETSC_FALSE;
770: break;
771: case MAT_NO_NEW_NONZERO_LOCATIONS:
772: a->nonew = 1;
773: break;
774: case MAT_NEW_NONZERO_LOCATION_ERR:
775: a->nonew = -1;
776: break;
777: case MAT_NEW_NONZERO_ALLOCATION_ERR:
778: a->nonew = -2;
779: break;
780: case MAT_YES_NEW_NONZERO_LOCATIONS:
781: a->nonew = 0;
782: break;
783: case MAT_IGNORE_ZERO_ENTRIES:
784: a->ignorezeroentries = PETSC_TRUE;
785: break;
786: case MAT_USE_COMPRESSEDROW:
787: a->compressedrow.use = PETSC_TRUE;
788: break;
789: case MAT_DO_NOT_USE_COMPRESSEDROW:
790: a->compressedrow.use = PETSC_FALSE;
791: break;
792: case MAT_ROWS_SORTED:
793: case MAT_ROWS_UNSORTED:
794: case MAT_YES_NEW_DIAGONALS:
795: case MAT_IGNORE_OFF_PROC_ENTRIES:
796: case MAT_USE_HASH_TABLE:
797: PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
798: break;
799: case MAT_NO_NEW_DIAGONALS:
800: SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
801: default:
802: break;
803: }
804: MatSetOption_Inode(A,op);
805: return(0);
806: }
810: PetscErrorCode MatGetDiagonal_SeqAIJ(Mat A,Vec v)
811: {
812: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
814: PetscInt i,j,n;
815: PetscScalar *x,zero = 0.0;
818: VecSet(v,zero);
819: VecGetArray(v,&x);
820: VecGetLocalSize(v,&n);
821: if (n != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
822: for (i=0; i<A->rmap.n; i++) {
823: for (j=a->i[i]; j<a->i[i+1]; j++) {
824: if (a->j[j] == i) {
825: x[i] = a->a[j];
826: break;
827: }
828: }
829: }
830: VecRestoreArray(v,&x);
831: return(0);
832: }
836: PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
837: {
838: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
839: PetscScalar *x,*y;
840: PetscErrorCode ierr;
841: PetscInt m = A->rmap.n;
842: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
843: PetscScalar *v,alpha;
844: PetscInt n,i,*idx,*ii,*ridx=PETSC_NULL;
845: Mat_CompressedRow cprow = a->compressedrow;
846: PetscTruth usecprow = cprow.use;
847: #endif
850: if (zz != yy) {VecCopy(zz,yy);}
851: VecGetArray(xx,&x);
852: VecGetArray(yy,&y);
854: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
855: fortranmulttransposeaddaij_(&m,x,a->i,a->j,a->a,y);
856: #else
857: if (usecprow){
858: m = cprow.nrows;
859: ii = cprow.i;
860: ridx = cprow.rindex;
861: } else {
862: ii = a->i;
863: }
864: for (i=0; i<m; i++) {
865: idx = a->j + ii[i] ;
866: v = a->a + ii[i] ;
867: n = ii[i+1] - ii[i];
868: if (usecprow){
869: alpha = x[ridx[i]];
870: } else {
871: alpha = x[i];
872: }
873: while (n-->0) {y[*idx++] += alpha * *v++;}
874: }
875: #endif
876: PetscLogFlops(2*a->nz);
877: VecRestoreArray(xx,&x);
878: VecRestoreArray(yy,&y);
879: return(0);
880: }
884: PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec xx,Vec yy)
885: {
886: PetscScalar zero = 0.0;
890: VecSet(yy,zero);
891: MatMultTransposeAdd_SeqAIJ(A,xx,yy,yy);
892: return(0);
893: }
898: PetscErrorCode MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
899: {
900: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
901: PetscScalar *x,*y,*aa;
903: PetscInt m=A->rmap.n,*aj,*ii;
904: PetscInt n,i,j,nonzerorow=0,*ridx=PETSC_NULL;
905: PetscScalar sum;
906: PetscTruth usecprow=a->compressedrow.use;
907: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
908: PetscInt jrow;
909: #endif
911: #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
912: #pragma disjoint(*x,*y,*aa)
913: #endif
916: VecGetArray(xx,&x);
917: VecGetArray(yy,&y);
918: aj = a->j;
919: aa = a->a;
920: ii = a->i;
921: if (usecprow){ /* use compressed row format */
922: m = a->compressedrow.nrows;
923: ii = a->compressedrow.i;
924: ridx = a->compressedrow.rindex;
925: for (i=0; i<m; i++){
926: n = ii[i+1] - ii[i];
927: aj = a->j + ii[i];
928: aa = a->a + ii[i];
929: sum = 0.0;
930: nonzerorow += (n>0);
931: for (j=0; j<n; j++) sum += (*aa++)*x[*aj++];
932: y[*ridx++] = sum;
933: }
934: } else { /* do not use compressed row format */
935: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
936: fortranmultaij_(&m,x,ii,aj,aa,y);
937: #else
938: for (i=0; i<m; i++) {
939: jrow = ii[i];
940: n = ii[i+1] - jrow;
941: sum = 0.0;
942: nonzerorow += (n>0);
943: for (j=0; j<n; j++) {
944: sum += aa[jrow]*x[aj[jrow]]; jrow++;
945: }
946: y[i] = sum;
947: }
948: #endif
949: }
950: PetscLogFlops(2*a->nz - nonzerorow);
951: VecRestoreArray(xx,&x);
952: VecRestoreArray(yy,&y);
953: return(0);
954: }
958: PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
959: {
960: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
961: PetscScalar *x,*y,*z,*aa;
963: PetscInt m = A->rmap.n,*aj,*ii;
964: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
965: PetscInt n,i,jrow,j,*ridx=PETSC_NULL;
966: PetscScalar sum;
967: PetscTruth usecprow=a->compressedrow.use;
968: #endif
971: VecGetArray(xx,&x);
972: VecGetArray(yy,&y);
973: if (zz != yy) {
974: VecGetArray(zz,&z);
975: } else {
976: z = y;
977: }
979: aj = a->j;
980: aa = a->a;
981: ii = a->i;
982: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
983: fortranmultaddaij_(&m,x,ii,aj,aa,y,z);
984: #else
985: if (usecprow){ /* use compressed row format */
986: if (zz != yy){
987: PetscMemcpy(z,y,m*sizeof(PetscScalar));
988: }
989: m = a->compressedrow.nrows;
990: ii = a->compressedrow.i;
991: ridx = a->compressedrow.rindex;
992: for (i=0; i<m; i++){
993: n = ii[i+1] - ii[i];
994: aj = a->j + ii[i];
995: aa = a->a + ii[i];
996: sum = y[*ridx];
997: for (j=0; j<n; j++) sum += (*aa++)*x[*aj++];
998: z[*ridx++] = sum;
999: }
1000: } else { /* do not use compressed row format */
1001: for (i=0; i<m; i++) {
1002: jrow = ii[i];
1003: n = ii[i+1] - jrow;
1004: sum = y[i];
1005: for (j=0; j<n; j++) {
1006: sum += aa[jrow]*x[aj[jrow]]; jrow++;
1007: }
1008: z[i] = sum;
1009: }
1010: }
1011: #endif
1012: PetscLogFlops(2*a->nz);
1013: VecRestoreArray(xx,&x);
1014: VecRestoreArray(yy,&y);
1015: if (zz != yy) {
1016: VecRestoreArray(zz,&z);
1017: }
1018: return(0);
1019: }
1021: /*
1022: Adds diagonal pointers to sparse matrix structure.
1023: */
1026: PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat A)
1027: {
1028: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1030: PetscInt i,j,m = A->rmap.n;
1033: if (!a->diag) {
1034: PetscMalloc(m*sizeof(PetscInt),&a->diag);
1035: PetscLogObjectMemory(A, m*sizeof(PetscInt));
1036: }
1037: for (i=0; i<A->rmap.n; i++) {
1038: a->diag[i] = a->i[i+1];
1039: for (j=a->i[i]; j<a->i[i+1]; j++) {
1040: if (a->j[j] == i) {
1041: a->diag[i] = j;
1042: break;
1043: }
1044: }
1045: }
1046: return(0);
1047: }
1049: /*
1050: Checks for missing diagonals
1051: */
1054: PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat A,PetscTruth *missing,PetscInt *d)
1055: {
1056: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1057: PetscInt *diag,*jj = a->j,i;
1060: *missing = PETSC_FALSE;
1061: if (A->rmap.n > 0 && !jj) {
1062: *missing = PETSC_TRUE;
1063: if (d) *d = 0;
1064: PetscInfo(A,"Matrix has no entries therefor is missing diagonal");
1065: } else {
1066: diag = a->diag;
1067: for (i=0; i<A->rmap.n; i++) {
1068: if (jj[diag[i]] != i) {
1069: *missing = PETSC_TRUE;
1070: if (d) *d = i;
1071: PetscInfo1(A,"Matrix is missing diagonal number %D",i);
1072: }
1073: }
1074: }
1075: return(0);
1076: }
1080: PetscErrorCode MatRelax_SeqAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
1081: {
1082: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1083: PetscScalar *x,d,*xs,sum,*t,scale,*idiag=0,*mdiag;
1084: const PetscScalar *v = a->a, *b, *bs,*xb, *ts;
1085: PetscErrorCode ierr;
1086: PetscInt n = A->cmap.n,m = A->rmap.n,i;
1087: const PetscInt *idx,*diag;
1090: its = its*lits;
1091: if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
1093: diag = a->diag;
1094: if (!a->idiag) {
1095: PetscMalloc(3*m*sizeof(PetscScalar),&a->idiag);
1096: PetscLogObjectMemory(A, 3*m*sizeof(PetscScalar));
1097: a->ssor = a->idiag + m;
1098: mdiag = a->ssor + m;
1099: v = a->a;
1101: /* this is wrong when fshift omega changes each iteration */
1102: if (omega == 1.0 && !fshift) {
1103: for (i=0; i<m; i++) {
1104: mdiag[i] = v[diag[i]];
1105: if (!PetscAbsScalar(mdiag[i])) SETERRQ1(PETSC_ERR_ARG_INCOMP,"Zero diagonal on row %D",i);
1106: a->idiag[i] = 1.0/v[diag[i]];
1107: }
1108: PetscLogFlops(m);
1109: } else {
1110: for (i=0; i<m; i++) {
1111: mdiag[i] = v[diag[i]];
1112: a->idiag[i] = omega/(fshift + v[diag[i]]);
1113: }
1114: PetscLogFlops(2*m);
1115: }
1116: }
1117: t = a->ssor;
1118: idiag = a->idiag;
1119: mdiag = a->idiag + 2*m;
1121: VecGetArray(xx,&x);
1122: if (xx != bb) {
1123: VecGetArray(bb,(PetscScalar**)&b);
1124: } else {
1125: b = x;
1126: }
1128: /* We count flops by assuming the upper triangular and lower triangular parts have the same number of nonzeros */
1129: xs = x;
1130: if (flag == SOR_APPLY_UPPER) {
1131: /* apply (U + D/omega) to the vector */
1132: bs = b;
1133: for (i=0; i<m; i++) {
1134: d = fshift + a->a[diag[i]];
1135: n = a->i[i+1] - diag[i] - 1;
1136: idx = a->j + diag[i] + 1;
1137: v = a->a + diag[i] + 1;
1138: sum = b[i]*d/omega;
1139: SPARSEDENSEDOT(sum,bs,v,idx,n);
1140: x[i] = sum;
1141: }
1142: VecRestoreArray(xx,&x);
1143: if (bb != xx) {VecRestoreArray(bb,(PetscScalar**)&b);}
1144: PetscLogFlops(a->nz);
1145: return(0);
1146: }
1149: /* Let A = L + U + D; where L is lower trianglar,
1150: U is upper triangular, E is diagonal; This routine applies
1152: (L + E)^{-1} A (U + E)^{-1}
1154: to a vector efficiently using Eisenstat's trick. This is for
1155: the case of SSOR preconditioner, so E is D/omega where omega
1156: is the relaxation factor.
1157: */
1159: if (flag == SOR_APPLY_LOWER) {
1160: SETERRQ(PETSC_ERR_SUP,"SOR_APPLY_LOWER is not implemented");
1161: } else if (flag & SOR_EISENSTAT) {
1162: /* Let A = L + U + D; where L is lower trianglar,
1163: U is upper triangular, E is diagonal; This routine applies
1165: (L + E)^{-1} A (U + E)^{-1}
1167: to a vector efficiently using Eisenstat's trick. This is for
1168: the case of SSOR preconditioner, so E is D/omega where omega
1169: is the relaxation factor.
1170: */
1171: scale = (2.0/omega) - 1.0;
1173: /* x = (E + U)^{-1} b */
1174: for (i=m-1; i>=0; i--) {
1175: n = a->i[i+1] - diag[i] - 1;
1176: idx = a->j + diag[i] + 1;
1177: v = a->a + diag[i] + 1;
1178: sum = b[i];
1179: SPARSEDENSEMDOT(sum,xs,v,idx,n);
1180: x[i] = sum*idiag[i];
1181: }
1183: /* t = b - (2*E - D)x */
1184: v = a->a;
1185: for (i=0; i<m; i++) { t[i] = b[i] - scale*(v[*diag++])*x[i]; }
1187: /* t = (E + L)^{-1}t */
1188: ts = t;
1189: diag = a->diag;
1190: for (i=0; i<m; i++) {
1191: n = diag[i] - a->i[i];
1192: idx = a->j + a->i[i];
1193: v = a->a + a->i[i];
1194: sum = t[i];
1195: SPARSEDENSEMDOT(sum,ts,v,idx,n);
1196: t[i] = sum*idiag[i];
1197: /* x = x + t */
1198: x[i] += t[i];
1199: }
1201: PetscLogFlops(6*m-1 + 2*a->nz);
1202: VecRestoreArray(xx,&x);
1203: if (bb != xx) {VecRestoreArray(bb,(PetscScalar**)&b);}
1204: return(0);
1205: }
1206: if (flag & SOR_ZERO_INITIAL_GUESS) {
1207: if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
1208: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1209: fortranrelaxaijforwardzero_(&m,&omega,x,a->i,a->j,(PetscInt*)diag,idiag,a->a,(void*)b);
1210: #else
1211: for (i=0; i<m; i++) {
1212: n = diag[i] - a->i[i];
1213: idx = a->j + a->i[i];
1214: v = a->a + a->i[i];
1215: sum = b[i];
1216: SPARSEDENSEMDOT(sum,xs,v,idx,n);
1217: x[i] = sum*idiag[i];
1218: }
1219: #endif
1220: xb = x;
1221: PetscLogFlops(a->nz);
1222: } else xb = b;
1223: if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
1224: (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
1225: for (i=0; i<m; i++) {
1226: x[i] *= mdiag[i];
1227: }
1228: PetscLogFlops(m);
1229: }
1230: if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
1231: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1232: fortranrelaxaijbackwardzero_(&m,&omega,x,a->i,a->j,(PetscInt*)diag,idiag,a->a,(void*)xb);
1233: #else
1234: for (i=m-1; i>=0; i--) {
1235: n = a->i[i+1] - diag[i] - 1;
1236: idx = a->j + diag[i] + 1;
1237: v = a->a + diag[i] + 1;
1238: sum = xb[i];
1239: SPARSEDENSEMDOT(sum,xs,v,idx,n);
1240: x[i] = sum*idiag[i];
1241: }
1242: #endif
1243: PetscLogFlops(a->nz);
1244: }
1245: its--;
1246: }
1247: while (its--) {
1248: if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
1249: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1250: fortranrelaxaijforward_(&m,&omega,x,a->i,a->j,(PetscInt*)diag,a->a,(void*)b);
1251: #else
1252: for (i=0; i<m; i++) {
1253: n = a->i[i+1] - a->i[i];
1254: idx = a->j + a->i[i];
1255: v = a->a + a->i[i];
1256: sum = b[i];
1257: SPARSEDENSEMDOT(sum,xs,v,idx,n);
1258: x[i] = (1. - omega)*x[i] + (sum + mdiag[i]*x[i])*idiag[i];
1259: }
1260: #endif
1261: PetscLogFlops(a->nz);
1262: }
1263: if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
1264: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1265: fortranrelaxaijbackward_(&m,&omega,x,a->i,a->j,(PetscInt*)diag,a->a,(void*)b);
1266: #else
1267: for (i=m-1; i>=0; i--) {
1268: n = a->i[i+1] - a->i[i];
1269: idx = a->j + a->i[i];
1270: v = a->a + a->i[i];
1271: sum = b[i];
1272: SPARSEDENSEMDOT(sum,xs,v,idx,n);
1273: x[i] = (1. - omega)*x[i] + (sum + mdiag[i]*x[i])*idiag[i];
1274: }
1275: #endif
1276: PetscLogFlops(a->nz);
1277: }
1278: }
1279: VecRestoreArray(xx,&x);
1280: if (bb != xx) {VecRestoreArray(bb,(PetscScalar**)&b);}
1281: return(0);
1282: }
1286: PetscErrorCode MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,MatInfo *info)
1287: {
1288: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1291: info->rows_global = (double)A->rmap.n;
1292: info->columns_global = (double)A->cmap.n;
1293: info->rows_local = (double)A->rmap.n;
1294: info->columns_local = (double)A->cmap.n;
1295: info->block_size = 1.0;
1296: info->nz_allocated = (double)a->maxnz;
1297: info->nz_used = (double)a->nz;
1298: info->nz_unneeded = (double)(a->maxnz - a->nz);
1299: info->assemblies = (double)A->num_ass;
1300: info->mallocs = (double)a->reallocs;
1301: info->memory = A->mem;
1302: if (A->factor) {
1303: info->fill_ratio_given = A->info.fill_ratio_given;
1304: info->fill_ratio_needed = A->info.fill_ratio_needed;
1305: info->factor_mallocs = A->info.factor_mallocs;
1306: } else {
1307: info->fill_ratio_given = 0;
1308: info->fill_ratio_needed = 0;
1309: info->factor_mallocs = 0;
1310: }
1311: return(0);
1312: }
1316: PetscErrorCode MatZeroRows_SeqAIJ(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag)
1317: {
1318: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1319: PetscInt i,m = A->rmap.n - 1,d;
1321: PetscTruth missing;
1324: if (a->keepzeroedrows) {
1325: for (i=0; i<N; i++) {
1326: if (rows[i] < 0 || rows[i] > m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range", rows[i]);
1327: PetscMemzero(&a->a[a->i[rows[i]]],a->ilen[rows[i]]*sizeof(PetscScalar));
1328: }
1329: if (diag != 0.0) {
1330: MatMissingDiagonal_SeqAIJ(A,&missing,&d);
1331: if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry in row %D",d);
1332: for (i=0; i<N; i++) {
1333: a->a[a->diag[rows[i]]] = diag;
1334: }
1335: }
1336: A->same_nonzero = PETSC_TRUE;
1337: } else {
1338: if (diag != 0.0) {
1339: for (i=0; i<N; i++) {
1340: if (rows[i] < 0 || rows[i] > m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range", rows[i]);
1341: if (a->ilen[rows[i]] > 0) {
1342: a->ilen[rows[i]] = 1;
1343: a->a[a->i[rows[i]]] = diag;
1344: a->j[a->i[rows[i]]] = rows[i];
1345: } else { /* in case row was completely empty */
1346: MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],&diag,INSERT_VALUES);
1347: }
1348: }
1349: } else {
1350: for (i=0; i<N; i++) {
1351: if (rows[i] < 0 || rows[i] > m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range", rows[i]);
1352: a->ilen[rows[i]] = 0;
1353: }
1354: }
1355: A->same_nonzero = PETSC_FALSE;
1356: }
1357: MatAssemblyEnd_SeqAIJ(A,MAT_FINAL_ASSEMBLY);
1358: return(0);
1359: }
1363: PetscErrorCode MatGetRow_SeqAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
1364: {
1365: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1366: PetscInt *itmp;
1369: if (row < 0 || row >= A->rmap.n) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Row %D out of range",row);
1371: *nz = a->i[row+1] - a->i[row];
1372: if (v) *v = a->a + a->i[row];
1373: if (idx) {
1374: itmp = a->j + a->i[row];
1375: if (*nz) {
1376: *idx = itmp;
1377: }
1378: else *idx = 0;
1379: }
1380: return(0);
1381: }
1383: /* remove this function? */
1386: PetscErrorCode MatRestoreRow_SeqAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
1387: {
1389: return(0);
1390: }
1394: PetscErrorCode MatNorm_SeqAIJ(Mat A,NormType type,PetscReal *nrm)
1395: {
1396: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1397: PetscScalar *v = a->a;
1398: PetscReal sum = 0.0;
1400: PetscInt i,j;
1403: if (type == NORM_FROBENIUS) {
1404: for (i=0; i<a->nz; i++) {
1405: #if defined(PETSC_USE_COMPLEX)
1406: sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1407: #else
1408: sum += (*v)*(*v); v++;
1409: #endif
1410: }
1411: *nrm = sqrt(sum);
1412: } else if (type == NORM_1) {
1413: PetscReal *tmp;
1414: PetscInt *jj = a->j;
1415: PetscMalloc((A->cmap.n+1)*sizeof(PetscReal),&tmp);
1416: PetscMemzero(tmp,A->cmap.n*sizeof(PetscReal));
1417: *nrm = 0.0;
1418: for (j=0; j<a->nz; j++) {
1419: tmp[*jj++] += PetscAbsScalar(*v); v++;
1420: }
1421: for (j=0; j<A->cmap.n; j++) {
1422: if (tmp[j] > *nrm) *nrm = tmp[j];
1423: }
1424: PetscFree(tmp);
1425: } else if (type == NORM_INFINITY) {
1426: *nrm = 0.0;
1427: for (j=0; j<A->rmap.n; j++) {
1428: v = a->a + a->i[j];
1429: sum = 0.0;
1430: for (i=0; i<a->i[j+1]-a->i[j]; i++) {
1431: sum += PetscAbsScalar(*v); v++;
1432: }
1433: if (sum > *nrm) *nrm = sum;
1434: }
1435: } else {
1436: SETERRQ(PETSC_ERR_SUP,"No support for two norm");
1437: }
1438: return(0);
1439: }
1443: PetscErrorCode MatTranspose_SeqAIJ(Mat A,Mat *B)
1444: {
1445: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1446: Mat C;
1448: PetscInt i,*aj = a->j,*ai = a->i,m = A->rmap.n,len,*col;
1449: PetscScalar *array = a->a;
1452: if (!B && m != A->cmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Square matrix only for in-place");
1453: PetscMalloc((1+A->cmap.n)*sizeof(PetscInt),&col);
1454: PetscMemzero(col,(1+A->cmap.n)*sizeof(PetscInt));
1455:
1456: for (i=0; i<ai[m]; i++) col[aj[i]] += 1;
1457: MatCreate(A->comm,&C);
1458: MatSetSizes(C,A->cmap.n,m,A->cmap.n,m);
1459: MatSetType(C,A->type_name);
1460: MatSeqAIJSetPreallocation_SeqAIJ(C,0,col);
1461: PetscFree(col);
1462: for (i=0; i<m; i++) {
1463: len = ai[i+1]-ai[i];
1464: MatSetValues_SeqAIJ(C,len,aj,1,&i,array,INSERT_VALUES);
1465: array += len;
1466: aj += len;
1467: }
1469: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1470: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
1472: if (B) {
1473: *B = C;
1474: } else {
1475: MatHeaderCopy(A,C);
1476: }
1477: return(0);
1478: }
1483: PetscErrorCode MatIsTranspose_SeqAIJ(Mat A,Mat B,PetscReal tol,PetscTruth *f)
1484: {
1485: Mat_SeqAIJ *aij = (Mat_SeqAIJ *) A->data,*bij = (Mat_SeqAIJ*) A->data;
1486: PetscInt *adx,*bdx,*aii,*bii,*aptr,*bptr; PetscScalar *va,*vb;
1488: PetscInt ma,na,mb,nb, i;
1491: bij = (Mat_SeqAIJ *) B->data;
1492:
1493: MatGetSize(A,&ma,&na);
1494: MatGetSize(B,&mb,&nb);
1495: if (ma!=nb || na!=mb){
1496: *f = PETSC_FALSE;
1497: return(0);
1498: }
1499: aii = aij->i; bii = bij->i;
1500: adx = aij->j; bdx = bij->j;
1501: va = aij->a; vb = bij->a;
1502: PetscMalloc(ma*sizeof(PetscInt),&aptr);
1503: PetscMalloc(mb*sizeof(PetscInt),&bptr);
1504: for (i=0; i<ma; i++) aptr[i] = aii[i];
1505: for (i=0; i<mb; i++) bptr[i] = bii[i];
1507: *f = PETSC_TRUE;
1508: for (i=0; i<ma; i++) {
1509: while (aptr[i]<aii[i+1]) {
1510: PetscInt idc,idr;
1511: PetscScalar vc,vr;
1512: /* column/row index/value */
1513: idc = adx[aptr[i]];
1514: idr = bdx[bptr[idc]];
1515: vc = va[aptr[i]];
1516: vr = vb[bptr[idc]];
1517: if (i!=idr || PetscAbsScalar(vc-vr) > tol) {
1518: *f = PETSC_FALSE;
1519: goto done;
1520: } else {
1521: aptr[i]++;
1522: if (B || i!=idc) bptr[idc]++;
1523: }
1524: }
1525: }
1526: done:
1527: PetscFree(aptr);
1528: if (B) {
1529: PetscFree(bptr);
1530: }
1531: return(0);
1532: }
1537: PetscErrorCode MatIsSymmetric_SeqAIJ(Mat A,PetscReal tol,PetscTruth *f)
1538: {
1541: MatIsTranspose_SeqAIJ(A,A,tol,f);
1542: return(0);
1543: }
1547: PetscErrorCode MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr)
1548: {
1549: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1550: PetscScalar *l,*r,x,*v;
1552: PetscInt i,j,m = A->rmap.n,n = A->cmap.n,M,nz = a->nz,*jj;
1555: if (ll) {
1556: /* The local size is used so that VecMPI can be passed to this routine
1557: by MatDiagonalScale_MPIAIJ */
1558: VecGetLocalSize(ll,&m);
1559: if (m != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length");
1560: VecGetArray(ll,&l);
1561: v = a->a;
1562: for (i=0; i<m; i++) {
1563: x = l[i];
1564: M = a->i[i+1] - a->i[i];
1565: for (j=0; j<M; j++) { (*v++) *= x;}
1566: }
1567: VecRestoreArray(ll,&l);
1568: PetscLogFlops(nz);
1569: }
1570: if (rr) {
1571: VecGetLocalSize(rr,&n);
1572: if (n != A->cmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vector wrong length");
1573: VecGetArray(rr,&r);
1574: v = a->a; jj = a->j;
1575: for (i=0; i<nz; i++) {
1576: (*v++) *= r[*jj++];
1577: }
1578: VecRestoreArray(rr,&r);
1579: PetscLogFlops(nz);
1580: }
1581: return(0);
1582: }
1586: PetscErrorCode MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,PetscInt csize,MatReuse scall,Mat *B)
1587: {
1588: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*c;
1590: PetscInt *smap,i,k,kstart,kend,oldcols = A->cmap.n,*lens;
1591: PetscInt row,mat_i,*mat_j,tcol,first,step,*mat_ilen,sum,lensi;
1592: PetscInt *irow,*icol,nrows,ncols;
1593: PetscInt *starts,*j_new,*i_new,*aj = a->j,*ai = a->i,ii,*ailen = a->ilen;
1594: PetscScalar *a_new,*mat_a;
1595: Mat C;
1596: PetscTruth stride;
1599: ISSorted(isrow,(PetscTruth*)&i);
1600: if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"ISrow is not sorted");
1601: ISSorted(iscol,(PetscTruth*)&i);
1602: if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IScol is not sorted");
1604: ISGetIndices(isrow,&irow);
1605: ISGetLocalSize(isrow,&nrows);
1606: ISGetLocalSize(iscol,&ncols);
1608: ISStrideGetInfo(iscol,&first,&step);
1609: ISStride(iscol,&stride);
1610: if (stride && step == 1) {
1611: /* special case of contiguous rows */
1612: PetscMalloc((2*nrows+1)*sizeof(PetscInt),&lens);
1613: starts = lens + nrows;
1614: /* loop over new rows determining lens and starting points */
1615: for (i=0; i<nrows; i++) {
1616: kstart = ai[irow[i]];
1617: kend = kstart + ailen[irow[i]];
1618: for (k=kstart; k<kend; k++) {
1619: if (aj[k] >= first) {
1620: starts[i] = k;
1621: break;
1622: }
1623: }
1624: sum = 0;
1625: while (k < kend) {
1626: if (aj[k++] >= first+ncols) break;
1627: sum++;
1628: }
1629: lens[i] = sum;
1630: }
1631: /* create submatrix */
1632: if (scall == MAT_REUSE_MATRIX) {
1633: PetscInt n_cols,n_rows;
1634: MatGetSize(*B,&n_rows,&n_cols);
1635: if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size");
1636: MatZeroEntries(*B);
1637: C = *B;
1638: } else {
1639: MatCreate(A->comm,&C);
1640: MatSetSizes(C,nrows,ncols,PETSC_DETERMINE,PETSC_DETERMINE);
1641: MatSetType(C,A->type_name);
1642: MatSeqAIJSetPreallocation_SeqAIJ(C,0,lens);
1643: }
1644: c = (Mat_SeqAIJ*)C->data;
1646: /* loop over rows inserting into submatrix */
1647: a_new = c->a;
1648: j_new = c->j;
1649: i_new = c->i;
1651: for (i=0; i<nrows; i++) {
1652: ii = starts[i];
1653: lensi = lens[i];
1654: for (k=0; k<lensi; k++) {
1655: *j_new++ = aj[ii+k] - first;
1656: }
1657: PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(PetscScalar));
1658: a_new += lensi;
1659: i_new[i+1] = i_new[i] + lensi;
1660: c->ilen[i] = lensi;
1661: }
1662: PetscFree(lens);
1663: } else {
1664: ISGetIndices(iscol,&icol);
1665: PetscMalloc((1+oldcols)*sizeof(PetscInt),&smap);
1666:
1667: PetscMalloc((1+nrows)*sizeof(PetscInt),&lens);
1668: PetscMemzero(smap,oldcols*sizeof(PetscInt));
1669: for (i=0; i<ncols; i++) smap[icol[i]] = i+1;
1670: /* determine lens of each row */
1671: for (i=0; i<nrows; i++) {
1672: kstart = ai[irow[i]];
1673: kend = kstart + a->ilen[irow[i]];
1674: lens[i] = 0;
1675: for (k=kstart; k<kend; k++) {
1676: if (smap[aj[k]]) {
1677: lens[i]++;
1678: }
1679: }
1680: }
1681: /* Create and fill new matrix */
1682: if (scall == MAT_REUSE_MATRIX) {
1683: PetscTruth equal;
1685: c = (Mat_SeqAIJ *)((*B)->data);
1686: if ((*B)->rmap.n != nrows || (*B)->cmap.n != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size");
1687: PetscMemcmp(c->ilen,lens,(*B)->rmap.n*sizeof(PetscInt),&equal);
1688: if (!equal) {
1689: SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros");
1690: }
1691: PetscMemzero(c->ilen,(*B)->rmap.n*sizeof(PetscInt));
1692: C = *B;
1693: } else {
1694: MatCreate(A->comm,&C);
1695: MatSetSizes(C,nrows,ncols,PETSC_DETERMINE,PETSC_DETERMINE);
1696: MatSetType(C,A->type_name);
1697: MatSeqAIJSetPreallocation_SeqAIJ(C,0,lens);
1698: }
1699: c = (Mat_SeqAIJ *)(C->data);
1700: for (i=0; i<nrows; i++) {
1701: row = irow[i];
1702: kstart = ai[row];
1703: kend = kstart + a->ilen[row];
1704: mat_i = c->i[i];
1705: mat_j = c->j + mat_i;
1706: mat_a = c->a + mat_i;
1707: mat_ilen = c->ilen + i;
1708: for (k=kstart; k<kend; k++) {
1709: if ((tcol=smap[a->j[k]])) {
1710: *mat_j++ = tcol - 1;
1711: *mat_a++ = a->a[k];
1712: (*mat_ilen)++;
1714: }
1715: }
1716: }
1717: /* Free work space */
1718: ISRestoreIndices(iscol,&icol);
1719: PetscFree(smap);
1720: PetscFree(lens);
1721: }
1722: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1723: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
1725: ISRestoreIndices(isrow,&irow);
1726: *B = C;
1727: return(0);
1728: }
1730: /*
1731: */
1734: PetscErrorCode MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,MatFactorInfo *info)
1735: {
1736: Mat_SeqAIJ *a = (Mat_SeqAIJ*)inA->data;
1738: Mat outA;
1739: PetscTruth row_identity,col_identity;
1742: if (info->levels != 0) SETERRQ(PETSC_ERR_SUP,"Only levels=0 supported for in-place ilu");
1743: ISIdentity(row,&row_identity);
1744: ISIdentity(col,&col_identity);
1746: outA = inA;
1747: inA->factor = FACTOR_LU;
1748: PetscObjectReference((PetscObject)row);
1749: if (a->row) { ISDestroy(a->row); }
1750: a->row = row;
1751: PetscObjectReference((PetscObject)col);
1752: if (a->col) { ISDestroy(a->col); }
1753: a->col = col;
1755: /* Create the inverse permutation so that it can be used in MatLUFactorNumeric() */
1756: if (a->icol) {ISDestroy(a->icol);} /* need to remove old one */
1757: ISInvertPermutation(col,PETSC_DECIDE,&a->icol);
1758: PetscLogObjectParent(inA,a->icol);
1760: if (!a->solve_work) { /* this matrix may have been factored before */
1761: PetscMalloc((inA->rmap.n+1)*sizeof(PetscScalar),&a->solve_work);
1762: PetscLogObjectMemory(inA, (inA->rmap.n+1)*sizeof(PetscScalar));
1763: }
1765: MatMarkDiagonal_SeqAIJ(inA);
1766: if (row_identity && col_identity) {
1767: MatLUFactorNumeric_SeqAIJ(inA,info,&outA);
1768: } else {
1769: MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(inA,info,&outA);
1770: }
1771: return(0);
1772: }
1774: #include petscblaslapack.h
1777: PetscErrorCode MatScale_SeqAIJ(Mat inA,PetscScalar alpha)
1778: {
1779: Mat_SeqAIJ *a = (Mat_SeqAIJ*)inA->data;
1780: PetscBLASInt bnz = (PetscBLASInt)a->nz,one = 1;
1781: PetscScalar oalpha = alpha;
1786: BLASscal_(&bnz,&oalpha,a->a,&one);
1787: PetscLogFlops(a->nz);
1788: return(0);
1789: }
1793: PetscErrorCode MatGetSubMatrices_SeqAIJ(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1794: {
1796: PetscInt i;
1799: if (scall == MAT_INITIAL_MATRIX) {
1800: PetscMalloc((n+1)*sizeof(Mat),B);
1801: }
1803: for (i=0; i<n; i++) {
1804: MatGetSubMatrix_SeqAIJ(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);
1805: }
1806: return(0);
1807: }
1811: PetscErrorCode MatIncreaseOverlap_SeqAIJ(Mat A,PetscInt is_max,IS is[],PetscInt ov)
1812: {
1813: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1815: PetscInt row,i,j,k,l,m,n,*idx,*nidx,isz,val;
1816: PetscInt start,end,*ai,*aj;
1817: PetscBT table;
1820: m = A->rmap.n;
1821: ai = a->i;
1822: aj = a->j;
1824: if (ov < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"illegal negative overlap value used");
1826: PetscMalloc((m+1)*sizeof(PetscInt),&nidx);
1827: PetscBTCreate(m,table);
1829: for (i=0; i<is_max; i++) {
1830: /* Initialize the two local arrays */
1831: isz = 0;
1832: PetscBTMemzero(m,table);
1833:
1834: /* Extract the indices, assume there can be duplicate entries */
1835: ISGetIndices(is[i],&idx);
1836: ISGetLocalSize(is[i],&n);
1837:
1838: /* Enter these into the temp arrays. I.e., mark table[row], enter row into new index */
1839: for (j=0; j<n ; ++j){
1840: if(!PetscBTLookupSet(table,idx[j])) { nidx[isz++] = idx[j];}
1841: }
1842: ISRestoreIndices(is[i],&idx);
1843: ISDestroy(is[i]);
1844:
1845: k = 0;
1846: for (j=0; j<ov; j++){ /* for each overlap */
1847: n = isz;
1848: for (; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */
1849: row = nidx[k];
1850: start = ai[row];
1851: end = ai[row+1];
1852: for (l = start; l<end ; l++){
1853: val = aj[l] ;
1854: if (!PetscBTLookupSet(table,val)) {nidx[isz++] = val;}
1855: }
1856: }
1857: }
1858: ISCreateGeneral(PETSC_COMM_SELF,isz,nidx,(is+i));
1859: }
1860: PetscBTDestroy(table);
1861: PetscFree(nidx);
1862: return(0);
1863: }
1865: /* -------------------------------------------------------------- */
1868: PetscErrorCode MatPermute_SeqAIJ(Mat A,IS rowp,IS colp,Mat *B)
1869: {
1870: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1872: PetscInt i,nz,m = A->rmap.n,n = A->cmap.n,*col;
1873: PetscInt *row,*cnew,j,*lens;
1874: IS icolp,irowp;
1875: PetscInt *cwork;
1876: PetscScalar *vwork;
1879: ISInvertPermutation(rowp,PETSC_DECIDE,&irowp);
1880: ISGetIndices(irowp,&row);
1881: ISInvertPermutation(colp,PETSC_DECIDE,&icolp);
1882: ISGetIndices(icolp,&col);
1883:
1884: /* determine lengths of permuted rows */
1885: PetscMalloc((m+1)*sizeof(PetscInt),&lens);
1886: for (i=0; i<m; i++) {
1887: lens[row[i]] = a->i[i+1] - a->i[i];
1888: }
1889: MatCreate(A->comm,B);
1890: MatSetSizes(*B,m,n,m,n);
1891: MatSetType(*B,A->type_name);
1892: MatSeqAIJSetPreallocation_SeqAIJ(*B,0,lens);
1893: PetscFree(lens);
1895: PetscMalloc(n*sizeof(PetscInt),&cnew);
1896: for (i=0; i<m; i++) {
1897: MatGetRow_SeqAIJ(A,i,&nz,&cwork,&vwork);
1898: for (j=0; j<nz; j++) { cnew[j] = col[cwork[j]];}
1899: MatSetValues_SeqAIJ(*B,1,&row[i],nz,cnew,vwork,INSERT_VALUES);
1900: MatRestoreRow_SeqAIJ(A,i,&nz,&cwork,&vwork);
1901: }
1902: PetscFree(cnew);
1903: (*B)->assembled = PETSC_FALSE;
1904: MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);
1905: MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);
1906: ISRestoreIndices(irowp,&row);
1907: ISRestoreIndices(icolp,&col);
1908: ISDestroy(irowp);
1909: ISDestroy(icolp);
1910: return(0);
1911: }
1915: PetscErrorCode MatCopy_SeqAIJ(Mat A,Mat B,MatStructure str)
1916: {
1920: /* If the two matrices have the same copy implementation, use fast copy. */
1921: if (str == SAME_NONZERO_PATTERN && (A->ops->copy == B->ops->copy)) {
1922: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1923: Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data;
1925: if (a->i[A->rmap.n] != b->i[B->rmap.n]) {
1926: SETERRQ(PETSC_ERR_ARG_INCOMP,"Number of nonzeros in two matrices are different");
1927: }
1928: PetscMemcpy(b->a,a->a,(a->i[A->rmap.n])*sizeof(PetscScalar));
1929: } else {
1930: MatCopy_Basic(A,B,str);
1931: }
1932: return(0);
1933: }
1937: PetscErrorCode MatSetUpPreallocation_SeqAIJ(Mat A)
1938: {
1942: MatSeqAIJSetPreallocation_SeqAIJ(A,PETSC_DEFAULT,0);
1943: return(0);
1944: }
1948: PetscErrorCode MatGetArray_SeqAIJ(Mat A,PetscScalar *array[])
1949: {
1950: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1952: *array = a->a;
1953: return(0);
1954: }
1958: PetscErrorCode MatRestoreArray_SeqAIJ(Mat A,PetscScalar *array[])
1959: {
1961: return(0);
1962: }
1966: PetscErrorCode MatFDColoringApply_SeqAIJ(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
1967: {
1968: PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void *))coloring->f;
1970: PetscInt k,N,start,end,l,row,col,srow,**vscaleforrow,m1,m2;
1971: PetscScalar dx,*y,*xx,*w3_array;
1972: PetscScalar *vscale_array;
1973: PetscReal epsilon = coloring->error_rel,umin = coloring->umin;
1974: Vec w1,w2,w3;
1975: void *fctx = coloring->fctx;
1976: PetscTruth flg;
1979: if (!coloring->w1) {
1980: VecDuplicate(x1,&coloring->w1);
1981: PetscLogObjectParent(coloring,coloring->w1);
1982: VecDuplicate(x1,&coloring->w2);
1983: PetscLogObjectParent(coloring,coloring->w2);
1984: VecDuplicate(x1,&coloring->w3);
1985: PetscLogObjectParent(coloring,coloring->w3);
1986: }
1987: w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
1989: MatSetUnfactored(J);
1990: PetscOptionsHasName(coloring->prefix,"-mat_fd_coloring_dont_rezero",&flg);
1991: if (flg) {
1992: PetscInfo(coloring,"Not calling MatZeroEntries()\n");
1993: } else {
1994: PetscTruth assembled;
1995: MatAssembled(J,&assembled);
1996: if (assembled) {
1997: MatZeroEntries(J);
1998: }
1999: }
2001: VecGetOwnershipRange(x1,&start,&end);
2002: VecGetSize(x1,&N);
2004: /*
2005: This is a horrible, horrible, hack. See DMMGComputeJacobian_Multigrid() it inproperly sets
2006: coloring->F for the coarser grids from the finest
2007: */
2008: if (coloring->F) {
2009: VecGetLocalSize(coloring->F,&m1);
2010: VecGetLocalSize(w1,&m2);
2011: if (m1 != m2) {
2012: coloring->F = 0;
2013: }
2014: }
2016: if (coloring->F) {
2017: w1 = coloring->F;
2018: coloring->F = 0;
2019: } else {
2021: (*f)(sctx,x1,w1,fctx);
2023: }
2025: /*
2026: Compute all the scale factors and share with other processors
2027: */
2028: VecGetArray(x1,&xx);xx = xx - start;
2029: VecGetArray(coloring->vscale,&vscale_array);vscale_array = vscale_array - start;
2030: for (k=0; k<coloring->ncolors; k++) {
2031: /*
2032: Loop over each column associated with color adding the
2033: perturbation to the vector w3.
2034: */
2035: for (l=0; l<coloring->ncolumns[k]; l++) {
2036: col = coloring->columns[k][l]; /* column of the matrix we are probing for */
2037: dx = xx[col];
2038: if (dx == 0.0) dx = 1.0;
2039: #if !defined(PETSC_USE_COMPLEX)
2040: if (dx < umin && dx >= 0.0) dx = umin;
2041: else if (dx < 0.0 && dx > -umin) dx = -umin;
2042: #else
2043: if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin;
2044: else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
2045: #endif
2046: dx *= epsilon;
2047: vscale_array[col] = 1.0/dx;
2048: }
2049: }
2050: vscale_array = vscale_array + start;VecRestoreArray(coloring->vscale,&vscale_array);
2051: VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);
2052: VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);
2054: /* VecView(coloring->vscale,PETSC_VIEWER_STDOUT_WORLD);
2055: VecView(x1,PETSC_VIEWER_STDOUT_WORLD);*/
2057: if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
2058: else vscaleforrow = coloring->columnsforrow;
2060: VecGetArray(coloring->vscale,&vscale_array);
2061: /*
2062: Loop over each color
2063: */
2064: for (k=0; k<coloring->ncolors; k++) {
2065: coloring->currentcolor = k;
2066: VecCopy(x1,w3);
2067: VecGetArray(w3,&w3_array);w3_array = w3_array - start;
2068: /*
2069: Loop over each column associated with color adding the
2070: perturbation to the vector w3.
2071: */
2072: for (l=0; l<coloring->ncolumns[k]; l++) {
2073: col = coloring->columns[k][l]; /* column of the matrix we are probing for */
2074: dx = xx[col];
2075: if (dx == 0.0) dx = 1.0;
2076: #if !defined(PETSC_USE_COMPLEX)
2077: if (dx < umin && dx >= 0.0) dx = umin;
2078: else if (dx < 0.0 && dx > -umin) dx = -umin;
2079: #else
2080: if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin;
2081: else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
2082: #endif
2083: dx *= epsilon;
2084: if (!PetscAbsScalar(dx)) SETERRQ(PETSC_ERR_PLIB,"Computed 0 differencing parameter");
2085: w3_array[col] += dx;
2086: }
2087: w3_array = w3_array + start; VecRestoreArray(w3,&w3_array);
2089: /*
2090: Evaluate function at x1 + dx (here dx is a vector of perturbations)
2091: */
2094: (*f)(sctx,w3,w2,fctx);
2096: VecAXPY(w2,-1.0,w1);
2098: /*
2099: Loop over rows of vector, putting results into Jacobian matrix
2100: */
2101: VecGetArray(w2,&y);
2102: for (l=0; l<coloring->nrows[k]; l++) {
2103: row = coloring->rows[k][l];
2104: col = coloring->columnsforrow[k][l];
2105: y[row] *= vscale_array[vscaleforrow[k][l]];
2106: srow = row + start;
2107: MatSetValues_SeqAIJ(J,1,&srow,1,&col,y+row,INSERT_VALUES);
2108: }
2109: VecRestoreArray(w2,&y);
2110: }
2111: coloring->currentcolor = k;
2112: VecRestoreArray(coloring->vscale,&vscale_array);
2113: xx = xx + start; VecRestoreArray(x1,&xx);
2114: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
2115: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
2116: return(0);
2117: }
2119: #include petscblaslapack.h
2122: PetscErrorCode MatAXPY_SeqAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str)
2123: {
2125: PetscInt i;
2126: Mat_SeqAIJ *x = (Mat_SeqAIJ *)X->data,*y = (Mat_SeqAIJ *)Y->data;
2127: PetscBLASInt one=1,bnz = (PetscBLASInt)x->nz;
2130: if (str == SAME_NONZERO_PATTERN) {
2131: PetscScalar alpha = a;
2132: BLASaxpy_(&bnz,&alpha,x->a,&one,y->a,&one);
2133: } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
2134: if (y->xtoy && y->XtoY != X) {
2135: PetscFree(y->xtoy);
2136: MatDestroy(y->XtoY);
2137: }
2138: if (!y->xtoy) { /* get xtoy */
2139: MatAXPYGetxtoy_Private(X->rmap.n,x->i,x->j,PETSC_NULL, y->i,y->j,PETSC_NULL, &y->xtoy);
2140: y->XtoY = X;
2141: PetscObjectReference((PetscObject)X);
2142: }
2143: for (i=0; i<x->nz; i++) y->a[y->xtoy[i]] += a*(x->a[i]);
2144: PetscInfo3(0,"ratio of nnz(X)/nnz(Y): %d/%d = %G\n",x->nz,y->nz,(PetscReal)(x->nz)/y->nz);
2145: } else {
2146: MatAXPY_Basic(Y,a,X,str);
2147: }
2148: return(0);
2149: }
2153: PetscErrorCode MatSetBlockSize_SeqAIJ(Mat A,PetscInt bs)
2154: {
2156: return(0);
2157: }
2161: PetscErrorCode MatConjugate_SeqAIJ(Mat mat)
2162: {
2163: #if defined(PETSC_USE_COMPLEX)
2164: Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
2165: PetscInt i,nz;
2166: PetscScalar *a;
2169: nz = aij->nz;
2170: a = aij->a;
2171: for (i=0; i<nz; i++) {
2172: a[i] = PetscConj(a[i]);
2173: }
2174: #else
2176: #endif
2177: return(0);
2178: }
2182: PetscErrorCode MatGetRowMaxAbs_SeqAIJ(Mat A,Vec v,PetscInt idx[])
2183: {
2184: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
2186: PetscInt i,j,m = A->rmap.n,*ai,*aj,ncols,n;
2187: PetscReal atmp;
2188: PetscScalar *x;
2189: MatScalar *aa;
2192: if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2193: aa = a->a;
2194: ai = a->i;
2195: aj = a->j;
2197: VecSet(v,0.0);
2198: VecGetArray(v,&x);
2199: VecGetLocalSize(v,&n);
2200: if (n != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
2201: for (i=0; i<m; i++) {
2202: ncols = ai[1] - ai[0]; ai++;
2203: x[i] = 0.0; if (idx) idx[i] = 0;
2204: for (j=0; j<ncols; j++){
2205: atmp = PetscAbsScalar(*aa);
2206: if (PetscAbsScalar(x[i]) < atmp) {x[i] = atmp; if (idx) idx[i] = *aj;}
2207: aa++; aj++;
2208: }
2209: }
2210: VecRestoreArray(v,&x);
2211: return(0);
2212: }
2216: PetscErrorCode MatGetRowMax_SeqAIJ(Mat A,Vec v,PetscInt idx[])
2217: {
2218: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
2220: PetscInt i,j,m = A->rmap.n,*ai,*aj,ncols,n;
2221: PetscScalar *x;
2222: MatScalar *aa;
2225: if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2226: aa = a->a;
2227: ai = a->i;
2228: aj = a->j;
2230: VecSet(v,0.0);
2231: VecGetArray(v,&x);
2232: VecGetLocalSize(v,&n);
2233: if (n != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
2234: for (i=0; i<m; i++) {
2235: ncols = ai[1] - ai[0]; ai++;
2236: if (ncols == A->cmap.n) { /* row is dense */
2237: x[i] = *aa; if (idx) idx[i] = 0;
2238: } else { /* row is sparse so already KNOW maximum is 0.0 or higher */
2239: x[i] = 0.0;
2240: if (idx) {
2241: idx[i] = 0; /* in case ncols is zero */
2242: for (j=0;j<ncols;j++) { /* find first implicit 0.0 in the row */
2243: if (aj[j] > j) {
2244: idx[i] = j;
2245: break;
2246: }
2247: }
2248: }
2249: }
2250: for (j=0; j<ncols; j++){
2251: if (PetscRealPart(x[i]) < PetscRealPart(*aa)) {x[i] = *aa; if (idx) idx[i] = *aj;}
2252: aa++; aj++;
2253: }
2254: }
2255: VecRestoreArray(v,&x);
2256: return(0);
2257: }
2261: PetscErrorCode MatGetRowMin_SeqAIJ(Mat A,Vec v,PetscInt idx[])
2262: {
2263: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
2265: PetscInt i,j,m = A->rmap.n,*ai,*aj,ncols,n;
2266: PetscScalar *x;
2267: MatScalar *aa;
2270: if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2271: aa = a->a;
2272: ai = a->i;
2273: aj = a->j;
2275: VecSet(v,0.0);
2276: VecGetArray(v,&x);
2277: VecGetLocalSize(v,&n);
2278: if (n != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
2279: for (i=0; i<m; i++) {
2280: ncols = ai[1] - ai[0]; ai++;
2281: if (ncols == A->cmap.n) { /* row is dense */
2282: x[i] = *aa; if (idx) idx[i] = 0;
2283: } else { /* row is sparse so already KNOW minimum is 0.0 or lower */
2284: x[i] = 0.0;
2285: if (idx) { /* find first implicit 0.0 in the row */
2286: idx[i] = 0; /* in case ncols is zero */
2287: for (j=0;j<ncols;j++) {
2288: if (aj[j] > j) {
2289: idx[i] = j;
2290: break;
2291: }
2292: }
2293: }
2294: }
2295: for (j=0; j<ncols; j++){
2296: if (PetscRealPart(x[i]) > PetscRealPart(*aa)) {x[i] = *aa; if (idx) idx[i] = *aj;}
2297: aa++; aj++;
2298: }
2299: }
2300: VecRestoreArray(v,&x);
2301: return(0);
2302: }
2304: /* -------------------------------------------------------------------*/
2305: static struct _MatOps MatOps_Values = {MatSetValues_SeqAIJ,
2306: MatGetRow_SeqAIJ,
2307: MatRestoreRow_SeqAIJ,
2308: MatMult_SeqAIJ,
2309: /* 4*/ MatMultAdd_SeqAIJ,
2310: MatMultTranspose_SeqAIJ,
2311: MatMultTransposeAdd_SeqAIJ,
2312: MatSolve_SeqAIJ,
2313: MatSolveAdd_SeqAIJ,
2314: MatSolveTranspose_SeqAIJ,
2315: /*10*/ MatSolveTransposeAdd_SeqAIJ,
2316: MatLUFactor_SeqAIJ,
2317: 0,
2318: MatRelax_SeqAIJ,
2319: MatTranspose_SeqAIJ,
2320: /*15*/ MatGetInfo_SeqAIJ,
2321: MatEqual_SeqAIJ,
2322: MatGetDiagonal_SeqAIJ,
2323: MatDiagonalScale_SeqAIJ,
2324: MatNorm_SeqAIJ,
2325: /*20*/ 0,
2326: MatAssemblyEnd_SeqAIJ,
2327: MatCompress_SeqAIJ,
2328: MatSetOption_SeqAIJ,
2329: MatZeroEntries_SeqAIJ,
2330: /*25*/ MatZeroRows_SeqAIJ,
2331: MatLUFactorSymbolic_SeqAIJ,
2332: MatLUFactorNumeric_SeqAIJ,
2333: MatCholeskyFactorSymbolic_SeqAIJ,
2334: MatCholeskyFactorNumeric_SeqAIJ,
2335: /*30*/ MatSetUpPreallocation_SeqAIJ,
2336: MatILUFactorSymbolic_SeqAIJ,
2337: MatICCFactorSymbolic_SeqAIJ,
2338: MatGetArray_SeqAIJ,
2339: MatRestoreArray_SeqAIJ,
2340: /*35*/ MatDuplicate_SeqAIJ,
2341: 0,
2342: 0,
2343: MatILUFactor_SeqAIJ,
2344: 0,
2345: /*40*/ MatAXPY_SeqAIJ,
2346: MatGetSubMatrices_SeqAIJ,
2347: MatIncreaseOverlap_SeqAIJ,
2348: MatGetValues_SeqAIJ,
2349: MatCopy_SeqAIJ,
2350: /*45*/ MatGetRowMax_SeqAIJ,
2351: MatScale_SeqAIJ,
2352: 0,
2353: MatDiagonalSet_SeqAIJ,
2354: MatILUDTFactor_SeqAIJ,
2355: /*50*/ MatSetBlockSize_SeqAIJ,
2356: MatGetRowIJ_SeqAIJ,
2357: MatRestoreRowIJ_SeqAIJ,
2358: MatGetColumnIJ_SeqAIJ,
2359: MatRestoreColumnIJ_SeqAIJ,
2360: /*55*/ MatFDColoringCreate_SeqAIJ,
2361: 0,
2362: 0,
2363: MatPermute_SeqAIJ,
2364: 0,
2365: /*60*/ 0,
2366: MatDestroy_SeqAIJ,
2367: MatView_SeqAIJ,
2368: 0,
2369: 0,
2370: /*65*/ 0,
2371: 0,
2372: 0,
2373: 0,
2374: 0,
2375: /*70*/ MatGetRowMaxAbs_SeqAIJ,
2376: 0,
2377: MatSetColoring_SeqAIJ,
2378: #if defined(PETSC_HAVE_ADIC)
2379: MatSetValuesAdic_SeqAIJ,
2380: #else
2381: 0,
2382: #endif
2383: MatSetValuesAdifor_SeqAIJ,
2384: /*75*/ 0,
2385: 0,
2386: 0,
2387: 0,
2388: 0,
2389: /*80*/ 0,
2390: 0,
2391: 0,
2392: 0,
2393: MatLoad_SeqAIJ,
2394: /*85*/ MatIsSymmetric_SeqAIJ,
2395: 0,
2396: 0,
2397: 0,
2398: 0,
2399: /*90*/ MatMatMult_SeqAIJ_SeqAIJ,
2400: MatMatMultSymbolic_SeqAIJ_SeqAIJ,
2401: MatMatMultNumeric_SeqAIJ_SeqAIJ,
2402: MatPtAP_Basic,
2403: MatPtAPSymbolic_SeqAIJ,
2404: /*95*/ MatPtAPNumeric_SeqAIJ,
2405: MatMatMultTranspose_SeqAIJ_SeqAIJ,
2406: MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ,
2407: MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ,
2408: MatPtAPSymbolic_SeqAIJ_SeqAIJ,
2409: /*100*/MatPtAPNumeric_SeqAIJ_SeqAIJ,
2410: 0,
2411: 0,
2412: MatConjugate_SeqAIJ,
2413: 0,
2414: /*105*/MatSetValuesRow_SeqAIJ,
2415: MatRealPart_SeqAIJ,
2416: MatImaginaryPart_SeqAIJ,
2417: 0,
2418: 0,
2419: /*110*/MatMatSolve_SeqAIJ,
2420: 0,
2421: MatGetRowMin_SeqAIJ
2422: };
2427: PetscErrorCode MatSeqAIJSetColumnIndices_SeqAIJ(Mat mat,PetscInt *indices)
2428: {
2429: Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
2430: PetscInt i,nz,n;
2434: nz = aij->maxnz;
2435: n = mat->rmap.n;
2436: for (i=0; i<nz; i++) {
2437: aij->j[i] = indices[i];
2438: }
2439: aij->nz = nz;
2440: for (i=0; i<n; i++) {
2441: aij->ilen[i] = aij->imax[i];
2442: }
2444: return(0);
2445: }
2450: /*@
2451: MatSeqAIJSetColumnIndices - Set the column indices for all the rows
2452: in the matrix.
2454: Input Parameters:
2455: + mat - the SeqAIJ matrix
2456: - indices - the column indices
2458: Level: advanced
2460: Notes:
2461: This can be called if you have precomputed the nonzero structure of the
2462: matrix and want to provide it to the matrix object to improve the performance
2463: of the MatSetValues() operation.
2465: You MUST have set the correct numbers of nonzeros per row in the call to
2466: MatCreateSeqAIJ(), and the columns indices MUST be sorted.
2468: MUST be called before any calls to MatSetValues();
2470: The indices should start with zero, not one.
2472: @*/
2473: PetscErrorCode MatSeqAIJSetColumnIndices(Mat mat,PetscInt *indices)
2474: {
2475: PetscErrorCode ierr,(*f)(Mat,PetscInt *);
2480: PetscObjectQueryFunction((PetscObject)mat,"MatSeqAIJSetColumnIndices_C",(void (**)(void))&f);
2481: if (f) {
2482: (*f)(mat,indices);
2483: } else {
2484: SETERRQ(PETSC_ERR_SUP,"Wrong type of matrix to set column indices");
2485: }
2486: return(0);
2487: }
2489: /* ----------------------------------------------------------------------------------------*/
2494: PetscErrorCode MatStoreValues_SeqAIJ(Mat mat)
2495: {
2496: Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
2498: size_t nz = aij->i[mat->rmap.n];
2501: if (aij->nonew != 1) {
2502: SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
2503: }
2505: /* allocate space for values if not already there */
2506: if (!aij->saved_values) {
2507: PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);
2508: PetscLogObjectMemory(mat,(nz+1)*sizeof(PetscScalar));
2509: }
2511: /* copy values over */
2512: PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));
2513: return(0);
2514: }
2519: /*@
2520: MatStoreValues - Stashes a copy of the matrix values; this allows, for
2521: example, reuse of the linear part of a Jacobian, while recomputing the
2522: nonlinear portion.
2524: Collect on Mat
2526: Input Parameters:
2527: . mat - the matrix (currently only AIJ matrices support this option)
2529: Level: advanced
2531: Common Usage, with SNESSolve():
2532: $ Create Jacobian matrix
2533: $ Set linear terms into matrix
2534: $ Apply boundary conditions to matrix, at this time matrix must have
2535: $ final nonzero structure (i.e. setting the nonlinear terms and applying
2536: $ boundary conditions again will not change the nonzero structure
2537: $ MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS);
2538: $ MatStoreValues(mat);
2539: $ Call SNESSetJacobian() with matrix
2540: $ In your Jacobian routine
2541: $ MatRetrieveValues(mat);
2542: $ Set nonlinear terms in matrix
2543:
2544: Common Usage without SNESSolve(), i.e. when you handle nonlinear solve yourself:
2545: $ // build linear portion of Jacobian
2546: $ MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS);
2547: $ MatStoreValues(mat);
2548: $ loop over nonlinear iterations
2549: $ MatRetrieveValues(mat);
2550: $ // call MatSetValues(mat,...) to set nonliner portion of Jacobian
2551: $ // call MatAssemblyBegin/End() on matrix
2552: $ Solve linear system with Jacobian
2553: $ endloop
2555: Notes:
2556: Matrix must already be assemblied before calling this routine
2557: Must set the matrix option MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS); before
2558: calling this routine.
2560: When this is called multiple times it overwrites the previous set of stored values
2561: and does not allocated additional space.
2563: .seealso: MatRetrieveValues()
2565: @*/
2566: PetscErrorCode MatStoreValues(Mat mat)
2567: {
2568: PetscErrorCode ierr,(*f)(Mat);
2572: if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2573: if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2575: PetscObjectQueryFunction((PetscObject)mat,"MatStoreValues_C",(void (**)(void))&f);
2576: if (f) {
2577: (*f)(mat);
2578: } else {
2579: SETERRQ(PETSC_ERR_SUP,"Wrong type of matrix to store values");
2580: }
2581: return(0);
2582: }
2587: PetscErrorCode MatRetrieveValues_SeqAIJ(Mat mat)
2588: {
2589: Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
2591: PetscInt nz = aij->i[mat->rmap.n];
2594: if (aij->nonew != 1) {
2595: SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
2596: }
2597: if (!aij->saved_values) {
2598: SETERRQ(PETSC_ERR_ORDER,"Must call MatStoreValues(A);first");
2599: }
2600: /* copy values over */
2601: PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));
2602: return(0);
2603: }
2608: /*@
2609: MatRetrieveValues - Retrieves the copy of the matrix values; this allows, for
2610: example, reuse of the linear part of a Jacobian, while recomputing the
2611: nonlinear portion.
2613: Collect on Mat
2615: Input Parameters:
2616: . mat - the matrix (currently on AIJ matrices support this option)
2618: Level: advanced
2620: .seealso: MatStoreValues()
2622: @*/
2623: PetscErrorCode MatRetrieveValues(Mat mat)
2624: {
2625: PetscErrorCode ierr,(*f)(Mat);
2629: if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2630: if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2632: PetscObjectQueryFunction((PetscObject)mat,"MatRetrieveValues_C",(void (**)(void))&f);
2633: if (f) {
2634: (*f)(mat);
2635: } else {
2636: SETERRQ(PETSC_ERR_SUP,"Wrong type of matrix to retrieve values");
2637: }
2638: return(0);
2639: }
2642: /* --------------------------------------------------------------------------------*/
2645: /*@C
2646: MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format
2647: (the default parallel PETSc format). For good matrix assembly performance
2648: the user should preallocate the matrix storage by setting the parameter nz
2649: (or the array nnz). By setting these parameters accurately, performance
2650: during matrix assembly can be increased by more than a factor of 50.
2652: Collective on MPI_Comm
2654: Input Parameters:
2655: + comm - MPI communicator, set to PETSC_COMM_SELF
2656: . m - number of rows
2657: . n - number of columns
2658: . nz - number of nonzeros per row (same for all rows)
2659: - nnz - array containing the number of nonzeros in the various rows
2660: (possibly different for each row) or PETSC_NULL
2662: Output Parameter:
2663: . A - the matrix
2665: Notes:
2666: If nnz is given then nz is ignored
2668: The AIJ format (also called the Yale sparse matrix format or
2669: compressed row storage), is fully compatible with standard Fortran 77
2670: storage. That is, the stored row and column indices can begin at
2671: either one (as in Fortran) or zero. See the users' manual for details.
2673: Specify the preallocated storage with either nz or nnz (not both).
2674: Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
2675: allocation. For large problems you MUST preallocate memory or you
2676: will get TERRIBLE performance, see the users' manual chapter on matrices.
2678: By default, this format uses inodes (identical nodes) when possible, to
2679: improve numerical efficiency of matrix-vector products and solves. We
2680: search for consecutive rows with the same nonzero structure, thereby
2681: reusing matrix information to achieve increased efficiency.
2683: Options Database Keys:
2684: + -mat_no_inode - Do not use inodes
2685: . -mat_inode_limit <limit> - Sets inode limit (max limit=5)
2686: - -mat_aij_oneindex - Internally use indexing starting at 1
2687: rather than 0. Note that when calling MatSetValues(),
2688: the user still MUST index entries starting at 0!
2690: Level: intermediate
2692: .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues(), MatSeqAIJSetColumnIndices(), MatCreateSeqAIJWithArrays()
2694: @*/
2695: PetscErrorCode MatCreateSeqAIJ(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A)
2696: {
2700: MatCreate(comm,A);
2701: MatSetSizes(*A,m,n,m,n);
2702: MatSetType(*A,MATSEQAIJ);
2703: MatSeqAIJSetPreallocation_SeqAIJ(*A,nz,(PetscInt*)nnz);
2704: return(0);
2705: }
2709: /*@C
2710: MatSeqAIJSetPreallocation - For good matrix assembly performance
2711: the user should preallocate the matrix storage by setting the parameter nz
2712: (or the array nnz). By setting these parameters accurately, performance
2713: during matrix assembly can be increased by more than a factor of 50.
2715: Collective on MPI_Comm
2717: Input Parameters:
2718: + B - The matrix
2719: . nz - number of nonzeros per row (same for all rows)
2720: - nnz - array containing the number of nonzeros in the various rows
2721: (possibly different for each row) or PETSC_NULL
2723: Notes:
2724: If nnz is given then nz is ignored
2726: The AIJ format (also called the Yale sparse matrix format or
2727: compressed row storage), is fully compatible with standard Fortran 77
2728: storage. That is, the stored row and column indices can begin at
2729: either one (as in Fortran) or zero. See the users' manual for details.
2731: Specify the preallocated storage with either nz or nnz (not both).
2732: Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
2733: allocation. For large problems you MUST preallocate memory or you
2734: will get TERRIBLE performance, see the users' manual chapter on matrices.
2736: Developers: Use nz of MAT_SKIP_ALLOCATION to not allocate any space for the matrix
2737: entries or columns indices
2739: By default, this format uses inodes (identical nodes) when possible, to
2740: improve numerical efficiency of matrix-vector products and solves. We
2741: search for consecutive rows with the same nonzero structure, thereby
2742: reusing matrix information to achieve increased efficiency.
2744: Options Database Keys:
2745: + -mat_no_inode - Do not use inodes
2746: . -mat_inode_limit <limit> - Sets inode limit (max limit=5)
2747: - -mat_aij_oneindex - Internally use indexing starting at 1
2748: rather than 0. Note that when calling MatSetValues(),
2749: the user still MUST index entries starting at 0!
2751: Level: intermediate
2753: .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues(), MatSeqAIJSetColumnIndices(), MatCreateSeqAIJWithArrays()
2755: @*/
2756: PetscErrorCode MatSeqAIJSetPreallocation(Mat B,PetscInt nz,const PetscInt nnz[])
2757: {
2758: PetscErrorCode ierr,(*f)(Mat,PetscInt,const PetscInt[]);
2761: PetscObjectQueryFunction((PetscObject)B,"MatSeqAIJSetPreallocation_C",(void (**)(void))&f);
2762: if (f) {
2763: (*f)(B,nz,nnz);
2764: }
2765: return(0);
2766: }
2771: PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat B,PetscInt nz,PetscInt *nnz)
2772: {
2773: Mat_SeqAIJ *b;
2774: PetscTruth skipallocation = PETSC_FALSE;
2776: PetscInt i;
2779:
2780: if (nz == MAT_SKIP_ALLOCATION) {
2781: skipallocation = PETSC_TRUE;
2782: nz = 0;
2783: }
2785: B->rmap.bs = B->cmap.bs = 1;
2786: PetscMapSetUp(&B->rmap);
2787: PetscMapSetUp(&B->cmap);
2789: if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
2790: if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz);
2791: if (nnz) {
2792: for (i=0; i<B->rmap.n; i++) {
2793: if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]);
2794: 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);
2795: }
2796: }
2798: B->preallocated = PETSC_TRUE;
2799: b = (Mat_SeqAIJ*)B->data;
2801: if (!skipallocation) {
2802: PetscMalloc2(B->rmap.n,PetscInt,&b->imax,B->rmap.n,PetscInt,&b->ilen);
2803: PetscLogObjectMemory(B,2*B->rmap.n*sizeof(PetscInt));
2804: if (!nnz) {
2805: if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 10;
2806: else if (nz <= 0) nz = 1;
2807: for (i=0; i<B->rmap.n; i++) b->imax[i] = nz;
2808: nz = nz*B->rmap.n;
2809: } else {
2810: nz = 0;
2811: for (i=0; i<B->rmap.n; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
2812: }
2814: /* b->ilen will count nonzeros in each row so far. */
2815: for (i=0; i<B->rmap.n; i++) { b->ilen[i] = 0;}
2817: /* allocate the matrix space */
2818: PetscMalloc3(nz,PetscScalar,&b->a,nz,PetscInt,&b->j,B->rmap.n+1,PetscInt,&b->i);
2819: PetscLogObjectMemory(B,(B->rmap.n+1)*sizeof(PetscInt)+nz*(sizeof(PetscScalar)+sizeof(PetscInt)));
2820: b->i[0] = 0;
2821: for (i=1; i<B->rmap.n+1; i++) {
2822: b->i[i] = b->i[i-1] + b->imax[i-1];
2823: }
2824: b->singlemalloc = PETSC_TRUE;
2825: b->free_a = PETSC_TRUE;
2826: b->free_ij = PETSC_TRUE;
2827: } else {
2828: b->free_a = PETSC_FALSE;
2829: b->free_ij = PETSC_FALSE;
2830: }
2832: b->nz = 0;
2833: b->maxnz = nz;
2834: B->info.nz_unneeded = (double)b->maxnz;
2835: return(0);
2836: }
2839: #undef __FUNCT__
2841: /*@C
2842: MatSeqAIJSetPreallocationCSR - Allocates memory for a sparse sequential matrix in AIJ format.
2844: Input Parameters:
2845: + B - the matrix
2846: . i - the indices into j for the start of each row (starts with zero)
2847: . j - the column indices for each row (starts with zero) these must be sorted for each row
2848: - v - optional values in the matrix
2850: Contributed by: Lisandro Dalchin
2852: Level: developer
2854: .keywords: matrix, aij, compressed row, sparse, sequential
2856: .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatSeqAIJSetPreallocation(), MatCreateSeqAIJ(), SeqAIJ
2857: @*/
2858: PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat B,const PetscInt i[],const PetscInt j[],const PetscScalar v[])
2859: {
2860: PetscErrorCode (*f)(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
2865: PetscObjectQueryFunction((PetscObject)B,"MatSeqAIJSetPreallocationCSR_C",(void (**)(void))&f);
2866: if (f) {
2867: (*f)(B,i,j,v);
2868: }
2869: return(0);
2870: }
2873: #undef __FUNCT__
2875: PetscErrorCode MatSeqAIJSetPreallocationCSR_SeqAIJ(Mat B,const PetscInt Ii[],const PetscInt J[],const PetscScalar v[])
2876: {
2877: PetscInt i;
2878: PetscInt m,n;
2879: PetscInt nz;
2880: PetscInt *nnz, nz_max = 0;
2881: PetscScalar *values;
2885: MatGetSize(B, &m, &n);
2887: if (Ii[0]) {
2888: SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Ii[0] must be 0 it is %D", Ii[0]);
2889: }
2890: PetscMalloc((m+1) * sizeof(PetscInt), &nnz);
2891: for(i = 0; i < m; i++) {
2892: nz = Ii[i+1]- Ii[i];
2893: nz_max = PetscMax(nz_max, nz);
2894: if (nz < 0) {
2895: SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE, "Local row %D has a negative number of columns %D", i, nnz);
2896: }
2897: nnz[i] = nz;
2898: }
2899: MatSeqAIJSetPreallocation(B, 0, nnz);
2900: PetscFree(nnz);
2902: if (v) {
2903: values = (PetscScalar*) v;
2904: } else {
2905: PetscMalloc((nz_max+1)*sizeof(PetscScalar), &values);
2906: PetscMemzero(values, nz_max*sizeof(PetscScalar));
2907: }
2909: MatSetOption(B,MAT_COLUMNS_SORTED);
2911: for(i = 0; i < m; i++) {
2912: nz = Ii[i+1] - Ii[i];
2913: MatSetValues_SeqAIJ(B, 1, &i, nz, J+Ii[i], values + (v ? Ii[i] : 0), INSERT_VALUES);
2914: }
2916: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
2917: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
2918: MatSetOption(B,MAT_COLUMNS_UNSORTED);
2920: if (!v) {
2921: PetscFree(values);
2922: }
2923: return(0);
2924: }
2927: /*MC
2928: MATSEQAIJ - MATSEQAIJ = "seqaij" - A matrix type to be used for sequential sparse matrices,
2929: based on compressed sparse row format.
2931: Options Database Keys:
2932: . -mat_type seqaij - sets the matrix type to "seqaij" during a call to MatSetFromOptions()
2934: Level: beginner
2936: .seealso: MatCreateSeqAIJ(), MatSetFromOptions(), MatSetType(), MatCreate(), MatType
2937: M*/
2946: PetscErrorCode MatCreate_SeqAIJ(Mat B)
2947: {
2948: Mat_SeqAIJ *b;
2950: PetscMPIInt size;
2953: MPI_Comm_size(B->comm,&size);
2954: if (size > 1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Comm must be of size 1");
2956: PetscNew(Mat_SeqAIJ,&b);
2957: B->data = (void*)b;
2958: PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
2959: B->factor = 0;
2960: B->mapping = 0;
2961: b->row = 0;
2962: b->col = 0;
2963: b->icol = 0;
2964: b->reallocs = 0;
2965: b->sorted = PETSC_FALSE;
2966: b->ignorezeroentries = PETSC_FALSE;
2967: b->roworiented = PETSC_TRUE;
2968: b->nonew = 0;
2969: b->diag = 0;
2970: b->solve_work = 0;
2971: B->spptr = 0;
2972: b->saved_values = 0;
2973: b->idiag = 0;
2974: b->ssor = 0;
2975: b->keepzeroedrows = PETSC_FALSE;
2976: b->xtoy = 0;
2977: b->XtoY = 0;
2978: b->compressedrow.use = PETSC_FALSE;
2979: b->compressedrow.nrows = B->rmap.n;
2980: b->compressedrow.i = PETSC_NULL;
2981: b->compressedrow.rindex = PETSC_NULL;
2982: b->compressedrow.checked = PETSC_FALSE;
2983: B->same_nonzero = PETSC_FALSE;
2985: PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);
2986: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqAIJSetColumnIndices_C",
2987: "MatSeqAIJSetColumnIndices_SeqAIJ",
2988: MatSeqAIJSetColumnIndices_SeqAIJ);
2989: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
2990: "MatStoreValues_SeqAIJ",
2991: MatStoreValues_SeqAIJ);
2992: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
2993: "MatRetrieveValues_SeqAIJ",
2994: MatRetrieveValues_SeqAIJ);
2995: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_seqsbaij_C",
2996: "MatConvert_SeqAIJ_SeqSBAIJ",
2997: MatConvert_SeqAIJ_SeqSBAIJ);
2998: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_seqbaij_C",
2999: "MatConvert_SeqAIJ_SeqBAIJ",
3000: MatConvert_SeqAIJ_SeqBAIJ);
3001: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_seqcsrperm_C",
3002: "MatConvert_SeqAIJ_SeqCSRPERM",
3003: MatConvert_SeqAIJ_SeqCSRPERM);
3004: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_seqcrl_C",
3005: "MatConvert_SeqAIJ_SeqCRL",
3006: MatConvert_SeqAIJ_SeqCRL);
3007: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatIsTranspose_C",
3008: "MatIsTranspose_SeqAIJ",
3009: MatIsTranspose_SeqAIJ);
3010: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqAIJSetPreallocation_C",
3011: "MatSeqAIJSetPreallocation_SeqAIJ",
3012: MatSeqAIJSetPreallocation_SeqAIJ);
3013: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqAIJSetPreallocationCSR_C",
3014: "MatSeqAIJSetPreallocationCSR_SeqAIJ",
3015: MatSeqAIJSetPreallocationCSR_SeqAIJ);
3016: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatReorderForNonzeroDiagonal_C",
3017: "MatReorderForNonzeroDiagonal_SeqAIJ",
3018: MatReorderForNonzeroDiagonal_SeqAIJ);
3019: MatCreate_Inode(B);
3020: PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);
3021: return(0);
3022: }
3027: PetscErrorCode MatDuplicate_SeqAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
3028: {
3029: Mat C;
3030: Mat_SeqAIJ *c,*a = (Mat_SeqAIJ*)A->data;
3032: PetscInt i,m = A->rmap.n;
3035: *B = 0;
3036: MatCreate(A->comm,&C);
3037: MatSetSizes(C,A->rmap.n,A->cmap.n,A->rmap.n,A->cmap.n);
3038: MatSetType(C,A->type_name);
3039: PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));
3040:
3041: PetscMapCopy(A->comm,&A->rmap,&C->rmap);
3042: PetscMapCopy(A->comm,&A->cmap,&C->cmap);
3044: c = (Mat_SeqAIJ*)C->data;
3046: C->factor = A->factor;
3048: c->row = 0;
3049: c->col = 0;
3050: c->icol = 0;
3051: c->reallocs = 0;
3053: C->assembled = PETSC_TRUE;
3054:
3055: PetscMalloc2(m,PetscInt,&c->imax,m,PetscInt,&c->ilen);
3056: PetscLogObjectMemory(C, 2*m*sizeof(PetscInt));
3057: for (i=0; i<m; i++) {
3058: c->imax[i] = a->imax[i];
3059: c->ilen[i] = a->ilen[i];
3060: }
3062: /* allocate the matrix space */
3063: PetscMalloc3(a->i[m],PetscScalar,&c->a,a->i[m],PetscInt,&c->j,m+1,PetscInt,&c->i);
3064: PetscLogObjectMemory(C, a->i[m]*(sizeof(PetscScalar)+sizeof(PetscInt))+(m+1)*sizeof(PetscInt));
3065: c->singlemalloc = PETSC_TRUE;
3066: PetscMemcpy(c->i,a->i,(m+1)*sizeof(PetscInt));
3067: if (m > 0) {
3068: PetscMemcpy(c->j,a->j,(a->i[m])*sizeof(PetscInt));
3069: if (cpvalues == MAT_COPY_VALUES) {
3070: PetscMemcpy(c->a,a->a,(a->i[m])*sizeof(PetscScalar));
3071: } else {
3072: PetscMemzero(c->a,(a->i[m])*sizeof(PetscScalar));
3073: }
3074: }
3076: c->sorted = a->sorted;
3077: c->ignorezeroentries = a->ignorezeroentries;
3078: c->roworiented = a->roworiented;
3079: c->nonew = a->nonew;
3080: if (a->diag) {
3081: PetscMalloc((m+1)*sizeof(PetscInt),&c->diag);
3082: PetscLogObjectMemory(C,(m+1)*sizeof(PetscInt));
3083: for (i=0; i<m; i++) {
3084: c->diag[i] = a->diag[i];
3085: }
3086: } else c->diag = 0;
3087: c->solve_work = 0;
3088: c->saved_values = 0;
3089: c->idiag = 0;
3090: c->ssor = 0;
3091: c->keepzeroedrows = a->keepzeroedrows;
3092: c->free_a = PETSC_TRUE;
3093: c->free_ij = PETSC_TRUE;
3094: c->xtoy = 0;
3095: c->XtoY = 0;
3097: c->nz = a->nz;
3098: c->maxnz = a->maxnz;
3099: C->preallocated = PETSC_TRUE;
3101: c->compressedrow.use = a->compressedrow.use;
3102: c->compressedrow.nrows = a->compressedrow.nrows;
3103: c->compressedrow.checked = a->compressedrow.checked;
3104: if ( a->compressedrow.checked && a->compressedrow.use){
3105: i = a->compressedrow.nrows;
3106: PetscMalloc((2*i+1)*sizeof(PetscInt),&c->compressedrow.i);
3107: c->compressedrow.rindex = c->compressedrow.i + i + 1;
3108: PetscMemcpy(c->compressedrow.i,a->compressedrow.i,(i+1)*sizeof(PetscInt));
3109: PetscMemcpy(c->compressedrow.rindex,a->compressedrow.rindex,i*sizeof(PetscInt));
3110: } else {
3111: c->compressedrow.use = PETSC_FALSE;
3112: c->compressedrow.i = PETSC_NULL;
3113: c->compressedrow.rindex = PETSC_NULL;
3114: }
3115: C->same_nonzero = A->same_nonzero;
3116: MatDuplicate_Inode(A,cpvalues,&C);
3118: *B = C;
3119: PetscFListDuplicate(A->qlist,&C->qlist);
3120: return(0);
3121: }
3125: PetscErrorCode MatLoad_SeqAIJ(PetscViewer viewer, MatType type,Mat *A)
3126: {
3127: Mat_SeqAIJ *a;
3128: Mat B;
3130: PetscInt i,sum,nz,header[4],*rowlengths = 0,M,N;
3131: int fd;
3132: PetscMPIInt size;
3133: MPI_Comm comm;
3134:
3136: PetscObjectGetComm((PetscObject)viewer,&comm);
3137: MPI_Comm_size(comm,&size);
3138: if (size > 1) SETERRQ(PETSC_ERR_ARG_SIZ,"view must have one processor");
3139: PetscViewerBinaryGetDescriptor(viewer,&fd);
3140: PetscBinaryRead(fd,header,4,PETSC_INT);
3141: if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object in file");
3142: M = header[1]; N = header[2]; nz = header[3];
3144: if (nz < 0) {
3145: SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format on disk,cannot load as SeqAIJ");
3146: }
3148: /* read in row lengths */
3149: PetscMalloc(M*sizeof(PetscInt),&rowlengths);
3150: PetscBinaryRead(fd,rowlengths,M,PETSC_INT);
3152: /* check if sum of rowlengths is same as nz */
3153: for (i=0,sum=0; i< M; i++) sum +=rowlengths[i];
3154: if (sum != nz) SETERRQ2(PETSC_ERR_FILE_READ,"Inconsistant matrix data in file. no-nonzeros = %d, sum-row-lengths = %d\n",nz,sum);
3156: /* create our matrix */
3157: MatCreate(comm,&B);
3158: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,M,N);
3159: MatSetType(B,type);
3160: MatSeqAIJSetPreallocation_SeqAIJ(B,0,rowlengths);
3161: a = (Mat_SeqAIJ*)B->data;
3163: PetscBinaryRead(fd,a->j,nz,PETSC_INT);
3165: /* read in nonzero values */
3166: PetscBinaryRead(fd,a->a,nz,PETSC_SCALAR);
3168: /* set matrix "i" values */
3169: a->i[0] = 0;
3170: for (i=1; i<= M; i++) {
3171: a->i[i] = a->i[i-1] + rowlengths[i-1];
3172: a->ilen[i-1] = rowlengths[i-1];
3173: }
3174: PetscFree(rowlengths);
3176: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
3177: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
3178: *A = B;
3179: return(0);
3180: }
3184: PetscErrorCode MatEqual_SeqAIJ(Mat A,Mat B,PetscTruth* flg)
3185: {
3186: Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data,*b = (Mat_SeqAIJ *)B->data;
3190: /* If the matrix dimensions are not equal,or no of nonzeros */
3191: if ((A->rmap.n != B->rmap.n) || (A->cmap.n != B->cmap.n) ||(a->nz != b->nz)) {
3192: *flg = PETSC_FALSE;
3193: return(0);
3194: }
3195:
3196: /* if the a->i are the same */
3197: PetscMemcmp(a->i,b->i,(A->rmap.n+1)*sizeof(PetscInt),flg);
3198: if (!*flg) return(0);
3199:
3200: /* if a->j are the same */
3201: PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),flg);
3202: if (!*flg) return(0);
3203:
3204: /* if a->a are the same */
3205: PetscMemcmp(a->a,b->a,(a->nz)*sizeof(PetscScalar),flg);
3207: return(0);
3208:
3209: }
3213: /*@
3214: MatCreateSeqAIJWithArrays - Creates an sequential AIJ matrix using matrix elements (in CSR format)
3215: provided by the user.
3217: Collective on MPI_Comm
3219: Input Parameters:
3220: + comm - must be an MPI communicator of size 1
3221: . m - number of rows
3222: . n - number of columns
3223: . i - row indices
3224: . j - column indices
3225: - a - matrix values
3227: Output Parameter:
3228: . mat - the matrix
3230: Level: intermediate
3232: Notes:
3233: The i, j, and a arrays are not copied by this routine, the user must free these arrays
3234: once the matrix is destroyed
3236: You cannot set new nonzero locations into this matrix, that will generate an error.
3238: The i and j indices are 0 based
3240: The format which is used for the sparse matrix input, is equivalent to a
3241: row-major ordering.. i.e for the following matrix, the input data expected is
3242: as shown:
3244: 1 0 0
3245: 2 0 3
3246: 4 5 6
3248: i = {0,1,3,6} [size = nrow+1 = 3+1]
3249: j = {0,0,2,0,1,2} [size = nz = 6]
3250: v = {1,2,3,4,5,6} [size = nz = 6]
3252: .seealso: MatCreate(), MatCreateMPIAIJ(), MatCreateSeqAIJ(), MatCreateMPIAIJWithArrays(), MatMPIAIJSetPreallocationCSR()
3254: @*/
3255: PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt* i,PetscInt*j,PetscScalar *a,Mat *mat)
3256: {
3258: PetscInt ii;
3259: Mat_SeqAIJ *aij;
3262: if (i[0]) {
3263: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"i (row indices) must start with 0");
3264: }
3265: MatCreate(comm,mat);
3266: MatSetSizes(*mat,m,n,m,n);
3267: MatSetType(*mat,MATSEQAIJ);
3268: MatSeqAIJSetPreallocation_SeqAIJ(*mat,MAT_SKIP_ALLOCATION,0);
3269: aij = (Mat_SeqAIJ*)(*mat)->data;
3270: PetscMalloc2(m,PetscInt,&aij->imax,m,PetscInt,&aij->ilen);
3272: aij->i = i;
3273: aij->j = j;
3274: aij->a = a;
3275: aij->singlemalloc = PETSC_FALSE;
3276: aij->nonew = -1; /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/
3277: aij->free_a = PETSC_FALSE;
3278: aij->free_ij = PETSC_FALSE;
3280: for (ii=0; ii<m; ii++) {
3281: aij->ilen[ii] = aij->imax[ii] = i[ii+1] - i[ii];
3282: #if defined(PETSC_USE_DEBUG)
3283: 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]);
3284: #endif
3285: }
3286: #if defined(PETSC_USE_DEBUG)
3287: for (ii=0; ii<aij->i[m]; ii++) {
3288: if (j[ii] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Negative column index at location = %d index = %d",ii,j[ii]);
3289: if (j[ii] > n - 1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column index to large at location = %d index = %d",ii,j[ii]);
3290: }
3291: #endif
3293: MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
3294: MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
3295: return(0);
3296: }
3300: PetscErrorCode MatSetColoring_SeqAIJ(Mat A,ISColoring coloring)
3301: {
3303: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
3306: if (coloring->ctype == IS_COLORING_GLOBAL) {
3307: ISColoringReference(coloring);
3308: a->coloring = coloring;
3309: } else if (coloring->ctype == IS_COLORING_GHOSTED) {
3310: PetscInt i,*larray;
3311: ISColoring ocoloring;
3312: ISColoringValue *colors;
3314: /* set coloring for diagonal portion */
3315: PetscMalloc((A->cmap.n+1)*sizeof(PetscInt),&larray);
3316: for (i=0; i<A->cmap.n; i++) {
3317: larray[i] = i;
3318: }
3319: ISGlobalToLocalMappingApply(A->mapping,IS_GTOLM_MASK,A->cmap.n,larray,PETSC_NULL,larray);
3320: PetscMalloc((A->cmap.n+1)*sizeof(ISColoringValue),&colors);
3321: for (i=0; i<A->cmap.n; i++) {
3322: colors[i] = coloring->colors[larray[i]];
3323: }
3324: PetscFree(larray);
3325: ISColoringCreate(PETSC_COMM_SELF,coloring->n,A->cmap.n,colors,&ocoloring);
3326: a->coloring = ocoloring;
3327: }
3328: return(0);
3329: }
3331: #if defined(PETSC_HAVE_ADIC)
3333: #include "adic/ad_utils.h"
3338: PetscErrorCode MatSetValuesAdic_SeqAIJ(Mat A,void *advalues)
3339: {
3340: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
3341: PetscInt m = A->rmap.n,*ii = a->i,*jj = a->j,nz,i,j,nlen;
3342: PetscScalar *v = a->a,*values = ((PetscScalar*)advalues)+1;
3343: ISColoringValue *color;
3346: if (!a->coloring) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Coloring not set for matrix");
3347: nlen = PetscADGetDerivTypeSize()/sizeof(PetscScalar);
3348: color = a->coloring->colors;
3349: /* loop over rows */
3350: for (i=0; i<m; i++) {
3351: nz = ii[i+1] - ii[i];
3352: /* loop over columns putting computed value into matrix */
3353: for (j=0; j<nz; j++) {
3354: *v++ = values[color[*jj++]];
3355: }
3356: values += nlen; /* jump to next row of derivatives */
3357: }
3358: return(0);
3359: }
3360: #endif
3364: PetscErrorCode MatSetValuesAdifor_SeqAIJ(Mat A,PetscInt nl,void *advalues)
3365: {
3366: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
3367: PetscInt m = A->rmap.n,*ii = a->i,*jj = a->j,nz,i,j;
3368: PetscScalar *v = a->a,*values = (PetscScalar *)advalues;
3369: ISColoringValue *color;
3372: if (!a->coloring) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Coloring not set for matrix");
3373: color = a->coloring->colors;
3374: /* loop over rows */
3375: for (i=0; i<m; i++) {
3376: nz = ii[i+1] - ii[i];
3377: /* loop over columns putting computed value into matrix */
3378: for (j=0; j<nz; j++) {
3379: *v++ = values[color[*jj++]];
3380: }
3381: values += nl; /* jump to next row of derivatives */
3382: }
3383: return(0);
3384: }
3386: /*
3387: Special version for direct calls from Fortran
3388: */
3389: #if defined(PETSC_HAVE_FORTRAN_CAPS)
3390: #define matsetvaluesseqaij_ MATSETVALUESSEQAIJ
3391: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
3392: #define matsetvaluesseqaij_ matsetvaluesseqaij
3393: #endif
3395: /* Change these macros so can be used in void function */
3396: #undef CHKERRQ
3397: #define CHKERRQ(ierr) CHKERRABORT(A->comm,ierr)
3398: #undef SETERRQ2
3399: #define SETERRQ2(ierr,b,c,d) CHKERRABORT(A->comm,ierr)
3404: void PETSC_STDCALL matsetvaluesseqaij_(Mat *AA,PetscInt *mm,const PetscInt im[],PetscInt *nn,const PetscInt in[],const PetscScalar v[],InsertMode *isis, PetscErrorCode *_ierr)
3405: {
3406: Mat A = *AA;
3407: PetscInt m = *mm, n = *nn;
3408: InsertMode is = *isis;
3409: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
3410: PetscInt *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N;
3411: PetscInt *imax,*ai,*ailen;
3413: PetscInt *aj,nonew = a->nonew,lastcol = -1;
3414: PetscScalar *ap,value,*aa;
3415: PetscTruth ignorezeroentries = a->ignorezeroentries;
3416: PetscTruth roworiented = a->roworiented;
3419: MatPreallocated(A);
3420: imax = a->imax;
3421: ai = a->i;
3422: ailen = a->ilen;
3423: aj = a->j;
3424: aa = a->a;
3426: for (k=0; k<m; k++) { /* loop over added rows */
3427: row = im[k];
3428: if (row < 0) continue;
3429: #if defined(PETSC_USE_DEBUG)
3430: if (row >= A->rmap.n) SETERRABORT(A->comm,PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
3431: #endif
3432: rp = aj + ai[row]; ap = aa + ai[row];
3433: rmax = imax[row]; nrow = ailen[row];
3434: low = 0;
3435: high = nrow;
3436: for (l=0; l<n; l++) { /* loop over added columns */
3437: if (in[l] < 0) continue;
3438: #if defined(PETSC_USE_DEBUG)
3439: if (in[l] >= A->cmap.n) SETERRABORT(A->comm,PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
3440: #endif
3441: col = in[l];
3442: if (roworiented) {
3443: value = v[l + k*n];
3444: } else {
3445: value = v[k + l*m];
3446: }
3447: if (value == 0.0 && ignorezeroentries && (is == ADD_VALUES)) continue;
3449: if (col <= lastcol) low = 0; else high = nrow;
3450: lastcol = col;
3451: while (high-low > 5) {
3452: t = (low+high)/2;
3453: if (rp[t] > col) high = t;
3454: else low = t;
3455: }
3456: for (i=low; i<high; i++) {
3457: if (rp[i] > col) break;
3458: if (rp[i] == col) {
3459: if (is == ADD_VALUES) ap[i] += value;
3460: else ap[i] = value;
3461: goto noinsert;
3462: }
3463: }
3464: if (value == 0.0 && ignorezeroentries) goto noinsert;
3465: if (nonew == 1) goto noinsert;
3466: if (nonew == -1) SETERRABORT(A->comm,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
3467: MatSeqXAIJReallocateAIJ(A,A->rmap.n,1,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar);
3468: N = nrow++ - 1; a->nz++; high++;
3469: /* shift up all the later entries in this row */
3470: for (ii=N; ii>=i; ii--) {
3471: rp[ii+1] = rp[ii];
3472: ap[ii+1] = ap[ii];
3473: }
3474: rp[i] = col;
3475: ap[i] = value;
3476: noinsert:;
3477: low = i + 1;
3478: }
3479: ailen[row] = nrow;
3480: }
3481: A->same_nonzero = PETSC_FALSE;
3482: PetscFunctionReturnVoid();
3483: }