Actual source code: aij.c

  1: #define PETSCMAT_DLL

  3: /*
  4:     Defines the basic matrix operations for the AIJ (compressed row)
  5:   matrix storage format.
  6: */

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

 15: PetscErrorCode  MatDiagonalSet_SeqAIJ(Mat Y,Vec D,InsertMode is)
 16: {
 18:   Mat_SeqAIJ     *aij = (Mat_SeqAIJ*) Y->data;
 19:   PetscInt       i,*diag, m = Y->rmap.n;
 20:   PetscScalar    *v,*aa = aij->a;
 21:   PetscTruth     missing;

 24:   if (Y->assembled) {
 25:     MatMissingDiagonal_SeqAIJ(Y,&missing,PETSC_NULL);
 26:     if (!missing) {
 27:       diag = aij->diag;
 28:       VecGetArray(D,&v);
 29:       if (is == INSERT_VALUES) {
 30:         for (i=0; i<m; i++) {
 31:           aa[diag[i]] = v[i];
 32:         }
 33:       } else {
 34:         for (i=0; i<m; i++) {
 35:           aa[diag[i]] += v[i];
 36:         }
 37:       }
 38:       VecRestoreArray(D,&v);
 39:       return(0);
 40:     }
 41:   }
 42:   MatDiagonalSet_Default(Y,D,is);
 43:   return(0);
 44: }

 48: PetscErrorCode MatGetRowIJ_SeqAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *m,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
 49: {
 50:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
 52:   PetscInt       i,ishift;
 53: 
 55:   *m     = A->rmap.n;
 56:   if (!ia) return(0);
 57:   ishift = 0;
 58:   if (symmetric && !A->structurally_symmetric) {
 59:     MatToSymmetricIJ_SeqAIJ(A->rmap.n,a->i,a->j,ishift,oshift,ia,ja);
 60:   } else if (oshift == 1) {
 61:     PetscInt nz = a->i[A->rmap.n];
 62:     /* malloc space and  add 1 to i and j indices */
 63:     PetscMalloc((A->rmap.n+1)*sizeof(PetscInt),ia);
 64:     PetscMalloc((nz+1)*sizeof(PetscInt),ja);
 65:     for (i=0; i<nz; i++) (*ja)[i] = a->j[i] + 1;
 66:     for (i=0; i<A->rmap.n+1; i++) (*ia)[i] = a->i[i] + 1;
 67:   } else {
 68:     *ia = a->i; *ja = a->j;
 69:   }
 70:   return(0);
 71: }

 75: PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *n,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
 76: {
 78: 
 80:   if (!ia) return(0);
 81:   if ((symmetric && !A->structurally_symmetric) || oshift == 1) {
 82:     PetscFree(*ia);
 83:     PetscFree(*ja);
 84:   }
 85:   return(0);
 86: }

 90: PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
 91: {
 92:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
 94:   PetscInt       i,*collengths,*cia,*cja,n = A->cmap.n,m = A->rmap.n;
 95:   PetscInt       nz = a->i[m],row,*jj,mr,col;
 96: 
 98:   *nn = n;
 99:   if (!ia) return(0);
100:   if (symmetric) {
101:     MatToSymmetricIJ_SeqAIJ(A->rmap.n,a->i,a->j,0,oshift,ia,ja);
102:   } else {
103:     PetscMalloc((n+1)*sizeof(PetscInt),&collengths);
104:     PetscMemzero(collengths,n*sizeof(PetscInt));
105:     PetscMalloc((n+1)*sizeof(PetscInt),&cia);
106:     PetscMalloc((nz+1)*sizeof(PetscInt),&cja);
107:     jj = a->j;
108:     for (i=0; i<nz; i++) {
109:       collengths[jj[i]]++;
110:     }
111:     cia[0] = oshift;
112:     for (i=0; i<n; i++) {
113:       cia[i+1] = cia[i] + collengths[i];
114:     }
115:     PetscMemzero(collengths,n*sizeof(PetscInt));
116:     jj   = a->j;
117:     for (row=0; row<m; row++) {
118:       mr = a->i[row+1] - a->i[row];
119:       for (i=0; i<mr; i++) {
120:         col = *jj++;
121:         cja[cia[col] + collengths[col]++ - oshift] = row + oshift;
122:       }
123:     }
124:     PetscFree(collengths);
125:     *ia = cia; *ja = cja;
126:   }
127:   return(0);
128: }

132: PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *n,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
133: {

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

139:   PetscFree(*ia);
140:   PetscFree(*ja);
141: 
142:   return(0);
143: }

147: PetscErrorCode MatSetValuesRow_SeqAIJ(Mat A,PetscInt row,const PetscScalar v[])
148: {
149:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
150:   PetscInt       *ai = a->i;

154:   PetscMemcpy(a->a+ai[row],v,(ai[row+1]-ai[row])*sizeof(PetscScalar));
155:   return(0);
156: }

158: #define CHUNKSIZE   15

162: PetscErrorCode MatSetValues_SeqAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
163: {
164:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
165:   PetscInt       *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N;
166:   PetscInt       *imax = a->imax,*ai = a->i,*ailen = a->ilen;
168:   PetscInt       *aj = a->j,nonew = a->nonew,lastcol = -1;
169:   PetscScalar    *ap,value,*aa = a->a;
170:   PetscTruth     ignorezeroentries = ((a->ignorezeroentries && is == ADD_VALUES) ? PETSC_TRUE:PETSC_FALSE);
171:   PetscTruth     roworiented = a->roworiented;

174:   for (k=0; k<m; k++) { /* loop over added rows */
175:     row  = im[k];
176:     if (row < 0) continue;
177: #if defined(PETSC_USE_DEBUG)  
178:     if (row >= A->rmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap.n-1);
179: #endif
180:     rp   = aj + ai[row]; ap = aa + ai[row];
181:     rmax = imax[row]; nrow = ailen[row];
182:     low  = 0;
183:     high = nrow;
184:     for (l=0; l<n; l++) { /* loop over added columns */
185:       if (in[l] < 0) continue;
186: #if defined(PETSC_USE_DEBUG)  
187:       if (in[l] >= A->cmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->cmap.n-1);
188: #endif
189:       col = in[l];
190:       if (roworiented) {
191:         value = v[l + k*n];
192:       } else {
193:         value = v[k + l*m];
194:       }
195:       if (value == 0.0 && ignorezeroentries && (is == ADD_VALUES)) continue;

197:       if (col <= lastcol) low = 0; else high = nrow;
198:       lastcol = col;
199:       while (high-low > 5) {
200:         t = (low+high)/2;
201:         if (rp[t] > col) high = t;
202:         else             low  = t;
203:       }
204:       for (i=low; i<high; i++) {
205:         if (rp[i] > col) break;
206:         if (rp[i] == col) {
207:           if (is == ADD_VALUES) ap[i] += value;
208:           else                  ap[i] = value;
209:           goto noinsert;
210:         }
211:       }
212:       if (value == 0.0 && ignorezeroentries) goto noinsert;
213:       if (nonew == 1) goto noinsert;
214:       if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero at (%D,%D) in the matrix",row,col);
215:       MatSeqXAIJReallocateAIJ(A,A->rmap.n,1,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar);
216:       N = nrow++ - 1; a->nz++; high++;
217:       /* shift up all the later entries in this row */
218:       for (ii=N; ii>=i; ii--) {
219:         rp[ii+1] = rp[ii];
220:         ap[ii+1] = ap[ii];
221:       }
222:       rp[i] = col;
223:       ap[i] = value;
224:       noinsert:;
225:       low = i + 1;
226:     }
227:     ailen[row] = nrow;
228:   }
229:   A->same_nonzero = PETSC_FALSE;
230:   return(0);
231: }


236: PetscErrorCode MatGetValues_SeqAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],PetscScalar v[])
237: {
238:   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data;
239:   PetscInt     *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
240:   PetscInt     *ai = a->i,*ailen = a->ilen;
241:   PetscScalar  *ap,*aa = a->a,zero = 0.0;

244:   for (k=0; k<m; k++) { /* loop over rows */
245:     row  = im[k];
246:     if (row < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",row);
247:     if (row >= A->rmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap.n-1);
248:     rp   = aj + ai[row]; ap = aa + ai[row];
249:     nrow = ailen[row];
250:     for (l=0; l<n; l++) { /* loop over columns */
251:       if (in[l] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %D",in[l]);
252:       if (in[l] >= A->cmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->cmap.n-1);
253:       col = in[l] ;
254:       high = nrow; low = 0; /* assume unsorted */
255:       while (high-low > 5) {
256:         t = (low+high)/2;
257:         if (rp[t] > col) high = t;
258:         else             low  = t;
259:       }
260:       for (i=low; i<high; i++) {
261:         if (rp[i] > col) break;
262:         if (rp[i] == col) {
263:           *v++ = ap[i];
264:           goto finished;
265:         }
266:       }
267:       *v++ = zero;
268:       finished:;
269:     }
270:   }
271:   return(0);
272: }


277: PetscErrorCode MatView_SeqAIJ_Binary(Mat A,PetscViewer viewer)
278: {
279:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
281:   PetscInt       i,*col_lens;
282:   int            fd;

285:   PetscViewerBinaryGetDescriptor(viewer,&fd);
286:   PetscMalloc((4+A->rmap.n)*sizeof(PetscInt),&col_lens);
287:   col_lens[0] = MAT_FILE_COOKIE;
288:   col_lens[1] = A->rmap.n;
289:   col_lens[2] = A->cmap.n;
290:   col_lens[3] = a->nz;

292:   /* store lengths of each row and write (including header) to file */
293:   for (i=0; i<A->rmap.n; i++) {
294:     col_lens[4+i] = a->i[i+1] - a->i[i];
295:   }
296:   PetscBinaryWrite(fd,col_lens,4+A->rmap.n,PETSC_INT,PETSC_TRUE);
297:   PetscFree(col_lens);

299:   /* store column indices (zero start index) */
300:   PetscBinaryWrite(fd,a->j,a->nz,PETSC_INT,PETSC_FALSE);

302:   /* store nonzero values */
303:   PetscBinaryWrite(fd,a->a,a->nz,PETSC_SCALAR,PETSC_FALSE);
304:   return(0);
305: }

307: EXTERN PetscErrorCode MatSeqAIJFactorInfo_Matlab(Mat,PetscViewer);

311: PetscErrorCode MatView_SeqAIJ_ASCII(Mat A,PetscViewer viewer)
312: {
313:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
314:   PetscErrorCode    ierr;
315:   PetscInt          i,j,m = A->rmap.n,shift=0;
316:   const char        *name;
317:   PetscViewerFormat format;

320:   PetscObjectGetName((PetscObject)A,&name);
321:   PetscViewerGetFormat(viewer,&format);
322:   if (format == PETSC_VIEWER_ASCII_MATLAB) {
323:     PetscInt nofinalvalue = 0;
324:     if ((a->i[m] == a->i[m-1]) || (a->j[a->nz-1] != A->cmap.n-!shift)) {
325:       nofinalvalue = 1;
326:     }
327:     PetscViewerASCIIUseTabs(viewer,PETSC_NO);
328:     PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",m,A->cmap.n);
329:     PetscViewerASCIIPrintf(viewer,"%% Nonzeros = %D \n",a->nz);
330:     PetscViewerASCIIPrintf(viewer,"zzz = zeros(%D,3);\n",a->nz+nofinalvalue);
331:     PetscViewerASCIIPrintf(viewer,"zzz = [\n");

333:     for (i=0; i<m; i++) {
334:       for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
335: #if defined(PETSC_USE_COMPLEX)
336:         PetscViewerASCIIPrintf(viewer,"%D %D  %18.16e + %18.16ei \n",i+1,a->j[j]+!shift,PetscRealPart(a->a[j]),PetscImaginaryPart(a->a[j]));
337: #else
338:         PetscViewerASCIIPrintf(viewer,"%D %D  %18.16e\n",i+1,a->j[j]+!shift,a->a[j]);
339: #endif
340:       }
341:     }
342:     if (nofinalvalue) {
343:       PetscViewerASCIIPrintf(viewer,"%D %D  %18.16e\n",m,A->cmap.n,0.0);
344:     }
345:     PetscViewerASCIIPrintf(viewer,"];\n %s = spconvert(zzz);\n",name);
346:     PetscViewerASCIIUseTabs(viewer,PETSC_YES);
347:   } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) {
348:      return(0);
349:   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
350:     PetscViewerASCIIUseTabs(viewer,PETSC_NO);
351:     for (i=0; i<m; i++) {
352:       PetscViewerASCIIPrintf(viewer,"row %D:",i);
353:       for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
354: #if defined(PETSC_USE_COMPLEX)
355:         if (PetscImaginaryPart(a->a[j]) > 0.0 && PetscRealPart(a->a[j]) != 0.0) {
356:           PetscViewerASCIIPrintf(viewer," (%D, %G + %G i)",a->j[j]+shift,PetscRealPart(a->a[j]),PetscImaginaryPart(a->a[j]));
357:         } else if (PetscImaginaryPart(a->a[j]) < 0.0 && PetscRealPart(a->a[j]) != 0.0) {
358:           PetscViewerASCIIPrintf(viewer," (%D, %G - %G i)",a->j[j]+shift,PetscRealPart(a->a[j]),-PetscImaginaryPart(a->a[j]));
359:         } else if (PetscRealPart(a->a[j]) != 0.0) {
360:           PetscViewerASCIIPrintf(viewer," (%D, %G) ",a->j[j]+shift,PetscRealPart(a->a[j]));
361:         }
362: #else
363:         if (a->a[j] != 0.0) {PetscViewerASCIIPrintf(viewer," (%D, %G) ",a->j[j]+shift,a->a[j]);}
364: #endif
365:       }
366:       PetscViewerASCIIPrintf(viewer,"\n");
367:     }
368:     PetscViewerASCIIUseTabs(viewer,PETSC_YES);
369:   } else if (format == PETSC_VIEWER_ASCII_SYMMODU) {
370:     PetscInt nzd=0,fshift=1,*sptr;
371:     PetscViewerASCIIUseTabs(viewer,PETSC_NO);
372:     PetscMalloc((m+1)*sizeof(PetscInt),&sptr);
373:     for (i=0; i<m; i++) {
374:       sptr[i] = nzd+1;
375:       for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
376:         if (a->j[j] >= i) {
377: #if defined(PETSC_USE_COMPLEX)
378:           if (PetscImaginaryPart(a->a[j]) != 0.0 || PetscRealPart(a->a[j]) != 0.0) nzd++;
379: #else
380:           if (a->a[j] != 0.0) nzd++;
381: #endif
382:         }
383:       }
384:     }
385:     sptr[m] = nzd+1;
386:     PetscViewerASCIIPrintf(viewer," %D %D\n\n",m,nzd);
387:     for (i=0; i<m+1; i+=6) {
388:       if (i+4<m) {PetscViewerASCIIPrintf(viewer," %D %D %D %D %D %D\n",sptr[i],sptr[i+1],sptr[i+2],sptr[i+3],sptr[i+4],sptr[i+5]);}
389:       else if (i+3<m) {PetscViewerASCIIPrintf(viewer," %D %D %D %D %D\n",sptr[i],sptr[i+1],sptr[i+2],sptr[i+3],sptr[i+4]);}
390:       else if (i+2<m) {PetscViewerASCIIPrintf(viewer," %D %D %D %D\n",sptr[i],sptr[i+1],sptr[i+2],sptr[i+3]);}
391:       else if (i+1<m) {PetscViewerASCIIPrintf(viewer," %D %D %D\n",sptr[i],sptr[i+1],sptr[i+2]);}
392:       else if (i<m)   {PetscViewerASCIIPrintf(viewer," %D %D\n",sptr[i],sptr[i+1]);}
393:       else            {PetscViewerASCIIPrintf(viewer," %D\n",sptr[i]);}
394:     }
395:     PetscViewerASCIIPrintf(viewer,"\n");
396:     PetscFree(sptr);
397:     for (i=0; i<m; i++) {
398:       for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
399:         if (a->j[j] >= i) {PetscViewerASCIIPrintf(viewer," %D ",a->j[j]+fshift);}
400:       }
401:       PetscViewerASCIIPrintf(viewer,"\n");
402:     }
403:     PetscViewerASCIIPrintf(viewer,"\n");
404:     for (i=0; i<m; i++) {
405:       for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
406:         if (a->j[j] >= i) {
407: #if defined(PETSC_USE_COMPLEX)
408:           if (PetscImaginaryPart(a->a[j]) != 0.0 || PetscRealPart(a->a[j]) != 0.0) {
409:             PetscViewerASCIIPrintf(viewer," %18.16e %18.16e ",PetscRealPart(a->a[j]),PetscImaginaryPart(a->a[j]));
410:           }
411: #else
412:           if (a->a[j] != 0.0) {PetscViewerASCIIPrintf(viewer," %18.16e ",a->a[j]);}
413: #endif
414:         }
415:       }
416:       PetscViewerASCIIPrintf(viewer,"\n");
417:     }
418:     PetscViewerASCIIUseTabs(viewer,PETSC_YES);
419:   } else if (format == PETSC_VIEWER_ASCII_DENSE) {
420:     PetscInt         cnt = 0,jcnt;
421:     PetscScalar value;

423:     PetscViewerASCIIUseTabs(viewer,PETSC_NO);
424:     for (i=0; i<m; i++) {
425:       jcnt = 0;
426:       for (j=0; j<A->cmap.n; j++) {
427:         if (jcnt < a->i[i+1]-a->i[i] && j == a->j[cnt]) {
428:           value = a->a[cnt++];
429:           jcnt++;
430:         } else {
431:           value = 0.0;
432:         }
433: #if defined(PETSC_USE_COMPLEX)
434:         PetscViewerASCIIPrintf(viewer," %7.5e+%7.5e i ",PetscRealPart(value),PetscImaginaryPart(value));
435: #else
436:         PetscViewerASCIIPrintf(viewer," %7.5e ",value);
437: #endif
438:       }
439:       PetscViewerASCIIPrintf(viewer,"\n");
440:     }
441:     PetscViewerASCIIUseTabs(viewer,PETSC_YES);
442:   } else {
443:     PetscViewerASCIIUseTabs(viewer,PETSC_NO);
444:     for (i=0; i<m; i++) {
445:       PetscViewerASCIIPrintf(viewer,"row %D:",i);
446:       for (j=a->i[i]+shift; j<a->i[i+1]+shift; j++) {
447: #if defined(PETSC_USE_COMPLEX)
448:         if (PetscImaginaryPart(a->a[j]) > 0.0) {
449:           PetscViewerASCIIPrintf(viewer," (%D, %G + %G i)",a->j[j]+shift,PetscRealPart(a->a[j]),PetscImaginaryPart(a->a[j]));
450:         } else if (PetscImaginaryPart(a->a[j]) < 0.0) {
451:           PetscViewerASCIIPrintf(viewer," (%D, %G - %G i)",a->j[j]+shift,PetscRealPart(a->a[j]),-PetscImaginaryPart(a->a[j]));
452:         } else {
453:           PetscViewerASCIIPrintf(viewer," (%D, %G) ",a->j[j]+shift,PetscRealPart(a->a[j]));
454:         }
455: #else
456:         PetscViewerASCIIPrintf(viewer," (%D, %G) ",a->j[j]+shift,a->a[j]);
457: #endif
458:       }
459:       PetscViewerASCIIPrintf(viewer,"\n");
460:     }
461:     PetscViewerASCIIUseTabs(viewer,PETSC_YES);
462:   }
463:   PetscViewerFlush(viewer);
464:   return(0);
465: }

469: PetscErrorCode MatView_SeqAIJ_Draw_Zoom(PetscDraw draw,void *Aa)
470: {
471:   Mat               A = (Mat) Aa;
472:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
473:   PetscErrorCode    ierr;
474:   PetscInt          i,j,m = A->rmap.n,color;
475:   PetscReal         xl,yl,xr,yr,x_l,x_r,y_l,y_r,maxv = 0.0;
476:   PetscViewer       viewer;
477:   PetscViewerFormat format;

480:   PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);
481:   PetscViewerGetFormat(viewer,&format);

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

486:   if (format != PETSC_VIEWER_DRAW_CONTOUR) {
487:     /* Blue for negative, Cyan for zero and  Red for positive */
488:     color = PETSC_DRAW_BLUE;
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]; j<a->i[i+1]; j++) {
492:         x_l = a->j[j] ; 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:     color = PETSC_DRAW_CYAN;
502:     for (i=0; i<m; i++) {
503:       y_l = m - i - 1.0; y_r = y_l + 1.0;
504:       for (j=a->i[i]; j<a->i[i+1]; j++) {
505:         x_l = a->j[j]; x_r = x_l + 1.0;
506:         if (a->a[j] !=  0.) continue;
507:         PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
508:       }
509:     }
510:     color = PETSC_DRAW_RED;
511:     for (i=0; i<m; i++) {
512:       y_l = m - i - 1.0; y_r = y_l + 1.0;
513:       for (j=a->i[i]; j<a->i[i+1]; j++) {
514:         x_l = a->j[j]; x_r = x_l + 1.0;
515: #if defined(PETSC_USE_COMPLEX)
516:         if (PetscRealPart(a->a[j]) <=  0.) continue;
517: #else
518:         if (a->a[j] <=  0.) continue;
519: #endif
520:         PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
521:       }
522:     }
523:   } else {
524:     /* use contour shading to indicate magnitude of values */
525:     /* first determine max of all nonzero values */
526:     PetscInt    nz = a->nz,count;
527:     PetscDraw   popup;
528:     PetscReal scale;

530:     for (i=0; i<nz; i++) {
531:       if (PetscAbsScalar(a->a[i]) > maxv) maxv = PetscAbsScalar(a->a[i]);
532:     }
533:     scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv;
534:     PetscDrawGetPopup(draw,&popup);
535:     if (popup) {PetscDrawScalePopup(popup,0.0,maxv);}
536:     count = 0;
537:     for (i=0; i<m; i++) {
538:       y_l = m - i - 1.0; y_r = y_l + 1.0;
539:       for (j=a->i[i]; j<a->i[i+1]; j++) {
540:         x_l = a->j[j]; x_r = x_l + 1.0;
541:         color = PETSC_DRAW_BASIC_COLORS + (PetscInt)(scale*PetscAbsScalar(a->a[count]));
542:         PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
543:         count++;
544:       }
545:     }
546:   }
547:   return(0);
548: }

552: PetscErrorCode MatView_SeqAIJ_Draw(Mat A,PetscViewer viewer)
553: {
555:   PetscDraw      draw;
556:   PetscReal      xr,yr,xl,yl,h,w;
557:   PetscTruth     isnull;

560:   PetscViewerDrawGetDraw(viewer,0,&draw);
561:   PetscDrawIsNull(draw,&isnull);
562:   if (isnull) return(0);

564:   PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);
565:   xr  = A->cmap.n; yr = A->rmap.n; h = yr/10.0; w = xr/10.0;
566:   xr += w;    yr += h;  xl = -w;     yl = -h;
567:   PetscDrawSetCoordinates(draw,xl,yl,xr,yr);
568:   PetscDrawZoom(draw,MatView_SeqAIJ_Draw_Zoom,A);
569:   PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);
570:   return(0);
571: }

575: PetscErrorCode MatView_SeqAIJ(Mat A,PetscViewer viewer)
576: {
578:   PetscTruth     iascii,isbinary,isdraw;

581:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
582:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
583:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
584:   if (iascii) {
585:     MatView_SeqAIJ_ASCII(A,viewer);
586:   } else if (isbinary) {
587:     MatView_SeqAIJ_Binary(A,viewer);
588:   } else if (isdraw) {
589:     MatView_SeqAIJ_Draw(A,viewer);
590:   } else {
591:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by SeqAIJ matrices",((PetscObject)viewer)->type_name);
592:   }
593:   MatView_Inode(A,viewer);
594:   return(0);
595: }

599: PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode)
600: {
601:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
603:   PetscInt       fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
604:   PetscInt       m = A->rmap.n,*ip,N,*ailen = a->ilen,rmax = 0;
605:   PetscScalar    *aa = a->a,*ap;
606:   PetscReal      ratio=0.6;

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

611:   if (m) rmax = ailen[0]; /* determine row with most nonzeros */
612:   for (i=1; i<m; i++) {
613:     /* move each row back by the amount of empty slots (fshift) before it*/
614:     fshift += imax[i-1] - ailen[i-1];
615:     rmax   = PetscMax(rmax,ailen[i]);
616:     if (fshift) {
617:       ip = aj + ai[i] ;
618:       ap = aa + ai[i] ;
619:       N  = ailen[i];
620:       for (j=0; j<N; j++) {
621:         ip[j-fshift] = ip[j];
622:         ap[j-fshift] = ap[j];
623:       }
624:     }
625:     ai[i] = ai[i-1] + ailen[i-1];
626:   }
627:   if (m) {
628:     fshift += imax[m-1] - ailen[m-1];
629:     ai[m]  = ai[m-1] + ailen[m-1];
630:   }
631:   /* reset ilen and imax for each row */
632:   for (i=0; i<m; i++) {
633:     ailen[i] = imax[i] = ai[i+1] - ai[i];
634:   }
635:   a->nz = ai[m];

637:   MatMarkDiagonal_SeqAIJ(A);
638:   PetscInfo4(A,"Matrix size: %D X %D; storage space: %D unneeded,%D used\n",m,A->cmap.n,fshift,a->nz);
639:   PetscInfo1(A,"Number of mallocs during MatSetValues() is %D\n",a->reallocs);
640:   PetscInfo1(A,"Maximum nonzeros in any row is %D\n",rmax);

642:   a->reallocs          = 0;
643:   A->info.nz_unneeded  = (double)fshift;
644:   a->rmax              = rmax;

646:   /* check for zero rows. If found a large number of zero rows, use CompressedRow functions */
647:   Mat_CheckCompressedRow(A,&a->compressedrow,a->i,m,ratio);
648:   A->same_nonzero = PETSC_TRUE;

650:   MatAssemblyEnd_Inode(A,mode);
651:   return(0);
652: }

656: PetscErrorCode MatRealPart_SeqAIJ(Mat A)
657: {
658:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
659:   PetscInt       i,nz = a->nz;
660:   PetscScalar    *aa = a->a;

663:   for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]);
664:   return(0);
665: }

669: PetscErrorCode MatImaginaryPart_SeqAIJ(Mat A)
670: {
671:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
672:   PetscInt       i,nz = a->nz;
673:   PetscScalar    *aa = a->a;

676:   for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]);
677:   return(0);
678: }

682: PetscErrorCode MatZeroEntries_SeqAIJ(Mat A)
683: {
684:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;

688:   PetscMemzero(a->a,(a->i[A->rmap.n])*sizeof(PetscScalar));
689:   return(0);
690: }

694: PetscErrorCode MatDestroy_SeqAIJ(Mat A)
695: {
696:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;

700: #if defined(PETSC_USE_LOG)
701:   PetscLogObjectState((PetscObject)A,"Rows=%D, Cols=%D, NZ=%D",A->rmap.n,A->cmap.n,a->nz);
702: #endif
703:   MatSeqXAIJFreeAIJ(A,&a->a,&a->j,&a->i);
704:   if (a->row) {
705:     ISDestroy(a->row);
706:   }
707:   if (a->col) {
708:     ISDestroy(a->col);
709:   }
710:   PetscFree(a->diag);
711:   PetscFree2(a->imax,a->ilen);
712:   PetscFree(a->idiag);
713:   PetscFree(a->solve_work);
714:   if (a->icol) {ISDestroy(a->icol);}
715:   PetscFree(a->saved_values);
716:   if (a->coloring) {ISColoringDestroy(a->coloring);}
717:   PetscFree(a->xtoy);
718:   if (a->compressedrow.use){PetscFree(a->compressedrow.i);}

720:   MatDestroy_Inode(A);

722:   PetscFree(a);

724:   PetscObjectChangeTypeName((PetscObject)A,0);
725:   PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSeqAIJSetColumnIndices_C","",PETSC_NULL);
726:   PetscObjectComposeFunctionDynamic((PetscObject)A,"MatStoreValues_C","",PETSC_NULL);
727:   PetscObjectComposeFunctionDynamic((PetscObject)A,"MatRetrieveValues_C","",PETSC_NULL);
728:   PetscObjectComposeFunctionDynamic((PetscObject)A,"MatConvert_seqaij_seqsbaij_C","",PETSC_NULL);
729:   PetscObjectComposeFunctionDynamic((PetscObject)A,"MatConvert_seqaij_seqbaij_C","",PETSC_NULL);
730:   PetscObjectComposeFunctionDynamic((PetscObject)A,"MatConvert_seqaij_seqcsrperm_C","",PETSC_NULL);
731:   PetscObjectComposeFunctionDynamic((PetscObject)A,"MatIsTranspose_C","",PETSC_NULL);
732:   PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSeqAIJSetPreallocation_C","",PETSC_NULL);
733:   PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSeqAIJSetPreallocationCSR_C","",PETSC_NULL);
734:   PetscObjectComposeFunctionDynamic((PetscObject)A,"MatReorderForNonzeroDiagonal_C","",PETSC_NULL);
735:   return(0);
736: }

740: PetscErrorCode MatCompress_SeqAIJ(Mat A)
741: {
743:   return(0);
744: }

748: PetscErrorCode MatSetOption_SeqAIJ(Mat A,MatOption op)
749: {
750:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;

754:   switch (op) {
755:     case MAT_ROW_ORIENTED:
756:       a->roworiented       = PETSC_TRUE;
757:       break;
758:     case MAT_KEEP_ZEROED_ROWS:
759:       a->keepzeroedrows    = PETSC_TRUE;
760:       break;
761:     case MAT_COLUMN_ORIENTED:
762:       a->roworiented       = PETSC_FALSE;
763:       break;
764:     case MAT_COLUMNS_SORTED:
765:       a->sorted            = PETSC_TRUE;
766:       break;
767:     case MAT_COLUMNS_UNSORTED:
768:       a->sorted            = PETSC_FALSE;
769:       break;
770:     case MAT_NO_NEW_NONZERO_LOCATIONS:
771:       a->nonew             = 1;
772:       break;
773:     case MAT_NEW_NONZERO_LOCATION_ERR:
774:       a->nonew             = -1;
775:       break;
776:     case MAT_NEW_NONZERO_ALLOCATION_ERR:
777:       a->nonew             = -2;
778:       break;
779:     case MAT_YES_NEW_NONZERO_LOCATIONS:
780:       a->nonew             = 0;
781:       break;
782:     case MAT_IGNORE_ZERO_ENTRIES:
783:       a->ignorezeroentries = PETSC_TRUE;
784:       break;
785:     case MAT_USE_COMPRESSEDROW:
786:       a->compressedrow.use = PETSC_TRUE;
787:       break;
788:     case MAT_DO_NOT_USE_COMPRESSEDROW:
789:       a->compressedrow.use = PETSC_FALSE;
790:       break;
791:     case MAT_ROWS_SORTED:
792:     case MAT_ROWS_UNSORTED:
793:     case MAT_YES_NEW_DIAGONALS:
794:     case MAT_IGNORE_OFF_PROC_ENTRIES:
795:     case MAT_USE_HASH_TABLE:
796:       PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
797:       break;
798:     case MAT_NO_NEW_DIAGONALS:
799:       SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
800:     default:
801:       break;
802:   }
803:   MatSetOption_Inode(A,op);
804:   return(0);
805: }

809: PetscErrorCode MatGetDiagonal_SeqAIJ(Mat A,Vec v)
810: {
811:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
813:   PetscInt       i,j,n;
814:   PetscScalar    *x,zero = 0.0;

817:   VecSet(v,zero);
818:   VecGetArray(v,&x);
819:   VecGetLocalSize(v,&n);
820:   if (n != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
821:   for (i=0; i<A->rmap.n; i++) {
822:     for (j=a->i[i]; j<a->i[i+1]; j++) {
823:       if (a->j[j] == i) {
824:         x[i] = a->a[j];
825:         break;
826:       }
827:     }
828:   }
829:   VecRestoreArray(v,&x);
830:   return(0);
831: }

835: PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
836: {
837:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
838:   PetscScalar       *x,*y;
839:   PetscErrorCode    ierr;
840:   PetscInt          m = A->rmap.n;
841: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
842:   PetscScalar       *v,alpha;
843:   PetscInt          n,i,*idx,*ii,*ridx=PETSC_NULL;
844:   Mat_CompressedRow cprow = a->compressedrow;
845:   PetscTruth        usecprow = cprow.use;
846: #endif

849:   if (zz != yy) {VecCopy(zz,yy);}
850:   VecGetArray(xx,&x);
851:   VecGetArray(yy,&y);

853: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
854:   fortranmulttransposeaddaij_(&m,x,a->i,a->j,a->a,y);
855: #else
856:   if (usecprow){
857:     m    = cprow.nrows;
858:     ii   = cprow.i;
859:     ridx = cprow.rindex;
860:   } else {
861:     ii = a->i;
862:   }
863:   for (i=0; i<m; i++) {
864:     idx   = a->j + ii[i] ;
865:     v     = a->a + ii[i] ;
866:     n     = ii[i+1] - ii[i];
867:     if (usecprow){
868:       alpha = x[ridx[i]];
869:     } else {
870:       alpha = x[i];
871:     }
872:     while (n-->0) {y[*idx++] += alpha * *v++;}
873:   }
874: #endif
875:   PetscLogFlops(2*a->nz);
876:   VecRestoreArray(xx,&x);
877:   VecRestoreArray(yy,&y);
878:   return(0);
879: }

883: PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec xx,Vec yy)
884: {
885:   PetscScalar    zero = 0.0;

889:   VecSet(yy,zero);
890:   MatMultTransposeAdd_SeqAIJ(A,xx,yy,yy);
891:   return(0);
892: }


897: PetscErrorCode MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
898: {
899:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
900:   PetscScalar    *x,*y,*aa;
902:   PetscInt       m=A->rmap.n,*aj,*ii;
903:   PetscInt       n,i,j,*ridx=PETSC_NULL;
904:   PetscScalar    sum;
905:   PetscTruth     usecprow=a->compressedrow.use;
906: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
907:   PetscInt       jrow;
908: #endif

910: #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
911: #pragma disjoint(*x,*y,*aa)
912: #endif

915:   VecGetArray(xx,&x);
916:   VecGetArray(yy,&y);
917:   aj  = a->j;
918:   aa  = a->a;
919:   ii  = a->i;
920:   if (usecprow){ /* use compressed row format */
921:     m    = a->compressedrow.nrows;
922:     ii   = a->compressedrow.i;
923:     ridx = a->compressedrow.rindex;
924:     for (i=0; i<m; i++){
925:       n   = ii[i+1] - ii[i];
926:       aj  = a->j + ii[i];
927:       aa  = a->a + ii[i];
928:       sum = 0.0;
929:       for (j=0; j<n; j++) sum += (*aa++)*x[*aj++];
930:       y[*ridx++] = sum;
931:     }
932:   } else { /* do not use compressed row format */
933: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
934:     fortranmultaij_(&m,x,ii,aj,aa,y);
935: #else
936:     for (i=0; i<m; i++) {
937:       jrow = ii[i];
938:       n    = ii[i+1] - jrow;
939:       sum  = 0.0;
940:       for (j=0; j<n; j++) {
941:         sum += aa[jrow]*x[aj[jrow]]; jrow++;
942:       }
943:       y[i] = sum;
944:     }
945: #endif
946:   }
947:   PetscLogFlops(2*a->nz - m);
948:   VecRestoreArray(xx,&x);
949:   VecRestoreArray(yy,&y);
950:   return(0);
951: }

955: PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
956: {
957:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
958:   PetscScalar    *x,*y,*z,*aa;
960:   PetscInt       m = A->rmap.n,*aj,*ii;
961: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
962:   PetscInt       n,i,jrow,j,*ridx=PETSC_NULL;
963:   PetscScalar    sum;
964:   PetscTruth     usecprow=a->compressedrow.use;
965: #endif

968:   VecGetArray(xx,&x);
969:   VecGetArray(yy,&y);
970:   if (zz != yy) {
971:     VecGetArray(zz,&z);
972:   } else {
973:     z = y;
974:   }

976:   aj  = a->j;
977:   aa  = a->a;
978:   ii  = a->i;
979: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
980:   fortranmultaddaij_(&m,x,ii,aj,aa,y,z);
981: #else
982:   if (usecprow){ /* use compressed row format */
983:     if (zz != yy){
984:       PetscMemcpy(z,y,m*sizeof(PetscScalar));
985:     }
986:     m    = a->compressedrow.nrows;
987:     ii   = a->compressedrow.i;
988:     ridx = a->compressedrow.rindex;
989:     for (i=0; i<m; i++){
990:       n  = ii[i+1] - ii[i];
991:       aj  = a->j + ii[i];
992:       aa  = a->a + ii[i];
993:       sum = y[*ridx];
994:       for (j=0; j<n; j++) sum += (*aa++)*x[*aj++];
995:       z[*ridx++] = sum;
996:     }
997:   } else { /* do not use compressed row format */
998:     for (i=0; i<m; i++) {
999:       jrow = ii[i];
1000:       n    = ii[i+1] - jrow;
1001:       sum  = y[i];
1002:       for (j=0; j<n; j++) {
1003:         sum += aa[jrow]*x[aj[jrow]]; jrow++;
1004:       }
1005:       z[i] = sum;
1006:     }
1007:   }
1008: #endif
1009:   PetscLogFlops(2*a->nz);
1010:   VecRestoreArray(xx,&x);
1011:   VecRestoreArray(yy,&y);
1012:   if (zz != yy) {
1013:     VecRestoreArray(zz,&z);
1014:   }
1015:   return(0);
1016: }

1018: /*
1019:      Adds diagonal pointers to sparse matrix structure.
1020: */
1023: PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat A)
1024: {
1025:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1027:   PetscInt       i,j,m = A->rmap.n;

1030:   if (!a->diag) {
1031:     PetscMalloc(m*sizeof(PetscInt),&a->diag);
1032:   }
1033:   for (i=0; i<A->rmap.n; i++) {
1034:     a->diag[i] = a->i[i+1];
1035:     for (j=a->i[i]; j<a->i[i+1]; j++) {
1036:       if (a->j[j] == i) {
1037:         a->diag[i] = j;
1038:         break;
1039:       }
1040:     }
1041:   }
1042:   return(0);
1043: }

1045: /*
1046:      Checks for missing diagonals
1047: */
1050: PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat A,PetscTruth *missing,PetscInt *d)
1051: {
1052:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1053:   PetscInt       *diag,*jj = a->j,i;

1056:   *missing = PETSC_FALSE;
1057:   if (A->rmap.n > 0 && !jj) {
1058:     *missing  = PETSC_TRUE;
1059:     if (d) *d = 0;
1060:     PetscInfo(A,"Matrix has no entries therefor is missing diagonal");
1061:   } else {
1062:     diag = a->diag;
1063:     for (i=0; i<A->rmap.n; i++) {
1064:       if (jj[diag[i]] != i) {
1065:         *missing = PETSC_TRUE;
1066:         if (d) *d = i;
1067:         PetscInfo1(A,"Matrix is missing diagonal number %D",i);
1068:       }
1069:     }
1070:   }
1071:   return(0);
1072: }

1076: PetscErrorCode MatRelax_SeqAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
1077: {
1078:   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data;
1079:   PetscScalar        *x,d,*xs,sum,*t,scale,*idiag=0,*mdiag;
1080:   const PetscScalar  *v = a->a, *b, *bs,*xb, *ts;
1081:   PetscErrorCode     ierr;
1082:   PetscInt           n = A->cmap.n,m = A->rmap.n,i;
1083:   const PetscInt     *idx,*diag;

1086:   its = its*lits;
1087:   if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);

1089:   diag = a->diag;
1090:   if (!a->idiag) {
1091:     PetscMalloc(3*m*sizeof(PetscScalar),&a->idiag);
1092:     a->ssor  = a->idiag + m;
1093:     mdiag    = a->ssor + m;
1094:     v        = a->a;

1096:     /* this is wrong when fshift omega changes each iteration */
1097:     if (omega == 1.0 && !fshift) {
1098:       for (i=0; i<m; i++) {
1099:         mdiag[i]    = v[diag[i]];
1100:         if (!PetscAbsScalar(mdiag[i])) SETERRQ1(PETSC_ERR_ARG_INCOMP,"Zero diagonal on row %D",i);
1101:         a->idiag[i] = 1.0/v[diag[i]];
1102:       }
1103:       PetscLogFlops(m);
1104:     } else {
1105:       for (i=0; i<m; i++) {
1106:         mdiag[i]    = v[diag[i]];
1107:         a->idiag[i] = omega/(fshift + v[diag[i]]);
1108:       }
1109:       PetscLogFlops(2*m);
1110:     }
1111:   }
1112:   t     = a->ssor;
1113:   idiag = a->idiag;
1114:   mdiag = a->idiag + 2*m;

1116:   VecGetArray(xx,&x);
1117:   if (xx != bb) {
1118:     VecGetArray(bb,(PetscScalar**)&b);
1119:   } else {
1120:     b = x;
1121:   }

1123:   /* We count flops by assuming the upper triangular and lower triangular parts have the same number of nonzeros */
1124:   xs   = x;
1125:   if (flag == SOR_APPLY_UPPER) {
1126:    /* apply (U + D/omega) to the vector */
1127:     bs = b;
1128:     for (i=0; i<m; i++) {
1129:         d    = fshift + a->a[diag[i]];
1130:         n    = a->i[i+1] - diag[i] - 1;
1131:         idx  = a->j + diag[i] + 1;
1132:         v    = a->a + diag[i] + 1;
1133:         sum  = b[i]*d/omega;
1134:         SPARSEDENSEDOT(sum,bs,v,idx,n);
1135:         x[i] = sum;
1136:     }
1137:     VecRestoreArray(xx,&x);
1138:     if (bb != xx) {VecRestoreArray(bb,(PetscScalar**)&b);}
1139:     PetscLogFlops(a->nz);
1140:     return(0);
1141:   }


1144:     /* Let  A = L + U + D; where L is lower trianglar,
1145:     U is upper triangular, E is diagonal; This routine applies

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

1149:     to a vector efficiently using Eisenstat's trick. This is for
1150:     the case of SSOR preconditioner, so E is D/omega where omega
1151:     is the relaxation factor.
1152:     */

1154:   if (flag == SOR_APPLY_LOWER) {
1155:     SETERRQ(PETSC_ERR_SUP,"SOR_APPLY_LOWER is not implemented");
1156:   } else if (flag & SOR_EISENSTAT) {
1157:     /* Let  A = L + U + D; where L is lower trianglar,
1158:     U is upper triangular, E is diagonal; This routine applies

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

1162:     to a vector efficiently using Eisenstat's trick. This is for
1163:     the case of SSOR preconditioner, so E is D/omega where omega
1164:     is the relaxation factor.
1165:     */
1166:     scale = (2.0/omega) - 1.0;

1168:     /*  x = (E + U)^{-1} b */
1169:     for (i=m-1; i>=0; i--) {
1170:       n    = a->i[i+1] - diag[i] - 1;
1171:       idx  = a->j + diag[i] + 1;
1172:       v    = a->a + diag[i] + 1;
1173:       sum  = b[i];
1174:       SPARSEDENSEMDOT(sum,xs,v,idx,n);
1175:       x[i] = sum*idiag[i];
1176:     }

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

1182:     /*  t = (E + L)^{-1}t */
1183:     ts = t;
1184:     diag = a->diag;
1185:     for (i=0; i<m; i++) {
1186:       n    = diag[i] - a->i[i];
1187:       idx  = a->j + a->i[i];
1188:       v    = a->a + a->i[i];
1189:       sum  = t[i];
1190:       SPARSEDENSEMDOT(sum,ts,v,idx,n);
1191:       t[i] = sum*idiag[i];
1192:       /*  x = x + t */
1193:       x[i] += t[i];
1194:     }

1196:     PetscLogFlops(6*m-1 + 2*a->nz);
1197:     VecRestoreArray(xx,&x);
1198:     if (bb != xx) {VecRestoreArray(bb,(PetscScalar**)&b);}
1199:     return(0);
1200:   }
1201:   if (flag & SOR_ZERO_INITIAL_GUESS) {
1202:     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
1203: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1204:       fortranrelaxaijforwardzero_(&m,&omega,x,a->i,a->j,(PetscInt*)diag,idiag,a->a,(void*)b);
1205: #else
1206:       for (i=0; i<m; i++) {
1207:         n    = diag[i] - a->i[i];
1208:         idx  = a->j + a->i[i];
1209:         v    = a->a + a->i[i];
1210:         sum  = b[i];
1211:         SPARSEDENSEMDOT(sum,xs,v,idx,n);
1212:         x[i] = sum*idiag[i];
1213:       }
1214: #endif
1215:       xb = x;
1216:       PetscLogFlops(a->nz);
1217:     } else xb = b;
1218:     if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
1219:         (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
1220:       for (i=0; i<m; i++) {
1221:         x[i] *= mdiag[i];
1222:       }
1223:       PetscLogFlops(m);
1224:     }
1225:     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
1226: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1227:       fortranrelaxaijbackwardzero_(&m,&omega,x,a->i,a->j,(PetscInt*)diag,idiag,a->a,(void*)xb);
1228: #else
1229:       for (i=m-1; i>=0; i--) {
1230:         n    = a->i[i+1] - diag[i] - 1;
1231:         idx  = a->j + diag[i] + 1;
1232:         v    = a->a + diag[i] + 1;
1233:         sum  = xb[i];
1234:         SPARSEDENSEMDOT(sum,xs,v,idx,n);
1235:         x[i] = sum*idiag[i];
1236:       }
1237: #endif
1238:       PetscLogFlops(a->nz);
1239:     }
1240:     its--;
1241:   }
1242:   while (its--) {
1243:     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
1244: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1245:       fortranrelaxaijforward_(&m,&omega,x,a->i,a->j,(PetscInt*)diag,a->a,(void*)b);
1246: #else
1247:       for (i=0; i<m; i++) {
1248:         n    = a->i[i+1] - a->i[i];
1249:         idx  = a->j + a->i[i];
1250:         v    = a->a + a->i[i];
1251:         sum  = b[i];
1252:         SPARSEDENSEMDOT(sum,xs,v,idx,n);
1253:         x[i] = (1. - omega)*x[i] + (sum + mdiag[i]*x[i])*idiag[i];
1254:       }
1255: #endif 
1256:       PetscLogFlops(a->nz);
1257:     }
1258:     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
1259: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1260:       fortranrelaxaijbackward_(&m,&omega,x,a->i,a->j,(PetscInt*)diag,a->a,(void*)b);
1261: #else
1262:       for (i=m-1; i>=0; i--) {
1263:         n    = a->i[i+1] - a->i[i];
1264:         idx  = a->j + a->i[i];
1265:         v    = a->a + a->i[i];
1266:         sum  = b[i];
1267:         SPARSEDENSEMDOT(sum,xs,v,idx,n);
1268:         x[i] = (1. - omega)*x[i] + (sum + mdiag[i]*x[i])*idiag[i];
1269:       }
1270: #endif
1271:       PetscLogFlops(a->nz);
1272:     }
1273:   }
1274:   VecRestoreArray(xx,&x);
1275:   if (bb != xx) {VecRestoreArray(bb,(PetscScalar**)&b);}
1276:   return(0);
1277: }

1281: PetscErrorCode MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,MatInfo *info)
1282: {
1283:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;

1286:   info->rows_global    = (double)A->rmap.n;
1287:   info->columns_global = (double)A->cmap.n;
1288:   info->rows_local     = (double)A->rmap.n;
1289:   info->columns_local  = (double)A->cmap.n;
1290:   info->block_size     = 1.0;
1291:   info->nz_allocated   = (double)a->maxnz;
1292:   info->nz_used        = (double)a->nz;
1293:   info->nz_unneeded    = (double)(a->maxnz - a->nz);
1294:   info->assemblies     = (double)A->num_ass;
1295:   info->mallocs        = (double)a->reallocs;
1296:   info->memory         = A->mem;
1297:   if (A->factor) {
1298:     info->fill_ratio_given  = A->info.fill_ratio_given;
1299:     info->fill_ratio_needed = A->info.fill_ratio_needed;
1300:     info->factor_mallocs    = A->info.factor_mallocs;
1301:   } else {
1302:     info->fill_ratio_given  = 0;
1303:     info->fill_ratio_needed = 0;
1304:     info->factor_mallocs    = 0;
1305:   }
1306:   return(0);
1307: }

1311: PetscErrorCode MatZeroRows_SeqAIJ(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag)
1312: {
1313:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1314:   PetscInt       i,m = A->rmap.n - 1,d;
1316:   PetscTruth     missing;

1319:   if (a->keepzeroedrows) {
1320:     for (i=0; i<N; i++) {
1321:       if (rows[i] < 0 || rows[i] > m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range", rows[i]);
1322:       PetscMemzero(&a->a[a->i[rows[i]]],a->ilen[rows[i]]*sizeof(PetscScalar));
1323:     }
1324:     if (diag != 0.0) {
1325:       MatMissingDiagonal_SeqAIJ(A,&missing,&d);
1326:       if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry in row %D",d);
1327:       for (i=0; i<N; i++) {
1328:         a->a[a->diag[rows[i]]] = diag;
1329:       }
1330:     }
1331:     A->same_nonzero = PETSC_TRUE;
1332:   } else {
1333:     if (diag != 0.0) {
1334:       for (i=0; i<N; i++) {
1335:         if (rows[i] < 0 || rows[i] > m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range", rows[i]);
1336:         if (a->ilen[rows[i]] > 0) {
1337:           a->ilen[rows[i]]          = 1;
1338:           a->a[a->i[rows[i]]] = diag;
1339:           a->j[a->i[rows[i]]] = rows[i];
1340:         } else { /* in case row was completely empty */
1341:           MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],&diag,INSERT_VALUES);
1342:         }
1343:       }
1344:     } else {
1345:       for (i=0; i<N; i++) {
1346:         if (rows[i] < 0 || rows[i] > m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range", rows[i]);
1347:         a->ilen[rows[i]] = 0;
1348:       }
1349:     }
1350:     A->same_nonzero = PETSC_FALSE;
1351:   }
1352:   MatAssemblyEnd_SeqAIJ(A,MAT_FINAL_ASSEMBLY);
1353:   return(0);
1354: }

1358: PetscErrorCode MatGetRow_SeqAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
1359: {
1360:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1361:   PetscInt   *itmp;

1364:   if (row < 0 || row >= A->rmap.n) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Row %D out of range",row);

1366:   *nz = a->i[row+1] - a->i[row];
1367:   if (v) *v = a->a + a->i[row];
1368:   if (idx) {
1369:     itmp = a->j + a->i[row];
1370:     if (*nz) {
1371:       *idx = itmp;
1372:     }
1373:     else *idx = 0;
1374:   }
1375:   return(0);
1376: }

1378: /* remove this function? */
1381: PetscErrorCode MatRestoreRow_SeqAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
1382: {
1384:   return(0);
1385: }

1389: PetscErrorCode MatNorm_SeqAIJ(Mat A,NormType type,PetscReal *nrm)
1390: {
1391:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1392:   PetscScalar    *v = a->a;
1393:   PetscReal      sum = 0.0;
1395:   PetscInt       i,j;

1398:   if (type == NORM_FROBENIUS) {
1399:     for (i=0; i<a->nz; i++) {
1400: #if defined(PETSC_USE_COMPLEX)
1401:       sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1402: #else
1403:       sum += (*v)*(*v); v++;
1404: #endif
1405:     }
1406:     *nrm = sqrt(sum);
1407:   } else if (type == NORM_1) {
1408:     PetscReal *tmp;
1409:     PetscInt    *jj = a->j;
1410:     PetscMalloc((A->cmap.n+1)*sizeof(PetscReal),&tmp);
1411:     PetscMemzero(tmp,A->cmap.n*sizeof(PetscReal));
1412:     *nrm = 0.0;
1413:     for (j=0; j<a->nz; j++) {
1414:         tmp[*jj++] += PetscAbsScalar(*v);  v++;
1415:     }
1416:     for (j=0; j<A->cmap.n; j++) {
1417:       if (tmp[j] > *nrm) *nrm = tmp[j];
1418:     }
1419:     PetscFree(tmp);
1420:   } else if (type == NORM_INFINITY) {
1421:     *nrm = 0.0;
1422:     for (j=0; j<A->rmap.n; j++) {
1423:       v = a->a + a->i[j];
1424:       sum = 0.0;
1425:       for (i=0; i<a->i[j+1]-a->i[j]; i++) {
1426:         sum += PetscAbsScalar(*v); v++;
1427:       }
1428:       if (sum > *nrm) *nrm = sum;
1429:     }
1430:   } else {
1431:     SETERRQ(PETSC_ERR_SUP,"No support for two norm");
1432:   }
1433:   return(0);
1434: }

1438: PetscErrorCode MatTranspose_SeqAIJ(Mat A,Mat *B)
1439: {
1440:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1441:   Mat            C;
1443:   PetscInt       i,*aj = a->j,*ai = a->i,m = A->rmap.n,len,*col;
1444:   PetscScalar    *array = a->a;

1447:   if (!B && m != A->cmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Square matrix only for in-place");
1448:   PetscMalloc((1+A->cmap.n)*sizeof(PetscInt),&col);
1449:   PetscMemzero(col,(1+A->cmap.n)*sizeof(PetscInt));
1450: 
1451:   for (i=0; i<ai[m]; i++) col[aj[i]] += 1;
1452:   MatCreate(A->comm,&C);
1453:   MatSetSizes(C,A->cmap.n,m,A->cmap.n,m);
1454:   MatSetType(C,A->type_name);
1455:   MatSeqAIJSetPreallocation_SeqAIJ(C,0,col);
1456:   PetscFree(col);
1457:   for (i=0; i<m; i++) {
1458:     len    = ai[i+1]-ai[i];
1459:     MatSetValues_SeqAIJ(C,len,aj,1,&i,array,INSERT_VALUES);
1460:     array += len;
1461:     aj    += len;
1462:   }

1464:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1465:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

1467:   if (B) {
1468:     *B = C;
1469:   } else {
1470:     MatHeaderCopy(A,C);
1471:   }
1472:   return(0);
1473: }

1478: PetscErrorCode  MatIsTranspose_SeqAIJ(Mat A,Mat B,PetscReal tol,PetscTruth *f)
1479: {
1480:   Mat_SeqAIJ     *aij = (Mat_SeqAIJ *) A->data,*bij = (Mat_SeqAIJ*) A->data;
1481:   PetscInt       *adx,*bdx,*aii,*bii,*aptr,*bptr; PetscScalar *va,*vb;
1483:   PetscInt       ma,na,mb,nb, i;

1486:   bij = (Mat_SeqAIJ *) B->data;
1487: 
1488:   MatGetSize(A,&ma,&na);
1489:   MatGetSize(B,&mb,&nb);
1490:   if (ma!=nb || na!=mb){
1491:     *f = PETSC_FALSE;
1492:     return(0);
1493:   }
1494:   aii = aij->i; bii = bij->i;
1495:   adx = aij->j; bdx = bij->j;
1496:   va  = aij->a; vb = bij->a;
1497:   PetscMalloc(ma*sizeof(PetscInt),&aptr);
1498:   PetscMalloc(mb*sizeof(PetscInt),&bptr);
1499:   for (i=0; i<ma; i++) aptr[i] = aii[i];
1500:   for (i=0; i<mb; i++) bptr[i] = bii[i];

1502:   *f = PETSC_TRUE;
1503:   for (i=0; i<ma; i++) {
1504:     while (aptr[i]<aii[i+1]) {
1505:       PetscInt         idc,idr;
1506:       PetscScalar vc,vr;
1507:       /* column/row index/value */
1508:       idc = adx[aptr[i]];
1509:       idr = bdx[bptr[idc]];
1510:       vc  = va[aptr[i]];
1511:       vr  = vb[bptr[idc]];
1512:       if (i!=idr || PetscAbsScalar(vc-vr) > tol) {
1513:         *f = PETSC_FALSE;
1514:         goto done;
1515:       } else {
1516:         aptr[i]++;
1517:         if (B || i!=idc) bptr[idc]++;
1518:       }
1519:     }
1520:   }
1521:  done:
1522:   PetscFree(aptr);
1523:   if (B) {
1524:     PetscFree(bptr);
1525:   }
1526:   return(0);
1527: }

1532: PetscErrorCode MatIsSymmetric_SeqAIJ(Mat A,PetscReal tol,PetscTruth *f)
1533: {
1536:   MatIsTranspose_SeqAIJ(A,A,tol,f);
1537:   return(0);
1538: }

1542: PetscErrorCode MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr)
1543: {
1544:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1545:   PetscScalar    *l,*r,x,*v;
1547:   PetscInt       i,j,m = A->rmap.n,n = A->cmap.n,M,nz = a->nz,*jj;

1550:   if (ll) {
1551:     /* The local size is used so that VecMPI can be passed to this routine
1552:        by MatDiagonalScale_MPIAIJ */
1553:     VecGetLocalSize(ll,&m);
1554:     if (m != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length");
1555:     VecGetArray(ll,&l);
1556:     v = a->a;
1557:     for (i=0; i<m; i++) {
1558:       x = l[i];
1559:       M = a->i[i+1] - a->i[i];
1560:       for (j=0; j<M; j++) { (*v++) *= x;}
1561:     }
1562:     VecRestoreArray(ll,&l);
1563:     PetscLogFlops(nz);
1564:   }
1565:   if (rr) {
1566:     VecGetLocalSize(rr,&n);
1567:     if (n != A->cmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vector wrong length");
1568:     VecGetArray(rr,&r);
1569:     v = a->a; jj = a->j;
1570:     for (i=0; i<nz; i++) {
1571:       (*v++) *= r[*jj++];
1572:     }
1573:     VecRestoreArray(rr,&r);
1574:     PetscLogFlops(nz);
1575:   }
1576:   return(0);
1577: }

1581: PetscErrorCode MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,PetscInt csize,MatReuse scall,Mat *B)
1582: {
1583:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data,*c;
1585:   PetscInt       *smap,i,k,kstart,kend,oldcols = A->cmap.n,*lens;
1586:   PetscInt       row,mat_i,*mat_j,tcol,first,step,*mat_ilen,sum,lensi;
1587:   PetscInt       *irow,*icol,nrows,ncols;
1588:   PetscInt       *starts,*j_new,*i_new,*aj = a->j,*ai = a->i,ii,*ailen = a->ilen;
1589:   PetscScalar    *a_new,*mat_a;
1590:   Mat            C;
1591:   PetscTruth     stride;

1594:   ISSorted(isrow,(PetscTruth*)&i);
1595:   if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"ISrow is not sorted");
1596:   ISSorted(iscol,(PetscTruth*)&i);
1597:   if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IScol is not sorted");

1599:   ISGetIndices(isrow,&irow);
1600:   ISGetLocalSize(isrow,&nrows);
1601:   ISGetLocalSize(iscol,&ncols);

1603:   ISStrideGetInfo(iscol,&first,&step);
1604:   ISStride(iscol,&stride);
1605:   if (stride && step == 1) {
1606:     /* special case of contiguous rows */
1607:     PetscMalloc((2*nrows+1)*sizeof(PetscInt),&lens);
1608:     starts = lens + nrows;
1609:     /* loop over new rows determining lens and starting points */
1610:     for (i=0; i<nrows; i++) {
1611:       kstart  = ai[irow[i]];
1612:       kend    = kstart + ailen[irow[i]];
1613:       for (k=kstart; k<kend; k++) {
1614:         if (aj[k] >= first) {
1615:           starts[i] = k;
1616:           break;
1617:         }
1618:       }
1619:       sum = 0;
1620:       while (k < kend) {
1621:         if (aj[k++] >= first+ncols) break;
1622:         sum++;
1623:       }
1624:       lens[i] = sum;
1625:     }
1626:     /* create submatrix */
1627:     if (scall == MAT_REUSE_MATRIX) {
1628:       PetscInt n_cols,n_rows;
1629:       MatGetSize(*B,&n_rows,&n_cols);
1630:       if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size");
1631:       MatZeroEntries(*B);
1632:       C = *B;
1633:     } else {
1634:       MatCreate(A->comm,&C);
1635:       MatSetSizes(C,nrows,ncols,PETSC_DETERMINE,PETSC_DETERMINE);
1636:       MatSetType(C,A->type_name);
1637:       MatSeqAIJSetPreallocation_SeqAIJ(C,0,lens);
1638:     }
1639:     c = (Mat_SeqAIJ*)C->data;

1641:     /* loop over rows inserting into submatrix */
1642:     a_new    = c->a;
1643:     j_new    = c->j;
1644:     i_new    = c->i;

1646:     for (i=0; i<nrows; i++) {
1647:       ii    = starts[i];
1648:       lensi = lens[i];
1649:       for (k=0; k<lensi; k++) {
1650:         *j_new++ = aj[ii+k] - first;
1651:       }
1652:       PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(PetscScalar));
1653:       a_new      += lensi;
1654:       i_new[i+1]  = i_new[i] + lensi;
1655:       c->ilen[i]  = lensi;
1656:     }
1657:     PetscFree(lens);
1658:   } else {
1659:     ISGetIndices(iscol,&icol);
1660:     PetscMalloc((1+oldcols)*sizeof(PetscInt),&smap);
1661: 
1662:     PetscMalloc((1+nrows)*sizeof(PetscInt),&lens);
1663:     PetscMemzero(smap,oldcols*sizeof(PetscInt));
1664:     for (i=0; i<ncols; i++) smap[icol[i]] = i+1;
1665:     /* determine lens of each row */
1666:     for (i=0; i<nrows; i++) {
1667:       kstart  = ai[irow[i]];
1668:       kend    = kstart + a->ilen[irow[i]];
1669:       lens[i] = 0;
1670:       for (k=kstart; k<kend; k++) {
1671:         if (smap[aj[k]]) {
1672:           lens[i]++;
1673:         }
1674:       }
1675:     }
1676:     /* Create and fill new matrix */
1677:     if (scall == MAT_REUSE_MATRIX) {
1678:       PetscTruth equal;

1680:       c = (Mat_SeqAIJ *)((*B)->data);
1681:       if ((*B)->rmap.n  != nrows || (*B)->cmap.n != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size");
1682:       PetscMemcmp(c->ilen,lens,(*B)->rmap.n*sizeof(PetscInt),&equal);
1683:       if (!equal) {
1684:         SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros");
1685:       }
1686:       PetscMemzero(c->ilen,(*B)->rmap.n*sizeof(PetscInt));
1687:       C = *B;
1688:     } else {
1689:       MatCreate(A->comm,&C);
1690:       MatSetSizes(C,nrows,ncols,PETSC_DETERMINE,PETSC_DETERMINE);
1691:       MatSetType(C,A->type_name);
1692:       MatSeqAIJSetPreallocation_SeqAIJ(C,0,lens);
1693:     }
1694:     c = (Mat_SeqAIJ *)(C->data);
1695:     for (i=0; i<nrows; i++) {
1696:       row    = irow[i];
1697:       kstart = ai[row];
1698:       kend   = kstart + a->ilen[row];
1699:       mat_i  = c->i[i];
1700:       mat_j  = c->j + mat_i;
1701:       mat_a  = c->a + mat_i;
1702:       mat_ilen = c->ilen + i;
1703:       for (k=kstart; k<kend; k++) {
1704:         if ((tcol=smap[a->j[k]])) {
1705:           *mat_j++ = tcol - 1;
1706:           *mat_a++ = a->a[k];
1707:           (*mat_ilen)++;

1709:         }
1710:       }
1711:     }
1712:     /* Free work space */
1713:     ISRestoreIndices(iscol,&icol);
1714:     PetscFree(smap);
1715:     PetscFree(lens);
1716:   }
1717:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1718:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

1720:   ISRestoreIndices(isrow,&irow);
1721:   *B = C;
1722:   return(0);
1723: }

1725: /*
1726: */
1729: PetscErrorCode MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,MatFactorInfo *info)
1730: {
1731:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)inA->data;
1733:   Mat            outA;
1734:   PetscTruth     row_identity,col_identity;

1737:   if (info->levels != 0) SETERRQ(PETSC_ERR_SUP,"Only levels=0 supported for in-place ilu");
1738:   ISIdentity(row,&row_identity);
1739:   ISIdentity(col,&col_identity);

1741:   outA          = inA;
1742:   inA->factor   = FACTOR_LU;
1743:   a->row        = row;
1744:   a->col        = col;
1745:   PetscObjectReference((PetscObject)row);
1746:   PetscObjectReference((PetscObject)col);

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

1753:   if (!a->solve_work) { /* this matrix may have been factored before */
1754:      PetscMalloc((inA->rmap.n+1)*sizeof(PetscScalar),&a->solve_work);
1755:   }

1757:   MatMarkDiagonal_SeqAIJ(inA);
1758:   if (row_identity && col_identity) {
1759:     MatLUFactorNumeric_SeqAIJ(inA,info,&outA);
1760:   } else {
1761:     MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(inA,info,&outA);
1762:   }
1763:   return(0);
1764: }

1766:  #include petscblaslapack.h
1769: PetscErrorCode MatScale_SeqAIJ(Mat inA,PetscScalar alpha)
1770: {
1771:   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)inA->data;
1772:   PetscBLASInt bnz = (PetscBLASInt)a->nz,one = 1;
1773:   PetscScalar oalpha = alpha;


1778:   BLASscal_(&bnz,&oalpha,a->a,&one);
1779:   PetscLogFlops(a->nz);
1780:   return(0);
1781: }

1785: PetscErrorCode MatGetSubMatrices_SeqAIJ(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1786: {
1788:   PetscInt       i;

1791:   if (scall == MAT_INITIAL_MATRIX) {
1792:     PetscMalloc((n+1)*sizeof(Mat),B);
1793:   }

1795:   for (i=0; i<n; i++) {
1796:     MatGetSubMatrix_SeqAIJ(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);
1797:   }
1798:   return(0);
1799: }

1803: PetscErrorCode MatIncreaseOverlap_SeqAIJ(Mat A,PetscInt is_max,IS is[],PetscInt ov)
1804: {
1805:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1807:   PetscInt       row,i,j,k,l,m,n,*idx,*nidx,isz,val;
1808:   PetscInt       start,end,*ai,*aj;
1809:   PetscBT        table;

1812:   m     = A->rmap.n;
1813:   ai    = a->i;
1814:   aj    = a->j;

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

1818:   PetscMalloc((m+1)*sizeof(PetscInt),&nidx);
1819:   PetscBTCreate(m,table);

1821:   for (i=0; i<is_max; i++) {
1822:     /* Initialize the two local arrays */
1823:     isz  = 0;
1824:     PetscBTMemzero(m,table);
1825: 
1826:     /* Extract the indices, assume there can be duplicate entries */
1827:     ISGetIndices(is[i],&idx);
1828:     ISGetLocalSize(is[i],&n);
1829: 
1830:     /* Enter these into the temp arrays. I.e., mark table[row], enter row into new index */
1831:     for (j=0; j<n ; ++j){
1832:       if(!PetscBTLookupSet(table,idx[j])) { nidx[isz++] = idx[j];}
1833:     }
1834:     ISRestoreIndices(is[i],&idx);
1835:     ISDestroy(is[i]);
1836: 
1837:     k = 0;
1838:     for (j=0; j<ov; j++){ /* for each overlap */
1839:       n = isz;
1840:       for (; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */
1841:         row   = nidx[k];
1842:         start = ai[row];
1843:         end   = ai[row+1];
1844:         for (l = start; l<end ; l++){
1845:           val = aj[l] ;
1846:           if (!PetscBTLookupSet(table,val)) {nidx[isz++] = val;}
1847:         }
1848:       }
1849:     }
1850:     ISCreateGeneral(PETSC_COMM_SELF,isz,nidx,(is+i));
1851:   }
1852:   PetscBTDestroy(table);
1853:   PetscFree(nidx);
1854:   return(0);
1855: }

1857: /* -------------------------------------------------------------- */
1860: PetscErrorCode MatPermute_SeqAIJ(Mat A,IS rowp,IS colp,Mat *B)
1861: {
1862:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1864:   PetscInt       i,nz,m = A->rmap.n,n = A->cmap.n,*col;
1865:   PetscInt       *row,*cnew,j,*lens;
1866:   IS             icolp,irowp;
1867:   PetscInt       *cwork;
1868:   PetscScalar    *vwork;

1871:   ISInvertPermutation(rowp,PETSC_DECIDE,&irowp);
1872:   ISGetIndices(irowp,&row);
1873:   ISInvertPermutation(colp,PETSC_DECIDE,&icolp);
1874:   ISGetIndices(icolp,&col);
1875: 
1876:   /* determine lengths of permuted rows */
1877:   PetscMalloc((m+1)*sizeof(PetscInt),&lens);
1878:   for (i=0; i<m; i++) {
1879:     lens[row[i]] = a->i[i+1] - a->i[i];
1880:   }
1881:   MatCreate(A->comm,B);
1882:   MatSetSizes(*B,m,n,m,n);
1883:   MatSetType(*B,A->type_name);
1884:   MatSeqAIJSetPreallocation_SeqAIJ(*B,0,lens);
1885:   PetscFree(lens);

1887:   PetscMalloc(n*sizeof(PetscInt),&cnew);
1888:   for (i=0; i<m; i++) {
1889:     MatGetRow_SeqAIJ(A,i,&nz,&cwork,&vwork);
1890:     for (j=0; j<nz; j++) { cnew[j] = col[cwork[j]];}
1891:     MatSetValues_SeqAIJ(*B,1,&row[i],nz,cnew,vwork,INSERT_VALUES);
1892:     MatRestoreRow_SeqAIJ(A,i,&nz,&cwork,&vwork);
1893:   }
1894:   PetscFree(cnew);
1895:   (*B)->assembled     = PETSC_FALSE;
1896:   MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);
1897:   MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);
1898:   ISRestoreIndices(irowp,&row);
1899:   ISRestoreIndices(icolp,&col);
1900:   ISDestroy(irowp);
1901:   ISDestroy(icolp);
1902:   return(0);
1903: }

1907: PetscErrorCode MatCopy_SeqAIJ(Mat A,Mat B,MatStructure str)
1908: {

1912:   /* If the two matrices have the same copy implementation, use fast copy. */
1913:   if (str == SAME_NONZERO_PATTERN && (A->ops->copy == B->ops->copy)) {
1914:     Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1915:     Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data;

1917:     if (a->i[A->rmap.n] != b->i[B->rmap.n]) {
1918:       SETERRQ(PETSC_ERR_ARG_INCOMP,"Number of nonzeros in two matrices are different");
1919:     }
1920:     PetscMemcpy(b->a,a->a,(a->i[A->rmap.n])*sizeof(PetscScalar));
1921:   } else {
1922:     MatCopy_Basic(A,B,str);
1923:   }
1924:   return(0);
1925: }

1929: PetscErrorCode MatSetUpPreallocation_SeqAIJ(Mat A)
1930: {

1934:    MatSeqAIJSetPreallocation_SeqAIJ(A,PETSC_DEFAULT,0);
1935:   return(0);
1936: }

1940: PetscErrorCode MatGetArray_SeqAIJ(Mat A,PetscScalar *array[])
1941: {
1942:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1944:   *array = a->a;
1945:   return(0);
1946: }

1950: PetscErrorCode MatRestoreArray_SeqAIJ(Mat A,PetscScalar *array[])
1951: {
1953:   return(0);
1954: }

1958: PetscErrorCode MatFDColoringApply_SeqAIJ(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
1959: {
1960:   PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void *))coloring->f;
1962:   PetscInt       k,N,start,end,l,row,col,srow,**vscaleforrow,m1,m2;
1963:   PetscScalar    dx,*y,*xx,*w3_array;
1964:   PetscScalar    *vscale_array;
1965:   PetscReal      epsilon = coloring->error_rel,umin = coloring->umin;
1966:   Vec            w1,w2,w3;
1967:   void           *fctx = coloring->fctx;
1968:   PetscTruth     flg;

1971:   if (!coloring->w1) {
1972:     VecDuplicate(x1,&coloring->w1);
1973:     PetscLogObjectParent(coloring,coloring->w1);
1974:     VecDuplicate(x1,&coloring->w2);
1975:     PetscLogObjectParent(coloring,coloring->w2);
1976:     VecDuplicate(x1,&coloring->w3);
1977:     PetscLogObjectParent(coloring,coloring->w3);
1978:   }
1979:   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;

1981:   MatSetUnfactored(J);
1982:   PetscOptionsHasName(coloring->prefix,"-mat_fd_coloring_dont_rezero",&flg);
1983:   if (flg) {
1984:     PetscInfo(coloring,"Not calling MatZeroEntries()\n");
1985:   } else {
1986:     PetscTruth assembled;
1987:     MatAssembled(J,&assembled);
1988:     if (assembled) {
1989:       MatZeroEntries(J);
1990:     }
1991:   }

1993:   VecGetOwnershipRange(x1,&start,&end);
1994:   VecGetSize(x1,&N);

1996:   /*
1997:        This is a horrible, horrible, hack. See DMMGComputeJacobian_Multigrid() it inproperly sets
1998:      coloring->F for the coarser grids from the finest
1999:   */
2000:   if (coloring->F) {
2001:     VecGetLocalSize(coloring->F,&m1);
2002:     VecGetLocalSize(w1,&m2);
2003:     if (m1 != m2) {
2004:       coloring->F = 0;
2005:     }
2006:   }

2008:   if (coloring->F) {
2009:     w1          = coloring->F;
2010:     coloring->F = 0;
2011:   } else {
2013:     (*f)(sctx,x1,w1,fctx);
2015:   }

2017:   /* 
2018:       Compute all the scale factors and share with other processors
2019:   */
2020:   VecGetArray(x1,&xx);xx = xx - start;
2021:   VecGetArray(coloring->vscale,&vscale_array);vscale_array = vscale_array - start;
2022:   for (k=0; k<coloring->ncolors; k++) {
2023:     /*
2024:        Loop over each column associated with color adding the 
2025:        perturbation to the vector w3.
2026:     */
2027:     for (l=0; l<coloring->ncolumns[k]; l++) {
2028:       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
2029:       dx  = xx[col];
2030:       if (dx == 0.0) dx = 1.0;
2031: #if !defined(PETSC_USE_COMPLEX)
2032:       if (dx < umin && dx >= 0.0)      dx = umin;
2033:       else if (dx < 0.0 && dx > -umin) dx = -umin;
2034: #else
2035:       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
2036:       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
2037: #endif
2038:       dx                *= epsilon;
2039:       vscale_array[col] = 1.0/dx;
2040:     }
2041:   }
2042:   vscale_array = vscale_array + start;VecRestoreArray(coloring->vscale,&vscale_array);
2043:   VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);
2044:   VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);

2046:   /*  VecView(coloring->vscale,PETSC_VIEWER_STDOUT_WORLD);
2047:       VecView(x1,PETSC_VIEWER_STDOUT_WORLD);*/

2049:   if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
2050:   else                        vscaleforrow = coloring->columnsforrow;

2052:   VecGetArray(coloring->vscale,&vscale_array);
2053:   /*
2054:       Loop over each color
2055:   */
2056:   for (k=0; k<coloring->ncolors; k++) {
2057:     coloring->currentcolor = k;
2058:     VecCopy(x1,w3);
2059:     VecGetArray(w3,&w3_array);w3_array = w3_array - start;
2060:     /*
2061:        Loop over each column associated with color adding the 
2062:        perturbation to the vector w3.
2063:     */
2064:     for (l=0; l<coloring->ncolumns[k]; l++) {
2065:       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
2066:       dx  = xx[col];
2067:       if (dx == 0.0) dx = 1.0;
2068: #if !defined(PETSC_USE_COMPLEX)
2069:       if (dx < umin && dx >= 0.0)      dx = umin;
2070:       else if (dx < 0.0 && dx > -umin) dx = -umin;
2071: #else
2072:       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
2073:       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
2074: #endif
2075:       dx            *= epsilon;
2076:       if (!PetscAbsScalar(dx)) SETERRQ(PETSC_ERR_PLIB,"Computed 0 differencing parameter");
2077:       w3_array[col] += dx;
2078:     }
2079:     w3_array = w3_array + start; VecRestoreArray(w3,&w3_array);

2081:     /*
2082:        Evaluate function at x1 + dx (here dx is a vector of perturbations)
2083:     */

2086:     (*f)(sctx,w3,w2,fctx);
2088:     VecAXPY(w2,-1.0,w1);

2090:     /*
2091:        Loop over rows of vector, putting results into Jacobian matrix
2092:     */
2093:     VecGetArray(w2,&y);
2094:     for (l=0; l<coloring->nrows[k]; l++) {
2095:       row    = coloring->rows[k][l];
2096:       col    = coloring->columnsforrow[k][l];
2097:       y[row] *= vscale_array[vscaleforrow[k][l]];
2098:       srow   = row + start;
2099:       MatSetValues_SeqAIJ(J,1,&srow,1,&col,y+row,INSERT_VALUES);
2100:     }
2101:     VecRestoreArray(w2,&y);
2102:   }
2103:   coloring->currentcolor = k;
2104:   VecRestoreArray(coloring->vscale,&vscale_array);
2105:   xx = xx + start; VecRestoreArray(x1,&xx);
2106:   MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
2107:   MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
2108:   return(0);
2109: }

2111:  #include petscblaslapack.h
2114: PetscErrorCode MatAXPY_SeqAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str)
2115: {
2117:   PetscInt       i;
2118:   Mat_SeqAIJ     *x  = (Mat_SeqAIJ *)X->data,*y = (Mat_SeqAIJ *)Y->data;
2119:   PetscBLASInt   one=1,bnz = (PetscBLASInt)x->nz;

2122:   if (str == SAME_NONZERO_PATTERN) {
2123:     PetscScalar alpha = a;
2124:     BLASaxpy_(&bnz,&alpha,x->a,&one,y->a,&one);
2125:   } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
2126:     if (y->xtoy && y->XtoY != X) {
2127:       PetscFree(y->xtoy);
2128:       MatDestroy(y->XtoY);
2129:     }
2130:     if (!y->xtoy) { /* get xtoy */
2131:       MatAXPYGetxtoy_Private(X->rmap.n,x->i,x->j,PETSC_NULL, y->i,y->j,PETSC_NULL, &y->xtoy);
2132:       y->XtoY = X;
2133:     }
2134:     for (i=0; i<x->nz; i++) y->a[y->xtoy[i]] += a*(x->a[i]);
2135:     PetscInfo3(0,"ratio of nnz(X)/nnz(Y): %d/%d = %G\n",x->nz,y->nz,(PetscReal)(x->nz)/y->nz);
2136:   } else {
2137:     MatAXPY_Basic(Y,a,X,str);
2138:   }
2139:   return(0);
2140: }

2144: PetscErrorCode MatSetBlockSize_SeqAIJ(Mat A,PetscInt bs)
2145: {
2147:   return(0);
2148: }

2152: PetscErrorCode  MatConjugate_SeqAIJ(Mat mat)
2153: {
2154: #if defined(PETSC_USE_COMPLEX)
2155:   Mat_SeqAIJ  *aij = (Mat_SeqAIJ *)mat->data;
2156:   PetscInt    i,nz;
2157:   PetscScalar *a;

2160:   nz = aij->nz;
2161:   a  = aij->a;
2162:   for (i=0; i<nz; i++) {
2163:     a[i] = PetscConj(a[i]);
2164:   }
2165: #else
2167: #endif
2168:   return(0);
2169: }

2173: PetscErrorCode MatGetRowMax_SeqAIJ(Mat A,Vec v)
2174: {
2175:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
2177:   PetscInt       i,j,m = A->rmap.n,*ai,*aj,ncols,n;
2178:   PetscReal      atmp;
2179:   PetscScalar    *x,zero = 0.0;
2180:   MatScalar      *aa;

2183:   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2184:   aa   = a->a;
2185:   ai   = a->i;
2186:   aj   = a->j;

2188:   VecSet(v,zero);
2189:   VecGetArray(v,&x);
2190:   VecGetLocalSize(v,&n);
2191:   if (n != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
2192:   for (i=0; i<m; i++) {
2193:     ncols = ai[1] - ai[0]; ai++;
2194:     for (j=0; j<ncols; j++){
2195:       atmp = PetscAbsScalar(*aa); aa++;
2196:       if (PetscAbsScalar(x[i]) < atmp) x[i] = atmp;
2197:       aj++;
2198:     }
2199:   }
2200:   VecRestoreArray(v,&x);
2201:   return(0);
2202: }

2204: /* -------------------------------------------------------------------*/
2205: static struct _MatOps MatOps_Values = {MatSetValues_SeqAIJ,
2206:        MatGetRow_SeqAIJ,
2207:        MatRestoreRow_SeqAIJ,
2208:        MatMult_SeqAIJ,
2209: /* 4*/ MatMultAdd_SeqAIJ,
2210:        MatMultTranspose_SeqAIJ,
2211:        MatMultTransposeAdd_SeqAIJ,
2212:        MatSolve_SeqAIJ,
2213:        MatSolveAdd_SeqAIJ,
2214:        MatSolveTranspose_SeqAIJ,
2215: /*10*/ MatSolveTransposeAdd_SeqAIJ,
2216:        MatLUFactor_SeqAIJ,
2217:        0,
2218:        MatRelax_SeqAIJ,
2219:        MatTranspose_SeqAIJ,
2220: /*15*/ MatGetInfo_SeqAIJ,
2221:        MatEqual_SeqAIJ,
2222:        MatGetDiagonal_SeqAIJ,
2223:        MatDiagonalScale_SeqAIJ,
2224:        MatNorm_SeqAIJ,
2225: /*20*/ 0,
2226:        MatAssemblyEnd_SeqAIJ,
2227:        MatCompress_SeqAIJ,
2228:        MatSetOption_SeqAIJ,
2229:        MatZeroEntries_SeqAIJ,
2230: /*25*/ MatZeroRows_SeqAIJ,
2231:        MatLUFactorSymbolic_SeqAIJ,
2232:        MatLUFactorNumeric_SeqAIJ,
2233:        MatCholeskyFactorSymbolic_SeqAIJ,
2234:        MatCholeskyFactorNumeric_SeqAIJ,
2235: /*30*/ MatSetUpPreallocation_SeqAIJ,
2236:        MatILUFactorSymbolic_SeqAIJ,
2237:        MatICCFactorSymbolic_SeqAIJ,
2238:        MatGetArray_SeqAIJ,
2239:        MatRestoreArray_SeqAIJ,
2240: /*35*/ MatDuplicate_SeqAIJ,
2241:        0,
2242:        0,
2243:        MatILUFactor_SeqAIJ,
2244:        0,
2245: /*40*/ MatAXPY_SeqAIJ,
2246:        MatGetSubMatrices_SeqAIJ,
2247:        MatIncreaseOverlap_SeqAIJ,
2248:        MatGetValues_SeqAIJ,
2249:        MatCopy_SeqAIJ,
2250: /*45*/ 0,
2251:        MatScale_SeqAIJ,
2252:        0,
2253:        MatDiagonalSet_SeqAIJ,
2254:        MatILUDTFactor_SeqAIJ,
2255: /*50*/ MatSetBlockSize_SeqAIJ,
2256:        MatGetRowIJ_SeqAIJ,
2257:        MatRestoreRowIJ_SeqAIJ,
2258:        MatGetColumnIJ_SeqAIJ,
2259:        MatRestoreColumnIJ_SeqAIJ,
2260: /*55*/ MatFDColoringCreate_SeqAIJ,
2261:        0,
2262:        0,
2263:        MatPermute_SeqAIJ,
2264:        0,
2265: /*60*/ 0,
2266:        MatDestroy_SeqAIJ,
2267:        MatView_SeqAIJ,
2268:        0,
2269:        0,
2270: /*65*/ 0,
2271:        0,
2272:        0,
2273:        0,
2274:        0,
2275: /*70*/ MatGetRowMax_SeqAIJ,
2276:        0,
2277:        MatSetColoring_SeqAIJ,
2278: #if defined(PETSC_HAVE_ADIC)
2279:        MatSetValuesAdic_SeqAIJ,
2280: #else
2281:        0,
2282: #endif
2283:        MatSetValuesAdifor_SeqAIJ,
2284: /*75*/ 0,
2285:        0,
2286:        0,
2287:        0,
2288:        0,
2289: /*80*/ 0,
2290:        0,
2291:        0,
2292:        0,
2293:        MatLoad_SeqAIJ,
2294: /*85*/ MatIsSymmetric_SeqAIJ,
2295:        0,
2296:        0,
2297:        0,
2298:        0,
2299: /*90*/ MatMatMult_SeqAIJ_SeqAIJ,
2300:        MatMatMultSymbolic_SeqAIJ_SeqAIJ,
2301:        MatMatMultNumeric_SeqAIJ_SeqAIJ,
2302:        MatPtAP_Basic,
2303:        MatPtAPSymbolic_SeqAIJ,
2304: /*95*/ MatPtAPNumeric_SeqAIJ,
2305:        MatMatMultTranspose_SeqAIJ_SeqAIJ,
2306:        MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ,
2307:        MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ,
2308:        MatPtAPSymbolic_SeqAIJ_SeqAIJ,
2309: /*100*/MatPtAPNumeric_SeqAIJ_SeqAIJ,
2310:        0,
2311:        0,
2312:        MatConjugate_SeqAIJ,
2313:        0,
2314: /*105*/MatSetValuesRow_SeqAIJ,
2315:        MatRealPart_SeqAIJ,
2316:        MatImaginaryPart_SeqAIJ,
2317:        0,
2318:        0,
2319: /*110*/MatMatSolve_SeqAIJ
2320: };

2325: PetscErrorCode  MatSeqAIJSetColumnIndices_SeqAIJ(Mat mat,PetscInt *indices)
2326: {
2327:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
2328:   PetscInt   i,nz,n;


2332:   nz = aij->maxnz;
2333:   n  = mat->cmap.n;
2334:   for (i=0; i<nz; i++) {
2335:     aij->j[i] = indices[i];
2336:   }
2337:   aij->nz = nz;
2338:   for (i=0; i<n; i++) {
2339:     aij->ilen[i] = aij->imax[i];
2340:   }

2342:   return(0);
2343: }

2348: /*@
2349:     MatSeqAIJSetColumnIndices - Set the column indices for all the rows
2350:        in the matrix.

2352:   Input Parameters:
2353: +  mat - the SeqAIJ matrix
2354: -  indices - the column indices

2356:   Level: advanced

2358:   Notes:
2359:     This can be called if you have precomputed the nonzero structure of the 
2360:   matrix and want to provide it to the matrix object to improve the performance
2361:   of the MatSetValues() operation.

2363:     You MUST have set the correct numbers of nonzeros per row in the call to 
2364:   MatCreateSeqAIJ(), and the columns indices MUST be sorted.

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

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

2370: @*/
2371: PetscErrorCode  MatSeqAIJSetColumnIndices(Mat mat,PetscInt *indices)
2372: {
2373:   PetscErrorCode ierr,(*f)(Mat,PetscInt *);

2378:   PetscObjectQueryFunction((PetscObject)mat,"MatSeqAIJSetColumnIndices_C",(void (**)(void))&f);
2379:   if (f) {
2380:     (*f)(mat,indices);
2381:   } else {
2382:     SETERRQ(PETSC_ERR_SUP,"Wrong type of matrix to set column indices");
2383:   }
2384:   return(0);
2385: }

2387: /* ----------------------------------------------------------------------------------------*/

2392: PetscErrorCode  MatStoreValues_SeqAIJ(Mat mat)
2393: {
2394:   Mat_SeqAIJ     *aij = (Mat_SeqAIJ *)mat->data;
2396:   size_t         nz = aij->i[mat->rmap.n];

2399:   if (aij->nonew != 1) {
2400:     SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
2401:   }

2403:   /* allocate space for values if not already there */
2404:   if (!aij->saved_values) {
2405:     PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);
2406:   }

2408:   /* copy values over */
2409:   PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));
2410:   return(0);
2411: }

2416: /*@
2417:     MatStoreValues - Stashes a copy of the matrix values; this allows, for 
2418:        example, reuse of the linear part of a Jacobian, while recomputing the 
2419:        nonlinear portion.

2421:    Collect on Mat

2423:   Input Parameters:
2424: .  mat - the matrix (currently only AIJ matrices support this option)

2426:   Level: advanced

2428:   Common Usage, with SNESSolve():
2429: $    Create Jacobian matrix
2430: $    Set linear terms into matrix
2431: $    Apply boundary conditions to matrix, at this time matrix must have 
2432: $      final nonzero structure (i.e. setting the nonlinear terms and applying 
2433: $      boundary conditions again will not change the nonzero structure
2434: $    MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS);
2435: $    MatStoreValues(mat);
2436: $    Call SNESSetJacobian() with matrix
2437: $    In your Jacobian routine
2438: $      MatRetrieveValues(mat);
2439: $      Set nonlinear terms in matrix
2440:  
2441:   Common Usage without SNESSolve(), i.e. when you handle nonlinear solve yourself:
2442: $    // build linear portion of Jacobian 
2443: $    MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS);
2444: $    MatStoreValues(mat);
2445: $    loop over nonlinear iterations
2446: $       MatRetrieveValues(mat);
2447: $       // call MatSetValues(mat,...) to set nonliner portion of Jacobian 
2448: $       // call MatAssemblyBegin/End() on matrix
2449: $       Solve linear system with Jacobian
2450: $    endloop 

2452:   Notes:
2453:     Matrix must already be assemblied before calling this routine
2454:     Must set the matrix option MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS); before 
2455:     calling this routine.

2457:     When this is called multiple times it overwrites the previous set of stored values
2458:     and does not allocated additional space.

2460: .seealso: MatRetrieveValues()

2462: @*/
2463: PetscErrorCode  MatStoreValues(Mat mat)
2464: {
2465:   PetscErrorCode ierr,(*f)(Mat);

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

2472:   PetscObjectQueryFunction((PetscObject)mat,"MatStoreValues_C",(void (**)(void))&f);
2473:   if (f) {
2474:     (*f)(mat);
2475:   } else {
2476:     SETERRQ(PETSC_ERR_SUP,"Wrong type of matrix to store values");
2477:   }
2478:   return(0);
2479: }

2484: PetscErrorCode  MatRetrieveValues_SeqAIJ(Mat mat)
2485: {
2486:   Mat_SeqAIJ     *aij = (Mat_SeqAIJ *)mat->data;
2488:   PetscInt       nz = aij->i[mat->rmap.n];

2491:   if (aij->nonew != 1) {
2492:     SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
2493:   }
2494:   if (!aij->saved_values) {
2495:     SETERRQ(PETSC_ERR_ORDER,"Must call MatStoreValues(A);first");
2496:   }
2497:   /* copy values over */
2498:   PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));
2499:   return(0);
2500: }

2505: /*@
2506:     MatRetrieveValues - Retrieves the copy of the matrix values; this allows, for 
2507:        example, reuse of the linear part of a Jacobian, while recomputing the 
2508:        nonlinear portion.

2510:    Collect on Mat

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

2515:   Level: advanced

2517: .seealso: MatStoreValues()

2519: @*/
2520: PetscErrorCode  MatRetrieveValues(Mat mat)
2521: {
2522:   PetscErrorCode ierr,(*f)(Mat);

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

2529:   PetscObjectQueryFunction((PetscObject)mat,"MatRetrieveValues_C",(void (**)(void))&f);
2530:   if (f) {
2531:     (*f)(mat);
2532:   } else {
2533:     SETERRQ(PETSC_ERR_SUP,"Wrong type of matrix to retrieve values");
2534:   }
2535:   return(0);
2536: }


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

2549:    Collective on MPI_Comm

2551:    Input Parameters:
2552: +  comm - MPI communicator, set to PETSC_COMM_SELF
2553: .  m - number of rows
2554: .  n - number of columns
2555: .  nz - number of nonzeros per row (same for all rows)
2556: -  nnz - array containing the number of nonzeros in the various rows 
2557:          (possibly different for each row) or PETSC_NULL

2559:    Output Parameter:
2560: .  A - the matrix 

2562:    Notes:
2563:    If nnz is given then nz is ignored

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

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

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

2580:    Options Database Keys:
2581: +  -mat_no_inode  - Do not use inodes
2582: .  -mat_inode_limit <limit> - Sets inode limit (max limit=5)
2583: -  -mat_aij_oneindex - Internally use indexing starting at 1
2584:         rather than 0.  Note that when calling MatSetValues(),
2585:         the user still MUST index entries starting at 0!

2587:    Level: intermediate

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

2591: @*/
2592: PetscErrorCode  MatCreateSeqAIJ(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A)
2593: {

2597:   MatCreate(comm,A);
2598:   MatSetSizes(*A,m,n,m,n);
2599:   MatSetType(*A,MATSEQAIJ);
2600:   MatSeqAIJSetPreallocation_SeqAIJ(*A,nz,(PetscInt*)nnz);
2601:   return(0);
2602: }

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

2612:    Collective on MPI_Comm

2614:    Input Parameters:
2615: +  B - The matrix
2616: .  nz - number of nonzeros per row (same for all rows)
2617: -  nnz - array containing the number of nonzeros in the various rows 
2618:          (possibly different for each row) or PETSC_NULL

2620:    Notes:
2621:      If nnz is given then nz is ignored

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

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

2633:    Developers: Use nz of MAT_SKIP_ALLOCATION to not allocate any space for the matrix
2634:    entries or columns indices

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

2641:    Options Database Keys:
2642: +  -mat_no_inode  - Do not use inodes
2643: .  -mat_inode_limit <limit> - Sets inode limit (max limit=5)
2644: -  -mat_aij_oneindex - Internally use indexing starting at 1
2645:         rather than 0.  Note that when calling MatSetValues(),
2646:         the user still MUST index entries starting at 0!

2648:    Level: intermediate

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

2652: @*/
2653: PetscErrorCode  MatSeqAIJSetPreallocation(Mat B,PetscInt nz,const PetscInt nnz[])
2654: {
2655:   PetscErrorCode ierr,(*f)(Mat,PetscInt,const PetscInt[]);

2658:   PetscObjectQueryFunction((PetscObject)B,"MatSeqAIJSetPreallocation_C",(void (**)(void))&f);
2659:   if (f) {
2660:     (*f)(B,nz,nnz);
2661:   }
2662:   return(0);
2663: }

2668: PetscErrorCode  MatSeqAIJSetPreallocation_SeqAIJ(Mat B,PetscInt nz,PetscInt *nnz)
2669: {
2670:   Mat_SeqAIJ     *b;
2671:   PetscTruth     skipallocation = PETSC_FALSE;
2673:   PetscInt       i;

2676: 
2677:   if (nz == MAT_SKIP_ALLOCATION) {
2678:     skipallocation = PETSC_TRUE;
2679:     nz             = 0;
2680:   }

2682:   B->rmap.bs = B->cmap.bs = 1;
2683:   PetscMapInitialize(B->comm,&B->rmap);
2684:   PetscMapInitialize(B->comm,&B->cmap);

2686:   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
2687:   if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz);
2688:   if (nnz) {
2689:     for (i=0; i<B->rmap.n; i++) {
2690:       if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]);
2691:       if (nnz[i] > B->cmap.n) SETERRQ3(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be greater than row length: local row %d value %d rowlength %d",i,nnz[i],B->cmap.n);
2692:     }
2693:   }

2695:   B->preallocated = PETSC_TRUE;
2696:   b = (Mat_SeqAIJ*)B->data;

2698:   if (!skipallocation) {
2699:     PetscMalloc2(B->rmap.n,PetscInt,&b->imax,B->rmap.n,PetscInt,&b->ilen);
2700:     if (!nnz) {
2701:       if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 10;
2702:       else if (nz <= 0)        nz = 1;
2703:       for (i=0; i<B->rmap.n; i++) b->imax[i] = nz;
2704:       nz = nz*B->rmap.n;
2705:     } else {
2706:       nz = 0;
2707:       for (i=0; i<B->rmap.n; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
2708:     }

2710:     /* b->ilen will count nonzeros in each row so far. */
2711:     for (i=0; i<B->rmap.n; i++) { b->ilen[i] = 0;}

2713:     /* allocate the matrix space */
2714:     PetscMalloc3(nz,PetscScalar,&b->a,nz,PetscInt,&b->j,B->rmap.n+1,PetscInt,&b->i);
2715:     b->i[0] = 0;
2716:     for (i=1; i<B->rmap.n+1; i++) {
2717:       b->i[i] = b->i[i-1] + b->imax[i-1];
2718:     }
2719:     b->singlemalloc = PETSC_TRUE;
2720:     b->free_a       = PETSC_TRUE;
2721:     b->free_ij      = PETSC_TRUE;
2722:   } else {
2723:     b->free_a       = PETSC_FALSE;
2724:     b->free_ij      = PETSC_FALSE;
2725:   }

2727:   b->nz                = 0;
2728:   b->maxnz             = nz;
2729:   B->info.nz_unneeded  = (double)b->maxnz;
2730:   return(0);
2731: }

2734: #undef  __FUNCT__
2736: /*@C
2737:    MatSeqAIJSetPreallocationCSR - Allocates memory for a sparse sequential matrix in AIJ format.  

2739:    Input Parameters:
2740: +  B - the matrix 
2741: .  i - the indices into j for the start of each row (starts with zero)
2742: .  j - the column indices for each row (starts with zero) these must be sorted for each row
2743: -  v - optional values in the matrix

2745:    Contributed by: Lisandro Dalchin

2747:    Level: developer

2749: .keywords: matrix, aij, compressed row, sparse, sequential

2751: .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatSeqAIJSetPreallocation(), MatCreateSeqAIJ(), SeqAIJ
2752: @*/
2753: PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat B,const PetscInt i[],const PetscInt j[],const PetscScalar v[])
2754: {
2755:   PetscErrorCode (*f)(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);

2760:   PetscObjectQueryFunction((PetscObject)B,"MatSeqAIJSetPreallocationCSR_C",(void (**)(void))&f);
2761:   if (f) {
2762:     (*f)(B,i,j,v);
2763:   }
2764:   return(0);
2765: }

2768: #undef  __FUNCT__
2770: PetscErrorCode  MatSeqAIJSetPreallocationCSR_SeqAIJ(Mat B,const PetscInt Ii[],const PetscInt J[],const PetscScalar v[])
2771: {
2772:   PetscInt       i;
2773:   PetscInt       m,n;
2774:   PetscInt       nz;
2775:   PetscInt       *nnz, nz_max = 0;
2776:   PetscScalar    *values;

2780:   MatGetSize(B, &m, &n);

2782:   if (Ii[0]) {
2783:     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Ii[0] must be 0 it is %D", Ii[0]);
2784:   }
2785:   PetscMalloc((m+1) * sizeof(PetscInt), &nnz);
2786:   for(i = 0; i < m; i++) {
2787:     nz     = Ii[i+1]- Ii[i];
2788:     nz_max = PetscMax(nz_max, nz);
2789:     if (nz < 0) {
2790:       SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE, "Local row %D has a negative number of columns %D", i, nnz);
2791:     }
2792:     nnz[i] = nz;
2793:   }
2794:   MatSeqAIJSetPreallocation(B, 0, nnz);
2795:   PetscFree(nnz);

2797:   if (v) {
2798:     values = (PetscScalar*) v;
2799:   } else {
2800:     PetscMalloc((nz_max+1)*sizeof(PetscScalar), &values);
2801:     PetscMemzero(values, nz_max*sizeof(PetscScalar));
2802:   }

2804:   MatSetOption(B,MAT_COLUMNS_SORTED);

2806:   for(i = 0; i < m; i++) {
2807:     nz  = Ii[i+1] - Ii[i];
2808:     MatSetValues_SeqAIJ(B, 1, &i, nz, J+Ii[i], values + (v ? Ii[i] : 0), INSERT_VALUES);
2809:   }

2811:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
2812:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
2813:   MatSetOption(B,MAT_COLUMNS_UNSORTED);

2815:   if (!v) {
2816:     PetscFree(values);
2817:   }
2818:   return(0);
2819: }

2822: /*MC
2823:    MATSEQAIJ - MATSEQAIJ = "seqaij" - A matrix type to be used for sequential sparse matrices, 
2824:    based on compressed sparse row format.

2826:    Options Database Keys:
2827: . -mat_type seqaij - sets the matrix type to "seqaij" during a call to MatSetFromOptions()

2829:   Level: beginner

2831: .seealso: MatCreateSeqAIJ(), MatSetFromOptions(), MatSetType(), MatCreate(), MatType
2832: M*/


2841: PetscErrorCode  MatCreate_SeqAIJ(Mat B)
2842: {
2843:   Mat_SeqAIJ     *b;
2845:   PetscMPIInt    size;

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

2851:   PetscNew(Mat_SeqAIJ,&b);
2852:   B->data             = (void*)b;
2853:   PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
2854:   B->factor           = 0;
2855:   B->mapping          = 0;
2856:   b->row              = 0;
2857:   b->col              = 0;
2858:   b->icol             = 0;
2859:   b->reallocs         = 0;
2860:   b->sorted            = PETSC_FALSE;
2861:   b->ignorezeroentries = PETSC_FALSE;
2862:   b->roworiented       = PETSC_TRUE;
2863:   b->nonew             = 0;
2864:   b->diag              = 0;
2865:   b->solve_work        = 0;
2866:   B->spptr             = 0;
2867:   b->saved_values      = 0;
2868:   b->idiag             = 0;
2869:   b->ssor              = 0;
2870:   b->keepzeroedrows    = PETSC_FALSE;
2871:   b->xtoy              = 0;
2872:   b->XtoY              = 0;
2873:   b->compressedrow.use     = PETSC_FALSE;
2874:   b->compressedrow.nrows   = B->rmap.n;
2875:   b->compressedrow.i       = PETSC_NULL;
2876:   b->compressedrow.rindex  = PETSC_NULL;
2877:   b->compressedrow.checked = PETSC_FALSE;
2878:   B->same_nonzero          = PETSC_FALSE;

2880:   PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);
2881:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqAIJSetColumnIndices_C",
2882:                                      "MatSeqAIJSetColumnIndices_SeqAIJ",
2883:                                      MatSeqAIJSetColumnIndices_SeqAIJ);
2884:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
2885:                                      "MatStoreValues_SeqAIJ",
2886:                                      MatStoreValues_SeqAIJ);
2887:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
2888:                                      "MatRetrieveValues_SeqAIJ",
2889:                                      MatRetrieveValues_SeqAIJ);
2890:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_seqsbaij_C",
2891:                                      "MatConvert_SeqAIJ_SeqSBAIJ",
2892:                                       MatConvert_SeqAIJ_SeqSBAIJ);
2893:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_seqbaij_C",
2894:                                      "MatConvert_SeqAIJ_SeqBAIJ",
2895:                                       MatConvert_SeqAIJ_SeqBAIJ);
2896:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_seqcsrperm_C",
2897:                                      "MatConvert_SeqAIJ_SeqCSRPERM",
2898:                                       MatConvert_SeqAIJ_SeqCSRPERM);
2899:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_seqcrl_C",
2900:                                      "MatConvert_SeqAIJ_SeqCRL",
2901:                                       MatConvert_SeqAIJ_SeqCRL);
2902:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatIsTranspose_C",
2903:                                      "MatIsTranspose_SeqAIJ",
2904:                                       MatIsTranspose_SeqAIJ);
2905:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqAIJSetPreallocation_C",
2906:                                      "MatSeqAIJSetPreallocation_SeqAIJ",
2907:                                       MatSeqAIJSetPreallocation_SeqAIJ);
2908:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqAIJSetPreallocationCSR_C",
2909:                                      "MatSeqAIJSetPreallocationCSR_SeqAIJ",
2910:                                       MatSeqAIJSetPreallocationCSR_SeqAIJ);
2911:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatReorderForNonzeroDiagonal_C",
2912:                                      "MatReorderForNonzeroDiagonal_SeqAIJ",
2913:                                       MatReorderForNonzeroDiagonal_SeqAIJ);
2914:   MatCreate_Inode(B);
2915:   PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);
2916:   return(0);
2917: }

2922: PetscErrorCode MatDuplicate_SeqAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
2923: {
2924:   Mat            C;
2925:   Mat_SeqAIJ     *c,*a = (Mat_SeqAIJ*)A->data;
2927:   PetscInt       i,m = A->rmap.n;

2930:   *B = 0;
2931:   MatCreate(A->comm,&C);
2932:   MatSetSizes(C,A->rmap.n,A->cmap.n,A->rmap.n,A->cmap.n);
2933:   MatSetType(C,A->type_name);
2934:   PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));
2935: 
2936:   PetscMapCopy(A->comm,&A->rmap,&C->rmap);
2937:   PetscMapCopy(A->comm,&A->cmap,&C->cmap);

2939:   c = (Mat_SeqAIJ*)C->data;

2941:   C->factor           = A->factor;

2943:   c->row            = 0;
2944:   c->col            = 0;
2945:   c->icol           = 0;
2946:   c->reallocs       = 0;

2948:   C->assembled      = PETSC_TRUE;
2949: 
2950:   PetscMalloc2(m,PetscInt,&c->imax,m,PetscInt,&c->ilen);
2951:   for (i=0; i<m; i++) {
2952:     c->imax[i] = a->imax[i];
2953:     c->ilen[i] = a->ilen[i];
2954:   }

2956:   /* allocate the matrix space */
2957:   PetscMalloc3(a->i[m],PetscScalar,&c->a,a->i[m],PetscInt,&c->j,m+1,PetscInt,&c->i);
2958:   c->singlemalloc = PETSC_TRUE;
2959:   PetscMemcpy(c->i,a->i,(m+1)*sizeof(PetscInt));
2960:   if (m > 0) {
2961:     PetscMemcpy(c->j,a->j,(a->i[m])*sizeof(PetscInt));
2962:     if (cpvalues == MAT_COPY_VALUES) {
2963:       PetscMemcpy(c->a,a->a,(a->i[m])*sizeof(PetscScalar));
2964:     } else {
2965:       PetscMemzero(c->a,(a->i[m])*sizeof(PetscScalar));
2966:     }
2967:   }

2969:   c->sorted            = a->sorted;
2970:   c->ignorezeroentries = a->ignorezeroentries;
2971:   c->roworiented       = a->roworiented;
2972:   c->nonew             = a->nonew;
2973:   if (a->diag) {
2974:     PetscMalloc((m+1)*sizeof(PetscInt),&c->diag);
2975:     PetscLogObjectMemory(C,(m+1)*sizeof(PetscInt));
2976:     for (i=0; i<m; i++) {
2977:       c->diag[i] = a->diag[i];
2978:     }
2979:   } else c->diag        = 0;
2980:   c->solve_work         = 0;
2981:   c->saved_values          = 0;
2982:   c->idiag                 = 0;
2983:   c->ssor                  = 0;
2984:   c->keepzeroedrows        = a->keepzeroedrows;
2985:   c->free_a                = PETSC_TRUE;
2986:   c->free_ij               = PETSC_TRUE;
2987:   c->xtoy                  = 0;
2988:   c->XtoY                  = 0;

2990:   c->nz                 = a->nz;
2991:   c->maxnz              = a->maxnz;
2992:   C->preallocated       = PETSC_TRUE;

2994:   c->compressedrow.use     = a->compressedrow.use;
2995:   c->compressedrow.nrows   = a->compressedrow.nrows;
2996:   c->compressedrow.checked = a->compressedrow.checked;
2997:   if ( a->compressedrow.checked && a->compressedrow.use){
2998:     i = a->compressedrow.nrows;
2999:     PetscMalloc((2*i+1)*sizeof(PetscInt),&c->compressedrow.i);
3000:     c->compressedrow.rindex = c->compressedrow.i + i + 1;
3001:     PetscMemcpy(c->compressedrow.i,a->compressedrow.i,(i+1)*sizeof(PetscInt));
3002:     PetscMemcpy(c->compressedrow.rindex,a->compressedrow.rindex,i*sizeof(PetscInt));
3003:   } else {
3004:     c->compressedrow.use    = PETSC_FALSE;
3005:     c->compressedrow.i      = PETSC_NULL;
3006:     c->compressedrow.rindex = PETSC_NULL;
3007:   }
3008:   C->same_nonzero = A->same_nonzero;
3009:   MatDuplicate_Inode(A,cpvalues,&C);

3011:   *B = C;
3012:   PetscFListDuplicate(A->qlist,&C->qlist);
3013:   return(0);
3014: }

3018: PetscErrorCode MatLoad_SeqAIJ(PetscViewer viewer, MatType type,Mat *A)
3019: {
3020:   Mat_SeqAIJ     *a;
3021:   Mat            B;
3023:   PetscInt       i,sum,nz,header[4],*rowlengths = 0,M,N;
3024:   int            fd;
3025:   PetscMPIInt    size;
3026:   MPI_Comm       comm;
3027: 
3029:   PetscObjectGetComm((PetscObject)viewer,&comm);
3030:   MPI_Comm_size(comm,&size);
3031:   if (size > 1) SETERRQ(PETSC_ERR_ARG_SIZ,"view must have one processor");
3032:   PetscViewerBinaryGetDescriptor(viewer,&fd);
3033:   PetscBinaryRead(fd,header,4,PETSC_INT);
3034:   if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object in file");
3035:   M = header[1]; N = header[2]; nz = header[3];

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

3041:   /* read in row lengths */
3042:   PetscMalloc(M*sizeof(PetscInt),&rowlengths);
3043:   PetscBinaryRead(fd,rowlengths,M,PETSC_INT);

3045:   /* check if sum of rowlengths is same as nz */
3046:   for (i=0,sum=0; i< M; i++) sum +=rowlengths[i];
3047:   if (sum != nz) SETERRQ2(PETSC_ERR_FILE_READ,"Inconsistant matrix data in file. no-nonzeros = %d, sum-row-lengths = %d\n",nz,sum);

3049:   /* create our matrix */
3050:   MatCreate(comm,&B);
3051:   MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,M,N);
3052:   MatSetType(B,type);
3053:   MatSeqAIJSetPreallocation_SeqAIJ(B,0,rowlengths);
3054:   a = (Mat_SeqAIJ*)B->data;

3056:   PetscBinaryRead(fd,a->j,nz,PETSC_INT);

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

3061:   /* set matrix "i" values */
3062:   a->i[0] = 0;
3063:   for (i=1; i<= M; i++) {
3064:     a->i[i]      = a->i[i-1] + rowlengths[i-1];
3065:     a->ilen[i-1] = rowlengths[i-1];
3066:   }
3067:   PetscFree(rowlengths);

3069:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
3070:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
3071:   *A = B;
3072:   return(0);
3073: }

3077: PetscErrorCode MatEqual_SeqAIJ(Mat A,Mat B,PetscTruth* flg)
3078: {
3079:   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data,*b = (Mat_SeqAIJ *)B->data;

3083:   /* If the  matrix dimensions are not equal,or no of nonzeros */
3084:   if ((A->rmap.n != B->rmap.n) || (A->cmap.n != B->cmap.n) ||(a->nz != b->nz)) {
3085:     *flg = PETSC_FALSE;
3086:     return(0);
3087:   }
3088: 
3089:   /* if the a->i are the same */
3090:   PetscMemcmp(a->i,b->i,(A->rmap.n+1)*sizeof(PetscInt),flg);
3091:   if (!*flg) return(0);
3092: 
3093:   /* if a->j are the same */
3094:   PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),flg);
3095:   if (!*flg) return(0);
3096: 
3097:   /* if a->a are the same */
3098:   PetscMemcmp(a->a,b->a,(a->nz)*sizeof(PetscScalar),flg);

3100:   return(0);
3101: 
3102: }

3106: /*@
3107:      MatCreateSeqAIJWithArrays - Creates an sequential AIJ matrix using matrix elements (in CSR format)
3108:               provided by the user.

3110:       Collective on MPI_Comm

3112:    Input Parameters:
3113: +   comm - must be an MPI communicator of size 1
3114: .   m - number of rows
3115: .   n - number of columns
3116: .   i - row indices
3117: .   j - column indices
3118: -   a - matrix values

3120:    Output Parameter:
3121: .   mat - the matrix

3123:    Level: intermediate

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

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

3131:        The i and j indices are 0 based

3133: .seealso: MatCreate(), MatCreateMPIAIJ(), MatCreateSeqAIJ(), MatCreateMPIAIJWithArrays(), MatMPIAIJSetPreallocationCSR()

3135: @*/
3136: PetscErrorCode  MatCreateSeqAIJWithArrays(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt* i,PetscInt*j,PetscScalar *a,Mat *mat)
3137: {
3139:   PetscInt       ii;
3140:   Mat_SeqAIJ     *aij;

3143:   if (i[0]) {
3144:     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"i (row indices) must start with 0");
3145:   }
3146:   MatCreate(comm,mat);
3147:   MatSetSizes(*mat,m,n,m,n);
3148:   MatSetType(*mat,MATSEQAIJ);
3149:   MatSeqAIJSetPreallocation_SeqAIJ(*mat,MAT_SKIP_ALLOCATION,0);
3150:   aij  = (Mat_SeqAIJ*)(*mat)->data;
3151:   PetscMalloc2(m,PetscInt,&aij->imax,m,PetscInt,&aij->ilen);

3153:   aij->i = i;
3154:   aij->j = j;
3155:   aij->a = a;
3156:   aij->singlemalloc = PETSC_FALSE;
3157:   aij->nonew        = -1;             /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/
3158:   aij->free_a       = PETSC_FALSE;
3159:   aij->free_ij      = PETSC_FALSE;

3161:   for (ii=0; ii<m; ii++) {
3162:     aij->ilen[ii] = aij->imax[ii] = i[ii+1] - i[ii];
3163: #if defined(PETSC_USE_DEBUG)
3164:     if (i[ii+1] - i[ii] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Negative row length in i (row indices) row = %d length = %d",ii,i[ii+1] - i[ii]);
3165: #endif    
3166:   }
3167: #if defined(PETSC_USE_DEBUG)
3168:   for (ii=0; ii<aij->i[m]; ii++) {
3169:     if (j[ii] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Negative column index at location = %d index = %d",ii,j[ii]);
3170:     if (j[ii] > n - 1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column index to large at location = %d index = %d",ii,j[ii]);
3171:   }
3172: #endif    

3174:   MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
3175:   MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
3176:   return(0);
3177: }

3181: PetscErrorCode MatSetColoring_SeqAIJ(Mat A,ISColoring coloring)
3182: {
3184:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;

3187:   if (coloring->ctype == IS_COLORING_GLOBAL) {
3188:     ISColoringReference(coloring);
3189:     a->coloring = coloring;
3190:   } else if (coloring->ctype == IS_COLORING_GHOSTED) {
3191:     PetscInt             i,*larray;
3192:     ISColoring      ocoloring;
3193:     ISColoringValue *colors;

3195:     /* set coloring for diagonal portion */
3196:     PetscMalloc((A->cmap.n+1)*sizeof(PetscInt),&larray);
3197:     for (i=0; i<A->cmap.n; i++) {
3198:       larray[i] = i;
3199:     }
3200:     ISGlobalToLocalMappingApply(A->mapping,IS_GTOLM_MASK,A->cmap.n,larray,PETSC_NULL,larray);
3201:     PetscMalloc((A->cmap.n+1)*sizeof(ISColoringValue),&colors);
3202:     for (i=0; i<A->cmap.n; i++) {
3203:       colors[i] = coloring->colors[larray[i]];
3204:     }
3205:     PetscFree(larray);
3206:     ISColoringCreate(PETSC_COMM_SELF,coloring->n,A->cmap.n,colors,&ocoloring);
3207:     a->coloring = ocoloring;
3208:   }
3209:   return(0);
3210: }

3212: #if defined(PETSC_HAVE_ADIC)
3214: #include "adic/ad_utils.h"

3219: PetscErrorCode MatSetValuesAdic_SeqAIJ(Mat A,void *advalues)
3220: {
3221:   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
3222:   PetscInt        m = A->rmap.n,*ii = a->i,*jj = a->j,nz,i,j,nlen;
3223:   PetscScalar     *v = a->a,*values = ((PetscScalar*)advalues)+1;
3224:   ISColoringValue *color;

3227:   if (!a->coloring) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Coloring not set for matrix");
3228:   nlen  = PetscADGetDerivTypeSize()/sizeof(PetscScalar);
3229:   color = a->coloring->colors;
3230:   /* loop over rows */
3231:   for (i=0; i<m; i++) {
3232:     nz = ii[i+1] - ii[i];
3233:     /* loop over columns putting computed value into matrix */
3234:     for (j=0; j<nz; j++) {
3235:       *v++ = values[color[*jj++]];
3236:     }
3237:     values += nlen; /* jump to next row of derivatives */
3238:   }
3239:   return(0);
3240: }
3241: #endif

3245: PetscErrorCode MatSetValuesAdifor_SeqAIJ(Mat A,PetscInt nl,void *advalues)
3246: {
3247:   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
3248:   PetscInt             m = A->rmap.n,*ii = a->i,*jj = a->j,nz,i,j;
3249:   PetscScalar     *v = a->a,*values = (PetscScalar *)advalues;
3250:   ISColoringValue *color;

3253:   if (!a->coloring) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Coloring not set for matrix");
3254:   color = a->coloring->colors;
3255:   /* loop over rows */
3256:   for (i=0; i<m; i++) {
3257:     nz = ii[i+1] - ii[i];
3258:     /* loop over columns putting computed value into matrix */
3259:     for (j=0; j<nz; j++) {
3260:       *v++ = values[color[*jj++]];
3261:     }
3262:     values += nl; /* jump to next row of derivatives */
3263:   }
3264:   return(0);
3265: }

3267: /*
3268:     Special version for direct calls from Fortran 
3269: */
3270: #if defined(PETSC_HAVE_FORTRAN_CAPS)
3271: #define matsetvaluesseqaij_ MATSETVALUESSEQAIJ
3272: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
3273: #define matsetvaluesseqaij_ matsetvaluesseqaij
3274: #endif

3276: /* Change these macros so can be used in void function */
3277: #undef CHKERRQ
3278: #define CHKERRQ(ierr) CHKERRABORT(A->comm,ierr) 
3279: #undef SETERRQ2
3280: #define SETERRQ2(ierr,b,c,d) CHKERRABORT(A->comm,ierr) 

3285: void PETSC_STDCALL matsetvaluesseqaij_(Mat *AA,PetscInt *mm,const PetscInt im[],PetscInt *nn,const PetscInt in[],const PetscScalar v[],InsertMode *isis, PetscErrorCode *_ierr)
3286: {
3287:   Mat            A = *AA;
3288:   PetscInt       m = *mm, n = *nn;
3289:   InsertMode     is = *isis;
3290:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
3291:   PetscInt       *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N;
3292:   PetscInt       *imax,*ai,*ailen;
3294:   PetscInt       *aj,nonew = a->nonew,lastcol = -1;
3295:   PetscScalar    *ap,value,*aa;
3296:   PetscTruth     ignorezeroentries = ((a->ignorezeroentries && is == ADD_VALUES) ? PETSC_TRUE:PETSC_FALSE);
3297:   PetscTruth     roworiented = a->roworiented;

3300:   MatPreallocated(A);
3301:   imax = a->imax;
3302:   ai = a->i;
3303:   ailen = a->ilen;
3304:   aj = a->j;
3305:   aa = a->a;

3307:   for (k=0; k<m; k++) { /* loop over added rows */
3308:     row  = im[k];
3309:     if (row < 0) continue;
3310: #if defined(PETSC_USE_DEBUG)  
3311:     if (row >= A->rmap.n) SETERRABORT(A->comm,PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
3312: #endif
3313:     rp   = aj + ai[row]; ap = aa + ai[row];
3314:     rmax = imax[row]; nrow = ailen[row];
3315:     low  = 0;
3316:     high = nrow;
3317:     for (l=0; l<n; l++) { /* loop over added columns */
3318:       if (in[l] < 0) continue;
3319: #if defined(PETSC_USE_DEBUG)  
3320:       if (in[l] >= A->cmap.n) SETERRABORT(A->comm,PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
3321: #endif
3322:       col = in[l];
3323:       if (roworiented) {
3324:         value = v[l + k*n];
3325:       } else {
3326:         value = v[k + l*m];
3327:       }
3328:       if (value == 0.0 && ignorezeroentries && (is == ADD_VALUES)) continue;

3330:       if (col <= lastcol) low = 0; else high = nrow;
3331:       lastcol = col;
3332:       while (high-low > 5) {
3333:         t = (low+high)/2;
3334:         if (rp[t] > col) high = t;
3335:         else             low  = t;
3336:       }
3337:       for (i=low; i<high; i++) {
3338:         if (rp[i] > col) break;
3339:         if (rp[i] == col) {
3340:           if (is == ADD_VALUES) ap[i] += value;
3341:           else                  ap[i] = value;
3342:           goto noinsert;
3343:         }
3344:       }
3345:       if (value == 0.0 && ignorezeroentries) goto noinsert;
3346:       if (nonew == 1) goto noinsert;
3347:       if (nonew == -1) SETERRABORT(A->comm,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
3348:       MatSeqXAIJReallocateAIJ(A,A->rmap.n,1,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar);
3349:       N = nrow++ - 1; a->nz++; high++;
3350:       /* shift up all the later entries in this row */
3351:       for (ii=N; ii>=i; ii--) {
3352:         rp[ii+1] = rp[ii];
3353:         ap[ii+1] = ap[ii];
3354:       }
3355:       rp[i] = col;
3356:       ap[i] = value;
3357:       noinsert:;
3358:       low = i + 1;
3359:     }
3360:     ailen[row] = nrow;
3361:   }
3362:   A->same_nonzero = PETSC_FALSE;
3363:   PetscFunctionReturnVoid();
3364: }