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