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