Actual source code: aij.c

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

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

 14: PetscErrorCode MatGetRowIJ_SeqAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *m,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
 15: {
 16:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
 18:   PetscInt       i,ishift;
 19: 
 21:   *m     = A->m;
 22:   if (!ia) return(0);
 23:   ishift = 0;
 24:   if (symmetric && !A->structurally_symmetric) {
 25:     MatToSymmetricIJ_SeqAIJ(A->m,a->i,a->j,ishift,oshift,ia,ja);
 26:   } else if (oshift == 1) {
 27:     PetscInt nz = a->i[A->m];
 28:     /* malloc space and  add 1 to i and j indices */
 29:     PetscMalloc((A->m+1)*sizeof(PetscInt),ia);
 30:     PetscMalloc((nz+1)*sizeof(PetscInt),ja);
 31:     for (i=0; i<nz; i++) (*ja)[i] = a->j[i] + 1;
 32:     for (i=0; i<A->m+1; i++) (*ia)[i] = a->i[i] + 1;
 33:   } else {
 34:     *ia = a->i; *ja = a->j;
 35:   }
 36:   return(0);
 37: }

 41: PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *n,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
 42: {
 44: 
 46:   if (!ia) return(0);
 47:   if ((symmetric && !A->structurally_symmetric) || oshift == 1) {
 48:     PetscFree(*ia);
 49:     PetscFree(*ja);
 50:   }
 51:   return(0);
 52: }

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

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

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

105:   PetscFree(*ia);
106:   PetscFree(*ja);
107: 
108:   return(0);
109: }

111: #define CHUNKSIZE   15

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

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

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

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

173:         /* malloc new storage space */
174:         len     = ((size_t) new_nz)*(sizeof(PetscInt)+sizeof(PetscScalar))+(A->m+1)*sizeof(PetscInt);
175:         PetscMalloc(len,&new_a);
176:         new_j   = (PetscInt*)(new_a + new_nz);
177:         new_i   = new_j + new_nz;

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

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

220: PetscErrorCode MatGetValues_SeqAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],PetscScalar v[])
221: {
222:   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)A->data;
223:   PetscInt     *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
224:   PetscInt     *ai = a->i,*ailen = a->ilen;
225:   PetscScalar  *ap,*aa = a->a,zero = 0.0;

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


261: PetscErrorCode MatView_SeqAIJ_Binary(Mat A,PetscViewer viewer)
262: {
263:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
265:   PetscInt       i,*col_lens;
266:   int            fd;

269:   PetscViewerBinaryGetDescriptor(viewer,&fd);
270:   PetscMalloc((4+A->m)*sizeof(PetscInt),&col_lens);
271:   col_lens[0] = MAT_FILE_COOKIE;
272:   col_lens[1] = A->m;
273:   col_lens[2] = A->n;
274:   col_lens[3] = a->nz;

276:   /* store lengths of each row and write (including header) to file */
277:   for (i=0; i<A->m; i++) {
278:     col_lens[4+i] = a->i[i+1] - a->i[i];
279:   }
280:   PetscBinaryWrite(fd,col_lens,4+A->m,PETSC_INT,PETSC_TRUE);
281:   PetscFree(col_lens);

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

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

291: EXTERN PetscErrorCode MatSeqAIJFactorInfo_Matlab(Mat,PetscViewer);

295: PetscErrorCode MatView_SeqAIJ_ASCII(Mat A,PetscViewer viewer)
296: {
297:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
298:   PetscErrorCode    ierr;
299:   PetscInt          i,j,m = A->m,shift=0;
300:   char              *name;
301:   PetscViewerFormat format;

304:   PetscObjectGetName((PetscObject)A,&name);
305:   PetscViewerGetFormat(viewer,&format);
306:   if (format == PETSC_VIEWER_ASCII_INFO_DETAIL || format == PETSC_VIEWER_ASCII_INFO) {
307:     if (a->inode.size) {
308:       PetscViewerASCIIPrintf(viewer,"using I-node routines: found %D nodes, limit used is %D\n",a->inode.node_count,a->inode.limit);
309:     } else {
310:       PetscViewerASCIIPrintf(viewer,"not using I-node routines\n");
311:     }
312:   } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
313:     PetscInt nofinalvalue = 0;
314:     if ((a->i[m] == a->i[m-1]) || (a->j[a->nz-1] != A->n-!shift)) {
315:       nofinalvalue = 1;
316:     }
317:     PetscViewerASCIIUseTabs(viewer,PETSC_NO);
318:     PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",m,A->n);
319:     PetscViewerASCIIPrintf(viewer,"%% Nonzeros = %D \n",a->nz);
320:     PetscViewerASCIIPrintf(viewer,"zzz = zeros(%D,3);\n",a->nz+nofinalvalue);
321:     PetscViewerASCIIPrintf(viewer,"zzz = [\n");

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

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

459: PetscErrorCode MatView_SeqAIJ_Draw_Zoom(PetscDraw draw,void *Aa)
460: {
461:   Mat               A = (Mat) Aa;
462:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
463:   PetscErrorCode    ierr;
464:   PetscInt          i,j,m = A->m,color;
465:   PetscReal         xl,yl,xr,yr,x_l,x_r,y_l,y_r,maxv = 0.0;
466:   PetscViewer       viewer;
467:   PetscViewerFormat format;

470:   PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);
471:   PetscViewerGetFormat(viewer,&format);

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

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

520:     for (i=0; i<nz; i++) {
521:       if (PetscAbsScalar(a->a[i]) > maxv) maxv = PetscAbsScalar(a->a[i]);
522:     }
523:     scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv;
524:     PetscDrawGetPopup(draw,&popup);
525:     if (popup) {PetscDrawScalePopup(popup,0.0,maxv);}
526:     count = 0;
527:     for (i=0; i<m; i++) {
528:       y_l = m - i - 1.0; y_r = y_l + 1.0;
529:       for (j=a->i[i]; j<a->i[i+1]; j++) {
530:         x_l = a->j[j]; x_r = x_l + 1.0;
531:         color = PETSC_DRAW_BASIC_COLORS + (PetscInt)(scale*PetscAbsScalar(a->a[count]));
532:         PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
533:         count++;
534:       }
535:     }
536:   }
537:   return(0);
538: }

542: PetscErrorCode MatView_SeqAIJ_Draw(Mat A,PetscViewer viewer)
543: {
545:   PetscDraw      draw;
546:   PetscReal      xr,yr,xl,yl,h,w;
547:   PetscTruth     isnull;

550:   PetscViewerDrawGetDraw(viewer,0,&draw);
551:   PetscDrawIsNull(draw,&isnull);
552:   if (isnull) return(0);

554:   PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);
555:   xr  = A->n; yr = A->m; h = yr/10.0; w = xr/10.0;
556:   xr += w;    yr += h;  xl = -w;     yl = -h;
557:   PetscDrawSetCoordinates(draw,xl,yl,xr,yr);
558:   PetscDrawZoom(draw,MatView_SeqAIJ_Draw_Zoom,A);
559:   PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);
560:   return(0);
561: }

565: PetscErrorCode MatView_SeqAIJ(Mat A,PetscViewer viewer)
566: {
567:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
569:   PetscTruth     issocket,iascii,isbinary,isdraw;

572:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);
573:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
574:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
575:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
576:   if (issocket) {
577:     PetscViewerSocketPutSparse_Private(viewer,A->m,A->n,a->nz,a->a,a->i,a->j);
578:   } else if (iascii) {
579:     MatView_SeqAIJ_ASCII(A,viewer);
580:   } else if (isbinary) {
581:     MatView_SeqAIJ_Binary(A,viewer);
582:   } else if (isdraw) {
583:     MatView_SeqAIJ_Draw(A,viewer);
584:   } else {
585:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by SeqAIJ matrices",((PetscObject)viewer)->type_name);
586:   }
587:   return(0);
588: }

592: PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode)
593: {
594:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
596:   PetscInt       fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
597:   PetscInt       m = A->m,*ip,N,*ailen = a->ilen,rmax = 0;
598:   PetscScalar    *aa = a->a,*ap;

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

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

629:   /* diagonals may have moved, so kill the diagonal pointers */
630:   if (fshift && a->diag) {
631:     PetscFree(a->diag);
632:     PetscLogObjectMemory(A,-(m+1)*sizeof(PetscInt));
633:     a->diag = 0;
634:   }
635:   PetscLogInfo(A,"MatAssemblyEnd_SeqAIJ:Matrix size: %D X %D; storage space: %D unneeded,%D used\n",m,A->n,fshift,a->nz);
636:   PetscLogInfo(A,"MatAssemblyEnd_SeqAIJ:Number of mallocs during MatSetValues() is %D\n",a->reallocs);
637:   PetscLogInfo(A,"MatAssemblyEnd_SeqAIJ:Maximum nonzeros in any row is %D\n",rmax);
638:   a->reallocs          = 0;
639:   A->info.nz_unneeded  = (double)fshift;
640:   a->rmax              = rmax;

642:   /* check out for identical nodes. If found, use inode functions */
643:   Mat_AIJ_CheckInode(A,(PetscTruth)(!fshift));

645:   return(0);
646: }

650: PetscErrorCode MatZeroEntries_SeqAIJ(Mat A)
651: {
652:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;

656:   PetscMemzero(a->a,(a->i[A->m])*sizeof(PetscScalar));
657:   return(0);
658: }

662: PetscErrorCode MatDestroy_SeqAIJ(Mat A)
663: {
664:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;

668: #if defined(PETSC_USE_LOG)
669:   PetscLogObjectState((PetscObject)A,"Rows=%D, Cols=%D, NZ=%D",A->m,A->n,a->nz);
670: #endif
671:   if (a->freedata) {
672:     PetscFree(a->a);
673:     if (!a->singlemalloc) {
674:       PetscFree(a->i);
675:       PetscFree(a->j);
676:     }
677:   }
678:   if (a->row) {
679:     ISDestroy(a->row);
680:   }
681:   if (a->col) {
682:     ISDestroy(a->col);
683:   }
684:   if (a->diag) {PetscFree(a->diag);}
685:   if (a->ilen) {PetscFree(a->ilen);}
686:   if (a->imax) {PetscFree(a->imax);}
687:   if (a->idiag) {PetscFree(a->idiag);}
688:   if (a->solve_work) {PetscFree(a->solve_work);}
689:   if (a->inode.size) {PetscFree(a->inode.size);}
690:   if (a->icol) {ISDestroy(a->icol);}
691:   if (a->saved_values) {PetscFree(a->saved_values);}
692:   if (a->coloring) {ISColoringDestroy(a->coloring);}
693:   if (a->xtoy) {PetscFree(a->xtoy);}
694: 
695:   PetscFree(a);

697:   PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSeqAIJSetColumnIndices_C","",PETSC_NULL);
698:   PetscObjectComposeFunctionDynamic((PetscObject)A,"MatStoreValues_C","",PETSC_NULL);
699:   PetscObjectComposeFunctionDynamic((PetscObject)A,"MatRetrieveValues_C","",PETSC_NULL);
700:   PetscObjectComposeFunctionDynamic((PetscObject)A,"MatConvert_seqaij_seqsbaij_C","",PETSC_NULL);
701:   PetscObjectComposeFunctionDynamic((PetscObject)A,"MatConvert_seqaij_seqbaij_C","",PETSC_NULL);
702:   PetscObjectComposeFunctionDynamic((PetscObject)A,"MatIsTranspose_C","",PETSC_NULL);
703:   PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSeqAIJSetPreallocation_C","",PETSC_NULL);
704:   PetscObjectComposeFunctionDynamic((PetscObject)A,"MatReorderForNonzeroDiagonal_C","",PETSC_NULL);
705:   PetscObjectComposeFunctionDynamic((PetscObject)A,"MatAdjustForInodes_C","",PETSC_NULL);
706:   PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSeqAIJGetInodeSizes_C","",PETSC_NULL);
707:   return(0);
708: }

712: PetscErrorCode MatCompress_SeqAIJ(Mat A)
713: {
715:   return(0);
716: }

720: PetscErrorCode MatSetOption_SeqAIJ(Mat A,MatOption op)
721: {
722:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;

725:   switch (op) {
726:     case MAT_ROW_ORIENTED:
727:       a->roworiented       = PETSC_TRUE;
728:       break;
729:     case MAT_KEEP_ZEROED_ROWS:
730:       a->keepzeroedrows    = PETSC_TRUE;
731:       break;
732:     case MAT_COLUMN_ORIENTED:
733:       a->roworiented       = PETSC_FALSE;
734:       break;
735:     case MAT_COLUMNS_SORTED:
736:       a->sorted            = PETSC_TRUE;
737:       break;
738:     case MAT_COLUMNS_UNSORTED:
739:       a->sorted            = PETSC_FALSE;
740:       break;
741:     case MAT_NO_NEW_NONZERO_LOCATIONS:
742:       a->nonew             = 1;
743:       break;
744:     case MAT_NEW_NONZERO_LOCATION_ERR:
745:       a->nonew             = -1;
746:       break;
747:     case MAT_NEW_NONZERO_ALLOCATION_ERR:
748:       a->nonew             = -2;
749:       break;
750:     case MAT_YES_NEW_NONZERO_LOCATIONS:
751:       a->nonew             = 0;
752:       break;
753:     case MAT_IGNORE_ZERO_ENTRIES:
754:       a->ignorezeroentries = PETSC_TRUE;
755:       break;
756:     case MAT_USE_INODES:
757:       a->inode.use         = PETSC_TRUE;
758:       break;
759:     case MAT_DO_NOT_USE_INODES:
760:       a->inode.use         = PETSC_FALSE;
761:       break;
762:     case MAT_ROWS_SORTED:
763:     case MAT_ROWS_UNSORTED:
764:     case MAT_YES_NEW_DIAGONALS:
765:     case MAT_IGNORE_OFF_PROC_ENTRIES:
766:     case MAT_USE_HASH_TABLE:
767:       PetscLogInfo(A,"MatSetOption_SeqAIJ:Option ignored\n");
768:       break;
769:     case MAT_NO_NEW_DIAGONALS:
770:       SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
771:     case MAT_INODE_LIMIT_1:
772:       a->inode.limit  = 1;
773:       break;
774:     case MAT_INODE_LIMIT_2:
775:       a->inode.limit  = 2;
776:       break;
777:     case MAT_INODE_LIMIT_3:
778:       a->inode.limit  = 3;
779:       break;
780:     case MAT_INODE_LIMIT_4:
781:       a->inode.limit  = 4;
782:       break;
783:     case MAT_INODE_LIMIT_5:
784:       a->inode.limit  = 5;
785:       break;
786:     case MAT_SYMMETRIC:
787:     case MAT_STRUCTURALLY_SYMMETRIC:
788:     case MAT_NOT_SYMMETRIC:
789:     case MAT_NOT_STRUCTURALLY_SYMMETRIC:
790:     case MAT_HERMITIAN:
791:     case MAT_NOT_HERMITIAN:
792:     case MAT_SYMMETRY_ETERNAL:
793:     case MAT_NOT_SYMMETRY_ETERNAL:
794:       break;
795:     default:
796:       SETERRQ(PETSC_ERR_SUP,"unknown option");
797:   }
798:   return(0);
799: }

803: PetscErrorCode MatGetDiagonal_SeqAIJ(Mat A,Vec v)
804: {
805:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
807:   PetscInt       i,j,n;
808:   PetscScalar    *x,zero = 0.0;

811:   VecSet(&zero,v);
812:   VecGetArray(v,&x);
813:   VecGetLocalSize(v,&n);
814:   if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
815:   for (i=0; i<A->m; i++) {
816:     for (j=a->i[i]; j<a->i[i+1]; j++) {
817:       if (a->j[j] == i) {
818:         x[i] = a->a[j];
819:         break;
820:       }
821:     }
822:   }
823:   VecRestoreArray(v,&x);
824:   return(0);
825: }


830: PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
831: {
832:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
833:   PetscScalar    *x,*y;
835:   PetscInt       m = A->m;
836: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
837:   PetscScalar    *v,alpha;
838:   PetscInt       n,i,*idx;
839: #endif

842:   if (zz != yy) {VecCopy(zz,yy);}
843:   VecGetArray(xx,&x);
844:   VecGetArray(yy,&y);

846: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
847:   fortranmulttransposeaddaij_(&m,x,a->i,a->j,a->a,y);
848: #else
849:   for (i=0; i<m; i++) {
850:     idx   = a->j + a->i[i] ;
851:     v     = a->a + a->i[i] ;
852:     n     = a->i[i+1] - a->i[i];
853:     alpha = x[i];
854:     while (n-->0) {y[*idx++] += alpha * *v++;}
855:   }
856: #endif
857:   PetscLogFlops(2*a->nz);
858:   VecRestoreArray(xx,&x);
859:   VecRestoreArray(yy,&y);
860:   return(0);
861: }

865: PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec xx,Vec yy)
866: {
867:   PetscScalar    zero = 0.0;

871:   VecSet(&zero,yy);
872:   MatMultTransposeAdd_SeqAIJ(A,xx,yy,yy);
873:   return(0);
874: }


879: PetscErrorCode MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
880: {
881:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
882:   PetscScalar    *x,*y,*v;
884:   PetscInt       m = A->m,*idx,*ii;
885: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
886:   PetscInt       n,i,jrow,j;
887:   PetscScalar    sum;
888: #endif

890: #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
891: #pragma disjoint(*x,*y,*v)
892: #endif

895:   VecGetArray(xx,&x);
896:   VecGetArray(yy,&y);
897:   idx  = a->j;
898:   v    = a->a;
899:   ii   = a->i;
900: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
901:   fortranmultaij_(&m,x,ii,idx,v,y);
902: #else
903:   for (i=0; i<m; i++) {
904:     jrow = ii[i];
905:     n    = ii[i+1] - jrow;
906:     sum  = 0.0;
907:     for (j=0; j<n; j++) {
908:       sum += v[jrow]*x[idx[jrow]]; jrow++;
909:      }
910:     y[i] = sum;
911:   }
912: #endif
913:   PetscLogFlops(2*a->nz - m);
914:   VecRestoreArray(xx,&x);
915:   VecRestoreArray(yy,&y);
916:   return(0);
917: }

921: PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
922: {
923:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
924:   PetscScalar    *x,*y,*z,*v;
926:   PetscInt       m = A->m,*idx,*ii;
927: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
928:   PetscInt       n,i,jrow,j;
929: PetscScalar      sum;
930: #endif

933:   VecGetArray(xx,&x);
934:   VecGetArray(yy,&y);
935:   if (zz != yy) {
936:     VecGetArray(zz,&z);
937:   } else {
938:     z = y;
939:   }
940: 
941:   idx  = a->j;
942:   v    = a->a;
943:   ii   = a->i;
944: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
945:   fortranmultaddaij_(&m,x,ii,idx,v,y,z);
946: #else
947:   for (i=0; i<m; i++) {
948:     jrow = ii[i];
949:     n    = ii[i+1] - jrow;
950:     sum  = y[i];
951:     for (j=0; j<n; j++) {
952:       sum += v[jrow]*x[idx[jrow]]; jrow++;
953:      }
954:     z[i] = sum;
955:   }
956: #endif
957:   PetscLogFlops(2*a->nz);
958:   VecRestoreArray(xx,&x);
959:   VecRestoreArray(yy,&y);
960:   if (zz != yy) {
961:     VecRestoreArray(zz,&z);
962:   }
963:   return(0);
964: }

966: /*
967:      Adds diagonal pointers to sparse matrix structure.
968: */
971: PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat A)
972: {
973:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
975:   PetscInt       i,j,*diag,m = A->m;

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

980:   PetscMalloc((m+1)*sizeof(PetscInt),&diag);
981:   PetscLogObjectMemory(A,(m+1)*sizeof(PetscInt));
982:   for (i=0; i<A->m; i++) {
983:     diag[i] = a->i[i+1];
984:     for (j=a->i[i]; j<a->i[i+1]; j++) {
985:       if (a->j[j] == i) {
986:         diag[i] = j;
987:         break;
988:       }
989:     }
990:   }
991:   a->diag = diag;
992:   return(0);
993: }

995: /*
996:      Checks for missing diagonals
997: */
1000: PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat A)
1001: {
1002:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1004:   PetscInt       *diag,*jj = a->j,i;

1007:   MatMarkDiagonal_SeqAIJ(A);
1008:   diag = a->diag;
1009:   for (i=0; i<A->m; i++) {
1010:     if (jj[diag[i]] != i) {
1011:       SETERRQ1(PETSC_ERR_PLIB,"Matrix is missing diagonal number %D",i);
1012:     }
1013:   }
1014:   return(0);
1015: }

1019: PetscErrorCode MatRelax_SeqAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
1020: {
1021:   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data;
1022:   PetscScalar        *x,d,*xs,sum,*t,scale,*idiag=0,*mdiag;
1023:   const PetscScalar  *v = a->a, *b, *bs,*xb, *ts;
1024:   PetscErrorCode     ierr;
1025:   PetscInt           n = A->n,m = A->m,i;
1026:   const PetscInt     *idx,*diag;

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

1032:   if (!a->diag) {MatMarkDiagonal_SeqAIJ(A);}
1033:   diag = a->diag;
1034:   if (!a->idiag) {
1035:     PetscMalloc(3*m*sizeof(PetscScalar),&a->idiag);
1036:     a->ssor  = a->idiag + m;
1037:     mdiag    = a->ssor + m;

1039:     v        = a->a;

1041:     /* this is wrong when fshift omega changes each iteration */
1042:     if (omega == 1.0 && !fshift) {
1043:       for (i=0; i<m; i++) {
1044:         mdiag[i]    = v[diag[i]];
1045:         a->idiag[i] = 1.0/v[diag[i]];
1046:       }
1047:       PetscLogFlops(m);
1048:     } else {
1049:       for (i=0; i<m; i++) {
1050:         mdiag[i]    = v[diag[i]];
1051:         a->idiag[i] = omega/(fshift + v[diag[i]]);
1052:       }
1053:       PetscLogFlops(2*m);
1054:     }
1055:   }
1056:   t     = a->ssor;
1057:   idiag = a->idiag;
1058:   mdiag = a->idiag + 2*m;

1060:   VecGetArray(xx,&x);
1061:   if (xx != bb) {
1062:     VecGetArray(bb,(PetscScalar**)&b);
1063:   } else {
1064:     b = x;
1065:   }

1067:   /* We count flops by assuming the upper triangular and lower triangular parts have the same number of nonzeros */
1068:   xs   = x;
1069:   if (flag == SOR_APPLY_UPPER) {
1070:    /* apply (U + D/omega) to the vector */
1071:     bs = b;
1072:     for (i=0; i<m; i++) {
1073:         d    = fshift + a->a[diag[i]];
1074:         n    = a->i[i+1] - diag[i] - 1;
1075:         idx  = a->j + diag[i] + 1;
1076:         v    = a->a + diag[i] + 1;
1077:         sum  = b[i]*d/omega;
1078:         SPARSEDENSEDOT(sum,bs,v,idx,n);
1079:         x[i] = sum;
1080:     }
1081:     VecRestoreArray(xx,&x);
1082:     if (bb != xx) {VecRestoreArray(bb,(PetscScalar**)&b);}
1083:     PetscLogFlops(a->nz);
1084:     return(0);
1085:   }


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

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

1093:     to a vector efficiently using Eisenstat's trick. This is for
1094:     the case of SSOR preconditioner, so E is D/omega where omega
1095:     is the relaxation factor.
1096:     */

1098:   if (flag == SOR_APPLY_LOWER) {
1099:     SETERRQ(PETSC_ERR_SUP,"SOR_APPLY_LOWER is not implemented");
1100:   } else if (flag & SOR_EISENSTAT) {
1101:     /* Let  A = L + U + D; where L is lower trianglar,
1102:     U is upper triangular, E is diagonal; This routine applies

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

1106:     to a vector efficiently using Eisenstat's trick. This is for
1107:     the case of SSOR preconditioner, so E is D/omega where omega
1108:     is the relaxation factor.
1109:     */
1110:     scale = (2.0/omega) - 1.0;

1112:     /*  x = (E + U)^{-1} b */
1113:     for (i=m-1; i>=0; i--) {
1114:       n    = a->i[i+1] - diag[i] - 1;
1115:       idx  = a->j + diag[i] + 1;
1116:       v    = a->a + diag[i] + 1;
1117:       sum  = b[i];
1118:       SPARSEDENSEMDOT(sum,xs,v,idx,n);
1119:       x[i] = sum*idiag[i];
1120:     }

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

1126:     /*  t = (E + L)^{-1}t */
1127:     ts = t;
1128:     diag = a->diag;
1129:     for (i=0; i<m; i++) {
1130:       n    = diag[i] - a->i[i];
1131:       idx  = a->j + a->i[i];
1132:       v    = a->a + a->i[i];
1133:       sum  = t[i];
1134:       SPARSEDENSEMDOT(sum,ts,v,idx,n);
1135:       t[i] = sum*idiag[i];
1136:       /*  x = x + t */
1137:       x[i] += t[i];
1138:     }

1140:     PetscLogFlops(6*m-1 + 2*a->nz);
1141:     VecRestoreArray(xx,&x);
1142:     if (bb != xx) {VecRestoreArray(bb,(PetscScalar**)&b);}
1143:     return(0);
1144:   }
1145:   if (flag & SOR_ZERO_INITIAL_GUESS) {
1146:     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
1147: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1148:       fortranrelaxaijforwardzero_(&m,&omega,x,a->i,a->j,(PetscInt*)diag,idiag,a->a,(void*)b);
1149: #else
1150:       for (i=0; i<m; i++) {
1151:         n    = diag[i] - a->i[i];
1152:         idx  = a->j + a->i[i];
1153:         v    = a->a + a->i[i];
1154:         sum  = b[i];
1155:         SPARSEDENSEMDOT(sum,xs,v,idx,n);
1156:         x[i] = sum*idiag[i];
1157:       }
1158: #endif
1159:       xb = x;
1160:       PetscLogFlops(a->nz);
1161:     } else xb = b;
1162:     if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) &&
1163:         (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) {
1164:       for (i=0; i<m; i++) {
1165:         x[i] *= mdiag[i];
1166:       }
1167:       PetscLogFlops(m);
1168:     }
1169:     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
1170: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1171:       fortranrelaxaijbackwardzero_(&m,&omega,x,a->i,a->j,(PetscInt*)diag,idiag,a->a,(void*)xb);
1172: #else
1173:       for (i=m-1; i>=0; i--) {
1174:         n    = a->i[i+1] - diag[i] - 1;
1175:         idx  = a->j + diag[i] + 1;
1176:         v    = a->a + diag[i] + 1;
1177:         sum  = xb[i];
1178:         SPARSEDENSEMDOT(sum,xs,v,idx,n);
1179:         x[i] = sum*idiag[i];
1180:       }
1181: #endif
1182:       PetscLogFlops(a->nz);
1183:     }
1184:     its--;
1185:   }
1186:   while (its--) {
1187:     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
1188: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1189:       fortranrelaxaijforward_(&m,&omega,x,a->i,a->j,(PetscInt*)diag,a->a,(void*)b);
1190: #else
1191:       for (i=0; i<m; i++) {
1192:         d    = fshift + a->a[diag[i]];
1193:         n    = a->i[i+1] - a->i[i];
1194:         idx  = a->j + a->i[i];
1195:         v    = a->a + a->i[i];
1196:         sum  = b[i];
1197:         SPARSEDENSEMDOT(sum,xs,v,idx,n);
1198:         x[i] = (1. - omega)*x[i] + (sum + mdiag[i]*x[i])*idiag[i];
1199:       }
1200: #endif 
1201:       PetscLogFlops(a->nz);
1202:     }
1203:     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
1204: #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ)
1205:       fortranrelaxaijbackward_(&m,&omega,x,a->i,a->j,(PetscInt*)diag,a->a,(void*)b);
1206: #else
1207:       for (i=m-1; i>=0; i--) {
1208:         d    = fshift + a->a[diag[i]];
1209:         n    = a->i[i+1] - a->i[i];
1210:         idx  = a->j + a->i[i];
1211:         v    = a->a + a->i[i];
1212:         sum  = b[i];
1213:         SPARSEDENSEMDOT(sum,xs,v,idx,n);
1214:         x[i] = (1. - omega)*x[i] + (sum + mdiag[i]*x[i])*idiag[i];
1215:       }
1216: #endif
1217:       PetscLogFlops(a->nz);
1218:     }
1219:   }
1220:   VecRestoreArray(xx,&x);
1221:   if (bb != xx) {VecRestoreArray(bb,(PetscScalar**)&b);}
1222:   return(0);
1223: }

1227: PetscErrorCode MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,MatInfo *info)
1228: {
1229:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;

1232:   info->rows_global    = (double)A->m;
1233:   info->columns_global = (double)A->n;
1234:   info->rows_local     = (double)A->m;
1235:   info->columns_local  = (double)A->n;
1236:   info->block_size     = 1.0;
1237:   info->nz_allocated   = (double)a->maxnz;
1238:   info->nz_used        = (double)a->nz;
1239:   info->nz_unneeded    = (double)(a->maxnz - a->nz);
1240:   info->assemblies     = (double)A->num_ass;
1241:   info->mallocs        = (double)a->reallocs;
1242:   info->memory         = A->mem;
1243:   if (A->factor) {
1244:     info->fill_ratio_given  = A->info.fill_ratio_given;
1245:     info->fill_ratio_needed = A->info.fill_ratio_needed;
1246:     info->factor_mallocs    = A->info.factor_mallocs;
1247:   } else {
1248:     info->fill_ratio_given  = 0;
1249:     info->fill_ratio_needed = 0;
1250:     info->factor_mallocs    = 0;
1251:   }
1252:   return(0);
1253: }

1257: PetscErrorCode MatZeroRows_SeqAIJ(Mat A,IS is,const PetscScalar *diag)
1258: {
1259:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1261:   PetscInt       i,N,*rows,m = A->m - 1;

1264:   ISGetLocalSize(is,&N);
1265:   ISGetIndices(is,&rows);
1266:   if (a->keepzeroedrows) {
1267:     for (i=0; i<N; i++) {
1268:       if (rows[i] < 0 || rows[i] > m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range", rows[i]);
1269:       PetscMemzero(&a->a[a->i[rows[i]]],a->ilen[rows[i]]*sizeof(PetscScalar));
1270:     }
1271:     if (diag) {
1272:       MatMissingDiagonal_SeqAIJ(A);
1273:       MatMarkDiagonal_SeqAIJ(A);
1274:       for (i=0; i<N; i++) {
1275:         a->a[a->diag[rows[i]]] = *diag;
1276:       }
1277:     }
1278:   } else {
1279:     if (diag) {
1280:       for (i=0; i<N; i++) {
1281:         if (rows[i] < 0 || rows[i] > m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range", rows[i]);
1282:         if (a->ilen[rows[i]] > 0) {
1283:           a->ilen[rows[i]]          = 1;
1284:           a->a[a->i[rows[i]]] = *diag;
1285:           a->j[a->i[rows[i]]] = rows[i];
1286:         } else { /* in case row was completely empty */
1287:           MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES);
1288:         }
1289:       }
1290:     } else {
1291:       for (i=0; i<N; i++) {
1292:         if (rows[i] < 0 || rows[i] > m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range", rows[i]);
1293:         a->ilen[rows[i]] = 0;
1294:       }
1295:     }
1296:   }
1297:   ISRestoreIndices(is,&rows);
1298:   MatAssemblyEnd_SeqAIJ(A,MAT_FINAL_ASSEMBLY);
1299:   return(0);
1300: }

1304: PetscErrorCode MatGetRow_SeqAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
1305: {
1306:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1307:   PetscInt   *itmp;

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

1312:   *nz = a->i[row+1] - a->i[row];
1313:   if (v) *v = a->a + a->i[row];
1314:   if (idx) {
1315:     itmp = a->j + a->i[row];
1316:     if (*nz) {
1317:       *idx = itmp;
1318:     }
1319:     else *idx = 0;
1320:   }
1321:   return(0);
1322: }

1324: /* remove this function? */
1327: PetscErrorCode MatRestoreRow_SeqAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
1328: {
1330:   return(0);
1331: }

1335: PetscErrorCode MatNorm_SeqAIJ(Mat A,NormType type,PetscReal *nrm)
1336: {
1337:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1338:   PetscScalar    *v = a->a;
1339:   PetscReal      sum = 0.0;
1341:   PetscInt       i,j;

1344:   if (type == NORM_FROBENIUS) {
1345:     for (i=0; i<a->nz; i++) {
1346: #if defined(PETSC_USE_COMPLEX)
1347:       sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
1348: #else
1349:       sum += (*v)*(*v); v++;
1350: #endif
1351:     }
1352:     *nrm = sqrt(sum);
1353:   } else if (type == NORM_1) {
1354:     PetscReal *tmp;
1355:     PetscInt    *jj = a->j;
1356:     PetscMalloc((A->n+1)*sizeof(PetscReal),&tmp);
1357:     PetscMemzero(tmp,A->n*sizeof(PetscReal));
1358:     *nrm = 0.0;
1359:     for (j=0; j<a->nz; j++) {
1360:         tmp[*jj++] += PetscAbsScalar(*v);  v++;
1361:     }
1362:     for (j=0; j<A->n; j++) {
1363:       if (tmp[j] > *nrm) *nrm = tmp[j];
1364:     }
1365:     PetscFree(tmp);
1366:   } else if (type == NORM_INFINITY) {
1367:     *nrm = 0.0;
1368:     for (j=0; j<A->m; j++) {
1369:       v = a->a + a->i[j];
1370:       sum = 0.0;
1371:       for (i=0; i<a->i[j+1]-a->i[j]; i++) {
1372:         sum += PetscAbsScalar(*v); v++;
1373:       }
1374:       if (sum > *nrm) *nrm = sum;
1375:     }
1376:   } else {
1377:     SETERRQ(PETSC_ERR_SUP,"No support for two norm");
1378:   }
1379:   return(0);
1380: }

1384: PetscErrorCode MatTranspose_SeqAIJ(Mat A,Mat *B)
1385: {
1386:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1387:   Mat            C;
1389:   PetscInt       i,*aj = a->j,*ai = a->i,m = A->m,len,*col;
1390:   PetscScalar    *array = a->a;

1393:   if (!B && m != A->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Square matrix only for in-place");
1394:   PetscMalloc((1+A->n)*sizeof(PetscInt),&col);
1395:   PetscMemzero(col,(1+A->n)*sizeof(PetscInt));
1396: 
1397:   for (i=0; i<ai[m]; i++) col[aj[i]] += 1;
1398:   MatCreate(A->comm,A->n,m,A->n,m,&C);
1399:   MatSetType(C,A->type_name);
1400:   MatSeqAIJSetPreallocation(C,0,col);
1401:   PetscFree(col);
1402:   for (i=0; i<m; i++) {
1403:     len    = ai[i+1]-ai[i];
1404:     MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES);
1405:     array += len;
1406:     aj    += len;
1407:   }

1409:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1410:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

1412:   if (B) {
1413:     *B = C;
1414:   } else {
1415:     MatHeaderCopy(A,C);
1416:   }
1417:   return(0);
1418: }

1423: PetscErrorCode MatIsTranspose_SeqAIJ(Mat A,Mat B,PetscReal tol,PetscTruth *f)
1424: {
1425:   Mat_SeqAIJ     *aij = (Mat_SeqAIJ *) A->data,*bij = (Mat_SeqAIJ*) A->data;
1426:   PetscInt       *adx,*bdx,*aii,*bii,*aptr,*bptr; PetscScalar *va,*vb;
1428:   PetscInt       ma,na,mb,nb, i;

1431:   bij = (Mat_SeqAIJ *) B->data;
1432: 
1433:   MatGetSize(A,&ma,&na);
1434:   MatGetSize(B,&mb,&nb);
1435:   if (ma!=nb || na!=mb){
1436:     *f = PETSC_FALSE;
1437:     return(0);
1438:   }
1439:   aii = aij->i; bii = bij->i;
1440:   adx = aij->j; bdx = bij->j;
1441:   va  = aij->a; vb = bij->a;
1442:   PetscMalloc(ma*sizeof(PetscInt),&aptr);
1443:   PetscMalloc(mb*sizeof(PetscInt),&bptr);
1444:   for (i=0; i<ma; i++) aptr[i] = aii[i];
1445:   for (i=0; i<mb; i++) bptr[i] = bii[i];

1447:   *f = PETSC_TRUE;
1448:   for (i=0; i<ma; i++) {
1449:     while (aptr[i]<aii[i+1]) {
1450:       PetscInt         idc,idr;
1451:       PetscScalar vc,vr;
1452:       /* column/row index/value */
1453:       idc = adx[aptr[i]];
1454:       idr = bdx[bptr[idc]];
1455:       vc  = va[aptr[i]];
1456:       vr  = vb[bptr[idc]];
1457:       if (i!=idr || PetscAbsScalar(vc-vr) > tol) {
1458:         *f = PETSC_FALSE;
1459:         goto done;
1460:       } else {
1461:         aptr[i]++;
1462:         if (B || i!=idc) bptr[idc]++;
1463:       }
1464:     }
1465:   }
1466:  done:
1467:   PetscFree(aptr);
1468:   if (B) {
1469:     PetscFree(bptr);
1470:   }
1471:   return(0);
1472: }

1477: PetscErrorCode MatIsSymmetric_SeqAIJ(Mat A,PetscReal tol,PetscTruth *f)
1478: {
1481:   MatIsTranspose_SeqAIJ(A,A,tol,f);
1482:   return(0);
1483: }

1487: PetscErrorCode MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr)
1488: {
1489:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1490:   PetscScalar    *l,*r,x,*v;
1492:   PetscInt       i,j,m = A->m,n = A->n,M,nz = a->nz,*jj;

1495:   if (ll) {
1496:     /* The local size is used so that VecMPI can be passed to this routine
1497:        by MatDiagonalScale_MPIAIJ */
1498:     VecGetLocalSize(ll,&m);
1499:     if (m != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length");
1500:     VecGetArray(ll,&l);
1501:     v = a->a;
1502:     for (i=0; i<m; i++) {
1503:       x = l[i];
1504:       M = a->i[i+1] - a->i[i];
1505:       for (j=0; j<M; j++) { (*v++) *= x;}
1506:     }
1507:     VecRestoreArray(ll,&l);
1508:     PetscLogFlops(nz);
1509:   }
1510:   if (rr) {
1511:     VecGetLocalSize(rr,&n);
1512:     if (n != A->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vector wrong length");
1513:     VecGetArray(rr,&r);
1514:     v = a->a; jj = a->j;
1515:     for (i=0; i<nz; i++) {
1516:       (*v++) *= r[*jj++];
1517:     }
1518:     VecRestoreArray(rr,&r);
1519:     PetscLogFlops(nz);
1520:   }
1521:   return(0);
1522: }

1526: PetscErrorCode MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,PetscInt csize,MatReuse scall,Mat *B)
1527: {
1528:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data,*c;
1530:   PetscInt       *smap,i,k,kstart,kend,oldcols = A->n,*lens;
1531:   PetscInt       row,mat_i,*mat_j,tcol,first,step,*mat_ilen,sum,lensi;
1532:   PetscInt       *irow,*icol,nrows,ncols;
1533:   PetscInt       *starts,*j_new,*i_new,*aj = a->j,*ai = a->i,ii,*ailen = a->ilen;
1534:   PetscScalar    *a_new,*mat_a;
1535:   Mat            C;
1536:   PetscTruth     stride;

1539:   ISSorted(isrow,(PetscTruth*)&i);
1540:   if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"ISrow is not sorted");
1541:   ISSorted(iscol,(PetscTruth*)&i);
1542:   if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IScol is not sorted");

1544:   ISGetIndices(isrow,&irow);
1545:   ISGetLocalSize(isrow,&nrows);
1546:   ISGetLocalSize(iscol,&ncols);

1548:   ISStrideGetInfo(iscol,&first,&step);
1549:   ISStride(iscol,&stride);
1550:   if (stride && step == 1) {
1551:     /* special case of contiguous rows */
1552:     PetscMalloc((2*nrows+1)*sizeof(PetscInt),&lens);
1553:     starts = lens + nrows;
1554:     /* loop over new rows determining lens and starting points */
1555:     for (i=0; i<nrows; i++) {
1556:       kstart  = ai[irow[i]];
1557:       kend    = kstart + ailen[irow[i]];
1558:       for (k=kstart; k<kend; k++) {
1559:         if (aj[k] >= first) {
1560:           starts[i] = k;
1561:           break;
1562:         }
1563:       }
1564:       sum = 0;
1565:       while (k < kend) {
1566:         if (aj[k++] >= first+ncols) break;
1567:         sum++;
1568:       }
1569:       lens[i] = sum;
1570:     }
1571:     /* create submatrix */
1572:     if (scall == MAT_REUSE_MATRIX) {
1573:       PetscInt n_cols,n_rows;
1574:       MatGetSize(*B,&n_rows,&n_cols);
1575:       if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size");
1576:       MatZeroEntries(*B);
1577:       C = *B;
1578:     } else {
1579:       MatCreate(A->comm,nrows,ncols,PETSC_DETERMINE,PETSC_DETERMINE,&C);
1580:       MatSetType(C,A->type_name);
1581:       MatSeqAIJSetPreallocation(C,0,lens);
1582:     }
1583:     c = (Mat_SeqAIJ*)C->data;

1585:     /* loop over rows inserting into submatrix */
1586:     a_new    = c->a;
1587:     j_new    = c->j;
1588:     i_new    = c->i;

1590:     for (i=0; i<nrows; i++) {
1591:       ii    = starts[i];
1592:       lensi = lens[i];
1593:       for (k=0; k<lensi; k++) {
1594:         *j_new++ = aj[ii+k] - first;
1595:       }
1596:       PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(PetscScalar));
1597:       a_new      += lensi;
1598:       i_new[i+1]  = i_new[i] + lensi;
1599:       c->ilen[i]  = lensi;
1600:     }
1601:     PetscFree(lens);
1602:   } else {
1603:     ISGetIndices(iscol,&icol);
1604:     PetscMalloc((1+oldcols)*sizeof(PetscInt),&smap);
1605: 
1606:     PetscMalloc((1+nrows)*sizeof(PetscInt),&lens);
1607:     PetscMemzero(smap,oldcols*sizeof(PetscInt));
1608:     for (i=0; i<ncols; i++) smap[icol[i]] = i+1;
1609:     /* determine lens of each row */
1610:     for (i=0; i<nrows; i++) {
1611:       kstart  = ai[irow[i]];
1612:       kend    = kstart + a->ilen[irow[i]];
1613:       lens[i] = 0;
1614:       for (k=kstart; k<kend; k++) {
1615:         if (smap[aj[k]]) {
1616:           lens[i]++;
1617:         }
1618:       }
1619:     }
1620:     /* Create and fill new matrix */
1621:     if (scall == MAT_REUSE_MATRIX) {
1622:       PetscTruth equal;

1624:       c = (Mat_SeqAIJ *)((*B)->data);
1625:       if ((*B)->m  != nrows || (*B)->n != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size");
1626:       PetscMemcmp(c->ilen,lens,(*B)->m*sizeof(PetscInt),&equal);
1627:       if (!equal) {
1628:         SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros");
1629:       }
1630:       PetscMemzero(c->ilen,(*B)->m*sizeof(PetscInt));
1631:       C = *B;
1632:     } else {
1633:       MatCreate(A->comm,nrows,ncols,PETSC_DETERMINE,PETSC_DETERMINE,&C);
1634:       MatSetType(C,A->type_name);
1635:       MatSeqAIJSetPreallocation(C,0,lens);
1636:     }
1637:     c = (Mat_SeqAIJ *)(C->data);
1638:     for (i=0; i<nrows; i++) {
1639:       row    = irow[i];
1640:       kstart = ai[row];
1641:       kend   = kstart + a->ilen[row];
1642:       mat_i  = c->i[i];
1643:       mat_j  = c->j + mat_i;
1644:       mat_a  = c->a + mat_i;
1645:       mat_ilen = c->ilen + i;
1646:       for (k=kstart; k<kend; k++) {
1647:         if ((tcol=smap[a->j[k]])) {
1648:           *mat_j++ = tcol - 1;
1649:           *mat_a++ = a->a[k];
1650:           (*mat_ilen)++;

1652:         }
1653:       }
1654:     }
1655:     /* Free work space */
1656:     ISRestoreIndices(iscol,&icol);
1657:     PetscFree(smap);
1658:     PetscFree(lens);
1659:   }
1660:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1661:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

1663:   ISRestoreIndices(isrow,&irow);
1664:   *B = C;
1665:   return(0);
1666: }

1668: /*
1669: */
1672: PetscErrorCode MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,MatFactorInfo *info)
1673: {
1674:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)inA->data;
1676:   Mat            outA;
1677:   PetscTruth     row_identity,col_identity;

1680:   if (info->levels != 0) SETERRQ(PETSC_ERR_SUP,"Only levels=0 supported for in-place ilu");
1681:   ISIdentity(row,&row_identity);
1682:   ISIdentity(col,&col_identity);
1683:   if (!row_identity || !col_identity) {
1684:     SETERRQ(PETSC_ERR_ARG_WRONG,"Row and column permutations must be identity for in-place ILU");
1685:   }

1687:   outA          = inA;
1688:   inA->factor   = FACTOR_LU;
1689:   a->row        = row;
1690:   a->col        = col;
1691:   PetscObjectReference((PetscObject)row);
1692:   PetscObjectReference((PetscObject)col);

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

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

1703:   if (!a->diag) {
1704:     MatMarkDiagonal_SeqAIJ(inA);
1705:   }
1706:   MatLUFactorNumeric_SeqAIJ(inA,&outA);
1707:   return(0);
1708: }

1710:  #include petscblaslapack.h
1713: PetscErrorCode MatScale_SeqAIJ(const PetscScalar *alpha,Mat inA)
1714: {
1715:   Mat_SeqAIJ   *a = (Mat_SeqAIJ*)inA->data;
1716:   PetscBLASInt bnz = (PetscBLASInt)a->nz,one = 1;

1719:   BLscal_(&bnz,(PetscScalar*)alpha,a->a,&one);
1720:   PetscLogFlops(a->nz);
1721:   return(0);
1722: }

1726: PetscErrorCode MatGetSubMatrices_SeqAIJ(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1727: {
1729:   PetscInt       i;

1732:   if (scall == MAT_INITIAL_MATRIX) {
1733:     PetscMalloc((n+1)*sizeof(Mat),B);
1734:   }

1736:   for (i=0; i<n; i++) {
1737:     MatGetSubMatrix_SeqAIJ(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);
1738:   }
1739:   return(0);
1740: }

1744: PetscErrorCode MatIncreaseOverlap_SeqAIJ(Mat A,PetscInt is_max,IS is[],PetscInt ov)
1745: {
1746:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1748:   PetscInt       row,i,j,k,l,m,n,*idx,*nidx,isz,val;
1749:   PetscInt       start,end,*ai,*aj;
1750:   PetscBT        table;

1753:   m     = A->m;
1754:   ai    = a->i;
1755:   aj    = a->j;

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

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

1762:   for (i=0; i<is_max; i++) {
1763:     /* Initialize the two local arrays */
1764:     isz  = 0;
1765:     PetscBTMemzero(m,table);
1766: 
1767:     /* Extract the indices, assume there can be duplicate entries */
1768:     ISGetIndices(is[i],&idx);
1769:     ISGetLocalSize(is[i],&n);
1770: 
1771:     /* Enter these into the temp arrays. I.e., mark table[row], enter row into new index */
1772:     for (j=0; j<n ; ++j){
1773:       if(!PetscBTLookupSet(table,idx[j])) { nidx[isz++] = idx[j];}
1774:     }
1775:     ISRestoreIndices(is[i],&idx);
1776:     ISDestroy(is[i]);
1777: 
1778:     k = 0;
1779:     for (j=0; j<ov; j++){ /* for each overlap */
1780:       n = isz;
1781:       for (; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */
1782:         row   = nidx[k];
1783:         start = ai[row];
1784:         end   = ai[row+1];
1785:         for (l = start; l<end ; l++){
1786:           val = aj[l] ;
1787:           if (!PetscBTLookupSet(table,val)) {nidx[isz++] = val;}
1788:         }
1789:       }
1790:     }
1791:     ISCreateGeneral(PETSC_COMM_SELF,isz,nidx,(is+i));
1792:   }
1793:   PetscBTDestroy(table);
1794:   PetscFree(nidx);
1795:   return(0);
1796: }

1798: /* -------------------------------------------------------------- */
1801: PetscErrorCode MatPermute_SeqAIJ(Mat A,IS rowp,IS colp,Mat *B)
1802: {
1803:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1805:   PetscInt       i,nz,m = A->m,n = A->n,*col;
1806:   PetscInt       *row,*cnew,j,*lens;
1807:   IS             icolp,irowp;
1808:   PetscInt       *cwork;
1809:   PetscScalar    *vwork;

1812:   ISInvertPermutation(rowp,PETSC_DECIDE,&irowp);
1813:   ISGetIndices(irowp,&row);
1814:   ISInvertPermutation(colp,PETSC_DECIDE,&icolp);
1815:   ISGetIndices(icolp,&col);
1816: 
1817:   /* determine lengths of permuted rows */
1818:   PetscMalloc((m+1)*sizeof(PetscInt),&lens);
1819:   for (i=0; i<m; i++) {
1820:     lens[row[i]] = a->i[i+1] - a->i[i];
1821:   }
1822:   MatCreate(A->comm,m,n,m,n,B);
1823:   MatSetType(*B,A->type_name);
1824:   MatSeqAIJSetPreallocation(*B,0,lens);
1825:   PetscFree(lens);

1827:   PetscMalloc(n*sizeof(PetscInt),&cnew);
1828:   for (i=0; i<m; i++) {
1829:     MatGetRow_SeqAIJ(A,i,&nz,&cwork,&vwork);
1830:     for (j=0; j<nz; j++) { cnew[j] = col[cwork[j]];}
1831:     MatSetValues_SeqAIJ(*B,1,&row[i],nz,cnew,vwork,INSERT_VALUES);
1832:     MatRestoreRow_SeqAIJ(A,i,&nz,&cwork,&vwork);
1833:   }
1834:   PetscFree(cnew);
1835:   (*B)->assembled     = PETSC_FALSE;
1836:   MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);
1837:   MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);
1838:   ISRestoreIndices(irowp,&row);
1839:   ISRestoreIndices(icolp,&col);
1840:   ISDestroy(irowp);
1841:   ISDestroy(icolp);
1842:   return(0);
1843: }

1847: PetscErrorCode MatPrintHelp_SeqAIJ(Mat A)
1848: {
1849:   static PetscTruth called = PETSC_FALSE;
1850:   MPI_Comm          comm = A->comm;
1851:   PetscErrorCode    ierr;

1854:   if (called) {return(0);} else called = PETSC_TRUE;
1855:   (*PetscHelpPrintf)(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n");
1856:   (*PetscHelpPrintf)(comm,"  -mat_lu_pivotthreshold <threshold>: Set pivoting threshold\n");
1857:   (*PetscHelpPrintf)(comm,"  -mat_aij_oneindex: internal indices begin at 1 instead of the default 0.\n");
1858:   (*PetscHelpPrintf)(comm,"  -mat_aij_no_inode: Do not use inodes\n");
1859:   (*PetscHelpPrintf)(comm,"  -mat_aij_inode_limit <limit>: Set inode limit (max limit=5)\n");
1860:   return(0);
1861: }

1865: PetscErrorCode MatCopy_SeqAIJ(Mat A,Mat B,MatStructure str)
1866: {

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

1875:     if (a->i[A->m] != b->i[B->m]) {
1876:       SETERRQ(PETSC_ERR_ARG_INCOMP,"Number of nonzeros in two matrices are different");
1877:     }
1878:     PetscMemcpy(b->a,a->a,(a->i[A->m])*sizeof(PetscScalar));
1879:   } else {
1880:     MatCopy_Basic(A,B,str);
1881:   }
1882:   return(0);
1883: }

1887: PetscErrorCode MatSetUpPreallocation_SeqAIJ(Mat A)
1888: {

1892:    MatSeqAIJSetPreallocation(A,PETSC_DEFAULT,0);
1893:   return(0);
1894: }

1898: PetscErrorCode MatGetArray_SeqAIJ(Mat A,PetscScalar *array[])
1899: {
1900:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1902:   *array = a->a;
1903:   return(0);
1904: }

1908: PetscErrorCode MatRestoreArray_SeqAIJ(Mat A,PetscScalar *array[])
1909: {
1911:   return(0);
1912: }

1916: PetscErrorCode MatFDColoringApply_SeqAIJ(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
1917: {
1918:   PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void *))coloring->f;
1920:   PetscInt       k,N,start,end,l,row,col,srow,**vscaleforrow,m1,m2;
1921:   PetscScalar    dx,mone = -1.0,*y,*xx,*w3_array;
1922:   PetscScalar    *vscale_array;
1923:   PetscReal      epsilon = coloring->error_rel,umin = coloring->umin;
1924:   Vec            w1,w2,w3;
1925:   void           *fctx = coloring->fctx;
1926:   PetscTruth     flg;

1929:   if (!coloring->w1) {
1930:     VecDuplicate(x1,&coloring->w1);
1931:     PetscLogObjectParent(coloring,coloring->w1);
1932:     VecDuplicate(x1,&coloring->w2);
1933:     PetscLogObjectParent(coloring,coloring->w2);
1934:     VecDuplicate(x1,&coloring->w3);
1935:     PetscLogObjectParent(coloring,coloring->w3);
1936:   }
1937:   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;

1939:   MatSetUnfactored(J);
1940:   PetscOptionsHasName(coloring->prefix,"-mat_fd_coloring_dont_rezero",&flg);
1941:   if (flg) {
1942:     PetscLogInfo(coloring,"MatFDColoringApply_SeqAIJ: Not calling MatZeroEntries()\n");
1943:   } else {
1944:     PetscTruth assembled;
1945:     MatAssembled(J,&assembled);
1946:     if (assembled) {
1947:       MatZeroEntries(J);
1948:     }
1949:   }

1951:   VecGetOwnershipRange(x1,&start,&end);
1952:   VecGetSize(x1,&N);

1954:   /*
1955:        This is a horrible, horrible, hack. See DMMGComputeJacobian_Multigrid() it inproperly sets
1956:      coloring->F for the coarser grids from the finest
1957:   */
1958:   if (coloring->F) {
1959:     VecGetLocalSize(coloring->F,&m1);
1960:     VecGetLocalSize(w1,&m2);
1961:     if (m1 != m2) {
1962:       coloring->F = 0;
1963:     }
1964:   }

1966:   if (coloring->F) {
1967:     w1          = coloring->F;
1968:     coloring->F = 0;
1969:   } else {
1970:     PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);
1971:     (*f)(sctx,x1,w1,fctx);
1972:     PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);
1973:   }

1975:   /* 
1976:       Compute all the scale factors and share with other processors
1977:   */
1978:   VecGetArray(x1,&xx);xx = xx - start;
1979:   VecGetArray(coloring->vscale,&vscale_array);vscale_array = vscale_array - start;
1980:   for (k=0; k<coloring->ncolors; k++) {
1981:     /*
1982:        Loop over each column associated with color adding the 
1983:        perturbation to the vector w3.
1984:     */
1985:     for (l=0; l<coloring->ncolumns[k]; l++) {
1986:       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
1987:       dx  = xx[col];
1988:       if (dx == 0.0) dx = 1.0;
1989: #if !defined(PETSC_USE_COMPLEX)
1990:       if (dx < umin && dx >= 0.0)      dx = umin;
1991:       else if (dx < 0.0 && dx > -umin) dx = -umin;
1992: #else
1993:       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
1994:       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
1995: #endif
1996:       dx                *= epsilon;
1997:       vscale_array[col] = 1.0/dx;
1998:     }
1999:   }
2000:   vscale_array = vscale_array + start;VecRestoreArray(coloring->vscale,&vscale_array);
2001:   VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);
2002:   VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);

2004:   /*  VecView(coloring->vscale,PETSC_VIEWER_STDOUT_WORLD);
2005:       VecView(x1,PETSC_VIEWER_STDOUT_WORLD);*/

2007:   if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
2008:   else                        vscaleforrow = coloring->columnsforrow;

2010:   VecGetArray(coloring->vscale,&vscale_array);
2011:   /*
2012:       Loop over each color
2013:   */
2014:   for (k=0; k<coloring->ncolors; k++) {
2015:     coloring->currentcolor = k;
2016:     VecCopy(x1,w3);
2017:     VecGetArray(w3,&w3_array);w3_array = w3_array - start;
2018:     /*
2019:        Loop over each column associated with color adding the 
2020:        perturbation to the vector w3.
2021:     */
2022:     for (l=0; l<coloring->ncolumns[k]; l++) {
2023:       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
2024:       dx  = xx[col];
2025:       if (dx == 0.0) dx = 1.0;
2026: #if !defined(PETSC_USE_COMPLEX)
2027:       if (dx < umin && dx >= 0.0)      dx = umin;
2028:       else if (dx < 0.0 && dx > -umin) dx = -umin;
2029: #else
2030:       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
2031:       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
2032: #endif
2033:       dx            *= epsilon;
2034:       if (!PetscAbsScalar(dx)) SETERRQ(PETSC_ERR_PLIB,"Computed 0 differencing parameter");
2035:       w3_array[col] += dx;
2036:     }
2037:     w3_array = w3_array + start; VecRestoreArray(w3,&w3_array);

2039:     /*
2040:        Evaluate function at x1 + dx (here dx is a vector of perturbations)
2041:     */

2043:     PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);
2044:     (*f)(sctx,w3,w2,fctx);
2045:     PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);
2046:     VecAXPY(&mone,w1,w2);

2048:     /*
2049:        Loop over rows of vector, putting results into Jacobian matrix
2050:     */
2051:     VecGetArray(w2,&y);
2052:     for (l=0; l<coloring->nrows[k]; l++) {
2053:       row    = coloring->rows[k][l];
2054:       col    = coloring->columnsforrow[k][l];
2055:       y[row] *= vscale_array[vscaleforrow[k][l]];
2056:       srow   = row + start;
2057:       MatSetValues_SeqAIJ(J,1,&srow,1,&col,y+row,INSERT_VALUES);
2058:     }
2059:     VecRestoreArray(w2,&y);
2060:   }
2061:   coloring->currentcolor = k;
2062:   VecRestoreArray(coloring->vscale,&vscale_array);
2063:   xx = xx + start; VecRestoreArray(x1,&xx);
2064:   MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
2065:   MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
2066:   return(0);
2067: }

2069:  #include petscblaslapack.h
2072: PetscErrorCode MatAXPY_SeqAIJ(const PetscScalar a[],Mat X,Mat Y,MatStructure str)
2073: {
2075:   PetscInt       i;
2076:   Mat_SeqAIJ     *x  = (Mat_SeqAIJ *)X->data,*y = (Mat_SeqAIJ *)Y->data;
2077:   PetscBLASInt   one=1,bnz = (PetscBLASInt)x->nz;

2080:   if (str == SAME_NONZERO_PATTERN) {
2081:     BLaxpy_(&bnz,(PetscScalar*)a,x->a,&one,y->a,&one);
2082:   } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
2083:     if (y->xtoy && y->XtoY != X) {
2084:       PetscFree(y->xtoy);
2085:       MatDestroy(y->XtoY);
2086:     }
2087:     if (!y->xtoy) { /* get xtoy */
2088:       MatAXPYGetxtoy_Private(X->m,x->i,x->j,PETSC_NULL, y->i,y->j,PETSC_NULL, &y->xtoy);
2089:       y->XtoY = X;
2090:     }
2091:     for (i=0; i<x->nz; i++) y->a[y->xtoy[i]] += (*a)*(x->a[i]);
2092:     PetscLogInfo(0,"MatAXPY_SeqAIJ: ratio of nnz(X)/nnz(Y): %d/%d = %g\n",x->nz,y->nz,(PetscReal)(x->nz)/y->nz);
2093:   } else {
2094:     MatAXPY_Basic(a,X,Y,str);
2095:   }
2096:   return(0);
2097: }

2101: PetscErrorCode MatSetBlockSize_SeqAIJ(Mat A,PetscInt bs)
2102: {
2104:   return(0);
2105: }

2107: /* -------------------------------------------------------------------*/
2108: static struct _MatOps MatOps_Values = {MatSetValues_SeqAIJ,
2109:        MatGetRow_SeqAIJ,
2110:        MatRestoreRow_SeqAIJ,
2111:        MatMult_SeqAIJ,
2112: /* 4*/ MatMultAdd_SeqAIJ,
2113:        MatMultTranspose_SeqAIJ,
2114:        MatMultTransposeAdd_SeqAIJ,
2115:        MatSolve_SeqAIJ,
2116:        MatSolveAdd_SeqAIJ,
2117:        MatSolveTranspose_SeqAIJ,
2118: /*10*/ MatSolveTransposeAdd_SeqAIJ,
2119:        MatLUFactor_SeqAIJ,
2120:        0,
2121:        MatRelax_SeqAIJ,
2122:        MatTranspose_SeqAIJ,
2123: /*15*/ MatGetInfo_SeqAIJ,
2124:        MatEqual_SeqAIJ,
2125:        MatGetDiagonal_SeqAIJ,
2126:        MatDiagonalScale_SeqAIJ,
2127:        MatNorm_SeqAIJ,
2128: /*20*/ 0,
2129:        MatAssemblyEnd_SeqAIJ,
2130:        MatCompress_SeqAIJ,
2131:        MatSetOption_SeqAIJ,
2132:        MatZeroEntries_SeqAIJ,
2133: /*25*/ MatZeroRows_SeqAIJ,
2134:        MatLUFactorSymbolic_SeqAIJ,
2135:        MatLUFactorNumeric_SeqAIJ,
2136:        MatCholeskyFactorSymbolic_SeqAIJ,
2137:        MatCholeskyFactorNumeric_SeqAIJ,
2138: /*30*/ MatSetUpPreallocation_SeqAIJ,
2139:        MatILUFactorSymbolic_SeqAIJ,
2140:        MatICCFactorSymbolic_SeqAIJ,
2141:        MatGetArray_SeqAIJ,
2142:        MatRestoreArray_SeqAIJ,
2143: /*35*/ MatDuplicate_SeqAIJ,
2144:        0,
2145:        0,
2146:        MatILUFactor_SeqAIJ,
2147:        0,
2148: /*40*/ MatAXPY_SeqAIJ,
2149:        MatGetSubMatrices_SeqAIJ,
2150:        MatIncreaseOverlap_SeqAIJ,
2151:        MatGetValues_SeqAIJ,
2152:        MatCopy_SeqAIJ,
2153: /*45*/ MatPrintHelp_SeqAIJ,
2154:        MatScale_SeqAIJ,
2155:        0,
2156:        0,
2157:        MatILUDTFactor_SeqAIJ,
2158: /*50*/ MatSetBlockSize_SeqAIJ,
2159:        MatGetRowIJ_SeqAIJ,
2160:        MatRestoreRowIJ_SeqAIJ,
2161:        MatGetColumnIJ_SeqAIJ,
2162:        MatRestoreColumnIJ_SeqAIJ,
2163: /*55*/ MatFDColoringCreate_SeqAIJ,
2164:        0,
2165:        0,
2166:        MatPermute_SeqAIJ,
2167:        0,
2168: /*60*/ 0,
2169:        MatDestroy_SeqAIJ,
2170:        MatView_SeqAIJ,
2171:        MatGetPetscMaps_Petsc,
2172:        0,
2173: /*65*/ 0,
2174:        0,
2175:        0,
2176:        0,
2177:        0,
2178: /*70*/ 0,
2179:        0,
2180:        MatSetColoring_SeqAIJ,
2181:        MatSetValuesAdic_SeqAIJ,
2182:        MatSetValuesAdifor_SeqAIJ,
2183: /*75*/ MatFDColoringApply_SeqAIJ,
2184:        0,
2185:        0,
2186:        0,
2187:        0,
2188: /*80*/ 0,
2189:        0,
2190:        0,
2191:        0,
2192:        MatLoad_SeqAIJ,
2193: /*85*/ MatIsSymmetric_SeqAIJ,
2194:        0,
2195:        0,
2196:        0,
2197:        0,
2198: /*90*/ MatMatMult_SeqAIJ_SeqAIJ,
2199:        MatMatMultSymbolic_SeqAIJ_SeqAIJ,
2200:        MatMatMultNumeric_SeqAIJ_SeqAIJ,
2201:        MatPtAP_SeqAIJ_SeqAIJ,
2202:        MatPtAPSymbolic_SeqAIJ_SeqAIJ,
2203: /*95*/ MatPtAPNumeric_SeqAIJ_SeqAIJ,
2204:        MatMatMultTranspose_SeqAIJ_SeqAIJ,
2205:        MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ,
2206:        MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ,
2207: };

2212: PetscErrorCode MatSeqAIJSetColumnIndices_SeqAIJ(Mat mat,PetscInt *indices)
2213: {
2214:   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data;
2215:   PetscInt   i,nz,n;


2219:   nz = aij->maxnz;
2220:   n  = mat->n;
2221:   for (i=0; i<nz; i++) {
2222:     aij->j[i] = indices[i];
2223:   }
2224:   aij->nz = nz;
2225:   for (i=0; i<n; i++) {
2226:     aij->ilen[i] = aij->imax[i];
2227:   }

2229:   return(0);
2230: }

2235: /*@
2236:     MatSeqAIJSetColumnIndices - Set the column indices for all the rows
2237:        in the matrix.

2239:   Input Parameters:
2240: +  mat - the SeqAIJ matrix
2241: -  indices - the column indices

2243:   Level: advanced

2245:   Notes:
2246:     This can be called if you have precomputed the nonzero structure of the 
2247:   matrix and want to provide it to the matrix object to improve the performance
2248:   of the MatSetValues() operation.

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

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

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

2257: @*/
2258: PetscErrorCode MatSeqAIJSetColumnIndices(Mat mat,PetscInt *indices)
2259: {
2260:   PetscErrorCode ierr,(*f)(Mat,PetscInt *);

2265:   PetscObjectQueryFunction((PetscObject)mat,"MatSeqAIJSetColumnIndices_C",(void (**)(void))&f);
2266:   if (f) {
2267:     (*f)(mat,indices);
2268:   } else {
2269:     SETERRQ(PETSC_ERR_SUP,"Wrong type of matrix to set column indices");
2270:   }
2271:   return(0);
2272: }

2274: /* ----------------------------------------------------------------------------------------*/

2279: PetscErrorCode MatStoreValues_SeqAIJ(Mat mat)
2280: {
2281:   Mat_SeqAIJ     *aij = (Mat_SeqAIJ *)mat->data;
2283:   size_t         nz = aij->i[mat->m];

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

2290:   /* allocate space for values if not already there */
2291:   if (!aij->saved_values) {
2292:     PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);
2293:   }

2295:   /* copy values over */
2296:   PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));
2297:   return(0);
2298: }

2303: /*@
2304:     MatStoreValues - Stashes a copy of the matrix values; this allows, for 
2305:        example, reuse of the linear part of a Jacobian, while recomputing the 
2306:        nonlinear portion.

2308:    Collect on Mat

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

2313:   Level: advanced

2315:   Common Usage, with SNESSolve():
2316: $    Create Jacobian matrix
2317: $    Set linear terms into matrix
2318: $    Apply boundary conditions to matrix, at this time matrix must have 
2319: $      final nonzero structure (i.e. setting the nonlinear terms and applying 
2320: $      boundary conditions again will not change the nonzero structure
2321: $    MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS);
2322: $    MatStoreValues(mat);
2323: $    Call SNESSetJacobian() with matrix
2324: $    In your Jacobian routine
2325: $      MatRetrieveValues(mat);
2326: $      Set nonlinear terms in matrix
2327:  
2328:   Common Usage without SNESSolve(), i.e. when you handle nonlinear solve yourself:
2329: $    // build linear portion of Jacobian 
2330: $    MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS);
2331: $    MatStoreValues(mat);
2332: $    loop over nonlinear iterations
2333: $       MatRetrieveValues(mat);
2334: $       // call MatSetValues(mat,...) to set nonliner portion of Jacobian 
2335: $       // call MatAssemblyBegin/End() on matrix
2336: $       Solve linear system with Jacobian
2337: $    endloop 

2339:   Notes:
2340:     Matrix must already be assemblied before calling this routine
2341:     Must set the matrix option MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS); before 
2342:     calling this routine.

2344: .seealso: MatRetrieveValues()

2346: @*/
2347: PetscErrorCode MatStoreValues(Mat mat)
2348: {
2349:   PetscErrorCode ierr,(*f)(Mat);

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

2356:   PetscObjectQueryFunction((PetscObject)mat,"MatStoreValues_C",(void (**)(void))&f);
2357:   if (f) {
2358:     (*f)(mat);
2359:   } else {
2360:     SETERRQ(PETSC_ERR_SUP,"Wrong type of matrix to store values");
2361:   }
2362:   return(0);
2363: }

2368: PetscErrorCode MatRetrieveValues_SeqAIJ(Mat mat)
2369: {
2370:   Mat_SeqAIJ     *aij = (Mat_SeqAIJ *)mat->data;
2372:   PetscInt       nz = aij->i[mat->m];

2375:   if (aij->nonew != 1) {
2376:     SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
2377:   }
2378:   if (!aij->saved_values) {
2379:     SETERRQ(PETSC_ERR_ORDER,"Must call MatStoreValues(A);first");
2380:   }
2381:   /* copy values over */
2382:   PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));
2383:   return(0);
2384: }

2389: /*@
2390:     MatRetrieveValues - Retrieves the copy of the matrix values; this allows, for 
2391:        example, reuse of the linear part of a Jacobian, while recomputing the 
2392:        nonlinear portion.

2394:    Collect on Mat

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

2399:   Level: advanced

2401: .seealso: MatStoreValues()

2403: @*/
2404: PetscErrorCode MatRetrieveValues(Mat mat)
2405: {
2406:   PetscErrorCode ierr,(*f)(Mat);

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

2413:   PetscObjectQueryFunction((PetscObject)mat,"MatRetrieveValues_C",(void (**)(void))&f);
2414:   if (f) {
2415:     (*f)(mat);
2416:   } else {
2417:     SETERRQ(PETSC_ERR_SUP,"Wrong type of matrix to retrieve values");
2418:   }
2419:   return(0);
2420: }


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

2433:    Collective on MPI_Comm

2435:    Input Parameters:
2436: +  comm - MPI communicator, set to PETSC_COMM_SELF
2437: .  m - number of rows
2438: .  n - number of columns
2439: .  nz - number of nonzeros per row (same for all rows)
2440: -  nnz - array containing the number of nonzeros in the various rows 
2441:          (possibly different for each row) or PETSC_NULL

2443:    Output Parameter:
2444: .  A - the matrix 

2446:    Notes:
2447:    If nnz is given then nz is ignored

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

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

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

2464:    Options Database Keys:
2465: +  -mat_aij_no_inode  - Do not use inodes
2466: .  -mat_aij_inode_limit <limit> - Sets inode limit (max limit=5)
2467: -  -mat_aij_oneindex - Internally use indexing starting at 1
2468:         rather than 0.  Note that when calling MatSetValues(),
2469:         the user still MUST index entries starting at 0!

2471:    Level: intermediate

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

2475: @*/
2476: PetscErrorCode MatCreateSeqAIJ(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A)
2477: {

2481:   MatCreate(comm,m,n,m,n,A);
2482:   MatSetType(*A,MATSEQAIJ);
2483:   MatSeqAIJSetPreallocation(*A,nz,nnz);
2484:   return(0);
2485: }

2487: #define SKIP_ALLOCATION -4

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

2497:    Collective on MPI_Comm

2499:    Input Parameters:
2500: +  comm - MPI communicator, set to PETSC_COMM_SELF
2501: .  m - number of rows
2502: .  n - number of columns
2503: .  nz - number of nonzeros per row (same for all rows)
2504: -  nnz - array containing the number of nonzeros in the various rows 
2505:          (possibly different for each row) or PETSC_NULL

2507:    Output Parameter:
2508: .  A - the matrix 

2510:    Notes:
2511:      If nnz is given then nz is ignored

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

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

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

2528:    Options Database Keys:
2529: +  -mat_aij_no_inode  - Do not use inodes
2530: .  -mat_aij_inode_limit <limit> - Sets inode limit (max limit=5)
2531: -  -mat_aij_oneindex - Internally use indexing starting at 1
2532:         rather than 0.  Note that when calling MatSetValues(),
2533:         the user still MUST index entries starting at 0!

2535:    Level: intermediate

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

2539: @*/
2540: PetscErrorCode MatSeqAIJSetPreallocation(Mat B,PetscInt nz,const PetscInt nnz[])
2541: {
2542:   PetscErrorCode ierr,(*f)(Mat,PetscInt,const PetscInt[]);

2545:   PetscObjectQueryFunction((PetscObject)B,"MatSeqAIJSetPreallocation_C",(void (**)(void))&f);
2546:   if (f) {
2547:     (*f)(B,nz,nnz);
2548:   }
2549:   return(0);
2550: }

2555: PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat B,PetscInt nz,PetscInt *nnz)
2556: {
2557:   Mat_SeqAIJ     *b;
2558:   size_t         len = 0;
2559:   PetscTruth     skipallocation = PETSC_FALSE;
2561:   PetscInt       i;

2564: 
2565:   if (nz == SKIP_ALLOCATION) {
2566:     skipallocation = PETSC_TRUE;
2567:     nz             = 0;
2568:   }

2570:   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
2571:   if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz);
2572:   if (nnz) {
2573:     for (i=0; i<B->m; i++) {
2574:       if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]);
2575:       if (nnz[i] > B->n) SETERRQ3(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be greater than row length: local row %d value %d rowlength %d",i,nnz[i],B->n);
2576:     }
2577:   }

2579:   B->preallocated = PETSC_TRUE;
2580:   b = (Mat_SeqAIJ*)B->data;

2582:   PetscMalloc((B->m+1)*sizeof(PetscInt),&b->imax);
2583:   if (!nnz) {
2584:     if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 10;
2585:     else if (nz <= 0)        nz = 1;
2586:     for (i=0; i<B->m; i++) b->imax[i] = nz;
2587:     nz = nz*B->m;
2588:   } else {
2589:     nz = 0;
2590:     for (i=0; i<B->m; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
2591:   }

2593:   if (!skipallocation) {
2594:     /* allocate the matrix space */
2595:     len             = ((size_t) nz)*(sizeof(PetscInt) + sizeof(PetscScalar)) + (B->m+1)*sizeof(PetscInt);
2596:     PetscMalloc(len,&b->a);
2597:     b->j            = (PetscInt*)(b->a + nz);
2598:     PetscMemzero(b->j,nz*sizeof(PetscInt));
2599:     b->i            = b->j + nz;
2600:     b->i[0] = 0;
2601:     for (i=1; i<B->m+1; i++) {
2602:       b->i[i] = b->i[i-1] + b->imax[i-1];
2603:     }
2604:     b->singlemalloc = PETSC_TRUE;
2605:     b->freedata     = PETSC_TRUE;
2606:   } else {
2607:     b->freedata     = PETSC_FALSE;
2608:   }

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

2615:   b->nz                = 0;
2616:   b->maxnz             = nz;
2617:   B->info.nz_unneeded  = (double)b->maxnz;
2618:   return(0);
2619: }

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

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

2629:   Level: beginner

2631: .seealso: MatCreateSeqAIJ(), MatSetFromOptions(), MatSetType(), MatCreate(), MatType
2632: M*/

2637: PetscErrorCode MatCreate_SeqAIJ(Mat B)
2638: {
2639:   Mat_SeqAIJ     *b;
2641:   PetscMPIInt    size;

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

2647:   B->m = B->M = PetscMax(B->m,B->M);
2648:   B->n = B->N = PetscMax(B->n,B->N);

2650:   PetscNew(Mat_SeqAIJ,&b);
2651:   B->data             = (void*)b;
2652:   PetscMemzero(b,sizeof(Mat_SeqAIJ));
2653:   PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
2654:   B->factor           = 0;
2655:   B->lupivotthreshold = 1.0;
2656:   B->mapping          = 0;
2657:   PetscOptionsGetReal(B->prefix,"-mat_lu_pivotthreshold",&B->lupivotthreshold,PETSC_NULL);
2658:   PetscOptionsHasName(B->prefix,"-pc_ilu_preserve_row_sums",&b->ilu_preserve_row_sums);
2659:   b->row              = 0;
2660:   b->col              = 0;
2661:   b->icol             = 0;
2662:   b->reallocs         = 0;
2663:   b->lu_shift         = PETSC_FALSE;
2664: 
2665:   PetscMapCreateMPI(B->comm,B->m,B->m,&B->rmap);
2666:   PetscMapCreateMPI(B->comm,B->n,B->n,&B->cmap);

2668:   b->sorted            = PETSC_FALSE;
2669:   b->ignorezeroentries = PETSC_FALSE;
2670:   b->roworiented       = PETSC_TRUE;
2671:   b->nonew             = 0;
2672:   b->diag              = 0;
2673:   b->solve_work        = 0;
2674:   B->spptr             = 0;
2675:   b->inode.use         = PETSC_TRUE;
2676:   b->inode.node_count  = 0;
2677:   b->inode.size        = 0;
2678:   b->inode.limit       = 5;
2679:   b->inode.max_limit   = 5;
2680:   b->saved_values      = 0;
2681:   b->idiag             = 0;
2682:   b->ssor              = 0;
2683:   b->keepzeroedrows    = PETSC_FALSE;
2684:   b->xtoy              = 0;
2685:   b->XtoY              = 0;

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

2689:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqAIJSetColumnIndices_C",
2690:                                      "MatSeqAIJSetColumnIndices_SeqAIJ",
2691:                                      MatSeqAIJSetColumnIndices_SeqAIJ);
2692:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
2693:                                      "MatStoreValues_SeqAIJ",
2694:                                      MatStoreValues_SeqAIJ);
2695:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
2696:                                      "MatRetrieveValues_SeqAIJ",
2697:                                      MatRetrieveValues_SeqAIJ);
2698: #if !defined(PETSC_USE_64BIT_INT)
2699:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_seqsbaij_C",
2700:                                      "MatConvert_SeqAIJ_SeqSBAIJ",
2701:                                       MatConvert_SeqAIJ_SeqSBAIJ);
2702: #endif
2703:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_seqbaij_C",
2704:                                      "MatConvert_SeqAIJ_SeqBAIJ",
2705:                                       MatConvert_SeqAIJ_SeqBAIJ);
2706:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatIsTranspose_C",
2707:                                      "MatIsTranspose_SeqAIJ",
2708:                                       MatIsTranspose_SeqAIJ);
2709:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqAIJSetPreallocation_C",
2710:                                      "MatSeqAIJSetPreallocation_SeqAIJ",
2711:                                       MatSeqAIJSetPreallocation_SeqAIJ);
2712:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatReorderForNonzeroDiagonal_C",
2713:                                      "MatReorderForNonzeroDiagonal_SeqAIJ",
2714:                                       MatReorderForNonzeroDiagonal_SeqAIJ);
2715:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatAdjustForInodes_C",
2716:                                      "MatAdjustForInodes_SeqAIJ",
2717:                                       MatAdjustForInodes_SeqAIJ);
2718:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqAIJGetInodeSizes_C",
2719:                                      "MatSeqAIJGetInodeSizes_SeqAIJ",
2720:                                       MatSeqAIJGetInodeSizes_SeqAIJ);
2721:   return(0);
2722: }

2727: PetscErrorCode MatDuplicate_SeqAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
2728: {
2729:   Mat            C;
2730:   Mat_SeqAIJ     *c,*a = (Mat_SeqAIJ*)A->data;
2732:   PetscInt       i,m = A->m;
2733:   size_t         len;

2736:   *B = 0;
2737:   MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);
2738:   MatSetType(C,A->type_name);
2739:   PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));
2740: 
2741:   c    = (Mat_SeqAIJ*)C->data;

2743:   C->factor         = A->factor;
2744:   c->row            = 0;
2745:   c->col            = 0;
2746:   c->icol           = 0;
2747:   c->keepzeroedrows = a->keepzeroedrows;
2748:   C->assembled      = PETSC_TRUE;

2750:   C->M          = A->m;
2751:   C->N          = A->n;
2752:   C->bs         = A->bs;

2754:   PetscMalloc((m+1)*sizeof(PetscInt),&c->imax);
2755:   PetscMalloc((m+1)*sizeof(PetscInt),&c->ilen);
2756:   for (i=0; i<m; i++) {
2757:     c->imax[i] = a->imax[i];
2758:     c->ilen[i] = a->ilen[i];
2759:   }

2761:   /* allocate the matrix space */
2762:   c->singlemalloc = PETSC_TRUE;
2763:   len   = ((size_t) (m+1))*sizeof(PetscInt)+(a->i[m])*(sizeof(PetscScalar)+sizeof(PetscInt));
2764:   PetscMalloc(len,&c->a);
2765:   c->j  = (PetscInt*)(c->a + a->i[m] );
2766:   c->i  = c->j + a->i[m];
2767:   PetscMemcpy(c->i,a->i,(m+1)*sizeof(PetscInt));
2768:   if (m > 0) {
2769:     PetscMemcpy(c->j,a->j,(a->i[m])*sizeof(PetscInt));
2770:     if (cpvalues == MAT_COPY_VALUES) {
2771:       PetscMemcpy(c->a,a->a,(a->i[m])*sizeof(PetscScalar));
2772:     } else {
2773:       PetscMemzero(c->a,(a->i[m])*sizeof(PetscScalar));
2774:     }
2775:   }

2777:   PetscLogObjectMemory(C,len+2*(m+1)*sizeof(PetscInt)+sizeof(struct _p_Mat)+sizeof(Mat_SeqAIJ));
2778:   c->sorted      = a->sorted;
2779:   c->roworiented = a->roworiented;
2780:   c->nonew       = a->nonew;
2781:   c->ilu_preserve_row_sums = a->ilu_preserve_row_sums;
2782:   c->saved_values = 0;
2783:   c->idiag        = 0;
2784:   c->ssor         = 0;
2785:   c->ignorezeroentries = a->ignorezeroentries;
2786:   c->freedata     = PETSC_TRUE;

2788:   if (a->diag) {
2789:     PetscMalloc((m+1)*sizeof(PetscInt),&c->diag);
2790:     PetscLogObjectMemory(C,(m+1)*sizeof(PetscInt));
2791:     for (i=0; i<m; i++) {
2792:       c->diag[i] = a->diag[i];
2793:     }
2794:   } else c->diag        = 0;
2795:   c->inode.use          = a->inode.use;
2796:   c->inode.limit        = a->inode.limit;
2797:   c->inode.max_limit    = a->inode.max_limit;
2798:   if (a->inode.size){
2799:     PetscMalloc((m+1)*sizeof(PetscInt),&c->inode.size);
2800:     c->inode.node_count = a->inode.node_count;
2801:     PetscMemcpy(c->inode.size,a->inode.size,(m+1)*sizeof(PetscInt));
2802:   } else {
2803:     c->inode.size       = 0;
2804:     c->inode.node_count = 0;
2805:   }
2806:   c->nz                 = a->nz;
2807:   c->maxnz              = a->maxnz;
2808:   c->solve_work         = 0;
2809:   C->preallocated       = PETSC_TRUE;

2811:   *B = C;
2812:   PetscFListDuplicate(A->qlist,&C->qlist);
2813:   return(0);
2814: }

2818: PetscErrorCode MatLoad_SeqAIJ(PetscViewer viewer,const MatType type,Mat *A)
2819: {
2820:   Mat_SeqAIJ     *a;
2821:   Mat            B;
2823:   PetscInt       i,nz,header[4],*rowlengths = 0,M,N;
2824:   int            fd;
2825:   PetscMPIInt    size;
2826:   MPI_Comm       comm;
2827: 
2829:   PetscObjectGetComm((PetscObject)viewer,&comm);
2830:   MPI_Comm_size(comm,&size);
2831:   if (size > 1) SETERRQ(PETSC_ERR_ARG_SIZ,"view must have one processor");
2832:   PetscViewerBinaryGetDescriptor(viewer,&fd);
2833:   PetscBinaryRead(fd,header,4,PETSC_INT);
2834:   if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object in file");
2835:   M = header[1]; N = header[2]; nz = header[3];

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

2841:   /* read in row lengths */
2842:   PetscMalloc(M*sizeof(PetscInt),&rowlengths);
2843:   PetscBinaryRead(fd,rowlengths,M,PETSC_INT);

2845:   /* create our matrix */
2846:   MatCreate(comm,PETSC_DECIDE,PETSC_DECIDE,M,N,&B);
2847:   MatSetType(B,type);
2848:   MatSeqAIJSetPreallocation(B,0,rowlengths);
2849:   a = (Mat_SeqAIJ*)B->data;

2851:   /* read in column indices and adjust for Fortran indexing*/
2852:   PetscBinaryRead(fd,a->j,nz,PETSC_INT);

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

2857:   /* set matrix "i" values */
2858:   a->i[0] = 0;
2859:   for (i=1; i<= M; i++) {
2860:     a->i[i]      = a->i[i-1] + rowlengths[i-1];
2861:     a->ilen[i-1] = rowlengths[i-1];
2862:   }
2863:   PetscFree(rowlengths);

2865:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
2866:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
2867:   *A = B;
2868:   return(0);
2869: }

2873: PetscErrorCode MatEqual_SeqAIJ(Mat A,Mat B,PetscTruth* flg)
2874: {
2875:   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data,*b = (Mat_SeqAIJ *)B->data;

2879:   /* If the  matrix dimensions are not equal,or no of nonzeros */
2880:   if ((A->m != B->m) || (A->n != B->n) ||(a->nz != b->nz)) {
2881:     *flg = PETSC_FALSE;
2882:     return(0);
2883:   }
2884: 
2885:   /* if the a->i are the same */
2886:   PetscMemcmp(a->i,b->i,(A->m+1)*sizeof(PetscInt),flg);
2887:   if (*flg == PETSC_FALSE) return(0);
2888: 
2889:   /* if a->j are the same */
2890:   PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),flg);
2891:   if (*flg == PETSC_FALSE) return(0);
2892: 
2893:   /* if a->a are the same */
2894:   PetscMemcmp(a->a,b->a,(a->nz)*sizeof(PetscScalar),flg);

2896:   return(0);
2897: 
2898: }

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

2906:       Coolective on MPI_Comm

2908:    Input Parameters:
2909: +   comm - must be an MPI communicator of size 1
2910: .   m - number of rows
2911: .   n - number of columns
2912: .   i - row indices
2913: .   j - column indices
2914: -   a - matrix values

2916:    Output Parameter:
2917: .   mat - the matrix

2919:    Level: intermediate

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

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

2927:        The i and j indices are 0 based

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

2931: @*/
2932: PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt* i,PetscInt*j,PetscScalar *a,Mat *mat)
2933: {
2935:   PetscInt       ii;
2936:   Mat_SeqAIJ     *aij;

2939:   MatCreate(comm,m,n,m,n,mat);
2940:   MatSetType(*mat,MATSEQAIJ);
2941:   MatSeqAIJSetPreallocation(*mat,SKIP_ALLOCATION,0);
2942:   aij  = (Mat_SeqAIJ*)(*mat)->data;

2944:   if (i[0] != 0) {
2945:     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"i (row indices) must start with 0");
2946:   }
2947:   aij->i = i;
2948:   aij->j = j;
2949:   aij->a = a;
2950:   aij->singlemalloc = PETSC_FALSE;
2951:   aij->nonew        = -1;             /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/
2952:   aij->freedata     = PETSC_FALSE;

2954:   for (ii=0; ii<m; ii++) {
2955:     aij->ilen[ii] = aij->imax[ii] = i[ii+1] - i[ii];
2956: #if defined(PETSC_USE_BOPT_g)
2957:     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]);
2958: #endif    
2959:   }
2960: #if defined(PETSC_USE_BOPT_g)
2961:   for (ii=0; ii<aij->i[m]; ii++) {
2962:     if (j[ii] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Negative column index at location = %d index = %d",ii,j[ii]);
2963:     if (j[ii] > n - 1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column index to large at location = %d index = %d",ii,j[ii]);
2964:   }
2965: #endif    

2967:   MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
2968:   MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
2969:   return(0);
2970: }

2974: PetscErrorCode MatSetColoring_SeqAIJ(Mat A,ISColoring coloring)
2975: {
2977:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;

2980:   if (coloring->ctype == IS_COLORING_LOCAL) {
2981:     ISColoringReference(coloring);
2982:     a->coloring = coloring;
2983:   } else if (coloring->ctype == IS_COLORING_GHOSTED) {
2984:     PetscInt             i,*larray;
2985:     ISColoring      ocoloring;
2986:     ISColoringValue *colors;

2988:     /* set coloring for diagonal portion */
2989:     PetscMalloc((A->n+1)*sizeof(PetscInt),&larray);
2990:     for (i=0; i<A->n; i++) {
2991:       larray[i] = i;
2992:     }
2993:     ISGlobalToLocalMappingApply(A->mapping,IS_GTOLM_MASK,A->n,larray,PETSC_NULL,larray);
2994:     PetscMalloc((A->n+1)*sizeof(ISColoringValue),&colors);
2995:     for (i=0; i<A->n; i++) {
2996:       colors[i] = coloring->colors[larray[i]];
2997:     }
2998:     PetscFree(larray);
2999:     ISColoringCreate(PETSC_COMM_SELF,A->n,colors,&ocoloring);
3000:     a->coloring = ocoloring;
3001:   }
3002:   return(0);
3003: }

3005: #if defined(PETSC_HAVE_ADIC) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
3007: #include "adic/ad_utils.h"

3012: PetscErrorCode MatSetValuesAdic_SeqAIJ(Mat A,void *advalues)
3013: {
3014:   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
3015:   PetscInt        m = A->m,*ii = a->i,*jj = a->j,nz,i,j,nlen;
3016:   PetscScalar     *v = a->a,*values = ((PetscScalar*)advalues)+1;
3017:   ISColoringValue *color;

3020:   if (!a->coloring) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Coloring not set for matrix");
3021:   nlen  = PetscADGetDerivTypeSize()/sizeof(PetscScalar);
3022:   color = a->coloring->colors;
3023:   /* loop over rows */
3024:   for (i=0; i<m; i++) {
3025:     nz = ii[i+1] - ii[i];
3026:     /* loop over columns putting computed value into matrix */
3027:     for (j=0; j<nz; j++) {
3028:       *v++ = values[color[*jj++]];
3029:     }
3030:     values += nlen; /* jump to next row of derivatives */
3031:   }
3032:   return(0);
3033: }

3035: #else

3039: PetscErrorCode MatSetValuesAdic_SeqAIJ(Mat A,void *advalues)
3040: {
3042:   SETERRQ(PETSC_ERR_SUP_SYS,"PETSc installed without ADIC");
3043: }

3045: #endif

3049: PetscErrorCode MatSetValuesAdifor_SeqAIJ(Mat A,PetscInt nl,void *advalues)
3050: {
3051:   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
3052:   PetscInt             m = A->m,*ii = a->i,*jj = a->j,nz,i,j;
3053:   PetscScalar     *v = a->a,*values = (PetscScalar *)advalues;
3054:   ISColoringValue *color;

3057:   if (!a->coloring) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Coloring not set for matrix");
3058:   color = a->coloring->colors;
3059:   /* loop over rows */
3060:   for (i=0; i<m; i++) {
3061:     nz = ii[i+1] - ii[i];
3062:     /* loop over columns putting computed value into matrix */
3063:     for (j=0; j<nz; j++) {
3064:       *v++ = values[color[*jj++]];
3065:     }
3066:     values += nl; /* jump to next row of derivatives */
3067:   }
3068:   return(0);
3069: }