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