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: for (i=0; i<A->rmap.n+1; i++) (*ia)[i] = a->i[i] + 1;
65: if (ja) {
66: PetscMalloc((nz+1)*sizeof(PetscInt),ja);
67: for (i=0; i<nz; i++) (*ja)[i] = a->j[i] + 1;
68: }
69: } else {
70: *ia = a->i;
71: if (ja) *ja = a->j;
72: }
73: return(0);
74: }
78: PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscTruth blockcompressed,PetscInt *n,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
79: {
81:
83: if (!ia) return(0);
84: if ((symmetric && !A->structurally_symmetric) || oshift == 1) {
85: PetscFree(*ia);
86: if (ja) {PetscFree(*ja);}
87: }
88: return(0);
89: }
93: PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscTruth blockcompressed,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
94: {
95: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
97: PetscInt i,*collengths,*cia,*cja,n = A->cmap.n,m = A->rmap.n;
98: PetscInt nz = a->i[m],row,*jj,mr,col;
99:
101: *nn = n;
102: if (!ia) return(0);
103: if (symmetric) {
104: MatToSymmetricIJ_SeqAIJ(A->rmap.n,a->i,a->j,0,oshift,ia,ja);
105: } else {
106: PetscMalloc((n+1)*sizeof(PetscInt),&collengths);
107: PetscMemzero(collengths,n*sizeof(PetscInt));
108: PetscMalloc((n+1)*sizeof(PetscInt),&cia);
109: PetscMalloc((nz+1)*sizeof(PetscInt),&cja);
110: jj = a->j;
111: for (i=0; i<nz; i++) {
112: collengths[jj[i]]++;
113: }
114: cia[0] = oshift;
115: for (i=0; i<n; i++) {
116: cia[i+1] = cia[i] + collengths[i];
117: }
118: PetscMemzero(collengths,n*sizeof(PetscInt));
119: jj = a->j;
120: for (row=0; row<m; row++) {
121: mr = a->i[row+1] - a->i[row];
122: for (i=0; i<mr; i++) {
123: col = *jj++;
124: cja[cia[col] + collengths[col]++ - oshift] = row + oshift;
125: }
126: }
127: PetscFree(collengths);
128: *ia = cia; *ja = cja;
129: }
130: return(0);
131: }
135: PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscTruth blockcompressed,PetscInt *n,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
136: {
140: if (!ia) return(0);
142: PetscFree(*ia);
143: PetscFree(*ja);
144:
145: return(0);
146: }
150: PetscErrorCode MatSetValuesRow_SeqAIJ(Mat A,PetscInt row,const PetscScalar v[])
151: {
152: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
153: PetscInt *ai = a->i;
157: PetscMemcpy(a->a+ai[row],v,(ai[row+1]-ai[row])*sizeof(PetscScalar));
158: return(0);
159: }
161: #define CHUNKSIZE 15
165: PetscErrorCode MatSetValues_SeqAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
166: {
167: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
168: PetscInt *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N;
169: PetscInt *imax = a->imax,*ai = a->i,*ailen = a->ilen;
171: PetscInt *aj = a->j,nonew = a->nonew,lastcol = -1;
172: PetscScalar *ap,value,*aa = a->a;
173: PetscTruth ignorezeroentries = a->ignorezeroentries;
174: PetscTruth roworiented = a->roworiented;
177: for (k=0; k<m; k++) { /* loop over added rows */
178: row = im[k];
179: if (row < 0) continue;
180: #if defined(PETSC_USE_DEBUG)
181: if (row >= A->rmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap.n-1);
182: #endif
183: rp = aj + ai[row]; ap = aa + ai[row];
184: rmax = imax[row]; nrow = ailen[row];
185: low = 0;
186: high = nrow;
187: for (l=0; l<n; l++) { /* loop over added columns */
188: if (in[l] < 0) continue;
189: #if defined(PETSC_USE_DEBUG)
190: if (in[l] >= A->cmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->cmap.n-1);
191: #endif
192: col = in[l];
193: if (roworiented) {
194: value = v[l + k*n];
195: } else {
196: value = v[k + l*m];
197: }
198: if (value == 0.0 && ignorezeroentries && (is == ADD_VALUES)) continue;
200: if (col <= lastcol) low = 0; else high = nrow;
201: lastcol = col;
202: while (high-low > 5) {
203: t = (low+high)/2;
204: if (rp[t] > col) high = t;
205: else low = t;
206: }
207: for (i=low; i<high; i++) {
208: if (rp[i] > col) break;
209: if (rp[i] == col) {
210: if (is == ADD_VALUES) ap[i] += value;
211: else ap[i] = value;
212: goto noinsert;
213: }
214: }
215: if (value == 0.0 && ignorezeroentries) goto noinsert;
216: if (nonew == 1) goto noinsert;
217: if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero at (%D,%D) in the matrix",row,col);
218: MatSeqXAIJReallocateAIJ(A,A->rmap.n,1,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar);
219: N = nrow++ - 1; a->nz++; high++;
220: /* shift up all the later entries in this row */
221: for (ii=N; ii>=i; ii--) {
222: rp[ii+1] = rp[ii];
223: ap[ii+1] = ap[ii];
224: }
225: rp[i] = col;
226: ap[i] = value;
227: noinsert:;
228: low = i + 1;
229: }
230: ailen[row] = nrow;
231: }
232: A->same_nonzero = PETSC_FALSE;
233: return(0);
234: }
239: PetscErrorCode MatGetValues_SeqAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],PetscScalar v[])
240: {
241: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
242: PetscInt *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
243: PetscInt *ai = a->i,*ailen = a->ilen;
244: PetscScalar *ap,*aa = a->a;
247: for (k=0; k<m; k++) { /* loop over rows */
248: row = im[k];
249: if (row < 0) {v += n; continue;} /* SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",row); */
250: if (row >= A->rmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap.n-1);
251: rp = aj + ai[row]; ap = aa + ai[row];
252: nrow = ailen[row];
253: for (l=0; l<n; l++) { /* loop over columns */
254: if (in[l] < 0) {v++; continue;} /* SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %D",in[l]); */
255: if (in[l] >= A->cmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->cmap.n-1);
256: col = in[l] ;
257: high = nrow; low = 0; /* assume unsorted */
258: while (high-low > 5) {
259: t = (low+high)/2;
260: if (rp[t] > col) high = t;
261: else low = t;
262: }
263: for (i=low; i<high; i++) {
264: if (rp[i] > col) break;
265: if (rp[i] == col) {
266: *v++ = ap[i];
267: goto finished;
268: }
269: }
270: *v++ = 0.0;
271: finished:;
272: }
273: }
274: return(0);
275: }
280: PetscErrorCode MatView_SeqAIJ_Binary(Mat A,PetscViewer viewer)
281: {
282: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
284: PetscInt i,*col_lens;
285: int fd;
288: PetscViewerBinaryGetDescriptor(viewer,&fd);
289: PetscMalloc((4+A->rmap.n)*sizeof(PetscInt),&col_lens);
290: col_lens[0] = MAT_FILE_COOKIE;
291: col_lens[1] = A->rmap.n;
292: col_lens[2] = A->cmap.n;
293: col_lens[3] = a->nz;
295: /* store lengths of each row and write (including header) to file */
296: for (i=0; i<A->rmap.n; i++) {
297: col_lens[4+i] = a->i[i+1] - a->i[i];
298: }
299: PetscBinaryWrite(fd,col_lens,4+A->rmap.n,PETSC_INT,PETSC_TRUE);
300: PetscFree(col_lens);
302: /* store column indices (zero start index) */
303: PetscBinaryWrite(fd,a->j,a->nz,PETSC_INT,PETSC_FALSE);
305: /* store nonzero values */
306: PetscBinaryWrite(fd,a->a,a->nz,PETSC_SCALAR,PETSC_FALSE);
307: return(0);
308: }
310: EXTERN PetscErrorCode MatSeqAIJFactorInfo_Matlab(Mat,PetscViewer);
314: PetscErrorCode MatView_SeqAIJ_ASCII(Mat A,PetscViewer viewer)
315: {
316: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
317: PetscErrorCode ierr;
318: PetscInt i,j,m = A->rmap.n,shift=0;
319: const char *name;
320: PetscViewerFormat format;
323: PetscObjectGetName((PetscObject)A,&name);
324: PetscViewerGetFormat(viewer,&format);
325: if (format == PETSC_VIEWER_ASCII_MATLAB) {
326: PetscInt nofinalvalue = 0;
327: if ((a->i[m] == a->i[m-1]) || (a->j[a->nz-1] != A->cmap.n-!shift)) {
328: nofinalvalue = 1;
329: }
330: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
331: PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",m,A->cmap.n);
332: PetscViewerASCIIPrintf(viewer,"%% Nonzeros = %D \n",a->nz);
333: PetscViewerASCIIPrintf(viewer,"zzz = zeros(%D,3);\n",a->nz+nofinalvalue);
334: PetscViewerASCIIPrintf(viewer,"zzz = [\n");
336: for (i=0; i<m; i++) {
337: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
338: #if defined(PETSC_USE_COMPLEX)
339: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e + %18.16ei \n",i+1,a->j[j]+!shift,PetscRealPart(a->a[j]),PetscImaginaryPart(a->a[j]));
340: #else
341: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",i+1,a->j[j]+!shift,a->a[j]);
342: #endif
343: }
344: }
345: if (nofinalvalue) {
346: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",m,A->cmap.n,0.0);
347: }
348: PetscViewerASCIIPrintf(viewer,"];\n %s = spconvert(zzz);\n",name);
349: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
350: } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) {
351: return(0);
352: } else if (format == PETSC_VIEWER_ASCII_COMMON) {
353: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
354: for (i=0; i<m; i++) {
355: PetscViewerASCIIPrintf(viewer,"row %D:",i);
356: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
357: #if defined(PETSC_USE_COMPLEX)
358: if (PetscImaginaryPart(a->a[j]) > 0.0 && PetscRealPart(a->a[j]) != 0.0) {
359: PetscViewerASCIIPrintf(viewer," (%D, %G + %G i)",a->j[j]+shift,PetscRealPart(a->a[j]),PetscImaginaryPart(a->a[j]));
360: } else if (PetscImaginaryPart(a->a[j]) < 0.0 && PetscRealPart(a->a[j]) != 0.0) {
361: PetscViewerASCIIPrintf(viewer," (%D, %G - %G i)",a->j[j]+shift,PetscRealPart(a->a[j]),-PetscImaginaryPart(a->a[j]));
362: } else if (PetscRealPart(a->a[j]) != 0.0) {
363: PetscViewerASCIIPrintf(viewer," (%D, %G) ",a->j[j]+shift,PetscRealPart(a->a[j]));
364: }
365: #else
366: if (a->a[j] != 0.0) {PetscViewerASCIIPrintf(viewer," (%D, %G) ",a->j[j]+shift,a->a[j]);}
367: #endif
368: }
369: PetscViewerASCIIPrintf(viewer,"\n");
370: }
371: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
372: } else if (format == PETSC_VIEWER_ASCII_SYMMODU) {
373: PetscInt nzd=0,fshift=1,*sptr;
374: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
375: PetscMalloc((m+1)*sizeof(PetscInt),&sptr);
376: for (i=0; i<m; i++) {
377: sptr[i] = nzd+1;
378: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
379: if (a->j[j] >= i) {
380: #if defined(PETSC_USE_COMPLEX)
381: if (PetscImaginaryPart(a->a[j]) != 0.0 || PetscRealPart(a->a[j]) != 0.0) nzd++;
382: #else
383: if (a->a[j] != 0.0) nzd++;
384: #endif
385: }
386: }
387: }
388: sptr[m] = nzd+1;
389: PetscViewerASCIIPrintf(viewer," %D %D\n\n",m,nzd);
390: for (i=0; i<m+1; i+=6) {
391: 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]);}
392: 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]);}
393: else if (i+2<m) {PetscViewerASCIIPrintf(viewer," %D %D %D %D\n",sptr[i],sptr[i+1],sptr[i+2],sptr[i+3]);}
394: else if (i+1<m) {PetscViewerASCIIPrintf(viewer," %D %D %D\n",sptr[i],sptr[i+1],sptr[i+2]);}
395: else if (i<m) {PetscViewerASCIIPrintf(viewer," %D %D\n",sptr[i],sptr[i+1]);}
396: else {PetscViewerASCIIPrintf(viewer," %D\n",sptr[i]);}
397: }
398: PetscViewerASCIIPrintf(viewer,"\n");
399: PetscFree(sptr);
400: for (i=0; i<m; i++) {
401: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
402: if (a->j[j] >= i) {PetscViewerASCIIPrintf(viewer," %D ",a->j[j]+fshift);}
403: }
404: PetscViewerASCIIPrintf(viewer,"\n");
405: }
406: PetscViewerASCIIPrintf(viewer,"\n");
407: for (i=0; i<m; i++) {
408: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
409: if (a->j[j] >= i) {
410: #if defined(PETSC_USE_COMPLEX)
411: if (PetscImaginaryPart(a->a[j]) != 0.0 || PetscRealPart(a->a[j]) != 0.0) {
412: PetscViewerASCIIPrintf(viewer," %18.16e %18.16e ",PetscRealPart(a->a[j]),PetscImaginaryPart(a->a[j]));
413: }
414: #else
415: if (a->a[j] != 0.0) {PetscViewerASCIIPrintf(viewer," %18.16e ",a->a[j]);}
416: #endif
417: }
418: }
419: PetscViewerASCIIPrintf(viewer,"\n");
420: }
421: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
422: } else if (format == PETSC_VIEWER_ASCII_DENSE) {
423: PetscInt cnt = 0,jcnt;
424: PetscScalar value;
426: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
427: for (i=0; i<m; i++) {
428: jcnt = 0;
429: for (j=0; j<A->cmap.n; j++) {
430: if (jcnt < a->i[i+1]-a->i[i] && j == a->j[cnt]) {
431: value = a->a[cnt++];
432: jcnt++;
433: } else {
434: value = 0.0;
435: }
436: #if defined(PETSC_USE_COMPLEX)
437: PetscViewerASCIIPrintf(viewer," %7.5e+%7.5e i ",PetscRealPart(value),PetscImaginaryPart(value));
438: #else
439: PetscViewerASCIIPrintf(viewer," %7.5e ",value);
440: #endif
441: }
442: PetscViewerASCIIPrintf(viewer,"\n");
443: }
444: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
445: } else {
446: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
447: for (i=0; i<m; i++) {
448: PetscViewerASCIIPrintf(viewer,"row %D:",i);
449: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
450: #if defined(PETSC_USE_COMPLEX)
451: if (PetscImaginaryPart(a->a[j]) > 0.0) {
452: PetscViewerASCIIPrintf(viewer," (%D, %G + %G i)",a->j[j]+shift,PetscRealPart(a->a[j]),PetscImaginaryPart(a->a[j]));
453: } else if (PetscImaginaryPart(a->a[j]) < 0.0) {
454: PetscViewerASCIIPrintf(viewer," (%D, %G - %G i)",a->j[j]+shift,PetscRealPart(a->a[j]),-PetscImaginaryPart(a->a[j]));
455: } else {
456: PetscViewerASCIIPrintf(viewer," (%D, %G) ",a->j[j]+shift,PetscRealPart(a->a[j]));
457: }
458: #else
459: PetscViewerASCIIPrintf(viewer," (%D, %G) ",a->j[j]+shift,a->a[j]);
460: #endif
461: }
462: PetscViewerASCIIPrintf(viewer,"\n");
463: }
464: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
465: }
466: PetscViewerFlush(viewer);
467: return(0);
468: }
472: PetscErrorCode MatView_SeqAIJ_Draw_Zoom(PetscDraw draw,void *Aa)
473: {
474: Mat A = (Mat) Aa;
475: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
476: PetscErrorCode ierr;
477: PetscInt i,j,m = A->rmap.n,color;
478: PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r,maxv = 0.0;
479: PetscViewer viewer;
480: PetscViewerFormat format;
483: PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);
484: PetscViewerGetFormat(viewer,&format);
486: PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);
487: /* loop over matrix elements drawing boxes */
489: if (format != PETSC_VIEWER_DRAW_CONTOUR) {
490: /* Blue for negative, Cyan for zero and Red for positive */
491: color = PETSC_DRAW_BLUE;
492: for (i=0; i<m; i++) {
493: y_l = m - i - 1.0; y_r = y_l + 1.0;
494: for (j=a->i[i]; j<a->i[i+1]; j++) {
495: x_l = a->j[j] ; x_r = x_l + 1.0;
496: #if defined(PETSC_USE_COMPLEX)
497: if (PetscRealPart(a->a[j]) >= 0.) continue;
498: #else
499: if (a->a[j] >= 0.) continue;
500: #endif
501: PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
502: }
503: }
504: color = PETSC_DRAW_CYAN;
505: for (i=0; i<m; i++) {
506: y_l = m - i - 1.0; y_r = y_l + 1.0;
507: for (j=a->i[i]; j<a->i[i+1]; j++) {
508: x_l = a->j[j]; x_r = x_l + 1.0;
509: if (a->a[j] != 0.) continue;
510: PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
511: }
512: }
513: color = PETSC_DRAW_RED;
514: for (i=0; i<m; i++) {
515: y_l = m - i - 1.0; y_r = y_l + 1.0;
516: for (j=a->i[i]; j<a->i[i+1]; j++) {
517: x_l = a->j[j]; x_r = x_l + 1.0;
518: #if defined(PETSC_USE_COMPLEX)
519: if (PetscRealPart(a->a[j]) <= 0.) continue;
520: #else
521: if (a->a[j] <= 0.) continue;
522: #endif
523: PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
524: }
525: }
526: } else {
527: /* use contour shading to indicate magnitude of values */
528: /* first determine max of all nonzero values */
529: PetscInt nz = a->nz,count;
530: PetscDraw popup;
531: PetscReal scale;
533: for (i=0; i<nz; i++) {
534: if (PetscAbsScalar(a->a[i]) > maxv) maxv = PetscAbsScalar(a->a[i]);
535: }
536: scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv;
537: PetscDrawGetPopup(draw,&popup);
538: if (popup) {PetscDrawScalePopup(popup,0.0,maxv);}
539: count = 0;
540: for (i=0; i<m; i++) {
541: y_l = m - i - 1.0; y_r = y_l + 1.0;
542: for (j=a->i[i]; j<a->i[i+1]; j++) {
543: x_l = a->j[j]; x_r = x_l + 1.0;
544: color = PETSC_DRAW_BASIC_COLORS + (PetscInt)(scale*PetscAbsScalar(a->a[count]));
545: PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
546: count++;
547: }
548: }
549: }
550: return(0);
551: }
555: PetscErrorCode MatView_SeqAIJ_Draw(Mat A,PetscViewer viewer)
556: {
558: PetscDraw draw;
559: PetscReal xr,yr,xl,yl,h,w;
560: PetscTruth isnull;
563: PetscViewerDrawGetDraw(viewer,0,&draw);
564: PetscDrawIsNull(draw,&isnull);
565: if (isnull) return(0);
567: PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);
568: xr = A->cmap.n; yr = A->rmap.n; h = yr/10.0; w = xr/10.0;
569: xr += w; yr += h; xl = -w; yl = -h;
570: PetscDrawSetCoordinates(draw,xl,yl,xr,yr);
571: PetscDrawZoom(draw,MatView_SeqAIJ_Draw_Zoom,A);
572: PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);
573: return(0);
574: }
578: PetscErrorCode MatView_SeqAIJ(Mat A,PetscViewer viewer)
579: {
581: PetscTruth iascii,isbinary,isdraw;
584: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
585: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
586: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
587: if (iascii) {
588: MatView_SeqAIJ_ASCII(A,viewer);
589: } else if (isbinary) {
590: MatView_SeqAIJ_Binary(A,viewer);
591: } else if (isdraw) {
592: MatView_SeqAIJ_Draw(A,viewer);
593: } else {
594: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by SeqAIJ matrices",((PetscObject)viewer)->type_name);
595: }
596: MatView_Inode(A,viewer);
597: return(0);
598: }
602: PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode)
603: {
604: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
606: PetscInt fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
607: PetscInt m = A->rmap.n,*ip,N,*ailen = a->ilen,rmax = 0;
608: PetscScalar *aa = a->a,*ap;
609: PetscReal ratio=0.6;
612: if (mode == MAT_FLUSH_ASSEMBLY) return(0);
614: if (m) rmax = ailen[0]; /* determine row with most nonzeros */
615: for (i=1; i<m; i++) {
616: /* move each row back by the amount of empty slots (fshift) before it*/
617: fshift += imax[i-1] - ailen[i-1];
618: rmax = PetscMax(rmax,ailen[i]);
619: if (fshift) {
620: ip = aj + ai[i] ;
621: ap = aa + ai[i] ;
622: N = ailen[i];
623: for (j=0; j<N; j++) {
624: ip[j-fshift] = ip[j];
625: ap[j-fshift] = ap[j];
626: }
627: }
628: ai[i] = ai[i-1] + ailen[i-1];
629: }
630: if (m) {
631: fshift += imax[m-1] - ailen[m-1];
632: ai[m] = ai[m-1] + ailen[m-1];
633: }
634: /* reset ilen and imax for each row */
635: for (i=0; i<m; i++) {
636: ailen[i] = imax[i] = ai[i+1] - ai[i];
637: }
638: a->nz = ai[m];
640: MatMarkDiagonal_SeqAIJ(A);
641: PetscInfo4(A,"Matrix size: %D X %D; storage space: %D unneeded,%D used\n",m,A->cmap.n,fshift,a->nz);
642: PetscInfo1(A,"Number of mallocs during MatSetValues() is %D\n",a->reallocs);
643: PetscInfo1(A,"Maximum nonzeros in any row is %D\n",rmax);
645: a->reallocs = 0;
646: A->info.nz_unneeded = (double)fshift;
647: a->rmax = rmax;
649: /* check for zero rows. If found a large number of zero rows, use CompressedRow functions */
650: Mat_CheckCompressedRow(A,&a->compressedrow,a->i,m,ratio);
651: A->same_nonzero = PETSC_TRUE;
653: MatAssemblyEnd_Inode(A,mode);
655: a->idiagvalid = PETSC_FALSE;
656: return(0);
657: }
661: PetscErrorCode MatRealPart_SeqAIJ(Mat A)
662: {
663: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
664: PetscInt i,nz = a->nz;
665: PetscScalar *aa = a->a;
668: for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]);
669: return(0);
670: }
674: PetscErrorCode MatImaginaryPart_SeqAIJ(Mat A)
675: {
676: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
677: PetscInt i,nz = a->nz;
678: PetscScalar *aa = a->a;
681: for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]);
682: return(0);
683: }
687: PetscErrorCode MatZeroEntries_SeqAIJ(Mat A)
688: {
689: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
693: PetscMemzero(a->a,(a->i[A->rmap.n])*sizeof(PetscScalar));
694: return(0);
695: }
699: PetscErrorCode MatDestroy_SeqAIJ(Mat A)
700: {
701: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
705: #if defined(PETSC_USE_LOG)
706: PetscLogObjectState((PetscObject)A,"Rows=%D, Cols=%D, NZ=%D",A->rmap.n,A->cmap.n,a->nz);
707: #endif
708: MatSeqXAIJFreeAIJ(A,&a->a,&a->j,&a->i);
709: if (a->row) {
710: ISDestroy(a->row);
711: }
712: if (a->col) {
713: ISDestroy(a->col);
714: }
715: PetscFree(a->diag);
716: PetscFree2(a->imax,a->ilen);
717: PetscFree3(a->idiag,a->mdiag,a->ssor_work);
718: PetscFree(a->solve_work);
719: if (a->icol) {ISDestroy(a->icol);}
720: PetscFree(a->saved_values);
721: if (a->coloring) {ISColoringDestroy(a->coloring);}
722: PetscFree(a->xtoy);
723: if (a->XtoY) {MatDestroy(a->XtoY);}
724: if (a->compressedrow.use){PetscFree(a->compressedrow.i);}
726: MatDestroy_Inode(A);
728: PetscFree(a);
730: PetscObjectChangeTypeName((PetscObject)A,0);
731: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSeqAIJSetColumnIndices_C","",PETSC_NULL);
732: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatStoreValues_C","",PETSC_NULL);
733: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatRetrieveValues_C","",PETSC_NULL);
734: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatConvert_seqaij_seqsbaij_C","",PETSC_NULL);
735: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatConvert_seqaij_seqbaij_C","",PETSC_NULL);
736: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatConvert_seqaij_seqcsrperm_C","",PETSC_NULL);
737: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatIsTranspose_C","",PETSC_NULL);
738: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSeqAIJSetPreallocation_C","",PETSC_NULL);
739: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSeqAIJSetPreallocationCSR_C","",PETSC_NULL);
740: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatReorderForNonzeroDiagonal_C","",PETSC_NULL);
741: return(0);
742: }
746: PetscErrorCode MatCompress_SeqAIJ(Mat A)
747: {
749: return(0);
750: }
754: PetscErrorCode MatSetOption_SeqAIJ(Mat A,MatOption op,PetscTruth flg)
755: {
756: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
760: switch (op) {
761: case MAT_ROW_ORIENTED:
762: a->roworiented = flg;
763: break;
764: case MAT_KEEP_ZEROED_ROWS:
765: a->keepzeroedrows = flg;
766: break;
767: case MAT_NEW_NONZERO_LOCATIONS:
768: a->nonew = (flg ? 0 : 1);
769: break;
770: case MAT_NEW_NONZERO_LOCATION_ERR:
771: a->nonew = (flg ? -1 : 0);
772: break;
773: case MAT_NEW_NONZERO_ALLOCATION_ERR:
774: a->nonew = (flg ? -2 : 0);
775: break;
776: case MAT_IGNORE_ZERO_ENTRIES:
777: a->ignorezeroentries = flg;
778: break;
779: case MAT_USE_COMPRESSEDROW:
780: a->compressedrow.use = flg;
781: break;
782: case MAT_NEW_DIAGONALS:
783: case MAT_IGNORE_OFF_PROC_ENTRIES:
784: case MAT_USE_HASH_TABLE:
785: PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
786: break;
787: default:
788: break;
789: }
790: MatSetOption_Inode(A,op,flg);
791: return(0);
792: }
796: PetscErrorCode MatGetDiagonal_SeqAIJ(Mat A,Vec v)
797: {
798: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
800: PetscInt i,j,n;
801: PetscScalar *x,zero = 0.0;
804: VecSet(v,zero);
805: VecGetArray(v,&x);
806: VecGetLocalSize(v,&n);
807: if (n != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
808: for (i=0; i<A->rmap.n; i++) {
809: for (j=a->i[i]; j<a->i[i+1]; j++) {
810: if (a->j[j] == i) {
811: x[i] = a->a[j];
812: break;
813: }
814: }
815: }
816: VecRestoreArray(v,&x);
817: return(0);
818: }
822: PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
823: {
824: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
825: PetscScalar *x,*y;
826: PetscErrorCode ierr;
827: PetscInt m = A->rmap.n;
828: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
829: PetscScalar *v,alpha;
830: PetscInt n,i,*idx,*ii,*ridx=PETSC_NULL;
831: Mat_CompressedRow cprow = a->compressedrow;
832: PetscTruth usecprow = cprow.use;
833: #endif
836: if (zz != yy) {VecCopy(zz,yy);}
837: VecGetArray(xx,&x);
838: VecGetArray(yy,&y);
840: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
841: fortranmulttransposeaddaij_(&m,x,a->i,a->j,a->a,y);
842: #else
843: if (usecprow){
844: m = cprow.nrows;
845: ii = cprow.i;
846: ridx = cprow.rindex;
847: } else {
848: ii = a->i;
849: }
850: for (i=0; i<m; i++) {
851: idx = a->j + ii[i] ;
852: v = a->a + ii[i] ;
853: n = ii[i+1] - ii[i];
854: if (usecprow){
855: alpha = x[ridx[i]];
856: } else {
857: alpha = x[i];
858: }
859: while (n-->0) {y[*idx++] += alpha * *v++;}
860: }
861: #endif
862: PetscLogFlops(2*a->nz);
863: VecRestoreArray(xx,&x);
864: VecRestoreArray(yy,&y);
865: return(0);
866: }
870: PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec xx,Vec yy)
871: {
872: PetscScalar zero = 0.0;
876: VecSet(yy,zero);
877: MatMultTransposeAdd_SeqAIJ(A,xx,yy,yy);
878: return(0);
879: }
884: PetscErrorCode MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
885: {
886: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
887: PetscScalar *y;
888: const PetscScalar *x,*aa;
889: PetscErrorCode ierr;
890: PetscInt m=A->rmap.n,*aj,*ii;
891: PetscInt n,i,j,nonzerorow=0,*ridx=PETSC_NULL;
892: PetscScalar sum;
893: PetscTruth usecprow=a->compressedrow.use;
894: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
895: PetscInt jrow;
896: #endif
898: #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
899: #pragma disjoint(*x,*y,*aa)
900: #endif
903: VecGetArray(xx,(PetscScalar**)&x);
904: VecGetArray(yy,&y);
905: aj = a->j;
906: aa = a->a;
907: ii = a->i;
908: if (usecprow){ /* use compressed row format */
909: m = a->compressedrow.nrows;
910: ii = a->compressedrow.i;
911: ridx = a->compressedrow.rindex;
912: for (i=0; i<m; i++){
913: n = ii[i+1] - ii[i];
914: aj = a->j + ii[i];
915: aa = a->a + ii[i];
916: sum = 0.0;
917: nonzerorow += (n>0);
918: for (j=0; j<n; j++) sum += (*aa++)*x[*aj++];
919: y[*ridx++] = sum;
920: }
921: } else { /* do not use compressed row format */
922: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
923: fortranmultaij_(&m,x,ii,aj,aa,y);
924: #else
925: for (i=0; i<m; i++) {
926: jrow = ii[i];
927: n = ii[i+1] - jrow;
928: sum = 0.0;
929: nonzerorow += (n>0);
930: for (j=0; j<n; j++) {
931: sum += aa[jrow]*x[aj[jrow]]; jrow++;
932: }
933: y[i] = sum;
934: }
935: #endif
936: }
937: PetscLogFlops(2*a->nz - nonzerorow);
938: VecRestoreArray(xx,(PetscScalar**)&x);
939: VecRestoreArray(yy,&y);
940: return(0);
941: }
945: PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
946: {
947: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
948: PetscScalar *x,*y,*z,*aa;
950: PetscInt m = A->rmap.n,*aj,*ii;
951: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
952: PetscInt n,i,jrow,j,*ridx=PETSC_NULL;
953: PetscScalar sum;
954: PetscTruth usecprow=a->compressedrow.use;
955: #endif
958: VecGetArray(xx,&x);
959: VecGetArray(yy,&y);
960: if (zz != yy) {
961: VecGetArray(zz,&z);
962: } else {
963: z = y;
964: }
966: aj = a->j;
967: aa = a->a;
968: ii = a->i;
969: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
970: fortranmultaddaij_(&m,x,ii,aj,aa,y,z);
971: #else
972: if (usecprow){ /* use compressed row format */
973: if (zz != yy){
974: PetscMemcpy(z,y,m*sizeof(PetscScalar));
975: }
976: m = a->compressedrow.nrows;
977: ii = a->compressedrow.i;
978: ridx = a->compressedrow.rindex;
979: for (i=0; i<m; i++){
980: n = ii[i+1] - ii[i];
981: aj = a->j + ii[i];
982: aa = a->a + ii[i];
983: sum = y[*ridx];
984: for (j=0; j<n; j++) sum += (*aa++)*x[*aj++];
985: z[*ridx++] = sum;
986: }
987: } else { /* do not use compressed row format */
988: for (i=0; i<m; i++) {
989: jrow = ii[i];
990: n = ii[i+1] - jrow;
991: sum = y[i];
992: for (j=0; j<n; j++) {
993: sum += aa[jrow]*x[aj[jrow]]; jrow++;
994: }
995: z[i] = sum;
996: }
997: }
998: #endif
999: PetscLogFlops(2*a->nz);
1000: VecRestoreArray(xx,&x);
1001: VecRestoreArray(yy,&y);
1002: if (zz != yy) {
1003: VecRestoreArray(zz,&z);
1004: }
1005: return(0);
1006: }
1008: /*
1009: Adds diagonal pointers to sparse matrix structure.
1010: */
1013: PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat A)
1014: {
1015: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1017: PetscInt i,j,m = A->rmap.n;
1020: if (!a->diag) {
1021: PetscMalloc(m*sizeof(PetscInt),&a->diag);
1022: PetscLogObjectMemory(A, m*sizeof(PetscInt));
1023: }
1024: for (i=0; i<A->rmap.n; i++) {
1025: a->diag[i] = a->i[i+1];
1026: for (j=a->i[i]; j<a->i[i+1]; j++) {
1027: if (a->j[j] == i) {
1028: a->diag[i] = j;
1029: break;
1030: }
1031: }
1032: }
1033: return(0);
1034: }
1036: /*
1037: Checks for missing diagonals
1038: */
1041: PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat A,PetscTruth *missing,PetscInt *d)
1042: {
1043: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1044: PetscInt *diag,*jj = a->j,i;
1047: *missing = PETSC_FALSE;
1048: if (A->rmap.n > 0 && !jj) {
1049: *missing = PETSC_TRUE;
1050: if (d) *d = 0;
1051: PetscInfo(A,"Matrix has no entries therefor is missing diagonal");
1052: } else {
1053: diag = a->diag;
1054: for (i=0; i<A->rmap.n; i++) {
1055: if (jj[diag[i]] != i) {
1056: *missing = PETSC_TRUE;
1057: if (d) *d = i;
1058: PetscInfo1(A,"Matrix is missing diagonal number %D",i);
1059: }
1060: }
1061: }
1062: return(0);
1063: }
1068: PetscErrorCode MatInvertDiagonal_SeqAIJ(Mat A,PetscScalar omega,PetscScalar fshift)
1069: {
1070: Mat_SeqAIJ *a = (Mat_SeqAIJ*) A->data;
1072: PetscInt i,*diag,m = A->rmap.n;
1073: PetscScalar *v = a->a,*idiag,*mdiag;
1076: if (a->idiagvalid) return(0);
1077: MatMarkDiagonal_SeqAIJ(A);
1078: diag = a->diag;
1079: if (!a->idiag) {
1080: PetscMalloc3(m,PetscScalar,&a->idiag,m,PetscScalar,&a->mdiag,m,PetscScalar,&a->ssor_work);
1081: PetscLogObjectMemory(A, 3*m*sizeof(PetscScalar));
1082: v = a->a;
1083: }
1084: mdiag = a->mdiag;
1085: idiag = a->idiag;
1086:
1087: if (omega == 1.0 && !PetscAbsScalar(fshift)) {
1088: for (i=0; i<m; i++) {
1089: mdiag[i] = v[diag[i]];
1090: if (!PetscAbsScalar(mdiag[i])) SETERRQ1(PETSC_ERR_ARG_INCOMP,"Zero diagonal on row %D",i);
1091: idiag[i] = 1.0/v[diag[i]];
1092: }
1093: PetscLogFlops(m);
1094: } else {
1095: for (i=0; i<m; i++) {
1096: mdiag[i] = v[diag[i]];
1097: idiag[i] = omega/(fshift + v[diag[i]]);
1098: }
1099: PetscLogFlops(2*m);
1100: }
1101: a->idiagvalid = PETSC_TRUE;
1102: return(0);
1103: }
1108: PetscErrorCode MatRelax_SeqAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
1109: {
1110: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1111: PetscScalar *x,d,*xs,sum,*t,scale,*idiag=0,*mdiag;
1112: const PetscScalar *v = a->a, *b, *bs,*xb, *ts;
1113: PetscErrorCode ierr;
1114: PetscInt n = A->cmap.n,m = A->rmap.n,i;
1115: const PetscInt *idx,*diag;
1118: its = its*lits;
1120: if (fshift != a->fshift || omega != a->omega) a->idiagvalid = PETSC_FALSE; /* must recompute idiag[] */
1121: if (!a->idiagvalid) {MatInvertDiagonal_SeqAIJ(A,omega,fshift);}
1122: a->fshift = fshift;
1123: a->omega = omega;
1125: diag = a->diag;
1126: t = a->ssor_work;
1127: idiag = a->idiag;
1128: mdiag = a->mdiag;
1130: VecGetArray(xx,&x);
1131: if (xx != bb) {
1132: VecGetArray(bb,(PetscScalar**)&b);
1133: } else {
1134: b = x;
1135: }
1136: CHKMEMQ;
1137: /* We count flops by assuming the upper triangular and lower triangular parts have the same number of nonzeros */
1138: xs = x;
1139: if (flag == SOR_APPLY_UPPER) {
1140: /* apply (U + D/omega) to the vector */
1141: bs = b;
1142: for (i=0; i<m; i++) {
1143: d = fshift + mdiag[i];
1144: n = a->i[i+1] - diag[i] - 1;
1145: idx = a->j + diag[i] + 1;
1146: v = a->a + diag[i] + 1;
1147: sum = b[i]*d/omega;
1148: SPARSEDENSEDOT(sum,bs,v,idx,n);
1149: x[i] = sum;
1150: }
1151: VecRestoreArray(xx,&x);
1152: if (bb != xx) {VecRestoreArray(bb,(PetscScalar**)&b);}
1153: PetscLogFlops(a->nz);
1154: return(0);
1155: }
1157: if (flag == SOR_APPLY_LOWER) {
1158: SETERRQ(PETSC_ERR_SUP,"SOR_APPLY_LOWER is not implemented");
1159: } else if (flag & SOR_EISENSTAT) {
1160: /* Let A = L + U + D; where L is lower trianglar,
1161: U is upper triangular, E is diagonal; This routine applies
1163: (L + E)^{-1} A (U + E)^{-1}
1165: to a vector efficiently using Eisenstat's trick. This is for
1166: the case of SSOR preconditioner, so E is D/omega where omega
1167: is the relaxation factor.
1168: */
1169: scale = (2.0/omega) - 1.0;
1171: /* x = (E + U)^{-1} b */
1172: for (i=m-1; i>=0; i--) {
1173: n = a->i[i+1] - diag[i] - 1;
1174: idx = a->j + diag[i] + 1;
1175: v = a->a + diag[i] + 1;
1176: sum = b[i];
1177: SPARSEDENSEMDOT(sum,xs,v,idx,n);
1178: x[i] = sum*idiag[i];
1179: }
1181: /* t = b - (2*E - D)x */
1182: v = a->a;
1183: for (i=0; i<m; i++) { t[i] = b[i] - scale*(v[*diag++])*x[i]; }
1185: /* t = (E + L)^{-1}t */
1186: ts = t;
1187: diag = a->diag;
1188: for (i=0; i<m; i++) {
1189: n = diag[i] - a->i[i];
1190: idx = a->j + a->i[i];
1191: v = a->a + a->i[i];
1192: sum = t[i];
1193: SPARSEDENSEMDOT(sum,ts,v,idx,n);
1194: t[i] = sum*idiag[i];
1195: /* x = x + t */
1196: x[i] += t[i];
1197: }
1199: PetscLogFlops(6*m-1 + 2*a->nz);
1200: VecRestoreArray(xx,&x);
1201: if (bb != xx) {VecRestoreArray(bb,(PetscScalar**)&b);}
1202: return(0);
1203: }
1204: if (flag & SOR_ZERO_INITIAL_GUESS) {
1205: if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
1206: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1207: fortranrelaxaijforwardzero_(&m,&omega,x,a->i,a->j,(PetscInt*)diag,idiag,a->a,(void*)b);
1208: #else
1209: for (i=0; i<m; i++) {
1210: n = diag[i] - a->i[i];
1211: idx = a->j + a->i[i];
1212: v = a->a + a->i[i];
1213: sum = b[i];
1214: SPARSEDENSEMDOT(sum,xs,v,idx,n);
1215: x[i] = sum*idiag[i];
1216: }
1217: #endif
1218: xb = x;
1219: PetscLogFlops(a->nz);
1220: } else xb = b;
1221: if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
1222: (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
1223: for (i=0; i<m; i++) {
1224: x[i] *= mdiag[i];
1225: }
1226: PetscLogFlops(m);
1227: }
1228: if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
1229: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1230: fortranrelaxaijbackwardzero_(&m,&omega,x,a->i,a->j,(PetscInt*)diag,idiag,a->a,(void*)xb);
1231: #else
1232: for (i=m-1; i>=0; i--) {
1233: n = a->i[i+1] - diag[i] - 1;
1234: idx = a->j + diag[i] + 1;
1235: v = a->a + diag[i] + 1;
1236: sum = xb[i];
1237: SPARSEDENSEMDOT(sum,xs,v,idx,n);
1238: x[i] = sum*idiag[i];
1239: }
1240: #endif
1241: PetscLogFlops(a->nz);
1242: }
1243: its--;
1244: }
1245: while (its--) {
1246: if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
1247: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1248: fortranrelaxaijforward_(&m,&omega,x,a->i,a->j,(PetscInt*)diag,a->a,(void*)b);
1249: #else
1250: for (i=0; i<m; i++) {
1251: n = a->i[i+1] - a->i[i];
1252: idx = a->j + a->i[i];
1253: v = a->a + a->i[i];
1254: sum = b[i];
1255: SPARSEDENSEMDOT(sum,xs,v,idx,n);
1256: x[i] = (1. - omega)*x[i] + (sum + mdiag[i]*x[i])*idiag[i];
1257: }
1258: #endif
1259: PetscLogFlops(a->nz);
1260: }
1261: if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
1262: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1263: fortranrelaxaijbackward_(&m,&omega,x,a->i,a->j,(PetscInt*)diag,a->a,(void*)b);
1264: #else
1265: for (i=m-1; i>=0; i--) {
1266: n = a->i[i+1] - a->i[i];
1267: idx = a->j + a->i[i];
1268: v = a->a + a->i[i];
1269: sum = b[i];
1270: SPARSEDENSEMDOT(sum,xs,v,idx,n);
1271: x[i] = (1. - omega)*x[i] + (sum + mdiag[i]*x[i])*idiag[i];
1272: }
1273: #endif
1274: PetscLogFlops(a->nz);
1275: }
1276: }
1277: VecRestoreArray(xx,&x);
1278: if (bb != xx) {VecRestoreArray(bb,(PetscScalar**)&b);}
1279: CHKMEMQ; return(0);
1280: }
1285: PetscErrorCode MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,MatInfo *info)
1286: {
1287: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1290: info->rows_global = (double)A->rmap.n;
1291: info->columns_global = (double)A->cmap.n;
1292: info->rows_local = (double)A->rmap.n;
1293: info->columns_local = (double)A->cmap.n;
1294: info->block_size = 1.0;
1295: info->nz_allocated = (double)a->maxnz;
1296: info->nz_used = (double)a->nz;
1297: info->nz_unneeded = (double)(a->maxnz - a->nz);
1298: info->assemblies = (double)A->num_ass;
1299: info->mallocs = (double)a->reallocs;
1300: info->memory = ((PetscObject)A)->mem;
1301: if (A->factor) {
1302: info->fill_ratio_given = A->info.fill_ratio_given;
1303: info->fill_ratio_needed = A->info.fill_ratio_needed;
1304: info->factor_mallocs = A->info.factor_mallocs;
1305: } else {
1306: info->fill_ratio_given = 0;
1307: info->fill_ratio_needed = 0;
1308: info->factor_mallocs = 0;
1309: }
1310: return(0);
1311: }
1315: PetscErrorCode MatZeroRows_SeqAIJ(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag)
1316: {
1317: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1318: PetscInt i,m = A->rmap.n - 1,d;
1320: PetscTruth missing;
1323: if (a->keepzeroedrows) {
1324: for (i=0; i<N; i++) {
1325: if (rows[i] < 0 || rows[i] > m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range", rows[i]);
1326: PetscMemzero(&a->a[a->i[rows[i]]],a->ilen[rows[i]]*sizeof(PetscScalar));
1327: }
1328: if (diag != 0.0) {
1329: MatMissingDiagonal_SeqAIJ(A,&missing,&d);
1330: if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry in row %D",d);
1331: for (i=0; i<N; i++) {
1332: a->a[a->diag[rows[i]]] = diag;
1333: }
1334: }
1335: A->same_nonzero = PETSC_TRUE;
1336: } else {
1337: if (diag != 0.0) {
1338: for (i=0; i<N; i++) {
1339: if (rows[i] < 0 || rows[i] > m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range", rows[i]);
1340: if (a->ilen[rows[i]] > 0) {
1341: a->ilen[rows[i]] = 1;
1342: a->a[a->i[rows[i]]] = diag;
1343: a->j[a->i[rows[i]]] = rows[i];
1344: } else { /* in case row was completely empty */
1345: MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],&diag,INSERT_VALUES);
1346: }
1347: }
1348: } else {
1349: for (i=0; i<N; i++) {
1350: if (rows[i] < 0 || rows[i] > m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range", rows[i]);
1351: a->ilen[rows[i]] = 0;
1352: }
1353: }
1354: A->same_nonzero = PETSC_FALSE;
1355: }
1356: MatAssemblyEnd_SeqAIJ(A,MAT_FINAL_ASSEMBLY);
1357: return(0);
1358: }
1362: PetscErrorCode MatGetRow_SeqAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
1363: {
1364: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1365: PetscInt *itmp;
1368: if (row < 0 || row >= A->rmap.n) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Row %D out of range",row);
1370: *nz = a->i[row+1] - a->i[row];
1371: if (v) *v = a->a + a->i[row];
1372: if (idx) {
1373: itmp = a->j + a->i[row];
1374: if (*nz) {
1375: *idx = itmp;
1376: }
1377: else *idx = 0;
1378: }
1379: return(0);
1380: }
1382: /* remove this function? */
1385: PetscErrorCode MatRestoreRow_SeqAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
1386: {
1388: return(0);
1389: }
1393: PetscErrorCode MatNorm_SeqAIJ(Mat A,NormType type,PetscReal *nrm)
1394: {
1395: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1396: PetscScalar *v = a->a;
1397: PetscReal sum = 0.0;
1399: PetscInt i,j;
1402: if (type == NORM_FROBENIUS) {
1403: for (i=0; i<a->nz; i++) {
1404: #if defined(PETSC_USE_COMPLEX)
1405: sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1406: #else
1407: sum += (*v)*(*v); v++;
1408: #endif
1409: }
1410: *nrm = sqrt(sum);
1411: } else if (type == NORM_1) {
1412: PetscReal *tmp;
1413: PetscInt *jj = a->j;
1414: PetscMalloc((A->cmap.n+1)*sizeof(PetscReal),&tmp);
1415: PetscMemzero(tmp,A->cmap.n*sizeof(PetscReal));
1416: *nrm = 0.0;
1417: for (j=0; j<a->nz; j++) {
1418: tmp[*jj++] += PetscAbsScalar(*v); v++;
1419: }
1420: for (j=0; j<A->cmap.n; j++) {
1421: if (tmp[j] > *nrm) *nrm = tmp[j];
1422: }
1423: PetscFree(tmp);
1424: } else if (type == NORM_INFINITY) {
1425: *nrm = 0.0;
1426: for (j=0; j<A->rmap.n; j++) {
1427: v = a->a + a->i[j];
1428: sum = 0.0;
1429: for (i=0; i<a->i[j+1]-a->i[j]; i++) {
1430: sum += PetscAbsScalar(*v); v++;
1431: }
1432: if (sum > *nrm) *nrm = sum;
1433: }
1434: } else {
1435: SETERRQ(PETSC_ERR_SUP,"No support for two norm");
1436: }
1437: return(0);
1438: }
1442: PetscErrorCode MatTranspose_SeqAIJ(Mat A,Mat *B)
1443: {
1444: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1445: Mat C;
1447: PetscInt i,*aj = a->j,*ai = a->i,m = A->rmap.n,len,*col;
1448: PetscScalar *array = a->a;
1451: if (!B && m != A->cmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Square matrix only for in-place");
1452: PetscMalloc((1+A->cmap.n)*sizeof(PetscInt),&col);
1453: PetscMemzero(col,(1+A->cmap.n)*sizeof(PetscInt));
1454:
1455: for (i=0; i<ai[m]; i++) col[aj[i]] += 1;
1456: MatCreate(((PetscObject)A)->comm,&C);
1457: MatSetSizes(C,A->cmap.n,m,A->cmap.n,m);
1458: MatSetType(C,((PetscObject)A)->type_name);
1459: MatSeqAIJSetPreallocation_SeqAIJ(C,0,col);
1460: PetscFree(col);
1461: for (i=0; i<m; i++) {
1462: len = ai[i+1]-ai[i];
1463: MatSetValues_SeqAIJ(C,len,aj,1,&i,array,INSERT_VALUES);
1464: array += len;
1465: aj += len;
1466: }
1468: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1469: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
1471: if (B) {
1472: *B = C;
1473: } else {
1474: MatHeaderCopy(A,C);
1475: }
1476: return(0);
1477: }
1482: PetscErrorCode MatIsTranspose_SeqAIJ(Mat A,Mat B,PetscReal tol,PetscTruth *f)
1483: {
1484: Mat_SeqAIJ *aij = (Mat_SeqAIJ *) A->data,*bij = (Mat_SeqAIJ*) A->data;
1485: PetscInt *adx,*bdx,*aii,*bii,*aptr,*bptr; PetscScalar *va,*vb;
1487: PetscInt ma,na,mb,nb, i;
1490: bij = (Mat_SeqAIJ *) B->data;
1491:
1492: MatGetSize(A,&ma,&na);
1493: MatGetSize(B,&mb,&nb);
1494: if (ma!=nb || na!=mb){
1495: *f = PETSC_FALSE;
1496: return(0);
1497: }
1498: aii = aij->i; bii = bij->i;
1499: adx = aij->j; bdx = bij->j;
1500: va = aij->a; vb = bij->a;
1501: PetscMalloc(ma*sizeof(PetscInt),&aptr);
1502: PetscMalloc(mb*sizeof(PetscInt),&bptr);
1503: for (i=0; i<ma; i++) aptr[i] = aii[i];
1504: for (i=0; i<mb; i++) bptr[i] = bii[i];
1506: *f = PETSC_TRUE;
1507: for (i=0; i<ma; i++) {
1508: while (aptr[i]<aii[i+1]) {
1509: PetscInt idc,idr;
1510: PetscScalar vc,vr;
1511: /* column/row index/value */
1512: idc = adx[aptr[i]];
1513: idr = bdx[bptr[idc]];
1514: vc = va[aptr[i]];
1515: vr = vb[bptr[idc]];
1516: if (i!=idr || PetscAbsScalar(vc-vr) > tol) {
1517: *f = PETSC_FALSE;
1518: goto done;
1519: } else {
1520: aptr[i]++;
1521: if (B || i!=idc) bptr[idc]++;
1522: }
1523: }
1524: }
1525: done:
1526: PetscFree(aptr);
1527: if (B) {
1528: PetscFree(bptr);
1529: }
1530: return(0);
1531: }
1537: PetscErrorCode MatIsHermitianTranspose_SeqAIJ(Mat A,Mat B,PetscReal tol,PetscTruth *f)
1538: {
1539: Mat_SeqAIJ *aij = (Mat_SeqAIJ *) A->data,*bij = (Mat_SeqAIJ*) A->data;
1540: PetscInt *adx,*bdx,*aii,*bii,*aptr,*bptr; PetscScalar *va,*vb;
1542: PetscInt ma,na,mb,nb, i;
1545: bij = (Mat_SeqAIJ *) B->data;
1546:
1547: MatGetSize(A,&ma,&na);
1548: MatGetSize(B,&mb,&nb);
1549: if (ma!=nb || na!=mb){
1550: *f = PETSC_FALSE;
1551: return(0);
1552: }
1553: aii = aij->i; bii = bij->i;
1554: adx = aij->j; bdx = bij->j;
1555: va = aij->a; vb = bij->a;
1556: PetscMalloc(ma*sizeof(PetscInt),&aptr);
1557: PetscMalloc(mb*sizeof(PetscInt),&bptr);
1558: for (i=0; i<ma; i++) aptr[i] = aii[i];
1559: for (i=0; i<mb; i++) bptr[i] = bii[i];
1561: *f = PETSC_TRUE;
1562: for (i=0; i<ma; i++) {
1563: while (aptr[i]<aii[i+1]) {
1564: PetscInt idc,idr;
1565: PetscScalar vc,vr;
1566: /* column/row index/value */
1567: idc = adx[aptr[i]];
1568: idr = bdx[bptr[idc]];
1569: vc = va[aptr[i]];
1570: vr = vb[bptr[idc]];
1571: if (i!=idr || PetscAbsScalar(vc-PetscConj(vr)) > tol) {
1572: *f = PETSC_FALSE;
1573: goto done;
1574: } else {
1575: aptr[i]++;
1576: if (B || i!=idc) bptr[idc]++;
1577: }
1578: }
1579: }
1580: done:
1581: PetscFree(aptr);
1582: if (B) {
1583: PetscFree(bptr);
1584: }
1585: return(0);
1586: }
1591: PetscErrorCode MatIsSymmetric_SeqAIJ(Mat A,PetscReal tol,PetscTruth *f)
1592: {
1595: MatIsTranspose_SeqAIJ(A,A,tol,f);
1596: return(0);
1597: }
1601: PetscErrorCode MatIsHermitian_SeqAIJ(Mat A,PetscReal tol,PetscTruth *f)
1602: {
1605: MatIsHermitianTranspose_SeqAIJ(A,A,tol,f);
1606: return(0);
1607: }
1611: PetscErrorCode MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr)
1612: {
1613: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1614: PetscScalar *l,*r,x,*v;
1616: PetscInt i,j,m = A->rmap.n,n = A->cmap.n,M,nz = a->nz,*jj;
1619: if (ll) {
1620: /* The local size is used so that VecMPI can be passed to this routine
1621: by MatDiagonalScale_MPIAIJ */
1622: VecGetLocalSize(ll,&m);
1623: if (m != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length");
1624: VecGetArray(ll,&l);
1625: v = a->a;
1626: for (i=0; i<m; i++) {
1627: x = l[i];
1628: M = a->i[i+1] - a->i[i];
1629: for (j=0; j<M; j++) { (*v++) *= x;}
1630: }
1631: VecRestoreArray(ll,&l);
1632: PetscLogFlops(nz);
1633: }
1634: if (rr) {
1635: VecGetLocalSize(rr,&n);
1636: if (n != A->cmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vector wrong length");
1637: VecGetArray(rr,&r);
1638: v = a->a; jj = a->j;
1639: for (i=0; i<nz; i++) {
1640: (*v++) *= r[*jj++];
1641: }
1642: VecRestoreArray(rr,&r);
1643: PetscLogFlops(nz);
1644: }
1645: return(0);
1646: }
1650: PetscErrorCode MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,PetscInt csize,MatReuse scall,Mat *B)
1651: {
1652: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*c;
1654: PetscInt *smap,i,k,kstart,kend,oldcols = A->cmap.n,*lens;
1655: PetscInt row,mat_i,*mat_j,tcol,first,step,*mat_ilen,sum,lensi;
1656: PetscInt *irow,*icol,nrows,ncols;
1657: PetscInt *starts,*j_new,*i_new,*aj = a->j,*ai = a->i,ii,*ailen = a->ilen;
1658: PetscScalar *a_new,*mat_a;
1659: Mat C;
1660: PetscTruth stride;
1663: ISSorted(isrow,(PetscTruth*)&i);
1664: if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"ISrow is not sorted");
1665: ISSorted(iscol,(PetscTruth*)&i);
1666: if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IScol is not sorted");
1668: ISGetIndices(isrow,&irow);
1669: ISGetLocalSize(isrow,&nrows);
1670: ISGetLocalSize(iscol,&ncols);
1672: ISStrideGetInfo(iscol,&first,&step);
1673: ISStride(iscol,&stride);
1674: if (stride && step == 1) {
1675: /* special case of contiguous rows */
1676: PetscMalloc((2*nrows+1)*sizeof(PetscInt),&lens);
1677: starts = lens + nrows;
1678: /* loop over new rows determining lens and starting points */
1679: for (i=0; i<nrows; i++) {
1680: kstart = ai[irow[i]];
1681: kend = kstart + ailen[irow[i]];
1682: for (k=kstart; k<kend; k++) {
1683: if (aj[k] >= first) {
1684: starts[i] = k;
1685: break;
1686: }
1687: }
1688: sum = 0;
1689: while (k < kend) {
1690: if (aj[k++] >= first+ncols) break;
1691: sum++;
1692: }
1693: lens[i] = sum;
1694: }
1695: /* create submatrix */
1696: if (scall == MAT_REUSE_MATRIX) {
1697: PetscInt n_cols,n_rows;
1698: MatGetSize(*B,&n_rows,&n_cols);
1699: if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size");
1700: MatZeroEntries(*B);
1701: C = *B;
1702: } else {
1703: MatCreate(((PetscObject)A)->comm,&C);
1704: MatSetSizes(C,nrows,ncols,PETSC_DETERMINE,PETSC_DETERMINE);
1705: MatSetType(C,((PetscObject)A)->type_name);
1706: MatSeqAIJSetPreallocation_SeqAIJ(C,0,lens);
1707: }
1708: c = (Mat_SeqAIJ*)C->data;
1710: /* loop over rows inserting into submatrix */
1711: a_new = c->a;
1712: j_new = c->j;
1713: i_new = c->i;
1715: for (i=0; i<nrows; i++) {
1716: ii = starts[i];
1717: lensi = lens[i];
1718: for (k=0; k<lensi; k++) {
1719: *j_new++ = aj[ii+k] - first;
1720: }
1721: PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(PetscScalar));
1722: a_new += lensi;
1723: i_new[i+1] = i_new[i] + lensi;
1724: c->ilen[i] = lensi;
1725: }
1726: PetscFree(lens);
1727: } else {
1728: ISGetIndices(iscol,&icol);
1729: PetscMalloc((1+oldcols)*sizeof(PetscInt),&smap);
1730:
1731: PetscMalloc((1+nrows)*sizeof(PetscInt),&lens);
1732: PetscMemzero(smap,oldcols*sizeof(PetscInt));
1733: for (i=0; i<ncols; i++) smap[icol[i]] = i+1;
1734: /* determine lens of each row */
1735: for (i=0; i<nrows; i++) {
1736: kstart = ai[irow[i]];
1737: kend = kstart + a->ilen[irow[i]];
1738: lens[i] = 0;
1739: for (k=kstart; k<kend; k++) {
1740: if (smap[aj[k]]) {
1741: lens[i]++;
1742: }
1743: }
1744: }
1745: /* Create and fill new matrix */
1746: if (scall == MAT_REUSE_MATRIX) {
1747: PetscTruth equal;
1749: c = (Mat_SeqAIJ *)((*B)->data);
1750: if ((*B)->rmap.n != nrows || (*B)->cmap.n != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size");
1751: PetscMemcmp(c->ilen,lens,(*B)->rmap.n*sizeof(PetscInt),&equal);
1752: if (!equal) {
1753: SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros");
1754: }
1755: PetscMemzero(c->ilen,(*B)->rmap.n*sizeof(PetscInt));
1756: C = *B;
1757: } else {
1758: MatCreate(((PetscObject)A)->comm,&C);
1759: MatSetSizes(C,nrows,ncols,PETSC_DETERMINE,PETSC_DETERMINE);
1760: MatSetType(C,((PetscObject)A)->type_name);
1761: MatSeqAIJSetPreallocation_SeqAIJ(C,0,lens);
1762: }
1763: c = (Mat_SeqAIJ *)(C->data);
1764: for (i=0; i<nrows; i++) {
1765: row = irow[i];
1766: kstart = ai[row];
1767: kend = kstart + a->ilen[row];
1768: mat_i = c->i[i];
1769: mat_j = c->j + mat_i;
1770: mat_a = c->a + mat_i;
1771: mat_ilen = c->ilen + i;
1772: for (k=kstart; k<kend; k++) {
1773: if ((tcol=smap[a->j[k]])) {
1774: *mat_j++ = tcol - 1;
1775: *mat_a++ = a->a[k];
1776: (*mat_ilen)++;
1778: }
1779: }
1780: }
1781: /* Free work space */
1782: ISRestoreIndices(iscol,&icol);
1783: PetscFree(smap);
1784: PetscFree(lens);
1785: }
1786: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1787: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
1789: ISRestoreIndices(isrow,&irow);
1790: *B = C;
1791: return(0);
1792: }
1794: /*
1795: */
1798: PetscErrorCode MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,MatFactorInfo *info)
1799: {
1800: Mat_SeqAIJ *a = (Mat_SeqAIJ*)inA->data;
1802: Mat outA;
1803: PetscTruth row_identity,col_identity;
1806: if (info->levels != 0) SETERRQ(PETSC_ERR_SUP,"Only levels=0 supported for in-place ilu");
1807: ISIdentity(row,&row_identity);
1808: ISIdentity(col,&col_identity);
1810: outA = inA;
1811: inA->factor = FACTOR_LU;
1812: PetscObjectReference((PetscObject)row);
1813: if (a->row) { ISDestroy(a->row); }
1814: a->row = row;
1815: PetscObjectReference((PetscObject)col);
1816: if (a->col) { ISDestroy(a->col); }
1817: a->col = col;
1819: /* Create the inverse permutation so that it can be used in MatLUFactorNumeric() */
1820: if (a->icol) {ISDestroy(a->icol);} /* need to remove old one */
1821: ISInvertPermutation(col,PETSC_DECIDE,&a->icol);
1822: PetscLogObjectParent(inA,a->icol);
1824: if (!a->solve_work) { /* this matrix may have been factored before */
1825: PetscMalloc((inA->rmap.n+1)*sizeof(PetscScalar),&a->solve_work);
1826: PetscLogObjectMemory(inA, (inA->rmap.n+1)*sizeof(PetscScalar));
1827: }
1829: MatMarkDiagonal_SeqAIJ(inA);
1830: if (row_identity && col_identity) {
1831: MatLUFactorNumeric_SeqAIJ(inA,info,&outA);
1832: } else {
1833: MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(inA,info,&outA);
1834: }
1835: return(0);
1836: }
1838: #include petscblaslapack.h
1841: PetscErrorCode MatScale_SeqAIJ(Mat inA,PetscScalar alpha)
1842: {
1843: Mat_SeqAIJ *a = (Mat_SeqAIJ*)inA->data;
1844: PetscBLASInt bnz = (PetscBLASInt)a->nz,one = 1;
1845: PetscScalar oalpha = alpha;
1850: BLASscal_(&bnz,&oalpha,a->a,&one);
1851: PetscLogFlops(a->nz);
1852: return(0);
1853: }
1857: PetscErrorCode MatGetSubMatrices_SeqAIJ(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1858: {
1860: PetscInt i;
1863: if (scall == MAT_INITIAL_MATRIX) {
1864: PetscMalloc((n+1)*sizeof(Mat),B);
1865: }
1867: for (i=0; i<n; i++) {
1868: MatGetSubMatrix_SeqAIJ(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);
1869: }
1870: return(0);
1871: }
1875: PetscErrorCode MatIncreaseOverlap_SeqAIJ(Mat A,PetscInt is_max,IS is[],PetscInt ov)
1876: {
1877: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1879: PetscInt row,i,j,k,l,m,n,*idx,*nidx,isz,val;
1880: PetscInt start,end,*ai,*aj;
1881: PetscBT table;
1884: m = A->rmap.n;
1885: ai = a->i;
1886: aj = a->j;
1888: if (ov < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"illegal negative overlap value used");
1890: PetscMalloc((m+1)*sizeof(PetscInt),&nidx);
1891: PetscBTCreate(m,table);
1893: for (i=0; i<is_max; i++) {
1894: /* Initialize the two local arrays */
1895: isz = 0;
1896: PetscBTMemzero(m,table);
1897:
1898: /* Extract the indices, assume there can be duplicate entries */
1899: ISGetIndices(is[i],&idx);
1900: ISGetLocalSize(is[i],&n);
1901:
1902: /* Enter these into the temp arrays. I.e., mark table[row], enter row into new index */
1903: for (j=0; j<n ; ++j){
1904: if(!PetscBTLookupSet(table,idx[j])) { nidx[isz++] = idx[j];}
1905: }
1906: ISRestoreIndices(is[i],&idx);
1907: ISDestroy(is[i]);
1908:
1909: k = 0;
1910: for (j=0; j<ov; j++){ /* for each overlap */
1911: n = isz;
1912: for (; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */
1913: row = nidx[k];
1914: start = ai[row];
1915: end = ai[row+1];
1916: for (l = start; l<end ; l++){
1917: val = aj[l] ;
1918: if (!PetscBTLookupSet(table,val)) {nidx[isz++] = val;}
1919: }
1920: }
1921: }
1922: ISCreateGeneral(PETSC_COMM_SELF,isz,nidx,(is+i));
1923: }
1924: PetscBTDestroy(table);
1925: PetscFree(nidx);
1926: return(0);
1927: }
1929: /* -------------------------------------------------------------- */
1932: PetscErrorCode MatPermute_SeqAIJ(Mat A,IS rowp,IS colp,Mat *B)
1933: {
1934: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1936: PetscInt i,nz,m = A->rmap.n,n = A->cmap.n,*col;
1937: PetscInt *row,*cnew,j,*lens;
1938: IS icolp,irowp;
1939: PetscInt *cwork;
1940: PetscScalar *vwork;
1943: ISInvertPermutation(rowp,PETSC_DECIDE,&irowp);
1944: ISGetIndices(irowp,&row);
1945: ISInvertPermutation(colp,PETSC_DECIDE,&icolp);
1946: ISGetIndices(icolp,&col);
1947:
1948: /* determine lengths of permuted rows */
1949: PetscMalloc((m+1)*sizeof(PetscInt),&lens);
1950: for (i=0; i<m; i++) {
1951: lens[row[i]] = a->i[i+1] - a->i[i];
1952: }
1953: MatCreate(((PetscObject)A)->comm,B);
1954: MatSetSizes(*B,m,n,m,n);
1955: MatSetType(*B,((PetscObject)A)->type_name);
1956: MatSeqAIJSetPreallocation_SeqAIJ(*B,0,lens);
1957: PetscFree(lens);
1959: PetscMalloc(n*sizeof(PetscInt),&cnew);
1960: for (i=0; i<m; i++) {
1961: MatGetRow_SeqAIJ(A,i,&nz,&cwork,&vwork);
1962: for (j=0; j<nz; j++) { cnew[j] = col[cwork[j]];}
1963: MatSetValues_SeqAIJ(*B,1,&row[i],nz,cnew,vwork,INSERT_VALUES);
1964: MatRestoreRow_SeqAIJ(A,i,&nz,&cwork,&vwork);
1965: }
1966: PetscFree(cnew);
1967: (*B)->assembled = PETSC_FALSE;
1968: MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);
1969: MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);
1970: ISRestoreIndices(irowp,&row);
1971: ISRestoreIndices(icolp,&col);
1972: ISDestroy(irowp);
1973: ISDestroy(icolp);
1974: return(0);
1975: }
1979: PetscErrorCode MatCopy_SeqAIJ(Mat A,Mat B,MatStructure str)
1980: {
1984: /* If the two matrices have the same copy implementation, use fast copy. */
1985: if (str == SAME_NONZERO_PATTERN && (A->ops->copy == B->ops->copy)) {
1986: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1987: Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data;
1989: if (a->i[A->rmap.n] != b->i[B->rmap.n]) {
1990: SETERRQ(PETSC_ERR_ARG_INCOMP,"Number of nonzeros in two matrices are different");
1991: }
1992: PetscMemcpy(b->a,a->a,(a->i[A->rmap.n])*sizeof(PetscScalar));
1993: } else {
1994: MatCopy_Basic(A,B,str);
1995: }
1996: return(0);
1997: }
2001: PetscErrorCode MatSetUpPreallocation_SeqAIJ(Mat A)
2002: {
2006: MatSeqAIJSetPreallocation_SeqAIJ(A,PETSC_DEFAULT,0);
2007: return(0);
2008: }
2012: PetscErrorCode MatGetArray_SeqAIJ(Mat A,PetscScalar *array[])
2013: {
2014: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
2016: *array = a->a;
2017: return(0);
2018: }
2022: PetscErrorCode MatRestoreArray_SeqAIJ(Mat A,PetscScalar *array[])
2023: {
2025: return(0);
2026: }
2030: PetscErrorCode MatFDColoringApply_SeqAIJ(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
2031: {
2032: PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void *))coloring->f;
2034: PetscInt k,N,start,end,l,row,col,srow,**vscaleforrow,m1,m2;
2035: PetscScalar dx,*y,*xx,*w3_array;
2036: PetscScalar *vscale_array;
2037: PetscReal epsilon = coloring->error_rel,umin = coloring->umin;
2038: Vec w1,w2,w3;
2039: void *fctx = coloring->fctx;
2040: PetscTruth flg;
2043: if (!coloring->w1) {
2044: VecDuplicate(x1,&coloring->w1);
2045: PetscLogObjectParent(coloring,coloring->w1);
2046: VecDuplicate(x1,&coloring->w2);
2047: PetscLogObjectParent(coloring,coloring->w2);
2048: VecDuplicate(x1,&coloring->w3);
2049: PetscLogObjectParent(coloring,coloring->w3);
2050: }
2051: w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
2053: MatSetUnfactored(J);
2054: PetscOptionsHasName(((PetscObject)coloring)->prefix,"-mat_fd_coloring_dont_rezero",&flg);
2055: if (flg) {
2056: PetscInfo(coloring,"Not calling MatZeroEntries()\n");
2057: } else {
2058: PetscTruth assembled;
2059: MatAssembled(J,&assembled);
2060: if (assembled) {
2061: MatZeroEntries(J);
2062: }
2063: }
2065: VecGetOwnershipRange(x1,&start,&end);
2066: VecGetSize(x1,&N);
2068: /*
2069: This is a horrible, horrible, hack. See DMMGComputeJacobian_Multigrid() it inproperly sets
2070: coloring->F for the coarser grids from the finest
2071: */
2072: if (coloring->F) {
2073: VecGetLocalSize(coloring->F,&m1);
2074: VecGetLocalSize(w1,&m2);
2075: if (m1 != m2) {
2076: coloring->F = 0;
2077: }
2078: }
2080: if (coloring->F) {
2081: w1 = coloring->F;
2082: coloring->F = 0;
2083: } else {
2085: (*f)(sctx,x1,w1,fctx);
2087: }
2089: /*
2090: Compute all the scale factors and share with other processors
2091: */
2092: VecGetArray(x1,&xx);xx = xx - start;
2093: VecGetArray(coloring->vscale,&vscale_array);vscale_array = vscale_array - start;
2094: for (k=0; k<coloring->ncolors; k++) {
2095: /*
2096: Loop over each column associated with color adding the
2097: perturbation to the vector w3.
2098: */
2099: for (l=0; l<coloring->ncolumns[k]; l++) {
2100: col = coloring->columns[k][l]; /* column of the matrix we are probing for */
2101: dx = xx[col];
2102: if (dx == 0.0) dx = 1.0;
2103: #if !defined(PETSC_USE_COMPLEX)
2104: if (dx < umin && dx >= 0.0) dx = umin;
2105: else if (dx < 0.0 && dx > -umin) dx = -umin;
2106: #else
2107: if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin;
2108: else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
2109: #endif
2110: dx *= epsilon;
2111: vscale_array[col] = 1.0/dx;
2112: }
2113: }
2114: vscale_array = vscale_array + start;VecRestoreArray(coloring->vscale,&vscale_array);
2115: VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);
2116: VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);
2118: /* VecView(coloring->vscale,PETSC_VIEWER_STDOUT_WORLD);
2119: VecView(x1,PETSC_VIEWER_STDOUT_WORLD);*/
2121: if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
2122: else vscaleforrow = coloring->columnsforrow;
2124: VecGetArray(coloring->vscale,&vscale_array);
2125: /*
2126: Loop over each color
2127: */
2128: for (k=0; k<coloring->ncolors; k++) {
2129: coloring->currentcolor = k;
2130: VecCopy(x1,w3);
2131: VecGetArray(w3,&w3_array);w3_array = w3_array - start;
2132: /*
2133: Loop over each column associated with color adding the
2134: perturbation to the vector w3.
2135: */
2136: for (l=0; l<coloring->ncolumns[k]; l++) {
2137: col = coloring->columns[k][l]; /* column of the matrix we are probing for */
2138: dx = xx[col];
2139: if (dx == 0.0) dx = 1.0;
2140: #if !defined(PETSC_USE_COMPLEX)
2141: if (dx < umin && dx >= 0.0) dx = umin;
2142: else if (dx < 0.0 && dx > -umin) dx = -umin;
2143: #else
2144: if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin;
2145: else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
2146: #endif
2147: dx *= epsilon;
2148: if (!PetscAbsScalar(dx)) SETERRQ(PETSC_ERR_PLIB,"Computed 0 differencing parameter");
2149: w3_array[col] += dx;
2150: }
2151: w3_array = w3_array + start; VecRestoreArray(w3,&w3_array);
2153: /*
2154: Evaluate function at x1 + dx (here dx is a vector of perturbations)
2155: */
2158: (*f)(sctx,w3,w2,fctx);
2160: VecAXPY(w2,-1.0,w1);
2162: /*
2163: Loop over rows of vector, putting results into Jacobian matrix
2164: */
2165: VecGetArray(w2,&y);
2166: for (l=0; l<coloring->nrows[k]; l++) {
2167: row = coloring->rows[k][l];
2168: col = coloring->columnsforrow[k][l];
2169: y[row] *= vscale_array[vscaleforrow[k][l]];
2170: srow = row + start;
2171: MatSetValues_SeqAIJ(J,1,&srow,1,&col,y+row,INSERT_VALUES);
2172: }
2173: VecRestoreArray(w2,&y);
2174: }
2175: coloring->currentcolor = k;
2176: VecRestoreArray(coloring->vscale,&vscale_array);
2177: xx = xx + start; VecRestoreArray(x1,&xx);
2178: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
2179: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
2180: return(0);
2181: }
2183: #include petscblaslapack.h
2186: PetscErrorCode MatAXPY_SeqAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str)
2187: {
2189: PetscInt i;
2190: Mat_SeqAIJ *x = (Mat_SeqAIJ *)X->data,*y = (Mat_SeqAIJ *)Y->data;
2191: PetscBLASInt one=1,bnz = (PetscBLASInt)x->nz;
2194: if (str == SAME_NONZERO_PATTERN) {
2195: PetscScalar alpha = a;
2196: BLASaxpy_(&bnz,&alpha,x->a,&one,y->a,&one);
2197: } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
2198: if (y->xtoy && y->XtoY != X) {
2199: PetscFree(y->xtoy);
2200: MatDestroy(y->XtoY);
2201: }
2202: if (!y->xtoy) { /* get xtoy */
2203: MatAXPYGetxtoy_Private(X->rmap.n,x->i,x->j,PETSC_NULL, y->i,y->j,PETSC_NULL, &y->xtoy);
2204: y->XtoY = X;
2205: PetscObjectReference((PetscObject)X);
2206: }
2207: for (i=0; i<x->nz; i++) y->a[y->xtoy[i]] += a*(x->a[i]);
2208: PetscInfo3(Y,"ratio of nnz(X)/nnz(Y): %d/%d = %G\n",x->nz,y->nz,(PetscReal)(x->nz)/y->nz);
2209: } else {
2210: MatAXPY_Basic(Y,a,X,str);
2211: }
2212: return(0);
2213: }
2217: PetscErrorCode MatSetBlockSize_SeqAIJ(Mat A,PetscInt bs)
2218: {
2220: return(0);
2221: }
2225: PetscErrorCode MatConjugate_SeqAIJ(Mat mat)
2226: {
2227: #if defined(PETSC_USE_COMPLEX)
2228: Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
2229: PetscInt i,nz;
2230: PetscScalar *a;
2233: nz = aij->nz;
2234: a = aij->a;
2235: for (i=0; i<nz; i++) {
2236: a[i] = PetscConj(a[i]);
2237: }
2238: #else
2240: #endif
2241: return(0);
2242: }
2246: PetscErrorCode MatGetRowMaxAbs_SeqAIJ(Mat A,Vec v,PetscInt idx[])
2247: {
2248: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
2250: PetscInt i,j,m = A->rmap.n,*ai,*aj,ncols,n;
2251: PetscReal atmp;
2252: PetscScalar *x;
2253: MatScalar *aa;
2256: if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2257: aa = a->a;
2258: ai = a->i;
2259: aj = a->j;
2261: VecSet(v,0.0);
2262: VecGetArray(v,&x);
2263: VecGetLocalSize(v,&n);
2264: if (n != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
2265: for (i=0; i<m; i++) {
2266: ncols = ai[1] - ai[0]; ai++;
2267: x[i] = 0.0; if (idx) idx[i] = 0;
2268: for (j=0; j<ncols; j++){
2269: atmp = PetscAbsScalar(*aa);
2270: if (PetscAbsScalar(x[i]) < atmp) {x[i] = atmp; if (idx) idx[i] = *aj;}
2271: aa++; aj++;
2272: }
2273: }
2274: VecRestoreArray(v,&x);
2275: return(0);
2276: }
2280: PetscErrorCode MatGetRowMax_SeqAIJ(Mat A,Vec v,PetscInt idx[])
2281: {
2282: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
2284: PetscInt i,j,m = A->rmap.n,*ai,*aj,ncols,n;
2285: PetscScalar *x;
2286: MatScalar *aa;
2289: if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2290: aa = a->a;
2291: ai = a->i;
2292: aj = a->j;
2294: VecSet(v,0.0);
2295: VecGetArray(v,&x);
2296: VecGetLocalSize(v,&n);
2297: if (n != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
2298: for (i=0; i<m; i++) {
2299: ncols = ai[1] - ai[0]; ai++;
2300: if (ncols == A->cmap.n) { /* row is dense */
2301: x[i] = *aa; if (idx) idx[i] = 0;
2302: } else { /* row is sparse so already KNOW maximum is 0.0 or higher */
2303: x[i] = 0.0;
2304: if (idx) {
2305: idx[i] = 0; /* in case ncols is zero */
2306: for (j=0;j<ncols;j++) { /* find first implicit 0.0 in the row */
2307: if (aj[j] > j) {
2308: idx[i] = j;
2309: break;
2310: }
2311: }
2312: }
2313: }
2314: for (j=0; j<ncols; j++){
2315: if (PetscRealPart(x[i]) < PetscRealPart(*aa)) {x[i] = *aa; if (idx) idx[i] = *aj;}
2316: aa++; aj++;
2317: }
2318: }
2319: VecRestoreArray(v,&x);
2320: return(0);
2321: }
2325: PetscErrorCode MatGetRowMin_SeqAIJ(Mat A,Vec v,PetscInt idx[])
2326: {
2327: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
2329: PetscInt i,j,m = A->rmap.n,*ai,*aj,ncols,n;
2330: PetscScalar *x;
2331: MatScalar *aa;
2334: if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2335: aa = a->a;
2336: ai = a->i;
2337: aj = a->j;
2339: VecSet(v,0.0);
2340: VecGetArray(v,&x);
2341: VecGetLocalSize(v,&n);
2342: if (n != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
2343: for (i=0; i<m; i++) {
2344: ncols = ai[1] - ai[0]; ai++;
2345: if (ncols == A->cmap.n) { /* row is dense */
2346: x[i] = *aa; if (idx) idx[i] = 0;
2347: } else { /* row is sparse so already KNOW minimum is 0.0 or lower */
2348: x[i] = 0.0;
2349: if (idx) { /* find first implicit 0.0 in the row */
2350: idx[i] = 0; /* in case ncols is zero */
2351: for (j=0;j<ncols;j++) {
2352: if (aj[j] > j) {
2353: idx[i] = j;
2354: break;
2355: }
2356: }
2357: }
2358: }
2359: for (j=0; j<ncols; j++){
2360: if (PetscRealPart(x[i]) > PetscRealPart(*aa)) {x[i] = *aa; if (idx) idx[i] = *aj;}
2361: aa++; aj++;
2362: }
2363: }
2364: VecRestoreArray(v,&x);
2365: return(0);
2366: }
2368: /* -------------------------------------------------------------------*/
2369: static struct _MatOps MatOps_Values = {MatSetValues_SeqAIJ,
2370: MatGetRow_SeqAIJ,
2371: MatRestoreRow_SeqAIJ,
2372: MatMult_SeqAIJ,
2373: /* 4*/ MatMultAdd_SeqAIJ,
2374: MatMultTranspose_SeqAIJ,
2375: MatMultTransposeAdd_SeqAIJ,
2376: MatSolve_SeqAIJ,
2377: MatSolveAdd_SeqAIJ,
2378: MatSolveTranspose_SeqAIJ,
2379: /*10*/ MatSolveTransposeAdd_SeqAIJ,
2380: MatLUFactor_SeqAIJ,
2381: 0,
2382: MatRelax_SeqAIJ,
2383: MatTranspose_SeqAIJ,
2384: /*15*/ MatGetInfo_SeqAIJ,
2385: MatEqual_SeqAIJ,
2386: MatGetDiagonal_SeqAIJ,
2387: MatDiagonalScale_SeqAIJ,
2388: MatNorm_SeqAIJ,
2389: /*20*/ 0,
2390: MatAssemblyEnd_SeqAIJ,
2391: MatCompress_SeqAIJ,
2392: MatSetOption_SeqAIJ,
2393: MatZeroEntries_SeqAIJ,
2394: /*25*/ MatZeroRows_SeqAIJ,
2395: MatLUFactorSymbolic_SeqAIJ,
2396: MatLUFactorNumeric_SeqAIJ,
2397: MatCholeskyFactorSymbolic_SeqAIJ,
2398: MatCholeskyFactorNumeric_SeqAIJ,
2399: /*30*/ MatSetUpPreallocation_SeqAIJ,
2400: MatILUFactorSymbolic_SeqAIJ,
2401: MatICCFactorSymbolic_SeqAIJ,
2402: MatGetArray_SeqAIJ,
2403: MatRestoreArray_SeqAIJ,
2404: /*35*/ MatDuplicate_SeqAIJ,
2405: 0,
2406: 0,
2407: MatILUFactor_SeqAIJ,
2408: 0,
2409: /*40*/ MatAXPY_SeqAIJ,
2410: MatGetSubMatrices_SeqAIJ,
2411: MatIncreaseOverlap_SeqAIJ,
2412: MatGetValues_SeqAIJ,
2413: MatCopy_SeqAIJ,
2414: /*45*/ MatGetRowMax_SeqAIJ,
2415: MatScale_SeqAIJ,
2416: 0,
2417: MatDiagonalSet_SeqAIJ,
2418: MatILUDTFactor_SeqAIJ,
2419: /*50*/ MatSetBlockSize_SeqAIJ,
2420: MatGetRowIJ_SeqAIJ,
2421: MatRestoreRowIJ_SeqAIJ,
2422: MatGetColumnIJ_SeqAIJ,
2423: MatRestoreColumnIJ_SeqAIJ,
2424: /*55*/ MatFDColoringCreate_SeqAIJ,
2425: 0,
2426: 0,
2427: MatPermute_SeqAIJ,
2428: 0,
2429: /*60*/ 0,
2430: MatDestroy_SeqAIJ,
2431: MatView_SeqAIJ,
2432: 0,
2433: 0,
2434: /*65*/ 0,
2435: 0,
2436: 0,
2437: 0,
2438: 0,
2439: /*70*/ MatGetRowMaxAbs_SeqAIJ,
2440: 0,
2441: MatSetColoring_SeqAIJ,
2442: #if defined(PETSC_HAVE_ADIC)
2443: MatSetValuesAdic_SeqAIJ,
2444: #else
2445: 0,
2446: #endif
2447: MatSetValuesAdifor_SeqAIJ,
2448: /*75*/ 0,
2449: 0,
2450: 0,
2451: 0,
2452: 0,
2453: /*80*/ 0,
2454: 0,
2455: 0,
2456: 0,
2457: MatLoad_SeqAIJ,
2458: /*85*/ MatIsSymmetric_SeqAIJ,
2459: MatIsHermitian_SeqAIJ,
2460: 0,
2461: 0,
2462: 0,
2463: /*90*/ MatMatMult_SeqAIJ_SeqAIJ,
2464: MatMatMultSymbolic_SeqAIJ_SeqAIJ,
2465: MatMatMultNumeric_SeqAIJ_SeqAIJ,
2466: MatPtAP_Basic,
2467: MatPtAPSymbolic_SeqAIJ,
2468: /*95*/ MatPtAPNumeric_SeqAIJ,
2469: MatMatMultTranspose_SeqAIJ_SeqAIJ,
2470: MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ,
2471: MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ,
2472: MatPtAPSymbolic_SeqAIJ_SeqAIJ,
2473: /*100*/MatPtAPNumeric_SeqAIJ_SeqAIJ,
2474: 0,
2475: 0,
2476: MatConjugate_SeqAIJ,
2477: 0,
2478: /*105*/MatSetValuesRow_SeqAIJ,
2479: MatRealPart_SeqAIJ,
2480: MatImaginaryPart_SeqAIJ,
2481: 0,
2482: 0,
2483: /*110*/MatMatSolve_SeqAIJ,
2484: 0,
2485: MatGetRowMin_SeqAIJ,
2486: 0,
2487: MatMissingDiagonal_SeqAIJ
2488: };
2493: PetscErrorCode MatSeqAIJSetColumnIndices_SeqAIJ(Mat mat,PetscInt *indices)
2494: {
2495: Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
2496: PetscInt i,nz,n;
2500: nz = aij->maxnz;
2501: n = mat->rmap.n;
2502: for (i=0; i<nz; i++) {
2503: aij->j[i] = indices[i];
2504: }
2505: aij->nz = nz;
2506: for (i=0; i<n; i++) {
2507: aij->ilen[i] = aij->imax[i];
2508: }
2510: return(0);
2511: }
2516: /*@
2517: MatSeqAIJSetColumnIndices - Set the column indices for all the rows
2518: in the matrix.
2520: Input Parameters:
2521: + mat - the SeqAIJ matrix
2522: - indices - the column indices
2524: Level: advanced
2526: Notes:
2527: This can be called if you have precomputed the nonzero structure of the
2528: matrix and want to provide it to the matrix object to improve the performance
2529: of the MatSetValues() operation.
2531: You MUST have set the correct numbers of nonzeros per row in the call to
2532: MatCreateSeqAIJ(), and the columns indices MUST be sorted.
2534: MUST be called before any calls to MatSetValues();
2536: The indices should start with zero, not one.
2538: @*/
2539: PetscErrorCode MatSeqAIJSetColumnIndices(Mat mat,PetscInt *indices)
2540: {
2541: PetscErrorCode ierr,(*f)(Mat,PetscInt *);
2546: PetscObjectQueryFunction((PetscObject)mat,"MatSeqAIJSetColumnIndices_C",(void (**)(void))&f);
2547: if (f) {
2548: (*f)(mat,indices);
2549: } else {
2550: SETERRQ(PETSC_ERR_SUP,"Wrong type of matrix to set column indices");
2551: }
2552: return(0);
2553: }
2555: /* ----------------------------------------------------------------------------------------*/
2560: PetscErrorCode MatStoreValues_SeqAIJ(Mat mat)
2561: {
2562: Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
2564: size_t nz = aij->i[mat->rmap.n];
2567: if (aij->nonew != 1) {
2568: SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first");
2569: }
2571: /* allocate space for values if not already there */
2572: if (!aij->saved_values) {
2573: PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);
2574: PetscLogObjectMemory(mat,(nz+1)*sizeof(PetscScalar));
2575: }
2577: /* copy values over */
2578: PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));
2579: return(0);
2580: }
2585: /*@
2586: MatStoreValues - Stashes a copy of the matrix values; this allows, for
2587: example, reuse of the linear part of a Jacobian, while recomputing the
2588: nonlinear portion.
2590: Collect on Mat
2592: Input Parameters:
2593: . mat - the matrix (currently only AIJ matrices support this option)
2595: Level: advanced
2597: Common Usage, with SNESSolve():
2598: $ Create Jacobian matrix
2599: $ Set linear terms into matrix
2600: $ Apply boundary conditions to matrix, at this time matrix must have
2601: $ final nonzero structure (i.e. setting the nonlinear terms and applying
2602: $ boundary conditions again will not change the nonzero structure
2603: $ MatSetOption(mat,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);
2604: $ MatStoreValues(mat);
2605: $ Call SNESSetJacobian() with matrix
2606: $ In your Jacobian routine
2607: $ MatRetrieveValues(mat);
2608: $ Set nonlinear terms in matrix
2609:
2610: Common Usage without SNESSolve(), i.e. when you handle nonlinear solve yourself:
2611: $ // build linear portion of Jacobian
2612: $ MatSetOption(mat,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);
2613: $ MatStoreValues(mat);
2614: $ loop over nonlinear iterations
2615: $ MatRetrieveValues(mat);
2616: $ // call MatSetValues(mat,...) to set nonliner portion of Jacobian
2617: $ // call MatAssemblyBegin/End() on matrix
2618: $ Solve linear system with Jacobian
2619: $ endloop
2621: Notes:
2622: Matrix must already be assemblied before calling this routine
2623: Must set the matrix option MatSetOption(mat,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE); before
2624: calling this routine.
2626: When this is called multiple times it overwrites the previous set of stored values
2627: and does not allocated additional space.
2629: .seealso: MatRetrieveValues()
2631: @*/
2632: PetscErrorCode MatStoreValues(Mat mat)
2633: {
2634: PetscErrorCode ierr,(*f)(Mat);
2638: if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2639: if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2641: PetscObjectQueryFunction((PetscObject)mat,"MatStoreValues_C",(void (**)(void))&f);
2642: if (f) {
2643: (*f)(mat);
2644: } else {
2645: SETERRQ(PETSC_ERR_SUP,"Wrong type of matrix to store values");
2646: }
2647: return(0);
2648: }
2653: PetscErrorCode MatRetrieveValues_SeqAIJ(Mat mat)
2654: {
2655: Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
2657: PetscInt nz = aij->i[mat->rmap.n];
2660: if (aij->nonew != 1) {
2661: SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first");
2662: }
2663: if (!aij->saved_values) {
2664: SETERRQ(PETSC_ERR_ORDER,"Must call MatStoreValues(A);first");
2665: }
2666: /* copy values over */
2667: PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));
2668: return(0);
2669: }
2674: /*@
2675: MatRetrieveValues - Retrieves the copy of the matrix values; this allows, for
2676: example, reuse of the linear part of a Jacobian, while recomputing the
2677: nonlinear portion.
2679: Collect on Mat
2681: Input Parameters:
2682: . mat - the matrix (currently on AIJ matrices support this option)
2684: Level: advanced
2686: .seealso: MatStoreValues()
2688: @*/
2689: PetscErrorCode MatRetrieveValues(Mat mat)
2690: {
2691: PetscErrorCode ierr,(*f)(Mat);
2695: if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2696: if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2698: PetscObjectQueryFunction((PetscObject)mat,"MatRetrieveValues_C",(void (**)(void))&f);
2699: if (f) {
2700: (*f)(mat);
2701: } else {
2702: SETERRQ(PETSC_ERR_SUP,"Wrong type of matrix to retrieve values");
2703: }
2704: return(0);
2705: }
2708: /* --------------------------------------------------------------------------------*/
2711: /*@C
2712: MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format
2713: (the default parallel PETSc format). For good matrix assembly performance
2714: the user should preallocate the matrix storage by setting the parameter nz
2715: (or the array nnz). By setting these parameters accurately, performance
2716: during matrix assembly can be increased by more than a factor of 50.
2718: Collective on MPI_Comm
2720: Input Parameters:
2721: + comm - MPI communicator, set to PETSC_COMM_SELF
2722: . m - number of rows
2723: . n - number of columns
2724: . nz - number of nonzeros per row (same for all rows)
2725: - nnz - array containing the number of nonzeros in the various rows
2726: (possibly different for each row) or PETSC_NULL
2728: Output Parameter:
2729: . A - the matrix
2731: It is recommended that one use the MatCreate(), MatSetType() and/or MatSetFromOptions(),
2732: MatXXXXSetPreallocation() paradgm instead of this routine directly. This is definitely
2733: true if you plan to use the external direct solvers such as SuperLU, MUMPS or Spooles.
2734: [MatXXXXSetPreallocation() is, for example, MatSeqAIJSetPreallocation]
2736: Notes:
2737: If nnz is given then nz is ignored
2739: The AIJ format (also called the Yale sparse matrix format or
2740: compressed row storage), is fully compatible with standard Fortran 77
2741: storage. That is, the stored row and column indices can begin at
2742: either one (as in Fortran) or zero. See the users' manual for details.
2744: Specify the preallocated storage with either nz or nnz (not both).
2745: Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
2746: allocation. For large problems you MUST preallocate memory or you
2747: will get TERRIBLE performance, see the users' manual chapter on matrices.
2749: By default, this format uses inodes (identical nodes) when possible, to
2750: improve numerical efficiency of matrix-vector products and solves. We
2751: search for consecutive rows with the same nonzero structure, thereby
2752: reusing matrix information to achieve increased efficiency.
2754: Options Database Keys:
2755: + -mat_no_inode - Do not use inodes
2756: . -mat_inode_limit <limit> - Sets inode limit (max limit=5)
2757: - -mat_aij_oneindex - Internally use indexing starting at 1
2758: rather than 0. Note that when calling MatSetValues(),
2759: the user still MUST index entries starting at 0!
2761: Level: intermediate
2763: .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues(), MatSeqAIJSetColumnIndices(), MatCreateSeqAIJWithArrays()
2765: @*/
2766: PetscErrorCode MatCreateSeqAIJ(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A)
2767: {
2771: MatCreate(comm,A);
2772: MatSetSizes(*A,m,n,m,n);
2773: MatSetType(*A,MATSEQAIJ);
2774: MatSeqAIJSetPreallocation_SeqAIJ(*A,nz,(PetscInt*)nnz);
2775: return(0);
2776: }
2780: /*@C
2781: MatSeqAIJSetPreallocation - For good matrix assembly performance
2782: the user should preallocate the matrix storage by setting the parameter nz
2783: (or the array nnz). By setting these parameters accurately, performance
2784: during matrix assembly can be increased by more than a factor of 50.
2786: Collective on MPI_Comm
2788: Input Parameters:
2789: + B - The matrix
2790: . nz - number of nonzeros per row (same for all rows)
2791: - nnz - array containing the number of nonzeros in the various rows
2792: (possibly different for each row) or PETSC_NULL
2794: Notes:
2795: If nnz is given then nz is ignored
2797: The AIJ format (also called the Yale sparse matrix format or
2798: compressed row storage), is fully compatible with standard Fortran 77
2799: storage. That is, the stored row and column indices can begin at
2800: either one (as in Fortran) or zero. See the users' manual for details.
2802: Specify the preallocated storage with either nz or nnz (not both).
2803: Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
2804: allocation. For large problems you MUST preallocate memory or you
2805: will get TERRIBLE performance, see the users' manual chapter on matrices.
2807: You can call MatGetInfo() to get information on how effective the preallocation was;
2808: for example the fields mallocs,nz_allocated,nz_used,nz_unneeded;
2809: You can also run with the option -info and look for messages with the string
2810: malloc in them to see if additional memory allocation was needed.
2812: Developers: Use nz of MAT_SKIP_ALLOCATION to not allocate any space for the matrix
2813: entries or columns indices
2815: By default, this format uses inodes (identical nodes) when possible, to
2816: improve numerical efficiency of matrix-vector products and solves. We
2817: search for consecutive rows with the same nonzero structure, thereby
2818: reusing matrix information to achieve increased efficiency.
2820: Options Database Keys:
2821: + -mat_no_inode - Do not use inodes
2822: . -mat_inode_limit <limit> - Sets inode limit (max limit=5)
2823: - -mat_aij_oneindex - Internally use indexing starting at 1
2824: rather than 0. Note that when calling MatSetValues(),
2825: the user still MUST index entries starting at 0!
2827: Level: intermediate
2829: .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues(), MatSeqAIJSetColumnIndices(), MatCreateSeqAIJWithArrays(), MatGetInfo()
2831: @*/
2832: PetscErrorCode MatSeqAIJSetPreallocation(Mat B,PetscInt nz,const PetscInt nnz[])
2833: {
2834: PetscErrorCode ierr,(*f)(Mat,PetscInt,const PetscInt[]);
2837: PetscObjectQueryFunction((PetscObject)B,"MatSeqAIJSetPreallocation_C",(void (**)(void))&f);
2838: if (f) {
2839: (*f)(B,nz,nnz);
2840: }
2841: return(0);
2842: }
2847: PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat B,PetscInt nz,PetscInt *nnz)
2848: {
2849: Mat_SeqAIJ *b;
2850: PetscTruth skipallocation = PETSC_FALSE;
2852: PetscInt i;
2855:
2856: if (nz == MAT_SKIP_ALLOCATION) {
2857: skipallocation = PETSC_TRUE;
2858: nz = 0;
2859: }
2861: B->rmap.bs = B->cmap.bs = 1;
2862: PetscMapSetUp(&B->rmap);
2863: PetscMapSetUp(&B->cmap);
2865: if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
2866: if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz);
2867: if (nnz) {
2868: for (i=0; i<B->rmap.n; i++) {
2869: if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]);
2870: 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);
2871: }
2872: }
2874: B->preallocated = PETSC_TRUE;
2875: b = (Mat_SeqAIJ*)B->data;
2877: if (!skipallocation) {
2878: if (!b->imax) {
2879: PetscMalloc2(B->rmap.n,PetscInt,&b->imax,B->rmap.n,PetscInt,&b->ilen);
2880: PetscLogObjectMemory(B,2*B->rmap.n*sizeof(PetscInt));
2881: }
2882: if (!nnz) {
2883: if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 10;
2884: else if (nz <= 0) nz = 1;
2885: for (i=0; i<B->rmap.n; i++) b->imax[i] = nz;
2886: nz = nz*B->rmap.n;
2887: } else {
2888: nz = 0;
2889: for (i=0; i<B->rmap.n; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
2890: }
2891: /* b->ilen will count nonzeros in each row so far. */
2892: for (i=0; i<B->rmap.n; i++) { b->ilen[i] = 0; }
2894: /* allocate the matrix space */
2895: MatSeqXAIJFreeAIJ(B,&b->a,&b->j,&b->i);
2896: PetscMalloc3(nz,PetscScalar,&b->a,nz,PetscInt,&b->j,B->rmap.n+1,PetscInt,&b->i);
2897: PetscLogObjectMemory(B,(B->rmap.n+1)*sizeof(PetscInt)+nz*(sizeof(PetscScalar)+sizeof(PetscInt)));
2898: b->i[0] = 0;
2899: for (i=1; i<B->rmap.n+1; i++) {
2900: b->i[i] = b->i[i-1] + b->imax[i-1];
2901: }
2902: b->singlemalloc = PETSC_TRUE;
2903: b->free_a = PETSC_TRUE;
2904: b->free_ij = PETSC_TRUE;
2905: } else {
2906: b->free_a = PETSC_FALSE;
2907: b->free_ij = PETSC_FALSE;
2908: }
2910: b->nz = 0;
2911: b->maxnz = nz;
2912: B->info.nz_unneeded = (double)b->maxnz;
2913: return(0);
2914: }
2917: #undef __FUNCT__
2919: /*@
2920: MatSeqAIJSetPreallocationCSR - Allocates memory for a sparse sequential matrix in AIJ format.
2922: Input Parameters:
2923: + B - the matrix
2924: . i - the indices into j for the start of each row (starts with zero)
2925: . j - the column indices for each row (starts with zero) these must be sorted for each row
2926: - v - optional values in the matrix
2928: Contributed by: Lisandro Dalchin
2930: Level: developer
2932: The i,j,v values are COPIED with this routine; to avoid the copy use MatCreateSeqAIJWithArrays()
2934: .keywords: matrix, aij, compressed row, sparse, sequential
2936: .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatSeqAIJSetPreallocation(), MatCreateSeqAIJ(), SeqAIJ
2937: @*/
2938: PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat B,const PetscInt i[],const PetscInt j[],const PetscScalar v[])
2939: {
2940: PetscErrorCode (*f)(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
2945: PetscObjectQueryFunction((PetscObject)B,"MatSeqAIJSetPreallocationCSR_C",(void (**)(void))&f);
2946: if (f) {
2947: (*f)(B,i,j,v);
2948: }
2949: return(0);
2950: }
2953: #undef __FUNCT__
2955: PetscErrorCode MatSeqAIJSetPreallocationCSR_SeqAIJ(Mat B,const PetscInt Ii[],const PetscInt J[],const PetscScalar v[])
2956: {
2957: PetscInt i;
2958: PetscInt m,n;
2959: PetscInt nz;
2960: PetscInt *nnz, nz_max = 0;
2961: PetscScalar *values;
2965: MatGetSize(B, &m, &n);
2967: if (Ii[0]) {
2968: SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Ii[0] must be 0 it is %D", Ii[0]);
2969: }
2970: PetscMalloc((m+1) * sizeof(PetscInt), &nnz);
2971: for(i = 0; i < m; i++) {
2972: nz = Ii[i+1]- Ii[i];
2973: nz_max = PetscMax(nz_max, nz);
2974: if (nz < 0) {
2975: SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE, "Local row %D has a negative number of columns %D", i, nnz);
2976: }
2977: nnz[i] = nz;
2978: }
2979: MatSeqAIJSetPreallocation(B, 0, nnz);
2980: PetscFree(nnz);
2982: if (v) {
2983: values = (PetscScalar*) v;
2984: } else {
2985: PetscMalloc((nz_max+1)*sizeof(PetscScalar), &values);
2986: PetscMemzero(values, nz_max*sizeof(PetscScalar));
2987: }
2989: for(i = 0; i < m; i++) {
2990: nz = Ii[i+1] - Ii[i];
2991: MatSetValues_SeqAIJ(B, 1, &i, nz, J+Ii[i], values + (v ? Ii[i] : 0), INSERT_VALUES);
2992: }
2994: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
2995: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
2997: if (!v) {
2998: PetscFree(values);
2999: }
3000: return(0);
3001: }
3004: /*MC
3005: MATSEQAIJ - MATSEQAIJ = "seqaij" - A matrix type to be used for sequential sparse matrices,
3006: based on compressed sparse row format.
3008: Options Database Keys:
3009: . -mat_type seqaij - sets the matrix type to "seqaij" during a call to MatSetFromOptions()
3011: Level: beginner
3013: .seealso: MatCreateSeqAIJ(), MatSetFromOptions(), MatSetType(), MatCreate(), MatType
3014: M*/
3023: PetscErrorCode MatCreate_SeqAIJ(Mat B)
3024: {
3025: Mat_SeqAIJ *b;
3027: PetscMPIInt size;
3030: MPI_Comm_size(((PetscObject)B)->comm,&size);
3031: if (size > 1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Comm must be of size 1");
3033: PetscNewLog(B,Mat_SeqAIJ,&b);
3034: B->data = (void*)b;
3035: PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
3036: B->factor = 0;
3037: B->mapping = 0;
3038: b->row = 0;
3039: b->col = 0;
3040: b->icol = 0;
3041: b->reallocs = 0;
3042: b->ignorezeroentries = PETSC_FALSE;
3043: b->roworiented = PETSC_TRUE;
3044: b->nonew = 0;
3045: b->diag = 0;
3046: b->solve_work = 0;
3047: B->spptr = 0;
3048: b->saved_values = 0;
3049: b->idiag = 0;
3050: b->mdiag = 0;
3051: b->ssor_work = 0;
3052: b->omega = 1.0;
3053: b->fshift = 0.0;
3054: b->idiagvalid = PETSC_FALSE;
3055: b->keepzeroedrows = PETSC_FALSE;
3056: b->xtoy = 0;
3057: b->XtoY = 0;
3058: b->compressedrow.use = PETSC_FALSE;
3059: b->compressedrow.nrows = B->rmap.n;
3060: b->compressedrow.i = PETSC_NULL;
3061: b->compressedrow.rindex = PETSC_NULL;
3062: b->compressedrow.checked = PETSC_FALSE;
3063: B->same_nonzero = PETSC_FALSE;
3065: PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);
3066: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqAIJSetColumnIndices_C",
3067: "MatSeqAIJSetColumnIndices_SeqAIJ",
3068: MatSeqAIJSetColumnIndices_SeqAIJ);
3069: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
3070: "MatStoreValues_SeqAIJ",
3071: MatStoreValues_SeqAIJ);
3072: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
3073: "MatRetrieveValues_SeqAIJ",
3074: MatRetrieveValues_SeqAIJ);
3075: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_seqsbaij_C",
3076: "MatConvert_SeqAIJ_SeqSBAIJ",
3077: MatConvert_SeqAIJ_SeqSBAIJ);
3078: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_seqbaij_C",
3079: "MatConvert_SeqAIJ_SeqBAIJ",
3080: MatConvert_SeqAIJ_SeqBAIJ);
3081: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_seqcsrperm_C",
3082: "MatConvert_SeqAIJ_SeqCSRPERM",
3083: MatConvert_SeqAIJ_SeqCSRPERM);
3084: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_seqcrl_C",
3085: "MatConvert_SeqAIJ_SeqCRL",
3086: MatConvert_SeqAIJ_SeqCRL);
3087: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatIsTranspose_C",
3088: "MatIsTranspose_SeqAIJ",
3089: MatIsTranspose_SeqAIJ);
3090: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatIsHermitianTranspose_C",
3091: "MatIsHermitianTranspose_SeqAIJ",
3092: MatIsTranspose_SeqAIJ);
3093: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqAIJSetPreallocation_C",
3094: "MatSeqAIJSetPreallocation_SeqAIJ",
3095: MatSeqAIJSetPreallocation_SeqAIJ);
3096: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqAIJSetPreallocationCSR_C",
3097: "MatSeqAIJSetPreallocationCSR_SeqAIJ",
3098: MatSeqAIJSetPreallocationCSR_SeqAIJ);
3099: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatReorderForNonzeroDiagonal_C",
3100: "MatReorderForNonzeroDiagonal_SeqAIJ",
3101: MatReorderForNonzeroDiagonal_SeqAIJ);
3102: MatCreate_Inode(B);
3103: PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);
3104: return(0);
3105: }
3110: PetscErrorCode MatDuplicate_SeqAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
3111: {
3112: Mat C;
3113: Mat_SeqAIJ *c,*a = (Mat_SeqAIJ*)A->data;
3115: PetscInt i,m = A->rmap.n;
3118: *B = 0;
3119: MatCreate(((PetscObject)A)->comm,&C);
3120: MatSetSizes(C,A->rmap.n,A->cmap.n,A->rmap.n,A->cmap.n);
3121: MatSetType(C,((PetscObject)A)->type_name);
3122: PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));
3123:
3124: PetscMapCopy(((PetscObject)A)->comm,&A->rmap,&C->rmap);
3125: PetscMapCopy(((PetscObject)A)->comm,&A->cmap,&C->cmap);
3127: c = (Mat_SeqAIJ*)C->data;
3129: C->factor = A->factor;
3131: c->row = 0;
3132: c->col = 0;
3133: c->icol = 0;
3134: c->reallocs = 0;
3136: C->assembled = PETSC_TRUE;
3137:
3138: PetscMalloc2(m,PetscInt,&c->imax,m,PetscInt,&c->ilen);
3139: PetscLogObjectMemory(C, 2*m*sizeof(PetscInt));
3140: for (i=0; i<m; i++) {
3141: c->imax[i] = a->imax[i];
3142: c->ilen[i] = a->ilen[i];
3143: }
3145: /* allocate the matrix space */
3146: PetscMalloc3(a->i[m],PetscScalar,&c->a,a->i[m],PetscInt,&c->j,m+1,PetscInt,&c->i);
3147: PetscLogObjectMemory(C, a->i[m]*(sizeof(PetscScalar)+sizeof(PetscInt))+(m+1)*sizeof(PetscInt));
3148: c->singlemalloc = PETSC_TRUE;
3149: PetscMemcpy(c->i,a->i,(m+1)*sizeof(PetscInt));
3150: if (m > 0) {
3151: PetscMemcpy(c->j,a->j,(a->i[m])*sizeof(PetscInt));
3152: if (cpvalues == MAT_COPY_VALUES) {
3153: PetscMemcpy(c->a,a->a,(a->i[m])*sizeof(PetscScalar));
3154: } else {
3155: PetscMemzero(c->a,(a->i[m])*sizeof(PetscScalar));
3156: }
3157: }
3159: c->ignorezeroentries = a->ignorezeroentries;
3160: c->roworiented = a->roworiented;
3161: c->nonew = a->nonew;
3162: if (a->diag) {
3163: PetscMalloc((m+1)*sizeof(PetscInt),&c->diag);
3164: PetscLogObjectMemory(C,(m+1)*sizeof(PetscInt));
3165: for (i=0; i<m; i++) {
3166: c->diag[i] = a->diag[i];
3167: }
3168: } else c->diag = 0;
3169: c->solve_work = 0;
3170: c->saved_values = 0;
3171: c->idiag = 0;
3172: c->ssor_work = 0;
3173: c->keepzeroedrows = a->keepzeroedrows;
3174: c->free_a = PETSC_TRUE;
3175: c->free_ij = PETSC_TRUE;
3176: c->xtoy = 0;
3177: c->XtoY = 0;
3179: c->nz = a->nz;
3180: c->maxnz = a->maxnz;
3181: C->preallocated = PETSC_TRUE;
3183: c->compressedrow.use = a->compressedrow.use;
3184: c->compressedrow.nrows = a->compressedrow.nrows;
3185: c->compressedrow.checked = a->compressedrow.checked;
3186: if ( a->compressedrow.checked && a->compressedrow.use){
3187: i = a->compressedrow.nrows;
3188: PetscMalloc((2*i+1)*sizeof(PetscInt),&c->compressedrow.i);
3189: c->compressedrow.rindex = c->compressedrow.i + i + 1;
3190: PetscMemcpy(c->compressedrow.i,a->compressedrow.i,(i+1)*sizeof(PetscInt));
3191: PetscMemcpy(c->compressedrow.rindex,a->compressedrow.rindex,i*sizeof(PetscInt));
3192: } else {
3193: c->compressedrow.use = PETSC_FALSE;
3194: c->compressedrow.i = PETSC_NULL;
3195: c->compressedrow.rindex = PETSC_NULL;
3196: }
3197: C->same_nonzero = A->same_nonzero;
3198: MatDuplicate_Inode(A,cpvalues,&C);
3200: *B = C;
3201: PetscFListDuplicate(((PetscObject)A)->qlist,&((PetscObject)C)->qlist);
3202: return(0);
3203: }
3207: PetscErrorCode MatLoad_SeqAIJ(PetscViewer viewer, MatType type,Mat *A)
3208: {
3209: Mat_SeqAIJ *a;
3210: Mat B;
3212: PetscInt i,sum,nz,header[4],*rowlengths = 0,M,N;
3213: int fd;
3214: PetscMPIInt size;
3215: MPI_Comm comm;
3216:
3218: PetscObjectGetComm((PetscObject)viewer,&comm);
3219: MPI_Comm_size(comm,&size);
3220: if (size > 1) SETERRQ(PETSC_ERR_ARG_SIZ,"view must have one processor");
3221: PetscViewerBinaryGetDescriptor(viewer,&fd);
3222: PetscBinaryRead(fd,header,4,PETSC_INT);
3223: if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object in file");
3224: M = header[1]; N = header[2]; nz = header[3];
3226: if (nz < 0) {
3227: SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format on disk,cannot load as SeqAIJ");
3228: }
3230: /* read in row lengths */
3231: PetscMalloc(M*sizeof(PetscInt),&rowlengths);
3232: PetscBinaryRead(fd,rowlengths,M,PETSC_INT);
3234: /* check if sum of rowlengths is same as nz */
3235: for (i=0,sum=0; i< M; i++) sum +=rowlengths[i];
3236: if (sum != nz) SETERRQ2(PETSC_ERR_FILE_READ,"Inconsistant matrix data in file. no-nonzeros = %d, sum-row-lengths = %d\n",nz,sum);
3238: /* create our matrix */
3239: MatCreate(comm,&B);
3240: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,M,N);
3241: MatSetType(B,type);
3242: MatSeqAIJSetPreallocation_SeqAIJ(B,0,rowlengths);
3243: a = (Mat_SeqAIJ*)B->data;
3245: PetscBinaryRead(fd,a->j,nz,PETSC_INT);
3247: /* read in nonzero values */
3248: PetscBinaryRead(fd,a->a,nz,PETSC_SCALAR);
3250: /* set matrix "i" values */
3251: a->i[0] = 0;
3252: for (i=1; i<= M; i++) {
3253: a->i[i] = a->i[i-1] + rowlengths[i-1];
3254: a->ilen[i-1] = rowlengths[i-1];
3255: }
3256: PetscFree(rowlengths);
3258: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
3259: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
3260: *A = B;
3261: return(0);
3262: }
3266: PetscErrorCode MatEqual_SeqAIJ(Mat A,Mat B,PetscTruth* flg)
3267: {
3268: Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data,*b = (Mat_SeqAIJ *)B->data;
3272: /* If the matrix dimensions are not equal,or no of nonzeros */
3273: if ((A->rmap.n != B->rmap.n) || (A->cmap.n != B->cmap.n) ||(a->nz != b->nz)) {
3274: *flg = PETSC_FALSE;
3275: return(0);
3276: }
3277:
3278: /* if the a->i are the same */
3279: PetscMemcmp(a->i,b->i,(A->rmap.n+1)*sizeof(PetscInt),flg);
3280: if (!*flg) return(0);
3281:
3282: /* if a->j are the same */
3283: PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),flg);
3284: if (!*flg) return(0);
3285:
3286: /* if a->a are the same */
3287: PetscMemcmp(a->a,b->a,(a->nz)*sizeof(PetscScalar),flg);
3289: return(0);
3290:
3291: }
3295: /*@
3296: MatCreateSeqAIJWithArrays - Creates an sequential AIJ matrix using matrix elements (in CSR format)
3297: provided by the user.
3299: Collective on MPI_Comm
3301: Input Parameters:
3302: + comm - must be an MPI communicator of size 1
3303: . m - number of rows
3304: . n - number of columns
3305: . i - row indices
3306: . j - column indices
3307: - a - matrix values
3309: Output Parameter:
3310: . mat - the matrix
3312: Level: intermediate
3314: Notes:
3315: The i, j, and a arrays are not copied by this routine, the user must free these arrays
3316: once the matrix is destroyed
3318: You cannot set new nonzero locations into this matrix, that will generate an error.
3320: The i and j indices are 0 based
3322: The format which is used for the sparse matrix input, is equivalent to a
3323: row-major ordering.. i.e for the following matrix, the input data expected is
3324: as shown:
3326: 1 0 0
3327: 2 0 3
3328: 4 5 6
3330: i = {0,1,3,6} [size = nrow+1 = 3+1]
3331: j = {0,0,2,0,1,2} [size = nz = 6]; values must be sorted for each row
3332: v = {1,2,3,4,5,6} [size = nz = 6]
3334:
3335: .seealso: MatCreate(), MatCreateMPIAIJ(), MatCreateSeqAIJ(), MatCreateMPIAIJWithArrays(), MatMPIAIJSetPreallocationCSR()
3337: @*/
3338: PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt* i,PetscInt*j,PetscScalar *a,Mat *mat)
3339: {
3341: PetscInt ii;
3342: Mat_SeqAIJ *aij;
3343: #if defined(PETSC_USE_DEBUG)
3344: PetscInt jj;
3345: #endif
3348: if (i[0]) {
3349: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"i (row indices) must start with 0");
3350: }
3351: MatCreate(comm,mat);
3352: MatSetSizes(*mat,m,n,m,n);
3353: MatSetType(*mat,MATSEQAIJ);
3354: MatSeqAIJSetPreallocation_SeqAIJ(*mat,MAT_SKIP_ALLOCATION,0);
3355: aij = (Mat_SeqAIJ*)(*mat)->data;
3356: PetscMalloc2(m,PetscInt,&aij->imax,m,PetscInt,&aij->ilen);
3358: aij->i = i;
3359: aij->j = j;
3360: aij->a = a;
3361: aij->singlemalloc = PETSC_FALSE;
3362: aij->nonew = -1; /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/
3363: aij->free_a = PETSC_FALSE;
3364: aij->free_ij = PETSC_FALSE;
3366: for (ii=0; ii<m; ii++) {
3367: aij->ilen[ii] = aij->imax[ii] = i[ii+1] - i[ii];
3368: #if defined(PETSC_USE_DEBUG)
3369: 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]);
3370: for (jj=i[ii]+1; jj<i[ii+1]; jj++) {
3371: if (j[jj] < j[jj-1]) SETERRQ3(PETSC_ERR_ARG_OUTOFRANGE,"Column entry number %D (actual colum %D) in row %D is not sorted",jj-i[ii],j[jj],ii);
3372: if (j[jj] == j[jj]-1) SETERRQ3(PETSC_ERR_ARG_OUTOFRANGE,"Column entry number %D (actual colum %D) in row %D is identical to previous entry",jj-i[ii],j[jj],ii);
3373: }
3374: #endif
3375: }
3376: #if defined(PETSC_USE_DEBUG)
3377: for (ii=0; ii<aij->i[m]; ii++) {
3378: if (j[ii] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Negative column index at location = %d index = %d",ii,j[ii]);
3379: if (j[ii] > n - 1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column index to large at location = %d index = %d",ii,j[ii]);
3380: }
3381: #endif
3383: MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
3384: MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
3385: return(0);
3386: }
3390: PetscErrorCode MatSetColoring_SeqAIJ(Mat A,ISColoring coloring)
3391: {
3393: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
3396: if (coloring->ctype == IS_COLORING_GLOBAL) {
3397: ISColoringReference(coloring);
3398: a->coloring = coloring;
3399: } else if (coloring->ctype == IS_COLORING_GHOSTED) {
3400: PetscInt i,*larray;
3401: ISColoring ocoloring;
3402: ISColoringValue *colors;
3404: /* set coloring for diagonal portion */
3405: PetscMalloc((A->cmap.n+1)*sizeof(PetscInt),&larray);
3406: for (i=0; i<A->cmap.n; i++) {
3407: larray[i] = i;
3408: }
3409: ISGlobalToLocalMappingApply(A->mapping,IS_GTOLM_MASK,A->cmap.n,larray,PETSC_NULL,larray);
3410: PetscMalloc((A->cmap.n+1)*sizeof(ISColoringValue),&colors);
3411: for (i=0; i<A->cmap.n; i++) {
3412: colors[i] = coloring->colors[larray[i]];
3413: }
3414: PetscFree(larray);
3415: ISColoringCreate(PETSC_COMM_SELF,coloring->n,A->cmap.n,colors,&ocoloring);
3416: a->coloring = ocoloring;
3417: }
3418: return(0);
3419: }
3421: #if defined(PETSC_HAVE_ADIC)
3423: #include "adic/ad_utils.h"
3428: PetscErrorCode MatSetValuesAdic_SeqAIJ(Mat A,void *advalues)
3429: {
3430: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
3431: PetscInt m = A->rmap.n,*ii = a->i,*jj = a->j,nz,i,j,nlen;
3432: PetscScalar *v = a->a,*values = ((PetscScalar*)advalues)+1;
3433: ISColoringValue *color;
3436: if (!a->coloring) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Coloring not set for matrix");
3437: nlen = PetscADGetDerivTypeSize()/sizeof(PetscScalar);
3438: color = a->coloring->colors;
3439: /* loop over rows */
3440: for (i=0; i<m; i++) {
3441: nz = ii[i+1] - ii[i];
3442: /* loop over columns putting computed value into matrix */
3443: for (j=0; j<nz; j++) {
3444: *v++ = values[color[*jj++]];
3445: }
3446: values += nlen; /* jump to next row of derivatives */
3447: }
3448: return(0);
3449: }
3450: #endif
3454: PetscErrorCode MatSetValuesAdifor_SeqAIJ(Mat A,PetscInt nl,void *advalues)
3455: {
3456: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
3457: PetscInt m = A->rmap.n,*ii = a->i,*jj = a->j,nz,i,j;
3458: PetscScalar *v = a->a,*values = (PetscScalar *)advalues;
3459: ISColoringValue *color;
3462: if (!a->coloring) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Coloring not set for matrix");
3463: color = a->coloring->colors;
3464: /* loop over rows */
3465: for (i=0; i<m; i++) {
3466: nz = ii[i+1] - ii[i];
3467: /* loop over columns putting computed value into matrix */
3468: for (j=0; j<nz; j++) {
3469: *v++ = values[color[*jj++]];
3470: }
3471: values += nl; /* jump to next row of derivatives */
3472: }
3473: return(0);
3474: }
3476: /*
3477: Special version for direct calls from Fortran
3478: */
3479: #if defined(PETSC_HAVE_FORTRAN_CAPS)
3480: #define matsetvaluesseqaij_ MATSETVALUESSEQAIJ
3481: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
3482: #define matsetvaluesseqaij_ matsetvaluesseqaij
3483: #endif
3485: /* Change these macros so can be used in void function */
3486: #undef CHKERRQ
3487: #define CHKERRQ(ierr) CHKERRABORT(((PetscObject)A)->comm,ierr)
3488: #undef SETERRQ2
3489: #define SETERRQ2(ierr,b,c,d) CHKERRABORT(((PetscObject)A)->comm,ierr)
3494: void PETSC_STDCALL matsetvaluesseqaij_(Mat *AA,PetscInt *mm,const PetscInt im[],PetscInt *nn,const PetscInt in[],const PetscScalar v[],InsertMode *isis, PetscErrorCode *_ierr)
3495: {
3496: Mat A = *AA;
3497: PetscInt m = *mm, n = *nn;
3498: InsertMode is = *isis;
3499: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
3500: PetscInt *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N;
3501: PetscInt *imax,*ai,*ailen;
3503: PetscInt *aj,nonew = a->nonew,lastcol = -1;
3504: PetscScalar *ap,value,*aa;
3505: PetscTruth ignorezeroentries = a->ignorezeroentries;
3506: PetscTruth roworiented = a->roworiented;
3509: MatPreallocated(A);
3510: imax = a->imax;
3511: ai = a->i;
3512: ailen = a->ilen;
3513: aj = a->j;
3514: aa = a->a;
3516: for (k=0; k<m; k++) { /* loop over added rows */
3517: row = im[k];
3518: if (row < 0) continue;
3519: #if defined(PETSC_USE_DEBUG)
3520: if (row >= A->rmap.n) SETERRABORT(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
3521: #endif
3522: rp = aj + ai[row]; ap = aa + ai[row];
3523: rmax = imax[row]; nrow = ailen[row];
3524: low = 0;
3525: high = nrow;
3526: for (l=0; l<n; l++) { /* loop over added columns */
3527: if (in[l] < 0) continue;
3528: #if defined(PETSC_USE_DEBUG)
3529: if (in[l] >= A->cmap.n) SETERRABORT(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
3530: #endif
3531: col = in[l];
3532: if (roworiented) {
3533: value = v[l + k*n];
3534: } else {
3535: value = v[k + l*m];
3536: }
3537: if (value == 0.0 && ignorezeroentries && (is == ADD_VALUES)) continue;
3539: if (col <= lastcol) low = 0; else high = nrow;
3540: lastcol = col;
3541: while (high-low > 5) {
3542: t = (low+high)/2;
3543: if (rp[t] > col) high = t;
3544: else low = t;
3545: }
3546: for (i=low; i<high; i++) {
3547: if (rp[i] > col) break;
3548: if (rp[i] == col) {
3549: if (is == ADD_VALUES) ap[i] += value;
3550: else ap[i] = value;
3551: goto noinsert;
3552: }
3553: }
3554: if (value == 0.0 && ignorezeroentries) goto noinsert;
3555: if (nonew == 1) goto noinsert;
3556: if (nonew == -1) SETERRABORT(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
3557: MatSeqXAIJReallocateAIJ(A,A->rmap.n,1,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar);
3558: N = nrow++ - 1; a->nz++; high++;
3559: /* shift up all the later entries in this row */
3560: for (ii=N; ii>=i; ii--) {
3561: rp[ii+1] = rp[ii];
3562: ap[ii+1] = ap[ii];
3563: }
3564: rp[i] = col;
3565: ap[i] = value;
3566: noinsert:;
3567: low = i + 1;
3568: }
3569: ailen[row] = nrow;
3570: }
3571: A->same_nonzero = PETSC_FALSE;
3572: PetscFunctionReturnVoid();
3573: }