Actual source code: aij.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  2: /*
  3:     Defines the basic matrix operations for the AIJ (compressed row)
  4:   matrix storage format.
  5: */


  8: #include <../src/mat/impls/aij/seq/aij.h>          /*I "petscmat.h" I*/
  9: #include <petscblaslapack.h>
 10: #include <petscbt.h>
 11: #include <petsc-private/kernels/blocktranspose.h>
 12: #if defined(PETSC_THREADCOMM_ACTIVE)
 13: #include <petscthreadcomm.h>
 14: #endif

 18: PetscErrorCode MatGetColumnNorms_SeqAIJ(Mat A,NormType type,PetscReal *norms)
 19: {
 21:   PetscInt       i,m,n;
 22:   Mat_SeqAIJ     *aij = (Mat_SeqAIJ*)A->data;

 25:   MatGetSize(A,&m,&n);
 26:   PetscMemzero(norms,n*sizeof(PetscReal));
 27:   if (type == NORM_2) {
 28:     for (i=0; i<aij->i[m]; i++) {
 29:       norms[aij->j[i]] += PetscAbsScalar(aij->a[i]*aij->a[i]);
 30:     }
 31:   } else if (type == NORM_1) {
 32:     for (i=0; i<aij->i[m]; i++) {
 33:       norms[aij->j[i]] += PetscAbsScalar(aij->a[i]);
 34:     }
 35:   } else if (type == NORM_INFINITY) {
 36:     for (i=0; i<aij->i[m]; i++) {
 37:       norms[aij->j[i]] = PetscMax(PetscAbsScalar(aij->a[i]),norms[aij->j[i]]);
 38:     }
 39:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown NormType");

 41:   if (type == NORM_2) {
 42:     for (i=0; i<n; i++) norms[i] = PetscSqrtReal(norms[i]);
 43:   }
 44:   return(0);
 45: }

 49: PetscErrorCode MatFindOffBlockDiagonalEntries_SeqAIJ(Mat A,IS *is)
 50: {
 51:   Mat_SeqAIJ      *a  = (Mat_SeqAIJ*)A->data;
 52:   PetscInt        i,m=A->rmap->n,cnt = 0, bs = A->rmap->bs;
 53:   const PetscInt  *jj = a->j,*ii = a->i;
 54:   PetscInt        *rows;
 55:   PetscErrorCode  ierr;

 58:   for (i=0; i<m; i++) {
 59:     if ((ii[i] != ii[i+1]) && ((jj[ii[i]] < bs*(i/bs)) || (jj[ii[i+1]-1] > bs*((i+bs)/bs)-1))) {
 60:       cnt++;
 61:     }
 62:   }
 63:   PetscMalloc1(cnt,&rows);
 64:   cnt  = 0;
 65:   for (i=0; i<m; i++) {
 66:     if ((ii[i] != ii[i+1]) && ((jj[ii[i]] < bs*(i/bs)) || (jj[ii[i+1]-1] > bs*((i+bs)/bs)-1))) {
 67:       rows[cnt] = i;
 68:       cnt++;
 69:     }
 70:   }
 71:   ISCreateGeneral(PETSC_COMM_SELF,cnt,rows,PETSC_OWN_POINTER,is);
 72:   return(0);
 73: }

 77: PetscErrorCode MatFindZeroDiagonals_SeqAIJ_Private(Mat A,PetscInt *nrows,PetscInt **zrows)
 78: {
 79:   Mat_SeqAIJ      *a  = (Mat_SeqAIJ*)A->data;
 80:   const MatScalar *aa = a->a;
 81:   PetscInt        i,m=A->rmap->n,cnt = 0;
 82:   const PetscInt  *jj = a->j,*diag;
 83:   PetscInt        *rows;
 84:   PetscErrorCode  ierr;

 87:   MatMarkDiagonal_SeqAIJ(A);
 88:   diag = a->diag;
 89:   for (i=0; i<m; i++) {
 90:     if ((jj[diag[i]] != i) || (aa[diag[i]] == 0.0)) {
 91:       cnt++;
 92:     }
 93:   }
 94:   PetscMalloc1(cnt,&rows);
 95:   cnt  = 0;
 96:   for (i=0; i<m; i++) {
 97:     if ((jj[diag[i]] != i) || (aa[diag[i]] == 0.0)) {
 98:       rows[cnt++] = i;
 99:     }
100:   }
101:   *nrows = cnt;
102:   *zrows = rows;
103:   return(0);
104: }

108: PetscErrorCode MatFindZeroDiagonals_SeqAIJ(Mat A,IS *zrows)
109: {
110:   PetscInt       nrows,*rows;

114:   *zrows = NULL;
115:   MatFindZeroDiagonals_SeqAIJ_Private(A,&nrows,&rows);
116:   ISCreateGeneral(PetscObjectComm((PetscObject)A),nrows,rows,PETSC_OWN_POINTER,zrows);
117:   return(0);
118: }

122: PetscErrorCode MatFindNonzeroRows_SeqAIJ(Mat A,IS *keptrows)
123: {
124:   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
125:   const MatScalar *aa;
126:   PetscInt        m=A->rmap->n,cnt = 0;
127:   const PetscInt  *ii;
128:   PetscInt        n,i,j,*rows;
129:   PetscErrorCode  ierr;

132:   *keptrows = 0;
133:   ii        = a->i;
134:   for (i=0; i<m; i++) {
135:     n = ii[i+1] - ii[i];
136:     if (!n) {
137:       cnt++;
138:       goto ok1;
139:     }
140:     aa = a->a + ii[i];
141:     for (j=0; j<n; j++) {
142:       if (aa[j] != 0.0) goto ok1;
143:     }
144:     cnt++;
145: ok1:;
146:   }
147:   if (!cnt) return(0);
148:   PetscMalloc1((A->rmap->n-cnt),&rows);
149:   cnt  = 0;
150:   for (i=0; i<m; i++) {
151:     n = ii[i+1] - ii[i];
152:     if (!n) continue;
153:     aa = a->a + ii[i];
154:     for (j=0; j<n; j++) {
155:       if (aa[j] != 0.0) {
156:         rows[cnt++] = i;
157:         break;
158:       }
159:     }
160:   }
161:   ISCreateGeneral(PETSC_COMM_SELF,cnt,rows,PETSC_OWN_POINTER,keptrows);
162:   return(0);
163: }

167: PetscErrorCode  MatDiagonalSet_SeqAIJ(Mat Y,Vec D,InsertMode is)
168: {
170:   Mat_SeqAIJ     *aij = (Mat_SeqAIJ*) Y->data;
171:   PetscInt       i,*diag, m = Y->rmap->n;
172:   MatScalar      *aa = aij->a;
173:   PetscScalar    *v;
174:   PetscBool      missing;

177:   if (Y->assembled) {
178:     MatMissingDiagonal_SeqAIJ(Y,&missing,NULL);
179:     if (!missing) {
180:       diag = aij->diag;
181:       VecGetArray(D,&v);
182:       if (is == INSERT_VALUES) {
183:         for (i=0; i<m; i++) {
184:           aa[diag[i]] = v[i];
185:         }
186:       } else {
187:         for (i=0; i<m; i++) {
188:           aa[diag[i]] += v[i];
189:         }
190:       }
191:       VecRestoreArray(D,&v);
192:       return(0);
193:     }
194:     MatSeqAIJInvalidateDiagonal(Y);
195:   }
196:   MatDiagonalSet_Default(Y,D,is);
197:   return(0);
198: }

202: PetscErrorCode MatGetRowIJ_SeqAIJ(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *m,const PetscInt *ia[],const PetscInt *ja[],PetscBool  *done)
203: {
204:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
206:   PetscInt       i,ishift;

209:   *m = A->rmap->n;
210:   if (!ia) return(0);
211:   ishift = 0;
212:   if (symmetric && !A->structurally_symmetric) {
213:     MatToSymmetricIJ_SeqAIJ(A->rmap->n,a->i,a->j,ishift,oshift,(PetscInt**)ia,(PetscInt**)ja);
214:   } else if (oshift == 1) {
215:     PetscInt *tia;
216:     PetscInt nz = a->i[A->rmap->n];
217:     /* malloc space and  add 1 to i and j indices */
218:     PetscMalloc1((A->rmap->n+1),&tia);
219:     for (i=0; i<A->rmap->n+1; i++) tia[i] = a->i[i] + 1;
220:     *ia = tia;
221:     if (ja) {
222:       PetscInt *tja;
223:       PetscMalloc1((nz+1),&tja);
224:       for (i=0; i<nz; i++) tja[i] = a->j[i] + 1;
225:       *ja = tja;
226:     }
227:   } else {
228:     *ia = a->i;
229:     if (ja) *ja = a->j;
230:   }
231:   return(0);
232: }

236: PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *n,const PetscInt *ia[],const PetscInt *ja[],PetscBool  *done)
237: {

241:   if (!ia) return(0);
242:   if ((symmetric && !A->structurally_symmetric) || oshift == 1) {
243:     PetscFree(*ia);
244:     if (ja) {PetscFree(*ja);}
245:   }
246:   return(0);
247: }

251: PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *nn,const PetscInt *ia[],const PetscInt *ja[],PetscBool  *done)
252: {
253:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
255:   PetscInt       i,*collengths,*cia,*cja,n = A->cmap->n,m = A->rmap->n;
256:   PetscInt       nz = a->i[m],row,*jj,mr,col;

259:   *nn = n;
260:   if (!ia) return(0);
261:   if (symmetric) {
262:     MatToSymmetricIJ_SeqAIJ(A->rmap->n,a->i,a->j,0,oshift,(PetscInt**)ia,(PetscInt**)ja);
263:   } else {
264:     PetscCalloc1(n+1,&collengths);
265:     PetscMalloc1((n+1),&cia);
266:     PetscMalloc1((nz+1),&cja);
267:     jj   = a->j;
268:     for (i=0; i<nz; i++) {
269:       collengths[jj[i]]++;
270:     }
271:     cia[0] = oshift;
272:     for (i=0; i<n; i++) {
273:       cia[i+1] = cia[i] + collengths[i];
274:     }
275:     PetscMemzero(collengths,n*sizeof(PetscInt));
276:     jj   = a->j;
277:     for (row=0; row<m; row++) {
278:       mr = a->i[row+1] - a->i[row];
279:       for (i=0; i<mr; i++) {
280:         col = *jj++;

282:         cja[cia[col] + collengths[col]++ - oshift] = row + oshift;
283:       }
284:     }
285:     PetscFree(collengths);
286:     *ia  = cia; *ja = cja;
287:   }
288:   return(0);
289: }

293: PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *n,const PetscInt *ia[],const PetscInt *ja[],PetscBool  *done)
294: {

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

300:   PetscFree(*ia);
301:   PetscFree(*ja);
302:   return(0);
303: }

305: /*
306:  MatGetColumnIJ_SeqAIJ_Color() and MatRestoreColumnIJ_SeqAIJ_Color() are customized from
307:  MatGetColumnIJ_SeqAIJ() and MatRestoreColumnIJ_SeqAIJ() by adding an output
308:  spidx[], index of a->a, to be used in MatTransposeColoringCreate_SeqAIJ() and MatFDColoringCreate_SeqXAIJ()
309: */
312: PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *nn,const PetscInt *ia[],const PetscInt *ja[],PetscInt *spidx[],PetscBool  *done)
313: {
314:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
316:   PetscInt       i,*collengths,*cia,*cja,n = A->cmap->n,m = A->rmap->n;
317:   PetscInt       nz = a->i[m],row,*jj,mr,col;
318:   PetscInt       *cspidx;

321:   *nn = n;
322:   if (!ia) return(0);

324:   PetscCalloc1(n+1,&collengths);
325:   PetscMalloc1((n+1),&cia);
326:   PetscMalloc1((nz+1),&cja);
327:   PetscMalloc1((nz+1),&cspidx);
328:   jj   = a->j;
329:   for (i=0; i<nz; i++) {
330:     collengths[jj[i]]++;
331:   }
332:   cia[0] = oshift;
333:   for (i=0; i<n; i++) {
334:     cia[i+1] = cia[i] + collengths[i];
335:   }
336:   PetscMemzero(collengths,n*sizeof(PetscInt));
337:   jj   = a->j;
338:   for (row=0; row<m; row++) {
339:     mr = a->i[row+1] - a->i[row];
340:     for (i=0; i<mr; i++) {
341:       col = *jj++;
342:       cspidx[cia[col] + collengths[col] - oshift] = a->i[row] + i; /* index of a->j */
343:       cja[cia[col] + collengths[col]++ - oshift]  = row + oshift;
344:     }
345:   }
346:   PetscFree(collengths);
347:   *ia    = cia; *ja = cja;
348:   *spidx = cspidx;
349:   return(0);
350: }

354: PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *n,const PetscInt *ia[],const PetscInt *ja[],PetscInt *spidx[],PetscBool  *done)
355: {

359:   MatRestoreColumnIJ_SeqAIJ(A,oshift,symmetric,inodecompressed,n,ia,ja,done);
360:   PetscFree(*spidx);
361:   return(0);
362: }

366: PetscErrorCode MatSetValuesRow_SeqAIJ(Mat A,PetscInt row,const PetscScalar v[])
367: {
368:   Mat_SeqAIJ     *a  = (Mat_SeqAIJ*)A->data;
369:   PetscInt       *ai = a->i;

373:   PetscMemcpy(a->a+ai[row],v,(ai[row+1]-ai[row])*sizeof(PetscScalar));
374:   return(0);
375: }

377: /*
378:     MatSeqAIJSetValuesLocalFast - An optimized version of MatSetValuesLocal() for SeqAIJ matrices with several assumptions

380:       -   a single row of values is set with each call
381:       -   no row or column indices are negative or (in error) larger than the number of rows or columns
382:       -   the values are always added to the matrix, not set
383:       -   no new locations are introduced in the nonzero structure of the matrix

385:      This does NOT assume the global column indices are sorted

387: */

389: #include <petsc-private/isimpl.h>
392: PetscErrorCode MatSeqAIJSetValuesLocalFast(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
393: {
394:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
395:   PetscInt       low,high,t,row,nrow,i,col,l;
396:   const PetscInt *rp,*ai = a->i,*ailen = a->ilen,*aj = a->j;
397:   PetscInt       lastcol = -1;
398:   MatScalar      *ap,value,*aa = a->a;
399:   const PetscInt *ridx = A->rmap->mapping->indices,*cidx = A->cmap->mapping->indices;

401:   row = ridx[im[0]];
402:   rp   = aj + ai[row];
403:   ap = aa + ai[row];
404:   nrow = ailen[row];
405:   low  = 0;
406:   high = nrow;
407:   for (l=0; l<n; l++) { /* loop over added columns */
408:     col = cidx[in[l]];
409:     value = v[l];

411:     if (col <= lastcol) low = 0;
412:     else high = nrow;
413:     lastcol = col;
414:     while (high-low > 5) {
415:       t = (low+high)/2;
416:       if (rp[t] > col) high = t;
417:       else low = t;
418:     }
419:     for (i=low; i<high; i++) {
420:       if (rp[i] == col) {
421:         ap[i] += value;
422:         low = i + 1;
423:         break;
424:       }
425:     }
426:   }
427:   return 0;
428: }

432: PetscErrorCode MatSetValues_SeqAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
433: {
434:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
435:   PetscInt       *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N;
436:   PetscInt       *imax = a->imax,*ai = a->i,*ailen = a->ilen;
438:   PetscInt       *aj = a->j,nonew = a->nonew,lastcol = -1;
439:   MatScalar      *ap,value,*aa = a->a;
440:   PetscBool      ignorezeroentries = a->ignorezeroentries;
441:   PetscBool      roworiented       = a->roworiented;

444:   for (k=0; k<m; k++) { /* loop over added rows */
445:     row = im[k];
446:     if (row < 0) continue;
447: #if defined(PETSC_USE_DEBUG)
448:     if (row >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap->n-1);
449: #endif
450:     rp   = aj + ai[row]; ap = aa + ai[row];
451:     rmax = imax[row]; nrow = ailen[row];
452:     low  = 0;
453:     high = nrow;
454:     for (l=0; l<n; l++) { /* loop over added columns */
455:       if (in[l] < 0) continue;
456: #if defined(PETSC_USE_DEBUG)
457:       if (in[l] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->cmap->n-1);
458: #endif
459:       col = in[l];
460:       if (roworiented) {
461:         value = v[l + k*n];
462:       } else {
463:         value = v[k + l*m];
464:       }
465:       if ((value == 0.0 && ignorezeroentries) && (is == ADD_VALUES)) continue;

467:       if (col <= lastcol) low = 0;
468:       else high = nrow;
469:       lastcol = col;
470:       while (high-low > 5) {
471:         t = (low+high)/2;
472:         if (rp[t] > col) high = t;
473:         else low = t;
474:       }
475:       for (i=low; i<high; i++) {
476:         if (rp[i] > col) break;
477:         if (rp[i] == col) {
478:           if (is == ADD_VALUES) ap[i] += value;
479:           else ap[i] = value;
480:           low = i + 1;
481:           goto noinsert;
482:         }
483:       }
484:       if (value == 0.0 && ignorezeroentries) goto noinsert;
485:       if (nonew == 1) goto noinsert;
486:       if (nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero at (%D,%D) in the matrix",row,col);
487:       MatSeqXAIJReallocateAIJ(A,A->rmap->n,1,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar);
488:       N = nrow++ - 1; a->nz++; high++;
489:       /* shift up all the later entries in this row */
490:       for (ii=N; ii>=i; ii--) {
491:         rp[ii+1] = rp[ii];
492:         ap[ii+1] = ap[ii];
493:       }
494:       rp[i] = col;
495:       ap[i] = value;
496:       low   = i + 1;
497:       A->nonzerostate++;
498: noinsert:;
499:     }
500:     ailen[row] = nrow;
501:   }
502:   return(0);
503: }


508: PetscErrorCode MatGetValues_SeqAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],PetscScalar v[])
509: {
510:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
511:   PetscInt   *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
512:   PetscInt   *ai = a->i,*ailen = a->ilen;
513:   MatScalar  *ap,*aa = a->a;

516:   for (k=0; k<m; k++) { /* loop over rows */
517:     row = im[k];
518:     if (row < 0) {v += n; continue;} /* SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",row); */
519:     if (row >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap->n-1);
520:     rp   = aj + ai[row]; ap = aa + ai[row];
521:     nrow = ailen[row];
522:     for (l=0; l<n; l++) { /* loop over columns */
523:       if (in[l] < 0) {v++; continue;} /* SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %D",in[l]); */
524:       if (in[l] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->cmap->n-1);
525:       col  = in[l];
526:       high = nrow; low = 0; /* assume unsorted */
527:       while (high-low > 5) {
528:         t = (low+high)/2;
529:         if (rp[t] > col) high = t;
530:         else low = t;
531:       }
532:       for (i=low; i<high; i++) {
533:         if (rp[i] > col) break;
534:         if (rp[i] == col) {
535:           *v++ = ap[i];
536:           goto finished;
537:         }
538:       }
539:       *v++ = 0.0;
540: finished:;
541:     }
542:   }
543:   return(0);
544: }


549: PetscErrorCode MatView_SeqAIJ_Binary(Mat A,PetscViewer viewer)
550: {
551:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
553:   PetscInt       i,*col_lens;
554:   int            fd;
555:   FILE           *file;

558:   PetscViewerBinaryGetDescriptor(viewer,&fd);
559:   PetscMalloc1((4+A->rmap->n),&col_lens);

561:   col_lens[0] = MAT_FILE_CLASSID;
562:   col_lens[1] = A->rmap->n;
563:   col_lens[2] = A->cmap->n;
564:   col_lens[3] = a->nz;

566:   /* store lengths of each row and write (including header) to file */
567:   for (i=0; i<A->rmap->n; i++) {
568:     col_lens[4+i] = a->i[i+1] - a->i[i];
569:   }
570:   PetscBinaryWrite(fd,col_lens,4+A->rmap->n,PETSC_INT,PETSC_TRUE);
571:   PetscFree(col_lens);

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

576:   /* store nonzero values */
577:   PetscBinaryWrite(fd,a->a,a->nz,PETSC_SCALAR,PETSC_FALSE);

579:   PetscViewerBinaryGetInfoPointer(viewer,&file);
580:   if (file) {
581:     fprintf(file,"-matload_block_size %d\n",(int)PetscAbs(A->rmap->bs));
582:   }
583:   return(0);
584: }

586: extern PetscErrorCode MatSeqAIJFactorInfo_Matlab(Mat,PetscViewer);

590: PetscErrorCode MatView_SeqAIJ_ASCII(Mat A,PetscViewer viewer)
591: {
592:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
593:   PetscErrorCode    ierr;
594:   PetscInt          i,j,m = A->rmap->n;
595:   const char        *name;
596:   PetscViewerFormat format;

599:   PetscViewerGetFormat(viewer,&format);
600:   if (format == PETSC_VIEWER_ASCII_MATLAB) {
601:     PetscInt nofinalvalue = 0;
602:     if (m && ((a->i[m] == a->i[m-1]) || (a->j[a->nz-1] != A->cmap->n-1))) {
603:       /* Need a dummy value to ensure the dimension of the matrix. */
604:       nofinalvalue = 1;
605:     }
606:     PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
607:     PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",m,A->cmap->n);
608:     PetscViewerASCIIPrintf(viewer,"%% Nonzeros = %D \n",a->nz);
609: #if defined(PETSC_USE_COMPLEX)
610:     PetscViewerASCIIPrintf(viewer,"zzz = zeros(%D,4);\n",a->nz+nofinalvalue);
611: #else
612:     PetscViewerASCIIPrintf(viewer,"zzz = zeros(%D,3);\n",a->nz+nofinalvalue);
613: #endif
614:     PetscViewerASCIIPrintf(viewer,"zzz = [\n");

616:     for (i=0; i<m; i++) {
617:       for (j=a->i[i]; j<a->i[i+1]; j++) {
618: #if defined(PETSC_USE_COMPLEX)
619:         PetscViewerASCIIPrintf(viewer,"%D %D  %18.16e %18.16e\n",i+1,a->j[j]+1,(double)PetscRealPart(a->a[j]),(double)PetscImaginaryPart(a->a[j]));
620: #else
621:         PetscViewerASCIIPrintf(viewer,"%D %D  %18.16e\n",i+1,a->j[j]+1,(double)a->a[j]);
622: #endif
623:       }
624:     }
625:     if (nofinalvalue) {
626: #if defined(PETSC_USE_COMPLEX)
627:       PetscViewerASCIIPrintf(viewer,"%D %D  %18.16e %18.16e\n",m,A->cmap->n,0.,0.);
628: #else
629:       PetscViewerASCIIPrintf(viewer,"%D %D  %18.16e\n",m,A->cmap->n,0.0);
630: #endif
631:     }
632:     PetscObjectGetName((PetscObject)A,&name);
633:     PetscViewerASCIIPrintf(viewer,"];\n %s = spconvert(zzz);\n",name);
634:     PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
635:   } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) {
636:     return(0);
637:   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
638:     PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
639:     for (i=0; i<m; i++) {
640:       PetscViewerASCIIPrintf(viewer,"row %D:",i);
641:       for (j=a->i[i]; j<a->i[i+1]; j++) {
642: #if defined(PETSC_USE_COMPLEX)
643:         if (PetscImaginaryPart(a->a[j]) > 0.0 && PetscRealPart(a->a[j]) != 0.0) {
644:           PetscViewerASCIIPrintf(viewer," (%D, %g + %g i)",a->j[j],(double)PetscRealPart(a->a[j]),(double)PetscImaginaryPart(a->a[j]));
645:         } else if (PetscImaginaryPart(a->a[j]) < 0.0 && PetscRealPart(a->a[j]) != 0.0) {
646:           PetscViewerASCIIPrintf(viewer," (%D, %g - %g i)",a->j[j],(double)PetscRealPart(a->a[j]),(double)-PetscImaginaryPart(a->a[j]));
647:         } else if (PetscRealPart(a->a[j]) != 0.0) {
648:           PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[j],(double)PetscRealPart(a->a[j]));
649:         }
650: #else
651:         if (a->a[j] != 0.0) {PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[j],(double)a->a[j]);}
652: #endif
653:       }
654:       PetscViewerASCIIPrintf(viewer,"\n");
655:     }
656:     PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
657:   } else if (format == PETSC_VIEWER_ASCII_SYMMODU) {
658:     PetscInt nzd=0,fshift=1,*sptr;
659:     PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
660:     PetscMalloc1((m+1),&sptr);
661:     for (i=0; i<m; i++) {
662:       sptr[i] = nzd+1;
663:       for (j=a->i[i]; j<a->i[i+1]; j++) {
664:         if (a->j[j] >= i) {
665: #if defined(PETSC_USE_COMPLEX)
666:           if (PetscImaginaryPart(a->a[j]) != 0.0 || PetscRealPart(a->a[j]) != 0.0) nzd++;
667: #else
668:           if (a->a[j] != 0.0) nzd++;
669: #endif
670:         }
671:       }
672:     }
673:     sptr[m] = nzd+1;
674:     PetscViewerASCIIPrintf(viewer," %D %D\n\n",m,nzd);
675:     for (i=0; i<m+1; i+=6) {
676:       if (i+4<m) {
677:         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]);
678:       } else if (i+3<m) {
679:         PetscViewerASCIIPrintf(viewer," %D %D %D %D %D\n",sptr[i],sptr[i+1],sptr[i+2],sptr[i+3],sptr[i+4]);
680:       } else if (i+2<m) {
681:         PetscViewerASCIIPrintf(viewer," %D %D %D %D\n",sptr[i],sptr[i+1],sptr[i+2],sptr[i+3]);
682:       } else if (i+1<m) {
683:         PetscViewerASCIIPrintf(viewer," %D %D %D\n",sptr[i],sptr[i+1],sptr[i+2]);
684:       } else if (i<m) {
685:         PetscViewerASCIIPrintf(viewer," %D %D\n",sptr[i],sptr[i+1]);
686:       } else {
687:         PetscViewerASCIIPrintf(viewer," %D\n",sptr[i]);
688:       }
689:     }
690:     PetscViewerASCIIPrintf(viewer,"\n");
691:     PetscFree(sptr);
692:     for (i=0; i<m; i++) {
693:       for (j=a->i[i]; j<a->i[i+1]; j++) {
694:         if (a->j[j] >= i) {PetscViewerASCIIPrintf(viewer," %D ",a->j[j]+fshift);}
695:       }
696:       PetscViewerASCIIPrintf(viewer,"\n");
697:     }
698:     PetscViewerASCIIPrintf(viewer,"\n");
699:     for (i=0; i<m; i++) {
700:       for (j=a->i[i]; j<a->i[i+1]; j++) {
701:         if (a->j[j] >= i) {
702: #if defined(PETSC_USE_COMPLEX)
703:           if (PetscImaginaryPart(a->a[j]) != 0.0 || PetscRealPart(a->a[j]) != 0.0) {
704:             PetscViewerASCIIPrintf(viewer," %18.16e %18.16e ",(double)PetscRealPart(a->a[j]),(double)PetscImaginaryPart(a->a[j]));
705:           }
706: #else
707:           if (a->a[j] != 0.0) {PetscViewerASCIIPrintf(viewer," %18.16e ",(double)a->a[j]);}
708: #endif
709:         }
710:       }
711:       PetscViewerASCIIPrintf(viewer,"\n");
712:     }
713:     PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
714:   } else if (format == PETSC_VIEWER_ASCII_DENSE) {
715:     PetscInt    cnt = 0,jcnt;
716:     PetscScalar value;
717: #if defined(PETSC_USE_COMPLEX)
718:     PetscBool   realonly = PETSC_TRUE;

720:     for (i=0; i<a->i[m]; i++) {
721:       if (PetscImaginaryPart(a->a[i]) != 0.0) {
722:         realonly = PETSC_FALSE;
723:         break;
724:       }
725:     }
726: #endif

728:     PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
729:     for (i=0; i<m; i++) {
730:       jcnt = 0;
731:       for (j=0; j<A->cmap->n; j++) {
732:         if (jcnt < a->i[i+1]-a->i[i] && j == a->j[cnt]) {
733:           value = a->a[cnt++];
734:           jcnt++;
735:         } else {
736:           value = 0.0;
737:         }
738: #if defined(PETSC_USE_COMPLEX)
739:         if (realonly) {
740:           PetscViewerASCIIPrintf(viewer," %7.5e ",(double)PetscRealPart(value));
741:         } else {
742:           PetscViewerASCIIPrintf(viewer," %7.5e+%7.5e i ",(double)PetscRealPart(value),(double)PetscImaginaryPart(value));
743:         }
744: #else
745:         PetscViewerASCIIPrintf(viewer," %7.5e ",(double)value);
746: #endif
747:       }
748:       PetscViewerASCIIPrintf(viewer,"\n");
749:     }
750:     PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
751:   } else if (format == PETSC_VIEWER_ASCII_MATRIXMARKET) {
752:     PetscInt fshift=1;
753:     PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
754: #if defined(PETSC_USE_COMPLEX)
755:     PetscViewerASCIIPrintf(viewer,"%%matrix complex general\n");
756: #else
757:     PetscViewerASCIIPrintf(viewer,"%%matrix real general\n");
758: #endif
759:     PetscViewerASCIIPrintf(viewer,"%D %D %D\n", m, A->cmap->n, a->nz);
760:     for (i=0; i<m; i++) {
761:       for (j=a->i[i]; j<a->i[i+1]; j++) {
762: #if defined(PETSC_USE_COMPLEX)
763:         if (PetscImaginaryPart(a->a[j]) > 0.0) {
764:           PetscViewerASCIIPrintf(viewer,"%D %D, %g %g\n", i+fshift,a->j[j]+fshift,(double)PetscRealPart(a->a[j]),(double)PetscImaginaryPart(a->a[j]));
765:         } else if (PetscImaginaryPart(a->a[j]) < 0.0) {
766:           PetscViewerASCIIPrintf(viewer,"%D %D, %g -%g\n", i+fshift,a->j[j]+fshift,(double)PetscRealPart(a->a[j]),(double)-PetscImaginaryPart(a->a[j]));
767:         } else {
768:           PetscViewerASCIIPrintf(viewer,"%D %D, %g\n", i+fshift,a->j[j]+fshift,(double)PetscRealPart(a->a[j]));
769:         }
770: #else
771:         PetscViewerASCIIPrintf(viewer,"%D %D %g\n", i+fshift, a->j[j]+fshift, (double)a->a[j]);
772: #endif
773:       }
774:     }
775:     PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
776:   } else {
777:     PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
778:     if (A->factortype) {
779:       for (i=0; i<m; i++) {
780:         PetscViewerASCIIPrintf(viewer,"row %D:",i);
781:         /* L part */
782:         for (j=a->i[i]; j<a->i[i+1]; j++) {
783: #if defined(PETSC_USE_COMPLEX)
784:           if (PetscImaginaryPart(a->a[j]) > 0.0) {
785:             PetscViewerASCIIPrintf(viewer," (%D, %g + %g i)",a->j[j],(double)PetscRealPart(a->a[j]),(double)PetscImaginaryPart(a->a[j]));
786:           } else if (PetscImaginaryPart(a->a[j]) < 0.0) {
787:             PetscViewerASCIIPrintf(viewer," (%D, %g - %g i)",a->j[j],(double)PetscRealPart(a->a[j]),(double)(-PetscImaginaryPart(a->a[j])));
788:           } else {
789:             PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[j],(double)PetscRealPart(a->a[j]));
790:           }
791: #else
792:           PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[j],(double)a->a[j]);
793: #endif
794:         }
795:         /* diagonal */
796:         j = a->diag[i];
797: #if defined(PETSC_USE_COMPLEX)
798:         if (PetscImaginaryPart(a->a[j]) > 0.0) {
799:           PetscViewerASCIIPrintf(viewer," (%D, %g + %g i)",a->j[j],(double)PetscRealPart(1.0/a->a[j]),(double)PetscImaginaryPart(1.0/a->a[j]));
800:         } else if (PetscImaginaryPart(a->a[j]) < 0.0) {
801:           PetscViewerASCIIPrintf(viewer," (%D, %g - %g i)",a->j[j],(double)PetscRealPart(1.0/a->a[j]),(double)(-PetscImaginaryPart(1.0/a->a[j])));
802:         } else {
803:           PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[j],(double)PetscRealPart(1.0/a->a[j]));
804:         }
805: #else
806:         PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[j],(double)(1.0/a->a[j]));
807: #endif

809:         /* U part */
810:         for (j=a->diag[i+1]+1; j<a->diag[i]; j++) {
811: #if defined(PETSC_USE_COMPLEX)
812:           if (PetscImaginaryPart(a->a[j]) > 0.0) {
813:             PetscViewerASCIIPrintf(viewer," (%D, %g + %g i)",a->j[j],(double)PetscRealPart(a->a[j]),(double)PetscImaginaryPart(a->a[j]));
814:           } else if (PetscImaginaryPart(a->a[j]) < 0.0) {
815:             PetscViewerASCIIPrintf(viewer," (%D, %g - %g i)",a->j[j],(double)PetscRealPart(a->a[j]),(double)(-PetscImaginaryPart(a->a[j])));
816:           } else {
817:             PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[j],(double)PetscRealPart(a->a[j]));
818:           }
819: #else
820:           PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[j],(double)a->a[j]);
821: #endif
822:         }
823:         PetscViewerASCIIPrintf(viewer,"\n");
824:       }
825:     } else {
826:       for (i=0; i<m; i++) {
827:         PetscViewerASCIIPrintf(viewer,"row %D:",i);
828:         for (j=a->i[i]; j<a->i[i+1]; j++) {
829: #if defined(PETSC_USE_COMPLEX)
830:           if (PetscImaginaryPart(a->a[j]) > 0.0) {
831:             PetscViewerASCIIPrintf(viewer," (%D, %g + %g i)",a->j[j],(double)PetscRealPart(a->a[j]),(double)PetscImaginaryPart(a->a[j]));
832:           } else if (PetscImaginaryPart(a->a[j]) < 0.0) {
833:             PetscViewerASCIIPrintf(viewer," (%D, %g - %g i)",a->j[j],(double)PetscRealPart(a->a[j]),(double)-PetscImaginaryPart(a->a[j]));
834:           } else {
835:             PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[j],(double)PetscRealPart(a->a[j]));
836:           }
837: #else
838:           PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[j],(double)a->a[j]);
839: #endif
840:         }
841:         PetscViewerASCIIPrintf(viewer,"\n");
842:       }
843:     }
844:     PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
845:   }
846:   PetscViewerFlush(viewer);
847:   return(0);
848: }

850: #include <petscdraw.h>
853: PetscErrorCode MatView_SeqAIJ_Draw_Zoom(PetscDraw draw,void *Aa)
854: {
855:   Mat               A  = (Mat) Aa;
856:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
857:   PetscErrorCode    ierr;
858:   PetscInt          i,j,m = A->rmap->n,color;
859:   PetscReal         xl,yl,xr,yr,x_l,x_r,y_l,y_r,maxv = 0.0;
860:   PetscViewer       viewer;
861:   PetscViewerFormat format;

864:   PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);
865:   PetscViewerGetFormat(viewer,&format);

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

870:   if (format != PETSC_VIEWER_DRAW_CONTOUR) {
871:     /* Blue for negative, Cyan for zero and  Red for positive */
872:     color = PETSC_DRAW_BLUE;
873:     for (i=0; i<m; i++) {
874:       y_l = m - i - 1.0; y_r = y_l + 1.0;
875:       for (j=a->i[i]; j<a->i[i+1]; j++) {
876:         x_l = a->j[j]; x_r = x_l + 1.0;
877:         if (PetscRealPart(a->a[j]) >=  0.) continue;
878:         PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
879:       }
880:     }
881:     color = PETSC_DRAW_CYAN;
882:     for (i=0; i<m; i++) {
883:       y_l = m - i - 1.0; y_r = y_l + 1.0;
884:       for (j=a->i[i]; j<a->i[i+1]; j++) {
885:         x_l = a->j[j]; x_r = x_l + 1.0;
886:         if (a->a[j] !=  0.) continue;
887:         PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
888:       }
889:     }
890:     color = PETSC_DRAW_RED;
891:     for (i=0; i<m; i++) {
892:       y_l = m - i - 1.0; y_r = y_l + 1.0;
893:       for (j=a->i[i]; j<a->i[i+1]; j++) {
894:         x_l = a->j[j]; x_r = x_l + 1.0;
895:         if (PetscRealPart(a->a[j]) <=  0.) continue;
896:         PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
897:       }
898:     }
899:   } else {
900:     /* use contour shading to indicate magnitude of values */
901:     /* first determine max of all nonzero values */
902:     PetscInt  nz = a->nz,count;
903:     PetscDraw popup;
904:     PetscReal scale;

906:     for (i=0; i<nz; i++) {
907:       if (PetscAbsScalar(a->a[i]) > maxv) maxv = PetscAbsScalar(a->a[i]);
908:     }
909:     scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv;
910:     PetscDrawGetPopup(draw,&popup);
911:     if (popup) {
912:       PetscDrawScalePopup(popup,0.0,maxv);
913:     }
914:     count = 0;
915:     for (i=0; i<m; i++) {
916:       y_l = m - i - 1.0; y_r = y_l + 1.0;
917:       for (j=a->i[i]; j<a->i[i+1]; j++) {
918:         x_l   = a->j[j]; x_r = x_l + 1.0;
919:         color = PETSC_DRAW_BASIC_COLORS + (PetscInt)(scale*PetscAbsScalar(a->a[count]));
920:         PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
921:         count++;
922:       }
923:     }
924:   }
925:   return(0);
926: }

928: #include <petscdraw.h>
931: PetscErrorCode MatView_SeqAIJ_Draw(Mat A,PetscViewer viewer)
932: {
934:   PetscDraw      draw;
935:   PetscReal      xr,yr,xl,yl,h,w;
936:   PetscBool      isnull;

939:   PetscViewerDrawGetDraw(viewer,0,&draw);
940:   PetscDrawIsNull(draw,&isnull);
941:   if (isnull) return(0);

943:   PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);
944:   xr   = A->cmap->n; yr = A->rmap->n; h = yr/10.0; w = xr/10.0;
945:   xr  += w;    yr += h;  xl = -w;     yl = -h;
946:   PetscDrawSetCoordinates(draw,xl,yl,xr,yr);
947:   PetscDrawZoom(draw,MatView_SeqAIJ_Draw_Zoom,A);
948:   PetscObjectCompose((PetscObject)A,"Zoomviewer",NULL);
949:   return(0);
950: }

954: PetscErrorCode MatView_SeqAIJ(Mat A,PetscViewer viewer)
955: {
957:   PetscBool      iascii,isbinary,isdraw;

960:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
961:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
962:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
963:   if (iascii) {
964:     MatView_SeqAIJ_ASCII(A,viewer);
965:   } else if (isbinary) {
966:     MatView_SeqAIJ_Binary(A,viewer);
967:   } else if (isdraw) {
968:     MatView_SeqAIJ_Draw(A,viewer);
969:   }
970:   MatView_SeqAIJ_Inode(A,viewer);
971:   return(0);
972: }

976: PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode)
977: {
978:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
980:   PetscInt       fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
981:   PetscInt       m      = A->rmap->n,*ip,N,*ailen = a->ilen,rmax = 0;
982:   MatScalar      *aa    = a->a,*ap;
983:   PetscReal      ratio  = 0.6;

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

988:   if (m) rmax = ailen[0]; /* determine row with most nonzeros */
989:   for (i=1; i<m; i++) {
990:     /* move each row back by the amount of empty slots (fshift) before it*/
991:     fshift += imax[i-1] - ailen[i-1];
992:     rmax    = PetscMax(rmax,ailen[i]);
993:     if (fshift) {
994:       ip = aj + ai[i];
995:       ap = aa + ai[i];
996:       N  = ailen[i];
997:       for (j=0; j<N; j++) {
998:         ip[j-fshift] = ip[j];
999:         ap[j-fshift] = ap[j];
1000:       }
1001:     }
1002:     ai[i] = ai[i-1] + ailen[i-1];
1003:   }
1004:   if (m) {
1005:     fshift += imax[m-1] - ailen[m-1];
1006:     ai[m]   = ai[m-1] + ailen[m-1];
1007:   }

1009:   /* reset ilen and imax for each row */
1010:   a->nonzerorowcnt = 0;
1011:   for (i=0; i<m; i++) {
1012:     ailen[i] = imax[i] = ai[i+1] - ai[i];
1013:     a->nonzerorowcnt += ((ai[i+1] - ai[i]) > 0);
1014:   }
1015:   a->nz = ai[m];
1016:   if (fshift && a->nounused == -1) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB, "Unused space detected in matrix: %D X %D, %D unneeded", m, A->cmap->n, fshift);

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

1023:   A->info.mallocs    += a->reallocs;
1024:   a->reallocs         = 0;
1025:   A->info.nz_unneeded = (PetscReal)fshift;
1026:   a->rmax             = rmax;

1028:   MatCheckCompressedRow(A,a->nonzerorowcnt,&a->compressedrow,a->i,m,ratio);
1029:   MatAssemblyEnd_SeqAIJ_Inode(A,mode);
1030:   MatSeqAIJInvalidateDiagonal(A);
1031:   return(0);
1032: }

1036: PetscErrorCode MatRealPart_SeqAIJ(Mat A)
1037: {
1038:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1039:   PetscInt       i,nz = a->nz;
1040:   MatScalar      *aa = a->a;

1044:   for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]);
1045:   MatSeqAIJInvalidateDiagonal(A);
1046:   return(0);
1047: }

1051: PetscErrorCode MatImaginaryPart_SeqAIJ(Mat A)
1052: {
1053:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1054:   PetscInt       i,nz = a->nz;
1055:   MatScalar      *aa = a->a;

1059:   for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]);
1060:   MatSeqAIJInvalidateDiagonal(A);
1061:   return(0);
1062: }

1064: #if defined(PETSC_THREADCOMM_ACTIVE)
1065: PetscErrorCode MatZeroEntries_SeqAIJ_Kernel(PetscInt thread_id,Mat A)
1066: {
1068:   PetscInt       *trstarts=A->rmap->trstarts;
1069:   PetscInt       n,start,end;
1070:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;

1072:   start = trstarts[thread_id];
1073:   end   = trstarts[thread_id+1];
1074:   n     = a->i[end] - a->i[start];
1075:   PetscMemzero(a->a+a->i[start],n*sizeof(PetscScalar));
1076:   return 0;
1077: }

1081: PetscErrorCode MatZeroEntries_SeqAIJ(Mat A)
1082: {

1086:   PetscThreadCommRunKernel(PetscObjectComm((PetscObject)A),(PetscThreadKernel)MatZeroEntries_SeqAIJ_Kernel,1,A);
1087:   MatSeqAIJInvalidateDiagonal(A);
1088:   return(0);
1089: }
1090: #else
1093: PetscErrorCode MatZeroEntries_SeqAIJ(Mat A)
1094: {
1095:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;

1099:   PetscMemzero(a->a,(a->i[A->rmap->n])*sizeof(PetscScalar));
1100:   MatSeqAIJInvalidateDiagonal(A);
1101:   return(0);
1102: }
1103: #endif

1105: extern PetscErrorCode MatDestroy_Redundant(Mat_Redundant **);

1109: PetscErrorCode MatDestroy_SeqAIJ(Mat A)
1110: {
1111:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;

1115: #if defined(PETSC_USE_LOG)
1116:   PetscLogObjectState((PetscObject)A,"Rows=%D, Cols=%D, NZ=%D",A->rmap->n,A->cmap->n,a->nz);
1117: #endif
1118:   MatDestroy_Redundant(&a->redundant);
1119:   MatSeqXAIJFreeAIJ(A,&a->a,&a->j,&a->i);
1120:   ISDestroy(&a->row);
1121:   ISDestroy(&a->col);
1122:   PetscFree(a->diag);
1123:   PetscFree(a->ibdiag);
1124:   PetscFree2(a->imax,a->ilen);
1125:   PetscFree3(a->idiag,a->mdiag,a->ssor_work);
1126:   PetscFree(a->solve_work);
1127:   ISDestroy(&a->icol);
1128:   PetscFree(a->saved_values);
1129:   ISColoringDestroy(&a->coloring);
1130:   PetscFree(a->xtoy);
1131:   MatDestroy(&a->XtoY);
1132:   PetscFree2(a->compressedrow.i,a->compressedrow.rindex);
1133:   PetscFree(a->matmult_abdense);

1135:   MatDestroy_SeqAIJ_Inode(A);
1136:   PetscFree(A->data);

1138:   PetscObjectChangeTypeName((PetscObject)A,0);
1139:   PetscObjectComposeFunction((PetscObject)A,"MatSeqAIJSetColumnIndices_C",NULL);
1140:   PetscObjectComposeFunction((PetscObject)A,"MatStoreValues_C",NULL);
1141:   PetscObjectComposeFunction((PetscObject)A,"MatRetrieveValues_C",NULL);
1142:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqaij_seqsbaij_C",NULL);
1143:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqaij_seqbaij_C",NULL);
1144:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqaij_seqaijperm_C",NULL);
1145:   PetscObjectComposeFunction((PetscObject)A,"MatIsTranspose_C",NULL);
1146:   PetscObjectComposeFunction((PetscObject)A,"MatSeqAIJSetPreallocation_C",NULL);
1147:   PetscObjectComposeFunction((PetscObject)A,"MatSeqAIJSetPreallocationCSR_C",NULL);
1148:   PetscObjectComposeFunction((PetscObject)A,"MatReorderForNonzeroDiagonal_C",NULL);
1149:   return(0);
1150: }

1154: PetscErrorCode MatSetOption_SeqAIJ(Mat A,MatOption op,PetscBool flg)
1155: {
1156:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;

1160:   switch (op) {
1161:   case MAT_ROW_ORIENTED:
1162:     a->roworiented = flg;
1163:     break;
1164:   case MAT_KEEP_NONZERO_PATTERN:
1165:     a->keepnonzeropattern = flg;
1166:     break;
1167:   case MAT_NEW_NONZERO_LOCATIONS:
1168:     a->nonew = (flg ? 0 : 1);
1169:     break;
1170:   case MAT_NEW_NONZERO_LOCATION_ERR:
1171:     a->nonew = (flg ? -1 : 0);
1172:     break;
1173:   case MAT_NEW_NONZERO_ALLOCATION_ERR:
1174:     a->nonew = (flg ? -2 : 0);
1175:     break;
1176:   case MAT_UNUSED_NONZERO_LOCATION_ERR:
1177:     a->nounused = (flg ? -1 : 0);
1178:     break;
1179:   case MAT_IGNORE_ZERO_ENTRIES:
1180:     a->ignorezeroentries = flg;
1181:     break;
1182:   case MAT_SPD:
1183:   case MAT_SYMMETRIC:
1184:   case MAT_STRUCTURALLY_SYMMETRIC:
1185:   case MAT_HERMITIAN:
1186:   case MAT_SYMMETRY_ETERNAL:
1187:     /* These options are handled directly by MatSetOption() */
1188:     break;
1189:   case MAT_NEW_DIAGONALS:
1190:   case MAT_IGNORE_OFF_PROC_ENTRIES:
1191:   case MAT_USE_HASH_TABLE:
1192:     PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
1193:     break;
1194:   case MAT_USE_INODES:
1195:     /* Not an error because MatSetOption_SeqAIJ_Inode handles this one */
1196:     break;
1197:   default:
1198:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %d",op);
1199:   }
1200:   MatSetOption_SeqAIJ_Inode(A,op,flg);
1201:   return(0);
1202: }

1206: PetscErrorCode MatGetDiagonal_SeqAIJ(Mat A,Vec v)
1207: {
1208:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1210:   PetscInt       i,j,n,*ai=a->i,*aj=a->j,nz;
1211:   PetscScalar    *aa=a->a,*x,zero=0.0;

1214:   VecGetLocalSize(v,&n);
1215:   if (n != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");

1217:   if (A->factortype == MAT_FACTOR_ILU || A->factortype == MAT_FACTOR_LU) {
1218:     PetscInt *diag=a->diag;
1219:     VecGetArray(v,&x);
1220:     for (i=0; i<n; i++) x[i] = 1.0/aa[diag[i]];
1221:     VecRestoreArray(v,&x);
1222:     return(0);
1223:   }

1225:   VecSet(v,zero);
1226:   VecGetArray(v,&x);
1227:   for (i=0; i<n; i++) {
1228:     nz = ai[i+1] - ai[i];
1229:     if (!nz) x[i] = 0.0;
1230:     for (j=ai[i]; j<ai[i+1]; j++) {
1231:       if (aj[j] == i) {
1232:         x[i] = aa[j];
1233:         break;
1234:       }
1235:     }
1236:   }
1237:   VecRestoreArray(v,&x);
1238:   return(0);
1239: }

1241: #include <../src/mat/impls/aij/seq/ftn-kernels/fmult.h>
1244: PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
1245: {
1246:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1247:   PetscScalar    *x,*y;
1249:   PetscInt       m = A->rmap->n;
1250: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
1251:   MatScalar         *v;
1252:   PetscScalar       alpha;
1253:   PetscInt          n,i,j,*idx,*ii,*ridx=NULL;
1254:   Mat_CompressedRow cprow    = a->compressedrow;
1255:   PetscBool         usecprow = cprow.use;
1256: #endif

1259:   if (zz != yy) {VecCopy(zz,yy);}
1260:   VecGetArray(xx,&x);
1261:   VecGetArray(yy,&y);

1263: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
1264:   fortranmulttransposeaddaij_(&m,x,a->i,a->j,a->a,y);
1265: #else
1266:   if (usecprow) {
1267:     m    = cprow.nrows;
1268:     ii   = cprow.i;
1269:     ridx = cprow.rindex;
1270:   } else {
1271:     ii = a->i;
1272:   }
1273:   for (i=0; i<m; i++) {
1274:     idx = a->j + ii[i];
1275:     v   = a->a + ii[i];
1276:     n   = ii[i+1] - ii[i];
1277:     if (usecprow) {
1278:       alpha = x[ridx[i]];
1279:     } else {
1280:       alpha = x[i];
1281:     }
1282:     for (j=0; j<n; j++) y[idx[j]] += alpha*v[j];
1283:   }
1284: #endif
1285:   PetscLogFlops(2.0*a->nz);
1286:   VecRestoreArray(xx,&x);
1287:   VecRestoreArray(yy,&y);
1288:   return(0);
1289: }

1293: PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec xx,Vec yy)
1294: {

1298:   VecSet(yy,0.0);
1299:   MatMultTransposeAdd_SeqAIJ(A,xx,yy,yy);
1300:   return(0);
1301: }

1303: #include <../src/mat/impls/aij/seq/ftn-kernels/fmult.h>
1304: #if defined(PETSC_THREADCOMM_ACTIVE)
1305: PetscErrorCode MatMult_SeqAIJ_Kernel(PetscInt thread_id,Mat A,Vec xx,Vec yy)
1306: {
1307:   PetscErrorCode    ierr;
1308:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
1309:   PetscScalar       *y;
1310:   const PetscScalar *x;
1311:   const MatScalar   *aa;
1312:   PetscInt          *trstarts=A->rmap->trstarts;
1313:   PetscInt          n,start,end,i;
1314:   const PetscInt    *aj,*ai;
1315:   PetscScalar       sum;

1317:   VecGetArrayRead(xx,&x);
1318:   VecGetArray(yy,&y);
1319:   start = trstarts[thread_id];
1320:   end   = trstarts[thread_id+1];
1321:   aj    = a->j;
1322:   aa    = a->a;
1323:   ai    = a->i;
1324:   for (i=start; i<end; i++) {
1325:     n   = ai[i+1] - ai[i];
1326:     aj  = a->j + ai[i];
1327:     aa  = a->a + ai[i];
1328:     sum = 0.0;
1329:     PetscSparseDensePlusDot(sum,x,aa,aj,n);
1330:     y[i] = sum;
1331:   }
1332:   VecRestoreArrayRead(xx,&x);
1333:   VecRestoreArray(yy,&y);
1334:   return 0;
1335: }

1339: PetscErrorCode MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
1340: {
1341:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
1342:   PetscScalar       *y;
1343:   const PetscScalar *x;
1344:   const MatScalar   *aa;
1345:   PetscErrorCode    ierr;
1346:   PetscInt          m=A->rmap->n;
1347:   const PetscInt    *aj,*ii,*ridx=NULL;
1348:   PetscInt          n,i;
1349:   PetscScalar       sum;
1350:   PetscBool         usecprow=a->compressedrow.use;

1352: #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
1353: #pragma disjoint(*x,*y,*aa)
1354: #endif

1357:   aj = a->j;
1358:   aa = a->a;
1359:   ii = a->i;
1360:   if (usecprow) { /* use compressed row format */
1361:     VecGetArrayRead(xx,&x);
1362:     VecGetArray(yy,&y);
1363:     PetscMemzero(y,m*sizeof(PetscScalar));
1364:     m    = a->compressedrow.nrows;
1365:     ii   = a->compressedrow.i;
1366:     ridx = a->compressedrow.rindex;
1367:     for (i=0; i<m; i++) {
1368:       n           = ii[i+1] - ii[i];
1369:       aj          = a->j + ii[i];
1370:       aa          = a->a + ii[i];
1371:       sum         = 0.0;
1372:       PetscSparseDensePlusDot(sum,x,aa,aj,n);
1373:       /* for (j=0; j<n; j++) sum += (*aa++)*x[*aj++]; */
1374:       y[*ridx++] = sum;
1375:     }
1376:     VecRestoreArrayRead(xx,&x);
1377:     VecRestoreArray(yy,&y);
1378:   } else { /* do not use compressed row format */
1379: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
1380:     fortranmultaij_(&m,x,ii,aj,aa,y);
1381: #else
1382:     PetscThreadCommRunKernel(PetscObjectComm((PetscObject)A),(PetscThreadKernel)MatMult_SeqAIJ_Kernel,3,A,xx,yy);
1383: #endif
1384:   }
1385:   PetscLogFlops(2.0*a->nz - a->nonzerorowcnt);
1386:   return(0);
1387: }
1388: #else
1391: PetscErrorCode MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
1392: {
1393:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
1394:   PetscScalar       *y;
1395:   const PetscScalar *x;
1396:   const MatScalar   *aa;
1397:   PetscErrorCode    ierr;
1398:   PetscInt          m=A->rmap->n;
1399:   const PetscInt    *aj,*ii,*ridx=NULL;
1400:   PetscInt          n,i;
1401:   PetscScalar       sum;
1402:   PetscBool         usecprow=a->compressedrow.use;

1404: #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
1405: #pragma disjoint(*x,*y,*aa)
1406: #endif

1409:   VecGetArrayRead(xx,&x);
1410:   VecGetArray(yy,&y);
1411:   aj   = a->j;
1412:   aa   = a->a;
1413:   ii   = a->i;
1414:   if (usecprow) { /* use compressed row format */
1415:     PetscMemzero(y,m*sizeof(PetscScalar));
1416:     m    = a->compressedrow.nrows;
1417:     ii   = a->compressedrow.i;
1418:     ridx = a->compressedrow.rindex;
1419:     for (i=0; i<m; i++) {
1420:       n           = ii[i+1] - ii[i];
1421:       aj          = a->j + ii[i];
1422:       aa          = a->a + ii[i];
1423:       sum         = 0.0;
1424:       PetscSparseDensePlusDot(sum,x,aa,aj,n);
1425:       /* for (j=0; j<n; j++) sum += (*aa++)*x[*aj++]; */
1426:       y[*ridx++] = sum;
1427:     }
1428:   } else { /* do not use compressed row format */
1429: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
1430:     fortranmultaij_(&m,x,ii,aj,aa,y);
1431: #else
1432: #if defined(PETSC_THREADCOMM_ACTIVE)
1433:     PetscThreadCommRunKernel(PetscObjectComm((PetscObject)A),(PetscThreadKernel)MatMult_SeqAIJ_Kernel,3,A,xx,yy);
1434: #else
1435:     for (i=0; i<m; i++) {
1436:       n           = ii[i+1] - ii[i];
1437:       aj          = a->j + ii[i];
1438:       aa          = a->a + ii[i];
1439:       sum         = 0.0;
1440:       PetscSparseDensePlusDot(sum,x,aa,aj,n);
1441:       y[i] = sum;
1442:     }
1443: #endif
1444: #endif
1445:   }
1446:   PetscLogFlops(2.0*a->nz - a->nonzerorowcnt);
1447:   VecRestoreArrayRead(xx,&x);
1448:   VecRestoreArray(yy,&y);
1449:   return(0);
1450: }
1451: #endif

1455: PetscErrorCode MatMultMax_SeqAIJ(Mat A,Vec xx,Vec yy)
1456: {
1457:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
1458:   PetscScalar       *y;
1459:   const PetscScalar *x;
1460:   const MatScalar   *aa;
1461:   PetscErrorCode    ierr;
1462:   PetscInt          m=A->rmap->n;
1463:   const PetscInt    *aj,*ii,*ridx=NULL;
1464:   PetscInt          n,i,nonzerorow=0;
1465:   PetscScalar       sum;
1466:   PetscBool         usecprow=a->compressedrow.use;

1468: #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
1469: #pragma disjoint(*x,*y,*aa)
1470: #endif

1473:   VecGetArrayRead(xx,&x);
1474:   VecGetArray(yy,&y);
1475:   aj   = a->j;
1476:   aa   = a->a;
1477:   ii   = a->i;
1478:   if (usecprow) { /* use compressed row format */
1479:     m    = a->compressedrow.nrows;
1480:     ii   = a->compressedrow.i;
1481:     ridx = a->compressedrow.rindex;
1482:     for (i=0; i<m; i++) {
1483:       n           = ii[i+1] - ii[i];
1484:       aj          = a->j + ii[i];
1485:       aa          = a->a + ii[i];
1486:       sum         = 0.0;
1487:       nonzerorow += (n>0);
1488:       PetscSparseDenseMaxDot(sum,x,aa,aj,n);
1489:       /* for (j=0; j<n; j++) sum += (*aa++)*x[*aj++]; */
1490:       y[*ridx++] = sum;
1491:     }
1492:   } else { /* do not use compressed row format */
1493:     for (i=0; i<m; i++) {
1494:       n           = ii[i+1] - ii[i];
1495:       aj          = a->j + ii[i];
1496:       aa          = a->a + ii[i];
1497:       sum         = 0.0;
1498:       nonzerorow += (n>0);
1499:       PetscSparseDenseMaxDot(sum,x,aa,aj,n);
1500:       y[i] = sum;
1501:     }
1502:   }
1503:   PetscLogFlops(2.0*a->nz - nonzerorow);
1504:   VecRestoreArrayRead(xx,&x);
1505:   VecRestoreArray(yy,&y);
1506:   return(0);
1507: }

1511: PetscErrorCode MatMultAddMax_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
1512: {
1513:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
1514:   PetscScalar       *y,*z;
1515:   const PetscScalar *x;
1516:   const MatScalar   *aa;
1517:   PetscErrorCode    ierr;
1518:   PetscInt          m = A->rmap->n,*aj,*ii;
1519:   PetscInt          n,i,*ridx=NULL;
1520:   PetscScalar       sum;
1521:   PetscBool         usecprow=a->compressedrow.use;

1524:   VecGetArrayRead(xx,&x);
1525:   VecGetArray(yy,&y);
1526:   if (zz != yy) {
1527:     VecGetArray(zz,&z);
1528:   } else {
1529:     z = y;
1530:   }

1532:   aj = a->j;
1533:   aa = a->a;
1534:   ii = a->i;
1535:   if (usecprow) { /* use compressed row format */
1536:     if (zz != yy) {
1537:       PetscMemcpy(z,y,m*sizeof(PetscScalar));
1538:     }
1539:     m    = a->compressedrow.nrows;
1540:     ii   = a->compressedrow.i;
1541:     ridx = a->compressedrow.rindex;
1542:     for (i=0; i<m; i++) {
1543:       n   = ii[i+1] - ii[i];
1544:       aj  = a->j + ii[i];
1545:       aa  = a->a + ii[i];
1546:       sum = y[*ridx];
1547:       PetscSparseDenseMaxDot(sum,x,aa,aj,n);
1548:       z[*ridx++] = sum;
1549:     }
1550:   } else { /* do not use compressed row format */
1551:     for (i=0; i<m; i++) {
1552:       n   = ii[i+1] - ii[i];
1553:       aj  = a->j + ii[i];
1554:       aa  = a->a + ii[i];
1555:       sum = y[i];
1556:       PetscSparseDenseMaxDot(sum,x,aa,aj,n);
1557:       z[i] = sum;
1558:     }
1559:   }
1560:   PetscLogFlops(2.0*a->nz);
1561:   VecRestoreArrayRead(xx,&x);
1562:   VecRestoreArray(yy,&y);
1563:   if (zz != yy) {
1564:     VecRestoreArray(zz,&z);
1565:   }
1566:   return(0);
1567: }

1569: #include <../src/mat/impls/aij/seq/ftn-kernels/fmultadd.h>
1572: PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
1573: {
1574:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
1575:   PetscScalar       *y,*z;
1576:   const PetscScalar *x;
1577:   const MatScalar   *aa;
1578:   PetscErrorCode    ierr;
1579:   PetscInt          m = A->rmap->n,*aj,*ii;
1580:   PetscInt          n,i,*ridx=NULL;
1581:   PetscScalar       sum;
1582:   PetscBool         usecprow=a->compressedrow.use;

1585:   VecGetArrayRead(xx,&x);
1586:   VecGetArray(yy,&y);
1587:   if (zz != yy) {
1588:     VecGetArray(zz,&z);
1589:   } else {
1590:     z = y;
1591:   }

1593:   aj = a->j;
1594:   aa = a->a;
1595:   ii = a->i;
1596:   if (usecprow) { /* use compressed row format */
1597:     if (zz != yy) {
1598:       PetscMemcpy(z,y,m*sizeof(PetscScalar));
1599:     }
1600:     m    = a->compressedrow.nrows;
1601:     ii   = a->compressedrow.i;
1602:     ridx = a->compressedrow.rindex;
1603:     for (i=0; i<m; i++) {
1604:       n   = ii[i+1] - ii[i];
1605:       aj  = a->j + ii[i];
1606:       aa  = a->a + ii[i];
1607:       sum = y[*ridx];
1608:       PetscSparseDensePlusDot(sum,x,aa,aj,n);
1609:       z[*ridx++] = sum;
1610:     }
1611:   } else { /* do not use compressed row format */
1612: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
1613:     fortranmultaddaij_(&m,x,ii,aj,aa,y,z);
1614: #else
1615:     for (i=0; i<m; i++) {
1616:       n   = ii[i+1] - ii[i];
1617:       aj  = a->j + ii[i];
1618:       aa  = a->a + ii[i];
1619:       sum = y[i];
1620:       PetscSparseDensePlusDot(sum,x,aa,aj,n);
1621:       z[i] = sum;
1622:     }
1623: #endif
1624:   }
1625:   PetscLogFlops(2.0*a->nz);
1626:   VecRestoreArrayRead(xx,&x);
1627:   VecRestoreArray(yy,&y);
1628:   if (zz != yy) {
1629:     VecRestoreArray(zz,&z);
1630:   }
1631: #if defined(PETSC_HAVE_CUSP)
1632:   /*
1633:   VecView(xx,0);
1634:   VecView(zz,0);
1635:   MatView(A,0);
1636:   */
1637: #endif
1638:   return(0);
1639: }

1641: /*
1642:      Adds diagonal pointers to sparse matrix structure.
1643: */
1646: PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat A)
1647: {
1648:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1650:   PetscInt       i,j,m = A->rmap->n;

1653:   if (!a->diag) {
1654:     PetscMalloc1(m,&a->diag);
1655:     PetscLogObjectMemory((PetscObject)A, m*sizeof(PetscInt));
1656:   }
1657:   for (i=0; i<A->rmap->n; i++) {
1658:     a->diag[i] = a->i[i+1];
1659:     for (j=a->i[i]; j<a->i[i+1]; j++) {
1660:       if (a->j[j] == i) {
1661:         a->diag[i] = j;
1662:         break;
1663:       }
1664:     }
1665:   }
1666:   return(0);
1667: }

1669: /*
1670:      Checks for missing diagonals
1671: */
1674: PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat A,PetscBool  *missing,PetscInt *d)
1675: {
1676:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1677:   PetscInt   *diag,*ii = a->i,i;

1680:   *missing = PETSC_FALSE;
1681:   if (A->rmap->n > 0 && !ii) {
1682:     *missing = PETSC_TRUE;
1683:     if (d) *d = 0;
1684:     PetscInfo(A,"Matrix has no entries therefore is missing diagonal");
1685:   } else {
1686:     diag = a->diag;
1687:     for (i=0; i<A->rmap->n; i++) {
1688:       if (diag[i] >= ii[i+1]) {
1689:         *missing = PETSC_TRUE;
1690:         if (d) *d = i;
1691:         PetscInfo1(A,"Matrix is missing diagonal number %D",i);
1692:         break;
1693:       }
1694:     }
1695:   }
1696:   return(0);
1697: }

1701: PetscErrorCode  MatInvertDiagonal_SeqAIJ(Mat A,PetscScalar omega,PetscScalar fshift)
1702: {
1703:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*) A->data;
1705:   PetscInt       i,*diag,m = A->rmap->n;
1706:   MatScalar      *v = a->a;
1707:   PetscScalar    *idiag,*mdiag;

1710:   if (a->idiagvalid) return(0);
1711:   MatMarkDiagonal_SeqAIJ(A);
1712:   diag = a->diag;
1713:   if (!a->idiag) {
1714:     PetscMalloc3(m,&a->idiag,m,&a->mdiag,m,&a->ssor_work);
1715:     PetscLogObjectMemory((PetscObject)A, 3*m*sizeof(PetscScalar));
1716:     v    = a->a;
1717:   }
1718:   mdiag = a->mdiag;
1719:   idiag = a->idiag;

1721:   if (omega == 1.0 && !PetscAbsScalar(fshift)) {
1722:     for (i=0; i<m; i++) {
1723:       mdiag[i] = v[diag[i]];
1724:       if (!PetscAbsScalar(mdiag[i])) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Zero diagonal on row %D",i);
1725:       idiag[i] = 1.0/v[diag[i]];
1726:     }
1727:     PetscLogFlops(m);
1728:   } else {
1729:     for (i=0; i<m; i++) {
1730:       mdiag[i] = v[diag[i]];
1731:       idiag[i] = omega/(fshift + v[diag[i]]);
1732:     }
1733:     PetscLogFlops(2.0*m);
1734:   }
1735:   a->idiagvalid = PETSC_TRUE;
1736:   return(0);
1737: }

1739: #include <../src/mat/impls/aij/seq/ftn-kernels/frelax.h>
1742: PetscErrorCode MatSOR_SeqAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
1743: {
1744:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
1745:   PetscScalar       *x,d,sum,*t,scale;
1746:   const MatScalar   *v = a->a,*idiag=0,*mdiag;
1747:   const PetscScalar *b, *bs,*xb, *ts;
1748:   PetscErrorCode    ierr;
1749:   PetscInt          n = A->cmap->n,m = A->rmap->n,i;
1750:   const PetscInt    *idx,*diag;

1753:   its = its*lits;

1755:   if (fshift != a->fshift || omega != a->omega) a->idiagvalid = PETSC_FALSE; /* must recompute idiag[] */
1756:   if (!a->idiagvalid) {MatInvertDiagonal_SeqAIJ(A,omega,fshift);}
1757:   a->fshift = fshift;
1758:   a->omega  = omega;

1760:   diag  = a->diag;
1761:   t     = a->ssor_work;
1762:   idiag = a->idiag;
1763:   mdiag = a->mdiag;

1765:   VecGetArray(xx,&x);
1766:   VecGetArrayRead(bb,&b);
1767:   /* We count flops by assuming the upper triangular and lower triangular parts have the same number of nonzeros */
1768:   if (flag == SOR_APPLY_UPPER) {
1769:     /* apply (U + D/omega) to the vector */
1770:     bs = b;
1771:     for (i=0; i<m; i++) {
1772:       d   = fshift + mdiag[i];
1773:       n   = a->i[i+1] - diag[i] - 1;
1774:       idx = a->j + diag[i] + 1;
1775:       v   = a->a + diag[i] + 1;
1776:       sum = b[i]*d/omega;
1777:       PetscSparseDensePlusDot(sum,bs,v,idx,n);
1778:       x[i] = sum;
1779:     }
1780:     VecRestoreArray(xx,&x);
1781:     VecRestoreArrayRead(bb,&b);
1782:     PetscLogFlops(a->nz);
1783:     return(0);
1784:   }

1786:   if (flag == SOR_APPLY_LOWER) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"SOR_APPLY_LOWER is not implemented");
1787:   else if (flag & SOR_EISENSTAT) {
1788:     /* Let  A = L + U + D; where L is lower trianglar,
1789:     U is upper triangular, E = D/omega; This routine applies

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

1793:     to a vector efficiently using Eisenstat's trick.
1794:     */
1795:     scale = (2.0/omega) - 1.0;

1797:     /*  x = (E + U)^{-1} b */
1798:     for (i=m-1; i>=0; i--) {
1799:       n   = a->i[i+1] - diag[i] - 1;
1800:       idx = a->j + diag[i] + 1;
1801:       v   = a->a + diag[i] + 1;
1802:       sum = b[i];
1803:       PetscSparseDenseMinusDot(sum,x,v,idx,n);
1804:       x[i] = sum*idiag[i];
1805:     }

1807:     /*  t = b - (2*E - D)x */
1808:     v = a->a;
1809:     for (i=0; i<m; i++) t[i] = b[i] - scale*(v[*diag++])*x[i];

1811:     /*  t = (E + L)^{-1}t */
1812:     ts   = t;
1813:     diag = a->diag;
1814:     for (i=0; i<m; i++) {
1815:       n   = diag[i] - a->i[i];
1816:       idx = a->j + a->i[i];
1817:       v   = a->a + a->i[i];
1818:       sum = t[i];
1819:       PetscSparseDenseMinusDot(sum,ts,v,idx,n);
1820:       t[i] = sum*idiag[i];
1821:       /*  x = x + t */
1822:       x[i] += t[i];
1823:     }

1825:     PetscLogFlops(6.0*m-1 + 2.0*a->nz);
1826:     VecRestoreArray(xx,&x);
1827:     VecRestoreArrayRead(bb,&b);
1828:     return(0);
1829:   }
1830:   if (flag & SOR_ZERO_INITIAL_GUESS) {
1831:     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
1832:       for (i=0; i<m; i++) {
1833:         n   = diag[i] - a->i[i];
1834:         idx = a->j + a->i[i];
1835:         v   = a->a + a->i[i];
1836:         sum = b[i];
1837:         PetscSparseDenseMinusDot(sum,x,v,idx,n);
1838:         t[i] = sum;
1839:         x[i] = sum*idiag[i];
1840:       }
1841:       xb   = t;
1842:       PetscLogFlops(a->nz);
1843:     } else xb = b;
1844:     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
1845:       for (i=m-1; i>=0; i--) {
1846:         n   = a->i[i+1] - diag[i] - 1;
1847:         idx = a->j + diag[i] + 1;
1848:         v   = a->a + diag[i] + 1;
1849:         sum = xb[i];
1850:         PetscSparseDenseMinusDot(sum,x,v,idx,n);
1851:         if (xb == b) {
1852:           x[i] = sum*idiag[i];
1853:         } else {
1854:           x[i] = (1-omega)*x[i] + sum*idiag[i];  /* omega in idiag */
1855:         }
1856:       }
1857:       PetscLogFlops(a->nz); /* assumes 1/2 in upper */
1858:     }
1859:     its--;
1860:   }
1861:   while (its--) {
1862:     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
1863:       for (i=0; i<m; i++) {
1864:         /* lower */
1865:         n   = diag[i] - a->i[i];
1866:         idx = a->j + a->i[i];
1867:         v   = a->a + a->i[i];
1868:         sum = b[i];
1869:         PetscSparseDenseMinusDot(sum,x,v,idx,n);
1870:         t[i] = sum;             /* save application of the lower-triangular part */
1871:         /* upper */
1872:         n   = a->i[i+1] - diag[i] - 1;
1873:         idx = a->j + diag[i] + 1;
1874:         v   = a->a + diag[i] + 1;
1875:         PetscSparseDenseMinusDot(sum,x,v,idx,n);
1876:         x[i] = (1. - omega)*x[i] + sum*idiag[i]; /* omega in idiag */
1877:       }
1878:       xb   = t;
1879:       PetscLogFlops(2.0*a->nz);
1880:     } else xb = b;
1881:     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
1882:       for (i=m-1; i>=0; i--) {
1883:         sum = xb[i];
1884:         if (xb == b) {
1885:           /* whole matrix (no checkpointing available) */
1886:           n   = a->i[i+1] - a->i[i];
1887:           idx = a->j + a->i[i];
1888:           v   = a->a + a->i[i];
1889:           PetscSparseDenseMinusDot(sum,x,v,idx,n);
1890:           x[i] = (1. - omega)*x[i] + (sum + mdiag[i]*x[i])*idiag[i];
1891:         } else { /* lower-triangular part has been saved, so only apply upper-triangular */
1892:           n   = a->i[i+1] - diag[i] - 1;
1893:           idx = a->j + diag[i] + 1;
1894:           v   = a->a + diag[i] + 1;
1895:           PetscSparseDenseMinusDot(sum,x,v,idx,n);
1896:           x[i] = (1. - omega)*x[i] + sum*idiag[i];  /* omega in idiag */
1897:         }
1898:       }
1899:       if (xb == b) {
1900:         PetscLogFlops(2.0*a->nz);
1901:       } else {
1902:         PetscLogFlops(a->nz); /* assumes 1/2 in upper */
1903:       }
1904:     }
1905:   }
1906:   VecRestoreArray(xx,&x);
1907:   VecRestoreArrayRead(bb,&b);
1908:   return(0);
1909: }


1914: PetscErrorCode MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,MatInfo *info)
1915: {
1916:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;

1919:   info->block_size   = 1.0;
1920:   info->nz_allocated = (double)a->maxnz;
1921:   info->nz_used      = (double)a->nz;
1922:   info->nz_unneeded  = (double)(a->maxnz - a->nz);
1923:   info->assemblies   = (double)A->num_ass;
1924:   info->mallocs      = (double)A->info.mallocs;
1925:   info->memory       = ((PetscObject)A)->mem;
1926:   if (A->factortype) {
1927:     info->fill_ratio_given  = A->info.fill_ratio_given;
1928:     info->fill_ratio_needed = A->info.fill_ratio_needed;
1929:     info->factor_mallocs    = A->info.factor_mallocs;
1930:   } else {
1931:     info->fill_ratio_given  = 0;
1932:     info->fill_ratio_needed = 0;
1933:     info->factor_mallocs    = 0;
1934:   }
1935:   return(0);
1936: }

1940: PetscErrorCode MatZeroRows_SeqAIJ(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
1941: {
1942:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
1943:   PetscInt          i,m = A->rmap->n - 1,d = 0;
1944:   PetscErrorCode    ierr;
1945:   const PetscScalar *xx;
1946:   PetscScalar       *bb;
1947:   PetscBool         missing;

1950:   if (x && b) {
1951:     VecGetArrayRead(x,&xx);
1952:     VecGetArray(b,&bb);
1953:     for (i=0; i<N; i++) {
1954:       if (rows[i] < 0 || rows[i] > m) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range", rows[i]);
1955:       bb[rows[i]] = diag*xx[rows[i]];
1956:     }
1957:     VecRestoreArrayRead(x,&xx);
1958:     VecRestoreArray(b,&bb);
1959:   }

1961:   if (a->keepnonzeropattern) {
1962:     for (i=0; i<N; i++) {
1963:       if (rows[i] < 0 || rows[i] > m) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range", rows[i]);
1964:       PetscMemzero(&a->a[a->i[rows[i]]],a->ilen[rows[i]]*sizeof(PetscScalar));
1965:     }
1966:     if (diag != 0.0) {
1967:       MatMissingDiagonal_SeqAIJ(A,&missing,&d);
1968:       if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry in row %D",d);
1969:       for (i=0; i<N; i++) {
1970:         a->a[a->diag[rows[i]]] = diag;
1971:       }
1972:     }
1973:   } else {
1974:     if (diag != 0.0) {
1975:       for (i=0; i<N; i++) {
1976:         if (rows[i] < 0 || rows[i] > m) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range", rows[i]);
1977:         if (a->ilen[rows[i]] > 0) {
1978:           a->ilen[rows[i]]    = 1;
1979:           a->a[a->i[rows[i]]] = diag;
1980:           a->j[a->i[rows[i]]] = rows[i];
1981:         } else { /* in case row was completely empty */
1982:           MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],&diag,INSERT_VALUES);
1983:         }
1984:       }
1985:     } else {
1986:       for (i=0; i<N; i++) {
1987:         if (rows[i] < 0 || rows[i] > m) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range", rows[i]);
1988:         a->ilen[rows[i]] = 0;
1989:       }
1990:     }
1991:     A->nonzerostate++;
1992:   }
1993:   MatAssemblyEnd_SeqAIJ(A,MAT_FINAL_ASSEMBLY);
1994:   return(0);
1995: }

1999: PetscErrorCode MatZeroRowsColumns_SeqAIJ(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
2000: {
2001:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
2002:   PetscInt          i,j,m = A->rmap->n - 1,d = 0;
2003:   PetscErrorCode    ierr;
2004:   PetscBool         missing,*zeroed,vecs = PETSC_FALSE;
2005:   const PetscScalar *xx;
2006:   PetscScalar       *bb;

2009:   if (x && b) {
2010:     VecGetArrayRead(x,&xx);
2011:     VecGetArray(b,&bb);
2012:     vecs = PETSC_TRUE;
2013:   }
2014:   PetscCalloc1(A->rmap->n,&zeroed);
2015:   for (i=0; i<N; i++) {
2016:     if (rows[i] < 0 || rows[i] > m) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range", rows[i]);
2017:     PetscMemzero(&a->a[a->i[rows[i]]],a->ilen[rows[i]]*sizeof(PetscScalar));

2019:     zeroed[rows[i]] = PETSC_TRUE;
2020:   }
2021:   for (i=0; i<A->rmap->n; i++) {
2022:     if (!zeroed[i]) {
2023:       for (j=a->i[i]; j<a->i[i+1]; j++) {
2024:         if (zeroed[a->j[j]]) {
2025:           if (vecs) bb[i] -= a->a[j]*xx[a->j[j]];
2026:           a->a[j] = 0.0;
2027:         }
2028:       }
2029:     } else if (vecs) bb[i] = diag*xx[i];
2030:   }
2031:   if (x && b) {
2032:     VecRestoreArrayRead(x,&xx);
2033:     VecRestoreArray(b,&bb);
2034:   }
2035:   PetscFree(zeroed);
2036:   if (diag != 0.0) {
2037:     MatMissingDiagonal_SeqAIJ(A,&missing,&d);
2038:     if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry in row %D",d);
2039:     for (i=0; i<N; i++) {
2040:       a->a[a->diag[rows[i]]] = diag;
2041:     }
2042:   }
2043:   MatAssemblyEnd_SeqAIJ(A,MAT_FINAL_ASSEMBLY);
2044:   return(0);
2045: }

2049: PetscErrorCode MatGetRow_SeqAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
2050: {
2051:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
2052:   PetscInt   *itmp;

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

2057:   *nz = a->i[row+1] - a->i[row];
2058:   if (v) *v = a->a + a->i[row];
2059:   if (idx) {
2060:     itmp = a->j + a->i[row];
2061:     if (*nz) *idx = itmp;
2062:     else *idx = 0;
2063:   }
2064:   return(0);
2065: }

2067: /* remove this function? */
2070: PetscErrorCode MatRestoreRow_SeqAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
2071: {
2073:   return(0);
2074: }

2078: PetscErrorCode MatNorm_SeqAIJ(Mat A,NormType type,PetscReal *nrm)
2079: {
2080:   Mat_SeqAIJ     *a  = (Mat_SeqAIJ*)A->data;
2081:   MatScalar      *v  = a->a;
2082:   PetscReal      sum = 0.0;
2084:   PetscInt       i,j;

2087:   if (type == NORM_FROBENIUS) {
2088:     for (i=0; i<a->nz; i++) {
2089:       sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
2090:     }
2091:     *nrm = PetscSqrtReal(sum);
2092:   } else if (type == NORM_1) {
2093:     PetscReal *tmp;
2094:     PetscInt  *jj = a->j;
2095:     PetscCalloc1(A->cmap->n+1,&tmp);
2096:     *nrm = 0.0;
2097:     for (j=0; j<a->nz; j++) {
2098:       tmp[*jj++] += PetscAbsScalar(*v);  v++;
2099:     }
2100:     for (j=0; j<A->cmap->n; j++) {
2101:       if (tmp[j] > *nrm) *nrm = tmp[j];
2102:     }
2103:     PetscFree(tmp);
2104:   } else if (type == NORM_INFINITY) {
2105:     *nrm = 0.0;
2106:     for (j=0; j<A->rmap->n; j++) {
2107:       v   = a->a + a->i[j];
2108:       sum = 0.0;
2109:       for (i=0; i<a->i[j+1]-a->i[j]; i++) {
2110:         sum += PetscAbsScalar(*v); v++;
2111:       }
2112:       if (sum > *nrm) *nrm = sum;
2113:     }
2114:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for two norm");
2115:   return(0);
2116: }

2118: /* Merged from MatGetSymbolicTranspose_SeqAIJ() - replace MatGetSymbolicTranspose_SeqAIJ()? */
2121: PetscErrorCode MatTransposeSymbolic_SeqAIJ(Mat A,Mat *B)
2122: {
2124:   PetscInt       i,j,anzj;
2125:   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b;
2126:   PetscInt       an=A->cmap->N,am=A->rmap->N;
2127:   PetscInt       *ati,*atj,*atfill,*ai=a->i,*aj=a->j;

2130:   /* Allocate space for symbolic transpose info and work array */
2131:   PetscCalloc1((an+1),&ati);
2132:   PetscMalloc1(ai[am],&atj);
2133:   PetscMalloc1(an,&atfill);

2135:   /* Walk through aj and count ## of non-zeros in each row of A^T. */
2136:   /* Note: offset by 1 for fast conversion into csr format. */
2137:   for (i=0;i<ai[am];i++) ati[aj[i]+1] += 1;
2138:   /* Form ati for csr format of A^T. */
2139:   for (i=0;i<an;i++) ati[i+1] += ati[i];

2141:   /* Copy ati into atfill so we have locations of the next free space in atj */
2142:   PetscMemcpy(atfill,ati,an*sizeof(PetscInt));

2144:   /* Walk through A row-wise and mark nonzero entries of A^T. */
2145:   for (i=0;i<am;i++) {
2146:     anzj = ai[i+1] - ai[i];
2147:     for (j=0;j<anzj;j++) {
2148:       atj[atfill[*aj]] = i;
2149:       atfill[*aj++]   += 1;
2150:     }
2151:   }

2153:   /* Clean up temporary space and complete requests. */
2154:   PetscFree(atfill);
2155:   MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),an,am,ati,atj,NULL,B);
2156:   MatSetBlockSizes(*B,PetscAbs(A->cmap->bs),PetscAbs(A->rmap->bs));

2158:   b          = (Mat_SeqAIJ*)((*B)->data);
2159:   b->free_a  = PETSC_FALSE;
2160:   b->free_ij = PETSC_TRUE;
2161:   b->nonew   = 0;
2162:   return(0);
2163: }

2167: PetscErrorCode MatTranspose_SeqAIJ(Mat A,MatReuse reuse,Mat *B)
2168: {
2169:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
2170:   Mat            C;
2172:   PetscInt       i,*aj = a->j,*ai = a->i,m = A->rmap->n,len,*col;
2173:   MatScalar      *array = a->a;

2176:   if (reuse == MAT_REUSE_MATRIX && A == *B && m != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Square matrix only for in-place");

2178:   if (reuse == MAT_INITIAL_MATRIX || *B == A) {
2179:     PetscCalloc1((1+A->cmap->n),&col);

2181:     for (i=0; i<ai[m]; i++) col[aj[i]] += 1;
2182:     MatCreate(PetscObjectComm((PetscObject)A),&C);
2183:     MatSetSizes(C,A->cmap->n,m,A->cmap->n,m);
2184:     MatSetBlockSizes(C,PetscAbs(A->cmap->bs),PetscAbs(A->rmap->bs));
2185:     MatSetType(C,((PetscObject)A)->type_name);
2186:     MatSeqAIJSetPreallocation_SeqAIJ(C,0,col);
2187:     PetscFree(col);
2188:   } else {
2189:     C = *B;
2190:   }

2192:   for (i=0; i<m; i++) {
2193:     len    = ai[i+1]-ai[i];
2194:     MatSetValues_SeqAIJ(C,len,aj,1,&i,array,INSERT_VALUES);
2195:     array += len;
2196:     aj    += len;
2197:   }
2198:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
2199:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

2201:   if (reuse == MAT_INITIAL_MATRIX || *B != A) {
2202:     *B = C;
2203:   } else {
2204:     MatHeaderMerge(A,C);
2205:   }
2206:   return(0);
2207: }

2211: PetscErrorCode  MatIsTranspose_SeqAIJ(Mat A,Mat B,PetscReal tol,PetscBool  *f)
2212: {
2213:   Mat_SeqAIJ     *aij = (Mat_SeqAIJ*) A->data,*bij = (Mat_SeqAIJ*) A->data;
2214:   PetscInt       *adx,*bdx,*aii,*bii,*aptr,*bptr;
2215:   MatScalar      *va,*vb;
2217:   PetscInt       ma,na,mb,nb, i;

2220:   bij = (Mat_SeqAIJ*) B->data;

2222:   MatGetSize(A,&ma,&na);
2223:   MatGetSize(B,&mb,&nb);
2224:   if (ma!=nb || na!=mb) {
2225:     *f = PETSC_FALSE;
2226:     return(0);
2227:   }
2228:   aii  = aij->i; bii = bij->i;
2229:   adx  = aij->j; bdx = bij->j;
2230:   va   = aij->a; vb = bij->a;
2231:   PetscMalloc1(ma,&aptr);
2232:   PetscMalloc1(mb,&bptr);
2233:   for (i=0; i<ma; i++) aptr[i] = aii[i];
2234:   for (i=0; i<mb; i++) bptr[i] = bii[i];

2236:   *f = PETSC_TRUE;
2237:   for (i=0; i<ma; i++) {
2238:     while (aptr[i]<aii[i+1]) {
2239:       PetscInt    idc,idr;
2240:       PetscScalar vc,vr;
2241:       /* column/row index/value */
2242:       idc = adx[aptr[i]];
2243:       idr = bdx[bptr[idc]];
2244:       vc  = va[aptr[i]];
2245:       vr  = vb[bptr[idc]];
2246:       if (i!=idr || PetscAbsScalar(vc-vr) > tol) {
2247:         *f = PETSC_FALSE;
2248:         goto done;
2249:       } else {
2250:         aptr[i]++;
2251:         if (B || i!=idc) bptr[idc]++;
2252:       }
2253:     }
2254:   }
2255: done:
2256:   PetscFree(aptr);
2257:   PetscFree(bptr);
2258:   return(0);
2259: }

2263: PetscErrorCode  MatIsHermitianTranspose_SeqAIJ(Mat A,Mat B,PetscReal tol,PetscBool  *f)
2264: {
2265:   Mat_SeqAIJ     *aij = (Mat_SeqAIJ*) A->data,*bij = (Mat_SeqAIJ*) A->data;
2266:   PetscInt       *adx,*bdx,*aii,*bii,*aptr,*bptr;
2267:   MatScalar      *va,*vb;
2269:   PetscInt       ma,na,mb,nb, i;

2272:   bij = (Mat_SeqAIJ*) B->data;

2274:   MatGetSize(A,&ma,&na);
2275:   MatGetSize(B,&mb,&nb);
2276:   if (ma!=nb || na!=mb) {
2277:     *f = PETSC_FALSE;
2278:     return(0);
2279:   }
2280:   aii  = aij->i; bii = bij->i;
2281:   adx  = aij->j; bdx = bij->j;
2282:   va   = aij->a; vb = bij->a;
2283:   PetscMalloc1(ma,&aptr);
2284:   PetscMalloc1(mb,&bptr);
2285:   for (i=0; i<ma; i++) aptr[i] = aii[i];
2286:   for (i=0; i<mb; i++) bptr[i] = bii[i];

2288:   *f = PETSC_TRUE;
2289:   for (i=0; i<ma; i++) {
2290:     while (aptr[i]<aii[i+1]) {
2291:       PetscInt    idc,idr;
2292:       PetscScalar vc,vr;
2293:       /* column/row index/value */
2294:       idc = adx[aptr[i]];
2295:       idr = bdx[bptr[idc]];
2296:       vc  = va[aptr[i]];
2297:       vr  = vb[bptr[idc]];
2298:       if (i!=idr || PetscAbsScalar(vc-PetscConj(vr)) > tol) {
2299:         *f = PETSC_FALSE;
2300:         goto done;
2301:       } else {
2302:         aptr[i]++;
2303:         if (B || i!=idc) bptr[idc]++;
2304:       }
2305:     }
2306:   }
2307: done:
2308:   PetscFree(aptr);
2309:   PetscFree(bptr);
2310:   return(0);
2311: }

2315: PetscErrorCode MatIsSymmetric_SeqAIJ(Mat A,PetscReal tol,PetscBool  *f)
2316: {

2320:   MatIsTranspose_SeqAIJ(A,A,tol,f);
2321:   return(0);
2322: }

2326: PetscErrorCode MatIsHermitian_SeqAIJ(Mat A,PetscReal tol,PetscBool  *f)
2327: {

2331:   MatIsHermitianTranspose_SeqAIJ(A,A,tol,f);
2332:   return(0);
2333: }

2337: PetscErrorCode MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr)
2338: {
2339:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
2340:   PetscScalar    *l,*r,x;
2341:   MatScalar      *v;
2343:   PetscInt       i,j,m = A->rmap->n,n = A->cmap->n,M,nz = a->nz,*jj;

2346:   if (ll) {
2347:     /* The local size is used so that VecMPI can be passed to this routine
2348:        by MatDiagonalScale_MPIAIJ */
2349:     VecGetLocalSize(ll,&m);
2350:     if (m != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length");
2351:     VecGetArray(ll,&l);
2352:     v    = a->a;
2353:     for (i=0; i<m; i++) {
2354:       x = l[i];
2355:       M = a->i[i+1] - a->i[i];
2356:       for (j=0; j<M; j++) (*v++) *= x;
2357:     }
2358:     VecRestoreArray(ll,&l);
2359:     PetscLogFlops(nz);
2360:   }
2361:   if (rr) {
2362:     VecGetLocalSize(rr,&n);
2363:     if (n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vector wrong length");
2364:     VecGetArray(rr,&r);
2365:     v    = a->a; jj = a->j;
2366:     for (i=0; i<nz; i++) (*v++) *= r[*jj++];
2367:     VecRestoreArray(rr,&r);
2368:     PetscLogFlops(nz);
2369:   }
2370:   MatSeqAIJInvalidateDiagonal(A);
2371:   return(0);
2372: }

2376: PetscErrorCode MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,PetscInt csize,MatReuse scall,Mat *B)
2377: {
2378:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data,*c;
2380:   PetscInt       *smap,i,k,kstart,kend,oldcols = A->cmap->n,*lens;
2381:   PetscInt       row,mat_i,*mat_j,tcol,first,step,*mat_ilen,sum,lensi;
2382:   const PetscInt *irow,*icol;
2383:   PetscInt       nrows,ncols;
2384:   PetscInt       *starts,*j_new,*i_new,*aj = a->j,*ai = a->i,ii,*ailen = a->ilen;
2385:   MatScalar      *a_new,*mat_a;
2386:   Mat            C;
2387:   PetscBool      stride,sorted;

2390:   ISSorted(isrow,&sorted);
2391:   if (!sorted) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"ISrow is not sorted");
2392:   ISSorted(iscol,&sorted);
2393:   if (!sorted) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"IScol is not sorted");

2395:   ISGetIndices(isrow,&irow);
2396:   ISGetLocalSize(isrow,&nrows);
2397:   ISGetLocalSize(iscol,&ncols);

2399:   PetscObjectTypeCompare((PetscObject)iscol,ISSTRIDE,&stride);
2400:   if (stride) {
2401:     ISStrideGetInfo(iscol,&first,&step);
2402:   } else {
2403:     first = 0;
2404:     step  = 0;
2405:   }
2406:   if (stride && step == 1) {
2407:     /* special case of contiguous rows */
2408:     PetscMalloc2(nrows,&lens,nrows,&starts);
2409:     /* loop over new rows determining lens and starting points */
2410:     for (i=0; i<nrows; i++) {
2411:       kstart = ai[irow[i]];
2412:       kend   = kstart + ailen[irow[i]];
2413:       starts[i] = kstart;
2414:       for (k=kstart; k<kend; k++) {
2415:         if (aj[k] >= first) {
2416:           starts[i] = k;
2417:           break;
2418:         }
2419:       }
2420:       sum = 0;
2421:       while (k < kend) {
2422:         if (aj[k++] >= first+ncols) break;
2423:         sum++;
2424:       }
2425:       lens[i] = sum;
2426:     }
2427:     /* create submatrix */
2428:     if (scall == MAT_REUSE_MATRIX) {
2429:       PetscInt n_cols,n_rows;
2430:       MatGetSize(*B,&n_rows,&n_cols);
2431:       if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size");
2432:       MatZeroEntries(*B);
2433:       C    = *B;
2434:     } else {
2435:       PetscInt rbs,cbs;
2436:       MatCreate(PetscObjectComm((PetscObject)A),&C);
2437:       MatSetSizes(C,nrows,ncols,PETSC_DETERMINE,PETSC_DETERMINE);
2438:       ISGetBlockSize(isrow,&rbs);
2439:       ISGetBlockSize(iscol,&cbs);
2440:       MatSetBlockSizes(C,rbs,cbs);
2441:       MatSetType(C,((PetscObject)A)->type_name);
2442:       MatSeqAIJSetPreallocation_SeqAIJ(C,0,lens);
2443:     }
2444:     c = (Mat_SeqAIJ*)C->data;

2446:     /* loop over rows inserting into submatrix */
2447:     a_new = c->a;
2448:     j_new = c->j;
2449:     i_new = c->i;

2451:     for (i=0; i<nrows; i++) {
2452:       ii    = starts[i];
2453:       lensi = lens[i];
2454:       for (k=0; k<lensi; k++) {
2455:         *j_new++ = aj[ii+k] - first;
2456:       }
2457:       PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(PetscScalar));
2458:       a_new     += lensi;
2459:       i_new[i+1] = i_new[i] + lensi;
2460:       c->ilen[i] = lensi;
2461:     }
2462:     PetscFree2(lens,starts);
2463:   } else {
2464:     ISGetIndices(iscol,&icol);
2465:     PetscCalloc1(oldcols,&smap);
2466:     PetscMalloc1((1+nrows),&lens);
2467:     for (i=0; i<ncols; i++) {
2468: #if defined(PETSC_USE_DEBUG)
2469:       if (icol[i] >= oldcols) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Requesting column beyond largest column icol[%D] %D <= A->cmap->n %D",i,icol[i],oldcols);
2470: #endif
2471:       smap[icol[i]] = i+1;
2472:     }

2474:     /* determine lens of each row */
2475:     for (i=0; i<nrows; i++) {
2476:       kstart  = ai[irow[i]];
2477:       kend    = kstart + a->ilen[irow[i]];
2478:       lens[i] = 0;
2479:       for (k=kstart; k<kend; k++) {
2480:         if (smap[aj[k]]) {
2481:           lens[i]++;
2482:         }
2483:       }
2484:     }
2485:     /* Create and fill new matrix */
2486:     if (scall == MAT_REUSE_MATRIX) {
2487:       PetscBool equal;

2489:       c = (Mat_SeqAIJ*)((*B)->data);
2490:       if ((*B)->rmap->n  != nrows || (*B)->cmap->n != ncols) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size");
2491:       PetscMemcmp(c->ilen,lens,(*B)->rmap->n*sizeof(PetscInt),&equal);
2492:       if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros");
2493:       PetscMemzero(c->ilen,(*B)->rmap->n*sizeof(PetscInt));
2494:       C    = *B;
2495:     } else {
2496:       PetscInt rbs,cbs;
2497:       MatCreate(PetscObjectComm((PetscObject)A),&C);
2498:       MatSetSizes(C,nrows,ncols,PETSC_DETERMINE,PETSC_DETERMINE);
2499:       ISGetBlockSize(isrow,&rbs);
2500:       ISGetBlockSize(iscol,&cbs);
2501:       MatSetBlockSizes(C,rbs,cbs);
2502:       MatSetType(C,((PetscObject)A)->type_name);
2503:       MatSeqAIJSetPreallocation_SeqAIJ(C,0,lens);
2504:     }
2505:     c = (Mat_SeqAIJ*)(C->data);
2506:     for (i=0; i<nrows; i++) {
2507:       row      = irow[i];
2508:       kstart   = ai[row];
2509:       kend     = kstart + a->ilen[row];
2510:       mat_i    = c->i[i];
2511:       mat_j    = c->j + mat_i;
2512:       mat_a    = c->a + mat_i;
2513:       mat_ilen = c->ilen + i;
2514:       for (k=kstart; k<kend; k++) {
2515:         if ((tcol=smap[a->j[k]])) {
2516:           *mat_j++ = tcol - 1;
2517:           *mat_a++ = a->a[k];
2518:           (*mat_ilen)++;

2520:         }
2521:       }
2522:     }
2523:     /* Free work space */
2524:     ISRestoreIndices(iscol,&icol);
2525:     PetscFree(smap);
2526:     PetscFree(lens);
2527:   }
2528:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
2529:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

2531:   ISRestoreIndices(isrow,&irow);
2532:   *B   = C;
2533:   return(0);
2534: }

2538: PetscErrorCode  MatGetMultiProcBlock_SeqAIJ(Mat mat,MPI_Comm subComm,MatReuse scall,Mat *subMat)
2539: {
2541:   Mat            B;

2544:   if (scall == MAT_INITIAL_MATRIX) {
2545:     MatCreate(subComm,&B);
2546:     MatSetSizes(B,mat->rmap->n,mat->cmap->n,mat->rmap->n,mat->cmap->n);
2547:     MatSetBlockSizesFromMats(B,mat,mat);
2548:     MatSetType(B,MATSEQAIJ);
2549:     MatDuplicateNoCreate_SeqAIJ(B,mat,MAT_COPY_VALUES,PETSC_TRUE);
2550:     *subMat = B;
2551:   } else {
2552:     MatCopy_SeqAIJ(mat,*subMat,SAME_NONZERO_PATTERN);
2553:   }
2554:   return(0);
2555: }

2559: PetscErrorCode MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,const MatFactorInfo *info)
2560: {
2561:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)inA->data;
2563:   Mat            outA;
2564:   PetscBool      row_identity,col_identity;

2567:   if (info->levels != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only levels=0 supported for in-place ilu");

2569:   ISIdentity(row,&row_identity);
2570:   ISIdentity(col,&col_identity);

2572:   outA             = inA;
2573:   outA->factortype = MAT_FACTOR_LU;

2575:   PetscObjectReference((PetscObject)row);
2576:   ISDestroy(&a->row);

2578:   a->row = row;

2580:   PetscObjectReference((PetscObject)col);
2581:   ISDestroy(&a->col);

2583:   a->col = col;

2585:   /* Create the inverse permutation so that it can be used in MatLUFactorNumeric() */
2586:   ISDestroy(&a->icol);
2587:   ISInvertPermutation(col,PETSC_DECIDE,&a->icol);
2588:   PetscLogObjectParent((PetscObject)inA,(PetscObject)a->icol);

2590:   if (!a->solve_work) { /* this matrix may have been factored before */
2591:     PetscMalloc1((inA->rmap->n+1),&a->solve_work);
2592:     PetscLogObjectMemory((PetscObject)inA, (inA->rmap->n+1)*sizeof(PetscScalar));
2593:   }

2595:   MatMarkDiagonal_SeqAIJ(inA);
2596:   if (row_identity && col_identity) {
2597:     MatLUFactorNumeric_SeqAIJ_inplace(outA,inA,info);
2598:   } else {
2599:     MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(outA,inA,info);
2600:   }
2601:   return(0);
2602: }

2606: PetscErrorCode MatScale_SeqAIJ(Mat inA,PetscScalar alpha)
2607: {
2608:   Mat_SeqAIJ     *a     = (Mat_SeqAIJ*)inA->data;
2609:   PetscScalar    oalpha = alpha;
2611:   PetscBLASInt   one = 1,bnz;

2614:   PetscBLASIntCast(a->nz,&bnz);
2615:   PetscStackCallBLAS("BLASscal",BLASscal_(&bnz,&oalpha,a->a,&one));
2616:   PetscLogFlops(a->nz);
2617:   MatSeqAIJInvalidateDiagonal(inA);
2618:   return(0);
2619: }

2623: PetscErrorCode MatGetSubMatrices_SeqAIJ(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
2624: {
2626:   PetscInt       i;

2629:   if (scall == MAT_INITIAL_MATRIX) {
2630:     PetscMalloc1((n+1),B);
2631:   }

2633:   for (i=0; i<n; i++) {
2634:     MatGetSubMatrix_SeqAIJ(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);
2635:   }
2636:   return(0);
2637: }

2641: PetscErrorCode MatIncreaseOverlap_SeqAIJ(Mat A,PetscInt is_max,IS is[],PetscInt ov)
2642: {
2643:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
2645:   PetscInt       row,i,j,k,l,m,n,*nidx,isz,val;
2646:   const PetscInt *idx;
2647:   PetscInt       start,end,*ai,*aj;
2648:   PetscBT        table;

2651:   m  = A->rmap->n;
2652:   ai = a->i;
2653:   aj = a->j;

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

2657:   PetscMalloc1((m+1),&nidx);
2658:   PetscBTCreate(m,&table);

2660:   for (i=0; i<is_max; i++) {
2661:     /* Initialize the two local arrays */
2662:     isz  = 0;
2663:     PetscBTMemzero(m,table);

2665:     /* Extract the indices, assume there can be duplicate entries */
2666:     ISGetIndices(is[i],&idx);
2667:     ISGetLocalSize(is[i],&n);

2669:     /* Enter these into the temp arrays. I.e., mark table[row], enter row into new index */
2670:     for (j=0; j<n; ++j) {
2671:       if (!PetscBTLookupSet(table,idx[j])) nidx[isz++] = idx[j];
2672:     }
2673:     ISRestoreIndices(is[i],&idx);
2674:     ISDestroy(&is[i]);

2676:     k = 0;
2677:     for (j=0; j<ov; j++) { /* for each overlap */
2678:       n = isz;
2679:       for (; k<n; k++) { /* do only those rows in nidx[k], which are not done yet */
2680:         row   = nidx[k];
2681:         start = ai[row];
2682:         end   = ai[row+1];
2683:         for (l = start; l<end; l++) {
2684:           val = aj[l];
2685:           if (!PetscBTLookupSet(table,val)) nidx[isz++] = val;
2686:         }
2687:       }
2688:     }
2689:     ISCreateGeneral(PETSC_COMM_SELF,isz,nidx,PETSC_COPY_VALUES,(is+i));
2690:   }
2691:   PetscBTDestroy(&table);
2692:   PetscFree(nidx);
2693:   return(0);
2694: }

2696: /* -------------------------------------------------------------- */
2699: PetscErrorCode MatPermute_SeqAIJ(Mat A,IS rowp,IS colp,Mat *B)
2700: {
2701:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
2703:   PetscInt       i,nz = 0,m = A->rmap->n,n = A->cmap->n;
2704:   const PetscInt *row,*col;
2705:   PetscInt       *cnew,j,*lens;
2706:   IS             icolp,irowp;
2707:   PetscInt       *cwork = NULL;
2708:   PetscScalar    *vwork = NULL;

2711:   ISInvertPermutation(rowp,PETSC_DECIDE,&irowp);
2712:   ISGetIndices(irowp,&row);
2713:   ISInvertPermutation(colp,PETSC_DECIDE,&icolp);
2714:   ISGetIndices(icolp,&col);

2716:   /* determine lengths of permuted rows */
2717:   PetscMalloc1((m+1),&lens);
2718:   for (i=0; i<m; i++) lens[row[i]] = a->i[i+1] - a->i[i];
2719:   MatCreate(PetscObjectComm((PetscObject)A),B);
2720:   MatSetSizes(*B,m,n,m,n);
2721:   MatSetBlockSizesFromMats(*B,A,A);
2722:   MatSetType(*B,((PetscObject)A)->type_name);
2723:   MatSeqAIJSetPreallocation_SeqAIJ(*B,0,lens);
2724:   PetscFree(lens);

2726:   PetscMalloc1(n,&cnew);
2727:   for (i=0; i<m; i++) {
2728:     MatGetRow_SeqAIJ(A,i,&nz,&cwork,&vwork);
2729:     for (j=0; j<nz; j++) cnew[j] = col[cwork[j]];
2730:     MatSetValues_SeqAIJ(*B,1,&row[i],nz,cnew,vwork,INSERT_VALUES);
2731:     MatRestoreRow_SeqAIJ(A,i,&nz,&cwork,&vwork);
2732:   }
2733:   PetscFree(cnew);

2735:   (*B)->assembled = PETSC_FALSE;

2737:   MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);
2738:   MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);
2739:   ISRestoreIndices(irowp,&row);
2740:   ISRestoreIndices(icolp,&col);
2741:   ISDestroy(&irowp);
2742:   ISDestroy(&icolp);
2743:   return(0);
2744: }

2748: PetscErrorCode MatCopy_SeqAIJ(Mat A,Mat B,MatStructure str)
2749: {

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

2758:     if (a->i[A->rmap->n] != b->i[B->rmap->n]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of nonzeros in two matrices are different");
2759:     PetscMemcpy(b->a,a->a,(a->i[A->rmap->n])*sizeof(PetscScalar));
2760:   } else {
2761:     MatCopy_Basic(A,B,str);
2762:   }
2763:   return(0);
2764: }

2768: PetscErrorCode MatSetUp_SeqAIJ(Mat A)
2769: {

2773:    MatSeqAIJSetPreallocation_SeqAIJ(A,PETSC_DEFAULT,0);
2774:   return(0);
2775: }

2779: PetscErrorCode MatSeqAIJGetArray_SeqAIJ(Mat A,PetscScalar *array[])
2780: {
2781:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;

2784:   *array = a->a;
2785:   return(0);
2786: }

2790: PetscErrorCode MatSeqAIJRestoreArray_SeqAIJ(Mat A,PetscScalar *array[])
2791: {
2793:   return(0);
2794: }

2796: /*
2797:    Computes the number of nonzeros per row needed for preallocation when X and Y
2798:    have different nonzero structure.
2799: */
2802: PetscErrorCode MatAXPYGetPreallocation_SeqX_private(PetscInt m,const PetscInt *xi,const PetscInt *xj,const PetscInt *yi,const PetscInt *yj,PetscInt *nnz)
2803: {
2804:   PetscInt       i,j,k,nzx,nzy;

2807:   /* Set the number of nonzeros in the new matrix */
2808:   for (i=0; i<m; i++) {
2809:     const PetscInt *xjj = xj+xi[i],*yjj = yj+yi[i];
2810:     nzx = xi[i+1] - xi[i];
2811:     nzy = yi[i+1] - yi[i];
2812:     nnz[i] = 0;
2813:     for (j=0,k=0; j<nzx; j++) {                   /* Point in X */
2814:       for (; k<nzy && yjj[k]<xjj[j]; k++) nnz[i]++; /* Catch up to X */
2815:       if (k<nzy && yjj[k]==xjj[j]) k++;             /* Skip duplicate */
2816:       nnz[i]++;
2817:     }
2818:     for (; k<nzy; k++) nnz[i]++;
2819:   }
2820:   return(0);
2821: }

2825: PetscErrorCode MatAXPYGetPreallocation_SeqAIJ(Mat Y,Mat X,PetscInt *nnz)
2826: {
2827:   PetscInt       m = Y->rmap->N;
2828:   Mat_SeqAIJ     *x = (Mat_SeqAIJ*)X->data;
2829:   Mat_SeqAIJ     *y = (Mat_SeqAIJ*)Y->data;

2833:   /* Set the number of nonzeros in the new matrix */
2834:   MatAXPYGetPreallocation_SeqX_private(m,x->i,x->j,y->i,y->j,nnz);
2835:   return(0);
2836: }

2840: PetscErrorCode MatAXPY_SeqAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str)
2841: {
2843:   Mat_SeqAIJ     *x = (Mat_SeqAIJ*)X->data,*y = (Mat_SeqAIJ*)Y->data;
2844:   PetscBLASInt   one=1,bnz;

2847:   PetscBLASIntCast(x->nz,&bnz);
2848:   if (str == SAME_NONZERO_PATTERN) {
2849:     PetscScalar alpha = a;
2850:     PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&bnz,&alpha,x->a,&one,y->a,&one));
2851:     MatSeqAIJInvalidateDiagonal(Y);
2852:     PetscObjectStateIncrease((PetscObject)Y);
2853:   } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
2854:     MatAXPY_Basic(Y,a,X,str);
2855:   } else {
2856:     Mat      B;
2857:     PetscInt *nnz;
2858:     PetscMalloc1(Y->rmap->N,&nnz);
2859:     MatCreate(PetscObjectComm((PetscObject)Y),&B);
2860:     PetscObjectSetName((PetscObject)B,((PetscObject)Y)->name);
2861:     MatSetSizes(B,Y->rmap->n,Y->cmap->n,Y->rmap->N,Y->cmap->N);
2862:     MatSetBlockSizesFromMats(B,Y,Y);
2863:     MatSetType(B,(MatType) ((PetscObject)Y)->type_name);
2864:     MatAXPYGetPreallocation_SeqAIJ(Y,X,nnz);
2865:     MatSeqAIJSetPreallocation(B,0,nnz);
2866:     MatAXPY_BasicWithPreallocation(B,Y,a,X,str);
2867:     MatHeaderReplace(Y,B);
2868:     PetscFree(nnz);
2869:   }
2870:   return(0);
2871: }

2875: PetscErrorCode  MatConjugate_SeqAIJ(Mat mat)
2876: {
2877: #if defined(PETSC_USE_COMPLEX)
2878:   Mat_SeqAIJ  *aij = (Mat_SeqAIJ*)mat->data;
2879:   PetscInt    i,nz;
2880:   PetscScalar *a;

2883:   nz = aij->nz;
2884:   a  = aij->a;
2885:   for (i=0; i<nz; i++) a[i] = PetscConj(a[i]);
2886: #else
2888: #endif
2889:   return(0);
2890: }

2894: PetscErrorCode MatGetRowMaxAbs_SeqAIJ(Mat A,Vec v,PetscInt idx[])
2895: {
2896:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
2898:   PetscInt       i,j,m = A->rmap->n,*ai,*aj,ncols,n;
2899:   PetscReal      atmp;
2900:   PetscScalar    *x;
2901:   MatScalar      *aa;

2904:   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2905:   aa = a->a;
2906:   ai = a->i;
2907:   aj = a->j;

2909:   VecSet(v,0.0);
2910:   VecGetArray(v,&x);
2911:   VecGetLocalSize(v,&n);
2912:   if (n != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
2913:   for (i=0; i<m; i++) {
2914:     ncols = ai[1] - ai[0]; ai++;
2915:     x[i]  = 0.0;
2916:     for (j=0; j<ncols; j++) {
2917:       atmp = PetscAbsScalar(*aa);
2918:       if (PetscAbsScalar(x[i]) < atmp) {x[i] = atmp; if (idx) idx[i] = *aj;}
2919:       aa++; aj++;
2920:     }
2921:   }
2922:   VecRestoreArray(v,&x);
2923:   return(0);
2924: }

2928: PetscErrorCode MatGetRowMax_SeqAIJ(Mat A,Vec v,PetscInt idx[])
2929: {
2930:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
2932:   PetscInt       i,j,m = A->rmap->n,*ai,*aj,ncols,n;
2933:   PetscScalar    *x;
2934:   MatScalar      *aa;

2937:   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2938:   aa = a->a;
2939:   ai = a->i;
2940:   aj = a->j;

2942:   VecSet(v,0.0);
2943:   VecGetArray(v,&x);
2944:   VecGetLocalSize(v,&n);
2945:   if (n != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
2946:   for (i=0; i<m; i++) {
2947:     ncols = ai[1] - ai[0]; ai++;
2948:     if (ncols == A->cmap->n) { /* row is dense */
2949:       x[i] = *aa; if (idx) idx[i] = 0;
2950:     } else {  /* row is sparse so already KNOW maximum is 0.0 or higher */
2951:       x[i] = 0.0;
2952:       if (idx) {
2953:         idx[i] = 0; /* in case ncols is zero */
2954:         for (j=0;j<ncols;j++) { /* find first implicit 0.0 in the row */
2955:           if (aj[j] > j) {
2956:             idx[i] = j;
2957:             break;
2958:           }
2959:         }
2960:       }
2961:     }
2962:     for (j=0; j<ncols; j++) {
2963:       if (PetscRealPart(x[i]) < PetscRealPart(*aa)) {x[i] = *aa; if (idx) idx[i] = *aj;}
2964:       aa++; aj++;
2965:     }
2966:   }
2967:   VecRestoreArray(v,&x);
2968:   return(0);
2969: }

2973: PetscErrorCode MatGetRowMinAbs_SeqAIJ(Mat A,Vec v,PetscInt idx[])
2974: {
2975:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
2977:   PetscInt       i,j,m = A->rmap->n,*ai,*aj,ncols,n;
2978:   PetscReal      atmp;
2979:   PetscScalar    *x;
2980:   MatScalar      *aa;

2983:   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2984:   aa = a->a;
2985:   ai = a->i;
2986:   aj = a->j;

2988:   VecSet(v,0.0);
2989:   VecGetArray(v,&x);
2990:   VecGetLocalSize(v,&n);
2991:   if (n != A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector, %D vs. %D rows", A->rmap->n, n);
2992:   for (i=0; i<m; i++) {
2993:     ncols = ai[1] - ai[0]; ai++;
2994:     if (ncols) {
2995:       /* Get first nonzero */
2996:       for (j = 0; j < ncols; j++) {
2997:         atmp = PetscAbsScalar(aa[j]);
2998:         if (atmp > 1.0e-12) {
2999:           x[i] = atmp;
3000:           if (idx) idx[i] = aj[j];
3001:           break;
3002:         }
3003:       }
3004:       if (j == ncols) {x[i] = PetscAbsScalar(*aa); if (idx) idx[i] = *aj;}
3005:     } else {
3006:       x[i] = 0.0; if (idx) idx[i] = 0;
3007:     }
3008:     for (j = 0; j < ncols; j++) {
3009:       atmp = PetscAbsScalar(*aa);
3010:       if (atmp > 1.0e-12 && PetscAbsScalar(x[i]) > atmp) {x[i] = atmp; if (idx) idx[i] = *aj;}
3011:       aa++; aj++;
3012:     }
3013:   }
3014:   VecRestoreArray(v,&x);
3015:   return(0);
3016: }

3020: PetscErrorCode MatGetRowMin_SeqAIJ(Mat A,Vec v,PetscInt idx[])
3021: {
3022:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
3024:   PetscInt       i,j,m = A->rmap->n,*ai,*aj,ncols,n;
3025:   PetscScalar    *x;
3026:   MatScalar      *aa;

3029:   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
3030:   aa = a->a;
3031:   ai = a->i;
3032:   aj = a->j;

3034:   VecSet(v,0.0);
3035:   VecGetArray(v,&x);
3036:   VecGetLocalSize(v,&n);
3037:   if (n != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
3038:   for (i=0; i<m; i++) {
3039:     ncols = ai[1] - ai[0]; ai++;
3040:     if (ncols == A->cmap->n) { /* row is dense */
3041:       x[i] = *aa; if (idx) idx[i] = 0;
3042:     } else {  /* row is sparse so already KNOW minimum is 0.0 or lower */
3043:       x[i] = 0.0;
3044:       if (idx) {   /* find first implicit 0.0 in the row */
3045:         idx[i] = 0; /* in case ncols is zero */
3046:         for (j=0; j<ncols; j++) {
3047:           if (aj[j] > j) {
3048:             idx[i] = j;
3049:             break;
3050:           }
3051:         }
3052:       }
3053:     }
3054:     for (j=0; j<ncols; j++) {
3055:       if (PetscRealPart(x[i]) > PetscRealPart(*aa)) {x[i] = *aa; if (idx) idx[i] = *aj;}
3056:       aa++; aj++;
3057:     }
3058:   }
3059:   VecRestoreArray(v,&x);
3060:   return(0);
3061: }

3063: #include <petscblaslapack.h>
3064: #include <petsc-private/kernels/blockinvert.h>

3068: PetscErrorCode  MatInvertBlockDiagonal_SeqAIJ(Mat A,const PetscScalar **values)
3069: {
3070:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*) A->data;
3072:   PetscInt       i,bs = PetscAbs(A->rmap->bs),mbs = A->rmap->n/bs,ipvt[5],bs2 = bs*bs,*v_pivots,ij[7],*IJ,j;
3073:   MatScalar      *diag,work[25],*v_work;
3074:   PetscReal      shift = 0.0;

3077:   if (a->ibdiagvalid) {
3078:     if (values) *values = a->ibdiag;
3079:     return(0);
3080:   }
3081:   MatMarkDiagonal_SeqAIJ(A);
3082:   if (!a->ibdiag) {
3083:     PetscMalloc1(bs2*mbs,&a->ibdiag);
3084:     PetscLogObjectMemory((PetscObject)A,bs2*mbs*sizeof(PetscScalar));
3085:   }
3086:   diag = a->ibdiag;
3087:   if (values) *values = a->ibdiag;
3088:   /* factor and invert each block */
3089:   switch (bs) {
3090:   case 1:
3091:     for (i=0; i<mbs; i++) {
3092:       MatGetValues(A,1,&i,1,&i,diag+i);
3093:       diag[i] = (PetscScalar)1.0 / (diag[i] + shift);
3094:     }
3095:     break;
3096:   case 2:
3097:     for (i=0; i<mbs; i++) {
3098:       ij[0] = 2*i; ij[1] = 2*i + 1;
3099:       MatGetValues(A,2,ij,2,ij,diag);
3100:       PetscKernel_A_gets_inverse_A_2(diag,shift);
3101:       PetscKernel_A_gets_transpose_A_2(diag);
3102:       diag += 4;
3103:     }
3104:     break;
3105:   case 3:
3106:     for (i=0; i<mbs; i++) {
3107:       ij[0] = 3*i; ij[1] = 3*i + 1; ij[2] = 3*i + 2;
3108:       MatGetValues(A,3,ij,3,ij,diag);
3109:       PetscKernel_A_gets_inverse_A_3(diag,shift);
3110:       PetscKernel_A_gets_transpose_A_3(diag);
3111:       diag += 9;
3112:     }
3113:     break;
3114:   case 4:
3115:     for (i=0; i<mbs; i++) {
3116:       ij[0] = 4*i; ij[1] = 4*i + 1; ij[2] = 4*i + 2; ij[3] = 4*i + 3;
3117:       MatGetValues(A,4,ij,4,ij,diag);
3118:       PetscKernel_A_gets_inverse_A_4(diag,shift);
3119:       PetscKernel_A_gets_transpose_A_4(diag);
3120:       diag += 16;
3121:     }
3122:     break;
3123:   case 5:
3124:     for (i=0; i<mbs; i++) {
3125:       ij[0] = 5*i; ij[1] = 5*i + 1; ij[2] = 5*i + 2; ij[3] = 5*i + 3; ij[4] = 5*i + 4;
3126:       MatGetValues(A,5,ij,5,ij,diag);
3127:       PetscKernel_A_gets_inverse_A_5(diag,ipvt,work,shift);
3128:       PetscKernel_A_gets_transpose_A_5(diag);
3129:       diag += 25;
3130:     }
3131:     break;
3132:   case 6:
3133:     for (i=0; i<mbs; i++) {
3134:       ij[0] = 6*i; ij[1] = 6*i + 1; ij[2] = 6*i + 2; ij[3] = 6*i + 3; ij[4] = 6*i + 4; ij[5] = 6*i + 5;
3135:       MatGetValues(A,6,ij,6,ij,diag);
3136:       PetscKernel_A_gets_inverse_A_6(diag,shift);
3137:       PetscKernel_A_gets_transpose_A_6(diag);
3138:       diag += 36;
3139:     }
3140:     break;
3141:   case 7:
3142:     for (i=0; i<mbs; i++) {
3143:       ij[0] = 7*i; ij[1] = 7*i + 1; ij[2] = 7*i + 2; ij[3] = 7*i + 3; ij[4] = 7*i + 4; ij[5] = 7*i + 5; ij[5] = 7*i + 6;
3144:       MatGetValues(A,7,ij,7,ij,diag);
3145:       PetscKernel_A_gets_inverse_A_7(diag,shift);
3146:       PetscKernel_A_gets_transpose_A_7(diag);
3147:       diag += 49;
3148:     }
3149:     break;
3150:   default:
3151:     PetscMalloc3(bs,&v_work,bs,&v_pivots,bs,&IJ);
3152:     for (i=0; i<mbs; i++) {
3153:       for (j=0; j<bs; j++) {
3154:         IJ[j] = bs*i + j;
3155:       }
3156:       MatGetValues(A,bs,IJ,bs,IJ,diag);
3157:       PetscKernel_A_gets_inverse_A(bs,diag,v_pivots,v_work);
3158:       PetscKernel_A_gets_transpose_A_N(diag,bs);
3159:       diag += bs2;
3160:     }
3161:     PetscFree3(v_work,v_pivots,IJ);
3162:   }
3163:   a->ibdiagvalid = PETSC_TRUE;
3164:   return(0);
3165: }

3169: static PetscErrorCode  MatSetRandom_SeqAIJ(Mat x,PetscRandom rctx)
3170: {
3172:   Mat_SeqAIJ     *aij = (Mat_SeqAIJ*)x->data;
3173:   PetscScalar    a;
3174:   PetscInt       m,n,i,j,col;

3177:   if (!x->assembled) {
3178:     MatGetSize(x,&m,&n);
3179:     for (i=0; i<m; i++) {
3180:       for (j=0; j<aij->imax[i]; j++) {
3181:         PetscRandomGetValue(rctx,&a);
3182:         col  = (PetscInt)(n*PetscRealPart(a));
3183:         MatSetValues(x,1,&i,1,&col,&a,ADD_VALUES);
3184:       }
3185:     }
3186:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet coded");
3187:   MatAssemblyBegin(x,MAT_FINAL_ASSEMBLY);
3188:   MatAssemblyEnd(x,MAT_FINAL_ASSEMBLY);
3189:   return(0);
3190: }

3192: /* -------------------------------------------------------------------*/
3193: static struct _MatOps MatOps_Values = { MatSetValues_SeqAIJ,
3194:                                         MatGetRow_SeqAIJ,
3195:                                         MatRestoreRow_SeqAIJ,
3196:                                         MatMult_SeqAIJ,
3197:                                 /*  4*/ MatMultAdd_SeqAIJ,
3198:                                         MatMultTranspose_SeqAIJ,
3199:                                         MatMultTransposeAdd_SeqAIJ,
3200:                                         0,
3201:                                         0,
3202:                                         0,
3203:                                 /* 10*/ 0,
3204:                                         MatLUFactor_SeqAIJ,
3205:                                         0,
3206:                                         MatSOR_SeqAIJ,
3207:                                         MatTranspose_SeqAIJ,
3208:                                 /*1 5*/ MatGetInfo_SeqAIJ,
3209:                                         MatEqual_SeqAIJ,
3210:                                         MatGetDiagonal_SeqAIJ,
3211:                                         MatDiagonalScale_SeqAIJ,
3212:                                         MatNorm_SeqAIJ,
3213:                                 /* 20*/ 0,
3214:                                         MatAssemblyEnd_SeqAIJ,
3215:                                         MatSetOption_SeqAIJ,
3216:                                         MatZeroEntries_SeqAIJ,
3217:                                 /* 24*/ MatZeroRows_SeqAIJ,
3218:                                         0,
3219:                                         0,
3220:                                         0,
3221:                                         0,
3222:                                 /* 29*/ MatSetUp_SeqAIJ,
3223:                                         0,
3224:                                         0,
3225:                                         0,
3226:                                         0,
3227:                                 /* 34*/ MatDuplicate_SeqAIJ,
3228:                                         0,
3229:                                         0,
3230:                                         MatILUFactor_SeqAIJ,
3231:                                         0,
3232:                                 /* 39*/ MatAXPY_SeqAIJ,
3233:                                         MatGetSubMatrices_SeqAIJ,
3234:                                         MatIncreaseOverlap_SeqAIJ,
3235:                                         MatGetValues_SeqAIJ,
3236:                                         MatCopy_SeqAIJ,
3237:                                 /* 44*/ MatGetRowMax_SeqAIJ,
3238:                                         MatScale_SeqAIJ,
3239:                                         0,
3240:                                         MatDiagonalSet_SeqAIJ,
3241:                                         MatZeroRowsColumns_SeqAIJ,
3242:                                 /* 49*/ MatSetRandom_SeqAIJ,
3243:                                         MatGetRowIJ_SeqAIJ,
3244:                                         MatRestoreRowIJ_SeqAIJ,
3245:                                         MatGetColumnIJ_SeqAIJ,
3246:                                         MatRestoreColumnIJ_SeqAIJ,
3247:                                 /* 54*/ MatFDColoringCreate_SeqXAIJ,
3248:                                         0,
3249:                                         0,
3250:                                         MatPermute_SeqAIJ,
3251:                                         0,
3252:                                 /* 59*/ 0,
3253:                                         MatDestroy_SeqAIJ,
3254:                                         MatView_SeqAIJ,
3255:                                         0,
3256:                                         MatMatMatMult_SeqAIJ_SeqAIJ_SeqAIJ,
3257:                                 /* 64*/ MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ,
3258:                                         MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ,
3259:                                         0,
3260:                                         0,
3261:                                         0,
3262:                                 /* 69*/ MatGetRowMaxAbs_SeqAIJ,
3263:                                         MatGetRowMinAbs_SeqAIJ,
3264:                                         0,
3265:                                         MatSetColoring_SeqAIJ,
3266:                                         0,
3267:                                 /* 74*/ MatSetValuesAdifor_SeqAIJ,
3268:                                         MatFDColoringApply_AIJ,
3269:                                         0,
3270:                                         0,
3271:                                         0,
3272:                                 /* 79*/ MatFindZeroDiagonals_SeqAIJ,
3273:                                         0,
3274:                                         0,
3275:                                         0,
3276:                                         MatLoad_SeqAIJ,
3277:                                 /* 84*/ MatIsSymmetric_SeqAIJ,
3278:                                         MatIsHermitian_SeqAIJ,
3279:                                         0,
3280:                                         0,
3281:                                         0,
3282:                                 /* 89*/ MatMatMult_SeqAIJ_SeqAIJ,
3283:                                         MatMatMultSymbolic_SeqAIJ_SeqAIJ,
3284:                                         MatMatMultNumeric_SeqAIJ_SeqAIJ,
3285:                                         MatPtAP_SeqAIJ_SeqAIJ,
3286:                                         MatPtAPSymbolic_SeqAIJ_SeqAIJ_DenseAxpy,
3287:                                 /* 94*/ MatPtAPNumeric_SeqAIJ_SeqAIJ,
3288:                                         MatMatTransposeMult_SeqAIJ_SeqAIJ,
3289:                                         MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ,
3290:                                         MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ,
3291:                                         0,
3292:                                 /* 99*/ 0,
3293:                                         0,
3294:                                         0,
3295:                                         MatConjugate_SeqAIJ,
3296:                                         0,
3297:                                 /*104*/ MatSetValuesRow_SeqAIJ,
3298:                                         MatRealPart_SeqAIJ,
3299:                                         MatImaginaryPart_SeqAIJ,
3300:                                         0,
3301:                                         0,
3302:                                 /*109*/ MatMatSolve_SeqAIJ,
3303:                                         0,
3304:                                         MatGetRowMin_SeqAIJ,
3305:                                         0,
3306:                                         MatMissingDiagonal_SeqAIJ,
3307:                                 /*114*/ 0,
3308:                                         0,
3309:                                         0,
3310:                                         0,
3311:                                         0,
3312:                                 /*119*/ 0,
3313:                                         0,
3314:                                         0,
3315:                                         0,
3316:                                         MatGetMultiProcBlock_SeqAIJ,
3317:                                 /*124*/ MatFindNonzeroRows_SeqAIJ,
3318:                                         MatGetColumnNorms_SeqAIJ,
3319:                                         MatInvertBlockDiagonal_SeqAIJ,
3320:                                         0,
3321:                                         0,
3322:                                 /*129*/ 0,
3323:                                         MatTransposeMatMult_SeqAIJ_SeqAIJ,
3324:                                         MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ,
3325:                                         MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ,
3326:                                         MatTransposeColoringCreate_SeqAIJ,
3327:                                 /*134*/ MatTransColoringApplySpToDen_SeqAIJ,
3328:                                         MatTransColoringApplyDenToSp_SeqAIJ,
3329:                                         MatRARt_SeqAIJ_SeqAIJ,
3330:                                         MatRARtSymbolic_SeqAIJ_SeqAIJ,
3331:                                         MatRARtNumeric_SeqAIJ_SeqAIJ,
3332:                                  /*139*/0,
3333:                                         0,
3334:                                         0,
3335:                                         MatFDColoringSetUp_SeqXAIJ,
3336:                                         MatFindOffBlockDiagonalEntries_SeqAIJ
3337: };

3341: PetscErrorCode  MatSeqAIJSetColumnIndices_SeqAIJ(Mat mat,PetscInt *indices)
3342: {
3343:   Mat_SeqAIJ *aij = (Mat_SeqAIJ*)mat->data;
3344:   PetscInt   i,nz,n;

3347:   nz = aij->maxnz;
3348:   n  = mat->rmap->n;
3349:   for (i=0; i<nz; i++) {
3350:     aij->j[i] = indices[i];
3351:   }
3352:   aij->nz = nz;
3353:   for (i=0; i<n; i++) {
3354:     aij->ilen[i] = aij->imax[i];
3355:   }
3356:   return(0);
3357: }

3361: /*@
3362:     MatSeqAIJSetColumnIndices - Set the column indices for all the rows
3363:        in the matrix.

3365:   Input Parameters:
3366: +  mat - the SeqAIJ matrix
3367: -  indices - the column indices

3369:   Level: advanced

3371:   Notes:
3372:     This can be called if you have precomputed the nonzero structure of the
3373:   matrix and want to provide it to the matrix object to improve the performance
3374:   of the MatSetValues() operation.

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

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

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

3383: @*/
3384: PetscErrorCode  MatSeqAIJSetColumnIndices(Mat mat,PetscInt *indices)
3385: {

3391:   PetscUseMethod(mat,"MatSeqAIJSetColumnIndices_C",(Mat,PetscInt*),(mat,indices));
3392:   return(0);
3393: }

3395: /* ----------------------------------------------------------------------------------------*/

3399: PetscErrorCode  MatStoreValues_SeqAIJ(Mat mat)
3400: {
3401:   Mat_SeqAIJ     *aij = (Mat_SeqAIJ*)mat->data;
3403:   size_t         nz = aij->i[mat->rmap->n];

3406:   if (!aij->nonew) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first");

3408:   /* allocate space for values if not already there */
3409:   if (!aij->saved_values) {
3410:     PetscMalloc1((nz+1),&aij->saved_values);
3411:     PetscLogObjectMemory((PetscObject)mat,(nz+1)*sizeof(PetscScalar));
3412:   }

3414:   /* copy values over */
3415:   PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));
3416:   return(0);
3417: }

3421: /*@
3422:     MatStoreValues - Stashes a copy of the matrix values; this allows, for
3423:        example, reuse of the linear part of a Jacobian, while recomputing the
3424:        nonlinear portion.

3426:    Collect on Mat

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

3431:   Level: advanced

3433:   Common Usage, with SNESSolve():
3434: $    Create Jacobian matrix
3435: $    Set linear terms into matrix
3436: $    Apply boundary conditions to matrix, at this time matrix must have
3437: $      final nonzero structure (i.e. setting the nonlinear terms and applying
3438: $      boundary conditions again will not change the nonzero structure
3439: $    MatSetOption(mat,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);
3440: $    MatStoreValues(mat);
3441: $    Call SNESSetJacobian() with matrix
3442: $    In your Jacobian routine
3443: $      MatRetrieveValues(mat);
3444: $      Set nonlinear terms in matrix

3446:   Common Usage without SNESSolve(), i.e. when you handle nonlinear solve yourself:
3447: $    // build linear portion of Jacobian
3448: $    MatSetOption(mat,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);
3449: $    MatStoreValues(mat);
3450: $    loop over nonlinear iterations
3451: $       MatRetrieveValues(mat);
3452: $       // call MatSetValues(mat,...) to set nonliner portion of Jacobian
3453: $       // call MatAssemblyBegin/End() on matrix
3454: $       Solve linear system with Jacobian
3455: $    endloop

3457:   Notes:
3458:     Matrix must already be assemblied before calling this routine
3459:     Must set the matrix option MatSetOption(mat,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE); before
3460:     calling this routine.

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

3465: .seealso: MatRetrieveValues()

3467: @*/
3468: PetscErrorCode  MatStoreValues(Mat mat)
3469: {

3474:   if (!mat->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
3475:   if (mat->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
3476:   PetscUseMethod(mat,"MatStoreValues_C",(Mat),(mat));
3477:   return(0);
3478: }

3482: PetscErrorCode  MatRetrieveValues_SeqAIJ(Mat mat)
3483: {
3484:   Mat_SeqAIJ     *aij = (Mat_SeqAIJ*)mat->data;
3486:   PetscInt       nz = aij->i[mat->rmap->n];

3489:   if (!aij->nonew) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first");
3490:   if (!aij->saved_values) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Must call MatStoreValues(A);first");
3491:   /* copy values over */
3492:   PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));
3493:   return(0);
3494: }

3498: /*@
3499:     MatRetrieveValues - Retrieves the copy of the matrix values; this allows, for
3500:        example, reuse of the linear part of a Jacobian, while recomputing the
3501:        nonlinear portion.

3503:    Collect on Mat

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

3508:   Level: advanced

3510: .seealso: MatStoreValues()

3512: @*/
3513: PetscErrorCode  MatRetrieveValues(Mat mat)
3514: {

3519:   if (!mat->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
3520:   if (mat->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
3521:   PetscUseMethod(mat,"MatRetrieveValues_C",(Mat),(mat));
3522:   return(0);
3523: }


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

3536:    Collective on MPI_Comm

3538:    Input Parameters:
3539: +  comm - MPI communicator, set to PETSC_COMM_SELF
3540: .  m - number of rows
3541: .  n - number of columns
3542: .  nz - number of nonzeros per row (same for all rows)
3543: -  nnz - array containing the number of nonzeros in the various rows
3544:          (possibly different for each row) or NULL

3546:    Output Parameter:
3547: .  A - the matrix

3549:    It is recommended that one use the MatCreate(), MatSetType() and/or MatSetFromOptions(),
3550:    MatXXXXSetPreallocation() paradgm instead of this routine directly.
3551:    [MatXXXXSetPreallocation() is, for example, MatSeqAIJSetPreallocation]

3553:    Notes:
3554:    If nnz is given then nz is ignored

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

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

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

3571:    Options Database Keys:
3572: +  -mat_no_inode  - Do not use inodes
3573: -  -mat_inode_limit <limit> - Sets inode limit (max limit=5)

3575:    Level: intermediate

3577: .seealso: MatCreate(), MatCreateAIJ(), MatSetValues(), MatSeqAIJSetColumnIndices(), MatCreateSeqAIJWithArrays()

3579: @*/
3580: PetscErrorCode  MatCreateSeqAIJ(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A)
3581: {

3585:   MatCreate(comm,A);
3586:   MatSetSizes(*A,m,n,m,n);
3587:   MatSetType(*A,MATSEQAIJ);
3588:   MatSeqAIJSetPreallocation_SeqAIJ(*A,nz,nnz);
3589:   return(0);
3590: }

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

3600:    Collective on MPI_Comm

3602:    Input Parameters:
3603: +  B - The matrix
3604: .  nz - number of nonzeros per row (same for all rows)
3605: -  nnz - array containing the number of nonzeros in the various rows
3606:          (possibly different for each row) or NULL

3608:    Notes:
3609:      If nnz is given then nz is ignored

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

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

3621:    You can call MatGetInfo() to get information on how effective the preallocation was;
3622:    for example the fields mallocs,nz_allocated,nz_used,nz_unneeded;
3623:    You can also run with the option -info and look for messages with the string
3624:    malloc in them to see if additional memory allocation was needed.

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

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

3634:    Options Database Keys:
3635: +  -mat_no_inode  - Do not use inodes
3636: .  -mat_inode_limit <limit> - Sets inode limit (max limit=5)
3637: -  -mat_aij_oneindex - Internally use indexing starting at 1
3638:         rather than 0.  Note that when calling MatSetValues(),
3639:         the user still MUST index entries starting at 0!

3641:    Level: intermediate

3643: .seealso: MatCreate(), MatCreateAIJ(), MatSetValues(), MatSeqAIJSetColumnIndices(), MatCreateSeqAIJWithArrays(), MatGetInfo()

3645: @*/
3646: PetscErrorCode  MatSeqAIJSetPreallocation(Mat B,PetscInt nz,const PetscInt nnz[])
3647: {

3653:   PetscTryMethod(B,"MatSeqAIJSetPreallocation_C",(Mat,PetscInt,const PetscInt[]),(B,nz,nnz));
3654:   return(0);
3655: }

3659: PetscErrorCode  MatSeqAIJSetPreallocation_SeqAIJ(Mat B,PetscInt nz,const PetscInt *nnz)
3660: {
3661:   Mat_SeqAIJ     *b;
3662:   PetscBool      skipallocation = PETSC_FALSE,realalloc = PETSC_FALSE;
3664:   PetscInt       i;

3667:   if (nz >= 0 || nnz) realalloc = PETSC_TRUE;
3668:   if (nz == MAT_SKIP_ALLOCATION) {
3669:     skipallocation = PETSC_TRUE;
3670:     nz             = 0;
3671:   }

3673:   PetscLayoutSetUp(B->rmap);
3674:   PetscLayoutSetUp(B->cmap);

3676:   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
3677:   if (nz < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %D",nz);
3678:   if (nnz) {
3679:     for (i=0; i<B->rmap->n; i++) {
3680:       if (nnz[i] < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %D value %D",i,nnz[i]);
3681:       if (nnz[i] > B->cmap->n) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be greater than row length: local row %D value %d rowlength %D",i,nnz[i],B->cmap->n);
3682:     }
3683:   }

3685:   B->preallocated = PETSC_TRUE;

3687:   b = (Mat_SeqAIJ*)B->data;

3689:   if (!skipallocation) {
3690:     if (!b->imax) {
3691:       PetscMalloc2(B->rmap->n,&b->imax,B->rmap->n,&b->ilen);
3692:       PetscLogObjectMemory((PetscObject)B,2*B->rmap->n*sizeof(PetscInt));
3693:     }
3694:     if (!nnz) {
3695:       if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 10;
3696:       else if (nz < 0) nz = 1;
3697:       for (i=0; i<B->rmap->n; i++) b->imax[i] = nz;
3698:       nz = nz*B->rmap->n;
3699:     } else {
3700:       nz = 0;
3701:       for (i=0; i<B->rmap->n; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
3702:     }
3703:     /* b->ilen will count nonzeros in each row so far. */
3704:     for (i=0; i<B->rmap->n; i++) b->ilen[i] = 0;

3706:     /* allocate the matrix space */
3707:     MatSeqXAIJFreeAIJ(B,&b->a,&b->j,&b->i);
3708:     PetscMalloc3(nz,&b->a,nz,&b->j,B->rmap->n+1,&b->i);
3709:     PetscLogObjectMemory((PetscObject)B,(B->rmap->n+1)*sizeof(PetscInt)+nz*(sizeof(PetscScalar)+sizeof(PetscInt)));
3710:     b->i[0] = 0;
3711:     for (i=1; i<B->rmap->n+1; i++) {
3712:       b->i[i] = b->i[i-1] + b->imax[i-1];
3713:     }
3714:     b->singlemalloc = PETSC_TRUE;
3715:     b->free_a       = PETSC_TRUE;
3716:     b->free_ij      = PETSC_TRUE;
3717: #if defined(PETSC_THREADCOMM_ACTIVE)
3718:     MatZeroEntries_SeqAIJ(B);
3719: #endif
3720:   } else {
3721:     b->free_a  = PETSC_FALSE;
3722:     b->free_ij = PETSC_FALSE;
3723:   }

3725:   b->nz               = 0;
3726:   b->maxnz            = nz;
3727:   B->info.nz_unneeded = (double)b->maxnz;
3728:   if (realalloc) {
3729:     MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);
3730:   }
3731:   return(0);
3732: }

3734: #undef  __FUNCT__
3736: /*@
3737:    MatSeqAIJSetPreallocationCSR - Allocates memory for a sparse sequential matrix in AIJ format.

3739:    Input Parameters:
3740: +  B - the matrix
3741: .  i - the indices into j for the start of each row (starts with zero)
3742: .  j - the column indices for each row (starts with zero) these must be sorted for each row
3743: -  v - optional values in the matrix

3745:    Level: developer

3747:    The i,j,v values are COPIED with this routine; to avoid the copy use MatCreateSeqAIJWithArrays()

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

3751: .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatSeqAIJSetPreallocation(), MatCreateSeqAIJ(), SeqAIJ
3752: @*/
3753: PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat B,const PetscInt i[],const PetscInt j[],const PetscScalar v[])
3754: {

3760:   PetscTryMethod(B,"MatSeqAIJSetPreallocationCSR_C",(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]),(B,i,j,v));
3761:   return(0);
3762: }

3764: #undef  __FUNCT__
3766: PetscErrorCode  MatSeqAIJSetPreallocationCSR_SeqAIJ(Mat B,const PetscInt Ii[],const PetscInt J[],const PetscScalar v[])
3767: {
3768:   PetscInt       i;
3769:   PetscInt       m,n;
3770:   PetscInt       nz;
3771:   PetscInt       *nnz, nz_max = 0;
3772:   PetscScalar    *values;

3776:   if (Ii[0]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Ii[0] must be 0 it is %D", Ii[0]);

3778:   PetscLayoutSetUp(B->rmap);
3779:   PetscLayoutSetUp(B->cmap);

3781:   MatGetSize(B, &m, &n);
3782:   PetscMalloc1((m+1), &nnz);
3783:   for (i = 0; i < m; i++) {
3784:     nz     = Ii[i+1]- Ii[i];
3785:     nz_max = PetscMax(nz_max, nz);
3786:     if (nz < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Local row %D has a negative number of columns %D", i, nnz);
3787:     nnz[i] = nz;
3788:   }
3789:   MatSeqAIJSetPreallocation(B, 0, nnz);
3790:   PetscFree(nnz);

3792:   if (v) {
3793:     values = (PetscScalar*) v;
3794:   } else {
3795:     PetscCalloc1(nz_max, &values);
3796:   }

3798:   for (i = 0; i < m; i++) {
3799:     nz   = Ii[i+1] - Ii[i];
3800:     MatSetValues_SeqAIJ(B, 1, &i, nz, J+Ii[i], values + (v ? Ii[i] : 0), INSERT_VALUES);
3801:   }

3803:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
3804:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

3806:   if (!v) {
3807:     PetscFree(values);
3808:   }
3809:   MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
3810:   return(0);
3811: }

3813: #include <../src/mat/impls/dense/seq/dense.h>
3814: #include <petsc-private/kernels/petscaxpy.h>

3818: /*
3819:     Computes (B'*A')' since computing B*A directly is untenable

3821:                n                       p                          p
3822:         (              )       (              )         (                  )
3823:       m (      A       )  *  n (       B      )   =   m (         C        )
3824:         (              )       (              )         (                  )

3826: */
3827: PetscErrorCode MatMatMultNumeric_SeqDense_SeqAIJ(Mat A,Mat B,Mat C)
3828: {
3829:   PetscErrorCode    ierr;
3830:   Mat_SeqDense      *sub_a = (Mat_SeqDense*)A->data;
3831:   Mat_SeqAIJ        *sub_b = (Mat_SeqAIJ*)B->data;
3832:   Mat_SeqDense      *sub_c = (Mat_SeqDense*)C->data;
3833:   PetscInt          i,n,m,q,p;
3834:   const PetscInt    *ii,*idx;
3835:   const PetscScalar *b,*a,*a_q;
3836:   PetscScalar       *c,*c_q;

3839:   m    = A->rmap->n;
3840:   n    = A->cmap->n;
3841:   p    = B->cmap->n;
3842:   a    = sub_a->v;
3843:   b    = sub_b->a;
3844:   c    = sub_c->v;
3845:   PetscMemzero(c,m*p*sizeof(PetscScalar));

3847:   ii  = sub_b->i;
3848:   idx = sub_b->j;
3849:   for (i=0; i<n; i++) {
3850:     q = ii[i+1] - ii[i];
3851:     while (q-->0) {
3852:       c_q = c + m*(*idx);
3853:       a_q = a + m*i;
3854:       PetscKernelAXPY(c_q,*b,a_q,m);
3855:       idx++;
3856:       b++;
3857:     }
3858:   }
3859:   return(0);
3860: }

3864: PetscErrorCode MatMatMultSymbolic_SeqDense_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
3865: {
3867:   PetscInt       m=A->rmap->n,n=B->cmap->n;
3868:   Mat            Cmat;

3871:   if (A->cmap->n != B->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"A->cmap->n %D != B->rmap->n %D\n",A->cmap->n,B->rmap->n);
3872:   MatCreate(PetscObjectComm((PetscObject)A),&Cmat);
3873:   MatSetSizes(Cmat,m,n,m,n);
3874:   MatSetBlockSizesFromMats(Cmat,A,B);
3875:   MatSetType(Cmat,MATSEQDENSE);
3876:   MatSeqDenseSetPreallocation(Cmat,NULL);

3878:   Cmat->ops->matmultnumeric = MatMatMultNumeric_SeqDense_SeqAIJ;

3880:   *C = Cmat;
3881:   return(0);
3882: }

3884: /* ----------------------------------------------------------------*/
3887: PetscErrorCode MatMatMult_SeqDense_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
3888: {

3892:   if (scall == MAT_INITIAL_MATRIX) {
3893:     PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);
3894:     MatMatMultSymbolic_SeqDense_SeqAIJ(A,B,fill,C);
3895:     PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);
3896:   }
3897:   PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);
3898:   MatMatMultNumeric_SeqDense_SeqAIJ(A,B,*C);
3899:   PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);
3900:   return(0);
3901: }


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

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

3911:   Level: beginner

3913: .seealso: MatCreateSeqAIJ(), MatSetFromOptions(), MatSetType(), MatCreate(), MatType
3914: M*/

3916: /*MC
3917:    MATAIJ - MATAIJ = "aij" - A matrix type to be used for sparse matrices.

3919:    This matrix type is identical to MATSEQAIJ when constructed with a single process communicator,
3920:    and MATMPIAIJ otherwise.  As a result, for single process communicators,
3921:   MatSeqAIJSetPreallocation is supported, and similarly MatMPIAIJSetPreallocation is supported
3922:   for communicators controlling multiple processes.  It is recommended that you call both of
3923:   the above preallocation routines for simplicity.

3925:    Options Database Keys:
3926: . -mat_type aij - sets the matrix type to "aij" during a call to MatSetFromOptions()

3928:   Developer Notes: Subclasses include MATAIJCUSP, MATAIJPERM, MATAIJCRL, and also automatically switches over to use inodes when
3929:    enough exist.

3931:   Level: beginner

3933: .seealso: MatCreateAIJ(), MatCreateSeqAIJ(), MATSEQAIJ,MATMPIAIJ
3934: M*/

3936: /*MC
3937:    MATAIJCRL - MATAIJCRL = "aijcrl" - A matrix type to be used for sparse matrices.

3939:    This matrix type is identical to MATSEQAIJCRL when constructed with a single process communicator,
3940:    and MATMPIAIJCRL otherwise.  As a result, for single process communicators,
3941:    MatSeqAIJSetPreallocation() is supported, and similarly MatMPIAIJSetPreallocation() is supported
3942:   for communicators controlling multiple processes.  It is recommended that you call both of
3943:   the above preallocation routines for simplicity.

3945:    Options Database Keys:
3946: . -mat_type aijcrl - sets the matrix type to "aijcrl" during a call to MatSetFromOptions()

3948:   Level: beginner

3950: .seealso: MatCreateMPIAIJCRL,MATSEQAIJCRL,MATMPIAIJCRL, MATSEQAIJCRL, MATMPIAIJCRL
3951: M*/

3953: #if defined(PETSC_HAVE_PASTIX)
3954: PETSC_EXTERN PetscErrorCode MatGetFactor_seqaij_pastix(Mat,MatFactorType,Mat*);
3955: #endif
3956: #if defined(PETSC_HAVE_ESSL) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_REAL_SINGLE) && !defined(PETSC_USE_REAL___FLOAT128)
3957: PETSC_EXTERN PetscErrorCode MatGetFactor_seqaij_essl(Mat,MatFactorType,Mat*);
3958: #endif
3959: PETSC_EXTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJCRL(Mat,MatType,MatReuse,Mat*);
3960: PETSC_EXTERN PetscErrorCode MatGetFactor_seqaij_petsc(Mat,MatFactorType,Mat*);
3961: PETSC_EXTERN PetscErrorCode MatGetFactor_seqaij_bas(Mat,MatFactorType,Mat*);
3962: extern PetscErrorCode  MatGetFactorAvailable_seqaij_petsc(Mat,MatFactorType,PetscBool*);
3963: #if defined(PETSC_HAVE_MUMPS)
3964: PETSC_EXTERN PetscErrorCode MatGetFactor_aij_mumps(Mat,MatFactorType,Mat*);
3965: #endif
3966: #if defined(PETSC_HAVE_SUPERLU)
3967: PETSC_EXTERN PetscErrorCode MatGetFactor_seqaij_superlu(Mat,MatFactorType,Mat*);
3968: #endif
3969: #if defined(PETSC_HAVE_MKL_PARDISO)
3970: PETSC_EXTERN PetscErrorCode MatGetFactor_aij_mkl_pardiso(Mat,MatFactorType,Mat*);
3971: #endif
3972: #if defined(PETSC_HAVE_SUPERLU_DIST)
3973: PETSC_EXTERN PetscErrorCode MatGetFactor_seqaij_superlu_dist(Mat,MatFactorType,Mat*);
3974: #endif
3975: #if defined(PETSC_HAVE_SUITESPARSE)
3976: PETSC_EXTERN PetscErrorCode MatGetFactor_seqaij_umfpack(Mat,MatFactorType,Mat*);
3977: #endif
3978: #if defined(PETSC_HAVE_SUITESPARSE)
3979: PETSC_EXTERN PetscErrorCode MatGetFactor_seqaij_cholmod(Mat,MatFactorType,Mat*);
3980: #endif
3981: #if defined(PETSC_HAVE_SUITESPARSE)
3982: PETSC_EXTERN PetscErrorCode MatGetFactor_seqaij_klu(Mat,MatFactorType,Mat*);
3983: #endif
3984: #if defined(PETSC_HAVE_LUSOL)
3985: PETSC_EXTERN PetscErrorCode MatGetFactor_seqaij_lusol(Mat,MatFactorType,Mat*);
3986: #endif
3987: #if defined(PETSC_HAVE_MATLAB_ENGINE)
3988: PETSC_EXTERN PetscErrorCode MatGetFactor_seqaij_matlab(Mat,MatFactorType,Mat*);
3989: extern PetscErrorCode  MatlabEnginePut_SeqAIJ(PetscObject,void*);
3990: extern PetscErrorCode  MatlabEngineGet_SeqAIJ(PetscObject,void*);
3991: #endif
3992: #if defined(PETSC_HAVE_CLIQUE)
3993: PETSC_EXTERN PetscErrorCode MatGetFactor_aij_clique(Mat,MatFactorType,Mat*);
3994: #endif


3999: /*@C
4000:    MatSeqAIJGetArray - gives access to the array where the data for a SeqSeqAIJ matrix is stored

4002:    Not Collective

4004:    Input Parameter:
4005: .  mat - a MATSEQAIJ matrix

4007:    Output Parameter:
4008: .   array - pointer to the data

4010:    Level: intermediate

4012: .seealso: MatSeqAIJRestoreArray(), MatSeqAIJGetArrayF90()
4013: @*/
4014: PetscErrorCode  MatSeqAIJGetArray(Mat A,PetscScalar **array)
4015: {

4019:   PetscUseMethod(A,"MatSeqAIJGetArray_C",(Mat,PetscScalar**),(A,array));
4020:   return(0);
4021: }

4025: /*@C
4026:    MatSeqAIJGetMaxRowNonzeros - returns the maximum number of nonzeros in any row

4028:    Not Collective

4030:    Input Parameter:
4031: .  mat - a MATSEQAIJ matrix

4033:    Output Parameter:
4034: .   nz - the maximum number of nonzeros in any row

4036:    Level: intermediate

4038: .seealso: MatSeqAIJRestoreArray(), MatSeqAIJGetArrayF90()
4039: @*/
4040: PetscErrorCode  MatSeqAIJGetMaxRowNonzeros(Mat A,PetscInt *nz)
4041: {
4042:   Mat_SeqAIJ     *aij = (Mat_SeqAIJ*)A->data;

4045:   *nz = aij->rmax;
4046:   return(0);
4047: }

4051: /*@C
4052:    MatSeqAIJRestoreArray - returns access to the array where the data for a MATSEQAIJ matrix is stored obtained by MatSeqAIJGetArray()

4054:    Not Collective

4056:    Input Parameters:
4057: .  mat - a MATSEQAIJ matrix
4058: .  array - pointer to the data

4060:    Level: intermediate

4062: .seealso: MatSeqAIJGetArray(), MatSeqAIJRestoreArrayF90()
4063: @*/
4064: PetscErrorCode  MatSeqAIJRestoreArray(Mat A,PetscScalar **array)
4065: {

4069:   PetscUseMethod(A,"MatSeqAIJRestoreArray_C",(Mat,PetscScalar**),(A,array));
4070:   return(0);
4071: }

4075: PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJ(Mat B)
4076: {
4077:   Mat_SeqAIJ     *b;
4079:   PetscMPIInt    size;

4082:   MPI_Comm_size(PetscObjectComm((PetscObject)B),&size);
4083:   if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Comm must be of size 1");

4085:   PetscNewLog(B,&b);

4087:   B->data = (void*)b;

4089:   PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));

4091:   b->row                = 0;
4092:   b->col                = 0;
4093:   b->icol               = 0;
4094:   b->reallocs           = 0;
4095:   b->ignorezeroentries  = PETSC_FALSE;
4096:   b->roworiented        = PETSC_TRUE;
4097:   b->nonew              = 0;
4098:   b->diag               = 0;
4099:   b->solve_work         = 0;
4100:   B->spptr              = 0;
4101:   b->saved_values       = 0;
4102:   b->idiag              = 0;
4103:   b->mdiag              = 0;
4104:   b->ssor_work          = 0;
4105:   b->omega              = 1.0;
4106:   b->fshift             = 0.0;
4107:   b->idiagvalid         = PETSC_FALSE;
4108:   b->ibdiagvalid        = PETSC_FALSE;
4109:   b->keepnonzeropattern = PETSC_FALSE;
4110:   b->xtoy               = 0;
4111:   b->XtoY               = 0;

4113:   PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);
4114:   PetscObjectComposeFunction((PetscObject)B,"MatSeqAIJGetArray_C",MatSeqAIJGetArray_SeqAIJ);
4115:   PetscObjectComposeFunction((PetscObject)B,"MatSeqAIJRestoreArray_C",MatSeqAIJRestoreArray_SeqAIJ);

4117: #if defined(PETSC_HAVE_MATLAB_ENGINE)
4118:   PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_matlab_C",MatGetFactor_seqaij_matlab);
4119:   PetscObjectComposeFunction((PetscObject)B,"PetscMatlabEnginePut_C",MatlabEnginePut_SeqAIJ);
4120:   PetscObjectComposeFunction((PetscObject)B,"PetscMatlabEngineGet_C",MatlabEngineGet_SeqAIJ);
4121: #endif
4122: #if defined(PETSC_HAVE_PASTIX)
4123:   PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_pastix_C",MatGetFactor_seqaij_pastix);
4124: #endif
4125: #if defined(PETSC_HAVE_ESSL) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_REAL_SINGLE) && !defined(PETSC_USE_REAL___FLOAT128)
4126:   PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_essl_C",MatGetFactor_seqaij_essl);
4127: #endif
4128: #if defined(PETSC_HAVE_SUPERLU)
4129:   PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_superlu_C",MatGetFactor_seqaij_superlu);
4130: #endif
4131: #if defined(PETSC_HAVE_MKL_PARDISO)
4132:   PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_mkl_pardiso_C",MatGetFactor_aij_mkl_pardiso);
4133: #endif
4134: #if defined(PETSC_HAVE_SUPERLU_DIST)
4135:   PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_superlu_dist_C",MatGetFactor_seqaij_superlu_dist);
4136: #endif
4137: #if defined(PETSC_HAVE_MUMPS)
4138:   PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_mumps_C",MatGetFactor_aij_mumps);
4139: #endif
4140: #if defined(PETSC_HAVE_SUITESPARSE)
4141:   PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_umfpack_C",MatGetFactor_seqaij_umfpack);
4142: #endif
4143: #if defined(PETSC_HAVE_SUITESPARSE)
4144:   PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_cholmod_C",MatGetFactor_seqaij_cholmod);
4145: #endif
4146: #if defined(PETSC_HAVE_SUITESPARSE)
4147:   PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_klu_C",MatGetFactor_seqaij_klu);
4148: #endif
4149: #if defined(PETSC_HAVE_LUSOL)
4150:   PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_lusol_C",MatGetFactor_seqaij_lusol);
4151: #endif
4152: #if defined(PETSC_HAVE_CLIQUE)
4153:   PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_clique_C",MatGetFactor_aij_clique);
4154: #endif

4156:   PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_petsc_C",MatGetFactor_seqaij_petsc);
4157:   PetscObjectComposeFunction((PetscObject)B,"MatGetFactorAvailable_petsc_C",MatGetFactorAvailable_seqaij_petsc);
4158:   PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_bas_C",MatGetFactor_seqaij_bas);
4159:   PetscObjectComposeFunction((PetscObject)B,"MatSeqAIJSetColumnIndices_C",MatSeqAIJSetColumnIndices_SeqAIJ);
4160:   PetscObjectComposeFunction((PetscObject)B,"MatStoreValues_C",MatStoreValues_SeqAIJ);
4161:   PetscObjectComposeFunction((PetscObject)B,"MatRetrieveValues_C",MatRetrieveValues_SeqAIJ);
4162:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_seqsbaij_C",MatConvert_SeqAIJ_SeqSBAIJ);
4163:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_seqbaij_C",MatConvert_SeqAIJ_SeqBAIJ);
4164:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_seqaijperm_C",MatConvert_SeqAIJ_SeqAIJPERM);
4165:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_seqaijcrl_C",MatConvert_SeqAIJ_SeqAIJCRL);
4166:   PetscObjectComposeFunction((PetscObject)B,"MatIsTranspose_C",MatIsTranspose_SeqAIJ);
4167:   PetscObjectComposeFunction((PetscObject)B,"MatIsHermitianTranspose_C",MatIsTranspose_SeqAIJ);
4168:   PetscObjectComposeFunction((PetscObject)B,"MatSeqAIJSetPreallocation_C",MatSeqAIJSetPreallocation_SeqAIJ);
4169:   PetscObjectComposeFunction((PetscObject)B,"MatSeqAIJSetPreallocationCSR_C",MatSeqAIJSetPreallocationCSR_SeqAIJ);
4170:   PetscObjectComposeFunction((PetscObject)B,"MatReorderForNonzeroDiagonal_C",MatReorderForNonzeroDiagonal_SeqAIJ);
4171:   PetscObjectComposeFunction((PetscObject)B,"MatMatMult_seqdense_seqaij_C",MatMatMult_SeqDense_SeqAIJ);
4172:   PetscObjectComposeFunction((PetscObject)B,"MatMatMultSymbolic_seqdense_seqaij_C",MatMatMultSymbolic_SeqDense_SeqAIJ);
4173:   PetscObjectComposeFunction((PetscObject)B,"MatMatMultNumeric_seqdense_seqaij_C",MatMatMultNumeric_SeqDense_SeqAIJ);
4174:   MatCreate_SeqAIJ_Inode(B);
4175:   PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);
4176:   return(0);
4177: }

4181: /*
4182:     Given a matrix generated with MatGetFactor() duplicates all the information in A into B
4183: */
4184: PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat C,Mat A,MatDuplicateOption cpvalues,PetscBool mallocmatspace)
4185: {
4186:   Mat_SeqAIJ     *c,*a = (Mat_SeqAIJ*)A->data;
4188:   PetscInt       i,m = A->rmap->n;

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

4193:   C->factortype = A->factortype;
4194:   c->row        = 0;
4195:   c->col        = 0;
4196:   c->icol       = 0;
4197:   c->reallocs   = 0;

4199:   C->assembled = PETSC_TRUE;

4201:   PetscLayoutReference(A->rmap,&C->rmap);
4202:   PetscLayoutReference(A->cmap,&C->cmap);

4204:   PetscMalloc2(m,&c->imax,m,&c->ilen);
4205:   PetscLogObjectMemory((PetscObject)C, 2*m*sizeof(PetscInt));
4206:   for (i=0; i<m; i++) {
4207:     c->imax[i] = a->imax[i];
4208:     c->ilen[i] = a->ilen[i];
4209:   }

4211:   /* allocate the matrix space */
4212:   if (mallocmatspace) {
4213:     PetscMalloc3(a->i[m],&c->a,a->i[m],&c->j,m+1,&c->i);
4214:     PetscLogObjectMemory((PetscObject)C, a->i[m]*(sizeof(PetscScalar)+sizeof(PetscInt))+(m+1)*sizeof(PetscInt));

4216:     c->singlemalloc = PETSC_TRUE;

4218:     PetscMemcpy(c->i,a->i,(m+1)*sizeof(PetscInt));
4219:     if (m > 0) {
4220:       PetscMemcpy(c->j,a->j,(a->i[m])*sizeof(PetscInt));
4221:       if (cpvalues == MAT_COPY_VALUES) {
4222:         PetscMemcpy(c->a,a->a,(a->i[m])*sizeof(PetscScalar));
4223:       } else {
4224:         PetscMemzero(c->a,(a->i[m])*sizeof(PetscScalar));
4225:       }
4226:     }
4227:   }

4229:   c->ignorezeroentries = a->ignorezeroentries;
4230:   c->roworiented       = a->roworiented;
4231:   c->nonew             = a->nonew;
4232:   if (a->diag) {
4233:     PetscMalloc1((m+1),&c->diag);
4234:     PetscLogObjectMemory((PetscObject)C,(m+1)*sizeof(PetscInt));
4235:     for (i=0; i<m; i++) {
4236:       c->diag[i] = a->diag[i];
4237:     }
4238:   } else c->diag = 0;

4240:   c->solve_work         = 0;
4241:   c->saved_values       = 0;
4242:   c->idiag              = 0;
4243:   c->ssor_work          = 0;
4244:   c->keepnonzeropattern = a->keepnonzeropattern;
4245:   c->free_a             = PETSC_TRUE;
4246:   c->free_ij            = PETSC_TRUE;
4247:   c->xtoy               = 0;
4248:   c->XtoY               = 0;

4250:   c->rmax         = a->rmax;
4251:   c->nz           = a->nz;
4252:   c->maxnz        = a->nz;       /* Since we allocate exactly the right amount */
4253:   C->preallocated = PETSC_TRUE;

4255:   c->compressedrow.use   = a->compressedrow.use;
4256:   c->compressedrow.nrows = a->compressedrow.nrows;
4257:   if (a->compressedrow.use) {
4258:     i    = a->compressedrow.nrows;
4259:     PetscMalloc2(i+1,&c->compressedrow.i,i,&c->compressedrow.rindex);
4260:     PetscMemcpy(c->compressedrow.i,a->compressedrow.i,(i+1)*sizeof(PetscInt));
4261:     PetscMemcpy(c->compressedrow.rindex,a->compressedrow.rindex,i*sizeof(PetscInt));
4262:   } else {
4263:     c->compressedrow.use    = PETSC_FALSE;
4264:     c->compressedrow.i      = NULL;
4265:     c->compressedrow.rindex = NULL;
4266:   }
4267:   C->nonzerostate = A->nonzerostate;

4269:   MatDuplicate_SeqAIJ_Inode(A,cpvalues,&C);
4270:   PetscFunctionListDuplicate(((PetscObject)A)->qlist,&((PetscObject)C)->qlist);
4271:   return(0);
4272: }

4276: PetscErrorCode MatDuplicate_SeqAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
4277: {

4281:   MatCreate(PetscObjectComm((PetscObject)A),B);
4282:   MatSetSizes(*B,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);
4283:   if (!(A->rmap->n % A->rmap->bs) && !(A->cmap->n % A->cmap->bs)) {
4284:     MatSetBlockSizesFromMats(*B,A,A);
4285:   }
4286:   MatSetType(*B,((PetscObject)A)->type_name);
4287:   MatDuplicateNoCreate_SeqAIJ(*B,A,cpvalues,PETSC_TRUE);
4288:   return(0);
4289: }

4293: PetscErrorCode MatLoad_SeqAIJ(Mat newMat, PetscViewer viewer)
4294: {
4295:   Mat_SeqAIJ     *a;
4297:   PetscInt       i,sum,nz,header[4],*rowlengths = 0,M,N,rows,cols;
4298:   int            fd;
4299:   PetscMPIInt    size;
4300:   MPI_Comm       comm;
4301:   PetscInt       bs = newMat->rmap->bs;

4304:   PetscObjectGetComm((PetscObject)viewer,&comm);
4305:   MPI_Comm_size(comm,&size);
4306:   if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"view must have one processor");

4308:   PetscOptionsBegin(comm,NULL,"Options for loading SEQAIJ matrix","Mat");
4309:   PetscOptionsInt("-matload_block_size","Set the blocksize used to store the matrix","MatLoad",bs,&bs,NULL);
4310:   PetscOptionsEnd();
4311:   if (bs < 0) bs = 1;
4312:   MatSetBlockSize(newMat,bs);

4314:   PetscViewerBinaryGetDescriptor(viewer,&fd);
4315:   PetscBinaryRead(fd,header,4,PETSC_INT);
4316:   if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"not matrix object in file");
4317:   M = header[1]; N = header[2]; nz = header[3];

4319:   if (nz < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format on disk,cannot load as SeqAIJ");

4321:   /* read in row lengths */
4322:   PetscMalloc1(M,&rowlengths);
4323:   PetscBinaryRead(fd,rowlengths,M,PETSC_INT);

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

4329:   /* set global size if not set already*/
4330:   if (newMat->rmap->n < 0 && newMat->rmap->N < 0 && newMat->cmap->n < 0 && newMat->cmap->N < 0) {
4331:     MatSetSizes(newMat,PETSC_DECIDE,PETSC_DECIDE,M,N);
4332:   } else {
4333:     /* if sizes and type are already set, check if the vector global sizes are correct */
4334:     MatGetSize(newMat,&rows,&cols);
4335:     if (rows < 0 && cols < 0) { /* user might provide local size instead of global size */
4336:       MatGetLocalSize(newMat,&rows,&cols);
4337:     }
4338:     if (M != rows ||  N != cols) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Matrix in file of different length (%D, %D) than the input matrix (%D, %D)",M,N,rows,cols);
4339:   }
4340:   MatSeqAIJSetPreallocation_SeqAIJ(newMat,0,rowlengths);
4341:   a    = (Mat_SeqAIJ*)newMat->data;

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

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

4348:   /* set matrix "i" values */
4349:   a->i[0] = 0;
4350:   for (i=1; i<= M; i++) {
4351:     a->i[i]      = a->i[i-1] + rowlengths[i-1];
4352:     a->ilen[i-1] = rowlengths[i-1];
4353:   }
4354:   PetscFree(rowlengths);

4356:   MatAssemblyBegin(newMat,MAT_FINAL_ASSEMBLY);
4357:   MatAssemblyEnd(newMat,MAT_FINAL_ASSEMBLY);
4358:   return(0);
4359: }

4363: PetscErrorCode MatEqual_SeqAIJ(Mat A,Mat B,PetscBool * flg)
4364: {
4365:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data,*b = (Mat_SeqAIJ*)B->data;
4367: #if defined(PETSC_USE_COMPLEX)
4368:   PetscInt k;
4369: #endif

4372:   /* If the  matrix dimensions are not equal,or no of nonzeros */
4373:   if ((A->rmap->n != B->rmap->n) || (A->cmap->n != B->cmap->n) ||(a->nz != b->nz)) {
4374:     *flg = PETSC_FALSE;
4375:     return(0);
4376:   }

4378:   /* if the a->i are the same */
4379:   PetscMemcmp(a->i,b->i,(A->rmap->n+1)*sizeof(PetscInt),flg);
4380:   if (!*flg) return(0);

4382:   /* if a->j are the same */
4383:   PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),flg);
4384:   if (!*flg) return(0);

4386:   /* if a->a are the same */
4387: #if defined(PETSC_USE_COMPLEX)
4388:   for (k=0; k<a->nz; k++) {
4389:     if (PetscRealPart(a->a[k]) != PetscRealPart(b->a[k]) || PetscImaginaryPart(a->a[k]) != PetscImaginaryPart(b->a[k])) {
4390:       *flg = PETSC_FALSE;
4391:       return(0);
4392:     }
4393:   }
4394: #else
4395:   PetscMemcmp(a->a,b->a,(a->nz)*sizeof(PetscScalar),flg);
4396: #endif
4397:   return(0);
4398: }

4402: /*@
4403:      MatCreateSeqAIJWithArrays - Creates an sequential AIJ matrix using matrix elements (in CSR format)
4404:               provided by the user.

4406:       Collective on MPI_Comm

4408:    Input Parameters:
4409: +   comm - must be an MPI communicator of size 1
4410: .   m - number of rows
4411: .   n - number of columns
4412: .   i - row indices
4413: .   j - column indices
4414: -   a - matrix values

4416:    Output Parameter:
4417: .   mat - the matrix

4419:    Level: intermediate

4421:    Notes:
4422:        The i, j, and a arrays are not copied by this routine, the user must free these arrays
4423:     once the matrix is destroyed and not before

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

4427:        The i and j indices are 0 based

4429:        The format which is used for the sparse matrix input, is equivalent to a
4430:     row-major ordering.. i.e for the following matrix, the input data expected is
4431:     as shown:

4433:         1 0 0
4434:         2 0 3
4435:         4 5 6

4437:         i =  {0,1,3,6}  [size = nrow+1  = 3+1]
4438:         j =  {0,0,2,0,1,2}  [size = nz = 6]; values must be sorted for each row
4439:         v =  {1,2,3,4,5,6}  [size = nz = 6]


4442: .seealso: MatCreate(), MatCreateAIJ(), MatCreateSeqAIJ(), MatCreateMPIAIJWithArrays(), MatMPIAIJSetPreallocationCSR()

4444: @*/
4445: PetscErrorCode  MatCreateSeqAIJWithArrays(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt *i,PetscInt *j,PetscScalar *a,Mat *mat)
4446: {
4448:   PetscInt       ii;
4449:   Mat_SeqAIJ     *aij;
4450: #if defined(PETSC_USE_DEBUG)
4451:   PetscInt jj;
4452: #endif

4455:   if (i[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"i (row indices) must start with 0");
4456:   MatCreate(comm,mat);
4457:   MatSetSizes(*mat,m,n,m,n);
4458:   /* MatSetBlockSizes(*mat,,); */
4459:   MatSetType(*mat,MATSEQAIJ);
4460:   MatSeqAIJSetPreallocation_SeqAIJ(*mat,MAT_SKIP_ALLOCATION,0);
4461:   aij  = (Mat_SeqAIJ*)(*mat)->data;
4462:   PetscMalloc2(m,&aij->imax,m,&aij->ilen);

4464:   aij->i            = i;
4465:   aij->j            = j;
4466:   aij->a            = a;
4467:   aij->singlemalloc = PETSC_FALSE;
4468:   aij->nonew        = -1;             /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/
4469:   aij->free_a       = PETSC_FALSE;
4470:   aij->free_ij      = PETSC_FALSE;

4472:   for (ii=0; ii<m; ii++) {
4473:     aij->ilen[ii] = aij->imax[ii] = i[ii+1] - i[ii];
4474: #if defined(PETSC_USE_DEBUG)
4475:     if (i[ii+1] - i[ii] < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row length in i (row indices) row = %D length = %D",ii,i[ii+1] - i[ii]);
4476:     for (jj=i[ii]+1; jj<i[ii+1]; jj++) {
4477:       if (j[jj] < j[jj-1]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column entry number %D (actual colum %D) in row %D is not sorted",jj-i[ii],j[jj],ii);
4478:       if (j[jj] == j[jj]-1) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column entry number %D (actual colum %D) in row %D is identical to previous entry",jj-i[ii],j[jj],ii);
4479:     }
4480: #endif
4481:   }
4482: #if defined(PETSC_USE_DEBUG)
4483:   for (ii=0; ii<aij->i[m]; ii++) {
4484:     if (j[ii] < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column index at location = %D index = %D",ii,j[ii]);
4485:     if (j[ii] > n - 1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column index to large at location = %D index = %D",ii,j[ii]);
4486:   }
4487: #endif

4489:   MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
4490:   MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
4491:   return(0);
4492: }
4495: /*@C
4496:      MatCreateSeqAIJFromTriple - Creates an sequential AIJ matrix using matrix elements (in COO format)
4497:               provided by the user.

4499:       Collective on MPI_Comm

4501:    Input Parameters:
4502: +   comm - must be an MPI communicator of size 1
4503: .   m   - number of rows
4504: .   n   - number of columns
4505: .   i   - row indices
4506: .   j   - column indices
4507: .   a   - matrix values
4508: .   nz  - number of nonzeros
4509: -   idx - 0 or 1 based

4511:    Output Parameter:
4512: .   mat - the matrix

4514:    Level: intermediate

4516:    Notes:
4517:        The i and j indices are 0 based

4519:        The format which is used for the sparse matrix input, is equivalent to a
4520:     row-major ordering.. i.e for the following matrix, the input data expected is
4521:     as shown:

4523:         1 0 0
4524:         2 0 3
4525:         4 5 6

4527:         i =  {0,1,1,2,2,2}
4528:         j =  {0,0,2,0,1,2}
4529:         v =  {1,2,3,4,5,6}


4532: .seealso: MatCreate(), MatCreateAIJ(), MatCreateSeqAIJ(), MatCreateSeqAIJWithArrays(), MatMPIAIJSetPreallocationCSR()

4534: @*/
4535: PetscErrorCode  MatCreateSeqAIJFromTriple(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt *i,PetscInt *j,PetscScalar *a,Mat *mat,PetscInt nz,PetscBool idx)
4536: {
4538:   PetscInt       ii, *nnz, one = 1,row,col;


4542:   PetscCalloc1(m,&nnz);
4543:   for (ii = 0; ii < nz; ii++) {
4544:     nnz[i[ii] - !!idx] += 1;
4545:   }
4546:   MatCreate(comm,mat);
4547:   MatSetSizes(*mat,m,n,m,n);
4548:   MatSetType(*mat,MATSEQAIJ);
4549:   MatSeqAIJSetPreallocation_SeqAIJ(*mat,0,nnz);
4550:   for (ii = 0; ii < nz; ii++) {
4551:     if (idx) {
4552:       row = i[ii] - 1;
4553:       col = j[ii] - 1;
4554:     } else {
4555:       row = i[ii];
4556:       col = j[ii];
4557:     }
4558:     MatSetValues(*mat,one,&row,one,&col,&a[ii],ADD_VALUES);
4559:   }
4560:   MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
4561:   MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
4562:   PetscFree(nnz);
4563:   return(0);
4564: }

4568: PetscErrorCode MatSetColoring_SeqAIJ(Mat A,ISColoring coloring)
4569: {
4571:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;

4574:   if (coloring->ctype == IS_COLORING_GLOBAL) {
4575:     ISColoringReference(coloring);
4576:     a->coloring = coloring;
4577:   } else if (coloring->ctype == IS_COLORING_GHOSTED) {
4578:     PetscInt        i,*larray;
4579:     ISColoring      ocoloring;
4580:     ISColoringValue *colors;

4582:     /* set coloring for diagonal portion */
4583:     PetscMalloc1(A->cmap->n,&larray);
4584:     for (i=0; i<A->cmap->n; i++) larray[i] = i;
4585:     ISGlobalToLocalMappingApply(A->cmap->mapping,IS_GTOLM_MASK,A->cmap->n,larray,NULL,larray);
4586:     PetscMalloc1(A->cmap->n,&colors);
4587:     for (i=0; i<A->cmap->n; i++) colors[i] = coloring->colors[larray[i]];
4588:     PetscFree(larray);
4589:     ISColoringCreate(PETSC_COMM_SELF,coloring->n,A->cmap->n,colors,&ocoloring);
4590:     a->coloring = ocoloring;
4591:   }
4592:   return(0);
4593: }

4597: PetscErrorCode MatSetValuesAdifor_SeqAIJ(Mat A,PetscInt nl,void *advalues)
4598: {
4599:   Mat_SeqAIJ      *a      = (Mat_SeqAIJ*)A->data;
4600:   PetscInt        m       = A->rmap->n,*ii = a->i,*jj = a->j,nz,i,j;
4601:   MatScalar       *v      = a->a;
4602:   PetscScalar     *values = (PetscScalar*)advalues;
4603:   ISColoringValue *color;

4606:   if (!a->coloring) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Coloring not set for matrix");
4607:   color = a->coloring->colors;
4608:   /* loop over rows */
4609:   for (i=0; i<m; i++) {
4610:     nz = ii[i+1] - ii[i];
4611:     /* loop over columns putting computed value into matrix */
4612:     for (j=0; j<nz; j++) *v++ = values[color[*jj++]];
4613:     values += nl; /* jump to next row of derivatives */
4614:   }
4615:   return(0);
4616: }

4620: PetscErrorCode MatSeqAIJInvalidateDiagonal(Mat A)
4621: {
4622:   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;

4626:   a->idiagvalid  = PETSC_FALSE;
4627:   a->ibdiagvalid = PETSC_FALSE;

4629:   MatSeqAIJInvalidateDiagonal_Inode(A);
4630:   return(0);
4631: }

4633: /*
4634:     Special version for direct calls from Fortran
4635: */
4636: #include <petsc-private/fortranimpl.h>
4637: #if defined(PETSC_HAVE_FORTRAN_CAPS)
4638: #define matsetvaluesseqaij_ MATSETVALUESSEQAIJ
4639: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
4640: #define matsetvaluesseqaij_ matsetvaluesseqaij
4641: #endif

4643: /* Change these macros so can be used in void function */
4644: #undef CHKERRQ
4645: #define CHKERRQ(ierr) CHKERRABORT(PetscObjectComm((PetscObject)A),ierr)
4646: #undef SETERRQ2
4647: #define SETERRQ2(comm,ierr,b,c,d) CHKERRABORT(comm,ierr)
4648: #undef SETERRQ3
4649: #define SETERRQ3(comm,ierr,b,c,d,e) CHKERRABORT(comm,ierr)

4653: PETSC_EXTERN void PETSC_STDCALL matsetvaluesseqaij_(Mat *AA,PetscInt *mm,const PetscInt im[],PetscInt *nn,const PetscInt in[],const PetscScalar v[],InsertMode *isis, PetscErrorCode *_ierr)
4654: {
4655:   Mat            A  = *AA;
4656:   PetscInt       m  = *mm, n = *nn;
4657:   InsertMode     is = *isis;
4658:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
4659:   PetscInt       *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N;
4660:   PetscInt       *imax,*ai,*ailen;
4662:   PetscInt       *aj,nonew = a->nonew,lastcol = -1;
4663:   MatScalar      *ap,value,*aa;
4664:   PetscBool      ignorezeroentries = a->ignorezeroentries;
4665:   PetscBool      roworiented       = a->roworiented;

4668:   MatCheckPreallocated(A,1);
4669:   imax  = a->imax;
4670:   ai    = a->i;
4671:   ailen = a->ilen;
4672:   aj    = a->j;
4673:   aa    = a->a;

4675:   for (k=0; k<m; k++) { /* loop over added rows */
4676:     row = im[k];
4677:     if (row < 0) continue;
4678: #if defined(PETSC_USE_DEBUG)
4679:     if (row >= A->rmap->n) SETERRABORT(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
4680: #endif
4681:     rp   = aj + ai[row]; ap = aa + ai[row];
4682:     rmax = imax[row]; nrow = ailen[row];
4683:     low  = 0;
4684:     high = nrow;
4685:     for (l=0; l<n; l++) { /* loop over added columns */
4686:       if (in[l] < 0) continue;
4687: #if defined(PETSC_USE_DEBUG)
4688:       if (in[l] >= A->cmap->n) SETERRABORT(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
4689: #endif
4690:       col = in[l];
4691:       if (roworiented) value = v[l + k*n];
4692:       else value = v[k + l*m];

4694:       if (value == 0.0 && ignorezeroentries && (is == ADD_VALUES)) continue;

4696:       if (col <= lastcol) low = 0;
4697:       else high = nrow;
4698:       lastcol = col;
4699:       while (high-low > 5) {
4700:         t = (low+high)/2;
4701:         if (rp[t] > col) high = t;
4702:         else             low  = t;
4703:       }
4704:       for (i=low; i<high; i++) {
4705:         if (rp[i] > col) break;
4706:         if (rp[i] == col) {
4707:           if (is == ADD_VALUES) ap[i] += value;
4708:           else                  ap[i] = value;
4709:           goto noinsert;
4710:         }
4711:       }
4712:       if (value == 0.0 && ignorezeroentries) goto noinsert;
4713:       if (nonew == 1) goto noinsert;
4714:       if (nonew == -1) SETERRABORT(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
4715:       MatSeqXAIJReallocateAIJ(A,A->rmap->n,1,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar);
4716:       N = nrow++ - 1; a->nz++; high++;
4717:       /* shift up all the later entries in this row */
4718:       for (ii=N; ii>=i; ii--) {
4719:         rp[ii+1] = rp[ii];
4720:         ap[ii+1] = ap[ii];
4721:       }
4722:       rp[i] = col;
4723:       ap[i] = value;
4724:       A->nonzerostate++;
4725: noinsert:;
4726:       low = i + 1;
4727:     }
4728:     ailen[row] = nrow;
4729:   }
4730:   PetscFunctionReturnVoid();
4731: }