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