Actual source code: aij.c

petsc-3.12.5 2020-03-29
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:   PetscArrayzero(norms,n);
 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;
177: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
178:   PetscBool         inserted = PETSC_FALSE;
179: #endif

182:   if (Y->assembled) {
183:     MatMissingDiagonal_SeqAIJ(Y,&missing,NULL);
184:     if (!missing) {
185:       diag = aij->diag;
186:       VecGetArrayRead(D,&v);
187:       if (is == INSERT_VALUES) {
188: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
189:         inserted = PETSC_TRUE;
190: #endif
191:         for (i=0; i<m; i++) {
192:           aa[diag[i]] = v[i];
193:         }
194:       } else {
195:         for (i=0; i<m; i++) {
196: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
197:           if (v[i] != 0.0) inserted = PETSC_TRUE;
198: #endif
199:           aa[diag[i]] += v[i];
200:         }
201:       }
202: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
203:       if (inserted) Y->offloadmask = PETSC_OFFLOAD_CPU;
204: #endif
205:       VecRestoreArrayRead(D,&v);
206:       return(0);
207:     }
208:     MatSeqAIJInvalidateDiagonal(Y);
209:   }
210:   MatDiagonalSet_Default(Y,D,is);
211:   return(0);
212: }

214: PetscErrorCode MatGetRowIJ_SeqAIJ(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *m,const PetscInt *ia[],const PetscInt *ja[],PetscBool  *done)
215: {
216:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
218:   PetscInt       i,ishift;

221:   *m = A->rmap->n;
222:   if (!ia) return(0);
223:   ishift = 0;
224:   if (symmetric && !A->structurally_symmetric) {
225:     MatToSymmetricIJ_SeqAIJ(A->rmap->n,a->i,a->j,PETSC_TRUE,ishift,oshift,(PetscInt**)ia,(PetscInt**)ja);
226:   } else if (oshift == 1) {
227:     PetscInt *tia;
228:     PetscInt nz = a->i[A->rmap->n];
229:     /* malloc space and  add 1 to i and j indices */
230:     PetscMalloc1(A->rmap->n+1,&tia);
231:     for (i=0; i<A->rmap->n+1; i++) tia[i] = a->i[i] + 1;
232:     *ia = tia;
233:     if (ja) {
234:       PetscInt *tja;
235:       PetscMalloc1(nz+1,&tja);
236:       for (i=0; i<nz; i++) tja[i] = a->j[i] + 1;
237:       *ja = tja;
238:     }
239:   } else {
240:     *ia = a->i;
241:     if (ja) *ja = a->j;
242:   }
243:   return(0);
244: }

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

251:   if (!ia) return(0);
252:   if ((symmetric && !A->structurally_symmetric) || oshift == 1) {
253:     PetscFree(*ia);
254:     if (ja) {PetscFree(*ja);}
255:   }
256:   return(0);
257: }

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

267:   *nn = n;
268:   if (!ia) return(0);
269:   if (symmetric) {
270:     MatToSymmetricIJ_SeqAIJ(A->rmap->n,a->i,a->j,PETSC_TRUE,0,oshift,(PetscInt**)ia,(PetscInt**)ja);
271:   } else {
272:     PetscCalloc1(n,&collengths);
273:     PetscMalloc1(n+1,&cia);
274:     PetscMalloc1(nz,&cja);
275:     jj   = a->j;
276:     for (i=0; i<nz; i++) {
277:       collengths[jj[i]]++;
278:     }
279:     cia[0] = oshift;
280:     for (i=0; i<n; i++) {
281:       cia[i+1] = cia[i] + collengths[i];
282:     }
283:     PetscArrayzero(collengths,n);
284:     jj   = a->j;
285:     for (row=0; row<m; row++) {
286:       mr = a->i[row+1] - a->i[row];
287:       for (i=0; i<mr; i++) {
288:         col = *jj++;

290:         cja[cia[col] + collengths[col]++ - oshift] = row + oshift;
291:       }
292:     }
293:     PetscFree(collengths);
294:     *ia  = cia; *ja = cja;
295:   }
296:   return(0);
297: }

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

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

306:   PetscFree(*ia);
307:   PetscFree(*ja);
308:   return(0);
309: }

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

326:   *nn = n;
327:   if (!ia) return(0);

329:   PetscCalloc1(n,&collengths);
330:   PetscMalloc1(n+1,&cia);
331:   PetscMalloc1(nz,&cja);
332:   PetscMalloc1(nz,&cspidx);
333:   jj   = a->j;
334:   for (i=0; i<nz; i++) {
335:     collengths[jj[i]]++;
336:   }
337:   cia[0] = oshift;
338:   for (i=0; i<n; i++) {
339:     cia[i+1] = cia[i] + collengths[i];
340:   }
341:   PetscArrayzero(collengths,n);
342:   jj   = a->j;
343:   for (row=0; row<m; row++) {
344:     mr = a->i[row+1] - a->i[row];
345:     for (i=0; i<mr; i++) {
346:       col         = *jj++;
347:       tmp         = cia[col] + collengths[col]++ - oshift;
348:       cspidx[tmp] = a->i[row] + i; /* index of a->j */
349:       cja[tmp]    = row + oshift;
350:     }
351:   }
352:   PetscFree(collengths);
353:   *ia    = cia;
354:   *ja    = cja;
355:   *spidx = cspidx;
356:   return(0);
357: }

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

364:   MatRestoreColumnIJ_SeqAIJ(A,oshift,symmetric,inodecompressed,n,ia,ja,done);
365:   PetscFree(*spidx);
366:   return(0);
367: }

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

376:   PetscArraycpy(a->a+ai[row],v,ai[row+1]-ai[row]);
377: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
378:   if (A->offloadmask != PETSC_OFFLOAD_UNALLOCATED && ai[row+1]-ai[row]) A->offloadmask = PETSC_OFFLOAD_CPU;
379: #endif
380:   return(0);
381: }

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

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

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

393: */

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

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

415:     if (col <= lastcol) low = 0;
416:     else high = nrow;
417:     lastcol = col;
418:     while (high-low > 5) {
419:       t = (low+high)/2;
420:       if (rp[t] > col) high = t;
421:       else low = t;
422:     }
423:     for (i=low; i<high; i++) {
424:       if (rp[i] == col) {
425:         ap[i] += value;
426:         low = i + 1;
427:         break;
428:       }
429:     }
430:   }
431: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
432:   if (A->offloadmask != PETSC_OFFLOAD_UNALLOCATED && m && n) A->offloadmask = PETSC_OFFLOAD_CPU;
433: #endif
434:   return 0;
435: }

437: PetscErrorCode MatSetValues_SeqAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
438: {
439:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
440:   PetscInt       *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N;
441:   PetscInt       *imax = a->imax,*ai = a->i,*ailen = a->ilen;
443:   PetscInt       *aj = a->j,nonew = a->nonew,lastcol = -1;
444:   MatScalar      *ap=NULL,value=0.0,*aa = a->a;
445:   PetscBool      ignorezeroentries = a->ignorezeroentries;
446:   PetscBool      roworiented       = a->roworiented;
447: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
448:   PetscBool      inserted          = PETSC_FALSE;
449: #endif

452:   for (k=0; k<m; k++) { /* loop over added rows */
453:     row = im[k];
454:     if (row < 0) continue;
455: #if defined(PETSC_USE_DEBUG)
456:     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);
457: #endif
458:     rp   = aj + ai[row];
459:     if (!A->structure_only) ap = aa + ai[row];
460:     rmax = imax[row]; nrow = ailen[row];
461:     low  = 0;
462:     high = nrow;
463:     for (l=0; l<n; l++) { /* loop over added columns */
464:       if (in[l] < 0) continue;
465: #if defined(PETSC_USE_DEBUG)
466:       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);
467: #endif
468:       col = in[l];
469:       if (v && !A->structure_only) value = roworiented ? v[l + k*n] : v[k + l*m];
470:       if (!A->structure_only && value == 0.0 && ignorezeroentries && is == ADD_VALUES && row != col) continue;

472:       if (col <= lastcol) low = 0;
473:       else high = nrow;
474:       lastcol = col;
475:       while (high-low > 5) {
476:         t = (low+high)/2;
477:         if (rp[t] > col) high = t;
478:         else low = t;
479:       }
480:       for (i=low; i<high; i++) {
481:         if (rp[i] > col) break;
482:         if (rp[i] == col) {
483:           if (!A->structure_only) {
484:             if (is == ADD_VALUES) {
485:               ap[i] += value;
486:               (void)PetscLogFlops(1.0);
487:             }
488:             else ap[i] = value;
489: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
490:             inserted = PETSC_TRUE;
491: #endif
492:           }
493:           low = i + 1;
494:           goto noinsert;
495:         }
496:       }
497:       if (value == 0.0 && ignorezeroentries && row != col) goto noinsert;
498:       if (nonew == 1) goto noinsert;
499:       if (nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero at (%D,%D) in the matrix",row,col);
500:       if (A->structure_only) {
501:         MatSeqXAIJReallocateAIJ_structure_only(A,A->rmap->n,1,nrow,row,col,rmax,ai,aj,rp,imax,nonew,MatScalar);
502:       } else {
503:         MatSeqXAIJReallocateAIJ(A,A->rmap->n,1,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar);
504:       }
505:       N = nrow++ - 1; a->nz++; high++;
506:       /* shift up all the later entries in this row */
507:       PetscArraymove(rp+i+1,rp+i,N-i+1);
508:       rp[i] = col;
509:       if (!A->structure_only){
510:         PetscArraymove(ap+i+1,ap+i,N-i+1);
511:         ap[i] = value;
512:       }
513:       low = i + 1;
514:       A->nonzerostate++;
515: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
516:       inserted = PETSC_TRUE;
517: #endif
518: noinsert:;
519:     }
520:     ailen[row] = nrow;
521:   }
522: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
523:   if (A->offloadmask != PETSC_OFFLOAD_UNALLOCATED && inserted) A->offloadmask = PETSC_OFFLOAD_CPU;
524: #endif
525:   return(0);
526: }

528: PetscErrorCode MatSetValues_SeqAIJ_SortedFull(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
529: {
530:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
531:   PetscInt       *rp,k,row;
532:   PetscInt       *ai = a->i,*ailen = a->ilen;
534:   PetscInt       *aj = a->j;
535:   MatScalar      *aa = a->a,*ap;

538:   for (k=0; k<m; k++) { /* loop over added rows */
539:     row  = im[k];
540:     rp   = aj + ai[row];
541:     ap   = aa + ai[row];
542:     if (!A->was_assembled) {
543:       PetscMemcpy(rp,in,n*sizeof(PetscInt));
544:     }
545:     if (!A->structure_only) {
546:       if (v) {
547:         PetscMemcpy(ap,v,n*sizeof(PetscScalar));
548:         v   += n;
549:       } else {
550:         PetscMemzero(ap,n*sizeof(PetscScalar));
551:       }
552:     }
553:     ailen[row] = n;
554:     a->nz      += n;
555:   }
556: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
557:   if (A->offloadmask != PETSC_OFFLOAD_UNALLOCATED && m && n) A->offloadmask = PETSC_OFFLOAD_CPU;
558: #endif
559:   return(0);
560: }


563: PetscErrorCode MatGetValues_SeqAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],PetscScalar v[])
564: {
565:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
566:   PetscInt   *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
567:   PetscInt   *ai = a->i,*ailen = a->ilen;
568:   MatScalar  *ap,*aa = a->a;

571:   for (k=0; k<m; k++) { /* loop over rows */
572:     row = im[k];
573:     if (row < 0) {v += n; continue;} /* SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",row); */
574:     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);
575:     rp   = aj + ai[row]; ap = aa + ai[row];
576:     nrow = ailen[row];
577:     for (l=0; l<n; l++) { /* loop over columns */
578:       if (in[l] < 0) {v++; continue;} /* SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %D",in[l]); */
579:       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);
580:       col  = in[l];
581:       high = nrow; low = 0; /* assume unsorted */
582:       while (high-low > 5) {
583:         t = (low+high)/2;
584:         if (rp[t] > col) high = t;
585:         else low = t;
586:       }
587:       for (i=low; i<high; i++) {
588:         if (rp[i] > col) break;
589:         if (rp[i] == col) {
590:           *v++ = ap[i];
591:           goto finished;
592:         }
593:       }
594:       *v++ = 0.0;
595: finished:;
596:     }
597:   }
598:   return(0);
599: }


602: PetscErrorCode MatView_SeqAIJ_Binary(Mat A,PetscViewer viewer)
603: {
604:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
606:   PetscInt       i,*col_lens;
607:   int            fd;
608:   FILE           *file;

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

614:   col_lens[0] = MAT_FILE_CLASSID;
615:   col_lens[1] = A->rmap->n;
616:   col_lens[2] = A->cmap->n;
617:   col_lens[3] = a->nz;

619:   /* store lengths of each row and write (including header) to file */
620:   for (i=0; i<A->rmap->n; i++) {
621:     col_lens[4+i] = a->i[i+1] - a->i[i];
622:   }
623:   PetscBinaryWrite(fd,col_lens,4+A->rmap->n,PETSC_INT,PETSC_TRUE);
624:   PetscFree(col_lens);

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

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

632:   PetscViewerBinaryGetInfoPointer(viewer,&file);
633:   if (file) {
634:     fprintf(file,"-matload_block_size %d\n",(int)PetscAbs(A->rmap->bs));
635:   }
636:   return(0);
637: }

639: static PetscErrorCode MatView_SeqAIJ_ASCII_structonly(Mat A,PetscViewer viewer)
640: {
642:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
643:   PetscInt       i,k,m=A->rmap->N;

646:   PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
647:   for (i=0; i<m; i++) {
648:     PetscViewerASCIIPrintf(viewer,"row %D:",i);
649:     for (k=a->i[i]; k<a->i[i+1]; k++) {
650:       PetscViewerASCIIPrintf(viewer," (%D) ",a->j[k]);
651:     }
652:     PetscViewerASCIIPrintf(viewer,"\n");
653:   }
654:   PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
655:   return(0);
656: }

658: extern PetscErrorCode MatSeqAIJFactorInfo_Matlab(Mat,PetscViewer);

660: PetscErrorCode MatView_SeqAIJ_ASCII(Mat A,PetscViewer viewer)
661: {
662:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
663:   PetscErrorCode    ierr;
664:   PetscInt          i,j,m = A->rmap->n;
665:   const char        *name;
666:   PetscViewerFormat format;

669:   if (A->structure_only) {
670:     MatView_SeqAIJ_ASCII_structonly(A,viewer);
671:     return(0);
672:   }

674:   PetscViewerGetFormat(viewer,&format);
675:   if (format == PETSC_VIEWER_ASCII_MATLAB) {
676:     PetscInt nofinalvalue = 0;
677:     if (m && ((a->i[m] == a->i[m-1]) || (a->j[a->nz-1] != A->cmap->n-1))) {
678:       /* Need a dummy value to ensure the dimension of the matrix. */
679:       nofinalvalue = 1;
680:     }
681:     PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
682:     PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",m,A->cmap->n);
683:     PetscViewerASCIIPrintf(viewer,"%% Nonzeros = %D \n",a->nz);
684: #if defined(PETSC_USE_COMPLEX)
685:     PetscViewerASCIIPrintf(viewer,"zzz = zeros(%D,4);\n",a->nz+nofinalvalue);
686: #else
687:     PetscViewerASCIIPrintf(viewer,"zzz = zeros(%D,3);\n",a->nz+nofinalvalue);
688: #endif
689:     PetscViewerASCIIPrintf(viewer,"zzz = [\n");

691:     for (i=0; i<m; i++) {
692:       for (j=a->i[i]; j<a->i[i+1]; j++) {
693: #if defined(PETSC_USE_COMPLEX)
694:         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]));
695: #else
696:         PetscViewerASCIIPrintf(viewer,"%D %D  %18.16e\n",i+1,a->j[j]+1,(double)a->a[j]);
697: #endif
698:       }
699:     }
700:     if (nofinalvalue) {
701: #if defined(PETSC_USE_COMPLEX)
702:       PetscViewerASCIIPrintf(viewer,"%D %D  %18.16e %18.16e\n",m,A->cmap->n,0.,0.);
703: #else
704:       PetscViewerASCIIPrintf(viewer,"%D %D  %18.16e\n",m,A->cmap->n,0.0);
705: #endif
706:     }
707:     PetscObjectGetName((PetscObject)A,&name);
708:     PetscViewerASCIIPrintf(viewer,"];\n %s = spconvert(zzz);\n",name);
709:     PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
710:   } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
711:     return(0);
712:   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
713:     PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
714:     for (i=0; i<m; i++) {
715:       PetscViewerASCIIPrintf(viewer,"row %D:",i);
716:       for (j=a->i[i]; j<a->i[i+1]; j++) {
717: #if defined(PETSC_USE_COMPLEX)
718:         if (PetscImaginaryPart(a->a[j]) > 0.0 && PetscRealPart(a->a[j]) != 0.0) {
719:           PetscViewerASCIIPrintf(viewer," (%D, %g + %g i)",a->j[j],(double)PetscRealPart(a->a[j]),(double)PetscImaginaryPart(a->a[j]));
720:         } else if (PetscImaginaryPart(a->a[j]) < 0.0 && PetscRealPart(a->a[j]) != 0.0) {
721:           PetscViewerASCIIPrintf(viewer," (%D, %g - %g i)",a->j[j],(double)PetscRealPart(a->a[j]),(double)-PetscImaginaryPart(a->a[j]));
722:         } else if (PetscRealPart(a->a[j]) != 0.0) {
723:           PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[j],(double)PetscRealPart(a->a[j]));
724:         }
725: #else
726:         if (a->a[j] != 0.0) {PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[j],(double)a->a[j]);}
727: #endif
728:       }
729:       PetscViewerASCIIPrintf(viewer,"\n");
730:     }
731:     PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
732:   } else if (format == PETSC_VIEWER_ASCII_SYMMODU) {
733:     PetscInt nzd=0,fshift=1,*sptr;
734:     PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
735:     PetscMalloc1(m+1,&sptr);
736:     for (i=0; i<m; i++) {
737:       sptr[i] = nzd+1;
738:       for (j=a->i[i]; j<a->i[i+1]; j++) {
739:         if (a->j[j] >= i) {
740: #if defined(PETSC_USE_COMPLEX)
741:           if (PetscImaginaryPart(a->a[j]) != 0.0 || PetscRealPart(a->a[j]) != 0.0) nzd++;
742: #else
743:           if (a->a[j] != 0.0) nzd++;
744: #endif
745:         }
746:       }
747:     }
748:     sptr[m] = nzd+1;
749:     PetscViewerASCIIPrintf(viewer," %D %D\n\n",m,nzd);
750:     for (i=0; i<m+1; i+=6) {
751:       if (i+4<m) {
752:         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]);
753:       } else if (i+3<m) {
754:         PetscViewerASCIIPrintf(viewer," %D %D %D %D %D\n",sptr[i],sptr[i+1],sptr[i+2],sptr[i+3],sptr[i+4]);
755:       } else if (i+2<m) {
756:         PetscViewerASCIIPrintf(viewer," %D %D %D %D\n",sptr[i],sptr[i+1],sptr[i+2],sptr[i+3]);
757:       } else if (i+1<m) {
758:         PetscViewerASCIIPrintf(viewer," %D %D %D\n",sptr[i],sptr[i+1],sptr[i+2]);
759:       } else if (i<m) {
760:         PetscViewerASCIIPrintf(viewer," %D %D\n",sptr[i],sptr[i+1]);
761:       } else {
762:         PetscViewerASCIIPrintf(viewer," %D\n",sptr[i]);
763:       }
764:     }
765:     PetscViewerASCIIPrintf(viewer,"\n");
766:     PetscFree(sptr);
767:     for (i=0; i<m; i++) {
768:       for (j=a->i[i]; j<a->i[i+1]; j++) {
769:         if (a->j[j] >= i) {PetscViewerASCIIPrintf(viewer," %D ",a->j[j]+fshift);}
770:       }
771:       PetscViewerASCIIPrintf(viewer,"\n");
772:     }
773:     PetscViewerASCIIPrintf(viewer,"\n");
774:     for (i=0; i<m; i++) {
775:       for (j=a->i[i]; j<a->i[i+1]; j++) {
776:         if (a->j[j] >= i) {
777: #if defined(PETSC_USE_COMPLEX)
778:           if (PetscImaginaryPart(a->a[j]) != 0.0 || PetscRealPart(a->a[j]) != 0.0) {
779:             PetscViewerASCIIPrintf(viewer," %18.16e %18.16e ",(double)PetscRealPart(a->a[j]),(double)PetscImaginaryPart(a->a[j]));
780:           }
781: #else
782:           if (a->a[j] != 0.0) {PetscViewerASCIIPrintf(viewer," %18.16e ",(double)a->a[j]);}
783: #endif
784:         }
785:       }
786:       PetscViewerASCIIPrintf(viewer,"\n");
787:     }
788:     PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
789:   } else if (format == PETSC_VIEWER_ASCII_DENSE) {
790:     PetscInt    cnt = 0,jcnt;
791:     PetscScalar value;
792: #if defined(PETSC_USE_COMPLEX)
793:     PetscBool   realonly = PETSC_TRUE;

795:     for (i=0; i<a->i[m]; i++) {
796:       if (PetscImaginaryPart(a->a[i]) != 0.0) {
797:         realonly = PETSC_FALSE;
798:         break;
799:       }
800:     }
801: #endif

803:     PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
804:     for (i=0; i<m; i++) {
805:       jcnt = 0;
806:       for (j=0; j<A->cmap->n; j++) {
807:         if (jcnt < a->i[i+1]-a->i[i] && j == a->j[cnt]) {
808:           value = a->a[cnt++];
809:           jcnt++;
810:         } else {
811:           value = 0.0;
812:         }
813: #if defined(PETSC_USE_COMPLEX)
814:         if (realonly) {
815:           PetscViewerASCIIPrintf(viewer," %7.5e ",(double)PetscRealPart(value));
816:         } else {
817:           PetscViewerASCIIPrintf(viewer," %7.5e+%7.5e i ",(double)PetscRealPart(value),(double)PetscImaginaryPart(value));
818:         }
819: #else
820:         PetscViewerASCIIPrintf(viewer," %7.5e ",(double)value);
821: #endif
822:       }
823:       PetscViewerASCIIPrintf(viewer,"\n");
824:     }
825:     PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
826:   } else if (format == PETSC_VIEWER_ASCII_MATRIXMARKET) {
827:     PetscInt fshift=1;
828:     PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
829: #if defined(PETSC_USE_COMPLEX)
830:     PetscViewerASCIIPrintf(viewer,"%%%%MatrixMarket matrix coordinate complex general\n");
831: #else
832:     PetscViewerASCIIPrintf(viewer,"%%%%MatrixMarket matrix coordinate real general\n");
833: #endif
834:     PetscViewerASCIIPrintf(viewer,"%D %D %D\n", m, A->cmap->n, a->nz);
835:     for (i=0; i<m; i++) {
836:       for (j=a->i[i]; j<a->i[i+1]; j++) {
837: #if defined(PETSC_USE_COMPLEX)
838:         PetscViewerASCIIPrintf(viewer,"%D %D %g %g\n", i+fshift,a->j[j]+fshift,(double)PetscRealPart(a->a[j]),(double)PetscImaginaryPart(a->a[j]));
839: #else
840:         PetscViewerASCIIPrintf(viewer,"%D %D %g\n", i+fshift, a->j[j]+fshift, (double)a->a[j]);
841: #endif
842:       }
843:     }
844:     PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
845:   } else {
846:     PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
847:     if (A->factortype) {
848:       for (i=0; i<m; i++) {
849:         PetscViewerASCIIPrintf(viewer,"row %D:",i);
850:         /* L part */
851:         for (j=a->i[i]; j<a->i[i+1]; j++) {
852: #if defined(PETSC_USE_COMPLEX)
853:           if (PetscImaginaryPart(a->a[j]) > 0.0) {
854:             PetscViewerASCIIPrintf(viewer," (%D, %g + %g i)",a->j[j],(double)PetscRealPart(a->a[j]),(double)PetscImaginaryPart(a->a[j]));
855:           } else if (PetscImaginaryPart(a->a[j]) < 0.0) {
856:             PetscViewerASCIIPrintf(viewer," (%D, %g - %g i)",a->j[j],(double)PetscRealPart(a->a[j]),(double)(-PetscImaginaryPart(a->a[j])));
857:           } else {
858:             PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[j],(double)PetscRealPart(a->a[j]));
859:           }
860: #else
861:           PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[j],(double)a->a[j]);
862: #endif
863:         }
864:         /* diagonal */
865:         j = a->diag[i];
866: #if defined(PETSC_USE_COMPLEX)
867:         if (PetscImaginaryPart(a->a[j]) > 0.0) {
868:           PetscViewerASCIIPrintf(viewer," (%D, %g + %g i)",a->j[j],(double)PetscRealPart(1.0/a->a[j]),(double)PetscImaginaryPart(1.0/a->a[j]));
869:         } else if (PetscImaginaryPart(a->a[j]) < 0.0) {
870:           PetscViewerASCIIPrintf(viewer," (%D, %g - %g i)",a->j[j],(double)PetscRealPart(1.0/a->a[j]),(double)(-PetscImaginaryPart(1.0/a->a[j])));
871:         } else {
872:           PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[j],(double)PetscRealPart(1.0/a->a[j]));
873:         }
874: #else
875:         PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[j],(double)(1.0/a->a[j]));
876: #endif

878:         /* U part */
879:         for (j=a->diag[i+1]+1; j<a->diag[i]; j++) {
880: #if defined(PETSC_USE_COMPLEX)
881:           if (PetscImaginaryPart(a->a[j]) > 0.0) {
882:             PetscViewerASCIIPrintf(viewer," (%D, %g + %g i)",a->j[j],(double)PetscRealPart(a->a[j]),(double)PetscImaginaryPart(a->a[j]));
883:           } else if (PetscImaginaryPart(a->a[j]) < 0.0) {
884:             PetscViewerASCIIPrintf(viewer," (%D, %g - %g i)",a->j[j],(double)PetscRealPart(a->a[j]),(double)(-PetscImaginaryPart(a->a[j])));
885:           } else {
886:             PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[j],(double)PetscRealPart(a->a[j]));
887:           }
888: #else
889:           PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[j],(double)a->a[j]);
890: #endif
891:         }
892:         PetscViewerASCIIPrintf(viewer,"\n");
893:       }
894:     } else {
895:       for (i=0; i<m; i++) {
896:         PetscViewerASCIIPrintf(viewer,"row %D:",i);
897:         for (j=a->i[i]; j<a->i[i+1]; j++) {
898: #if defined(PETSC_USE_COMPLEX)
899:           if (PetscImaginaryPart(a->a[j]) > 0.0) {
900:             PetscViewerASCIIPrintf(viewer," (%D, %g + %g i)",a->j[j],(double)PetscRealPart(a->a[j]),(double)PetscImaginaryPart(a->a[j]));
901:           } else if (PetscImaginaryPart(a->a[j]) < 0.0) {
902:             PetscViewerASCIIPrintf(viewer," (%D, %g - %g i)",a->j[j],(double)PetscRealPart(a->a[j]),(double)-PetscImaginaryPart(a->a[j]));
903:           } else {
904:             PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[j],(double)PetscRealPart(a->a[j]));
905:           }
906: #else
907:           PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[j],(double)a->a[j]);
908: #endif
909:         }
910:         PetscViewerASCIIPrintf(viewer,"\n");
911:       }
912:     }
913:     PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
914:   }
915:   PetscViewerFlush(viewer);
916:   return(0);
917: }

919:  #include <petscdraw.h>
920: PetscErrorCode MatView_SeqAIJ_Draw_Zoom(PetscDraw draw,void *Aa)
921: {
922:   Mat               A  = (Mat) Aa;
923:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
924:   PetscErrorCode    ierr;
925:   PetscInt          i,j,m = A->rmap->n;
926:   int               color;
927:   PetscReal         xl,yl,xr,yr,x_l,x_r,y_l,y_r;
928:   PetscViewer       viewer;
929:   PetscViewerFormat format;

932:   PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);
933:   PetscViewerGetFormat(viewer,&format);
934:   PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);

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

938:   if (format != PETSC_VIEWER_DRAW_CONTOUR) {
939:     PetscDrawCollectiveBegin(draw);
940:     /* Blue for negative, Cyan for zero and  Red for positive */
941:     color = PETSC_DRAW_BLUE;
942:     for (i=0; i<m; i++) {
943:       y_l = m - i - 1.0; y_r = y_l + 1.0;
944:       for (j=a->i[i]; j<a->i[i+1]; j++) {
945:         x_l = a->j[j]; x_r = x_l + 1.0;
946:         if (PetscRealPart(a->a[j]) >=  0.) continue;
947:         PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
948:       }
949:     }
950:     color = PETSC_DRAW_CYAN;
951:     for (i=0; i<m; i++) {
952:       y_l = m - i - 1.0; y_r = y_l + 1.0;
953:       for (j=a->i[i]; j<a->i[i+1]; j++) {
954:         x_l = a->j[j]; x_r = x_l + 1.0;
955:         if (a->a[j] !=  0.) continue;
956:         PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
957:       }
958:     }
959:     color = PETSC_DRAW_RED;
960:     for (i=0; i<m; i++) {
961:       y_l = m - i - 1.0; y_r = y_l + 1.0;
962:       for (j=a->i[i]; j<a->i[i+1]; j++) {
963:         x_l = a->j[j]; x_r = x_l + 1.0;
964:         if (PetscRealPart(a->a[j]) <=  0.) continue;
965:         PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
966:       }
967:     }
968:     PetscDrawCollectiveEnd(draw);
969:   } else {
970:     /* use contour shading to indicate magnitude of values */
971:     /* first determine max of all nonzero values */
972:     PetscReal minv = 0.0, maxv = 0.0;
973:     PetscInt  nz = a->nz, count = 0;
974:     PetscDraw popup;

976:     for (i=0; i<nz; i++) {
977:       if (PetscAbsScalar(a->a[i]) > maxv) maxv = PetscAbsScalar(a->a[i]);
978:     }
979:     if (minv >= maxv) maxv = minv + PETSC_SMALL;
980:     PetscDrawGetPopup(draw,&popup);
981:     PetscDrawScalePopup(popup,minv,maxv);

983:     PetscDrawCollectiveBegin(draw);
984:     for (i=0; i<m; i++) {
985:       y_l = m - i - 1.0;
986:       y_r = y_l + 1.0;
987:       for (j=a->i[i]; j<a->i[i+1]; j++) {
988:         x_l = a->j[j];
989:         x_r = x_l + 1.0;
990:         color = PetscDrawRealToColor(PetscAbsScalar(a->a[count]),minv,maxv);
991:         PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);
992:         count++;
993:       }
994:     }
995:     PetscDrawCollectiveEnd(draw);
996:   }
997:   return(0);
998: }

1000:  #include <petscdraw.h>
1001: PetscErrorCode MatView_SeqAIJ_Draw(Mat A,PetscViewer viewer)
1002: {
1004:   PetscDraw      draw;
1005:   PetscReal      xr,yr,xl,yl,h,w;
1006:   PetscBool      isnull;

1009:   PetscViewerDrawGetDraw(viewer,0,&draw);
1010:   PetscDrawIsNull(draw,&isnull);
1011:   if (isnull) return(0);

1013:   xr   = A->cmap->n; yr  = A->rmap->n; h = yr/10.0; w = xr/10.0;
1014:   xr  += w;          yr += h;         xl = -w;     yl = -h;
1015:   PetscDrawSetCoordinates(draw,xl,yl,xr,yr);
1016:   PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);
1017:   PetscDrawZoom(draw,MatView_SeqAIJ_Draw_Zoom,A);
1018:   PetscObjectCompose((PetscObject)A,"Zoomviewer",NULL);
1019:   PetscDrawSave(draw);
1020:   return(0);
1021: }

1023: PetscErrorCode MatView_SeqAIJ(Mat A,PetscViewer viewer)
1024: {
1026:   PetscBool      iascii,isbinary,isdraw;

1029:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1030:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
1031:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
1032:   if (iascii) {
1033:     MatView_SeqAIJ_ASCII(A,viewer);
1034:   } else if (isbinary) {
1035:     MatView_SeqAIJ_Binary(A,viewer);
1036:   } else if (isdraw) {
1037:     MatView_SeqAIJ_Draw(A,viewer);
1038:   }
1039:   MatView_SeqAIJ_Inode(A,viewer);
1040:   return(0);
1041: }

1043: PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode)
1044: {
1045:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1047:   PetscInt       fshift = 0,i,*ai = a->i,*aj = a->j,*imax = a->imax;
1048:   PetscInt       m      = A->rmap->n,*ip,N,*ailen = a->ilen,rmax = 0;
1049:   MatScalar      *aa    = a->a,*ap;
1050:   PetscReal      ratio  = 0.6;

1053:   if (mode == MAT_FLUSH_ASSEMBLY) return(0);
1054:   MatSeqAIJInvalidateDiagonal(A);
1055:   if (A->was_assembled && A->ass_nonzerostate == A->nonzerostate) return(0);

1057:   if (m) rmax = ailen[0]; /* determine row with most nonzeros */
1058:   for (i=1; i<m; i++) {
1059:     /* move each row back by the amount of empty slots (fshift) before it*/
1060:     fshift += imax[i-1] - ailen[i-1];
1061:     rmax    = PetscMax(rmax,ailen[i]);
1062:     if (fshift) {
1063:       ip = aj + ai[i];
1064:       ap = aa + ai[i];
1065:       N  = ailen[i];
1066:       PetscArraymove(ip-fshift,ip,N);
1067:       if (!A->structure_only) {
1068:         PetscArraymove(ap-fshift,ap,N);
1069:       }
1070:     }
1071:     ai[i] = ai[i-1] + ailen[i-1];
1072:   }
1073:   if (m) {
1074:     fshift += imax[m-1] - ailen[m-1];
1075:     ai[m]   = ai[m-1] + ailen[m-1];
1076:   }

1078:   /* reset ilen and imax for each row */
1079:   a->nonzerorowcnt = 0;
1080:   if (A->structure_only) {
1081:     PetscFree(a->imax);
1082:     PetscFree(a->ilen);
1083:   } else { /* !A->structure_only */
1084:     for (i=0; i<m; i++) {
1085:       ailen[i] = imax[i] = ai[i+1] - ai[i];
1086:       a->nonzerorowcnt += ((ai[i+1] - ai[i]) > 0);
1087:     }
1088:   }
1089:   a->nz = ai[m];
1090:   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);

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

1097:   A->info.mallocs    += a->reallocs;
1098:   a->reallocs         = 0;
1099:   A->info.nz_unneeded = (PetscReal)fshift;
1100:   a->rmax             = rmax;

1102:   if (!A->structure_only) {
1103:     MatCheckCompressedRow(A,a->nonzerorowcnt,&a->compressedrow,a->i,m,ratio);
1104:   }
1105:   MatAssemblyEnd_SeqAIJ_Inode(A,mode);
1106:   return(0);
1107: }

1109: PetscErrorCode MatRealPart_SeqAIJ(Mat A)
1110: {
1111:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1112:   PetscInt       i,nz = a->nz;
1113:   MatScalar      *aa = a->a;

1117:   for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]);
1118:   MatSeqAIJInvalidateDiagonal(A);
1119: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
1120:   if (A->offloadmask != PETSC_OFFLOAD_UNALLOCATED) A->offloadmask = PETSC_OFFLOAD_CPU;
1121: #endif
1122:   return(0);
1123: }

1125: PetscErrorCode MatImaginaryPart_SeqAIJ(Mat A)
1126: {
1127:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1128:   PetscInt       i,nz = a->nz;
1129:   MatScalar      *aa = a->a;

1133:   for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]);
1134:   MatSeqAIJInvalidateDiagonal(A);
1135: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
1136:   if (A->offloadmask != PETSC_OFFLOAD_UNALLOCATED) A->offloadmask = PETSC_OFFLOAD_CPU;
1137: #endif
1138:   return(0);
1139: }

1141: PetscErrorCode MatZeroEntries_SeqAIJ(Mat A)
1142: {
1143:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;

1147:   PetscArrayzero(a->a,a->i[A->rmap->n]);
1148:   MatSeqAIJInvalidateDiagonal(A);
1149: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
1150:   if (A->offloadmask != PETSC_OFFLOAD_UNALLOCATED) A->offloadmask = PETSC_OFFLOAD_CPU;
1151: #endif
1152:   return(0);
1153: }

1155: PetscErrorCode MatDestroy_SeqAIJ(Mat A)
1156: {
1157:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;

1161: #if defined(PETSC_USE_LOG)
1162:   PetscLogObjectState((PetscObject)A,"Rows=%D, Cols=%D, NZ=%D",A->rmap->n,A->cmap->n,a->nz);
1163: #endif
1164:   MatSeqXAIJFreeAIJ(A,&a->a,&a->j,&a->i);
1165:   ISDestroy(&a->row);
1166:   ISDestroy(&a->col);
1167:   PetscFree(a->diag);
1168:   PetscFree(a->ibdiag);
1169:   PetscFree(a->imax);
1170:   PetscFree(a->ilen);
1171:   PetscFree(a->ipre);
1172:   PetscFree3(a->idiag,a->mdiag,a->ssor_work);
1173:   PetscFree(a->solve_work);
1174:   ISDestroy(&a->icol);
1175:   PetscFree(a->saved_values);
1176:   ISColoringDestroy(&a->coloring);
1177:   PetscFree2(a->compressedrow.i,a->compressedrow.rindex);
1178:   PetscFree(a->matmult_abdense);

1180:   MatDestroy_SeqAIJ_Inode(A);
1181:   PetscFree(A->data);

1183:   PetscObjectChangeTypeName((PetscObject)A,0);
1184:   PetscObjectComposeFunction((PetscObject)A,"MatSeqAIJSetColumnIndices_C",NULL);
1185:   PetscObjectComposeFunction((PetscObject)A,"MatStoreValues_C",NULL);
1186:   PetscObjectComposeFunction((PetscObject)A,"MatRetrieveValues_C",NULL);
1187:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqaij_seqsbaij_C",NULL);
1188:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqaij_seqbaij_C",NULL);
1189:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqaij_seqaijperm_C",NULL);
1190: #if defined(PETSC_HAVE_ELEMENTAL)
1191:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqaij_elemental_C",NULL);
1192: #endif
1193: #if defined(PETSC_HAVE_HYPRE)
1194:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqaij_hypre_C",NULL);
1195:   PetscObjectComposeFunction((PetscObject)A,"MatMatMatMult_transpose_seqaij_seqaij_C",NULL);
1196: #endif
1197:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqaij_seqdense_C",NULL);
1198:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqaij_seqsell_C",NULL);
1199:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqaij_is_C",NULL);
1200:   PetscObjectComposeFunction((PetscObject)A,"MatIsTranspose_C",NULL);
1201:   PetscObjectComposeFunction((PetscObject)A,"MatSeqAIJSetPreallocation_C",NULL);
1202:   PetscObjectComposeFunction((PetscObject)A,"MatResetPreallocation_C",NULL);
1203:   PetscObjectComposeFunction((PetscObject)A,"MatSeqAIJSetPreallocationCSR_C",NULL);
1204:   PetscObjectComposeFunction((PetscObject)A,"MatReorderForNonzeroDiagonal_C",NULL);
1205:   PetscObjectComposeFunction((PetscObject)A,"MatPtAP_is_seqaij_C",NULL);
1206:   return(0);
1207: }

1209: PetscErrorCode MatSetOption_SeqAIJ(Mat A,MatOption op,PetscBool flg)
1210: {
1211:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;

1215:   switch (op) {
1216:   case MAT_ROW_ORIENTED:
1217:     a->roworiented = flg;
1218:     break;
1219:   case MAT_KEEP_NONZERO_PATTERN:
1220:     a->keepnonzeropattern = flg;
1221:     break;
1222:   case MAT_NEW_NONZERO_LOCATIONS:
1223:     a->nonew = (flg ? 0 : 1);
1224:     break;
1225:   case MAT_NEW_NONZERO_LOCATION_ERR:
1226:     a->nonew = (flg ? -1 : 0);
1227:     break;
1228:   case MAT_NEW_NONZERO_ALLOCATION_ERR:
1229:     a->nonew = (flg ? -2 : 0);
1230:     break;
1231:   case MAT_UNUSED_NONZERO_LOCATION_ERR:
1232:     a->nounused = (flg ? -1 : 0);
1233:     break;
1234:   case MAT_IGNORE_ZERO_ENTRIES:
1235:     a->ignorezeroentries = flg;
1236:     break;
1237:   case MAT_SPD:
1238:   case MAT_SYMMETRIC:
1239:   case MAT_STRUCTURALLY_SYMMETRIC:
1240:   case MAT_HERMITIAN:
1241:   case MAT_SYMMETRY_ETERNAL:
1242:   case MAT_STRUCTURE_ONLY:
1243:     /* These options are handled directly by MatSetOption() */
1244:     break;
1245:   case MAT_NEW_DIAGONALS:
1246:   case MAT_IGNORE_OFF_PROC_ENTRIES:
1247:   case MAT_USE_HASH_TABLE:
1248:     PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
1249:     break;
1250:   case MAT_USE_INODES:
1251:     /* Not an error because MatSetOption_SeqAIJ_Inode handles this one */
1252:     break;
1253:   case MAT_SUBMAT_SINGLEIS:
1254:     A->submat_singleis = flg;
1255:     break;
1256:   case MAT_SORTED_FULL:
1257:     if (flg) A->ops->setvalues = MatSetValues_SeqAIJ_SortedFull;
1258:     else     A->ops->setvalues = MatSetValues_SeqAIJ;
1259:     break;
1260:   default:
1261:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %d",op);
1262:   }
1263:   MatSetOption_SeqAIJ_Inode(A,op,flg);
1264:   return(0);
1265: }

1267: PetscErrorCode MatGetDiagonal_SeqAIJ(Mat A,Vec v)
1268: {
1269:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1271:   PetscInt       i,j,n,*ai=a->i,*aj=a->j;
1272:   PetscScalar    *aa=a->a,*x;

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

1278:   if (A->factortype == MAT_FACTOR_ILU || A->factortype == MAT_FACTOR_LU) {
1279:     PetscInt *diag=a->diag;
1280:     VecGetArrayWrite(v,&x);
1281:     for (i=0; i<n; i++) x[i] = 1.0/aa[diag[i]];
1282:     VecRestoreArrayWrite(v,&x);
1283:     return(0);
1284:   }

1286:   VecGetArrayWrite(v,&x);
1287:   for (i=0; i<n; i++) {
1288:     x[i] = 0.0;
1289:     for (j=ai[i]; j<ai[i+1]; j++) {
1290:       if (aj[j] == i) {
1291:         x[i] = aa[j];
1292:         break;
1293:       }
1294:     }
1295:   }
1296:   VecRestoreArrayWrite(v,&x);
1297:   return(0);
1298: }

1300: #include <../src/mat/impls/aij/seq/ftn-kernels/fmult.h>
1301: PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy)
1302: {
1303:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
1304:   PetscScalar       *y;
1305:   const PetscScalar *x;
1306:   PetscErrorCode    ierr;
1307:   PetscInt          m = A->rmap->n;
1308: #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
1309:   const MatScalar   *v;
1310:   PetscScalar       alpha;
1311:   PetscInt          n,i,j;
1312:   const PetscInt    *idx,*ii,*ridx=NULL;
1313:   Mat_CompressedRow cprow    = a->compressedrow;
1314:   PetscBool         usecprow = cprow.use;
1315: #endif

1318:   if (zz != yy) {VecCopy(zz,yy);}
1319:   VecGetArrayRead(xx,&x);
1320:   VecGetArray(yy,&y);

1322: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ)
1323:   fortranmulttransposeaddaij_(&m,x,a->i,a->j,a->a,y);
1324: #else
1325:   if (usecprow) {
1326:     m    = cprow.nrows;
1327:     ii   = cprow.i;
1328:     ridx = cprow.rindex;
1329:   } else {
1330:     ii = a->i;
1331:   }
1332:   for (i=0; i<m; i++) {
1333:     idx = a->j + ii[i];
1334:     v   = a->a + ii[i];
1335:     n   = ii[i+1] - ii[i];
1336:     if (usecprow) {
1337:       alpha = x[ridx[i]];
1338:     } else {
1339:       alpha = x[i];
1340:     }
1341:     for (j=0; j<n; j++) y[idx[j]] += alpha*v[j];
1342:   }
1343: #endif
1344:   PetscLogFlops(2.0*a->nz);
1345:   VecRestoreArrayRead(xx,&x);
1346:   VecRestoreArray(yy,&y);
1347:   return(0);
1348: }

1350: PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec xx,Vec yy)
1351: {

1355:   VecSet(yy,0.0);
1356:   MatMultTransposeAdd_SeqAIJ(A,xx,yy,yy);
1357:   return(0);
1358: }

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

1362: PetscErrorCode MatMult_SeqAIJ(Mat A,Vec xx,Vec yy)
1363: {
1364:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
1365:   PetscScalar       *y;
1366:   const PetscScalar *x;
1367:   const MatScalar   *aa;
1368:   PetscErrorCode    ierr;
1369:   PetscInt          m=A->rmap->n;
1370:   const PetscInt    *aj,*ii,*ridx=NULL;
1371:   PetscInt          n,i;
1372:   PetscScalar       sum;
1373:   PetscBool         usecprow=a->compressedrow.use;

1375: #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
1376: #pragma disjoint(*x,*y,*aa)
1377: #endif

1380:   VecGetArrayRead(xx,&x);
1381:   VecGetArray(yy,&y);
1382:   ii   = a->i;
1383:   if (usecprow) { /* use compressed row format */
1384:     PetscArrayzero(y,m);
1385:     m    = a->compressedrow.nrows;
1386:     ii   = a->compressedrow.i;
1387:     ridx = a->compressedrow.rindex;
1388:     for (i=0; i<m; i++) {
1389:       n           = ii[i+1] - ii[i];
1390:       aj          = a->j + ii[i];
1391:       aa          = a->a + ii[i];
1392:       sum         = 0.0;
1393:       PetscSparseDensePlusDot(sum,x,aa,aj,n);
1394:       /* for (j=0; j<n; j++) sum += (*aa++)*x[*aj++]; */
1395:       y[*ridx++] = sum;
1396:     }
1397:   } else { /* do not use compressed row format */
1398: #if defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ)
1399:     aj   = a->j;
1400:     aa   = a->a;
1401:     fortranmultaij_(&m,x,ii,aj,aa,y);
1402: #else
1403:     for (i=0; i<m; i++) {
1404:       n           = ii[i+1] - ii[i];
1405:       aj          = a->j + ii[i];
1406:       aa          = a->a + ii[i];
1407:       sum         = 0.0;
1408:       PetscSparseDensePlusDot(sum,x,aa,aj,n);
1409:       y[i] = sum;
1410:     }
1411: #endif
1412:   }
1413:   PetscLogFlops(2.0*a->nz - a->nonzerorowcnt);
1414:   VecRestoreArrayRead(xx,&x);
1415:   VecRestoreArray(yy,&y);
1416:   return(0);
1417: }

1419: PetscErrorCode MatMultMax_SeqAIJ(Mat A,Vec xx,Vec yy)
1420: {
1421:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
1422:   PetscScalar       *y;
1423:   const PetscScalar *x;
1424:   const MatScalar   *aa;
1425:   PetscErrorCode    ierr;
1426:   PetscInt          m=A->rmap->n;
1427:   const PetscInt    *aj,*ii,*ridx=NULL;
1428:   PetscInt          n,i,nonzerorow=0;
1429:   PetscScalar       sum;
1430:   PetscBool         usecprow=a->compressedrow.use;

1432: #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
1433: #pragma disjoint(*x,*y,*aa)
1434: #endif

1437:   VecGetArrayRead(xx,&x);
1438:   VecGetArray(yy,&y);
1439:   if (usecprow) { /* use compressed row format */
1440:     m    = a->compressedrow.nrows;
1441:     ii   = a->compressedrow.i;
1442:     ridx = a->compressedrow.rindex;
1443:     for (i=0; i<m; i++) {
1444:       n           = ii[i+1] - ii[i];
1445:       aj          = a->j + ii[i];
1446:       aa          = a->a + ii[i];
1447:       sum         = 0.0;
1448:       nonzerorow += (n>0);
1449:       PetscSparseDenseMaxDot(sum,x,aa,aj,n);
1450:       /* for (j=0; j<n; j++) sum += (*aa++)*x[*aj++]; */
1451:       y[*ridx++] = sum;
1452:     }
1453:   } else { /* do not use compressed row format */
1454:     ii = a->i;
1455:     for (i=0; i<m; i++) {
1456:       n           = ii[i+1] - ii[i];
1457:       aj          = a->j + ii[i];
1458:       aa          = a->a + ii[i];
1459:       sum         = 0.0;
1460:       nonzerorow += (n>0);
1461:       PetscSparseDenseMaxDot(sum,x,aa,aj,n);
1462:       y[i] = sum;
1463:     }
1464:   }
1465:   PetscLogFlops(2.0*a->nz - nonzerorow);
1466:   VecRestoreArrayRead(xx,&x);
1467:   VecRestoreArray(yy,&y);
1468:   return(0);
1469: }

1471: PetscErrorCode MatMultAddMax_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
1472: {
1473:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
1474:   PetscScalar       *y,*z;
1475:   const PetscScalar *x;
1476:   const MatScalar   *aa;
1477:   PetscErrorCode    ierr;
1478:   PetscInt          m = A->rmap->n,*aj,*ii;
1479:   PetscInt          n,i,*ridx=NULL;
1480:   PetscScalar       sum;
1481:   PetscBool         usecprow=a->compressedrow.use;

1484:   VecGetArrayRead(xx,&x);
1485:   VecGetArrayPair(yy,zz,&y,&z);
1486:   if (usecprow) { /* use compressed row format */
1487:     if (zz != yy) {
1488:       PetscArraycpy(z,y,m);
1489:     }
1490:     m    = a->compressedrow.nrows;
1491:     ii   = a->compressedrow.i;
1492:     ridx = a->compressedrow.rindex;
1493:     for (i=0; i<m; i++) {
1494:       n   = ii[i+1] - ii[i];
1495:       aj  = a->j + ii[i];
1496:       aa  = a->a + ii[i];
1497:       sum = y[*ridx];
1498:       PetscSparseDenseMaxDot(sum,x,aa,aj,n);
1499:       z[*ridx++] = sum;
1500:     }
1501:   } else { /* do not use compressed row format */
1502:     ii = a->i;
1503:     for (i=0; i<m; i++) {
1504:       n   = ii[i+1] - ii[i];
1505:       aj  = a->j + ii[i];
1506:       aa  = a->a + ii[i];
1507:       sum = y[i];
1508:       PetscSparseDenseMaxDot(sum,x,aa,aj,n);
1509:       z[i] = sum;
1510:     }
1511:   }
1512:   PetscLogFlops(2.0*a->nz);
1513:   VecRestoreArrayRead(xx,&x);
1514:   VecRestoreArrayPair(yy,zz,&y,&z);
1515:   return(0);
1516: }

1518: #include <../src/mat/impls/aij/seq/ftn-kernels/fmultadd.h>
1519: PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz)
1520: {
1521:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
1522:   PetscScalar       *y,*z;
1523:   const PetscScalar *x;
1524:   const MatScalar   *aa;
1525:   PetscErrorCode    ierr;
1526:   const PetscInt    *aj,*ii,*ridx=NULL;
1527:   PetscInt          m = A->rmap->n,n,i;
1528:   PetscScalar       sum;
1529:   PetscBool         usecprow=a->compressedrow.use;

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

1572: /*
1573:      Adds diagonal pointers to sparse matrix structure.
1574: */
1575: PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat A)
1576: {
1577:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1579:   PetscInt       i,j,m = A->rmap->n;

1582:   if (!a->diag) {
1583:     PetscMalloc1(m,&a->diag);
1584:     PetscLogObjectMemory((PetscObject)A, m*sizeof(PetscInt));
1585:   }
1586:   for (i=0; i<A->rmap->n; i++) {
1587:     a->diag[i] = a->i[i+1];
1588:     for (j=a->i[i]; j<a->i[i+1]; j++) {
1589:       if (a->j[j] == i) {
1590:         a->diag[i] = j;
1591:         break;
1592:       }
1593:     }
1594:   }
1595:   return(0);
1596: }

1598: PetscErrorCode MatShift_SeqAIJ(Mat A,PetscScalar v)
1599: {
1600:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
1601:   const PetscInt    *diag = (const PetscInt*)a->diag;
1602:   const PetscInt    *ii = (const PetscInt*) a->i;
1603:   PetscInt          i,*mdiag = NULL;
1604:   PetscErrorCode    ierr;
1605:   PetscInt          cnt = 0; /* how many diagonals are missing */

1608:   if (!A->preallocated || !a->nz) {
1609:     MatSeqAIJSetPreallocation(A,1,NULL);
1610:     MatShift_Basic(A,v);
1611:     return(0);
1612:   }

1614:   if (a->diagonaldense) {
1615:     cnt = 0;
1616:   } else {
1617:     PetscCalloc1(A->rmap->n,&mdiag);
1618:     for (i=0; i<A->rmap->n; i++) {
1619:       if (diag[i] >= ii[i+1]) {
1620:         cnt++;
1621:         mdiag[i] = 1;
1622:       }
1623:     }
1624:   }
1625:   if (!cnt) {
1626:     MatShift_Basic(A,v);
1627:   } else {
1628:     PetscScalar *olda = a->a;  /* preserve pointers to current matrix nonzeros structure and values */
1629:     PetscInt    *oldj = a->j, *oldi = a->i;
1630:     PetscBool   singlemalloc = a->singlemalloc,free_a = a->free_a,free_ij = a->free_ij;

1632:     a->a = NULL;
1633:     a->j = NULL;
1634:     a->i = NULL;
1635:     /* increase the values in imax for each row where a diagonal is being inserted then reallocate the matrix data structures */
1636:     for (i=0; i<A->rmap->n; i++) {
1637:       a->imax[i] += mdiag[i];
1638:       a->imax[i] = PetscMin(a->imax[i],A->cmap->n);
1639:     }
1640:     MatSeqAIJSetPreallocation_SeqAIJ(A,0,a->imax);

1642:     /* copy old values into new matrix data structure */
1643:     for (i=0; i<A->rmap->n; i++) {
1644:       MatSetValues(A,1,&i,a->imax[i] - mdiag[i],&oldj[oldi[i]],&olda[oldi[i]],ADD_VALUES);
1645:       if (i < A->cmap->n) {
1646:         MatSetValue(A,i,i,v,ADD_VALUES);
1647:       }
1648:     }
1649:     MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
1650:     MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
1651:     if (singlemalloc) {
1652:       PetscFree3(olda,oldj,oldi);
1653:     } else {
1654:       if (free_a)  {PetscFree(olda);}
1655:       if (free_ij) {PetscFree(oldj);}
1656:       if (free_ij) {PetscFree(oldi);}
1657:     }
1658:   }
1659:   PetscFree(mdiag);
1660:   a->diagonaldense = PETSC_TRUE;
1661:   return(0);
1662: }

1664: /*
1665:      Checks for missing diagonals
1666: */
1667: PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat A,PetscBool  *missing,PetscInt *d)
1668: {
1669:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
1670:   PetscInt       *diag,*ii = a->i,i;

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

1695:  #include <petscblaslapack.h>
1696:  #include <petsc/private/kernels/blockinvert.h>

1698: /*
1699:     Note that values is allocated externally by the PC and then passed into this routine
1700: */
1701: PetscErrorCode MatInvertVariableBlockDiagonal_SeqAIJ(Mat A,PetscInt nblocks,const PetscInt *bsizes,PetscScalar *diag)
1702: {
1703:   PetscErrorCode  ierr;
1704:   PetscInt        n = A->rmap->n, i, ncnt = 0, *indx,j,bsizemax = 0,*v_pivots;
1705:   PetscBool       allowzeropivot,zeropivotdetected=PETSC_FALSE;
1706:   const PetscReal shift = 0.0;
1707:   PetscInt        ipvt[5];
1708:   PetscScalar     work[25],*v_work;

1711:   allowzeropivot = PetscNot(A->erroriffailure);
1712:   for (i=0; i<nblocks; i++) ncnt += bsizes[i];
1713:   if (ncnt != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Total blocksizes %D doesn't match number matrix rows %D",ncnt,n);
1714:   for (i=0; i<nblocks; i++) {
1715:     bsizemax = PetscMax(bsizemax,bsizes[i]);
1716:   }
1717:   PetscMalloc1(bsizemax,&indx);
1718:   if (bsizemax > 7) {
1719:     PetscMalloc2(bsizemax,&v_work,bsizemax,&v_pivots);
1720:   }
1721:   ncnt = 0;
1722:   for (i=0; i<nblocks; i++) {
1723:     for (j=0; j<bsizes[i]; j++) indx[j] = ncnt+j;
1724:     MatGetValues(A,bsizes[i],indx,bsizes[i],indx,diag);
1725:     switch (bsizes[i]) {
1726:     case 1:
1727:       *diag = 1.0/(*diag);
1728:       break;
1729:     case 2:
1730:       PetscKernel_A_gets_inverse_A_2(diag,shift,allowzeropivot,&zeropivotdetected);
1731:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1732:       PetscKernel_A_gets_transpose_A_2(diag);
1733:       break;
1734:     case 3:
1735:       PetscKernel_A_gets_inverse_A_3(diag,shift,allowzeropivot,&zeropivotdetected);
1736:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1737:       PetscKernel_A_gets_transpose_A_3(diag);
1738:       break;
1739:     case 4:
1740:       PetscKernel_A_gets_inverse_A_4(diag,shift,allowzeropivot,&zeropivotdetected);
1741:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1742:       PetscKernel_A_gets_transpose_A_4(diag);
1743:       break;
1744:     case 5:
1745:       PetscKernel_A_gets_inverse_A_5(diag,ipvt,work,shift,allowzeropivot,&zeropivotdetected);
1746:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1747:       PetscKernel_A_gets_transpose_A_5(diag);
1748:       break;
1749:     case 6:
1750:       PetscKernel_A_gets_inverse_A_6(diag,shift,allowzeropivot,&zeropivotdetected);
1751:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1752:       PetscKernel_A_gets_transpose_A_6(diag);
1753:       break;
1754:     case 7:
1755:       PetscKernel_A_gets_inverse_A_7(diag,shift,allowzeropivot,&zeropivotdetected);
1756:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1757:       PetscKernel_A_gets_transpose_A_7(diag);
1758:       break;
1759:     default:
1760:       PetscKernel_A_gets_inverse_A(bsizes[i],diag,v_pivots,v_work,allowzeropivot,&zeropivotdetected);
1761:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1762:       PetscKernel_A_gets_transpose_A_N(diag,bsizes[i]);
1763:     }
1764:     ncnt   += bsizes[i];
1765:     diag += bsizes[i]*bsizes[i];
1766:   }
1767:   if (bsizemax > 7) {
1768:     PetscFree2(v_work,v_pivots);
1769:   }
1770:   PetscFree(indx);
1771:   return(0);
1772: }

1774: /*
1775:    Negative shift indicates do not generate an error if there is a zero diagonal, just invert it anyways
1776: */
1777: PetscErrorCode  MatInvertDiagonal_SeqAIJ(Mat A,PetscScalar omega,PetscScalar fshift)
1778: {
1779:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*) A->data;
1781:   PetscInt       i,*diag,m = A->rmap->n;
1782:   MatScalar      *v = a->a;
1783:   PetscScalar    *idiag,*mdiag;

1786:   if (a->idiagvalid) return(0);
1787:   MatMarkDiagonal_SeqAIJ(A);
1788:   diag = a->diag;
1789:   if (!a->idiag) {
1790:     PetscMalloc3(m,&a->idiag,m,&a->mdiag,m,&a->ssor_work);
1791:     PetscLogObjectMemory((PetscObject)A, 3*m*sizeof(PetscScalar));
1792:     v    = a->a;
1793:   }
1794:   mdiag = a->mdiag;
1795:   idiag = a->idiag;

1797:   if (omega == 1.0 && PetscRealPart(fshift) <= 0.0) {
1798:     for (i=0; i<m; i++) {
1799:       mdiag[i] = v[diag[i]];
1800:       if (!PetscAbsScalar(mdiag[i])) { /* zero diagonal */
1801:         if (PetscRealPart(fshift)) {
1802:           PetscInfo1(A,"Zero diagonal on row %D\n",i);
1803:           A->factorerrortype             = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1804:           A->factorerror_zeropivot_value = 0.0;
1805:           A->factorerror_zeropivot_row   = i;
1806:         } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Zero diagonal on row %D",i);
1807:       }
1808:       idiag[i] = 1.0/v[diag[i]];
1809:     }
1810:     PetscLogFlops(m);
1811:   } else {
1812:     for (i=0; i<m; i++) {
1813:       mdiag[i] = v[diag[i]];
1814:       idiag[i] = omega/(fshift + v[diag[i]]);
1815:     }
1816:     PetscLogFlops(2.0*m);
1817:   }
1818:   a->idiagvalid = PETSC_TRUE;
1819:   return(0);
1820: }

1822: #include <../src/mat/impls/aij/seq/ftn-kernels/frelax.h>
1823: PetscErrorCode MatSOR_SeqAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
1824: {
1825:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
1826:   PetscScalar       *x,d,sum,*t,scale;
1827:   const MatScalar   *v,*idiag=0,*mdiag;
1828:   const PetscScalar *b, *bs,*xb, *ts;
1829:   PetscErrorCode    ierr;
1830:   PetscInt          n,m = A->rmap->n,i;
1831:   const PetscInt    *idx,*diag;

1834:   its = its*lits;

1836:   if (fshift != a->fshift || omega != a->omega) a->idiagvalid = PETSC_FALSE; /* must recompute idiag[] */
1837:   if (!a->idiagvalid) {MatInvertDiagonal_SeqAIJ(A,omega,fshift);}
1838:   a->fshift = fshift;
1839:   a->omega  = omega;

1841:   diag  = a->diag;
1842:   t     = a->ssor_work;
1843:   idiag = a->idiag;
1844:   mdiag = a->mdiag;

1846:   VecGetArray(xx,&x);
1847:   VecGetArrayRead(bb,&b);
1848:   /* We count flops by assuming the upper triangular and lower triangular parts have the same number of nonzeros */
1849:   if (flag == SOR_APPLY_UPPER) {
1850:     /* apply (U + D/omega) to the vector */
1851:     bs = b;
1852:     for (i=0; i<m; i++) {
1853:       d   = fshift + mdiag[i];
1854:       n   = a->i[i+1] - diag[i] - 1;
1855:       idx = a->j + diag[i] + 1;
1856:       v   = a->a + diag[i] + 1;
1857:       sum = b[i]*d/omega;
1858:       PetscSparseDensePlusDot(sum,bs,v,idx,n);
1859:       x[i] = sum;
1860:     }
1861:     VecRestoreArray(xx,&x);
1862:     VecRestoreArrayRead(bb,&b);
1863:     PetscLogFlops(a->nz);
1864:     return(0);
1865:   }

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

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

1874:     to a vector efficiently using Eisenstat's trick.
1875:     */
1876:     scale = (2.0/omega) - 1.0;

1878:     /*  x = (E + U)^{-1} b */
1879:     for (i=m-1; i>=0; i--) {
1880:       n   = a->i[i+1] - diag[i] - 1;
1881:       idx = a->j + diag[i] + 1;
1882:       v   = a->a + diag[i] + 1;
1883:       sum = b[i];
1884:       PetscSparseDenseMinusDot(sum,x,v,idx,n);
1885:       x[i] = sum*idiag[i];
1886:     }

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

1892:     /*  t = (E + L)^{-1}t */
1893:     ts   = t;
1894:     diag = a->diag;
1895:     for (i=0; i<m; i++) {
1896:       n   = diag[i] - a->i[i];
1897:       idx = a->j + a->i[i];
1898:       v   = a->a + a->i[i];
1899:       sum = t[i];
1900:       PetscSparseDenseMinusDot(sum,ts,v,idx,n);
1901:       t[i] = sum*idiag[i];
1902:       /*  x = x + t */
1903:       x[i] += t[i];
1904:     }

1906:     PetscLogFlops(6.0*m-1 + 2.0*a->nz);
1907:     VecRestoreArray(xx,&x);
1908:     VecRestoreArrayRead(bb,&b);
1909:     return(0);
1910:   }
1911:   if (flag & SOR_ZERO_INITIAL_GUESS) {
1912:     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
1913:       for (i=0; i<m; i++) {
1914:         n   = diag[i] - a->i[i];
1915:         idx = a->j + a->i[i];
1916:         v   = a->a + a->i[i];
1917:         sum = b[i];
1918:         PetscSparseDenseMinusDot(sum,x,v,idx,n);
1919:         t[i] = sum;
1920:         x[i] = sum*idiag[i];
1921:       }
1922:       xb   = t;
1923:       PetscLogFlops(a->nz);
1924:     } else xb = b;
1925:     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
1926:       for (i=m-1; i>=0; i--) {
1927:         n   = a->i[i+1] - diag[i] - 1;
1928:         idx = a->j + diag[i] + 1;
1929:         v   = a->a + diag[i] + 1;
1930:         sum = xb[i];
1931:         PetscSparseDenseMinusDot(sum,x,v,idx,n);
1932:         if (xb == b) {
1933:           x[i] = sum*idiag[i];
1934:         } else {
1935:           x[i] = (1-omega)*x[i] + sum*idiag[i];  /* omega in idiag */
1936:         }
1937:       }
1938:       PetscLogFlops(a->nz); /* assumes 1/2 in upper */
1939:     }
1940:     its--;
1941:   }
1942:   while (its--) {
1943:     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
1944:       for (i=0; i<m; i++) {
1945:         /* lower */
1946:         n   = diag[i] - a->i[i];
1947:         idx = a->j + a->i[i];
1948:         v   = a->a + a->i[i];
1949:         sum = b[i];
1950:         PetscSparseDenseMinusDot(sum,x,v,idx,n);
1951:         t[i] = sum;             /* save application of the lower-triangular part */
1952:         /* upper */
1953:         n   = a->i[i+1] - diag[i] - 1;
1954:         idx = a->j + diag[i] + 1;
1955:         v   = a->a + diag[i] + 1;
1956:         PetscSparseDenseMinusDot(sum,x,v,idx,n);
1957:         x[i] = (1. - omega)*x[i] + sum*idiag[i]; /* omega in idiag */
1958:       }
1959:       xb   = t;
1960:       PetscLogFlops(2.0*a->nz);
1961:     } else xb = b;
1962:     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
1963:       for (i=m-1; i>=0; i--) {
1964:         sum = xb[i];
1965:         if (xb == b) {
1966:           /* whole matrix (no checkpointing available) */
1967:           n   = a->i[i+1] - a->i[i];
1968:           idx = a->j + a->i[i];
1969:           v   = a->a + a->i[i];
1970:           PetscSparseDenseMinusDot(sum,x,v,idx,n);
1971:           x[i] = (1. - omega)*x[i] + (sum + mdiag[i]*x[i])*idiag[i];
1972:         } else { /* lower-triangular part has been saved, so only apply upper-triangular */
1973:           n   = a->i[i+1] - diag[i] - 1;
1974:           idx = a->j + diag[i] + 1;
1975:           v   = a->a + diag[i] + 1;
1976:           PetscSparseDenseMinusDot(sum,x,v,idx,n);
1977:           x[i] = (1. - omega)*x[i] + sum*idiag[i];  /* omega in idiag */
1978:         }
1979:       }
1980:       if (xb == b) {
1981:         PetscLogFlops(2.0*a->nz);
1982:       } else {
1983:         PetscLogFlops(a->nz); /* assumes 1/2 in upper */
1984:       }
1985:     }
1986:   }
1987:   VecRestoreArray(xx,&x);
1988:   VecRestoreArrayRead(bb,&b);
1989:   return(0);
1990: }


1993: PetscErrorCode MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,MatInfo *info)
1994: {
1995:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;

1998:   info->block_size   = 1.0;
1999:   info->nz_allocated = a->maxnz;
2000:   info->nz_used      = a->nz;
2001:   info->nz_unneeded  = (a->maxnz - a->nz);
2002:   info->assemblies   = A->num_ass;
2003:   info->mallocs      = A->info.mallocs;
2004:   info->memory       = ((PetscObject)A)->mem;
2005:   if (A->factortype) {
2006:     info->fill_ratio_given  = A->info.fill_ratio_given;
2007:     info->fill_ratio_needed = A->info.fill_ratio_needed;
2008:     info->factor_mallocs    = A->info.factor_mallocs;
2009:   } else {
2010:     info->fill_ratio_given  = 0;
2011:     info->fill_ratio_needed = 0;
2012:     info->factor_mallocs    = 0;
2013:   }
2014:   return(0);
2015: }

2017: PetscErrorCode MatZeroRows_SeqAIJ(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
2018: {
2019:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
2020:   PetscInt          i,m = A->rmap->n - 1;
2021:   PetscErrorCode    ierr;
2022:   const PetscScalar *xx;
2023:   PetscScalar       *bb;
2024:   PetscInt          d = 0;

2027:   if (x && b) {
2028:     VecGetArrayRead(x,&xx);
2029:     VecGetArray(b,&bb);
2030:     for (i=0; i<N; i++) {
2031:       if (rows[i] < 0 || rows[i] > m) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range", rows[i]);
2032:       if (rows[i] >= A->cmap->n) continue;
2033:       bb[rows[i]] = diag*xx[rows[i]];
2034:     }
2035:     VecRestoreArrayRead(x,&xx);
2036:     VecRestoreArray(b,&bb);
2037:   }

2039:   if (a->keepnonzeropattern) {
2040:     for (i=0; i<N; i++) {
2041:       if (rows[i] < 0 || rows[i] > m) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range", rows[i]);
2042:       PetscArrayzero(&a->a[a->i[rows[i]]],a->ilen[rows[i]]);
2043:     }
2044:     if (diag != 0.0) {
2045:       for (i=0; i<N; i++) {
2046:         d = rows[i];
2047:         if (rows[i] >= A->cmap->n) continue;
2048:         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);
2049:       }
2050:       for (i=0; i<N; i++) {
2051:         if (rows[i] >= A->cmap->n) continue;
2052:         a->a[a->diag[rows[i]]] = diag;
2053:       }
2054:     }
2055:   } else {
2056:     if (diag != 0.0) {
2057:       for (i=0; i<N; i++) {
2058:         if (rows[i] < 0 || rows[i] > m) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range", rows[i]);
2059:         if (a->ilen[rows[i]] > 0) {
2060:           if (rows[i] >= A->cmap->n) {
2061:             a->ilen[rows[i]] = 0;
2062:           } else {
2063:             a->ilen[rows[i]]    = 1;
2064:             a->a[a->i[rows[i]]] = diag;
2065:             a->j[a->i[rows[i]]] = rows[i];
2066:           }
2067:         } else if (rows[i] < A->cmap->n) { /* in case row was completely empty */
2068:           MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],&diag,INSERT_VALUES);
2069:         }
2070:       }
2071:     } else {
2072:       for (i=0; i<N; i++) {
2073:         if (rows[i] < 0 || rows[i] > m) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range", rows[i]);
2074:         a->ilen[rows[i]] = 0;
2075:       }
2076:     }
2077:     A->nonzerostate++;
2078:   }
2079: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
2080:   if (A->offloadmask != PETSC_OFFLOAD_UNALLOCATED) A->offloadmask = PETSC_OFFLOAD_CPU;
2081: #endif
2082:   (*A->ops->assemblyend)(A,MAT_FINAL_ASSEMBLY);
2083:   return(0);
2084: }

2086: PetscErrorCode MatZeroRowsColumns_SeqAIJ(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
2087: {
2088:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
2089:   PetscInt          i,j,m = A->rmap->n - 1,d = 0;
2090:   PetscErrorCode    ierr;
2091:   PetscBool         missing,*zeroed,vecs = PETSC_FALSE;
2092:   const PetscScalar *xx;
2093:   PetscScalar       *bb;

2096:   if (x && b) {
2097:     VecGetArrayRead(x,&xx);
2098:     VecGetArray(b,&bb);
2099:     vecs = PETSC_TRUE;
2100:   }
2101:   PetscCalloc1(A->rmap->n,&zeroed);
2102:   for (i=0; i<N; i++) {
2103:     if (rows[i] < 0 || rows[i] > m) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range", rows[i]);
2104:     PetscArrayzero(&a->a[a->i[rows[i]]],a->ilen[rows[i]]);

2106:     zeroed[rows[i]] = PETSC_TRUE;
2107:   }
2108:   for (i=0; i<A->rmap->n; i++) {
2109:     if (!zeroed[i]) {
2110:       for (j=a->i[i]; j<a->i[i+1]; j++) {
2111:         if (a->j[j] < A->rmap->n && zeroed[a->j[j]]) {
2112:           if (vecs) bb[i] -= a->a[j]*xx[a->j[j]];
2113:           a->a[j] = 0.0;
2114:         }
2115:       }
2116:     } else if (vecs && i < A->cmap->N) bb[i] = diag*xx[i];
2117:   }
2118:   if (x && b) {
2119:     VecRestoreArrayRead(x,&xx);
2120:     VecRestoreArray(b,&bb);
2121:   }
2122:   PetscFree(zeroed);
2123:   if (diag != 0.0) {
2124:     MatMissingDiagonal_SeqAIJ(A,&missing,&d);
2125:     if (missing) {
2126:       for (i=0; i<N; i++) {
2127:         if (rows[i] >= A->cmap->N) continue;
2128:         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]);
2129:         MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],&diag,INSERT_VALUES);
2130:       }
2131:     } else {
2132:       for (i=0; i<N; i++) {
2133:         a->a[a->diag[rows[i]]] = diag;
2134:       }
2135:     }
2136:   }
2137: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
2138:   if (A->offloadmask != PETSC_OFFLOAD_UNALLOCATED) A->offloadmask = PETSC_OFFLOAD_CPU;
2139: #endif
2140:   (*A->ops->assemblyend)(A,MAT_FINAL_ASSEMBLY);
2141:   return(0);
2142: }

2144: PetscErrorCode MatGetRow_SeqAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
2145: {
2146:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
2147:   PetscInt   *itmp;

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

2152:   *nz = a->i[row+1] - a->i[row];
2153:   if (v) *v = a->a + a->i[row];
2154:   if (idx) {
2155:     itmp = a->j + a->i[row];
2156:     if (*nz) *idx = itmp;
2157:     else *idx = 0;
2158:   }
2159:   return(0);
2160: }

2162: /* remove this function? */
2163: PetscErrorCode MatRestoreRow_SeqAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
2164: {
2166:   return(0);
2167: }

2169: PetscErrorCode MatNorm_SeqAIJ(Mat A,NormType type,PetscReal *nrm)
2170: {
2171:   Mat_SeqAIJ     *a  = (Mat_SeqAIJ*)A->data;
2172:   MatScalar      *v  = a->a;
2173:   PetscReal      sum = 0.0;
2175:   PetscInt       i,j;

2178:   if (type == NORM_FROBENIUS) {
2179: #if defined(PETSC_USE_REAL___FP16)
2180:     PetscBLASInt one = 1,nz = a->nz;
2181:     *nrm = BLASnrm2_(&nz,v,&one);
2182: #else
2183:     for (i=0; i<a->nz; i++) {
2184:       sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
2185:     }
2186:     *nrm = PetscSqrtReal(sum);
2187: #endif
2188:     PetscLogFlops(2*a->nz);
2189:   } else if (type == NORM_1) {
2190:     PetscReal *tmp;
2191:     PetscInt  *jj = a->j;
2192:     PetscCalloc1(A->cmap->n+1,&tmp);
2193:     *nrm = 0.0;
2194:     for (j=0; j<a->nz; j++) {
2195:       tmp[*jj++] += PetscAbsScalar(*v);  v++;
2196:     }
2197:     for (j=0; j<A->cmap->n; j++) {
2198:       if (tmp[j] > *nrm) *nrm = tmp[j];
2199:     }
2200:     PetscFree(tmp);
2201:     PetscLogFlops(PetscMax(a->nz-1,0));
2202:   } else if (type == NORM_INFINITY) {
2203:     *nrm = 0.0;
2204:     for (j=0; j<A->rmap->n; j++) {
2205:       v   = a->a + a->i[j];
2206:       sum = 0.0;
2207:       for (i=0; i<a->i[j+1]-a->i[j]; i++) {
2208:         sum += PetscAbsScalar(*v); v++;
2209:       }
2210:       if (sum > *nrm) *nrm = sum;
2211:     }
2212:     PetscLogFlops(PetscMax(a->nz-1,0));
2213:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for two norm");
2214:   return(0);
2215: }

2217: /* Merged from MatGetSymbolicTranspose_SeqAIJ() - replace MatGetSymbolicTranspose_SeqAIJ()? */
2218: PetscErrorCode MatTransposeSymbolic_SeqAIJ(Mat A,Mat *B)
2219: {
2221:   PetscInt       i,j,anzj;
2222:   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b;
2223:   PetscInt       an=A->cmap->N,am=A->rmap->N;
2224:   PetscInt       *ati,*atj,*atfill,*ai=a->i,*aj=a->j;

2227:   /* Allocate space for symbolic transpose info and work array */
2228:   PetscCalloc1(an+1,&ati);
2229:   PetscMalloc1(ai[am],&atj);
2230:   PetscMalloc1(an,&atfill);

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

2238:   /* Copy ati into atfill so we have locations of the next free space in atj */
2239:   PetscArraycpy(atfill,ati,an);

2241:   /* Walk through A row-wise and mark nonzero entries of A^T. */
2242:   for (i=0;i<am;i++) {
2243:     anzj = ai[i+1] - ai[i];
2244:     for (j=0;j<anzj;j++) {
2245:       atj[atfill[*aj]] = i;
2246:       atfill[*aj++]   += 1;
2247:     }
2248:   }

2250:   /* Clean up temporary space and complete requests. */
2251:   PetscFree(atfill);
2252:   MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),an,am,ati,atj,NULL,B);
2253:   MatSetBlockSizes(*B,PetscAbs(A->cmap->bs),PetscAbs(A->rmap->bs));
2254:   MatSetType(*B,((PetscObject)A)->type_name);

2256:   b          = (Mat_SeqAIJ*)((*B)->data);
2257:   b->free_a  = PETSC_FALSE;
2258:   b->free_ij = PETSC_TRUE;
2259:   b->nonew   = 0;
2260:   return(0);
2261: }

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

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

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

2311: PetscErrorCode  MatIsHermitianTranspose_SeqAIJ(Mat A,Mat B,PetscReal tol,PetscBool  *f)
2312: {
2313:   Mat_SeqAIJ     *aij = (Mat_SeqAIJ*) A->data,*bij = (Mat_SeqAIJ*) B->data;
2314:   PetscInt       *adx,*bdx,*aii,*bii,*aptr,*bptr;
2315:   MatScalar      *va,*vb;
2317:   PetscInt       ma,na,mb,nb, i;

2320:   MatGetSize(A,&ma,&na);
2321:   MatGetSize(B,&mb,&nb);
2322:   if (ma!=nb || na!=mb) {
2323:     *f = PETSC_FALSE;
2324:     return(0);
2325:   }
2326:   aii  = aij->i; bii = bij->i;
2327:   adx  = aij->j; bdx = bij->j;
2328:   va   = aij->a; vb = bij->a;
2329:   PetscMalloc1(ma,&aptr);
2330:   PetscMalloc1(mb,&bptr);
2331:   for (i=0; i<ma; i++) aptr[i] = aii[i];
2332:   for (i=0; i<mb; i++) bptr[i] = bii[i];

2334:   *f = PETSC_TRUE;
2335:   for (i=0; i<ma; i++) {
2336:     while (aptr[i]<aii[i+1]) {
2337:       PetscInt    idc,idr;
2338:       PetscScalar vc,vr;
2339:       /* column/row index/value */
2340:       idc = adx[aptr[i]];
2341:       idr = bdx[bptr[idc]];
2342:       vc  = va[aptr[i]];
2343:       vr  = vb[bptr[idc]];
2344:       if (i!=idr || PetscAbsScalar(vc-PetscConj(vr)) > tol) {
2345:         *f = PETSC_FALSE;
2346:         goto done;
2347:       } else {
2348:         aptr[i]++;
2349:         if (B || i!=idc) bptr[idc]++;
2350:       }
2351:     }
2352:   }
2353: done:
2354:   PetscFree(aptr);
2355:   PetscFree(bptr);
2356:   return(0);
2357: }

2359: PetscErrorCode MatIsSymmetric_SeqAIJ(Mat A,PetscReal tol,PetscBool  *f)
2360: {

2364:   MatIsTranspose_SeqAIJ(A,A,tol,f);
2365:   return(0);
2366: }

2368: PetscErrorCode MatIsHermitian_SeqAIJ(Mat A,PetscReal tol,PetscBool  *f)
2369: {

2373:   MatIsHermitianTranspose_SeqAIJ(A,A,tol,f);
2374:   return(0);
2375: }

2377: PetscErrorCode MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr)
2378: {
2379:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
2380:   const PetscScalar *l,*r;
2381:   PetscScalar       x;
2382:   MatScalar         *v;
2383:   PetscErrorCode    ierr;
2384:   PetscInt          i,j,m = A->rmap->n,n = A->cmap->n,M,nz = a->nz;
2385:   const PetscInt    *jj;

2388:   if (ll) {
2389:     /* The local size is used so that VecMPI can be passed to this routine
2390:        by MatDiagonalScale_MPIAIJ */
2391:     VecGetLocalSize(ll,&m);
2392:     if (m != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length");
2393:     VecGetArrayRead(ll,&l);
2394:     v    = a->a;
2395:     for (i=0; i<m; i++) {
2396:       x = l[i];
2397:       M = a->i[i+1] - a->i[i];
2398:       for (j=0; j<M; j++) (*v++) *= x;
2399:     }
2400:     VecRestoreArrayRead(ll,&l);
2401:     PetscLogFlops(nz);
2402:   }
2403:   if (rr) {
2404:     VecGetLocalSize(rr,&n);
2405:     if (n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vector wrong length");
2406:     VecGetArrayRead(rr,&r);
2407:     v    = a->a; jj = a->j;
2408:     for (i=0; i<nz; i++) (*v++) *= r[*jj++];
2409:     VecRestoreArrayRead(rr,&r);
2410:     PetscLogFlops(nz);
2411:   }
2412:   MatSeqAIJInvalidateDiagonal(A);
2413: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
2414:   if (A->offloadmask != PETSC_OFFLOAD_UNALLOCATED) A->offloadmask = PETSC_OFFLOAD_CPU;
2415: #endif
2416:   return(0);
2417: }

2419: PetscErrorCode MatCreateSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,PetscInt csize,MatReuse scall,Mat *B)
2420: {
2421:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data,*c;
2423:   PetscInt       *smap,i,k,kstart,kend,oldcols = A->cmap->n,*lens;
2424:   PetscInt       row,mat_i,*mat_j,tcol,first,step,*mat_ilen,sum,lensi;
2425:   const PetscInt *irow,*icol;
2426:   PetscInt       nrows,ncols;
2427:   PetscInt       *starts,*j_new,*i_new,*aj = a->j,*ai = a->i,ii,*ailen = a->ilen;
2428:   MatScalar      *a_new,*mat_a;
2429:   Mat            C;
2430:   PetscBool      stride;


2434:   ISGetIndices(isrow,&irow);
2435:   ISGetLocalSize(isrow,&nrows);
2436:   ISGetLocalSize(iscol,&ncols);

2438:   PetscObjectTypeCompare((PetscObject)iscol,ISSTRIDE,&stride);
2439:   if (stride) {
2440:     ISStrideGetInfo(iscol,&first,&step);
2441:   } else {
2442:     first = 0;
2443:     step  = 0;
2444:   }
2445:   if (stride && step == 1) {
2446:     /* special case of contiguous rows */
2447:     PetscMalloc2(nrows,&lens,nrows,&starts);
2448:     /* loop over new rows determining lens and starting points */
2449:     for (i=0; i<nrows; i++) {
2450:       kstart = ai[irow[i]];
2451:       kend   = kstart + ailen[irow[i]];
2452:       starts[i] = kstart;
2453:       for (k=kstart; k<kend; k++) {
2454:         if (aj[k] >= first) {
2455:           starts[i] = k;
2456:           break;
2457:         }
2458:       }
2459:       sum = 0;
2460:       while (k < kend) {
2461:         if (aj[k++] >= first+ncols) break;
2462:         sum++;
2463:       }
2464:       lens[i] = sum;
2465:     }
2466:     /* create submatrix */
2467:     if (scall == MAT_REUSE_MATRIX) {
2468:       PetscInt n_cols,n_rows;
2469:       MatGetSize(*B,&n_rows,&n_cols);
2470:       if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size");
2471:       MatZeroEntries(*B);
2472:       C    = *B;
2473:     } else {
2474:       PetscInt rbs,cbs;
2475:       MatCreate(PetscObjectComm((PetscObject)A),&C);
2476:       MatSetSizes(C,nrows,ncols,PETSC_DETERMINE,PETSC_DETERMINE);
2477:       ISGetBlockSize(isrow,&rbs);
2478:       ISGetBlockSize(iscol,&cbs);
2479:       MatSetBlockSizes(C,rbs,cbs);
2480:       MatSetType(C,((PetscObject)A)->type_name);
2481:       MatSeqAIJSetPreallocation_SeqAIJ(C,0,lens);
2482:     }
2483:     c = (Mat_SeqAIJ*)C->data;

2485:     /* loop over rows inserting into submatrix */
2486:     a_new = c->a;
2487:     j_new = c->j;
2488:     i_new = c->i;

2490:     for (i=0; i<nrows; i++) {
2491:       ii    = starts[i];
2492:       lensi = lens[i];
2493:       for (k=0; k<lensi; k++) {
2494:         *j_new++ = aj[ii+k] - first;
2495:       }
2496:       PetscArraycpy(a_new,a->a + starts[i],lensi);
2497:       a_new     += lensi;
2498:       i_new[i+1] = i_new[i] + lensi;
2499:       c->ilen[i] = lensi;
2500:     }
2501:     PetscFree2(lens,starts);
2502:   } else {
2503:     ISGetIndices(iscol,&icol);
2504:     PetscCalloc1(oldcols,&smap);
2505:     PetscMalloc1(1+nrows,&lens);
2506:     for (i=0; i<ncols; i++) {
2507: #if defined(PETSC_USE_DEBUG)
2508:       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);
2509: #endif
2510:       smap[icol[i]] = i+1;
2511:     }

2513:     /* determine lens of each row */
2514:     for (i=0; i<nrows; i++) {
2515:       kstart  = ai[irow[i]];
2516:       kend    = kstart + a->ilen[irow[i]];
2517:       lens[i] = 0;
2518:       for (k=kstart; k<kend; k++) {
2519:         if (smap[aj[k]]) {
2520:           lens[i]++;
2521:         }
2522:       }
2523:     }
2524:     /* Create and fill new matrix */
2525:     if (scall == MAT_REUSE_MATRIX) {
2526:       PetscBool equal;

2528:       c = (Mat_SeqAIJ*)((*B)->data);
2529:       if ((*B)->rmap->n  != nrows || (*B)->cmap->n != ncols) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size");
2530:       PetscArraycmp(c->ilen,lens,(*B)->rmap->n,&equal);
2531:       if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros");
2532:       PetscArrayzero(c->ilen,(*B)->rmap->n);
2533:       C    = *B;
2534:     } else {
2535:       PetscInt rbs,cbs;
2536:       MatCreate(PetscObjectComm((PetscObject)A),&C);
2537:       MatSetSizes(C,nrows,ncols,PETSC_DETERMINE,PETSC_DETERMINE);
2538:       ISGetBlockSize(isrow,&rbs);
2539:       ISGetBlockSize(iscol,&cbs);
2540:       MatSetBlockSizes(C,rbs,cbs);
2541:       MatSetType(C,((PetscObject)A)->type_name);
2542:       MatSeqAIJSetPreallocation_SeqAIJ(C,0,lens);
2543:     }
2544:     c = (Mat_SeqAIJ*)(C->data);
2545:     for (i=0; i<nrows; i++) {
2546:       row      = irow[i];
2547:       kstart   = ai[row];
2548:       kend     = kstart + a->ilen[row];
2549:       mat_i    = c->i[i];
2550:       mat_j    = c->j + mat_i;
2551:       mat_a    = c->a + mat_i;
2552:       mat_ilen = c->ilen + i;
2553:       for (k=kstart; k<kend; k++) {
2554:         if ((tcol=smap[a->j[k]])) {
2555:           *mat_j++ = tcol - 1;
2556:           *mat_a++ = a->a[k];
2557:           (*mat_ilen)++;

2559:         }
2560:       }
2561:     }
2562:     /* Free work space */
2563:     ISRestoreIndices(iscol,&icol);
2564:     PetscFree(smap);
2565:     PetscFree(lens);
2566:     /* sort */
2567:     for (i = 0; i < nrows; i++) {
2568:       PetscInt ilen;

2570:       mat_i = c->i[i];
2571:       mat_j = c->j + mat_i;
2572:       mat_a = c->a + mat_i;
2573:       ilen  = c->ilen[i];
2574:       PetscSortIntWithScalarArray(ilen,mat_j,mat_a);
2575:     }
2576:   }
2577: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
2578:   MatPinToCPU(C,A->pinnedtocpu);
2579: #endif
2580:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
2581:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

2583:   ISRestoreIndices(isrow,&irow);
2584:   *B   = C;
2585:   return(0);
2586: }

2588: PetscErrorCode  MatGetMultiProcBlock_SeqAIJ(Mat mat,MPI_Comm subComm,MatReuse scall,Mat *subMat)
2589: {
2591:   Mat            B;

2594:   if (scall == MAT_INITIAL_MATRIX) {
2595:     MatCreate(subComm,&B);
2596:     MatSetSizes(B,mat->rmap->n,mat->cmap->n,mat->rmap->n,mat->cmap->n);
2597:     MatSetBlockSizesFromMats(B,mat,mat);
2598:     MatSetType(B,MATSEQAIJ);
2599:     MatDuplicateNoCreate_SeqAIJ(B,mat,MAT_COPY_VALUES,PETSC_TRUE);
2600:     *subMat = B;
2601:   } else {
2602:     MatCopy_SeqAIJ(mat,*subMat,SAME_NONZERO_PATTERN);
2603:   }
2604:   return(0);
2605: }

2607: PetscErrorCode MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,const MatFactorInfo *info)
2608: {
2609:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)inA->data;
2611:   Mat            outA;
2612:   PetscBool      row_identity,col_identity;

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

2617:   ISIdentity(row,&row_identity);
2618:   ISIdentity(col,&col_identity);

2620:   outA             = inA;
2621:   outA->factortype = MAT_FACTOR_LU;
2622:   PetscFree(inA->solvertype);
2623:   PetscStrallocpy(MATSOLVERPETSC,&inA->solvertype);

2625:   PetscObjectReference((PetscObject)row);
2626:   ISDestroy(&a->row);

2628:   a->row = row;

2630:   PetscObjectReference((PetscObject)col);
2631:   ISDestroy(&a->col);

2633:   a->col = col;

2635:   /* Create the inverse permutation so that it can be used in MatLUFactorNumeric() */
2636:   ISDestroy(&a->icol);
2637:   ISInvertPermutation(col,PETSC_DECIDE,&a->icol);
2638:   PetscLogObjectParent((PetscObject)inA,(PetscObject)a->icol);

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

2645:   MatMarkDiagonal_SeqAIJ(inA);
2646:   if (row_identity && col_identity) {
2647:     MatLUFactorNumeric_SeqAIJ_inplace(outA,inA,info);
2648:   } else {
2649:     MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(outA,inA,info);
2650:   }
2651:   return(0);
2652: }

2654: PetscErrorCode MatScale_SeqAIJ(Mat inA,PetscScalar alpha)
2655: {
2656:   Mat_SeqAIJ     *a     = (Mat_SeqAIJ*)inA->data;
2657:   PetscScalar    oalpha = alpha;
2659:   PetscBLASInt   one = 1,bnz;

2662:   PetscBLASIntCast(a->nz,&bnz);
2663:   PetscStackCallBLAS("BLASscal",BLASscal_(&bnz,&oalpha,a->a,&one));
2664:   PetscLogFlops(a->nz);
2665:   MatSeqAIJInvalidateDiagonal(inA);
2666: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
2667:   if (inA->offloadmask != PETSC_OFFLOAD_UNALLOCATED) inA->offloadmask = PETSC_OFFLOAD_CPU;
2668: #endif
2669:   return(0);
2670: }

2672: PetscErrorCode MatDestroySubMatrix_Private(Mat_SubSppt *submatj)
2673: {
2675:   PetscInt       i;

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

2681:     for (i=0; i<submatj->nrqr; ++i) {
2682:       PetscFree(submatj->sbuf2[i]);
2683:     }
2684:     PetscFree3(submatj->sbuf2,submatj->req_size,submatj->req_source1);

2686:     if (submatj->rbuf1) {
2687:       PetscFree(submatj->rbuf1[0]);
2688:       PetscFree(submatj->rbuf1);
2689:     }

2691:     for (i=0; i<submatj->nrqs; ++i) {
2692:       PetscFree(submatj->rbuf3[i]);
2693:     }
2694:     PetscFree3(submatj->req_source2,submatj->rbuf2,submatj->rbuf3);
2695:     PetscFree(submatj->pa);
2696:   }

2698: #if defined(PETSC_USE_CTABLE)
2699:   PetscTableDestroy((PetscTable*)&submatj->rmap);
2700:   if (submatj->cmap_loc) {PetscFree(submatj->cmap_loc);}
2701:   PetscFree(submatj->rmap_loc);
2702: #else
2703:   PetscFree(submatj->rmap);
2704: #endif

2706:   if (!submatj->allcolumns) {
2707: #if defined(PETSC_USE_CTABLE)
2708:     PetscTableDestroy((PetscTable*)&submatj->cmap);
2709: #else
2710:     PetscFree(submatj->cmap);
2711: #endif
2712:   }
2713:   PetscFree(submatj->row2proc);

2715:   PetscFree(submatj);
2716:   return(0);
2717: }

2719: PetscErrorCode MatDestroySubMatrix_SeqAIJ(Mat C)
2720: {
2722:   Mat_SeqAIJ     *c = (Mat_SeqAIJ*)C->data;
2723:   Mat_SubSppt    *submatj = c->submatis1;

2726:   (*submatj->destroy)(C);
2727:   MatDestroySubMatrix_Private(submatj);
2728:   return(0);
2729: }

2731: PetscErrorCode MatDestroySubMatrices_SeqAIJ(PetscInt n,Mat *mat[])
2732: {
2734:   PetscInt       i;
2735:   Mat            C;
2736:   Mat_SeqAIJ     *c;
2737:   Mat_SubSppt    *submatj;

2740:   for (i=0; i<n; i++) {
2741:     C       = (*mat)[i];
2742:     c       = (Mat_SeqAIJ*)C->data;
2743:     submatj = c->submatis1;
2744:     if (submatj) {
2745:       if (--((PetscObject)C)->refct <= 0) {
2746:         (*submatj->destroy)(C);
2747:         MatDestroySubMatrix_Private(submatj);
2748:         PetscFree(C->defaultvectype);
2749:         PetscLayoutDestroy(&C->rmap);
2750:         PetscLayoutDestroy(&C->cmap);
2751:         PetscHeaderDestroy(&C);
2752:       }
2753:     } else {
2754:       MatDestroy(&C);
2755:     }
2756:   }

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

2761:   PetscFree(*mat);
2762:   return(0);
2763: }

2765: PetscErrorCode MatCreateSubMatrices_SeqAIJ(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
2766: {
2768:   PetscInt       i;

2771:   if (scall == MAT_INITIAL_MATRIX) {
2772:     PetscCalloc1(n+1,B);
2773:   }

2775:   for (i=0; i<n; i++) {
2776:     MatCreateSubMatrix_SeqAIJ(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);
2777:   }
2778:   return(0);
2779: }

2781: PetscErrorCode MatIncreaseOverlap_SeqAIJ(Mat A,PetscInt is_max,IS is[],PetscInt ov)
2782: {
2783:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
2785:   PetscInt       row,i,j,k,l,m,n,*nidx,isz,val;
2786:   const PetscInt *idx;
2787:   PetscInt       start,end,*ai,*aj;
2788:   PetscBT        table;

2791:   m  = A->rmap->n;
2792:   ai = a->i;
2793:   aj = a->j;

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

2797:   PetscMalloc1(m+1,&nidx);
2798:   PetscBTCreate(m,&table);

2800:   for (i=0; i<is_max; i++) {
2801:     /* Initialize the two local arrays */
2802:     isz  = 0;
2803:     PetscBTMemzero(m,table);

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

2809:     /* Enter these into the temp arrays. I.e., mark table[row], enter row into new index */
2810:     for (j=0; j<n; ++j) {
2811:       if (!PetscBTLookupSet(table,idx[j])) nidx[isz++] = idx[j];
2812:     }
2813:     ISRestoreIndices(is[i],&idx);
2814:     ISDestroy(&is[i]);

2816:     k = 0;
2817:     for (j=0; j<ov; j++) { /* for each overlap */
2818:       n = isz;
2819:       for (; k<n; k++) { /* do only those rows in nidx[k], which are not done yet */
2820:         row   = nidx[k];
2821:         start = ai[row];
2822:         end   = ai[row+1];
2823:         for (l = start; l<end; l++) {
2824:           val = aj[l];
2825:           if (!PetscBTLookupSet(table,val)) nidx[isz++] = val;
2826:         }
2827:       }
2828:     }
2829:     ISCreateGeneral(PETSC_COMM_SELF,isz,nidx,PETSC_COPY_VALUES,(is+i));
2830:   }
2831:   PetscBTDestroy(&table);
2832:   PetscFree(nidx);
2833:   return(0);
2834: }

2836: /* -------------------------------------------------------------- */
2837: PetscErrorCode MatPermute_SeqAIJ(Mat A,IS rowp,IS colp,Mat *B)
2838: {
2839:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
2841:   PetscInt       i,nz = 0,m = A->rmap->n,n = A->cmap->n;
2842:   const PetscInt *row,*col;
2843:   PetscInt       *cnew,j,*lens;
2844:   IS             icolp,irowp;
2845:   PetscInt       *cwork = NULL;
2846:   PetscScalar    *vwork = NULL;

2849:   ISInvertPermutation(rowp,PETSC_DECIDE,&irowp);
2850:   ISGetIndices(irowp,&row);
2851:   ISInvertPermutation(colp,PETSC_DECIDE,&icolp);
2852:   ISGetIndices(icolp,&col);

2854:   /* determine lengths of permuted rows */
2855:   PetscMalloc1(m+1,&lens);
2856:   for (i=0; i<m; i++) lens[row[i]] = a->i[i+1] - a->i[i];
2857:   MatCreate(PetscObjectComm((PetscObject)A),B);
2858:   MatSetSizes(*B,m,n,m,n);
2859:   MatSetBlockSizesFromMats(*B,A,A);
2860:   MatSetType(*B,((PetscObject)A)->type_name);
2861:   MatSeqAIJSetPreallocation_SeqAIJ(*B,0,lens);
2862:   PetscFree(lens);

2864:   PetscMalloc1(n,&cnew);
2865:   for (i=0; i<m; i++) {
2866:     MatGetRow_SeqAIJ(A,i,&nz,&cwork,&vwork);
2867:     for (j=0; j<nz; j++) cnew[j] = col[cwork[j]];
2868:     MatSetValues_SeqAIJ(*B,1,&row[i],nz,cnew,vwork,INSERT_VALUES);
2869:     MatRestoreRow_SeqAIJ(A,i,&nz,&cwork,&vwork);
2870:   }
2871:   PetscFree(cnew);

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

2875: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
2876:   MatPinToCPU(*B,A->pinnedtocpu);
2877: #endif
2878:   MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);
2879:   MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);
2880:   ISRestoreIndices(irowp,&row);
2881:   ISRestoreIndices(icolp,&col);
2882:   ISDestroy(&irowp);
2883:   ISDestroy(&icolp);
2884:   if (rowp == colp) {
2885:     if (A->symmetric) {
2886:       MatSetOption(*B,MAT_SYMMETRIC,PETSC_TRUE);
2887:     }
2888:     if (A->hermitian) {
2889:       MatSetOption(*B,MAT_HERMITIAN,PETSC_TRUE);
2890:     }
2891:   }
2892:   return(0);
2893: }

2895: PetscErrorCode MatCopy_SeqAIJ(Mat A,Mat B,MatStructure str)
2896: {

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

2905:     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");
2906:     PetscArraycpy(b->a,a->a,a->i[A->rmap->n]);
2907:     PetscObjectStateIncrease((PetscObject)B);
2908:   } else {
2909:     MatCopy_Basic(A,B,str);
2910:   }
2911:   return(0);
2912: }

2914: PetscErrorCode MatSetUp_SeqAIJ(Mat A)
2915: {

2919:    MatSeqAIJSetPreallocation_SeqAIJ(A,PETSC_DEFAULT,0);
2920:   return(0);
2921: }

2923: PETSC_INTERN PetscErrorCode MatSeqAIJGetArray_SeqAIJ(Mat A,PetscScalar *array[])
2924: {
2925:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;

2928:   *array = a->a;
2929:   return(0);
2930: }

2932: PETSC_INTERN PetscErrorCode MatSeqAIJRestoreArray_SeqAIJ(Mat A,PetscScalar *array[])
2933: {
2935:   *array = NULL;
2936:   return(0);
2937: }

2939: /*
2940:    Computes the number of nonzeros per row needed for preallocation when X and Y
2941:    have different nonzero structure.
2942: */
2943: PetscErrorCode MatAXPYGetPreallocation_SeqX_private(PetscInt m,const PetscInt *xi,const PetscInt *xj,const PetscInt *yi,const PetscInt *yj,PetscInt *nnz)
2944: {
2945:   PetscInt       i,j,k,nzx,nzy;

2948:   /* Set the number of nonzeros in the new matrix */
2949:   for (i=0; i<m; i++) {
2950:     const PetscInt *xjj = xj+xi[i],*yjj = yj+yi[i];
2951:     nzx = xi[i+1] - xi[i];
2952:     nzy = yi[i+1] - yi[i];
2953:     nnz[i] = 0;
2954:     for (j=0,k=0; j<nzx; j++) {                   /* Point in X */
2955:       for (; k<nzy && yjj[k]<xjj[j]; k++) nnz[i]++; /* Catch up to X */
2956:       if (k<nzy && yjj[k]==xjj[j]) k++;             /* Skip duplicate */
2957:       nnz[i]++;
2958:     }
2959:     for (; k<nzy; k++) nnz[i]++;
2960:   }
2961:   return(0);
2962: }

2964: PetscErrorCode MatAXPYGetPreallocation_SeqAIJ(Mat Y,Mat X,PetscInt *nnz)
2965: {
2966:   PetscInt       m = Y->rmap->N;
2967:   Mat_SeqAIJ     *x = (Mat_SeqAIJ*)X->data;
2968:   Mat_SeqAIJ     *y = (Mat_SeqAIJ*)Y->data;

2972:   /* Set the number of nonzeros in the new matrix */
2973:   MatAXPYGetPreallocation_SeqX_private(m,x->i,x->j,y->i,y->j,nnz);
2974:   return(0);
2975: }

2977: PetscErrorCode MatAXPY_SeqAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str)
2978: {
2980:   Mat_SeqAIJ     *x = (Mat_SeqAIJ*)X->data,*y = (Mat_SeqAIJ*)Y->data;
2981:   PetscBLASInt   one=1,bnz;

2984:   PetscBLASIntCast(x->nz,&bnz);
2985:   if (str == SAME_NONZERO_PATTERN) {
2986:     PetscScalar alpha = a;
2987:     PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&bnz,&alpha,x->a,&one,y->a,&one));
2988:     MatSeqAIJInvalidateDiagonal(Y);
2989:     PetscObjectStateIncrease((PetscObject)Y);
2990:     /* the MatAXPY_Basic* subroutines calls MatAssembly, so the matrix on the GPU
2991:        will be updated */
2992: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
2993:     if (Y->offloadmask != PETSC_OFFLOAD_UNALLOCATED) {
2994:       Y->offloadmask = PETSC_OFFLOAD_CPU;
2995:     }
2996: #endif
2997:   } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
2998:     MatAXPY_Basic(Y,a,X,str);
2999:   } else {
3000:     Mat      B;
3001:     PetscInt *nnz;
3002:     PetscMalloc1(Y->rmap->N,&nnz);
3003:     MatCreate(PetscObjectComm((PetscObject)Y),&B);
3004:     PetscObjectSetName((PetscObject)B,((PetscObject)Y)->name);
3005:     MatSetSizes(B,Y->rmap->n,Y->cmap->n,Y->rmap->N,Y->cmap->N);
3006:     MatSetBlockSizesFromMats(B,Y,Y);
3007:     MatSetType(B,(MatType) ((PetscObject)Y)->type_name);
3008:     MatAXPYGetPreallocation_SeqAIJ(Y,X,nnz);
3009:     MatSeqAIJSetPreallocation(B,0,nnz);
3010:     MatAXPY_BasicWithPreallocation(B,Y,a,X,str);
3011:     MatHeaderReplace(Y,&B);
3012:     PetscFree(nnz);
3013:   }
3014:   return(0);
3015: }

3017: PetscErrorCode  MatConjugate_SeqAIJ(Mat mat)
3018: {
3019: #if defined(PETSC_USE_COMPLEX)
3020:   Mat_SeqAIJ  *aij = (Mat_SeqAIJ*)mat->data;
3021:   PetscInt    i,nz;
3022:   PetscScalar *a;

3025:   nz = aij->nz;
3026:   a  = aij->a;
3027:   for (i=0; i<nz; i++) a[i] = PetscConj(a[i]);
3028: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
3029:   if (mat->offloadmask != PETSC_OFFLOAD_UNALLOCATED) mat->offloadmask = PETSC_OFFLOAD_CPU;
3030: #endif
3031: #else
3033: #endif
3034:   return(0);
3035: }

3037: PetscErrorCode MatGetRowMaxAbs_SeqAIJ(Mat A,Vec v,PetscInt idx[])
3038: {
3039:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
3041:   PetscInt       i,j,m = A->rmap->n,*ai,*aj,ncols,n;
3042:   PetscReal      atmp;
3043:   PetscScalar    *x;
3044:   MatScalar      *aa;

3047:   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
3048:   aa = a->a;
3049:   ai = a->i;
3050:   aj = a->j;

3052:   VecSet(v,0.0);
3053:   VecGetArray(v,&x);
3054:   VecGetLocalSize(v,&n);
3055:   if (n != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
3056:   for (i=0; i<m; i++) {
3057:     ncols = ai[1] - ai[0]; ai++;
3058:     x[i]  = 0.0;
3059:     for (j=0; j<ncols; j++) {
3060:       atmp = PetscAbsScalar(*aa);
3061:       if (PetscAbsScalar(x[i]) < atmp) {x[i] = atmp; if (idx) idx[i] = *aj;}
3062:       aa++; aj++;
3063:     }
3064:   }
3065:   VecRestoreArray(v,&x);
3066:   return(0);
3067: }

3069: PetscErrorCode MatGetRowMax_SeqAIJ(Mat A,Vec v,PetscInt idx[])
3070: {
3071:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
3073:   PetscInt       i,j,m = A->rmap->n,*ai,*aj,ncols,n;
3074:   PetscScalar    *x;
3075:   MatScalar      *aa;

3078:   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
3079:   aa = a->a;
3080:   ai = a->i;
3081:   aj = a->j;

3083:   VecSet(v,0.0);
3084:   VecGetArray(v,&x);
3085:   VecGetLocalSize(v,&n);
3086:   if (n != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
3087:   for (i=0; i<m; i++) {
3088:     ncols = ai[1] - ai[0]; ai++;
3089:     if (ncols == A->cmap->n) { /* row is dense */
3090:       x[i] = *aa; if (idx) idx[i] = 0;
3091:     } else {  /* row is sparse so already KNOW maximum is 0.0 or higher */
3092:       x[i] = 0.0;
3093:       if (idx) {
3094:         idx[i] = 0; /* in case ncols is zero */
3095:         for (j=0;j<ncols;j++) { /* find first implicit 0.0 in the row */
3096:           if (aj[j] > j) {
3097:             idx[i] = j;
3098:             break;
3099:           }
3100:         }
3101:       }
3102:     }
3103:     for (j=0; j<ncols; j++) {
3104:       if (PetscRealPart(x[i]) < PetscRealPart(*aa)) {x[i] = *aa; if (idx) idx[i] = *aj;}
3105:       aa++; aj++;
3106:     }
3107:   }
3108:   VecRestoreArray(v,&x);
3109:   return(0);
3110: }

3112: PetscErrorCode MatGetRowMinAbs_SeqAIJ(Mat A,Vec v,PetscInt idx[])
3113: {
3114:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
3116:   PetscInt       i,j,m = A->rmap->n,*ai,*aj,ncols,n;
3117:   PetscReal      atmp;
3118:   PetscScalar    *x;
3119:   MatScalar      *aa;

3122:   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
3123:   aa = a->a;
3124:   ai = a->i;
3125:   aj = a->j;

3127:   VecSet(v,0.0);
3128:   VecGetArray(v,&x);
3129:   VecGetLocalSize(v,&n);
3130:   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);
3131:   for (i=0; i<m; i++) {
3132:     ncols = ai[1] - ai[0]; ai++;
3133:     if (ncols) {
3134:       /* Get first nonzero */
3135:       for (j = 0; j < ncols; j++) {
3136:         atmp = PetscAbsScalar(aa[j]);
3137:         if (atmp > 1.0e-12) {
3138:           x[i] = atmp;
3139:           if (idx) idx[i] = aj[j];
3140:           break;
3141:         }
3142:       }
3143:       if (j == ncols) {x[i] = PetscAbsScalar(*aa); if (idx) idx[i] = *aj;}
3144:     } else {
3145:       x[i] = 0.0; if (idx) idx[i] = 0;
3146:     }
3147:     for (j = 0; j < ncols; j++) {
3148:       atmp = PetscAbsScalar(*aa);
3149:       if (atmp > 1.0e-12 && PetscAbsScalar(x[i]) > atmp) {x[i] = atmp; if (idx) idx[i] = *aj;}
3150:       aa++; aj++;
3151:     }
3152:   }
3153:   VecRestoreArray(v,&x);
3154:   return(0);
3155: }

3157: PetscErrorCode MatGetRowMin_SeqAIJ(Mat A,Vec v,PetscInt idx[])
3158: {
3159:   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
3160:   PetscErrorCode  ierr;
3161:   PetscInt        i,j,m = A->rmap->n,ncols,n;
3162:   const PetscInt  *ai,*aj;
3163:   PetscScalar     *x;
3164:   const MatScalar *aa;

3167:   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
3168:   aa = a->a;
3169:   ai = a->i;
3170:   aj = a->j;

3172:   VecSet(v,0.0);
3173:   VecGetArray(v,&x);
3174:   VecGetLocalSize(v,&n);
3175:   if (n != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
3176:   for (i=0; i<m; i++) {
3177:     ncols = ai[1] - ai[0]; ai++;
3178:     if (ncols == A->cmap->n) { /* row is dense */
3179:       x[i] = *aa; if (idx) idx[i] = 0;
3180:     } else {  /* row is sparse so already KNOW minimum is 0.0 or lower */
3181:       x[i] = 0.0;
3182:       if (idx) {   /* find first implicit 0.0 in the row */
3183:         idx[i] = 0; /* in case ncols is zero */
3184:         for (j=0; j<ncols; j++) {
3185:           if (aj[j] > j) {
3186:             idx[i] = j;
3187:             break;
3188:           }
3189:         }
3190:       }
3191:     }
3192:     for (j=0; j<ncols; j++) {
3193:       if (PetscRealPart(x[i]) > PetscRealPart(*aa)) {x[i] = *aa; if (idx) idx[i] = *aj;}
3194:       aa++; aj++;
3195:     }
3196:   }
3197:   VecRestoreArray(v,&x);
3198:   return(0);
3199: }

3201: PetscErrorCode MatInvertBlockDiagonal_SeqAIJ(Mat A,const PetscScalar **values)
3202: {
3203:   Mat_SeqAIJ      *a = (Mat_SeqAIJ*) A->data;
3204:   PetscErrorCode  ierr;
3205:   PetscInt        i,bs = PetscAbs(A->rmap->bs),mbs = A->rmap->n/bs,ipvt[5],bs2 = bs*bs,*v_pivots,ij[7],*IJ,j;
3206:   MatScalar       *diag,work[25],*v_work;
3207:   const PetscReal shift = 0.0;
3208:   PetscBool       allowzeropivot,zeropivotdetected=PETSC_FALSE;

3211:   allowzeropivot = PetscNot(A->erroriffailure);
3212:   if (a->ibdiagvalid) {
3213:     if (values) *values = a->ibdiag;
3214:     return(0);
3215:   }
3216:   MatMarkDiagonal_SeqAIJ(A);
3217:   if (!a->ibdiag) {
3218:     PetscMalloc1(bs2*mbs,&a->ibdiag);
3219:     PetscLogObjectMemory((PetscObject)A,bs2*mbs*sizeof(PetscScalar));
3220:   }
3221:   diag = a->ibdiag;
3222:   if (values) *values = a->ibdiag;
3223:   /* factor and invert each block */
3224:   switch (bs) {
3225:   case 1:
3226:     for (i=0; i<mbs; i++) {
3227:       MatGetValues(A,1,&i,1,&i,diag+i);
3228:       if (PetscAbsScalar(diag[i] + shift) < PETSC_MACHINE_EPSILON) {
3229:         if (allowzeropivot) {
3230:           A->factorerrortype             = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3231:           A->factorerror_zeropivot_value = PetscAbsScalar(diag[i]);
3232:           A->factorerror_zeropivot_row   = i;
3233:           PetscInfo3(A,"Zero pivot, row %D pivot %g tolerance %g\n",i,(double)PetscAbsScalar(diag[i]),(double)PETSC_MACHINE_EPSILON);
3234:         } 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);
3235:       }
3236:       diag[i] = (PetscScalar)1.0 / (diag[i] + shift);
3237:     }
3238:     break;
3239:   case 2:
3240:     for (i=0; i<mbs; i++) {
3241:       ij[0] = 2*i; ij[1] = 2*i + 1;
3242:       MatGetValues(A,2,ij,2,ij,diag);
3243:       PetscKernel_A_gets_inverse_A_2(diag,shift,allowzeropivot,&zeropivotdetected);
3244:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3245:       PetscKernel_A_gets_transpose_A_2(diag);
3246:       diag += 4;
3247:     }
3248:     break;
3249:   case 3:
3250:     for (i=0; i<mbs; i++) {
3251:       ij[0] = 3*i; ij[1] = 3*i + 1; ij[2] = 3*i + 2;
3252:       MatGetValues(A,3,ij,3,ij,diag);
3253:       PetscKernel_A_gets_inverse_A_3(diag,shift,allowzeropivot,&zeropivotdetected);
3254:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3255:       PetscKernel_A_gets_transpose_A_3(diag);
3256:       diag += 9;
3257:     }
3258:     break;
3259:   case 4:
3260:     for (i=0; i<mbs; i++) {
3261:       ij[0] = 4*i; ij[1] = 4*i + 1; ij[2] = 4*i + 2; ij[3] = 4*i + 3;
3262:       MatGetValues(A,4,ij,4,ij,diag);
3263:       PetscKernel_A_gets_inverse_A_4(diag,shift,allowzeropivot,&zeropivotdetected);
3264:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3265:       PetscKernel_A_gets_transpose_A_4(diag);
3266:       diag += 16;
3267:     }
3268:     break;
3269:   case 5:
3270:     for (i=0; i<mbs; i++) {
3271:       ij[0] = 5*i; ij[1] = 5*i + 1; ij[2] = 5*i + 2; ij[3] = 5*i + 3; ij[4] = 5*i + 4;
3272:       MatGetValues(A,5,ij,5,ij,diag);
3273:       PetscKernel_A_gets_inverse_A_5(diag,ipvt,work,shift,allowzeropivot,&zeropivotdetected);
3274:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3275:       PetscKernel_A_gets_transpose_A_5(diag);
3276:       diag += 25;
3277:     }
3278:     break;
3279:   case 6:
3280:     for (i=0; i<mbs; i++) {
3281:       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;
3282:       MatGetValues(A,6,ij,6,ij,diag);
3283:       PetscKernel_A_gets_inverse_A_6(diag,shift,allowzeropivot,&zeropivotdetected);
3284:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3285:       PetscKernel_A_gets_transpose_A_6(diag);
3286:       diag += 36;
3287:     }
3288:     break;
3289:   case 7:
3290:     for (i=0; i<mbs; i++) {
3291:       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;
3292:       MatGetValues(A,7,ij,7,ij,diag);
3293:       PetscKernel_A_gets_inverse_A_7(diag,shift,allowzeropivot,&zeropivotdetected);
3294:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3295:       PetscKernel_A_gets_transpose_A_7(diag);
3296:       diag += 49;
3297:     }
3298:     break;
3299:   default:
3300:     PetscMalloc3(bs,&v_work,bs,&v_pivots,bs,&IJ);
3301:     for (i=0; i<mbs; i++) {
3302:       for (j=0; j<bs; j++) {
3303:         IJ[j] = bs*i + j;
3304:       }
3305:       MatGetValues(A,bs,IJ,bs,IJ,diag);
3306:       PetscKernel_A_gets_inverse_A(bs,diag,v_pivots,v_work,allowzeropivot,&zeropivotdetected);
3307:       if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3308:       PetscKernel_A_gets_transpose_A_N(diag,bs);
3309:       diag += bs2;
3310:     }
3311:     PetscFree3(v_work,v_pivots,IJ);
3312:   }
3313:   a->ibdiagvalid = PETSC_TRUE;
3314:   return(0);
3315: }

3317: static PetscErrorCode  MatSetRandom_SeqAIJ(Mat x,PetscRandom rctx)
3318: {
3320:   Mat_SeqAIJ     *aij = (Mat_SeqAIJ*)x->data;
3321:   PetscScalar    a;
3322:   PetscInt       m,n,i,j,col;

3325:   if (!x->assembled) {
3326:     MatGetSize(x,&m,&n);
3327:     for (i=0; i<m; i++) {
3328:       for (j=0; j<aij->imax[i]; j++) {
3329:         PetscRandomGetValue(rctx,&a);
3330:         col  = (PetscInt)(n*PetscRealPart(a));
3331:         MatSetValues(x,1,&i,1,&col,&a,ADD_VALUES);
3332:       }
3333:     }
3334:   } else {
3335:     for (i=0; i<aij->nz; i++) {PetscRandomGetValue(rctx,aij->a+i);}
3336:   }
3337:   MatAssemblyBegin(x,MAT_FINAL_ASSEMBLY);
3338:   MatAssemblyEnd(x,MAT_FINAL_ASSEMBLY);
3339:   return(0);
3340: }

3342: /* Like MatSetRandom_SeqAIJ, but do not set values on columns in range of [low, high) */
3343: PetscErrorCode  MatSetRandomSkipColumnRange_SeqAIJ_Private(Mat x,PetscInt low,PetscInt high,PetscRandom rctx)
3344: {
3346:   Mat_SeqAIJ     *aij = (Mat_SeqAIJ*)x->data;
3347:   PetscScalar    a;
3348:   PetscInt       m,n,i,j,col,nskip;

3351:   nskip = high - low;
3352:   MatGetSize(x,&m,&n);
3353:   n    -= nskip; /* shrink number of columns where nonzeros can be set */
3354:   for (i=0; i<m; i++) {
3355:     for (j=0; j<aij->imax[i]; j++) {
3356:       PetscRandomGetValue(rctx,&a);
3357:       col  = (PetscInt)(n*PetscRealPart(a));
3358:       if (col >= low) col += nskip; /* shift col rightward to skip the hole */
3359:       MatSetValues(x,1,&i,1,&col,&a,ADD_VALUES);
3360:     }
3361:   }
3362:   MatAssemblyBegin(x,MAT_FINAL_ASSEMBLY);
3363:   MatAssemblyEnd(x,MAT_FINAL_ASSEMBLY);
3364:   return(0);
3365: }


3368: /* -------------------------------------------------------------------*/
3369: static struct _MatOps MatOps_Values = { MatSetValues_SeqAIJ,
3370:                                         MatGetRow_SeqAIJ,
3371:                                         MatRestoreRow_SeqAIJ,
3372:                                         MatMult_SeqAIJ,
3373:                                 /*  4*/ MatMultAdd_SeqAIJ,
3374:                                         MatMultTranspose_SeqAIJ,
3375:                                         MatMultTransposeAdd_SeqAIJ,
3376:                                         0,
3377:                                         0,
3378:                                         0,
3379:                                 /* 10*/ 0,
3380:                                         MatLUFactor_SeqAIJ,
3381:                                         0,
3382:                                         MatSOR_SeqAIJ,
3383:                                         MatTranspose_SeqAIJ,
3384:                                 /*1 5*/ MatGetInfo_SeqAIJ,
3385:                                         MatEqual_SeqAIJ,
3386:                                         MatGetDiagonal_SeqAIJ,
3387:                                         MatDiagonalScale_SeqAIJ,
3388:                                         MatNorm_SeqAIJ,
3389:                                 /* 20*/ 0,
3390:                                         MatAssemblyEnd_SeqAIJ,
3391:                                         MatSetOption_SeqAIJ,
3392:                                         MatZeroEntries_SeqAIJ,
3393:                                 /* 24*/ MatZeroRows_SeqAIJ,
3394:                                         0,
3395:                                         0,
3396:                                         0,
3397:                                         0,
3398:                                 /* 29*/ MatSetUp_SeqAIJ,
3399:                                         0,
3400:                                         0,
3401:                                         0,
3402:                                         0,
3403:                                 /* 34*/ MatDuplicate_SeqAIJ,
3404:                                         0,
3405:                                         0,
3406:                                         MatILUFactor_SeqAIJ,
3407:                                         0,
3408:                                 /* 39*/ MatAXPY_SeqAIJ,
3409:                                         MatCreateSubMatrices_SeqAIJ,
3410:                                         MatIncreaseOverlap_SeqAIJ,
3411:                                         MatGetValues_SeqAIJ,
3412:                                         MatCopy_SeqAIJ,
3413:                                 /* 44*/ MatGetRowMax_SeqAIJ,
3414:                                         MatScale_SeqAIJ,
3415:                                         MatShift_SeqAIJ,
3416:                                         MatDiagonalSet_SeqAIJ,
3417:                                         MatZeroRowsColumns_SeqAIJ,
3418:                                 /* 49*/ MatSetRandom_SeqAIJ,
3419:                                         MatGetRowIJ_SeqAIJ,
3420:                                         MatRestoreRowIJ_SeqAIJ,
3421:                                         MatGetColumnIJ_SeqAIJ,
3422:                                         MatRestoreColumnIJ_SeqAIJ,
3423:                                 /* 54*/ MatFDColoringCreate_SeqXAIJ,
3424:                                         0,
3425:                                         0,
3426:                                         MatPermute_SeqAIJ,
3427:                                         0,
3428:                                 /* 59*/ 0,
3429:                                         MatDestroy_SeqAIJ,
3430:                                         MatView_SeqAIJ,
3431:                                         0,
3432:                                         MatMatMatMult_SeqAIJ_SeqAIJ_SeqAIJ,
3433:                                 /* 64*/ MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ,
3434:                                         MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ,
3435:                                         0,
3436:                                         0,
3437:                                         0,
3438:                                 /* 69*/ MatGetRowMaxAbs_SeqAIJ,
3439:                                         MatGetRowMinAbs_SeqAIJ,
3440:                                         0,
3441:                                         0,
3442:                                         0,
3443:                                 /* 74*/ 0,
3444:                                         MatFDColoringApply_AIJ,
3445:                                         0,
3446:                                         0,
3447:                                         0,
3448:                                 /* 79*/ MatFindZeroDiagonals_SeqAIJ,
3449:                                         0,
3450:                                         0,
3451:                                         0,
3452:                                         MatLoad_SeqAIJ,
3453:                                 /* 84*/ MatIsSymmetric_SeqAIJ,
3454:                                         MatIsHermitian_SeqAIJ,
3455:                                         0,
3456:                                         0,
3457:                                         0,
3458:                                 /* 89*/ MatMatMult_SeqAIJ_SeqAIJ,
3459:                                         MatMatMultSymbolic_SeqAIJ_SeqAIJ,
3460:                                         MatMatMultNumeric_SeqAIJ_SeqAIJ,
3461:                                         MatPtAP_SeqAIJ_SeqAIJ,
3462:                                         MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy,
3463:                                 /* 94*/ MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy,
3464:                                         MatMatTransposeMult_SeqAIJ_SeqAIJ,
3465:                                         MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ,
3466:                                         MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ,
3467:                                         0,
3468:                                 /* 99*/ 0,
3469:                                         0,
3470:                                         0,
3471:                                         MatConjugate_SeqAIJ,
3472:                                         0,
3473:                                 /*104*/ MatSetValuesRow_SeqAIJ,
3474:                                         MatRealPart_SeqAIJ,
3475:                                         MatImaginaryPart_SeqAIJ,
3476:                                         0,
3477:                                         0,
3478:                                 /*109*/ MatMatSolve_SeqAIJ,
3479:                                         0,
3480:                                         MatGetRowMin_SeqAIJ,
3481:                                         0,
3482:                                         MatMissingDiagonal_SeqAIJ,
3483:                                 /*114*/ 0,
3484:                                         0,
3485:                                         0,
3486:                                         0,
3487:                                         0,
3488:                                 /*119*/ 0,
3489:                                         0,
3490:                                         0,
3491:                                         0,
3492:                                         MatGetMultiProcBlock_SeqAIJ,
3493:                                 /*124*/ MatFindNonzeroRows_SeqAIJ,
3494:                                         MatGetColumnNorms_SeqAIJ,
3495:                                         MatInvertBlockDiagonal_SeqAIJ,
3496:                                         MatInvertVariableBlockDiagonal_SeqAIJ,
3497:                                         0,
3498:                                 /*129*/ 0,
3499:                                         MatTransposeMatMult_SeqAIJ_SeqAIJ,
3500:                                         MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ,
3501:                                         MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ,
3502:                                         MatTransposeColoringCreate_SeqAIJ,
3503:                                 /*134*/ MatTransColoringApplySpToDen_SeqAIJ,
3504:                                         MatTransColoringApplyDenToSp_SeqAIJ,
3505:                                         MatRARt_SeqAIJ_SeqAIJ,
3506:                                         MatRARtSymbolic_SeqAIJ_SeqAIJ,
3507:                                         MatRARtNumeric_SeqAIJ_SeqAIJ,
3508:                                  /*139*/0,
3509:                                         0,
3510:                                         0,
3511:                                         MatFDColoringSetUp_SeqXAIJ,
3512:                                         MatFindOffBlockDiagonalEntries_SeqAIJ,
3513:                                  /*144*/MatCreateMPIMatConcatenateSeqMat_SeqAIJ,
3514:                                         MatDestroySubMatrices_SeqAIJ
3515: };

3517: PetscErrorCode  MatSeqAIJSetColumnIndices_SeqAIJ(Mat mat,PetscInt *indices)
3518: {
3519:   Mat_SeqAIJ *aij = (Mat_SeqAIJ*)mat->data;
3520:   PetscInt   i,nz,n;

3523:   nz = aij->maxnz;
3524:   n  = mat->rmap->n;
3525:   for (i=0; i<nz; i++) {
3526:     aij->j[i] = indices[i];
3527:   }
3528:   aij->nz = nz;
3529:   for (i=0; i<n; i++) {
3530:     aij->ilen[i] = aij->imax[i];
3531:   }
3532:   return(0);
3533: }

3535: /*
3536:  * When a sparse matrix has many zero columns, we should compact them out to save the space
3537:  * This happens in MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable()
3538:  * */
3539: PetscErrorCode  MatSeqAIJCompactOutExtraColumns_SeqAIJ(Mat mat, ISLocalToGlobalMapping *mapping)
3540: {
3541:   Mat_SeqAIJ *aij = (Mat_SeqAIJ*)mat->data;
3542:   PetscTable         gid1_lid1;
3543:   PetscTablePosition tpos;
3544:   PetscInt           gid,lid,i,j,ncols,ec;
3545:   PetscInt           *garray;
3546:   PetscErrorCode  ierr;

3551:   /* use a table */
3552:   PetscTableCreate(mat->rmap->n,mat->cmap->N+1,&gid1_lid1);
3553:   ec = 0;
3554:   for (i=0; i<mat->rmap->n; i++) {
3555:     ncols = aij->i[i+1] - aij->i[i];
3556:     for (j=0; j<ncols; j++) {
3557:       PetscInt data,gid1 = aij->j[aij->i[i] + j] + 1;
3558:       PetscTableFind(gid1_lid1,gid1,&data);
3559:       if (!data) {
3560:         /* one based table */
3561:         PetscTableAdd(gid1_lid1,gid1,++ec,INSERT_VALUES);
3562:       }
3563:     }
3564:   }
3565:   /* form array of columns we need */
3566:   PetscMalloc1(ec+1,&garray);
3567:   PetscTableGetHeadPosition(gid1_lid1,&tpos);
3568:   while (tpos) {
3569:     PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);
3570:     gid--;
3571:     lid--;
3572:     garray[lid] = gid;
3573:   }
3574:   PetscSortInt(ec,garray); /* sort, and rebuild */
3575:   PetscTableRemoveAll(gid1_lid1);
3576:   for (i=0; i<ec; i++) {
3577:     PetscTableAdd(gid1_lid1,garray[i]+1,i+1,INSERT_VALUES);
3578:   }
3579:   /* compact out the extra columns in B */
3580:   for (i=0; i<mat->rmap->n; i++) {
3581:         ncols = aij->i[i+1] - aij->i[i];
3582:     for (j=0; j<ncols; j++) {
3583:       PetscInt gid1 = aij->j[aij->i[i] + j] + 1;
3584:       PetscTableFind(gid1_lid1,gid1,&lid);
3585:       lid--;
3586:       aij->j[aij->i[i] + j] = lid;
3587:     }
3588:   }
3589:   PetscLayoutDestroy(&mat->cmap);
3590:   PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)mat),ec,ec,1,&mat->cmap);
3591:   PetscTableDestroy(&gid1_lid1);
3592:   ISLocalToGlobalMappingCreate(PETSC_COMM_SELF,mat->cmap->bs,mat->cmap->n,garray,PETSC_OWN_POINTER,mapping);
3593:   ISLocalToGlobalMappingSetType(*mapping,ISLOCALTOGLOBALMAPPINGHASH);
3594:   return(0);
3595: }

3597: /*@
3598:     MatSeqAIJSetColumnIndices - Set the column indices for all the rows
3599:        in the matrix.

3601:   Input Parameters:
3602: +  mat - the SeqAIJ matrix
3603: -  indices - the column indices

3605:   Level: advanced

3607:   Notes:
3608:     This can be called if you have precomputed the nonzero structure of the
3609:   matrix and want to provide it to the matrix object to improve the performance
3610:   of the MatSetValues() operation.

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

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

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

3619: @*/
3620: PetscErrorCode  MatSeqAIJSetColumnIndices(Mat mat,PetscInt *indices)
3621: {

3627:   PetscUseMethod(mat,"MatSeqAIJSetColumnIndices_C",(Mat,PetscInt*),(mat,indices));
3628:   return(0);
3629: }

3631: /* ----------------------------------------------------------------------------------------*/

3633: PetscErrorCode  MatStoreValues_SeqAIJ(Mat mat)
3634: {
3635:   Mat_SeqAIJ     *aij = (Mat_SeqAIJ*)mat->data;
3637:   size_t         nz = aij->i[mat->rmap->n];

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

3642:   /* allocate space for values if not already there */
3643:   if (!aij->saved_values) {
3644:     PetscMalloc1(nz+1,&aij->saved_values);
3645:     PetscLogObjectMemory((PetscObject)mat,(nz+1)*sizeof(PetscScalar));
3646:   }

3648:   /* copy values over */
3649:   PetscArraycpy(aij->saved_values,aij->a,nz);
3650:   return(0);
3651: }

3653: /*@
3654:     MatStoreValues - Stashes a copy of the matrix values; this allows, for
3655:        example, reuse of the linear part of a Jacobian, while recomputing the
3656:        nonlinear portion.

3658:    Collect on Mat

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

3663:   Level: advanced

3665:   Common Usage, with SNESSolve():
3666: $    Create Jacobian matrix
3667: $    Set linear terms into matrix
3668: $    Apply boundary conditions to matrix, at this time matrix must have
3669: $      final nonzero structure (i.e. setting the nonlinear terms and applying
3670: $      boundary conditions again will not change the nonzero structure
3671: $    MatSetOption(mat,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);
3672: $    MatStoreValues(mat);
3673: $    Call SNESSetJacobian() with matrix
3674: $    In your Jacobian routine
3675: $      MatRetrieveValues(mat);
3676: $      Set nonlinear terms in matrix

3678:   Common Usage without SNESSolve(), i.e. when you handle nonlinear solve yourself:
3679: $    // build linear portion of Jacobian
3680: $    MatSetOption(mat,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);
3681: $    MatStoreValues(mat);
3682: $    loop over nonlinear iterations
3683: $       MatRetrieveValues(mat);
3684: $       // call MatSetValues(mat,...) to set nonliner portion of Jacobian
3685: $       // call MatAssemblyBegin/End() on matrix
3686: $       Solve linear system with Jacobian
3687: $    endloop

3689:   Notes:
3690:     Matrix must already be assemblied before calling this routine
3691:     Must set the matrix option MatSetOption(mat,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE); before
3692:     calling this routine.

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

3697: .seealso: MatRetrieveValues()

3699: @*/
3700: PetscErrorCode  MatStoreValues(Mat mat)
3701: {

3706:   if (!mat->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
3707:   if (mat->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
3708:   PetscUseMethod(mat,"MatStoreValues_C",(Mat),(mat));
3709:   return(0);
3710: }

3712: PetscErrorCode  MatRetrieveValues_SeqAIJ(Mat mat)
3713: {
3714:   Mat_SeqAIJ     *aij = (Mat_SeqAIJ*)mat->data;
3716:   PetscInt       nz = aij->i[mat->rmap->n];

3719:   if (!aij->nonew) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first");
3720:   if (!aij->saved_values) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Must call MatStoreValues(A);first");
3721:   /* copy values over */
3722:   PetscArraycpy(aij->a,aij->saved_values,nz);
3723:   return(0);
3724: }

3726: /*@
3727:     MatRetrieveValues - Retrieves the copy of the matrix values; this allows, for
3728:        example, reuse of the linear part of a Jacobian, while recomputing the
3729:        nonlinear portion.

3731:    Collect on Mat

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

3736:   Level: advanced

3738: .seealso: MatStoreValues()

3740: @*/
3741: PetscErrorCode  MatRetrieveValues(Mat mat)
3742: {

3747:   if (!mat->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
3748:   if (mat->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
3749:   PetscUseMethod(mat,"MatRetrieveValues_C",(Mat),(mat));
3750:   return(0);
3751: }


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

3762:    Collective

3764:    Input Parameters:
3765: +  comm - MPI communicator, set to PETSC_COMM_SELF
3766: .  m - number of rows
3767: .  n - number of columns
3768: .  nz - number of nonzeros per row (same for all rows)
3769: -  nnz - array containing the number of nonzeros in the various rows
3770:          (possibly different for each row) or NULL

3772:    Output Parameter:
3773: .  A - the matrix

3775:    It is recommended that one use the MatCreate(), MatSetType() and/or MatSetFromOptions(),
3776:    MatXXXXSetPreallocation() paradigm instead of this routine directly.
3777:    [MatXXXXSetPreallocation() is, for example, MatSeqAIJSetPreallocation]

3779:    Notes:
3780:    If nnz is given then nz is ignored

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

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

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

3797:    Options Database Keys:
3798: +  -mat_no_inode  - Do not use inodes
3799: -  -mat_inode_limit <limit> - Sets inode limit (max limit=5)

3801:    Level: intermediate

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

3805: @*/
3806: PetscErrorCode  MatCreateSeqAIJ(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A)
3807: {

3811:   MatCreate(comm,A);
3812:   MatSetSizes(*A,m,n,m,n);
3813:   MatSetType(*A,MATSEQAIJ);
3814:   MatSeqAIJSetPreallocation_SeqAIJ(*A,nz,nnz);
3815:   return(0);
3816: }

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

3824:    Collective

3826:    Input Parameters:
3827: +  B - The matrix
3828: .  nz - number of nonzeros per row (same for all rows)
3829: -  nnz - array containing the number of nonzeros in the various rows
3830:          (possibly different for each row) or NULL

3832:    Notes:
3833:      If nnz is given then nz is ignored

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

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

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

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

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

3858:    Options Database Keys:
3859: +  -mat_no_inode  - Do not use inodes
3860: -  -mat_inode_limit <limit> - Sets inode limit (max limit=5)

3862:    Level: intermediate

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

3866: @*/
3867: PetscErrorCode  MatSeqAIJSetPreallocation(Mat B,PetscInt nz,const PetscInt nnz[])
3868: {

3874:   PetscTryMethod(B,"MatSeqAIJSetPreallocation_C",(Mat,PetscInt,const PetscInt[]),(B,nz,nnz));
3875:   return(0);
3876: }

3878: PetscErrorCode  MatSeqAIJSetPreallocation_SeqAIJ(Mat B,PetscInt nz,const PetscInt *nnz)
3879: {
3880:   Mat_SeqAIJ     *b;
3881:   PetscBool      skipallocation = PETSC_FALSE,realalloc = PETSC_FALSE;
3883:   PetscInt       i;

3886:   if (nz >= 0 || nnz) realalloc = PETSC_TRUE;
3887:   if (nz == MAT_SKIP_ALLOCATION) {
3888:     skipallocation = PETSC_TRUE;
3889:     nz             = 0;
3890:   }
3891:   PetscLayoutSetUp(B->rmap);
3892:   PetscLayoutSetUp(B->cmap);

3894:   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
3895:   if (nz < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %D",nz);
3896: #if defined(PETSC_USE_DEBUG)
3897:   if (nnz) {
3898:     for (i=0; i<B->rmap->n; i++) {
3899:       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]);
3900:       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);
3901:     }
3902:   }
3903: #endif

3905:   B->preallocated = PETSC_TRUE;

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

3909:   if (!skipallocation) {
3910:     if (!b->imax) {
3911:       PetscMalloc1(B->rmap->n,&b->imax);
3912:       PetscLogObjectMemory((PetscObject)B,B->rmap->n*sizeof(PetscInt));
3913:     }
3914:     if (!b->ilen) {
3915:       /* b->ilen will count nonzeros in each row so far. */
3916:       PetscCalloc1(B->rmap->n,&b->ilen);
3917:       PetscLogObjectMemory((PetscObject)B,B->rmap->n*sizeof(PetscInt));
3918:     } else {
3919:       PetscMemzero(b->ilen,B->rmap->n*sizeof(PetscInt));
3920:     }
3921:     if (!b->ipre) {
3922:       PetscMalloc1(B->rmap->n,&b->ipre);
3923:       PetscLogObjectMemory((PetscObject)B,B->rmap->n*sizeof(PetscInt));
3924:     }
3925:     if (!nnz) {
3926:       if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 10;
3927:       else if (nz < 0) nz = 1;
3928:       nz = PetscMin(nz,B->cmap->n);
3929:       for (i=0; i<B->rmap->n; i++) b->imax[i] = nz;
3930:       nz = nz*B->rmap->n;
3931:     } else {
3932:       PetscInt64 nz64 = 0;
3933:       for (i=0; i<B->rmap->n; i++) {b->imax[i] = nnz[i]; nz64 += nnz[i];}
3934:       PetscIntCast(nz64,&nz);
3935:     }

3937:     /* allocate the matrix space */
3938:     /* FIXME: should B's old memory be unlogged? */
3939:     MatSeqXAIJFreeAIJ(B,&b->a,&b->j,&b->i);
3940:     if (B->structure_only) {
3941:       PetscMalloc1(nz,&b->j);
3942:       PetscMalloc1(B->rmap->n+1,&b->i);
3943:       PetscLogObjectMemory((PetscObject)B,(B->rmap->n+1)*sizeof(PetscInt)+nz*sizeof(PetscInt));
3944:     } else {
3945:       PetscMalloc3(nz,&b->a,nz,&b->j,B->rmap->n+1,&b->i);
3946:       PetscLogObjectMemory((PetscObject)B,(B->rmap->n+1)*sizeof(PetscInt)+nz*(sizeof(PetscScalar)+sizeof(PetscInt)));
3947:     }
3948:     b->i[0] = 0;
3949:     for (i=1; i<B->rmap->n+1; i++) {
3950:       b->i[i] = b->i[i-1] + b->imax[i-1];
3951:     }
3952:     if (B->structure_only) {
3953:       b->singlemalloc = PETSC_FALSE;
3954:       b->free_a       = PETSC_FALSE;
3955:     } else {
3956:       b->singlemalloc = PETSC_TRUE;
3957:       b->free_a       = PETSC_TRUE;
3958:     }
3959:     b->free_ij      = PETSC_TRUE;
3960:   } else {
3961:     b->free_a  = PETSC_FALSE;
3962:     b->free_ij = PETSC_FALSE;
3963:   }

3965:   if (b->ipre && nnz != b->ipre  && b->imax) {
3966:     /* reserve user-requested sparsity */
3967:     PetscArraycpy(b->ipre,b->imax,B->rmap->n);
3968:   }


3971:   b->nz               = 0;
3972:   b->maxnz            = nz;
3973:   B->info.nz_unneeded = (double)b->maxnz;
3974:   if (realalloc) {
3975:     MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);
3976:   }
3977:   B->was_assembled = PETSC_FALSE;
3978:   B->assembled     = PETSC_FALSE;
3979:   return(0);
3980: }


3983: PetscErrorCode MatResetPreallocation_SeqAIJ(Mat A)
3984: {
3985:   Mat_SeqAIJ     *a;
3986:   PetscInt       i;


3992:   /* Check local size. If zero, then return */
3993:   if (!A->rmap->n) return(0);

3995:   a = (Mat_SeqAIJ*)A->data;
3996:   /* if no saved info, we error out */
3997:   if (!a->ipre) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"No saved preallocation info \n");

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

4001:   PetscArraycpy(a->imax,a->ipre,A->rmap->n);
4002:   PetscArrayzero(a->ilen,A->rmap->n);
4003:   a->i[0] = 0;
4004:   for (i=1; i<A->rmap->n+1; i++) {
4005:     a->i[i] = a->i[i-1] + a->imax[i-1];
4006:   }
4007:   A->preallocated     = PETSC_TRUE;
4008:   a->nz               = 0;
4009:   a->maxnz            = a->i[A->rmap->n];
4010:   A->info.nz_unneeded = (double)a->maxnz;
4011:   A->was_assembled    = PETSC_FALSE;
4012:   A->assembled        = PETSC_FALSE;
4013:   return(0);
4014: }

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

4019:    Input Parameters:
4020: +  B - the matrix
4021: .  i - the indices into j for the start of each row (starts with zero)
4022: .  j - the column indices for each row (starts with zero) these must be sorted for each row
4023: -  v - optional values in the matrix

4025:    Level: developer

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

4029: .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatSeqAIJSetPreallocation(), MatCreateSeqAIJ(), MATSEQAIJ
4030: @*/
4031: PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat B,const PetscInt i[],const PetscInt j[],const PetscScalar v[])
4032: {

4038:   PetscTryMethod(B,"MatSeqAIJSetPreallocationCSR_C",(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]),(B,i,j,v));
4039:   return(0);
4040: }

4042: PetscErrorCode  MatSeqAIJSetPreallocationCSR_SeqAIJ(Mat B,const PetscInt Ii[],const PetscInt J[],const PetscScalar v[])
4043: {
4044:   PetscInt       i;
4045:   PetscInt       m,n;
4046:   PetscInt       nz;
4047:   PetscInt       *nnz, nz_max = 0;

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

4053:   PetscLayoutSetUp(B->rmap);
4054:   PetscLayoutSetUp(B->cmap);

4056:   MatGetSize(B, &m, &n);
4057:   PetscMalloc1(m+1, &nnz);
4058:   for (i = 0; i < m; i++) {
4059:     nz     = Ii[i+1]- Ii[i];
4060:     nz_max = PetscMax(nz_max, nz);
4061:     if (nz < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Local row %D has a negative number of columns %D", i, nnz);
4062:     nnz[i] = nz;
4063:   }
4064:   MatSeqAIJSetPreallocation(B, 0, nnz);
4065:   PetscFree(nnz);

4067:   for (i = 0; i < m; i++) {
4068:     MatSetValues_SeqAIJ(B, 1, &i, Ii[i+1] - Ii[i], J+Ii[i], v ? v + Ii[i] : NULL, INSERT_VALUES);
4069:   }

4071:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
4072:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

4074:   MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
4075:   return(0);
4076: }

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

4081: /*
4082:     Computes (B'*A')' since computing B*A directly is untenable

4084:                n                       p                          p
4085:         (              )       (              )         (                  )
4086:       m (      A       )  *  n (       B      )   =   m (         C        )
4087:         (              )       (              )         (                  )

4089: */
4090: PetscErrorCode MatMatMultNumeric_SeqDense_SeqAIJ(Mat A,Mat B,Mat C)
4091: {
4092:   PetscErrorCode    ierr;
4093:   Mat_SeqDense      *sub_a = (Mat_SeqDense*)A->data;
4094:   Mat_SeqAIJ        *sub_b = (Mat_SeqAIJ*)B->data;
4095:   Mat_SeqDense      *sub_c = (Mat_SeqDense*)C->data;
4096:   PetscInt          i,n,m,q,p;
4097:   const PetscInt    *ii,*idx;
4098:   const PetscScalar *b,*a,*a_q;
4099:   PetscScalar       *c,*c_q;

4102:   m    = A->rmap->n;
4103:   n    = A->cmap->n;
4104:   p    = B->cmap->n;
4105:   a    = sub_a->v;
4106:   b    = sub_b->a;
4107:   c    = sub_c->v;
4108:   PetscArrayzero(c,m*p);

4110:   ii  = sub_b->i;
4111:   idx = sub_b->j;
4112:   for (i=0; i<n; i++) {
4113:     q = ii[i+1] - ii[i];
4114:     while (q-->0) {
4115:       c_q = c + m*(*idx);
4116:       a_q = a + m*i;
4117:       PetscKernelAXPY(c_q,*b,a_q,m);
4118:       idx++;
4119:       b++;
4120:     }
4121:   }
4122:   return(0);
4123: }

4125: PetscErrorCode MatMatMultSymbolic_SeqDense_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
4126: {
4128:   PetscInt       m=A->rmap->n,n=B->cmap->n;
4129:   Mat            Cmat;

4132:   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);
4133:   MatCreate(PetscObjectComm((PetscObject)A),&Cmat);
4134:   MatSetSizes(Cmat,m,n,m,n);
4135:   MatSetBlockSizesFromMats(Cmat,A,B);
4136:   MatSetType(Cmat,MATSEQDENSE);
4137:   MatSeqDenseSetPreallocation(Cmat,NULL);

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

4141:   *C = Cmat;
4142:   return(0);
4143: }

4145: /* ----------------------------------------------------------------*/
4146: PETSC_INTERN PetscErrorCode MatMatMult_SeqDense_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
4147: {

4151:   if (scall == MAT_INITIAL_MATRIX) {
4152:     PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);
4153:     MatMatMultSymbolic_SeqDense_SeqAIJ(A,B,fill,C);
4154:     PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);
4155:   }
4156:   PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);
4157:   MatMatMultNumeric_SeqDense_SeqAIJ(A,B,*C);
4158:   PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);
4159:   return(0);
4160: }


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

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

4170:    Level: beginner

4172:    Notes:
4173:     MatSetValues() may be called for this matrix type with a NULL argument for the numerical values,
4174:     in this case the values associated with the rows and columns one passes in are set to zero
4175:     in the matrix

4177:     MatSetOptions(,MAT_STRUCTURE_ONLY,PETSC_TRUE) may be called for this matrix type. In this no
4178:     space is allocated for the nonzero entries and any entries passed with MatSetValues() are ignored

4180:   Developer Notes:
4181:     It would be nice if all matrix formats supported passing NULL in for the numerical values

4183: .seealso: MatCreateSeqAIJ(), MatSetFromOptions(), MatSetType(), MatCreate(), MatType
4184: M*/

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

4189:    This matrix type is identical to MATSEQAIJ when constructed with a single process communicator,
4190:    and MATMPIAIJ otherwise.  As a result, for single process communicators,
4191:   MatSeqAIJSetPreallocation is supported, and similarly MatMPIAIJSetPreallocation() is supported
4192:   for communicators controlling multiple processes.  It is recommended that you call both of
4193:   the above preallocation routines for simplicity.

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

4198:   Developer Notes:
4199:     Subclasses include MATAIJCUSPARSE, MATAIJPERM, MATAIJSELL, MATAIJMKL, MATAIJCRL, and also automatically switches over to use inodes when
4200:    enough exist.

4202:   Level: beginner

4204: .seealso: MatCreateAIJ(), MatCreateSeqAIJ(), MATSEQAIJ,MATMPIAIJ
4205: M*/

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

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

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

4219:   Level: beginner

4221: .seealso: MatCreateMPIAIJCRL,MATSEQAIJCRL,MATMPIAIJCRL, MATSEQAIJCRL, MATMPIAIJCRL
4222: M*/

4224: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJCRL(Mat,MatType,MatReuse,Mat*);
4225: #if defined(PETSC_HAVE_ELEMENTAL)
4226: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_Elemental(Mat,MatType,MatReuse,Mat*);
4227: #endif
4228: #if defined(PETSC_HAVE_HYPRE)
4229: PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat A,MatType,MatReuse,Mat*);
4230: PETSC_INTERN PetscErrorCode MatMatMatMult_Transpose_AIJ_AIJ(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
4231: #endif
4232: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqDense(Mat,MatType,MatReuse,Mat*);

4234: PETSC_EXTERN PetscErrorCode MatConvert_SeqAIJ_SeqSELL(Mat,MatType,MatReuse,Mat*);
4235: PETSC_INTERN PetscErrorCode MatConvert_XAIJ_IS(Mat,MatType,MatReuse,Mat*);
4236: PETSC_INTERN PetscErrorCode MatPtAP_IS_XAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);

4238: /*@C
4239:    MatSeqAIJGetArray - gives read/write access to the array where the data for a MATSEQAIJ matrix is stored

4241:    Not Collective

4243:    Input Parameter:
4244: .  mat - a MATSEQAIJ matrix

4246:    Output Parameter:
4247: .   array - pointer to the data

4249:    Level: intermediate

4251: .seealso: MatSeqAIJRestoreArray(), MatSeqAIJGetArrayF90()
4252: @*/
4253: PetscErrorCode  MatSeqAIJGetArray(Mat A,PetscScalar **array)
4254: {

4258:   PetscUseMethod(A,"MatSeqAIJGetArray_C",(Mat,PetscScalar**),(A,array));
4259:   return(0);
4260: }

4262: /*@C
4263:    MatSeqAIJGetArrayRead - gives read-only access to the array where the data for a MATSEQAIJ matrix is stored

4265:    Not Collective

4267:    Input Parameter:
4268: .  mat - a MATSEQAIJ matrix

4270:    Output Parameter:
4271: .   array - pointer to the data

4273:    Level: intermediate

4275: .seealso: MatSeqAIJGetArray(), MatSeqAIJRestoreArrayRead()
4276: @*/
4277: PetscErrorCode  MatSeqAIJGetArrayRead(Mat A,const PetscScalar **array)
4278: {
4279: #if defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_VIENNACL)
4280:   PetscOffloadMask oval;
4281: #endif

4285: #if defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_VIENNACL)
4286:   oval = A->offloadmask;
4287: #endif
4288:   MatSeqAIJGetArray(A,(PetscScalar**)array);
4289: #if defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_VIENNACL)
4290:   if (oval == PETSC_OFFLOAD_GPU || oval == PETSC_OFFLOAD_BOTH) A->offloadmask = PETSC_OFFLOAD_BOTH;
4291: #endif
4292:   return(0);
4293: }

4295: /*@C
4296:    MatSeqAIJRestoreArrayRead - restore the read-only access array obtained from MatSeqAIJGetArrayRead

4298:    Not Collective

4300:    Input Parameter:
4301: .  mat - a MATSEQAIJ matrix

4303:    Output Parameter:
4304: .   array - pointer to the data

4306:    Level: intermediate

4308: .seealso: MatSeqAIJGetArray(), MatSeqAIJGetArrayRead()
4309: @*/
4310: PetscErrorCode  MatSeqAIJRestoreArrayRead(Mat A,const PetscScalar **array)
4311: {
4312: #if defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_VIENNACL)
4313:   PetscOffloadMask oval;
4314: #endif

4318: #if defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_VIENNACL)
4319:   oval = A->offloadmask;
4320: #endif
4321:   MatSeqAIJRestoreArray(A,(PetscScalar**)array);
4322: #if defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_VIENNACL)
4323:   A->offloadmask = oval;
4324: #endif
4325:   return(0);
4326: }

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

4331:    Not Collective

4333:    Input Parameter:
4334: .  mat - a MATSEQAIJ matrix

4336:    Output Parameter:
4337: .   nz - the maximum number of nonzeros in any row

4339:    Level: intermediate

4341: .seealso: MatSeqAIJRestoreArray(), MatSeqAIJGetArrayF90()
4342: @*/
4343: PetscErrorCode  MatSeqAIJGetMaxRowNonzeros(Mat A,PetscInt *nz)
4344: {
4345:   Mat_SeqAIJ     *aij = (Mat_SeqAIJ*)A->data;

4348:   *nz = aij->rmax;
4349:   return(0);
4350: }

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

4355:    Not Collective

4357:    Input Parameters:
4358: +  mat - a MATSEQAIJ matrix
4359: -  array - pointer to the data

4361:    Level: intermediate

4363: .seealso: MatSeqAIJGetArray(), MatSeqAIJRestoreArrayF90()
4364: @*/
4365: PetscErrorCode  MatSeqAIJRestoreArray(Mat A,PetscScalar **array)
4366: {

4370:   PetscUseMethod(A,"MatSeqAIJRestoreArray_C",(Mat,PetscScalar**),(A,array));
4371:   return(0);
4372: }

4374: #if defined(PETSC_HAVE_CUDA)
4375: PETSC_EXTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJCUSPARSE(Mat);
4376: #endif

4378: PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJ(Mat B)
4379: {
4380:   Mat_SeqAIJ     *b;
4382:   PetscMPIInt    size;

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

4388:   PetscNewLog(B,&b);

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

4392:   PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
4393:   if (B->sortedfull) B->ops->setvalues = MatSetValues_SeqAIJ_SortedFull;

4395:   b->row                = 0;
4396:   b->col                = 0;
4397:   b->icol               = 0;
4398:   b->reallocs           = 0;
4399:   b->ignorezeroentries  = PETSC_FALSE;
4400:   b->roworiented        = PETSC_TRUE;
4401:   b->nonew              = 0;
4402:   b->diag               = 0;
4403:   b->solve_work         = 0;
4404:   B->spptr              = 0;
4405:   b->saved_values       = 0;
4406:   b->idiag              = 0;
4407:   b->mdiag              = 0;
4408:   b->ssor_work          = 0;
4409:   b->omega              = 1.0;
4410:   b->fshift             = 0.0;
4411:   b->idiagvalid         = PETSC_FALSE;
4412:   b->ibdiagvalid        = PETSC_FALSE;
4413:   b->keepnonzeropattern = PETSC_FALSE;

4415:   PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);
4416:   PetscObjectComposeFunction((PetscObject)B,"MatSeqAIJGetArray_C",MatSeqAIJGetArray_SeqAIJ);
4417:   PetscObjectComposeFunction((PetscObject)B,"MatSeqAIJRestoreArray_C",MatSeqAIJRestoreArray_SeqAIJ);

4419: #if defined(PETSC_HAVE_MATLAB_ENGINE)
4420:   PetscObjectComposeFunction((PetscObject)B,"PetscMatlabEnginePut_C",MatlabEnginePut_SeqAIJ);
4421:   PetscObjectComposeFunction((PetscObject)B,"PetscMatlabEngineGet_C",MatlabEngineGet_SeqAIJ);
4422: #endif

4424:   PetscObjectComposeFunction((PetscObject)B,"MatSeqAIJSetColumnIndices_C",MatSeqAIJSetColumnIndices_SeqAIJ);
4425:   PetscObjectComposeFunction((PetscObject)B,"MatStoreValues_C",MatStoreValues_SeqAIJ);
4426:   PetscObjectComposeFunction((PetscObject)B,"MatRetrieveValues_C",MatRetrieveValues_SeqAIJ);
4427:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_seqsbaij_C",MatConvert_SeqAIJ_SeqSBAIJ);
4428:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_seqbaij_C",MatConvert_SeqAIJ_SeqBAIJ);
4429:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_seqaijperm_C",MatConvert_SeqAIJ_SeqAIJPERM);
4430:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_seqaijsell_C",MatConvert_SeqAIJ_SeqAIJSELL);
4431: #if defined(PETSC_HAVE_MKL_SPARSE)
4432:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_seqaijmkl_C",MatConvert_SeqAIJ_SeqAIJMKL);
4433: #endif
4434: #if defined(PETSC_HAVE_CUDA)
4435:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_seqaijcusparse_C",MatConvert_SeqAIJ_SeqAIJCUSPARSE);
4436: #endif
4437:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_seqaijcrl_C",MatConvert_SeqAIJ_SeqAIJCRL);
4438: #if defined(PETSC_HAVE_ELEMENTAL)
4439:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_elemental_C",MatConvert_SeqAIJ_Elemental);
4440: #endif
4441: #if defined(PETSC_HAVE_HYPRE)
4442:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_hypre_C",MatConvert_AIJ_HYPRE);
4443:   PetscObjectComposeFunction((PetscObject)B,"MatMatMatMult_transpose_seqaij_seqaij_C",MatMatMatMult_Transpose_AIJ_AIJ);
4444: #endif
4445:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_seqdense_C",MatConvert_SeqAIJ_SeqDense);
4446:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_seqsell_C",MatConvert_SeqAIJ_SeqSELL);
4447:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_is_C",MatConvert_XAIJ_IS);
4448:   PetscObjectComposeFunction((PetscObject)B,"MatIsTranspose_C",MatIsTranspose_SeqAIJ);
4449:   PetscObjectComposeFunction((PetscObject)B,"MatIsHermitianTranspose_C",MatIsTranspose_SeqAIJ);
4450:   PetscObjectComposeFunction((PetscObject)B,"MatSeqAIJSetPreallocation_C",MatSeqAIJSetPreallocation_SeqAIJ);
4451:   PetscObjectComposeFunction((PetscObject)B,"MatResetPreallocation_C",MatResetPreallocation_SeqAIJ);
4452:   PetscObjectComposeFunction((PetscObject)B,"MatSeqAIJSetPreallocationCSR_C",MatSeqAIJSetPreallocationCSR_SeqAIJ);
4453:   PetscObjectComposeFunction((PetscObject)B,"MatReorderForNonzeroDiagonal_C",MatReorderForNonzeroDiagonal_SeqAIJ);
4454:   PetscObjectComposeFunction((PetscObject)B,"MatMatMult_seqdense_seqaij_C",MatMatMult_SeqDense_SeqAIJ);
4455:   PetscObjectComposeFunction((PetscObject)B,"MatMatMultSymbolic_seqdense_seqaij_C",MatMatMultSymbolic_SeqDense_SeqAIJ);
4456:   PetscObjectComposeFunction((PetscObject)B,"MatMatMultNumeric_seqdense_seqaij_C",MatMatMultNumeric_SeqDense_SeqAIJ);
4457:   PetscObjectComposeFunction((PetscObject)B,"MatPtAP_is_seqaij_C",MatPtAP_IS_XAIJ);
4458:   MatCreate_SeqAIJ_Inode(B);
4459:   PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);
4460:   MatSeqAIJSetTypeFromOptions(B);  /* this allows changing the matrix subtype to say MATSEQAIJPERM */
4461:   return(0);
4462: }

4464: /*
4465:     Given a matrix generated with MatGetFactor() duplicates all the information in A into B
4466: */
4467: PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat C,Mat A,MatDuplicateOption cpvalues,PetscBool mallocmatspace)
4468: {
4469:   Mat_SeqAIJ     *c,*a = (Mat_SeqAIJ*)A->data;
4471:   PetscInt       m = A->rmap->n,i;

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

4476:   C->factortype = A->factortype;
4477:   c->row        = 0;
4478:   c->col        = 0;
4479:   c->icol       = 0;
4480:   c->reallocs   = 0;

4482:   C->assembled = PETSC_TRUE;

4484:   PetscLayoutReference(A->rmap,&C->rmap);
4485:   PetscLayoutReference(A->cmap,&C->cmap);

4487:   PetscMalloc1(m,&c->imax);
4488:   PetscMemcpy(c->imax,a->imax,m*sizeof(PetscInt));
4489:   PetscMalloc1(m,&c->ilen);
4490:   PetscMemcpy(c->ilen,a->ilen,m*sizeof(PetscInt));
4491:   PetscLogObjectMemory((PetscObject)C, 2*m*sizeof(PetscInt));

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

4498:     c->singlemalloc = PETSC_TRUE;

4500:     PetscArraycpy(c->i,a->i,m+1);
4501:     if (m > 0) {
4502:       PetscArraycpy(c->j,a->j,a->i[m]);
4503:       if (cpvalues == MAT_COPY_VALUES) {
4504:         PetscArraycpy(c->a,a->a,a->i[m]);
4505:       } else {
4506:         PetscArrayzero(c->a,a->i[m]);
4507:       }
4508:     }
4509:   }

4511:   c->ignorezeroentries = a->ignorezeroentries;
4512:   c->roworiented       = a->roworiented;
4513:   c->nonew             = a->nonew;
4514:   if (a->diag) {
4515:     PetscMalloc1(m+1,&c->diag);
4516:     PetscMemcpy(c->diag,a->diag,m*sizeof(PetscInt));
4517:     PetscLogObjectMemory((PetscObject)C,(m+1)*sizeof(PetscInt));
4518:   } else c->diag = NULL;

4520:   c->solve_work         = 0;
4521:   c->saved_values       = 0;
4522:   c->idiag              = 0;
4523:   c->ssor_work          = 0;
4524:   c->keepnonzeropattern = a->keepnonzeropattern;
4525:   c->free_a             = PETSC_TRUE;
4526:   c->free_ij            = PETSC_TRUE;

4528:   c->rmax         = a->rmax;
4529:   c->nz           = a->nz;
4530:   c->maxnz        = a->nz;       /* Since we allocate exactly the right amount */
4531:   C->preallocated = PETSC_TRUE;

4533:   c->compressedrow.use   = a->compressedrow.use;
4534:   c->compressedrow.nrows = a->compressedrow.nrows;
4535:   if (a->compressedrow.use) {
4536:     i    = a->compressedrow.nrows;
4537:     PetscMalloc2(i+1,&c->compressedrow.i,i,&c->compressedrow.rindex);
4538:     PetscArraycpy(c->compressedrow.i,a->compressedrow.i,i+1);
4539:     PetscArraycpy(c->compressedrow.rindex,a->compressedrow.rindex,i);
4540:   } else {
4541:     c->compressedrow.use    = PETSC_FALSE;
4542:     c->compressedrow.i      = NULL;
4543:     c->compressedrow.rindex = NULL;
4544:   }
4545:   c->nonzerorowcnt = a->nonzerorowcnt;
4546:   C->nonzerostate  = A->nonzerostate;

4548:   MatDuplicate_SeqAIJ_Inode(A,cpvalues,&C);
4549:   PetscFunctionListDuplicate(((PetscObject)A)->qlist,&((PetscObject)C)->qlist);
4550:   return(0);
4551: }

4553: PetscErrorCode MatDuplicate_SeqAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
4554: {

4558:   MatCreate(PetscObjectComm((PetscObject)A),B);
4559:   MatSetSizes(*B,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);
4560:   if (!(A->rmap->n % A->rmap->bs) && !(A->cmap->n % A->cmap->bs)) {
4561:     MatSetBlockSizesFromMats(*B,A,A);
4562:   }
4563:   MatSetType(*B,((PetscObject)A)->type_name);
4564:   MatDuplicateNoCreate_SeqAIJ(*B,A,cpvalues,PETSC_TRUE);
4565:   return(0);
4566: }

4568: PetscErrorCode MatLoad_SeqAIJ(Mat newMat, PetscViewer viewer)
4569: {
4570:   PetscBool      isbinary, ishdf5;

4576:   /* force binary viewer to load .info file if it has not yet done so */
4577:   PetscViewerSetUp(viewer);
4578:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
4579:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,  &ishdf5);
4580:   if (isbinary) {
4581:     MatLoad_SeqAIJ_Binary(newMat,viewer);
4582:   } else if (ishdf5) {
4583: #if defined(PETSC_HAVE_HDF5)
4584:     MatLoad_AIJ_HDF5(newMat,viewer);
4585: #else
4586:     SETERRQ(PetscObjectComm((PetscObject)newMat),PETSC_ERR_SUP,"HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
4587: #endif
4588:   } else {
4589:     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);
4590:   }
4591:   return(0);
4592: }

4594: PetscErrorCode MatLoad_SeqAIJ_Binary(Mat newMat, PetscViewer viewer)
4595: {
4596:   Mat_SeqAIJ     *a;
4598:   PetscInt       i,sum,nz,header[4],*rowlengths = 0,M,N,rows,cols;
4599:   int            fd;
4600:   PetscMPIInt    size;
4601:   MPI_Comm       comm;
4602:   PetscInt       bs = newMat->rmap->bs;

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

4609:   PetscOptionsBegin(comm,NULL,"Options for loading SEQAIJ matrix","Mat");
4610:   PetscOptionsInt("-matload_block_size","Set the blocksize used to store the matrix","MatLoad",bs,&bs,NULL);
4611:   PetscOptionsEnd();
4612:   if (bs < 0) bs = 1;
4613:   MatSetBlockSize(newMat,bs);

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

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

4622:   /* read in row lengths */
4623:   PetscMalloc1(M,&rowlengths);
4624:   PetscBinaryRead(fd,rowlengths,M,NULL,PETSC_INT);

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

4630:   /* set global size if not set already*/
4631:   if (newMat->rmap->n < 0 && newMat->rmap->N < 0 && newMat->cmap->n < 0 && newMat->cmap->N < 0) {
4632:     MatSetSizes(newMat,PETSC_DECIDE,PETSC_DECIDE,M,N);
4633:   } else {
4634:     /* if sizes and type are already set, check if the matrix  global sizes are correct */
4635:     MatGetSize(newMat,&rows,&cols);
4636:     if (rows < 0 && cols < 0) { /* user might provide local size instead of global size */
4637:       MatGetLocalSize(newMat,&rows,&cols);
4638:     }
4639:     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);
4640:   }
4641:   MatSeqAIJSetPreallocation_SeqAIJ(newMat,0,rowlengths);
4642:   a    = (Mat_SeqAIJ*)newMat->data;

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

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

4649:   /* set matrix "i" values */
4650:   a->i[0] = 0;
4651:   for (i=1; i<= M; i++) {
4652:     a->i[i]      = a->i[i-1] + rowlengths[i-1];
4653:     a->ilen[i-1] = rowlengths[i-1];
4654:   }
4655:   PetscFree(rowlengths);

4657:   MatAssemblyBegin(newMat,MAT_FINAL_ASSEMBLY);
4658:   MatAssemblyEnd(newMat,MAT_FINAL_ASSEMBLY);
4659:   return(0);
4660: }

4662: PetscErrorCode MatEqual_SeqAIJ(Mat A,Mat B,PetscBool * flg)
4663: {
4664:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data,*b = (Mat_SeqAIJ*)B->data;
4666: #if defined(PETSC_USE_COMPLEX)
4667:   PetscInt k;
4668: #endif

4671:   /* If the  matrix dimensions are not equal,or no of nonzeros */
4672:   if ((A->rmap->n != B->rmap->n) || (A->cmap->n != B->cmap->n) ||(a->nz != b->nz)) {
4673:     *flg = PETSC_FALSE;
4674:     return(0);
4675:   }

4677:   /* if the a->i are the same */
4678:   PetscArraycmp(a->i,b->i,A->rmap->n+1,flg);
4679:   if (!*flg) return(0);

4681:   /* if a->j are the same */
4682:   PetscArraycmp(a->j,b->j,a->nz,flg);
4683:   if (!*flg) return(0);

4685:   /* if a->a are the same */
4686: #if defined(PETSC_USE_COMPLEX)
4687:   for (k=0; k<a->nz; k++) {
4688:     if (PetscRealPart(a->a[k]) != PetscRealPart(b->a[k]) || PetscImaginaryPart(a->a[k]) != PetscImaginaryPart(b->a[k])) {
4689:       *flg = PETSC_FALSE;
4690:       return(0);
4691:     }
4692:   }
4693: #else
4694:   PetscArraycmp(a->a,b->a,a->nz,flg);
4695: #endif
4696:   return(0);
4697: }

4699: /*@
4700:      MatCreateSeqAIJWithArrays - Creates an sequential AIJ matrix using matrix elements (in CSR format)
4701:               provided by the user.

4703:       Collective

4705:    Input Parameters:
4706: +   comm - must be an MPI communicator of size 1
4707: .   m - number of rows
4708: .   n - number of columns
4709: .   i - row indices; that is i[0] = 0, i[row] = i[row-1] + number of elements in that row of the matrix
4710: .   j - column indices
4711: -   a - matrix values

4713:    Output Parameter:
4714: .   mat - the matrix

4716:    Level: intermediate

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

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

4724:        The i and j indices are 0 based

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

4730: $        1 0 0
4731: $        2 0 3
4732: $        4 5 6
4733: $
4734: $        i =  {0,1,3,6}  [size = nrow+1  = 3+1]
4735: $        j =  {0,0,2,0,1,2}  [size = 6]; values must be sorted for each row
4736: $        v =  {1,2,3,4,5,6}  [size = 6]


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

4741: @*/
4742: PetscErrorCode  MatCreateSeqAIJWithArrays(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt i[],PetscInt j[],PetscScalar a[],Mat *mat)
4743: {
4745:   PetscInt       ii;
4746:   Mat_SeqAIJ     *aij;
4747: #if defined(PETSC_USE_DEBUG)
4748:   PetscInt jj;
4749: #endif

4752:   if (m > 0 && i[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"i (row indices) must start with 0");
4753:   MatCreate(comm,mat);
4754:   MatSetSizes(*mat,m,n,m,n);
4755:   /* MatSetBlockSizes(*mat,,); */
4756:   MatSetType(*mat,MATSEQAIJ);
4757:   MatSeqAIJSetPreallocation_SeqAIJ(*mat,MAT_SKIP_ALLOCATION,0);
4758:   aij  = (Mat_SeqAIJ*)(*mat)->data;
4759:   PetscMalloc1(m,&aij->imax);
4760:   PetscMalloc1(m,&aij->ilen);

4762:   aij->i            = i;
4763:   aij->j            = j;
4764:   aij->a            = a;
4765:   aij->singlemalloc = PETSC_FALSE;
4766:   aij->nonew        = -1;             /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/
4767:   aij->free_a       = PETSC_FALSE;
4768:   aij->free_ij      = PETSC_FALSE;

4770:   for (ii=0; ii<m; ii++) {
4771:     aij->ilen[ii] = aij->imax[ii] = i[ii+1] - i[ii];
4772: #if defined(PETSC_USE_DEBUG)
4773:     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]);
4774:     for (jj=i[ii]+1; jj<i[ii+1]; jj++) {
4775:       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);
4776:       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);
4777:     }
4778: #endif
4779:   }
4780: #if defined(PETSC_USE_DEBUG)
4781:   for (ii=0; ii<aij->i[m]; ii++) {
4782:     if (j[ii] < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column index at location = %D index = %D",ii,j[ii]);
4783:     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]);
4784:   }
4785: #endif

4787:   MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
4788:   MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
4789:   return(0);
4790: }
4791: /*@C
4792:      MatCreateSeqAIJFromTriple - Creates an sequential AIJ matrix using matrix elements (in COO format)
4793:               provided by the user.

4795:       Collective

4797:    Input Parameters:
4798: +   comm - must be an MPI communicator of size 1
4799: .   m   - number of rows
4800: .   n   - number of columns
4801: .   i   - row indices
4802: .   j   - column indices
4803: .   a   - matrix values
4804: .   nz  - number of nonzeros
4805: -   idx - 0 or 1 based

4807:    Output Parameter:
4808: .   mat - the matrix

4810:    Level: intermediate

4812:    Notes:
4813:        The i and j indices are 0 based

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

4819:         1 0 0
4820:         2 0 3
4821:         4 5 6

4823:         i =  {0,1,1,2,2,2}
4824:         j =  {0,0,2,0,1,2}
4825:         v =  {1,2,3,4,5,6}


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

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


4838:   PetscCalloc1(m,&nnz);
4839:   for (ii = 0; ii < nz; ii++) {
4840:     nnz[i[ii] - !!idx] += 1;
4841:   }
4842:   MatCreate(comm,mat);
4843:   MatSetSizes(*mat,m,n,m,n);
4844:   MatSetType(*mat,MATSEQAIJ);
4845:   MatSeqAIJSetPreallocation_SeqAIJ(*mat,0,nnz);
4846:   for (ii = 0; ii < nz; ii++) {
4847:     if (idx) {
4848:       row = i[ii] - 1;
4849:       col = j[ii] - 1;
4850:     } else {
4851:       row = i[ii];
4852:       col = j[ii];
4853:     }
4854:     MatSetValues(*mat,one,&row,one,&col,&a[ii],ADD_VALUES);
4855:   }
4856:   MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
4857:   MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
4858:   PetscFree(nnz);
4859:   return(0);
4860: }

4862: PetscErrorCode MatSeqAIJInvalidateDiagonal(Mat A)
4863: {
4864:   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;

4868:   a->idiagvalid  = PETSC_FALSE;
4869:   a->ibdiagvalid = PETSC_FALSE;

4871:   MatSeqAIJInvalidateDiagonal_Inode(A);
4872:   return(0);
4873: }

4875: PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqAIJ(MPI_Comm comm,Mat inmat,PetscInt n,MatReuse scall,Mat *outmat)
4876: {
4878:   PetscMPIInt    size;

4881:   MPI_Comm_size(comm,&size);
4882:   if (size == 1) {
4883:     if (scall == MAT_INITIAL_MATRIX) {
4884:       MatDuplicate(inmat,MAT_COPY_VALUES,outmat);
4885:     } else {
4886:       MatCopy(inmat,*outmat,SAME_NONZERO_PATTERN);
4887:     }
4888:   } else {
4889:     MatCreateMPIMatConcatenateSeqMat_MPIAIJ(comm,inmat,n,scall,outmat);
4890:   }
4891:   return(0);
4892: }

4894: /*
4895:  Permute A into C's *local* index space using rowemb,colemb.
4896:  The embedding are supposed to be injections and the above implies that the range of rowemb is a subset
4897:  of [0,m), colemb is in [0,n).
4898:  If pattern == DIFFERENT_NONZERO_PATTERN, C is preallocated according to A.
4899:  */
4900: PetscErrorCode MatSetSeqMat_SeqAIJ(Mat C,IS rowemb,IS colemb,MatStructure pattern,Mat B)
4901: {
4902:   /* If making this function public, change the error returned in this function away from _PLIB. */
4904:   Mat_SeqAIJ     *Baij;
4905:   PetscBool      seqaij;
4906:   PetscInt       m,n,*nz,i,j,count;
4907:   PetscScalar    v;
4908:   const PetscInt *rowindices,*colindices;

4911:   if (!B) return(0);
4912:   /* Check to make sure the target matrix (and embeddings) are compatible with C and each other. */
4913:   PetscObjectBaseTypeCompare((PetscObject)B,MATSEQAIJ,&seqaij);
4914:   if (!seqaij) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Input matrix is of wrong type");
4915:   if (rowemb) {
4916:     ISGetLocalSize(rowemb,&m);
4917:     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);
4918:   } else {
4919:     if (C->rmap->n != B->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Input matrix is row-incompatible with the target matrix");
4920:   }
4921:   if (colemb) {
4922:     ISGetLocalSize(colemb,&n);
4923:     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);
4924:   } else {
4925:     if (C->cmap->n != B->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Input matrix is col-incompatible with the target matrix");
4926:   }

4928:   Baij = (Mat_SeqAIJ*)(B->data);
4929:   if (pattern == DIFFERENT_NONZERO_PATTERN) {
4930:     PetscMalloc1(B->rmap->n,&nz);
4931:     for (i=0; i<B->rmap->n; i++) {
4932:       nz[i] = Baij->i[i+1] - Baij->i[i];
4933:     }
4934:     MatSeqAIJSetPreallocation(C,0,nz);
4935:     PetscFree(nz);
4936:   }
4937:   if (pattern == SUBSET_NONZERO_PATTERN) {
4938:     MatZeroEntries(C);
4939:   }
4940:   count = 0;
4941:   rowindices = NULL;
4942:   colindices = NULL;
4943:   if (rowemb) {
4944:     ISGetIndices(rowemb,&rowindices);
4945:   }
4946:   if (colemb) {
4947:     ISGetIndices(colemb,&colindices);
4948:   }
4949:   for (i=0; i<B->rmap->n; i++) {
4950:     PetscInt row;
4951:     row = i;
4952:     if (rowindices) row = rowindices[i];
4953:     for (j=Baij->i[i]; j<Baij->i[i+1]; j++) {
4954:       PetscInt col;
4955:       col  = Baij->j[count];
4956:       if (colindices) col = colindices[col];
4957:       v    = Baij->a[count];
4958:       MatSetValues(C,1,&row,1,&col,&v,INSERT_VALUES);
4959:       ++count;
4960:     }
4961:   }
4962:   /* FIXME: set C's nonzerostate correctly. */
4963:   /* Assembly for C is necessary. */
4964:   C->preallocated = PETSC_TRUE;
4965:   C->assembled     = PETSC_TRUE;
4966:   C->was_assembled = PETSC_FALSE;
4967:   return(0);
4968: }

4970: PetscFunctionList MatSeqAIJList = NULL;

4972: /*@C
4973:    MatSeqAIJSetType - Converts a MATSEQAIJ matrix to a subtype

4975:    Collective on Mat

4977:    Input Parameters:
4978: +  mat      - the matrix object
4979: -  matype   - matrix type

4981:    Options Database Key:
4982: .  -mat_seqai_type  <method> - for example seqaijcrl


4985:   Level: intermediate

4987: .seealso: PCSetType(), VecSetType(), MatCreate(), MatType, Mat
4988: @*/
4989: PetscErrorCode  MatSeqAIJSetType(Mat mat, MatType matype)
4990: {
4991:   PetscErrorCode ierr,(*r)(Mat,MatType,MatReuse,Mat*);
4992:   PetscBool      sametype;

4996:   PetscObjectTypeCompare((PetscObject)mat,matype,&sametype);
4997:   if (sametype) return(0);

4999:    PetscFunctionListFind(MatSeqAIJList,matype,&r);
5000:   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown Mat type given: %s",matype);
5001:   (*r)(mat,matype,MAT_INPLACE_MATRIX,&mat);
5002:   return(0);
5003: }


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

5009:    Not Collective

5011:    Input Parameters:
5012: +  name - name of a new user-defined matrix type, for example MATSEQAIJCRL
5013: -  function - routine to convert to subtype

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


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

5022:    Level: advanced

5024: .seealso: MatSeqAIJRegisterAll()


5027:   Level: advanced
5028: @*/
5029: PetscErrorCode  MatSeqAIJRegister(const char sname[],PetscErrorCode (*function)(Mat,MatType,MatReuse,Mat *))
5030: {

5034:   MatInitializePackage();
5035:   PetscFunctionListAdd(&MatSeqAIJList,sname,function);
5036:   return(0);
5037: }

5039: PetscBool MatSeqAIJRegisterAllCalled = PETSC_FALSE;

5041: /*@C
5042:   MatSeqAIJRegisterAll - Registers all of the matrix subtypes of SeqAIJ

5044:   Not Collective

5046:   Level: advanced

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

5050: .seealso:  MatRegisterAll(), MatSeqAIJRegister()
5051: @*/
5052: PetscErrorCode  MatSeqAIJRegisterAll(void)
5053: {

5057:   if (MatSeqAIJRegisterAllCalled) return(0);
5058:   MatSeqAIJRegisterAllCalled = PETSC_TRUE;

5060:   MatSeqAIJRegister(MATSEQAIJCRL,      MatConvert_SeqAIJ_SeqAIJCRL);
5061:   MatSeqAIJRegister(MATSEQAIJPERM,     MatConvert_SeqAIJ_SeqAIJPERM);
5062:   MatSeqAIJRegister(MATSEQAIJSELL,     MatConvert_SeqAIJ_SeqAIJSELL);
5063: #if defined(PETSC_HAVE_MKL_SPARSE)
5064:   MatSeqAIJRegister(MATSEQAIJMKL,      MatConvert_SeqAIJ_SeqAIJMKL);
5065: #endif
5066: #if defined(PETSC_HAVE_VIENNACL) && defined(PETSC_HAVE_VIENNACL_NO_CUDA)
5067:   MatSeqAIJRegister(MATMPIAIJVIENNACL, MatConvert_SeqAIJ_SeqAIJViennaCL);
5068: #endif
5069:   return(0);
5070: }

5072: /*
5073:     Special version for direct calls from Fortran
5074: */
5075:  #include <petsc/private/fortranimpl.h>
5076: #if defined(PETSC_HAVE_FORTRAN_CAPS)
5077: #define matsetvaluesseqaij_ MATSETVALUESSEQAIJ
5078: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
5079: #define matsetvaluesseqaij_ matsetvaluesseqaij
5080: #endif

5082: /* Change these macros so can be used in void function */
5083: #undef CHKERRQ
5084: #define CHKERRQ(ierr) CHKERRABORT(PetscObjectComm((PetscObject)A),ierr)
5085: #undef SETERRQ2
5086: #define SETERRQ2(comm,ierr,b,c,d) CHKERRABORT(comm,ierr)
5087: #undef SETERRQ3
5088: #define SETERRQ3(comm,ierr,b,c,d,e) CHKERRABORT(comm,ierr)

5090: 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)
5091: {
5092:   Mat            A  = *AA;
5093:   PetscInt       m  = *mm, n = *nn;
5094:   InsertMode     is = *isis;
5095:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
5096:   PetscInt       *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N;
5097:   PetscInt       *imax,*ai,*ailen;
5099:   PetscInt       *aj,nonew = a->nonew,lastcol = -1;
5100:   MatScalar      *ap,value,*aa;
5101:   PetscBool      ignorezeroentries = a->ignorezeroentries;
5102:   PetscBool      roworiented       = a->roworiented;

5105:   MatCheckPreallocated(A,1);
5106:   imax  = a->imax;
5107:   ai    = a->i;
5108:   ailen = a->ilen;
5109:   aj    = a->j;
5110:   aa    = a->a;

5112:   for (k=0; k<m; k++) { /* loop over added rows */
5113:     row = im[k];
5114:     if (row < 0) continue;
5115: #if defined(PETSC_USE_DEBUG)
5116:     if (row >= A->rmap->n) SETERRABORT(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
5117: #endif
5118:     rp   = aj + ai[row]; ap = aa + ai[row];
5119:     rmax = imax[row]; nrow = ailen[row];
5120:     low  = 0;
5121:     high = nrow;
5122:     for (l=0; l<n; l++) { /* loop over added columns */
5123:       if (in[l] < 0) continue;
5124: #if defined(PETSC_USE_DEBUG)
5125:       if (in[l] >= A->cmap->n) SETERRABORT(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
5126: #endif
5127:       col = in[l];
5128:       if (roworiented) value = v[l + k*n];
5129:       else value = v[k + l*m];

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

5133:       if (col <= lastcol) low = 0;
5134:       else high = nrow;
5135:       lastcol = col;
5136:       while (high-low > 5) {
5137:         t = (low+high)/2;
5138:         if (rp[t] > col) high = t;
5139:         else             low  = t;
5140:       }
5141:       for (i=low; i<high; i++) {
5142:         if (rp[i] > col) break;
5143:         if (rp[i] == col) {
5144:           if (is == ADD_VALUES) ap[i] += value;
5145:           else                  ap[i] = value;
5146:           goto noinsert;
5147:         }
5148:       }
5149:       if (value == 0.0 && ignorezeroentries) goto noinsert;
5150:       if (nonew == 1) goto noinsert;
5151:       if (nonew == -1) SETERRABORT(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
5152:       MatSeqXAIJReallocateAIJ(A,A->rmap->n,1,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar);
5153:       N = nrow++ - 1; a->nz++; high++;
5154:       /* shift up all the later entries in this row */
5155:       for (ii=N; ii>=i; ii--) {
5156:         rp[ii+1] = rp[ii];
5157:         ap[ii+1] = ap[ii];
5158:       }
5159:       rp[i] = col;
5160:       ap[i] = value;
5161:       A->nonzerostate++;
5162: noinsert:;
5163:       low = i + 1;
5164:     }
5165:     ailen[row] = nrow;
5166:   }
5167:   PetscFunctionReturnVoid();
5168: }