Actual source code: aij.c

  1: /*$Id: aij.c,v 1.370 2001/04/10 19:35:19 bsmith Exp $*/
  2: /*
  3:     Defines the basic matrix operations for the AIJ (compressed row)
  4:   matrix storage format.
  5: */

  7: #include "petscsys.h"
  8: #include "src/mat/impls/aij/seq/aij.h"
  9: #include "src/vec/vecimpl.h"
 10: #include "src/inline/spops.h"
 11: #include "src/inline/dot.h"
 12: #include "petscbt.h"


 15: EXTERN int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**);

 17: int MatGetRowIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *m,int **ia,int **ja,PetscTruth *done)
 18: {
 19:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
 20:   int        ierr,i,ishift;
 21: 
 23:   *m     = A->m;
 24:   if (!ia) return(0);
 25:   ishift = a->indexshift;
 26:   if (symmetric && !A->structurally_symmetric) {
 27:     MatToSymmetricIJ_SeqAIJ(A->m,a->i,a->j,ishift,oshift,ia,ja);
 28:   } else if (oshift == 0 && ishift == -1) {
 29:     int nz = a->i[A->m] - 1;
 30:     /* malloc space and  subtract 1 from i and j indices */
 31:     PetscMalloc((A->m+1)*sizeof(int),ia);
 32:     PetscMalloc((nz+1)*sizeof(int),ja);
 33:     for (i=0; i<nz; i++) (*ja)[i] = a->j[i] - 1;
 34:     for (i=0; i<A->m+1; i++) (*ia)[i] = a->i[i] - 1;
 35:   } else if (oshift == 1 && ishift == 0) {
 36:     int nz = a->i[A->m];
 37:     /* malloc space and  add 1 to i and j indices */
 38:     PetscMalloc((A->m+1)*sizeof(int),ia);
 39:     PetscMalloc((nz+1)*sizeof(int),ja);
 40:     for (i=0; i<nz; i++) (*ja)[i] = a->j[i] + 1;
 41:     for (i=0; i<A->m+1; i++) (*ia)[i] = a->i[i] + 1;
 42:   } else {
 43:     *ia = a->i; *ja = a->j;
 44:   }
 45:   return(0);
 46: }

 48: int MatRestoreRowIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *n,int **ia,int **ja,PetscTruth *done)
 49: {
 50:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
 51:   int        ishift = a->indexshift,ierr;
 52: 
 54:   if (!ia) return(0);
 55:   if ((symmetric && !A->structurally_symmetric) || (oshift == 0 && ishift == -1) || (oshift == 1 && ishift == 0)) {
 56:     PetscFree(*ia);
 57:     PetscFree(*ja);
 58:   }
 59:   return(0);
 60: }

 62: int MatGetColumnIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done)
 63: {
 64:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
 65:   int        ierr,i,ishift = a->indexshift,*collengths,*cia,*cja,n = A->n,m = A->m;
 66:   int        nz = a->i[m]+ishift,row,*jj,mr,col;
 67: 
 69:   *nn     = A->n;
 70:   if (!ia) return(0);
 71:   if (symmetric) {
 72:     MatToSymmetricIJ_SeqAIJ(A->m,a->i,a->j,ishift,oshift,ia,ja);
 73:   } else {
 74:     PetscMalloc((n+1)*sizeof(int),&collengths);
 75:     PetscMemzero(collengths,n*sizeof(int));
 76:     PetscMalloc((n+1)*sizeof(int),&cia);
 77:     PetscMalloc((nz+1)*sizeof(int),&cja);
 78:     jj = a->j;
 79:     for (i=0; i<nz; i++) {
 80:       collengths[jj[i] + ishift]++;
 81:     }
 82:     cia[0] = oshift;
 83:     for (i=0; i<n; i++) {
 84:       cia[i+1] = cia[i] + collengths[i];
 85:     }
 86:     PetscMemzero(collengths,n*sizeof(int));
 87:     jj   = a->j;
 88:     for (row=0; row<m; row++) {
 89:       mr = a->i[row+1] - a->i[row];
 90:       for (i=0; i<mr; i++) {
 91:         col = *jj++ + ishift;
 92:         cja[cia[col] + collengths[col]++ - oshift] = row + oshift;
 93:       }
 94:     }
 95:     PetscFree(collengths);
 96:     *ia = cia; *ja = cja;
 97:   }
 98:   return(0);
 99: }

101: int MatRestoreColumnIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *n,int **ia,int **ja,PetscTruth *done)
102: {

106:   if (!ia) return(0);

108:   PetscFree(*ia);
109:   PetscFree(*ja);
110: 
111:   return(0);
112: }

114: #define CHUNKSIZE   15

116: int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
117: {
118:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
119:   int        *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted = a->sorted;
120:   int        *imax = a->imax,*ai = a->i,*ailen = a->ilen;
121:   int        *aj = a->j,nonew = a->nonew,shift = a->indexshift,ierr;
122:   Scalar     *ap,value,*aa = a->a;
123:   PetscTruth ignorezeroentries = ((a->ignorezeroentries && is == ADD_VALUES) ? PETSC_TRUE:PETSC_FALSE);
124:   PetscTruth roworiented = a->roworiented;

127:   for (k=0; k<m; k++) { /* loop over added rows */
128:     row  = im[k];
129:     if (row < 0) continue;
130: #if defined(PETSC_USE_BOPT_g)  
131:     if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,A->m);
132: #endif
133:     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
134:     rmax = imax[row]; nrow = ailen[row];
135:     low = 0;
136:     for (l=0; l<n; l++) { /* loop over added columns */
137:       if (in[l] < 0) continue;
138: #if defined(PETSC_USE_BOPT_g)  
139:       if (in[l] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],A->n);
140: #endif
141:       col = in[l] - shift;
142:       if (roworiented) {
143:         value = v[l + k*n];
144:       } else {
145:         value = v[k + l*m];
146:       }
147:       if (value == 0.0 && ignorezeroentries) continue;

149:       if (!sorted) low = 0; high = nrow;
150:       while (high-low > 5) {
151:         t = (low+high)/2;
152:         if (rp[t] > col) high = t;
153:         else             low  = t;
154:       }
155:       for (i=low; i<high; i++) {
156:         if (rp[i] > col) break;
157:         if (rp[i] == col) {
158:           if (is == ADD_VALUES) ap[i] += value;
159:           else                  ap[i] = value;
160:           goto noinsert;
161:         }
162:       }
163:       if (nonew == 1) goto noinsert;
164:       else if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero at (%d,%d) in the matrix",row,col);
165:       if (nrow >= rmax) {
166:         /* there is no extra room in row, therefore enlarge */
167:         int    new_nz = ai[A->m] + CHUNKSIZE,len,*new_i,*new_j;
168:         Scalar *new_a;

170:         if (nonew == -2) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero at (%d,%d) in the matrix requiring new malloc()",row,col);

172:         /* malloc new storage space */
173:         len     = new_nz*(sizeof(int)+sizeof(Scalar))+(A->m+1)*sizeof(int);
174:         ierr    = PetscMalloc(len,&new_a);
175:         new_j   = (int*)(new_a + new_nz);
176:         new_i   = new_j + new_nz;

178:         /* copy over old data into new slots */
179:         for (ii=0; ii<row+1; ii++) {new_i[ii] = ai[ii];}
180:         for (ii=row+1; ii<A->m+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;}
181:         PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int));
182:         len  = (new_nz - CHUNKSIZE - ai[row] - nrow - shift);
183:         PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow,len*sizeof(int));
184:         PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar));
185:         PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow,len*sizeof(Scalar));
186:         /* free up old matrix storage */
187:         PetscFree(a->a);
188:         if (!a->singlemalloc) {
189:           PetscFree(a->i);
190:           PetscFree(a->j);
191:         }
192:         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
193:         a->singlemalloc = PETSC_TRUE;

195:         rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
196:         rmax = imax[row] = imax[row] + CHUNKSIZE;
197:         PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar)));
198:         a->maxnz += CHUNKSIZE;
199:         a->reallocs++;
200:       }
201:       N = nrow++ - 1; a->nz++;
202:       /* shift up all the later entries in this row */
203:       for (ii=N; ii>=i; ii--) {
204:         rp[ii+1] = rp[ii];
205:         ap[ii+1] = ap[ii];
206:       }
207:       rp[i] = col;
208:       ap[i] = value;
209:       noinsert:;
210:       low = i + 1;
211:     }
212:     ailen[row] = nrow;
213:   }
214:   return(0);
215: }

217: int MatGetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v)
218: {
219:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
220:   int        *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
221:   int        *ai = a->i,*ailen = a->ilen,shift = a->indexshift;
222:   Scalar     *ap,*aa = a->a,zero = 0.0;

225:   for (k=0; k<m; k++) { /* loop over rows */
226:     row  = im[k];
227:     if (row < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %d",row);
228:     if (row >= A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: %d",row);
229:     rp   = aj + ai[row] + shift; ap = aa + ai[row] + shift;
230:     nrow = ailen[row];
231:     for (l=0; l<n; l++) { /* loop over columns */
232:       if (in[l] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %d",in[l]);
233:       if (in[l] >= A->n) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: %d",in[l]);
234:       col = in[l] - shift;
235:       high = nrow; low = 0; /* assume unsorted */
236:       while (high-low > 5) {
237:         t = (low+high)/2;
238:         if (rp[t] > col) high = t;
239:         else             low  = t;
240:       }
241:       for (i=low; i<high; i++) {
242:         if (rp[i] > col) break;
243:         if (rp[i] == col) {
244:           *v++ = ap[i];
245:           goto finished;
246:         }
247:       }
248:       *v++ = zero;
249:       finished:;
250:     }
251:   }
252:   return(0);
253: }


256: int MatView_SeqAIJ_Binary(Mat A,PetscViewer viewer)
257: {
258:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
259:   int        i,fd,*col_lens,ierr;

262:   PetscViewerBinaryGetDescriptor(viewer,&fd);
263:   PetscMalloc((4+A->m)*sizeof(int),&col_lens);
264:   col_lens[0] = MAT_COOKIE;
265:   col_lens[1] = A->m;
266:   col_lens[2] = A->n;
267:   col_lens[3] = a->nz;

269:   /* store lengths of each row and write (including header) to file */
270:   for (i=0; i<A->m; i++) {
271:     col_lens[4+i] = a->i[i+1] - a->i[i];
272:   }
273:   PetscBinaryWrite(fd,col_lens,4+A->m,PETSC_INT,1);
274:   PetscFree(col_lens);

276:   /* store column indices (zero start index) */
277:   if (a->indexshift) {
278:     for (i=0; i<a->nz; i++) a->j[i]--;
279:   }
280:   PetscBinaryWrite(fd,a->j,a->nz,PETSC_INT,0);
281:   if (a->indexshift) {
282:     for (i=0; i<a->nz; i++) a->j[i]++;
283:   }

285:   /* store nonzero values */
286:   PetscBinaryWrite(fd,a->a,a->nz,PETSC_SCALAR,0);
287:   return(0);
288: }

290: int MatView_SeqAIJ_ASCII(Mat A,PetscViewer viewer)
291: {
292:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
293:   int               ierr,i,j,m = A->m,shift = a->indexshift;
294:   char              *name;
295:   PetscViewerFormat format;

298:   PetscObjectGetName((PetscObject)A,&name);
299:   PetscViewerGetFormat(viewer,&format);
300:   if (format == PETSC_VIEWER_ASCII_INFO_LONG || format == PETSC_VIEWER_ASCII_INFO) {
301:     if (a->inode.size) {
302:       PetscViewerASCIIPrintf(viewer,"using I-node routines: found %d nodes, limit used is %dn",a->inode.node_count,a->inode.limit);
303:     } else {
304:       PetscViewerASCIIPrintf(viewer,"not using I-node routinesn");
305:     }
306:   } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
307:     int nofinalvalue = 0;
308:     if ((a->i[m] == a->i[m-1]) || (a->j[a->nz-1] != A->n-!shift)) {
309:       nofinalvalue = 1;
310:     }
311:     PetscViewerASCIIUseTabs(viewer,PETSC_NO);
312:     PetscViewerASCIIPrintf(viewer,"%% Size = %d %d n",m,A->n);
313:     PetscViewerASCIIPrintf(viewer,"%% Nonzeros = %d n",a->nz);
314:     PetscViewerASCIIPrintf(viewer,"zzz = zeros(%d,3);n",a->nz+nofinalvalue);
315:     PetscViewerASCIIPrintf(viewer,"zzz = [n");

317:     for (i=0; i<m; i++) {
318:       for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
319: #if defined(PETSC_USE_COMPLEX)
320:         PetscViewerASCIIPrintf(viewer,"%d %d  %18.16e + %18.16ei n",i+1,a->j[j]+!shift,PetscRealPart(a->a[j]),PetscImaginaryPart(a->a[j]));
321: #else
322:         PetscViewerASCIIPrintf(viewer,"%d %d  %18.16en",i+1,a->j[j]+!shift,a->a[j]);
323: #endif
324:       }
325:     }
326:     if (nofinalvalue) {
327:       PetscViewerASCIIPrintf(viewer,"%d %d  %18.16en",m,A->n,0.0);
328:     }
329:     PetscViewerASCIIPrintf(viewer,"];n %s = spconvert(zzz);n",name);
330:     PetscViewerASCIIUseTabs(viewer,PETSC_YES);
331:   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
332:     PetscViewerASCIIUseTabs(viewer,PETSC_NO);
333:     for (i=0; i<m; i++) {
334:       PetscViewerASCIIPrintf(viewer,"row %d:",i);
335:       for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
336: #if defined(PETSC_USE_COMPLEX)
337:         if (PetscImaginaryPart(a->a[j]) > 0.0 && PetscRealPart(a->a[j]) != 0.0) {
338:           PetscViewerASCIIPrintf(viewer," %d %g + %g i",a->j[j]+shift,PetscRealPart(a->a[j]),PetscImaginaryPart(a->a[j]));
339:         } else if (PetscImaginaryPart(a->a[j]) < 0.0 && PetscRealPart(a->a[j]) != 0.0) {
340:           PetscViewerASCIIPrintf(viewer," %d %g - %g i",a->j[j]+shift,PetscRealPart(a->a[j]),-PetscImaginaryPart(a->a[j]));
341:         } else if (PetscRealPart(a->a[j]) != 0.0) {
342:           PetscViewerASCIIPrintf(viewer," %d %g ",a->j[j]+shift,PetscRealPart(a->a[j]));
343:         }
344: #else
345:         if (a->a[j] != 0.0) {PetscViewerASCIIPrintf(viewer," %d %g ",a->j[j]+shift,a->a[j]);}
346: #endif
347:       }
348:       PetscViewerASCIIPrintf(viewer,"n");
349:     }
350:     PetscViewerASCIIUseTabs(viewer,PETSC_YES);
351:   } else if (format == PETSC_VIEWER_ASCII_SYMMODU) {
352:     int nzd=0,fshift=1,*sptr;
353:     PetscViewerASCIIUseTabs(viewer,PETSC_NO);
354:     PetscMalloc((m+1)*sizeof(int),&sptr);
355:     for (i=0; i<m; i++) {
356:       sptr[i] = nzd+1;
357:       for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
358:         if (a->j[j] >= i) {
359: #if defined(PETSC_USE_COMPLEX)
360:           if (PetscImaginaryPart(a->a[j]) != 0.0 || PetscRealPart(a->a[j]) != 0.0) nzd++;
361: #else
362:           if (a->a[j] != 0.0) nzd++;
363: #endif
364:         }
365:       }
366:     }
367:     sptr[m] = nzd+1;
368:     PetscViewerASCIIPrintf(viewer," %d %dnn",m,nzd);
369:     for (i=0; i<m+1; i+=6) {
370:       if (i+4<m) {PetscViewerASCIIPrintf(viewer," %d %d %d %d %d %dn",sptr[i],sptr[i+1],sptr[i+2],sptr[i+3],sptr[i+4],sptr[i+5]);}
371:       else if (i+3<m) {PetscViewerASCIIPrintf(viewer," %d %d %d %d %dn",sptr[i],sptr[i+1],sptr[i+2],sptr[i+3],sptr[i+4]);}
372:       else if (i+2<m) {PetscViewerASCIIPrintf(viewer," %d %d %d %dn",sptr[i],sptr[i+1],sptr[i+2],sptr[i+3]);}
373:       else if (i+1<m) {PetscViewerASCIIPrintf(viewer," %d %d %dn",sptr[i],sptr[i+1],sptr[i+2]);}
374:       else if (i<m)   {PetscViewerASCIIPrintf(viewer," %d %dn",sptr[i],sptr[i+1]);}
375:       else            {PetscViewerASCIIPrintf(viewer," %dn",sptr[i]);}
376:     }
377:     PetscViewerASCIIPrintf(viewer,"n");
378:     PetscFree(sptr);
379:     for (i=0; i<m; i++) {
380:       for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
381:         if (a->j[j] >= i) {PetscViewerASCIIPrintf(viewer," %d ",a->j[j]+fshift);}
382:       }
383:       PetscViewerASCIIPrintf(viewer,"n");
384:     }
385:     PetscViewerASCIIPrintf(viewer,"n");
386:     for (i=0; i<m; i++) {
387:       for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
388:         if (a->j[j] >= i) {
389: #if defined(PETSC_USE_COMPLEX)
390:           if (PetscImaginaryPart(a->a[j]) != 0.0 || PetscRealPart(a->a[j]) != 0.0) {
391:             PetscViewerASCIIPrintf(viewer," %18.16e %18.16e ",PetscRealPart(a->a[j]),PetscImaginaryPart(a->a[j]));
392:           }
393: #else
394:           if (a->a[j] != 0.0) {PetscViewerASCIIPrintf(viewer," %18.16e ",a->a[j]);}
395: #endif
396:         }
397:       }
398:       PetscViewerASCIIPrintf(viewer,"n");
399:     }
400:     PetscViewerASCIIUseTabs(viewer,PETSC_YES);
401:   } else if (format == PETSC_VIEWER_ASCII_DENSE) {
402:     int    cnt = 0,jcnt;
403:     Scalar value;

405:     PetscViewerASCIIUseTabs(viewer,PETSC_NO);
406:     for (i=0; i<m; i++) {
407:       jcnt = 0;
408:       for (j=0; j<A->n; j++) {
409:         if (jcnt < a->i[i+1]-a->i[i] && j == a->j[cnt]) {
410:           value = a->a[cnt++];
411:           jcnt++;
412:         } else {
413:           value = 0.0;
414:         }
415: #if defined(PETSC_USE_COMPLEX)
416:         PetscViewerASCIIPrintf(viewer," %7.5e+%7.5e i ",PetscRealPart(value),PetscImaginaryPart(value));
417: #else
418:         PetscViewerASCIIPrintf(viewer," %7.5e ",value);
419: #endif
420:       }
421:       PetscViewerASCIIPrintf(viewer,"n");
422:     }
423:     PetscViewerASCIIUseTabs(viewer,PETSC_YES);
424:   } else {
425:     PetscViewerASCIIUseTabs(viewer,PETSC_NO);
426:     for (i=0; i<m; i++) {
427:       PetscViewerASCIIPrintf(viewer,"row %d:",i);
428:       for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
429: #if defined(PETSC_USE_COMPLEX)
430:         if (PetscImaginaryPart(a->a[j]) > 0.0) {
431:           PetscViewerASCIIPrintf(viewer," %d %g + %g i",a->j[j]+shift,PetscRealPart(a->a[j]),PetscImaginaryPart(a->a[j]));
432:         } else if (PetscImaginaryPart(a->a[j]) < 0.0) {
433:           PetscViewerASCIIPrintf(viewer," %d %g - %g i",a->j[j]+shift,PetscRealPart(a->a[j]),-PetscImaginaryPart(a->a[j]));
434:         } else {
435:           PetscViewerASCIIPrintf(viewer," %d %g ",a->j[j]+shift,PetscRealPart(a->a[j]));
436:         }
437: #else
438:         PetscViewerASCIIPrintf(viewer," %d %g ",a->j[j]+shift,a->a[j]);
439: #endif
440:       }
441:       PetscViewerASCIIPrintf(viewer,"n");
442:     }
443:     PetscViewerASCIIUseTabs(viewer,PETSC_YES);
444:   }
445:   PetscViewerFlush(viewer);
446:   return(0);
447: }

449: int MatView_SeqAIJ_Draw_Zoom(PetscDraw draw,void *Aa)
450: {
451:   Mat               A = (Mat) Aa;
452:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
453:   int               ierr,i,j,m = A->m,shift = a->indexshift,color;
454:   PetscReal         xl,yl,xr,yr,x_l,x_r,y_l,y_r,maxv = 0.0;
455:   PetscViewer       viewer;
456:   PetscViewerFormat format;

459:   PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);
460:   PetscViewerGetFormat(viewer,&format);

462:   PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);
463:   /* loop over matrix elements drawing boxes */

465:   if (format != PETSC_VIEWER_DRAW_CONTOUR) {
466:     /* Blue for negative, Cyan for zero and  Red for positive */
467:     color = PETSC_DRAW_BLUE;
468:     for (i=0; i<m; i++) {
469:       y_l = m - i - 1.0; y_r = y_l + 1.0;
470:       for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
471:         x_l = a->j[j] + shift; x_r = x_l + 1.0;
472: #if defined(PETSC_USE_COMPLEX)
473:         if (PetscRealPart(a->a[j]) >=  0.) continue;
474: #else
475:         if (a->a[j] >=  0.) continue;
476: #endif
477:         PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
478:       }
479:     }
480:     color = PETSC_DRAW_CYAN;
481:     for (i=0; i<m; i++) {
482:       y_l = m - i - 1.0; y_r = y_l + 1.0;
483:       for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
484:         x_l = a->j[j] + shift; x_r = x_l + 1.0;
485:         if (a->a[j] !=  0.) continue;
486:         PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
487:       }
488:     }
489:     color = PETSC_DRAW_RED;
490:     for (i=0; i<m; i++) {
491:       y_l = m - i - 1.0; y_r = y_l + 1.0;
492:       for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
493:         x_l = a->j[j] + shift; x_r = x_l + 1.0;
494: #if defined(PETSC_USE_COMPLEX)
495:         if (PetscRealPart(a->a[j]) <=  0.) continue;
496: #else
497:         if (a->a[j] <=  0.) continue;
498: #endif
499:         PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
500:       }
501:     }
502:   } else {
503:     /* use contour shading to indicate magnitude of values */
504:     /* first determine max of all nonzero values */
505:     int    nz = a->nz,count;
506:     PetscDraw   popup;
507:     PetscReal scale;

509:     for (i=0; i<nz; i++) {
510:       if (PetscAbsScalar(a->a[i]) > maxv) maxv = PetscAbsScalar(a->a[i]);
511:     }
512:     scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv;
513:     ierr  = PetscDrawGetPopup(draw,&popup);
514:     if (popup) {ierr  = PetscDrawScalePopup(popup,0.0,maxv);}
515:     count = 0;
516:     for (i=0; i<m; i++) {
517:       y_l = m - i - 1.0; y_r = y_l + 1.0;
518:       for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
519:         x_l = a->j[j] + shift; x_r = x_l + 1.0;
520:         color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(a->a[count]));
521:         ierr  = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
522:         count++;
523:       }
524:     }
525:   }
526:   return(0);
527: }

529: int MatView_SeqAIJ_Draw(Mat A,PetscViewer viewer)
530: {
531:   int        ierr;
532:   PetscDraw  draw;
533:   PetscReal  xr,yr,xl,yl,h,w;
534:   PetscTruth isnull;

537:   PetscViewerDrawGetDraw(viewer,0,&draw);
538:   PetscDrawIsNull(draw,&isnull);
539:   if (isnull) return(0);

541:   PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);
542:   xr  = A->n; yr = A->m; h = yr/10.0; w = xr/10.0;
543:   xr += w;    yr += h;  xl = -w;     yl = -h;
544:   PetscDrawSetCoordinates(draw,xl,yl,xr,yr);
545:   PetscDrawZoom(draw,MatView_SeqAIJ_Draw_Zoom,A);
546:   PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);
547:   return(0);
548: }

550: int MatView_SeqAIJ(Mat A,PetscViewer viewer)
551: {
552:   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)A->data;
553:   int         ierr;
554:   PetscTruth  issocket,isascii,isbinary,isdraw;

557:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);
558:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
559:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
560:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
561:   if (issocket) {
562:     if (a->indexshift) {
563:       SETERRQ(1,"Can only socket send sparse matrix with 0 based indexing");
564:     }
565:     PetscViewerSocketPutSparse_Private(viewer,A->m,A->n,a->nz,a->a,a->i,a->j);
566:   } else if (isascii) {
567:     MatView_SeqAIJ_ASCII(A,viewer);
568:   } else if (isbinary) {
569:     MatView_SeqAIJ_Binary(A,viewer);
570:   } else if (isdraw) {
571:     MatView_SeqAIJ_Draw(A,viewer);
572:   } else {
573:     SETERRQ1(1,"Viewer type %s not supported by SeqAIJ matrices",((PetscObject)viewer)->type_name);
574:   }
575:   return(0);
576: }

578: EXTERN int Mat_AIJ_CheckInode(Mat);
579: int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode)
580: {
581:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
582:   int        fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax,ierr;
583:   int        m = A->m,*ip,N,*ailen = a->ilen,shift = a->indexshift,rmax = 0;
584:   Scalar     *aa = a->a,*ap;

587:   if (mode == MAT_FLUSH_ASSEMBLY) return(0);

589:   if (m) rmax = ailen[0]; /* determine row with most nonzeros */
590:   for (i=1; i<m; i++) {
591:     /* move each row back by the amount of empty slots (fshift) before it*/
592:     fshift += imax[i-1] - ailen[i-1];
593:     rmax   = PetscMax(rmax,ailen[i]);
594:     if (fshift) {
595:       ip = aj + ai[i] + shift;
596:       ap = aa + ai[i] + shift;
597:       N  = ailen[i];
598:       for (j=0; j<N; j++) {
599:         ip[j-fshift] = ip[j];
600:         ap[j-fshift] = ap[j];
601:       }
602:     }
603:     ai[i] = ai[i-1] + ailen[i-1];
604:   }
605:   if (m) {
606:     fshift += imax[m-1] - ailen[m-1];
607:     ai[m]  = ai[m-1] + ailen[m-1];
608:   }
609:   /* reset ilen and imax for each row */
610:   for (i=0; i<m; i++) {
611:     ailen[i] = imax[i] = ai[i+1] - ai[i];
612:   }
613:   a->nz = ai[m] + shift;

615:   /* diagonals may have moved, so kill the diagonal pointers */
616:   if (fshift && a->diag) {
617:     PetscFree(a->diag);
618:     PetscLogObjectMemory(A,-(m+1)*sizeof(int));
619:     a->diag = 0;
620:   }
621:   PetscLogInfo(A,"MatAssemblyEnd_SeqAIJ:Matrix size: %d X %d; storage space: %d unneeded,%d usedn",m,A->n,fshift,a->nz);
622:   PetscLogInfo(A,"MatAssemblyEnd_SeqAIJ:Number of mallocs during MatSetValues() is %dn",a->reallocs);
623:   PetscLogInfo(A,"MatAssemblyEnd_SeqAIJ:Most nonzeros in any row is %dn",rmax);
624:   a->reallocs          = 0;
625:   A->info.nz_unneeded  = (double)fshift;
626:   a->rmax              = rmax;

628:   /* check out for identical nodes. If found, use inode functions */
629:   Mat_AIJ_CheckInode(A);
630:   return(0);
631: }

633: int MatZeroEntries_SeqAIJ(Mat A)
634: {
635:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
636:   int        ierr;

639:   PetscMemzero(a->a,(a->i[A->m]+a->indexshift)*sizeof(Scalar));
640:   return(0);
641: }

643: int MatDestroy_SeqAIJ(Mat A)
644: {
645:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
646:   int        ierr;

649: #if defined(PETSC_USE_LOG)
650:   PetscLogObjectState((PetscObject)A,"Rows=%d, Cols=%d, NZ=%d",A->m,A->n,a->nz);
651: #endif
652:   if (a->freedata) {
653:     PetscFree(a->a);
654:     if (!a->singlemalloc) {
655:       PetscFree(a->i);
656:       PetscFree(a->j);
657:     }
658:   }
659:   if (a->row) {
660:     ISDestroy(a->row);
661:   }
662:   if (a->col) {
663:     ISDestroy(a->col);
664:   }
665:   if (a->diag) {PetscFree(a->diag);}
666:   if (a->ilen) {PetscFree(a->ilen);}
667:   if (a->imax) {PetscFree(a->imax);}
668:   if (a->idiag) {PetscFree(a->idiag);}
669:   if (a->solve_work) {PetscFree(a->solve_work);}
670:   if (a->inode.size) {PetscFree(a->inode.size);}
671:   if (a->icol) {ISDestroy(a->icol);}
672:   if (a->saved_values) {PetscFree(a->saved_values);}
673:   PetscFree(a);
674:   return(0);
675: }

677: int MatCompress_SeqAIJ(Mat A)
678: {
680:   return(0);
681: }

683: int MatSetOption_SeqAIJ(Mat A,MatOption op)
684: {
685:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;

688:   if      (op == MAT_ROW_ORIENTED)                 a->roworiented       = PETSC_TRUE;
689:   else if (op == MAT_KEEP_ZEROED_ROWS)             a->keepzeroedrows    = PETSC_TRUE;
690:   else if (op == MAT_COLUMN_ORIENTED)              a->roworiented       = PETSC_FALSE;
691:   else if (op == MAT_COLUMNS_SORTED)               a->sorted            = PETSC_TRUE;
692:   else if (op == MAT_COLUMNS_UNSORTED)             a->sorted            = PETSC_FALSE;
693:   else if (op == MAT_NO_NEW_NONZERO_LOCATIONS)     a->nonew             = 1;
694:   else if (op == MAT_NEW_NONZERO_LOCATION_ERR)     a->nonew             = -1;
695:   else if (op == MAT_NEW_NONZERO_ALLOCATION_ERR)   a->nonew             = -2;
696:   else if (op == MAT_YES_NEW_NONZERO_LOCATIONS)    a->nonew             = 0;
697:   else if (op == MAT_IGNORE_ZERO_ENTRIES)          a->ignorezeroentries = PETSC_TRUE;
698:   else if (op == MAT_USE_INODES)                   a->inode.use         = PETSC_TRUE;
699:   else if (op == MAT_DO_NOT_USE_INODES)            a->inode.use         = PETSC_FALSE;
700:   else if (op == MAT_ROWS_SORTED ||
701:            op == MAT_ROWS_UNSORTED ||
702:            op == MAT_YES_NEW_DIAGONALS ||
703:            op == MAT_IGNORE_OFF_PROC_ENTRIES||
704:            op == MAT_USE_HASH_TABLE)
705:     PetscLogInfo(A,"MatSetOption_SeqAIJ:Option ignoredn");
706:   else if (op == MAT_NO_NEW_DIAGONALS) {
707:     SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
708:   } else if (op == MAT_INODE_LIMIT_1)          a->inode.limit  = 1;
709:   else if (op == MAT_INODE_LIMIT_2)            a->inode.limit  = 2;
710:   else if (op == MAT_INODE_LIMIT_3)            a->inode.limit  = 3;
711:   else if (op == MAT_INODE_LIMIT_4)            a->inode.limit  = 4;
712:   else if (op == MAT_INODE_LIMIT_5)            a->inode.limit  = 5;
713:   else SETERRQ(PETSC_ERR_SUP,"unknown option");
714:   return(0);
715: }

717: int MatGetDiagonal_SeqAIJ(Mat A,Vec v)
718: {
719:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
720:   int        i,j,n,shift = a->indexshift,ierr;
721:   Scalar     *x,zero = 0.0;

724:   VecSet(&zero,v);
725:   VecGetArray(v,&x);
726:   VecGetLocalSize(v,&n);
727:   if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
728:   for (i=0; i<A->m; i++) {
729:     for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
730:       if (a->j[j]+shift == i) {
731:         x[i] = a->a[j];
732:         break;
733:       }
734:     }
735:   }
736:   VecRestoreArray(v,&x);
737:   return(0);
738: }

740: int MatMultTranspose_SeqAIJ(Mat A,Vec xx,Vec yy)
741: {
742:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
743:   Scalar     *x,*y,*v,alpha,zero = 0.0;
744:   int        ierr,m = A->m,n,i,*idx,shift = a->indexshift;

747:   VecSet(&zero,yy);
748:   VecGetArray(xx,&x);
749:   VecGetArray(yy,&y);
750:   y = y + shift; /* shift for Fortran start by 1 indexing */
751:   for (i=0; i<m; i++) {
752:     idx   = a->j + a->i[i] + shift;
753:     v     = a->a + a->i[i] + shift;
754:     n     = a->i[i+1] - a->i[i];
755:     alpha = x[i];
756:     while (n-->0) {y[*idx++] += alpha * *v++;}
757:   }
758:   PetscLogFlops(2*a->nz - A->n);
759:   VecRestoreArray(xx,&x);
760:   VecRestoreArray(yy,&y);
761:   return(0);
762: }

764: int MatMultTransposeAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
765: {
766:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
767:   Scalar     *x,*y,*v,alpha;
768:   int        ierr,m = A->m,n,i,*idx,shift = a->indexshift;

771:   if (zz != yy) {VecCopy(zz,yy);}
772:   VecGetArray(xx,&x);
773:   VecGetArray(yy,&y);
774:   y = y + shift; /* shift for Fortran start by 1 indexing */
775:   for (i=0; i<m; i++) {
776:     idx   = a->j + a->i[i] + shift;
777:     v     = a->a + a->i[i] + shift;
778:     n     = a->i[i+1] - a->i[i];
779:     alpha = x[i];
780:     while (n-->0) {y[*idx++] += alpha * *v++;}
781:   }
782:   PetscLogFlops(2*a->nz);
783:   VecRestoreArray(xx,&x);
784:   VecRestoreArray(yy,&y);
785:   return(0);
786: }

788: int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
789: {
790:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
791:   Scalar     *x,*y,*v,sum;
792:   int        ierr,m = A->m,*idx,shift = a->indexshift,*ii;
793: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
794:   int        n,i,jrow,j;
795: #endif

797: #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
798: #pragma disjoint(*x,*y,*v)
799: #endif

802:   VecGetArray(xx,&x);
803:   VecGetArray(yy,&y);
804:   x    = x + shift;    /* shift for Fortran start by 1 indexing */
805:   idx  = a->j;
806:   v    = a->a;
807:   ii   = a->i;
808: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
809:   fortranmultaij_(&m,x,ii,idx+shift,v+shift,y);
810: #else
811:   v    += shift; /* shift for Fortran start by 1 indexing */
812:   idx  += shift;
813:   for (i=0; i<m; i++) {
814:     jrow = ii[i];
815:     n    = ii[i+1] - jrow;
816:     sum  = 0.0;
817:     for (j=0; j<n; j++) {
818:       sum += v[jrow]*x[idx[jrow]]; jrow++;
819:      }
820:     y[i] = sum;
821:   }
822: #endif
823:   PetscLogFlops(2*a->nz - m);
824:   VecRestoreArray(xx,&x);
825:   VecRestoreArray(yy,&y);
826:   return(0);
827: }

829: int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
830: {
831:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
832:   Scalar     *x,*y,*z,*v,sum;
833:   int        ierr,m = A->m,*idx,shift = a->indexshift,*ii;
834: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
835:   int        n,i,jrow,j;
836: #endif

839:   VecGetArray(xx,&x);
840:   VecGetArray(yy,&y);
841:   if (zz != yy) {
842:     VecGetArray(zz,&z);
843:   } else {
844:     z = y;
845:   }
846:   x    = x + shift; /* shift for Fortran start by 1 indexing */
847:   idx  = a->j;
848:   v    = a->a;
849:   ii   = a->i;
850: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
851:   fortranmultaddaij_(&m,x,ii,idx+shift,v+shift,y,z);
852: #else
853:   v   += shift; /* shift for Fortran start by 1 indexing */
854:   idx += shift;
855:   for (i=0; i<m; i++) {
856:     jrow = ii[i];
857:     n    = ii[i+1] - jrow;
858:     sum  = y[i];
859:     for (j=0; j<n; j++) {
860:       sum += v[jrow]*x[idx[jrow]]; jrow++;
861:      }
862:     z[i] = sum;
863:   }
864: #endif
865:   PetscLogFlops(2*a->nz);
866:   VecRestoreArray(xx,&x);
867:   VecRestoreArray(yy,&y);
868:   if (zz != yy) {
869:     VecRestoreArray(zz,&z);
870:   }
871:   return(0);
872: }

874: /*
875:      Adds diagonal pointers to sparse matrix structure.
876: */
877: int MatMarkDiagonal_SeqAIJ(Mat A)
878: {
879:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
880:   int        i,j,*diag,m = A->m,shift = a->indexshift,ierr;

883:   if (a->diag) return(0);

885:   PetscMalloc((m+1)*sizeof(int),&diag);
886:   PetscLogObjectMemory(A,(m+1)*sizeof(int));
887:   for (i=0; i<A->m; i++) {
888:     diag[i] = a->i[i+1];
889:     for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
890:       if (a->j[j]+shift == i) {
891:         diag[i] = j - shift;
892:         break;
893:       }
894:     }
895:   }
896:   a->diag = diag;
897:   return(0);
898: }

900: /*
901:      Checks for missing diagonals
902: */
903: int MatMissingDiagonal_SeqAIJ(Mat A)
904: {
905:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
906:   int        *diag,*jj = a->j,i,shift = a->indexshift,ierr;

909:   MatMarkDiagonal_SeqAIJ(A);
910:   diag = a->diag;
911:   for (i=0; i<A->m; i++) {
912:     if (jj[diag[i]+shift] != i-shift) {
913:       SETERRQ1(1,"Matrix is missing diagonal number %d",i);
914:     }
915:   }
916:   return(0);
917: }

919: int MatRelax_SeqAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,int its,Vec xx)
920: {
921:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
922:   Scalar     *x,*b,*bs, d,*xs,sum,*v = a->a,*t=0,scale,*ts,*xb,*idiag=0;
923:   int        ierr,*idx,*diag,n = A->n,m = A->m,i,shift = a->indexshift;

926:   VecGetArray(xx,&x);
927:   if (xx != bb) {
928:     VecGetArray(bb,&b);
929:   } else {
930:     b = x;
931:   }

933:   if (!a->diag) {MatMarkDiagonal_SeqAIJ(A);}
934:   diag = a->diag;
935:   xs   = x + shift; /* shifted by one for index start of a or a->j*/
936:   if (flag == SOR_APPLY_UPPER) {
937:    /* apply (U + D/omega) to the vector */
938:     bs = b + shift;
939:     for (i=0; i<m; i++) {
940:         d    = fshift + a->a[diag[i] + shift];
941:         n    = a->i[i+1] - diag[i] - 1;
942:         PetscLogFlops(2*n-1);
943:         idx  = a->j + diag[i] + (!shift);
944:         v    = a->a + diag[i] + (!shift);
945:         sum  = b[i]*d/omega;
946:         SPARSEDENSEDOT(sum,bs,v,idx,n);
947:         x[i] = sum;
948:     }
949:     VecRestoreArray(xx,&x);
950:     if (bb != xx) {VecRestoreArray(bb,&b);}
951:     return(0);
952:   }

954:   /* setup workspace for Eisenstat */
955:   if (flag & SOR_EISENSTAT) {
956:     if (!a->idiag) {
957:       ierr     = PetscMalloc(2*m*sizeof(Scalar),&a->idiag);
958:       a->ssor  = a->idiag + m;
959:       v        = a->a;
960:       for (i=0; i<m; i++) { a->idiag[i] = 1.0/v[diag[i]];}
961:     }
962:     t     = a->ssor;
963:     idiag = a->idiag;
964:   }
965:     /* Let  A = L + U + D; where L is lower trianglar,
966:     U is upper triangular, E is diagonal; This routine applies

968:             (L + E)^{-1} A (U + E)^{-1}

970:     to a vector efficiently using Eisenstat's trick. This is for
971:     the case of SSOR preconditioner, so E is D/omega where omega
972:     is the relaxation factor.
973:     */

975:   if (flag == SOR_APPLY_LOWER) {
976:     SETERRQ(PETSC_ERR_SUP,"SOR_APPLY_LOWER is not implemented");
977:   } else if ((flag & SOR_EISENSTAT) && omega == 1.0 && shift == 0 && fshift == 0.0) {
978:     /* special case for omega = 1.0 saves flops and some integer ops */
979:     Scalar *v2;
980: 
981:     v2    = a->a;
982:     /*  x = (E + U)^{-1} b */
983:     for (i=m-1; i>=0; i--) {
984:       n    = a->i[i+1] - diag[i] - 1;
985:       idx  = a->j + diag[i] + 1;
986:       v    = a->a + diag[i] + 1;
987:       sum  = b[i];
988:       SPARSEDENSEMDOT(sum,xs,v,idx,n);
989:       x[i] = sum*idiag[i];

991:       /*  t = b - (2*E - D)x */
992:       t[i] = b[i] - (v2[diag[i]])*x[i];
993:     }

995:     /*  t = (E + L)^{-1}t */
996:     diag = a->diag;
997:     for (i=0; i<m; i++) {
998:       n    = diag[i] - a->i[i];
999:       idx  = a->j + a->i[i];
1000:       v    = a->a + a->i[i];
1001:       sum  = t[i];
1002:       SPARSEDENSEMDOT(sum,t,v,idx,n);
1003:       t[i]  = sum*idiag[i];

1005:       /*  x = x + t */
1006:       x[i] += t[i];
1007:     }

1009:     PetscLogFlops(3*m-1 + 2*a->nz);
1010:     VecRestoreArray(xx,&x);
1011:     if (bb != xx) {VecRestoreArray(bb,&b);}
1012:     return(0);
1013:   } else if (flag & SOR_EISENSTAT) {
1014:     /* Let  A = L + U + D; where L is lower trianglar,
1015:     U is upper triangular, E is diagonal; This routine applies

1017:             (L + E)^{-1} A (U + E)^{-1}

1019:     to a vector efficiently using Eisenstat's trick. This is for
1020:     the case of SSOR preconditioner, so E is D/omega where omega
1021:     is the relaxation factor.
1022:     */
1023:     scale = (2.0/omega) - 1.0;

1025:     /*  x = (E + U)^{-1} b */
1026:     for (i=m-1; i>=0; i--) {
1027:       d    = fshift + a->a[diag[i] + shift];
1028:       n    = a->i[i+1] - diag[i] - 1;
1029:       idx  = a->j + diag[i] + (!shift);
1030:       v    = a->a + diag[i] + (!shift);
1031:       sum  = b[i];
1032:       SPARSEDENSEMDOT(sum,xs,v,idx,n);
1033:       x[i] = omega*(sum/d);
1034:     }

1036:     /*  t = b - (2*E - D)x */
1037:     v = a->a;
1038:     for (i=0; i<m; i++) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; }

1040:     /*  t = (E + L)^{-1}t */
1041:     ts = t + shift; /* shifted by one for index start of a or a->j*/
1042:     diag = a->diag;
1043:     for (i=0; i<m; i++) {
1044:       d    = fshift + a->a[diag[i]+shift];
1045:       n    = diag[i] - a->i[i];
1046:       idx  = a->j + a->i[i] + shift;
1047:       v    = a->a + a->i[i] + shift;
1048:       sum  = t[i];
1049:       SPARSEDENSEMDOT(sum,ts,v,idx,n);
1050:       t[i] = omega*(sum/d);
1051:       /*  x = x + t */
1052:       x[i] += t[i];
1053:     }

1055:     PetscLogFlops(6*m-1 + 2*a->nz);
1056:     VecRestoreArray(xx,&x);
1057:     if (bb != xx) {VecRestoreArray(bb,&b);}
1058:     return(0);
1059:   }
1060:   if (flag & SOR_ZERO_INITIAL_GUESS) {
1061:     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
1062:       for (i=0; i<m; i++) {
1063:         d    = fshift + a->a[diag[i]+shift];
1064:         n    = diag[i] - a->i[i];
1065:         PetscLogFlops(2*n-1);
1066:         idx  = a->j + a->i[i] + shift;
1067:         v    = a->a + a->i[i] + shift;
1068:         sum  = b[i];
1069:         SPARSEDENSEMDOT(sum,xs,v,idx,n);
1070:         x[i] = omega*(sum/d);
1071:       }
1072:       xb = x;
1073:     } else xb = b;
1074:     if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
1075:         (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
1076:       for (i=0; i<m; i++) {
1077:         x[i] *= a->a[diag[i]+shift];
1078:       }
1079:       PetscLogFlops(m);
1080:     }
1081:     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
1082:       for (i=m-1; i>=0; i--) {
1083:         d    = fshift + a->a[diag[i] + shift];
1084:         n    = a->i[i+1] - diag[i] - 1;
1085:         PetscLogFlops(2*n-1);
1086:         idx  = a->j + diag[i] + (!shift);
1087:         v    = a->a + diag[i] + (!shift);
1088:         sum  = xb[i];
1089:         SPARSEDENSEMDOT(sum,xs,v,idx,n);
1090:         x[i] = omega*(sum/d);
1091:       }
1092:     }
1093:     its--;
1094:   }
1095:   while (its--) {
1096:     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
1097:       for (i=0; i<m; i++) {
1098:         d    = fshift + a->a[diag[i]+shift];
1099:         n    = a->i[i+1] - a->i[i];
1100:         PetscLogFlops(2*n-1);
1101:         idx  = a->j + a->i[i] + shift;
1102:         v    = a->a + a->i[i] + shift;
1103:         sum  = b[i];
1104:         SPARSEDENSEMDOT(sum,xs,v,idx,n);
1105:         x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d;
1106:       }
1107:     }
1108:     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
1109:       for (i=m-1; i>=0; i--) {
1110:         d    = fshift + a->a[diag[i] + shift];
1111:         n    = a->i[i+1] - a->i[i];
1112:         PetscLogFlops(2*n-1);
1113:         idx  = a->j + a->i[i] + shift;
1114:         v    = a->a + a->i[i] + shift;
1115:         sum  = b[i];
1116:         SPARSEDENSEMDOT(sum,xs,v,idx,n);
1117:         x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d;
1118:       }
1119:     }
1120:   }
1121:   VecRestoreArray(xx,&x);
1122:   if (bb != xx) {VecRestoreArray(bb,&b);}
1123:   return(0);
1124: }

1126: int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,MatInfo *info)
1127: {
1128:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;

1131:   info->rows_global    = (double)A->m;
1132:   info->columns_global = (double)A->n;
1133:   info->rows_local     = (double)A->m;
1134:   info->columns_local  = (double)A->n;
1135:   info->block_size     = 1.0;
1136:   info->nz_allocated   = (double)a->maxnz;
1137:   info->nz_used        = (double)a->nz;
1138:   info->nz_unneeded    = (double)(a->maxnz - a->nz);
1139:   info->assemblies     = (double)A->num_ass;
1140:   info->mallocs        = (double)a->reallocs;
1141:   info->memory         = A->mem;
1142:   if (A->factor) {
1143:     info->fill_ratio_given  = A->info.fill_ratio_given;
1144:     info->fill_ratio_needed = A->info.fill_ratio_needed;
1145:     info->factor_mallocs    = A->info.factor_mallocs;
1146:   } else {
1147:     info->fill_ratio_given  = 0;
1148:     info->fill_ratio_needed = 0;
1149:     info->factor_mallocs    = 0;
1150:   }
1151:   return(0);
1152: }

1154: EXTERN int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,MatLUInfo*,Mat*);
1155: EXTERN int MatLUFactorNumeric_SeqAIJ(Mat,Mat*);
1156: EXTERN int MatLUFactor_SeqAIJ(Mat,IS,IS,MatLUInfo*);
1157: EXTERN int MatSolve_SeqAIJ(Mat,Vec,Vec);
1158: EXTERN int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
1159: EXTERN int MatSolveTranspose_SeqAIJ(Mat,Vec,Vec);
1160: EXTERN int MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec);

1162: int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag)
1163: {
1164:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1165:   int         i,ierr,N,*rows,m = A->m - 1,shift = a->indexshift;

1168:   ISGetLocalSize(is,&N);
1169:   ISGetIndices(is,&rows);
1170:   if (a->keepzeroedrows) {
1171:     for (i=0; i<N; i++) {
1172:       if (rows[i] < 0 || rows[i] > m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"row out of range");
1173:       PetscMemzero(&a->a[a->i[rows[i]]+shift],a->ilen[rows[i]]*sizeof(Scalar));
1174:     }
1175:     if (diag) {
1176:       MatMissingDiagonal_SeqAIJ(A);
1177:       MatMarkDiagonal_SeqAIJ(A);
1178:       for (i=0; i<N; i++) {
1179:         a->a[a->diag[rows[i]]] = *diag;
1180:       }
1181:     }
1182:   } else {
1183:     if (diag) {
1184:       for (i=0; i<N; i++) {
1185:         if (rows[i] < 0 || rows[i] > m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"row out of range");
1186:         if (a->ilen[rows[i]] > 0) {
1187:           a->ilen[rows[i]]          = 1;
1188:           a->a[a->i[rows[i]]+shift] = *diag;
1189:           a->j[a->i[rows[i]]+shift] = rows[i]+shift;
1190:         } else { /* in case row was completely empty */
1191:           MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES);
1192:         }
1193:       }
1194:     } else {
1195:       for (i=0; i<N; i++) {
1196:         if (rows[i] < 0 || rows[i] > m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"row out of range");
1197:         a->ilen[rows[i]] = 0;
1198:       }
1199:     }
1200:   }
1201:   ISRestoreIndices(is,&rows);
1202:   MatAssemblyEnd_SeqAIJ(A,MAT_FINAL_ASSEMBLY);
1203:   return(0);
1204: }

1206: int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n)
1207: {
1209:   if (m) *m = 0;
1210:   if (n) *n = A->m;
1211:   return(0);
1212: }

1214: int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
1215: {
1216:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1217:   int        *itmp,i,shift = a->indexshift,ierr;

1220:   if (row < 0 || row >= A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Row %d out of range",row);

1222:   *nz = a->i[row+1] - a->i[row];
1223:   if (v) *v = a->a + a->i[row] + shift;
1224:   if (idx) {
1225:     itmp = a->j + a->i[row] + shift;
1226:     if (*nz && shift) {
1227:       PetscMalloc((*nz)*sizeof(int),idx);
1228:       for (i=0; i<(*nz); i++) {(*idx)[i] = itmp[i] + shift;}
1229:     } else if (*nz) {
1230:       *idx = itmp;
1231:     }
1232:     else *idx = 0;
1233:   }
1234:   return(0);
1235: }

1237: int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
1238: {
1239:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;

1243:   if (idx) {if (*idx && a->indexshift) {PetscFree(*idx);}}
1244:   return(0);
1245: }

1247: int MatNorm_SeqAIJ(Mat A,NormType type,PetscReal *norm)
1248: {
1249:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1250:   Scalar     *v = a->a;
1251:   PetscReal  sum = 0.0;
1252:   int        i,j,shift = a->indexshift,ierr;

1255:   if (type == NORM_FROBENIUS) {
1256:     for (i=0; i<a->nz; i++) {
1257: #if defined(PETSC_USE_COMPLEX)
1258:       sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1259: #else
1260:       sum += (*v)*(*v); v++;
1261: #endif
1262:     }
1263:     *norm = sqrt(sum);
1264:   } else if (type == NORM_1) {
1265:     PetscReal *tmp;
1266:     int    *jj = a->j;
1267:     PetscMalloc((A->n+1)*sizeof(PetscReal),&tmp);
1268:     PetscMemzero(tmp,A->n*sizeof(PetscReal));
1269:     *norm = 0.0;
1270:     for (j=0; j<a->nz; j++) {
1271:         tmp[*jj++ + shift] += PetscAbsScalar(*v);  v++;
1272:     }
1273:     for (j=0; j<A->n; j++) {
1274:       if (tmp[j] > *norm) *norm = tmp[j];
1275:     }
1276:     PetscFree(tmp);
1277:   } else if (type == NORM_INFINITY) {
1278:     *norm = 0.0;
1279:     for (j=0; j<A->m; j++) {
1280:       v = a->a + a->i[j] + shift;
1281:       sum = 0.0;
1282:       for (i=0; i<a->i[j+1]-a->i[j]; i++) {
1283:         sum += PetscAbsScalar(*v); v++;
1284:       }
1285:       if (sum > *norm) *norm = sum;
1286:     }
1287:   } else {
1288:     SETERRQ(PETSC_ERR_SUP,"No support for two norm");
1289:   }
1290:   return(0);
1291: }

1293: int MatTranspose_SeqAIJ(Mat A,Mat *B)
1294: {
1295:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1296:   Mat        C;
1297:   int        i,ierr,*aj = a->j,*ai = a->i,m = A->m,len,*col;
1298:   int        shift = a->indexshift;
1299:   Scalar     *array = a->a;

1302:   if (!B && m != A->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Square matrix only for in-place");
1303:   PetscMalloc((1+A->n)*sizeof(int),&col);
1304:   PetscMemzero(col,(1+A->n)*sizeof(int));
1305:   if (shift) {
1306:     for (i=0; i<ai[m]-1; i++) aj[i] -= 1;
1307:   }
1308:   for (i=0; i<ai[m]+shift; i++) col[aj[i]] += 1;
1309:   MatCreateSeqAIJ(A->comm,A->n,m,0,col,&C);
1310:   PetscFree(col);
1311:   for (i=0; i<m; i++) {
1312:     len    = ai[i+1]-ai[i];
1313:     ierr   = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES);
1314:     array += len;
1315:     aj    += len;
1316:   }
1317:   if (shift) {
1318:     for (i=0; i<ai[m]-1; i++) aj[i] += 1;
1319:   }

1321:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1322:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

1324:   if (B) {
1325:     *B = C;
1326:   } else {
1327:     MatHeaderCopy(A,C);
1328:   }
1329:   return(0);
1330: }

1332: int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr)
1333: {
1334:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1335:   Scalar     *l,*r,x,*v;
1336:   int        ierr,i,j,m = A->m,n = A->n,M,nz = a->nz,*jj,shift = a->indexshift;

1339:   if (ll) {
1340:     /* The local size is used so that VecMPI can be passed to this routine
1341:        by MatDiagonalScale_MPIAIJ */
1342:     VecGetLocalSize(ll,&m);
1343:     if (m != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length");
1344:     VecGetArray(ll,&l);
1345:     v = a->a;
1346:     for (i=0; i<m; i++) {
1347:       x = l[i];
1348:       M = a->i[i+1] - a->i[i];
1349:       for (j=0; j<M; j++) { (*v++) *= x;}
1350:     }
1351:     VecRestoreArray(ll,&l);
1352:     PetscLogFlops(nz);
1353:   }
1354:   if (rr) {
1355:     VecGetLocalSize(rr,&n);
1356:     if (n != A->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vector wrong length");
1357:     VecGetArray(rr,&r);
1358:     v = a->a; jj = a->j;
1359:     for (i=0; i<nz; i++) {
1360:       (*v++) *= r[*jj++ + shift];
1361:     }
1362:     VecRestoreArray(rr,&r);
1363:     PetscLogFlops(nz);
1364:   }
1365:   return(0);
1366: }

1368: int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,int csize,MatReuse scall,Mat *B)
1369: {
1370:   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data,*c;
1371:   int          *smap,i,k,kstart,kend,ierr,oldcols = A->n,*lens;
1372:   int          row,mat_i,*mat_j,tcol,first,step,*mat_ilen,sum,lensi;
1373:   int          *irow,*icol,nrows,ncols,shift = a->indexshift,*ssmap;
1374:   int          *starts,*j_new,*i_new,*aj = a->j,*ai = a->i,ii,*ailen = a->ilen;
1375:   Scalar       *a_new,*mat_a;
1376:   Mat          C;
1377:   PetscTruth   stride;

1380:   ISSorted(isrow,(PetscTruth*)&i);
1381:   if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"ISrow is not sorted");
1382:   ISSorted(iscol,(PetscTruth*)&i);
1383:   if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IScol is not sorted");

1385:   ISGetIndices(isrow,&irow);
1386:   ISGetLocalSize(isrow,&nrows);
1387:   ISGetLocalSize(iscol,&ncols);

1389:   ISStrideGetInfo(iscol,&first,&step);
1390:   ISStride(iscol,&stride);
1391:   if (stride && step == 1) {
1392:     /* special case of contiguous rows */
1393:     ierr   = PetscMalloc((ncols+nrows+1)*sizeof(int),&lens);
1394:     starts = lens + ncols;
1395:     /* loop over new rows determining lens and starting points */
1396:     for (i=0; i<nrows; i++) {
1397:       kstart  = ai[irow[i]]+shift;
1398:       kend    = kstart + ailen[irow[i]];
1399:       for (k=kstart; k<kend; k++) {
1400:         if (aj[k]+shift >= first) {
1401:           starts[i] = k;
1402:           break;
1403:         }
1404:       }
1405:       sum = 0;
1406:       while (k < kend) {
1407:         if (aj[k++]+shift >= first+ncols) break;
1408:         sum++;
1409:       }
1410:       lens[i] = sum;
1411:     }
1412:     /* create submatrix */
1413:     if (scall == MAT_REUSE_MATRIX) {
1414:       int n_cols,n_rows;
1415:       MatGetSize(*B,&n_rows,&n_cols);
1416:       if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size");
1417:       MatZeroEntries(*B);
1418:       C = *B;
1419:     } else {
1420:       MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);
1421:     }
1422:     c = (Mat_SeqAIJ*)C->data;

1424:     /* loop over rows inserting into submatrix */
1425:     a_new    = c->a;
1426:     j_new    = c->j;
1427:     i_new    = c->i;
1428:     i_new[0] = -shift;
1429:     for (i=0; i<nrows; i++) {
1430:       ii    = starts[i];
1431:       lensi = lens[i];
1432:       for (k=0; k<lensi; k++) {
1433:         *j_new++ = aj[ii+k] - first;
1434:       }
1435:       PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar));
1436:       a_new      += lensi;
1437:       i_new[i+1]  = i_new[i] + lensi;
1438:       c->ilen[i]  = lensi;
1439:     }
1440:     PetscFree(lens);
1441:   } else {
1442:     ierr  = ISGetIndices(iscol,&icol);
1443:     ierr  = PetscMalloc((1+oldcols)*sizeof(int),&smap);
1444:     ssmap = smap + shift;
1445:     ierr  = PetscMalloc((1+nrows)*sizeof(int),&lens);
1446:     ierr  = PetscMemzero(smap,oldcols*sizeof(int));
1447:     for (i=0; i<ncols; i++) smap[icol[i]] = i+1;
1448:     /* determine lens of each row */
1449:     for (i=0; i<nrows; i++) {
1450:       kstart  = ai[irow[i]]+shift;
1451:       kend    = kstart + a->ilen[irow[i]];
1452:       lens[i] = 0;
1453:       for (k=kstart; k<kend; k++) {
1454:         if (ssmap[aj[k]]) {
1455:           lens[i]++;
1456:         }
1457:       }
1458:     }
1459:     /* Create and fill new matrix */
1460:     if (scall == MAT_REUSE_MATRIX) {
1461:       PetscTruth equal;

1463:       c = (Mat_SeqAIJ *)((*B)->data);
1464:       if ((*B)->m  != nrows || (*B)->n != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size");
1465:       PetscMemcmp(c->ilen,lens,(*B)->m*sizeof(int),&equal);
1466:       if (!equal) {
1467:         SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros");
1468:       }
1469:       PetscMemzero(c->ilen,(*B)->m*sizeof(int));
1470:       C = *B;
1471:     } else {
1472:       MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);
1473:     }
1474:     c = (Mat_SeqAIJ *)(C->data);
1475:     for (i=0; i<nrows; i++) {
1476:       row    = irow[i];
1477:       kstart = ai[row]+shift;
1478:       kend   = kstart + a->ilen[row];
1479:       mat_i  = c->i[i]+shift;
1480:       mat_j  = c->j + mat_i;
1481:       mat_a  = c->a + mat_i;
1482:       mat_ilen = c->ilen + i;
1483:       for (k=kstart; k<kend; k++) {
1484:         if ((tcol=ssmap[a->j[k]])) {
1485:           *mat_j++ = tcol - (!shift);
1486:           *mat_a++ = a->a[k];
1487:           (*mat_ilen)++;

1489:         }
1490:       }
1491:     }
1492:     /* Free work space */
1493:     ISRestoreIndices(iscol,&icol);
1494:     PetscFree(smap);
1495:     PetscFree(lens);
1496:   }
1497:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1498:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

1500:   ISRestoreIndices(isrow,&irow);
1501:   *B = C;
1502:   return(0);
1503: }

1505: /*
1506: */
1507: int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,MatILUInfo *info)
1508: {
1509:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)inA->data;
1510:   int        ierr;
1511:   Mat        outA;
1512:   PetscTruth row_identity,col_identity;

1515:   if (info && info->levels != 0) SETERRQ(PETSC_ERR_SUP,"Only levels=0 supported for in-place ilu");
1516:   ISIdentity(row,&row_identity);
1517:   ISIdentity(col,&col_identity);
1518:   if (!row_identity || !col_identity) {
1519:     SETERRQ(1,"Row and column permutations must be identity for in-place ILU");
1520:   }

1522:   outA          = inA;
1523:   inA->factor   = FACTOR_LU;
1524:   a->row        = row;
1525:   a->col        = col;
1526:   PetscObjectReference((PetscObject)row);
1527:   PetscObjectReference((PetscObject)col);

1529:   /* Create the inverse permutation so that it can be used in MatLUFactorNumeric() */
1530:   if (a->icol) {ISDestroy(a->icol);} /* need to remove old one */
1531:   ISInvertPermutation(col,PETSC_DECIDE,&a->icol);
1532:   PetscLogObjectParent(inA,a->icol);

1534:   if (!a->solve_work) { /* this matrix may have been factored before */
1535:      PetscMalloc((inA->m+1)*sizeof(Scalar),&a->solve_work);
1536:   }

1538:   if (!a->diag) {
1539:     MatMarkDiagonal_SeqAIJ(inA);
1540:   }
1541:   MatLUFactorNumeric_SeqAIJ(inA,&outA);
1542:   return(0);
1543: }

1545: #include "petscblaslapack.h"
1546: int MatScale_SeqAIJ(Scalar *alpha,Mat inA)
1547: {
1548:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)inA->data;
1549:   int        one = 1;

1552:   BLscal_(&a->nz,alpha,a->a,&one);
1553:   PetscLogFlops(a->nz);
1554:   return(0);
1555: }

1557: int MatGetSubMatrices_SeqAIJ(Mat A,int n,IS *irow,IS *icol,MatReuse scall,Mat **B)
1558: {
1559:   int ierr,i;

1562:   if (scall == MAT_INITIAL_MATRIX) {
1563:     PetscMalloc((n+1)*sizeof(Mat),B);
1564:   }

1566:   for (i=0; i<n; i++) {
1567:     MatGetSubMatrix_SeqAIJ(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);
1568:   }
1569:   return(0);
1570: }

1572: int MatGetBlockSize_SeqAIJ(Mat A,int *bs)
1573: {
1575:   *bs = 1;
1576:   return(0);
1577: }

1579: int MatIncreaseOverlap_SeqAIJ(Mat A,int is_max,IS *is,int ov)
1580: {
1581:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1582:   int        shift,row,i,j,k,l,m,n,*idx,ierr,*nidx,isz,val;
1583:   int        start,end,*ai,*aj;
1584:   PetscBT    table;

1587:   shift = a->indexshift;
1588:   m     = A->m;
1589:   ai    = a->i;
1590:   aj    = a->j+shift;

1592:   if (ov < 0)  SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"illegal overlap value used");

1594:   PetscMalloc((m+1)*sizeof(int),&nidx);
1595:   PetscBTCreate(m,table);

1597:   for (i=0; i<is_max; i++) {
1598:     /* Initialize the two local arrays */
1599:     isz  = 0;
1600:     PetscBTMemzero(m,table);
1601: 
1602:     /* Extract the indices, assume there can be duplicate entries */
1603:     ISGetIndices(is[i],&idx);
1604:     ISGetLocalSize(is[i],&n);
1605: 
1606:     /* Enter these into the temp arrays. I.e., mark table[row], enter row into new index */
1607:     for (j=0; j<n ; ++j){
1608:       if(!PetscBTLookupSet(table,idx[j])) { nidx[isz++] = idx[j];}
1609:     }
1610:     ISRestoreIndices(is[i],&idx);
1611:     ISDestroy(is[i]);
1612: 
1613:     k = 0;
1614:     for (j=0; j<ov; j++){ /* for each overlap */
1615:       n = isz;
1616:       for (; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */
1617:         row   = nidx[k];
1618:         start = ai[row];
1619:         end   = ai[row+1];
1620:         for (l = start; l<end ; l++){
1621:           val = aj[l] + shift;
1622:           if (!PetscBTLookupSet(table,val)) {nidx[isz++] = val;}
1623:         }
1624:       }
1625:     }
1626:     ISCreateGeneral(PETSC_COMM_SELF,isz,nidx,(is+i));
1627:   }
1628:   PetscBTDestroy(table);
1629:   PetscFree(nidx);
1630:   return(0);
1631: }

1633: /* -------------------------------------------------------------- */
1634: int MatPermute_SeqAIJ(Mat A,IS rowp,IS colp,Mat *B)
1635: {
1636:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1637:   Scalar     *vwork;
1638:   int        i,ierr,nz,m = A->m,n = A->n,*cwork;
1639:   int        *row,*col,*cnew,j,*lens;
1640:   IS         icolp,irowp;

1643:   ISInvertPermutation(rowp,PETSC_DECIDE,&irowp);
1644:   ISGetIndices(irowp,&row);
1645:   ISInvertPermutation(colp,PETSC_DECIDE,&icolp);
1646:   ISGetIndices(icolp,&col);
1647: 
1648:   /* determine lengths of permuted rows */
1649:   PetscMalloc((m+1)*sizeof(int),&lens);
1650:   for (i=0; i<m; i++) {
1651:     lens[row[i]] = a->i[i+1] - a->i[i];
1652:   }
1653:   MatCreateSeqAIJ(A->comm,m,n,0,lens,B);
1654:   PetscFree(lens);

1656:   PetscMalloc(n*sizeof(int),&cnew);
1657:   for (i=0; i<m; i++) {
1658:     MatGetRow(A,i,&nz,&cwork,&vwork);
1659:     for (j=0; j<nz; j++) { cnew[j] = col[cwork[j]];}
1660:     MatSetValues(*B,1,&row[i],nz,cnew,vwork,INSERT_VALUES);
1661:     MatRestoreRow(A,i,&nz,&cwork,&vwork);
1662:   }
1663:   PetscFree(cnew);
1664:   MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);
1665:   MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);
1666:   ISRestoreIndices(irowp,&row);
1667:   ISRestoreIndices(icolp,&col);
1668:   ISDestroy(irowp);
1669:   ISDestroy(icolp);
1670:   return(0);
1671: }

1673: int MatPrintHelp_SeqAIJ(Mat A)
1674: {
1675:   static PetscTruth called = PETSC_FALSE;
1676:   MPI_Comm          comm = A->comm;
1677:   int               ierr;

1680:   if (called) {return(0);} else called = PETSC_TRUE;
1681:   (*PetscHelpPrintf)(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):n");
1682:   (*PetscHelpPrintf)(comm,"  -mat_lu_pivotthreshold <threshold>: Set pivoting thresholdn");
1683:   (*PetscHelpPrintf)(comm,"  -mat_aij_oneindex: internal indices begin at 1 instead of the default 0.n");
1684:   (*PetscHelpPrintf)(comm,"  -mat_aij_no_inode: Do not use inodesn");
1685:   (*PetscHelpPrintf)(comm,"  -mat_aij_inode_limit <limit>: Set inode limit (max limit=5)n");
1686: #if defined(PETSC_HAVE_ESSL)
1687:   (*PetscHelpPrintf)(comm,"  -mat_aij_essl: Use IBM sparse LU factorization and solve.n");
1688: #endif
1689: #if defined(PETSC_HAVE_LUSOL)
1690:   (*PetscHelpPrintf)(comm,"  -mat_aij_lusol: Use the Stanford LUSOL sparse factorization and solve.n");
1691: #endif
1692: #if defined(PETSC_HAVE_MATLAB_ENGINE) && !defined(PETSC_USE_COMPLEX)
1693:   (*PetscHelpPrintf)(comm,"  -mat_aij_matlab: Use Matlab engine sparse LU factorization and solve.n");
1694: #endif
1695:   return(0);
1696: }
1697: EXTERN int MatEqual_SeqAIJ(Mat A,Mat B,PetscTruth* flg);
1698: EXTERN int MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring);
1699: EXTERN int MatILUDTFactor_SeqAIJ(Mat,MatILUInfo*,IS,IS,Mat*);
1700: int MatCopy_SeqAIJ(Mat A,Mat B,MatStructure str)
1701: {
1702:   int        ierr;
1703:   PetscTruth flg;

1706:   PetscTypeCompare((PetscObject)B,MATSEQAIJ,&flg);
1707:   if (str == SAME_NONZERO_PATTERN && flg) {
1708:     Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1709:     Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data;

1711:     if (a->i[A->m]+a->indexshift != b->i[B->m]+a->indexshift) {
1712:       SETERRQ(1,"Number of nonzeros in two matrices are different");
1713:     }
1714:     PetscMemcpy(b->a,a->a,(a->i[A->m]+a->indexshift)*sizeof(Scalar));
1715:   } else {
1716:     MatCopy_Basic(A,B,str);
1717:   }
1718:   return(0);
1719: }

1721: int MatSetUpPreallocation_SeqAIJ(Mat A)
1722: {
1723:   int        ierr;

1726:    MatSeqAIJSetPreallocation(A,PETSC_DEFAULT,0);
1727:   return(0);
1728: }

1730: int MatGetArray_SeqAIJ(Mat A,Scalar **array)
1731: {
1732:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1734:   *array = a->a;
1735:   return(0);
1736: }

1738: int MatRestoreArray_SeqAIJ(Mat A,Scalar **array)
1739: {
1741:   return(0);
1742: }

1744: /* -------------------------------------------------------------------*/
1745: static struct _MatOps MatOps_Values = {MatSetValues_SeqAIJ,
1746:        MatGetRow_SeqAIJ,
1747:        MatRestoreRow_SeqAIJ,
1748:        MatMult_SeqAIJ,
1749:        MatMultAdd_SeqAIJ,
1750:        MatMultTranspose_SeqAIJ,
1751:        MatMultTransposeAdd_SeqAIJ,
1752:        MatSolve_SeqAIJ,
1753:        MatSolveAdd_SeqAIJ,
1754:        MatSolveTranspose_SeqAIJ,
1755:        MatSolveTransposeAdd_SeqAIJ,
1756:        MatLUFactor_SeqAIJ,
1757:        0,
1758:        MatRelax_SeqAIJ,
1759:        MatTranspose_SeqAIJ,
1760:        MatGetInfo_SeqAIJ,
1761:        MatEqual_SeqAIJ,
1762:        MatGetDiagonal_SeqAIJ,
1763:        MatDiagonalScale_SeqAIJ,
1764:        MatNorm_SeqAIJ,
1765:        0,
1766:        MatAssemblyEnd_SeqAIJ,
1767:        MatCompress_SeqAIJ,
1768:        MatSetOption_SeqAIJ,
1769:        MatZeroEntries_SeqAIJ,
1770:        MatZeroRows_SeqAIJ,
1771:        MatLUFactorSymbolic_SeqAIJ,
1772:        MatLUFactorNumeric_SeqAIJ,
1773:        0,
1774:        0,
1775:        MatSetUpPreallocation_SeqAIJ,
1776:        0,
1777:        MatGetOwnershipRange_SeqAIJ,
1778:        MatILUFactorSymbolic_SeqAIJ,
1779:        0,
1780:        MatGetArray_SeqAIJ,
1781:        MatRestoreArray_SeqAIJ,
1782:        MatDuplicate_SeqAIJ,
1783:        0,
1784:        0,
1785:        MatILUFactor_SeqAIJ,
1786:        0,
1787:        0,
1788:        MatGetSubMatrices_SeqAIJ,
1789:        MatIncreaseOverlap_SeqAIJ,
1790:        MatGetValues_SeqAIJ,
1791:        MatCopy_SeqAIJ,
1792:        MatPrintHelp_SeqAIJ,
1793:        MatScale_SeqAIJ,
1794:        0,
1795:        0,
1796:        MatILUDTFactor_SeqAIJ,
1797:        MatGetBlockSize_SeqAIJ,
1798:        MatGetRowIJ_SeqAIJ,
1799:        MatRestoreRowIJ_SeqAIJ,
1800:        MatGetColumnIJ_SeqAIJ,
1801:        MatRestoreColumnIJ_SeqAIJ,
1802:        MatFDColoringCreate_SeqAIJ,
1803:        0,
1804:        0,
1805:        MatPermute_SeqAIJ,
1806:        0,
1807:        0,
1808:        MatDestroy_SeqAIJ,
1809:        MatView_SeqAIJ,
1810:        MatGetMaps_Petsc};

1812: EXTERN int MatUseSuperLU_SeqAIJ(Mat);
1813: EXTERN int MatUseEssl_SeqAIJ(Mat);
1814: EXTERN int MatUseLUSOL_SeqAIJ(Mat);
1815: EXTERN int MatUseMatlab_SeqAIJ(Mat);
1816: EXTERN int MatUseDXML_SeqAIJ(Mat);

1818: EXTERN_C_BEGIN

1820: int MatSeqAIJSetColumnIndices_SeqAIJ(Mat mat,int *indices)
1821: {
1822:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
1823:   int        i,nz,n;


1827:   nz = aij->maxnz;
1828:   n  = mat->n;
1829:   for (i=0; i<nz; i++) {
1830:     aij->j[i] = indices[i];
1831:   }
1832:   aij->nz = nz;
1833:   for (i=0; i<n; i++) {
1834:     aij->ilen[i] = aij->imax[i];
1835:   }

1837:   return(0);
1838: }
1839: EXTERN_C_END

1841: /*@
1842:     MatSeqAIJSetColumnIndices - Set the column indices for all the rows
1843:        in the matrix.

1845:   Input Parameters:
1846: +  mat - the SeqAIJ matrix
1847: -  indices - the column indices

1849:   Level: advanced

1851:   Notes:
1852:     This can be called if you have precomputed the nonzero structure of the 
1853:   matrix and want to provide it to the matrix object to improve the performance
1854:   of the MatSetValues() operation.

1856:     You MUST have set the correct numbers of nonzeros per row in the call to 
1857:   MatCreateSeqAIJ().

1859:     MUST be called before any calls to MatSetValues();

1861:     The indices should start with zero, not one.

1863: @*/
1864: int MatSeqAIJSetColumnIndices(Mat mat,int *indices)
1865: {
1866:   int ierr,(*f)(Mat,int *);

1870:   PetscObjectQueryFunction((PetscObject)mat,"MatSeqAIJSetColumnIndices_C",(void (**)())&f);
1871:   if (f) {
1872:     (*f)(mat,indices);
1873:   } else {
1874:     SETERRQ(1,"Wrong type of matrix to set column indices");
1875:   }
1876:   return(0);
1877: }

1879: /* ----------------------------------------------------------------------------------------*/

1881: EXTERN_C_BEGIN
1882: int MatStoreValues_SeqAIJ(Mat mat)
1883: {
1884:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
1885:   int        nz = aij->i[mat->m]+aij->indexshift,ierr;

1888:   if (aij->nonew != 1) {
1889:     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
1890:   }

1892:   /* allocate space for values if not already there */
1893:   if (!aij->saved_values) {
1894:     PetscMalloc((nz+1)*sizeof(Scalar),&aij->saved_values);
1895:   }

1897:   /* copy values over */
1898:   PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(Scalar));
1899:   return(0);
1900: }
1901: EXTERN_C_END

1903: /*@
1904:     MatStoreValues - Stashes a copy of the matrix values; this allows, for 
1905:        example, reuse of the linear part of a Jacobian, while recomputing the 
1906:        nonlinear portion.

1908:    Collect on Mat

1910:   Input Parameters:
1911: .  mat - the matrix (currently on AIJ matrices support this option)

1913:   Level: advanced

1915:   Common Usage, with SNESSolve():
1916: $    Create Jacobian matrix
1917: $    Set linear terms into matrix
1918: $    Apply boundary conditions to matrix, at this time matrix must have 
1919: $      final nonzero structure (i.e. setting the nonlinear terms and applying 
1920: $      boundary conditions again will not change the nonzero structure
1921: $    MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS);
1922: $    MatStoreValues(mat);
1923: $    Call SNESSetJacobian() with matrix
1924: $    In your Jacobian routine
1925: $      MatRetrieveValues(mat);
1926: $      Set nonlinear terms in matrix
1927:  
1928:   Common Usage without SNESSolve(), i.e. when you handle nonlinear solve yourself:
1929: $    // build linear portion of Jacobian 
1930: $    MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS);
1931: $    MatStoreValues(mat);
1932: $    loop over nonlinear iterations
1933: $       MatRetrieveValues(mat);
1934: $       // call MatSetValues(mat,...) to set nonliner portion of Jacobian 
1935: $       // call MatAssemblyBegin/End() on matrix
1936: $       Solve linear system with Jacobian
1937: $    endloop 

1939:   Notes:
1940:     Matrix must already be assemblied before calling this routine
1941:     Must set the matrix option MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS); before 
1942:     calling this routine.

1944: .seealso: MatRetrieveValues()

1946: @*/
1947: int MatStoreValues(Mat mat)
1948: {
1949:   int ierr,(*f)(Mat);

1953:   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
1954:   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");

1956:   PetscObjectQueryFunction((PetscObject)mat,"MatStoreValues_C",(void (**)())&f);
1957:   if (f) {
1958:     (*f)(mat);
1959:   } else {
1960:     SETERRQ(1,"Wrong type of matrix to store values");
1961:   }
1962:   return(0);
1963: }

1965: EXTERN_C_BEGIN
1966: int MatRetrieveValues_SeqAIJ(Mat mat)
1967: {
1968:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
1969:   int        nz = aij->i[mat->m]+aij->indexshift,ierr;

1972:   if (aij->nonew != 1) {
1973:     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
1974:   }
1975:   if (!aij->saved_values) {
1976:     SETERRQ(1,"Must call MatStoreValues(A);first");
1977:   }

1979:   /* copy values over */
1980:   PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(Scalar));
1981:   return(0);
1982: }
1983: EXTERN_C_END

1985: /*@
1986:     MatRetrieveValues - Retrieves the copy of the matrix values; this allows, for 
1987:        example, reuse of the linear part of a Jacobian, while recomputing the 
1988:        nonlinear portion.

1990:    Collect on Mat

1992:   Input Parameters:
1993: .  mat - the matrix (currently on AIJ matrices support this option)

1995:   Level: advanced

1997: .seealso: MatStoreValues()

1999: @*/
2000: int MatRetrieveValues(Mat mat)
2001: {
2002:   int ierr,(*f)(Mat);

2006:   if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
2007:   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");

2009:   PetscObjectQueryFunction((PetscObject)mat,"MatRetrieveValues_C",(void (**)())&f);
2010:   if (f) {
2011:     (*f)(mat);
2012:   } else {
2013:     SETERRQ(1,"Wrong type of matrix to retrieve values");
2014:   }
2015:   return(0);
2016: }

2018: /*
2019:    This allows SeqAIJ matrices to be passed to the matlab engine
2020: */
2021: #if defined(PETSC_HAVE_MATLAB_ENGINE) && !defined(PETSC_USE_COMPLEX)
2022: #include "engine.h"   /* Matlab include file */
2023: #include "mex.h"      /* Matlab include file */
2024: EXTERN_C_BEGIN
2025: int MatMatlabEnginePut_SeqAIJ(PetscObject obj,void *engine)
2026: {
2027:   int        ierr,i,*ai,*aj;
2028:   Mat        B = (Mat)obj;
2029:   Scalar     *array;
2030:   mxArray    *mat;
2031:   Mat_SeqAIJ *aij = (Mat_SeqAIJ*)B->data;

2034:   mat  = mxCreateSparse(B->n,B->m,aij->nz,mxREAL);
2035:   PetscMemcpy(mxGetPr(mat),aij->a,aij->nz*sizeof(Scalar));
2036:   /* Matlab stores by column, not row so we pass in the transpose of the matrix */
2037:   PetscMemcpy(mxGetIr(mat),aij->j,aij->nz*sizeof(int));
2038:   PetscMemcpy(mxGetJc(mat),aij->i,(B->m+1)*sizeof(int));

2040:   /* Matlab indices start at 0 for sparse (what a surprise) */
2041:   if (aij->indexshift) {
2042:     for (i=0; i<B->m+1; i++) {
2043:       ai[i]--;
2044:     }
2045:     for (i=0; i<aij->nz; i++) {
2046:       aj[i]--;
2047:     }
2048:   }
2049:   PetscObjectName(obj);
2050:   mxSetName(mat,obj->name);
2051:   engPutArray((Engine *)engine,mat);
2052:   return(0);
2053: }
2054: EXTERN_C_END

2056: EXTERN_C_BEGIN
2057: int MatMatlabEngineGet_SeqAIJ(PetscObject obj,void *engine)
2058: {
2059:   int        ierr,ii;
2060:   Mat        mat = (Mat)obj;
2061:   Mat_SeqAIJ *aij = (Mat_SeqAIJ*)mat->data;
2062:   mxArray    *mmat;

2065:   PetscFree(aij->a);
2066:   aij->indexshift = 0;

2068:   mmat = engGetArray((Engine *)engine,obj->name);

2070:   aij->nz           = (mxGetJc(mmat))[mat->m];
2071:   ierr              = PetscMalloc(aij->nz*(sizeof(int)+sizeof(Scalar))+(mat->m+1)*sizeof(int),&aij->a);
2072:   aij->j            = (int*)(aij->a + aij->nz);
2073:   aij->i            = aij->j + aij->nz;
2074:   aij->singlemalloc = PETSC_TRUE;
2075:   aij->freedata     = PETSC_TRUE;

2077:   PetscMemcpy(aij->a,mxGetPr(mmat),aij->nz*sizeof(Scalar));
2078:   /* Matlab stores by column, not row so we pass in the transpose of the matrix */
2079:   PetscMemcpy(aij->j,mxGetIr(mmat),aij->nz*sizeof(int));
2080:   PetscMemcpy(aij->i,mxGetJc(mmat),(mat->m+1)*sizeof(int));

2082:   for (ii=0; ii<mat->m; ii++) {
2083:     aij->ilen[ii] = aij->imax[ii] = aij->i[ii+1] - aij->i[ii];
2084:   }

2086:   MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
2087:   MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);

2089:   return(0);
2090: }
2091: EXTERN_C_END
2092: #endif

2094: /* --------------------------------------------------------------------------------*/

2096: /*@C
2097:    MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format
2098:    (the default parallel PETSc format).  For good matrix assembly performance
2099:    the user should preallocate the matrix storage by setting the parameter nz
2100:    (or the array nnz).  By setting these parameters accurately, performance
2101:    during matrix assembly can be increased by more than a factor of 50.

2103:    Collective on MPI_Comm

2105:    Input Parameters:
2106: +  comm - MPI communicator, set to PETSC_COMM_SELF
2107: .  m - number of rows
2108: .  n - number of columns
2109: .  nz - number of nonzeros per row (same for all rows)
2110: -  nnz - array containing the number of nonzeros in the various rows 
2111:          (possibly different for each row) or PETSC_NULL

2113:    Output Parameter:
2114: .  A - the matrix 

2116:    Notes:
2117:    The AIJ format (also called the Yale sparse matrix format or
2118:    compressed row storage), is fully compatible with standard Fortran 77
2119:    storage.  That is, the stored row and column indices can begin at
2120:    either one (as in Fortran) or zero.  See the users' manual for details.

2122:    Specify the preallocated storage with either nz or nnz (not both).
2123:    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 
2124:    allocation.  For large problems you MUST preallocate memory or you 
2125:    will get TERRIBLE performance, see the users' manual chapter on matrices.

2127:    By default, this format uses inodes (identical nodes) when possible, to 
2128:    improve numerical efficiency of matrix-vector products and solves. We 
2129:    search for consecutive rows with the same nonzero structure, thereby
2130:    reusing matrix information to achieve increased efficiency.

2132:    Options Database Keys:
2133: +  -mat_aij_no_inode  - Do not use inodes
2134: .  -mat_aij_inode_limit <limit> - Sets inode limit (max limit=5)
2135: -  -mat_aij_oneindex - Internally use indexing starting at 1
2136:         rather than 0.  Note that when calling MatSetValues(),
2137:         the user still MUST index entries starting at 0!

2139:    Level: intermediate

2141: .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues(), MatSeqAIJSetColumnIndices(), MatCreateSeqAIJWithArrays()

2143: @*/
2144: int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz,Mat *A)
2145: {

2149:   MatCreate(comm,m,n,m,n,A);
2150:   MatSetType(*A,MATSEQAIJ);
2151:   MatSeqAIJSetPreallocation(*A,nz,nnz);
2152:   return(0);
2153: }

2155: /*@C
2156:    MatSeqAIJSetPreallocation - For good matrix assembly performance
2157:    the user should preallocate the matrix storage by setting the parameter nz
2158:    (or the array nnz).  By setting these parameters accurately, performance
2159:    during matrix assembly can be increased by more than a factor of 50.

2161:    Collective on MPI_Comm

2163:    Input Parameters:
2164: +  comm - MPI communicator, set to PETSC_COMM_SELF
2165: .  m - number of rows
2166: .  n - number of columns
2167: .  nz - number of nonzeros per row (same for all rows)
2168: -  nnz - array containing the number of nonzeros in the various rows 
2169:          (possibly different for each row) or PETSC_NULL

2171:    Output Parameter:
2172: .  A - the matrix 

2174:    Notes:
2175:    The AIJ format (also called the Yale sparse matrix format or
2176:    compressed row storage), is fully compatible with standard Fortran 77
2177:    storage.  That is, the stored row and column indices can begin at
2178:    either one (as in Fortran) or zero.  See the users' manual for details.

2180:    Specify the preallocated storage with either nz or nnz (not both).
2181:    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 
2182:    allocation.  For large problems you MUST preallocate memory or you 
2183:    will get TERRIBLE performance, see the users' manual chapter on matrices.

2185:    By default, this format uses inodes (identical nodes) when possible, to 
2186:    improve numerical efficiency of matrix-vector products and solves. We 
2187:    search for consecutive rows with the same nonzero structure, thereby
2188:    reusing matrix information to achieve increased efficiency.

2190:    Options Database Keys:
2191: +  -mat_aij_no_inode  - Do not use inodes
2192: .  -mat_aij_inode_limit <limit> - Sets inode limit (max limit=5)
2193: -  -mat_aij_oneindex - Internally use indexing starting at 1
2194:         rather than 0.  Note that when calling MatSetValues(),
2195:         the user still MUST index entries starting at 0!

2197:    Level: intermediate

2199: .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues(), MatSeqAIJSetColumnIndices(), MatCreateSeqAIJWithArrays()

2201: @*/
2202: int MatSeqAIJSetPreallocation(Mat B,int nz,int *nnz)
2203: {
2204:   Mat_SeqAIJ *b;
2205:   int        i,len,ierr;
2206:   PetscTruth flg2;

2209:   PetscTypeCompare((PetscObject)B,MATSEQAIJ,&flg2);
2210:   if (!flg2) return(0);
2211: 
2212:   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
2213:   if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz);
2214:   if (nnz) {
2215:     for (i=0; i<B->m; i++) {
2216:       if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]);
2217:     }
2218:   }

2220:   B->preallocated = PETSC_TRUE;
2221:   b = (Mat_SeqAIJ*)B->data;

2223:   PetscMalloc((B->m+1)*sizeof(int),&b->imax);
2224:   if (!nnz) {
2225:     if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 10;
2226:     else if (nz <= 0)        nz = 1;
2227:     for (i=0; i<B->m; i++) b->imax[i] = nz;
2228:     nz = nz*B->m;
2229:   } else {
2230:     nz = 0;
2231:     for (i=0; i<B->m; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
2232:   }

2234:   /* allocate the matrix space */
2235:   len             = nz*(sizeof(int) + sizeof(Scalar)) + (B->m+1)*sizeof(int);
2236:   ierr            = PetscMalloc(len,&b->a);
2237:   b->j            = (int*)(b->a + nz);
2238:   ierr            = PetscMemzero(b->j,nz*sizeof(int));
2239:   b->i            = b->j + nz;
2240:   b->singlemalloc = PETSC_TRUE;
2241:   b->freedata     = PETSC_TRUE;

2243:   b->i[0] = -b->indexshift;
2244:   for (i=1; i<B->m+1; i++) {
2245:     b->i[i] = b->i[i-1] + b->imax[i-1];
2246:   }

2248:   /* b->ilen will count nonzeros in each row so far. */
2249:   PetscMalloc((B->m+1)*sizeof(int),&b->ilen);
2250:   PetscLogObjectMemory(B,len+2*(B->m+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqAIJ));
2251:   for (i=0; i<B->m; i++) { b->ilen[i] = 0;}

2253:   b->nz                = 0;
2254:   b->maxnz             = nz;
2255:   B->info.nz_unneeded  = (double)b->maxnz;
2256:   return(0);
2257: }

2259: EXTERN_C_BEGIN
2260: int MatCreate_SeqAIJ(Mat B)
2261: {
2262:   Mat_SeqAIJ *b;
2263:   int        ierr,size;
2264:   PetscTruth flg;

2267:   MPI_Comm_size(B->comm,&size);
2268:   if (size > 1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Comm must be of size 1");

2270:   B->m = B->M = PetscMax(B->m,B->M);
2271:   B->n = B->N = PetscMax(B->n,B->N);

2273:   PetscNew(Mat_SeqAIJ,&b);
2274:   B->data             = (void*)b;
2275:   PetscMemzero(b,sizeof(Mat_SeqAIJ));
2276:   PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
2277:   B->factor           = 0;
2278:   B->lupivotthreshold = 1.0;
2279:   B->mapping          = 0;
2280:   PetscOptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold,PETSC_NULL);
2281:   PetscOptionsHasName(PETSC_NULL,"-pc_ilu_preserve_row_sums",&b->ilu_preserve_row_sums);
2282:   b->row              = 0;
2283:   b->col              = 0;
2284:   b->icol             = 0;
2285:   b->indexshift       = 0;
2286:   b->reallocs         = 0;
2287:   PetscOptionsHasName(PETSC_NULL,"-mat_aij_oneindex",&flg);
2288:   if (flg) b->indexshift = -1;
2289: 
2290:   MapCreateMPI(B->comm,B->m,B->m,&B->rmap);
2291:   MapCreateMPI(B->comm,B->n,B->n,&B->cmap);

2293:   b->sorted            = PETSC_FALSE;
2294:   b->ignorezeroentries = PETSC_FALSE;
2295:   b->roworiented       = PETSC_TRUE;
2296:   b->nonew             = 0;
2297:   b->diag              = 0;
2298:   b->solve_work        = 0;
2299:   b->spptr             = 0;
2300:   b->inode.use         = PETSC_TRUE;
2301:   b->inode.node_count  = 0;
2302:   b->inode.size        = 0;
2303:   b->inode.limit       = 5;
2304:   b->inode.max_limit   = 5;
2305:   b->saved_values      = 0;
2306:   b->idiag             = 0;
2307:   b->ssor              = 0;
2308:   b->keepzeroedrows    = PETSC_FALSE;

2310:   PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);

2312:   /*  SuperLU is not currently supported through PETSc */
2313: #if defined(PETSC_HAVE_SUPERLU)
2314:   PetscOptionsHasName(PETSC_NULL,"-mat_aij_superlu",&flg);
2315:   if (flg) { MatUseSuperLU_SeqAIJ(B); }
2316: #endif
2317:   PetscOptionsHasName(PETSC_NULL,"-mat_aij_essl",&flg);
2318:   if (flg) { MatUseEssl_SeqAIJ(B); }
2319:   PetscOptionsHasName(PETSC_NULL,"-mat_aij_lusol",&flg);
2320:   if (flg) { MatUseLUSOL_SeqAIJ(B); }
2321:   PetscOptionsHasName(PETSC_NULL,"-mat_aij_matlab",&flg);
2322:   if (flg) {MatUseMatlab_SeqAIJ(B);}
2323:   PetscOptionsHasName(PETSC_NULL,"-mat_aij_dxml",&flg);
2324:   if (flg) {
2325:     if (!b->indexshift) SETERRQ(PETSC_ERR_LIB,"need -mat_aij_oneindex with -mat_aij_dxml");
2326:     MatUseDXML_SeqAIJ(B);
2327:   }
2328:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqAIJSetColumnIndices_C",
2329:                                      "MatSeqAIJSetColumnIndices_SeqAIJ",
2330:                                      MatSeqAIJSetColumnIndices_SeqAIJ);
2331:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
2332:                                      "MatStoreValues_SeqAIJ",
2333:                                      MatStoreValues_SeqAIJ);
2334:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
2335:                                      "MatRetrieveValues_SeqAIJ",
2336:                                      MatRetrieveValues_SeqAIJ);
2337: #if defined(PETSC_HAVE_MATLAB_ENGINE) && !defined(PETSC_USE_COMPLEX)
2338:   PetscObjectComposeFunctionDynamic((PetscObject)B,"PetscMatlabEnginePut_C","MatMatlabEnginePut_SeqAIJ",MatMatlabEnginePut_SeqAIJ);
2339:   PetscObjectComposeFunctionDynamic((PetscObject)B,"PetscMatlabEngineGet_C","MatMatlabEngineGet_SeqAIJ",MatMatlabEngineGet_SeqAIJ);
2340: #endif
2341:   return(0);
2342: }
2343: EXTERN_C_END

2345: int MatDuplicate_SeqAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
2346: {
2347:   Mat        C;
2348:   Mat_SeqAIJ *c,*a = (Mat_SeqAIJ*)A->data;
2349:   int        i,len,m = A->m,shift = a->indexshift,ierr;

2352:   *B = 0;
2353:   MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);
2354:   MatSetType(C,MATSEQAIJ);
2355:   c    = (Mat_SeqAIJ*)C->data;

2357:   C->factor         = A->factor;
2358:   c->row            = 0;
2359:   c->col            = 0;
2360:   c->icol           = 0;
2361:   c->indexshift     = shift;
2362:   c->keepzeroedrows = a->keepzeroedrows;
2363:   C->assembled      = PETSC_TRUE;

2365:   C->M          = A->m;
2366:   C->N          = A->n;

2368:   PetscMalloc((m+1)*sizeof(int),&c->imax);
2369:   PetscMalloc((m+1)*sizeof(int),&c->ilen);
2370:   for (i=0; i<m; i++) {
2371:     c->imax[i] = a->imax[i];
2372:     c->ilen[i] = a->ilen[i];
2373:   }

2375:   /* allocate the matrix space */
2376:   c->singlemalloc = PETSC_TRUE;
2377:   len   = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int));
2378:   ierr  = PetscMalloc(len,&c->a);
2379:   c->j  = (int*)(c->a + a->i[m] + shift);
2380:   c->i  = c->j + a->i[m] + shift;
2381:   PetscMemcpy(c->i,a->i,(m+1)*sizeof(int));
2382:   if (m > 0) {
2383:     PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int));
2384:     if (cpvalues == MAT_COPY_VALUES) {
2385:       PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar));
2386:     } else {
2387:       PetscMemzero(c->a,(a->i[m]+shift)*sizeof(Scalar));
2388:     }
2389:   }

2391:   PetscLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqAIJ));
2392:   c->sorted      = a->sorted;
2393:   c->roworiented = a->roworiented;
2394:   c->nonew       = a->nonew;
2395:   c->ilu_preserve_row_sums = a->ilu_preserve_row_sums;
2396:   c->saved_values = 0;
2397:   c->idiag        = 0;
2398:   c->ssor         = 0;
2399:   c->ignorezeroentries = a->ignorezeroentries;
2400:   c->freedata     = PETSC_TRUE;

2402:   if (a->diag) {
2403:     PetscMalloc((m+1)*sizeof(int),&c->diag);
2404:     PetscLogObjectMemory(C,(m+1)*sizeof(int));
2405:     for (i=0; i<m; i++) {
2406:       c->diag[i] = a->diag[i];
2407:     }
2408:   } else c->diag        = 0;
2409:   c->inode.use          = a->inode.use;
2410:   c->inode.limit        = a->inode.limit;
2411:   c->inode.max_limit    = a->inode.max_limit;
2412:   if (a->inode.size){
2413:     ierr                = PetscMalloc((m+1)*sizeof(int),&c->inode.size);
2414:     c->inode.node_count = a->inode.node_count;
2415:     ierr                = PetscMemcpy(c->inode.size,a->inode.size,(m+1)*sizeof(int));
2416:   } else {
2417:     c->inode.size       = 0;
2418:     c->inode.node_count = 0;
2419:   }
2420:   c->nz                 = a->nz;
2421:   c->maxnz              = a->maxnz;
2422:   c->solve_work         = 0;
2423:   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
2424:   C->preallocated       = PETSC_TRUE;

2426:   *B = C;
2427:   PetscFListDuplicate(A->qlist,&C->qlist);
2428:   return(0);
2429: }

2431: EXTERN_C_BEGIN
2432: int MatLoad_SeqAIJ(PetscViewer viewer,MatType type,Mat *A)
2433: {
2434:   Mat_SeqAIJ   *a;
2435:   Mat          B;
2436:   int          i,nz,ierr,fd,header[4],size,*rowlengths = 0,M,N,shift;
2437:   MPI_Comm     comm;
2438: 
2440:   PetscObjectGetComm((PetscObject)viewer,&comm);
2441:   MPI_Comm_size(comm,&size);
2442:   if (size > 1) SETERRQ(PETSC_ERR_ARG_SIZ,"view must have one processor");
2443:   PetscViewerBinaryGetDescriptor(viewer,&fd);
2444:   PetscBinaryRead(fd,header,4,PETSC_INT);
2445:   if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object in file");
2446:   M = header[1]; N = header[2]; nz = header[3];

2448:   if (nz < 0) {
2449:     SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format on disk,cannot load as SeqAIJ");
2450:   }

2452:   /* read in row lengths */
2453:   PetscMalloc(M*sizeof(int),&rowlengths);
2454:   PetscBinaryRead(fd,rowlengths,M,PETSC_INT);

2456:   /* create our matrix */
2457:   MatCreateSeqAIJ(comm,M,N,0,rowlengths,A);
2458:   B = *A;
2459:   a = (Mat_SeqAIJ*)B->data;
2460:   shift = a->indexshift;

2462:   /* read in column indices and adjust for Fortran indexing*/
2463:   PetscBinaryRead(fd,a->j,nz,PETSC_INT);
2464:   if (shift) {
2465:     for (i=0; i<nz; i++) {
2466:       a->j[i] += 1;
2467:     }
2468:   }

2470:   /* read in nonzero values */
2471:   PetscBinaryRead(fd,a->a,nz,PETSC_SCALAR);

2473:   /* set matrix "i" values */
2474:   a->i[0] = -shift;
2475:   for (i=1; i<= M; i++) {
2476:     a->i[i]      = a->i[i-1] + rowlengths[i-1];
2477:     a->ilen[i-1] = rowlengths[i-1];
2478:   }
2479:   PetscFree(rowlengths);

2481:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
2482:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
2483:   return(0);
2484: }
2485: EXTERN_C_END

2487: int MatEqual_SeqAIJ(Mat A,Mat B,PetscTruth* flg)
2488: {
2489:   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data,*b = (Mat_SeqAIJ *)B->data;
2490:   int        ierr;
2491:   PetscTruth flag;

2494:   PetscTypeCompare((PetscObject)B,MATSEQAIJ,&flag);
2495:   if (!flag) SETERRQ(PETSC_ERR_ARG_INCOMP,"Matrices must be same type");

2497:   /* If the  matrix dimensions are not equal,or no of nonzeros or shift */
2498:   if ((A->m != B->m) || (A->n != B->n) ||(a->nz != b->nz)|| (a->indexshift != b->indexshift)) {
2499:     *flg = PETSC_FALSE;
2500:     return(0);
2501:   }
2502: 
2503:   /* if the a->i are the same */
2504:   PetscMemcmp(a->i,b->i,(A->m+1)*sizeof(int),flg);
2505:   if (*flg == PETSC_FALSE) return(0);
2506: 
2507:   /* if a->j are the same */
2508:   PetscMemcmp(a->j,b->j,(a->nz)*sizeof(int),flg);
2509:   if (*flg == PETSC_FALSE) return(0);
2510: 
2511:   /* if a->a are the same */
2512:   PetscMemcmp(a->a,b->a,(a->nz)*sizeof(Scalar),flg);

2514:   return(0);
2515: 
2516: }

2518: /*@C
2519:      MatCreateSeqAIJWithArrays - Creates an sequential AIJ matrix using matrix elements (in CSR format)
2520:               provided by the user.

2522:       Coolective on MPI_Comm

2524:    Input Parameters:
2525: +   comm - must be an MPI communicator of size 1
2526: .   m - number of rows
2527: .   n - number of columns
2528: .   i - row indices
2529: .   j - column indices
2530: -   a - matrix values

2532:    Output Parameter:
2533: .   mat - the matrix

2535:    Level: intermediate

2537:    Notes:
2538:        The i, j, and a arrays are not copied by this routine, the user must free these arrays
2539:     once the matrix is destroyed

2541:        You cannot set new nonzero locations into this matrix, that will generate an error.

2543:        The i and j indices can be either 0- or 1 based

2545: .seealso: MatCreate(), MatCreateMPIAIJ(), MatCreateSeqAIJ()

2547: @*/
2548: int MatCreateSeqAIJWithArrays(MPI_Comm comm,int m,int n,int* i,int*j,Scalar *a,Mat *mat)
2549: {
2550:   int        ierr,ii;
2551:   Mat_SeqAIJ *aij;

2554:   MatCreateSeqAIJ(comm,m,n,0,0,mat);
2555:   aij  = (Mat_SeqAIJ*)(*mat)->data;
2556:   PetscFree(aij->a);

2558:   if (i[0] == 1) {
2559:     aij->indexshift = -1;
2560:   } else if (i[0]) {
2561:     SETERRQ(1,"i (row indices) do not start with 0 or 1");
2562:   }
2563:   aij->i = i;
2564:   aij->j = j;
2565:   aij->a = a;
2566:   aij->singlemalloc = PETSC_FALSE;
2567:   aij->nonew        = -1;             /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/
2568:   aij->freedata     = PETSC_FALSE;

2570:   for (ii=0; ii<m; ii++) {
2571:     aij->ilen[ii] = aij->imax[ii] = i[ii+1] - i[ii];
2572: #if defined(PETSC_USE_BOPT_g)
2573:     if (i[ii+1] - i[ii] < 0) SETERRQ2(1,"Negative row length in i (row indices) row = %d length = %d",ii,i[ii+1] - i[ii]);
2574: #endif    
2575:   }
2576: #if defined(PETSC_USE_BOPT_g)
2577:   for (ii=0; ii<aij->i[m]; ii++) {
2578:     if (j[ii] < -aij->indexshift) SETERRQ2(1,"Negative column index at location = %d index = %d",ii,j[ii]);
2579:     if (j[ii] > n - 1 -aij->indexshift) SETERRQ2(1,"Column index to large at location = %d index = %d",ii,j[ii]);
2580:   }
2581: #endif    

2583:   /* changes indices to start at 0 */
2584:   if (i[0]) {
2585:     aij->indexshift = 0;
2586:     for (ii=0; ii<m; ii++) {
2587:       i[ii]--;
2588:     }
2589:     for (ii=0; ii<i[m]; ii++) {
2590:       j[ii]--;
2591:     }
2592:   }

2594:   MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
2595:   MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
2596:   return(0);
2597: }