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