Actual source code: aij.c

petsc-master 2019-05-20
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>
  9:  #include <petscblaslapack.h>
 10:  #include <petscbt.h>
 11:  #include <petsc/private/kernels/blocktranspose.h>

 13: PetscErrorCode MatSeqAIJSetTypeFromOptions(Mat A)
 14: {
 15:   PetscErrorCode       ierr;
 16:   PetscBool            flg;
 17:   char                 type[256];

 20:   PetscObjectOptionsBegin((PetscObject)A);
 21:   PetscOptionsFList("-mat_seqaij_type","Matrix SeqAIJ type","MatSeqAIJSetType",MatSeqAIJList,"seqaij",type,256,&flg);
 22:   if (flg) {
 23:     MatSeqAIJSetType(A,type);
 24:   }
 25:   PetscOptionsEnd();
 26:   return(0);
 27: }

 29: PetscErrorCode MatGetColumnNorms_SeqAIJ(Mat A,NormType type,PetscReal *norms)
 30: {
 32:   PetscInt       i,m,n;
 33:   Mat_SeqAIJ     *aij = (Mat_SeqAIJ*)A->data;

 36:   MatGetSize(A,&m,&n);
 37:   PetscMemzero(norms,n*sizeof(PetscReal));
 38:   if (type == NORM_2) {
 39:     for (i=0; i<aij->i[m]; i++) {
 40:       norms[aij->j[i]] += PetscAbsScalar(aij->a[i]*aij->a[i]);
 41:     }
 42:   } else if (type == NORM_1) {
 43:     for (i=0; i<aij->i[m]; i++) {
 44:       norms[aij->j[i]] += PetscAbsScalar(aij->a[i]);
 45:     }
 46:   } else if (type == NORM_INFINITY) {
 47:     for (i=0; i<aij->i[m]; i++) {
 48:       norms[aij->j[i]] = PetscMax(PetscAbsScalar(aij->a[i]),norms[aij->j[i]]);
 49:     }
 50:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown NormType");

 52:   if (type == NORM_2) {
 53:     for (i=0; i<n; i++) norms[i] = PetscSqrtReal(norms[i]);
 54:   }
 55:   return(0);
 56: }

 58: PetscErrorCode MatFindOffBlockDiagonalEntries_SeqAIJ(Mat A,IS *is)
 59: {
 60:   Mat_SeqAIJ      *a  = (Mat_SeqAIJ*)A->data;
 61:   PetscInt        i,m=A->rmap->n,cnt = 0, bs = A->rmap->bs;
 62:   const PetscInt  *jj = a->j,*ii = a->i;
 63:   PetscInt        *rows;
 64:   PetscErrorCode  ierr;

 67:   for (i=0; i<m; i++) {
 68:     if ((ii[i] != ii[i+1]) && ((jj[ii[i]] < bs*(i/bs)) || (jj[ii[i+1]-1] > bs*((i+bs)/bs)-1))) {
 69:       cnt++;
 70:     }
 71:   }
 72:   PetscMalloc1(cnt,&rows);
 73:   cnt  = 0;
 74:   for (i=0; i<m; i++) {
 75:     if ((ii[i] != ii[i+1]) && ((jj[ii[i]] < bs*(i/bs)) || (jj[ii[i+1]-1] > bs*((i+bs)/bs)-1))) {
 76:       rows[cnt] = i;
 77:       cnt++;
 78:     }
 79:   }
 80:   ISCreateGeneral(PETSC_COMM_SELF,cnt,rows,PETSC_OWN_POINTER,is);
 81:   return(0);
 82: }

 84: PetscErrorCode MatFindZeroDiagonals_SeqAIJ_Private(Mat A,PetscInt *nrows,PetscInt **zrows)
 85: {
 86:   Mat_SeqAIJ      *a  = (Mat_SeqAIJ*)A->data;
 87:   const MatScalar *aa = a->a;
 88:   PetscInt        i,m=A->rmap->n,cnt = 0;
 89:   const PetscInt  *ii = a->i,*jj = a->j,*diag;
 90:   PetscInt        *rows;
 91:   PetscErrorCode  ierr;

 94:   MatMarkDiagonal_SeqAIJ(A);
 95:   diag = a->diag;
 96:   for (i=0; i<m; i++) {
 97:     if ((diag[i] >= ii[i+1]) || (jj[diag[i]] != i) || (aa[diag[i]] == 0.0)) {
 98:       cnt++;
 99:     }
100:   }
101:   PetscMalloc1(cnt,&rows);
102:   cnt  = 0;
103:   for (i=0; i<m; i++) {
104:     if ((diag[i] >= ii[i+1]) || (jj[diag[i]] != i) || (aa[diag[i]] == 0.0)) {
105:       rows[cnt++] = i;
106:     }
107:   }
108:   *nrows = cnt;
109:   *zrows = rows;
110:   return(0);
111: }

113: PetscErrorCode MatFindZeroDiagonals_SeqAIJ(Mat A,IS *zrows)
114: {
115:   PetscInt       nrows,*rows;

119:   *zrows = NULL;
120:   MatFindZeroDiagonals_SeqAIJ_Private(A,&nrows,&rows);
121:   ISCreateGeneral(PetscObjectComm((PetscObject)A),nrows,rows,PETSC_OWN_POINTER,zrows);
122:   return(0);
123: }

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

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

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

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

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,PETSC_TRUE,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: }

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

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

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

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

278:         cja[cia[col] + collengths[col]++ - oshift] = row + oshift;
279:       }
280:     }
281:     PetscFree(collengths);
282:     *ia  = cia; *ja = cja;
283:   }
284:   return(0);
285: }

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

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

294:   PetscFree(*ia);
295:   PetscFree(*ja);
296:   return(0);
297: }

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

313:   *nn = n;
314:   if (!ia) return(0);

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

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

349:   MatRestoreColumnIJ_SeqAIJ(A,oshift,symmetric,inodecompressed,n,ia,ja,done);
350:   PetscFree(*spidx);
351:   return(0);
352: }

354: PetscErrorCode MatSetValuesRow_SeqAIJ(Mat A,PetscInt row,const PetscScalar v[])
355: {
356:   Mat_SeqAIJ     *a  = (Mat_SeqAIJ*)A->data;
357:   PetscInt       *ai = a->i;

361:   PetscMemcpy(a->a+ai[row],v,(ai[row+1]-ai[row])*sizeof(PetscScalar));
362:   return(0);
363: }

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

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

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

375: */

377:  #include <petsc/private/isimpl.h>
378: PetscErrorCode MatSeqAIJSetValuesLocalFast(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
379: {
380:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
381:   PetscInt       low,high,t,row,nrow,i,col,l;
382:   const PetscInt *rp,*ai = a->i,*ailen = a->ilen,*aj = a->j;
383:   PetscInt       lastcol = -1;
384:   MatScalar      *ap,value,*aa = a->a;
385:   const PetscInt *ridx = A->rmap->mapping->indices,*cidx = A->cmap->mapping->indices;

387:   row = ridx[im[0]];
388:   rp   = aj + ai[row];
389:   ap = aa + ai[row];
390:   nrow = ailen[row];
391:   low  = 0;
392:   high = nrow;
393:   for (l=0; l<n; l++) { /* loop over added columns */
394:     col = cidx[in[l]];
395:     value = v[l];

397:     if (col <= lastcol) low = 0;
398:     else high = nrow;
399:     lastcol = col;
400:     while (high-low > 5) {
401:       t = (low+high)/2;
402:       if (rp[t] > col) high = t;
403:       else low = t;
404:     }
405:     for (i=low; i<high; i++) {
406:       if (rp[i] == col) {
407:         ap[i] += value;
408:         low = i + 1;
409:         break;
410:       }
411:     }
412:   }
413:   return 0;
414: }

416: PetscErrorCode MatSetValues_SeqAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
417: {
418:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
419:   PetscInt       *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N;
420:   PetscInt       *imax = a->imax,*ai = a->i,*ailen = a->ilen;
422:   PetscInt       *aj = a->j,nonew = a->nonew,lastcol = -1;
423:   MatScalar      *ap=NULL,value=0.0,*aa = a->a;
424:   PetscBool      ignorezeroentries = a->ignorezeroentries;
425:   PetscBool      roworiented       = a->roworiented;

428:   for (k=0; k<m; k++) { /* loop over added rows */
429:     row = im[k];
430:     if (row < 0) continue;
431: #if defined(PETSC_USE_DEBUG)
432:     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);
433: #endif
434:     rp   = aj + ai[row];
435:     if (!A->structure_only) ap = aa + ai[row];
436:     rmax = imax[row]; nrow = ailen[row];
437:     low  = 0;
438:     high = nrow;
439:     for (l=0; l<n; l++) { /* loop over added columns */
440:       if (in[l] < 0) continue;
441: #if defined(PETSC_USE_DEBUG)
442:       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);
443: #endif
444:       col = in[l];
445:       if (!A->structure_only) {
446:         if (roworiented) {
447:           value = v[l + k*n];
448:         } else {
449:           value = v[k + l*m];
450:         }
451:       } else { /* A->structure_only */
452:         value = 1; /* avoid 'continue' below?  */
453:       }
454:       if ((value == 0.0 && ignorezeroentries) && (is == ADD_VALUES) && row != col) continue;

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


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

512:   for (k=0; k<m; k++) { /* loop over rows */
513:     row = im[k];
514:     if (row < 0) {v += n; continue;} /* SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",row); */
515:     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);
516:     rp   = aj + ai[row]; ap = aa + ai[row];
517:     nrow = ailen[row];
518:     for (l=0; l<n; l++) { /* loop over columns */
519:       if (in[l] < 0) {v++; continue;} /* SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %D",in[l]); */
520:       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);
521:       col  = in[l];
522:       high = nrow; low = 0; /* assume unsorted */
523:       while (high-low > 5) {
524:         t = (low+high)/2;
525:         if (rp[t] > col) high = t;
526:         else low = t;
527:       }
528:       for (i=low; i<high; i++) {
529:         if (rp[i] > col) break;
530:         if (rp[i] == col) {
531:           *v++ = ap[i];
532:           goto finished;
533:         }
534:       }
535:       *v++ = 0.0;
536: finished:;
537:     }
538:   }
539:   return(0);
540: }


543: PetscErrorCode MatView_SeqAIJ_Binary(Mat A,PetscViewer viewer)
544: {
545:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
547:   PetscInt       i,*col_lens;
548:   int            fd;
549:   FILE           *file;

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

555:   col_lens[0] = MAT_FILE_CLASSID;
556:   col_lens[1] = A->rmap->n;
557:   col_lens[2] = A->cmap->n;
558:   col_lens[3] = a->nz;

560:   /* store lengths of each row and write (including header) to file */
561:   for (i=0; i<A->rmap->n; i++) {
562:     col_lens[4+i] = a->i[i+1] - a->i[i];
563:   }
564:   PetscBinaryWrite(fd,col_lens,4+A->rmap->n,PETSC_INT,PETSC_TRUE);
565:   PetscFree(col_lens);

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

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

573:   PetscViewerBinaryGetInfoPointer(viewer,&file);
574:   if (file) {
575:     fprintf(file,"-matload_block_size %d\n",(int)PetscAbs(A->rmap->bs));
576:   }
577:   return(0);
578: }

580: static PetscErrorCode MatView_SeqAIJ_ASCII_structonly(Mat A,PetscViewer viewer)
581: {
583:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
584:   PetscInt       i,k,m=A->rmap->N;

587:   PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
588:   for (i=0; i<m; i++) {
589:     PetscViewerASCIIPrintf(viewer,"row %D:",i);
590:     for (k=a->i[i]; k<a->i[i+1]; k++) {
591:       PetscViewerASCIIPrintf(viewer," (%D) ",a->j[k]);
592:     }
593:     PetscViewerASCIIPrintf(viewer,"\n");
594:   }
595:   PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
596:   return(0);
597: }

599: extern PetscErrorCode MatSeqAIJFactorInfo_Matlab(Mat,PetscViewer);

601: PetscErrorCode MatView_SeqAIJ_ASCII(Mat A,PetscViewer viewer)
602: {
603:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
604:   PetscErrorCode    ierr;
605:   PetscInt          i,j,m = A->rmap->n;
606:   const char        *name;
607:   PetscViewerFormat format;

610:   if (A->structure_only) {
611:     MatView_SeqAIJ_ASCII_structonly(A,viewer);
612:     return(0);
613:   }

615:   PetscViewerGetFormat(viewer,&format);
616:   if (format == PETSC_VIEWER_ASCII_MATLAB) {
617:     PetscInt nofinalvalue = 0;
618:     if (m && ((a->i[m] == a->i[m-1]) || (a->j[a->nz-1] != A->cmap->n-1))) {
619:       /* Need a dummy value to ensure the dimension of the matrix. */
620:       nofinalvalue = 1;
621:     }
622:     PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
623:     PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",m,A->cmap->n);
624:     PetscViewerASCIIPrintf(viewer,"%% Nonzeros = %D \n",a->nz);
625: #if defined(PETSC_USE_COMPLEX)
626:     PetscViewerASCIIPrintf(viewer,"zzz = zeros(%D,4);\n",a->nz+nofinalvalue);
627: #else
628:     PetscViewerASCIIPrintf(viewer,"zzz = zeros(%D,3);\n",a->nz+nofinalvalue);
629: #endif
630:     PetscViewerASCIIPrintf(viewer,"zzz = [\n");

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

736:     for (i=0; i<a->i[m]; i++) {
737:       if (PetscImaginaryPart(a->a[i]) != 0.0) {
738:         realonly = PETSC_FALSE;
739:         break;
740:       }
741:     }
742: #endif

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

819:         /* U part */
820:         for (j=a->diag[i+1]+1; j<a->diag[i]; j++) {
821: #if defined(PETSC_USE_COMPLEX)
822:           if (PetscImaginaryPart(a->a[j]) > 0.0) {
823:             PetscViewerASCIIPrintf(viewer," (%D, %g + %g i)",a->j[j],(double)PetscRealPart(a->a[j]),(double)PetscImaginaryPart(a->a[j]));
824:           } else if (PetscImaginaryPart(a->a[j]) < 0.0) {
825:             PetscViewerASCIIPrintf(viewer," (%D, %g - %g i)",a->j[j],(double)PetscRealPart(a->a[j]),(double)(-PetscImaginaryPart(a->a[j])));
826:           } else {
827:             PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[j],(double)PetscRealPart(a->a[j]));
828:           }
829: #else
830:           PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[j],(double)a->a[j]);
831: #endif
832:         }
833:         PetscViewerASCIIPrintf(viewer,"\n");
834:       }
835:     } else {
836:       for (i=0; i<m; i++) {
837:         PetscViewerASCIIPrintf(viewer,"row %D:",i);
838:         for (j=a->i[i]; j<a->i[i+1]; j++) {
839: #if defined(PETSC_USE_COMPLEX)
840:           if (PetscImaginaryPart(a->a[j]) > 0.0) {
841:             PetscViewerASCIIPrintf(viewer," (%D, %g + %g i)",a->j[j],(double)PetscRealPart(a->a[j]),(double)PetscImaginaryPart(a->a[j]));
842:           } else if (PetscImaginaryPart(a->a[j]) < 0.0) {
843:             PetscViewerASCIIPrintf(viewer," (%D, %g - %g i)",a->j[j],(double)PetscRealPart(a->a[j]),(double)-PetscImaginaryPart(a->a[j]));
844:           } else {
845:             PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[j],(double)PetscRealPart(a->a[j]));
846:           }
847: #else
848:           PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[j],(double)a->a[j]);
849: #endif
850:         }
851:         PetscViewerASCIIPrintf(viewer,"\n");
852:       }
853:     }
854:     PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
855:   }
856:   PetscViewerFlush(viewer);
857:   return(0);
858: }

860:  #include <petscdraw.h>
861: PetscErrorCode MatView_SeqAIJ_Draw_Zoom(PetscDraw draw,void *Aa)
862: {
863:   Mat               A  = (Mat) Aa;
864:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
865:   PetscErrorCode    ierr;
866:   PetscInt          i,j,m = A->rmap->n;
867:   int               color;
868:   PetscReal         xl,yl,xr,yr,x_l,x_r,y_l,y_r;
869:   PetscViewer       viewer;
870:   PetscViewerFormat format;

873:   PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);
874:   PetscViewerGetFormat(viewer,&format);
875:   PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);

877:   /* loop over matrix elements drawing boxes */

879:   if (format != PETSC_VIEWER_DRAW_CONTOUR) {
880:     PetscDrawCollectiveBegin(draw);
881:     /* Blue for negative, Cyan for zero and  Red for positive */
882:     color = PETSC_DRAW_BLUE;
883:     for (i=0; i<m; i++) {
884:       y_l = m - i - 1.0; y_r = y_l + 1.0;
885:       for (j=a->i[i]; j<a->i[i+1]; j++) {
886:         x_l = a->j[j]; x_r = x_l + 1.0;
887:         if (PetscRealPart(a->a[j]) >=  0.) continue;
888:         PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
889:       }
890:     }
891:     color = PETSC_DRAW_CYAN;
892:     for (i=0; i<m; i++) {
893:       y_l = m - i - 1.0; y_r = y_l + 1.0;
894:       for (j=a->i[i]; j<a->i[i+1]; j++) {
895:         x_l = a->j[j]; x_r = x_l + 1.0;
896:         if (a->a[j] !=  0.) continue;
897:         PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
898:       }
899:     }
900:     color = PETSC_DRAW_RED;
901:     for (i=0; i<m; i++) {
902:       y_l = m - i - 1.0; y_r = y_l + 1.0;
903:       for (j=a->i[i]; j<a->i[i+1]; j++) {
904:         x_l = a->j[j]; x_r = x_l + 1.0;
905:         if (PetscRealPart(a->a[j]) <=  0.) continue;
906:         PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
907:       }
908:     }
909:     PetscDrawCollectiveEnd(draw);
910:   } else {
911:     /* use contour shading to indicate magnitude of values */
912:     /* first determine max of all nonzero values */
913:     PetscReal minv = 0.0, maxv = 0.0;
914:     PetscInt  nz = a->nz, count = 0;
915:     PetscDraw popup;

917:     for (i=0; i<nz; i++) {
918:       if (PetscAbsScalar(a->a[i]) > maxv) maxv = PetscAbsScalar(a->a[i]);
919:     }
920:     if (minv >= maxv) maxv = minv + PETSC_SMALL;
921:     PetscDrawGetPopup(draw,&popup);
922:     PetscDrawScalePopup(popup,minv,maxv);

924:     PetscDrawCollectiveBegin(draw);
925:     for (i=0; i<m; i++) {
926:       y_l = m - i - 1.0;
927:       y_r = y_l + 1.0;
928:       for (j=a->i[i]; j<a->i[i+1]; j++) {
929:         x_l = a->j[j];
930:         x_r = x_l + 1.0;
931:         color = PetscDrawRealToColor(PetscAbsScalar(a->a[count]),minv,maxv);
932:         PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
933:         count++;
934:       }
935:     }
936:     PetscDrawCollectiveEnd(draw);
937:   }
938:   return(0);
939: }

941:  #include <petscdraw.h>
942: PetscErrorCode MatView_SeqAIJ_Draw(Mat A,PetscViewer viewer)
943: {
945:   PetscDraw      draw;
946:   PetscReal      xr,yr,xl,yl,h,w;
947:   PetscBool      isnull;

950:   PetscViewerDrawGetDraw(viewer,0,&draw);
951:   PetscDrawIsNull(draw,&isnull);
952:   if (isnull) return(0);

954:   xr   = A->cmap->n; yr  = A->rmap->n; h = yr/10.0; w = xr/10.0;
955:   xr  += w;          yr += h;         xl = -w;     yl = -h;
956:   PetscDrawSetCoordinates(draw,xl,yl,xr,yr);
957:   PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);
958:   PetscDrawZoom(draw,MatView_SeqAIJ_Draw_Zoom,A);
959:   PetscObjectCompose((PetscObject)A,"Zoomviewer",NULL);
960:   PetscDrawSave(draw);
961:   return(0);
962: }

964: PetscErrorCode MatView_SeqAIJ(Mat A,PetscViewer viewer)
965: {
967:   PetscBool      iascii,isbinary,isdraw;

970:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
971:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
972:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
973:   if (iascii) {
974:     MatView_SeqAIJ_ASCII(A,viewer);
975:   } else if (isbinary) {
976:     MatView_SeqAIJ_Binary(A,viewer);
977:   } else if (isdraw) {
978:     MatView_SeqAIJ_Draw(A,viewer);
979:   }
980:   MatView_SeqAIJ_Inode(A,viewer);
981:   return(0);
982: }

984: PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode)
985: {
986:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
988:   PetscInt       fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
989:   PetscInt       m      = A->rmap->n,*ip,N,*ailen = a->ilen,rmax = 0;
990:   MatScalar      *aa    = a->a,*ap;
991:   PetscReal      ratio  = 0.6;

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

996:   if (m) rmax = ailen[0]; /* determine row with most nonzeros */
997:   for (i=1; i<m; i++) {
998:     /* move each row back by the amount of empty slots (fshift) before it*/
999:     fshift += imax[i-1] - ailen[i-1];
1000:     rmax    = PetscMax(rmax,ailen[i]);
1001:     if (fshift) {
1002:       ip = aj + ai[i];
1003:       ap = aa + ai[i];
1004:       N  = ailen[i];
1005:       for (j=0; j<N; j++) {
1006:         ip[j-fshift] = ip[j];
1007:         if (!A->structure_only) ap[j-fshift] = ap[j];
1008:       }
1009:     }
1010:     ai[i] = ai[i-1] + ailen[i-1];
1011:   }
1012:   if (m) {
1013:     fshift += imax[m-1] - ailen[m-1];
1014:     ai[m]   = ai[m-1] + ailen[m-1];
1015:   }

1017:   /* reset ilen and imax for each row */
1018:   a->nonzerorowcnt = 0;
1019:   if (A->structure_only) {
1020:     PetscFree2(a->imax,a->ilen);
1021:   } else { /* !A->structure_only */
1022:     for (i=0; i<m; i++) {
1023:       ailen[i] = imax[i] = ai[i+1] - ai[i];
1024:       a->nonzerorowcnt += ((ai[i+1] - ai[i]) > 0);
1025:     }
1026:   }
1027:   a->nz = ai[m];
1028:   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);

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

1035:   A->info.mallocs    += a->reallocs;
1036:   a->reallocs         = 0;
1037:   A->info.nz_unneeded = (PetscReal)fshift;
1038:   a->rmax             = rmax;

1040:   if (!A->structure_only) {
1041:     MatCheckCompressedRow(A,a->nonzerorowcnt,&a->compressedrow,a->i,m,ratio);
1042:   }
1043:   MatAssemblyEnd_SeqAIJ_Inode(A,mode);
1044:   MatSeqAIJInvalidateDiagonal(A);
1045:   return(0);
1046: }

1048: PetscErrorCode MatRealPart_SeqAIJ(Mat A)
1049: {
1050:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1051:   PetscInt       i,nz = a->nz;
1052:   MatScalar      *aa = a->a;

1056:   for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]);
1057:   MatSeqAIJInvalidateDiagonal(A);
1058:   return(0);
1059: }

1061: PetscErrorCode MatImaginaryPart_SeqAIJ(Mat A)
1062: {
1063:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1064:   PetscInt       i,nz = a->nz;
1065:   MatScalar      *aa = a->a;

1069:   for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]);
1070:   MatSeqAIJInvalidateDiagonal(A);
1071:   return(0);
1072: }

1074: PetscErrorCode MatZeroEntries_SeqAIJ(Mat A)
1075: {
1076:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;

1080:   PetscMemzero(a->a,(a->i[A->rmap->n])*sizeof(PetscScalar));
1081:   MatSeqAIJInvalidateDiagonal(A);
1082:   return(0);
1083: }

1085: PetscErrorCode MatDestroy_SeqAIJ(Mat A)
1086: {
1087:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;

1091: #if defined(PETSC_USE_LOG)
1092:   PetscLogObjectState((PetscObject)A,"Rows=%D, Cols=%D, NZ=%D",A->rmap->n,A->cmap->n,a->nz);
1093: #endif
1094:   MatSeqXAIJFreeAIJ(A,&a->a,&a->j,&a->i);
1095:   ISDestroy(&a->row);
1096:   ISDestroy(&a->col);
1097:   PetscFree(a->diag);
1098:   PetscFree(a->ibdiag);
1099:   PetscFree2(a->imax,a->ilen);
1100:   PetscFree(a->ipre);
1101:   PetscFree3(a->idiag,a->mdiag,a->ssor_work);
1102:   PetscFree(a->solve_work);
1103:   ISDestroy(&a->icol);
1104:   PetscFree(a->saved_values);
1105:   ISColoringDestroy(&a->coloring);
1106:   PetscFree2(a->compressedrow.i,a->compressedrow.rindex);
1107:   PetscFree(a->matmult_abdense);

1109:   MatDestroy_SeqAIJ_Inode(A);
1110:   PetscFree(A->data);

1112:   PetscObjectChangeTypeName((PetscObject)A,0);
1113:   PetscObjectComposeFunction((PetscObject)A,"MatSeqAIJSetColumnIndices_C",NULL);
1114:   PetscObjectComposeFunction((PetscObject)A,"MatStoreValues_C",NULL);
1115:   PetscObjectComposeFunction((PetscObject)A,"MatRetrieveValues_C",NULL);
1116:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqaij_seqsbaij_C",NULL);
1117:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqaij_seqbaij_C",NULL);
1118:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqaij_seqaijperm_C",NULL);
1119: #if defined(PETSC_HAVE_ELEMENTAL)
1120:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqaij_elemental_C",NULL);
1121: #endif
1122: #if defined(PETSC_HAVE_HYPRE)
1123:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqaij_hypre_C",NULL);
1124:   PetscObjectComposeFunction((PetscObject)A,"MatMatMatMult_transpose_seqaij_seqaij_C",NULL);
1125: #endif
1126:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqaij_seqdense_C",NULL);
1127:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqaij_seqsell_C",NULL);
1128:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqaij_is_C",NULL);
1129:   PetscObjectComposeFunction((PetscObject)A,"MatIsTranspose_C",NULL);
1130:   PetscObjectComposeFunction((PetscObject)A,"MatSeqAIJSetPreallocation_C",NULL);
1131:   PetscObjectComposeFunction((PetscObject)A,"MatResetPreallocation_C",NULL);
1132:   PetscObjectComposeFunction((PetscObject)A,"MatSeqAIJSetPreallocationCSR_C",NULL);
1133:   PetscObjectComposeFunction((PetscObject)A,"MatReorderForNonzeroDiagonal_C",NULL);
1134:   PetscObjectComposeFunction((PetscObject)A,"MatPtAP_is_seqaij_C",NULL);
1135:   return(0);
1136: }

1138: PetscErrorCode MatSetOption_SeqAIJ(Mat A,MatOption op,PetscBool flg)
1139: {
1140:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;

1144:   switch (op) {
1145:   case MAT_ROW_ORIENTED:
1146:     a->roworiented = flg;
1147:     break;
1148:   case MAT_KEEP_NONZERO_PATTERN:
1149:     a->keepnonzeropattern = flg;
1150:     break;
1151:   case MAT_NEW_NONZERO_LOCATIONS:
1152:     a->nonew = (flg ? 0 : 1);
1153:     break;
1154:   case MAT_NEW_NONZERO_LOCATION_ERR:
1155:     a->nonew = (flg ? -1 : 0);
1156:     break;
1157:   case MAT_NEW_NONZERO_ALLOCATION_ERR:
1158:     a->nonew = (flg ? -2 : 0);
1159:     break;
1160:   case MAT_UNUSED_NONZERO_LOCATION_ERR:
1161:     a->nounused = (flg ? -1 : 0);
1162:     break;
1163:   case MAT_IGNORE_ZERO_ENTRIES:
1164:     a->ignorezeroentries = flg;
1165:     break;
1166:   case MAT_SPD:
1167:   case MAT_SYMMETRIC:
1168:   case MAT_STRUCTURALLY_SYMMETRIC:
1169:   case MAT_HERMITIAN:
1170:   case MAT_SYMMETRY_ETERNAL:
1171:   case MAT_STRUCTURE_ONLY:
1172:     /* These options are handled directly by MatSetOption() */
1173:     break;
1174:   case MAT_NEW_DIAGONALS:
1175:   case MAT_IGNORE_OFF_PROC_ENTRIES:
1176:   case MAT_USE_HASH_TABLE:
1177:     PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
1178:     break;
1179:   case MAT_USE_INODES:
1180:     /* Not an error because MatSetOption_SeqAIJ_Inode handles this one */
1181:     break;
1182:   case MAT_SUBMAT_SINGLEIS:
1183:     A->submat_singleis = flg;
1184:     break;
1185:   default:
1186:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %d",op);
1187:   }
1188:   MatSetOption_SeqAIJ_Inode(A,op,flg);
1189:   return(0);
1190: }

1192: PetscErrorCode MatGetDiagonal_SeqAIJ(Mat A,Vec v)
1193: {
1194:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1196:   PetscInt       i,j,n,*ai=a->i,*aj=a->j,nz;
1197:   PetscScalar    *aa=a->a,*x,zero=0.0;

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

1203:   if (A->factortype == MAT_FACTOR_ILU || A->factortype == MAT_FACTOR_LU) {
1204:     PetscInt *diag=a->diag;
1205:     VecGetArray(v,&x);
1206:     for (i=0; i<n; i++) x[i] = 1.0/aa[diag[i]];
1207:     VecRestoreArray(v,&x);
1208:     return(0);
1209:   }

1211:   VecSet(v,zero);
1212:   VecGetArray(v,&x);
1213:   for (i=0; i<n; i++) {
1214:     nz = ai[i+1] - ai[i];
1215:     if (!nz) x[i] = 0.0;
1216:     for (j=ai[i]; j<ai[i+1]; j++) {
1217:       if (aj[j] == i) {
1218:         x[i] = aa[j];
1219:         break;
1220:       }
1221:     }
1222:   }
1223:   VecRestoreArray(v,&x);
1224:   return(0);
1225: }

1227: #include <../src/mat/impls/aij/seq/ftn-kernels/fmult.h>
1228: PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
1229: {
1230:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
1231:   PetscScalar       *y;
1232:   const PetscScalar *x;
1233:   PetscErrorCode    ierr;
1234:   PetscInt          m = A->rmap->n;
1235: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
1236:   const MatScalar   *v;
1237:   PetscScalar       alpha;
1238:   PetscInt          n,i,j;
1239:   const PetscInt    *idx,*ii,*ridx=NULL;
1240:   Mat_CompressedRow cprow    = a->compressedrow;
1241:   PetscBool         usecprow = cprow.use;
1242: #endif

1245:   if (zz != yy) {VecCopy(zz,yy);}
1246:   VecGetArrayRead(xx,&x);
1247:   VecGetArray(yy,&y);

1249: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
1250:   fortranmulttransposeaddaij_(&m,x,a->i,a->j,a->a,y);
1251: #else
1252:   if (usecprow) {
1253:     m    = cprow.nrows;
1254:     ii   = cprow.i;
1255:     ridx = cprow.rindex;
1256:   } else {
1257:     ii = a->i;
1258:   }
1259:   for (i=0; i<m; i++) {
1260:     idx = a->j + ii[i];
1261:     v   = a->a + ii[i];
1262:     n   = ii[i+1] - ii[i];
1263:     if (usecprow) {
1264:       alpha = x[ridx[i]];
1265:     } else {
1266:       alpha = x[i];
1267:     }
1268:     for (j=0; j<n; j++) y[idx[j]] += alpha*v[j];
1269:   }
1270: #endif
1271:   PetscLogFlops(2.0*a->nz);
1272:   VecRestoreArrayRead(xx,&x);
1273:   VecRestoreArray(yy,&y);
1274:   return(0);
1275: }

1277: PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec xx,Vec yy)
1278: {

1282:   VecSet(yy,0.0);
1283:   MatMultTransposeAdd_SeqAIJ(A,xx,yy,yy);
1284:   return(0);
1285: }

1287: #include <../src/mat/impls/aij/seq/ftn-kernels/fmult.h>

1289: PetscErrorCode MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
1290: {
1291:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
1292:   PetscScalar       *y;
1293:   const PetscScalar *x;
1294:   const MatScalar   *aa;
1295:   PetscErrorCode    ierr;
1296:   PetscInt          m=A->rmap->n;
1297:   const PetscInt    *aj,*ii,*ridx=NULL;
1298:   PetscInt          n,i;
1299:   PetscScalar       sum;
1300:   PetscBool         usecprow=a->compressedrow.use;

1302: #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
1303: #pragma disjoint(*x,*y,*aa)
1304: #endif

1307:   VecGetArrayRead(xx,&x);
1308:   VecGetArray(yy,&y);
1309:   ii   = a->i;
1310:   if (usecprow) { /* use compressed row format */
1311:     PetscMemzero(y,m*sizeof(PetscScalar));
1312:     m    = a->compressedrow.nrows;
1313:     ii   = a->compressedrow.i;
1314:     ridx = a->compressedrow.rindex;
1315:     for (i=0; i<m; i++) {
1316:       n           = ii[i+1] - ii[i];
1317:       aj          = a->j + ii[i];
1318:       aa          = a->a + ii[i];
1319:       sum         = 0.0;
1320:       PetscSparseDensePlusDot(sum,x,aa,aj,n);
1321:       /* for (j=0; j<n; j++) sum += (*aa++)*x[*aj++]; */
1322:       y[*ridx++] = sum;
1323:     }
1324:   } else { /* do not use compressed row format */
1325: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
1326:     aj   = a->j;
1327:     aa   = a->a;
1328:     fortranmultaij_(&m,x,ii,aj,aa,y);
1329: #else
1330:     for (i=0; i<m; i++) {
1331:       n           = ii[i+1] - ii[i];
1332:       aj          = a->j + ii[i];
1333:       aa          = a->a + ii[i];
1334:       sum         = 0.0;
1335:       PetscSparseDensePlusDot(sum,x,aa,aj,n);
1336:       y[i] = sum;
1337:     }
1338: #endif
1339:   }
1340:   PetscLogFlops(2.0*a->nz - a->nonzerorowcnt);
1341:   VecRestoreArrayRead(xx,&x);
1342:   VecRestoreArray(yy,&y);
1343:   return(0);
1344: }

1346: PetscErrorCode MatMultMax_SeqAIJ(Mat A,Vec xx,Vec yy)
1347: {
1348:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
1349:   PetscScalar       *y;
1350:   const PetscScalar *x;
1351:   const MatScalar   *aa;
1352:   PetscErrorCode    ierr;
1353:   PetscInt          m=A->rmap->n;
1354:   const PetscInt    *aj,*ii,*ridx=NULL;
1355:   PetscInt          n,i,nonzerorow=0;
1356:   PetscScalar       sum;
1357:   PetscBool         usecprow=a->compressedrow.use;

1359: #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
1360: #pragma disjoint(*x,*y,*aa)
1361: #endif

1364:   VecGetArrayRead(xx,&x);
1365:   VecGetArray(yy,&y);
1366:   if (usecprow) { /* use compressed row format */
1367:     m    = a->compressedrow.nrows;
1368:     ii   = a->compressedrow.i;
1369:     ridx = a->compressedrow.rindex;
1370:     for (i=0; i<m; i++) {
1371:       n           = ii[i+1] - ii[i];
1372:       aj          = a->j + ii[i];
1373:       aa          = a->a + ii[i];
1374:       sum         = 0.0;
1375:       nonzerorow += (n>0);
1376:       PetscSparseDenseMaxDot(sum,x,aa,aj,n);
1377:       /* for (j=0; j<n; j++) sum += (*aa++)*x[*aj++]; */
1378:       y[*ridx++] = sum;
1379:     }
1380:   } else { /* do not use compressed row format */
1381:     ii = a->i;
1382:     for (i=0; i<m; i++) {
1383:       n           = ii[i+1] - ii[i];
1384:       aj          = a->j + ii[i];
1385:       aa          = a->a + ii[i];
1386:       sum         = 0.0;
1387:       nonzerorow += (n>0);
1388:       PetscSparseDenseMaxDot(sum,x,aa,aj,n);
1389:       y[i] = sum;
1390:     }
1391:   }
1392:   PetscLogFlops(2.0*a->nz - nonzerorow);
1393:   VecRestoreArrayRead(xx,&x);
1394:   VecRestoreArray(yy,&y);
1395:   return(0);
1396: }

1398: PetscErrorCode MatMultAddMax_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
1399: {
1400:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
1401:   PetscScalar       *y,*z;
1402:   const PetscScalar *x;
1403:   const MatScalar   *aa;
1404:   PetscErrorCode    ierr;
1405:   PetscInt          m = A->rmap->n,*aj,*ii;
1406:   PetscInt          n,i,*ridx=NULL;
1407:   PetscScalar       sum;
1408:   PetscBool         usecprow=a->compressedrow.use;

1411:   VecGetArrayRead(xx,&x);
1412:   VecGetArrayPair(yy,zz,&y,&z);
1413:   if (usecprow) { /* use compressed row format */
1414:     if (zz != yy) {
1415:       PetscMemcpy(z,y,m*sizeof(PetscScalar));
1416:     }
1417:     m    = a->compressedrow.nrows;
1418:     ii   = a->compressedrow.i;
1419:     ridx = a->compressedrow.rindex;
1420:     for (i=0; i<m; i++) {
1421:       n   = ii[i+1] - ii[i];
1422:       aj  = a->j + ii[i];
1423:       aa  = a->a + ii[i];
1424:       sum = y[*ridx];
1425:       PetscSparseDenseMaxDot(sum,x,aa,aj,n);
1426:       z[*ridx++] = sum;
1427:     }
1428:   } else { /* do not use compressed row format */
1429:     ii = a->i;
1430:     for (i=0; i<m; i++) {
1431:       n   = ii[i+1] - ii[i];
1432:       aj  = a->j + ii[i];
1433:       aa  = a->a + ii[i];
1434:       sum = y[i];
1435:       PetscSparseDenseMaxDot(sum,x,aa,aj,n);
1436:       z[i] = sum;
1437:     }
1438:   }
1439:   PetscLogFlops(2.0*a->nz);
1440:   VecRestoreArrayRead(xx,&x);
1441:   VecRestoreArrayPair(yy,zz,&y,&z);
1442:   return(0);
1443: }

1445: #include <../src/mat/impls/aij/seq/ftn-kernels/fmultadd.h>
1446: PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
1447: {
1448:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
1449:   PetscScalar       *y,*z;
1450:   const PetscScalar *x;
1451:   const MatScalar   *aa;
1452:   PetscErrorCode    ierr;
1453:   const PetscInt    *aj,*ii,*ridx=NULL;
1454:   PetscInt          m = A->rmap->n,n,i;
1455:   PetscScalar       sum;
1456:   PetscBool         usecprow=a->compressedrow.use;

1459:   VecGetArrayRead(xx,&x);
1460:   VecGetArrayPair(yy,zz,&y,&z);
1461:   if (usecprow) { /* use compressed row format */
1462:     if (zz != yy) {
1463:       PetscMemcpy(z,y,m*sizeof(PetscScalar));
1464:     }
1465:     m    = a->compressedrow.nrows;
1466:     ii   = a->compressedrow.i;
1467:     ridx = a->compressedrow.rindex;
1468:     for (i=0; i<m; i++) {
1469:       n   = ii[i+1] - ii[i];
1470:       aj  = a->j + ii[i];
1471:       aa  = a->a + ii[i];
1472:       sum = y[*ridx];
1473:       PetscSparseDensePlusDot(sum,x,aa,aj,n);
1474:       z[*ridx++] = sum;
1475:     }
1476:   } else { /* do not use compressed row format */
1477:     ii = a->i;
1478: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ)
1479:     aj = a->j;
1480:     aa = a->a;
1481:     fortranmultaddaij_(&m,x,ii,aj,aa,y,z);
1482: #else
1483:     for (i=0; i<m; i++) {
1484:       n   = ii[i+1] - ii[i];
1485:       aj  = a->j + ii[i];
1486:       aa  = a->a + ii[i];
1487:       sum = y[i];
1488:       PetscSparseDensePlusDot(sum,x,aa,aj,n);
1489:       z[i] = sum;
1490:     }
1491: #endif
1492:   }
1493:   PetscLogFlops(2.0*a->nz);
1494:   VecRestoreArrayRead(xx,&x);
1495:   VecRestoreArrayPair(yy,zz,&y,&z);
1496:   return(0);
1497: }

1499: /*
1500:      Adds diagonal pointers to sparse matrix structure.
1501: */
1502: PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat A)
1503: {
1504:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1506:   PetscInt       i,j,m = A->rmap->n;

1509:   if (!a->diag) {
1510:     PetscMalloc1(m,&a->diag);
1511:     PetscLogObjectMemory((PetscObject)A, m*sizeof(PetscInt));
1512:   }
1513:   for (i=0; i<A->rmap->n; i++) {
1514:     a->diag[i] = a->i[i+1];
1515:     for (j=a->i[i]; j<a->i[i+1]; j++) {
1516:       if (a->j[j] == i) {
1517:         a->diag[i] = j;
1518:         break;
1519:       }
1520:     }
1521:   }
1522:   return(0);
1523: }

1525: PetscErrorCode MatShift_SeqAIJ(Mat A,PetscScalar v)
1526: {
1527:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
1528:   const PetscInt    *diag = (const PetscInt*)a->diag;
1529:   const PetscInt    *ii = (const PetscInt*) a->i;
1530:   PetscInt          i,*mdiag = NULL;
1531:   PetscErrorCode    ierr;
1532:   PetscInt          cnt = 0; /* how many diagonals are missing */

1535:   if (!A->preallocated || !a->nz) {
1536:     MatSeqAIJSetPreallocation(A,1,NULL);
1537:     MatShift_Basic(A,v);
1538:     return(0);
1539:   }

1541:   if (a->diagonaldense) {
1542:     cnt = 0;
1543:   } else {
1544:     PetscCalloc1(A->rmap->n,&mdiag);
1545:     for (i=0; i<A->rmap->n; i++) {
1546:       if (diag[i] >= ii[i+1]) {
1547:         cnt++;
1548:         mdiag[i] = 1;
1549:       }
1550:     }
1551:   }
1552:   if (!cnt) {
1553:     MatShift_Basic(A,v);
1554:   } else {
1555:     PetscScalar *olda = a->a;  /* preserve pointers to current matrix nonzeros structure and values */
1556:     PetscInt    *oldj = a->j, *oldi = a->i;
1557:     PetscBool   singlemalloc = a->singlemalloc,free_a = a->free_a,free_ij = a->free_ij;

1559:     a->a = NULL;
1560:     a->j = NULL;
1561:     a->i = NULL;
1562:     /* increase the values in imax for each row where a diagonal is being inserted then reallocate the matrix data structures */
1563:     for (i=0; i<A->rmap->n; i++) {
1564:       a->imax[i] += mdiag[i];
1565:       a->imax[i] = PetscMin(a->imax[i],A->cmap->n);
1566:     }
1567:     MatSeqAIJSetPreallocation_SeqAIJ(A,0,a->imax);

1569:     /* copy old values into new matrix data structure */
1570:     for (i=0; i<A->rmap->n; i++) {
1571:       MatSetValues(A,1,&i,a->imax[i] - mdiag[i],&oldj[oldi[i]],&olda[oldi[i]],ADD_VALUES);
1572:       if (i < A->cmap->n) {
1573:         MatSetValue(A,i,i,v,ADD_VALUES);
1574:       }
1575:     }
1576:     MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
1577:     MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
1578:     if (singlemalloc) {
1579:       PetscFree3(olda,oldj,oldi);
1580:     } else {
1581:       if (free_a)  {PetscFree(olda);}
1582:       if (free_ij) {PetscFree(oldj);}
1583:       if (free_ij) {PetscFree(oldi);}
1584:     }
1585:   }
1586:   PetscFree(mdiag);
1587:   a->diagonaldense = PETSC_TRUE;
1588:   return(0);
1589: }

1591: /*
1592:      Checks for missing diagonals
1593: */
1594: PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat A,PetscBool  *missing,PetscInt *d)
1595: {
1596:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1597:   PetscInt       *diag,*ii = a->i,i;

1601:   *missing = PETSC_FALSE;
1602:   if (A->rmap->n > 0 && !ii) {
1603:     *missing = PETSC_TRUE;
1604:     if (d) *d = 0;
1605:     PetscInfo(A,"Matrix has no entries therefore is missing diagonal\n");
1606:   } else {
1607:     diag = a->diag;
1608:     for (i=0; i<A->rmap->n; i++) {
1609:       if (diag[i] >= ii[i+1]) {
1610:         *missing = PETSC_TRUE;
1611:         if (d) *d = i;
1612:         PetscInfo1(A,"Matrix is missing diagonal number %D\n",i);
1613:         break;
1614:       }
1615:     }
1616:   }
1617:   return(0);
1618: }

1620:  #include <petscblaslapack.h>
1621:  #include <petsc/private/kernels/blockinvert.h>

1623: /*
1624:     Note that values is allocated externally by the PC and then passed into this routine
1625: */
1626: PetscErrorCode MatInvertVariableBlockDiagonal_SeqAIJ(Mat A,PetscInt nblocks,const PetscInt *bsizes,PetscScalar *diag)
1627: {
1628:   PetscErrorCode  ierr;
1629:   PetscInt        n = A->rmap->n, i, ncnt = 0, *indx,j,bsizemax = 0,*v_pivots;
1630:   PetscBool       allowzeropivot,zeropivotdetected=PETSC_FALSE;
1631:   const PetscReal shift = 0.0;
1632:   PetscInt        ipvt[5];
1633:   PetscScalar     work[25],*v_work;

1636:   allowzeropivot = PetscNot(A->erroriffailure);
1637:   for (i=0; i<nblocks; i++) ncnt += bsizes[i];
1638:   if (ncnt != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Total blocksizes %D doesn't match number matrix rows %D",ncnt,n);
1639:   for (i=0; i<nblocks; i++) {
1640:     bsizemax = PetscMax(bsizemax,bsizes[i]);
1641:   }
1642:   PetscMalloc1(bsizemax,&indx);
1643:   if (bsizemax > 7) {
1644:     PetscMalloc2(bsizemax,&v_work,bsizemax,&v_pivots);
1645:   }
1646:   ncnt = 0;
1647:   for (i=0; i<nblocks; i++) {
1648:     for (j=0; j<bsizes[i]; j++) indx[j] = ncnt+j;
1649:     MatGetValues(A,bsizes[i],indx,bsizes[i],indx,diag);
1650:     switch (bsizes[i]) {
1651:     case 1:
1652:       *diag = 1.0/(*diag);
1653:       break;
1654:     case 2:
1655:       PetscKernel_A_gets_inverse_A_2(diag,shift,allowzeropivot,&zeropivotdetected);
1656:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1657:       PetscKernel_A_gets_transpose_A_2(diag);
1658:       break;
1659:     case 3:
1660:       PetscKernel_A_gets_inverse_A_3(diag,shift,allowzeropivot,&zeropivotdetected);
1661:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1662:       PetscKernel_A_gets_transpose_A_3(diag);
1663:       break;
1664:     case 4:
1665:       PetscKernel_A_gets_inverse_A_4(diag,shift,allowzeropivot,&zeropivotdetected);
1666:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1667:       PetscKernel_A_gets_transpose_A_4(diag);
1668:       break;
1669:     case 5:
1670:       PetscKernel_A_gets_inverse_A_5(diag,ipvt,work,shift,allowzeropivot,&zeropivotdetected);
1671:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1672:       PetscKernel_A_gets_transpose_A_5(diag);
1673:       break;
1674:     case 6:
1675:       PetscKernel_A_gets_inverse_A_6(diag,shift,allowzeropivot,&zeropivotdetected);
1676:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1677:       PetscKernel_A_gets_transpose_A_6(diag);
1678:       break;
1679:     case 7:
1680:       PetscKernel_A_gets_inverse_A_7(diag,shift,allowzeropivot,&zeropivotdetected);
1681:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1682:       PetscKernel_A_gets_transpose_A_7(diag);
1683:       break;
1684:     default:
1685:       PetscKernel_A_gets_inverse_A(bsizes[i],diag,v_pivots,v_work,allowzeropivot,&zeropivotdetected);
1686:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1687:       PetscKernel_A_gets_transpose_A_N(diag,bsizes[i]);
1688:     }
1689:     ncnt   += bsizes[i];
1690:     diag += bsizes[i]*bsizes[i];
1691:   }
1692:   if (bsizemax > 7) {
1693:     PetscFree2(v_work,v_pivots);
1694:   }
1695:   PetscFree(indx);
1696:   return(0);
1697: }

1699: /*
1700:    Negative shift indicates do not generate an error if there is a zero diagonal, just invert it anyways
1701: */
1702: PetscErrorCode  MatInvertDiagonal_SeqAIJ(Mat A,PetscScalar omega,PetscScalar fshift)
1703: {
1704:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*) A->data;
1706:   PetscInt       i,*diag,m = A->rmap->n;
1707:   MatScalar      *v = a->a;
1708:   PetscScalar    *idiag,*mdiag;

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

1722:   if (omega == 1.0 && PetscRealPart(fshift) <= 0.0) {
1723:     for (i=0; i<m; i++) {
1724:       mdiag[i] = v[diag[i]];
1725:       if (!PetscAbsScalar(mdiag[i])) { /* zero diagonal */
1726:         if (PetscRealPart(fshift)) {
1727:           PetscInfo1(A,"Zero diagonal on row %D\n",i);
1728:           A->factorerrortype             = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1729:           A->factorerror_zeropivot_value = 0.0;
1730:           A->factorerror_zeropivot_row   = i;
1731:         } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Zero diagonal on row %D",i);
1732:       }
1733:       idiag[i] = 1.0/v[diag[i]];
1734:     }
1735:     PetscLogFlops(m);
1736:   } else {
1737:     for (i=0; i<m; i++) {
1738:       mdiag[i] = v[diag[i]];
1739:       idiag[i] = omega/(fshift + v[diag[i]]);
1740:     }
1741:     PetscLogFlops(2.0*m);
1742:   }
1743:   a->idiagvalid = PETSC_TRUE;
1744:   return(0);
1745: }

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

1759:   its = its*lits;

1761:   if (fshift != a->fshift || omega != a->omega) a->idiagvalid = PETSC_FALSE; /* must recompute idiag[] */
1762:   if (!a->idiagvalid) {MatInvertDiagonal_SeqAIJ(A,omega,fshift);}
1763:   a->fshift = fshift;
1764:   a->omega  = omega;

1766:   diag  = a->diag;
1767:   t     = a->ssor_work;
1768:   idiag = a->idiag;
1769:   mdiag = a->mdiag;

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

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

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

1799:     to a vector efficiently using Eisenstat's trick.
1800:     */
1801:     scale = (2.0/omega) - 1.0;

1803:     /*  x = (E + U)^{-1} b */
1804:     for (i=m-1; i>=0; i--) {
1805:       n   = a->i[i+1] - diag[i] - 1;
1806:       idx = a->j + diag[i] + 1;
1807:       v   = a->a + diag[i] + 1;
1808:       sum = b[i];
1809:       PetscSparseDenseMinusDot(sum,x,v,idx,n);
1810:       x[i] = sum*idiag[i];
1811:     }

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

1817:     /*  t = (E + L)^{-1}t */
1818:     ts   = t;
1819:     diag = a->diag;
1820:     for (i=0; i<m; i++) {
1821:       n   = diag[i] - a->i[i];
1822:       idx = a->j + a->i[i];
1823:       v   = a->a + a->i[i];
1824:       sum = t[i];
1825:       PetscSparseDenseMinusDot(sum,ts,v,idx,n);
1826:       t[i] = sum*idiag[i];
1827:       /*  x = x + t */
1828:       x[i] += t[i];
1829:     }

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


1918: PetscErrorCode MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,MatInfo *info)
1919: {
1920:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;

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

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

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

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

2008: PetscErrorCode MatZeroRowsColumns_SeqAIJ(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
2009: {
2010:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
2011:   PetscInt          i,j,m = A->rmap->n - 1,d = 0;
2012:   PetscErrorCode    ierr;
2013:   PetscBool         missing,*zeroed,vecs = PETSC_FALSE;
2014:   const PetscScalar *xx;
2015:   PetscScalar       *bb;

2018:   if (x && b) {
2019:     VecGetArrayRead(x,&xx);
2020:     VecGetArray(b,&bb);
2021:     vecs = PETSC_TRUE;
2022:   }
2023:   PetscCalloc1(A->rmap->n,&zeroed);
2024:   for (i=0; i<N; i++) {
2025:     if (rows[i] < 0 || rows[i] > m) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range", rows[i]);
2026:     PetscMemzero(&a->a[a->i[rows[i]]],a->ilen[rows[i]]*sizeof(PetscScalar));

2028:     zeroed[rows[i]] = PETSC_TRUE;
2029:   }
2030:   for (i=0; i<A->rmap->n; i++) {
2031:     if (!zeroed[i]) {
2032:       for (j=a->i[i]; j<a->i[i+1]; j++) {
2033:         if (a->j[j] < A->rmap->n && zeroed[a->j[j]]) {
2034:           if (vecs) bb[i] -= a->a[j]*xx[a->j[j]];
2035:           a->a[j] = 0.0;
2036:         }
2037:       }
2038:     } else if (vecs && i < A->cmap->N) bb[i] = diag*xx[i];
2039:   }
2040:   if (x && b) {
2041:     VecRestoreArrayRead(x,&xx);
2042:     VecRestoreArray(b,&bb);
2043:   }
2044:   PetscFree(zeroed);
2045:   if (diag != 0.0) {
2046:     MatMissingDiagonal_SeqAIJ(A,&missing,&d);
2047:     if (missing) {
2048:       for (i=0; i<N; i++) {
2049:         if (rows[i] >= A->cmap->N) continue;
2050:         if (a->nonew && rows[i] >= d) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry in row %D (%D)",d,rows[i]);
2051:         MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],&diag,INSERT_VALUES);
2052:       }
2053:     } else {
2054:       for (i=0; i<N; i++) {
2055:         a->a[a->diag[rows[i]]] = diag;
2056:       }
2057:     }
2058:   }
2059:   (*A->ops->assemblyend)(A,MAT_FINAL_ASSEMBLY);
2060:   return(0);
2061: }

2063: PetscErrorCode MatGetRow_SeqAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
2064: {
2065:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
2066:   PetscInt   *itmp;

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

2071:   *nz = a->i[row+1] - a->i[row];
2072:   if (v) *v = a->a + a->i[row];
2073:   if (idx) {
2074:     itmp = a->j + a->i[row];
2075:     if (*nz) *idx = itmp;
2076:     else *idx = 0;
2077:   }
2078:   return(0);
2079: }

2081: /* remove this function? */
2082: PetscErrorCode MatRestoreRow_SeqAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
2083: {
2085:   return(0);
2086: }

2088: PetscErrorCode MatNorm_SeqAIJ(Mat A,NormType type,PetscReal *nrm)
2089: {
2090:   Mat_SeqAIJ     *a  = (Mat_SeqAIJ*)A->data;
2091:   MatScalar      *v  = a->a;
2092:   PetscReal      sum = 0.0;
2094:   PetscInt       i,j;

2097:   if (type == NORM_FROBENIUS) {
2098: #if defined(PETSC_USE_REAL___FP16)
2099:     PetscBLASInt one = 1,nz = a->nz;
2100:     *nrm = BLASnrm2_(&nz,v,&one);
2101: #else
2102:     for (i=0; i<a->nz; i++) {
2103:       sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
2104:     }
2105:     *nrm = PetscSqrtReal(sum);
2106: #endif
2107:     PetscLogFlops(2*a->nz);
2108:   } else if (type == NORM_1) {
2109:     PetscReal *tmp;
2110:     PetscInt  *jj = a->j;
2111:     PetscCalloc1(A->cmap->n+1,&tmp);
2112:     *nrm = 0.0;
2113:     for (j=0; j<a->nz; j++) {
2114:       tmp[*jj++] += PetscAbsScalar(*v);  v++;
2115:     }
2116:     for (j=0; j<A->cmap->n; j++) {
2117:       if (tmp[j] > *nrm) *nrm = tmp[j];
2118:     }
2119:     PetscFree(tmp);
2120:     PetscLogFlops(PetscMax(a->nz-1,0));
2121:   } else if (type == NORM_INFINITY) {
2122:     *nrm = 0.0;
2123:     for (j=0; j<A->rmap->n; j++) {
2124:       v   = a->a + a->i[j];
2125:       sum = 0.0;
2126:       for (i=0; i<a->i[j+1]-a->i[j]; i++) {
2127:         sum += PetscAbsScalar(*v); v++;
2128:       }
2129:       if (sum > *nrm) *nrm = sum;
2130:     }
2131:     PetscLogFlops(PetscMax(a->nz-1,0));
2132:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for two norm");
2133:   return(0);
2134: }

2136: /* Merged from MatGetSymbolicTranspose_SeqAIJ() - replace MatGetSymbolicTranspose_SeqAIJ()? */
2137: PetscErrorCode MatTransposeSymbolic_SeqAIJ(Mat A,Mat *B)
2138: {
2140:   PetscInt       i,j,anzj;
2141:   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b;
2142:   PetscInt       an=A->cmap->N,am=A->rmap->N;
2143:   PetscInt       *ati,*atj,*atfill,*ai=a->i,*aj=a->j;

2146:   /* Allocate space for symbolic transpose info and work array */
2147:   PetscCalloc1(an+1,&ati);
2148:   PetscMalloc1(ai[am],&atj);
2149:   PetscMalloc1(an,&atfill);

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

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

2160:   /* Walk through A row-wise and mark nonzero entries of A^T. */
2161:   for (i=0;i<am;i++) {
2162:     anzj = ai[i+1] - ai[i];
2163:     for (j=0;j<anzj;j++) {
2164:       atj[atfill[*aj]] = i;
2165:       atfill[*aj++]   += 1;
2166:     }
2167:   }

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

2174:   b          = (Mat_SeqAIJ*)((*B)->data);
2175:   b->free_a  = PETSC_FALSE;
2176:   b->free_ij = PETSC_TRUE;
2177:   b->nonew   = 0;
2178:   return(0);
2179: }

2181: PetscErrorCode  MatIsTranspose_SeqAIJ(Mat A,Mat B,PetscReal tol,PetscBool  *f)
2182: {
2183:   Mat_SeqAIJ     *aij = (Mat_SeqAIJ*) A->data,*bij = (Mat_SeqAIJ*) B->data;
2184:   PetscInt       *adx,*bdx,*aii,*bii,*aptr,*bptr;
2185:   MatScalar      *va,*vb;
2187:   PetscInt       ma,na,mb,nb, i;

2190:   MatGetSize(A,&ma,&na);
2191:   MatGetSize(B,&mb,&nb);
2192:   if (ma!=nb || na!=mb) {
2193:     *f = PETSC_FALSE;
2194:     return(0);
2195:   }
2196:   aii  = aij->i; bii = bij->i;
2197:   adx  = aij->j; bdx = bij->j;
2198:   va   = aij->a; vb = bij->a;
2199:   PetscMalloc1(ma,&aptr);
2200:   PetscMalloc1(mb,&bptr);
2201:   for (i=0; i<ma; i++) aptr[i] = aii[i];
2202:   for (i=0; i<mb; i++) bptr[i] = bii[i];

2204:   *f = PETSC_TRUE;
2205:   for (i=0; i<ma; i++) {
2206:     while (aptr[i]<aii[i+1]) {
2207:       PetscInt    idc,idr;
2208:       PetscScalar vc,vr;
2209:       /* column/row index/value */
2210:       idc = adx[aptr[i]];
2211:       idr = bdx[bptr[idc]];
2212:       vc  = va[aptr[i]];
2213:       vr  = vb[bptr[idc]];
2214:       if (i!=idr || PetscAbsScalar(vc-vr) > tol) {
2215:         *f = PETSC_FALSE;
2216:         goto done;
2217:       } else {
2218:         aptr[i]++;
2219:         if (B || i!=idc) bptr[idc]++;
2220:       }
2221:     }
2222:   }
2223: done:
2224:   PetscFree(aptr);
2225:   PetscFree(bptr);
2226:   return(0);
2227: }

2229: PetscErrorCode  MatIsHermitianTranspose_SeqAIJ(Mat A,Mat B,PetscReal tol,PetscBool  *f)
2230: {
2231:   Mat_SeqAIJ     *aij = (Mat_SeqAIJ*) A->data,*bij = (Mat_SeqAIJ*) B->data;
2232:   PetscInt       *adx,*bdx,*aii,*bii,*aptr,*bptr;
2233:   MatScalar      *va,*vb;
2235:   PetscInt       ma,na,mb,nb, i;

2238:   MatGetSize(A,&ma,&na);
2239:   MatGetSize(B,&mb,&nb);
2240:   if (ma!=nb || na!=mb) {
2241:     *f = PETSC_FALSE;
2242:     return(0);
2243:   }
2244:   aii  = aij->i; bii = bij->i;
2245:   adx  = aij->j; bdx = bij->j;
2246:   va   = aij->a; vb = bij->a;
2247:   PetscMalloc1(ma,&aptr);
2248:   PetscMalloc1(mb,&bptr);
2249:   for (i=0; i<ma; i++) aptr[i] = aii[i];
2250:   for (i=0; i<mb; i++) bptr[i] = bii[i];

2252:   *f = PETSC_TRUE;
2253:   for (i=0; i<ma; i++) {
2254:     while (aptr[i]<aii[i+1]) {
2255:       PetscInt    idc,idr;
2256:       PetscScalar vc,vr;
2257:       /* column/row index/value */
2258:       idc = adx[aptr[i]];
2259:       idr = bdx[bptr[idc]];
2260:       vc  = va[aptr[i]];
2261:       vr  = vb[bptr[idc]];
2262:       if (i!=idr || PetscAbsScalar(vc-PetscConj(vr)) > tol) {
2263:         *f = PETSC_FALSE;
2264:         goto done;
2265:       } else {
2266:         aptr[i]++;
2267:         if (B || i!=idc) bptr[idc]++;
2268:       }
2269:     }
2270:   }
2271: done:
2272:   PetscFree(aptr);
2273:   PetscFree(bptr);
2274:   return(0);
2275: }

2277: PetscErrorCode MatIsSymmetric_SeqAIJ(Mat A,PetscReal tol,PetscBool  *f)
2278: {

2282:   MatIsTranspose_SeqAIJ(A,A,tol,f);
2283:   return(0);
2284: }

2286: PetscErrorCode MatIsHermitian_SeqAIJ(Mat A,PetscReal tol,PetscBool  *f)
2287: {

2291:   MatIsHermitianTranspose_SeqAIJ(A,A,tol,f);
2292:   return(0);
2293: }

2295: PetscErrorCode MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr)
2296: {
2297:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
2298:   const PetscScalar *l,*r;
2299:   PetscScalar       x;
2300:   MatScalar         *v;
2301:   PetscErrorCode    ierr;
2302:   PetscInt          i,j,m = A->rmap->n,n = A->cmap->n,M,nz = a->nz;
2303:   const PetscInt    *jj;

2306:   if (ll) {
2307:     /* The local size is used so that VecMPI can be passed to this routine
2308:        by MatDiagonalScale_MPIAIJ */
2309:     VecGetLocalSize(ll,&m);
2310:     if (m != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length");
2311:     VecGetArrayRead(ll,&l);
2312:     v    = a->a;
2313:     for (i=0; i<m; i++) {
2314:       x = l[i];
2315:       M = a->i[i+1] - a->i[i];
2316:       for (j=0; j<M; j++) (*v++) *= x;
2317:     }
2318:     VecRestoreArrayRead(ll,&l);
2319:     PetscLogFlops(nz);
2320:   }
2321:   if (rr) {
2322:     VecGetLocalSize(rr,&n);
2323:     if (n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vector wrong length");
2324:     VecGetArrayRead(rr,&r);
2325:     v    = a->a; jj = a->j;
2326:     for (i=0; i<nz; i++) (*v++) *= r[*jj++];
2327:     VecRestoreArrayRead(rr,&r);
2328:     PetscLogFlops(nz);
2329:   }
2330:   MatSeqAIJInvalidateDiagonal(A);
2331:   return(0);
2332: }

2334: PetscErrorCode MatCreateSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,PetscInt csize,MatReuse scall,Mat *B)
2335: {
2336:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data,*c;
2338:   PetscInt       *smap,i,k,kstart,kend,oldcols = A->cmap->n,*lens;
2339:   PetscInt       row,mat_i,*mat_j,tcol,first,step,*mat_ilen,sum,lensi;
2340:   const PetscInt *irow,*icol;
2341:   PetscInt       nrows,ncols;
2342:   PetscInt       *starts,*j_new,*i_new,*aj = a->j,*ai = a->i,ii,*ailen = a->ilen;
2343:   MatScalar      *a_new,*mat_a;
2344:   Mat            C;
2345:   PetscBool      stride;


2349:   ISGetIndices(isrow,&irow);
2350:   ISGetLocalSize(isrow,&nrows);
2351:   ISGetLocalSize(iscol,&ncols);

2353:   PetscObjectTypeCompare((PetscObject)iscol,ISSTRIDE,&stride);
2354:   if (stride) {
2355:     ISStrideGetInfo(iscol,&first,&step);
2356:   } else {
2357:     first = 0;
2358:     step  = 0;
2359:   }
2360:   if (stride && step == 1) {
2361:     /* special case of contiguous rows */
2362:     PetscMalloc2(nrows,&lens,nrows,&starts);
2363:     /* loop over new rows determining lens and starting points */
2364:     for (i=0; i<nrows; i++) {
2365:       kstart = ai[irow[i]];
2366:       kend   = kstart + ailen[irow[i]];
2367:       starts[i] = kstart;
2368:       for (k=kstart; k<kend; k++) {
2369:         if (aj[k] >= first) {
2370:           starts[i] = k;
2371:           break;
2372:         }
2373:       }
2374:       sum = 0;
2375:       while (k < kend) {
2376:         if (aj[k++] >= first+ncols) break;
2377:         sum++;
2378:       }
2379:       lens[i] = sum;
2380:     }
2381:     /* create submatrix */
2382:     if (scall == MAT_REUSE_MATRIX) {
2383:       PetscInt n_cols,n_rows;
2384:       MatGetSize(*B,&n_rows,&n_cols);
2385:       if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size");
2386:       MatZeroEntries(*B);
2387:       C    = *B;
2388:     } else {
2389:       PetscInt rbs,cbs;
2390:       MatCreate(PetscObjectComm((PetscObject)A),&C);
2391:       MatSetSizes(C,nrows,ncols,PETSC_DETERMINE,PETSC_DETERMINE);
2392:       ISGetBlockSize(isrow,&rbs);
2393:       ISGetBlockSize(iscol,&cbs);
2394:       MatSetBlockSizes(C,rbs,cbs);
2395:       MatSetType(C,((PetscObject)A)->type_name);
2396:       MatSeqAIJSetPreallocation_SeqAIJ(C,0,lens);
2397:     }
2398:     c = (Mat_SeqAIJ*)C->data;

2400:     /* loop over rows inserting into submatrix */
2401:     a_new = c->a;
2402:     j_new = c->j;
2403:     i_new = c->i;

2405:     for (i=0; i<nrows; i++) {
2406:       ii    = starts[i];
2407:       lensi = lens[i];
2408:       for (k=0; k<lensi; k++) {
2409:         *j_new++ = aj[ii+k] - first;
2410:       }
2411:       PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(PetscScalar));
2412:       a_new     += lensi;
2413:       i_new[i+1] = i_new[i] + lensi;
2414:       c->ilen[i] = lensi;
2415:     }
2416:     PetscFree2(lens,starts);
2417:   } else {
2418:     ISGetIndices(iscol,&icol);
2419:     PetscCalloc1(oldcols,&smap);
2420:     PetscMalloc1(1+nrows,&lens);
2421:     for (i=0; i<ncols; i++) {
2422: #if defined(PETSC_USE_DEBUG)
2423:       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);
2424: #endif
2425:       smap[icol[i]] = i+1;
2426:     }

2428:     /* determine lens of each row */
2429:     for (i=0; i<nrows; i++) {
2430:       kstart  = ai[irow[i]];
2431:       kend    = kstart + a->ilen[irow[i]];
2432:       lens[i] = 0;
2433:       for (k=kstart; k<kend; k++) {
2434:         if (smap[aj[k]]) {
2435:           lens[i]++;
2436:         }
2437:       }
2438:     }
2439:     /* Create and fill new matrix */
2440:     if (scall == MAT_REUSE_MATRIX) {
2441:       PetscBool equal;

2443:       c = (Mat_SeqAIJ*)((*B)->data);
2444:       if ((*B)->rmap->n  != nrows || (*B)->cmap->n != ncols) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size");
2445:       PetscMemcmp(c->ilen,lens,(*B)->rmap->n*sizeof(PetscInt),&equal);
2446:       if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros");
2447:       PetscMemzero(c->ilen,(*B)->rmap->n*sizeof(PetscInt));
2448:       C    = *B;
2449:     } else {
2450:       PetscInt rbs,cbs;
2451:       MatCreate(PetscObjectComm((PetscObject)A),&C);
2452:       MatSetSizes(C,nrows,ncols,PETSC_DETERMINE,PETSC_DETERMINE);
2453:       ISGetBlockSize(isrow,&rbs);
2454:       ISGetBlockSize(iscol,&cbs);
2455:       MatSetBlockSizes(C,rbs,cbs);
2456:       MatSetType(C,((PetscObject)A)->type_name);
2457:       MatSeqAIJSetPreallocation_SeqAIJ(C,0,lens);
2458:     }
2459:     c = (Mat_SeqAIJ*)(C->data);
2460:     for (i=0; i<nrows; i++) {
2461:       row      = irow[i];
2462:       kstart   = ai[row];
2463:       kend     = kstart + a->ilen[row];
2464:       mat_i    = c->i[i];
2465:       mat_j    = c->j + mat_i;
2466:       mat_a    = c->a + mat_i;
2467:       mat_ilen = c->ilen + i;
2468:       for (k=kstart; k<kend; k++) {
2469:         if ((tcol=smap[a->j[k]])) {
2470:           *mat_j++ = tcol - 1;
2471:           *mat_a++ = a->a[k];
2472:           (*mat_ilen)++;

2474:         }
2475:       }
2476:     }
2477:     /* Free work space */
2478:     ISRestoreIndices(iscol,&icol);
2479:     PetscFree(smap);
2480:     PetscFree(lens);
2481:     /* sort */
2482:     for (i = 0; i < nrows; i++) {
2483:       PetscInt ilen;

2485:       mat_i = c->i[i];
2486:       mat_j = c->j + mat_i;
2487:       mat_a = c->a + mat_i;
2488:       ilen  = c->ilen[i];
2489:       PetscSortIntWithScalarArray(ilen,mat_j,mat_a);
2490:     }
2491:   }
2492:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
2493:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

2495:   ISRestoreIndices(isrow,&irow);
2496:   *B   = C;
2497:   return(0);
2498: }

2500: PetscErrorCode  MatGetMultiProcBlock_SeqAIJ(Mat mat,MPI_Comm subComm,MatReuse scall,Mat *subMat)
2501: {
2503:   Mat            B;

2506:   if (scall == MAT_INITIAL_MATRIX) {
2507:     MatCreate(subComm,&B);
2508:     MatSetSizes(B,mat->rmap->n,mat->cmap->n,mat->rmap->n,mat->cmap->n);
2509:     MatSetBlockSizesFromMats(B,mat,mat);
2510:     MatSetType(B,MATSEQAIJ);
2511:     MatDuplicateNoCreate_SeqAIJ(B,mat,MAT_COPY_VALUES,PETSC_TRUE);
2512:     *subMat = B;
2513:   } else {
2514:     MatCopy_SeqAIJ(mat,*subMat,SAME_NONZERO_PATTERN);
2515:   }
2516:   return(0);
2517: }

2519: PetscErrorCode MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,const MatFactorInfo *info)
2520: {
2521:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)inA->data;
2523:   Mat            outA;
2524:   PetscBool      row_identity,col_identity;

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

2529:   ISIdentity(row,&row_identity);
2530:   ISIdentity(col,&col_identity);

2532:   outA             = inA;
2533:   outA->factortype = MAT_FACTOR_LU;
2534:   PetscFree(inA->solvertype);
2535:   PetscStrallocpy(MATSOLVERPETSC,&inA->solvertype);

2537:   PetscObjectReference((PetscObject)row);
2538:   ISDestroy(&a->row);

2540:   a->row = row;

2542:   PetscObjectReference((PetscObject)col);
2543:   ISDestroy(&a->col);

2545:   a->col = col;

2547:   /* Create the inverse permutation so that it can be used in MatLUFactorNumeric() */
2548:   ISDestroy(&a->icol);
2549:   ISInvertPermutation(col,PETSC_DECIDE,&a->icol);
2550:   PetscLogObjectParent((PetscObject)inA,(PetscObject)a->icol);

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

2557:   MatMarkDiagonal_SeqAIJ(inA);
2558:   if (row_identity && col_identity) {
2559:     MatLUFactorNumeric_SeqAIJ_inplace(outA,inA,info);
2560:   } else {
2561:     MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(outA,inA,info);
2562:   }
2563:   return(0);
2564: }

2566: PetscErrorCode MatScale_SeqAIJ(Mat inA,PetscScalar alpha)
2567: {
2568:   Mat_SeqAIJ     *a     = (Mat_SeqAIJ*)inA->data;
2569:   PetscScalar    oalpha = alpha;
2571:   PetscBLASInt   one = 1,bnz;

2574:   PetscBLASIntCast(a->nz,&bnz);
2575:   PetscStackCallBLAS("BLASscal",BLASscal_(&bnz,&oalpha,a->a,&one));
2576:   PetscLogFlops(a->nz);
2577:   MatSeqAIJInvalidateDiagonal(inA);
2578:   return(0);
2579: }

2581: PetscErrorCode MatDestroySubMatrix_Private(Mat_SubSppt *submatj)
2582: {
2584:   PetscInt       i;

2587:   if (!submatj->id) { /* delete data that are linked only to submats[id=0] */
2588:     PetscFree4(submatj->sbuf1,submatj->ptr,submatj->tmp,submatj->ctr);

2590:     for (i=0; i<submatj->nrqr; ++i) {
2591:       PetscFree(submatj->sbuf2[i]);
2592:     }
2593:     PetscFree3(submatj->sbuf2,submatj->req_size,submatj->req_source1);

2595:     if (submatj->rbuf1) {
2596:       PetscFree(submatj->rbuf1[0]);
2597:       PetscFree(submatj->rbuf1);
2598:     }

2600:     for (i=0; i<submatj->nrqs; ++i) {
2601:       PetscFree(submatj->rbuf3[i]);
2602:     }
2603:     PetscFree3(submatj->req_source2,submatj->rbuf2,submatj->rbuf3);
2604:     PetscFree(submatj->pa);
2605:   }

2607: #if defined(PETSC_USE_CTABLE)
2608:   PetscTableDestroy((PetscTable*)&submatj->rmap);
2609:   if (submatj->cmap_loc) {PetscFree(submatj->cmap_loc);}
2610:   PetscFree(submatj->rmap_loc);
2611: #else
2612:   PetscFree(submatj->rmap);
2613: #endif

2615:   if (!submatj->allcolumns) {
2616: #if defined(PETSC_USE_CTABLE)
2617:     PetscTableDestroy((PetscTable*)&submatj->cmap);
2618: #else
2619:     PetscFree(submatj->cmap);
2620: #endif
2621:   }
2622:   PetscFree(submatj->row2proc);

2624:   PetscFree(submatj);
2625:   return(0);
2626: }

2628: PetscErrorCode MatDestroySubMatrix_SeqAIJ(Mat C)
2629: {
2631:   Mat_SeqAIJ     *c = (Mat_SeqAIJ*)C->data;
2632:   Mat_SubSppt    *submatj = c->submatis1;

2635:   (*submatj->destroy)(C);
2636:   MatDestroySubMatrix_Private(submatj);
2637:   return(0);
2638: }

2640: PetscErrorCode MatDestroySubMatrices_SeqAIJ(PetscInt n,Mat *mat[])
2641: {
2643:   PetscInt       i;
2644:   Mat            C;
2645:   Mat_SeqAIJ     *c;
2646:   Mat_SubSppt    *submatj;

2649:   for (i=0; i<n; i++) {
2650:     C       = (*mat)[i];
2651:     c       = (Mat_SeqAIJ*)C->data;
2652:     submatj = c->submatis1;
2653:     if (submatj) {
2654:       if (--((PetscObject)C)->refct <= 0) {
2655:         (*submatj->destroy)(C);
2656:         MatDestroySubMatrix_Private(submatj);
2657:         PetscFree(C->defaultvectype);
2658:         PetscLayoutDestroy(&C->rmap);
2659:         PetscLayoutDestroy(&C->cmap);
2660:         PetscHeaderDestroy(&C);
2661:       }
2662:     } else {
2663:       MatDestroy(&C);
2664:     }
2665:   }

2667:   /* Destroy Dummy submatrices created for reuse */
2668:   MatDestroySubMatrices_Dummy(n,mat);

2670:   PetscFree(*mat);
2671:   return(0);
2672: }

2674: PetscErrorCode MatCreateSubMatrices_SeqAIJ(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
2675: {
2677:   PetscInt       i;

2680:   if (scall == MAT_INITIAL_MATRIX) {
2681:     PetscCalloc1(n+1,B);
2682:   }

2684:   for (i=0; i<n; i++) {
2685:     MatCreateSubMatrix_SeqAIJ(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);
2686:   }
2687:   return(0);
2688: }

2690: PetscErrorCode MatIncreaseOverlap_SeqAIJ(Mat A,PetscInt is_max,IS is[],PetscInt ov)
2691: {
2692:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
2694:   PetscInt       row,i,j,k,l,m,n,*nidx,isz,val;
2695:   const PetscInt *idx;
2696:   PetscInt       start,end,*ai,*aj;
2697:   PetscBT        table;

2700:   m  = A->rmap->n;
2701:   ai = a->i;
2702:   aj = a->j;

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

2706:   PetscMalloc1(m+1,&nidx);
2707:   PetscBTCreate(m,&table);

2709:   for (i=0; i<is_max; i++) {
2710:     /* Initialize the two local arrays */
2711:     isz  = 0;
2712:     PetscBTMemzero(m,table);

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

2718:     /* Enter these into the temp arrays. I.e., mark table[row], enter row into new index */
2719:     for (j=0; j<n; ++j) {
2720:       if (!PetscBTLookupSet(table,idx[j])) nidx[isz++] = idx[j];
2721:     }
2722:     ISRestoreIndices(is[i],&idx);
2723:     ISDestroy(&is[i]);

2725:     k = 0;
2726:     for (j=0; j<ov; j++) { /* for each overlap */
2727:       n = isz;
2728:       for (; k<n; k++) { /* do only those rows in nidx[k], which are not done yet */
2729:         row   = nidx[k];
2730:         start = ai[row];
2731:         end   = ai[row+1];
2732:         for (l = start; l<end; l++) {
2733:           val = aj[l];
2734:           if (!PetscBTLookupSet(table,val)) nidx[isz++] = val;
2735:         }
2736:       }
2737:     }
2738:     ISCreateGeneral(PETSC_COMM_SELF,isz,nidx,PETSC_COPY_VALUES,(is+i));
2739:   }
2740:   PetscBTDestroy(&table);
2741:   PetscFree(nidx);
2742:   return(0);
2743: }

2745: /* -------------------------------------------------------------- */
2746: PetscErrorCode MatPermute_SeqAIJ(Mat A,IS rowp,IS colp,Mat *B)
2747: {
2748:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
2750:   PetscInt       i,nz = 0,m = A->rmap->n,n = A->cmap->n;
2751:   const PetscInt *row,*col;
2752:   PetscInt       *cnew,j,*lens;
2753:   IS             icolp,irowp;
2754:   PetscInt       *cwork = NULL;
2755:   PetscScalar    *vwork = NULL;

2758:   ISInvertPermutation(rowp,PETSC_DECIDE,&irowp);
2759:   ISGetIndices(irowp,&row);
2760:   ISInvertPermutation(colp,PETSC_DECIDE,&icolp);
2761:   ISGetIndices(icolp,&col);

2763:   /* determine lengths of permuted rows */
2764:   PetscMalloc1(m+1,&lens);
2765:   for (i=0; i<m; i++) lens[row[i]] = a->i[i+1] - a->i[i];
2766:   MatCreate(PetscObjectComm((PetscObject)A),B);
2767:   MatSetSizes(*B,m,n,m,n);
2768:   MatSetBlockSizesFromMats(*B,A,A);
2769:   MatSetType(*B,((PetscObject)A)->type_name);
2770:   MatSeqAIJSetPreallocation_SeqAIJ(*B,0,lens);
2771:   PetscFree(lens);

2773:   PetscMalloc1(n,&cnew);
2774:   for (i=0; i<m; i++) {
2775:     MatGetRow_SeqAIJ(A,i,&nz,&cwork,&vwork);
2776:     for (j=0; j<nz; j++) cnew[j] = col[cwork[j]];
2777:     MatSetValues_SeqAIJ(*B,1,&row[i],nz,cnew,vwork,INSERT_VALUES);
2778:     MatRestoreRow_SeqAIJ(A,i,&nz,&cwork,&vwork);
2779:   }
2780:   PetscFree(cnew);

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

2784:   MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);
2785:   MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);
2786:   ISRestoreIndices(irowp,&row);
2787:   ISRestoreIndices(icolp,&col);
2788:   ISDestroy(&irowp);
2789:   ISDestroy(&icolp);
2790:   return(0);
2791: }

2793: PetscErrorCode MatCopy_SeqAIJ(Mat A,Mat B,MatStructure str)
2794: {

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

2803:     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");
2804:     PetscMemcpy(b->a,a->a,(a->i[A->rmap->n])*sizeof(PetscScalar));
2805:     PetscObjectStateIncrease((PetscObject)B);
2806:   } else {
2807:     MatCopy_Basic(A,B,str);
2808:   }
2809:   return(0);
2810: }

2812: PetscErrorCode MatSetUp_SeqAIJ(Mat A)
2813: {

2817:    MatSeqAIJSetPreallocation_SeqAIJ(A,PETSC_DEFAULT,0);
2818:   return(0);
2819: }

2821: PetscErrorCode MatSeqAIJGetArray_SeqAIJ(Mat A,PetscScalar *array[])
2822: {
2823:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;

2826:   *array = a->a;
2827:   return(0);
2828: }

2830: PetscErrorCode MatSeqAIJRestoreArray_SeqAIJ(Mat A,PetscScalar *array[])
2831: {
2833:   return(0);
2834: }

2836: /*
2837:    Computes the number of nonzeros per row needed for preallocation when X and Y
2838:    have different nonzero structure.
2839: */
2840: PetscErrorCode MatAXPYGetPreallocation_SeqX_private(PetscInt m,const PetscInt *xi,const PetscInt *xj,const PetscInt *yi,const PetscInt *yj,PetscInt *nnz)
2841: {
2842:   PetscInt       i,j,k,nzx,nzy;

2845:   /* Set the number of nonzeros in the new matrix */
2846:   for (i=0; i<m; i++) {
2847:     const PetscInt *xjj = xj+xi[i],*yjj = yj+yi[i];
2848:     nzx = xi[i+1] - xi[i];
2849:     nzy = yi[i+1] - yi[i];
2850:     nnz[i] = 0;
2851:     for (j=0,k=0; j<nzx; j++) {                   /* Point in X */
2852:       for (; k<nzy && yjj[k]<xjj[j]; k++) nnz[i]++; /* Catch up to X */
2853:       if (k<nzy && yjj[k]==xjj[j]) k++;             /* Skip duplicate */
2854:       nnz[i]++;
2855:     }
2856:     for (; k<nzy; k++) nnz[i]++;
2857:   }
2858:   return(0);
2859: }

2861: PetscErrorCode MatAXPYGetPreallocation_SeqAIJ(Mat Y,Mat X,PetscInt *nnz)
2862: {
2863:   PetscInt       m = Y->rmap->N;
2864:   Mat_SeqAIJ     *x = (Mat_SeqAIJ*)X->data;
2865:   Mat_SeqAIJ     *y = (Mat_SeqAIJ*)Y->data;

2869:   /* Set the number of nonzeros in the new matrix */
2870:   MatAXPYGetPreallocation_SeqX_private(m,x->i,x->j,y->i,y->j,nnz);
2871:   return(0);
2872: }

2874: PetscErrorCode MatAXPY_SeqAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str)
2875: {
2877:   Mat_SeqAIJ     *x = (Mat_SeqAIJ*)X->data,*y = (Mat_SeqAIJ*)Y->data;
2878:   PetscBLASInt   one=1,bnz;

2881:   PetscBLASIntCast(x->nz,&bnz);
2882:   if (str == SAME_NONZERO_PATTERN) {
2883:     PetscScalar alpha = a;
2884:     PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&bnz,&alpha,x->a,&one,y->a,&one));
2885:     MatSeqAIJInvalidateDiagonal(Y);
2886:     PetscObjectStateIncrease((PetscObject)Y);
2887:   } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
2888:     MatAXPY_Basic(Y,a,X,str);
2889:   } else {
2890:     Mat      B;
2891:     PetscInt *nnz;
2892:     PetscMalloc1(Y->rmap->N,&nnz);
2893:     MatCreate(PetscObjectComm((PetscObject)Y),&B);
2894:     PetscObjectSetName((PetscObject)B,((PetscObject)Y)->name);
2895:     MatSetSizes(B,Y->rmap->n,Y->cmap->n,Y->rmap->N,Y->cmap->N);
2896:     MatSetBlockSizesFromMats(B,Y,Y);
2897:     MatSetType(B,(MatType) ((PetscObject)Y)->type_name);
2898:     MatAXPYGetPreallocation_SeqAIJ(Y,X,nnz);
2899:     MatSeqAIJSetPreallocation(B,0,nnz);
2900:     MatAXPY_BasicWithPreallocation(B,Y,a,X,str);
2901:     MatHeaderReplace(Y,&B);
2902:     PetscFree(nnz);
2903:   }
2904:   return(0);
2905: }

2907: PetscErrorCode  MatConjugate_SeqAIJ(Mat mat)
2908: {
2909: #if defined(PETSC_USE_COMPLEX)
2910:   Mat_SeqAIJ  *aij = (Mat_SeqAIJ*)mat->data;
2911:   PetscInt    i,nz;
2912:   PetscScalar *a;

2915:   nz = aij->nz;
2916:   a  = aij->a;
2917:   for (i=0; i<nz; i++) a[i] = PetscConj(a[i]);
2918: #else
2920: #endif
2921:   return(0);
2922: }

2924: PetscErrorCode MatGetRowMaxAbs_SeqAIJ(Mat A,Vec v,PetscInt idx[])
2925: {
2926:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
2928:   PetscInt       i,j,m = A->rmap->n,*ai,*aj,ncols,n;
2929:   PetscReal      atmp;
2930:   PetscScalar    *x;
2931:   MatScalar      *aa;

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

2939:   VecSet(v,0.0);
2940:   VecGetArray(v,&x);
2941:   VecGetLocalSize(v,&n);
2942:   if (n != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
2943:   for (i=0; i<m; i++) {
2944:     ncols = ai[1] - ai[0]; ai++;
2945:     x[i]  = 0.0;
2946:     for (j=0; j<ncols; j++) {
2947:       atmp = PetscAbsScalar(*aa);
2948:       if (PetscAbsScalar(x[i]) < atmp) {x[i] = atmp; if (idx) idx[i] = *aj;}
2949:       aa++; aj++;
2950:     }
2951:   }
2952:   VecRestoreArray(v,&x);
2953:   return(0);
2954: }

2956: PetscErrorCode MatGetRowMax_SeqAIJ(Mat A,Vec v,PetscInt idx[])
2957: {
2958:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
2960:   PetscInt       i,j,m = A->rmap->n,*ai,*aj,ncols,n;
2961:   PetscScalar    *x;
2962:   MatScalar      *aa;

2965:   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2966:   aa = a->a;
2967:   ai = a->i;
2968:   aj = a->j;

2970:   VecSet(v,0.0);
2971:   VecGetArray(v,&x);
2972:   VecGetLocalSize(v,&n);
2973:   if (n != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
2974:   for (i=0; i<m; i++) {
2975:     ncols = ai[1] - ai[0]; ai++;
2976:     if (ncols == A->cmap->n) { /* row is dense */
2977:       x[i] = *aa; if (idx) idx[i] = 0;
2978:     } else {  /* row is sparse so already KNOW maximum is 0.0 or higher */
2979:       x[i] = 0.0;
2980:       if (idx) {
2981:         idx[i] = 0; /* in case ncols is zero */
2982:         for (j=0;j<ncols;j++) { /* find first implicit 0.0 in the row */
2983:           if (aj[j] > j) {
2984:             idx[i] = j;
2985:             break;
2986:           }
2987:         }
2988:       }
2989:     }
2990:     for (j=0; j<ncols; j++) {
2991:       if (PetscRealPart(x[i]) < PetscRealPart(*aa)) {x[i] = *aa; if (idx) idx[i] = *aj;}
2992:       aa++; aj++;
2993:     }
2994:   }
2995:   VecRestoreArray(v,&x);
2996:   return(0);
2997: }

2999: PetscErrorCode MatGetRowMinAbs_SeqAIJ(Mat A,Vec v,PetscInt idx[])
3000: {
3001:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
3003:   PetscInt       i,j,m = A->rmap->n,*ai,*aj,ncols,n;
3004:   PetscReal      atmp;
3005:   PetscScalar    *x;
3006:   MatScalar      *aa;

3009:   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
3010:   aa = a->a;
3011:   ai = a->i;
3012:   aj = a->j;

3014:   VecSet(v,0.0);
3015:   VecGetArray(v,&x);
3016:   VecGetLocalSize(v,&n);
3017:   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);
3018:   for (i=0; i<m; i++) {
3019:     ncols = ai[1] - ai[0]; ai++;
3020:     if (ncols) {
3021:       /* Get first nonzero */
3022:       for (j = 0; j < ncols; j++) {
3023:         atmp = PetscAbsScalar(aa[j]);
3024:         if (atmp > 1.0e-12) {
3025:           x[i] = atmp;
3026:           if (idx) idx[i] = aj[j];
3027:           break;
3028:         }
3029:       }
3030:       if (j == ncols) {x[i] = PetscAbsScalar(*aa); if (idx) idx[i] = *aj;}
3031:     } else {
3032:       x[i] = 0.0; if (idx) idx[i] = 0;
3033:     }
3034:     for (j = 0; j < ncols; j++) {
3035:       atmp = PetscAbsScalar(*aa);
3036:       if (atmp > 1.0e-12 && PetscAbsScalar(x[i]) > atmp) {x[i] = atmp; if (idx) idx[i] = *aj;}
3037:       aa++; aj++;
3038:     }
3039:   }
3040:   VecRestoreArray(v,&x);
3041:   return(0);
3042: }

3044: PetscErrorCode MatGetRowMin_SeqAIJ(Mat A,Vec v,PetscInt idx[])
3045: {
3046:   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
3047:   PetscErrorCode  ierr;
3048:   PetscInt        i,j,m = A->rmap->n,ncols,n;
3049:   const PetscInt  *ai,*aj;
3050:   PetscScalar     *x;
3051:   const MatScalar *aa;

3054:   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
3055:   aa = a->a;
3056:   ai = a->i;
3057:   aj = a->j;

3059:   VecSet(v,0.0);
3060:   VecGetArray(v,&x);
3061:   VecGetLocalSize(v,&n);
3062:   if (n != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
3063:   for (i=0; i<m; i++) {
3064:     ncols = ai[1] - ai[0]; ai++;
3065:     if (ncols == A->cmap->n) { /* row is dense */
3066:       x[i] = *aa; if (idx) idx[i] = 0;
3067:     } else {  /* row is sparse so already KNOW minimum is 0.0 or lower */
3068:       x[i] = 0.0;
3069:       if (idx) {   /* find first implicit 0.0 in the row */
3070:         idx[i] = 0; /* in case ncols is zero */
3071:         for (j=0; j<ncols; j++) {
3072:           if (aj[j] > j) {
3073:             idx[i] = j;
3074:             break;
3075:           }
3076:         }
3077:       }
3078:     }
3079:     for (j=0; j<ncols; j++) {
3080:       if (PetscRealPart(x[i]) > PetscRealPart(*aa)) {x[i] = *aa; if (idx) idx[i] = *aj;}
3081:       aa++; aj++;
3082:     }
3083:   }
3084:   VecRestoreArray(v,&x);
3085:   return(0);
3086: }

3088: PetscErrorCode MatInvertBlockDiagonal_SeqAIJ(Mat A,const PetscScalar **values)
3089: {
3090:   Mat_SeqAIJ      *a = (Mat_SeqAIJ*) A->data;
3091:   PetscErrorCode  ierr;
3092:   PetscInt        i,bs = PetscAbs(A->rmap->bs),mbs = A->rmap->n/bs,ipvt[5],bs2 = bs*bs,*v_pivots,ij[7],*IJ,j;
3093:   MatScalar       *diag,work[25],*v_work;
3094:   const PetscReal shift = 0.0;
3095:   PetscBool       allowzeropivot,zeropivotdetected=PETSC_FALSE;

3098:   allowzeropivot = PetscNot(A->erroriffailure);
3099:   if (a->ibdiagvalid) {
3100:     if (values) *values = a->ibdiag;
3101:     return(0);
3102:   }
3103:   MatMarkDiagonal_SeqAIJ(A);
3104:   if (!a->ibdiag) {
3105:     PetscMalloc1(bs2*mbs,&a->ibdiag);
3106:     PetscLogObjectMemory((PetscObject)A,bs2*mbs*sizeof(PetscScalar));
3107:   }
3108:   diag = a->ibdiag;
3109:   if (values) *values = a->ibdiag;
3110:   /* factor and invert each block */
3111:   switch (bs) {
3112:   case 1:
3113:     for (i=0; i<mbs; i++) {
3114:       MatGetValues(A,1,&i,1,&i,diag+i);
3115:       if (PetscAbsScalar(diag[i] + shift) < PETSC_MACHINE_EPSILON) {
3116:         if (allowzeropivot) {
3117:           A->factorerrortype             = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3118:           A->factorerror_zeropivot_value = PetscAbsScalar(diag[i]);
3119:           A->factorerror_zeropivot_row   = i;
3120:           PetscInfo3(A,"Zero pivot, row %D pivot %g tolerance %g\n",i,(double)PetscAbsScalar(diag[i]),(double)PETSC_MACHINE_EPSILON);
3121:         } else SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot, row %D pivot %g tolerance %g",i,(double)PetscAbsScalar(diag[i]),(double)PETSC_MACHINE_EPSILON);
3122:       }
3123:       diag[i] = (PetscScalar)1.0 / (diag[i] + shift);
3124:     }
3125:     break;
3126:   case 2:
3127:     for (i=0; i<mbs; i++) {
3128:       ij[0] = 2*i; ij[1] = 2*i + 1;
3129:       MatGetValues(A,2,ij,2,ij,diag);
3130:       PetscKernel_A_gets_inverse_A_2(diag,shift,allowzeropivot,&zeropivotdetected);
3131:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3132:       PetscKernel_A_gets_transpose_A_2(diag);
3133:       diag += 4;
3134:     }
3135:     break;
3136:   case 3:
3137:     for (i=0; i<mbs; i++) {
3138:       ij[0] = 3*i; ij[1] = 3*i + 1; ij[2] = 3*i + 2;
3139:       MatGetValues(A,3,ij,3,ij,diag);
3140:       PetscKernel_A_gets_inverse_A_3(diag,shift,allowzeropivot,&zeropivotdetected);
3141:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3142:       PetscKernel_A_gets_transpose_A_3(diag);
3143:       diag += 9;
3144:     }
3145:     break;
3146:   case 4:
3147:     for (i=0; i<mbs; i++) {
3148:       ij[0] = 4*i; ij[1] = 4*i + 1; ij[2] = 4*i + 2; ij[3] = 4*i + 3;
3149:       MatGetValues(A,4,ij,4,ij,diag);
3150:       PetscKernel_A_gets_inverse_A_4(diag,shift,allowzeropivot,&zeropivotdetected);
3151:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3152:       PetscKernel_A_gets_transpose_A_4(diag);
3153:       diag += 16;
3154:     }
3155:     break;
3156:   case 5:
3157:     for (i=0; i<mbs; i++) {
3158:       ij[0] = 5*i; ij[1] = 5*i + 1; ij[2] = 5*i + 2; ij[3] = 5*i + 3; ij[4] = 5*i + 4;
3159:       MatGetValues(A,5,ij,5,ij,diag);
3160:       PetscKernel_A_gets_inverse_A_5(diag,ipvt,work,shift,allowzeropivot,&zeropivotdetected);
3161:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3162:       PetscKernel_A_gets_transpose_A_5(diag);
3163:       diag += 25;
3164:     }
3165:     break;
3166:   case 6:
3167:     for (i=0; i<mbs; i++) {
3168:       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;
3169:       MatGetValues(A,6,ij,6,ij,diag);
3170:       PetscKernel_A_gets_inverse_A_6(diag,shift,allowzeropivot,&zeropivotdetected);
3171:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3172:       PetscKernel_A_gets_transpose_A_6(diag);
3173:       diag += 36;
3174:     }
3175:     break;
3176:   case 7:
3177:     for (i=0; i<mbs; i++) {
3178:       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;
3179:       MatGetValues(A,7,ij,7,ij,diag);
3180:       PetscKernel_A_gets_inverse_A_7(diag,shift,allowzeropivot,&zeropivotdetected);
3181:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3182:       PetscKernel_A_gets_transpose_A_7(diag);
3183:       diag += 49;
3184:     }
3185:     break;
3186:   default:
3187:     PetscMalloc3(bs,&v_work,bs,&v_pivots,bs,&IJ);
3188:     for (i=0; i<mbs; i++) {
3189:       for (j=0; j<bs; j++) {
3190:         IJ[j] = bs*i + j;
3191:       }
3192:       MatGetValues(A,bs,IJ,bs,IJ,diag);
3193:       PetscKernel_A_gets_inverse_A(bs,diag,v_pivots,v_work,allowzeropivot,&zeropivotdetected);
3194:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3195:       PetscKernel_A_gets_transpose_A_N(diag,bs);
3196:       diag += bs2;
3197:     }
3198:     PetscFree3(v_work,v_pivots,IJ);
3199:   }
3200:   a->ibdiagvalid = PETSC_TRUE;
3201:   return(0);
3202: }

3204: static PetscErrorCode  MatSetRandom_SeqAIJ(Mat x,PetscRandom rctx)
3205: {
3207:   Mat_SeqAIJ     *aij = (Mat_SeqAIJ*)x->data;
3208:   PetscScalar    a;
3209:   PetscInt       m,n,i,j,col;

3212:   if (!x->assembled) {
3213:     MatGetSize(x,&m,&n);
3214:     for (i=0; i<m; i++) {
3215:       for (j=0; j<aij->imax[i]; j++) {
3216:         PetscRandomGetValue(rctx,&a);
3217:         col  = (PetscInt)(n*PetscRealPart(a));
3218:         MatSetValues(x,1,&i,1,&col,&a,ADD_VALUES);
3219:       }
3220:     }
3221:   } else {
3222:     for (i=0; i<aij->nz; i++) {PetscRandomGetValue(rctx,aij->a+i);}
3223:   }
3224:   MatAssemblyBegin(x,MAT_FINAL_ASSEMBLY);
3225:   MatAssemblyEnd(x,MAT_FINAL_ASSEMBLY);
3226:   return(0);
3227: }

3229: /* Like MatSetRandom_SeqAIJ, but do not set values on columns in range of [low, high) */
3230: PetscErrorCode  MatSetRandomSkipColumnRange_SeqAIJ_Private(Mat x,PetscInt low,PetscInt high,PetscRandom rctx)
3231: {
3233:   Mat_SeqAIJ     *aij = (Mat_SeqAIJ*)x->data;
3234:   PetscScalar    a;
3235:   PetscInt       m,n,i,j,col,nskip;

3238:   nskip = high - low;
3239:   MatGetSize(x,&m,&n);
3240:   n    -= nskip; /* shrink number of columns where nonzeros can be set */
3241:   for (i=0; i<m; i++) {
3242:     for (j=0; j<aij->imax[i]; j++) {
3243:       PetscRandomGetValue(rctx,&a);
3244:       col  = (PetscInt)(n*PetscRealPart(a));
3245:       if (col >= low) col += nskip; /* shift col rightward to skip the hole */
3246:       MatSetValues(x,1,&i,1,&col,&a,ADD_VALUES);
3247:     }
3248:   }
3249:   MatAssemblyBegin(x,MAT_FINAL_ASSEMBLY);
3250:   MatAssemblyEnd(x,MAT_FINAL_ASSEMBLY);
3251:   return(0);
3252: }


3255: /* -------------------------------------------------------------------*/
3256: static struct _MatOps MatOps_Values = { MatSetValues_SeqAIJ,
3257:                                         MatGetRow_SeqAIJ,
3258:                                         MatRestoreRow_SeqAIJ,
3259:                                         MatMult_SeqAIJ,
3260:                                 /*  4*/ MatMultAdd_SeqAIJ,
3261:                                         MatMultTranspose_SeqAIJ,
3262:                                         MatMultTransposeAdd_SeqAIJ,
3263:                                         0,
3264:                                         0,
3265:                                         0,
3266:                                 /* 10*/ 0,
3267:                                         MatLUFactor_SeqAIJ,
3268:                                         0,
3269:                                         MatSOR_SeqAIJ,
3270:                                         MatTranspose_SeqAIJ,
3271:                                 /*1 5*/ MatGetInfo_SeqAIJ,
3272:                                         MatEqual_SeqAIJ,
3273:                                         MatGetDiagonal_SeqAIJ,
3274:                                         MatDiagonalScale_SeqAIJ,
3275:                                         MatNorm_SeqAIJ,
3276:                                 /* 20*/ 0,
3277:                                         MatAssemblyEnd_SeqAIJ,
3278:                                         MatSetOption_SeqAIJ,
3279:                                         MatZeroEntries_SeqAIJ,
3280:                                 /* 24*/ MatZeroRows_SeqAIJ,
3281:                                         0,
3282:                                         0,
3283:                                         0,
3284:                                         0,
3285:                                 /* 29*/ MatSetUp_SeqAIJ,
3286:                                         0,
3287:                                         0,
3288:                                         0,
3289:                                         0,
3290:                                 /* 34*/ MatDuplicate_SeqAIJ,
3291:                                         0,
3292:                                         0,
3293:                                         MatILUFactor_SeqAIJ,
3294:                                         0,
3295:                                 /* 39*/ MatAXPY_SeqAIJ,
3296:                                         MatCreateSubMatrices_SeqAIJ,
3297:                                         MatIncreaseOverlap_SeqAIJ,
3298:                                         MatGetValues_SeqAIJ,
3299:                                         MatCopy_SeqAIJ,
3300:                                 /* 44*/ MatGetRowMax_SeqAIJ,
3301:                                         MatScale_SeqAIJ,
3302:                                         MatShift_SeqAIJ,
3303:                                         MatDiagonalSet_SeqAIJ,
3304:                                         MatZeroRowsColumns_SeqAIJ,
3305:                                 /* 49*/ MatSetRandom_SeqAIJ,
3306:                                         MatGetRowIJ_SeqAIJ,
3307:                                         MatRestoreRowIJ_SeqAIJ,
3308:                                         MatGetColumnIJ_SeqAIJ,
3309:                                         MatRestoreColumnIJ_SeqAIJ,
3310:                                 /* 54*/ MatFDColoringCreate_SeqXAIJ,
3311:                                         0,
3312:                                         0,
3313:                                         MatPermute_SeqAIJ,
3314:                                         0,
3315:                                 /* 59*/ 0,
3316:                                         MatDestroy_SeqAIJ,
3317:                                         MatView_SeqAIJ,
3318:                                         0,
3319:                                         MatMatMatMult_SeqAIJ_SeqAIJ_SeqAIJ,
3320:                                 /* 64*/ MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ,
3321:                                         MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ,
3322:                                         0,
3323:                                         0,
3324:                                         0,
3325:                                 /* 69*/ MatGetRowMaxAbs_SeqAIJ,
3326:                                         MatGetRowMinAbs_SeqAIJ,
3327:                                         0,
3328:                                         0,
3329:                                         0,
3330:                                 /* 74*/ 0,
3331:                                         MatFDColoringApply_AIJ,
3332:                                         0,
3333:                                         0,
3334:                                         0,
3335:                                 /* 79*/ MatFindZeroDiagonals_SeqAIJ,
3336:                                         0,
3337:                                         0,
3338:                                         0,
3339:                                         MatLoad_SeqAIJ,
3340:                                 /* 84*/ MatIsSymmetric_SeqAIJ,
3341:                                         MatIsHermitian_SeqAIJ,
3342:                                         0,
3343:                                         0,
3344:                                         0,
3345:                                 /* 89*/ MatMatMult_SeqAIJ_SeqAIJ,
3346:                                         MatMatMultSymbolic_SeqAIJ_SeqAIJ,
3347:                                         MatMatMultNumeric_SeqAIJ_SeqAIJ,
3348:                                         MatPtAP_SeqAIJ_SeqAIJ,
3349:                                         MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy,
3350:                                 /* 94*/ MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy,
3351:                                         MatMatTransposeMult_SeqAIJ_SeqAIJ,
3352:                                         MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ,
3353:                                         MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ,
3354:                                         0,
3355:                                 /* 99*/ 0,
3356:                                         0,
3357:                                         0,
3358:                                         MatConjugate_SeqAIJ,
3359:                                         0,
3360:                                 /*104*/ MatSetValuesRow_SeqAIJ,
3361:                                         MatRealPart_SeqAIJ,
3362:                                         MatImaginaryPart_SeqAIJ,
3363:                                         0,
3364:                                         0,
3365:                                 /*109*/ MatMatSolve_SeqAIJ,
3366:                                         0,
3367:                                         MatGetRowMin_SeqAIJ,
3368:                                         0,
3369:                                         MatMissingDiagonal_SeqAIJ,
3370:                                 /*114*/ 0,
3371:                                         0,
3372:                                         0,
3373:                                         0,
3374:                                         0,
3375:                                 /*119*/ 0,
3376:                                         0,
3377:                                         0,
3378:                                         0,
3379:                                         MatGetMultiProcBlock_SeqAIJ,
3380:                                 /*124*/ MatFindNonzeroRows_SeqAIJ,
3381:                                         MatGetColumnNorms_SeqAIJ,
3382:                                         MatInvertBlockDiagonal_SeqAIJ,
3383:                                         MatInvertVariableBlockDiagonal_SeqAIJ,
3384:                                         0,
3385:                                 /*129*/ 0,
3386:                                         MatTransposeMatMult_SeqAIJ_SeqAIJ,
3387:                                         MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ,
3388:                                         MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ,
3389:                                         MatTransposeColoringCreate_SeqAIJ,
3390:                                 /*134*/ MatTransColoringApplySpToDen_SeqAIJ,
3391:                                         MatTransColoringApplyDenToSp_SeqAIJ,
3392:                                         MatRARt_SeqAIJ_SeqAIJ,
3393:                                         MatRARtSymbolic_SeqAIJ_SeqAIJ,
3394:                                         MatRARtNumeric_SeqAIJ_SeqAIJ,
3395:                                  /*139*/0,
3396:                                         0,
3397:                                         0,
3398:                                         MatFDColoringSetUp_SeqXAIJ,
3399:                                         MatFindOffBlockDiagonalEntries_SeqAIJ,
3400:                                  /*144*/MatCreateMPIMatConcatenateSeqMat_SeqAIJ,
3401:                                         MatDestroySubMatrices_SeqAIJ
3402: };

3404: PetscErrorCode  MatSeqAIJSetColumnIndices_SeqAIJ(Mat mat,PetscInt *indices)
3405: {
3406:   Mat_SeqAIJ *aij = (Mat_SeqAIJ*)mat->data;
3407:   PetscInt   i,nz,n;

3410:   nz = aij->maxnz;
3411:   n  = mat->rmap->n;
3412:   for (i=0; i<nz; i++) {
3413:     aij->j[i] = indices[i];
3414:   }
3415:   aij->nz = nz;
3416:   for (i=0; i<n; i++) {
3417:     aij->ilen[i] = aij->imax[i];
3418:   }
3419:   return(0);
3420: }

3422: /*
3423:  * When a sparse matrix has many zero columns, we should compact them out to save the space
3424:  * This happens in MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable()
3425:  * */
3426: PetscErrorCode  MatSeqAIJCompactOutExtraColumns_SeqAIJ(Mat mat, ISLocalToGlobalMapping *mapping)
3427: {
3428:   Mat_SeqAIJ *aij = (Mat_SeqAIJ*)mat->data;
3429:   PetscTable         gid1_lid1;
3430:   PetscTablePosition tpos;
3431:   PetscInt           gid,lid,i,j,ncols,ec;
3432:   PetscInt           *garray;
3433:   PetscErrorCode  ierr;

3438:   /* use a table */
3439:   PetscTableCreate(mat->rmap->n,mat->cmap->N+1,&gid1_lid1);
3440:   ec = 0;
3441:   for (i=0; i<mat->rmap->n; i++) {
3442:     ncols = aij->i[i+1] - aij->i[i];
3443:     for (j=0; j<ncols; j++) {
3444:       PetscInt data,gid1 = aij->j[aij->i[i] + j] + 1;
3445:       PetscTableFind(gid1_lid1,gid1,&data);
3446:       if (!data) {
3447:         /* one based table */
3448:         PetscTableAdd(gid1_lid1,gid1,++ec,INSERT_VALUES);
3449:       }
3450:     }
3451:   }
3452:   /* form array of columns we need */
3453:   PetscMalloc1(ec+1,&garray);
3454:   PetscTableGetHeadPosition(gid1_lid1,&tpos);
3455:   while (tpos) {
3456:     PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);
3457:     gid--;
3458:     lid--;
3459:     garray[lid] = gid;
3460:   }
3461:   PetscSortInt(ec,garray); /* sort, and rebuild */
3462:   PetscTableRemoveAll(gid1_lid1);
3463:   for (i=0; i<ec; i++) {
3464:     PetscTableAdd(gid1_lid1,garray[i]+1,i+1,INSERT_VALUES);
3465:   }
3466:   /* compact out the extra columns in B */
3467:   for (i=0; i<mat->rmap->n; i++) {
3468:         ncols = aij->i[i+1] - aij->i[i];
3469:     for (j=0; j<ncols; j++) {
3470:       PetscInt gid1 = aij->j[aij->i[i] + j] + 1;
3471:       PetscTableFind(gid1_lid1,gid1,&lid);
3472:       lid--;
3473:       aij->j[aij->i[i] + j] = lid;
3474:     }
3475:   }
3476:   mat->cmap->n = mat->cmap->N = ec;
3477:   mat->cmap->bs = 1;

3479:   PetscTableDestroy(&gid1_lid1);
3480:   PetscLayoutSetUp((mat->cmap));
3481:   ISLocalToGlobalMappingCreate(PETSC_COMM_SELF,mat->cmap->bs,mat->cmap->n,garray,PETSC_OWN_POINTER,mapping);
3482:   ISLocalToGlobalMappingSetType(*mapping,ISLOCALTOGLOBALMAPPINGHASH);
3483:   return(0);
3484: }

3486: /*@
3487:     MatSeqAIJSetColumnIndices - Set the column indices for all the rows
3488:        in the matrix.

3490:   Input Parameters:
3491: +  mat - the SeqAIJ matrix
3492: -  indices - the column indices

3494:   Level: advanced

3496:   Notes:
3497:     This can be called if you have precomputed the nonzero structure of the
3498:   matrix and want to provide it to the matrix object to improve the performance
3499:   of the MatSetValues() operation.

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

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

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

3508: @*/
3509: PetscErrorCode  MatSeqAIJSetColumnIndices(Mat mat,PetscInt *indices)
3510: {

3516:   PetscUseMethod(mat,"MatSeqAIJSetColumnIndices_C",(Mat,PetscInt*),(mat,indices));
3517:   return(0);
3518: }

3520: /* ----------------------------------------------------------------------------------------*/

3522: PetscErrorCode  MatStoreValues_SeqAIJ(Mat mat)
3523: {
3524:   Mat_SeqAIJ     *aij = (Mat_SeqAIJ*)mat->data;
3526:   size_t         nz = aij->i[mat->rmap->n];

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

3531:   /* allocate space for values if not already there */
3532:   if (!aij->saved_values) {
3533:     PetscMalloc1(nz+1,&aij->saved_values);
3534:     PetscLogObjectMemory((PetscObject)mat,(nz+1)*sizeof(PetscScalar));
3535:   }

3537:   /* copy values over */
3538:   PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));
3539:   return(0);
3540: }

3542: /*@
3543:     MatStoreValues - Stashes a copy of the matrix values; this allows, for
3544:        example, reuse of the linear part of a Jacobian, while recomputing the
3545:        nonlinear portion.

3547:    Collect on Mat

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

3552:   Level: advanced

3554:   Common Usage, with SNESSolve():
3555: $    Create Jacobian matrix
3556: $    Set linear terms into matrix
3557: $    Apply boundary conditions to matrix, at this time matrix must have
3558: $      final nonzero structure (i.e. setting the nonlinear terms and applying
3559: $      boundary conditions again will not change the nonzero structure
3560: $    MatSetOption(mat,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);
3561: $    MatStoreValues(mat);
3562: $    Call SNESSetJacobian() with matrix
3563: $    In your Jacobian routine
3564: $      MatRetrieveValues(mat);
3565: $      Set nonlinear terms in matrix

3567:   Common Usage without SNESSolve(), i.e. when you handle nonlinear solve yourself:
3568: $    // build linear portion of Jacobian
3569: $    MatSetOption(mat,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);
3570: $    MatStoreValues(mat);
3571: $    loop over nonlinear iterations
3572: $       MatRetrieveValues(mat);
3573: $       // call MatSetValues(mat,...) to set nonliner portion of Jacobian
3574: $       // call MatAssemblyBegin/End() on matrix
3575: $       Solve linear system with Jacobian
3576: $    endloop

3578:   Notes:
3579:     Matrix must already be assemblied before calling this routine
3580:     Must set the matrix option MatSetOption(mat,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE); before
3581:     calling this routine.

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

3586: .seealso: MatRetrieveValues()

3588: @*/
3589: PetscErrorCode  MatStoreValues(Mat mat)
3590: {

3595:   if (!mat->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
3596:   if (mat->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
3597:   PetscUseMethod(mat,"MatStoreValues_C",(Mat),(mat));
3598:   return(0);
3599: }

3601: PetscErrorCode  MatRetrieveValues_SeqAIJ(Mat mat)
3602: {
3603:   Mat_SeqAIJ     *aij = (Mat_SeqAIJ*)mat->data;
3605:   PetscInt       nz = aij->i[mat->rmap->n];

3608:   if (!aij->nonew) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first");
3609:   if (!aij->saved_values) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Must call MatStoreValues(A);first");
3610:   /* copy values over */
3611:   PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));
3612:   return(0);
3613: }

3615: /*@
3616:     MatRetrieveValues - Retrieves the copy of the matrix values; this allows, for
3617:        example, reuse of the linear part of a Jacobian, while recomputing the
3618:        nonlinear portion.

3620:    Collect on Mat

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

3625:   Level: advanced

3627: .seealso: MatStoreValues()

3629: @*/
3630: PetscErrorCode  MatRetrieveValues(Mat mat)
3631: {

3636:   if (!mat->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
3637:   if (mat->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
3638:   PetscUseMethod(mat,"MatRetrieveValues_C",(Mat),(mat));
3639:   return(0);
3640: }


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

3651:    Collective on MPI_Comm

3653:    Input Parameters:
3654: +  comm - MPI communicator, set to PETSC_COMM_SELF
3655: .  m - number of rows
3656: .  n - number of columns
3657: .  nz - number of nonzeros per row (same for all rows)
3658: -  nnz - array containing the number of nonzeros in the various rows
3659:          (possibly different for each row) or NULL

3661:    Output Parameter:
3662: .  A - the matrix

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

3668:    Notes:
3669:    If nnz is given then nz is ignored

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

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

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

3686:    Options Database Keys:
3687: +  -mat_no_inode  - Do not use inodes
3688: -  -mat_inode_limit <limit> - Sets inode limit (max limit=5)

3690:    Level: intermediate

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

3694: @*/
3695: PetscErrorCode  MatCreateSeqAIJ(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A)
3696: {

3700:   MatCreate(comm,A);
3701:   MatSetSizes(*A,m,n,m,n);
3702:   MatSetType(*A,MATSEQAIJ);
3703:   MatSeqAIJSetPreallocation_SeqAIJ(*A,nz,nnz);
3704:   return(0);
3705: }

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

3713:    Collective on MPI_Comm

3715:    Input Parameters:
3716: +  B - The matrix
3717: .  nz - number of nonzeros per row (same for all rows)
3718: -  nnz - array containing the number of nonzeros in the various rows
3719:          (possibly different for each row) or NULL

3721:    Notes:
3722:      If nnz is given then nz is ignored

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

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

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

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

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

3747:    Options Database Keys:
3748: +  -mat_no_inode  - Do not use inodes
3749: -  -mat_inode_limit <limit> - Sets inode limit (max limit=5)

3751:    Level: intermediate

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

3755: @*/
3756: PetscErrorCode  MatSeqAIJSetPreallocation(Mat B,PetscInt nz,const PetscInt nnz[])
3757: {

3763:   PetscTryMethod(B,"MatSeqAIJSetPreallocation_C",(Mat,PetscInt,const PetscInt[]),(B,nz,nnz));
3764:   return(0);
3765: }

3767: PetscErrorCode  MatSeqAIJSetPreallocation_SeqAIJ(Mat B,PetscInt nz,const PetscInt *nnz)
3768: {
3769:   Mat_SeqAIJ     *b;
3770:   PetscBool      skipallocation = PETSC_FALSE,realalloc = PETSC_FALSE;
3772:   PetscInt       i;

3775:   if (nz >= 0 || nnz) realalloc = PETSC_TRUE;
3776:   if (nz == MAT_SKIP_ALLOCATION) {
3777:     skipallocation = PETSC_TRUE;
3778:     nz             = 0;
3779:   }
3780:   PetscLayoutSetUp(B->rmap);
3781:   PetscLayoutSetUp(B->cmap);

3783:   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
3784:   if (nz < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %D",nz);
3785:   if (nnz) {
3786:     for (i=0; i<B->rmap->n; i++) {
3787:       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]);
3788:       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);
3789:     }
3790:   }

3792:   B->preallocated = PETSC_TRUE;

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

3796:   if (!skipallocation) {
3797:     if (!b->imax) {
3798:       PetscMalloc2(B->rmap->n,&b->imax,B->rmap->n,&b->ilen);
3799:       PetscLogObjectMemory((PetscObject)B,2*B->rmap->n*sizeof(PetscInt));
3800:     }
3801:     if (!b->ipre) {
3802:       PetscMalloc1(B->rmap->n,&b->ipre);
3803:       PetscLogObjectMemory((PetscObject)B,B->rmap->n*sizeof(PetscInt));
3804:     }
3805:     if (!nnz) {
3806:       if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 10;
3807:       else if (nz < 0) nz = 1;
3808:       nz = PetscMin(nz,B->cmap->n);
3809:       for (i=0; i<B->rmap->n; i++) b->imax[i] = nz;
3810:       nz = nz*B->rmap->n;
3811:     } else {
3812:       nz = 0;
3813:       for (i=0; i<B->rmap->n; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
3814:     }
3815:     /* b->ilen will count nonzeros in each row so far. */
3816:     for (i=0; i<B->rmap->n; i++) b->ilen[i] = 0;

3818:     /* allocate the matrix space */
3819:     /* FIXME: should B's old memory be unlogged? */
3820:     MatSeqXAIJFreeAIJ(B,&b->a,&b->j,&b->i);
3821:     if (B->structure_only) {
3822:       PetscMalloc1(nz,&b->j);
3823:       PetscMalloc1(B->rmap->n+1,&b->i);
3824:       PetscLogObjectMemory((PetscObject)B,(B->rmap->n+1)*sizeof(PetscInt)+nz*sizeof(PetscInt));
3825:     } else {
3826:       PetscMalloc3(nz,&b->a,nz,&b->j,B->rmap->n+1,&b->i);
3827:       PetscLogObjectMemory((PetscObject)B,(B->rmap->n+1)*sizeof(PetscInt)+nz*(sizeof(PetscScalar)+sizeof(PetscInt)));
3828:     }
3829:     b->i[0] = 0;
3830:     for (i=1; i<B->rmap->n+1; i++) {
3831:       b->i[i] = b->i[i-1] + b->imax[i-1];
3832:     }
3833:     if (B->structure_only) {
3834:       b->singlemalloc = PETSC_FALSE;
3835:       b->free_a       = PETSC_FALSE;
3836:     } else {
3837:       b->singlemalloc = PETSC_TRUE;
3838:       b->free_a       = PETSC_TRUE;
3839:     }
3840:     b->free_ij      = PETSC_TRUE;
3841:   } else {
3842:     b->free_a  = PETSC_FALSE;
3843:     b->free_ij = PETSC_FALSE;
3844:   }

3846:   if (b->ipre && nnz != b->ipre  && b->imax) {
3847:     /* reserve user-requested sparsity */
3848:     PetscMemcpy(b->ipre,b->imax,B->rmap->n*sizeof(PetscInt));
3849:   }


3852:   b->nz               = 0;
3853:   b->maxnz            = nz;
3854:   B->info.nz_unneeded = (double)b->maxnz;
3855:   if (realalloc) {
3856:     MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);
3857:   }
3858:   B->was_assembled = PETSC_FALSE;
3859:   B->assembled     = PETSC_FALSE;
3860:   return(0);
3861: }


3864: PetscErrorCode MatResetPreallocation_SeqAIJ(Mat A)
3865: {
3866:   Mat_SeqAIJ     *a;
3867:   PetscInt       i;

3872:   a = (Mat_SeqAIJ*)A->data;
3873:   /* if no saved info, we error out */
3874:   if (!a->ipre) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_NULL,"No saved preallocation info \n");

3876:   if (!a->i || !a->j || !a->a || !a->imax || !a->ilen) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_NULL,"Memory info is incomplete, and can not reset preallocation \n");

3878:   PetscMemcpy(a->imax,a->ipre,A->rmap->n*sizeof(PetscInt));
3879:   PetscMemzero(a->ilen,A->rmap->n*sizeof(PetscInt));
3880:   a->i[0] = 0;
3881:   for (i=1; i<A->rmap->n+1; i++) {
3882:     a->i[i] = a->i[i-1] + a->imax[i-1];
3883:   }
3884:   A->preallocated     = PETSC_TRUE;
3885:   a->nz               = 0;
3886:   a->maxnz            = a->i[A->rmap->n];
3887:   A->info.nz_unneeded = (double)a->maxnz;
3888:   A->was_assembled    = PETSC_FALSE;
3889:   A->assembled        = PETSC_FALSE;
3890:   return(0);
3891: }

3893: /*@
3894:    MatSeqAIJSetPreallocationCSR - Allocates memory for a sparse sequential matrix in AIJ format.

3896:    Input Parameters:
3897: +  B - the matrix
3898: .  i - the indices into j for the start of each row (starts with zero)
3899: .  j - the column indices for each row (starts with zero) these must be sorted for each row
3900: -  v - optional values in the matrix

3902:    Level: developer

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

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

3908: .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatSeqAIJSetPreallocation(), MatCreateSeqAIJ(), MATSEQAIJ
3909: @*/
3910: PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat B,const PetscInt i[],const PetscInt j[],const PetscScalar v[])
3911: {

3917:   PetscTryMethod(B,"MatSeqAIJSetPreallocationCSR_C",(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]),(B,i,j,v));
3918:   return(0);
3919: }

3921: PetscErrorCode  MatSeqAIJSetPreallocationCSR_SeqAIJ(Mat B,const PetscInt Ii[],const PetscInt J[],const PetscScalar v[])
3922: {
3923:   PetscInt       i;
3924:   PetscInt       m,n;
3925:   PetscInt       nz;
3926:   PetscInt       *nnz, nz_max = 0;
3927:   PetscScalar    *values;

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

3933:   PetscLayoutSetUp(B->rmap);
3934:   PetscLayoutSetUp(B->cmap);

3936:   MatGetSize(B, &m, &n);
3937:   PetscMalloc1(m+1, &nnz);
3938:   for (i = 0; i < m; i++) {
3939:     nz     = Ii[i+1]- Ii[i];
3940:     nz_max = PetscMax(nz_max, nz);
3941:     if (nz < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Local row %D has a negative number of columns %D", i, nnz);
3942:     nnz[i] = nz;
3943:   }
3944:   MatSeqAIJSetPreallocation(B, 0, nnz);
3945:   PetscFree(nnz);

3947:   if (v) {
3948:     values = (PetscScalar*) v;
3949:   } else {
3950:     PetscCalloc1(nz_max, &values);
3951:   }

3953:   for (i = 0; i < m; i++) {
3954:     nz   = Ii[i+1] - Ii[i];
3955:     MatSetValues_SeqAIJ(B, 1, &i, nz, J+Ii[i], values + (v ? Ii[i] : 0), INSERT_VALUES);
3956:   }

3958:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
3959:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

3961:   if (!v) {
3962:     PetscFree(values);
3963:   }
3964:   MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
3965:   return(0);
3966: }

3968:  #include <../src/mat/impls/dense/seq/dense.h>
3969:  #include <petsc/private/kernels/petscaxpy.h>

3971: /*
3972:     Computes (B'*A')' since computing B*A directly is untenable

3974:                n                       p                          p
3975:         (              )       (              )         (                  )
3976:       m (      A       )  *  n (       B      )   =   m (         C        )
3977:         (              )       (              )         (                  )

3979: */
3980: PetscErrorCode MatMatMultNumeric_SeqDense_SeqAIJ(Mat A,Mat B,Mat C)
3981: {
3982:   PetscErrorCode    ierr;
3983:   Mat_SeqDense      *sub_a = (Mat_SeqDense*)A->data;
3984:   Mat_SeqAIJ        *sub_b = (Mat_SeqAIJ*)B->data;
3985:   Mat_SeqDense      *sub_c = (Mat_SeqDense*)C->data;
3986:   PetscInt          i,n,m,q,p;
3987:   const PetscInt    *ii,*idx;
3988:   const PetscScalar *b,*a,*a_q;
3989:   PetscScalar       *c,*c_q;

3992:   m    = A->rmap->n;
3993:   n    = A->cmap->n;
3994:   p    = B->cmap->n;
3995:   a    = sub_a->v;
3996:   b    = sub_b->a;
3997:   c    = sub_c->v;
3998:   PetscMemzero(c,m*p*sizeof(PetscScalar));

4000:   ii  = sub_b->i;
4001:   idx = sub_b->j;
4002:   for (i=0; i<n; i++) {
4003:     q = ii[i+1] - ii[i];
4004:     while (q-->0) {
4005:       c_q = c + m*(*idx);
4006:       a_q = a + m*i;
4007:       PetscKernelAXPY(c_q,*b,a_q,m);
4008:       idx++;
4009:       b++;
4010:     }
4011:   }
4012:   return(0);
4013: }

4015: PetscErrorCode MatMatMultSymbolic_SeqDense_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
4016: {
4018:   PetscInt       m=A->rmap->n,n=B->cmap->n;
4019:   Mat            Cmat;

4022:   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);
4023:   MatCreate(PetscObjectComm((PetscObject)A),&Cmat);
4024:   MatSetSizes(Cmat,m,n,m,n);
4025:   MatSetBlockSizesFromMats(Cmat,A,B);
4026:   MatSetType(Cmat,MATSEQDENSE);
4027:   MatSeqDenseSetPreallocation(Cmat,NULL);

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

4031:   *C = Cmat;
4032:   return(0);
4033: }

4035: /* ----------------------------------------------------------------*/
4036: PETSC_INTERN PetscErrorCode MatMatMult_SeqDense_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
4037: {

4041:   if (scall == MAT_INITIAL_MATRIX) {
4042:     PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);
4043:     MatMatMultSymbolic_SeqDense_SeqAIJ(A,B,fill,C);
4044:     PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);
4045:   }
4046:   PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);
4047:   MatMatMultNumeric_SeqDense_SeqAIJ(A,B,*C);
4048:   PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);
4049:   return(0);
4050: }


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

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

4060:   Level: beginner

4062: .seealso: MatCreateSeqAIJ(), MatSetFromOptions(), MatSetType(), MatCreate(), MatType
4063: M*/

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

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

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

4077:   Developer Notes:
4078:     Subclasses include MATAIJCUSPARSE, MATAIJPERM, MATAIJSELL, MATAIJMKL, MATAIJCRL, and also automatically switches over to use inodes when
4079:    enough exist.

4081:   Level: beginner

4083: .seealso: MatCreateAIJ(), MatCreateSeqAIJ(), MATSEQAIJ,MATMPIAIJ
4084: M*/

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

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

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

4098:   Level: beginner

4100: .seealso: MatCreateMPIAIJCRL,MATSEQAIJCRL,MATMPIAIJCRL, MATSEQAIJCRL, MATMPIAIJCRL
4101: M*/

4103: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJCRL(Mat,MatType,MatReuse,Mat*);
4104: #if defined(PETSC_HAVE_ELEMENTAL)
4105: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_Elemental(Mat,MatType,MatReuse,Mat*);
4106: #endif
4107: #if defined(PETSC_HAVE_HYPRE)
4108: PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat A,MatType,MatReuse,Mat*);
4109: PETSC_INTERN PetscErrorCode MatMatMatMult_Transpose_AIJ_AIJ(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
4110: #endif
4111: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqDense(Mat,MatType,MatReuse,Mat*);

4113: PETSC_EXTERN PetscErrorCode MatConvert_SeqAIJ_SeqSELL(Mat,MatType,MatReuse,Mat*);
4114: PETSC_INTERN PetscErrorCode MatConvert_XAIJ_IS(Mat,MatType,MatReuse,Mat*);
4115: PETSC_INTERN PetscErrorCode MatPtAP_IS_XAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);

4117: /*@C
4118:    MatSeqAIJGetArray - gives access to the array where the data for a MATSEQAIJ matrix is stored

4120:    Not Collective

4122:    Input Parameter:
4123: .  mat - a MATSEQAIJ matrix

4125:    Output Parameter:
4126: .   array - pointer to the data

4128:    Level: intermediate

4130: .seealso: MatSeqAIJRestoreArray(), MatSeqAIJGetArrayF90()
4131: @*/
4132: PetscErrorCode  MatSeqAIJGetArray(Mat A,PetscScalar **array)
4133: {

4137:   PetscUseMethod(A,"MatSeqAIJGetArray_C",(Mat,PetscScalar**),(A,array));
4138:   return(0);
4139: }

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

4144:    Not Collective

4146:    Input Parameter:
4147: .  mat - a MATSEQAIJ matrix

4149:    Output Parameter:
4150: .   nz - the maximum number of nonzeros in any row

4152:    Level: intermediate

4154: .seealso: MatSeqAIJRestoreArray(), MatSeqAIJGetArrayF90()
4155: @*/
4156: PetscErrorCode  MatSeqAIJGetMaxRowNonzeros(Mat A,PetscInt *nz)
4157: {
4158:   Mat_SeqAIJ     *aij = (Mat_SeqAIJ*)A->data;

4161:   *nz = aij->rmax;
4162:   return(0);
4163: }

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

4168:    Not Collective

4170:    Input Parameters:
4171: .  mat - a MATSEQAIJ matrix
4172: .  array - pointer to the data

4174:    Level: intermediate

4176: .seealso: MatSeqAIJGetArray(), MatSeqAIJRestoreArrayF90()
4177: @*/
4178: PetscErrorCode  MatSeqAIJRestoreArray(Mat A,PetscScalar **array)
4179: {

4183:   PetscUseMethod(A,"MatSeqAIJRestoreArray_C",(Mat,PetscScalar**),(A,array));
4184:   return(0);
4185: }

4187: #if defined(PETSC_HAVE_CUDA)
4188: PETSC_EXTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJCUSPARSE(Mat);
4189: #endif

4191: PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJ(Mat B)
4192: {
4193:   Mat_SeqAIJ     *b;
4195:   PetscMPIInt    size;

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

4201:   PetscNewLog(B,&b);

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

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

4207:   b->row                = 0;
4208:   b->col                = 0;
4209:   b->icol               = 0;
4210:   b->reallocs           = 0;
4211:   b->ignorezeroentries  = PETSC_FALSE;
4212:   b->roworiented        = PETSC_TRUE;
4213:   b->nonew              = 0;
4214:   b->diag               = 0;
4215:   b->solve_work         = 0;
4216:   B->spptr              = 0;
4217:   b->saved_values       = 0;
4218:   b->idiag              = 0;
4219:   b->mdiag              = 0;
4220:   b->ssor_work          = 0;
4221:   b->omega              = 1.0;
4222:   b->fshift             = 0.0;
4223:   b->idiagvalid         = PETSC_FALSE;
4224:   b->ibdiagvalid        = PETSC_FALSE;
4225:   b->keepnonzeropattern = PETSC_FALSE;

4227:   PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);
4228:   PetscObjectComposeFunction((PetscObject)B,"MatSeqAIJGetArray_C",MatSeqAIJGetArray_SeqAIJ);
4229:   PetscObjectComposeFunction((PetscObject)B,"MatSeqAIJRestoreArray_C",MatSeqAIJRestoreArray_SeqAIJ);

4231: #if defined(PETSC_HAVE_MATLAB_ENGINE)
4232:   PetscObjectComposeFunction((PetscObject)B,"PetscMatlabEnginePut_C",MatlabEnginePut_SeqAIJ);
4233:   PetscObjectComposeFunction((PetscObject)B,"PetscMatlabEngineGet_C",MatlabEngineGet_SeqAIJ);
4234: #endif

4236:   PetscObjectComposeFunction((PetscObject)B,"MatSeqAIJSetColumnIndices_C",MatSeqAIJSetColumnIndices_SeqAIJ);
4237:   PetscObjectComposeFunction((PetscObject)B,"MatStoreValues_C",MatStoreValues_SeqAIJ);
4238:   PetscObjectComposeFunction((PetscObject)B,"MatRetrieveValues_C",MatRetrieveValues_SeqAIJ);
4239:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_seqsbaij_C",MatConvert_SeqAIJ_SeqSBAIJ);
4240:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_seqbaij_C",MatConvert_SeqAIJ_SeqBAIJ);
4241:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_seqaijperm_C",MatConvert_SeqAIJ_SeqAIJPERM);
4242:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_seqaijsell_C",MatConvert_SeqAIJ_SeqAIJSELL);
4243: #if defined(PETSC_HAVE_MKL_SPARSE)
4244:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_seqaijmkl_C",MatConvert_SeqAIJ_SeqAIJMKL);
4245: #endif
4246: #if defined(PETSC_HAVE_CUDA)
4247:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_seqaijcusparse_C",MatConvert_SeqAIJ_SeqAIJCUSPARSE);
4248: #endif
4249:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_seqaijcrl_C",MatConvert_SeqAIJ_SeqAIJCRL);
4250: #if defined(PETSC_HAVE_ELEMENTAL)
4251:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_elemental_C",MatConvert_SeqAIJ_Elemental);
4252: #endif
4253: #if defined(PETSC_HAVE_HYPRE)
4254:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_hypre_C",MatConvert_AIJ_HYPRE);
4255:   PetscObjectComposeFunction((PetscObject)B,"MatMatMatMult_transpose_seqaij_seqaij_C",MatMatMatMult_Transpose_AIJ_AIJ);
4256: #endif
4257:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_seqdense_C",MatConvert_SeqAIJ_SeqDense);
4258:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_seqsell_C",MatConvert_SeqAIJ_SeqSELL);
4259:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_is_C",MatConvert_XAIJ_IS);
4260:   PetscObjectComposeFunction((PetscObject)B,"MatIsTranspose_C",MatIsTranspose_SeqAIJ);
4261:   PetscObjectComposeFunction((PetscObject)B,"MatIsHermitianTranspose_C",MatIsTranspose_SeqAIJ);
4262:   PetscObjectComposeFunction((PetscObject)B,"MatSeqAIJSetPreallocation_C",MatSeqAIJSetPreallocation_SeqAIJ);
4263:   PetscObjectComposeFunction((PetscObject)B,"MatResetPreallocation_C",MatResetPreallocation_SeqAIJ);
4264:   PetscObjectComposeFunction((PetscObject)B,"MatSeqAIJSetPreallocationCSR_C",MatSeqAIJSetPreallocationCSR_SeqAIJ);
4265:   PetscObjectComposeFunction((PetscObject)B,"MatReorderForNonzeroDiagonal_C",MatReorderForNonzeroDiagonal_SeqAIJ);
4266:   PetscObjectComposeFunction((PetscObject)B,"MatMatMult_seqdense_seqaij_C",MatMatMult_SeqDense_SeqAIJ);
4267:   PetscObjectComposeFunction((PetscObject)B,"MatMatMultSymbolic_seqdense_seqaij_C",MatMatMultSymbolic_SeqDense_SeqAIJ);
4268:   PetscObjectComposeFunction((PetscObject)B,"MatMatMultNumeric_seqdense_seqaij_C",MatMatMultNumeric_SeqDense_SeqAIJ);
4269:   PetscObjectComposeFunction((PetscObject)B,"MatPtAP_is_seqaij_C",MatPtAP_IS_XAIJ);
4270:   MatCreate_SeqAIJ_Inode(B);
4271:   PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);
4272:   MatSeqAIJSetTypeFromOptions(B);  /* this allows changing the matrix subtype to say MATSEQAIJPERM */
4273:   return(0);
4274: }

4276: /*
4277:     Given a matrix generated with MatGetFactor() duplicates all the information in A into B
4278: */
4279: PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat C,Mat A,MatDuplicateOption cpvalues,PetscBool mallocmatspace)
4280: {
4281:   Mat_SeqAIJ     *c,*a = (Mat_SeqAIJ*)A->data;
4283:   PetscInt       i,m = A->rmap->n;

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

4288:   C->factortype = A->factortype;
4289:   c->row        = 0;
4290:   c->col        = 0;
4291:   c->icol       = 0;
4292:   c->reallocs   = 0;

4294:   C->assembled = PETSC_TRUE;

4296:   PetscLayoutReference(A->rmap,&C->rmap);
4297:   PetscLayoutReference(A->cmap,&C->cmap);

4299:   PetscMalloc2(m,&c->imax,m,&c->ilen);
4300:   PetscLogObjectMemory((PetscObject)C, 2*m*sizeof(PetscInt));
4301:   for (i=0; i<m; i++) {
4302:     c->imax[i] = a->imax[i];
4303:     c->ilen[i] = a->ilen[i];
4304:   }

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

4311:     c->singlemalloc = PETSC_TRUE;

4313:     PetscMemcpy(c->i,a->i,(m+1)*sizeof(PetscInt));
4314:     if (m > 0) {
4315:       PetscMemcpy(c->j,a->j,(a->i[m])*sizeof(PetscInt));
4316:       if (cpvalues == MAT_COPY_VALUES) {
4317:         PetscMemcpy(c->a,a->a,(a->i[m])*sizeof(PetscScalar));
4318:       } else {
4319:         PetscMemzero(c->a,(a->i[m])*sizeof(PetscScalar));
4320:       }
4321:     }
4322:   }

4324:   c->ignorezeroentries = a->ignorezeroentries;
4325:   c->roworiented       = a->roworiented;
4326:   c->nonew             = a->nonew;
4327:   if (a->diag) {
4328:     PetscMalloc1(m+1,&c->diag);
4329:     PetscLogObjectMemory((PetscObject)C,(m+1)*sizeof(PetscInt));
4330:     for (i=0; i<m; i++) {
4331:       c->diag[i] = a->diag[i];
4332:     }
4333:   } else c->diag = 0;

4335:   c->solve_work         = 0;
4336:   c->saved_values       = 0;
4337:   c->idiag              = 0;
4338:   c->ssor_work          = 0;
4339:   c->keepnonzeropattern = a->keepnonzeropattern;
4340:   c->free_a             = PETSC_TRUE;
4341:   c->free_ij            = PETSC_TRUE;

4343:   c->rmax         = a->rmax;
4344:   c->nz           = a->nz;
4345:   c->maxnz        = a->nz;       /* Since we allocate exactly the right amount */
4346:   C->preallocated = PETSC_TRUE;

4348:   c->compressedrow.use   = a->compressedrow.use;
4349:   c->compressedrow.nrows = a->compressedrow.nrows;
4350:   if (a->compressedrow.use) {
4351:     i    = a->compressedrow.nrows;
4352:     PetscMalloc2(i+1,&c->compressedrow.i,i,&c->compressedrow.rindex);
4353:     PetscMemcpy(c->compressedrow.i,a->compressedrow.i,(i+1)*sizeof(PetscInt));
4354:     PetscMemcpy(c->compressedrow.rindex,a->compressedrow.rindex,i*sizeof(PetscInt));
4355:   } else {
4356:     c->compressedrow.use    = PETSC_FALSE;
4357:     c->compressedrow.i      = NULL;
4358:     c->compressedrow.rindex = NULL;
4359:   }
4360:   c->nonzerorowcnt = a->nonzerorowcnt;
4361:   C->nonzerostate  = A->nonzerostate;

4363:   MatDuplicate_SeqAIJ_Inode(A,cpvalues,&C);
4364:   PetscFunctionListDuplicate(((PetscObject)A)->qlist,&((PetscObject)C)->qlist);
4365:   return(0);
4366: }

4368: PetscErrorCode MatDuplicate_SeqAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
4369: {

4373:   MatCreate(PetscObjectComm((PetscObject)A),B);
4374:   MatSetSizes(*B,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);
4375:   if (!(A->rmap->n % A->rmap->bs) && !(A->cmap->n % A->cmap->bs)) {
4376:     MatSetBlockSizesFromMats(*B,A,A);
4377:   }
4378:   MatSetType(*B,((PetscObject)A)->type_name);
4379:   MatDuplicateNoCreate_SeqAIJ(*B,A,cpvalues,PETSC_TRUE);
4380:   return(0);
4381: }

4383: PetscErrorCode MatLoad_SeqAIJ(Mat newMat, PetscViewer viewer)
4384: {
4385:   PetscBool      isbinary, ishdf5;

4391:   /* force binary viewer to load .info file if it has not yet done so */
4392:   PetscViewerSetUp(viewer);
4393:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
4394:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,  &ishdf5);
4395:   if (isbinary) {
4396:     MatLoad_SeqAIJ_Binary(newMat,viewer);
4397:   } else if (ishdf5) {
4398: #if defined(PETSC_HAVE_HDF5)
4399:     MatLoad_AIJ_HDF5(newMat,viewer);
4400: #else
4401:     SETERRQ(PetscObjectComm((PetscObject)newMat),PETSC_ERR_SUP,"HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
4402: #endif
4403:   } else {
4404:     SETERRQ2(PetscObjectComm((PetscObject)newMat),PETSC_ERR_SUP,"Viewer type %s not yet supported for reading %s matrices",((PetscObject)viewer)->type_name,((PetscObject)newMat)->type_name);
4405:   }
4406:   return(0);
4407: }

4409: PetscErrorCode MatLoad_SeqAIJ_Binary(Mat newMat, PetscViewer viewer)
4410: {
4411:   Mat_SeqAIJ     *a;
4413:   PetscInt       i,sum,nz,header[4],*rowlengths = 0,M,N,rows,cols;
4414:   int            fd;
4415:   PetscMPIInt    size;
4416:   MPI_Comm       comm;
4417:   PetscInt       bs = newMat->rmap->bs;

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

4424:   PetscOptionsBegin(comm,NULL,"Options for loading SEQAIJ matrix","Mat");
4425:   PetscOptionsInt("-matload_block_size","Set the blocksize used to store the matrix","MatLoad",bs,&bs,NULL);
4426:   PetscOptionsEnd();
4427:   if (bs < 0) bs = 1;
4428:   MatSetBlockSize(newMat,bs);

4430:   PetscViewerBinaryGetDescriptor(viewer,&fd);
4431:   PetscBinaryRead(fd,header,4,NULL,PETSC_INT);
4432:   if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"not matrix object in file");
4433:   M = header[1]; N = header[2]; nz = header[3];

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

4437:   /* read in row lengths */
4438:   PetscMalloc1(M,&rowlengths);
4439:   PetscBinaryRead(fd,rowlengths,M,NULL,PETSC_INT);

4441:   /* check if sum of rowlengths is same as nz */
4442:   for (i=0,sum=0; i< M; i++) sum +=rowlengths[i];
4443:   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);

4445:   /* set global size if not set already*/
4446:   if (newMat->rmap->n < 0 && newMat->rmap->N < 0 && newMat->cmap->n < 0 && newMat->cmap->N < 0) {
4447:     MatSetSizes(newMat,PETSC_DECIDE,PETSC_DECIDE,M,N);
4448:   } else {
4449:     /* if sizes and type are already set, check if the matrix  global sizes are correct */
4450:     MatGetSize(newMat,&rows,&cols);
4451:     if (rows < 0 && cols < 0) { /* user might provide local size instead of global size */
4452:       MatGetLocalSize(newMat,&rows,&cols);
4453:     }
4454:     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);
4455:   }
4456:   MatSeqAIJSetPreallocation_SeqAIJ(newMat,0,rowlengths);
4457:   a    = (Mat_SeqAIJ*)newMat->data;

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

4461:   /* read in nonzero values */
4462:   PetscBinaryRead(fd,a->a,nz,NULL,PETSC_SCALAR);

4464:   /* set matrix "i" values */
4465:   a->i[0] = 0;
4466:   for (i=1; i<= M; i++) {
4467:     a->i[i]      = a->i[i-1] + rowlengths[i-1];
4468:     a->ilen[i-1] = rowlengths[i-1];
4469:   }
4470:   PetscFree(rowlengths);

4472:   MatAssemblyBegin(newMat,MAT_FINAL_ASSEMBLY);
4473:   MatAssemblyEnd(newMat,MAT_FINAL_ASSEMBLY);
4474:   return(0);
4475: }

4477: PetscErrorCode MatEqual_SeqAIJ(Mat A,Mat B,PetscBool * flg)
4478: {
4479:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data,*b = (Mat_SeqAIJ*)B->data;
4481: #if defined(PETSC_USE_COMPLEX)
4482:   PetscInt k;
4483: #endif

4486:   /* If the  matrix dimensions are not equal,or no of nonzeros */
4487:   if ((A->rmap->n != B->rmap->n) || (A->cmap->n != B->cmap->n) ||(a->nz != b->nz)) {
4488:     *flg = PETSC_FALSE;
4489:     return(0);
4490:   }

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

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

4500:   /* if a->a are the same */
4501: #if defined(PETSC_USE_COMPLEX)
4502:   for (k=0; k<a->nz; k++) {
4503:     if (PetscRealPart(a->a[k]) != PetscRealPart(b->a[k]) || PetscImaginaryPart(a->a[k]) != PetscImaginaryPart(b->a[k])) {
4504:       *flg = PETSC_FALSE;
4505:       return(0);
4506:     }
4507:   }
4508: #else
4509:   PetscMemcmp(a->a,b->a,(a->nz)*sizeof(PetscScalar),flg);
4510: #endif
4511:   return(0);
4512: }

4514: /*@
4515:      MatCreateSeqAIJWithArrays - Creates an sequential AIJ matrix using matrix elements (in CSR format)
4516:               provided by the user.

4518:       Collective on MPI_Comm

4520:    Input Parameters:
4521: +   comm - must be an MPI communicator of size 1
4522: .   m - number of rows
4523: .   n - number of columns
4524: .   i - row indices; that is i[0] = 0, i[row] = i[row-1] + number of elements in that row of the matrix
4525: .   j - column indices
4526: -   a - matrix values

4528:    Output Parameter:
4529: .   mat - the matrix

4531:    Level: intermediate

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

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

4539:        The i and j indices are 0 based

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

4545: $        1 0 0
4546: $        2 0 3
4547: $        4 5 6
4548: $
4549: $        i =  {0,1,3,6}  [size = nrow+1  = 3+1]
4550: $        j =  {0,0,2,0,1,2}  [size = 6]; values must be sorted for each row
4551: $        v =  {1,2,3,4,5,6}  [size = 6]


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

4556: @*/
4557: PetscErrorCode  MatCreateSeqAIJWithArrays(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt i[],PetscInt j[],PetscScalar a[],Mat *mat)
4558: {
4560:   PetscInt       ii;
4561:   Mat_SeqAIJ     *aij;
4562: #if defined(PETSC_USE_DEBUG)
4563:   PetscInt jj;
4564: #endif

4567:   if (m > 0 && i[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"i (row indices) must start with 0");
4568:   MatCreate(comm,mat);
4569:   MatSetSizes(*mat,m,n,m,n);
4570:   /* MatSetBlockSizes(*mat,,); */
4571:   MatSetType(*mat,MATSEQAIJ);
4572:   MatSeqAIJSetPreallocation_SeqAIJ(*mat,MAT_SKIP_ALLOCATION,0);
4573:   aij  = (Mat_SeqAIJ*)(*mat)->data;
4574:   PetscMalloc2(m,&aij->imax,m,&aij->ilen);

4576:   aij->i            = i;
4577:   aij->j            = j;
4578:   aij->a            = a;
4579:   aij->singlemalloc = PETSC_FALSE;
4580:   aij->nonew        = -1;             /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/
4581:   aij->free_a       = PETSC_FALSE;
4582:   aij->free_ij      = PETSC_FALSE;

4584:   for (ii=0; ii<m; ii++) {
4585:     aij->ilen[ii] = aij->imax[ii] = i[ii+1] - i[ii];
4586: #if defined(PETSC_USE_DEBUG)
4587:     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]);
4588:     for (jj=i[ii]+1; jj<i[ii+1]; jj++) {
4589:       if (j[jj] < j[jj-1]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column entry number %D (actual column %D) in row %D is not sorted",jj-i[ii],j[jj],ii);
4590:       if (j[jj] == j[jj-1]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column entry number %D (actual column %D) in row %D is identical to previous entry",jj-i[ii],j[jj],ii);
4591:     }
4592: #endif
4593:   }
4594: #if defined(PETSC_USE_DEBUG)
4595:   for (ii=0; ii<aij->i[m]; ii++) {
4596:     if (j[ii] < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column index at location = %D index = %D",ii,j[ii]);
4597:     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]);
4598:   }
4599: #endif

4601:   MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
4602:   MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
4603:   return(0);
4604: }
4605: /*@C
4606:      MatCreateSeqAIJFromTriple - Creates an sequential AIJ matrix using matrix elements (in COO format)
4607:               provided by the user.

4609:       Collective on MPI_Comm

4611:    Input Parameters:
4612: +   comm - must be an MPI communicator of size 1
4613: .   m   - number of rows
4614: .   n   - number of columns
4615: .   i   - row indices
4616: .   j   - column indices
4617: .   a   - matrix values
4618: .   nz  - number of nonzeros
4619: -   idx - 0 or 1 based

4621:    Output Parameter:
4622: .   mat - the matrix

4624:    Level: intermediate

4626:    Notes:
4627:        The i and j indices are 0 based

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

4633:         1 0 0
4634:         2 0 3
4635:         4 5 6

4637:         i =  {0,1,1,2,2,2}
4638:         j =  {0,0,2,0,1,2}
4639:         v =  {1,2,3,4,5,6}


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

4644: @*/
4645: PetscErrorCode  MatCreateSeqAIJFromTriple(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt i[],PetscInt j[],PetscScalar a[],Mat *mat,PetscInt nz,PetscBool idx)
4646: {
4648:   PetscInt       ii, *nnz, one = 1,row,col;


4652:   PetscCalloc1(m,&nnz);
4653:   for (ii = 0; ii < nz; ii++) {
4654:     nnz[i[ii] - !!idx] += 1;
4655:   }
4656:   MatCreate(comm,mat);
4657:   MatSetSizes(*mat,m,n,m,n);
4658:   MatSetType(*mat,MATSEQAIJ);
4659:   MatSeqAIJSetPreallocation_SeqAIJ(*mat,0,nnz);
4660:   for (ii = 0; ii < nz; ii++) {
4661:     if (idx) {
4662:       row = i[ii] - 1;
4663:       col = j[ii] - 1;
4664:     } else {
4665:       row = i[ii];
4666:       col = j[ii];
4667:     }
4668:     MatSetValues(*mat,one,&row,one,&col,&a[ii],ADD_VALUES);
4669:   }
4670:   MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
4671:   MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
4672:   PetscFree(nnz);
4673:   return(0);
4674: }

4676: PetscErrorCode MatSeqAIJInvalidateDiagonal(Mat A)
4677: {
4678:   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;

4682:   a->idiagvalid  = PETSC_FALSE;
4683:   a->ibdiagvalid = PETSC_FALSE;

4685:   MatSeqAIJInvalidateDiagonal_Inode(A);
4686:   return(0);
4687: }

4689: PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqAIJ(MPI_Comm comm,Mat inmat,PetscInt n,MatReuse scall,Mat *outmat)
4690: {
4692:   PetscMPIInt    size;

4695:   MPI_Comm_size(comm,&size);
4696:   if (size == 1) {
4697:     if (scall == MAT_INITIAL_MATRIX) {
4698:       MatDuplicate(inmat,MAT_COPY_VALUES,outmat);
4699:     } else {
4700:       MatCopy(inmat,*outmat,SAME_NONZERO_PATTERN);
4701:     }
4702:   } else {
4703:     MatCreateMPIMatConcatenateSeqMat_MPIAIJ(comm,inmat,n,scall,outmat);
4704:   }
4705:   return(0);
4706: }

4708: /*
4709:  Permute A into C's *local* index space using rowemb,colemb.
4710:  The embedding are supposed to be injections and the above implies that the range of rowemb is a subset
4711:  of [0,m), colemb is in [0,n).
4712:  If pattern == DIFFERENT_NONZERO_PATTERN, C is preallocated according to A.
4713:  */
4714: PetscErrorCode MatSetSeqMat_SeqAIJ(Mat C,IS rowemb,IS colemb,MatStructure pattern,Mat B)
4715: {
4716:   /* If making this function public, change the error returned in this function away from _PLIB. */
4718:   Mat_SeqAIJ     *Baij;
4719:   PetscBool      seqaij;
4720:   PetscInt       m,n,*nz,i,j,count;
4721:   PetscScalar    v;
4722:   const PetscInt *rowindices,*colindices;

4725:   if (!B) return(0);
4726:   /* Check to make sure the target matrix (and embeddings) are compatible with C and each other. */
4727:   PetscObjectBaseTypeCompare((PetscObject)B,MATSEQAIJ,&seqaij);
4728:   if (!seqaij) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Input matrix is of wrong type");
4729:   if (rowemb) {
4730:     ISGetLocalSize(rowemb,&m);
4731:     if (m != B->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Row IS of size %D is incompatible with matrix row size %D",m,B->rmap->n);
4732:   } else {
4733:     if (C->rmap->n != B->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Input matrix is row-incompatible with the target matrix");
4734:   }
4735:   if (colemb) {
4736:     ISGetLocalSize(colemb,&n);
4737:     if (n != B->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Diag col IS of size %D is incompatible with input matrix col size %D",n,B->cmap->n);
4738:   } else {
4739:     if (C->cmap->n != B->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Input matrix is col-incompatible with the target matrix");
4740:   }

4742:   Baij = (Mat_SeqAIJ*)(B->data);
4743:   if (pattern == DIFFERENT_NONZERO_PATTERN) {
4744:     PetscMalloc1(B->rmap->n,&nz);
4745:     for (i=0; i<B->rmap->n; i++) {
4746:       nz[i] = Baij->i[i+1] - Baij->i[i];
4747:     }
4748:     MatSeqAIJSetPreallocation(C,0,nz);
4749:     PetscFree(nz);
4750:   }
4751:   if (pattern == SUBSET_NONZERO_PATTERN) {
4752:     MatZeroEntries(C);
4753:   }
4754:   count = 0;
4755:   rowindices = NULL;
4756:   colindices = NULL;
4757:   if (rowemb) {
4758:     ISGetIndices(rowemb,&rowindices);
4759:   }
4760:   if (colemb) {
4761:     ISGetIndices(colemb,&colindices);
4762:   }
4763:   for (i=0; i<B->rmap->n; i++) {
4764:     PetscInt row;
4765:     row = i;
4766:     if (rowindices) row = rowindices[i];
4767:     for (j=Baij->i[i]; j<Baij->i[i+1]; j++) {
4768:       PetscInt col;
4769:       col  = Baij->j[count];
4770:       if (colindices) col = colindices[col];
4771:       v    = Baij->a[count];
4772:       MatSetValues(C,1,&row,1,&col,&v,INSERT_VALUES);
4773:       ++count;
4774:     }
4775:   }
4776:   /* FIXME: set C's nonzerostate correctly. */
4777:   /* Assembly for C is necessary. */
4778:   C->preallocated = PETSC_TRUE;
4779:   C->assembled     = PETSC_TRUE;
4780:   C->was_assembled = PETSC_FALSE;
4781:   return(0);
4782: }

4784: PetscFunctionList MatSeqAIJList = NULL;

4786: /*@C
4787:    MatSeqAIJSetType - Converts a MATSEQAIJ matrix to a subtype

4789:    Collective on Mat

4791:    Input Parameters:
4792: +  mat      - the matrix object
4793: -  matype   - matrix type

4795:    Options Database Key:
4796: .  -mat_seqai_type  <method> - for example seqaijcrl


4799:   Level: intermediate

4801: .keywords: Mat, MatType, set, method

4803: .seealso: PCSetType(), VecSetType(), MatCreate(), MatType, Mat
4804: @*/
4805: PetscErrorCode  MatSeqAIJSetType(Mat mat, MatType matype)
4806: {
4807:   PetscErrorCode ierr,(*r)(Mat,MatType,MatReuse,Mat*);
4808:   PetscBool      sametype;

4812:   PetscObjectTypeCompare((PetscObject)mat,matype,&sametype);
4813:   if (sametype) return(0);

4815:    PetscFunctionListFind(MatSeqAIJList,matype,&r);
4816:   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown Mat type given: %s",matype);
4817:   (*r)(mat,matype,MAT_INPLACE_MATRIX,&mat);
4818:   return(0);
4819: }


4822: /*@C
4823:   MatSeqAIJRegister -  - Adds a new sub-matrix type for sequential AIJ matrices

4825:    Not Collective

4827:    Input Parameters:
4828: +  name - name of a new user-defined matrix type, for example MATSEQAIJCRL
4829: -  function - routine to convert to subtype

4831:    Notes:
4832:    MatSeqAIJRegister() may be called multiple times to add several user-defined solvers.


4835:    Then, your matrix can be chosen with the procedural interface at runtime via the option
4836: $     -mat_seqaij_type my_mat

4838:    Level: advanced

4840: .keywords: Mat, register

4842: .seealso: MatSeqAIJRegisterAll()


4845:   Level: advanced
4846: @*/
4847: PetscErrorCode  MatSeqAIJRegister(const char sname[],PetscErrorCode (*function)(Mat,MatType,MatReuse,Mat *))
4848: {

4852:   MatInitializePackage();
4853:   PetscFunctionListAdd(&MatSeqAIJList,sname,function);
4854:   return(0);
4855: }

4857: PetscBool MatSeqAIJRegisterAllCalled = PETSC_FALSE;

4859: /*@C
4860:   MatSeqAIJRegisterAll - Registers all of the matrix subtypes of SeqAIJ

4862:   Not Collective

4864:   Level: advanced

4866:   Developers Note: CUSP and CUSPARSE do not yet support the  MatConvert_SeqAIJ..() paradigm and thus cannot be registered here

4868: .keywords: KSP, register, all

4870: .seealso:  MatRegisterAll(), MatSeqAIJRegister()
4871: @*/
4872: PetscErrorCode  MatSeqAIJRegisterAll(void)
4873: {

4877:   if (MatSeqAIJRegisterAllCalled) return(0);
4878:   MatSeqAIJRegisterAllCalled = PETSC_TRUE;

4880:   MatSeqAIJRegister(MATSEQAIJCRL,      MatConvert_SeqAIJ_SeqAIJCRL);
4881:   MatSeqAIJRegister(MATSEQAIJPERM,     MatConvert_SeqAIJ_SeqAIJPERM);
4882:   MatSeqAIJRegister(MATSEQAIJSELL,     MatConvert_SeqAIJ_SeqAIJSELL);
4883: #if defined(PETSC_HAVE_MKL_SPARSE)
4884:   MatSeqAIJRegister(MATSEQAIJMKL,      MatConvert_SeqAIJ_SeqAIJMKL);
4885: #endif
4886: #if defined(PETSC_HAVE_VIENNACL) && defined(PETSC_HAVE_VIENNACL_NO_CUDA)
4887:   MatSeqAIJRegister(MATMPIAIJVIENNACL, MatConvert_SeqAIJ_SeqAIJViennaCL);
4888: #endif
4889:   return(0);
4890: }

4892: /*
4893:     Special version for direct calls from Fortran
4894: */
4895:  #include <petsc/private/fortranimpl.h>
4896: #if defined(PETSC_HAVE_FORTRAN_CAPS)
4897: #define matsetvaluesseqaij_ MATSETVALUESSEQAIJ
4898: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
4899: #define matsetvaluesseqaij_ matsetvaluesseqaij
4900: #endif

4902: /* Change these macros so can be used in void function */
4903: #undef CHKERRQ
4904: #define CHKERRQ(ierr) CHKERRABORT(PetscObjectComm((PetscObject)A),ierr)
4905: #undef SETERRQ2
4906: #define SETERRQ2(comm,ierr,b,c,d) CHKERRABORT(comm,ierr)
4907: #undef SETERRQ3
4908: #define SETERRQ3(comm,ierr,b,c,d,e) CHKERRABORT(comm,ierr)

4910: 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)
4911: {
4912:   Mat            A  = *AA;
4913:   PetscInt       m  = *mm, n = *nn;
4914:   InsertMode     is = *isis;
4915:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
4916:   PetscInt       *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N;
4917:   PetscInt       *imax,*ai,*ailen;
4919:   PetscInt       *aj,nonew = a->nonew,lastcol = -1;
4920:   MatScalar      *ap,value,*aa;
4921:   PetscBool      ignorezeroentries = a->ignorezeroentries;
4922:   PetscBool      roworiented       = a->roworiented;

4925:   MatCheckPreallocated(A,1);
4926:   imax  = a->imax;
4927:   ai    = a->i;
4928:   ailen = a->ilen;
4929:   aj    = a->j;
4930:   aa    = a->a;

4932:   for (k=0; k<m; k++) { /* loop over added rows */
4933:     row = im[k];
4934:     if (row < 0) continue;
4935: #if defined(PETSC_USE_DEBUG)
4936:     if (row >= A->rmap->n) SETERRABORT(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
4937: #endif
4938:     rp   = aj + ai[row]; ap = aa + ai[row];
4939:     rmax = imax[row]; nrow = ailen[row];
4940:     low  = 0;
4941:     high = nrow;
4942:     for (l=0; l<n; l++) { /* loop over added columns */
4943:       if (in[l] < 0) continue;
4944: #if defined(PETSC_USE_DEBUG)
4945:       if (in[l] >= A->cmap->n) SETERRABORT(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
4946: #endif
4947:       col = in[l];
4948:       if (roworiented) value = v[l + k*n];
4949:       else value = v[k + l*m];

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

4953:       if (col <= lastcol) low = 0;
4954:       else high = nrow;
4955:       lastcol = col;
4956:       while (high-low > 5) {
4957:         t = (low+high)/2;
4958:         if (rp[t] > col) high = t;
4959:         else             low  = t;
4960:       }
4961:       for (i=low; i<high; i++) {
4962:         if (rp[i] > col) break;
4963:         if (rp[i] == col) {
4964:           if (is == ADD_VALUES) ap[i] += value;
4965:           else                  ap[i] = value;
4966:           goto noinsert;
4967:         }
4968:       }
4969:       if (value == 0.0 && ignorezeroentries) goto noinsert;
4970:       if (nonew == 1) goto noinsert;
4971:       if (nonew == -1) SETERRABORT(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
4972:       MatSeqXAIJReallocateAIJ(A,A->rmap->n,1,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar);
4973:       N = nrow++ - 1; a->nz++; high++;
4974:       /* shift up all the later entries in this row */
4975:       for (ii=N; ii>=i; ii--) {
4976:         rp[ii+1] = rp[ii];
4977:         ap[ii+1] = ap[ii];
4978:       }
4979:       rp[i] = col;
4980:       ap[i] = value;
4981:       A->nonzerostate++;
4982: noinsert:;
4983:       low = i + 1;
4984:     }
4985:     ailen[row] = nrow;
4986:   }
4987:   PetscFunctionReturnVoid();
4988: }