xref: /petsc/src/mat/impls/aij/seq/aij.c (revision 76dd722b9f5340f7d2cb207df73dc90d5e00a5a5)
1 #ifndef lint
2 static char vcid[] = "$Id: aij.c,v 1.111 1995/11/06 19:49:26 balay Exp balay $";
3 #endif
4 
5 /*
6     Defines the basic matrix operations for the AIJ (compressed row)
7   matrix storage format.
8 */
9 #include "aij.h"
10 #include "vec/vecimpl.h"
11 #include "inline/spops.h"
12 
13 extern int MatToSymmetricIJ_SeqAIJ(Mat_SeqAIJ*,int**,int**);
14 
15 static int MatGetReordering_SeqAIJ(Mat A,MatOrdering type,IS *rperm, IS *cperm)
16 {
17   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
18   int        ierr, *ia, *ja,n,*idx,i;
19 
20   if (!a->assembled) SETERRQ(1,"MatGetReordering_SeqAIJ:Not for unassembled matrix");
21 
22   /*
23      this is tacky: In the future when we have written special factorization
24      and solve routines for the identity permutation we should use a
25      stride index set instead of the general one.
26   */
27   if (type  == ORDER_NATURAL) {
28     n = a->n;
29     idx = (int *) PetscMalloc( n*sizeof(int) ); CHKPTRQ(idx);
30     for ( i=0; i<n; i++ ) idx[i] = i;
31     ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,rperm); CHKERRQ(ierr);
32     ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,cperm); CHKERRQ(ierr);
33     PetscFree(idx);
34     ISSetPermutation(*rperm);
35     ISSetPermutation(*cperm);
36     ISSetIdentity(*rperm);
37     ISSetIdentity(*cperm);
38     return 0;
39   }
40 
41   ierr = MatToSymmetricIJ_SeqAIJ( a, &ia, &ja ); CHKERRQ(ierr);
42   ierr = MatGetReordering_IJ(a->n,ia,ja,type,rperm,cperm); CHKERRQ(ierr);
43   PetscFree(ia); PetscFree(ja);
44   return 0;
45 }
46 
47 #define CHUNKSIZE   10
48 
49 /* This version has row oriented v  */
50 static int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
51 {
52   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
53   int        *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax, N, sorted = a->sorted;
54   int        *imax = a->imax, *ai = a->i, *ailen = a->ilen;
55   int        *aj = a->j, nonew = a->nonew,shift = a->indexshift;
56   Scalar     *ap,value, *aa = a->a;
57 
58   for ( k=0; k<m; k++ ) { /* loop over added rows */
59     row  = im[k];
60     if (row < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative row");
61     if (row >= a->m) SETERRQ(1,"MatSetValues_SeqAIJ:Row too large");
62     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
63     rmax = imax[row]; nrow = ailen[row];
64     low = 0;
65     for ( l=0; l<n; l++ ) { /* loop over added columns */
66       if (in[l] < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative column");
67       if (in[l] >= a->n) SETERRQ(1,"MatSetValues_SeqAIJ:Column too large");
68       col = in[l] - shift; value = *v++;
69       if (!sorted) low = 0; high = nrow;
70       while (high-low > 5) {
71         t = (low+high)/2;
72         if (rp[t] > col) high = t;
73         else             low  = t;
74       }
75       for ( i=low; i<high; i++ ) {
76         if (rp[i] > col) break;
77         if (rp[i] == col) {
78           if (is == ADD_VALUES) ap[i] += value;
79           else                  ap[i] = value;
80           goto noinsert;
81         }
82       }
83       if (nonew) goto noinsert;
84       if (nrow >= rmax) {
85         /* there is no extra room in row, therefore enlarge */
86         int    new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j;
87         Scalar *new_a;
88 
89         /* malloc new storage space */
90         len     = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int);
91         new_a   = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a);
92         new_j   = (int *) (new_a + new_nz);
93         new_i   = new_j + new_nz;
94 
95         /* copy over old data into new slots */
96         for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];}
97         for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
98         PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int));
99         len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift);
100         PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow,
101                                                            len*sizeof(int));
102         PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar));
103         PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow,
104                                                            len*sizeof(Scalar));
105         /* free up old matrix storage */
106         PetscFree(a->a);
107         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);}
108         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
109         a->singlemalloc = 1;
110 
111         rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
112         rmax = imax[row] = imax[row] + CHUNKSIZE;
113         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar)));
114         a->maxnz += CHUNKSIZE;
115       }
116       N = nrow++ - 1; a->nz++;
117       /* shift up all the later entries in this row */
118       for ( ii=N; ii>=i; ii-- ) {
119         rp[ii+1] = rp[ii];
120         ap[ii+1] = ap[ii];
121       }
122       rp[i] = col;
123       ap[i] = value;
124       noinsert:;
125       low = i + 1;
126     }
127     ailen[row] = nrow;
128   }
129   return 0;
130 }
131 
132 #include "draw.h"
133 #include "pinclude/pviewer.h"
134 #include "sysio.h"
135 
136 static int MatView_SeqAIJ_Binary(Mat A,Viewer viewer)
137 {
138   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
139   int        i, fd, *col_lens, ierr;
140 
141   ierr = ViewerFileGetDescriptor_Private(viewer,&fd); CHKERRQ(ierr);
142   col_lens = (int *) PetscMalloc( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens);
143   col_lens[0] = MAT_COOKIE;
144   col_lens[1] = a->m;
145   col_lens[2] = a->n;
146   col_lens[3] = a->nz;
147 
148   /* store lengths of each row and write (including header) to file */
149   for ( i=0; i<a->m; i++ ) {
150     col_lens[4+i] = a->i[i+1] - a->i[i];
151   }
152   ierr = SYWrite(fd,col_lens,4+a->m,SYINT,1); CHKERRQ(ierr);
153   PetscFree(col_lens);
154 
155   /* store column indices (zero start index) */
156   if (a->indexshift) {
157     for ( i=0; i<a->nz; i++ ) a->j[i]--;
158   }
159   ierr = SYWrite(fd,a->j,a->nz,SYINT,0); CHKERRQ(ierr);
160   if (a->indexshift) {
161     for ( i=0; i<a->nz; i++ ) a->j[i]++;
162   }
163 
164   /* store nonzero values */
165   ierr = SYWrite(fd,a->a,a->nz,SYSCALAR,0); CHKERRQ(ierr);
166   return 0;
167 }
168 
169 static int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer)
170 {
171   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
172   int         ierr, i,j, m = a->m, shift = a->indexshift,format;
173   FILE        *fd;
174   char        *outputname;
175 
176   ierr = ViewerFileGetPointer_Private(viewer,&fd); CHKERRQ(ierr);
177   ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr);
178   ierr = ViewerFileGetFormat_Private(viewer,&format);
179   if (format == FILE_FORMAT_INFO) {
180     /* no need to print additional information */ ;
181   }
182   else if (format == FILE_FORMAT_MATLAB) {
183     int nz, nzalloc, mem;
184     MatGetInfo(A,MAT_LOCAL,&nz,&nzalloc,&mem);
185     fprintf(fd,"%% Size = %d %d \n",m,a->n);
186     fprintf(fd,"%% Nonzeros = %d \n",nz);
187     fprintf(fd,"zzz = zeros(%d,3);\n",nz);
188     fprintf(fd,"zzz = [\n");
189 
190     for (i=0; i<m; i++) {
191       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
192 #if defined(PETSC_COMPLEX)
193         fprintf(fd,"%d %d  %18.16e  %18.16e \n",i+1,a->j[j],real(a->a[j]),
194                    imag(a->a[j]));
195 #else
196         fprintf(fd,"%d %d  %18.16e\n", i+1, a->j[j], a->a[j]);
197 #endif
198       }
199     }
200     fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname);
201   }
202   else {
203     for ( i=0; i<m; i++ ) {
204       fprintf(fd,"row %d:",i);
205       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
206 #if defined(PETSC_COMPLEX)
207         if (imag(a->a[j]) != 0.0) {
208           fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j]));
209         }
210         else {
211           fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j]));
212         }
213 #else
214         fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]);
215 #endif
216       }
217       fprintf(fd,"\n");
218     }
219   }
220   fflush(fd);
221   return 0;
222 }
223 
224 static int MatView_SeqAIJ_Draw(Mat A,Viewer viewer)
225 {
226   Mat_SeqAIJ  *a = (Mat_SeqAIJ *) A->data;
227   int         ierr, i,j, m = a->m, shift = a->indexshift,pause,color;
228   double      xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r;
229   DrawCtx     draw = (DrawCtx) viewer;
230   DrawButton  button;
231 
232   xr  = a->n; yr = a->m; h = yr/10.0; w = xr/10.0;
233   xr += w;    yr += h;  xl = -w;     yl = -h;
234   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
235   /* loop over matrix elements drawing boxes */
236   color = DRAW_BLUE;
237   for ( i=0; i<m; i++ ) {
238     y_l = m - i - 1.0; y_r = y_l + 1.0;
239     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
240       x_l = a->j[j] + shift; x_r = x_l + 1.0;
241 #if defined(PETSC_COMPLEX)
242       if (real(a->a[j]) >=  0.) continue;
243 #else
244       if (a->a[j] >=  0.) continue;
245 #endif
246       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
247     }
248   }
249   color = DRAW_CYAN;
250   for ( i=0; i<m; i++ ) {
251     y_l = m - i - 1.0; y_r = y_l + 1.0;
252     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
253       x_l = a->j[j] + shift; x_r = x_l + 1.0;
254       if (a->a[j] !=  0.) continue;
255       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
256     }
257   }
258   color = DRAW_RED;
259   for ( i=0; i<m; i++ ) {
260     y_l = m - i - 1.0; y_r = y_l + 1.0;
261     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
262       x_l = a->j[j] + shift; x_r = x_l + 1.0;
263 #if defined(PETSC_COMPLEX)
264       if (real(a->a[j]) <=  0.) continue;
265 #else
266       if (a->a[j] <=  0.) continue;
267 #endif
268       DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
269     }
270   }
271   DrawFlush(draw);
272   DrawGetPause(draw,&pause);
273   if (pause >= 0) { PetscSleep(pause); return 0;}
274 
275   /* allow the matrix to zoom or shrink */
276   ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
277   while (button != BUTTON_RIGHT) {
278     DrawClear(draw);
279     if (button == BUTTON_LEFT) scale = .5;
280     else if (button == BUTTON_CENTER) scale = 2.;
281     xl = scale*(xl + w - xc) + xc - w*scale;
282     xr = scale*(xr - w - xc) + xc + w*scale;
283     yl = scale*(yl + h - yc) + yc - h*scale;
284     yr = scale*(yr - h - yc) + yc + h*scale;
285     w *= scale; h *= scale;
286     ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
287     color = DRAW_BLUE;
288     for ( i=0; i<m; i++ ) {
289       y_l = m - i - 1.0; y_r = y_l + 1.0;
290       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
291         x_l = a->j[j] + shift; x_r = x_l + 1.0;
292 #if defined(PETSC_COMPLEX)
293         if (real(a->a[j]) >=  0.) continue;
294 #else
295         if (a->a[j] >=  0.) continue;
296 #endif
297         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
298       }
299     }
300     color = DRAW_CYAN;
301     for ( i=0; i<m; i++ ) {
302       y_l = m - i - 1.0; y_r = y_l + 1.0;
303       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
304         x_l = a->j[j] + shift; x_r = x_l + 1.0;
305         if (a->a[j] !=  0.) continue;
306         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
307       }
308     }
309     color = DRAW_RED;
310     for ( i=0; i<m; i++ ) {
311       y_l = m - i - 1.0; y_r = y_l + 1.0;
312       for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
313         x_l = a->j[j] + shift; x_r = x_l + 1.0;
314 #if defined(PETSC_COMPLEX)
315         if (real(a->a[j]) <=  0.) continue;
316 #else
317         if (a->a[j] <=  0.) continue;
318 #endif
319         DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
320       }
321     }
322     ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
323   }
324   return 0;
325 }
326 
327 static int MatView_SeqAIJ(PetscObject obj,Viewer viewer)
328 {
329   Mat         A = (Mat) obj;
330   Mat_SeqAIJ  *a = (Mat_SeqAIJ*) A->data;
331   PetscObject vobj = (PetscObject) viewer;
332 
333   if (!a->assembled) SETERRQ(1,"MatView_SeqAIJ:Not for unassembled matrix");
334   if (!viewer) {
335     viewer = STDOUT_VIEWER_SELF; vobj = (PetscObject) viewer;
336   }
337   if (vobj->cookie == VIEWER_COOKIE) {
338     if (vobj->type == MATLAB_VIEWER) {
339       return ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j);
340     }
341     else if (vobj->type == ASCII_FILE_VIEWER || vobj->type == ASCII_FILES_VIEWER){
342       return MatView_SeqAIJ_ASCII(A,viewer);
343     }
344     else if (vobj->type == BINARY_FILE_VIEWER) {
345       return MatView_SeqAIJ_Binary(A,viewer);
346     }
347   }
348   else if (vobj->cookie == DRAW_COOKIE) {
349     if (vobj->type == NULLWINDOW) return 0;
350     else return MatView_SeqAIJ_Draw(A,viewer);
351   }
352   return 0;
353 }
354 
355 static int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode)
356 {
357   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
358   int        fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax;
359   int        m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift;
360   Scalar     *aa = a->a, *ap;
361 
362   if (mode == FLUSH_ASSEMBLY) return 0;
363 
364   for ( i=1; i<m; i++ ) {
365     /* move each row back by the amount of empty slots (fshift) before it*/
366     fshift += imax[i-1] - ailen[i-1];
367     if (fshift) {
368       ip = aj + ai[i] + shift; ap = aa + ai[i] + shift;
369       N = ailen[i];
370       for ( j=0; j<N; j++ ) {
371         ip[j-fshift] = ip[j];
372         ap[j-fshift] = ap[j];
373       }
374     }
375     ai[i] = ai[i-1] + ailen[i-1];
376   }
377   if (m) {
378     fshift += imax[m-1] - ailen[m-1];
379     ai[m] = ai[m-1] + ailen[m-1];
380   }
381   /* reset ilen and imax for each row */
382   for ( i=0; i<m; i++ ) {
383     ailen[i] = imax[i] = ai[i+1] - ai[i];
384   }
385   a->nz = ai[m] + shift;
386 
387   /* diagonals may have moved, so kill the diagonal pointers */
388   if (fshift && a->diag) {
389     PetscFree(a->diag);
390     PLogObjectMemory(A,-(m+1)*sizeof(int));
391     a->diag = 0;
392   }
393   /* check out for identical nodes. If found,use inode functions*/
394   ierr = Mat_SeqAIJ_CheckInode(A); CHKERRQ(ierr);
395   a->assembled = 1;
396   return 0;
397 }
398 
399 static int MatZeroEntries_SeqAIJ(Mat A)
400 {
401   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
402   PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar));
403   return 0;
404 }
405 
406 int MatDestroy_SeqAIJ(PetscObject obj)
407 {
408   Mat        A  = (Mat) obj;
409   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
410 
411 #if defined(PETSC_LOG)
412   PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
413 #endif
414   PetscFree(a->a);
415   if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);}
416   if (a->diag) PetscFree(a->diag);
417   if (a->ilen) PetscFree(a->ilen);
418   if (a->imax) PetscFree(a->imax);
419   if (a->solve_work) PetscFree(a->solve_work);
420   if (a->inode.size) PetscFree(a->inode.size);
421   PetscFree(a);
422   PLogObjectDestroy(A);
423   PetscHeaderDestroy(A);
424   return 0;
425 }
426 
427 static int MatCompress_SeqAIJ(Mat A)
428 {
429   return 0;
430 }
431 
432 static int MatSetOption_SeqAIJ(Mat A,MatOption op)
433 {
434   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
435   if      (op == ROW_ORIENTED)              a->roworiented = 1;
436   else if (op == COLUMNS_SORTED)            a->sorted      = 1;
437   else if (op == NO_NEW_NONZERO_LOCATIONS)  a->nonew       = 1;
438   else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew       = 0;
439   else if (op == ROWS_SORTED ||
440            op == SYMMETRIC_MATRIX ||
441            op == STRUCTURALLY_SYMMETRIC_MATRIX ||
442            op == YES_NEW_DIAGONALS)
443     PLogInfo((PetscObject)A,"Info:MatSetOption_SeqAIJ:Option ignored\n");
444   else if (op == COLUMN_ORIENTED)
445     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:COLUMN_ORIENTED");}
446   else if (op == NO_NEW_DIAGONALS)
447     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:NO_NEW_DIAGONALS");}
448   else
449     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:unknown option");}
450   return 0;
451 }
452 
453 static int MatGetDiagonal_SeqAIJ(Mat A,Vec v)
454 {
455   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
456   int        i,j, n,shift = a->indexshift;
457   Scalar     *x, zero = 0.0;
458 
459   if (!a->assembled) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Not for unassembled matrix");
460   VecSet(&zero,v);
461   VecGetArray(v,&x); VecGetLocalSize(v,&n);
462   if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Nonconforming matrix and vector");
463   for ( i=0; i<a->m; i++ ) {
464     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
465       if (a->j[j]+shift == i) {
466         x[i] = a->a[j];
467         break;
468       }
469     }
470   }
471   return 0;
472 }
473 
474 /* -------------------------------------------------------*/
475 /* Should check that shapes of vectors and matrices match */
476 /* -------------------------------------------------------*/
477 static int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy)
478 {
479   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
480   Scalar     *x, *y, *v, alpha;
481   int        m = a->m, n, i, *idx, shift = a->indexshift;
482 
483   if (!a->assembled) SETERRQ(1,"MatMultTrans_SeqAIJ:Not for unassembled matrix");
484   VecGetArray(xx,&x); VecGetArray(yy,&y);
485   PetscMemzero(y,a->n*sizeof(Scalar));
486   y = y + shift; /* shift for Fortran start by 1 indexing */
487   for ( i=0; i<m; i++ ) {
488     idx   = a->j + a->i[i] + shift;
489     v     = a->a + a->i[i] + shift;
490     n     = a->i[i+1] - a->i[i];
491     alpha = x[i];
492     while (n-->0) {y[*idx++] += alpha * *v++;}
493   }
494   PLogFlops(2*a->nz - a->n);
495   return 0;
496 }
497 
498 static int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
499 {
500   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
501   Scalar     *x, *y, *v, alpha;
502   int        m = a->m, n, i, *idx,shift = a->indexshift;
503 
504   if (!a->assembled) SETERRQ(1,"MatMultTransAdd_SeqAIJ:Not for unassembled matrix");
505   VecGetArray(xx,&x); VecGetArray(yy,&y);
506   if (zz != yy) VecCopy(zz,yy);
507   y = y + shift; /* shift for Fortran start by 1 indexing */
508   for ( i=0; i<m; i++ ) {
509     idx   = a->j + a->i[i] + shift;
510     v     = a->a + a->i[i] + shift;
511     n     = a->i[i+1] - a->i[i];
512     alpha = x[i];
513     while (n-->0) {y[*idx++] += alpha * *v++;}
514   }
515   return 0;
516 }
517 
518 static int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
519 {
520   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
521   Scalar     *x, *y, *v, sum;
522   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii;
523 
524   if (!a->assembled) SETERRQ(1,"MatMult_SeqAIJ:Not for unassembled matrix");
525   VecGetArray(xx,&x); VecGetArray(yy,&y);
526   x    = x + shift; /* shift for Fortran start by 1 indexing */
527   idx  = a->j;
528   v    = a->a;
529   ii   = a->i;
530   for ( i=0; i<m; i++ ) {
531     n    = ii[1] - ii[0]; ii++;
532     sum  = 0.0;
533     /* SPARSEDENSEDOT(sum,x,v,idx,n);  */
534     /* for ( j=n-1; j>-1; j--) sum += v[j]*x[idx[j]]; */
535     while (n--) sum += *v++ * x[*idx++];
536     y[i] = sum;
537   }
538   PLogFlops(2*a->nz - m);
539   return 0;
540 }
541 
542 static int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
543 {
544   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
545   Scalar     *x, *y, *z, *v, sum;
546   int        m = a->m, n, i, *idx, shift = a->indexshift,*ii;
547 
548   if (!a->assembled) SETERRQ(1,"MatMultAdd_SeqAIJ:Not for unassembled matrix");
549   VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z);
550   x    = x + shift; /* shift for Fortran start by 1 indexing */
551   idx  = a->j;
552   v    = a->a;
553   ii   = a->i;
554   for ( i=0; i<m; i++ ) {
555     n    = ii[1] - ii[0]; ii++;
556     sum  = y[i];
557     /* SPARSEDENSEDOT(sum,x,v,idx,n);  */
558     while (n--) sum += *v++ * x[*idx++];
559     z[i] = sum;
560   }
561   PLogFlops(2*a->nz);
562   return 0;
563 }
564 
565 /*
566      Adds diagonal pointers to sparse matrix structure.
567 */
568 
569 int MatMarkDiag_SeqAIJ(Mat A)
570 {
571   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
572   int        i,j, *diag, m = a->m,shift = a->indexshift;
573 
574   if (!a->assembled) SETERRQ(1,"MatMarkDiag_SeqAIJ:unassembled matrix");
575   diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag);
576   PLogObjectMemory(A,(m+1)*sizeof(int));
577   for ( i=0; i<a->m; i++ ) {
578     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
579       if (a->j[j]+shift == i) {
580         diag[i] = j - shift;
581         break;
582       }
583     }
584   }
585   a->diag = diag;
586   return 0;
587 }
588 
589 static int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag,
590                            double fshift,int its,Vec xx)
591 {
592   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
593   Scalar     *x, *b, *bs,  d, *xs, sum, *v = a->a,*t,scale,*ts, *xb;
594   int        ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift;
595 
596   VecGetArray(xx,&x); VecGetArray(bb,&b);
597   if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;}
598   diag = a->diag;
599   xs   = x + shift; /* shifted by one for index start of a or a->j*/
600   if (flag == SOR_APPLY_UPPER) {
601    /* apply ( U + D/omega) to the vector */
602     bs = b + shift;
603     for ( i=0; i<m; i++ ) {
604         d    = fshift + a->a[diag[i] + shift];
605         n    = a->i[i+1] - diag[i] - 1;
606         idx  = a->j + diag[i] + (!shift);
607         v    = a->a + diag[i] + (!shift);
608         sum  = b[i]*d/omega;
609         SPARSEDENSEDOT(sum,bs,v,idx,n);
610         x[i] = sum;
611     }
612     return 0;
613   }
614   if (flag == SOR_APPLY_LOWER) {
615     SETERRQ(1,"MatRelax_SeqAIJ:SOR_APPLY_LOWER is not done");
616   }
617   else if (flag & SOR_EISENSTAT) {
618     /* Let  A = L + U + D; where L is lower trianglar,
619     U is upper triangular, E is diagonal; This routine applies
620 
621             (L + E)^{-1} A (U + E)^{-1}
622 
623     to a vector efficiently using Eisenstat's trick. This is for
624     the case of SSOR preconditioner, so E is D/omega where omega
625     is the relaxation factor.
626     */
627     t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t);
628     scale = (2.0/omega) - 1.0;
629 
630     /*  x = (E + U)^{-1} b */
631     for ( i=m-1; i>=0; i-- ) {
632       d    = fshift + a->a[diag[i] + shift];
633       n    = a->i[i+1] - diag[i] - 1;
634       idx  = a->j + diag[i] + (!shift);
635       v    = a->a + diag[i] + (!shift);
636       sum  = b[i];
637       SPARSEDENSEMDOT(sum,xs,v,idx,n);
638       x[i] = omega*(sum/d);
639     }
640 
641     /*  t = b - (2*E - D)x */
642     v = a->a;
643     for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; }
644 
645     /*  t = (E + L)^{-1}t */
646     ts = t + shift; /* shifted by one for index start of a or a->j*/
647     diag = a->diag;
648     for ( i=0; i<m; i++ ) {
649       d    = fshift + a->a[diag[i]+shift];
650       n    = diag[i] - a->i[i];
651       idx  = a->j + a->i[i] + shift;
652       v    = a->a + a->i[i] + shift;
653       sum  = t[i];
654       SPARSEDENSEMDOT(sum,ts,v,idx,n);
655       t[i] = omega*(sum/d);
656     }
657 
658     /*  x = x + t */
659     for ( i=0; i<m; i++ ) { x[i] += t[i]; }
660     PetscFree(t);
661     return 0;
662   }
663   if (flag & SOR_ZERO_INITIAL_GUESS) {
664     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
665       for ( i=0; i<m; i++ ) {
666         d    = fshift + a->a[diag[i]+shift];
667         n    = diag[i] - a->i[i];
668         idx  = a->j + a->i[i] + shift;
669         v    = a->a + a->i[i] + shift;
670         sum  = b[i];
671         SPARSEDENSEMDOT(sum,xs,v,idx,n);
672         x[i] = omega*(sum/d);
673       }
674       xb = x;
675     }
676     else xb = b;
677     if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
678         (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
679       for ( i=0; i<m; i++ ) {
680         x[i] *= a->a[diag[i]+shift];
681       }
682     }
683     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
684       for ( i=m-1; i>=0; i-- ) {
685         d    = fshift + a->a[diag[i] + shift];
686         n    = a->i[i+1] - diag[i] - 1;
687         idx  = a->j + diag[i] + (!shift);
688         v    = a->a + diag[i] + (!shift);
689         sum  = xb[i];
690         SPARSEDENSEMDOT(sum,xs,v,idx,n);
691         x[i] = omega*(sum/d);
692       }
693     }
694     its--;
695   }
696   while (its--) {
697     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
698       for ( i=0; i<m; i++ ) {
699         d    = fshift + a->a[diag[i]+shift];
700         n    = a->i[i+1] - a->i[i];
701         idx  = a->j + a->i[i] + shift;
702         v    = a->a + a->i[i] + shift;
703         sum  = b[i];
704         SPARSEDENSEMDOT(sum,xs,v,idx,n);
705         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
706       }
707     }
708     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
709       for ( i=m-1; i>=0; i-- ) {
710         d    = fshift + a->a[diag[i] + shift];
711         n    = a->i[i+1] - a->i[i];
712         idx  = a->j + a->i[i] + shift;
713         v    = a->a + a->i[i] + shift;
714         sum  = b[i];
715         SPARSEDENSEMDOT(sum,xs,v,idx,n);
716         x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]);
717       }
718     }
719   }
720   return 0;
721 }
722 
723 static int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,int *nz,int *nzalloc,int *mem)
724 {
725   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
726   *nz      = a->nz;
727   *nzalloc = a->maxnz;
728   *mem     = (int)A->mem;
729   return 0;
730 }
731 
732 extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*);
733 extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*);
734 extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double);
735 extern int MatSolve_SeqAIJ(Mat,Vec,Vec);
736 extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
737 extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec);
738 extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec);
739 
740 static int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag)
741 {
742   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
743   int         i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift;
744 
745   ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr);
746   ierr = ISGetIndices(is,&rows); CHKERRQ(ierr);
747   if (diag) {
748     for ( i=0; i<N; i++ ) {
749       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
750       if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */
751         a->ilen[rows[i]] = 1;
752         a->a[a->i[rows[i]]+shift] = *diag;
753         a->j[a->i[rows[i]]+shift] = rows[i]+shift;
754       }
755       else {
756         ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES);
757         CHKERRQ(ierr);
758       }
759     }
760   }
761   else {
762     for ( i=0; i<N; i++ ) {
763       if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range");
764       a->ilen[rows[i]] = 0;
765     }
766   }
767   ISRestoreIndices(is,&rows);
768   ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
769   ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
770   return 0;
771 }
772 
773 static int MatGetSize_SeqAIJ(Mat A,int *m,int *n)
774 {
775   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
776   *m = a->m; *n = a->n;
777   return 0;
778 }
779 
780 static int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n)
781 {
782   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
783   *m = 0; *n = a->m;
784   return 0;
785 }
786 static int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
787 {
788   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
789   int        *itmp,i,ierr,shift = a->indexshift;
790 
791   if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqAIJ:Row out of range");
792 
793   if (!a->assembled) {
794     ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
795     ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr);
796   }
797   *nz = a->i[row+1] - a->i[row];
798   if (v) *v = a->a + a->i[row] + shift;
799   if (idx) {
800     if (*nz) {
801       itmp = a->j + a->i[row] + shift;
802       *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
803       for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;}
804     }
805     else *idx = 0;
806   }
807   return 0;
808 }
809 
810 static int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
811 {
812   if (idx) {if (*idx) PetscFree(*idx);}
813   return 0;
814 }
815 
816 static int MatNorm_SeqAIJ(Mat A,NormType type,double *norm)
817 {
818   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
819   Scalar     *v = a->a;
820   double     sum = 0.0;
821   int        i, j,shift = a->indexshift;
822 
823   if (!a->assembled) SETERRQ(1,"MatNorm_SeqAIJ:Not for unassembled matrix");
824   if (type == NORM_FROBENIUS) {
825     for (i=0; i<a->nz; i++ ) {
826 #if defined(PETSC_COMPLEX)
827       sum += real(conj(*v)*(*v)); v++;
828 #else
829       sum += (*v)*(*v); v++;
830 #endif
831     }
832     *norm = sqrt(sum);
833   }
834   else if (type == NORM_1) {
835     double *tmp;
836     int    *jj = a->j;
837     tmp = (double *) PetscMalloc( a->n*sizeof(double) ); CHKPTRQ(tmp);
838     PetscMemzero(tmp,a->n*sizeof(double));
839     *norm = 0.0;
840     for ( j=0; j<a->nz; j++ ) {
841         tmp[*jj++ + shift] += PetscAbsScalar(*v);  v++;
842     }
843     for ( j=0; j<a->n; j++ ) {
844       if (tmp[j] > *norm) *norm = tmp[j];
845     }
846     PetscFree(tmp);
847   }
848   else if (type == NORM_INFINITY) {
849     *norm = 0.0;
850     for ( j=0; j<a->m; j++ ) {
851       v = a->a + a->i[j] + shift;
852       sum = 0.0;
853       for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) {
854         sum += PetscAbsScalar(*v); v++;
855       }
856       if (sum > *norm) *norm = sum;
857     }
858   }
859   else {
860     SETERRQ(1,"MatNorm_SeqAIJ:No support for two norm yet");
861   }
862   return 0;
863 }
864 
865 static int MatTranspose_SeqAIJ(Mat A,Mat *B)
866 {
867   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
868   Mat        C;
869   int        i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col;
870   int        shift = a->indexshift;
871   Scalar     *array = a->a;
872 
873   if (!B && m != a->n) SETERRQ(1,"MatTranspose_SeqAIJ:Not for rectangular mat in place");
874   col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col);
875   PetscMemzero(col,(1+a->n)*sizeof(int));
876   if (shift) {
877     for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1;
878   }
879   for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1;
880   ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr);
881   PetscFree(col);
882   for ( i=0; i<m; i++ ) {
883     len = ai[i+1]-ai[i];
884     ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr);
885     array += len; aj += len;
886   }
887   if (shift) {
888     for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1;
889   }
890 
891   ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
892   ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
893 
894   if (B) {
895     *B = C;
896   } else {
897     /* This isn't really an in-place transpose */
898     PetscFree(a->a);
899     if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);}
900     if (a->diag) PetscFree(a->diag);
901     if (a->ilen) PetscFree(a->ilen);
902     if (a->imax) PetscFree(a->imax);
903     if (a->solve_work) PetscFree(a->solve_work);
904     PetscFree(a);
905     PetscMemcpy(A,C,sizeof(struct _Mat));
906     PetscHeaderDestroy(C);
907   }
908   return 0;
909 }
910 
911 static int MatScale_SeqAIJ(Mat A,Vec ll,Vec rr)
912 {
913   Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
914   Scalar     *l,*r,x,*v;
915   int        i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift;
916 
917   if (!a->assembled) SETERRQ(1,"MatScale_SeqAIJ:Not for unassembled matrix");
918   if (ll) {
919     VecGetArray(ll,&l); VecGetSize(ll,&m);
920     if (m != a->m) SETERRQ(1,"MatScale_SeqAIJ:Left scaling vector wrong length");
921     v = a->a;
922     for ( i=0; i<m; i++ ) {
923       x = l[i];
924       M = a->i[i+1] - a->i[i];
925       for ( j=0; j<M; j++ ) { (*v++) *= x;}
926     }
927   }
928   if (rr) {
929     VecGetArray(rr,&r); VecGetSize(rr,&n);
930     if (n != a->n) SETERRQ(1,"MatScale_SeqAIJ:Right scaling vector wrong length");
931     v = a->a; jj = a->j;
932     for ( i=0; i<nz; i++ ) {
933       (*v++) *= r[*jj++ + shift];
934     }
935   }
936   return 0;
937 }
938 
939 static int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B)
940 {
941   Mat_SeqAIJ   *a = (Mat_SeqAIJ *) A->data,*c;
942   int          nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens;
943   register int sum,lensi;
944   int          *irow, *icol, nrows, ncols, *cwork, shift = a->indexshift,*ssmap;
945   int          first,step,*starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen;
946   Scalar       *vwork,*a_new;
947   Mat          C;
948 
949   if (!a->assembled) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:Not for unassembled matrix");
950   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
951   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
952   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
953 
954   if (ISStrideGetInfo(iscol,&first,&step) && step == 1) {
955     /* special case of contiguous rows */
956     lens   = (int *) PetscMalloc((2*ncols+1)*sizeof(int)); CHKPTRQ(lens);
957     starts = lens + ncols;
958     /* loop over new rows determining lens and starting points */
959     for (i=0; i<nrows; i++) {
960       kstart  = ai[irow[i]]+shift;
961       kend    = kstart + ailen[irow[i]];
962       for ( k=kstart; k<kend; k++ ) {
963         if (aj[k] >= first) {
964           starts[i] = k;
965           break;
966 	}
967       }
968       sum = 0;
969       while (k < kend) {
970         if (aj[k++] >= first+ncols) break;
971         sum++;
972       }
973       lens[i] = sum;
974     }
975     /* create submatrix */
976     if (scall == MAT_REUSE_MATRIX) {
977       int n_cols,n_rows;
978       ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr);
979       if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:");
980       MatZeroEntries(*B);
981       C = *B;
982     }
983     else {
984       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
985     }
986     c = (Mat_SeqAIJ*) C->data;
987 
988     /* loop over rows inserting into submatrix */
989     a_new    = c->a;
990     j_new    = c->j;
991     i_new    = c->i;
992     i_new[0] = -shift;
993     for (i=0; i<nrows; i++) {
994       ii    = starts[i];
995       lensi = lens[i];
996       for ( k=0; k<lensi; k++ ) {
997         *j_new++ = aj[ii+k] - first;
998       }
999       PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar));
1000       a_new      += lensi;
1001       i_new[i+1]  = i_new[i] + lensi;
1002       c->ilen[i]  = lensi;
1003     }
1004     PetscFree(lens);
1005   }
1006   else {
1007     ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
1008     smap  = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap);
1009     ssmap = smap + shift;
1010     cwork = (int *) PetscMalloc((1+2*ncols)*sizeof(int)); CHKPTRQ(cwork);
1011     lens  = cwork + ncols;
1012     vwork = (Scalar *) PetscMalloc((1+ncols)*sizeof(Scalar)); CHKPTRQ(vwork);
1013     PetscMemzero(smap,oldcols*sizeof(int));
1014     for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
1015     /* determine lens of each row */
1016     for (i=0; i<nrows; i++) {
1017       kstart  = a->i[irow[i]]+shift;
1018       kend    = kstart + a->ilen[irow[i]];
1019       lens[i] = 0;
1020       for ( k=kstart; k<kend; k++ ) {
1021         if (ssmap[a->j[k]]) {
1022           lens[i]++;
1023         }
1024       }
1025     }
1026     /* Create and fill new matrix */
1027     if (scall == MAT_REUSE_MATRIX) {
1028       int n_cols,n_rows;
1029       ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr);
1030       if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:");
1031       MatZeroEntries(*B);
1032       C = *B;
1033     }
1034     else {
1035       ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr);
1036     }
1037     for (i=0; i<nrows; i++) {
1038       nznew  = 0;
1039       kstart = a->i[irow[i]]+shift;
1040       kend   = kstart + a->ilen[irow[i]];
1041       for ( k=kstart; k<kend; k++ ) {
1042         if (ssmap[a->j[k]]) {
1043           cwork[nznew]   = ssmap[a->j[k]] - 1;
1044           vwork[nznew++] = a->a[k];
1045         }
1046       }
1047       ierr = MatSetValues(C,1,&i,nznew,cwork,vwork,INSERT_VALUES); CHKERRQ(ierr);
1048     }
1049     /* Free work space */
1050     ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
1051     PetscFree(smap); PetscFree(cwork); PetscFree(vwork);
1052   }
1053   ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
1054   ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
1055 
1056   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
1057   *B = C;
1058   return 0;
1059 }
1060 
1061 /*
1062      note: This can only work for identity for row and col. It would
1063    be good to check this and otherwise generate an error.
1064 */
1065 static int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill)
1066 {
1067   Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data;
1068   int        ierr;
1069   Mat        outA;
1070 
1071   if (fill != 0) SETERRQ(1,"MatILUFactor_SeqAIJ:Only fill=0 supported");
1072 
1073   outA          = inA;
1074   inA->factor   = FACTOR_LU;
1075   a->row        = row;
1076   a->col        = col;
1077 
1078   a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work);
1079 
1080   if (!a->diag) {
1081     ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr);
1082   }
1083   ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr);
1084   return 0;
1085 }
1086 
1087 static int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,
1088                                     Mat **B)
1089 {
1090   int ierr,i;
1091 
1092   if (scall == MAT_INITIAL_MATRIX) {
1093     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
1094   }
1095 
1096   for ( i=0; i<n; i++ ) {
1097     ierr = MatGetSubMatrix(A,irow[i],icol[i],scall,&(*B)[i]); CHKERRQ(ierr);
1098   }
1099   return 0;
1100 }
1101 
1102 /* -------------------------------------------------------------------*/
1103 
1104 static struct _MatOps MatOps = {MatSetValues_SeqAIJ,
1105        MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ,
1106        MatMult_SeqAIJ,MatMultAdd_SeqAIJ,
1107        MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ,
1108        MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ,
1109        MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ,
1110        MatLUFactor_SeqAIJ,0,
1111        MatRelax_SeqAIJ,
1112        MatTranspose_SeqAIJ,
1113        MatGetInfo_SeqAIJ,0,
1114        MatGetDiagonal_SeqAIJ,MatScale_SeqAIJ,MatNorm_SeqAIJ,
1115        0,MatAssemblyEnd_SeqAIJ,
1116        MatCompress_SeqAIJ,
1117        MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ,
1118        MatGetReordering_SeqAIJ,
1119        MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0,
1120        MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ,
1121        MatILUFactorSymbolic_SeqAIJ,0,
1122        0,0,MatConvert_SeqAIJ,
1123        MatGetSubMatrix_SeqAIJ,0,
1124        MatCopyPrivate_SeqAIJ,0,0,
1125        MatILUFactor_SeqAIJ,0,0,
1126        MatGetSubMatrices_SeqAIJ};
1127 
1128 extern int MatUseSuperLU_SeqAIJ(Mat);
1129 extern int MatUseEssl_SeqAIJ(Mat);
1130 extern int MatUseDXML_SeqAIJ(Mat);
1131 
1132 /*@C
1133    MatCreateSeqAIJ - Creates a sparse matrix in AIJ format
1134    (the default uniprocessor PETSc format).
1135 
1136    Input Parameters:
1137 .  comm - MPI communicator, set to MPI_COMM_SELF
1138 .  m - number of rows
1139 .  n - number of columns
1140 .  nz - number of nonzeros per row (same for all rows)
1141 .  nzz - number of nonzeros per row or null (possibly different for each row)
1142 
1143    Output Parameter:
1144 .  A - the matrix
1145 
1146    Notes:
1147    The AIJ format (also called the Yale sparse matrix format or
1148    compressed row storage), is fully compatible with standard Fortran 77
1149    storage.  That is, the stored row and column indices can begin at
1150    either one (as in Fortran) or zero.  See the users manual for details.
1151 
1152    Specify the preallocated storage with either nz or nnz (not both).
1153    Set both nz and nnz to zero for PETSc to control dynamic memory
1154    allocation.
1155 
1156 .keywords: matrix, aij, compressed row, sparse
1157 
1158 .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues()
1159 @*/
1160 int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A)
1161 {
1162   Mat        B;
1163   Mat_SeqAIJ *b;
1164   int        i,len,ierr;
1165 
1166   *A      = 0;
1167   PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQAIJ,comm);
1168   PLogObjectCreate(B);
1169   B->data               = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b);
1170   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
1171   B->destroy          = MatDestroy_SeqAIJ;
1172   B->view             = MatView_SeqAIJ;
1173   B->factor           = 0;
1174   B->lupivotthreshold = 1.0;
1175   OptionsGetDouble(0,"-mat_lu_pivotthreshold",&B->lupivotthreshold);
1176   b->row              = 0;
1177   b->col              = 0;
1178   b->indexshift       = 0;
1179   if (OptionsHasName(0,"-mat_aij_oneindex")) b->indexshift = -1;
1180 
1181   b->m       = m;
1182   b->n       = n;
1183   b->imax    = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax);
1184   if (!nnz) {
1185     if (nz <= 0) nz = 1;
1186     for ( i=0; i<m; i++ ) b->imax[i] = nz;
1187     nz = nz*m;
1188   }
1189   else {
1190     nz = 0;
1191     for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
1192   }
1193 
1194   /* allocate the matrix space */
1195   len     = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int);
1196   b->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a);
1197   b->j  = (int *) (b->a + nz);
1198   PetscMemzero(b->j,nz*sizeof(int));
1199   b->i  = b->j + nz;
1200   b->singlemalloc = 1;
1201 
1202   b->i[0] = -b->indexshift;
1203   for (i=1; i<m+1; i++) {
1204     b->i[i] = b->i[i-1] + b->imax[i-1];
1205   }
1206 
1207   /* b->ilen will count nonzeros in each row so far. */
1208   b->ilen = (int *) PetscMalloc((m+1)*sizeof(int));
1209   PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1210   for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;}
1211 
1212   b->nz               = 0;
1213   b->maxnz            = nz;
1214   b->sorted           = 0;
1215   b->roworiented      = 1;
1216   b->nonew            = 0;
1217   b->diag             = 0;
1218   b->assembled        = 0;
1219   b->solve_work       = 0;
1220   b->spptr            = 0;
1221   b->inode.node_count = 0;
1222   b->inode.size       = 0;
1223 
1224   *A = B;
1225   if (OptionsHasName(0,"-mat_aij_superlu")) {
1226     ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr);
1227   }
1228   if (OptionsHasName(0,"-mat_aij_essl")) {
1229     ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr);
1230   }
1231   if (OptionsHasName(0,"-mat_aij_dxml")) {
1232     if (!b->indexshift) SETERRQ(1,"MatCreateSeqAIJ:need -mat_aij_oneindex with -mat_aij_dxml");
1233     ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr);
1234   }
1235 
1236   return 0;
1237 }
1238 
1239 int MatCopyPrivate_SeqAIJ(Mat A,Mat *B,int cpvalues)
1240 {
1241   Mat        C;
1242   Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data;
1243   int        i,len, m = a->m,shift = a->indexshift;
1244 
1245   *B = 0;
1246   if (!a->assembled) SETERRQ(1,"MatCopyPrivate_SeqAIJ:Cannot copy unassembled matrix");
1247   PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQAIJ,A->comm);
1248   PLogObjectCreate(C);
1249   C->data       = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c);
1250   PetscMemcpy(&C->ops,&MatOps,sizeof(struct _MatOps));
1251   C->destroy    = MatDestroy_SeqAIJ;
1252   C->view       = MatView_SeqAIJ;
1253   C->factor     = A->factor;
1254   c->row        = 0;
1255   c->col        = 0;
1256   c->indexshift = shift;
1257 
1258   c->m          = a->m;
1259   c->n          = a->n;
1260 
1261   c->imax       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax);
1262   c->ilen       = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen);
1263   for ( i=0; i<m; i++ ) {
1264     c->imax[i] = a->imax[i];
1265     c->ilen[i] = a->ilen[i];
1266   }
1267 
1268   /* allocate the matrix space */
1269   c->singlemalloc = 1;
1270   len     = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int));
1271   c->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a);
1272   c->j  = (int *) (c->a + a->i[m] + shift);
1273   c->i  = c->j + a->i[m] + shift;
1274   PetscMemcpy(c->i,a->i,(m+1)*sizeof(int));
1275   if (m > 0) {
1276     PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int));
1277     if (cpvalues == COPY_VALUES) {
1278       PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar));
1279     }
1280   }
1281 
1282   PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ));
1283   c->sorted      = a->sorted;
1284   c->roworiented = a->roworiented;
1285   c->nonew       = a->nonew;
1286 
1287   if (a->diag) {
1288     c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag);
1289     PLogObjectMemory(C,(m+1)*sizeof(int));
1290     for ( i=0; i<m; i++ ) {
1291       c->diag[i] = a->diag[i];
1292     }
1293   }
1294   else c->diag        = 0;
1295   if( a->inode.size){
1296     c->inode.size       = (int *) PetscMalloc( m *sizeof(int) ); CHKPTRQ(c->inode.size);
1297     c->inode.node_count = a->inode.node_count;
1298     PetscMemcpy( c->inode.size, a->inode.size, m*sizeof(int));
1299   } else {
1300     c->inode.size       = 0;
1301     c->inode.node_count = 0;
1302   }
1303   c->assembled        = 1;
1304   c->nz               = a->nz;
1305   c->maxnz            = a->maxnz;
1306   c->solve_work       = 0;
1307   c->spptr            = 0;      /* Dangerous -I'm throwing away a->spptr */
1308   b->inode.node_count = 0;
1309   b->inode.size       = 0;
1310 
1311   *B = C;
1312   return 0;
1313 }
1314 
1315 int MatLoad_SeqAIJ(Viewer bview,MatType type,Mat *A)
1316 {
1317   Mat_SeqAIJ   *a;
1318   Mat          B;
1319   int          i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift;
1320   PetscObject  vobj = (PetscObject) bview;
1321   MPI_Comm     comm = vobj->comm;
1322 
1323   MPI_Comm_size(comm,&size);
1324   if (size > 1) SETERRQ(1,"MatLoad_SeqAIJ:view must have one processor");
1325   ierr = ViewerFileGetDescriptor_Private(bview,&fd); CHKERRQ(ierr);
1326   ierr = SYRead(fd,header,4,SYINT); CHKERRQ(ierr);
1327   if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqAIJ:not matrix object in file");
1328   M = header[1]; N = header[2]; nz = header[3];
1329 
1330   /* read in row lengths */
1331   rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths);
1332   ierr = SYRead(fd,rowlengths,M,SYINT); CHKERRQ(ierr);
1333 
1334   /* create our matrix */
1335   ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr);
1336   B = *A;
1337   a = (Mat_SeqAIJ *) B->data;
1338   shift = a->indexshift;
1339 
1340   /* read in column indices and adjust for Fortran indexing*/
1341   ierr = SYRead(fd,a->j,nz,SYINT); CHKERRQ(ierr);
1342   if (shift) {
1343     for ( i=0; i<nz; i++ ) {
1344       a->j[i] += 1;
1345     }
1346   }
1347 
1348   /* read in nonzero values */
1349   ierr = SYRead(fd,a->a,nz,SYSCALAR); CHKERRQ(ierr);
1350 
1351   /* set matrix "i" values */
1352   a->i[0] = -shift;
1353   for ( i=1; i<= M; i++ ) {
1354     a->i[i]      = a->i[i-1] + rowlengths[i-1];
1355     a->ilen[i-1] = rowlengths[i-1];
1356   }
1357   PetscFree(rowlengths);
1358 
1359   ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
1360   ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr);
1361   return 0;
1362 }
1363 
1364 
1365 
1366