Actual source code: aij.c
1: /*$Id: aij.c,v 1.385 2001/09/07 20:09:22 bsmith Exp $*/
2: /*
3: Defines the basic matrix operations for the AIJ (compressed row)
4: matrix storage format.
5: */
7: #include src/mat/impls/aij/seq/aij.h
8: #include src/vec/vecimpl.h
9: #include src/inline/spops.h
10: #include src/inline/dot.h
11: #include petscbt.h
14: EXTERN int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**);
16: int MatGetRowIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *m,int **ia,int **ja,PetscTruth *done)
17: {
18: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
19: int ierr,i,ishift;
20:
22: *m = A->m;
23: if (!ia) return(0);
24: ishift = a->indexshift;
25: if (symmetric && !A->structurally_symmetric) {
26: MatToSymmetricIJ_SeqAIJ(A->m,a->i,a->j,ishift,oshift,ia,ja);
27: } else if (oshift == 0 && ishift == -1) {
28: int nz = a->i[A->m] - 1;
29: /* malloc space and subtract 1 from i and j indices */
30: PetscMalloc((A->m+1)*sizeof(int),ia);
31: PetscMalloc((nz+1)*sizeof(int),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 if (oshift == 1 && ishift == 0) {
35: int nz = a->i[A->m];
36: /* malloc space and add 1 to i and j indices */
37: PetscMalloc((A->m+1)*sizeof(int),ia);
38: PetscMalloc((nz+1)*sizeof(int),ja);
39: for (i=0; i<nz; i++) (*ja)[i] = a->j[i] + 1;
40: for (i=0; i<A->m+1; i++) (*ia)[i] = a->i[i] + 1;
41: } else {
42: *ia = a->i; *ja = a->j;
43: }
44: return(0);
45: }
47: int MatRestoreRowIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *n,int **ia,int **ja,PetscTruth *done)
48: {
49: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
50: int ishift = a->indexshift,ierr;
51:
53: if (!ia) return(0);
54: if ((symmetric && !A->structurally_symmetric) || (oshift == 0 && ishift == -1) || (oshift == 1 && ishift == 0)) {
55: PetscFree(*ia);
56: PetscFree(*ja);
57: }
58: return(0);
59: }
61: int MatGetColumnIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done)
62: {
63: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
64: int ierr,i,ishift = a->indexshift,*collengths,*cia,*cja,n = A->n,m = A->m;
65: int nz = a->i[m]+ishift,row,*jj,mr,col;
66:
68: *nn = A->n;
69: if (!ia) return(0);
70: if (symmetric) {
71: MatToSymmetricIJ_SeqAIJ(A->m,a->i,a->j,ishift,oshift,ia,ja);
72: } else {
73: PetscMalloc((n+1)*sizeof(int),&collengths);
74: PetscMemzero(collengths,n*sizeof(int));
75: PetscMalloc((n+1)*sizeof(int),&cia);
76: PetscMalloc((nz+1)*sizeof(int),&cja);
77: jj = a->j;
78: for (i=0; i<nz; i++) {
79: collengths[jj[i] + ishift]++;
80: }
81: cia[0] = oshift;
82: for (i=0; i<n; i++) {
83: cia[i+1] = cia[i] + collengths[i];
84: }
85: PetscMemzero(collengths,n*sizeof(int));
86: jj = a->j;
87: for (row=0; row<m; row++) {
88: mr = a->i[row+1] - a->i[row];
89: for (i=0; i<mr; i++) {
90: col = *jj++ + ishift;
91: cja[cia[col] + collengths[col]++ - oshift] = row + oshift;
92: }
93: }
94: PetscFree(collengths);
95: *ia = cia; *ja = cja;
96: }
97: return(0);
98: }
100: int MatRestoreColumnIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *n,int **ia,int **ja,PetscTruth *done)
101: {
105: if (!ia) return(0);
107: PetscFree(*ia);
108: PetscFree(*ja);
109:
110: return(0);
111: }
113: #define CHUNKSIZE 15
115: int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,PetscScalar *v,InsertMode is)
116: {
117: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
118: int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted = a->sorted;
119: int *imax = a->imax,*ai = a->i,*ailen = a->ilen;
120: int *aj = a->j,nonew = a->nonew,shift = a->indexshift,ierr;
121: PetscScalar *ap,value,*aa = a->a;
122: PetscTruth ignorezeroentries = ((a->ignorezeroentries && is == ADD_VALUES) ? PETSC_TRUE:PETSC_FALSE);
123: PetscTruth roworiented = a->roworiented;
126: for (k=0; k<m; k++) { /* loop over added rows */
127: row = im[k];
128: if (row < 0) continue;
129: #if defined(PETSC_USE_BOPT_g)
130: if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,A->m);
131: #endif
132: rp = aj + ai[row] + shift; ap = aa + ai[row] + shift;
133: rmax = imax[row]; nrow = ailen[row];
134: low = 0;
135: for (l=0; l<n; l++) { /* loop over added columns */
136: if (in[l] < 0) continue;
137: #if defined(PETSC_USE_BOPT_g)
138: if (in[l] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],A->n);
139: #endif
140: col = in[l] - shift;
141: if (roworiented) {
142: value = v[l + k*n];
143: } else {
144: value = v[k + l*m];
145: }
146: if (value == 0.0 && ignorezeroentries) continue;
148: if (!sorted) low = 0; high = nrow;
149: while (high-low > 5) {
150: t = (low+high)/2;
151: if (rp[t] > col) high = t;
152: else low = t;
153: }
154: for (i=low; i<high; i++) {
155: if (rp[i] > col) break;
156: if (rp[i] == col) {
157: if (is == ADD_VALUES) ap[i] += value;
158: else ap[i] = value;
159: goto noinsert;
160: }
161: }
162: if (nonew == 1) goto noinsert;
163: else if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero at (%d,%d) in the matrix",row,col);
164: if (nrow >= rmax) {
165: /* there is no extra room in row, therefore enlarge */
166: int new_nz = ai[A->m] + CHUNKSIZE,len,*new_i,*new_j;
167: PetscScalar *new_a;
169: if (nonew == -2) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero at (%d,%d) in the matrix requiring new malloc()",row,col);
171: /* malloc new storage space */
172: len = new_nz*(sizeof(int)+sizeof(PetscScalar))+(A->m+1)*sizeof(int);
173: ierr = PetscMalloc(len,&new_a);
174: new_j = (int*)(new_a + new_nz);
175: new_i = new_j + new_nz;
177: /* copy over old data into new slots */
178: for (ii=0; ii<row+1; ii++) {new_i[ii] = ai[ii];}
179: for (ii=row+1; ii<A->m+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;}
180: PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int));
181: len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift);
182: PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow,len*sizeof(int));
183: PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(PetscScalar));
184: PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow,len*sizeof(PetscScalar));
185: /* free up old matrix storage */
186: PetscFree(a->a);
187: if (!a->singlemalloc) {
188: PetscFree(a->i);
189: PetscFree(a->j);
190: }
191: aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
192: a->singlemalloc = PETSC_TRUE;
194: rp = aj + ai[row] + shift; ap = aa + ai[row] + shift;
195: rmax = imax[row] = imax[row] + CHUNKSIZE;
196: PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(PetscScalar)));
197: a->maxnz += CHUNKSIZE;
198: a->reallocs++;
199: }
200: N = nrow++ - 1; a->nz++;
201: /* shift up all the later entries in this row */
202: for (ii=N; ii>=i; ii--) {
203: rp[ii+1] = rp[ii];
204: ap[ii+1] = ap[ii];
205: }
206: rp[i] = col;
207: ap[i] = value;
208: noinsert:;
209: low = i + 1;
210: }
211: ailen[row] = nrow;
212: }
213: return(0);
214: }
216: int MatGetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,PetscScalar *v)
217: {
218: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
219: int *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
220: int *ai = a->i,*ailen = a->ilen,shift = a->indexshift;
221: PetscScalar *ap,*aa = a->a,zero = 0.0;
224: for (k=0; k<m; k++) { /* loop over rows */
225: row = im[k];
226: if (row < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %d",row);
227: if (row >= A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: %d",row);
228: rp = aj + ai[row] + shift; ap = aa + ai[row] + shift;
229: nrow = ailen[row];
230: for (l=0; l<n; l++) { /* loop over columns */
231: if (in[l] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %d",in[l]);
232: if (in[l] >= A->n) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: %d",in[l]);
233: col = in[l] - shift;
234: high = nrow; low = 0; /* assume unsorted */
235: while (high-low > 5) {
236: t = (low+high)/2;
237: if (rp[t] > col) high = t;
238: else low = t;
239: }
240: for (i=low; i<high; i++) {
241: if (rp[i] > col) break;
242: if (rp[i] == col) {
243: *v++ = ap[i];
244: goto finished;
245: }
246: }
247: *v++ = zero;
248: finished:;
249: }
250: }
251: return(0);
252: }
255: int MatView_SeqAIJ_Binary(Mat A,PetscViewer viewer)
256: {
257: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
258: int i,fd,*col_lens,ierr;
261: PetscViewerBinaryGetDescriptor(viewer,&fd);
262: PetscMalloc((4+A->m)*sizeof(int),&col_lens);
263: col_lens[0] = MAT_FILE_COOKIE;
264: col_lens[1] = A->m;
265: col_lens[2] = A->n;
266: col_lens[3] = a->nz;
268: /* store lengths of each row and write (including header) to file */
269: for (i=0; i<A->m; i++) {
270: col_lens[4+i] = a->i[i+1] - a->i[i];
271: }
272: PetscBinaryWrite(fd,col_lens,4+A->m,PETSC_INT,1);
273: PetscFree(col_lens);
275: /* store column indices (zero start index) */
276: if (a->indexshift) {
277: for (i=0; i<a->nz; i++) a->j[i]--;
278: }
279: PetscBinaryWrite(fd,a->j,a->nz,PETSC_INT,0);
280: if (a->indexshift) {
281: for (i=0; i<a->nz; i++) a->j[i]++;
282: }
284: /* store nonzero values */
285: PetscBinaryWrite(fd,a->a,a->nz,PETSC_SCALAR,0);
286: return(0);
287: }
289: extern int MatMPIAIJFactorInfo_SuperLu(Mat,PetscViewer);
290: extern int MatFactorInfo_Spooles(Mat,PetscViewer);
292: int MatView_SeqAIJ_ASCII(Mat A,PetscViewer viewer)
293: {
294: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
295: int ierr,i,j,m = A->m,shift = a->indexshift;
296: char *name;
297: PetscViewerFormat format;
300: PetscObjectGetName((PetscObject)A,&name);
301: PetscViewerGetFormat(viewer,&format);
302: if (format == PETSC_VIEWER_ASCII_INFO_LONG || format == PETSC_VIEWER_ASCII_INFO) {
303: if (a->inode.size) {
304: PetscViewerASCIIPrintf(viewer,"using I-node routines: found %d nodes, limit used is %dn",a->inode.node_count,a->inode.limit);
305: } else {
306: PetscViewerASCIIPrintf(viewer,"not using I-node routinesn");
307: }
308: } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
309: int nofinalvalue = 0;
310: if ((a->i[m] == a->i[m-1]) || (a->j[a->nz-1] != A->n-!shift)) {
311: nofinalvalue = 1;
312: }
313: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
314: PetscViewerASCIIPrintf(viewer,"%% Size = %d %d n",m,A->n);
315: PetscViewerASCIIPrintf(viewer,"%% Nonzeros = %d n",a->nz);
316: PetscViewerASCIIPrintf(viewer,"zzz = zeros(%d,3);n",a->nz+nofinalvalue);
317: PetscViewerASCIIPrintf(viewer,"zzz = [n");
319: for (i=0; i<m; i++) {
320: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
321: #if defined(PETSC_USE_COMPLEX)
322: PetscViewerASCIIPrintf(viewer,"%d %d %18.16e + %18.16ei n",i+1,a->j[j]+!shift,PetscRealPart(a->a[j]),PetscImaginaryPart(a->a[j]));
323: #else
324: PetscViewerASCIIPrintf(viewer,"%d %d %18.16en",i+1,a->j[j]+!shift,a->a[j]);
325: #endif
326: }
327: }
328: if (nofinalvalue) {
329: PetscViewerASCIIPrintf(viewer,"%d %d %18.16en",m,A->n,0.0);
330: }
331: PetscViewerASCIIPrintf(viewer,"];n %s = spconvert(zzz);n",name);
332: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
333: } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
334: #if defined(PETSC_HAVE_SUPERLUDIST) && !defined(PETSC_USE_SINGLE) && !defined(PETSC_USE_COMPLEX)
335: MatMPIAIJFactorInfo_SuperLu(A,viewer);
336: #endif
337: #if defined(PETSC_HAVE_SPOOLES) && !defined(PETSC_USE_SINGLE) && !defined(PETSC_USE_COMPLEX)
338: MatFactorInfo_Spooles(A,viewer);
339: #endif
340: return(0);
341: } else if (format == PETSC_VIEWER_ASCII_COMMON) {
342: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
343: for (i=0; i<m; i++) {
344: PetscViewerASCIIPrintf(viewer,"row %d:",i);
345: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
346: #if defined(PETSC_USE_COMPLEX)
347: if (PetscImaginaryPart(a->a[j]) > 0.0 && PetscRealPart(a->a[j]) != 0.0) {
348: PetscViewerASCIIPrintf(viewer," (%d, %g + %g i)",a->j[j]+shift,PetscRealPart(a->a[j]),PetscImaginaryPart(a->a[j]));
349: } else if (PetscImaginaryPart(a->a[j]) < 0.0 && PetscRealPart(a->a[j]) != 0.0) {
350: PetscViewerASCIIPrintf(viewer," (%d, %g - %g i)",a->j[j]+shift,PetscRealPart(a->a[j]),-PetscImaginaryPart(a->a[j]));
351: } else if (PetscRealPart(a->a[j]) != 0.0) {
352: PetscViewerASCIIPrintf(viewer," (%d, %g) ",a->j[j]+shift,PetscRealPart(a->a[j]));
353: }
354: #else
355: if (a->a[j] != 0.0) {PetscViewerASCIIPrintf(viewer," (%d, %g) ",a->j[j]+shift,a->a[j]);}
356: #endif
357: }
358: PetscViewerASCIIPrintf(viewer,"n");
359: }
360: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
361: } else if (format == PETSC_VIEWER_ASCII_SYMMODU) {
362: int nzd=0,fshift=1,*sptr;
363: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
364: PetscMalloc((m+1)*sizeof(int),&sptr);
365: for (i=0; i<m; i++) {
366: sptr[i] = nzd+1;
367: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
368: if (a->j[j] >= i) {
369: #if defined(PETSC_USE_COMPLEX)
370: if (PetscImaginaryPart(a->a[j]) != 0.0 || PetscRealPart(a->a[j]) != 0.0) nzd++;
371: #else
372: if (a->a[j] != 0.0) nzd++;
373: #endif
374: }
375: }
376: }
377: sptr[m] = nzd+1;
378: PetscViewerASCIIPrintf(viewer," %d %dnn",m,nzd);
379: for (i=0; i<m+1; i+=6) {
380: if (i+4<m) {PetscViewerASCIIPrintf(viewer," %d %d %d %d %d %dn",sptr[i],sptr[i+1],sptr[i+2],sptr[i+3],sptr[i+4],sptr[i+5]);}
381: else if (i+3<m) {PetscViewerASCIIPrintf(viewer," %d %d %d %d %dn",sptr[i],sptr[i+1],sptr[i+2],sptr[i+3],sptr[i+4]);}
382: else if (i+2<m) {PetscViewerASCIIPrintf(viewer," %d %d %d %dn",sptr[i],sptr[i+1],sptr[i+2],sptr[i+3]);}
383: else if (i+1<m) {PetscViewerASCIIPrintf(viewer," %d %d %dn",sptr[i],sptr[i+1],sptr[i+2]);}
384: else if (i<m) {PetscViewerASCIIPrintf(viewer," %d %dn",sptr[i],sptr[i+1]);}
385: else {PetscViewerASCIIPrintf(viewer," %dn",sptr[i]);}
386: }
387: PetscViewerASCIIPrintf(viewer,"n");
388: PetscFree(sptr);
389: for (i=0; i<m; i++) {
390: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
391: if (a->j[j] >= i) {PetscViewerASCIIPrintf(viewer," %d ",a->j[j]+fshift);}
392: }
393: PetscViewerASCIIPrintf(viewer,"n");
394: }
395: PetscViewerASCIIPrintf(viewer,"n");
396: for (i=0; i<m; i++) {
397: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
398: if (a->j[j] >= i) {
399: #if defined(PETSC_USE_COMPLEX)
400: if (PetscImaginaryPart(a->a[j]) != 0.0 || PetscRealPart(a->a[j]) != 0.0) {
401: PetscViewerASCIIPrintf(viewer," %18.16e %18.16e ",PetscRealPart(a->a[j]),PetscImaginaryPart(a->a[j]));
402: }
403: #else
404: if (a->a[j] != 0.0) {PetscViewerASCIIPrintf(viewer," %18.16e ",a->a[j]);}
405: #endif
406: }
407: }
408: PetscViewerASCIIPrintf(viewer,"n");
409: }
410: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
411: } else if (format == PETSC_VIEWER_ASCII_DENSE) {
412: int cnt = 0,jcnt;
413: PetscScalar value;
415: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
416: for (i=0; i<m; i++) {
417: jcnt = 0;
418: for (j=0; j<A->n; j++) {
419: if (jcnt < a->i[i+1]-a->i[i] && j == a->j[cnt]) {
420: value = a->a[cnt++];
421: jcnt++;
422: } else {
423: value = 0.0;
424: }
425: #if defined(PETSC_USE_COMPLEX)
426: PetscViewerASCIIPrintf(viewer," %7.5e+%7.5e i ",PetscRealPart(value),PetscImaginaryPart(value));
427: #else
428: PetscViewerASCIIPrintf(viewer," %7.5e ",value);
429: #endif
430: }
431: PetscViewerASCIIPrintf(viewer,"n");
432: }
433: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
434: } else {
435: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
436: for (i=0; i<m; i++) {
437: PetscViewerASCIIPrintf(viewer,"row %d:",i);
438: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
439: #if defined(PETSC_USE_COMPLEX)
440: if (PetscImaginaryPart(a->a[j]) > 0.0) {
441: PetscViewerASCIIPrintf(viewer," (%d, %g + %g i)",a->j[j]+shift,PetscRealPart(a->a[j]),PetscImaginaryPart(a->a[j]));
442: } else if (PetscImaginaryPart(a->a[j]) < 0.0) {
443: PetscViewerASCIIPrintf(viewer," (%d, %g - %g i)",a->j[j]+shift,PetscRealPart(a->a[j]),-PetscImaginaryPart(a->a[j]));
444: } else {
445: PetscViewerASCIIPrintf(viewer," (%d, %g) ",a->j[j]+shift,PetscRealPart(a->a[j]));
446: }
447: #else
448: PetscViewerASCIIPrintf(viewer," (%d, %g) ",a->j[j]+shift,a->a[j]);
449: #endif
450: }
451: PetscViewerASCIIPrintf(viewer,"n");
452: }
453: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
454: }
455: PetscViewerFlush(viewer);
456: return(0);
457: }
459: int MatView_SeqAIJ_Draw_Zoom(PetscDraw draw,void *Aa)
460: {
461: Mat A = (Mat) Aa;
462: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
463: int ierr,i,j,m = A->m,shift = a->indexshift,color;
464: PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r,maxv = 0.0;
465: PetscViewer viewer;
466: PetscViewerFormat format;
469: PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);
470: PetscViewerGetFormat(viewer,&format);
472: PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);
473: /* loop over matrix elements drawing boxes */
475: if (format != PETSC_VIEWER_DRAW_CONTOUR) {
476: /* Blue for negative, Cyan for zero and Red for positive */
477: color = PETSC_DRAW_BLUE;
478: for (i=0; i<m; i++) {
479: y_l = m - i - 1.0; y_r = y_l + 1.0;
480: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
481: x_l = a->j[j] + shift; x_r = x_l + 1.0;
482: #if defined(PETSC_USE_COMPLEX)
483: if (PetscRealPart(a->a[j]) >= 0.) continue;
484: #else
485: if (a->a[j] >= 0.) continue;
486: #endif
487: PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
488: }
489: }
490: color = PETSC_DRAW_CYAN;
491: for (i=0; i<m; i++) {
492: y_l = m - i - 1.0; y_r = y_l + 1.0;
493: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
494: x_l = a->j[j] + shift; x_r = x_l + 1.0;
495: if (a->a[j] != 0.) continue;
496: PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
497: }
498: }
499: color = PETSC_DRAW_RED;
500: for (i=0; i<m; i++) {
501: y_l = m - i - 1.0; y_r = y_l + 1.0;
502: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
503: x_l = a->j[j] + shift; x_r = x_l + 1.0;
504: #if defined(PETSC_USE_COMPLEX)
505: if (PetscRealPart(a->a[j]) <= 0.) continue;
506: #else
507: if (a->a[j] <= 0.) continue;
508: #endif
509: PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
510: }
511: }
512: } else {
513: /* use contour shading to indicate magnitude of values */
514: /* first determine max of all nonzero values */
515: int nz = a->nz,count;
516: PetscDraw popup;
517: PetscReal scale;
519: for (i=0; i<nz; i++) {
520: if (PetscAbsScalar(a->a[i]) > maxv) maxv = PetscAbsScalar(a->a[i]);
521: }
522: scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv;
523: ierr = PetscDrawGetPopup(draw,&popup);
524: if (popup) {ierr = PetscDrawScalePopup(popup,0.0,maxv);}
525: count = 0;
526: for (i=0; i<m; i++) {
527: y_l = m - i - 1.0; y_r = y_l + 1.0;
528: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
529: x_l = a->j[j] + shift; x_r = x_l + 1.0;
530: color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(a->a[count]));
531: ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
532: count++;
533: }
534: }
535: }
536: return(0);
537: }
539: int MatView_SeqAIJ_Draw(Mat A,PetscViewer viewer)
540: {
541: int ierr;
542: PetscDraw draw;
543: PetscReal xr,yr,xl,yl,h,w;
544: PetscTruth isnull;
547: PetscViewerDrawGetDraw(viewer,0,&draw);
548: PetscDrawIsNull(draw,&isnull);
549: if (isnull) return(0);
551: PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);
552: xr = A->n; yr = A->m; h = yr/10.0; w = xr/10.0;
553: xr += w; yr += h; xl = -w; yl = -h;
554: PetscDrawSetCoordinates(draw,xl,yl,xr,yr);
555: PetscDrawZoom(draw,MatView_SeqAIJ_Draw_Zoom,A);
556: PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);
557: return(0);
558: }
560: int MatView_SeqAIJ(Mat A,PetscViewer viewer)
561: {
562: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
563: int ierr;
564: PetscTruth issocket,isascii,isbinary,isdraw;
567: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);
568: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
569: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
570: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
571: if (issocket) {
572: if (a->indexshift) {
573: SETERRQ(1,"Can only socket send sparse matrix with 0 based indexing");
574: }
575: PetscViewerSocketPutSparse_Private(viewer,A->m,A->n,a->nz,a->a,a->i,a->j);
576: } else if (isascii) {
577: MatView_SeqAIJ_ASCII(A,viewer);
578: } else if (isbinary) {
579: MatView_SeqAIJ_Binary(A,viewer);
580: } else if (isdraw) {
581: MatView_SeqAIJ_Draw(A,viewer);
582: } else {
583: SETERRQ1(1,"Viewer type %s not supported by SeqAIJ matrices",((PetscObject)viewer)->type_name);
584: }
585: return(0);
586: }
588: EXTERN int Mat_AIJ_CheckInode(Mat,PetscTruth);
589: EXTERN int MatUseSuperLU_DIST_MPIAIJ(Mat);
590: EXTERN int MatUseSpooles_SeqAIJ(Mat);
591: int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode)
592: {
593: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
594: int fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax,ierr;
595: int m = A->m,*ip,N,*ailen = a->ilen,shift = a->indexshift,rmax = 0;
596: PetscScalar *aa = a->a,*ap;
597: #if defined(PETSC_HAVE_SUPERLUDIST) || defined(PETSC_HAVE_SPOOLES)
598: PetscTruth flag;
599: #endif
602: if (mode == MAT_FLUSH_ASSEMBLY) return(0);
604: if (m) rmax = ailen[0]; /* determine row with most nonzeros */
605: for (i=1; i<m; i++) {
606: /* move each row back by the amount of empty slots (fshift) before it*/
607: fshift += imax[i-1] - ailen[i-1];
608: rmax = PetscMax(rmax,ailen[i]);
609: if (fshift) {
610: ip = aj + ai[i] + shift;
611: ap = aa + ai[i] + shift;
612: N = ailen[i];
613: for (j=0; j<N; j++) {
614: ip[j-fshift] = ip[j];
615: ap[j-fshift] = ap[j];
616: }
617: }
618: ai[i] = ai[i-1] + ailen[i-1];
619: }
620: if (m) {
621: fshift += imax[m-1] - ailen[m-1];
622: ai[m] = ai[m-1] + ailen[m-1];
623: }
624: /* reset ilen and imax for each row */
625: for (i=0; i<m; i++) {
626: ailen[i] = imax[i] = ai[i+1] - ai[i];
627: }
628: a->nz = ai[m] + shift;
630: /* diagonals may have moved, so kill the diagonal pointers */
631: if (fshift && a->diag) {
632: PetscFree(a->diag);
633: PetscLogObjectMemory(A,-(m+1)*sizeof(int));
634: a->diag = 0;
635: }
636: PetscLogInfo(A,"MatAssemblyEnd_SeqAIJ:Matrix size: %d X %d; storage space: %d unneeded,%d usedn",m,A->n,fshift,a->nz);
637: PetscLogInfo(A,"MatAssemblyEnd_SeqAIJ:Number of mallocs during MatSetValues() is %dn",a->reallocs);
638: PetscLogInfo(A,"MatAssemblyEnd_SeqAIJ:Most nonzeros in any row is %dn",rmax);
639: a->reallocs = 0;
640: A->info.nz_unneeded = (double)fshift;
641: a->rmax = rmax;
643: /* check out for identical nodes. If found, use inode functions */
644: Mat_AIJ_CheckInode(A,(PetscTruth)(!fshift));
646: #if defined(PETSC_HAVE_SUPERLUDIST)
647: PetscOptionsHasName(PETSC_NULL,"-mat_aij_superlu_dist",&flag);
648: if (flag) { MatUseSuperLU_DIST_MPIAIJ(A); }
649: #endif
651: #if defined(PETSC_HAVE_SPOOLES)
652: PetscOptionsHasName(PETSC_NULL,"-mat_aij_spooles",&flag);
653: if (flag) { MatUseSpooles_SeqAIJ(A); }
654: #endif
656: return(0);
657: }
659: int MatZeroEntries_SeqAIJ(Mat A)
660: {
661: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
662: int ierr;
665: PetscMemzero(a->a,(a->i[A->m]+a->indexshift)*sizeof(PetscScalar));
666: return(0);
667: }
669: int MatDestroy_SeqAIJ(Mat A)
670: {
671: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
672: int ierr;
675: #if defined(PETSC_USE_LOG)
676: PetscLogObjectState((PetscObject)A,"Rows=%d, Cols=%d, NZ=%d",A->m,A->n,a->nz);
677: #endif
678: if (a->freedata) {
679: PetscFree(a->a);
680: if (!a->singlemalloc) {
681: PetscFree(a->i);
682: PetscFree(a->j);
683: }
684: }
685: if (a->row) {
686: ISDestroy(a->row);
687: }
688: if (a->col) {
689: ISDestroy(a->col);
690: }
691: if (a->diag) {PetscFree(a->diag);}
692: if (a->ilen) {PetscFree(a->ilen);}
693: if (a->imax) {PetscFree(a->imax);}
694: if (a->idiag) {PetscFree(a->idiag);}
695: if (a->solve_work) {PetscFree(a->solve_work);}
696: if (a->inode.size) {PetscFree(a->inode.size);}
697: if (a->icol) {ISDestroy(a->icol);}
698: if (a->saved_values) {PetscFree(a->saved_values);}
699: if (a->coloring) {ISColoringDestroy(a->coloring);}
700: PetscFree(a);
701: return(0);
702: }
704: int MatCompress_SeqAIJ(Mat A)
705: {
707: return(0);
708: }
710: int MatSetOption_SeqAIJ(Mat A,MatOption op)
711: {
712: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
715: switch (op) {
716: case MAT_ROW_ORIENTED:
717: a->roworiented = PETSC_TRUE;
718: break;
719: case MAT_KEEP_ZEROED_ROWS:
720: a->keepzeroedrows = PETSC_TRUE;
721: break;
722: case MAT_COLUMN_ORIENTED:
723: a->roworiented = PETSC_FALSE;
724: break;
725: case MAT_COLUMNS_SORTED:
726: a->sorted = PETSC_TRUE;
727: break;
728: case MAT_COLUMNS_UNSORTED:
729: a->sorted = PETSC_FALSE;
730: break;
731: case MAT_NO_NEW_NONZERO_LOCATIONS:
732: a->nonew = 1;
733: break;
734: case MAT_NEW_NONZERO_LOCATION_ERR:
735: a->nonew = -1;
736: break;
737: case MAT_NEW_NONZERO_ALLOCATION_ERR:
738: a->nonew = -2;
739: break;
740: case MAT_YES_NEW_NONZERO_LOCATIONS:
741: a->nonew = 0;
742: break;
743: case MAT_IGNORE_ZERO_ENTRIES:
744: a->ignorezeroentries = PETSC_TRUE;
745: break;
746: case MAT_USE_INODES:
747: a->inode.use = PETSC_TRUE;
748: break;
749: case MAT_DO_NOT_USE_INODES:
750: a->inode.use = PETSC_FALSE;
751: break;
752: case MAT_ROWS_SORTED:
753: case MAT_ROWS_UNSORTED:
754: case MAT_YES_NEW_DIAGONALS:
755: case MAT_IGNORE_OFF_PROC_ENTRIES:
756: case MAT_USE_HASH_TABLE:
757: case MAT_USE_SINGLE_PRECISION_SOLVES:
758: PetscLogInfo(A,"MatSetOption_SeqAIJ:Option ignoredn");
759: break;
760: case MAT_NO_NEW_DIAGONALS:
761: SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
762: case MAT_INODE_LIMIT_1:
763: a->inode.limit = 1;
764: break;
765: case MAT_INODE_LIMIT_2:
766: a->inode.limit = 2;
767: break;
768: case MAT_INODE_LIMIT_3:
769: a->inode.limit = 3;
770: break;
771: case MAT_INODE_LIMIT_4:
772: a->inode.limit = 4;
773: break;
774: case MAT_INODE_LIMIT_5:
775: a->inode.limit = 5;
776: break;
777: default:
778: SETERRQ(PETSC_ERR_SUP,"unknown option");
779: }
780: return(0);
781: }
783: int MatGetDiagonal_SeqAIJ(Mat A,Vec v)
784: {
785: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
786: int i,j,n,shift = a->indexshift,ierr;
787: PetscScalar *x,zero = 0.0;
790: VecSet(&zero,v);
791: VecGetArray(v,&x);
792: VecGetLocalSize(v,&n);
793: if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
794: for (i=0; i<A->m; i++) {
795: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
796: if (a->j[j]+shift == i) {
797: x[i] = a->a[j];
798: break;
799: }
800: }
801: }
802: VecRestoreArray(v,&x);
803: return(0);
804: }
807: int MatMultTransposeAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
808: {
809: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
810: PetscScalar *x,*y;
811: int ierr,m = A->m,shift = a->indexshift;
812: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
813: PetscScalar *v,alpha;
814: int n,i,*idx;
815: #endif
818: if (zz != yy) {VecCopy(zz,yy);}
819: VecGetArray(xx,&x);
820: VecGetArray(yy,&y);
821: y = y + shift; /* shift for Fortran start by 1 indexing */
823: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
824: fortranmulttransposeaddaij_(&m,x,a->i,a->j+shift,a->a+shift,y);
825: #else
826: for (i=0; i<m; i++) {
827: idx = a->j + a->i[i] + shift;
828: v = a->a + a->i[i] + shift;
829: n = a->i[i+1] - a->i[i];
830: alpha = x[i];
831: while (n-->0) {y[*idx++] += alpha * *v++;}
832: }
833: #endif
834: PetscLogFlops(2*a->nz);
835: VecRestoreArray(xx,&x);
836: VecRestoreArray(yy,&y);
837: return(0);
838: }
840: int MatMultTranspose_SeqAIJ(Mat A,Vec xx,Vec yy)
841: {
842: PetscScalar zero = 0.0;
843: int ierr;
846: VecSet(&zero,yy);
847: MatMultTransposeAdd_SeqAIJ(A,xx,yy,yy);
848: return(0);
849: }
852: int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
853: {
854: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
855: PetscScalar *x,*y,*v,sum;
856: int ierr,m = A->m,*idx,shift = a->indexshift,*ii;
857: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
858: int n,i,jrow,j;
859: #endif
861: #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
862: #pragma disjoint(*x,*y,*v)
863: #endif
866: VecGetArray(xx,&x);
867: VecGetArray(yy,&y);
868: x = x + shift; /* shift for Fortran start by 1 indexing */
869: idx = a->j;
870: v = a->a;
871: ii = a->i;
872: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
873: fortranmultaij_(&m,x,ii,idx+shift,v+shift,y);
874: #else
875: v += shift; /* shift for Fortran start by 1 indexing */
876: idx += shift;
877: for (i=0; i<m; i++) {
878: jrow = ii[i];
879: n = ii[i+1] - jrow;
880: sum = 0.0;
881: for (j=0; j<n; j++) {
882: sum += v[jrow]*x[idx[jrow]]; jrow++;
883: }
884: y[i] = sum;
885: }
886: #endif
887: PetscLogFlops(2*a->nz - m);
888: VecRestoreArray(xx,&x);
889: VecRestoreArray(yy,&y);
890: return(0);
891: }
893: int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
894: {
895: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
896: PetscScalar *x,*y,*z,*v,sum;
897: int ierr,m = A->m,*idx,shift = a->indexshift,*ii;
898: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
899: int n,i,jrow,j;
900: #endif
903: VecGetArray(xx,&x);
904: VecGetArray(yy,&y);
905: if (zz != yy) {
906: VecGetArray(zz,&z);
907: } else {
908: z = y;
909: }
910: x = x + shift; /* shift for Fortran start by 1 indexing */
911: idx = a->j;
912: v = a->a;
913: ii = a->i;
914: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
915: fortranmultaddaij_(&m,x,ii,idx+shift,v+shift,y,z);
916: #else
917: v += shift; /* shift for Fortran start by 1 indexing */
918: idx += shift;
919: for (i=0; i<m; i++) {
920: jrow = ii[i];
921: n = ii[i+1] - jrow;
922: sum = y[i];
923: for (j=0; j<n; j++) {
924: sum += v[jrow]*x[idx[jrow]]; jrow++;
925: }
926: z[i] = sum;
927: }
928: #endif
929: PetscLogFlops(2*a->nz);
930: VecRestoreArray(xx,&x);
931: VecRestoreArray(yy,&y);
932: if (zz != yy) {
933: VecRestoreArray(zz,&z);
934: }
935: return(0);
936: }
938: /*
939: Adds diagonal pointers to sparse matrix structure.
940: */
941: int MatMarkDiagonal_SeqAIJ(Mat A)
942: {
943: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
944: int i,j,*diag,m = A->m,shift = a->indexshift,ierr;
947: if (a->diag) return(0);
949: PetscMalloc((m+1)*sizeof(int),&diag);
950: PetscLogObjectMemory(A,(m+1)*sizeof(int));
951: for (i=0; i<A->m; i++) {
952: diag[i] = a->i[i+1];
953: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
954: if (a->j[j]+shift == i) {
955: diag[i] = j - shift;
956: break;
957: }
958: }
959: }
960: a->diag = diag;
961: return(0);
962: }
964: /*
965: Checks for missing diagonals
966: */
967: int MatMissingDiagonal_SeqAIJ(Mat A)
968: {
969: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
970: int *diag,*jj = a->j,i,shift = a->indexshift,ierr;
973: MatMarkDiagonal_SeqAIJ(A);
974: diag = a->diag;
975: for (i=0; i<A->m; i++) {
976: if (jj[diag[i]+shift] != i-shift) {
977: SETERRQ1(1,"Matrix is missing diagonal number %d",i);
978: }
979: }
980: return(0);
981: }
983: int MatRelax_SeqAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,int its,int lits,Vec xx)
984: {
985: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
986: PetscScalar *x,*b,*bs, d,*xs,sum,*v = a->a,*t=0,scale,*ts,*xb,*idiag=0;
987: int ierr,*idx,*diag,n = A->n,m = A->m,i,shift = a->indexshift;
990: its = its*lits;
991: if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %d and local its %d both positive",its,lits);
993: VecGetArray(xx,&x);
994: if (xx != bb) {
995: VecGetArray(bb,&b);
996: } else {
997: b = x;
998: }
1000: if (!a->diag) {MatMarkDiagonal_SeqAIJ(A);}
1001: diag = a->diag;
1002: xs = x + shift; /* shifted by one for index start of a or a->j*/
1003: if (flag == SOR_APPLY_UPPER) {
1004: /* apply (U + D/omega) to the vector */
1005: bs = b + shift;
1006: for (i=0; i<m; i++) {
1007: d = fshift + a->a[diag[i] + shift];
1008: n = a->i[i+1] - diag[i] - 1;
1009: PetscLogFlops(2*n-1);
1010: idx = a->j + diag[i] + (!shift);
1011: v = a->a + diag[i] + (!shift);
1012: sum = b[i]*d/omega;
1013: SPARSEDENSEDOT(sum,bs,v,idx,n);
1014: x[i] = sum;
1015: }
1016: VecRestoreArray(xx,&x);
1017: if (bb != xx) {VecRestoreArray(bb,&b);}
1018: return(0);
1019: }
1021: /* setup workspace for Eisenstat */
1022: if (flag & SOR_EISENSTAT) {
1023: if (!a->idiag) {
1024: ierr = PetscMalloc(2*m*sizeof(PetscScalar),&a->idiag);
1025: a->ssor = a->idiag + m;
1026: v = a->a;
1027: for (i=0; i<m; i++) { a->idiag[i] = 1.0/v[diag[i]];}
1028: }
1029: t = a->ssor;
1030: idiag = a->idiag;
1031: }
1032: /* Let A = L + U + D; where L is lower trianglar,
1033: U is upper triangular, E is diagonal; This routine applies
1035: (L + E)^{-1} A (U + E)^{-1}
1037: to a vector efficiently using Eisenstat's trick. This is for
1038: the case of SSOR preconditioner, so E is D/omega where omega
1039: is the relaxation factor.
1040: */
1042: if (flag == SOR_APPLY_LOWER) {
1043: SETERRQ(PETSC_ERR_SUP,"SOR_APPLY_LOWER is not implemented");
1044: } else if ((flag & SOR_EISENSTAT) && omega == 1.0 && shift == 0 && fshift == 0.0) {
1045: /* special case for omega = 1.0 saves flops and some integer ops */
1046: PetscScalar *v2;
1047:
1048: v2 = a->a;
1049: /* x = (E + U)^{-1} b */
1050: for (i=m-1; i>=0; i--) {
1051: n = a->i[i+1] - diag[i] - 1;
1052: idx = a->j + diag[i] + 1;
1053: v = a->a + diag[i] + 1;
1054: sum = b[i];
1055: SPARSEDENSEMDOT(sum,xs,v,idx,n);
1056: x[i] = sum*idiag[i];
1058: /* t = b - (2*E - D)x */
1059: t[i] = b[i] - (v2[diag[i]])*x[i];
1060: }
1062: /* t = (E + L)^{-1}t */
1063: diag = a->diag;
1064: for (i=0; i<m; i++) {
1065: n = diag[i] - a->i[i];
1066: idx = a->j + a->i[i];
1067: v = a->a + a->i[i];
1068: sum = t[i];
1069: SPARSEDENSEMDOT(sum,t,v,idx,n);
1070: t[i] = sum*idiag[i];
1072: /* x = x + t */
1073: x[i] += t[i];
1074: }
1076: PetscLogFlops(3*m-1 + 2*a->nz);
1077: VecRestoreArray(xx,&x);
1078: if (bb != xx) {VecRestoreArray(bb,&b);}
1079: return(0);
1080: } else if (flag & SOR_EISENSTAT) {
1081: /* Let A = L + U + D; where L is lower trianglar,
1082: U is upper triangular, E is diagonal; This routine applies
1084: (L + E)^{-1} A (U + E)^{-1}
1086: to a vector efficiently using Eisenstat's trick. This is for
1087: the case of SSOR preconditioner, so E is D/omega where omega
1088: is the relaxation factor.
1089: */
1090: scale = (2.0/omega) - 1.0;
1092: /* x = (E + U)^{-1} b */
1093: for (i=m-1; i>=0; i--) {
1094: d = fshift + a->a[diag[i] + shift];
1095: n = a->i[i+1] - diag[i] - 1;
1096: idx = a->j + diag[i] + (!shift);
1097: v = a->a + diag[i] + (!shift);
1098: sum = b[i];
1099: SPARSEDENSEMDOT(sum,xs,v,idx,n);
1100: x[i] = omega*(sum/d);
1101: }
1103: /* t = b - (2*E - D)x */
1104: v = a->a;
1105: for (i=0; i<m; i++) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; }
1107: /* t = (E + L)^{-1}t */
1108: ts = t + shift; /* shifted by one for index start of a or a->j*/
1109: diag = a->diag;
1110: for (i=0; i<m; i++) {
1111: d = fshift + a->a[diag[i]+shift];
1112: n = diag[i] - a->i[i];
1113: idx = a->j + a->i[i] + shift;
1114: v = a->a + a->i[i] + shift;
1115: sum = t[i];
1116: SPARSEDENSEMDOT(sum,ts,v,idx,n);
1117: t[i] = omega*(sum/d);
1118: /* x = x + t */
1119: x[i] += t[i];
1120: }
1122: PetscLogFlops(6*m-1 + 2*a->nz);
1123: VecRestoreArray(xx,&x);
1124: if (bb != xx) {VecRestoreArray(bb,&b);}
1125: return(0);
1126: }
1127: if (flag & SOR_ZERO_INITIAL_GUESS) {
1128: if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
1129: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1130: fortranrelaxaijforwardzero_(&m,&omega,x,a->i,a->j,diag,a->a,b);
1131: #else
1132: for (i=0; i<m; i++) {
1133: d = fshift + a->a[diag[i]+shift];
1134: n = diag[i] - a->i[i];
1135: PetscLogFlops(2*n-1);
1136: idx = a->j + a->i[i] + shift;
1137: v = a->a + a->i[i] + shift;
1138: sum = b[i];
1139: SPARSEDENSEMDOT(sum,xs,v,idx,n);
1140: x[i] = omega*(sum/d);
1141: }
1142: #endif
1143: xb = x;
1144: } else xb = b;
1145: if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
1146: (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
1147: for (i=0; i<m; i++) {
1148: x[i] *= a->a[diag[i]+shift];
1149: }
1150: PetscLogFlops(m);
1151: }
1152: if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
1153: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1154: fortranrelaxaijbackwardzero_(&m,&omega,x,a->i,a->j,diag,a->a,xb);
1155: #else
1156: for (i=m-1; i>=0; i--) {
1157: d = fshift + a->a[diag[i] + shift];
1158: n = a->i[i+1] - diag[i] - 1;
1159: PetscLogFlops(2*n-1);
1160: idx = a->j + diag[i] + (!shift);
1161: v = a->a + diag[i] + (!shift);
1162: sum = xb[i];
1163: SPARSEDENSEMDOT(sum,xs,v,idx,n);
1164: x[i] = omega*(sum/d);
1165: }
1166: #endif
1167: }
1168: its--;
1169: }
1170: while (its--) {
1171: if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
1172: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1173: fortranrelaxaijforward_(&m,&omega,x,a->i,a->j,diag,a->a,b);
1174: #else
1175: for (i=0; i<m; i++) {
1176: d = fshift + a->a[diag[i]+shift];
1177: n = a->i[i+1] - a->i[i];
1178: PetscLogFlops(2*n-1);
1179: idx = a->j + a->i[i] + shift;
1180: v = a->a + a->i[i] + shift;
1181: sum = b[i];
1182: SPARSEDENSEMDOT(sum,xs,v,idx,n);
1183: x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d;
1184: }
1185: #endif
1186: }
1187: if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
1188: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1189: fortranrelaxaijbackward_(&m,&omega,x,a->i,a->j,diag,a->a,b);
1190: #else
1191: for (i=m-1; i>=0; i--) {
1192: d = fshift + a->a[diag[i] + shift];
1193: n = a->i[i+1] - a->i[i];
1194: PetscLogFlops(2*n-1);
1195: idx = a->j + a->i[i] + shift;
1196: v = a->a + a->i[i] + shift;
1197: sum = b[i];
1198: SPARSEDENSEMDOT(sum,xs,v,idx,n);
1199: x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d;
1200: }
1201: #endif
1202: }
1203: }
1204: VecRestoreArray(xx,&x);
1205: if (bb != xx) {VecRestoreArray(bb,&b);}
1206: return(0);
1207: }
1209: int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,MatInfo *info)
1210: {
1211: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1214: info->rows_global = (double)A->m;
1215: info->columns_global = (double)A->n;
1216: info->rows_local = (double)A->m;
1217: info->columns_local = (double)A->n;
1218: info->block_size = 1.0;
1219: info->nz_allocated = (double)a->maxnz;
1220: info->nz_used = (double)a->nz;
1221: info->nz_unneeded = (double)(a->maxnz - a->nz);
1222: info->assemblies = (double)A->num_ass;
1223: info->mallocs = (double)a->reallocs;
1224: info->memory = A->mem;
1225: if (A->factor) {
1226: info->fill_ratio_given = A->info.fill_ratio_given;
1227: info->fill_ratio_needed = A->info.fill_ratio_needed;
1228: info->factor_mallocs = A->info.factor_mallocs;
1229: } else {
1230: info->fill_ratio_given = 0;
1231: info->fill_ratio_needed = 0;
1232: info->factor_mallocs = 0;
1233: }
1234: return(0);
1235: }
1237: EXTERN int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,MatLUInfo*,Mat*);
1238: EXTERN int MatLUFactorNumeric_SeqAIJ(Mat,Mat*);
1239: EXTERN int MatLUFactor_SeqAIJ(Mat,IS,IS,MatLUInfo*);
1240: EXTERN int MatSolve_SeqAIJ(Mat,Vec,Vec);
1241: EXTERN int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
1242: EXTERN int MatSolveTranspose_SeqAIJ(Mat,Vec,Vec);
1243: EXTERN int MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec);
1245: int MatZeroRows_SeqAIJ(Mat A,IS is,PetscScalar *diag)
1246: {
1247: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1248: int i,ierr,N,*rows,m = A->m - 1,shift = a->indexshift;
1251: ISGetLocalSize(is,&N);
1252: ISGetIndices(is,&rows);
1253: if (a->keepzeroedrows) {
1254: for (i=0; i<N; i++) {
1255: if (rows[i] < 0 || rows[i] > m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"row out of range");
1256: PetscMemzero(&a->a[a->i[rows[i]]+shift],a->ilen[rows[i]]*sizeof(PetscScalar));
1257: }
1258: if (diag) {
1259: MatMissingDiagonal_SeqAIJ(A);
1260: MatMarkDiagonal_SeqAIJ(A);
1261: for (i=0; i<N; i++) {
1262: a->a[a->diag[rows[i]]] = *diag;
1263: }
1264: }
1265: } else {
1266: if (diag) {
1267: for (i=0; i<N; i++) {
1268: if (rows[i] < 0 || rows[i] > m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"row out of range");
1269: if (a->ilen[rows[i]] > 0) {
1270: a->ilen[rows[i]] = 1;
1271: a->a[a->i[rows[i]]+shift] = *diag;
1272: a->j[a->i[rows[i]]+shift] = rows[i]+shift;
1273: } else { /* in case row was completely empty */
1274: MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES);
1275: }
1276: }
1277: } else {
1278: for (i=0; i<N; i++) {
1279: if (rows[i] < 0 || rows[i] > m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"row out of range");
1280: a->ilen[rows[i]] = 0;
1281: }
1282: }
1283: }
1284: ISRestoreIndices(is,&rows);
1285: MatAssemblyEnd_SeqAIJ(A,MAT_FINAL_ASSEMBLY);
1286: return(0);
1287: }
1289: int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,PetscScalar **v)
1290: {
1291: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1292: int *itmp,i,shift = a->indexshift,ierr;
1295: if (row < 0 || row >= A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Row %d out of range",row);
1297: *nz = a->i[row+1] - a->i[row];
1298: if (v) *v = a->a + a->i[row] + shift;
1299: if (idx) {
1300: itmp = a->j + a->i[row] + shift;
1301: if (*nz && shift) {
1302: PetscMalloc((*nz)*sizeof(int),idx);
1303: for (i=0; i<(*nz); i++) {(*idx)[i] = itmp[i] + shift;}
1304: } else if (*nz) {
1305: *idx = itmp;
1306: }
1307: else *idx = 0;
1308: }
1309: return(0);
1310: }
1312: int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,PetscScalar **v)
1313: {
1314: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1318: if (idx) {if (*idx && a->indexshift) {PetscFree(*idx);}}
1319: return(0);
1320: }
1322: int MatNorm_SeqAIJ(Mat A,NormType type,PetscReal *nrm)
1323: {
1324: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1325: PetscScalar *v = a->a;
1326: PetscReal sum = 0.0;
1327: int i,j,shift = a->indexshift,ierr;
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: int *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++ + shift] += 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] + shift;
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: }
1368: int MatTranspose_SeqAIJ(Mat A,Mat *B)
1369: {
1370: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1371: Mat C;
1372: int i,ierr,*aj = a->j,*ai = a->i,m = A->m,len,*col;
1373: int shift = a->indexshift;
1374: PetscScalar *array = a->a;
1377: if (!B && m != A->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Square matrix only for in-place");
1378: PetscMalloc((1+A->n)*sizeof(int),&col);
1379: PetscMemzero(col,(1+A->n)*sizeof(int));
1380: if (shift) {
1381: for (i=0; i<ai[m]-1; i++) aj[i] -= 1;
1382: }
1383: for (i=0; i<ai[m]+shift; i++) col[aj[i]] += 1;
1384: MatCreateSeqAIJ(A->comm,A->n,m,0,col,&C);
1385: PetscFree(col);
1386: for (i=0; i<m; i++) {
1387: len = ai[i+1]-ai[i];
1388: ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES);
1389: array += len;
1390: aj += len;
1391: }
1392: if (shift) {
1393: for (i=0; i<ai[m]-1; i++) aj[i] += 1;
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: }
1407: int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr)
1408: {
1409: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1410: PetscScalar *l,*r,x,*v;
1411: int ierr,i,j,m = A->m,n = A->n,M,nz = a->nz,*jj,shift = a->indexshift;
1414: if (ll) {
1415: /* The local size is used so that VecMPI can be passed to this routine
1416: by MatDiagonalScale_MPIAIJ */
1417: VecGetLocalSize(ll,&m);
1418: if (m != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length");
1419: VecGetArray(ll,&l);
1420: v = a->a;
1421: for (i=0; i<m; i++) {
1422: x = l[i];
1423: M = a->i[i+1] - a->i[i];
1424: for (j=0; j<M; j++) { (*v++) *= x;}
1425: }
1426: VecRestoreArray(ll,&l);
1427: PetscLogFlops(nz);
1428: }
1429: if (rr) {
1430: VecGetLocalSize(rr,&n);
1431: if (n != A->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vector wrong length");
1432: VecGetArray(rr,&r);
1433: v = a->a; jj = a->j;
1434: for (i=0; i<nz; i++) {
1435: (*v++) *= r[*jj++ + shift];
1436: }
1437: VecRestoreArray(rr,&r);
1438: PetscLogFlops(nz);
1439: }
1440: return(0);
1441: }
1443: int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,int csize,MatReuse scall,Mat *B)
1444: {
1445: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*c;
1446: int *smap,i,k,kstart,kend,ierr,oldcols = A->n,*lens;
1447: int row,mat_i,*mat_j,tcol,first,step,*mat_ilen,sum,lensi;
1448: int *irow,*icol,nrows,ncols,shift = a->indexshift,*ssmap;
1449: int *starts,*j_new,*i_new,*aj = a->j,*ai = a->i,ii,*ailen = a->ilen;
1450: PetscScalar *a_new,*mat_a;
1451: Mat C;
1452: PetscTruth stride;
1455: ISSorted(isrow,(PetscTruth*)&i);
1456: if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"ISrow is not sorted");
1457: ISSorted(iscol,(PetscTruth*)&i);
1458: if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IScol is not sorted");
1460: ISGetIndices(isrow,&irow);
1461: ISGetLocalSize(isrow,&nrows);
1462: ISGetLocalSize(iscol,&ncols);
1464: ISStrideGetInfo(iscol,&first,&step);
1465: ISStride(iscol,&stride);
1466: if (stride && step == 1) {
1467: /* special case of contiguous rows */
1468: ierr = PetscMalloc((2*nrows+1)*sizeof(int),&lens);
1469: starts = lens + nrows;
1470: /* loop over new rows determining lens and starting points */
1471: for (i=0; i<nrows; i++) {
1472: kstart = ai[irow[i]]+shift;
1473: kend = kstart + ailen[irow[i]];
1474: for (k=kstart; k<kend; k++) {
1475: if (aj[k]+shift >= first) {
1476: starts[i] = k;
1477: break;
1478: }
1479: }
1480: sum = 0;
1481: while (k < kend) {
1482: if (aj[k++]+shift >= first+ncols) break;
1483: sum++;
1484: }
1485: lens[i] = sum;
1486: }
1487: /* create submatrix */
1488: if (scall == MAT_REUSE_MATRIX) {
1489: int n_cols,n_rows;
1490: MatGetSize(*B,&n_rows,&n_cols);
1491: if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size");
1492: MatZeroEntries(*B);
1493: C = *B;
1494: } else {
1495: MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);
1496: }
1497: c = (Mat_SeqAIJ*)C->data;
1499: /* loop over rows inserting into submatrix */
1500: a_new = c->a;
1501: j_new = c->j;
1502: i_new = c->i;
1503: i_new[0] = -shift;
1504: for (i=0; i<nrows; i++) {
1505: ii = starts[i];
1506: lensi = lens[i];
1507: for (k=0; k<lensi; k++) {
1508: *j_new++ = aj[ii+k] - first;
1509: }
1510: PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(PetscScalar));
1511: a_new += lensi;
1512: i_new[i+1] = i_new[i] + lensi;
1513: c->ilen[i] = lensi;
1514: }
1515: PetscFree(lens);
1516: } else {
1517: ierr = ISGetIndices(iscol,&icol);
1518: ierr = PetscMalloc((1+oldcols)*sizeof(int),&smap);
1519: ssmap = smap + shift;
1520: ierr = PetscMalloc((1+nrows)*sizeof(int),&lens);
1521: ierr = PetscMemzero(smap,oldcols*sizeof(int));
1522: for (i=0; i<ncols; i++) smap[icol[i]] = i+1;
1523: /* determine lens of each row */
1524: for (i=0; i<nrows; i++) {
1525: kstart = ai[irow[i]]+shift;
1526: kend = kstart + a->ilen[irow[i]];
1527: lens[i] = 0;
1528: for (k=kstart; k<kend; k++) {
1529: if (ssmap[aj[k]]) {
1530: lens[i]++;
1531: }
1532: }
1533: }
1534: /* Create and fill new matrix */
1535: if (scall == MAT_REUSE_MATRIX) {
1536: PetscTruth equal;
1538: c = (Mat_SeqAIJ *)((*B)->data);
1539: if ((*B)->m != nrows || (*B)->n != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size");
1540: PetscMemcmp(c->ilen,lens,(*B)->m*sizeof(int),&equal);
1541: if (!equal) {
1542: SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros");
1543: }
1544: PetscMemzero(c->ilen,(*B)->m*sizeof(int));
1545: C = *B;
1546: } else {
1547: MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);
1548: }
1549: c = (Mat_SeqAIJ *)(C->data);
1550: for (i=0; i<nrows; i++) {
1551: row = irow[i];
1552: kstart = ai[row]+shift;
1553: kend = kstart + a->ilen[row];
1554: mat_i = c->i[i]+shift;
1555: mat_j = c->j + mat_i;
1556: mat_a = c->a + mat_i;
1557: mat_ilen = c->ilen + i;
1558: for (k=kstart; k<kend; k++) {
1559: if ((tcol=ssmap[a->j[k]])) {
1560: *mat_j++ = tcol - (!shift);
1561: *mat_a++ = a->a[k];
1562: (*mat_ilen)++;
1564: }
1565: }
1566: }
1567: /* Free work space */
1568: ISRestoreIndices(iscol,&icol);
1569: PetscFree(smap);
1570: PetscFree(lens);
1571: }
1572: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1573: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
1575: ISRestoreIndices(isrow,&irow);
1576: *B = C;
1577: return(0);
1578: }
1580: /*
1581: */
1582: int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,MatILUInfo *info)
1583: {
1584: Mat_SeqAIJ *a = (Mat_SeqAIJ*)inA->data;
1585: int ierr;
1586: Mat outA;
1587: PetscTruth row_identity,col_identity;
1590: if (info && info->levels != 0) SETERRQ(PETSC_ERR_SUP,"Only levels=0 supported for in-place ilu");
1591: ISIdentity(row,&row_identity);
1592: ISIdentity(col,&col_identity);
1593: if (!row_identity || !col_identity) {
1594: SETERRQ(1,"Row and column permutations must be identity for in-place ILU");
1595: }
1597: outA = inA;
1598: inA->factor = FACTOR_LU;
1599: a->row = row;
1600: a->col = col;
1601: PetscObjectReference((PetscObject)row);
1602: PetscObjectReference((PetscObject)col);
1604: /* Create the inverse permutation so that it can be used in MatLUFactorNumeric() */
1605: if (a->icol) {ISDestroy(a->icol);} /* need to remove old one */
1606: ISInvertPermutation(col,PETSC_DECIDE,&a->icol);
1607: PetscLogObjectParent(inA,a->icol);
1609: if (!a->solve_work) { /* this matrix may have been factored before */
1610: PetscMalloc((inA->m+1)*sizeof(PetscScalar),&a->solve_work);
1611: }
1613: if (!a->diag) {
1614: MatMarkDiagonal_SeqAIJ(inA);
1615: }
1616: MatLUFactorNumeric_SeqAIJ(inA,&outA);
1617: return(0);
1618: }
1620: #include petscblaslapack.h
1621: int MatScale_SeqAIJ(PetscScalar *alpha,Mat inA)
1622: {
1623: Mat_SeqAIJ *a = (Mat_SeqAIJ*)inA->data;
1624: int one = 1;
1627: BLscal_(&a->nz,alpha,a->a,&one);
1628: PetscLogFlops(a->nz);
1629: return(0);
1630: }
1632: int MatGetSubMatrices_SeqAIJ(Mat A,int n,IS *irow,IS *icol,MatReuse scall,Mat **B)
1633: {
1634: int ierr,i;
1637: if (scall == MAT_INITIAL_MATRIX) {
1638: PetscMalloc((n+1)*sizeof(Mat),B);
1639: }
1641: for (i=0; i<n; i++) {
1642: MatGetSubMatrix_SeqAIJ(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);
1643: }
1644: return(0);
1645: }
1647: int MatGetBlockSize_SeqAIJ(Mat A,int *bs)
1648: {
1650: *bs = 1;
1651: return(0);
1652: }
1654: int MatIncreaseOverlap_SeqAIJ(Mat A,int is_max,IS *is,int ov)
1655: {
1656: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1657: int shift,row,i,j,k,l,m,n,*idx,ierr,*nidx,isz,val;
1658: int start,end,*ai,*aj;
1659: PetscBT table;
1662: shift = a->indexshift;
1663: m = A->m;
1664: ai = a->i;
1665: aj = a->j+shift;
1667: if (ov < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"illegal overlap value used");
1669: PetscMalloc((m+1)*sizeof(int),&nidx);
1670: PetscBTCreate(m,table);
1672: for (i=0; i<is_max; i++) {
1673: /* Initialize the two local arrays */
1674: isz = 0;
1675: PetscBTMemzero(m,table);
1676:
1677: /* Extract the indices, assume there can be duplicate entries */
1678: ISGetIndices(is[i],&idx);
1679: ISGetLocalSize(is[i],&n);
1680:
1681: /* Enter these into the temp arrays. I.e., mark table[row], enter row into new index */
1682: for (j=0; j<n ; ++j){
1683: if(!PetscBTLookupSet(table,idx[j])) { nidx[isz++] = idx[j];}
1684: }
1685: ISRestoreIndices(is[i],&idx);
1686: ISDestroy(is[i]);
1687:
1688: k = 0;
1689: for (j=0; j<ov; j++){ /* for each overlap */
1690: n = isz;
1691: for (; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */
1692: row = nidx[k];
1693: start = ai[row];
1694: end = ai[row+1];
1695: for (l = start; l<end ; l++){
1696: val = aj[l] + shift;
1697: if (!PetscBTLookupSet(table,val)) {nidx[isz++] = val;}
1698: }
1699: }
1700: }
1701: ISCreateGeneral(PETSC_COMM_SELF,isz,nidx,(is+i));
1702: }
1703: PetscBTDestroy(table);
1704: PetscFree(nidx);
1705: return(0);
1706: }
1708: /* -------------------------------------------------------------- */
1709: int MatPermute_SeqAIJ(Mat A,IS rowp,IS colp,Mat *B)
1710: {
1711: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1712: PetscScalar *vwork;
1713: int i,ierr,nz,m = A->m,n = A->n,*cwork;
1714: int *row,*col,*cnew,j,*lens;
1715: IS icolp,irowp;
1718: ISInvertPermutation(rowp,PETSC_DECIDE,&irowp);
1719: ISGetIndices(irowp,&row);
1720: ISInvertPermutation(colp,PETSC_DECIDE,&icolp);
1721: ISGetIndices(icolp,&col);
1722:
1723: /* determine lengths of permuted rows */
1724: PetscMalloc((m+1)*sizeof(int),&lens);
1725: for (i=0; i<m; i++) {
1726: lens[row[i]] = a->i[i+1] - a->i[i];
1727: }
1728: MatCreateSeqAIJ(A->comm,m,n,0,lens,B);
1729: PetscFree(lens);
1731: PetscMalloc(n*sizeof(int),&cnew);
1732: for (i=0; i<m; i++) {
1733: MatGetRow(A,i,&nz,&cwork,&vwork);
1734: for (j=0; j<nz; j++) { cnew[j] = col[cwork[j]];}
1735: MatSetValues(*B,1,&row[i],nz,cnew,vwork,INSERT_VALUES);
1736: MatRestoreRow(A,i,&nz,&cwork,&vwork);
1737: }
1738: PetscFree(cnew);
1739: MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);
1740: MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);
1741: ISRestoreIndices(irowp,&row);
1742: ISRestoreIndices(icolp,&col);
1743: ISDestroy(irowp);
1744: ISDestroy(icolp);
1745: return(0);
1746: }
1748: int MatPrintHelp_SeqAIJ(Mat A)
1749: {
1750: static PetscTruth called = PETSC_FALSE;
1751: MPI_Comm comm = A->comm;
1752: int ierr;
1755: if (called) {return(0);} else called = PETSC_TRUE;
1756: (*PetscHelpPrintf)(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):n");
1757: (*PetscHelpPrintf)(comm," -mat_lu_pivotthreshold <threshold>: Set pivoting thresholdn");
1758: (*PetscHelpPrintf)(comm," -mat_aij_oneindex: internal indices begin at 1 instead of the default 0.n");
1759: (*PetscHelpPrintf)(comm," -mat_aij_no_inode: Do not use inodesn");
1760: (*PetscHelpPrintf)(comm," -mat_aij_inode_limit <limit>: Set inode limit (max limit=5)n");
1761: #if defined(PETSC_HAVE_ESSL)
1762: (*PetscHelpPrintf)(comm," -mat_aij_essl: Use IBM sparse LU factorization and solve.n");
1763: #endif
1764: #if defined(PETSC_HAVE_LUSOL)
1765: (*PetscHelpPrintf)(comm," -mat_aij_lusol: Use the Stanford LUSOL sparse factorization and solve.n");
1766: #endif
1767: #if defined(PETSC_HAVE_MATLAB_ENGINE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
1768: (*PetscHelpPrintf)(comm," -mat_aij_matlab: Use Matlab engine sparse LU factorization and solve.n");
1769: #endif
1770: return(0);
1771: }
1772: EXTERN int MatEqual_SeqAIJ(Mat A,Mat B,PetscTruth* flg);
1773: EXTERN int MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring);
1774: EXTERN int MatILUDTFactor_SeqAIJ(Mat,MatILUInfo*,IS,IS,Mat*);
1775: int MatCopy_SeqAIJ(Mat A,Mat B,MatStructure str)
1776: {
1777: int ierr;
1778: PetscTruth flg;
1781: PetscTypeCompare((PetscObject)B,MATSEQAIJ,&flg);
1782: if (str == SAME_NONZERO_PATTERN && flg) {
1783: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1784: Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data;
1786: if (a->i[A->m]+a->indexshift != b->i[B->m]+a->indexshift) {
1787: SETERRQ(1,"Number of nonzeros in two matrices are different");
1788: }
1789: PetscMemcpy(b->a,a->a,(a->i[A->m]+a->indexshift)*sizeof(PetscScalar));
1790: } else {
1791: MatCopy_Basic(A,B,str);
1792: }
1793: return(0);
1794: }
1796: int MatSetUpPreallocation_SeqAIJ(Mat A)
1797: {
1798: int ierr;
1801: MatSeqAIJSetPreallocation(A,PETSC_DEFAULT,0);
1802: return(0);
1803: }
1805: int MatGetArray_SeqAIJ(Mat A,PetscScalar **array)
1806: {
1807: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1809: *array = a->a;
1810: return(0);
1811: }
1813: int MatRestoreArray_SeqAIJ(Mat A,PetscScalar **array)
1814: {
1816: return(0);
1817: }
1819: int MatFDColoringApply_SeqAIJ(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
1820: {
1821: int (*f)(void *,Vec,Vec,void*) = (int (*)(void *,Vec,Vec,void *))coloring->f;
1822: int k,ierr,N,start,end,l,row,col,srow,**vscaleforrow,m1,m2;
1823: PetscScalar dx,mone = -1.0,*y,*xx,*w3_array;
1824: PetscScalar *vscale_array;
1825: PetscReal epsilon = coloring->error_rel,umin = coloring->umin;
1826: Vec w1,w2,w3;
1827: void *fctx = coloring->fctx;
1828: PetscTruth flg;
1831: if (!coloring->w1) {
1832: VecDuplicate(x1,&coloring->w1);
1833: PetscLogObjectParent(coloring,coloring->w1);
1834: VecDuplicate(x1,&coloring->w2);
1835: PetscLogObjectParent(coloring,coloring->w2);
1836: VecDuplicate(x1,&coloring->w3);
1837: PetscLogObjectParent(coloring,coloring->w3);
1838: }
1839: w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
1841: MatSetUnfactored(J);
1842: PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);
1843: if (flg) {
1844: PetscLogInfo(coloring,"MatFDColoringApply_SeqAIJ: Not calling MatZeroEntries()n");
1845: } else {
1846: MatZeroEntries(J);
1847: }
1849: VecGetOwnershipRange(x1,&start,&end);
1850: VecGetSize(x1,&N);
1852: /*
1853: This is a horrible, horrible, hack. See DMMGComputeJacobian_Multigrid() it inproperly sets
1854: coloring->F for the coarser grids from the finest
1855: */
1856: if (coloring->F) {
1857: VecGetLocalSize(coloring->F,&m1);
1858: VecGetLocalSize(w1,&m2);
1859: if (m1 != m2) {
1860: coloring->F = 0;
1861: }
1862: }
1864: if (coloring->F) {
1865: w1 = coloring->F;
1866: coloring->F = 0;
1867: } else {
1868: (*f)(sctx,x1,w1,fctx);
1869: }
1871: /*
1872: Compute all the scale factors and share with other processors
1873: */
1874: VecGetArray(x1,&xx);xx = xx - start;
1875: VecGetArray(coloring->vscale,&vscale_array);vscale_array = vscale_array - start;
1876: for (k=0; k<coloring->ncolors; k++) {
1877: /*
1878: Loop over each column associated with color adding the
1879: perturbation to the vector w3.
1880: */
1881: for (l=0; l<coloring->ncolumns[k]; l++) {
1882: col = coloring->columns[k][l]; /* column of the matrix we are probing for */
1883: dx = xx[col];
1884: if (dx == 0.0) dx = 1.0;
1885: #if !defined(PETSC_USE_COMPLEX)
1886: if (dx < umin && dx >= 0.0) dx = umin;
1887: else if (dx < 0.0 && dx > -umin) dx = -umin;
1888: #else
1889: if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin;
1890: else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
1891: #endif
1892: dx *= epsilon;
1893: vscale_array[col] = 1.0/dx;
1894: }
1895: }
1896: vscale_array = vscale_array + start;VecRestoreArray(coloring->vscale,&vscale_array);
1897: VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);
1898: VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);
1900: /* VecView(coloring->vscale,PETSC_VIEWER_STDOUT_WORLD);
1901: VecView(x1,PETSC_VIEWER_STDOUT_WORLD);*/
1903: if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
1904: else vscaleforrow = coloring->columnsforrow;
1906: VecGetArray(coloring->vscale,&vscale_array);
1907: /*
1908: Loop over each color
1909: */
1910: for (k=0; k<coloring->ncolors; k++) {
1911: VecCopy(x1,w3);
1912: VecGetArray(w3,&w3_array);w3_array = w3_array - start;
1913: /*
1914: Loop over each column associated with color adding the
1915: perturbation to the vector w3.
1916: */
1917: for (l=0; l<coloring->ncolumns[k]; l++) {
1918: col = coloring->columns[k][l]; /* column of the matrix we are probing for */
1919: dx = xx[col];
1920: if (dx == 0.0) dx = 1.0;
1921: #if !defined(PETSC_USE_COMPLEX)
1922: if (dx < umin && dx >= 0.0) dx = umin;
1923: else if (dx < 0.0 && dx > -umin) dx = -umin;
1924: #else
1925: if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin;
1926: else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
1927: #endif
1928: dx *= epsilon;
1929: if (!PetscAbsScalar(dx)) SETERRQ(1,"Computed 0 differencing parameter");
1930: w3_array[col] += dx;
1931: }
1932: w3_array = w3_array + start; VecRestoreArray(w3,&w3_array);
1934: /*
1935: Evaluate function at x1 + dx (here dx is a vector of perturbations)
1936: */
1938: (*f)(sctx,w3,w2,fctx);
1939: VecAXPY(&mone,w1,w2);
1941: /*
1942: Loop over rows of vector, putting results into Jacobian matrix
1943: */
1944: VecGetArray(w2,&y);
1945: for (l=0; l<coloring->nrows[k]; l++) {
1946: row = coloring->rows[k][l];
1947: col = coloring->columnsforrow[k][l];
1948: y[row] *= vscale_array[vscaleforrow[k][l]];
1949: srow = row + start;
1950: ierr = MatSetValues_SeqAIJ(J,1,&srow,1,&col,y+row,INSERT_VALUES);
1951: }
1952: VecRestoreArray(w2,&y);
1953: }
1954: VecRestoreArray(coloring->vscale,&vscale_array);
1955: xx = xx + start; ierr = VecRestoreArray(x1,&xx);
1956: ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
1957: ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
1958: return(0);
1959: }
1961: #include petscblaslapack.h
1963: int MatAXPY_SeqAIJ(PetscScalar *a,Mat X,Mat Y,MatStructure str)
1964: {
1965: int ierr,one=1;
1966: Mat_SeqAIJ *x = (Mat_SeqAIJ *)X->data,*y = (Mat_SeqAIJ *)Y->data;
1969: if (str == SAME_NONZERO_PATTERN) {
1970: BLaxpy_(&x->nz,a,x->a,&one,y->a,&one);
1971: } else {
1972: MatAXPY_Basic(a,X,Y,str);
1973: }
1974: return(0);
1975: }
1978: /* -------------------------------------------------------------------*/
1979: static struct _MatOps MatOps_Values = {MatSetValues_SeqAIJ,
1980: MatGetRow_SeqAIJ,
1981: MatRestoreRow_SeqAIJ,
1982: MatMult_SeqAIJ,
1983: MatMultAdd_SeqAIJ,
1984: MatMultTranspose_SeqAIJ,
1985: MatMultTransposeAdd_SeqAIJ,
1986: MatSolve_SeqAIJ,
1987: MatSolveAdd_SeqAIJ,
1988: MatSolveTranspose_SeqAIJ,
1989: MatSolveTransposeAdd_SeqAIJ,
1990: MatLUFactor_SeqAIJ,
1991: 0,
1992: MatRelax_SeqAIJ,
1993: MatTranspose_SeqAIJ,
1994: MatGetInfo_SeqAIJ,
1995: MatEqual_SeqAIJ,
1996: MatGetDiagonal_SeqAIJ,
1997: MatDiagonalScale_SeqAIJ,
1998: MatNorm_SeqAIJ,
1999: 0,
2000: MatAssemblyEnd_SeqAIJ,
2001: MatCompress_SeqAIJ,
2002: MatSetOption_SeqAIJ,
2003: MatZeroEntries_SeqAIJ,
2004: MatZeroRows_SeqAIJ,
2005: MatLUFactorSymbolic_SeqAIJ,
2006: MatLUFactorNumeric_SeqAIJ,
2007: 0,
2008: 0,
2009: MatSetUpPreallocation_SeqAIJ,
2010: MatILUFactorSymbolic_SeqAIJ,
2011: 0,
2012: MatGetArray_SeqAIJ,
2013: MatRestoreArray_SeqAIJ,
2014: MatDuplicate_SeqAIJ,
2015: 0,
2016: 0,
2017: MatILUFactor_SeqAIJ,
2018: 0,
2019: MatAXPY_SeqAIJ,
2020: MatGetSubMatrices_SeqAIJ,
2021: MatIncreaseOverlap_SeqAIJ,
2022: MatGetValues_SeqAIJ,
2023: MatCopy_SeqAIJ,
2024: MatPrintHelp_SeqAIJ,
2025: MatScale_SeqAIJ,
2026: 0,
2027: 0,
2028: MatILUDTFactor_SeqAIJ,
2029: MatGetBlockSize_SeqAIJ,
2030: MatGetRowIJ_SeqAIJ,
2031: MatRestoreRowIJ_SeqAIJ,
2032: MatGetColumnIJ_SeqAIJ,
2033: MatRestoreColumnIJ_SeqAIJ,
2034: MatFDColoringCreate_SeqAIJ,
2035: 0,
2036: 0,
2037: MatPermute_SeqAIJ,
2038: 0,
2039: 0,
2040: MatDestroy_SeqAIJ,
2041: MatView_SeqAIJ,
2042: MatGetPetscMaps_Petsc,
2043: 0,
2044: 0,
2045: 0,
2046: 0,
2047: 0,
2048: 0,
2049: 0,
2050: 0,
2051: MatSetColoring_SeqAIJ,
2052: MatSetValuesAdic_SeqAIJ,
2053: MatSetValuesAdifor_SeqAIJ,
2054: MatFDColoringApply_SeqAIJ};
2056: EXTERN int MatUseSuperLU_SeqAIJ(Mat);
2057: EXTERN int MatUseEssl_SeqAIJ(Mat);
2058: EXTERN int MatUseLUSOL_SeqAIJ(Mat);
2059: EXTERN int MatUseMatlab_SeqAIJ(Mat);
2060: EXTERN int MatUseDXML_SeqAIJ(Mat);
2062: EXTERN_C_BEGIN
2064: int MatSeqAIJSetColumnIndices_SeqAIJ(Mat mat,int *indices)
2065: {
2066: Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
2067: int i,nz,n;
2071: nz = aij->maxnz;
2072: n = mat->n;
2073: for (i=0; i<nz; i++) {
2074: aij->j[i] = indices[i];
2075: }
2076: aij->nz = nz;
2077: for (i=0; i<n; i++) {
2078: aij->ilen[i] = aij->imax[i];
2079: }
2081: return(0);
2082: }
2083: EXTERN_C_END
2085: /*@
2086: MatSeqAIJSetColumnIndices - Set the column indices for all the rows
2087: in the matrix.
2089: Input Parameters:
2090: + mat - the SeqAIJ matrix
2091: - indices - the column indices
2093: Level: advanced
2095: Notes:
2096: This can be called if you have precomputed the nonzero structure of the
2097: matrix and want to provide it to the matrix object to improve the performance
2098: of the MatSetValues() operation.
2100: You MUST have set the correct numbers of nonzeros per row in the call to
2101: MatCreateSeqAIJ().
2103: MUST be called before any calls to MatSetValues();
2105: The indices should start with zero, not one.
2107: @*/
2108: int MatSeqAIJSetColumnIndices(Mat mat,int *indices)
2109: {
2110: int ierr,(*f)(Mat,int *);
2114: PetscObjectQueryFunction((PetscObject)mat,"MatSeqAIJSetColumnIndices_C",(void (**)(void))&f);
2115: if (f) {
2116: (*f)(mat,indices);
2117: } else {
2118: SETERRQ(1,"Wrong type of matrix to set column indices");
2119: }
2120: return(0);
2121: }
2123: /* ----------------------------------------------------------------------------------------*/
2125: EXTERN_C_BEGIN
2126: int MatStoreValues_SeqAIJ(Mat mat)
2127: {
2128: Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
2129: int nz = aij->i[mat->m]+aij->indexshift,ierr;
2132: if (aij->nonew != 1) {
2133: SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
2134: }
2136: /* allocate space for values if not already there */
2137: if (!aij->saved_values) {
2138: PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);
2139: }
2141: /* copy values over */
2142: PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));
2143: return(0);
2144: }
2145: EXTERN_C_END
2147: /*@
2148: MatStoreValues - Stashes a copy of the matrix values; this allows, for
2149: example, reuse of the linear part of a Jacobian, while recomputing the
2150: nonlinear portion.
2152: Collect on Mat
2154: Input Parameters:
2155: . mat - the matrix (currently on AIJ matrices support this option)
2157: Level: advanced
2159: Common Usage, with SNESSolve():
2160: $ Create Jacobian matrix
2161: $ Set linear terms into matrix
2162: $ Apply boundary conditions to matrix, at this time matrix must have
2163: $ final nonzero structure (i.e. setting the nonlinear terms and applying
2164: $ boundary conditions again will not change the nonzero structure
2165: $ MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS);
2166: $ MatStoreValues(mat);
2167: $ Call SNESSetJacobian() with matrix
2168: $ In your Jacobian routine
2169: $ MatRetrieveValues(mat);
2170: $ Set nonlinear terms in matrix
2171:
2172: Common Usage without SNESSolve(), i.e. when you handle nonlinear solve yourself:
2173: $ // build linear portion of Jacobian
2174: $ MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS);
2175: $ MatStoreValues(mat);
2176: $ loop over nonlinear iterations
2177: $ MatRetrieveValues(mat);
2178: $ // call MatSetValues(mat,...) to set nonliner portion of Jacobian
2179: $ // call MatAssemblyBegin/End() on matrix
2180: $ Solve linear system with Jacobian
2181: $ endloop
2183: Notes:
2184: Matrix must already be assemblied before calling this routine
2185: Must set the matrix option MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS); before
2186: calling this routine.
2188: .seealso: MatRetrieveValues()
2190: @*/
2191: int MatStoreValues(Mat mat)
2192: {
2193: int ierr,(*f)(Mat);
2197: if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2198: if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2200: PetscObjectQueryFunction((PetscObject)mat,"MatStoreValues_C",(void (**)(void))&f);
2201: if (f) {
2202: (*f)(mat);
2203: } else {
2204: SETERRQ(1,"Wrong type of matrix to store values");
2205: }
2206: return(0);
2207: }
2209: EXTERN_C_BEGIN
2210: int MatRetrieveValues_SeqAIJ(Mat mat)
2211: {
2212: Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
2213: int nz = aij->i[mat->m]+aij->indexshift,ierr;
2216: if (aij->nonew != 1) {
2217: SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
2218: }
2219: if (!aij->saved_values) {
2220: SETERRQ(1,"Must call MatStoreValues(A);first");
2221: }
2223: /* copy values over */
2224: PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));
2225: return(0);
2226: }
2227: EXTERN_C_END
2229: /*@
2230: MatRetrieveValues - Retrieves the copy of the matrix values; this allows, for
2231: example, reuse of the linear part of a Jacobian, while recomputing the
2232: nonlinear portion.
2234: Collect on Mat
2236: Input Parameters:
2237: . mat - the matrix (currently on AIJ matrices support this option)
2239: Level: advanced
2241: .seealso: MatStoreValues()
2243: @*/
2244: int MatRetrieveValues(Mat mat)
2245: {
2246: int ierr,(*f)(Mat);
2250: if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2251: if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2253: PetscObjectQueryFunction((PetscObject)mat,"MatRetrieveValues_C",(void (**)(void))&f);
2254: if (f) {
2255: (*f)(mat);
2256: } else {
2257: SETERRQ(1,"Wrong type of matrix to retrieve values");
2258: }
2259: return(0);
2260: }
2262: /*
2263: This allows SeqAIJ matrices to be passed to the matlab engine
2264: */
2265: #if defined(PETSC_HAVE_MATLAB_ENGINE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
2266: #include "engine.h" /* Matlab include file */
2267: #include "mex.h" /* Matlab include file */
2268: EXTERN_C_BEGIN
2269: int MatMatlabEnginePut_SeqAIJ(PetscObject obj,void *engine)
2270: {
2271: int ierr,i,*ai,*aj;
2272: Mat B = (Mat)obj;
2273: PetscScalar *array;
2274: mxArray *mat;
2275: Mat_SeqAIJ *aij = (Mat_SeqAIJ*)B->data;
2278: mat = mxCreateSparse(B->n,B->m,aij->nz,mxREAL);
2279: PetscMemcpy(mxGetPr(mat),aij->a,aij->nz*sizeof(PetscScalar));
2280: /* Matlab stores by column, not row so we pass in the transpose of the matrix */
2281: PetscMemcpy(mxGetIr(mat),aij->j,aij->nz*sizeof(int));
2282: PetscMemcpy(mxGetJc(mat),aij->i,(B->m+1)*sizeof(int));
2284: /* Matlab indices start at 0 for sparse (what a surprise) */
2285: if (aij->indexshift) {
2286: for (i=0; i<B->m+1; i++) {
2287: ai[i]--;
2288: }
2289: for (i=0; i<aij->nz; i++) {
2290: aj[i]--;
2291: }
2292: }
2293: PetscObjectName(obj);
2294: mxSetName(mat,obj->name);
2295: engPutArray((Engine *)engine,mat);
2296: return(0);
2297: }
2298: EXTERN_C_END
2300: EXTERN_C_BEGIN
2301: int MatMatlabEngineGet_SeqAIJ(PetscObject obj,void *engine)
2302: {
2303: int ierr,ii;
2304: Mat mat = (Mat)obj;
2305: Mat_SeqAIJ *aij = (Mat_SeqAIJ*)mat->data;
2306: mxArray *mmat;
2309: PetscFree(aij->a);
2310: aij->indexshift = 0;
2312: mmat = engGetArray((Engine *)engine,obj->name);
2314: aij->nz = (mxGetJc(mmat))[mat->m];
2315: ierr = PetscMalloc(aij->nz*(sizeof(int)+sizeof(PetscScalar))+(mat->m+1)*sizeof(int),&aij->a);
2316: aij->j = (int*)(aij->a + aij->nz);
2317: aij->i = aij->j + aij->nz;
2318: aij->singlemalloc = PETSC_TRUE;
2319: aij->freedata = PETSC_TRUE;
2321: PetscMemcpy(aij->a,mxGetPr(mmat),aij->nz*sizeof(PetscScalar));
2322: /* Matlab stores by column, not row so we pass in the transpose of the matrix */
2323: PetscMemcpy(aij->j,mxGetIr(mmat),aij->nz*sizeof(int));
2324: PetscMemcpy(aij->i,mxGetJc(mmat),(mat->m+1)*sizeof(int));
2326: for (ii=0; ii<mat->m; ii++) {
2327: aij->ilen[ii] = aij->imax[ii] = aij->i[ii+1] - aij->i[ii];
2328: }
2330: MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
2331: MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);
2333: return(0);
2334: }
2335: EXTERN_C_END
2336: #endif
2338: /* --------------------------------------------------------------------------------*/
2339: /*@C
2340: MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format
2341: (the default parallel PETSc format). For good matrix assembly performance
2342: the user should preallocate the matrix storage by setting the parameter nz
2343: (or the array nnz). By setting these parameters accurately, performance
2344: during matrix assembly can be increased by more than a factor of 50.
2346: Collective on MPI_Comm
2348: Input Parameters:
2349: + comm - MPI communicator, set to PETSC_COMM_SELF
2350: . m - number of rows
2351: . n - number of columns
2352: . nz - number of nonzeros per row (same for all rows)
2353: - nnz - array containing the number of nonzeros in the various rows
2354: (possibly different for each row) or PETSC_NULL
2356: Output Parameter:
2357: . A - the matrix
2359: Notes:
2360: The AIJ format (also called the Yale sparse matrix format or
2361: compressed row storage), is fully compatible with standard Fortran 77
2362: storage. That is, the stored row and column indices can begin at
2363: either one (as in Fortran) or zero. See the users' manual for details.
2365: Specify the preallocated storage with either nz or nnz (not both).
2366: Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
2367: allocation. For large problems you MUST preallocate memory or you
2368: will get TERRIBLE performance, see the users' manual chapter on matrices.
2370: By default, this format uses inodes (identical nodes) when possible, to
2371: improve numerical efficiency of matrix-vector products and solves. We
2372: search for consecutive rows with the same nonzero structure, thereby
2373: reusing matrix information to achieve increased efficiency.
2375: Options Database Keys:
2376: + -mat_aij_no_inode - Do not use inodes
2377: . -mat_aij_inode_limit <limit> - Sets inode limit (max limit=5)
2378: - -mat_aij_oneindex - Internally use indexing starting at 1
2379: rather than 0. Note that when calling MatSetValues(),
2380: the user still MUST index entries starting at 0!
2382: Level: intermediate
2384: .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues(), MatSeqAIJSetColumnIndices(), MatCreateSeqAIJWithArrays()
2386: @*/
2387: int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz,Mat *A)
2388: {
2392: MatCreate(comm,m,n,m,n,A);
2393: MatSetType(*A,MATSEQAIJ);
2394: MatSeqAIJSetPreallocation(*A,nz,nnz);
2395: return(0);
2396: }
2398: #define SKIP_ALLOCATION -4
2400: /*@C
2401: MatSeqAIJSetPreallocation - For good matrix assembly performance
2402: the user should preallocate the matrix storage by setting the parameter nz
2403: (or the array nnz). By setting these parameters accurately, performance
2404: during matrix assembly can be increased by more than a factor of 50.
2406: Collective on MPI_Comm
2408: Input Parameters:
2409: + comm - MPI communicator, set to PETSC_COMM_SELF
2410: . m - number of rows
2411: . n - number of columns
2412: . nz - number of nonzeros per row (same for all rows)
2413: - nnz - array containing the number of nonzeros in the various rows
2414: (possibly different for each row) or PETSC_NULL
2416: Output Parameter:
2417: . A - the matrix
2419: Notes:
2420: The AIJ format (also called the Yale sparse matrix format or
2421: compressed row storage), is fully compatible with standard Fortran 77
2422: storage. That is, the stored row and column indices can begin at
2423: either one (as in Fortran) or zero. See the users' manual for details.
2425: Specify the preallocated storage with either nz or nnz (not both).
2426: Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
2427: allocation. For large problems you MUST preallocate memory or you
2428: will get TERRIBLE performance, see the users' manual chapter on matrices.
2430: By default, this format uses inodes (identical nodes) when possible, to
2431: improve numerical efficiency of matrix-vector products and solves. We
2432: search for consecutive rows with the same nonzero structure, thereby
2433: reusing matrix information to achieve increased efficiency.
2435: Options Database Keys:
2436: + -mat_aij_no_inode - Do not use inodes
2437: . -mat_aij_inode_limit <limit> - Sets inode limit (max limit=5)
2438: - -mat_aij_oneindex - Internally use indexing starting at 1
2439: rather than 0. Note that when calling MatSetValues(),
2440: the user still MUST index entries starting at 0!
2442: Level: intermediate
2444: .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues(), MatSeqAIJSetColumnIndices(), MatCreateSeqAIJWithArrays()
2446: @*/
2447: int MatSeqAIJSetPreallocation(Mat B,int nz,int *nnz)
2448: {
2449: Mat_SeqAIJ *b;
2450: int i,len=0,ierr;
2451: PetscTruth flg2,skipallocation = PETSC_FALSE;
2454: PetscTypeCompare((PetscObject)B,MATSEQAIJ,&flg2);
2455: if (!flg2) return(0);
2456:
2457: if (nz == SKIP_ALLOCATION) {
2458: skipallocation = PETSC_TRUE;
2459: nz = 0;
2460: }
2462: if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
2463: if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz);
2464: if (nnz) {
2465: for (i=0; i<B->m; i++) {
2466: if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]);
2467: 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);
2468: }
2469: }
2471: B->preallocated = PETSC_TRUE;
2472: b = (Mat_SeqAIJ*)B->data;
2474: PetscMalloc((B->m+1)*sizeof(int),&b->imax);
2475: if (!nnz) {
2476: if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 10;
2477: else if (nz <= 0) nz = 1;
2478: for (i=0; i<B->m; i++) b->imax[i] = nz;
2479: nz = nz*B->m;
2480: } else {
2481: nz = 0;
2482: for (i=0; i<B->m; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
2483: }
2485: if (!skipallocation) {
2486: /* allocate the matrix space */
2487: len = nz*(sizeof(int) + sizeof(PetscScalar)) + (B->m+1)*sizeof(int);
2488: ierr = PetscMalloc(len,&b->a);
2489: b->j = (int*)(b->a + nz);
2490: ierr = PetscMemzero(b->j,nz*sizeof(int));
2491: b->i = b->j + nz;
2492: b->i[0] = -b->indexshift;
2493: for (i=1; i<B->m+1; i++) {
2494: b->i[i] = b->i[i-1] + b->imax[i-1];
2495: }
2496: b->singlemalloc = PETSC_TRUE;
2497: b->freedata = PETSC_TRUE;
2498: } else {
2499: b->freedata = PETSC_FALSE;
2500: }
2502: /* b->ilen will count nonzeros in each row so far. */
2503: PetscMalloc((B->m+1)*sizeof(int),&b->ilen);
2504: PetscLogObjectMemory(B,len+2*(B->m+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqAIJ));
2505: for (i=0; i<B->m; i++) { b->ilen[i] = 0;}
2507: b->nz = 0;
2508: b->maxnz = nz;
2509: B->info.nz_unneeded = (double)b->maxnz;
2510: return(0);
2511: }
2513: EXTERN int RegisterApplyPtAPRoutines_Private(Mat);
2515: EXTERN_C_BEGIN
2516: int MatCreate_SeqAIJ(Mat B)
2517: {
2518: Mat_SeqAIJ *b;
2519: int ierr,size;
2520: PetscTruth flg;
2523: MPI_Comm_size(B->comm,&size);
2524: if (size > 1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Comm must be of size 1");
2526: B->m = B->M = PetscMax(B->m,B->M);
2527: B->n = B->N = PetscMax(B->n,B->N);
2529: PetscNew(Mat_SeqAIJ,&b);
2530: B->data = (void*)b;
2531: PetscMemzero(b,sizeof(Mat_SeqAIJ));
2532: PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
2533: B->factor = 0;
2534: B->lupivotthreshold = 1.0;
2535: B->mapping = 0;
2536: PetscOptionsGetReal(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold,PETSC_NULL);
2537: PetscOptionsHasName(PETSC_NULL,"-pc_ilu_preserve_row_sums",&b->ilu_preserve_row_sums);
2538: b->row = 0;
2539: b->col = 0;
2540: b->icol = 0;
2541: b->indexshift = 0;
2542: b->reallocs = 0;
2543: PetscOptionsHasName(PETSC_NULL,"-mat_aij_oneindex",&flg);
2544: if (flg) b->indexshift = -1;
2545:
2546: PetscMapCreateMPI(B->comm,B->m,B->m,&B->rmap);
2547: PetscMapCreateMPI(B->comm,B->n,B->n,&B->cmap);
2549: b->sorted = PETSC_FALSE;
2550: b->ignorezeroentries = PETSC_FALSE;
2551: b->roworiented = PETSC_TRUE;
2552: b->nonew = 0;
2553: b->diag = 0;
2554: b->solve_work = 0;
2555: B->spptr = 0;
2556: b->inode.use = PETSC_TRUE;
2557: b->inode.node_count = 0;
2558: b->inode.size = 0;
2559: b->inode.limit = 5;
2560: b->inode.max_limit = 5;
2561: b->saved_values = 0;
2562: b->idiag = 0;
2563: b->ssor = 0;
2564: b->keepzeroedrows = PETSC_FALSE;
2566: PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);
2568: #if defined(PETSC_HAVE_SUPERLU)
2569: PetscOptionsHasName(PETSC_NULL,"-mat_aij_superlu",&flg);
2570: if (flg) { MatUseSuperLU_SeqAIJ(B); }
2571: #endif
2572: PetscOptionsHasName(PETSC_NULL,"-mat_aij_essl",&flg);
2573: if (flg) { MatUseEssl_SeqAIJ(B); }
2574: PetscOptionsHasName(PETSC_NULL,"-mat_aij_lusol",&flg);
2575: if (flg) { MatUseLUSOL_SeqAIJ(B); }
2576: PetscOptionsHasName(PETSC_NULL,"-mat_aij_matlab",&flg);
2577: if (flg) {MatUseMatlab_SeqAIJ(B);}
2578: PetscOptionsHasName(PETSC_NULL,"-mat_aij_dxml",&flg);
2579: if (flg) {
2580: if (!b->indexshift) SETERRQ(PETSC_ERR_LIB,"need -mat_aij_oneindex with -mat_aij_dxml");
2581: MatUseDXML_SeqAIJ(B);
2582: }
2583: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqAIJSetColumnIndices_C",
2584: "MatSeqAIJSetColumnIndices_SeqAIJ",
2585: MatSeqAIJSetColumnIndices_SeqAIJ);
2586: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
2587: "MatStoreValues_SeqAIJ",
2588: MatStoreValues_SeqAIJ);
2589: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
2590: "MatRetrieveValues_SeqAIJ",
2591: MatRetrieveValues_SeqAIJ);
2592: #if defined(PETSC_HAVE_MATLAB_ENGINE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
2593: PetscObjectComposeFunctionDynamic((PetscObject)B,"PetscMatlabEnginePut_C","MatMatlabEnginePut_SeqAIJ",MatMatlabEnginePut_SeqAIJ);
2594: PetscObjectComposeFunctionDynamic((PetscObject)B,"PetscMatlabEngineGet_C","MatMatlabEngineGet_SeqAIJ",MatMatlabEngineGet_SeqAIJ);
2595: #endif
2596: RegisterApplyPtAPRoutines_Private(B);
2597: return(0);
2598: }
2599: EXTERN_C_END
2601: int MatDuplicate_SeqAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
2602: {
2603: Mat C;
2604: Mat_SeqAIJ *c,*a = (Mat_SeqAIJ*)A->data;
2605: int i,len,m = A->m,shift = a->indexshift,ierr;
2608: *B = 0;
2609: MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);
2610: MatSetType(C,MATSEQAIJ);
2611: c = (Mat_SeqAIJ*)C->data;
2613: C->factor = A->factor;
2614: c->row = 0;
2615: c->col = 0;
2616: c->icol = 0;
2617: c->indexshift = shift;
2618: c->keepzeroedrows = a->keepzeroedrows;
2619: C->assembled = PETSC_TRUE;
2621: C->M = A->m;
2622: C->N = A->n;
2624: PetscMalloc((m+1)*sizeof(int),&c->imax);
2625: PetscMalloc((m+1)*sizeof(int),&c->ilen);
2626: for (i=0; i<m; i++) {
2627: c->imax[i] = a->imax[i];
2628: c->ilen[i] = a->ilen[i];
2629: }
2631: /* allocate the matrix space */
2632: c->singlemalloc = PETSC_TRUE;
2633: len = (m+1)*sizeof(int)+(a->i[m])*(sizeof(PetscScalar)+sizeof(int));
2634: ierr = PetscMalloc(len,&c->a);
2635: c->j = (int*)(c->a + a->i[m] + shift);
2636: c->i = c->j + a->i[m] + shift;
2637: PetscMemcpy(c->i,a->i,(m+1)*sizeof(int));
2638: if (m > 0) {
2639: PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int));
2640: if (cpvalues == MAT_COPY_VALUES) {
2641: PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(PetscScalar));
2642: } else {
2643: PetscMemzero(c->a,(a->i[m]+shift)*sizeof(PetscScalar));
2644: }
2645: }
2647: PetscLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqAIJ));
2648: c->sorted = a->sorted;
2649: c->roworiented = a->roworiented;
2650: c->nonew = a->nonew;
2651: c->ilu_preserve_row_sums = a->ilu_preserve_row_sums;
2652: c->saved_values = 0;
2653: c->idiag = 0;
2654: c->ssor = 0;
2655: c->ignorezeroentries = a->ignorezeroentries;
2656: c->freedata = PETSC_TRUE;
2658: if (a->diag) {
2659: PetscMalloc((m+1)*sizeof(int),&c->diag);
2660: PetscLogObjectMemory(C,(m+1)*sizeof(int));
2661: for (i=0; i<m; i++) {
2662: c->diag[i] = a->diag[i];
2663: }
2664: } else c->diag = 0;
2665: c->inode.use = a->inode.use;
2666: c->inode.limit = a->inode.limit;
2667: c->inode.max_limit = a->inode.max_limit;
2668: if (a->inode.size){
2669: ierr = PetscMalloc((m+1)*sizeof(int),&c->inode.size);
2670: c->inode.node_count = a->inode.node_count;
2671: ierr = PetscMemcpy(c->inode.size,a->inode.size,(m+1)*sizeof(int));
2672: } else {
2673: c->inode.size = 0;
2674: c->inode.node_count = 0;
2675: }
2676: c->nz = a->nz;
2677: c->maxnz = a->maxnz;
2678: c->solve_work = 0;
2679: C->spptr = 0; /* Dangerous -I'm throwing away a->spptr */
2680: C->preallocated = PETSC_TRUE;
2682: *B = C;
2683: PetscFListDuplicate(A->qlist,&C->qlist);
2684: return(0);
2685: }
2687: EXTERN_C_BEGIN
2688: int MatLoad_SeqAIJ(PetscViewer viewer,MatType type,Mat *A)
2689: {
2690: Mat_SeqAIJ *a;
2691: Mat B;
2692: int i,nz,ierr,fd,header[4],size,*rowlengths = 0,M,N,shift;
2693: MPI_Comm comm;
2694:
2696: PetscObjectGetComm((PetscObject)viewer,&comm);
2697: MPI_Comm_size(comm,&size);
2698: if (size > 1) SETERRQ(PETSC_ERR_ARG_SIZ,"view must have one processor");
2699: PetscViewerBinaryGetDescriptor(viewer,&fd);
2700: PetscBinaryRead(fd,header,4,PETSC_INT);
2701: if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object in file");
2702: M = header[1]; N = header[2]; nz = header[3];
2704: if (nz < 0) {
2705: SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format on disk,cannot load as SeqAIJ");
2706: }
2708: /* read in row lengths */
2709: PetscMalloc(M*sizeof(int),&rowlengths);
2710: PetscBinaryRead(fd,rowlengths,M,PETSC_INT);
2712: /* create our matrix */
2713: MatCreateSeqAIJ(comm,M,N,0,rowlengths,A);
2714: B = *A;
2715: a = (Mat_SeqAIJ*)B->data;
2716: shift = a->indexshift;
2718: /* read in column indices and adjust for Fortran indexing*/
2719: PetscBinaryRead(fd,a->j,nz,PETSC_INT);
2720: if (shift) {
2721: for (i=0; i<nz; i++) {
2722: a->j[i] += 1;
2723: }
2724: }
2726: /* read in nonzero values */
2727: PetscBinaryRead(fd,a->a,nz,PETSC_SCALAR);
2729: /* set matrix "i" values */
2730: a->i[0] = -shift;
2731: for (i=1; i<= M; i++) {
2732: a->i[i] = a->i[i-1] + rowlengths[i-1];
2733: a->ilen[i-1] = rowlengths[i-1];
2734: }
2735: PetscFree(rowlengths);
2737: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
2738: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
2739: return(0);
2740: }
2741: EXTERN_C_END
2743: int MatEqual_SeqAIJ(Mat A,Mat B,PetscTruth* flg)
2744: {
2745: Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data,*b = (Mat_SeqAIJ *)B->data;
2746: int ierr;
2747: PetscTruth flag;
2750: PetscTypeCompare((PetscObject)B,MATSEQAIJ,&flag);
2751: if (!flag) SETERRQ(PETSC_ERR_ARG_INCOMP,"Matrices must be same type");
2753: /* If the matrix dimensions are not equal,or no of nonzeros or shift */
2754: if ((A->m != B->m) || (A->n != B->n) ||(a->nz != b->nz)|| (a->indexshift != b->indexshift)) {
2755: *flg = PETSC_FALSE;
2756: return(0);
2757: }
2758:
2759: /* if the a->i are the same */
2760: PetscMemcmp(a->i,b->i,(A->m+1)*sizeof(int),flg);
2761: if (*flg == PETSC_FALSE) return(0);
2762:
2763: /* if a->j are the same */
2764: PetscMemcmp(a->j,b->j,(a->nz)*sizeof(int),flg);
2765: if (*flg == PETSC_FALSE) return(0);
2766:
2767: /* if a->a are the same */
2768: PetscMemcmp(a->a,b->a,(a->nz)*sizeof(PetscScalar),flg);
2770: return(0);
2771:
2772: }
2774: /*@C
2775: MatCreateSeqAIJWithArrays - Creates an sequential AIJ matrix using matrix elements (in CSR format)
2776: provided by the user.
2778: Coolective on MPI_Comm
2780: Input Parameters:
2781: + comm - must be an MPI communicator of size 1
2782: . m - number of rows
2783: . n - number of columns
2784: . i - row indices
2785: . j - column indices
2786: - a - matrix values
2788: Output Parameter:
2789: . mat - the matrix
2791: Level: intermediate
2793: Notes:
2794: The i, j, and a arrays are not copied by this routine, the user must free these arrays
2795: once the matrix is destroyed
2797: You cannot set new nonzero locations into this matrix, that will generate an error.
2799: The i and j indices can be either 0- or 1 based
2801: .seealso: MatCreate(), MatCreateMPIAIJ(), MatCreateSeqAIJ()
2803: @*/
2804: int MatCreateSeqAIJWithArrays(MPI_Comm comm,int m,int n,int* i,int*j,PetscScalar *a,Mat *mat)
2805: {
2806: int ierr,ii;
2807: Mat_SeqAIJ *aij;
2810: MatCreateSeqAIJ(comm,m,n,SKIP_ALLOCATION,0,mat);
2811: aij = (Mat_SeqAIJ*)(*mat)->data;
2813: if (i[0] == 1) {
2814: aij->indexshift = -1;
2815: } else if (i[0]) {
2816: SETERRQ(1,"i (row indices) do not start with 0 or 1");
2817: }
2818: aij->i = i;
2819: aij->j = j;
2820: aij->a = a;
2821: aij->singlemalloc = PETSC_FALSE;
2822: aij->nonew = -1; /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/
2823: aij->freedata = PETSC_FALSE;
2825: for (ii=0; ii<m; ii++) {
2826: aij->ilen[ii] = aij->imax[ii] = i[ii+1] - i[ii];
2827: #if defined(PETSC_USE_BOPT_g)
2828: if (i[ii+1] - i[ii] < 0) SETERRQ2(1,"Negative row length in i (row indices) row = %d length = %d",ii,i[ii+1] - i[ii]);
2829: #endif
2830: }
2831: #if defined(PETSC_USE_BOPT_g)
2832: for (ii=0; ii<aij->i[m]; ii++) {
2833: if (j[ii] < -aij->indexshift) SETERRQ2(1,"Negative column index at location = %d index = %d",ii,j[ii]);
2834: if (j[ii] > n - 1 -aij->indexshift) SETERRQ2(1,"Column index to large at location = %d index = %d",ii,j[ii]);
2835: }
2836: #endif
2838: /* changes indices to start at 0 */
2839: if (i[0]) {
2840: aij->indexshift = 0;
2841: for (ii=0; ii<m; ii++) {
2842: i[ii]--;
2843: }
2844: for (ii=0; ii<i[m]; ii++) {
2845: j[ii]--;
2846: }
2847: }
2849: MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
2850: MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
2851: return(0);
2852: }
2854: int MatSetColoring_SeqAIJ(Mat A,ISColoring coloring)
2855: {
2856: int ierr;
2857: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
2860: if (coloring->ctype == IS_COLORING_LOCAL) {
2861: ierr = ISColoringReference(coloring);
2862: a->coloring = coloring;
2863: } else if (coloring->ctype == IS_COLORING_GHOSTED) {
2864: int *colors,i,*larray;
2865: ISColoring ocoloring;
2867: /* set coloring for diagonal portion */
2868: PetscMalloc((A->n+1)*sizeof(int),&larray);
2869: for (i=0; i<A->n; i++) {
2870: larray[i] = i;
2871: }
2872: ISGlobalToLocalMappingApply(A->mapping,IS_GTOLM_MASK,A->n,larray,PETSC_NULL,larray);
2873: PetscMalloc((A->n+1)*sizeof(int),&colors);
2874: for (i=0; i<A->n; i++) {
2875: colors[i] = coloring->colors[larray[i]];
2876: }
2877: PetscFree(larray);
2878: ISColoringCreate(PETSC_COMM_SELF,A->n,colors,&ocoloring);
2879: a->coloring = ocoloring;
2880: }
2881: return(0);
2882: }
2884: #if defined(PETSC_HAVE_ADIC) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
2885: EXTERN_C_BEGIN
2886: #include "adic/ad_utils.h"
2887: EXTERN_C_END
2889: int MatSetValuesAdic_SeqAIJ(Mat A,void *advalues)
2890: {
2891: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
2892: int m = A->m,*ii = a->i,*jj = a->j,nz,i,*color,j,nlen;
2893: PetscScalar *v = a->a,*values;
2894: char *cadvalues = (char *)advalues;
2897: if (!a->coloring) SETERRQ(1,"Coloring not set for matrix");
2898: nlen = PetscADGetDerivTypeSize();
2899: color = a->coloring->colors;
2900: /* loop over rows */
2901: for (i=0; i<m; i++) {
2902: nz = ii[i+1] - ii[i];
2903: /* loop over columns putting computed value into matrix */
2904: values = PetscADGetGradArray(cadvalues);
2905: for (j=0; j<nz; j++) {
2906: *v++ = values[color[*jj++]];
2907: }
2908: cadvalues += nlen; /* jump to next row of derivatives */
2909: }
2910: return(0);
2911: }
2913: #else
2915: int MatSetValuesAdic_SeqAIJ(Mat A,void *advalues)
2916: {
2918: SETERRQ(1,"PETSc installed without ADIC");
2919: }
2921: #endif
2923: int MatSetValuesAdifor_SeqAIJ(Mat A,int nl,void *advalues)
2924: {
2925: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
2926: int m = A->m,*ii = a->i,*jj = a->j,nz,i,*color,j;
2927: PetscScalar *v = a->a,*values = (PetscScalar *)advalues;
2930: if (!a->coloring) SETERRQ(1,"Coloring not set for matrix");
2931: color = a->coloring->colors;
2932: /* loop over rows */
2933: for (i=0; i<m; i++) {
2934: nz = ii[i+1] - ii[i];
2935: /* loop over columns putting computed value into matrix */
2936: for (j=0; j<nz; j++) {
2937: *v++ = values[color[*jj++]];
2938: }
2939: values += nl; /* jump to next row of derivatives */
2940: }
2941: return(0);
2942: }