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;
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: PetscScalar sum;
860: #endif
862: #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
863: #pragma disjoint(*x,*y,*v)
864: #endif
867: VecGetArray(xx,&x);
868: VecGetArray(yy,&y);
869: x = x + shift; /* shift for Fortran start by 1 indexing */
870: idx = a->j;
871: v = a->a;
872: ii = a->i;
873: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
874: fortranmultaij_(&m,x,ii,idx+shift,v+shift,y);
875: #else
876: v += shift; /* shift for Fortran start by 1 indexing */
877: idx += shift;
878: for (i=0; i<m; i++) {
879: jrow = ii[i];
880: n = ii[i+1] - jrow;
881: sum = 0.0;
882: for (j=0; j<n; j++) {
883: sum += v[jrow]*x[idx[jrow]]; jrow++;
884: }
885: y[i] = sum;
886: }
887: #endif
888: PetscLogFlops(2*a->nz - m);
889: VecRestoreArray(xx,&x);
890: VecRestoreArray(yy,&y);
891: return(0);
892: }
894: int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
895: {
896: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
897: PetscScalar *x,*y,*z,*v;
898: int ierr,m = A->m,*idx,shift = a->indexshift,*ii;
899: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
900: int n,i,jrow,j;
901: PetscScalar sum;
902: #endif
905: VecGetArray(xx,&x);
906: VecGetArray(yy,&y);
907: if (zz != yy) {
908: VecGetArray(zz,&z);
909: } else {
910: z = y;
911: }
912: x = x + shift; /* shift for Fortran start by 1 indexing */
913: idx = a->j;
914: v = a->a;
915: ii = a->i;
916: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
917: fortranmultaddaij_(&m,x,ii,idx+shift,v+shift,y,z);
918: #else
919: v += shift; /* shift for Fortran start by 1 indexing */
920: idx += shift;
921: for (i=0; i<m; i++) {
922: jrow = ii[i];
923: n = ii[i+1] - jrow;
924: sum = y[i];
925: for (j=0; j<n; j++) {
926: sum += v[jrow]*x[idx[jrow]]; jrow++;
927: }
928: z[i] = sum;
929: }
930: #endif
931: PetscLogFlops(2*a->nz);
932: VecRestoreArray(xx,&x);
933: VecRestoreArray(yy,&y);
934: if (zz != yy) {
935: VecRestoreArray(zz,&z);
936: }
937: return(0);
938: }
940: /*
941: Adds diagonal pointers to sparse matrix structure.
942: */
943: int MatMarkDiagonal_SeqAIJ(Mat A)
944: {
945: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
946: int i,j,*diag,m = A->m,shift = a->indexshift,ierr;
949: if (a->diag) return(0);
951: PetscMalloc((m+1)*sizeof(int),&diag);
952: PetscLogObjectMemory(A,(m+1)*sizeof(int));
953: for (i=0; i<A->m; i++) {
954: diag[i] = a->i[i+1];
955: for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
956: if (a->j[j]+shift == i) {
957: diag[i] = j - shift;
958: break;
959: }
960: }
961: }
962: a->diag = diag;
963: return(0);
964: }
966: /*
967: Checks for missing diagonals
968: */
969: int MatMissingDiagonal_SeqAIJ(Mat A)
970: {
971: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
972: int *diag,*jj = a->j,i,shift = a->indexshift,ierr;
975: MatMarkDiagonal_SeqAIJ(A);
976: diag = a->diag;
977: for (i=0; i<A->m; i++) {
978: if (jj[diag[i]+shift] != i-shift) {
979: SETERRQ1(1,"Matrix is missing diagonal number %d",i);
980: }
981: }
982: return(0);
983: }
985: int MatRelax_SeqAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,int its,int lits,Vec xx)
986: {
987: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
988: PetscScalar *x,*b,*bs, d,*xs,sum,*v = a->a,*t=0,scale,*ts,*xb,*idiag=0;
989: int ierr,*idx,*diag,n = A->n,m = A->m,i,shift = a->indexshift;
992: its = its*lits;
993: if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %d and local its %d both positive",its,lits);
995: VecGetArray(xx,&x);
996: if (xx != bb) {
997: VecGetArray(bb,&b);
998: } else {
999: b = x;
1000: }
1002: if (!a->diag) {MatMarkDiagonal_SeqAIJ(A);}
1003: diag = a->diag;
1004: xs = x + shift; /* shifted by one for index start of a or a->j*/
1005: if (flag == SOR_APPLY_UPPER) {
1006: /* apply (U + D/omega) to the vector */
1007: bs = b + shift;
1008: for (i=0; i<m; i++) {
1009: d = fshift + a->a[diag[i] + shift];
1010: n = a->i[i+1] - diag[i] - 1;
1011: PetscLogFlops(2*n-1);
1012: idx = a->j + diag[i] + (!shift);
1013: v = a->a + diag[i] + (!shift);
1014: sum = b[i]*d/omega;
1015: SPARSEDENSEDOT(sum,bs,v,idx,n);
1016: x[i] = sum;
1017: }
1018: VecRestoreArray(xx,&x);
1019: if (bb != xx) {VecRestoreArray(bb,&b);}
1020: return(0);
1021: }
1023: /* setup workspace for Eisenstat */
1024: if (flag & SOR_EISENSTAT) {
1025: if (!a->idiag) {
1026: ierr = PetscMalloc(2*m*sizeof(PetscScalar),&a->idiag);
1027: a->ssor = a->idiag + m;
1028: v = a->a;
1029: for (i=0; i<m; i++) { a->idiag[i] = 1.0/v[diag[i]];}
1030: }
1031: t = a->ssor;
1032: idiag = a->idiag;
1033: }
1034: /* Let A = L + U + D; where L is lower trianglar,
1035: U is upper triangular, E is diagonal; This routine applies
1037: (L + E)^{-1} A (U + E)^{-1}
1039: to a vector efficiently using Eisenstat's trick. This is for
1040: the case of SSOR preconditioner, so E is D/omega where omega
1041: is the relaxation factor.
1042: */
1044: if (flag == SOR_APPLY_LOWER) {
1045: SETERRQ(PETSC_ERR_SUP,"SOR_APPLY_LOWER is not implemented");
1046: } else if ((flag & SOR_EISENSTAT) && omega == 1.0 && shift == 0 && fshift == 0.0) {
1047: /* special case for omega = 1.0 saves flops and some integer ops */
1048: PetscScalar *v2;
1049:
1050: v2 = a->a;
1051: /* x = (E + U)^{-1} b */
1052: for (i=m-1; i>=0; i--) {
1053: n = a->i[i+1] - diag[i] - 1;
1054: idx = a->j + diag[i] + 1;
1055: v = a->a + diag[i] + 1;
1056: sum = b[i];
1057: SPARSEDENSEMDOT(sum,xs,v,idx,n);
1058: x[i] = sum*idiag[i];
1060: /* t = b - (2*E - D)x */
1061: t[i] = b[i] - (v2[diag[i]])*x[i];
1062: }
1064: /* t = (E + L)^{-1}t */
1065: diag = a->diag;
1066: for (i=0; i<m; i++) {
1067: n = diag[i] - a->i[i];
1068: idx = a->j + a->i[i];
1069: v = a->a + a->i[i];
1070: sum = t[i];
1071: SPARSEDENSEMDOT(sum,t,v,idx,n);
1072: t[i] = sum*idiag[i];
1074: /* x = x + t */
1075: x[i] += t[i];
1076: }
1078: PetscLogFlops(3*m-1 + 2*a->nz);
1079: VecRestoreArray(xx,&x);
1080: if (bb != xx) {VecRestoreArray(bb,&b);}
1081: return(0);
1082: } else if (flag & SOR_EISENSTAT) {
1083: /* Let A = L + U + D; where L is lower trianglar,
1084: U is upper triangular, E is diagonal; This routine applies
1086: (L + E)^{-1} A (U + E)^{-1}
1088: to a vector efficiently using Eisenstat's trick. This is for
1089: the case of SSOR preconditioner, so E is D/omega where omega
1090: is the relaxation factor.
1091: */
1092: scale = (2.0/omega) - 1.0;
1094: /* x = (E + U)^{-1} b */
1095: for (i=m-1; i>=0; i--) {
1096: d = fshift + a->a[diag[i] + shift];
1097: n = a->i[i+1] - diag[i] - 1;
1098: idx = a->j + diag[i] + (!shift);
1099: v = a->a + diag[i] + (!shift);
1100: sum = b[i];
1101: SPARSEDENSEMDOT(sum,xs,v,idx,n);
1102: x[i] = omega*(sum/d);
1103: }
1105: /* t = b - (2*E - D)x */
1106: v = a->a;
1107: for (i=0; i<m; i++) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; }
1109: /* t = (E + L)^{-1}t */
1110: ts = t + shift; /* shifted by one for index start of a or a->j*/
1111: diag = a->diag;
1112: for (i=0; i<m; i++) {
1113: d = fshift + a->a[diag[i]+shift];
1114: n = diag[i] - a->i[i];
1115: idx = a->j + a->i[i] + shift;
1116: v = a->a + a->i[i] + shift;
1117: sum = t[i];
1118: SPARSEDENSEMDOT(sum,ts,v,idx,n);
1119: t[i] = omega*(sum/d);
1120: /* x = x + t */
1121: x[i] += t[i];
1122: }
1124: PetscLogFlops(6*m-1 + 2*a->nz);
1125: VecRestoreArray(xx,&x);
1126: if (bb != xx) {VecRestoreArray(bb,&b);}
1127: return(0);
1128: }
1129: if (flag & SOR_ZERO_INITIAL_GUESS) {
1130: if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
1131: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1132: fortranrelaxaijforwardzero_(&m,&omega,x,a->i,a->j,diag,a->a,b);
1133: #else
1134: for (i=0; i<m; i++) {
1135: d = fshift + a->a[diag[i]+shift];
1136: n = diag[i] - a->i[i];
1137: PetscLogFlops(2*n-1);
1138: idx = a->j + a->i[i] + shift;
1139: v = a->a + a->i[i] + shift;
1140: sum = b[i];
1141: SPARSEDENSEMDOT(sum,xs,v,idx,n);
1142: x[i] = omega*(sum/d);
1143: }
1144: #endif
1145: xb = x;
1146: } else xb = b;
1147: if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
1148: (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
1149: for (i=0; i<m; i++) {
1150: x[i] *= a->a[diag[i]+shift];
1151: }
1152: PetscLogFlops(m);
1153: }
1154: if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
1155: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1156: fortranrelaxaijbackwardzero_(&m,&omega,x,a->i,a->j,diag,a->a,xb);
1157: #else
1158: for (i=m-1; i>=0; i--) {
1159: d = fshift + a->a[diag[i] + shift];
1160: n = a->i[i+1] - diag[i] - 1;
1161: PetscLogFlops(2*n-1);
1162: idx = a->j + diag[i] + (!shift);
1163: v = a->a + diag[i] + (!shift);
1164: sum = xb[i];
1165: SPARSEDENSEMDOT(sum,xs,v,idx,n);
1166: x[i] = omega*(sum/d);
1167: }
1168: #endif
1169: }
1170: its--;
1171: }
1172: while (its--) {
1173: if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
1174: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1175: fortranrelaxaijforward_(&m,&omega,x,a->i,a->j,diag,a->a,b);
1176: #else
1177: for (i=0; i<m; i++) {
1178: d = fshift + a->a[diag[i]+shift];
1179: n = a->i[i+1] - a->i[i];
1180: PetscLogFlops(2*n-1);
1181: idx = a->j + a->i[i] + shift;
1182: v = a->a + a->i[i] + shift;
1183: sum = b[i];
1184: SPARSEDENSEMDOT(sum,xs,v,idx,n);
1185: x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d;
1186: }
1187: #endif
1188: }
1189: if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
1190: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1191: fortranrelaxaijbackward_(&m,&omega,x,a->i,a->j,diag,a->a,b);
1192: #else
1193: for (i=m-1; i>=0; i--) {
1194: d = fshift + a->a[diag[i] + shift];
1195: n = a->i[i+1] - a->i[i];
1196: PetscLogFlops(2*n-1);
1197: idx = a->j + a->i[i] + shift;
1198: v = a->a + a->i[i] + shift;
1199: sum = b[i];
1200: SPARSEDENSEMDOT(sum,xs,v,idx,n);
1201: x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d;
1202: }
1203: #endif
1204: }
1205: }
1206: VecRestoreArray(xx,&x);
1207: if (bb != xx) {VecRestoreArray(bb,&b);}
1208: return(0);
1209: }
1211: int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,MatInfo *info)
1212: {
1213: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1216: info->rows_global = (double)A->m;
1217: info->columns_global = (double)A->n;
1218: info->rows_local = (double)A->m;
1219: info->columns_local = (double)A->n;
1220: info->block_size = 1.0;
1221: info->nz_allocated = (double)a->maxnz;
1222: info->nz_used = (double)a->nz;
1223: info->nz_unneeded = (double)(a->maxnz - a->nz);
1224: info->assemblies = (double)A->num_ass;
1225: info->mallocs = (double)a->reallocs;
1226: info->memory = A->mem;
1227: if (A->factor) {
1228: info->fill_ratio_given = A->info.fill_ratio_given;
1229: info->fill_ratio_needed = A->info.fill_ratio_needed;
1230: info->factor_mallocs = A->info.factor_mallocs;
1231: } else {
1232: info->fill_ratio_given = 0;
1233: info->fill_ratio_needed = 0;
1234: info->factor_mallocs = 0;
1235: }
1236: return(0);
1237: }
1239: EXTERN int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,MatLUInfo*,Mat*);
1240: EXTERN int MatLUFactorNumeric_SeqAIJ(Mat,Mat*);
1241: EXTERN int MatLUFactor_SeqAIJ(Mat,IS,IS,MatLUInfo*);
1242: EXTERN int MatSolve_SeqAIJ(Mat,Vec,Vec);
1243: EXTERN int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
1244: EXTERN int MatSolveTranspose_SeqAIJ(Mat,Vec,Vec);
1245: EXTERN int MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec);
1247: int MatZeroRows_SeqAIJ(Mat A,IS is,PetscScalar *diag)
1248: {
1249: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1250: int i,ierr,N,*rows,m = A->m - 1,shift = a->indexshift;
1253: ISGetLocalSize(is,&N);
1254: ISGetIndices(is,&rows);
1255: if (a->keepzeroedrows) {
1256: for (i=0; i<N; i++) {
1257: if (rows[i] < 0 || rows[i] > m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"row out of range");
1258: PetscMemzero(&a->a[a->i[rows[i]]+shift],a->ilen[rows[i]]*sizeof(PetscScalar));
1259: }
1260: if (diag) {
1261: MatMissingDiagonal_SeqAIJ(A);
1262: MatMarkDiagonal_SeqAIJ(A);
1263: for (i=0; i<N; i++) {
1264: a->a[a->diag[rows[i]]] = *diag;
1265: }
1266: }
1267: } else {
1268: if (diag) {
1269: for (i=0; i<N; i++) {
1270: if (rows[i] < 0 || rows[i] > m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"row out of range");
1271: if (a->ilen[rows[i]] > 0) {
1272: a->ilen[rows[i]] = 1;
1273: a->a[a->i[rows[i]]+shift] = *diag;
1274: a->j[a->i[rows[i]]+shift] = rows[i]+shift;
1275: } else { /* in case row was completely empty */
1276: MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES);
1277: }
1278: }
1279: } else {
1280: for (i=0; i<N; i++) {
1281: if (rows[i] < 0 || rows[i] > m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"row out of range");
1282: a->ilen[rows[i]] = 0;
1283: }
1284: }
1285: }
1286: ISRestoreIndices(is,&rows);
1287: MatAssemblyEnd_SeqAIJ(A,MAT_FINAL_ASSEMBLY);
1288: return(0);
1289: }
1291: int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,PetscScalar **v)
1292: {
1293: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1294: int *itmp,i,shift = a->indexshift,ierr;
1297: if (row < 0 || row >= A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Row %d out of range",row);
1299: *nz = a->i[row+1] - a->i[row];
1300: if (v) *v = a->a + a->i[row] + shift;
1301: if (idx) {
1302: itmp = a->j + a->i[row] + shift;
1303: if (*nz && shift) {
1304: PetscMalloc((*nz)*sizeof(int),idx);
1305: for (i=0; i<(*nz); i++) {(*idx)[i] = itmp[i] + shift;}
1306: } else if (*nz) {
1307: *idx = itmp;
1308: }
1309: else *idx = 0;
1310: }
1311: return(0);
1312: }
1314: int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,PetscScalar **v)
1315: {
1316: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1320: if (idx) {if (*idx && a->indexshift) {PetscFree(*idx);}}
1321: return(0);
1322: }
1324: int MatNorm_SeqAIJ(Mat A,NormType type,PetscReal *nrm)
1325: {
1326: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1327: PetscScalar *v = a->a;
1328: PetscReal sum = 0.0;
1329: int i,j,shift = a->indexshift,ierr;
1332: if (type == NORM_FROBENIUS) {
1333: for (i=0; i<a->nz; i++) {
1334: #if defined(PETSC_USE_COMPLEX)
1335: sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1336: #else
1337: sum += (*v)*(*v); v++;
1338: #endif
1339: }
1340: *nrm = sqrt(sum);
1341: } else if (type == NORM_1) {
1342: PetscReal *tmp;
1343: int *jj = a->j;
1344: PetscMalloc((A->n+1)*sizeof(PetscReal),&tmp);
1345: PetscMemzero(tmp,A->n*sizeof(PetscReal));
1346: *nrm = 0.0;
1347: for (j=0; j<a->nz; j++) {
1348: tmp[*jj++ + shift] += PetscAbsScalar(*v); v++;
1349: }
1350: for (j=0; j<A->n; j++) {
1351: if (tmp[j] > *nrm) *nrm = tmp[j];
1352: }
1353: PetscFree(tmp);
1354: } else if (type == NORM_INFINITY) {
1355: *nrm = 0.0;
1356: for (j=0; j<A->m; j++) {
1357: v = a->a + a->i[j] + shift;
1358: sum = 0.0;
1359: for (i=0; i<a->i[j+1]-a->i[j]; i++) {
1360: sum += PetscAbsScalar(*v); v++;
1361: }
1362: if (sum > *nrm) *nrm = sum;
1363: }
1364: } else {
1365: SETERRQ(PETSC_ERR_SUP,"No support for two norm");
1366: }
1367: return(0);
1368: }
1370: int MatTranspose_SeqAIJ(Mat A,Mat *B)
1371: {
1372: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1373: Mat C;
1374: int i,ierr,*aj = a->j,*ai = a->i,m = A->m,len,*col;
1375: int shift = a->indexshift;
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(int),&col);
1381: PetscMemzero(col,(1+A->n)*sizeof(int));
1382: if (shift) {
1383: for (i=0; i<ai[m]-1; i++) aj[i] -= 1;
1384: }
1385: for (i=0; i<ai[m]+shift; i++) col[aj[i]] += 1;
1386: MatCreateSeqAIJ(A->comm,A->n,m,0,col,&C);
1387: PetscFree(col);
1388: for (i=0; i<m; i++) {
1389: len = ai[i+1]-ai[i];
1390: ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES);
1391: array += len;
1392: aj += len;
1393: }
1394: if (shift) {
1395: for (i=0; i<ai[m]-1; i++) aj[i] += 1;
1396: }
1398: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1399: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
1401: if (B) {
1402: *B = C;
1403: } else {
1404: MatHeaderCopy(A,C);
1405: }
1406: return(0);
1407: }
1409: int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr)
1410: {
1411: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1412: PetscScalar *l,*r,x,*v;
1413: int ierr,i,j,m = A->m,n = A->n,M,nz = a->nz,*jj,shift = a->indexshift;
1416: if (ll) {
1417: /* The local size is used so that VecMPI can be passed to this routine
1418: by MatDiagonalScale_MPIAIJ */
1419: VecGetLocalSize(ll,&m);
1420: if (m != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length");
1421: VecGetArray(ll,&l);
1422: v = a->a;
1423: for (i=0; i<m; i++) {
1424: x = l[i];
1425: M = a->i[i+1] - a->i[i];
1426: for (j=0; j<M; j++) { (*v++) *= x;}
1427: }
1428: VecRestoreArray(ll,&l);
1429: PetscLogFlops(nz);
1430: }
1431: if (rr) {
1432: VecGetLocalSize(rr,&n);
1433: if (n != A->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vector wrong length");
1434: VecGetArray(rr,&r);
1435: v = a->a; jj = a->j;
1436: for (i=0; i<nz; i++) {
1437: (*v++) *= r[*jj++ + shift];
1438: }
1439: VecRestoreArray(rr,&r);
1440: PetscLogFlops(nz);
1441: }
1442: return(0);
1443: }
1445: int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,int csize,MatReuse scall,Mat *B)
1446: {
1447: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*c;
1448: int *smap,i,k,kstart,kend,ierr,oldcols = A->n,*lens;
1449: int row,mat_i,*mat_j,tcol,first,step,*mat_ilen,sum,lensi;
1450: int *irow,*icol,nrows,ncols,shift = a->indexshift,*ssmap;
1451: int *starts,*j_new,*i_new,*aj = a->j,*ai = a->i,ii,*ailen = a->ilen;
1452: PetscScalar *a_new,*mat_a;
1453: Mat C;
1454: PetscTruth stride;
1457: ISSorted(isrow,(PetscTruth*)&i);
1458: if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"ISrow is not sorted");
1459: ISSorted(iscol,(PetscTruth*)&i);
1460: if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IScol is not sorted");
1462: ISGetIndices(isrow,&irow);
1463: ISGetLocalSize(isrow,&nrows);
1464: ISGetLocalSize(iscol,&ncols);
1466: ISStrideGetInfo(iscol,&first,&step);
1467: ISStride(iscol,&stride);
1468: if (stride && step == 1) {
1469: /* special case of contiguous rows */
1470: ierr = PetscMalloc((2*nrows+1)*sizeof(int),&lens);
1471: starts = lens + nrows;
1472: /* loop over new rows determining lens and starting points */
1473: for (i=0; i<nrows; i++) {
1474: kstart = ai[irow[i]]+shift;
1475: kend = kstart + ailen[irow[i]];
1476: for (k=kstart; k<kend; k++) {
1477: if (aj[k]+shift >= first) {
1478: starts[i] = k;
1479: break;
1480: }
1481: }
1482: sum = 0;
1483: while (k < kend) {
1484: if (aj[k++]+shift >= first+ncols) break;
1485: sum++;
1486: }
1487: lens[i] = sum;
1488: }
1489: /* create submatrix */
1490: if (scall == MAT_REUSE_MATRIX) {
1491: int n_cols,n_rows;
1492: MatGetSize(*B,&n_rows,&n_cols);
1493: if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size");
1494: MatZeroEntries(*B);
1495: C = *B;
1496: } else {
1497: MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);
1498: }
1499: c = (Mat_SeqAIJ*)C->data;
1501: /* loop over rows inserting into submatrix */
1502: a_new = c->a;
1503: j_new = c->j;
1504: i_new = c->i;
1505: i_new[0] = -shift;
1506: for (i=0; i<nrows; i++) {
1507: ii = starts[i];
1508: lensi = lens[i];
1509: for (k=0; k<lensi; k++) {
1510: *j_new++ = aj[ii+k] - first;
1511: }
1512: PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(PetscScalar));
1513: a_new += lensi;
1514: i_new[i+1] = i_new[i] + lensi;
1515: c->ilen[i] = lensi;
1516: }
1517: PetscFree(lens);
1518: } else {
1519: ierr = ISGetIndices(iscol,&icol);
1520: ierr = PetscMalloc((1+oldcols)*sizeof(int),&smap);
1521: ssmap = smap + shift;
1522: ierr = PetscMalloc((1+nrows)*sizeof(int),&lens);
1523: ierr = PetscMemzero(smap,oldcols*sizeof(int));
1524: for (i=0; i<ncols; i++) smap[icol[i]] = i+1;
1525: /* determine lens of each row */
1526: for (i=0; i<nrows; i++) {
1527: kstart = ai[irow[i]]+shift;
1528: kend = kstart + a->ilen[irow[i]];
1529: lens[i] = 0;
1530: for (k=kstart; k<kend; k++) {
1531: if (ssmap[aj[k]]) {
1532: lens[i]++;
1533: }
1534: }
1535: }
1536: /* Create and fill new matrix */
1537: if (scall == MAT_REUSE_MATRIX) {
1538: PetscTruth equal;
1540: c = (Mat_SeqAIJ *)((*B)->data);
1541: if ((*B)->m != nrows || (*B)->n != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size");
1542: PetscMemcmp(c->ilen,lens,(*B)->m*sizeof(int),&equal);
1543: if (!equal) {
1544: SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros");
1545: }
1546: PetscMemzero(c->ilen,(*B)->m*sizeof(int));
1547: C = *B;
1548: } else {
1549: MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);
1550: }
1551: c = (Mat_SeqAIJ *)(C->data);
1552: for (i=0; i<nrows; i++) {
1553: row = irow[i];
1554: kstart = ai[row]+shift;
1555: kend = kstart + a->ilen[row];
1556: mat_i = c->i[i]+shift;
1557: mat_j = c->j + mat_i;
1558: mat_a = c->a + mat_i;
1559: mat_ilen = c->ilen + i;
1560: for (k=kstart; k<kend; k++) {
1561: if ((tcol=ssmap[a->j[k]])) {
1562: *mat_j++ = tcol - (!shift);
1563: *mat_a++ = a->a[k];
1564: (*mat_ilen)++;
1566: }
1567: }
1568: }
1569: /* Free work space */
1570: ISRestoreIndices(iscol,&icol);
1571: PetscFree(smap);
1572: PetscFree(lens);
1573: }
1574: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1575: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
1577: ISRestoreIndices(isrow,&irow);
1578: *B = C;
1579: return(0);
1580: }
1582: /*
1583: */
1584: int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,MatILUInfo *info)
1585: {
1586: Mat_SeqAIJ *a = (Mat_SeqAIJ*)inA->data;
1587: int ierr;
1588: Mat outA;
1589: PetscTruth row_identity,col_identity;
1592: if (info && info->levels != 0) SETERRQ(PETSC_ERR_SUP,"Only levels=0 supported for in-place ilu");
1593: ISIdentity(row,&row_identity);
1594: ISIdentity(col,&col_identity);
1595: if (!row_identity || !col_identity) {
1596: SETERRQ(1,"Row and column permutations must be identity for in-place ILU");
1597: }
1599: outA = inA;
1600: inA->factor = FACTOR_LU;
1601: a->row = row;
1602: a->col = col;
1603: PetscObjectReference((PetscObject)row);
1604: PetscObjectReference((PetscObject)col);
1606: /* Create the inverse permutation so that it can be used in MatLUFactorNumeric() */
1607: if (a->icol) {ISDestroy(a->icol);} /* need to remove old one */
1608: ISInvertPermutation(col,PETSC_DECIDE,&a->icol);
1609: PetscLogObjectParent(inA,a->icol);
1611: if (!a->solve_work) { /* this matrix may have been factored before */
1612: PetscMalloc((inA->m+1)*sizeof(PetscScalar),&a->solve_work);
1613: }
1615: if (!a->diag) {
1616: MatMarkDiagonal_SeqAIJ(inA);
1617: }
1618: MatLUFactorNumeric_SeqAIJ(inA,&outA);
1619: return(0);
1620: }
1622: #include petscblaslapack.h
1623: int MatScale_SeqAIJ(PetscScalar *alpha,Mat inA)
1624: {
1625: Mat_SeqAIJ *a = (Mat_SeqAIJ*)inA->data;
1626: int one = 1;
1629: BLscal_(&a->nz,alpha,a->a,&one);
1630: PetscLogFlops(a->nz);
1631: return(0);
1632: }
1634: int MatGetSubMatrices_SeqAIJ(Mat A,int n,IS *irow,IS *icol,MatReuse scall,Mat **B)
1635: {
1636: int ierr,i;
1639: if (scall == MAT_INITIAL_MATRIX) {
1640: PetscMalloc((n+1)*sizeof(Mat),B);
1641: }
1643: for (i=0; i<n; i++) {
1644: MatGetSubMatrix_SeqAIJ(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);
1645: }
1646: return(0);
1647: }
1649: int MatGetBlockSize_SeqAIJ(Mat A,int *bs)
1650: {
1652: *bs = 1;
1653: return(0);
1654: }
1656: int MatIncreaseOverlap_SeqAIJ(Mat A,int is_max,IS *is,int ov)
1657: {
1658: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1659: int shift,row,i,j,k,l,m,n,*idx,ierr,*nidx,isz,val;
1660: int start,end,*ai,*aj;
1661: PetscBT table;
1664: shift = a->indexshift;
1665: m = A->m;
1666: ai = a->i;
1667: aj = a->j+shift;
1669: if (ov < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"illegal overlap value used");
1671: PetscMalloc((m+1)*sizeof(int),&nidx);
1672: PetscBTCreate(m,table);
1674: for (i=0; i<is_max; i++) {
1675: /* Initialize the two local arrays */
1676: isz = 0;
1677: PetscBTMemzero(m,table);
1678:
1679: /* Extract the indices, assume there can be duplicate entries */
1680: ISGetIndices(is[i],&idx);
1681: ISGetLocalSize(is[i],&n);
1682:
1683: /* Enter these into the temp arrays. I.e., mark table[row], enter row into new index */
1684: for (j=0; j<n ; ++j){
1685: if(!PetscBTLookupSet(table,idx[j])) { nidx[isz++] = idx[j];}
1686: }
1687: ISRestoreIndices(is[i],&idx);
1688: ISDestroy(is[i]);
1689:
1690: k = 0;
1691: for (j=0; j<ov; j++){ /* for each overlap */
1692: n = isz;
1693: for (; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */
1694: row = nidx[k];
1695: start = ai[row];
1696: end = ai[row+1];
1697: for (l = start; l<end ; l++){
1698: val = aj[l] + shift;
1699: if (!PetscBTLookupSet(table,val)) {nidx[isz++] = val;}
1700: }
1701: }
1702: }
1703: ISCreateGeneral(PETSC_COMM_SELF,isz,nidx,(is+i));
1704: }
1705: PetscBTDestroy(table);
1706: PetscFree(nidx);
1707: return(0);
1708: }
1710: /* -------------------------------------------------------------- */
1711: int MatPermute_SeqAIJ(Mat A,IS rowp,IS colp,Mat *B)
1712: {
1713: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1714: PetscScalar *vwork;
1715: int i,ierr,nz,m = A->m,n = A->n,*cwork;
1716: int *row,*col,*cnew,j,*lens;
1717: IS icolp,irowp;
1720: ISInvertPermutation(rowp,PETSC_DECIDE,&irowp);
1721: ISGetIndices(irowp,&row);
1722: ISInvertPermutation(colp,PETSC_DECIDE,&icolp);
1723: ISGetIndices(icolp,&col);
1724:
1725: /* determine lengths of permuted rows */
1726: PetscMalloc((m+1)*sizeof(int),&lens);
1727: for (i=0; i<m; i++) {
1728: lens[row[i]] = a->i[i+1] - a->i[i];
1729: }
1730: MatCreateSeqAIJ(A->comm,m,n,0,lens,B);
1731: PetscFree(lens);
1733: PetscMalloc(n*sizeof(int),&cnew);
1734: for (i=0; i<m; i++) {
1735: MatGetRow(A,i,&nz,&cwork,&vwork);
1736: for (j=0; j<nz; j++) { cnew[j] = col[cwork[j]];}
1737: MatSetValues(*B,1,&row[i],nz,cnew,vwork,INSERT_VALUES);
1738: MatRestoreRow(A,i,&nz,&cwork,&vwork);
1739: }
1740: PetscFree(cnew);
1741: MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);
1742: MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);
1743: ISRestoreIndices(irowp,&row);
1744: ISRestoreIndices(icolp,&col);
1745: ISDestroy(irowp);
1746: ISDestroy(icolp);
1747: return(0);
1748: }
1750: int MatPrintHelp_SeqAIJ(Mat A)
1751: {
1752: static PetscTruth called = PETSC_FALSE;
1753: MPI_Comm comm = A->comm;
1754: int ierr;
1757: if (called) {return(0);} else called = PETSC_TRUE;
1758: (*PetscHelpPrintf)(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):n");
1759: (*PetscHelpPrintf)(comm," -mat_lu_pivotthreshold <threshold>: Set pivoting thresholdn");
1760: (*PetscHelpPrintf)(comm," -mat_aij_oneindex: internal indices begin at 1 instead of the default 0.n");
1761: (*PetscHelpPrintf)(comm," -mat_aij_no_inode: Do not use inodesn");
1762: (*PetscHelpPrintf)(comm," -mat_aij_inode_limit <limit>: Set inode limit (max limit=5)n");
1763: #if defined(PETSC_HAVE_ESSL)
1764: (*PetscHelpPrintf)(comm," -mat_aij_essl: Use IBM sparse LU factorization and solve.n");
1765: #endif
1766: #if defined(PETSC_HAVE_LUSOL)
1767: (*PetscHelpPrintf)(comm," -mat_aij_lusol: Use the Stanford LUSOL sparse factorization and solve.n");
1768: #endif
1769: #if defined(PETSC_HAVE_MATLAB_ENGINE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
1770: (*PetscHelpPrintf)(comm," -mat_aij_matlab: Use Matlab engine sparse LU factorization and solve.n");
1771: #endif
1772: return(0);
1773: }
1774: EXTERN int MatEqual_SeqAIJ(Mat A,Mat B,PetscTruth* flg);
1775: EXTERN int MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring);
1776: EXTERN int MatILUDTFactor_SeqAIJ(Mat,MatILUInfo*,IS,IS,Mat*);
1777: int MatCopy_SeqAIJ(Mat A,Mat B,MatStructure str)
1778: {
1779: int ierr;
1780: PetscTruth flg;
1783: PetscTypeCompare((PetscObject)B,MATSEQAIJ,&flg);
1784: if (str == SAME_NONZERO_PATTERN && flg) {
1785: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1786: Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data;
1788: if (a->i[A->m]+a->indexshift != b->i[B->m]+a->indexshift) {
1789: SETERRQ(1,"Number of nonzeros in two matrices are different");
1790: }
1791: PetscMemcpy(b->a,a->a,(a->i[A->m]+a->indexshift)*sizeof(PetscScalar));
1792: } else {
1793: MatCopy_Basic(A,B,str);
1794: }
1795: return(0);
1796: }
1798: int MatSetUpPreallocation_SeqAIJ(Mat A)
1799: {
1800: int ierr;
1803: MatSeqAIJSetPreallocation(A,PETSC_DEFAULT,0);
1804: return(0);
1805: }
1807: int MatGetArray_SeqAIJ(Mat A,PetscScalar **array)
1808: {
1809: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1811: *array = a->a;
1812: return(0);
1813: }
1815: int MatRestoreArray_SeqAIJ(Mat A,PetscScalar **array)
1816: {
1818: return(0);
1819: }
1821: int MatFDColoringApply_SeqAIJ(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
1822: {
1823: int (*f)(void *,Vec,Vec,void*) = (int (*)(void *,Vec,Vec,void *))coloring->f;
1824: int k,ierr,N,start,end,l,row,col,srow,**vscaleforrow,m1,m2;
1825: PetscScalar dx,mone = -1.0,*y,*xx,*w3_array;
1826: PetscScalar *vscale_array;
1827: PetscReal epsilon = coloring->error_rel,umin = coloring->umin;
1828: Vec w1,w2,w3;
1829: void *fctx = coloring->fctx;
1830: PetscTruth flg;
1833: if (!coloring->w1) {
1834: VecDuplicate(x1,&coloring->w1);
1835: PetscLogObjectParent(coloring,coloring->w1);
1836: VecDuplicate(x1,&coloring->w2);
1837: PetscLogObjectParent(coloring,coloring->w2);
1838: VecDuplicate(x1,&coloring->w3);
1839: PetscLogObjectParent(coloring,coloring->w3);
1840: }
1841: w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
1843: MatSetUnfactored(J);
1844: PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);
1845: if (flg) {
1846: PetscLogInfo(coloring,"MatFDColoringApply_SeqAIJ: Not calling MatZeroEntries()n");
1847: } else {
1848: MatZeroEntries(J);
1849: }
1851: VecGetOwnershipRange(x1,&start,&end);
1852: VecGetSize(x1,&N);
1854: /*
1855: This is a horrible, horrible, hack. See DMMGComputeJacobian_Multigrid() it inproperly sets
1856: coloring->F for the coarser grids from the finest
1857: */
1858: if (coloring->F) {
1859: VecGetLocalSize(coloring->F,&m1);
1860: VecGetLocalSize(w1,&m2);
1861: if (m1 != m2) {
1862: coloring->F = 0;
1863: }
1864: }
1866: if (coloring->F) {
1867: w1 = coloring->F;
1868: coloring->F = 0;
1869: } else {
1870: (*f)(sctx,x1,w1,fctx);
1871: }
1873: /*
1874: Compute all the scale factors and share with other processors
1875: */
1876: VecGetArray(x1,&xx);xx = xx - start;
1877: VecGetArray(coloring->vscale,&vscale_array);vscale_array = vscale_array - start;
1878: for (k=0; k<coloring->ncolors; k++) {
1879: /*
1880: Loop over each column associated with color adding the
1881: perturbation to the vector w3.
1882: */
1883: for (l=0; l<coloring->ncolumns[k]; l++) {
1884: col = coloring->columns[k][l]; /* column of the matrix we are probing for */
1885: dx = xx[col];
1886: if (dx == 0.0) dx = 1.0;
1887: #if !defined(PETSC_USE_COMPLEX)
1888: if (dx < umin && dx >= 0.0) dx = umin;
1889: else if (dx < 0.0 && dx > -umin) dx = -umin;
1890: #else
1891: if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin;
1892: else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
1893: #endif
1894: dx *= epsilon;
1895: vscale_array[col] = 1.0/dx;
1896: }
1897: }
1898: vscale_array = vscale_array + start;VecRestoreArray(coloring->vscale,&vscale_array);
1899: VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);
1900: VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);
1902: /* VecView(coloring->vscale,PETSC_VIEWER_STDOUT_WORLD);
1903: VecView(x1,PETSC_VIEWER_STDOUT_WORLD);*/
1905: if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
1906: else vscaleforrow = coloring->columnsforrow;
1908: VecGetArray(coloring->vscale,&vscale_array);
1909: /*
1910: Loop over each color
1911: */
1912: for (k=0; k<coloring->ncolors; k++) {
1913: VecCopy(x1,w3);
1914: VecGetArray(w3,&w3_array);w3_array = w3_array - start;
1915: /*
1916: Loop over each column associated with color adding the
1917: perturbation to the vector w3.
1918: */
1919: for (l=0; l<coloring->ncolumns[k]; l++) {
1920: col = coloring->columns[k][l]; /* column of the matrix we are probing for */
1921: dx = xx[col];
1922: if (dx == 0.0) dx = 1.0;
1923: #if !defined(PETSC_USE_COMPLEX)
1924: if (dx < umin && dx >= 0.0) dx = umin;
1925: else if (dx < 0.0 && dx > -umin) dx = -umin;
1926: #else
1927: if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin;
1928: else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
1929: #endif
1930: dx *= epsilon;
1931: if (!PetscAbsScalar(dx)) SETERRQ(1,"Computed 0 differencing parameter");
1932: w3_array[col] += dx;
1933: }
1934: w3_array = w3_array + start; VecRestoreArray(w3,&w3_array);
1936: /*
1937: Evaluate function at x1 + dx (here dx is a vector of perturbations)
1938: */
1940: (*f)(sctx,w3,w2,fctx);
1941: VecAXPY(&mone,w1,w2);
1943: /*
1944: Loop over rows of vector, putting results into Jacobian matrix
1945: */
1946: VecGetArray(w2,&y);
1947: for (l=0; l<coloring->nrows[k]; l++) {
1948: row = coloring->rows[k][l];
1949: col = coloring->columnsforrow[k][l];
1950: y[row] *= vscale_array[vscaleforrow[k][l]];
1951: srow = row + start;
1952: ierr = MatSetValues_SeqAIJ(J,1,&srow,1,&col,y+row,INSERT_VALUES);
1953: }
1954: VecRestoreArray(w2,&y);
1955: }
1956: VecRestoreArray(coloring->vscale,&vscale_array);
1957: xx = xx + start; ierr = VecRestoreArray(x1,&xx);
1958: ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
1959: ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
1960: return(0);
1961: }
1963: #include petscblaslapack.h
1965: int MatAXPY_SeqAIJ(PetscScalar *a,Mat X,Mat Y,MatStructure str)
1966: {
1967: int ierr,one=1;
1968: Mat_SeqAIJ *x = (Mat_SeqAIJ *)X->data,*y = (Mat_SeqAIJ *)Y->data;
1971: if (str == SAME_NONZERO_PATTERN) {
1972: BLaxpy_(&x->nz,a,x->a,&one,y->a,&one);
1973: } else {
1974: MatAXPY_Basic(a,X,Y,str);
1975: }
1976: return(0);
1977: }
1980: /* -------------------------------------------------------------------*/
1981: static struct _MatOps MatOps_Values = {MatSetValues_SeqAIJ,
1982: MatGetRow_SeqAIJ,
1983: MatRestoreRow_SeqAIJ,
1984: MatMult_SeqAIJ,
1985: MatMultAdd_SeqAIJ,
1986: MatMultTranspose_SeqAIJ,
1987: MatMultTransposeAdd_SeqAIJ,
1988: MatSolve_SeqAIJ,
1989: MatSolveAdd_SeqAIJ,
1990: MatSolveTranspose_SeqAIJ,
1991: MatSolveTransposeAdd_SeqAIJ,
1992: MatLUFactor_SeqAIJ,
1993: 0,
1994: MatRelax_SeqAIJ,
1995: MatTranspose_SeqAIJ,
1996: MatGetInfo_SeqAIJ,
1997: MatEqual_SeqAIJ,
1998: MatGetDiagonal_SeqAIJ,
1999: MatDiagonalScale_SeqAIJ,
2000: MatNorm_SeqAIJ,
2001: 0,
2002: MatAssemblyEnd_SeqAIJ,
2003: MatCompress_SeqAIJ,
2004: MatSetOption_SeqAIJ,
2005: MatZeroEntries_SeqAIJ,
2006: MatZeroRows_SeqAIJ,
2007: MatLUFactorSymbolic_SeqAIJ,
2008: MatLUFactorNumeric_SeqAIJ,
2009: 0,
2010: 0,
2011: MatSetUpPreallocation_SeqAIJ,
2012: MatILUFactorSymbolic_SeqAIJ,
2013: 0,
2014: MatGetArray_SeqAIJ,
2015: MatRestoreArray_SeqAIJ,
2016: MatDuplicate_SeqAIJ,
2017: 0,
2018: 0,
2019: MatILUFactor_SeqAIJ,
2020: 0,
2021: MatAXPY_SeqAIJ,
2022: MatGetSubMatrices_SeqAIJ,
2023: MatIncreaseOverlap_SeqAIJ,
2024: MatGetValues_SeqAIJ,
2025: MatCopy_SeqAIJ,
2026: MatPrintHelp_SeqAIJ,
2027: MatScale_SeqAIJ,
2028: 0,
2029: 0,
2030: MatILUDTFactor_SeqAIJ,
2031: MatGetBlockSize_SeqAIJ,
2032: MatGetRowIJ_SeqAIJ,
2033: MatRestoreRowIJ_SeqAIJ,
2034: MatGetColumnIJ_SeqAIJ,
2035: MatRestoreColumnIJ_SeqAIJ,
2036: MatFDColoringCreate_SeqAIJ,
2037: 0,
2038: 0,
2039: MatPermute_SeqAIJ,
2040: 0,
2041: 0,
2042: MatDestroy_SeqAIJ,
2043: MatView_SeqAIJ,
2044: MatGetPetscMaps_Petsc,
2045: 0,
2046: 0,
2047: 0,
2048: 0,
2049: 0,
2050: 0,
2051: 0,
2052: 0,
2053: MatSetColoring_SeqAIJ,
2054: MatSetValuesAdic_SeqAIJ,
2055: MatSetValuesAdifor_SeqAIJ,
2056: MatFDColoringApply_SeqAIJ};
2058: EXTERN int MatUseSuperLU_SeqAIJ(Mat);
2059: EXTERN int MatUseEssl_SeqAIJ(Mat);
2060: EXTERN int MatUseLUSOL_SeqAIJ(Mat);
2061: EXTERN int MatUseMatlab_SeqAIJ(Mat);
2062: EXTERN int MatUseDXML_SeqAIJ(Mat);
2064: EXTERN_C_BEGIN
2066: int MatSeqAIJSetColumnIndices_SeqAIJ(Mat mat,int *indices)
2067: {
2068: Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
2069: int i,nz,n;
2073: nz = aij->maxnz;
2074: n = mat->n;
2075: for (i=0; i<nz; i++) {
2076: aij->j[i] = indices[i];
2077: }
2078: aij->nz = nz;
2079: for (i=0; i<n; i++) {
2080: aij->ilen[i] = aij->imax[i];
2081: }
2083: return(0);
2084: }
2085: EXTERN_C_END
2087: /*@
2088: MatSeqAIJSetColumnIndices - Set the column indices for all the rows
2089: in the matrix.
2091: Input Parameters:
2092: + mat - the SeqAIJ matrix
2093: - indices - the column indices
2095: Level: advanced
2097: Notes:
2098: This can be called if you have precomputed the nonzero structure of the
2099: matrix and want to provide it to the matrix object to improve the performance
2100: of the MatSetValues() operation.
2102: You MUST have set the correct numbers of nonzeros per row in the call to
2103: MatCreateSeqAIJ().
2105: MUST be called before any calls to MatSetValues();
2107: The indices should start with zero, not one.
2109: @*/
2110: int MatSeqAIJSetColumnIndices(Mat mat,int *indices)
2111: {
2112: int ierr,(*f)(Mat,int *);
2116: PetscObjectQueryFunction((PetscObject)mat,"MatSeqAIJSetColumnIndices_C",(void (**)(void))&f);
2117: if (f) {
2118: (*f)(mat,indices);
2119: } else {
2120: SETERRQ(1,"Wrong type of matrix to set column indices");
2121: }
2122: return(0);
2123: }
2125: /* ----------------------------------------------------------------------------------------*/
2127: EXTERN_C_BEGIN
2128: int MatStoreValues_SeqAIJ(Mat mat)
2129: {
2130: Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
2131: int nz = aij->i[mat->m]+aij->indexshift,ierr;
2134: if (aij->nonew != 1) {
2135: SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
2136: }
2138: /* allocate space for values if not already there */
2139: if (!aij->saved_values) {
2140: PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);
2141: }
2143: /* copy values over */
2144: PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));
2145: return(0);
2146: }
2147: EXTERN_C_END
2149: /*@
2150: MatStoreValues - Stashes a copy of the matrix values; this allows, for
2151: example, reuse of the linear part of a Jacobian, while recomputing the
2152: nonlinear portion.
2154: Collect on Mat
2156: Input Parameters:
2157: . mat - the matrix (currently on AIJ matrices support this option)
2159: Level: advanced
2161: Common Usage, with SNESSolve():
2162: $ Create Jacobian matrix
2163: $ Set linear terms into matrix
2164: $ Apply boundary conditions to matrix, at this time matrix must have
2165: $ final nonzero structure (i.e. setting the nonlinear terms and applying
2166: $ boundary conditions again will not change the nonzero structure
2167: $ MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS);
2168: $ MatStoreValues(mat);
2169: $ Call SNESSetJacobian() with matrix
2170: $ In your Jacobian routine
2171: $ MatRetrieveValues(mat);
2172: $ Set nonlinear terms in matrix
2173:
2174: Common Usage without SNESSolve(), i.e. when you handle nonlinear solve yourself:
2175: $ // build linear portion of Jacobian
2176: $ MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS);
2177: $ MatStoreValues(mat);
2178: $ loop over nonlinear iterations
2179: $ MatRetrieveValues(mat);
2180: $ // call MatSetValues(mat,...) to set nonliner portion of Jacobian
2181: $ // call MatAssemblyBegin/End() on matrix
2182: $ Solve linear system with Jacobian
2183: $ endloop
2185: Notes:
2186: Matrix must already be assemblied before calling this routine
2187: Must set the matrix option MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS); before
2188: calling this routine.
2190: .seealso: MatRetrieveValues()
2192: @*/
2193: int MatStoreValues(Mat mat)
2194: {
2195: int ierr,(*f)(Mat);
2199: if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2200: if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2202: PetscObjectQueryFunction((PetscObject)mat,"MatStoreValues_C",(void (**)(void))&f);
2203: if (f) {
2204: (*f)(mat);
2205: } else {
2206: SETERRQ(1,"Wrong type of matrix to store values");
2207: }
2208: return(0);
2209: }
2211: EXTERN_C_BEGIN
2212: int MatRetrieveValues_SeqAIJ(Mat mat)
2213: {
2214: Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
2215: int nz = aij->i[mat->m]+aij->indexshift,ierr;
2218: if (aij->nonew != 1) {
2219: SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
2220: }
2221: if (!aij->saved_values) {
2222: SETERRQ(1,"Must call MatStoreValues(A);first");
2223: }
2225: /* copy values over */
2226: PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));
2227: return(0);
2228: }
2229: EXTERN_C_END
2231: /*@
2232: MatRetrieveValues - Retrieves the copy of the matrix values; this allows, for
2233: example, reuse of the linear part of a Jacobian, while recomputing the
2234: nonlinear portion.
2236: Collect on Mat
2238: Input Parameters:
2239: . mat - the matrix (currently on AIJ matrices support this option)
2241: Level: advanced
2243: .seealso: MatStoreValues()
2245: @*/
2246: int MatRetrieveValues(Mat mat)
2247: {
2248: int ierr,(*f)(Mat);
2252: if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2253: if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2255: PetscObjectQueryFunction((PetscObject)mat,"MatRetrieveValues_C",(void (**)(void))&f);
2256: if (f) {
2257: (*f)(mat);
2258: } else {
2259: SETERRQ(1,"Wrong type of matrix to retrieve values");
2260: }
2261: return(0);
2262: }
2264: /*
2265: This allows SeqAIJ matrices to be passed to the matlab engine
2266: */
2267: #if defined(PETSC_HAVE_MATLAB_ENGINE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
2268: #include "engine.h" /* Matlab include file */
2269: #include "mex.h" /* Matlab include file */
2270: EXTERN_C_BEGIN
2271: int MatMatlabEnginePut_SeqAIJ(PetscObject obj,void *engine)
2272: {
2273: int ierr,i,*ai,*aj;
2274: Mat B = (Mat)obj;
2275: PetscScalar *array;
2276: mxArray *mat;
2277: Mat_SeqAIJ *aij = (Mat_SeqAIJ*)B->data;
2280: mat = mxCreateSparse(B->n,B->m,aij->nz,mxREAL);
2281: PetscMemcpy(mxGetPr(mat),aij->a,aij->nz*sizeof(PetscScalar));
2282: /* Matlab stores by column, not row so we pass in the transpose of the matrix */
2283: PetscMemcpy(mxGetIr(mat),aij->j,aij->nz*sizeof(int));
2284: PetscMemcpy(mxGetJc(mat),aij->i,(B->m+1)*sizeof(int));
2286: /* Matlab indices start at 0 for sparse (what a surprise) */
2287: if (aij->indexshift) {
2288: for (i=0; i<B->m+1; i++) {
2289: ai[i]--;
2290: }
2291: for (i=0; i<aij->nz; i++) {
2292: aj[i]--;
2293: }
2294: }
2295: PetscObjectName(obj);
2296: mxSetName(mat,obj->name);
2297: engPutArray((Engine *)engine,mat);
2298: return(0);
2299: }
2300: EXTERN_C_END
2302: EXTERN_C_BEGIN
2303: int MatMatlabEngineGet_SeqAIJ(PetscObject obj,void *engine)
2304: {
2305: int ierr,ii;
2306: Mat mat = (Mat)obj;
2307: Mat_SeqAIJ *aij = (Mat_SeqAIJ*)mat->data;
2308: mxArray *mmat;
2311: PetscFree(aij->a);
2312: aij->indexshift = 0;
2314: mmat = engGetArray((Engine *)engine,obj->name);
2316: aij->nz = (mxGetJc(mmat))[mat->m];
2317: ierr = PetscMalloc(aij->nz*(sizeof(int)+sizeof(PetscScalar))+(mat->m+1)*sizeof(int),&aij->a);
2318: aij->j = (int*)(aij->a + aij->nz);
2319: aij->i = aij->j + aij->nz;
2320: aij->singlemalloc = PETSC_TRUE;
2321: aij->freedata = PETSC_TRUE;
2323: PetscMemcpy(aij->a,mxGetPr(mmat),aij->nz*sizeof(PetscScalar));
2324: /* Matlab stores by column, not row so we pass in the transpose of the matrix */
2325: PetscMemcpy(aij->j,mxGetIr(mmat),aij->nz*sizeof(int));
2326: PetscMemcpy(aij->i,mxGetJc(mmat),(mat->m+1)*sizeof(int));
2328: for (ii=0; ii<mat->m; ii++) {
2329: aij->ilen[ii] = aij->imax[ii] = aij->i[ii+1] - aij->i[ii];
2330: }
2332: MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
2333: MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);
2335: return(0);
2336: }
2337: EXTERN_C_END
2338: #endif
2340: /* --------------------------------------------------------------------------------*/
2341: /*@C
2342: MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format
2343: (the default parallel PETSc format). For good matrix assembly performance
2344: the user should preallocate the matrix storage by setting the parameter nz
2345: (or the array nnz). By setting these parameters accurately, performance
2346: during matrix assembly can be increased by more than a factor of 50.
2348: Collective on MPI_Comm
2350: Input Parameters:
2351: + comm - MPI communicator, set to PETSC_COMM_SELF
2352: . m - number of rows
2353: . n - number of columns
2354: . nz - number of nonzeros per row (same for all rows)
2355: - nnz - array containing the number of nonzeros in the various rows
2356: (possibly different for each row) or PETSC_NULL
2358: Output Parameter:
2359: . A - the matrix
2361: Notes:
2362: The AIJ format (also called the Yale sparse matrix format or
2363: compressed row storage), is fully compatible with standard Fortran 77
2364: storage. That is, the stored row and column indices can begin at
2365: either one (as in Fortran) or zero. See the users' manual for details.
2367: Specify the preallocated storage with either nz or nnz (not both).
2368: Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
2369: allocation. For large problems you MUST preallocate memory or you
2370: will get TERRIBLE performance, see the users' manual chapter on matrices.
2372: By default, this format uses inodes (identical nodes) when possible, to
2373: improve numerical efficiency of matrix-vector products and solves. We
2374: search for consecutive rows with the same nonzero structure, thereby
2375: reusing matrix information to achieve increased efficiency.
2377: Options Database Keys:
2378: + -mat_aij_no_inode - Do not use inodes
2379: . -mat_aij_inode_limit <limit> - Sets inode limit (max limit=5)
2380: - -mat_aij_oneindex - Internally use indexing starting at 1
2381: rather than 0. Note that when calling MatSetValues(),
2382: the user still MUST index entries starting at 0!
2384: Level: intermediate
2386: .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues(), MatSeqAIJSetColumnIndices(), MatCreateSeqAIJWithArrays()
2388: @*/
2389: int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz,Mat *A)
2390: {
2394: MatCreate(comm,m,n,m,n,A);
2395: MatSetType(*A,MATSEQAIJ);
2396: MatSeqAIJSetPreallocation(*A,nz,nnz);
2397: return(0);
2398: }
2400: #define SKIP_ALLOCATION -4
2402: /*@C
2403: MatSeqAIJSetPreallocation - For good matrix assembly performance
2404: the user should preallocate the matrix storage by setting the parameter nz
2405: (or the array nnz). By setting these parameters accurately, performance
2406: during matrix assembly can be increased by more than a factor of 50.
2408: Collective on MPI_Comm
2410: Input Parameters:
2411: + comm - MPI communicator, set to PETSC_COMM_SELF
2412: . m - number of rows
2413: . n - number of columns
2414: . nz - number of nonzeros per row (same for all rows)
2415: - nnz - array containing the number of nonzeros in the various rows
2416: (possibly different for each row) or PETSC_NULL
2418: Output Parameter:
2419: . A - the matrix
2421: Notes:
2422: The AIJ format (also called the Yale sparse matrix format or
2423: compressed row storage), is fully compatible with standard Fortran 77
2424: storage. That is, the stored row and column indices can begin at
2425: either one (as in Fortran) or zero. See the users' manual for details.
2427: Specify the preallocated storage with either nz or nnz (not both).
2428: Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
2429: allocation. For large problems you MUST preallocate memory or you
2430: will get TERRIBLE performance, see the users' manual chapter on matrices.
2432: By default, this format uses inodes (identical nodes) when possible, to
2433: improve numerical efficiency of matrix-vector products and solves. We
2434: search for consecutive rows with the same nonzero structure, thereby
2435: reusing matrix information to achieve increased efficiency.
2437: Options Database Keys:
2438: + -mat_aij_no_inode - Do not use inodes
2439: . -mat_aij_inode_limit <limit> - Sets inode limit (max limit=5)
2440: - -mat_aij_oneindex - Internally use indexing starting at 1
2441: rather than 0. Note that when calling MatSetValues(),
2442: the user still MUST index entries starting at 0!
2444: Level: intermediate
2446: .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues(), MatSeqAIJSetColumnIndices(), MatCreateSeqAIJWithArrays()
2448: @*/
2449: int MatSeqAIJSetPreallocation(Mat B,int nz,int *nnz)
2450: {
2451: Mat_SeqAIJ *b;
2452: int i,len=0,ierr;
2453: PetscTruth flg2,skipallocation = PETSC_FALSE;
2456: PetscTypeCompare((PetscObject)B,MATSEQAIJ,&flg2);
2457: if (!flg2) return(0);
2458:
2459: if (nz == SKIP_ALLOCATION) {
2460: skipallocation = PETSC_TRUE;
2461: nz = 0;
2462: }
2464: if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
2465: if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz);
2466: if (nnz) {
2467: for (i=0; i<B->m; i++) {
2468: if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]);
2469: 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);
2470: }
2471: }
2473: B->preallocated = PETSC_TRUE;
2474: b = (Mat_SeqAIJ*)B->data;
2476: PetscMalloc((B->m+1)*sizeof(int),&b->imax);
2477: if (!nnz) {
2478: if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 10;
2479: else if (nz <= 0) nz = 1;
2480: for (i=0; i<B->m; i++) b->imax[i] = nz;
2481: nz = nz*B->m;
2482: } else {
2483: nz = 0;
2484: for (i=0; i<B->m; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
2485: }
2487: if (!skipallocation) {
2488: /* allocate the matrix space */
2489: len = nz*(sizeof(int) + sizeof(PetscScalar)) + (B->m+1)*sizeof(int);
2490: ierr = PetscMalloc(len,&b->a);
2491: b->j = (int*)(b->a + nz);
2492: ierr = PetscMemzero(b->j,nz*sizeof(int));
2493: b->i = b->j + nz;
2494: b->i[0] = -b->indexshift;
2495: for (i=1; i<B->m+1; i++) {
2496: b->i[i] = b->i[i-1] + b->imax[i-1];
2497: }
2498: b->singlemalloc = PETSC_TRUE;
2499: b->freedata = PETSC_TRUE;
2500: } else {
2501: b->freedata = PETSC_FALSE;
2502: }
2504: /* b->ilen will count nonzeros in each row so far. */
2505: PetscMalloc((B->m+1)*sizeof(int),&b->ilen);
2506: PetscLogObjectMemory(B,len+2*(B->m+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqAIJ));
2507: for (i=0; i<B->m; i++) { b->ilen[i] = 0;}
2509: b->nz = 0;
2510: b->maxnz = nz;
2511: B->info.nz_unneeded = (double)b->maxnz;
2512: return(0);
2513: }
2515: EXTERN int RegisterApplyPtAPRoutines_Private(Mat);
2517: EXTERN_C_BEGIN
2518: int MatCreate_SeqAIJ(Mat B)
2519: {
2520: Mat_SeqAIJ *b;
2521: int ierr,size;
2522: PetscTruth flg;
2525: MPI_Comm_size(B->comm,&size);
2526: if (size > 1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Comm must be of size 1");
2528: B->m = B->M = PetscMax(B->m,B->M);
2529: B->n = B->N = PetscMax(B->n,B->N);
2531: PetscNew(Mat_SeqAIJ,&b);
2532: B->data = (void*)b;
2533: PetscMemzero(b,sizeof(Mat_SeqAIJ));
2534: PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
2535: B->factor = 0;
2536: B->lupivotthreshold = 1.0;
2537: B->mapping = 0;
2538: PetscOptionsGetReal(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold,PETSC_NULL);
2539: PetscOptionsHasName(PETSC_NULL,"-pc_ilu_preserve_row_sums",&b->ilu_preserve_row_sums);
2540: b->row = 0;
2541: b->col = 0;
2542: b->icol = 0;
2543: b->indexshift = 0;
2544: b->reallocs = 0;
2545: PetscOptionsHasName(PETSC_NULL,"-mat_aij_oneindex",&flg);
2546: if (flg) b->indexshift = -1;
2547:
2548: PetscMapCreateMPI(B->comm,B->m,B->m,&B->rmap);
2549: PetscMapCreateMPI(B->comm,B->n,B->n,&B->cmap);
2551: b->sorted = PETSC_FALSE;
2552: b->ignorezeroentries = PETSC_FALSE;
2553: b->roworiented = PETSC_TRUE;
2554: b->nonew = 0;
2555: b->diag = 0;
2556: b->solve_work = 0;
2557: B->spptr = 0;
2558: b->inode.use = PETSC_TRUE;
2559: b->inode.node_count = 0;
2560: b->inode.size = 0;
2561: b->inode.limit = 5;
2562: b->inode.max_limit = 5;
2563: b->saved_values = 0;
2564: b->idiag = 0;
2565: b->ssor = 0;
2566: b->keepzeroedrows = PETSC_FALSE;
2568: PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);
2570: #if defined(PETSC_HAVE_SUPERLU)
2571: PetscOptionsHasName(PETSC_NULL,"-mat_aij_superlu",&flg);
2572: if (flg) { MatUseSuperLU_SeqAIJ(B); }
2573: #endif
2574: PetscOptionsHasName(PETSC_NULL,"-mat_aij_essl",&flg);
2575: if (flg) { MatUseEssl_SeqAIJ(B); }
2576: PetscOptionsHasName(PETSC_NULL,"-mat_aij_lusol",&flg);
2577: if (flg) { MatUseLUSOL_SeqAIJ(B); }
2578: PetscOptionsHasName(PETSC_NULL,"-mat_aij_matlab",&flg);
2579: if (flg) {MatUseMatlab_SeqAIJ(B);}
2580: PetscOptionsHasName(PETSC_NULL,"-mat_aij_dxml",&flg);
2581: if (flg) {
2582: if (!b->indexshift) SETERRQ(PETSC_ERR_LIB,"need -mat_aij_oneindex with -mat_aij_dxml");
2583: MatUseDXML_SeqAIJ(B);
2584: }
2585: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqAIJSetColumnIndices_C",
2586: "MatSeqAIJSetColumnIndices_SeqAIJ",
2587: MatSeqAIJSetColumnIndices_SeqAIJ);
2588: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
2589: "MatStoreValues_SeqAIJ",
2590: MatStoreValues_SeqAIJ);
2591: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
2592: "MatRetrieveValues_SeqAIJ",
2593: MatRetrieveValues_SeqAIJ);
2594: #if defined(PETSC_HAVE_MATLAB_ENGINE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
2595: PetscObjectComposeFunctionDynamic((PetscObject)B,"PetscMatlabEnginePut_C","MatMatlabEnginePut_SeqAIJ",MatMatlabEnginePut_SeqAIJ);
2596: PetscObjectComposeFunctionDynamic((PetscObject)B,"PetscMatlabEngineGet_C","MatMatlabEngineGet_SeqAIJ",MatMatlabEngineGet_SeqAIJ);
2597: #endif
2598: RegisterApplyPtAPRoutines_Private(B);
2599: return(0);
2600: }
2601: EXTERN_C_END
2603: int MatDuplicate_SeqAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
2604: {
2605: Mat C;
2606: Mat_SeqAIJ *c,*a = (Mat_SeqAIJ*)A->data;
2607: int i,len,m = A->m,shift = a->indexshift,ierr;
2610: *B = 0;
2611: MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);
2612: MatSetType(C,MATSEQAIJ);
2613: c = (Mat_SeqAIJ*)C->data;
2615: C->factor = A->factor;
2616: c->row = 0;
2617: c->col = 0;
2618: c->icol = 0;
2619: c->indexshift = shift;
2620: c->keepzeroedrows = a->keepzeroedrows;
2621: C->assembled = PETSC_TRUE;
2623: C->M = A->m;
2624: C->N = A->n;
2626: PetscMalloc((m+1)*sizeof(int),&c->imax);
2627: PetscMalloc((m+1)*sizeof(int),&c->ilen);
2628: for (i=0; i<m; i++) {
2629: c->imax[i] = a->imax[i];
2630: c->ilen[i] = a->ilen[i];
2631: }
2633: /* allocate the matrix space */
2634: c->singlemalloc = PETSC_TRUE;
2635: len = (m+1)*sizeof(int)+(a->i[m])*(sizeof(PetscScalar)+sizeof(int));
2636: ierr = PetscMalloc(len,&c->a);
2637: c->j = (int*)(c->a + a->i[m] + shift);
2638: c->i = c->j + a->i[m] + shift;
2639: PetscMemcpy(c->i,a->i,(m+1)*sizeof(int));
2640: if (m > 0) {
2641: PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int));
2642: if (cpvalues == MAT_COPY_VALUES) {
2643: PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(PetscScalar));
2644: } else {
2645: PetscMemzero(c->a,(a->i[m]+shift)*sizeof(PetscScalar));
2646: }
2647: }
2649: PetscLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqAIJ));
2650: c->sorted = a->sorted;
2651: c->roworiented = a->roworiented;
2652: c->nonew = a->nonew;
2653: c->ilu_preserve_row_sums = a->ilu_preserve_row_sums;
2654: c->saved_values = 0;
2655: c->idiag = 0;
2656: c->ssor = 0;
2657: c->ignorezeroentries = a->ignorezeroentries;
2658: c->freedata = PETSC_TRUE;
2660: if (a->diag) {
2661: PetscMalloc((m+1)*sizeof(int),&c->diag);
2662: PetscLogObjectMemory(C,(m+1)*sizeof(int));
2663: for (i=0; i<m; i++) {
2664: c->diag[i] = a->diag[i];
2665: }
2666: } else c->diag = 0;
2667: c->inode.use = a->inode.use;
2668: c->inode.limit = a->inode.limit;
2669: c->inode.max_limit = a->inode.max_limit;
2670: if (a->inode.size){
2671: ierr = PetscMalloc((m+1)*sizeof(int),&c->inode.size);
2672: c->inode.node_count = a->inode.node_count;
2673: ierr = PetscMemcpy(c->inode.size,a->inode.size,(m+1)*sizeof(int));
2674: } else {
2675: c->inode.size = 0;
2676: c->inode.node_count = 0;
2677: }
2678: c->nz = a->nz;
2679: c->maxnz = a->maxnz;
2680: c->solve_work = 0;
2681: C->spptr = 0; /* Dangerous -I'm throwing away a->spptr */
2682: C->preallocated = PETSC_TRUE;
2684: *B = C;
2685: PetscFListDuplicate(A->qlist,&C->qlist);
2686: return(0);
2687: }
2689: EXTERN_C_BEGIN
2690: int MatLoad_SeqAIJ(PetscViewer viewer,MatType type,Mat *A)
2691: {
2692: Mat_SeqAIJ *a;
2693: Mat B;
2694: int i,nz,ierr,fd,header[4],size,*rowlengths = 0,M,N,shift;
2695: MPI_Comm comm;
2696:
2698: PetscObjectGetComm((PetscObject)viewer,&comm);
2699: MPI_Comm_size(comm,&size);
2700: if (size > 1) SETERRQ(PETSC_ERR_ARG_SIZ,"view must have one processor");
2701: PetscViewerBinaryGetDescriptor(viewer,&fd);
2702: PetscBinaryRead(fd,header,4,PETSC_INT);
2703: if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object in file");
2704: M = header[1]; N = header[2]; nz = header[3];
2706: if (nz < 0) {
2707: SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format on disk,cannot load as SeqAIJ");
2708: }
2710: /* read in row lengths */
2711: PetscMalloc(M*sizeof(int),&rowlengths);
2712: PetscBinaryRead(fd,rowlengths,M,PETSC_INT);
2714: /* create our matrix */
2715: MatCreateSeqAIJ(comm,M,N,0,rowlengths,A);
2716: B = *A;
2717: a = (Mat_SeqAIJ*)B->data;
2718: shift = a->indexshift;
2720: /* read in column indices and adjust for Fortran indexing*/
2721: PetscBinaryRead(fd,a->j,nz,PETSC_INT);
2722: if (shift) {
2723: for (i=0; i<nz; i++) {
2724: a->j[i] += 1;
2725: }
2726: }
2728: /* read in nonzero values */
2729: PetscBinaryRead(fd,a->a,nz,PETSC_SCALAR);
2731: /* set matrix "i" values */
2732: a->i[0] = -shift;
2733: for (i=1; i<= M; i++) {
2734: a->i[i] = a->i[i-1] + rowlengths[i-1];
2735: a->ilen[i-1] = rowlengths[i-1];
2736: }
2737: PetscFree(rowlengths);
2739: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
2740: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
2741: return(0);
2742: }
2743: EXTERN_C_END
2745: int MatEqual_SeqAIJ(Mat A,Mat B,PetscTruth* flg)
2746: {
2747: Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data,*b = (Mat_SeqAIJ *)B->data;
2748: int ierr;
2749: PetscTruth flag;
2752: PetscTypeCompare((PetscObject)B,MATSEQAIJ,&flag);
2753: if (!flag) SETERRQ(PETSC_ERR_ARG_INCOMP,"Matrices must be same type");
2755: /* If the matrix dimensions are not equal,or no of nonzeros or shift */
2756: if ((A->m != B->m) || (A->n != B->n) ||(a->nz != b->nz)|| (a->indexshift != b->indexshift)) {
2757: *flg = PETSC_FALSE;
2758: return(0);
2759: }
2760:
2761: /* if the a->i are the same */
2762: PetscMemcmp(a->i,b->i,(A->m+1)*sizeof(int),flg);
2763: if (*flg == PETSC_FALSE) return(0);
2764:
2765: /* if a->j are the same */
2766: PetscMemcmp(a->j,b->j,(a->nz)*sizeof(int),flg);
2767: if (*flg == PETSC_FALSE) return(0);
2768:
2769: /* if a->a are the same */
2770: PetscMemcmp(a->a,b->a,(a->nz)*sizeof(PetscScalar),flg);
2772: return(0);
2773:
2774: }
2776: /*@C
2777: MatCreateSeqAIJWithArrays - Creates an sequential AIJ matrix using matrix elements (in CSR format)
2778: provided by the user.
2780: Coolective on MPI_Comm
2782: Input Parameters:
2783: + comm - must be an MPI communicator of size 1
2784: . m - number of rows
2785: . n - number of columns
2786: . i - row indices
2787: . j - column indices
2788: - a - matrix values
2790: Output Parameter:
2791: . mat - the matrix
2793: Level: intermediate
2795: Notes:
2796: The i, j, and a arrays are not copied by this routine, the user must free these arrays
2797: once the matrix is destroyed
2799: You cannot set new nonzero locations into this matrix, that will generate an error.
2801: The i and j indices can be either 0- or 1 based
2803: .seealso: MatCreate(), MatCreateMPIAIJ(), MatCreateSeqAIJ()
2805: @*/
2806: int MatCreateSeqAIJWithArrays(MPI_Comm comm,int m,int n,int* i,int*j,PetscScalar *a,Mat *mat)
2807: {
2808: int ierr,ii;
2809: Mat_SeqAIJ *aij;
2812: MatCreateSeqAIJ(comm,m,n,SKIP_ALLOCATION,0,mat);
2813: aij = (Mat_SeqAIJ*)(*mat)->data;
2815: if (i[0] == 1) {
2816: aij->indexshift = -1;
2817: } else if (i[0]) {
2818: SETERRQ(1,"i (row indices) do not start with 0 or 1");
2819: }
2820: aij->i = i;
2821: aij->j = j;
2822: aij->a = a;
2823: aij->singlemalloc = PETSC_FALSE;
2824: aij->nonew = -1; /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/
2825: aij->freedata = PETSC_FALSE;
2827: for (ii=0; ii<m; ii++) {
2828: aij->ilen[ii] = aij->imax[ii] = i[ii+1] - i[ii];
2829: #if defined(PETSC_USE_BOPT_g)
2830: 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]);
2831: #endif
2832: }
2833: #if defined(PETSC_USE_BOPT_g)
2834: for (ii=0; ii<aij->i[m]; ii++) {
2835: if (j[ii] < -aij->indexshift) SETERRQ2(1,"Negative column index at location = %d index = %d",ii,j[ii]);
2836: if (j[ii] > n - 1 -aij->indexshift) SETERRQ2(1,"Column index to large at location = %d index = %d",ii,j[ii]);
2837: }
2838: #endif
2840: /* changes indices to start at 0 */
2841: if (i[0]) {
2842: aij->indexshift = 0;
2843: for (ii=0; ii<m; ii++) {
2844: i[ii]--;
2845: }
2846: for (ii=0; ii<i[m]; ii++) {
2847: j[ii]--;
2848: }
2849: }
2851: MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
2852: MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
2853: return(0);
2854: }
2856: int MatSetColoring_SeqAIJ(Mat A,ISColoring coloring)
2857: {
2858: int ierr;
2859: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
2862: if (coloring->ctype == IS_COLORING_LOCAL) {
2863: ierr = ISColoringReference(coloring);
2864: a->coloring = coloring;
2865: } else if (coloring->ctype == IS_COLORING_GHOSTED) {
2866: int *colors,i,*larray;
2867: ISColoring ocoloring;
2869: /* set coloring for diagonal portion */
2870: PetscMalloc((A->n+1)*sizeof(int),&larray);
2871: for (i=0; i<A->n; i++) {
2872: larray[i] = i;
2873: }
2874: ISGlobalToLocalMappingApply(A->mapping,IS_GTOLM_MASK,A->n,larray,PETSC_NULL,larray);
2875: PetscMalloc((A->n+1)*sizeof(int),&colors);
2876: for (i=0; i<A->n; i++) {
2877: colors[i] = coloring->colors[larray[i]];
2878: }
2879: PetscFree(larray);
2880: ISColoringCreate(PETSC_COMM_SELF,A->n,colors,&ocoloring);
2881: a->coloring = ocoloring;
2882: }
2883: return(0);
2884: }
2886: #if defined(PETSC_HAVE_ADIC) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
2887: EXTERN_C_BEGIN
2888: #include "adic/ad_utils.h"
2889: EXTERN_C_END
2891: int MatSetValuesAdic_SeqAIJ(Mat A,void *advalues)
2892: {
2893: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
2894: int m = A->m,*ii = a->i,*jj = a->j,nz,i,*color,j,nlen;
2895: PetscScalar *v = a->a,*values;
2896: char *cadvalues = (char *)advalues;
2899: if (!a->coloring) SETERRQ(1,"Coloring not set for matrix");
2900: nlen = PetscADGetDerivTypeSize();
2901: color = a->coloring->colors;
2902: /* loop over rows */
2903: for (i=0; i<m; i++) {
2904: nz = ii[i+1] - ii[i];
2905: /* loop over columns putting computed value into matrix */
2906: values = PetscADGetGradArray(cadvalues);
2907: for (j=0; j<nz; j++) {
2908: *v++ = values[color[*jj++]];
2909: }
2910: cadvalues += nlen; /* jump to next row of derivatives */
2911: }
2912: return(0);
2913: }
2915: #else
2917: int MatSetValuesAdic_SeqAIJ(Mat A,void *advalues)
2918: {
2920: SETERRQ(1,"PETSc installed without ADIC");
2921: }
2923: #endif
2925: int MatSetValuesAdifor_SeqAIJ(Mat A,int nl,void *advalues)
2926: {
2927: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
2928: int m = A->m,*ii = a->i,*jj = a->j,nz,i,*color,j;
2929: PetscScalar *v = a->a,*values = (PetscScalar *)advalues;
2932: if (!a->coloring) SETERRQ(1,"Coloring not set for matrix");
2933: color = a->coloring->colors;
2934: /* loop over rows */
2935: for (i=0; i<m; i++) {
2936: nz = ii[i+1] - ii[i];
2937: /* loop over columns putting computed value into matrix */
2938: for (j=0; j<nz; j++) {
2939: *v++ = values[color[*jj++]];
2940: }
2941: values += nl; /* jump to next row of derivatives */
2942: }
2943: return(0);
2944: }