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