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 MatGetRowIJ_SeqAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *m,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
16: {
17: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
19: PetscInt i,ishift;
20:
22: *m = A->m;
23: if (!ia) return(0);
24: ishift = 0;
25: if (symmetric && !A->structurally_symmetric) {
26: MatToSymmetricIJ_SeqAIJ(A->m,a->i,a->j,ishift,oshift,ia,ja);
27: } else if (oshift == 1) {
28: PetscInt nz = a->i[A->m];
29: /* malloc space and add 1 to i and j indices */
30: PetscMalloc((A->m+1)*sizeof(PetscInt),ia);
31: PetscMalloc((nz+1)*sizeof(PetscInt),ja);
32: for (i=0; i<nz; i++) (*ja)[i] = a->j[i] + 1;
33: for (i=0; i<A->m+1; i++) (*ia)[i] = a->i[i] + 1;
34: } else {
35: *ia = a->i; *ja = a->j;
36: }
37: return(0);
38: }
42: PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *n,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
43: {
45:
47: if (!ia) return(0);
48: if ((symmetric && !A->structurally_symmetric) || oshift == 1) {
49: PetscFree(*ia);
50: PetscFree(*ja);
51: }
52: return(0);
53: }
57: PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
58: {
59: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
61: PetscInt i,*collengths,*cia,*cja,n = A->n,m = A->m;
62: PetscInt nz = a->i[m],row,*jj,mr,col;
63:
65: *nn = A->n;
66: if (!ia) return(0);
67: if (symmetric) {
68: MatToSymmetricIJ_SeqAIJ(A->m,a->i,a->j,0,oshift,ia,ja);
69: } else {
70: PetscMalloc((n+1)*sizeof(PetscInt),&collengths);
71: PetscMemzero(collengths,n*sizeof(PetscInt));
72: PetscMalloc((n+1)*sizeof(PetscInt),&cia);
73: PetscMalloc((nz+1)*sizeof(PetscInt),&cja);
74: jj = a->j;
75: for (i=0; i<nz; i++) {
76: collengths[jj[i]]++;
77: }
78: cia[0] = oshift;
79: for (i=0; i<n; i++) {
80: cia[i+1] = cia[i] + collengths[i];
81: }
82: PetscMemzero(collengths,n*sizeof(PetscInt));
83: jj = a->j;
84: for (row=0; row<m; row++) {
85: mr = a->i[row+1] - a->i[row];
86: for (i=0; i<mr; i++) {
87: col = *jj++;
88: cja[cia[col] + collengths[col]++ - oshift] = row + oshift;
89: }
90: }
91: PetscFree(collengths);
92: *ia = cia; *ja = cja;
93: }
94: return(0);
95: }
99: PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *n,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
100: {
104: if (!ia) return(0);
106: PetscFree(*ia);
107: PetscFree(*ja);
108:
109: return(0);
110: }
112: #define CHUNKSIZE 15
116: PetscErrorCode MatSetValues_SeqAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
117: {
118: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
119: PetscInt *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N;
120: PetscInt *imax = a->imax,*ai = a->i,*ailen = a->ilen;
122: PetscInt *aj = a->j,nonew = a->nonew,lastcol = -1;
123: PetscScalar *ap,value,*aa = a->a;
124: PetscTruth ignorezeroentries = ((a->ignorezeroentries && is == ADD_VALUES) ? PETSC_TRUE:PETSC_FALSE);
125: PetscTruth roworiented = a->roworiented;
128: for (k=0; k<m; k++) { /* loop over added rows */
129: row = im[k];
130: if (row < 0) continue;
131: #if defined(PETSC_USE_DEBUG)
132: if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->m-1);
133: #endif
134: rp = aj + ai[row]; ap = aa + ai[row];
135: rmax = imax[row]; nrow = ailen[row];
136: low = 0;
137: high = nrow;
138: for (l=0; l<n; l++) { /* loop over added columns */
139: if (in[l] < 0) continue;
140: #if defined(PETSC_USE_DEBUG)
141: if (in[l] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->n-1);
142: #endif
143: col = in[l];
144: if (roworiented) {
145: value = v[l + k*n];
146: } else {
147: value = v[k + l*m];
148: }
149: if (value == 0.0 && ignorezeroentries && (is == ADD_VALUES)) continue;
151: if (col <= lastcol) low = 0; else high = nrow;
152: lastcol = col;
153: while (high-low > 5) {
154: t = (low+high)/2;
155: if (rp[t] > col) high = t;
156: else low = t;
157: }
158: for (i=low; i<high; i++) {
159: if (rp[i] > col) break;
160: if (rp[i] == col) {
161: if (is == ADD_VALUES) ap[i] += value;
162: else ap[i] = value;
163: goto noinsert;
164: }
165: }
166: if (value == 0.0 && ignorezeroentries) goto noinsert;
167: if (nonew == 1) goto noinsert;
168: if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero at (%D,%D) in the matrix",row,col);
169: MatSeqXAIJReallocateAIJ(a,1,nrow,row,col,rmax,aa,ai,aj,A->m,rp,ap,imax,nonew);
170: N = nrow++ - 1; a->nz++;
171: /* shift up all the later entries in this row */
172: for (ii=N; ii>=i; ii--) {
173: rp[ii+1] = rp[ii];
174: ap[ii+1] = ap[ii];
175: }
176: rp[i] = col;
177: ap[i] = value;
178: noinsert:;
179: low = i + 1;
180: }
181: ailen[row] = nrow;
182: }
183: A->same_nonzero = PETSC_FALSE;
184: return(0);
185: }
189: PetscErrorCode MatGetValues_SeqAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],PetscScalar v[])
190: {
191: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
192: PetscInt *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
193: PetscInt *ai = a->i,*ailen = a->ilen;
194: PetscScalar *ap,*aa = a->a,zero = 0.0;
197: for (k=0; k<m; k++) { /* loop over rows */
198: row = im[k];
199: if (row < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",row);
200: if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->m-1);
201: rp = aj + ai[row]; ap = aa + ai[row];
202: nrow = ailen[row];
203: for (l=0; l<n; l++) { /* loop over columns */
204: if (in[l] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %D",in[l]);
205: if (in[l] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->n-1);
206: col = in[l] ;
207: high = nrow; low = 0; /* assume unsorted */
208: while (high-low > 5) {
209: t = (low+high)/2;
210: if (rp[t] > col) high = t;
211: else low = t;
212: }
213: for (i=low; i<high; i++) {
214: if (rp[i] > col) break;
215: if (rp[i] == col) {
216: *v++ = ap[i];
217: goto finished;
218: }
219: }
220: *v++ = zero;
221: finished:;
222: }
223: }
224: return(0);
225: }
230: PetscErrorCode MatView_SeqAIJ_Binary(Mat A,PetscViewer viewer)
231: {
232: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
234: PetscInt i,*col_lens;
235: int fd;
238: PetscViewerBinaryGetDescriptor(viewer,&fd);
239: PetscMalloc((4+A->m)*sizeof(PetscInt),&col_lens);
240: col_lens[0] = MAT_FILE_COOKIE;
241: col_lens[1] = A->m;
242: col_lens[2] = A->n;
243: col_lens[3] = a->nz;
245: /* store lengths of each row and write (including header) to file */
246: for (i=0; i<A->m; i++) {
247: col_lens[4+i] = a->i[i+1] - a->i[i];
248: }
249: PetscBinaryWrite(fd,col_lens,4+A->m,PETSC_INT,PETSC_TRUE);
250: PetscFree(col_lens);
252: /* store column indices (zero start index) */
253: PetscBinaryWrite(fd,a->j,a->nz,PETSC_INT,PETSC_FALSE);
255: /* store nonzero values */
256: PetscBinaryWrite(fd,a->a,a->nz,PETSC_SCALAR,PETSC_FALSE);
257: return(0);
258: }
260: EXTERN PetscErrorCode MatSeqAIJFactorInfo_Matlab(Mat,PetscViewer);
264: PetscErrorCode MatView_SeqAIJ_ASCII(Mat A,PetscViewer viewer)
265: {
266: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
267: PetscErrorCode ierr;
268: PetscInt i,j,m = A->m,shift=0;
269: const char *name;
270: PetscViewerFormat format;
273: PetscObjectGetName((PetscObject)A,&name);
274: PetscViewerGetFormat(viewer,&format);
275: if (format == PETSC_VIEWER_ASCII_MATLAB) {
276: PetscInt nofinalvalue = 0;
277: if ((a->i[m] == a->i[m-1]) || (a->j[a->nz-1] != A->n-!shift)) {
278: nofinalvalue = 1;
279: }
280: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
281: PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",m,A->n);
282: PetscViewerASCIIPrintf(viewer,"%% Nonzeros = %D \n",a->nz);
283: PetscViewerASCIIPrintf(viewer,"zzz = zeros(%D,3);\n",a->nz+nofinalvalue);
284: PetscViewerASCIIPrintf(viewer,"zzz = [\n");
286: for (i=0; i<m; i++) {
287: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
288: #if defined(PETSC_USE_COMPLEX)
289: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e + %18.16ei \n",i+1,a->j[j]+!shift,PetscRealPart(a->a[j]),PetscImaginaryPart(a->a[j]));
290: #else
291: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",i+1,a->j[j]+!shift,a->a[j]);
292: #endif
293: }
294: }
295: if (nofinalvalue) {
296: PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",m,A->n,0.0);
297: }
298: PetscViewerASCIIPrintf(viewer,"];\n %s = spconvert(zzz);\n",name);
299: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
300: } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) {
301: return(0);
302: } else if (format == PETSC_VIEWER_ASCII_COMMON) {
303: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
304: for (i=0; i<m; i++) {
305: PetscViewerASCIIPrintf(viewer,"row %D:",i);
306: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
307: #if defined(PETSC_USE_COMPLEX)
308: if (PetscImaginaryPart(a->a[j]) > 0.0 && PetscRealPart(a->a[j]) != 0.0) {
309: PetscViewerASCIIPrintf(viewer," (%D, %g + %g i)",a->j[j]+shift,PetscRealPart(a->a[j]),PetscImaginaryPart(a->a[j]));
310: } else if (PetscImaginaryPart(a->a[j]) < 0.0 && PetscRealPart(a->a[j]) != 0.0) {
311: PetscViewerASCIIPrintf(viewer," (%D, %g - %g i)",a->j[j]+shift,PetscRealPart(a->a[j]),-PetscImaginaryPart(a->a[j]));
312: } else if (PetscRealPart(a->a[j]) != 0.0) {
313: PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[j]+shift,PetscRealPart(a->a[j]));
314: }
315: #else
316: if (a->a[j] != 0.0) {PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[j]+shift,a->a[j]);}
317: #endif
318: }
319: PetscViewerASCIIPrintf(viewer,"\n");
320: }
321: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
322: } else if (format == PETSC_VIEWER_ASCII_SYMMODU) {
323: PetscInt nzd=0,fshift=1,*sptr;
324: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
325: PetscMalloc((m+1)*sizeof(PetscInt),&sptr);
326: for (i=0; i<m; i++) {
327: sptr[i] = nzd+1;
328: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
329: if (a->j[j] >= i) {
330: #if defined(PETSC_USE_COMPLEX)
331: if (PetscImaginaryPart(a->a[j]) != 0.0 || PetscRealPart(a->a[j]) != 0.0) nzd++;
332: #else
333: if (a->a[j] != 0.0) nzd++;
334: #endif
335: }
336: }
337: }
338: sptr[m] = nzd+1;
339: PetscViewerASCIIPrintf(viewer," %D %D\n\n",m,nzd);
340: for (i=0; i<m+1; i+=6) {
341: 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]);}
342: 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]);}
343: else if (i+2<m) {PetscViewerASCIIPrintf(viewer," %D %D %D %D\n",sptr[i],sptr[i+1],sptr[i+2],sptr[i+3]);}
344: else if (i+1<m) {PetscViewerASCIIPrintf(viewer," %D %D %D\n",sptr[i],sptr[i+1],sptr[i+2]);}
345: else if (i<m) {PetscViewerASCIIPrintf(viewer," %D %D\n",sptr[i],sptr[i+1]);}
346: else {PetscViewerASCIIPrintf(viewer," %D\n",sptr[i]);}
347: }
348: PetscViewerASCIIPrintf(viewer,"\n");
349: PetscFree(sptr);
350: for (i=0; i<m; i++) {
351: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
352: if (a->j[j] >= i) {PetscViewerASCIIPrintf(viewer," %D ",a->j[j]+fshift);}
353: }
354: PetscViewerASCIIPrintf(viewer,"\n");
355: }
356: PetscViewerASCIIPrintf(viewer,"\n");
357: for (i=0; i<m; i++) {
358: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
359: if (a->j[j] >= i) {
360: #if defined(PETSC_USE_COMPLEX)
361: if (PetscImaginaryPart(a->a[j]) != 0.0 || PetscRealPart(a->a[j]) != 0.0) {
362: PetscViewerASCIIPrintf(viewer," %18.16e %18.16e ",PetscRealPart(a->a[j]),PetscImaginaryPart(a->a[j]));
363: }
364: #else
365: if (a->a[j] != 0.0) {PetscViewerASCIIPrintf(viewer," %18.16e ",a->a[j]);}
366: #endif
367: }
368: }
369: PetscViewerASCIIPrintf(viewer,"\n");
370: }
371: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
372: } else if (format == PETSC_VIEWER_ASCII_DENSE) {
373: PetscInt cnt = 0,jcnt;
374: PetscScalar value;
376: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
377: for (i=0; i<m; i++) {
378: jcnt = 0;
379: for (j=0; j<A->n; j++) {
380: if (jcnt < a->i[i+1]-a->i[i] && j == a->j[cnt]) {
381: value = a->a[cnt++];
382: jcnt++;
383: } else {
384: value = 0.0;
385: }
386: #if defined(PETSC_USE_COMPLEX)
387: PetscViewerASCIIPrintf(viewer," %7.5e+%7.5e i ",PetscRealPart(value),PetscImaginaryPart(value));
388: #else
389: PetscViewerASCIIPrintf(viewer," %7.5e ",value);
390: #endif
391: }
392: PetscViewerASCIIPrintf(viewer,"\n");
393: }
394: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
395: } else {
396: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
397: for (i=0; i<m; i++) {
398: PetscViewerASCIIPrintf(viewer,"row %D:",i);
399: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
400: #if defined(PETSC_USE_COMPLEX)
401: if (PetscImaginaryPart(a->a[j]) > 0.0) {
402: PetscViewerASCIIPrintf(viewer," (%D, %g + %g i)",a->j[j]+shift,PetscRealPart(a->a[j]),PetscImaginaryPart(a->a[j]));
403: } else if (PetscImaginaryPart(a->a[j]) < 0.0) {
404: PetscViewerASCIIPrintf(viewer," (%D, %g - %g i)",a->j[j]+shift,PetscRealPart(a->a[j]),-PetscImaginaryPart(a->a[j]));
405: } else {
406: PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[j]+shift,PetscRealPart(a->a[j]));
407: }
408: #else
409: PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[j]+shift,a->a[j]);
410: #endif
411: }
412: PetscViewerASCIIPrintf(viewer,"\n");
413: }
414: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
415: }
416: PetscViewerFlush(viewer);
417: return(0);
418: }
422: PetscErrorCode MatView_SeqAIJ_Draw_Zoom(PetscDraw draw,void *Aa)
423: {
424: Mat A = (Mat) Aa;
425: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
426: PetscErrorCode ierr;
427: PetscInt i,j,m = A->m,color;
428: PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r,maxv = 0.0;
429: PetscViewer viewer;
430: PetscViewerFormat format;
433: PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);
434: PetscViewerGetFormat(viewer,&format);
436: PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);
437: /* loop over matrix elements drawing boxes */
439: if (format != PETSC_VIEWER_DRAW_CONTOUR) {
440: /* Blue for negative, Cyan for zero and Red for positive */
441: color = PETSC_DRAW_BLUE;
442: for (i=0; i<m; i++) {
443: y_l = m - i - 1.0; y_r = y_l + 1.0;
444: for (j=a->i[i]; j<a->i[i+1]; j++) {
445: x_l = a->j[j] ; x_r = x_l + 1.0;
446: #if defined(PETSC_USE_COMPLEX)
447: if (PetscRealPart(a->a[j]) >= 0.) continue;
448: #else
449: if (a->a[j] >= 0.) continue;
450: #endif
451: PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
452: }
453: }
454: color = PETSC_DRAW_CYAN;
455: for (i=0; i<m; i++) {
456: y_l = m - i - 1.0; y_r = y_l + 1.0;
457: for (j=a->i[i]; j<a->i[i+1]; j++) {
458: x_l = a->j[j]; x_r = x_l + 1.0;
459: if (a->a[j] != 0.) continue;
460: PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
461: }
462: }
463: color = PETSC_DRAW_RED;
464: for (i=0; i<m; i++) {
465: y_l = m - i - 1.0; y_r = y_l + 1.0;
466: for (j=a->i[i]; j<a->i[i+1]; j++) {
467: x_l = a->j[j]; x_r = x_l + 1.0;
468: #if defined(PETSC_USE_COMPLEX)
469: if (PetscRealPart(a->a[j]) <= 0.) continue;
470: #else
471: if (a->a[j] <= 0.) continue;
472: #endif
473: PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
474: }
475: }
476: } else {
477: /* use contour shading to indicate magnitude of values */
478: /* first determine max of all nonzero values */
479: PetscInt nz = a->nz,count;
480: PetscDraw popup;
481: PetscReal scale;
483: for (i=0; i<nz; i++) {
484: if (PetscAbsScalar(a->a[i]) > maxv) maxv = PetscAbsScalar(a->a[i]);
485: }
486: scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv;
487: PetscDrawGetPopup(draw,&popup);
488: if (popup) {PetscDrawScalePopup(popup,0.0,maxv);}
489: count = 0;
490: for (i=0; i<m; i++) {
491: y_l = m - i - 1.0; y_r = y_l + 1.0;
492: for (j=a->i[i]; j<a->i[i+1]; j++) {
493: x_l = a->j[j]; x_r = x_l + 1.0;
494: color = PETSC_DRAW_BASIC_COLORS + (PetscInt)(scale*PetscAbsScalar(a->a[count]));
495: PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
496: count++;
497: }
498: }
499: }
500: return(0);
501: }
505: PetscErrorCode MatView_SeqAIJ_Draw(Mat A,PetscViewer viewer)
506: {
508: PetscDraw draw;
509: PetscReal xr,yr,xl,yl,h,w;
510: PetscTruth isnull;
513: PetscViewerDrawGetDraw(viewer,0,&draw);
514: PetscDrawIsNull(draw,&isnull);
515: if (isnull) return(0);
517: PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);
518: xr = A->n; yr = A->m; h = yr/10.0; w = xr/10.0;
519: xr += w; yr += h; xl = -w; yl = -h;
520: PetscDrawSetCoordinates(draw,xl,yl,xr,yr);
521: PetscDrawZoom(draw,MatView_SeqAIJ_Draw_Zoom,A);
522: PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);
523: return(0);
524: }
528: PetscErrorCode MatView_SeqAIJ(Mat A,PetscViewer viewer)
529: {
531: PetscTruth issocket,iascii,isbinary,isdraw;
534: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);
535: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
536: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
537: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
538: if (iascii) {
539: MatView_SeqAIJ_ASCII(A,viewer);
540: #if defined(PETSC_USE_SOCKET_VIEWER)
541: } else if (issocket) {
542: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
543: PetscViewerSocketPutSparse_Private(viewer,A->m,A->n,a->nz,a->a,a->i,a->j);
544: #endif
545: } else if (isbinary) {
546: MatView_SeqAIJ_Binary(A,viewer);
547: } else if (isdraw) {
548: MatView_SeqAIJ_Draw(A,viewer);
549: } else {
550: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by SeqAIJ matrices",((PetscObject)viewer)->type_name);
551: }
552: MatView_Inode(A,viewer);
553: return(0);
554: }
558: PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode)
559: {
560: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
562: PetscInt fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
563: PetscInt m = A->m,*ip,N,*ailen = a->ilen,rmax = 0;
564: PetscScalar *aa = a->a,*ap;
565: PetscReal ratio=0.6;
568: if (mode == MAT_FLUSH_ASSEMBLY) return(0);
570: if (m) rmax = ailen[0]; /* determine row with most nonzeros */
571: for (i=1; i<m; i++) {
572: /* move each row back by the amount of empty slots (fshift) before it*/
573: fshift += imax[i-1] - ailen[i-1];
574: rmax = PetscMax(rmax,ailen[i]);
575: if (fshift) {
576: ip = aj + ai[i] ;
577: ap = aa + ai[i] ;
578: N = ailen[i];
579: for (j=0; j<N; j++) {
580: ip[j-fshift] = ip[j];
581: ap[j-fshift] = ap[j];
582: }
583: }
584: ai[i] = ai[i-1] + ailen[i-1];
585: }
586: if (m) {
587: fshift += imax[m-1] - ailen[m-1];
588: ai[m] = ai[m-1] + ailen[m-1];
589: }
590: /* reset ilen and imax for each row */
591: for (i=0; i<m; i++) {
592: ailen[i] = imax[i] = ai[i+1] - ai[i];
593: }
594: a->nz = ai[m];
596: /* diagonals may have moved, so kill the diagonal pointers */
597: if (fshift && a->diag) {
598: PetscFree(a->diag);
599: PetscLogObjectMemory(A,-(m+1)*sizeof(PetscInt));
600: a->diag = 0;
601: }
602: PetscLogInfo((A,"MatAssemblyEnd_SeqAIJ:Matrix size: %D X %D; storage space: %D unneeded,%D used\n",m,A->n,fshift,a->nz));
603: PetscLogInfo((A,"MatAssemblyEnd_SeqAIJ:Number of mallocs during MatSetValues() is %D\n",a->reallocs));
604: PetscLogInfo((A,"MatAssemblyEnd_SeqAIJ:Maximum nonzeros in any row is %D\n",rmax));
605: a->reallocs = 0;
606: A->info.nz_unneeded = (double)fshift;
607: a->rmax = rmax;
609: /* check for zero rows. If found a large number of zero rows, use CompressedRow functions */
610: Mat_CheckCompressedRow(A,&a->compressedrow,a->i,m,ratio);
611: A->same_nonzero = PETSC_TRUE;
613: MatAssemblyEnd_Inode(A,mode);
614: return(0);
615: }
619: PetscErrorCode MatZeroEntries_SeqAIJ(Mat A)
620: {
621: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
625: PetscMemzero(a->a,(a->i[A->m])*sizeof(PetscScalar));
626: return(0);
627: }
631: PetscErrorCode MatDestroy_SeqAIJ(Mat A)
632: {
633: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
637: #if defined(PETSC_USE_LOG)
638: PetscLogObjectState((PetscObject)A,"Rows=%D, Cols=%D, NZ=%D",A->m,A->n,a->nz);
639: #endif
640: if (a->freedata){
641: MatSeqXAIJFreeAIJ(a->singlemalloc,&a->a,&a->j,&a->i);
642: }
643: if (a->row) {
644: ISDestroy(a->row);
645: }
646: if (a->col) {
647: ISDestroy(a->col);
648: }
649: if (a->diag) {PetscFree(a->diag);}
650: if (a->ilen) {PetscFree2(a->imax,a->ilen);}
651: if (a->idiag) {PetscFree(a->idiag);}
652: if (a->solve_work) {PetscFree(a->solve_work);}
653: if (a->icol) {ISDestroy(a->icol);}
654: if (a->saved_values) {PetscFree(a->saved_values);}
655: if (a->coloring) {ISColoringDestroy(a->coloring);}
656: if (a->xtoy) {PetscFree(a->xtoy);}
657: if (a->compressedrow.use){PetscFree(a->compressedrow.i);}
659: MatDestroy_Inode(A);
661: PetscFree(a);
663: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSeqAIJSetColumnIndices_C","",PETSC_NULL);
664: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatStoreValues_C","",PETSC_NULL);
665: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatRetrieveValues_C","",PETSC_NULL);
666: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatConvert_seqaij_seqsbaij_C","",PETSC_NULL);
667: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatConvert_seqaij_seqbaij_C","",PETSC_NULL);
668: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatIsTranspose_C","",PETSC_NULL);
669: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSeqAIJSetPreallocation_C","",PETSC_NULL);
670: PetscObjectComposeFunctionDynamic((PetscObject)A,"MatReorderForNonzeroDiagonal_C","",PETSC_NULL);
671: return(0);
672: }
676: PetscErrorCode MatCompress_SeqAIJ(Mat A)
677: {
679: return(0);
680: }
684: PetscErrorCode MatSetOption_SeqAIJ(Mat A,MatOption op)
685: {
686: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
690: switch (op) {
691: case MAT_ROW_ORIENTED:
692: a->roworiented = PETSC_TRUE;
693: break;
694: case MAT_KEEP_ZEROED_ROWS:
695: a->keepzeroedrows = PETSC_TRUE;
696: break;
697: case MAT_COLUMN_ORIENTED:
698: a->roworiented = PETSC_FALSE;
699: break;
700: case MAT_COLUMNS_SORTED:
701: a->sorted = PETSC_TRUE;
702: break;
703: case MAT_COLUMNS_UNSORTED:
704: a->sorted = PETSC_FALSE;
705: break;
706: case MAT_NO_NEW_NONZERO_LOCATIONS:
707: a->nonew = 1;
708: break;
709: case MAT_NEW_NONZERO_LOCATION_ERR:
710: a->nonew = -1;
711: break;
712: case MAT_NEW_NONZERO_ALLOCATION_ERR:
713: a->nonew = -2;
714: break;
715: case MAT_YES_NEW_NONZERO_LOCATIONS:
716: a->nonew = 0;
717: break;
718: case MAT_IGNORE_ZERO_ENTRIES:
719: a->ignorezeroentries = PETSC_TRUE;
720: break;
721: case MAT_USE_COMPRESSEDROW:
722: a->compressedrow.use = PETSC_TRUE;
723: break;
724: case MAT_DO_NOT_USE_COMPRESSEDROW:
725: a->compressedrow.use = PETSC_FALSE;
726: break;
727: case MAT_ROWS_SORTED:
728: case MAT_ROWS_UNSORTED:
729: case MAT_YES_NEW_DIAGONALS:
730: case MAT_IGNORE_OFF_PROC_ENTRIES:
731: case MAT_USE_HASH_TABLE:
732: PetscLogInfo((A,"MatSetOption_SeqAIJ:Option ignored\n"));
733: break;
734: case MAT_NO_NEW_DIAGONALS:
735: SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
736: default:
737: break;
738: }
739: MatSetOption_Inode(A,op);
740: return(0);
741: }
745: PetscErrorCode MatGetDiagonal_SeqAIJ(Mat A,Vec v)
746: {
747: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
749: PetscInt i,j,n;
750: PetscScalar *x,zero = 0.0;
753: VecSet(v,zero);
754: VecGetArray(v,&x);
755: VecGetLocalSize(v,&n);
756: if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
757: for (i=0; i<A->m; i++) {
758: for (j=a->i[i]; j<a->i[i+1]; j++) {
759: if (a->j[j] == i) {
760: x[i] = a->a[j];
761: break;
762: }
763: }
764: }
765: VecRestoreArray(v,&x);
766: return(0);
767: }
771: PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
772: {
773: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
774: PetscScalar *x,*y;
775: PetscErrorCode ierr;
776: PetscInt m = A->m;
777: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
778: PetscScalar *v,alpha;
779: PetscInt n,i,*idx,*ii,*ridx=PETSC_NULL;
780: Mat_CompressedRow cprow = a->compressedrow;
781: PetscTruth usecprow = cprow.use;
782: #endif
785: if (zz != yy) {VecCopy(zz,yy);}
786: VecGetArray(xx,&x);
787: VecGetArray(yy,&y);
789: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
790: fortranmulttransposeaddaij_(&m,x,a->i,a->j,a->a,y);
791: #else
792: if (usecprow){
793: m = cprow.nrows;
794: ii = cprow.i;
795: ridx = cprow.rindex;
796: } else {
797: ii = a->i;
798: }
799: for (i=0; i<m; i++) {
800: idx = a->j + ii[i] ;
801: v = a->a + ii[i] ;
802: n = ii[i+1] - ii[i];
803: if (usecprow){
804: alpha = x[ridx[i]];
805: } else {
806: alpha = x[i];
807: }
808: while (n-->0) {y[*idx++] += alpha * *v++;}
809: }
810: #endif
811: PetscLogFlops(2*a->nz);
812: VecRestoreArray(xx,&x);
813: VecRestoreArray(yy,&y);
814: return(0);
815: }
819: PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec xx,Vec yy)
820: {
821: PetscScalar zero = 0.0;
825: VecSet(yy,zero);
826: MatMultTransposeAdd_SeqAIJ(A,xx,yy,yy);
827: return(0);
828: }
833: PetscErrorCode MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
834: {
835: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
836: PetscScalar *x,*y,*aa;
838: PetscInt m=A->m,*aj,*ii;
839: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
840: PetscInt n,i,jrow,j,*ridx=PETSC_NULL;
841: PetscScalar sum;
842: PetscTruth usecprow=a->compressedrow.use;
843: #endif
845: #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
846: #pragma disjoint(*x,*y,*aa)
847: #endif
850: VecGetArray(xx,&x);
851: VecGetArray(yy,&y);
852: aj = a->j;
853: aa = a->a;
854: ii = a->i;
855: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
856: fortranmultaij_(&m,x,ii,aj,aa,y);
857: #else
858: if (usecprow){ /* use compressed row format */
859: m = a->compressedrow.nrows;
860: ii = a->compressedrow.i;
861: ridx = a->compressedrow.rindex;
862: for (i=0; i<m; i++){
863: n = ii[i+1] - ii[i];
864: aj = a->j + ii[i];
865: aa = a->a + ii[i];
866: sum = 0.0;
867: for (j=0; j<n; j++) sum += (*aa++)*x[*aj++];
868: y[*ridx++] = sum;
869: }
870: } else { /* do not use compressed row format */
871: for (i=0; i<m; i++) {
872: jrow = ii[i];
873: n = ii[i+1] - jrow;
874: sum = 0.0;
875: for (j=0; j<n; j++) {
876: sum += aa[jrow]*x[aj[jrow]]; jrow++;
877: }
878: y[i] = sum;
879: }
880: }
881: #endif
882: PetscLogFlops(2*a->nz - m);
883: VecRestoreArray(xx,&x);
884: VecRestoreArray(yy,&y);
885: return(0);
886: }
890: PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
891: {
892: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
893: PetscScalar *x,*y,*z,*aa;
895: PetscInt m = A->m,*aj,*ii;
896: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
897: PetscInt n,i,jrow,j,*ridx=PETSC_NULL;
898: PetscScalar sum;
899: PetscTruth usecprow=a->compressedrow.use;
900: #endif
903: VecGetArray(xx,&x);
904: VecGetArray(yy,&y);
905: if (zz != yy) {
906: VecGetArray(zz,&z);
907: } else {
908: z = y;
909: }
911: aj = a->j;
912: aa = a->a;
913: ii = a->i;
914: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
915: fortranmultaddaij_(&m,x,ii,aj,aa,y,z);
916: #else
917: if (usecprow){ /* use compressed row format */
918: if (zz != yy){
919: PetscMemcpy(z,y,m*sizeof(PetscScalar));
920: }
921: m = a->compressedrow.nrows;
922: ii = a->compressedrow.i;
923: ridx = a->compressedrow.rindex;
924: for (i=0; i<m; i++){
925: n = ii[i+1] - ii[i];
926: aj = a->j + ii[i];
927: aa = a->a + ii[i];
928: sum = y[*ridx];
929: for (j=0; j<n; j++) sum += (*aa++)*x[*aj++];
930: z[*ridx++] = sum;
931: }
932: } else { /* do not use compressed row format */
933: for (i=0; i<m; i++) {
934: jrow = ii[i];
935: n = ii[i+1] - jrow;
936: sum = y[i];
937: for (j=0; j<n; j++) {
938: sum += aa[jrow]*x[aj[jrow]]; jrow++;
939: }
940: z[i] = sum;
941: }
942: }
943: #endif
944: PetscLogFlops(2*a->nz);
945: VecRestoreArray(xx,&x);
946: VecRestoreArray(yy,&y);
947: if (zz != yy) {
948: VecRestoreArray(zz,&z);
949: }
950: return(0);
951: }
953: /*
954: Adds diagonal pointers to sparse matrix structure.
955: */
958: PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat A)
959: {
960: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
962: PetscInt i,j,*diag,m = A->m;
965: if (a->diag) return(0);
967: PetscMalloc((m+1)*sizeof(PetscInt),&diag);
968: PetscLogObjectMemory(A,(m+1)*sizeof(PetscInt));
969: for (i=0; i<A->m; i++) {
970: diag[i] = a->i[i+1];
971: for (j=a->i[i]; j<a->i[i+1]; j++) {
972: if (a->j[j] == i) {
973: diag[i] = j;
974: break;
975: }
976: }
977: }
978: a->diag = diag;
979: return(0);
980: }
982: /*
983: Checks for missing diagonals
984: */
987: PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat A)
988: {
989: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
991: PetscInt *diag,*jj = a->j,i;
994: MatMarkDiagonal_SeqAIJ(A);
995: diag = a->diag;
996: for (i=0; i<A->m; i++) {
997: if (jj[diag[i]] != i) {
998: SETERRQ1(PETSC_ERR_PLIB,"Matrix is missing diagonal number %D",i);
999: }
1000: }
1001: return(0);
1002: }
1006: PetscErrorCode MatRelax_SeqAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
1007: {
1008: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1009: PetscScalar *x,d,*xs,sum,*t,scale,*idiag=0,*mdiag;
1010: const PetscScalar *v = a->a, *b, *bs,*xb, *ts;
1011: PetscErrorCode ierr;
1012: PetscInt n = A->n,m = A->m,i;
1013: const PetscInt *idx,*diag;
1016: its = its*lits;
1017: if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
1019: if (!a->diag) {MatMarkDiagonal_SeqAIJ(A);}
1020: diag = a->diag;
1021: if (!a->idiag) {
1022: PetscMalloc(3*m*sizeof(PetscScalar),&a->idiag);
1023: a->ssor = a->idiag + m;
1024: mdiag = a->ssor + m;
1026: v = a->a;
1028: /* this is wrong when fshift omega changes each iteration */
1029: if (omega == 1.0 && !fshift) {
1030: for (i=0; i<m; i++) {
1031: mdiag[i] = v[diag[i]];
1032: a->idiag[i] = 1.0/v[diag[i]];
1033: }
1034: PetscLogFlops(m);
1035: } else {
1036: for (i=0; i<m; i++) {
1037: mdiag[i] = v[diag[i]];
1038: a->idiag[i] = omega/(fshift + v[diag[i]]);
1039: }
1040: PetscLogFlops(2*m);
1041: }
1042: }
1043: t = a->ssor;
1044: idiag = a->idiag;
1045: mdiag = a->idiag + 2*m;
1047: VecGetArray(xx,&x);
1048: if (xx != bb) {
1049: VecGetArray(bb,(PetscScalar**)&b);
1050: } else {
1051: b = x;
1052: }
1054: /* We count flops by assuming the upper triangular and lower triangular parts have the same number of nonzeros */
1055: xs = x;
1056: if (flag == SOR_APPLY_UPPER) {
1057: /* apply (U + D/omega) to the vector */
1058: bs = b;
1059: for (i=0; i<m; i++) {
1060: d = fshift + a->a[diag[i]];
1061: n = a->i[i+1] - diag[i] - 1;
1062: idx = a->j + diag[i] + 1;
1063: v = a->a + diag[i] + 1;
1064: sum = b[i]*d/omega;
1065: SPARSEDENSEDOT(sum,bs,v,idx,n);
1066: x[i] = sum;
1067: }
1068: VecRestoreArray(xx,&x);
1069: if (bb != xx) {VecRestoreArray(bb,(PetscScalar**)&b);}
1070: PetscLogFlops(a->nz);
1071: return(0);
1072: }
1075: /* Let A = L + U + D; where L is lower trianglar,
1076: U is upper triangular, E is diagonal; This routine applies
1078: (L + E)^{-1} A (U + E)^{-1}
1080: to a vector efficiently using Eisenstat's trick. This is for
1081: the case of SSOR preconditioner, so E is D/omega where omega
1082: is the relaxation factor.
1083: */
1085: if (flag == SOR_APPLY_LOWER) {
1086: SETERRQ(PETSC_ERR_SUP,"SOR_APPLY_LOWER is not implemented");
1087: } else if (flag & SOR_EISENSTAT) {
1088: /* Let A = L + U + D; where L is lower trianglar,
1089: U is upper triangular, E is diagonal; This routine applies
1091: (L + E)^{-1} A (U + E)^{-1}
1093: to a vector efficiently using Eisenstat's trick. This is for
1094: the case of SSOR preconditioner, so E is D/omega where omega
1095: is the relaxation factor.
1096: */
1097: scale = (2.0/omega) - 1.0;
1099: /* x = (E + U)^{-1} b */
1100: for (i=m-1; i>=0; i--) {
1101: n = a->i[i+1] - diag[i] - 1;
1102: idx = a->j + diag[i] + 1;
1103: v = a->a + diag[i] + 1;
1104: sum = b[i];
1105: SPARSEDENSEMDOT(sum,xs,v,idx,n);
1106: x[i] = sum*idiag[i];
1107: }
1109: /* t = b - (2*E - D)x */
1110: v = a->a;
1111: for (i=0; i<m; i++) { t[i] = b[i] - scale*(v[*diag++])*x[i]; }
1113: /* t = (E + L)^{-1}t */
1114: ts = t;
1115: diag = a->diag;
1116: for (i=0; i<m; i++) {
1117: n = diag[i] - a->i[i];
1118: idx = a->j + a->i[i];
1119: v = a->a + a->i[i];
1120: sum = t[i];
1121: SPARSEDENSEMDOT(sum,ts,v,idx,n);
1122: t[i] = sum*idiag[i];
1123: /* x = x + t */
1124: x[i] += t[i];
1125: }
1127: PetscLogFlops(6*m-1 + 2*a->nz);
1128: VecRestoreArray(xx,&x);
1129: if (bb != xx) {VecRestoreArray(bb,(PetscScalar**)&b);}
1130: return(0);
1131: }
1132: if (flag & SOR_ZERO_INITIAL_GUESS) {
1133: if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
1134: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1135: fortranrelaxaijforwardzero_(&m,&omega,x,a->i,a->j,(PetscInt*)diag,idiag,a->a,(void*)b);
1136: #else
1137: for (i=0; i<m; i++) {
1138: n = diag[i] - a->i[i];
1139: idx = a->j + a->i[i];
1140: v = a->a + a->i[i];
1141: sum = b[i];
1142: SPARSEDENSEMDOT(sum,xs,v,idx,n);
1143: x[i] = sum*idiag[i];
1144: }
1145: #endif
1146: xb = x;
1147: PetscLogFlops(a->nz);
1148: } else xb = b;
1149: if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
1150: (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
1151: for (i=0; i<m; i++) {
1152: x[i] *= mdiag[i];
1153: }
1154: PetscLogFlops(m);
1155: }
1156: if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
1157: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1158: fortranrelaxaijbackwardzero_(&m,&omega,x,a->i,a->j,(PetscInt*)diag,idiag,a->a,(void*)xb);
1159: #else
1160: for (i=m-1; i>=0; i--) {
1161: n = a->i[i+1] - diag[i] - 1;
1162: idx = a->j + diag[i] + 1;
1163: v = a->a + diag[i] + 1;
1164: sum = xb[i];
1165: SPARSEDENSEMDOT(sum,xs,v,idx,n);
1166: x[i] = sum*idiag[i];
1167: }
1168: #endif
1169: PetscLogFlops(a->nz);
1170: }
1171: its--;
1172: }
1173: while (its--) {
1174: if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
1175: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1176: fortranrelaxaijforward_(&m,&omega,x,a->i,a->j,(PetscInt*)diag,a->a,(void*)b);
1177: #else
1178: for (i=0; i<m; i++) {
1179: d = fshift + a->a[diag[i]];
1180: n = a->i[i+1] - a->i[i];
1181: idx = a->j + a->i[i];
1182: v = a->a + a->i[i];
1183: sum = b[i];
1184: SPARSEDENSEMDOT(sum,xs,v,idx,n);
1185: x[i] = (1. - omega)*x[i] + (sum + mdiag[i]*x[i])*idiag[i];
1186: }
1187: #endif
1188: PetscLogFlops(a->nz);
1189: }
1190: if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
1191: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1192: fortranrelaxaijbackward_(&m,&omega,x,a->i,a->j,(PetscInt*)diag,a->a,(void*)b);
1193: #else
1194: for (i=m-1; i>=0; i--) {
1195: d = fshift + a->a[diag[i]];
1196: n = a->i[i+1] - a->i[i];
1197: idx = a->j + a->i[i];
1198: v = a->a + a->i[i];
1199: sum = b[i];
1200: SPARSEDENSEMDOT(sum,xs,v,idx,n);
1201: x[i] = (1. - omega)*x[i] + (sum + mdiag[i]*x[i])*idiag[i];
1202: }
1203: #endif
1204: PetscLogFlops(a->nz);
1205: }
1206: }
1207: VecRestoreArray(xx,&x);
1208: if (bb != xx) {VecRestoreArray(bb,(PetscScalar**)&b);}
1209: return(0);
1210: }
1214: PetscErrorCode MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,MatInfo *info)
1215: {
1216: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1219: info->rows_global = (double)A->m;
1220: info->columns_global = (double)A->n;
1221: info->rows_local = (double)A->m;
1222: info->columns_local = (double)A->n;
1223: info->block_size = 1.0;
1224: info->nz_allocated = (double)a->maxnz;
1225: info->nz_used = (double)a->nz;
1226: info->nz_unneeded = (double)(a->maxnz - a->nz);
1227: info->assemblies = (double)A->num_ass;
1228: info->mallocs = (double)a->reallocs;
1229: info->memory = A->mem;
1230: if (A->factor) {
1231: info->fill_ratio_given = A->info.fill_ratio_given;
1232: info->fill_ratio_needed = A->info.fill_ratio_needed;
1233: info->factor_mallocs = A->info.factor_mallocs;
1234: } else {
1235: info->fill_ratio_given = 0;
1236: info->fill_ratio_needed = 0;
1237: info->factor_mallocs = 0;
1238: }
1239: return(0);
1240: }
1244: PetscErrorCode MatZeroRows_SeqAIJ(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag)
1245: {
1246: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1247: PetscInt i,m = A->m - 1;
1251: if (a->keepzeroedrows) {
1252: for (i=0; i<N; i++) {
1253: if (rows[i] < 0 || rows[i] > m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range", rows[i]);
1254: PetscMemzero(&a->a[a->i[rows[i]]],a->ilen[rows[i]]*sizeof(PetscScalar));
1255: }
1256: if (diag != 0.0) {
1257: MatMissingDiagonal_SeqAIJ(A);
1258: MatMarkDiagonal_SeqAIJ(A);
1259: for (i=0; i<N; i++) {
1260: a->a[a->diag[rows[i]]] = diag;
1261: }
1262: }
1263: A->same_nonzero = PETSC_TRUE;
1264: } else {
1265: if (diag != 0.0) {
1266: for (i=0; i<N; i++) {
1267: if (rows[i] < 0 || rows[i] > m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range", rows[i]);
1268: if (a->ilen[rows[i]] > 0) {
1269: a->ilen[rows[i]] = 1;
1270: a->a[a->i[rows[i]]] = diag;
1271: a->j[a->i[rows[i]]] = rows[i];
1272: } else { /* in case row was completely empty */
1273: MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],&diag,INSERT_VALUES);
1274: }
1275: }
1276: } else {
1277: for (i=0; i<N; i++) {
1278: if (rows[i] < 0 || rows[i] > m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range", rows[i]);
1279: a->ilen[rows[i]] = 0;
1280: }
1281: }
1282: A->same_nonzero = PETSC_FALSE;
1283: }
1284: MatAssemblyEnd_SeqAIJ(A,MAT_FINAL_ASSEMBLY);
1285: return(0);
1286: }
1290: PetscErrorCode MatGetRow_SeqAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
1291: {
1292: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1293: PetscInt *itmp;
1296: if (row < 0 || row >= A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Row %D out of range",row);
1298: *nz = a->i[row+1] - a->i[row];
1299: if (v) *v = a->a + a->i[row];
1300: if (idx) {
1301: itmp = a->j + a->i[row];
1302: if (*nz) {
1303: *idx = itmp;
1304: }
1305: else *idx = 0;
1306: }
1307: return(0);
1308: }
1310: /* remove this function? */
1313: PetscErrorCode MatRestoreRow_SeqAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
1314: {
1316: return(0);
1317: }
1321: PetscErrorCode MatNorm_SeqAIJ(Mat A,NormType type,PetscReal *nrm)
1322: {
1323: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1324: PetscScalar *v = a->a;
1325: PetscReal sum = 0.0;
1327: PetscInt i,j;
1330: if (type == NORM_FROBENIUS) {
1331: for (i=0; i<a->nz; i++) {
1332: #if defined(PETSC_USE_COMPLEX)
1333: sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1334: #else
1335: sum += (*v)*(*v); v++;
1336: #endif
1337: }
1338: *nrm = sqrt(sum);
1339: } else if (type == NORM_1) {
1340: PetscReal *tmp;
1341: PetscInt *jj = a->j;
1342: PetscMalloc((A->n+1)*sizeof(PetscReal),&tmp);
1343: PetscMemzero(tmp,A->n*sizeof(PetscReal));
1344: *nrm = 0.0;
1345: for (j=0; j<a->nz; j++) {
1346: tmp[*jj++] += PetscAbsScalar(*v); v++;
1347: }
1348: for (j=0; j<A->n; j++) {
1349: if (tmp[j] > *nrm) *nrm = tmp[j];
1350: }
1351: PetscFree(tmp);
1352: } else if (type == NORM_INFINITY) {
1353: *nrm = 0.0;
1354: for (j=0; j<A->m; j++) {
1355: v = a->a + a->i[j];
1356: sum = 0.0;
1357: for (i=0; i<a->i[j+1]-a->i[j]; i++) {
1358: sum += PetscAbsScalar(*v); v++;
1359: }
1360: if (sum > *nrm) *nrm = sum;
1361: }
1362: } else {
1363: SETERRQ(PETSC_ERR_SUP,"No support for two norm");
1364: }
1365: return(0);
1366: }
1370: PetscErrorCode MatTranspose_SeqAIJ(Mat A,Mat *B)
1371: {
1372: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1373: Mat C;
1375: PetscInt i,*aj = a->j,*ai = a->i,m = A->m,len,*col;
1376: PetscScalar *array = a->a;
1379: if (!B && m != A->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Square matrix only for in-place");
1380: PetscMalloc((1+A->n)*sizeof(PetscInt),&col);
1381: PetscMemzero(col,(1+A->n)*sizeof(PetscInt));
1382:
1383: for (i=0; i<ai[m]; i++) col[aj[i]] += 1;
1384: MatCreate(A->comm,&C);
1385: MatSetSizes(C,A->n,m,A->n,m);
1386: MatSetType(C,A->type_name);
1387: MatSeqAIJSetPreallocation_SeqAIJ(C,0,col);
1388: PetscFree(col);
1389: for (i=0; i<m; i++) {
1390: len = ai[i+1]-ai[i];
1391: MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES);
1392: array += len;
1393: aj += len;
1394: }
1396: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1397: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
1399: if (B) {
1400: *B = C;
1401: } else {
1402: MatHeaderCopy(A,C);
1403: }
1404: return(0);
1405: }
1410: PetscErrorCode PETSCMAT_DLLEXPORT MatIsTranspose_SeqAIJ(Mat A,Mat B,PetscReal tol,PetscTruth *f)
1411: {
1412: Mat_SeqAIJ *aij = (Mat_SeqAIJ *) A->data,*bij = (Mat_SeqAIJ*) A->data;
1413: PetscInt *adx,*bdx,*aii,*bii,*aptr,*bptr; PetscScalar *va,*vb;
1415: PetscInt ma,na,mb,nb, i;
1418: bij = (Mat_SeqAIJ *) B->data;
1419:
1420: MatGetSize(A,&ma,&na);
1421: MatGetSize(B,&mb,&nb);
1422: if (ma!=nb || na!=mb){
1423: *f = PETSC_FALSE;
1424: return(0);
1425: }
1426: aii = aij->i; bii = bij->i;
1427: adx = aij->j; bdx = bij->j;
1428: va = aij->a; vb = bij->a;
1429: PetscMalloc(ma*sizeof(PetscInt),&aptr);
1430: PetscMalloc(mb*sizeof(PetscInt),&bptr);
1431: for (i=0; i<ma; i++) aptr[i] = aii[i];
1432: for (i=0; i<mb; i++) bptr[i] = bii[i];
1434: *f = PETSC_TRUE;
1435: for (i=0; i<ma; i++) {
1436: while (aptr[i]<aii[i+1]) {
1437: PetscInt idc,idr;
1438: PetscScalar vc,vr;
1439: /* column/row index/value */
1440: idc = adx[aptr[i]];
1441: idr = bdx[bptr[idc]];
1442: vc = va[aptr[i]];
1443: vr = vb[bptr[idc]];
1444: if (i!=idr || PetscAbsScalar(vc-vr) > tol) {
1445: *f = PETSC_FALSE;
1446: goto done;
1447: } else {
1448: aptr[i]++;
1449: if (B || i!=idc) bptr[idc]++;
1450: }
1451: }
1452: }
1453: done:
1454: PetscFree(aptr);
1455: if (B) {
1456: PetscFree(bptr);
1457: }
1458: return(0);
1459: }
1464: PetscErrorCode MatIsSymmetric_SeqAIJ(Mat A,PetscReal tol,PetscTruth *f)
1465: {
1468: MatIsTranspose_SeqAIJ(A,A,tol,f);
1469: return(0);
1470: }
1474: PetscErrorCode MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr)
1475: {
1476: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1477: PetscScalar *l,*r,x,*v;
1479: PetscInt i,j,m = A->m,n = A->n,M,nz = a->nz,*jj;
1482: if (ll) {
1483: /* The local size is used so that VecMPI can be passed to this routine
1484: by MatDiagonalScale_MPIAIJ */
1485: VecGetLocalSize(ll,&m);
1486: if (m != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length");
1487: VecGetArray(ll,&l);
1488: v = a->a;
1489: for (i=0; i<m; i++) {
1490: x = l[i];
1491: M = a->i[i+1] - a->i[i];
1492: for (j=0; j<M; j++) { (*v++) *= x;}
1493: }
1494: VecRestoreArray(ll,&l);
1495: PetscLogFlops(nz);
1496: }
1497: if (rr) {
1498: VecGetLocalSize(rr,&n);
1499: if (n != A->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vector wrong length");
1500: VecGetArray(rr,&r);
1501: v = a->a; jj = a->j;
1502: for (i=0; i<nz; i++) {
1503: (*v++) *= r[*jj++];
1504: }
1505: VecRestoreArray(rr,&r);
1506: PetscLogFlops(nz);
1507: }
1508: return(0);
1509: }
1513: PetscErrorCode MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,PetscInt csize,MatReuse scall,Mat *B)
1514: {
1515: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*c;
1517: PetscInt *smap,i,k,kstart,kend,oldcols = A->n,*lens;
1518: PetscInt row,mat_i,*mat_j,tcol,first,step,*mat_ilen,sum,lensi;
1519: PetscInt *irow,*icol,nrows,ncols;
1520: PetscInt *starts,*j_new,*i_new,*aj = a->j,*ai = a->i,ii,*ailen = a->ilen;
1521: PetscScalar *a_new,*mat_a;
1522: Mat C;
1523: PetscTruth stride;
1526: ISSorted(isrow,(PetscTruth*)&i);
1527: if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"ISrow is not sorted");
1528: ISSorted(iscol,(PetscTruth*)&i);
1529: if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IScol is not sorted");
1531: ISGetIndices(isrow,&irow);
1532: ISGetLocalSize(isrow,&nrows);
1533: ISGetLocalSize(iscol,&ncols);
1535: ISStrideGetInfo(iscol,&first,&step);
1536: ISStride(iscol,&stride);
1537: if (stride && step == 1) {
1538: /* special case of contiguous rows */
1539: PetscMalloc((2*nrows+1)*sizeof(PetscInt),&lens);
1540: starts = lens + nrows;
1541: /* loop over new rows determining lens and starting points */
1542: for (i=0; i<nrows; i++) {
1543: kstart = ai[irow[i]];
1544: kend = kstart + ailen[irow[i]];
1545: for (k=kstart; k<kend; k++) {
1546: if (aj[k] >= first) {
1547: starts[i] = k;
1548: break;
1549: }
1550: }
1551: sum = 0;
1552: while (k < kend) {
1553: if (aj[k++] >= first+ncols) break;
1554: sum++;
1555: }
1556: lens[i] = sum;
1557: }
1558: /* create submatrix */
1559: if (scall == MAT_REUSE_MATRIX) {
1560: PetscInt n_cols,n_rows;
1561: MatGetSize(*B,&n_rows,&n_cols);
1562: if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size");
1563: MatZeroEntries(*B);
1564: C = *B;
1565: } else {
1566: MatCreate(A->comm,&C);
1567: MatSetSizes(C,nrows,ncols,PETSC_DETERMINE,PETSC_DETERMINE);
1568: MatSetType(C,A->type_name);
1569: MatSeqAIJSetPreallocation_SeqAIJ(C,0,lens);
1570: }
1571: c = (Mat_SeqAIJ*)C->data;
1573: /* loop over rows inserting into submatrix */
1574: a_new = c->a;
1575: j_new = c->j;
1576: i_new = c->i;
1578: for (i=0; i<nrows; i++) {
1579: ii = starts[i];
1580: lensi = lens[i];
1581: for (k=0; k<lensi; k++) {
1582: *j_new++ = aj[ii+k] - first;
1583: }
1584: PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(PetscScalar));
1585: a_new += lensi;
1586: i_new[i+1] = i_new[i] + lensi;
1587: c->ilen[i] = lensi;
1588: }
1589: PetscFree(lens);
1590: } else {
1591: ISGetIndices(iscol,&icol);
1592: PetscMalloc((1+oldcols)*sizeof(PetscInt),&smap);
1593:
1594: PetscMalloc((1+nrows)*sizeof(PetscInt),&lens);
1595: PetscMemzero(smap,oldcols*sizeof(PetscInt));
1596: for (i=0; i<ncols; i++) smap[icol[i]] = i+1;
1597: /* determine lens of each row */
1598: for (i=0; i<nrows; i++) {
1599: kstart = ai[irow[i]];
1600: kend = kstart + a->ilen[irow[i]];
1601: lens[i] = 0;
1602: for (k=kstart; k<kend; k++) {
1603: if (smap[aj[k]]) {
1604: lens[i]++;
1605: }
1606: }
1607: }
1608: /* Create and fill new matrix */
1609: if (scall == MAT_REUSE_MATRIX) {
1610: PetscTruth equal;
1612: c = (Mat_SeqAIJ *)((*B)->data);
1613: if ((*B)->m != nrows || (*B)->n != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size");
1614: PetscMemcmp(c->ilen,lens,(*B)->m*sizeof(PetscInt),&equal);
1615: if (!equal) {
1616: SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros");
1617: }
1618: PetscMemzero(c->ilen,(*B)->m*sizeof(PetscInt));
1619: C = *B;
1620: } else {
1621: MatCreate(A->comm,&C);
1622: MatSetSizes(C,nrows,ncols,PETSC_DETERMINE,PETSC_DETERMINE);
1623: MatSetType(C,A->type_name);
1624: MatSeqAIJSetPreallocation_SeqAIJ(C,0,lens);
1625: }
1626: c = (Mat_SeqAIJ *)(C->data);
1627: for (i=0; i<nrows; i++) {
1628: row = irow[i];
1629: kstart = ai[row];
1630: kend = kstart + a->ilen[row];
1631: mat_i = c->i[i];
1632: mat_j = c->j + mat_i;
1633: mat_a = c->a + mat_i;
1634: mat_ilen = c->ilen + i;
1635: for (k=kstart; k<kend; k++) {
1636: if ((tcol=smap[a->j[k]])) {
1637: *mat_j++ = tcol - 1;
1638: *mat_a++ = a->a[k];
1639: (*mat_ilen)++;
1641: }
1642: }
1643: }
1644: /* Free work space */
1645: ISRestoreIndices(iscol,&icol);
1646: PetscFree(smap);
1647: PetscFree(lens);
1648: }
1649: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1650: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
1652: ISRestoreIndices(isrow,&irow);
1653: *B = C;
1654: return(0);
1655: }
1657: /*
1658: */
1661: PetscErrorCode MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,MatFactorInfo *info)
1662: {
1663: Mat_SeqAIJ *a = (Mat_SeqAIJ*)inA->data;
1665: Mat outA;
1666: PetscTruth row_identity,col_identity;
1669: if (info->levels != 0) SETERRQ(PETSC_ERR_SUP,"Only levels=0 supported for in-place ilu");
1670: ISIdentity(row,&row_identity);
1671: ISIdentity(col,&col_identity);
1672: if (!row_identity || !col_identity) {
1673: SETERRQ(PETSC_ERR_ARG_WRONG,"Row and column permutations must be identity for in-place ILU");
1674: }
1676: outA = inA;
1677: inA->factor = FACTOR_LU;
1678: a->row = row;
1679: a->col = col;
1680: PetscObjectReference((PetscObject)row);
1681: PetscObjectReference((PetscObject)col);
1683: /* Create the inverse permutation so that it can be used in MatLUFactorNumeric() */
1684: if (a->icol) {ISDestroy(a->icol);} /* need to remove old one */
1685: ISInvertPermutation(col,PETSC_DECIDE,&a->icol);
1686: PetscLogObjectParent(inA,a->icol);
1688: if (!a->solve_work) { /* this matrix may have been factored before */
1689: PetscMalloc((inA->m+1)*sizeof(PetscScalar),&a->solve_work);
1690: }
1692: if (!a->diag) {
1693: MatMarkDiagonal_SeqAIJ(inA);
1694: }
1695: MatLUFactorNumeric_SeqAIJ(inA,info,&outA);
1696: return(0);
1697: }
1699: #include petscblaslapack.h
1702: PetscErrorCode MatScale_SeqAIJ(Mat inA,PetscScalar alpha)
1703: {
1704: Mat_SeqAIJ *a = (Mat_SeqAIJ*)inA->data;
1705: PetscBLASInt bnz = (PetscBLASInt)a->nz,one = 1;
1706: PetscScalar oalpha = alpha;
1711: BLASscal_(&bnz,&oalpha,a->a,&one);
1712: PetscLogFlops(a->nz);
1713: return(0);
1714: }
1718: PetscErrorCode MatGetSubMatrices_SeqAIJ(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1719: {
1721: PetscInt i;
1724: if (scall == MAT_INITIAL_MATRIX) {
1725: PetscMalloc((n+1)*sizeof(Mat),B);
1726: }
1728: for (i=0; i<n; i++) {
1729: MatGetSubMatrix_SeqAIJ(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);
1730: }
1731: return(0);
1732: }
1736: PetscErrorCode MatIncreaseOverlap_SeqAIJ(Mat A,PetscInt is_max,IS is[],PetscInt ov)
1737: {
1738: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1740: PetscInt row,i,j,k,l,m,n,*idx,*nidx,isz,val;
1741: PetscInt start,end,*ai,*aj;
1742: PetscBT table;
1745: m = A->m;
1746: ai = a->i;
1747: aj = a->j;
1749: if (ov < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"illegal negative overlap value used");
1751: PetscMalloc((m+1)*sizeof(PetscInt),&nidx);
1752: PetscBTCreate(m,table);
1754: for (i=0; i<is_max; i++) {
1755: /* Initialize the two local arrays */
1756: isz = 0;
1757: PetscBTMemzero(m,table);
1758:
1759: /* Extract the indices, assume there can be duplicate entries */
1760: ISGetIndices(is[i],&idx);
1761: ISGetLocalSize(is[i],&n);
1762:
1763: /* Enter these into the temp arrays. I.e., mark table[row], enter row into new index */
1764: for (j=0; j<n ; ++j){
1765: if(!PetscBTLookupSet(table,idx[j])) { nidx[isz++] = idx[j];}
1766: }
1767: ISRestoreIndices(is[i],&idx);
1768: ISDestroy(is[i]);
1769:
1770: k = 0;
1771: for (j=0; j<ov; j++){ /* for each overlap */
1772: n = isz;
1773: for (; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */
1774: row = nidx[k];
1775: start = ai[row];
1776: end = ai[row+1];
1777: for (l = start; l<end ; l++){
1778: val = aj[l] ;
1779: if (!PetscBTLookupSet(table,val)) {nidx[isz++] = val;}
1780: }
1781: }
1782: }
1783: ISCreateGeneral(PETSC_COMM_SELF,isz,nidx,(is+i));
1784: }
1785: PetscBTDestroy(table);
1786: PetscFree(nidx);
1787: return(0);
1788: }
1790: /* -------------------------------------------------------------- */
1793: PetscErrorCode MatPermute_SeqAIJ(Mat A,IS rowp,IS colp,Mat *B)
1794: {
1795: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1797: PetscInt i,nz,m = A->m,n = A->n,*col;
1798: PetscInt *row,*cnew,j,*lens;
1799: IS icolp,irowp;
1800: PetscInt *cwork;
1801: PetscScalar *vwork;
1804: ISInvertPermutation(rowp,PETSC_DECIDE,&irowp);
1805: ISGetIndices(irowp,&row);
1806: ISInvertPermutation(colp,PETSC_DECIDE,&icolp);
1807: ISGetIndices(icolp,&col);
1808:
1809: /* determine lengths of permuted rows */
1810: PetscMalloc((m+1)*sizeof(PetscInt),&lens);
1811: for (i=0; i<m; i++) {
1812: lens[row[i]] = a->i[i+1] - a->i[i];
1813: }
1814: MatCreate(A->comm,B);
1815: MatSetSizes(*B,m,n,m,n);
1816: MatSetType(*B,A->type_name);
1817: MatSeqAIJSetPreallocation_SeqAIJ(*B,0,lens);
1818: PetscFree(lens);
1820: PetscMalloc(n*sizeof(PetscInt),&cnew);
1821: for (i=0; i<m; i++) {
1822: MatGetRow_SeqAIJ(A,i,&nz,&cwork,&vwork);
1823: for (j=0; j<nz; j++) { cnew[j] = col[cwork[j]];}
1824: MatSetValues_SeqAIJ(*B,1,&row[i],nz,cnew,vwork,INSERT_VALUES);
1825: MatRestoreRow_SeqAIJ(A,i,&nz,&cwork,&vwork);
1826: }
1827: PetscFree(cnew);
1828: (*B)->assembled = PETSC_FALSE;
1829: MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);
1830: MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);
1831: ISRestoreIndices(irowp,&row);
1832: ISRestoreIndices(icolp,&col);
1833: ISDestroy(irowp);
1834: ISDestroy(icolp);
1835: return(0);
1836: }
1840: PetscErrorCode MatPrintHelp_SeqAIJ(Mat A)
1841: {
1842: static PetscTruth called = PETSC_FALSE;
1843: MPI_Comm comm = A->comm;
1844: PetscErrorCode ierr;
1847: MatPrintHelp_Inode(A);
1848: if (called) {return(0);} else called = PETSC_TRUE;
1849: (*PetscHelpPrintf)(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n");
1850: (*PetscHelpPrintf)(comm," -mat_lu_pivotthreshold <threshold>: Set pivoting threshold\n");
1851: (*PetscHelpPrintf)(comm," -mat_aij_oneindex: internal indices begin at 1 instead of the default 0.\n");
1852: (*PetscHelpPrintf)(comm," -mat_no_compressedrow: Do not use compressedrow\n");
1853: return(0);
1854: }
1858: PetscErrorCode MatCopy_SeqAIJ(Mat A,Mat B,MatStructure str)
1859: {
1863: /* If the two matrices have the same copy implementation, use fast copy. */
1864: if (str == SAME_NONZERO_PATTERN && (A->ops->copy == B->ops->copy)) {
1865: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1866: Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data;
1868: if (a->i[A->m] != b->i[B->m]) {
1869: SETERRQ(PETSC_ERR_ARG_INCOMP,"Number of nonzeros in two matrices are different");
1870: }
1871: PetscMemcpy(b->a,a->a,(a->i[A->m])*sizeof(PetscScalar));
1872: } else {
1873: MatCopy_Basic(A,B,str);
1874: }
1875: return(0);
1876: }
1880: PetscErrorCode MatSetUpPreallocation_SeqAIJ(Mat A)
1881: {
1885: MatSeqAIJSetPreallocation_SeqAIJ(A,PETSC_DEFAULT,0);
1886: return(0);
1887: }
1891: PetscErrorCode MatGetArray_SeqAIJ(Mat A,PetscScalar *array[])
1892: {
1893: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1895: *array = a->a;
1896: return(0);
1897: }
1901: PetscErrorCode MatRestoreArray_SeqAIJ(Mat A,PetscScalar *array[])
1902: {
1904: return(0);
1905: }
1909: PetscErrorCode MatFDColoringApply_SeqAIJ(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
1910: {
1911: PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void *))coloring->f;
1913: PetscInt k,N,start,end,l,row,col,srow,**vscaleforrow,m1,m2;
1914: PetscScalar dx,mone = -1.0,*y,*xx,*w3_array;
1915: PetscScalar *vscale_array;
1916: PetscReal epsilon = coloring->error_rel,umin = coloring->umin;
1917: Vec w1,w2,w3;
1918: void *fctx = coloring->fctx;
1919: PetscTruth flg;
1922: if (!coloring->w1) {
1923: VecDuplicate(x1,&coloring->w1);
1924: PetscLogObjectParent(coloring,coloring->w1);
1925: VecDuplicate(x1,&coloring->w2);
1926: PetscLogObjectParent(coloring,coloring->w2);
1927: VecDuplicate(x1,&coloring->w3);
1928: PetscLogObjectParent(coloring,coloring->w3);
1929: }
1930: w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
1932: MatSetUnfactored(J);
1933: PetscOptionsHasName(coloring->prefix,"-mat_fd_coloring_dont_rezero",&flg);
1934: if (flg) {
1935: PetscLogInfo((coloring,"MatFDColoringApply_SeqAIJ: Not calling MatZeroEntries()\n"));
1936: } else {
1937: PetscTruth assembled;
1938: MatAssembled(J,&assembled);
1939: if (assembled) {
1940: MatZeroEntries(J);
1941: }
1942: }
1944: VecGetOwnershipRange(x1,&start,&end);
1945: VecGetSize(x1,&N);
1947: /*
1948: This is a horrible, horrible, hack. See DMMGComputeJacobian_Multigrid() it inproperly sets
1949: coloring->F for the coarser grids from the finest
1950: */
1951: if (coloring->F) {
1952: VecGetLocalSize(coloring->F,&m1);
1953: VecGetLocalSize(w1,&m2);
1954: if (m1 != m2) {
1955: coloring->F = 0;
1956: }
1957: }
1959: if (coloring->F) {
1960: w1 = coloring->F;
1961: coloring->F = 0;
1962: } else {
1963: PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);
1964: (*f)(sctx,x1,w1,fctx);
1965: PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);
1966: }
1968: /*
1969: Compute all the scale factors and share with other processors
1970: */
1971: VecGetArray(x1,&xx);xx = xx - start;
1972: VecGetArray(coloring->vscale,&vscale_array);vscale_array = vscale_array - start;
1973: for (k=0; k<coloring->ncolors; k++) {
1974: /*
1975: Loop over each column associated with color adding the
1976: perturbation to the vector w3.
1977: */
1978: for (l=0; l<coloring->ncolumns[k]; l++) {
1979: col = coloring->columns[k][l]; /* column of the matrix we are probing for */
1980: dx = xx[col];
1981: if (dx == 0.0) dx = 1.0;
1982: #if !defined(PETSC_USE_COMPLEX)
1983: if (dx < umin && dx >= 0.0) dx = umin;
1984: else if (dx < 0.0 && dx > -umin) dx = -umin;
1985: #else
1986: if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin;
1987: else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
1988: #endif
1989: dx *= epsilon;
1990: vscale_array[col] = 1.0/dx;
1991: }
1992: }
1993: vscale_array = vscale_array + start;VecRestoreArray(coloring->vscale,&vscale_array);
1994: VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);
1995: VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);
1997: /* VecView(coloring->vscale,PETSC_VIEWER_STDOUT_WORLD);
1998: VecView(x1,PETSC_VIEWER_STDOUT_WORLD);*/
2000: if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
2001: else vscaleforrow = coloring->columnsforrow;
2003: VecGetArray(coloring->vscale,&vscale_array);
2004: /*
2005: Loop over each color
2006: */
2007: for (k=0; k<coloring->ncolors; k++) {
2008: coloring->currentcolor = k;
2009: VecCopy(x1,w3);
2010: VecGetArray(w3,&w3_array);w3_array = w3_array - start;
2011: /*
2012: Loop over each column associated with color adding the
2013: perturbation to the vector w3.
2014: */
2015: for (l=0; l<coloring->ncolumns[k]; l++) {
2016: col = coloring->columns[k][l]; /* column of the matrix we are probing for */
2017: dx = xx[col];
2018: if (dx == 0.0) dx = 1.0;
2019: #if !defined(PETSC_USE_COMPLEX)
2020: if (dx < umin && dx >= 0.0) dx = umin;
2021: else if (dx < 0.0 && dx > -umin) dx = -umin;
2022: #else
2023: if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin;
2024: else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
2025: #endif
2026: dx *= epsilon;
2027: if (!PetscAbsScalar(dx)) SETERRQ(PETSC_ERR_PLIB,"Computed 0 differencing parameter");
2028: w3_array[col] += dx;
2029: }
2030: w3_array = w3_array + start; VecRestoreArray(w3,&w3_array);
2032: /*
2033: Evaluate function at x1 + dx (here dx is a vector of perturbations)
2034: */
2036: PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);
2037: (*f)(sctx,w3,w2,fctx);
2038: PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);
2039: VecAXPY(w2,mone,w1);
2041: /*
2042: Loop over rows of vector, putting results into Jacobian matrix
2043: */
2044: VecGetArray(w2,&y);
2045: for (l=0; l<coloring->nrows[k]; l++) {
2046: row = coloring->rows[k][l];
2047: col = coloring->columnsforrow[k][l];
2048: y[row] *= vscale_array[vscaleforrow[k][l]];
2049: srow = row + start;
2050: MatSetValues_SeqAIJ(J,1,&srow,1,&col,y+row,INSERT_VALUES);
2051: }
2052: VecRestoreArray(w2,&y);
2053: }
2054: coloring->currentcolor = k;
2055: VecRestoreArray(coloring->vscale,&vscale_array);
2056: xx = xx + start; VecRestoreArray(x1,&xx);
2057: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
2058: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
2059: return(0);
2060: }
2062: #include petscblaslapack.h
2065: PetscErrorCode MatAXPY_SeqAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str)
2066: {
2068: PetscInt i;
2069: Mat_SeqAIJ *x = (Mat_SeqAIJ *)X->data,*y = (Mat_SeqAIJ *)Y->data;
2070: PetscBLASInt one=1,bnz = (PetscBLASInt)x->nz;
2073: if (str == SAME_NONZERO_PATTERN) {
2074: PetscScalar alpha = a;
2075: BLASaxpy_(&bnz,&alpha,x->a,&one,y->a,&one);
2076: } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
2077: if (y->xtoy && y->XtoY != X) {
2078: PetscFree(y->xtoy);
2079: MatDestroy(y->XtoY);
2080: }
2081: if (!y->xtoy) { /* get xtoy */
2082: MatAXPYGetxtoy_Private(X->m,x->i,x->j,PETSC_NULL, y->i,y->j,PETSC_NULL, &y->xtoy);
2083: y->XtoY = X;
2084: }
2085: for (i=0; i<x->nz; i++) y->a[y->xtoy[i]] += a*(x->a[i]);
2086: PetscLogInfo((0,"MatAXPY_SeqAIJ: ratio of nnz(X)/nnz(Y): %d/%d = %g\n",x->nz,y->nz,(PetscReal)(x->nz)/y->nz));
2087: } else {
2088: MatAXPY_Basic(Y,a,X,str);
2089: }
2090: return(0);
2091: }
2095: PetscErrorCode MatSetBlockSize_SeqAIJ(Mat A,PetscInt bs)
2096: {
2098: return(0);
2099: }
2103: PetscErrorCode PETSCMAT_DLLEXPORT MatConjugate_SeqAIJ(Mat mat)
2104: {
2105: #if defined(PETSC_USE_COMPLEX)
2106: Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
2107: PetscInt i,nz;
2108: PetscScalar *a;
2111: nz = aij->nz;
2112: a = aij->a;
2113: for (i=0; i<nz; i++) {
2114: a[i] = PetscConj(a[i]);
2115: }
2116: #else
2118: #endif
2119: return(0);
2120: }
2124: PetscErrorCode MatGetRowMax_SeqAIJ(Mat A,Vec v)
2125: {
2126: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
2128: PetscInt i,j,m = A->m,*ai,*aj,ncols,n;
2129: PetscReal atmp;
2130: PetscScalar *x,zero = 0.0;
2131: MatScalar *aa;
2134: if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2135: aa = a->a;
2136: ai = a->i;
2137: aj = a->j;
2139: VecSet(v,zero);
2140: VecGetArray(v,&x);
2141: VecGetLocalSize(v,&n);
2142: if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
2143: for (i=0; i<m; i++) {
2144: ncols = ai[1] - ai[0]; ai++;
2145: for (j=0; j<ncols; j++){
2146: atmp = PetscAbsScalar(*aa); aa++;
2147: if (PetscAbsScalar(x[i]) < atmp) x[i] = atmp;
2148: aj++;
2149: }
2150: }
2151: VecRestoreArray(v,&x);
2152: return(0);
2153: }
2155: /* -------------------------------------------------------------------*/
2156: static struct _MatOps MatOps_Values = {MatSetValues_SeqAIJ,
2157: MatGetRow_SeqAIJ,
2158: MatRestoreRow_SeqAIJ,
2159: MatMult_SeqAIJ,
2160: /* 4*/ MatMultAdd_SeqAIJ,
2161: MatMultTranspose_SeqAIJ,
2162: MatMultTransposeAdd_SeqAIJ,
2163: MatSolve_SeqAIJ,
2164: MatSolveAdd_SeqAIJ,
2165: MatSolveTranspose_SeqAIJ,
2166: /*10*/ MatSolveTransposeAdd_SeqAIJ,
2167: MatLUFactor_SeqAIJ,
2168: 0,
2169: MatRelax_SeqAIJ,
2170: MatTranspose_SeqAIJ,
2171: /*15*/ MatGetInfo_SeqAIJ,
2172: MatEqual_SeqAIJ,
2173: MatGetDiagonal_SeqAIJ,
2174: MatDiagonalScale_SeqAIJ,
2175: MatNorm_SeqAIJ,
2176: /*20*/ 0,
2177: MatAssemblyEnd_SeqAIJ,
2178: MatCompress_SeqAIJ,
2179: MatSetOption_SeqAIJ,
2180: MatZeroEntries_SeqAIJ,
2181: /*25*/ MatZeroRows_SeqAIJ,
2182: MatLUFactorSymbolic_SeqAIJ,
2183: MatLUFactorNumeric_SeqAIJ,
2184: MatCholeskyFactorSymbolic_SeqAIJ,
2185: MatCholeskyFactorNumeric_SeqAIJ,
2186: /*30*/ MatSetUpPreallocation_SeqAIJ,
2187: MatILUFactorSymbolic_SeqAIJ,
2188: MatICCFactorSymbolic_SeqAIJ,
2189: MatGetArray_SeqAIJ,
2190: MatRestoreArray_SeqAIJ,
2191: /*35*/ MatDuplicate_SeqAIJ,
2192: 0,
2193: 0,
2194: MatILUFactor_SeqAIJ,
2195: 0,
2196: /*40*/ MatAXPY_SeqAIJ,
2197: MatGetSubMatrices_SeqAIJ,
2198: MatIncreaseOverlap_SeqAIJ,
2199: MatGetValues_SeqAIJ,
2200: MatCopy_SeqAIJ,
2201: /*45*/ MatPrintHelp_SeqAIJ,
2202: MatScale_SeqAIJ,
2203: 0,
2204: 0,
2205: MatILUDTFactor_SeqAIJ,
2206: /*50*/ MatSetBlockSize_SeqAIJ,
2207: MatGetRowIJ_SeqAIJ,
2208: MatRestoreRowIJ_SeqAIJ,
2209: MatGetColumnIJ_SeqAIJ,
2210: MatRestoreColumnIJ_SeqAIJ,
2211: /*55*/ MatFDColoringCreate_SeqAIJ,
2212: 0,
2213: 0,
2214: MatPermute_SeqAIJ,
2215: 0,
2216: /*60*/ 0,
2217: MatDestroy_SeqAIJ,
2218: MatView_SeqAIJ,
2219: MatGetPetscMaps_Petsc,
2220: 0,
2221: /*65*/ 0,
2222: 0,
2223: 0,
2224: 0,
2225: 0,
2226: /*70*/ MatGetRowMax_SeqAIJ,
2227: 0,
2228: MatSetColoring_SeqAIJ,
2229: #if defined(PETSC_HAVE_ADIC)
2230: MatSetValuesAdic_SeqAIJ,
2231: #else
2232: 0,
2233: #endif
2234: MatSetValuesAdifor_SeqAIJ,
2235: /*75*/ MatFDColoringApply_SeqAIJ,
2236: 0,
2237: 0,
2238: 0,
2239: 0,
2240: /*80*/ 0,
2241: 0,
2242: 0,
2243: 0,
2244: MatLoad_SeqAIJ,
2245: /*85*/ MatIsSymmetric_SeqAIJ,
2246: 0,
2247: 0,
2248: 0,
2249: 0,
2250: /*90*/ MatMatMult_SeqAIJ_SeqAIJ,
2251: MatMatMultSymbolic_SeqAIJ_SeqAIJ,
2252: MatMatMultNumeric_SeqAIJ_SeqAIJ,
2253: MatPtAP_Basic,
2254: MatPtAPSymbolic_SeqAIJ,
2255: /*95*/ MatPtAPNumeric_SeqAIJ,
2256: MatMatMultTranspose_SeqAIJ_SeqAIJ,
2257: MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ,
2258: MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ,
2259: MatPtAPSymbolic_SeqAIJ_SeqAIJ,
2260: /*100*/MatPtAPNumeric_SeqAIJ_SeqAIJ,
2261: 0,
2262: 0,
2263: MatConjugate_SeqAIJ
2264: };
2269: PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetColumnIndices_SeqAIJ(Mat mat,PetscInt *indices)
2270: {
2271: Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
2272: PetscInt i,nz,n;
2276: nz = aij->maxnz;
2277: n = mat->n;
2278: for (i=0; i<nz; i++) {
2279: aij->j[i] = indices[i];
2280: }
2281: aij->nz = nz;
2282: for (i=0; i<n; i++) {
2283: aij->ilen[i] = aij->imax[i];
2284: }
2286: return(0);
2287: }
2292: /*@
2293: MatSeqAIJSetColumnIndices - Set the column indices for all the rows
2294: in the matrix.
2296: Input Parameters:
2297: + mat - the SeqAIJ matrix
2298: - indices - the column indices
2300: Level: advanced
2302: Notes:
2303: This can be called if you have precomputed the nonzero structure of the
2304: matrix and want to provide it to the matrix object to improve the performance
2305: of the MatSetValues() operation.
2307: You MUST have set the correct numbers of nonzeros per row in the call to
2308: MatCreateSeqAIJ(), and the columns indices MUST be sorted.
2310: MUST be called before any calls to MatSetValues();
2312: The indices should start with zero, not one.
2314: @*/
2315: PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetColumnIndices(Mat mat,PetscInt *indices)
2316: {
2317: PetscErrorCode ierr,(*f)(Mat,PetscInt *);
2322: PetscObjectQueryFunction((PetscObject)mat,"MatSeqAIJSetColumnIndices_C",(void (**)(void))&f);
2323: if (f) {
2324: (*f)(mat,indices);
2325: } else {
2326: SETERRQ(PETSC_ERR_SUP,"Wrong type of matrix to set column indices");
2327: }
2328: return(0);
2329: }
2331: /* ----------------------------------------------------------------------------------------*/
2336: PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues_SeqAIJ(Mat mat)
2337: {
2338: Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
2340: size_t nz = aij->i[mat->m];
2343: if (aij->nonew != 1) {
2344: SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
2345: }
2347: /* allocate space for values if not already there */
2348: if (!aij->saved_values) {
2349: PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);
2350: }
2352: /* copy values over */
2353: PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));
2354: return(0);
2355: }
2360: /*@
2361: MatStoreValues - Stashes a copy of the matrix values; this allows, for
2362: example, reuse of the linear part of a Jacobian, while recomputing the
2363: nonlinear portion.
2365: Collect on Mat
2367: Input Parameters:
2368: . mat - the matrix (currently only AIJ matrices support this option)
2370: Level: advanced
2372: Common Usage, with SNESSolve():
2373: $ Create Jacobian matrix
2374: $ Set linear terms into matrix
2375: $ Apply boundary conditions to matrix, at this time matrix must have
2376: $ final nonzero structure (i.e. setting the nonlinear terms and applying
2377: $ boundary conditions again will not change the nonzero structure
2378: $ MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS);
2379: $ MatStoreValues(mat);
2380: $ Call SNESSetJacobian() with matrix
2381: $ In your Jacobian routine
2382: $ MatRetrieveValues(mat);
2383: $ Set nonlinear terms in matrix
2384:
2385: Common Usage without SNESSolve(), i.e. when you handle nonlinear solve yourself:
2386: $ // build linear portion of Jacobian
2387: $ MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS);
2388: $ MatStoreValues(mat);
2389: $ loop over nonlinear iterations
2390: $ MatRetrieveValues(mat);
2391: $ // call MatSetValues(mat,...) to set nonliner portion of Jacobian
2392: $ // call MatAssemblyBegin/End() on matrix
2393: $ Solve linear system with Jacobian
2394: $ endloop
2396: Notes:
2397: Matrix must already be assemblied before calling this routine
2398: Must set the matrix option MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS); before
2399: calling this routine.
2401: When this is called multiple times it overwrites the previous set of stored values
2402: and does not allocated additional space.
2404: .seealso: MatRetrieveValues()
2406: @*/
2407: PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues(Mat mat)
2408: {
2409: PetscErrorCode ierr,(*f)(Mat);
2413: if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2414: if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2416: PetscObjectQueryFunction((PetscObject)mat,"MatStoreValues_C",(void (**)(void))&f);
2417: if (f) {
2418: (*f)(mat);
2419: } else {
2420: SETERRQ(PETSC_ERR_SUP,"Wrong type of matrix to store values");
2421: }
2422: return(0);
2423: }
2428: PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues_SeqAIJ(Mat mat)
2429: {
2430: Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
2432: PetscInt nz = aij->i[mat->m];
2435: if (aij->nonew != 1) {
2436: SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
2437: }
2438: if (!aij->saved_values) {
2439: SETERRQ(PETSC_ERR_ORDER,"Must call MatStoreValues(A);first");
2440: }
2441: /* copy values over */
2442: PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));
2443: return(0);
2444: }
2449: /*@
2450: MatRetrieveValues - Retrieves the copy of the matrix values; this allows, for
2451: example, reuse of the linear part of a Jacobian, while recomputing the
2452: nonlinear portion.
2454: Collect on Mat
2456: Input Parameters:
2457: . mat - the matrix (currently on AIJ matrices support this option)
2459: Level: advanced
2461: .seealso: MatStoreValues()
2463: @*/
2464: PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues(Mat mat)
2465: {
2466: PetscErrorCode ierr,(*f)(Mat);
2470: if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2471: if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2473: PetscObjectQueryFunction((PetscObject)mat,"MatRetrieveValues_C",(void (**)(void))&f);
2474: if (f) {
2475: (*f)(mat);
2476: } else {
2477: SETERRQ(PETSC_ERR_SUP,"Wrong type of matrix to retrieve values");
2478: }
2479: return(0);
2480: }
2483: /* --------------------------------------------------------------------------------*/
2486: /*@C
2487: MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format
2488: (the default parallel PETSc format). For good matrix assembly performance
2489: the user should preallocate the matrix storage by setting the parameter nz
2490: (or the array nnz). By setting these parameters accurately, performance
2491: during matrix assembly can be increased by more than a factor of 50.
2493: Collective on MPI_Comm
2495: Input Parameters:
2496: + comm - MPI communicator, set to PETSC_COMM_SELF
2497: . m - number of rows
2498: . n - number of columns
2499: . nz - number of nonzeros per row (same for all rows)
2500: - nnz - array containing the number of nonzeros in the various rows
2501: (possibly different for each row) or PETSC_NULL
2503: Output Parameter:
2504: . A - the matrix
2506: Notes:
2507: If nnz is given then nz is ignored
2509: The AIJ format (also called the Yale sparse matrix format or
2510: compressed row storage), is fully compatible with standard Fortran 77
2511: storage. That is, the stored row and column indices can begin at
2512: either one (as in Fortran) or zero. See the users' manual for details.
2514: Specify the preallocated storage with either nz or nnz (not both).
2515: Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
2516: allocation. For large problems you MUST preallocate memory or you
2517: will get TERRIBLE performance, see the users' manual chapter on matrices.
2519: By default, this format uses inodes (identical nodes) when possible, to
2520: improve numerical efficiency of matrix-vector products and solves. We
2521: search for consecutive rows with the same nonzero structure, thereby
2522: reusing matrix information to achieve increased efficiency.
2524: Options Database Keys:
2525: + -mat_no_inode - Do not use inodes
2526: . -mat_inode_limit <limit> - Sets inode limit (max limit=5)
2527: - -mat_aij_oneindex - Internally use indexing starting at 1
2528: rather than 0. Note that when calling MatSetValues(),
2529: the user still MUST index entries starting at 0!
2531: Level: intermediate
2533: .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues(), MatSeqAIJSetColumnIndices(), MatCreateSeqAIJWithArrays()
2535: @*/
2536: PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJ(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A)
2537: {
2541: MatCreate(comm,A);
2542: MatSetSizes(*A,m,n,m,n);
2543: MatSetType(*A,MATSEQAIJ);
2544: MatSeqAIJSetPreallocation_SeqAIJ(*A,nz,(PetscInt*)nnz);
2545: return(0);
2546: }
2550: /*@C
2551: MatSeqAIJSetPreallocation - For good matrix assembly performance
2552: the user should preallocate the matrix storage by setting the parameter nz
2553: (or the array nnz). By setting these parameters accurately, performance
2554: during matrix assembly can be increased by more than a factor of 50.
2556: Collective on MPI_Comm
2558: Input Parameters:
2559: + B - The matrix
2560: . nz - number of nonzeros per row (same for all rows)
2561: - nnz - array containing the number of nonzeros in the various rows
2562: (possibly different for each row) or PETSC_NULL
2564: Notes:
2565: If nnz is given then nz is ignored
2567: The AIJ format (also called the Yale sparse matrix format or
2568: compressed row storage), is fully compatible with standard Fortran 77
2569: storage. That is, the stored row and column indices can begin at
2570: either one (as in Fortran) or zero. See the users' manual for details.
2572: Specify the preallocated storage with either nz or nnz (not both).
2573: Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
2574: allocation. For large problems you MUST preallocate memory or you
2575: will get TERRIBLE performance, see the users' manual chapter on matrices.
2577: Developers: Use nz of MAT_SKIP_ALLOCATION to not allocate any space for the matrix
2578: entries or columns indices
2580: By default, this format uses inodes (identical nodes) when possible, to
2581: improve numerical efficiency of matrix-vector products and solves. We
2582: search for consecutive rows with the same nonzero structure, thereby
2583: reusing matrix information to achieve increased efficiency.
2585: Options Database Keys:
2586: + -mat_no_inode - Do not use inodes
2587: . -mat_inode_limit <limit> - Sets inode limit (max limit=5)
2588: - -mat_aij_oneindex - Internally use indexing starting at 1
2589: rather than 0. Note that when calling MatSetValues(),
2590: the user still MUST index entries starting at 0!
2592: Level: intermediate
2594: .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues(), MatSeqAIJSetColumnIndices(), MatCreateSeqAIJWithArrays()
2596: @*/
2597: PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocation(Mat B,PetscInt nz,const PetscInt nnz[])
2598: {
2599: PetscErrorCode ierr,(*f)(Mat,PetscInt,const PetscInt[]);
2602: PetscObjectQueryFunction((PetscObject)B,"MatSeqAIJSetPreallocation_C",(void (**)(void))&f);
2603: if (f) {
2604: (*f)(B,nz,nnz);
2605: }
2606: return(0);
2607: }
2612: PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocation_SeqAIJ(Mat B,PetscInt nz,PetscInt *nnz)
2613: {
2614: Mat_SeqAIJ *b;
2615: PetscTruth skipallocation = PETSC_FALSE;
2617: PetscInt i;
2620:
2621: if (nz == MAT_SKIP_ALLOCATION) {
2622: skipallocation = PETSC_TRUE;
2623: nz = 0;
2624: }
2626: if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
2627: if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz);
2628: if (nnz) {
2629: for (i=0; i<B->m; i++) {
2630: if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]);
2631: if (nnz[i] > B->n) SETERRQ3(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be greater than row length: local row %d value %d rowlength %d",i,nnz[i],B->n);
2632: }
2633: }
2635: B->preallocated = PETSC_TRUE;
2636: b = (Mat_SeqAIJ*)B->data;
2638: if (!skipallocation) {
2639: PetscMalloc2(B->m,PetscInt,&b->imax,B->m,PetscInt,&b->ilen);
2640: if (!nnz) {
2641: if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 10;
2642: else if (nz <= 0) nz = 1;
2643: for (i=0; i<B->m; i++) b->imax[i] = nz;
2644: nz = nz*B->m;
2645: } else {
2646: nz = 0;
2647: for (i=0; i<B->m; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
2648: }
2650: /* b->ilen will count nonzeros in each row so far. */
2651: for (i=0; i<B->m; i++) { b->ilen[i] = 0;}
2653: /* allocate the matrix space */
2654: PetscMalloc3(nz,PetscScalar,&b->a,nz,PetscInt,&b->j,B->m+1,PetscInt,&b->i);
2655: b->i[0] = 0;
2656: for (i=1; i<B->m+1; i++) {
2657: b->i[i] = b->i[i-1] + b->imax[i-1];
2658: }
2659: b->singlemalloc = PETSC_TRUE;
2660: b->freedata = PETSC_TRUE;
2661: } else {
2662: b->freedata = PETSC_FALSE;
2663: }
2665: b->nz = 0;
2666: b->maxnz = nz;
2667: B->info.nz_unneeded = (double)b->maxnz;
2668: return(0);
2669: }
2672: /*MC
2673: MATSEQAIJ - MATSEQAIJ = "seqaij" - A matrix type to be used for sequential sparse matrices,
2674: based on compressed sparse row format.
2676: Options Database Keys:
2677: . -mat_type seqaij - sets the matrix type to "seqaij" during a call to MatSetFromOptions()
2679: Level: beginner
2681: .seealso: MatCreateSeqAIJ(), MatSetFromOptions(), MatSetType(), MatCreate(), MatType
2682: M*/
2687: PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqAIJ(Mat B)
2688: {
2689: Mat_SeqAIJ *b;
2691: PetscMPIInt size;
2694: MPI_Comm_size(B->comm,&size);
2695: if (size > 1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Comm must be of size 1");
2697: B->m = B->M = PetscMax(B->m,B->M);
2698: B->n = B->N = PetscMax(B->n,B->N);
2700: PetscNew(Mat_SeqAIJ,&b);
2701: B->data = (void*)b;
2702: PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
2703: B->factor = 0;
2704: B->mapping = 0;
2705: b->row = 0;
2706: b->col = 0;
2707: b->icol = 0;
2708: b->reallocs = 0;
2709:
2710: PetscMapCreateMPI(B->comm,B->m,B->m,&B->rmap);
2711: PetscMapCreateMPI(B->comm,B->n,B->n,&B->cmap);
2713: b->sorted = PETSC_FALSE;
2714: b->ignorezeroentries = PETSC_FALSE;
2715: b->roworiented = PETSC_TRUE;
2716: b->nonew = 0;
2717: b->diag = 0;
2718: b->solve_work = 0;
2719: B->spptr = 0;
2720: b->saved_values = 0;
2721: b->idiag = 0;
2722: b->ssor = 0;
2723: b->keepzeroedrows = PETSC_FALSE;
2724: b->xtoy = 0;
2725: b->XtoY = 0;
2726: b->compressedrow.use = PETSC_FALSE;
2727: b->compressedrow.nrows = B->m;
2728: b->compressedrow.i = PETSC_NULL;
2729: b->compressedrow.rindex = PETSC_NULL;
2730: b->compressedrow.checked = PETSC_FALSE;
2731: B->same_nonzero = PETSC_FALSE;
2733: PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);
2735: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqAIJSetColumnIndices_C",
2736: "MatSeqAIJSetColumnIndices_SeqAIJ",
2737: MatSeqAIJSetColumnIndices_SeqAIJ);
2738: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
2739: "MatStoreValues_SeqAIJ",
2740: MatStoreValues_SeqAIJ);
2741: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
2742: "MatRetrieveValues_SeqAIJ",
2743: MatRetrieveValues_SeqAIJ);
2744: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_seqsbaij_C",
2745: "MatConvert_SeqAIJ_SeqSBAIJ",
2746: MatConvert_SeqAIJ_SeqSBAIJ);
2747: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_seqbaij_C",
2748: "MatConvert_SeqAIJ_SeqBAIJ",
2749: MatConvert_SeqAIJ_SeqBAIJ);
2750: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatIsTranspose_C",
2751: "MatIsTranspose_SeqAIJ",
2752: MatIsTranspose_SeqAIJ);
2753: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqAIJSetPreallocation_C",
2754: "MatSeqAIJSetPreallocation_SeqAIJ",
2755: MatSeqAIJSetPreallocation_SeqAIJ);
2756: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatReorderForNonzeroDiagonal_C",
2757: "MatReorderForNonzeroDiagonal_SeqAIJ",
2758: MatReorderForNonzeroDiagonal_SeqAIJ);
2759: MatCreate_Inode(B);
2760: return(0);
2761: }
2766: PetscErrorCode MatDuplicate_SeqAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
2767: {
2768: Mat C;
2769: Mat_SeqAIJ *c,*a = (Mat_SeqAIJ*)A->data;
2771: PetscInt i,m = A->m;
2774: *B = 0;
2775: MatCreate(A->comm,&C);
2776: MatSetSizes(C,A->m,A->n,A->m,A->n);
2777: MatSetType(C,A->type_name);
2778: PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));
2779:
2780: c = (Mat_SeqAIJ*)C->data;
2782: C->factor = A->factor;
2784: c->row = 0;
2785: c->col = 0;
2786: c->icol = 0;
2787: c->reallocs = 0;
2789: C->assembled = PETSC_TRUE;
2790:
2791: PetscMalloc2(m,PetscInt,&c->imax,m,PetscInt,&c->ilen);
2792: for (i=0; i<m; i++) {
2793: c->imax[i] = a->imax[i];
2794: c->ilen[i] = a->ilen[i];
2795: }
2797: /* allocate the matrix space */
2798: PetscMalloc3(a->i[m],PetscScalar,&c->a,a->i[m],PetscInt,&c->j,m+1,PetscInt,&c->i);
2799: c->singlemalloc = PETSC_TRUE;
2800: PetscMemcpy(c->i,a->i,(m+1)*sizeof(PetscInt));
2801: if (m > 0) {
2802: PetscMemcpy(c->j,a->j,(a->i[m])*sizeof(PetscInt));
2803: if (cpvalues == MAT_COPY_VALUES) {
2804: PetscMemcpy(c->a,a->a,(a->i[m])*sizeof(PetscScalar));
2805: } else {
2806: PetscMemzero(c->a,(a->i[m])*sizeof(PetscScalar));
2807: }
2808: }
2810: c->sorted = a->sorted;
2811: c->ignorezeroentries = a->ignorezeroentries;
2812: c->roworiented = a->roworiented;
2813: c->nonew = a->nonew;
2814: if (a->diag) {
2815: PetscMalloc((m+1)*sizeof(PetscInt),&c->diag);
2816: PetscLogObjectMemory(C,(m+1)*sizeof(PetscInt));
2817: for (i=0; i<m; i++) {
2818: c->diag[i] = a->diag[i];
2819: }
2820: } else c->diag = 0;
2821: c->solve_work = 0;
2822: c->saved_values = 0;
2823: c->idiag = 0;
2824: c->ssor = 0;
2825: c->keepzeroedrows = a->keepzeroedrows;
2826: c->freedata = PETSC_TRUE;
2827: c->xtoy = 0;
2828: c->XtoY = 0;
2830: c->nz = a->nz;
2831: c->maxnz = a->maxnz;
2832: C->preallocated = PETSC_TRUE;
2834: c->compressedrow.use = a->compressedrow.use;
2835: c->compressedrow.nrows = a->compressedrow.nrows;
2836: c->compressedrow.checked = a->compressedrow.checked;
2837: if ( a->compressedrow.checked && a->compressedrow.use){
2838: i = a->compressedrow.nrows;
2839: PetscMalloc((2*i+1)*sizeof(PetscInt),&c->compressedrow.i);
2840: c->compressedrow.rindex = c->compressedrow.i + i + 1;
2841: PetscMemcpy(c->compressedrow.i,a->compressedrow.i,(i+1)*sizeof(PetscInt));
2842: PetscMemcpy(c->compressedrow.rindex,a->compressedrow.rindex,i*sizeof(PetscInt));
2843: } else {
2844: c->compressedrow.use = PETSC_FALSE;
2845: c->compressedrow.i = PETSC_NULL;
2846: c->compressedrow.rindex = PETSC_NULL;
2847: }
2848: C->same_nonzero = A->same_nonzero;
2850: MatDuplicate_Inode(A,cpvalues,&C);
2852: *B = C;
2853: PetscFListDuplicate(A->qlist,&C->qlist);
2854: return(0);
2855: }
2859: PetscErrorCode MatLoad_SeqAIJ(PetscViewer viewer, MatType type,Mat *A)
2860: {
2861: Mat_SeqAIJ *a;
2862: Mat B;
2864: PetscInt i,sum,nz,header[4],*rowlengths = 0,M,N;
2865: int fd;
2866: PetscMPIInt size;
2867: MPI_Comm comm;
2868:
2870: PetscObjectGetComm((PetscObject)viewer,&comm);
2871: MPI_Comm_size(comm,&size);
2872: if (size > 1) SETERRQ(PETSC_ERR_ARG_SIZ,"view must have one processor");
2873: PetscViewerBinaryGetDescriptor(viewer,&fd);
2874: PetscBinaryRead(fd,header,4,PETSC_INT);
2875: if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object in file");
2876: M = header[1]; N = header[2]; nz = header[3];
2878: if (nz < 0) {
2879: SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format on disk,cannot load as SeqAIJ");
2880: }
2882: /* read in row lengths */
2883: PetscMalloc(M*sizeof(PetscInt),&rowlengths);
2884: PetscBinaryRead(fd,rowlengths,M,PETSC_INT);
2886: /* check if sum of rowlengths is same as nz */
2887: for (i=0,sum=0; i< M; i++) sum +=rowlengths[i];
2888: if (sum != nz) SETERRQ2(PETSC_ERR_FILE_READ,"Inconsistant matrix data in file. no-nonzeros = %d, sum-row-lengths = %d\n",nz,sum);
2890: /* create our matrix */
2891: MatCreate(comm,&B);
2892: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,M,N);
2893: MatSetType(B,type);
2894: MatSeqAIJSetPreallocation_SeqAIJ(B,0,rowlengths);
2895: a = (Mat_SeqAIJ*)B->data;
2897: /* read in column indices and adjust for Fortran indexing*/
2898: PetscBinaryRead(fd,a->j,nz,PETSC_INT);
2900: /* read in nonzero values */
2901: PetscBinaryRead(fd,a->a,nz,PETSC_SCALAR);
2903: /* set matrix "i" values */
2904: a->i[0] = 0;
2905: for (i=1; i<= M; i++) {
2906: a->i[i] = a->i[i-1] + rowlengths[i-1];
2907: a->ilen[i-1] = rowlengths[i-1];
2908: }
2909: PetscFree(rowlengths);
2911: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
2912: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
2913: *A = B;
2914: return(0);
2915: }
2919: PetscErrorCode MatEqual_SeqAIJ(Mat A,Mat B,PetscTruth* flg)
2920: {
2921: Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data,*b = (Mat_SeqAIJ *)B->data;
2925: /* If the matrix dimensions are not equal,or no of nonzeros */
2926: if ((A->m != B->m) || (A->n != B->n) ||(a->nz != b->nz)) {
2927: *flg = PETSC_FALSE;
2928: return(0);
2929: }
2930:
2931: /* if the a->i are the same */
2932: PetscMemcmp(a->i,b->i,(A->m+1)*sizeof(PetscInt),flg);
2933: if (!*flg) return(0);
2934:
2935: /* if a->j are the same */
2936: PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),flg);
2937: if (!*flg) return(0);
2938:
2939: /* if a->a are the same */
2940: PetscMemcmp(a->a,b->a,(a->nz)*sizeof(PetscScalar),flg);
2942: return(0);
2943:
2944: }
2948: /*@C
2949: MatCreateSeqAIJWithArrays - Creates an sequential AIJ matrix using matrix elements (in CSR format)
2950: provided by the user.
2952: Coolective on MPI_Comm
2954: Input Parameters:
2955: + comm - must be an MPI communicator of size 1
2956: . m - number of rows
2957: . n - number of columns
2958: . i - row indices
2959: . j - column indices
2960: - a - matrix values
2962: Output Parameter:
2963: . mat - the matrix
2965: Level: intermediate
2967: Notes:
2968: The i, j, and a arrays are not copied by this routine, the user must free these arrays
2969: once the matrix is destroyed
2971: You cannot set new nonzero locations into this matrix, that will generate an error.
2973: The i and j indices are 0 based
2975: .seealso: MatCreate(), MatCreateMPIAIJ(), MatCreateSeqAIJ()
2977: @*/
2978: PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJWithArrays(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt* i,PetscInt*j,PetscScalar *a,Mat *mat)
2979: {
2981: PetscInt ii;
2982: Mat_SeqAIJ *aij;
2985: if (i[0]) {
2986: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"i (row indices) must start with 0");
2987: }
2988: MatCreate(comm,mat);
2989: MatSetSizes(*mat,m,n,m,n);
2990: MatSetType(*mat,MATSEQAIJ);
2991: MatSeqAIJSetPreallocation_SeqAIJ(*mat,MAT_SKIP_ALLOCATION,0);
2992: aij = (Mat_SeqAIJ*)(*mat)->data;
2993: PetscMalloc2(m,PetscInt,&aij->imax,m,PetscInt,&aij->ilen);
2995: aij->i = i;
2996: aij->j = j;
2997: aij->a = a;
2998: aij->singlemalloc = PETSC_FALSE;
2999: aij->nonew = -1; /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/
3000: aij->freedata = PETSC_FALSE;
3002: for (ii=0; ii<m; ii++) {
3003: aij->ilen[ii] = aij->imax[ii] = i[ii+1] - i[ii];
3004: #if defined(PETSC_USE_DEBUG)
3005: 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]);
3006: #endif
3007: }
3008: #if defined(PETSC_USE_DEBUG)
3009: for (ii=0; ii<aij->i[m]; ii++) {
3010: if (j[ii] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Negative column index at location = %d index = %d",ii,j[ii]);
3011: if (j[ii] > n - 1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column index to large at location = %d index = %d",ii,j[ii]);
3012: }
3013: #endif
3015: MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
3016: MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
3017: return(0);
3018: }
3022: PetscErrorCode MatSetColoring_SeqAIJ(Mat A,ISColoring coloring)
3023: {
3025: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
3028: if (coloring->ctype == IS_COLORING_LOCAL) {
3029: ISColoringReference(coloring);
3030: a->coloring = coloring;
3031: } else if (coloring->ctype == IS_COLORING_GHOSTED) {
3032: PetscInt i,*larray;
3033: ISColoring ocoloring;
3034: ISColoringValue *colors;
3036: /* set coloring for diagonal portion */
3037: PetscMalloc((A->n+1)*sizeof(PetscInt),&larray);
3038: for (i=0; i<A->n; i++) {
3039: larray[i] = i;
3040: }
3041: ISGlobalToLocalMappingApply(A->mapping,IS_GTOLM_MASK,A->n,larray,PETSC_NULL,larray);
3042: PetscMalloc((A->n+1)*sizeof(ISColoringValue),&colors);
3043: for (i=0; i<A->n; i++) {
3044: colors[i] = coloring->colors[larray[i]];
3045: }
3046: PetscFree(larray);
3047: ISColoringCreate(PETSC_COMM_SELF,A->n,colors,&ocoloring);
3048: a->coloring = ocoloring;
3049: }
3050: return(0);
3051: }
3053: #if defined(PETSC_HAVE_ADIC)
3055: #include "adic/ad_utils.h"
3060: PetscErrorCode MatSetValuesAdic_SeqAIJ(Mat A,void *advalues)
3061: {
3062: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
3063: PetscInt m = A->m,*ii = a->i,*jj = a->j,nz,i,j,nlen;
3064: PetscScalar *v = a->a,*values = ((PetscScalar*)advalues)+1;
3065: ISColoringValue *color;
3068: if (!a->coloring) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Coloring not set for matrix");
3069: nlen = PetscADGetDerivTypeSize()/sizeof(PetscScalar);
3070: color = a->coloring->colors;
3071: /* loop over rows */
3072: for (i=0; i<m; i++) {
3073: nz = ii[i+1] - ii[i];
3074: /* loop over columns putting computed value into matrix */
3075: for (j=0; j<nz; j++) {
3076: *v++ = values[color[*jj++]];
3077: }
3078: values += nlen; /* jump to next row of derivatives */
3079: }
3080: return(0);
3081: }
3082: #endif
3086: PetscErrorCode MatSetValuesAdifor_SeqAIJ(Mat A,PetscInt nl,void *advalues)
3087: {
3088: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
3089: PetscInt m = A->m,*ii = a->i,*jj = a->j,nz,i,j;
3090: PetscScalar *v = a->a,*values = (PetscScalar *)advalues;
3091: ISColoringValue *color;
3094: if (!a->coloring) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Coloring not set for matrix");
3095: color = a->coloring->colors;
3096: /* loop over rows */
3097: for (i=0; i<m; i++) {
3098: nz = ii[i+1] - ii[i];
3099: /* loop over columns putting computed value into matrix */
3100: for (j=0; j<nz; j++) {
3101: *v++ = values[color[*jj++]];
3102: }
3103: values += nl; /* jump to next row of derivatives */
3104: }
3105: return(0);
3106: }