Actual source code: aij.c

  1: /*$Id: aij.c,v 1.385 2001/09/07 20:09:22 bsmith Exp $*/
  2: /*
  3:     Defines the basic matrix operations for the AIJ (compressed row)
  4:   matrix storage format.
  5: */

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


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

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

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

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

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

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

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

113: #define CHUNKSIZE   15

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

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

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

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

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

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

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

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

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


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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

577: EXTERN int Mat_AIJ_CheckInode(Mat,PetscTruth);
578: EXTERN int MatUseSuperLU_DIST_MPIAIJ(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:   PetscScalar  *aa = a->a,*ap;
585: #if defined(PETSC_HAVE_SUPERLUDIST) 
586:   PetscTruth   flag;
587: #endif

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

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

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

631:   /* check out for identical nodes. If found, use inode functions */
632:   Mat_AIJ_CheckInode(A,(PetscTruth)(!fshift));
633: #if defined(PETSC_HAVE_SUPERLUDIST) 
634:   PetscOptionsHasName(PETSC_NULL,"-mat_aij_superlu_dist",&flag);
635:   if (flag) { MatUseSuperLU_DIST_MPIAIJ(A); }
636: #endif 

638:   return(0);
639: }

641: int MatZeroEntries_SeqAIJ(Mat A)
642: {
643:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
644:   int        ierr;

647:   PetscMemzero(a->a,(a->i[A->m]+a->indexshift)*sizeof(PetscScalar));
648:   return(0);
649: }

651: int MatDestroy_SeqAIJ(Mat A)
652: {
653:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
654:   int        ierr;

657: #if defined(PETSC_USE_LOG)
658:   PetscLogObjectState((PetscObject)A,"Rows=%d, Cols=%d, NZ=%d",A->m,A->n,a->nz);
659: #endif
660:   if (a->freedata) {
661:     PetscFree(a->a);
662:     if (!a->singlemalloc) {
663:       PetscFree(a->i);
664:       PetscFree(a->j);
665:     }
666:   }
667:   if (a->row) {
668:     ISDestroy(a->row);
669:   }
670:   if (a->col) {
671:     ISDestroy(a->col);
672:   }
673:   if (a->diag) {PetscFree(a->diag);}
674:   if (a->ilen) {PetscFree(a->ilen);}
675:   if (a->imax) {PetscFree(a->imax);}
676:   if (a->idiag) {PetscFree(a->idiag);}
677:   if (a->solve_work) {PetscFree(a->solve_work);}
678:   if (a->inode.size) {PetscFree(a->inode.size);}
679:   if (a->icol) {ISDestroy(a->icol);}
680:   if (a->saved_values) {PetscFree(a->saved_values);}
681:   if (a->coloring) {ISColoringDestroy(a->coloring);}
682:   PetscFree(a);
683:   return(0);
684: }

686: int MatCompress_SeqAIJ(Mat A)
687: {
689:   return(0);
690: }

692: int MatSetOption_SeqAIJ(Mat A,MatOption op)
693: {
694:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;

697:   switch (op) {
698:     case MAT_ROW_ORIENTED:
699:       a->roworiented       = PETSC_TRUE;
700:       break;
701:     case MAT_KEEP_ZEROED_ROWS:
702:       a->keepzeroedrows    = PETSC_TRUE;
703:       break;
704:     case MAT_COLUMN_ORIENTED:
705:       a->roworiented       = PETSC_FALSE;
706:       break;
707:     case MAT_COLUMNS_SORTED:
708:       a->sorted            = PETSC_TRUE;
709:       break;
710:     case MAT_COLUMNS_UNSORTED:
711:       a->sorted            = PETSC_FALSE;
712:       break;
713:     case MAT_NO_NEW_NONZERO_LOCATIONS:
714:       a->nonew             = 1;
715:       break;
716:     case MAT_NEW_NONZERO_LOCATION_ERR:
717:       a->nonew             = -1;
718:       break;
719:     case MAT_NEW_NONZERO_ALLOCATION_ERR:
720:       a->nonew             = -2;
721:       break;
722:     case MAT_YES_NEW_NONZERO_LOCATIONS:
723:       a->nonew             = 0;
724:       break;
725:     case MAT_IGNORE_ZERO_ENTRIES:
726:       a->ignorezeroentries = PETSC_TRUE;
727:       break;
728:     case MAT_USE_INODES:
729:       a->inode.use         = PETSC_TRUE;
730:       break;
731:     case MAT_DO_NOT_USE_INODES:
732:       a->inode.use         = PETSC_FALSE;
733:       break;
734:     case MAT_ROWS_SORTED:
735:     case MAT_ROWS_UNSORTED:
736:     case MAT_YES_NEW_DIAGONALS:
737:     case MAT_IGNORE_OFF_PROC_ENTRIES:
738:     case MAT_USE_HASH_TABLE:
739:     case MAT_USE_SINGLE_PRECISION_SOLVES:
740:       PetscLogInfo(A,"MatSetOption_SeqAIJ:Option ignoredn");
741:       break;
742:     case MAT_NO_NEW_DIAGONALS:
743:       SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
744:     case MAT_INODE_LIMIT_1:
745:       a->inode.limit  = 1;
746:       break;
747:     case MAT_INODE_LIMIT_2:
748:       a->inode.limit  = 2;
749:       break;
750:     case MAT_INODE_LIMIT_3:
751:       a->inode.limit  = 3;
752:       break;
753:     case MAT_INODE_LIMIT_4:
754:       a->inode.limit  = 4;
755:       break;
756:     case MAT_INODE_LIMIT_5:
757:       a->inode.limit  = 5;
758:       break;
759:     default:
760:       SETERRQ(PETSC_ERR_SUP,"unknown option");
761:   }
762:   return(0);
763: }

765: int MatGetDiagonal_SeqAIJ(Mat A,Vec v)
766: {
767:   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data;
768:   int          i,j,n,shift = a->indexshift,ierr;
769:   PetscScalar  *x,zero = 0.0;

772:   VecSet(&zero,v);
773:   VecGetArray(v,&x);
774:   VecGetLocalSize(v,&n);
775:   if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
776:   for (i=0; i<A->m; i++) {
777:     for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
778:       if (a->j[j]+shift == i) {
779:         x[i] = a->a[j];
780:         break;
781:       }
782:     }
783:   }
784:   VecRestoreArray(v,&x);
785:   return(0);
786: }


789: int MatMultTransposeAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
790: {
791:   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data;
792:   PetscScalar  *x,*y;
793:   int          ierr,m = A->m,shift = a->indexshift;
794: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
795:   PetscScalar  *v,alpha;
796:   int          n,i,*idx;
797: #endif

800:   if (zz != yy) {VecCopy(zz,yy);}
801:   VecGetArray(xx,&x);
802:   VecGetArray(yy,&y);
803:   y = y + shift; /* shift for Fortran start by 1 indexing */

805: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
806:   fortranmulttransposeaddaij_(&m,x,a->i,a->j+shift,a->a+shift,y);
807: #else
808:   for (i=0; i<m; i++) {
809:     idx   = a->j + a->i[i] + shift;
810:     v     = a->a + a->i[i] + shift;
811:     n     = a->i[i+1] - a->i[i];
812:     alpha = x[i];
813:     while (n-->0) {y[*idx++] += alpha * *v++;}
814:   }
815: #endif
816:   PetscLogFlops(2*a->nz);
817:   VecRestoreArray(xx,&x);
818:   VecRestoreArray(yy,&y);
819:   return(0);
820: }

822: int MatMultTranspose_SeqAIJ(Mat A,Vec xx,Vec yy)
823: {
824:   PetscScalar  zero = 0.0;
825:   int          ierr;

828:   VecSet(&zero,yy);
829:   MatMultTransposeAdd_SeqAIJ(A,xx,yy,yy);
830:   return(0);
831: }


834: int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
835: {
836:   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data;
837:   PetscScalar  *x,*y,*v,sum;
838:   int          ierr,m = A->m,*idx,shift = a->indexshift,*ii;
839: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
840:   int          n,i,jrow,j;
841: #endif

843: #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
844: #pragma disjoint(*x,*y,*v)
845: #endif

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

875: int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
876: {
877:   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data;
878:   PetscScalar  *x,*y,*z,*v,sum;
879:   int          ierr,m = A->m,*idx,shift = a->indexshift,*ii;
880: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
881:   int          n,i,jrow,j;
882: #endif

885:   VecGetArray(xx,&x);
886:   VecGetArray(yy,&y);
887:   if (zz != yy) {
888:     VecGetArray(zz,&z);
889:   } else {
890:     z = y;
891:   }
892:   x    = x + shift; /* shift for Fortran start by 1 indexing */
893:   idx  = a->j;
894:   v    = a->a;
895:   ii   = a->i;
896: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
897:   fortranmultaddaij_(&m,x,ii,idx+shift,v+shift,y,z);
898: #else
899:   v   += shift; /* shift for Fortran start by 1 indexing */
900:   idx += shift;
901:   for (i=0; i<m; i++) {
902:     jrow = ii[i];
903:     n    = ii[i+1] - jrow;
904:     sum  = y[i];
905:     for (j=0; j<n; j++) {
906:       sum += v[jrow]*x[idx[jrow]]; jrow++;
907:      }
908:     z[i] = sum;
909:   }
910: #endif
911:   PetscLogFlops(2*a->nz);
912:   VecRestoreArray(xx,&x);
913:   VecRestoreArray(yy,&y);
914:   if (zz != yy) {
915:     VecRestoreArray(zz,&z);
916:   }
917:   return(0);
918: }

920: /*
921:      Adds diagonal pointers to sparse matrix structure.
922: */
923: int MatMarkDiagonal_SeqAIJ(Mat A)
924: {
925:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
926:   int        i,j,*diag,m = A->m,shift = a->indexshift,ierr;

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

931:   PetscMalloc((m+1)*sizeof(int),&diag);
932:   PetscLogObjectMemory(A,(m+1)*sizeof(int));
933:   for (i=0; i<A->m; i++) {
934:     diag[i] = a->i[i+1];
935:     for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
936:       if (a->j[j]+shift == i) {
937:         diag[i] = j - shift;
938:         break;
939:       }
940:     }
941:   }
942:   a->diag = diag;
943:   return(0);
944: }

946: /*
947:      Checks for missing diagonals
948: */
949: int MatMissingDiagonal_SeqAIJ(Mat A)
950: {
951:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
952:   int        *diag,*jj = a->j,i,shift = a->indexshift,ierr;

955:   MatMarkDiagonal_SeqAIJ(A);
956:   diag = a->diag;
957:   for (i=0; i<A->m; i++) {
958:     if (jj[diag[i]+shift] != i-shift) {
959:       SETERRQ1(1,"Matrix is missing diagonal number %d",i);
960:     }
961:   }
962:   return(0);
963: }

965: int MatRelax_SeqAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,int its,int lits,Vec xx)
966: {
967:   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data;
968:   PetscScalar  *x,*b,*bs, d,*xs,sum,*v = a->a,*t=0,scale,*ts,*xb,*idiag=0;
969:   int          ierr,*idx,*diag,n = A->n,m = A->m,i,shift = a->indexshift;

972:   its = its*lits;
973:   if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %d and local its %d both positive",its,lits);

975:   VecGetArray(xx,&x);
976:   if (xx != bb) {
977:     VecGetArray(bb,&b);
978:   } else {
979:     b = x;
980:   }

982:   if (!a->diag) {MatMarkDiagonal_SeqAIJ(A);}
983:   diag = a->diag;
984:   xs   = x + shift; /* shifted by one for index start of a or a->j*/
985:   if (flag == SOR_APPLY_UPPER) {
986:    /* apply (U + D/omega) to the vector */
987:     bs = b + shift;
988:     for (i=0; i<m; i++) {
989:         d    = fshift + a->a[diag[i] + shift];
990:         n    = a->i[i+1] - diag[i] - 1;
991:         PetscLogFlops(2*n-1);
992:         idx  = a->j + diag[i] + (!shift);
993:         v    = a->a + diag[i] + (!shift);
994:         sum  = b[i]*d/omega;
995:         SPARSEDENSEDOT(sum,bs,v,idx,n);
996:         x[i] = sum;
997:     }
998:     VecRestoreArray(xx,&x);
999:     if (bb != xx) {VecRestoreArray(bb,&b);}
1000:     return(0);
1001:   }

1003:   /* setup workspace for Eisenstat */
1004:   if (flag & SOR_EISENSTAT) {
1005:     if (!a->idiag) {
1006:       ierr     = PetscMalloc(2*m*sizeof(PetscScalar),&a->idiag);
1007:       a->ssor  = a->idiag + m;
1008:       v        = a->a;
1009:       for (i=0; i<m; i++) { a->idiag[i] = 1.0/v[diag[i]];}
1010:     }
1011:     t     = a->ssor;
1012:     idiag = a->idiag;
1013:   }
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:     */

1024:   if (flag == SOR_APPLY_LOWER) {
1025:     SETERRQ(PETSC_ERR_SUP,"SOR_APPLY_LOWER is not implemented");
1026:   } else if ((flag & SOR_EISENSTAT) && omega == 1.0 && shift == 0 && fshift == 0.0) {
1027:     /* special case for omega = 1.0 saves flops and some integer ops */
1028:     PetscScalar *v2;
1029: 
1030:     v2    = a->a;
1031:     /*  x = (E + U)^{-1} b */
1032:     for (i=m-1; i>=0; i--) {
1033:       n    = a->i[i+1] - diag[i] - 1;
1034:       idx  = a->j + diag[i] + 1;
1035:       v    = a->a + diag[i] + 1;
1036:       sum  = b[i];
1037:       SPARSEDENSEMDOT(sum,xs,v,idx,n);
1038:       x[i] = sum*idiag[i];

1040:       /*  t = b - (2*E - D)x */
1041:       t[i] = b[i] - (v2[diag[i]])*x[i];
1042:     }

1044:     /*  t = (E + L)^{-1}t */
1045:     diag = a->diag;
1046:     for (i=0; i<m; i++) {
1047:       n    = diag[i] - a->i[i];
1048:       idx  = a->j + a->i[i];
1049:       v    = a->a + a->i[i];
1050:       sum  = t[i];
1051:       SPARSEDENSEMDOT(sum,t,v,idx,n);
1052:       t[i]  = sum*idiag[i];

1054:       /*  x = x + t */
1055:       x[i] += t[i];
1056:     }

1058:     PetscLogFlops(3*m-1 + 2*a->nz);
1059:     VecRestoreArray(xx,&x);
1060:     if (bb != xx) {VecRestoreArray(bb,&b);}
1061:     return(0);
1062:   } else if (flag & SOR_EISENSTAT) {
1063:     /* Let  A = L + U + D; where L is lower trianglar,
1064:     U is upper triangular, E is diagonal; This routine applies

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

1068:     to a vector efficiently using Eisenstat's trick. This is for
1069:     the case of SSOR preconditioner, so E is D/omega where omega
1070:     is the relaxation factor.
1071:     */
1072:     scale = (2.0/omega) - 1.0;

1074:     /*  x = (E + U)^{-1} b */
1075:     for (i=m-1; i>=0; i--) {
1076:       d    = fshift + a->a[diag[i] + shift];
1077:       n    = a->i[i+1] - diag[i] - 1;
1078:       idx  = a->j + diag[i] + (!shift);
1079:       v    = a->a + diag[i] + (!shift);
1080:       sum  = b[i];
1081:       SPARSEDENSEMDOT(sum,xs,v,idx,n);
1082:       x[i] = omega*(sum/d);
1083:     }

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

1089:     /*  t = (E + L)^{-1}t */
1090:     ts = t + shift; /* shifted by one for index start of a or a->j*/
1091:     diag = a->diag;
1092:     for (i=0; i<m; i++) {
1093:       d    = fshift + a->a[diag[i]+shift];
1094:       n    = diag[i] - a->i[i];
1095:       idx  = a->j + a->i[i] + shift;
1096:       v    = a->a + a->i[i] + shift;
1097:       sum  = t[i];
1098:       SPARSEDENSEMDOT(sum,ts,v,idx,n);
1099:       t[i] = omega*(sum/d);
1100:       /*  x = x + t */
1101:       x[i] += t[i];
1102:     }

1104:     PetscLogFlops(6*m-1 + 2*a->nz);
1105:     VecRestoreArray(xx,&x);
1106:     if (bb != xx) {VecRestoreArray(bb,&b);}
1107:     return(0);
1108:   }
1109:   if (flag & SOR_ZERO_INITIAL_GUESS) {
1110:     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
1111: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1112:       fortranrelaxaijforwardzero_(&m,&omega,x,a->i,a->j,diag,a->a,b);
1113: #else
1114:       for (i=0; i<m; i++) {
1115:         d    = fshift + a->a[diag[i]+shift];
1116:         n    = diag[i] - a->i[i];
1117:         PetscLogFlops(2*n-1);
1118:         idx  = a->j + a->i[i] + shift;
1119:         v    = a->a + a->i[i] + shift;
1120:         sum  = b[i];
1121:         SPARSEDENSEMDOT(sum,xs,v,idx,n);
1122:         x[i] = omega*(sum/d);
1123:       }
1124: #endif
1125:       xb = x;
1126:     } else xb = b;
1127:     if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
1128:         (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
1129:       for (i=0; i<m; i++) {
1130:         x[i] *= a->a[diag[i]+shift];
1131:       }
1132:       PetscLogFlops(m);
1133:     }
1134:     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
1135: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1136:       fortranrelaxaijbackwardzero_(&m,&omega,x,a->i,a->j,diag,a->a,xb);
1137: #else
1138:       for (i=m-1; i>=0; i--) {
1139:         d    = fshift + a->a[diag[i] + shift];
1140:         n    = a->i[i+1] - diag[i] - 1;
1141:         PetscLogFlops(2*n-1);
1142:         idx  = a->j + diag[i] + (!shift);
1143:         v    = a->a + diag[i] + (!shift);
1144:         sum  = xb[i];
1145:         SPARSEDENSEMDOT(sum,xs,v,idx,n);
1146:         x[i] = omega*(sum/d);
1147:       }
1148: #endif
1149:     }
1150:     its--;
1151:   }
1152:   while (its--) {
1153:     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
1154: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1155:       fortranrelaxaijforward_(&m,&omega,x,a->i,a->j,diag,a->a,b);
1156: #else
1157:       for (i=0; i<m; i++) {
1158:         d    = fshift + a->a[diag[i]+shift];
1159:         n    = a->i[i+1] - a->i[i];
1160:         PetscLogFlops(2*n-1);
1161:         idx  = a->j + a->i[i] + shift;
1162:         v    = a->a + a->i[i] + shift;
1163:         sum  = b[i];
1164:         SPARSEDENSEMDOT(sum,xs,v,idx,n);
1165:         x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d;
1166:       }
1167: #endif
1168:     }
1169:     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
1170: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1171:       fortranrelaxaijbackward_(&m,&omega,x,a->i,a->j,diag,a->a,b);
1172: #else
1173:       for (i=m-1; i>=0; i--) {
1174:         d    = fshift + a->a[diag[i] + shift];
1175:         n    = a->i[i+1] - a->i[i];
1176:         PetscLogFlops(2*n-1);
1177:         idx  = a->j + a->i[i] + shift;
1178:         v    = a->a + a->i[i] + shift;
1179:         sum  = b[i];
1180:         SPARSEDENSEMDOT(sum,xs,v,idx,n);
1181:         x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d;
1182:       }
1183: #endif
1184:     }
1185:   }
1186:   VecRestoreArray(xx,&x);
1187:   if (bb != xx) {VecRestoreArray(bb,&b);}
1188:   return(0);
1189: }

1191: int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,MatInfo *info)
1192: {
1193:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;

1196:   info->rows_global    = (double)A->m;
1197:   info->columns_global = (double)A->n;
1198:   info->rows_local     = (double)A->m;
1199:   info->columns_local  = (double)A->n;
1200:   info->block_size     = 1.0;
1201:   info->nz_allocated   = (double)a->maxnz;
1202:   info->nz_used        = (double)a->nz;
1203:   info->nz_unneeded    = (double)(a->maxnz - a->nz);
1204:   info->assemblies     = (double)A->num_ass;
1205:   info->mallocs        = (double)a->reallocs;
1206:   info->memory         = A->mem;
1207:   if (A->factor) {
1208:     info->fill_ratio_given  = A->info.fill_ratio_given;
1209:     info->fill_ratio_needed = A->info.fill_ratio_needed;
1210:     info->factor_mallocs    = A->info.factor_mallocs;
1211:   } else {
1212:     info->fill_ratio_given  = 0;
1213:     info->fill_ratio_needed = 0;
1214:     info->factor_mallocs    = 0;
1215:   }
1216:   return(0);
1217: }

1219: EXTERN int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,MatLUInfo*,Mat*);
1220: EXTERN int MatLUFactorNumeric_SeqAIJ(Mat,Mat*);
1221: EXTERN int MatLUFactor_SeqAIJ(Mat,IS,IS,MatLUInfo*);
1222: EXTERN int MatSolve_SeqAIJ(Mat,Vec,Vec);
1223: EXTERN int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
1224: EXTERN int MatSolveTranspose_SeqAIJ(Mat,Vec,Vec);
1225: EXTERN int MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec);

1227: int MatZeroRows_SeqAIJ(Mat A,IS is,PetscScalar *diag)
1228: {
1229:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1230:   int         i,ierr,N,*rows,m = A->m - 1,shift = a->indexshift;

1233:   ISGetLocalSize(is,&N);
1234:   ISGetIndices(is,&rows);
1235:   if (a->keepzeroedrows) {
1236:     for (i=0; i<N; i++) {
1237:       if (rows[i] < 0 || rows[i] > m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"row out of range");
1238:       PetscMemzero(&a->a[a->i[rows[i]]+shift],a->ilen[rows[i]]*sizeof(PetscScalar));
1239:     }
1240:     if (diag) {
1241:       MatMissingDiagonal_SeqAIJ(A);
1242:       MatMarkDiagonal_SeqAIJ(A);
1243:       for (i=0; i<N; i++) {
1244:         a->a[a->diag[rows[i]]] = *diag;
1245:       }
1246:     }
1247:   } else {
1248:     if (diag) {
1249:       for (i=0; i<N; i++) {
1250:         if (rows[i] < 0 || rows[i] > m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"row out of range");
1251:         if (a->ilen[rows[i]] > 0) {
1252:           a->ilen[rows[i]]          = 1;
1253:           a->a[a->i[rows[i]]+shift] = *diag;
1254:           a->j[a->i[rows[i]]+shift] = rows[i]+shift;
1255:         } else { /* in case row was completely empty */
1256:           MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES);
1257:         }
1258:       }
1259:     } else {
1260:       for (i=0; i<N; i++) {
1261:         if (rows[i] < 0 || rows[i] > m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"row out of range");
1262:         a->ilen[rows[i]] = 0;
1263:       }
1264:     }
1265:   }
1266:   ISRestoreIndices(is,&rows);
1267:   MatAssemblyEnd_SeqAIJ(A,MAT_FINAL_ASSEMBLY);
1268:   return(0);
1269: }

1271: int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,PetscScalar **v)
1272: {
1273:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1274:   int        *itmp,i,shift = a->indexshift,ierr;

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

1279:   *nz = a->i[row+1] - a->i[row];
1280:   if (v) *v = a->a + a->i[row] + shift;
1281:   if (idx) {
1282:     itmp = a->j + a->i[row] + shift;
1283:     if (*nz && shift) {
1284:       PetscMalloc((*nz)*sizeof(int),idx);
1285:       for (i=0; i<(*nz); i++) {(*idx)[i] = itmp[i] + shift;}
1286:     } else if (*nz) {
1287:       *idx = itmp;
1288:     }
1289:     else *idx = 0;
1290:   }
1291:   return(0);
1292: }

1294: int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,PetscScalar **v)
1295: {
1296:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;

1300:   if (idx) {if (*idx && a->indexshift) {PetscFree(*idx);}}
1301:   return(0);
1302: }

1304: int MatNorm_SeqAIJ(Mat A,NormType type,PetscReal *nrm)
1305: {
1306:   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data;
1307:   PetscScalar  *v = a->a;
1308:   PetscReal    sum = 0.0;
1309:   int          i,j,shift = a->indexshift,ierr;

1312:   if (type == NORM_FROBENIUS) {
1313:     for (i=0; i<a->nz; i++) {
1314: #if defined(PETSC_USE_COMPLEX)
1315:       sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1316: #else
1317:       sum += (*v)*(*v); v++;
1318: #endif
1319:     }
1320:     *nrm = sqrt(sum);
1321:   } else if (type == NORM_1) {
1322:     PetscReal *tmp;
1323:     int    *jj = a->j;
1324:     PetscMalloc((A->n+1)*sizeof(PetscReal),&tmp);
1325:     PetscMemzero(tmp,A->n*sizeof(PetscReal));
1326:     *nrm = 0.0;
1327:     for (j=0; j<a->nz; j++) {
1328:         tmp[*jj++ + shift] += PetscAbsScalar(*v);  v++;
1329:     }
1330:     for (j=0; j<A->n; j++) {
1331:       if (tmp[j] > *nrm) *nrm = tmp[j];
1332:     }
1333:     PetscFree(tmp);
1334:   } else if (type == NORM_INFINITY) {
1335:     *nrm = 0.0;
1336:     for (j=0; j<A->m; j++) {
1337:       v = a->a + a->i[j] + shift;
1338:       sum = 0.0;
1339:       for (i=0; i<a->i[j+1]-a->i[j]; i++) {
1340:         sum += PetscAbsScalar(*v); v++;
1341:       }
1342:       if (sum > *nrm) *nrm = sum;
1343:     }
1344:   } else {
1345:     SETERRQ(PETSC_ERR_SUP,"No support for two norm");
1346:   }
1347:   return(0);
1348: }

1350: int MatTranspose_SeqAIJ(Mat A,Mat *B)
1351: {
1352:   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data;
1353:   Mat          C;
1354:   int          i,ierr,*aj = a->j,*ai = a->i,m = A->m,len,*col;
1355:   int          shift = a->indexshift;
1356:   PetscScalar  *array = a->a;

1359:   if (!B && m != A->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Square matrix only for in-place");
1360:   PetscMalloc((1+A->n)*sizeof(int),&col);
1361:   PetscMemzero(col,(1+A->n)*sizeof(int));
1362:   if (shift) {
1363:     for (i=0; i<ai[m]-1; i++) aj[i] -= 1;
1364:   }
1365:   for (i=0; i<ai[m]+shift; i++) col[aj[i]] += 1;
1366:   MatCreateSeqAIJ(A->comm,A->n,m,0,col,&C);
1367:   PetscFree(col);
1368:   for (i=0; i<m; i++) {
1369:     len    = ai[i+1]-ai[i];
1370:     ierr   = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES);
1371:     array += len;
1372:     aj    += len;
1373:   }
1374:   if (shift) {
1375:     for (i=0; i<ai[m]-1; i++) aj[i] += 1;
1376:   }

1378:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1379:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

1381:   if (B) {
1382:     *B = C;
1383:   } else {
1384:     MatHeaderCopy(A,C);
1385:   }
1386:   return(0);
1387: }

1389: int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr)
1390: {
1391:   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data;
1392:   PetscScalar  *l,*r,x,*v;
1393:   int          ierr,i,j,m = A->m,n = A->n,M,nz = a->nz,*jj,shift = a->indexshift;

1396:   if (ll) {
1397:     /* The local size is used so that VecMPI can be passed to this routine
1398:        by MatDiagonalScale_MPIAIJ */
1399:     VecGetLocalSize(ll,&m);
1400:     if (m != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length");
1401:     VecGetArray(ll,&l);
1402:     v = a->a;
1403:     for (i=0; i<m; i++) {
1404:       x = l[i];
1405:       M = a->i[i+1] - a->i[i];
1406:       for (j=0; j<M; j++) { (*v++) *= x;}
1407:     }
1408:     VecRestoreArray(ll,&l);
1409:     PetscLogFlops(nz);
1410:   }
1411:   if (rr) {
1412:     VecGetLocalSize(rr,&n);
1413:     if (n != A->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vector wrong length");
1414:     VecGetArray(rr,&r);
1415:     v = a->a; jj = a->j;
1416:     for (i=0; i<nz; i++) {
1417:       (*v++) *= r[*jj++ + shift];
1418:     }
1419:     VecRestoreArray(rr,&r);
1420:     PetscLogFlops(nz);
1421:   }
1422:   return(0);
1423: }

1425: int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,int csize,MatReuse scall,Mat *B)
1426: {
1427:   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data,*c;
1428:   int          *smap,i,k,kstart,kend,ierr,oldcols = A->n,*lens;
1429:   int          row,mat_i,*mat_j,tcol,first,step,*mat_ilen,sum,lensi;
1430:   int          *irow,*icol,nrows,ncols,shift = a->indexshift,*ssmap;
1431:   int          *starts,*j_new,*i_new,*aj = a->j,*ai = a->i,ii,*ailen = a->ilen;
1432:   PetscScalar  *a_new,*mat_a;
1433:   Mat          C;
1434:   PetscTruth   stride;

1437:   ISSorted(isrow,(PetscTruth*)&i);
1438:   if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"ISrow is not sorted");
1439:   ISSorted(iscol,(PetscTruth*)&i);
1440:   if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IScol is not sorted");

1442:   ISGetIndices(isrow,&irow);
1443:   ISGetLocalSize(isrow,&nrows);
1444:   ISGetLocalSize(iscol,&ncols);

1446:   ISStrideGetInfo(iscol,&first,&step);
1447:   ISStride(iscol,&stride);
1448:   if (stride && step == 1) {
1449:     /* special case of contiguous rows */
1450:     ierr   = PetscMalloc((2*nrows+1)*sizeof(int),&lens);
1451:     starts = lens + nrows;
1452:     /* loop over new rows determining lens and starting points */
1453:     for (i=0; i<nrows; i++) {
1454:       kstart  = ai[irow[i]]+shift;
1455:       kend    = kstart + ailen[irow[i]];
1456:       for (k=kstart; k<kend; k++) {
1457:         if (aj[k]+shift >= first) {
1458:           starts[i] = k;
1459:           break;
1460:         }
1461:       }
1462:       sum = 0;
1463:       while (k < kend) {
1464:         if (aj[k++]+shift >= first+ncols) break;
1465:         sum++;
1466:       }
1467:       lens[i] = sum;
1468:     }
1469:     /* create submatrix */
1470:     if (scall == MAT_REUSE_MATRIX) {
1471:       int n_cols,n_rows;
1472:       MatGetSize(*B,&n_rows,&n_cols);
1473:       if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size");
1474:       MatZeroEntries(*B);
1475:       C = *B;
1476:     } else {
1477:       MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);
1478:     }
1479:     c = (Mat_SeqAIJ*)C->data;

1481:     /* loop over rows inserting into submatrix */
1482:     a_new    = c->a;
1483:     j_new    = c->j;
1484:     i_new    = c->i;
1485:     i_new[0] = -shift;
1486:     for (i=0; i<nrows; i++) {
1487:       ii    = starts[i];
1488:       lensi = lens[i];
1489:       for (k=0; k<lensi; k++) {
1490:         *j_new++ = aj[ii+k] - first;
1491:       }
1492:       PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(PetscScalar));
1493:       a_new      += lensi;
1494:       i_new[i+1]  = i_new[i] + lensi;
1495:       c->ilen[i]  = lensi;
1496:     }
1497:     PetscFree(lens);
1498:   } else {
1499:     ierr  = ISGetIndices(iscol,&icol);
1500:     ierr  = PetscMalloc((1+oldcols)*sizeof(int),&smap);
1501:     ssmap = smap + shift;
1502:     ierr  = PetscMalloc((1+nrows)*sizeof(int),&lens);
1503:     ierr  = PetscMemzero(smap,oldcols*sizeof(int));
1504:     for (i=0; i<ncols; i++) smap[icol[i]] = i+1;
1505:     /* determine lens of each row */
1506:     for (i=0; i<nrows; i++) {
1507:       kstart  = ai[irow[i]]+shift;
1508:       kend    = kstart + a->ilen[irow[i]];
1509:       lens[i] = 0;
1510:       for (k=kstart; k<kend; k++) {
1511:         if (ssmap[aj[k]]) {
1512:           lens[i]++;
1513:         }
1514:       }
1515:     }
1516:     /* Create and fill new matrix */
1517:     if (scall == MAT_REUSE_MATRIX) {
1518:       PetscTruth equal;

1520:       c = (Mat_SeqAIJ *)((*B)->data);
1521:       if ((*B)->m  != nrows || (*B)->n != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size");
1522:       PetscMemcmp(c->ilen,lens,(*B)->m*sizeof(int),&equal);
1523:       if (!equal) {
1524:         SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros");
1525:       }
1526:       PetscMemzero(c->ilen,(*B)->m*sizeof(int));
1527:       C = *B;
1528:     } else {
1529:       MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);
1530:     }
1531:     c = (Mat_SeqAIJ *)(C->data);
1532:     for (i=0; i<nrows; i++) {
1533:       row    = irow[i];
1534:       kstart = ai[row]+shift;
1535:       kend   = kstart + a->ilen[row];
1536:       mat_i  = c->i[i]+shift;
1537:       mat_j  = c->j + mat_i;
1538:       mat_a  = c->a + mat_i;
1539:       mat_ilen = c->ilen + i;
1540:       for (k=kstart; k<kend; k++) {
1541:         if ((tcol=ssmap[a->j[k]])) {
1542:           *mat_j++ = tcol - (!shift);
1543:           *mat_a++ = a->a[k];
1544:           (*mat_ilen)++;

1546:         }
1547:       }
1548:     }
1549:     /* Free work space */
1550:     ISRestoreIndices(iscol,&icol);
1551:     PetscFree(smap);
1552:     PetscFree(lens);
1553:   }
1554:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1555:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

1557:   ISRestoreIndices(isrow,&irow);
1558:   *B = C;
1559:   return(0);
1560: }

1562: /*
1563: */
1564: int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,MatILUInfo *info)
1565: {
1566:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)inA->data;
1567:   int        ierr;
1568:   Mat        outA;
1569:   PetscTruth row_identity,col_identity;

1572:   if (info && info->levels != 0) SETERRQ(PETSC_ERR_SUP,"Only levels=0 supported for in-place ilu");
1573:   ISIdentity(row,&row_identity);
1574:   ISIdentity(col,&col_identity);
1575:   if (!row_identity || !col_identity) {
1576:     SETERRQ(1,"Row and column permutations must be identity for in-place ILU");
1577:   }

1579:   outA          = inA;
1580:   inA->factor   = FACTOR_LU;
1581:   a->row        = row;
1582:   a->col        = col;
1583:   PetscObjectReference((PetscObject)row);
1584:   PetscObjectReference((PetscObject)col);

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

1591:   if (!a->solve_work) { /* this matrix may have been factored before */
1592:      PetscMalloc((inA->m+1)*sizeof(PetscScalar),&a->solve_work);
1593:   }

1595:   if (!a->diag) {
1596:     MatMarkDiagonal_SeqAIJ(inA);
1597:   }
1598:   MatLUFactorNumeric_SeqAIJ(inA,&outA);
1599:   return(0);
1600: }

1602:  #include petscblaslapack.h
1603: int MatScale_SeqAIJ(PetscScalar *alpha,Mat inA)
1604: {
1605:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)inA->data;
1606:   int        one = 1;

1609:   BLscal_(&a->nz,alpha,a->a,&one);
1610:   PetscLogFlops(a->nz);
1611:   return(0);
1612: }

1614: int MatGetSubMatrices_SeqAIJ(Mat A,int n,IS *irow,IS *icol,MatReuse scall,Mat **B)
1615: {
1616:   int ierr,i;

1619:   if (scall == MAT_INITIAL_MATRIX) {
1620:     PetscMalloc((n+1)*sizeof(Mat),B);
1621:   }

1623:   for (i=0; i<n; i++) {
1624:     MatGetSubMatrix_SeqAIJ(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);
1625:   }
1626:   return(0);
1627: }

1629: int MatGetBlockSize_SeqAIJ(Mat A,int *bs)
1630: {
1632:   *bs = 1;
1633:   return(0);
1634: }

1636: int MatIncreaseOverlap_SeqAIJ(Mat A,int is_max,IS *is,int ov)
1637: {
1638:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1639:   int        shift,row,i,j,k,l,m,n,*idx,ierr,*nidx,isz,val;
1640:   int        start,end,*ai,*aj;
1641:   PetscBT    table;

1644:   shift = a->indexshift;
1645:   m     = A->m;
1646:   ai    = a->i;
1647:   aj    = a->j+shift;

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

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

1654:   for (i=0; i<is_max; i++) {
1655:     /* Initialize the two local arrays */
1656:     isz  = 0;
1657:     PetscBTMemzero(m,table);
1658: 
1659:     /* Extract the indices, assume there can be duplicate entries */
1660:     ISGetIndices(is[i],&idx);
1661:     ISGetLocalSize(is[i],&n);
1662: 
1663:     /* Enter these into the temp arrays. I.e., mark table[row], enter row into new index */
1664:     for (j=0; j<n ; ++j){
1665:       if(!PetscBTLookupSet(table,idx[j])) { nidx[isz++] = idx[j];}
1666:     }
1667:     ISRestoreIndices(is[i],&idx);
1668:     ISDestroy(is[i]);
1669: 
1670:     k = 0;
1671:     for (j=0; j<ov; j++){ /* for each overlap */
1672:       n = isz;
1673:       for (; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */
1674:         row   = nidx[k];
1675:         start = ai[row];
1676:         end   = ai[row+1];
1677:         for (l = start; l<end ; l++){
1678:           val = aj[l] + shift;
1679:           if (!PetscBTLookupSet(table,val)) {nidx[isz++] = val;}
1680:         }
1681:       }
1682:     }
1683:     ISCreateGeneral(PETSC_COMM_SELF,isz,nidx,(is+i));
1684:   }
1685:   PetscBTDestroy(table);
1686:   PetscFree(nidx);
1687:   return(0);
1688: }

1690: /* -------------------------------------------------------------- */
1691: int MatPermute_SeqAIJ(Mat A,IS rowp,IS colp,Mat *B)
1692: {
1693:   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data;
1694:   PetscScalar  *vwork;
1695:   int          i,ierr,nz,m = A->m,n = A->n,*cwork;
1696:   int          *row,*col,*cnew,j,*lens;
1697:   IS           icolp,irowp;

1700:   ISInvertPermutation(rowp,PETSC_DECIDE,&irowp);
1701:   ISGetIndices(irowp,&row);
1702:   ISInvertPermutation(colp,PETSC_DECIDE,&icolp);
1703:   ISGetIndices(icolp,&col);
1704: 
1705:   /* determine lengths of permuted rows */
1706:   PetscMalloc((m+1)*sizeof(int),&lens);
1707:   for (i=0; i<m; i++) {
1708:     lens[row[i]] = a->i[i+1] - a->i[i];
1709:   }
1710:   MatCreateSeqAIJ(A->comm,m,n,0,lens,B);
1711:   PetscFree(lens);

1713:   PetscMalloc(n*sizeof(int),&cnew);
1714:   for (i=0; i<m; i++) {
1715:     MatGetRow(A,i,&nz,&cwork,&vwork);
1716:     for (j=0; j<nz; j++) { cnew[j] = col[cwork[j]];}
1717:     MatSetValues(*B,1,&row[i],nz,cnew,vwork,INSERT_VALUES);
1718:     MatRestoreRow(A,i,&nz,&cwork,&vwork);
1719:   }
1720:   PetscFree(cnew);
1721:   MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);
1722:   MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);
1723:   ISRestoreIndices(irowp,&row);
1724:   ISRestoreIndices(icolp,&col);
1725:   ISDestroy(irowp);
1726:   ISDestroy(icolp);
1727:   return(0);
1728: }

1730: int MatPrintHelp_SeqAIJ(Mat A)
1731: {
1732:   static PetscTruth called = PETSC_FALSE;
1733:   MPI_Comm          comm = A->comm;
1734:   int               ierr;

1737:   if (called) {return(0);} else called = PETSC_TRUE;
1738:   (*PetscHelpPrintf)(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):n");
1739:   (*PetscHelpPrintf)(comm,"  -mat_lu_pivotthreshold <threshold>: Set pivoting thresholdn");
1740:   (*PetscHelpPrintf)(comm,"  -mat_aij_oneindex: internal indices begin at 1 instead of the default 0.n");
1741:   (*PetscHelpPrintf)(comm,"  -mat_aij_no_inode: Do not use inodesn");
1742:   (*PetscHelpPrintf)(comm,"  -mat_aij_inode_limit <limit>: Set inode limit (max limit=5)n");
1743: #if defined(PETSC_HAVE_ESSL)
1744:   (*PetscHelpPrintf)(comm,"  -mat_aij_essl: Use IBM sparse LU factorization and solve.n");
1745: #endif
1746: #if defined(PETSC_HAVE_LUSOL)
1747:   (*PetscHelpPrintf)(comm,"  -mat_aij_lusol: Use the Stanford LUSOL sparse factorization and solve.n");
1748: #endif
1749: #if defined(PETSC_HAVE_MATLAB_ENGINE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
1750:   (*PetscHelpPrintf)(comm,"  -mat_aij_matlab: Use Matlab engine sparse LU factorization and solve.n");
1751: #endif
1752:   return(0);
1753: }
1754: EXTERN int MatEqual_SeqAIJ(Mat A,Mat B,PetscTruth* flg);
1755: EXTERN int MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring);
1756: EXTERN int MatILUDTFactor_SeqAIJ(Mat,MatILUInfo*,IS,IS,Mat*);
1757: int MatCopy_SeqAIJ(Mat A,Mat B,MatStructure str)
1758: {
1759:   int        ierr;
1760:   PetscTruth flg;

1763:   PetscTypeCompare((PetscObject)B,MATSEQAIJ,&flg);
1764:   if (str == SAME_NONZERO_PATTERN && flg) {
1765:     Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1766:     Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data;

1768:     if (a->i[A->m]+a->indexshift != b->i[B->m]+a->indexshift) {
1769:       SETERRQ(1,"Number of nonzeros in two matrices are different");
1770:     }
1771:     PetscMemcpy(b->a,a->a,(a->i[A->m]+a->indexshift)*sizeof(PetscScalar));
1772:   } else {
1773:     MatCopy_Basic(A,B,str);
1774:   }
1775:   return(0);
1776: }

1778: int MatSetUpPreallocation_SeqAIJ(Mat A)
1779: {
1780:   int        ierr;

1783:    MatSeqAIJSetPreallocation(A,PETSC_DEFAULT,0);
1784:   return(0);
1785: }

1787: int MatGetArray_SeqAIJ(Mat A,PetscScalar **array)
1788: {
1789:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1791:   *array = a->a;
1792:   return(0);
1793: }

1795: int MatRestoreArray_SeqAIJ(Mat A,PetscScalar **array)
1796: {
1798:   return(0);
1799: }

1801: int MatFDColoringApply_SeqAIJ(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
1802: {
1803:   int           (*f)(void *,Vec,Vec,void*) = (int (*)(void *,Vec,Vec,void *))coloring->f;
1804:   int           k,ierr,N,start,end,l,row,col,srow,**vscaleforrow,m1,m2;
1805:   PetscScalar   dx,mone = -1.0,*y,*xx,*w3_array;
1806:   PetscScalar   *vscale_array;
1807:   PetscReal     epsilon = coloring->error_rel,umin = coloring->umin;
1808:   Vec           w1,w2,w3;
1809:   void          *fctx = coloring->fctx;
1810:   PetscTruth    flg;

1813:   if (!coloring->w1) {
1814:     VecDuplicate(x1,&coloring->w1);
1815:     PetscLogObjectParent(coloring,coloring->w1);
1816:     VecDuplicate(x1,&coloring->w2);
1817:     PetscLogObjectParent(coloring,coloring->w2);
1818:     VecDuplicate(x1,&coloring->w3);
1819:     PetscLogObjectParent(coloring,coloring->w3);
1820:   }
1821:   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;

1823:   MatSetUnfactored(J);
1824:   PetscOptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);
1825:   if (flg) {
1826:     PetscLogInfo(coloring,"MatFDColoringApply_SeqAIJ: Not calling MatZeroEntries()n");
1827:   } else {
1828:     MatZeroEntries(J);
1829:   }

1831:   VecGetOwnershipRange(x1,&start,&end);
1832:   VecGetSize(x1,&N);

1834:   /*
1835:        This is a horrible, horrible, hack. See DMMGComputeJacobian_Multigrid() it inproperly sets
1836:      coloring->F for the coarser grids from the finest
1837:   */
1838:   if (coloring->F) {
1839:     VecGetLocalSize(coloring->F,&m1);
1840:     VecGetLocalSize(w1,&m2);
1841:     if (m1 != m2) {
1842:       coloring->F = 0;
1843:     }
1844:   }

1846:   if (coloring->F) {
1847:     w1          = coloring->F;
1848:     coloring->F = 0;
1849:   } else {
1850:     (*f)(sctx,x1,w1,fctx);
1851:   }

1853:   /* 
1854:       Compute all the scale factors and share with other processors
1855:   */
1856:   VecGetArray(x1,&xx);xx = xx - start;
1857:   VecGetArray(coloring->vscale,&vscale_array);vscale_array = vscale_array - start;
1858:   for (k=0; k<coloring->ncolors; k++) {
1859:     /*
1860:        Loop over each column associated with color adding the 
1861:        perturbation to the vector w3.
1862:     */
1863:     for (l=0; l<coloring->ncolumns[k]; l++) {
1864:       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
1865:       dx  = xx[col];
1866:       if (dx == 0.0) dx = 1.0;
1867: #if !defined(PETSC_USE_COMPLEX)
1868:       if (dx < umin && dx >= 0.0)      dx = umin;
1869:       else if (dx < 0.0 && dx > -umin) dx = -umin;
1870: #else
1871:       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
1872:       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
1873: #endif
1874:       dx                *= epsilon;
1875:       vscale_array[col] = 1.0/dx;
1876:     }
1877:   }
1878:   vscale_array = vscale_array + start;VecRestoreArray(coloring->vscale,&vscale_array);
1879:   VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);
1880:   VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);

1882:   /*  VecView(coloring->vscale,PETSC_VIEWER_STDOUT_WORLD);
1883:       VecView(x1,PETSC_VIEWER_STDOUT_WORLD);*/

1885:   if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
1886:   else                        vscaleforrow = coloring->columnsforrow;

1888:   VecGetArray(coloring->vscale,&vscale_array);
1889:   /*
1890:       Loop over each color
1891:   */
1892:   for (k=0; k<coloring->ncolors; k++) {
1893:     VecCopy(x1,w3);
1894:     VecGetArray(w3,&w3_array);w3_array = w3_array - start;
1895:     /*
1896:        Loop over each column associated with color adding the 
1897:        perturbation to the vector w3.
1898:     */
1899:     for (l=0; l<coloring->ncolumns[k]; l++) {
1900:       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
1901:       dx  = xx[col];
1902:       if (dx == 0.0) dx = 1.0;
1903: #if !defined(PETSC_USE_COMPLEX)
1904:       if (dx < umin && dx >= 0.0)      dx = umin;
1905:       else if (dx < 0.0 && dx > -umin) dx = -umin;
1906: #else
1907:       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
1908:       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
1909: #endif
1910:       dx            *= epsilon;
1911:       if (!PetscAbsScalar(dx)) SETERRQ(1,"Computed 0 differencing parameter");
1912:       w3_array[col] += dx;
1913:     }
1914:     w3_array = w3_array + start; VecRestoreArray(w3,&w3_array);

1916:     /*
1917:        Evaluate function at x1 + dx (here dx is a vector of perturbations)
1918:     */

1920:     (*f)(sctx,w3,w2,fctx);
1921:     VecAXPY(&mone,w1,w2);

1923:     /*
1924:        Loop over rows of vector, putting results into Jacobian matrix
1925:     */
1926:     VecGetArray(w2,&y);
1927:     for (l=0; l<coloring->nrows[k]; l++) {
1928:       row    = coloring->rows[k][l];
1929:       col    = coloring->columnsforrow[k][l];
1930:       y[row] *= vscale_array[vscaleforrow[k][l]];
1931:       srow   = row + start;
1932:       ierr   = MatSetValues_SeqAIJ(J,1,&srow,1,&col,y+row,INSERT_VALUES);
1933:     }
1934:     VecRestoreArray(w2,&y);
1935:   }
1936:   VecRestoreArray(coloring->vscale,&vscale_array);
1937:   xx = xx + start; ierr  = VecRestoreArray(x1,&xx);
1938:   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
1939:   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
1940:   return(0);
1941: }

1943:  #include petscblaslapack.h

1945: int MatAXPY_SeqAIJ(PetscScalar *a,Mat X,Mat Y,MatStructure str)
1946: {
1947:   int        ierr,one;
1948:   Mat_SeqAIJ *x  = (Mat_SeqAIJ *)X->data,*y = (Mat_SeqAIJ *)Y->data;

1951:   if (str == SAME_NONZERO_PATTERN) {
1952:     BLaxpy_(&x->nz,a,x->a,&one,y->a,&one);
1953:   } else {
1954:     MatAXPY_Basic(a,X,Y,str);
1955:   }
1956:   return(0);
1957: }


1960: /* -------------------------------------------------------------------*/
1961: static struct _MatOps MatOps_Values = {MatSetValues_SeqAIJ,
1962:        MatGetRow_SeqAIJ,
1963:        MatRestoreRow_SeqAIJ,
1964:        MatMult_SeqAIJ,
1965:        MatMultAdd_SeqAIJ,
1966:        MatMultTranspose_SeqAIJ,
1967:        MatMultTransposeAdd_SeqAIJ,
1968:        MatSolve_SeqAIJ,
1969:        MatSolveAdd_SeqAIJ,
1970:        MatSolveTranspose_SeqAIJ,
1971:        MatSolveTransposeAdd_SeqAIJ,
1972:        MatLUFactor_SeqAIJ,
1973:        0,
1974:        MatRelax_SeqAIJ,
1975:        MatTranspose_SeqAIJ,
1976:        MatGetInfo_SeqAIJ,
1977:        MatEqual_SeqAIJ,
1978:        MatGetDiagonal_SeqAIJ,
1979:        MatDiagonalScale_SeqAIJ,
1980:        MatNorm_SeqAIJ,
1981:        0,
1982:        MatAssemblyEnd_SeqAIJ,
1983:        MatCompress_SeqAIJ,
1984:        MatSetOption_SeqAIJ,
1985:        MatZeroEntries_SeqAIJ,
1986:        MatZeroRows_SeqAIJ,
1987:        MatLUFactorSymbolic_SeqAIJ,
1988:        MatLUFactorNumeric_SeqAIJ,
1989:        0,
1990:        0,
1991:        MatSetUpPreallocation_SeqAIJ,
1992:        MatILUFactorSymbolic_SeqAIJ,
1993:        0,
1994:        MatGetArray_SeqAIJ,
1995:        MatRestoreArray_SeqAIJ,
1996:        MatDuplicate_SeqAIJ,
1997:        0,
1998:        0,
1999:        MatILUFactor_SeqAIJ,
2000:        0,
2001:        MatAXPY_SeqAIJ,
2002:        MatGetSubMatrices_SeqAIJ,
2003:        MatIncreaseOverlap_SeqAIJ,
2004:        MatGetValues_SeqAIJ,
2005:        MatCopy_SeqAIJ,
2006:        MatPrintHelp_SeqAIJ,
2007:        MatScale_SeqAIJ,
2008:        0,
2009:        0,
2010:        MatILUDTFactor_SeqAIJ,
2011:        MatGetBlockSize_SeqAIJ,
2012:        MatGetRowIJ_SeqAIJ,
2013:        MatRestoreRowIJ_SeqAIJ,
2014:        MatGetColumnIJ_SeqAIJ,
2015:        MatRestoreColumnIJ_SeqAIJ,
2016:        MatFDColoringCreate_SeqAIJ,
2017:        0,
2018:        0,
2019:        MatPermute_SeqAIJ,
2020:        0,
2021:        0,
2022:        MatDestroy_SeqAIJ,
2023:        MatView_SeqAIJ,
2024:        MatGetPetscMaps_Petsc,
2025:        0,
2026:        0,
2027:        0,
2028:        0,
2029:        0,
2030:        0,
2031:        0,
2032:        0,
2033:        MatSetColoring_SeqAIJ,
2034:        MatSetValuesAdic_SeqAIJ,
2035:        MatSetValuesAdifor_SeqAIJ,
2036:        MatFDColoringApply_SeqAIJ};

2038: EXTERN int MatUseSuperLU_SeqAIJ(Mat);
2039: EXTERN int MatUseEssl_SeqAIJ(Mat);
2040: EXTERN int MatUseLUSOL_SeqAIJ(Mat);
2041: EXTERN int MatUseMatlab_SeqAIJ(Mat);
2042: EXTERN int MatUseDXML_SeqAIJ(Mat);

2044: EXTERN_C_BEGIN

2046: int MatSeqAIJSetColumnIndices_SeqAIJ(Mat mat,int *indices)
2047: {
2048:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
2049:   int        i,nz,n;


2053:   nz = aij->maxnz;
2054:   n  = mat->n;
2055:   for (i=0; i<nz; i++) {
2056:     aij->j[i] = indices[i];
2057:   }
2058:   aij->nz = nz;
2059:   for (i=0; i<n; i++) {
2060:     aij->ilen[i] = aij->imax[i];
2061:   }

2063:   return(0);
2064: }
2065: EXTERN_C_END

2067: /*@
2068:     MatSeqAIJSetColumnIndices - Set the column indices for all the rows
2069:        in the matrix.

2071:   Input Parameters:
2072: +  mat - the SeqAIJ matrix
2073: -  indices - the column indices

2075:   Level: advanced

2077:   Notes:
2078:     This can be called if you have precomputed the nonzero structure of the 
2079:   matrix and want to provide it to the matrix object to improve the performance
2080:   of the MatSetValues() operation.

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

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

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

2089: @*/
2090: int MatSeqAIJSetColumnIndices(Mat mat,int *indices)
2091: {
2092:   int ierr,(*f)(Mat,int *);

2096:   PetscObjectQueryFunction((PetscObject)mat,"MatSeqAIJSetColumnIndices_C",(void (**)(void))&f);
2097:   if (f) {
2098:     (*f)(mat,indices);
2099:   } else {
2100:     SETERRQ(1,"Wrong type of matrix to set column indices");
2101:   }
2102:   return(0);
2103: }

2105: /* ----------------------------------------------------------------------------------------*/

2107: EXTERN_C_BEGIN
2108: int MatStoreValues_SeqAIJ(Mat mat)
2109: {
2110:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
2111:   int        nz = aij->i[mat->m]+aij->indexshift,ierr;

2114:   if (aij->nonew != 1) {
2115:     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
2116:   }

2118:   /* allocate space for values if not already there */
2119:   if (!aij->saved_values) {
2120:     PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);
2121:   }

2123:   /* copy values over */
2124:   PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));
2125:   return(0);
2126: }
2127: EXTERN_C_END

2129: /*@
2130:     MatStoreValues - Stashes a copy of the matrix values; this allows, for 
2131:        example, reuse of the linear part of a Jacobian, while recomputing the 
2132:        nonlinear portion.

2134:    Collect on Mat

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

2139:   Level: advanced

2141:   Common Usage, with SNESSolve():
2142: $    Create Jacobian matrix
2143: $    Set linear terms into matrix
2144: $    Apply boundary conditions to matrix, at this time matrix must have 
2145: $      final nonzero structure (i.e. setting the nonlinear terms and applying 
2146: $      boundary conditions again will not change the nonzero structure
2147: $    MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS);
2148: $    MatStoreValues(mat);
2149: $    Call SNESSetJacobian() with matrix
2150: $    In your Jacobian routine
2151: $      MatRetrieveValues(mat);
2152: $      Set nonlinear terms in matrix
2153:  
2154:   Common Usage without SNESSolve(), i.e. when you handle nonlinear solve yourself:
2155: $    // build linear portion of Jacobian 
2156: $    MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS);
2157: $    MatStoreValues(mat);
2158: $    loop over nonlinear iterations
2159: $       MatRetrieveValues(mat);
2160: $       // call MatSetValues(mat,...) to set nonliner portion of Jacobian 
2161: $       // call MatAssemblyBegin/End() on matrix
2162: $       Solve linear system with Jacobian
2163: $    endloop 

2165:   Notes:
2166:     Matrix must already be assemblied before calling this routine
2167:     Must set the matrix option MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS); before 
2168:     calling this routine.

2170: .seealso: MatRetrieveValues()

2172: @*/
2173: int MatStoreValues(Mat mat)
2174: {
2175:   int ierr,(*f)(Mat);

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

2182:   PetscObjectQueryFunction((PetscObject)mat,"MatStoreValues_C",(void (**)(void))&f);
2183:   if (f) {
2184:     (*f)(mat);
2185:   } else {
2186:     SETERRQ(1,"Wrong type of matrix to store values");
2187:   }
2188:   return(0);
2189: }

2191: EXTERN_C_BEGIN
2192: int MatRetrieveValues_SeqAIJ(Mat mat)
2193: {
2194:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
2195:   int        nz = aij->i[mat->m]+aij->indexshift,ierr;

2198:   if (aij->nonew != 1) {
2199:     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
2200:   }
2201:   if (!aij->saved_values) {
2202:     SETERRQ(1,"Must call MatStoreValues(A);first");
2203:   }

2205:   /* copy values over */
2206:   PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));
2207:   return(0);
2208: }
2209: EXTERN_C_END

2211: /*@
2212:     MatRetrieveValues - Retrieves the copy of the matrix values; this allows, for 
2213:        example, reuse of the linear part of a Jacobian, while recomputing the 
2214:        nonlinear portion.

2216:    Collect on Mat

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

2221:   Level: advanced

2223: .seealso: MatStoreValues()

2225: @*/
2226: int MatRetrieveValues(Mat mat)
2227: {
2228:   int ierr,(*f)(Mat);

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

2235:   PetscObjectQueryFunction((PetscObject)mat,"MatRetrieveValues_C",(void (**)(void))&f);
2236:   if (f) {
2237:     (*f)(mat);
2238:   } else {
2239:     SETERRQ(1,"Wrong type of matrix to retrieve values");
2240:   }
2241:   return(0);
2242: }

2244: /*
2245:    This allows SeqAIJ matrices to be passed to the matlab engine
2246: */
2247: #if defined(PETSC_HAVE_MATLAB_ENGINE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
2248: #include "engine.h"   /* Matlab include file */
2249: #include "mex.h"      /* Matlab include file */
2250: EXTERN_C_BEGIN
2251: int MatMatlabEnginePut_SeqAIJ(PetscObject obj,void *engine)
2252: {
2253:   int         ierr,i,*ai,*aj;
2254:   Mat         B = (Mat)obj;
2255:   PetscScalar *array;
2256:   mxArray     *mat;
2257:   Mat_SeqAIJ  *aij = (Mat_SeqAIJ*)B->data;

2260:   mat  = mxCreateSparse(B->n,B->m,aij->nz,mxREAL);
2261:   PetscMemcpy(mxGetPr(mat),aij->a,aij->nz*sizeof(PetscScalar));
2262:   /* Matlab stores by column, not row so we pass in the transpose of the matrix */
2263:   PetscMemcpy(mxGetIr(mat),aij->j,aij->nz*sizeof(int));
2264:   PetscMemcpy(mxGetJc(mat),aij->i,(B->m+1)*sizeof(int));

2266:   /* Matlab indices start at 0 for sparse (what a surprise) */
2267:   if (aij->indexshift) {
2268:     for (i=0; i<B->m+1; i++) {
2269:       ai[i]--;
2270:     }
2271:     for (i=0; i<aij->nz; i++) {
2272:       aj[i]--;
2273:     }
2274:   }
2275:   PetscObjectName(obj);
2276:   mxSetName(mat,obj->name);
2277:   engPutArray((Engine *)engine,mat);
2278:   return(0);
2279: }
2280: EXTERN_C_END

2282: EXTERN_C_BEGIN
2283: int MatMatlabEngineGet_SeqAIJ(PetscObject obj,void *engine)
2284: {
2285:   int        ierr,ii;
2286:   Mat        mat = (Mat)obj;
2287:   Mat_SeqAIJ *aij = (Mat_SeqAIJ*)mat->data;
2288:   mxArray    *mmat;

2291:   PetscFree(aij->a);
2292:   aij->indexshift = 0;

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

2296:   aij->nz           = (mxGetJc(mmat))[mat->m];
2297:   ierr              = PetscMalloc(aij->nz*(sizeof(int)+sizeof(PetscScalar))+(mat->m+1)*sizeof(int),&aij->a);
2298:   aij->j            = (int*)(aij->a + aij->nz);
2299:   aij->i            = aij->j + aij->nz;
2300:   aij->singlemalloc = PETSC_TRUE;
2301:   aij->freedata     = PETSC_TRUE;

2303:   PetscMemcpy(aij->a,mxGetPr(mmat),aij->nz*sizeof(PetscScalar));
2304:   /* Matlab stores by column, not row so we pass in the transpose of the matrix */
2305:   PetscMemcpy(aij->j,mxGetIr(mmat),aij->nz*sizeof(int));
2306:   PetscMemcpy(aij->i,mxGetJc(mmat),(mat->m+1)*sizeof(int));

2308:   for (ii=0; ii<mat->m; ii++) {
2309:     aij->ilen[ii] = aij->imax[ii] = aij->i[ii+1] - aij->i[ii];
2310:   }

2312:   MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
2313:   MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);

2315:   return(0);
2316: }
2317: EXTERN_C_END
2318: #endif

2320: /* --------------------------------------------------------------------------------*/

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

2329:    Collective on MPI_Comm

2331:    Input Parameters:
2332: +  comm - MPI communicator, set to PETSC_COMM_SELF
2333: .  m - number of rows
2334: .  n - number of columns
2335: .  nz - number of nonzeros per row (same for all rows)
2336: -  nnz - array containing the number of nonzeros in the various rows 
2337:          (possibly different for each row) or PETSC_NULL

2339:    Output Parameter:
2340: .  A - the matrix 

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

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

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

2358:    Options Database Keys:
2359: +  -mat_aij_no_inode  - Do not use inodes
2360: .  -mat_aij_inode_limit <limit> - Sets inode limit (max limit=5)
2361: -  -mat_aij_oneindex - Internally use indexing starting at 1
2362:         rather than 0.  Note that when calling MatSetValues(),
2363:         the user still MUST index entries starting at 0!

2365:    Level: intermediate

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

2369: @*/
2370: int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz,Mat *A)
2371: {

2375:   MatCreate(comm,m,n,m,n,A);
2376:   MatSetType(*A,MATSEQAIJ);
2377:   MatSeqAIJSetPreallocation(*A,nz,nnz);
2378:   return(0);
2379: }

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

2387:    Collective on MPI_Comm

2389:    Input Parameters:
2390: +  comm - MPI communicator, set to PETSC_COMM_SELF
2391: .  m - number of rows
2392: .  n - number of columns
2393: .  nz - number of nonzeros per row (same for all rows)
2394: -  nnz - array containing the number of nonzeros in the various rows 
2395:          (possibly different for each row) or PETSC_NULL

2397:    Output Parameter:
2398: .  A - the matrix 

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

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

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

2416:    Options Database Keys:
2417: +  -mat_aij_no_inode  - Do not use inodes
2418: .  -mat_aij_inode_limit <limit> - Sets inode limit (max limit=5)
2419: -  -mat_aij_oneindex - Internally use indexing starting at 1
2420:         rather than 0.  Note that when calling MatSetValues(),
2421:         the user still MUST index entries starting at 0!

2423:    Level: intermediate

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

2427: @*/
2428: int MatSeqAIJSetPreallocation(Mat B,int nz,int *nnz)
2429: {
2430:   Mat_SeqAIJ *b;
2431:   int        i,len,ierr;
2432:   PetscTruth flg2;

2435:   PetscTypeCompare((PetscObject)B,MATSEQAIJ,&flg2);
2436:   if (!flg2) return(0);
2437: 
2438:   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
2439:   if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz);
2440:   if (nnz) {
2441:     for (i=0; i<B->m; i++) {
2442:       if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]);
2443:       if (nnz[i] > B->n) SETERRQ3(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be greater than row length: local row %d value %d rowlength %d",i,nnz[i],B->n);
2444:     }
2445:   }

2447:   B->preallocated = PETSC_TRUE;
2448:   b = (Mat_SeqAIJ*)B->data;

2450:   PetscMalloc((B->m+1)*sizeof(int),&b->imax);
2451:   if (!nnz) {
2452:     if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 10;
2453:     else if (nz <= 0)        nz = 1;
2454:     for (i=0; i<B->m; i++) b->imax[i] = nz;
2455:     nz = nz*B->m;
2456:   } else {
2457:     nz = 0;
2458:     for (i=0; i<B->m; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
2459:   }

2461:   /* allocate the matrix space */
2462:   len             = nz*(sizeof(int) + sizeof(PetscScalar)) + (B->m+1)*sizeof(int);
2463:   ierr            = PetscMalloc(len,&b->a);
2464:   b->j            = (int*)(b->a + nz);
2465:   ierr            = PetscMemzero(b->j,nz*sizeof(int));
2466:   b->i            = b->j + nz;
2467:   b->singlemalloc = PETSC_TRUE;
2468:   b->freedata     = PETSC_TRUE;

2470:   b->i[0] = -b->indexshift;
2471:   for (i=1; i<B->m+1; i++) {
2472:     b->i[i] = b->i[i-1] + b->imax[i-1];
2473:   }

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

2480:   b->nz                = 0;
2481:   b->maxnz             = nz;
2482:   B->info.nz_unneeded  = (double)b->maxnz;
2483:   return(0);
2484: }

2486: EXTERN_C_BEGIN
2487: int MatCreate_SeqAIJ(Mat B)
2488: {
2489:   Mat_SeqAIJ *b;
2490:   int        ierr,size;
2491:   PetscTruth flg;

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

2497:   B->m = B->M = PetscMax(B->m,B->M);
2498:   B->n = B->N = PetscMax(B->n,B->N);

2500:   PetscNew(Mat_SeqAIJ,&b);
2501:   B->data             = (void*)b;
2502:   PetscMemzero(b,sizeof(Mat_SeqAIJ));
2503:   PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
2504:   B->factor           = 0;
2505:   B->lupivotthreshold = 1.0;
2506:   B->mapping          = 0;
2507:   PetscOptionsGetReal(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold,PETSC_NULL);
2508:   PetscOptionsHasName(PETSC_NULL,"-pc_ilu_preserve_row_sums",&b->ilu_preserve_row_sums);
2509:   b->row              = 0;
2510:   b->col              = 0;
2511:   b->icol             = 0;
2512:   b->indexshift       = 0;
2513:   b->reallocs         = 0;
2514:   PetscOptionsHasName(PETSC_NULL,"-mat_aij_oneindex",&flg);
2515:   if (flg) b->indexshift = -1;
2516: 
2517:   PetscMapCreateMPI(B->comm,B->m,B->m,&B->rmap);
2518:   PetscMapCreateMPI(B->comm,B->n,B->n,&B->cmap);

2520:   b->sorted            = PETSC_FALSE;
2521:   b->ignorezeroentries = PETSC_FALSE;
2522:   b->roworiented       = PETSC_TRUE;
2523:   b->nonew             = 0;
2524:   b->diag              = 0;
2525:   b->solve_work        = 0;
2526:   b->spptr             = 0;
2527:   b->inode.use         = PETSC_TRUE;
2528:   b->inode.node_count  = 0;
2529:   b->inode.size        = 0;
2530:   b->inode.limit       = 5;
2531:   b->inode.max_limit   = 5;
2532:   b->saved_values      = 0;
2533:   b->idiag             = 0;
2534:   b->ssor              = 0;
2535:   b->keepzeroedrows    = PETSC_FALSE;

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

2539: #if defined(PETSC_HAVE_SUPERLU)
2540:   PetscOptionsHasName(PETSC_NULL,"-mat_aij_superlu",&flg);
2541:   if (flg) { MatUseSuperLU_SeqAIJ(B); }
2542: #endif
2543:   PetscOptionsHasName(PETSC_NULL,"-mat_aij_essl",&flg);
2544:   if (flg) { MatUseEssl_SeqAIJ(B); }
2545:   PetscOptionsHasName(PETSC_NULL,"-mat_aij_lusol",&flg);
2546:   if (flg) { MatUseLUSOL_SeqAIJ(B); }
2547:   PetscOptionsHasName(PETSC_NULL,"-mat_aij_matlab",&flg);
2548:   if (flg) {MatUseMatlab_SeqAIJ(B);}
2549:   PetscOptionsHasName(PETSC_NULL,"-mat_aij_dxml",&flg);
2550:   if (flg) {
2551:     if (!b->indexshift) SETERRQ(PETSC_ERR_LIB,"need -mat_aij_oneindex with -mat_aij_dxml");
2552:     MatUseDXML_SeqAIJ(B);
2553:   }
2554:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqAIJSetColumnIndices_C",
2555:                                      "MatSeqAIJSetColumnIndices_SeqAIJ",
2556:                                      MatSeqAIJSetColumnIndices_SeqAIJ);
2557:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
2558:                                      "MatStoreValues_SeqAIJ",
2559:                                      MatStoreValues_SeqAIJ);
2560:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
2561:                                      "MatRetrieveValues_SeqAIJ",
2562:                                      MatRetrieveValues_SeqAIJ);
2563: #if defined(PETSC_HAVE_MATLAB_ENGINE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
2564:   PetscObjectComposeFunctionDynamic((PetscObject)B,"PetscMatlabEnginePut_C","MatMatlabEnginePut_SeqAIJ",MatMatlabEnginePut_SeqAIJ);
2565:   PetscObjectComposeFunctionDynamic((PetscObject)B,"PetscMatlabEngineGet_C","MatMatlabEngineGet_SeqAIJ",MatMatlabEngineGet_SeqAIJ);
2566: #endif
2567:   return(0);
2568: }
2569: EXTERN_C_END

2571: int MatDuplicate_SeqAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
2572: {
2573:   Mat        C;
2574:   Mat_SeqAIJ *c,*a = (Mat_SeqAIJ*)A->data;
2575:   int        i,len,m = A->m,shift = a->indexshift,ierr;

2578:   *B = 0;
2579:   MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);
2580:   MatSetType(C,MATSEQAIJ);
2581:   c    = (Mat_SeqAIJ*)C->data;

2583:   C->factor         = A->factor;
2584:   c->row            = 0;
2585:   c->col            = 0;
2586:   c->icol           = 0;
2587:   c->indexshift     = shift;
2588:   c->keepzeroedrows = a->keepzeroedrows;
2589:   C->assembled      = PETSC_TRUE;

2591:   C->M          = A->m;
2592:   C->N          = A->n;

2594:   PetscMalloc((m+1)*sizeof(int),&c->imax);
2595:   PetscMalloc((m+1)*sizeof(int),&c->ilen);
2596:   for (i=0; i<m; i++) {
2597:     c->imax[i] = a->imax[i];
2598:     c->ilen[i] = a->ilen[i];
2599:   }

2601:   /* allocate the matrix space */
2602:   c->singlemalloc = PETSC_TRUE;
2603:   len   = (m+1)*sizeof(int)+(a->i[m])*(sizeof(PetscScalar)+sizeof(int));
2604:   ierr  = PetscMalloc(len,&c->a);
2605:   c->j  = (int*)(c->a + a->i[m] + shift);
2606:   c->i  = c->j + a->i[m] + shift;
2607:   PetscMemcpy(c->i,a->i,(m+1)*sizeof(int));
2608:   if (m > 0) {
2609:     PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int));
2610:     if (cpvalues == MAT_COPY_VALUES) {
2611:       PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(PetscScalar));
2612:     } else {
2613:       PetscMemzero(c->a,(a->i[m]+shift)*sizeof(PetscScalar));
2614:     }
2615:   }

2617:   PetscLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqAIJ));
2618:   c->sorted      = a->sorted;
2619:   c->roworiented = a->roworiented;
2620:   c->nonew       = a->nonew;
2621:   c->ilu_preserve_row_sums = a->ilu_preserve_row_sums;
2622:   c->saved_values = 0;
2623:   c->idiag        = 0;
2624:   c->ssor         = 0;
2625:   c->ignorezeroentries = a->ignorezeroentries;
2626:   c->freedata     = PETSC_TRUE;

2628:   if (a->diag) {
2629:     PetscMalloc((m+1)*sizeof(int),&c->diag);
2630:     PetscLogObjectMemory(C,(m+1)*sizeof(int));
2631:     for (i=0; i<m; i++) {
2632:       c->diag[i] = a->diag[i];
2633:     }
2634:   } else c->diag        = 0;
2635:   c->inode.use          = a->inode.use;
2636:   c->inode.limit        = a->inode.limit;
2637:   c->inode.max_limit    = a->inode.max_limit;
2638:   if (a->inode.size){
2639:     ierr                = PetscMalloc((m+1)*sizeof(int),&c->inode.size);
2640:     c->inode.node_count = a->inode.node_count;
2641:     ierr                = PetscMemcpy(c->inode.size,a->inode.size,(m+1)*sizeof(int));
2642:   } else {
2643:     c->inode.size       = 0;
2644:     c->inode.node_count = 0;
2645:   }
2646:   c->nz                 = a->nz;
2647:   c->maxnz              = a->maxnz;
2648:   c->solve_work         = 0;
2649:   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
2650:   C->preallocated       = PETSC_TRUE;

2652:   *B = C;
2653:   PetscFListDuplicate(A->qlist,&C->qlist);
2654:   return(0);
2655: }

2657: EXTERN_C_BEGIN
2658: int MatLoad_SeqAIJ(PetscViewer viewer,MatType type,Mat *A)
2659: {
2660:   Mat_SeqAIJ   *a;
2661:   Mat          B;
2662:   int          i,nz,ierr,fd,header[4],size,*rowlengths = 0,M,N,shift;
2663:   MPI_Comm     comm;
2664: 
2666:   PetscObjectGetComm((PetscObject)viewer,&comm);
2667:   MPI_Comm_size(comm,&size);
2668:   if (size > 1) SETERRQ(PETSC_ERR_ARG_SIZ,"view must have one processor");
2669:   PetscViewerBinaryGetDescriptor(viewer,&fd);
2670:   PetscBinaryRead(fd,header,4,PETSC_INT);
2671:   if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object in file");
2672:   M = header[1]; N = header[2]; nz = header[3];

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

2678:   /* read in row lengths */
2679:   PetscMalloc(M*sizeof(int),&rowlengths);
2680:   PetscBinaryRead(fd,rowlengths,M,PETSC_INT);

2682:   /* create our matrix */
2683:   MatCreateSeqAIJ(comm,M,N,0,rowlengths,A);
2684:   B = *A;
2685:   a = (Mat_SeqAIJ*)B->data;
2686:   shift = a->indexshift;

2688:   /* read in column indices and adjust for Fortran indexing*/
2689:   PetscBinaryRead(fd,a->j,nz,PETSC_INT);
2690:   if (shift) {
2691:     for (i=0; i<nz; i++) {
2692:       a->j[i] += 1;
2693:     }
2694:   }

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

2699:   /* set matrix "i" values */
2700:   a->i[0] = -shift;
2701:   for (i=1; i<= M; i++) {
2702:     a->i[i]      = a->i[i-1] + rowlengths[i-1];
2703:     a->ilen[i-1] = rowlengths[i-1];
2704:   }
2705:   PetscFree(rowlengths);

2707:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
2708:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
2709:   return(0);
2710: }
2711: EXTERN_C_END

2713: int MatEqual_SeqAIJ(Mat A,Mat B,PetscTruth* flg)
2714: {
2715:   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data,*b = (Mat_SeqAIJ *)B->data;
2716:   int        ierr;
2717:   PetscTruth flag;

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

2723:   /* If the  matrix dimensions are not equal,or no of nonzeros or shift */
2724:   if ((A->m != B->m) || (A->n != B->n) ||(a->nz != b->nz)|| (a->indexshift != b->indexshift)) {
2725:     *flg = PETSC_FALSE;
2726:     return(0);
2727:   }
2728: 
2729:   /* if the a->i are the same */
2730:   PetscMemcmp(a->i,b->i,(A->m+1)*sizeof(int),flg);
2731:   if (*flg == PETSC_FALSE) return(0);
2732: 
2733:   /* if a->j are the same */
2734:   PetscMemcmp(a->j,b->j,(a->nz)*sizeof(int),flg);
2735:   if (*flg == PETSC_FALSE) return(0);
2736: 
2737:   /* if a->a are the same */
2738:   PetscMemcmp(a->a,b->a,(a->nz)*sizeof(PetscScalar),flg);

2740:   return(0);
2741: 
2742: }

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

2748:       Coolective on MPI_Comm

2750:    Input Parameters:
2751: +   comm - must be an MPI communicator of size 1
2752: .   m - number of rows
2753: .   n - number of columns
2754: .   i - row indices
2755: .   j - column indices
2756: -   a - matrix values

2758:    Output Parameter:
2759: .   mat - the matrix

2761:    Level: intermediate

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

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

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

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

2773: @*/
2774: int MatCreateSeqAIJWithArrays(MPI_Comm comm,int m,int n,int* i,int*j,PetscScalar *a,Mat *mat)
2775: {
2776:   int        ierr,ii;
2777:   Mat_SeqAIJ *aij;

2780:   MatCreateSeqAIJ(comm,m,n,0,0,mat);
2781:   aij  = (Mat_SeqAIJ*)(*mat)->data;
2782:   PetscFree(aij->a);

2784:   if (i[0] == 1) {
2785:     aij->indexshift = -1;
2786:   } else if (i[0]) {
2787:     SETERRQ(1,"i (row indices) do not start with 0 or 1");
2788:   }
2789:   aij->i = i;
2790:   aij->j = j;
2791:   aij->a = a;
2792:   aij->singlemalloc = PETSC_FALSE;
2793:   aij->nonew        = -1;             /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/
2794:   aij->freedata     = PETSC_FALSE;

2796:   for (ii=0; ii<m; ii++) {
2797:     aij->ilen[ii] = aij->imax[ii] = i[ii+1] - i[ii];
2798: #if defined(PETSC_USE_BOPT_g)
2799:     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]);
2800: #endif    
2801:   }
2802: #if defined(PETSC_USE_BOPT_g)
2803:   for (ii=0; ii<aij->i[m]; ii++) {
2804:     if (j[ii] < -aij->indexshift) SETERRQ2(1,"Negative column index at location = %d index = %d",ii,j[ii]);
2805:     if (j[ii] > n - 1 -aij->indexshift) SETERRQ2(1,"Column index to large at location = %d index = %d",ii,j[ii]);
2806:   }
2807: #endif    

2809:   /* changes indices to start at 0 */
2810:   if (i[0]) {
2811:     aij->indexshift = 0;
2812:     for (ii=0; ii<m; ii++) {
2813:       i[ii]--;
2814:     }
2815:     for (ii=0; ii<i[m]; ii++) {
2816:       j[ii]--;
2817:     }
2818:   }

2820:   MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
2821:   MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
2822:   return(0);
2823: }

2825: int MatSetColoring_SeqAIJ(Mat A,ISColoring coloring)
2826: {
2827:   int        ierr;
2828:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;

2831:   ierr        = ISColoringReference(coloring);
2832:   a->coloring = coloring;
2833:   return(0);
2834: }

2836: #if defined(PETSC_HAVE_ADIC) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
2837: EXTERN_C_BEGIN
2838: #include "adic/ad_utils.h"
2839: EXTERN_C_END

2841: int MatSetValuesAdic_SeqAIJ(Mat A,void *advalues)
2842: {
2843:   Mat_SeqAIJ  *a = (Mat_SeqAIJ*)A->data;
2844:   int         m = A->m,*ii = a->i,*jj = a->j,nz,i,*color,j,nlen;
2845:   PetscScalar *v = a->a,*values;
2846:   char        *cadvalues = (char *)advalues;

2849:   if (!a->coloring) SETERRQ(1,"Coloring not set for matrix");
2850:   nlen  = PetscADGetDerivTypeSize();
2851:   color = a->coloring->colors;
2852:   /* loop over rows */
2853:   for (i=0; i<m; i++) {
2854:     nz = ii[i+1] - ii[i];
2855:     /* loop over columns putting computed value into matrix */
2856:     values = PetscADGetGradArray(cadvalues);
2857:     for (j=0; j<nz; j++) {
2858:       *v++ = values[color[*jj++]];
2859:     }
2860:     cadvalues += nlen; /* jump to next row of derivatives */
2861:   }
2862:   return(0);
2863: }

2865: #else

2867: int MatSetValuesAdic_SeqAIJ(Mat A,void *advalues)
2868: {
2870:   SETERRQ(1,"PETSc installed without ADIC");
2871: }

2873: #endif

2875: int MatSetValuesAdifor_SeqAIJ(Mat A,int nl,void *advalues)
2876: {
2877:   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data;
2878:   int          m = A->m,*ii = a->i,*jj = a->j,nz,i,*color,j;
2879:   PetscScalar  *v = a->a,*values = (PetscScalar *)advalues;

2882:   if (!a->coloring) SETERRQ(1,"Coloring not set for matrix");
2883:   color = a->coloring->colors;
2884:   /* loop over rows */
2885:   for (i=0; i<m; i++) {
2886:     nz = ii[i+1] - ii[i];
2887:     /* loop over columns putting computed value into matrix */
2888:     for (j=0; j<nz; j++) {
2889:       *v++ = values[color[*jj++]];
2890:     }
2891:     values += nl; /* jump to next row of derivatives */
2892:   }
2893:   return(0);
2894: }