Actual source code: sbaij.c

petsc-3.13.6 2020-09-29
Report Typos and Errors

  2: /*
  3:     Defines the basic matrix operations for the SBAIJ (compressed row)
  4:   matrix storage format.
  5: */
  6:  #include <../src/mat/impls/baij/seq/baij.h>
  7:  #include <../src/mat/impls/sbaij/seq/sbaij.h>
  8:  #include <petscblaslapack.h>

 10:  #include <../src/mat/impls/sbaij/seq/relax.h>
 11: #define USESHORT
 12:  #include <../src/mat/impls/sbaij/seq/relax.h>

 14: #if defined(PETSC_HAVE_ELEMENTAL)
 15: PETSC_INTERN PetscErrorCode MatConvert_SeqSBAIJ_Elemental(Mat,MatType,MatReuse,Mat*);
 16: #endif
 17: PETSC_INTERN PetscErrorCode MatConvert_MPISBAIJ_Basic(Mat,MatType,MatReuse,Mat*);

 19: /*
 20:      Checks for missing diagonals
 21: */
 22: PetscErrorCode MatMissingDiagonal_SeqSBAIJ(Mat A,PetscBool  *missing,PetscInt *dd)
 23: {
 24:   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
 26:   PetscInt       *diag,*ii = a->i,i;

 29:   MatMarkDiagonal_SeqSBAIJ(A);
 30:   *missing = PETSC_FALSE;
 31:   if (A->rmap->n > 0 && !ii) {
 32:     *missing = PETSC_TRUE;
 33:     if (dd) *dd = 0;
 34:     PetscInfo(A,"Matrix has no entries therefore is missing diagonal\n");
 35:   } else {
 36:     diag = a->diag;
 37:     for (i=0; i<a->mbs; i++) {
 38:       if (diag[i] >= ii[i+1]) {
 39:         *missing = PETSC_TRUE;
 40:         if (dd) *dd = i;
 41:         break;
 42:       }
 43:     }
 44:   }
 45:   return(0);
 46: }

 48: PetscErrorCode MatMarkDiagonal_SeqSBAIJ(Mat A)
 49: {
 50:   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
 52:   PetscInt       i,j;

 55:   if (!a->diag) {
 56:     PetscMalloc1(a->mbs,&a->diag);
 57:     PetscLogObjectMemory((PetscObject)A,a->mbs*sizeof(PetscInt));
 58:     a->free_diag = PETSC_TRUE;
 59:   }
 60:   for (i=0; i<a->mbs; i++) {
 61:     a->diag[i] = a->i[i+1];
 62:     for (j=a->i[i]; j<a->i[i+1]; j++) {
 63:       if (a->j[j] == i) {
 64:         a->diag[i] = j;
 65:         break;
 66:       }
 67:     }
 68:   }
 69:   return(0);
 70: }

 72: static PetscErrorCode MatGetRowIJ_SeqSBAIJ(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *nn,const PetscInt *inia[],const PetscInt *inja[],PetscBool  *done)
 73: {
 74:   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
 76:   PetscInt       i,j,n = a->mbs,nz = a->i[n],*tia,*tja,bs = A->rmap->bs,k,l,cnt;
 77:   PetscInt       **ia = (PetscInt**)inia,**ja = (PetscInt**)inja;

 80:   *nn = n;
 81:   if (!ia) return(0);
 82:   if (symmetric) {
 83:     MatToSymmetricIJ_SeqAIJ(n,a->i,a->j,PETSC_FALSE,0,0,&tia,&tja);
 84:     nz   = tia[n];
 85:   } else {
 86:     tia = a->i; tja = a->j;
 87:   }

 89:   if (!blockcompressed && bs > 1) {
 90:     (*nn) *= bs;
 91:     /* malloc & create the natural set of indices */
 92:     PetscMalloc1((n+1)*bs,ia);
 93:     if (n) {
 94:       (*ia)[0] = oshift;
 95:       for (j=1; j<bs; j++) {
 96:         (*ia)[j] = (tia[1]-tia[0])*bs+(*ia)[j-1];
 97:       }
 98:     }

100:     for (i=1; i<n; i++) {
101:       (*ia)[i*bs] = (tia[i]-tia[i-1])*bs + (*ia)[i*bs-1];
102:       for (j=1; j<bs; j++) {
103:         (*ia)[i*bs+j] = (tia[i+1]-tia[i])*bs + (*ia)[i*bs+j-1];
104:       }
105:     }
106:     if (n) {
107:       (*ia)[n*bs] = (tia[n]-tia[n-1])*bs + (*ia)[n*bs-1];
108:     }

110:     if (inja) {
111:       PetscMalloc1(nz*bs*bs,ja);
112:       cnt = 0;
113:       for (i=0; i<n; i++) {
114:         for (j=0; j<bs; j++) {
115:           for (k=tia[i]; k<tia[i+1]; k++) {
116:             for (l=0; l<bs; l++) {
117:               (*ja)[cnt++] = bs*tja[k] + l;
118:             }
119:           }
120:         }
121:       }
122:     }

124:     if (symmetric) { /* deallocate memory allocated in MatToSymmetricIJ_SeqAIJ() */
125:       PetscFree(tia);
126:       PetscFree(tja);
127:     }
128:   } else if (oshift == 1) {
129:     if (symmetric) {
130:       nz = tia[A->rmap->n/bs];
131:       /*  add 1 to i and j indices */
132:       for (i=0; i<A->rmap->n/bs+1; i++) tia[i] = tia[i] + 1;
133:       *ia = tia;
134:       if (ja) {
135:         for (i=0; i<nz; i++) tja[i] = tja[i] + 1;
136:         *ja = tja;
137:       }
138:     } else {
139:       nz = a->i[A->rmap->n/bs];
140:       /* malloc space and  add 1 to i and j indices */
141:       PetscMalloc1(A->rmap->n/bs+1,ia);
142:       for (i=0; i<A->rmap->n/bs+1; i++) (*ia)[i] = a->i[i] + 1;
143:       if (ja) {
144:         PetscMalloc1(nz,ja);
145:         for (i=0; i<nz; i++) (*ja)[i] = a->j[i] + 1;
146:       }
147:     }
148:   } else {
149:     *ia = tia;
150:     if (ja) *ja = tja;
151:   }
152:   return(0);
153: }

155: static PetscErrorCode MatRestoreRowIJ_SeqSBAIJ(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *nn,const PetscInt *ia[],const PetscInt *ja[],PetscBool  *done)
156: {

160:   if (!ia) return(0);
161:   if ((!blockcompressed && A->rmap->bs > 1) || (symmetric || oshift == 1)) {
162:     PetscFree(*ia);
163:     if (ja) {PetscFree(*ja);}
164:   }
165:   return(0);
166: }

168: PetscErrorCode MatDestroy_SeqSBAIJ(Mat A)
169: {
170:   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;

174: #if defined(PETSC_USE_LOG)
175:   PetscLogObjectState((PetscObject)A,"Rows=%D, NZ=%D",A->rmap->N,a->nz);
176: #endif
177:   MatSeqXAIJFreeAIJ(A,&a->a,&a->j,&a->i);
178:   if (a->free_diag) {PetscFree(a->diag);}
179:   ISDestroy(&a->row);
180:   ISDestroy(&a->col);
181:   ISDestroy(&a->icol);
182:   PetscFree(a->idiag);
183:   PetscFree(a->inode.size);
184:   if (a->free_imax_ilen) {PetscFree2(a->imax,a->ilen);}
185:   PetscFree(a->solve_work);
186:   PetscFree(a->sor_work);
187:   PetscFree(a->solves_work);
188:   PetscFree(a->mult_work);
189:   PetscFree(a->saved_values);
190:   if (a->free_jshort) {PetscFree(a->jshort);}
191:   PetscFree(a->inew);
192:   MatDestroy(&a->parent);
193:   PetscFree(A->data);

195:   PetscObjectChangeTypeName((PetscObject)A,0);
196:   PetscObjectComposeFunction((PetscObject)A,"MatStoreValues_C",NULL);
197:   PetscObjectComposeFunction((PetscObject)A,"MatRetrieveValues_C",NULL);
198:   PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetColumnIndices_C",NULL);
199:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqaij_C",NULL);
200:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqbaij_C",NULL);
201:   PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetPreallocation_C",NULL);
202:   PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetPreallocationCSR_C",NULL);
203: #if defined(PETSC_HAVE_ELEMENTAL)
204:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_elemental_C",NULL);
205: #endif
206:   return(0);
207: }

209: PetscErrorCode MatSetOption_SeqSBAIJ(Mat A,MatOption op,PetscBool flg)
210: {
211:   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
212: #if defined(PETSC_USE_COMPLEX)
213:   PetscInt       bs;
214: #endif

218: #if defined(PETSC_USE_COMPLEX)
219:   MatGetBlockSize(A,&bs);
220: #endif
221:   switch (op) {
222:   case MAT_ROW_ORIENTED:
223:     a->roworiented = flg;
224:     break;
225:   case MAT_KEEP_NONZERO_PATTERN:
226:     a->keepnonzeropattern = flg;
227:     break;
228:   case MAT_NEW_NONZERO_LOCATIONS:
229:     a->nonew = (flg ? 0 : 1);
230:     break;
231:   case MAT_NEW_NONZERO_LOCATION_ERR:
232:     a->nonew = (flg ? -1 : 0);
233:     break;
234:   case MAT_NEW_NONZERO_ALLOCATION_ERR:
235:     a->nonew = (flg ? -2 : 0);
236:     break;
237:   case MAT_UNUSED_NONZERO_LOCATION_ERR:
238:     a->nounused = (flg ? -1 : 0);
239:     break;
240:   case MAT_NEW_DIAGONALS:
241:   case MAT_IGNORE_OFF_PROC_ENTRIES:
242:   case MAT_USE_HASH_TABLE:
243:   case MAT_SORTED_FULL:
244:     PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
245:     break;
246:   case MAT_HERMITIAN:
247: #if defined(PETSC_USE_COMPLEX)
248:     if (flg) { /* disable transpose ops */
249:       if (bs > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for Hermitian with block size greater than 1");
250:       A->ops->multtranspose    = NULL;
251:       A->ops->multtransposeadd = NULL;
252:       A->symmetric             = PETSC_FALSE;
253:     }
254: #endif
255:     break;
256:   case MAT_SYMMETRIC:
257:   case MAT_SPD:
258: #if defined(PETSC_USE_COMPLEX)
259:     if (flg) { /* An hermitian and symmetric matrix has zero imaginary part (restore back transpose ops) */
260:       A->ops->multtranspose    = A->ops->mult;
261:       A->ops->multtransposeadd = A->ops->multadd;
262:     }
263: #endif
264:     break;
265:     /* These options are handled directly by MatSetOption() */
266:   case MAT_STRUCTURALLY_SYMMETRIC:
267:   case MAT_SYMMETRY_ETERNAL:
268:   case MAT_STRUCTURE_ONLY:
269:     /* These options are handled directly by MatSetOption() */
270:     break;
271:   case MAT_IGNORE_LOWER_TRIANGULAR:
272:     a->ignore_ltriangular = flg;
273:     break;
274:   case MAT_ERROR_LOWER_TRIANGULAR:
275:     a->ignore_ltriangular = flg;
276:     break;
277:   case MAT_GETROW_UPPERTRIANGULAR:
278:     a->getrow_utriangular = flg;
279:     break;
280:   case MAT_SUBMAT_SINGLEIS:
281:     break;
282:   default:
283:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %d",op);
284:   }
285:   return(0);
286: }

288: PetscErrorCode MatGetRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
289: {
290:   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;

294:   if (A && !a->getrow_utriangular) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatGetRow is not supported for SBAIJ matrix format. Getting the upper triangular part of row, run with -mat_getrow_uppertriangular, call MatSetOption(mat,MAT_GETROW_UPPERTRIANGULAR,PETSC_TRUE) or MatGetRowUpperTriangular()");

296:   /* Get the upper triangular part of the row */
297:   MatGetRow_SeqBAIJ_private(A,row,nz,idx,v,a->i,a->j,a->a);
298:   return(0);
299: }

301: PetscErrorCode MatRestoreRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
302: {

306:   if (idx) {PetscFree(*idx);}
307:   if (v)   {PetscFree(*v);}
308:   return(0);
309: }

311: PetscErrorCode MatGetRowUpperTriangular_SeqSBAIJ(Mat A)
312: {
313:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;

316:   a->getrow_utriangular = PETSC_TRUE;
317:   return(0);
318: }

320: PetscErrorCode MatRestoreRowUpperTriangular_SeqSBAIJ(Mat A)
321: {
322:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;

325:   a->getrow_utriangular = PETSC_FALSE;
326:   return(0);
327: }

329: PetscErrorCode MatTranspose_SeqSBAIJ(Mat A,MatReuse reuse,Mat *B)
330: {

334:   if (reuse == MAT_INITIAL_MATRIX) {
335:     MatDuplicate(A,MAT_COPY_VALUES,B);
336:   } else if (reuse == MAT_REUSE_MATRIX) {
337:     MatCopy(A,*B,SAME_NONZERO_PATTERN);
338:   }
339:   return(0);
340: }

342: PetscErrorCode MatView_SeqSBAIJ_ASCII(Mat A,PetscViewer viewer)
343: {
344:   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
345:   PetscErrorCode    ierr;
346:   PetscInt          i,j,bs = A->rmap->bs,k,l,bs2=a->bs2;
347:   PetscViewerFormat format;
348:   PetscInt          *diag;

351:   PetscViewerGetFormat(viewer,&format);
352:   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
353:     PetscViewerASCIIPrintf(viewer,"  block size is %D\n",bs);
354:   } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
355:     Mat        aij;
356:     const char *matname;

358:     if (A->factortype && bs>1) {
359:       PetscPrintf(PETSC_COMM_SELF,"Warning: matrix is factored with bs>1. MatView() with PETSC_VIEWER_ASCII_MATLAB is not supported and ignored!\n");
360:       return(0);
361:     }
362:     MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&aij);
363:     PetscObjectGetName((PetscObject)A,&matname);
364:     PetscObjectSetName((PetscObject)aij,matname);
365:     MatView(aij,viewer);
366:     MatDestroy(&aij);
367:   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
368:     PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
369:     for (i=0; i<a->mbs; i++) {
370:       for (j=0; j<bs; j++) {
371:         PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);
372:         for (k=a->i[i]; k<a->i[i+1]; k++) {
373:           for (l=0; l<bs; l++) {
374: #if defined(PETSC_USE_COMPLEX)
375:             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
376:               PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",bs*a->j[k]+l,
377:                                             (double)PetscRealPart(a->a[bs2*k + l*bs + j]),(double)PetscImaginaryPart(a->a[bs2*k + l*bs + j]));
378:             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
379:               PetscViewerASCIIPrintf(viewer," (%D, %g - %g i) ",bs*a->j[k]+l,
380:                                             (double)PetscRealPart(a->a[bs2*k + l*bs + j]),-(double)PetscImaginaryPart(a->a[bs2*k + l*bs + j]));
381:             } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
382:               PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,(double)PetscRealPart(a->a[bs2*k + l*bs + j]));
383:             }
384: #else
385:             if (a->a[bs2*k + l*bs + j] != 0.0) {
386:               PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,(double)a->a[bs2*k + l*bs + j]);
387:             }
388: #endif
389:           }
390:         }
391:         PetscViewerASCIIPrintf(viewer,"\n");
392:       }
393:     }
394:     PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
395:   } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
396:     return(0);
397:   } else {
398:     PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
399:     if (A->factortype) { /* for factored matrix */
400:       if (bs>1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"matrix is factored with bs>1. Not implemented yet");

402:       diag=a->diag;
403:       for (i=0; i<a->mbs; i++) { /* for row block i */
404:         PetscViewerASCIIPrintf(viewer,"row %D:",i);
405:         /* diagonal entry */
406: #if defined(PETSC_USE_COMPLEX)
407:         if (PetscImaginaryPart(a->a[diag[i]]) > 0.0) {
408:           PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",a->j[diag[i]],(double)PetscRealPart(1.0/a->a[diag[i]]),(double)PetscImaginaryPart(1.0/a->a[diag[i]]));
409:         } else if (PetscImaginaryPart(a->a[diag[i]]) < 0.0) {
410:           PetscViewerASCIIPrintf(viewer," (%D, %g - %g i) ",a->j[diag[i]],(double)PetscRealPart(1.0/a->a[diag[i]]),-(double)PetscImaginaryPart(1.0/a->a[diag[i]]));
411:         } else {
412:           PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[diag[i]],(double)PetscRealPart(1.0/a->a[diag[i]]));
413:         }
414: #else
415:         PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[diag[i]],(double)(1.0/a->a[diag[i]]));
416: #endif
417:         /* off-diagonal entries */
418:         for (k=a->i[i]; k<a->i[i+1]-1; k++) {
419: #if defined(PETSC_USE_COMPLEX)
420:           if (PetscImaginaryPart(a->a[k]) > 0.0) {
421:             PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",bs*a->j[k],(double)PetscRealPart(a->a[k]),(double)PetscImaginaryPart(a->a[k]));
422:           } else if (PetscImaginaryPart(a->a[k]) < 0.0) {
423:             PetscViewerASCIIPrintf(viewer," (%D, %g - %g i) ",bs*a->j[k],(double)PetscRealPart(a->a[k]),-(double)PetscImaginaryPart(a->a[k]));
424:           } else {
425:             PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k],(double)PetscRealPart(a->a[k]));
426:           }
427: #else
428:           PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[k],(double)a->a[k]);
429: #endif
430:         }
431:         PetscViewerASCIIPrintf(viewer,"\n");
432:       }

434:     } else { /* for non-factored matrix */
435:       for (i=0; i<a->mbs; i++) { /* for row block i */
436:         for (j=0; j<bs; j++) {   /* for row bs*i + j */
437:           PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);
438:           for (k=a->i[i]; k<a->i[i+1]; k++) { /* for column block */
439:             for (l=0; l<bs; l++) {            /* for column */
440: #if defined(PETSC_USE_COMPLEX)
441:               if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) {
442:                 PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",bs*a->j[k]+l,
443:                                               (double)PetscRealPart(a->a[bs2*k + l*bs + j]),(double)PetscImaginaryPart(a->a[bs2*k + l*bs + j]));
444:               } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) {
445:                 PetscViewerASCIIPrintf(viewer," (%D, %g - %g i) ",bs*a->j[k]+l,
446:                                               (double)PetscRealPart(a->a[bs2*k + l*bs + j]),-(double)PetscImaginaryPart(a->a[bs2*k + l*bs + j]));
447:               } else {
448:                 PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,(double)PetscRealPart(a->a[bs2*k + l*bs + j]));
449:               }
450: #else
451:               PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,(double)a->a[bs2*k + l*bs + j]);
452: #endif
453:             }
454:           }
455:           PetscViewerASCIIPrintf(viewer,"\n");
456:         }
457:       }
458:     }
459:     PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
460:   }
461:   PetscViewerFlush(viewer);
462:   return(0);
463: }

465:  #include <petscdraw.h>
466: static PetscErrorCode MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw,void *Aa)
467: {
468:   Mat            A = (Mat) Aa;
469:   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data;
471:   PetscInt       row,i,j,k,l,mbs=a->mbs,color,bs=A->rmap->bs,bs2=a->bs2;
472:   PetscReal      xl,yl,xr,yr,x_l,x_r,y_l,y_r;
473:   MatScalar      *aa;
474:   PetscViewer    viewer;

477:   PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);
478:   PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);

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

482:   PetscDrawCollectiveBegin(draw);
483:   PetscDrawString(draw, .3*(xl+xr), .3*(yl+yr), PETSC_DRAW_BLACK, "symmetric");
484:   /* Blue for negative, Cyan for zero and  Red for positive */
485:   color = PETSC_DRAW_BLUE;
486:   for (i=0,row=0; i<mbs; i++,row+=bs) {
487:     for (j=a->i[i]; j<a->i[i+1]; j++) {
488:       y_l = A->rmap->N - row - 1.0; y_r = y_l + 1.0;
489:       x_l = a->j[j]*bs; x_r = x_l + 1.0;
490:       aa  = a->a + j*bs2;
491:       for (k=0; k<bs; k++) {
492:         for (l=0; l<bs; l++) {
493:           if (PetscRealPart(*aa++) >=  0.) continue;
494:           PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
495:         }
496:       }
497:     }
498:   }
499:   color = PETSC_DRAW_CYAN;
500:   for (i=0,row=0; i<mbs; i++,row+=bs) {
501:     for (j=a->i[i]; j<a->i[i+1]; j++) {
502:       y_l = A->rmap->N - row - 1.0; y_r = y_l + 1.0;
503:       x_l = a->j[j]*bs; x_r = x_l + 1.0;
504:       aa = a->a + j*bs2;
505:       for (k=0; k<bs; k++) {
506:         for (l=0; l<bs; l++) {
507:           if (PetscRealPart(*aa++) != 0.) continue;
508:           PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
509:         }
510:       }
511:     }
512:   }
513:   color = PETSC_DRAW_RED;
514:   for (i=0,row=0; i<mbs; i++,row+=bs) {
515:     for (j=a->i[i]; j<a->i[i+1]; j++) {
516:       y_l = A->rmap->N - row - 1.0; y_r = y_l + 1.0;
517:       x_l = a->j[j]*bs; x_r = x_l + 1.0;
518:       aa = a->a + j*bs2;
519:       for (k=0; k<bs; k++) {
520:         for (l=0; l<bs; l++) {
521:           if (PetscRealPart(*aa++) <= 0.) continue;
522:           PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
523:         }
524:       }
525:     }
526:   }
527:   PetscDrawCollectiveEnd(draw);
528:   return(0);
529: }

531: static PetscErrorCode MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer)
532: {
534:   PetscReal      xl,yl,xr,yr,w,h;
535:   PetscDraw      draw;
536:   PetscBool      isnull;

539:   PetscViewerDrawGetDraw(viewer,0,&draw);
540:   PetscDrawIsNull(draw,&isnull);
541:   if (isnull) return(0);

543:   xr   = A->rmap->N; yr = A->rmap->N; h = yr/10.0; w = xr/10.0;
544:   xr  += w;          yr += h;        xl = -w;     yl = -h;
545:   PetscDrawSetCoordinates(draw,xl,yl,xr,yr);
546:   PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);
547:   PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);
548:   PetscObjectCompose((PetscObject)A,"Zoomviewer",NULL);
549:   PetscDrawSave(draw);
550:   return(0);
551: }

553: /* Used for both MPIBAIJ and MPISBAIJ matrices */
554: #define MatView_SeqSBAIJ_Binary MatView_SeqBAIJ_Binary

556: PetscErrorCode MatView_SeqSBAIJ(Mat A,PetscViewer viewer)
557: {
559:   PetscBool      iascii,isbinary,isdraw;

562:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
563:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
564:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
565:   if (iascii) {
566:     MatView_SeqSBAIJ_ASCII(A,viewer);
567:   } else if (isbinary) {
568:     MatView_SeqSBAIJ_Binary(A,viewer);
569:   } else if (isdraw) {
570:     MatView_SeqSBAIJ_Draw(A,viewer);
571:   } else {
572:     Mat        B;
573:     const char *matname;
574:     MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&B);
575:     PetscObjectGetName((PetscObject)A,&matname);
576:     PetscObjectSetName((PetscObject)B,matname);
577:     MatView(B,viewer);
578:     MatDestroy(&B);
579:   }
580:   return(0);
581: }


584: PetscErrorCode MatGetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],PetscScalar v[])
585: {
586:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
587:   PetscInt     *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
588:   PetscInt     *ai = a->i,*ailen = a->ilen;
589:   PetscInt     brow,bcol,ridx,cidx,bs=A->rmap->bs,bs2=a->bs2;
590:   MatScalar    *ap,*aa = a->a;

593:   for (k=0; k<m; k++) { /* loop over rows */
594:     row = im[k]; brow = row/bs;
595:     if (row < 0) {v += n; continue;} /* SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",row); */
596:     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);
597:     rp   = aj + ai[brow]; ap = aa + bs2*ai[brow];
598:     nrow = ailen[brow];
599:     for (l=0; l<n; l++) { /* loop over columns */
600:       if (in[l] < 0) {v++; continue;} /* SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %D",in[l]); */
601:       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);
602:       col  = in[l];
603:       bcol = col/bs;
604:       cidx = col%bs;
605:       ridx = row%bs;
606:       high = nrow;
607:       low  = 0; /* assume unsorted */
608:       while (high-low > 5) {
609:         t = (low+high)/2;
610:         if (rp[t] > bcol) high = t;
611:         else              low  = t;
612:       }
613:       for (i=low; i<high; i++) {
614:         if (rp[i] > bcol) break;
615:         if (rp[i] == bcol) {
616:           *v++ = ap[bs2*i+bs*cidx+ridx];
617:           goto finished;
618:         }
619:       }
620:       *v++ = 0.0;
621: finished:;
622:     }
623:   }
624:   return(0);
625: }


628: PetscErrorCode MatSetValuesBlocked_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
629: {
630:   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
631:   PetscErrorCode    ierr;
632:   PetscInt          *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,lastcol = -1;
633:   PetscInt          *imax      =a->imax,*ai=a->i,*ailen=a->ilen;
634:   PetscInt          *aj        =a->j,nonew=a->nonew,bs2=a->bs2,bs=A->rmap->bs,stepval;
635:   PetscBool         roworiented=a->roworiented;
636:   const PetscScalar *value     = v;
637:   MatScalar         *ap,*aa = a->a,*bap;

640:   if (roworiented) stepval = (n-1)*bs;
641:   else stepval = (m-1)*bs;

643:   for (k=0; k<m; k++) { /* loop over added rows */
644:     row = im[k];
645:     if (row < 0) continue;
646: #if defined(PETSC_USE_DEBUG)
647:     if (row >= a->mbs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Block index row too large %D max %D",row,a->mbs-1);
648: #endif
649:     rp   = aj + ai[row];
650:     ap   = aa + bs2*ai[row];
651:     rmax = imax[row];
652:     nrow = ailen[row];
653:     low  = 0;
654:     high = nrow;
655:     for (l=0; l<n; l++) { /* loop over added columns */
656:       if (in[l] < 0) continue;
657:       col = in[l];
658: #if defined(PETSC_USE_DEBUG)
659:       if (col >= a->nbs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Block index column too large %D max %D",col,a->nbs-1);
660: #endif
661:       if (col < row) {
662:         if (a->ignore_ltriangular) continue; /* ignore lower triangular block */
663:         else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Lower triangular value cannot be set for sbaij format. Ignoring these values, run with -mat_ignore_lower_triangular or call MatSetOption(mat,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE)");
664:       }
665:       if (roworiented) value = v + k*(stepval+bs)*bs + l*bs;
666:       else value = v + l*(stepval+bs)*bs + k*bs;

668:       if (col <= lastcol) low = 0;
669:       else high = nrow;

671:       lastcol = col;
672:       while (high-low > 7) {
673:         t = (low+high)/2;
674:         if (rp[t] > col) high = t;
675:         else             low  = t;
676:       }
677:       for (i=low; i<high; i++) {
678:         if (rp[i] > col) break;
679:         if (rp[i] == col) {
680:           bap = ap +  bs2*i;
681:           if (roworiented) {
682:             if (is == ADD_VALUES) {
683:               for (ii=0; ii<bs; ii++,value+=stepval) {
684:                 for (jj=ii; jj<bs2; jj+=bs) {
685:                   bap[jj] += *value++;
686:                 }
687:               }
688:             } else {
689:               for (ii=0; ii<bs; ii++,value+=stepval) {
690:                 for (jj=ii; jj<bs2; jj+=bs) {
691:                   bap[jj] = *value++;
692:                 }
693:                }
694:             }
695:           } else {
696:             if (is == ADD_VALUES) {
697:               for (ii=0; ii<bs; ii++,value+=stepval) {
698:                 for (jj=0; jj<bs; jj++) {
699:                   *bap++ += *value++;
700:                 }
701:               }
702:             } else {
703:               for (ii=0; ii<bs; ii++,value+=stepval) {
704:                 for (jj=0; jj<bs; jj++) {
705:                   *bap++  = *value++;
706:                 }
707:               }
708:             }
709:           }
710:           goto noinsert2;
711:         }
712:       }
713:       if (nonew == 1) goto noinsert2;
714:       if (nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new block index nonzero block (%D, %D) in the matrix", row, col);
715:       MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar);
716:       N = nrow++ - 1; high++;
717:       /* shift up all the later entries in this row */
718:       PetscArraymove(rp+i+1,rp+i,N-i+1);
719:       PetscArraymove(ap+bs2*(i+1),ap+bs2*i,bs2*(N-i+1));
720:       PetscArrayzero(ap+bs2*i,bs2);
721:       rp[i] = col;
722:       bap   = ap +  bs2*i;
723:       if (roworiented) {
724:         for (ii=0; ii<bs; ii++,value+=stepval) {
725:           for (jj=ii; jj<bs2; jj+=bs) {
726:             bap[jj] = *value++;
727:           }
728:         }
729:       } else {
730:         for (ii=0; ii<bs; ii++,value+=stepval) {
731:           for (jj=0; jj<bs; jj++) {
732:             *bap++ = *value++;
733:           }
734:         }
735:        }
736:     noinsert2:;
737:       low = i;
738:     }
739:     ailen[row] = nrow;
740:   }
741:   return(0);
742: }

744: /*
745:     This is not yet used
746: */
747: PetscErrorCode MatAssemblyEnd_SeqSBAIJ_SeqAIJ_Inode(Mat A)
748: {
749:   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
751:   const PetscInt *ai = a->i, *aj = a->j,*cols;
752:   PetscInt       i   = 0,j,blk_size,m = A->rmap->n,node_count = 0,nzx,nzy,*ns,row,nz,cnt,cnt2,*counts;
753:   PetscBool      flag;

756:   PetscMalloc1(m,&ns);
757:   while (i < m) {
758:     nzx = ai[i+1] - ai[i];       /* Number of nonzeros */
759:     /* Limits the number of elements in a node to 'a->inode.limit' */
760:     for (j=i+1,blk_size=1; j<m && blk_size <a->inode.limit; ++j,++blk_size) {
761:       nzy = ai[j+1] - ai[j];
762:       if (nzy != (nzx - j + i)) break;
763:       PetscArraycmp(aj + ai[i] + j - i,aj + ai[j],nzy,&flag);
764:       if (!flag) break;
765:     }
766:     ns[node_count++] = blk_size;

768:     i = j;
769:   }
770:   if (!a->inode.size && m && node_count > .9*m) {
771:     PetscFree(ns);
772:     PetscInfo2(A,"Found %D nodes out of %D rows. Not using Inode routines\n",node_count,m);
773:   } else {
774:     a->inode.node_count = node_count;

776:     PetscMalloc1(node_count,&a->inode.size);
777:     PetscLogObjectMemory((PetscObject)A,node_count*sizeof(PetscInt));
778:     PetscArraycpy(a->inode.size,ns,node_count);
779:     PetscFree(ns);
780:     PetscInfo3(A,"Found %D nodes of %D. Limit used: %D. Using Inode routines\n",node_count,m,a->inode.limit);

782:     /* count collections of adjacent columns in each inode */
783:     row = 0;
784:     cnt = 0;
785:     for (i=0; i<node_count; i++) {
786:       cols = aj + ai[row] + a->inode.size[i];
787:       nz   = ai[row+1] - ai[row] - a->inode.size[i];
788:       for (j=1; j<nz; j++) {
789:         if (cols[j] != cols[j-1]+1) cnt++;
790:       }
791:       cnt++;
792:       row += a->inode.size[i];
793:     }
794:     PetscMalloc1(2*cnt,&counts);
795:     cnt  = 0;
796:     row  = 0;
797:     for (i=0; i<node_count; i++) {
798:       cols = aj + ai[row] + a->inode.size[i];
799:       counts[2*cnt] = cols[0];
800:       nz   = ai[row+1] - ai[row] - a->inode.size[i];
801:       cnt2 = 1;
802:       for (j=1; j<nz; j++) {
803:         if (cols[j] != cols[j-1]+1) {
804:           counts[2*(cnt++)+1] = cnt2;
805:           counts[2*cnt]       = cols[j];
806:           cnt2 = 1;
807:         } else cnt2++;
808:       }
809:       counts[2*(cnt++)+1] = cnt2;
810:       row += a->inode.size[i];
811:     }
812:     PetscIntView(2*cnt,counts,0);
813:   }
814:   return(0);
815: }

817: PetscErrorCode MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode)
818: {
819:   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
821:   PetscInt       fshift = 0,i,*ai = a->i,*aj = a->j,*imax = a->imax;
822:   PetscInt       m      = A->rmap->N,*ip,N,*ailen = a->ilen;
823:   PetscInt       mbs    = a->mbs,bs2 = a->bs2,rmax = 0;
824:   MatScalar      *aa    = a->a,*ap;

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

829:   if (m) rmax = ailen[0];
830:   for (i=1; i<mbs; i++) {
831:     /* move each row back by the amount of empty slots (fshift) before it*/
832:     fshift += imax[i-1] - ailen[i-1];
833:     rmax    = PetscMax(rmax,ailen[i]);
834:     if (fshift) {
835:       ip = aj + ai[i];
836:       ap = aa + bs2*ai[i];
837:       N  = ailen[i];
838:       PetscArraymove(ip-fshift,ip,N);
839:       PetscArraymove(ap-bs2*fshift,ap,bs2*N);
840:     }
841:     ai[i] = ai[i-1] + ailen[i-1];
842:   }
843:   if (mbs) {
844:     fshift += imax[mbs-1] - ailen[mbs-1];
845:     ai[mbs] = ai[mbs-1] + ailen[mbs-1];
846:   }
847:   /* reset ilen and imax for each row */
848:   for (i=0; i<mbs; i++) {
849:     ailen[i] = imax[i] = ai[i+1] - ai[i];
850:   }
851:   a->nz = ai[mbs];

853:   /* diagonals may have moved, reset it */
854:   if (a->diag) {
855:     PetscArraycpy(a->diag,ai,mbs);
856:   }
857:   if (fshift && a->nounused == -1) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_PLIB, "Unused space detected in matrix: %D X %D block size %D, %D unneeded", m, A->cmap->n, A->rmap->bs, fshift*bs2);

859:   PetscInfo5(A,"Matrix size: %D X %D, block size %D; storage space: %D unneeded, %D used\n",m,A->rmap->N,A->rmap->bs,fshift*bs2,a->nz*bs2);
860:   PetscInfo1(A,"Number of mallocs during MatSetValues is %D\n",a->reallocs);
861:   PetscInfo1(A,"Most nonzeros blocks in any row is %D\n",rmax);

863:   A->info.mallocs    += a->reallocs;
864:   a->reallocs         = 0;
865:   A->info.nz_unneeded = (PetscReal)fshift*bs2;
866:   a->idiagvalid       = PETSC_FALSE;
867:   a->rmax             = rmax;

869:   if (A->cmap->n < 65536 && A->cmap->bs == 1) {
870:     if (a->jshort && a->free_jshort) {
871:       /* when matrix data structure is changed, previous jshort must be replaced */
872:       PetscFree(a->jshort);
873:     }
874:     PetscMalloc1(a->i[A->rmap->n],&a->jshort);
875:     PetscLogObjectMemory((PetscObject)A,a->i[A->rmap->n]*sizeof(unsigned short));
876:     for (i=0; i<a->i[A->rmap->n]; i++) a->jshort[i] = a->j[i];
877:     A->ops->mult   = MatMult_SeqSBAIJ_1_ushort;
878:     A->ops->sor    = MatSOR_SeqSBAIJ_ushort;
879:     a->free_jshort = PETSC_TRUE;
880:   }
881:   return(0);
882: }

884: /*
885:    This function returns an array of flags which indicate the locations of contiguous
886:    blocks that should be zeroed. for eg: if bs = 3  and is = [0,1,2,3,5,6,7,8,9]
887:    then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)]
888:    Assume: sizes should be long enough to hold all the values.
889: */
890: PetscErrorCode MatZeroRows_SeqSBAIJ_Check_Blocks(PetscInt idx[],PetscInt n,PetscInt bs,PetscInt sizes[], PetscInt *bs_max)
891: {
892:   PetscInt  i,j,k,row;
893:   PetscBool flg;

896:   for (i=0,j=0; i<n; j++) {
897:     row = idx[i];
898:     if (row%bs!=0) { /* Not the begining of a block */
899:       sizes[j] = 1;
900:       i++;
901:     } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */
902:       sizes[j] = 1;         /* Also makes sure atleast 'bs' values exist for next else */
903:       i++;
904:     } else { /* Begining of the block, so check if the complete block exists */
905:       flg = PETSC_TRUE;
906:       for (k=1; k<bs; k++) {
907:         if (row+k != idx[i+k]) { /* break in the block */
908:           flg = PETSC_FALSE;
909:           break;
910:         }
911:       }
912:       if (flg) { /* No break in the bs */
913:         sizes[j] = bs;
914:         i       += bs;
915:       } else {
916:         sizes[j] = 1;
917:         i++;
918:       }
919:     }
920:   }
921:   *bs_max = j;
922:   return(0);
923: }


926: /* Only add/insert a(i,j) with i<=j (blocks).
927:    Any a(i,j) with i>j input by user is ingored.
928: */

930: PetscErrorCode MatSetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
931: {
932:   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
934:   PetscInt       *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,lastcol = -1;
935:   PetscInt       *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented;
936:   PetscInt       *aj  =a->j,nonew=a->nonew,bs=A->rmap->bs,brow,bcol;
937:   PetscInt       ridx,cidx,bs2=a->bs2;
938:   MatScalar      *ap,value,*aa=a->a,*bap;

941:   for (k=0; k<m; k++) { /* loop over added rows */
942:     row  = im[k];       /* row number */
943:     brow = row/bs;      /* block row number */
944:     if (row < 0) continue;
945: #if defined(PETSC_USE_DEBUG)
946:     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);
947: #endif
948:     rp   = aj + ai[brow]; /*ptr to beginning of column value of the row block*/
949:     ap   = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/
950:     rmax = imax[brow];  /* maximum space allocated for this row */
951:     nrow = ailen[brow]; /* actual length of this row */
952:     low  = 0;
953:     high = nrow;
954:     for (l=0; l<n; l++) { /* loop over added columns */
955:       if (in[l] < 0) continue;
956: #if defined(PETSC_USE_DEBUG)
957:       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);
958: #endif
959:       col  = in[l];
960:       bcol = col/bs;              /* block col number */

962:       if (brow > bcol) {
963:         if (a->ignore_ltriangular) continue; /* ignore lower triangular values */
964:         else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Lower triangular value cannot be set for sbaij format. Ignoring these values, run with -mat_ignore_lower_triangular or call MatSetOption(mat,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE)");
965:       }

967:       ridx = row % bs; cidx = col % bs; /*row and col index inside the block */
968:       if ((brow==bcol && ridx<=cidx) || (brow<bcol)) {
969:         /* element value a(k,l) */
970:         if (roworiented) value = v[l + k*n];
971:         else value = v[k + l*m];

973:         /* move pointer bap to a(k,l) quickly and add/insert value */
974:         if (col <= lastcol) low = 0;
975:         else high = nrow;

977:         lastcol = col;
978:         while (high-low > 7) {
979:           t = (low+high)/2;
980:           if (rp[t] > bcol) high = t;
981:           else              low  = t;
982:         }
983:         for (i=low; i<high; i++) {
984:           if (rp[i] > bcol) break;
985:           if (rp[i] == bcol) {
986:             bap = ap +  bs2*i + bs*cidx + ridx;
987:             if (is == ADD_VALUES) *bap += value;
988:             else                  *bap  = value;
989:             /* for diag block, add/insert its symmetric element a(cidx,ridx) */
990:             if (brow == bcol && ridx < cidx) {
991:               bap = ap +  bs2*i + bs*ridx + cidx;
992:               if (is == ADD_VALUES) *bap += value;
993:               else                  *bap  = value;
994:             }
995:             goto noinsert1;
996:           }
997:         }

999:         if (nonew == 1) goto noinsert1;
1000:         if (nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col);
1001:         MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar);

1003:         N = nrow++ - 1; high++;
1004:         /* shift up all the later entries in this row */
1005:         PetscArraymove(rp+i+1,rp+i,N-i+1);
1006:         PetscArraymove(ap+bs2*(i+1),ap+bs2*i,bs2*(N-i+1));
1007:         PetscArrayzero(ap+bs2*i,bs2);
1008:         rp[i]                      = bcol;
1009:         ap[bs2*i + bs*cidx + ridx] = value;
1010:         /* for diag block, add/insert its symmetric element a(cidx,ridx) */
1011:         if (brow == bcol && ridx < cidx) {
1012:           ap[bs2*i + bs*ridx + cidx] = value;
1013:         }
1014:         A->nonzerostate++;
1015: noinsert1:;
1016:         low = i;
1017:       }
1018:     }   /* end of loop over added columns */
1019:     ailen[brow] = nrow;
1020:   }   /* end of loop over added rows */
1021:   return(0);
1022: }

1024: PetscErrorCode MatICCFactor_SeqSBAIJ(Mat inA,IS row,const MatFactorInfo *info)
1025: {
1026:   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)inA->data;
1027:   Mat            outA;
1029:   PetscBool      row_identity;

1032:   if (info->levels != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only levels=0 is supported for in-place icc");
1033:   ISIdentity(row,&row_identity);
1034:   if (!row_identity) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix reordering is not supported");
1035:   if (inA->rmap->bs != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix block size %D is not supported",inA->rmap->bs); /* Need to replace MatCholeskyFactorSymbolic_SeqSBAIJ_MSR()! */

1037:   outA            = inA;
1038:   inA->factortype = MAT_FACTOR_ICC;
1039:   PetscFree(inA->solvertype);
1040:   PetscStrallocpy(MATSOLVERPETSC,&inA->solvertype);

1042:   MatMarkDiagonal_SeqSBAIJ(inA);
1043:   MatSeqSBAIJSetNumericFactorization_inplace(inA,row_identity);

1045:   PetscObjectReference((PetscObject)row);
1046:   ISDestroy(&a->row);
1047:   a->row = row;
1048:   PetscObjectReference((PetscObject)row);
1049:   ISDestroy(&a->col);
1050:   a->col = row;

1052:   /* Create the invert permutation so that it can be used in MatCholeskyFactorNumeric() */
1053:   if (a->icol) {ISInvertPermutation(row,PETSC_DECIDE, &a->icol);}
1054:   PetscLogObjectParent((PetscObject)inA,(PetscObject)a->icol);

1056:   if (!a->solve_work) {
1057:     PetscMalloc1(inA->rmap->N+inA->rmap->bs,&a->solve_work);
1058:     PetscLogObjectMemory((PetscObject)inA,(inA->rmap->N+inA->rmap->bs)*sizeof(PetscScalar));
1059:   }

1061:   MatCholeskyFactorNumeric(outA,inA,info);
1062:   return(0);
1063: }

1065: PetscErrorCode  MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,PetscInt *indices)
1066: {
1067:   Mat_SeqSBAIJ   *baij = (Mat_SeqSBAIJ*)mat->data;
1068:   PetscInt       i,nz,n;

1072:   nz = baij->maxnz;
1073:   n  = mat->cmap->n;
1074:   for (i=0; i<nz; i++) baij->j[i] = indices[i];

1076:   baij->nz = nz;
1077:   for (i=0; i<n; i++) baij->ilen[i] = baij->imax[i];

1079:   MatSetOption(mat,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
1080:   return(0);
1081: }

1083: /*@
1084:   MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows
1085:   in the matrix.

1087:   Input Parameters:
1088:   +  mat     - the SeqSBAIJ matrix
1089:   -  indices - the column indices

1091:   Level: advanced

1093:   Notes:
1094:   This can be called if you have precomputed the nonzero structure of the
1095:   matrix and want to provide it to the matrix object to improve the performance
1096:   of the MatSetValues() operation.

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

1101:   MUST be called before any calls to MatSetValues()

1103:   .seealso: MatCreateSeqSBAIJ
1104: @*/
1105: PetscErrorCode  MatSeqSBAIJSetColumnIndices(Mat mat,PetscInt *indices)
1106: {

1112:   PetscUseMethod(mat,"MatSeqSBAIJSetColumnIndices_C",(Mat,PetscInt*),(mat,indices));
1113:   return(0);
1114: }

1116: PetscErrorCode MatCopy_SeqSBAIJ(Mat A,Mat B,MatStructure str)
1117: {
1119:   PetscBool      isbaij;

1122:   PetscObjectTypeCompareAny((PetscObject)B,&isbaij,MATSEQSBAIJ,MATMPISBAIJ,"");
1123:   if (!isbaij) SETERRQ1(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"Not for matrix type %s",((PetscObject)B)->type_name);
1124:   /* If the two matrices have the same copy implementation and nonzero pattern, use fast copy. */
1125:   if (str == SAME_NONZERO_PATTERN && A->ops->copy == B->ops->copy) {
1126:     Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1127:     Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data;

1129:     if (a->i[a->mbs] != b->i[b->mbs]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of nonzeros in two matrices are different");
1130:     if (a->mbs != b->mbs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of rows in two matrices are different");
1131:     if (a->bs2 != b->bs2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Different block size");
1132:     PetscArraycpy(b->a,a->a,a->bs2*a->i[a->mbs]);
1133:     PetscObjectStateIncrease((PetscObject)B);
1134:   } else {
1135:     MatGetRowUpperTriangular(A);
1136:     MatCopy_Basic(A,B,str);
1137:     MatRestoreRowUpperTriangular(A);
1138:   }
1139:   return(0);
1140: }

1142: PetscErrorCode MatSetUp_SeqSBAIJ(Mat A)
1143: {

1147:   MatSeqSBAIJSetPreallocation(A,A->rmap->bs,PETSC_DEFAULT,0);
1148:   return(0);
1149: }

1151: static PetscErrorCode MatSeqSBAIJGetArray_SeqSBAIJ(Mat A,PetscScalar *array[])
1152: {
1153:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;

1156:   *array = a->a;
1157:   return(0);
1158: }

1160: static PetscErrorCode MatSeqSBAIJRestoreArray_SeqSBAIJ(Mat A,PetscScalar *array[])
1161: {
1163:   *array = NULL;
1164:   return(0);
1165: }

1167: PetscErrorCode MatAXPYGetPreallocation_SeqSBAIJ(Mat Y,Mat X,PetscInt *nnz)
1168: {
1169:   PetscInt       bs = Y->rmap->bs,mbs = Y->rmap->N/bs;
1170:   Mat_SeqSBAIJ   *x = (Mat_SeqSBAIJ*)X->data;
1171:   Mat_SeqSBAIJ   *y = (Mat_SeqSBAIJ*)Y->data;

1175:   /* Set the number of nonzeros in the new matrix */
1176:   MatAXPYGetPreallocation_SeqX_private(mbs,x->i,x->j,y->i,y->j,nnz);
1177:   return(0);
1178: }

1180: PetscErrorCode MatAXPY_SeqSBAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str)
1181: {
1182:   Mat_SeqSBAIJ   *x=(Mat_SeqSBAIJ*)X->data, *y=(Mat_SeqSBAIJ*)Y->data;
1184:   PetscInt       bs=Y->rmap->bs,bs2=bs*bs;
1185:   PetscBLASInt   one = 1;

1188:   if (str == SAME_NONZERO_PATTERN) {
1189:     PetscScalar  alpha = a;
1190:     PetscBLASInt bnz;
1191:     PetscBLASIntCast(x->nz*bs2,&bnz);
1192:     PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&bnz,&alpha,x->a,&one,y->a,&one));
1193:     PetscObjectStateIncrease((PetscObject)Y);
1194:   } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
1195:     MatSetOption(X,MAT_GETROW_UPPERTRIANGULAR,PETSC_TRUE);
1196:     MatAXPY_Basic(Y,a,X,str);
1197:     MatSetOption(X,MAT_GETROW_UPPERTRIANGULAR,PETSC_FALSE);
1198:   } else {
1199:     Mat      B;
1200:     PetscInt *nnz;
1201:     if (bs != X->rmap->bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrices must have same block size");
1202:     MatGetRowUpperTriangular(X);
1203:     MatGetRowUpperTriangular(Y);
1204:     PetscMalloc1(Y->rmap->N,&nnz);
1205:     MatCreate(PetscObjectComm((PetscObject)Y),&B);
1206:     PetscObjectSetName((PetscObject)B,((PetscObject)Y)->name);
1207:     MatSetSizes(B,Y->rmap->n,Y->cmap->n,Y->rmap->N,Y->cmap->N);
1208:     MatSetBlockSizesFromMats(B,Y,Y);
1209:     MatSetType(B,((PetscObject)Y)->type_name);
1210:     MatAXPYGetPreallocation_SeqSBAIJ(Y,X,nnz);
1211:     MatSeqSBAIJSetPreallocation(B,bs,0,nnz);

1213:     MatAXPY_BasicWithPreallocation(B,Y,a,X,str);

1215:     MatHeaderReplace(Y,&B);
1216:     PetscFree(nnz);
1217:     MatRestoreRowUpperTriangular(X);
1218:     MatRestoreRowUpperTriangular(Y);
1219:   }
1220:   return(0);
1221: }

1223: PetscErrorCode MatIsSymmetric_SeqSBAIJ(Mat A,PetscReal tol,PetscBool  *flg)
1224: {
1226:   *flg = PETSC_TRUE;
1227:   return(0);
1228: }

1230: PetscErrorCode MatIsStructurallySymmetric_SeqSBAIJ(Mat A,PetscBool  *flg)
1231: {
1233:   *flg = PETSC_TRUE;
1234:   return(0);
1235: }

1237: PetscErrorCode MatIsHermitian_SeqSBAIJ(Mat A,PetscReal tol,PetscBool  *flg)
1238: {
1240:   *flg = PETSC_FALSE;
1241:   return(0);
1242: }

1244: PetscErrorCode MatRealPart_SeqSBAIJ(Mat A)
1245: {
1246:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1247:   PetscInt     i,nz = a->bs2*a->i[a->mbs];
1248:   MatScalar    *aa = a->a;

1251:   for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]);
1252:   return(0);
1253: }

1255: PetscErrorCode MatImaginaryPart_SeqSBAIJ(Mat A)
1256: {
1257:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1258:   PetscInt     i,nz = a->bs2*a->i[a->mbs];
1259:   MatScalar    *aa = a->a;

1262:   for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]);
1263:   return(0);
1264: }

1266: PetscErrorCode MatZeroRowsColumns_SeqSBAIJ(Mat A,PetscInt is_n,const PetscInt is_idx[],PetscScalar diag,Vec x, Vec b)
1267: {
1268:   Mat_SeqSBAIJ      *baij=(Mat_SeqSBAIJ*)A->data;
1269:   PetscErrorCode    ierr;
1270:   PetscInt          i,j,k,count;
1271:   PetscInt          bs   =A->rmap->bs,bs2=baij->bs2,row,col;
1272:   PetscScalar       zero = 0.0;
1273:   MatScalar         *aa;
1274:   const PetscScalar *xx;
1275:   PetscScalar       *bb;
1276:   PetscBool         *zeroed,vecs = PETSC_FALSE;

1279:   /* fix right hand side if needed */
1280:   if (x && b) {
1281:     VecGetArrayRead(x,&xx);
1282:     VecGetArray(b,&bb);
1283:     vecs = PETSC_TRUE;
1284:   }

1286:   /* zero the columns */
1287:   PetscCalloc1(A->rmap->n,&zeroed);
1288:   for (i=0; i<is_n; i++) {
1289:     if (is_idx[i] < 0 || is_idx[i] >= A->rmap->N) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range",is_idx[i]);
1290:     zeroed[is_idx[i]] = PETSC_TRUE;
1291:   }
1292:   if (vecs) {
1293:     for (i=0; i<A->rmap->N; i++) {
1294:       row = i/bs;
1295:       for (j=baij->i[row]; j<baij->i[row+1]; j++) {
1296:         for (k=0; k<bs; k++) {
1297:           col = bs*baij->j[j] + k;
1298:           if (col <= i) continue;
1299:           aa = ((MatScalar*)(baij->a)) + j*bs2 + (i%bs) + bs*k;
1300:           if (!zeroed[i] && zeroed[col]) bb[i]   -= aa[0]*xx[col];
1301:           if (zeroed[i] && !zeroed[col]) bb[col] -= aa[0]*xx[i];
1302:         }
1303:       }
1304:     }
1305:     for (i=0; i<is_n; i++) bb[is_idx[i]] = diag*xx[is_idx[i]];
1306:   }

1308:   for (i=0; i<A->rmap->N; i++) {
1309:     if (!zeroed[i]) {
1310:       row = i/bs;
1311:       for (j=baij->i[row]; j<baij->i[row+1]; j++) {
1312:         for (k=0; k<bs; k++) {
1313:           col = bs*baij->j[j] + k;
1314:           if (zeroed[col]) {
1315:             aa = ((MatScalar*)(baij->a)) + j*bs2 + (i%bs) + bs*k;
1316:             aa[0] = 0.0;
1317:           }
1318:         }
1319:       }
1320:     }
1321:   }
1322:   PetscFree(zeroed);
1323:   if (vecs) {
1324:     VecRestoreArrayRead(x,&xx);
1325:     VecRestoreArray(b,&bb);
1326:   }

1328:   /* zero the rows */
1329:   for (i=0; i<is_n; i++) {
1330:     row   = is_idx[i];
1331:     count = (baij->i[row/bs +1] - baij->i[row/bs])*bs;
1332:     aa    = ((MatScalar*)(baij->a)) + baij->i[row/bs]*bs2 + (row%bs);
1333:     for (k=0; k<count; k++) {
1334:       aa[0] =  zero;
1335:       aa   += bs;
1336:     }
1337:     if (diag != 0.0) {
1338:       (*A->ops->setvalues)(A,1,&row,1,&row,&diag,INSERT_VALUES);
1339:     }
1340:   }
1341:   MatAssemblyEnd_SeqSBAIJ(A,MAT_FINAL_ASSEMBLY);
1342:   return(0);
1343: }

1345: PetscErrorCode MatShift_SeqSBAIJ(Mat Y,PetscScalar a)
1346: {
1348:   Mat_SeqSBAIJ    *aij = (Mat_SeqSBAIJ*)Y->data;

1351:   if (!Y->preallocated || !aij->nz) {
1352:     MatSeqSBAIJSetPreallocation(Y,Y->rmap->bs,1,NULL);
1353:   }
1354:   MatShift_Basic(Y,a);
1355:   return(0);
1356: }

1358: /* -------------------------------------------------------------------*/
1359: static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ,
1360:                                        MatGetRow_SeqSBAIJ,
1361:                                        MatRestoreRow_SeqSBAIJ,
1362:                                        MatMult_SeqSBAIJ_N,
1363:                                /*  4*/ MatMultAdd_SeqSBAIJ_N,
1364:                                        MatMult_SeqSBAIJ_N,       /* transpose versions are same as non-transpose versions */
1365:                                        MatMultAdd_SeqSBAIJ_N,
1366:                                        0,
1367:                                        0,
1368:                                        0,
1369:                                /* 10*/ 0,
1370:                                        0,
1371:                                        MatCholeskyFactor_SeqSBAIJ,
1372:                                        MatSOR_SeqSBAIJ,
1373:                                        MatTranspose_SeqSBAIJ,
1374:                                /* 15*/ MatGetInfo_SeqSBAIJ,
1375:                                        MatEqual_SeqSBAIJ,
1376:                                        MatGetDiagonal_SeqSBAIJ,
1377:                                        MatDiagonalScale_SeqSBAIJ,
1378:                                        MatNorm_SeqSBAIJ,
1379:                                /* 20*/ 0,
1380:                                        MatAssemblyEnd_SeqSBAIJ,
1381:                                        MatSetOption_SeqSBAIJ,
1382:                                        MatZeroEntries_SeqSBAIJ,
1383:                                /* 24*/ 0,
1384:                                        0,
1385:                                        0,
1386:                                        0,
1387:                                        0,
1388:                                /* 29*/ MatSetUp_SeqSBAIJ,
1389:                                        0,
1390:                                        0,
1391:                                        0,
1392:                                        0,
1393:                                /* 34*/ MatDuplicate_SeqSBAIJ,
1394:                                        0,
1395:                                        0,
1396:                                        0,
1397:                                        MatICCFactor_SeqSBAIJ,
1398:                                /* 39*/ MatAXPY_SeqSBAIJ,
1399:                                        MatCreateSubMatrices_SeqSBAIJ,
1400:                                        MatIncreaseOverlap_SeqSBAIJ,
1401:                                        MatGetValues_SeqSBAIJ,
1402:                                        MatCopy_SeqSBAIJ,
1403:                                /* 44*/ 0,
1404:                                        MatScale_SeqSBAIJ,
1405:                                        MatShift_SeqSBAIJ,
1406:                                        0,
1407:                                        MatZeroRowsColumns_SeqSBAIJ,
1408:                                /* 49*/ 0,
1409:                                        MatGetRowIJ_SeqSBAIJ,
1410:                                        MatRestoreRowIJ_SeqSBAIJ,
1411:                                        0,
1412:                                        0,
1413:                                /* 54*/ 0,
1414:                                        0,
1415:                                        0,
1416:                                        0,
1417:                                        MatSetValuesBlocked_SeqSBAIJ,
1418:                                /* 59*/ MatCreateSubMatrix_SeqSBAIJ,
1419:                                        0,
1420:                                        0,
1421:                                        0,
1422:                                        0,
1423:                                /* 64*/ 0,
1424:                                        0,
1425:                                        0,
1426:                                        0,
1427:                                        0,
1428:                                /* 69*/ MatGetRowMaxAbs_SeqSBAIJ,
1429:                                        0,
1430:                                        MatConvert_MPISBAIJ_Basic,
1431:                                        0,
1432:                                        0,
1433:                                /* 74*/ 0,
1434:                                        0,
1435:                                        0,
1436:                                        0,
1437:                                        0,
1438:                                /* 79*/ 0,
1439:                                        0,
1440:                                        0,
1441:                                        MatGetInertia_SeqSBAIJ,
1442:                                        MatLoad_SeqSBAIJ,
1443:                                /* 84*/ MatIsSymmetric_SeqSBAIJ,
1444:                                        MatIsHermitian_SeqSBAIJ,
1445:                                        MatIsStructurallySymmetric_SeqSBAIJ,
1446:                                        0,
1447:                                        0,
1448:                                /* 89*/ 0,
1449:                                        0,
1450:                                        0,
1451:                                        0,
1452:                                        0,
1453:                                /* 94*/ 0,
1454:                                        0,
1455:                                        0,
1456:                                        0,
1457:                                        0,
1458:                                /* 99*/ 0,
1459:                                        0,
1460:                                        0,
1461:                                        0,
1462:                                        0,
1463:                                /*104*/ 0,
1464:                                        MatRealPart_SeqSBAIJ,
1465:                                        MatImaginaryPart_SeqSBAIJ,
1466:                                        MatGetRowUpperTriangular_SeqSBAIJ,
1467:                                        MatRestoreRowUpperTriangular_SeqSBAIJ,
1468:                                /*109*/ 0,
1469:                                        0,
1470:                                        0,
1471:                                        0,
1472:                                        MatMissingDiagonal_SeqSBAIJ,
1473:                                /*114*/ 0,
1474:                                        0,
1475:                                        0,
1476:                                        0,
1477:                                        0,
1478:                                /*119*/ 0,
1479:                                        0,
1480:                                        0,
1481:                                        0,
1482:                                        0,
1483:                                /*124*/ 0,
1484:                                        0,
1485:                                        0,
1486:                                        0,
1487:                                        0,
1488:                                /*129*/ 0,
1489:                                        0,
1490:                                        0,
1491:                                        0,
1492:                                        0,
1493:                                /*134*/ 0,
1494:                                        0,
1495:                                        0,
1496:                                        0,
1497:                                        0,
1498:                                /*139*/ MatSetBlockSizes_Default,
1499:                                        0,
1500:                                        0,
1501:                                        0,
1502:                                        0,
1503:                                 /*144*/MatCreateMPIMatConcatenateSeqMat_SeqSBAIJ
1504: };

1506: PetscErrorCode  MatStoreValues_SeqSBAIJ(Mat mat)
1507: {
1508:   Mat_SeqSBAIJ   *aij = (Mat_SeqSBAIJ*)mat->data;
1509:   PetscInt       nz   = aij->i[mat->rmap->N]*mat->rmap->bs*aij->bs2;

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

1515:   /* allocate space for values if not already there */
1516:   if (!aij->saved_values) {
1517:     PetscMalloc1(nz+1,&aij->saved_values);
1518:   }

1520:   /* copy values over */
1521:   PetscArraycpy(aij->saved_values,aij->a,nz);
1522:   return(0);
1523: }

1525: PetscErrorCode  MatRetrieveValues_SeqSBAIJ(Mat mat)
1526: {
1527:   Mat_SeqSBAIJ   *aij = (Mat_SeqSBAIJ*)mat->data;
1529:   PetscInt       nz = aij->i[mat->rmap->N]*mat->rmap->bs*aij->bs2;

1532:   if (aij->nonew != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first");
1533:   if (!aij->saved_values) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Must call MatStoreValues(A);first");

1535:   /* copy values over */
1536:   PetscArraycpy(aij->a,aij->saved_values,nz);
1537:   return(0);
1538: }

1540: static PetscErrorCode  MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat B,PetscInt bs,PetscInt nz,PetscInt *nnz)
1541: {
1542:   Mat_SeqSBAIJ   *b = (Mat_SeqSBAIJ*)B->data;
1544:   PetscInt       i,mbs,nbs,bs2;
1545:   PetscBool      skipallocation = PETSC_FALSE,flg = PETSC_FALSE,realalloc = PETSC_FALSE;

1548:   if (nz >= 0 || nnz) realalloc = PETSC_TRUE;

1550:   MatSetBlockSize(B,PetscAbs(bs));
1551:   PetscLayoutSetUp(B->rmap);
1552:   PetscLayoutSetUp(B->cmap);
1553:   if (B->rmap->N > B->cmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"SEQSBAIJ matrix cannot have more rows %D than columns %D",B->rmap->N,B->cmap->N);
1554:   PetscLayoutGetBlockSize(B->rmap,&bs);

1556:   B->preallocated = PETSC_TRUE;

1558:   mbs = B->rmap->N/bs;
1559:   nbs = B->cmap->n/bs;
1560:   bs2 = bs*bs;

1562:   if (mbs*bs != B->rmap->N || nbs*bs!=B->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize");

1564:   if (nz == MAT_SKIP_ALLOCATION) {
1565:     skipallocation = PETSC_TRUE;
1566:     nz             = 0;
1567:   }

1569:   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3;
1570:   if (nz < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %D",nz);
1571:   if (nnz) {
1572:     for (i=0; i<mbs; i++) {
1573:       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]);
1574:       if (nnz[i] > nbs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be greater than block row length: local row %D value %D block rowlength %D",i,nnz[i],nbs);
1575:     }
1576:   }

1578:   B->ops->mult             = MatMult_SeqSBAIJ_N;
1579:   B->ops->multadd          = MatMultAdd_SeqSBAIJ_N;
1580:   B->ops->multtranspose    = MatMult_SeqSBAIJ_N;
1581:   B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_N;

1583:   PetscOptionsGetBool(((PetscObject)B)->options,((PetscObject)B)->prefix,"-mat_no_unroll",&flg,NULL);
1584:   if (!flg) {
1585:     switch (bs) {
1586:     case 1:
1587:       B->ops->mult             = MatMult_SeqSBAIJ_1;
1588:       B->ops->multadd          = MatMultAdd_SeqSBAIJ_1;
1589:       B->ops->multtranspose    = MatMult_SeqSBAIJ_1;
1590:       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_1;
1591:       break;
1592:     case 2:
1593:       B->ops->mult             = MatMult_SeqSBAIJ_2;
1594:       B->ops->multadd          = MatMultAdd_SeqSBAIJ_2;
1595:       B->ops->multtranspose    = MatMult_SeqSBAIJ_2;
1596:       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_2;
1597:       break;
1598:     case 3:
1599:       B->ops->mult             = MatMult_SeqSBAIJ_3;
1600:       B->ops->multadd          = MatMultAdd_SeqSBAIJ_3;
1601:       B->ops->multtranspose    = MatMult_SeqSBAIJ_3;
1602:       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_3;
1603:       break;
1604:     case 4:
1605:       B->ops->mult             = MatMult_SeqSBAIJ_4;
1606:       B->ops->multadd          = MatMultAdd_SeqSBAIJ_4;
1607:       B->ops->multtranspose    = MatMult_SeqSBAIJ_4;
1608:       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_4;
1609:       break;
1610:     case 5:
1611:       B->ops->mult             = MatMult_SeqSBAIJ_5;
1612:       B->ops->multadd          = MatMultAdd_SeqSBAIJ_5;
1613:       B->ops->multtranspose    = MatMult_SeqSBAIJ_5;
1614:       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_5;
1615:       break;
1616:     case 6:
1617:       B->ops->mult             = MatMult_SeqSBAIJ_6;
1618:       B->ops->multadd          = MatMultAdd_SeqSBAIJ_6;
1619:       B->ops->multtranspose    = MatMult_SeqSBAIJ_6;
1620:       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_6;
1621:       break;
1622:     case 7:
1623:       B->ops->mult             = MatMult_SeqSBAIJ_7;
1624:       B->ops->multadd          = MatMultAdd_SeqSBAIJ_7;
1625:       B->ops->multtranspose    = MatMult_SeqSBAIJ_7;
1626:       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_7;
1627:       break;
1628:     }
1629:   }

1631:   b->mbs = mbs;
1632:   b->nbs = nbs;
1633:   if (!skipallocation) {
1634:     if (!b->imax) {
1635:       PetscMalloc2(mbs,&b->imax,mbs,&b->ilen);

1637:       b->free_imax_ilen = PETSC_TRUE;

1639:       PetscLogObjectMemory((PetscObject)B,2*mbs*sizeof(PetscInt));
1640:     }
1641:     if (!nnz) {
1642:       if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
1643:       else if (nz <= 0) nz = 1;
1644:       nz = PetscMin(nbs,nz);
1645:       for (i=0; i<mbs; i++) b->imax[i] = nz;
1646:       nz = nz*mbs; /* total nz */
1647:     } else {
1648:       PetscInt64 nz64 = 0;
1649:       for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz64 += nnz[i];}
1650:       PetscIntCast(nz64,&nz);
1651:     }
1652:     /* b->ilen will count nonzeros in each block row so far. */
1653:     for (i=0; i<mbs; i++) b->ilen[i] = 0;
1654:     /* nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */

1656:     /* allocate the matrix space */
1657:     MatSeqXAIJFreeAIJ(B,&b->a,&b->j,&b->i);
1658:     PetscMalloc3(bs2*nz,&b->a,nz,&b->j,B->rmap->N+1,&b->i);
1659:     PetscLogObjectMemory((PetscObject)B,(B->rmap->N+1)*sizeof(PetscInt)+nz*(bs2*sizeof(PetscScalar)+sizeof(PetscInt)));
1660:     PetscArrayzero(b->a,nz*bs2);
1661:     PetscArrayzero(b->j,nz);

1663:     b->singlemalloc = PETSC_TRUE;

1665:     /* pointer to beginning of each row */
1666:     b->i[0] = 0;
1667:     for (i=1; i<mbs+1; i++) b->i[i] = b->i[i-1] + b->imax[i-1];

1669:     b->free_a  = PETSC_TRUE;
1670:     b->free_ij = PETSC_TRUE;
1671:   } else {
1672:     b->free_a  = PETSC_FALSE;
1673:     b->free_ij = PETSC_FALSE;
1674:   }

1676:   b->bs2     = bs2;
1677:   b->nz      = 0;
1678:   b->maxnz   = nz;
1679:   b->inew    = 0;
1680:   b->jnew    = 0;
1681:   b->anew    = 0;
1682:   b->a2anew  = 0;
1683:   b->permute = PETSC_FALSE;

1685:   B->was_assembled = PETSC_FALSE;
1686:   B->assembled     = PETSC_FALSE;
1687:   if (realalloc) {MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);}
1688:   return(0);
1689: }

1691: PetscErrorCode MatSeqSBAIJSetPreallocationCSR_SeqSBAIJ(Mat B,PetscInt bs,const PetscInt ii[],const PetscInt jj[], const PetscScalar V[])
1692: {
1693:   PetscInt       i,j,m,nz,anz, nz_max=0,*nnz;
1694:   PetscScalar    *values=0;
1695:   PetscBool      roworiented = ((Mat_SeqSBAIJ*)B->data)->roworiented;

1699:   if (bs < 1) SETERRQ1(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_OUTOFRANGE,"Invalid block size specified, must be positive but it is %D",bs);
1700:   PetscLayoutSetBlockSize(B->rmap,bs);
1701:   PetscLayoutSetBlockSize(B->cmap,bs);
1702:   PetscLayoutSetUp(B->rmap);
1703:   PetscLayoutSetUp(B->cmap);
1704:   PetscLayoutGetBlockSize(B->rmap,&bs);
1705:   m      = B->rmap->n/bs;

1707:   if (ii[0]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"ii[0] must be 0 but it is %D",ii[0]);
1708:   PetscMalloc1(m+1,&nnz);
1709:   for (i=0; i<m; i++) {
1710:     nz = ii[i+1] - ii[i];
1711:     if (nz < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row %D has a negative number of columns %D",i,nz);
1712:     anz = 0;
1713:     for (j=0; j<nz; j++) {
1714:       /* count only values on the diagonal or above */
1715:       if (jj[ii[i] + j] >= i) {
1716:         anz = nz - j;
1717:         break;
1718:       }
1719:     }
1720:     nz_max = PetscMax(nz_max,anz);
1721:     nnz[i] = anz;
1722:   }
1723:   MatSeqSBAIJSetPreallocation(B,bs,0,nnz);
1724:   PetscFree(nnz);

1726:   values = (PetscScalar*)V;
1727:   if (!values) {
1728:     PetscCalloc1(bs*bs*nz_max,&values);
1729:   }
1730:   for (i=0; i<m; i++) {
1731:     PetscInt          ncols  = ii[i+1] - ii[i];
1732:     const PetscInt    *icols = jj + ii[i];
1733:     if (!roworiented || bs == 1) {
1734:       const PetscScalar *svals = values + (V ? (bs*bs*ii[i]) : 0);
1735:       MatSetValuesBlocked_SeqSBAIJ(B,1,&i,ncols,icols,svals,INSERT_VALUES);
1736:     } else {
1737:       for (j=0; j<ncols; j++) {
1738:         const PetscScalar *svals = values + (V ? (bs*bs*(ii[i]+j)) : 0);
1739:         MatSetValuesBlocked_SeqSBAIJ(B,1,&i,1,&icols[j],svals,INSERT_VALUES);
1740:       }
1741:     }
1742:   }
1743:   if (!V) { PetscFree(values); }
1744:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
1745:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
1746:   MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
1747:   return(0);
1748: }

1750: /*
1751:    This is used to set the numeric factorization for both Cholesky and ICC symbolic factorization
1752: */
1753: PetscErrorCode MatSeqSBAIJSetNumericFactorization_inplace(Mat B,PetscBool natural)
1754: {
1756:   PetscBool      flg = PETSC_FALSE;
1757:   PetscInt       bs  = B->rmap->bs;

1760:   PetscOptionsGetBool(((PetscObject)B)->options,((PetscObject)B)->prefix,"-mat_no_unroll",&flg,NULL);
1761:   if (flg) bs = 8;

1763:   if (!natural) {
1764:     switch (bs) {
1765:     case 1:
1766:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_inplace;
1767:       break;
1768:     case 2:
1769:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2;
1770:       break;
1771:     case 3:
1772:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3;
1773:       break;
1774:     case 4:
1775:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4;
1776:       break;
1777:     case 5:
1778:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5;
1779:       break;
1780:     case 6:
1781:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6;
1782:       break;
1783:     case 7:
1784:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7;
1785:       break;
1786:     default:
1787:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N;
1788:       break;
1789:     }
1790:   } else {
1791:     switch (bs) {
1792:     case 1:
1793:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering_inplace;
1794:       break;
1795:     case 2:
1796:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering;
1797:       break;
1798:     case 3:
1799:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering;
1800:       break;
1801:     case 4:
1802:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering;
1803:       break;
1804:     case 5:
1805:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering;
1806:       break;
1807:     case 6:
1808:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering;
1809:       break;
1810:     case 7:
1811:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering;
1812:       break;
1813:     default:
1814:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering;
1815:       break;
1816:     }
1817:   }
1818:   return(0);
1819: }

1821: PETSC_INTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqAIJ(Mat, MatType,MatReuse,Mat*);
1822: PETSC_INTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqBAIJ(Mat, MatType,MatReuse,Mat*);

1824: PETSC_INTERN PetscErrorCode MatGetFactor_seqsbaij_petsc(Mat A,MatFactorType ftype,Mat *B)
1825: {
1826:   PetscInt       n = A->rmap->n;

1830: #if defined(PETSC_USE_COMPLEX)
1831:   if (A->hermitian && !A->symmetric && (ftype == MAT_FACTOR_CHOLESKY||ftype == MAT_FACTOR_ICC)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Hermitian CHOLESKY or ICC Factor is not supported");
1832: #endif

1834:   MatCreate(PetscObjectComm((PetscObject)A),B);
1835:   MatSetSizes(*B,n,n,n,n);
1836:   if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) {
1837:     MatSetType(*B,MATSEQSBAIJ);
1838:     MatSeqSBAIJSetPreallocation(*B,A->rmap->bs,MAT_SKIP_ALLOCATION,NULL);

1840:     (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqSBAIJ;
1841:     (*B)->ops->iccfactorsymbolic      = MatICCFactorSymbolic_SeqSBAIJ;
1842:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");

1844:   (*B)->factortype = ftype;
1845:   PetscFree((*B)->solvertype);
1846:   PetscStrallocpy(MATSOLVERPETSC,&(*B)->solvertype);
1847:   return(0);
1848: }

1850: /*@C
1851:    MatSeqSBAIJGetArray - gives access to the array where the data for a MATSEQSBAIJ matrix is stored

1853:    Not Collective

1855:    Input Parameter:
1856: .  mat - a MATSEQSBAIJ matrix

1858:    Output Parameter:
1859: .   array - pointer to the data

1861:    Level: intermediate

1863: .seealso: MatSeqSBAIJRestoreArray(), MatSeqAIJGetArray(), MatSeqAIJRestoreArray()
1864: @*/
1865: PetscErrorCode  MatSeqSBAIJGetArray(Mat A,PetscScalar **array)
1866: {

1870:   PetscUseMethod(A,"MatSeqSBAIJGetArray_C",(Mat,PetscScalar**),(A,array));
1871:   return(0);
1872: }

1874: /*@C
1875:    MatSeqSBAIJRestoreArray - returns access to the array where the data for a MATSEQSBAIJ matrix is stored obtained by MatSeqSBAIJGetArray()

1877:    Not Collective

1879:    Input Parameters:
1880: +  mat - a MATSEQSBAIJ matrix
1881: -  array - pointer to the data

1883:    Level: intermediate

1885: .seealso: MatSeqSBAIJGetArray(), MatSeqAIJGetArray(), MatSeqAIJRestoreArray()
1886: @*/
1887: PetscErrorCode  MatSeqSBAIJRestoreArray(Mat A,PetscScalar **array)
1888: {

1892:   PetscUseMethod(A,"MatSeqSBAIJRestoreArray_C",(Mat,PetscScalar**),(A,array));
1893:   return(0);
1894: }

1896: /*MC
1897:   MATSEQSBAIJ - MATSEQSBAIJ = "seqsbaij" - A matrix type to be used for sequential symmetric block sparse matrices,
1898:   based on block compressed sparse row format.  Only the upper triangular portion of the matrix is stored.

1900:   For complex numbers by default this matrix is symmetric, NOT Hermitian symmetric. To make it Hermitian symmetric you
1901:   can call MatSetOption(Mat, MAT_HERMITIAN).

1903:   Options Database Keys:
1904:   . -mat_type seqsbaij - sets the matrix type to "seqsbaij" during a call to MatSetFromOptions()

1906:   Notes:
1907:     By default if you insert values into the lower triangular part of the matrix they are simply ignored (since they are not
1908:      stored and it is assumed they symmetric to the upper triangular). If you call MatSetOption(Mat,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_FALSE) or use
1909:      the options database -mat_ignore_lower_triangular false it will generate an error if you try to set a value in the lower triangular portion.

1911:     The number of rows in the matrix must be less than or equal to the number of columns

1913:   Level: beginner

1915:   .seealso: MatCreateSeqSBAIJ(), MatType, MATMPISBAIJ
1916: M*/
1917: PETSC_EXTERN PetscErrorCode MatCreate_SeqSBAIJ(Mat B)
1918: {
1919:   Mat_SeqSBAIJ   *b;
1921:   PetscMPIInt    size;
1922:   PetscBool      no_unroll = PETSC_FALSE,no_inode = PETSC_FALSE;

1925:   MPI_Comm_size(PetscObjectComm((PetscObject)B),&size);
1926:   if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Comm must be of size 1");

1928:   PetscNewLog(B,&b);
1929:   B->data = (void*)b;
1930:   PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));

1932:   B->ops->destroy    = MatDestroy_SeqSBAIJ;
1933:   B->ops->view       = MatView_SeqSBAIJ;
1934:   b->row             = 0;
1935:   b->icol            = 0;
1936:   b->reallocs        = 0;
1937:   b->saved_values    = 0;
1938:   b->inode.limit     = 5;
1939:   b->inode.max_limit = 5;

1941:   b->roworiented        = PETSC_TRUE;
1942:   b->nonew              = 0;
1943:   b->diag               = 0;
1944:   b->solve_work         = 0;
1945:   b->mult_work          = 0;
1946:   B->spptr              = 0;
1947:   B->info.nz_unneeded   = (PetscReal)b->maxnz*b->bs2;
1948:   b->keepnonzeropattern = PETSC_FALSE;

1950:   b->inew    = 0;
1951:   b->jnew    = 0;
1952:   b->anew    = 0;
1953:   b->a2anew  = 0;
1954:   b->permute = PETSC_FALSE;

1956:   b->ignore_ltriangular = PETSC_TRUE;

1958:   PetscOptionsGetBool(((PetscObject)B)->options,((PetscObject)B)->prefix,"-mat_ignore_lower_triangular",&b->ignore_ltriangular,NULL);

1960:   b->getrow_utriangular = PETSC_FALSE;

1962:   PetscOptionsGetBool(((PetscObject)B)->options,((PetscObject)B)->prefix,"-mat_getrow_uppertriangular",&b->getrow_utriangular,NULL);

1964:   PetscObjectComposeFunction((PetscObject)B,"MatSeqSBAIJGetArray_C",MatSeqSBAIJGetArray_SeqSBAIJ);
1965:   PetscObjectComposeFunction((PetscObject)B,"MatSeqSBAIJRestoreArray_C",MatSeqSBAIJRestoreArray_SeqSBAIJ);
1966:   PetscObjectComposeFunction((PetscObject)B,"MatStoreValues_C",MatStoreValues_SeqSBAIJ);
1967:   PetscObjectComposeFunction((PetscObject)B,"MatRetrieveValues_C",MatRetrieveValues_SeqSBAIJ);
1968:   PetscObjectComposeFunction((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C",MatSeqSBAIJSetColumnIndices_SeqSBAIJ);
1969:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_seqaij_C",MatConvert_SeqSBAIJ_SeqAIJ);
1970:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_seqbaij_C",MatConvert_SeqSBAIJ_SeqBAIJ);
1971:   PetscObjectComposeFunction((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",MatSeqSBAIJSetPreallocation_SeqSBAIJ);
1972:   PetscObjectComposeFunction((PetscObject)B,"MatSeqSBAIJSetPreallocationCSR_C",MatSeqSBAIJSetPreallocationCSR_SeqSBAIJ);
1973: #if defined(PETSC_HAVE_ELEMENTAL)
1974:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_elemental_C",MatConvert_SeqSBAIJ_Elemental);
1975: #endif

1977:   B->symmetric                  = PETSC_TRUE;
1978:   B->structurally_symmetric     = PETSC_TRUE;
1979:   B->symmetric_set              = PETSC_TRUE;
1980:   B->structurally_symmetric_set = PETSC_TRUE;
1981:   B->symmetric_eternal          = PETSC_TRUE;
1982: #if defined(PETSC_USE_COMPLEX)
1983:   B->hermitian                  = PETSC_FALSE;
1984:   B->hermitian_set              = PETSC_FALSE;
1985: #else
1986:   B->hermitian                  = PETSC_TRUE;
1987:   B->hermitian_set              = PETSC_TRUE;
1988: #endif

1990:   PetscObjectChangeTypeName((PetscObject)B,MATSEQSBAIJ);

1992:   PetscOptionsBegin(PetscObjectComm((PetscObject)B),((PetscObject)B)->prefix,"Options for SEQSBAIJ matrix","Mat");
1993:   PetscOptionsBool("-mat_no_unroll","Do not optimize for inodes (slower)",NULL,no_unroll,&no_unroll,NULL);
1994:   if (no_unroll) {
1995:     PetscInfo(B,"Not using Inode routines due to -mat_no_unroll\n");
1996:   }
1997:   PetscOptionsBool("-mat_no_inode","Do not optimize for inodes (slower)",NULL,no_inode,&no_inode,NULL);
1998:   if (no_inode) {
1999:     PetscInfo(B,"Not using Inode routines due to -mat_no_inode\n");
2000:   }
2001:   PetscOptionsInt("-mat_inode_limit","Do not use inodes larger then this value",NULL,b->inode.limit,&b->inode.limit,NULL);
2002:   PetscOptionsEnd();
2003:   b->inode.use = (PetscBool)(!(no_unroll || no_inode));
2004:   if (b->inode.limit > b->inode.max_limit) b->inode.limit = b->inode.max_limit;
2005:   return(0);
2006: }

2008: /*@C
2009:    MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block
2010:    compressed row) format.  For good matrix assembly performance the
2011:    user should preallocate the matrix storage by setting the parameter nz
2012:    (or the array nnz).  By setting these parameters accurately, performance
2013:    during matrix assembly can be increased by more than a factor of 50.

2015:    Collective on Mat

2017:    Input Parameters:
2018: +  B - the symmetric matrix
2019: .  bs - size of block, the blocks are ALWAYS square. One can use MatSetBlockSizes() to set a different row and column blocksize but the row
2020:           blocksize always defines the size of the blocks. The column blocksize sets the blocksize of the vectors obtained with MatCreateVecs()
2021: .  nz - number of block nonzeros per block row (same for all rows)
2022: -  nnz - array containing the number of block nonzeros in the upper triangular plus
2023:          diagonal portion of each block (possibly different for each block row) or NULL

2025:    Options Database Keys:
2026: +   -mat_no_unroll - uses code that does not unroll the loops in the
2027:                      block calculations (much slower)
2028: -   -mat_block_size - size of the blocks to use (only works if a negative bs is passed in

2030:    Level: intermediate

2032:    Notes:
2033:    Specify the preallocated storage with either nz or nnz (not both).
2034:    Set nz=PETSC_DEFAULT and nnz=NULL for PETSc to control dynamic memory
2035:    allocation.  See Users-Manual: ch_mat for details.

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

2042:    If the nnz parameter is given then the nz parameter is ignored


2045: .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateSBAIJ()
2046: @*/
2047: PetscErrorCode  MatSeqSBAIJSetPreallocation(Mat B,PetscInt bs,PetscInt nz,const PetscInt nnz[])
2048: {

2055:   PetscTryMethod(B,"MatSeqSBAIJSetPreallocation_C",(Mat,PetscInt,PetscInt,const PetscInt[]),(B,bs,nz,nnz));
2056:   return(0);
2057: }

2059: /*@C
2060:    MatSeqSBAIJSetPreallocationCSR - Creates a sparse parallel matrix in SBAIJ format using the given nonzero structure and (optional) numerical values

2062:    Input Parameters:
2063: +  B - the matrix
2064: .  bs - size of block, the blocks are ALWAYS square.
2065: .  i - the indices into j for the start of each local row (starts with zero)
2066: .  j - the column indices for each local row (starts with zero) these must be sorted for each row
2067: -  v - optional values in the matrix

2069:    Level: advanced

2071:    Notes:
2072:    The order of the entries in values is specified by the MatOption MAT_ROW_ORIENTED.  For example, C programs
2073:    may want to use the default MAT_ROW_ORIENTED=PETSC_TRUE and use an array v[nnz][bs][bs] where the second index is
2074:    over rows within a block and the last index is over columns within a block row.  Fortran programs will likely set
2075:    MAT_ROW_ORIENTED=PETSC_FALSE and use a Fortran array v(bs,bs,nnz) in which the first index is over rows within a
2076:    block column and the second index is over columns within a block.

2078:    Any entries below the diagonal are ignored

2080:    Though this routine has Preallocation() in the name it also sets the exact nonzero locations of the matrix entries
2081:    and usually the numerical values as well

2083: .seealso: MatCreate(), MatCreateSeqSBAIJ(), MatSetValuesBlocked(), MatSeqSBAIJSetPreallocation(), MATSEQSBAIJ
2084: @*/
2085: PetscErrorCode MatSeqSBAIJSetPreallocationCSR(Mat B,PetscInt bs,const PetscInt i[],const PetscInt j[], const PetscScalar v[])
2086: {

2093:   PetscTryMethod(B,"MatSeqSBAIJSetPreallocationCSR_C",(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]),(B,bs,i,j,v));
2094:   return(0);
2095: }

2097: /*@C
2098:    MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block
2099:    compressed row) format.  For good matrix assembly performance the
2100:    user should preallocate the matrix storage by setting the parameter nz
2101:    (or the array nnz).  By setting these parameters accurately, performance
2102:    during matrix assembly can be increased by more than a factor of 50.

2104:    Collective

2106:    Input Parameters:
2107: +  comm - MPI communicator, set to PETSC_COMM_SELF
2108: .  bs - size of block, the blocks are ALWAYS square. One can use MatSetBlockSizes() to set a different row and column blocksize but the row
2109:           blocksize always defines the size of the blocks. The column blocksize sets the blocksize of the vectors obtained with MatCreateVecs()
2110: .  m - number of rows, or number of columns
2111: .  nz - number of block nonzeros per block row (same for all rows)
2112: -  nnz - array containing the number of block nonzeros in the upper triangular plus
2113:          diagonal portion of each block (possibly different for each block row) or NULL

2115:    Output Parameter:
2116: .  A - the symmetric matrix

2118:    Options Database Keys:
2119: +   -mat_no_unroll - uses code that does not unroll the loops in the
2120:                      block calculations (much slower)
2121: -   -mat_block_size - size of the blocks to use

2123:    Level: intermediate

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

2129:    Notes:
2130:    The number of rows and columns must be divisible by blocksize.
2131:    This matrix type does not support complex Hermitian operation.

2133:    Specify the preallocated storage with either nz or nnz (not both).
2134:    Set nz=PETSC_DEFAULT and nnz=NULL for PETSc to control dynamic memory
2135:    allocation.  See Users-Manual: ch_mat for details.

2137:    If the nnz parameter is given then the nz parameter is ignored

2139: .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateSBAIJ()
2140: @*/
2141: PetscErrorCode  MatCreateSeqSBAIJ(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A)
2142: {

2146:   MatCreate(comm,A);
2147:   MatSetSizes(*A,m,n,m,n);
2148:   MatSetType(*A,MATSEQSBAIJ);
2149:   MatSeqSBAIJSetPreallocation(*A,bs,nz,(PetscInt*)nnz);
2150:   return(0);
2151: }

2153: PetscErrorCode MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
2154: {
2155:   Mat            C;
2156:   Mat_SeqSBAIJ   *c,*a = (Mat_SeqSBAIJ*)A->data;
2158:   PetscInt       i,mbs = a->mbs,nz = a->nz,bs2 =a->bs2;

2161:   if (a->i[mbs] != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Corrupt matrix");

2163:   *B   = 0;
2164:   MatCreate(PetscObjectComm((PetscObject)A),&C);
2165:   MatSetSizes(C,A->rmap->N,A->cmap->n,A->rmap->N,A->cmap->n);
2166:   MatSetBlockSizesFromMats(C,A,A);
2167:   MatSetType(C,MATSEQSBAIJ);
2168:   c    = (Mat_SeqSBAIJ*)C->data;

2170:   C->preallocated       = PETSC_TRUE;
2171:   C->factortype         = A->factortype;
2172:   c->row                = 0;
2173:   c->icol               = 0;
2174:   c->saved_values       = 0;
2175:   c->keepnonzeropattern = a->keepnonzeropattern;
2176:   C->assembled          = PETSC_TRUE;

2178:   PetscLayoutReference(A->rmap,&C->rmap);
2179:   PetscLayoutReference(A->cmap,&C->cmap);
2180:   c->bs2 = a->bs2;
2181:   c->mbs = a->mbs;
2182:   c->nbs = a->nbs;

2184:   if (cpvalues == MAT_SHARE_NONZERO_PATTERN) {
2185:     c->imax           = a->imax;
2186:     c->ilen           = a->ilen;
2187:     c->free_imax_ilen = PETSC_FALSE;
2188:   } else {
2189:     PetscMalloc2((mbs+1),&c->imax,(mbs+1),&c->ilen);
2190:     PetscLogObjectMemory((PetscObject)C,2*(mbs+1)*sizeof(PetscInt));
2191:     for (i=0; i<mbs; i++) {
2192:       c->imax[i] = a->imax[i];
2193:       c->ilen[i] = a->ilen[i];
2194:     }
2195:     c->free_imax_ilen = PETSC_TRUE;
2196:   }

2198:   /* allocate the matrix space */
2199:   if (cpvalues == MAT_SHARE_NONZERO_PATTERN) {
2200:     PetscMalloc1(bs2*nz,&c->a);
2201:     PetscLogObjectMemory((PetscObject)C,nz*bs2*sizeof(MatScalar));
2202:     c->i            = a->i;
2203:     c->j            = a->j;
2204:     c->singlemalloc = PETSC_FALSE;
2205:     c->free_a       = PETSC_TRUE;
2206:     c->free_ij      = PETSC_FALSE;
2207:     c->parent       = A;
2208:     PetscObjectReference((PetscObject)A);
2209:     MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
2210:     MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
2211:   } else {
2212:     PetscMalloc3(bs2*nz,&c->a,nz,&c->j,mbs+1,&c->i);
2213:     PetscArraycpy(c->i,a->i,mbs+1);
2214:     PetscLogObjectMemory((PetscObject)C,(mbs+1)*sizeof(PetscInt) + nz*(bs2*sizeof(MatScalar) + sizeof(PetscInt)));
2215:     c->singlemalloc = PETSC_TRUE;
2216:     c->free_a       = PETSC_TRUE;
2217:     c->free_ij      = PETSC_TRUE;
2218:   }
2219:   if (mbs > 0) {
2220:     if (cpvalues != MAT_SHARE_NONZERO_PATTERN) {
2221:       PetscArraycpy(c->j,a->j,nz);
2222:     }
2223:     if (cpvalues == MAT_COPY_VALUES) {
2224:       PetscArraycpy(c->a,a->a,bs2*nz);
2225:     } else {
2226:       PetscArrayzero(c->a,bs2*nz);
2227:     }
2228:     if (a->jshort) {
2229:       /* cannot share jshort, it is reallocated in MatAssemblyEnd_SeqSBAIJ() */
2230:       /* if the parent matrix is reassembled, this child matrix will never notice */
2231:       PetscMalloc1(nz,&c->jshort);
2232:       PetscLogObjectMemory((PetscObject)C,nz*sizeof(unsigned short));
2233:       PetscArraycpy(c->jshort,a->jshort,nz);

2235:       c->free_jshort = PETSC_TRUE;
2236:     }
2237:   }

2239:   c->roworiented = a->roworiented;
2240:   c->nonew       = a->nonew;

2242:   if (a->diag) {
2243:     if (cpvalues == MAT_SHARE_NONZERO_PATTERN) {
2244:       c->diag      = a->diag;
2245:       c->free_diag = PETSC_FALSE;
2246:     } else {
2247:       PetscMalloc1(mbs,&c->diag);
2248:       PetscLogObjectMemory((PetscObject)C,mbs*sizeof(PetscInt));
2249:       for (i=0; i<mbs; i++) c->diag[i] = a->diag[i];
2250:       c->free_diag = PETSC_TRUE;
2251:     }
2252:   }
2253:   c->nz         = a->nz;
2254:   c->maxnz      = a->nz; /* Since we allocate exactly the right amount */
2255:   c->solve_work = 0;
2256:   c->mult_work  = 0;

2258:   *B   = C;
2259:   PetscFunctionListDuplicate(((PetscObject)A)->qlist,&((PetscObject)C)->qlist);
2260:   return(0);
2261: }

2263: /* Used for both SeqBAIJ and SeqSBAIJ matrices */
2264: #define MatLoad_SeqSBAIJ_Binary MatLoad_SeqBAIJ_Binary

2266: PetscErrorCode MatLoad_SeqSBAIJ(Mat mat,PetscViewer viewer)
2267: {
2269:   PetscBool      isbinary;

2272:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
2273:   if (!isbinary) SETERRQ2(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer type %s not yet supported for reading %s matrices",((PetscObject)viewer)->type_name,((PetscObject)mat)->type_name);
2274:   MatLoad_SeqSBAIJ_Binary(mat,viewer);
2275:   return(0);
2276: }

2278: /*@
2279:      MatCreateSeqSBAIJWithArrays - Creates an sequential SBAIJ matrix using matrix elements
2280:               (upper triangular entries in CSR format) provided by the user.

2282:      Collective

2284:    Input Parameters:
2285: +  comm - must be an MPI communicator of size 1
2286: .  bs - size of block
2287: .  m - number of rows
2288: .  n - number of columns
2289: .  i - row indices; that is i[0] = 0, i[row] = i[row-1] + number of block elements in that row block row of the matrix
2290: .  j - column indices
2291: -  a - matrix values

2293:    Output Parameter:
2294: .  mat - the matrix

2296:    Level: advanced

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

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

2304:        The i and j indices are 0 based

2306:        When block size is greater than 1 the matrix values must be stored using the SBAIJ storage format (see the SBAIJ code to determine this). For block size of 1
2307:        it is the regular CSR format excluding the lower triangular elements.

2309: .seealso: MatCreate(), MatCreateSBAIJ(), MatCreateSeqSBAIJ()

2311: @*/
2312: PetscErrorCode  MatCreateSeqSBAIJWithArrays(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt i[],PetscInt j[],PetscScalar a[],Mat *mat)
2313: {
2315:   PetscInt       ii;
2316:   Mat_SeqSBAIJ   *sbaij;

2319:   if (bs != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"block size %D > 1 is not supported yet",bs);
2320:   if (m > 0 && i[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"i (row indices) must start with 0");

2322:   MatCreate(comm,mat);
2323:   MatSetSizes(*mat,m,n,m,n);
2324:   MatSetType(*mat,MATSEQSBAIJ);
2325:   MatSeqSBAIJSetPreallocation(*mat,bs,MAT_SKIP_ALLOCATION,0);
2326:   sbaij = (Mat_SeqSBAIJ*)(*mat)->data;
2327:   PetscMalloc2(m,&sbaij->imax,m,&sbaij->ilen);
2328:   PetscLogObjectMemory((PetscObject)*mat,2*m*sizeof(PetscInt));

2330:   sbaij->i = i;
2331:   sbaij->j = j;
2332:   sbaij->a = a;

2334:   sbaij->singlemalloc   = PETSC_FALSE;
2335:   sbaij->nonew          = -1;             /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/
2336:   sbaij->free_a         = PETSC_FALSE;
2337:   sbaij->free_ij        = PETSC_FALSE;
2338:   sbaij->free_imax_ilen = PETSC_TRUE;

2340:   for (ii=0; ii<m; ii++) {
2341:     sbaij->ilen[ii] = sbaij->imax[ii] = i[ii+1] - i[ii];
2342: #if defined(PETSC_USE_DEBUG)
2343:     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]);
2344: #endif
2345:   }
2346: #if defined(PETSC_USE_DEBUG)
2347:   for (ii=0; ii<sbaij->i[m]; ii++) {
2348:     if (j[ii] < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column index at location = %d index = %d",ii,j[ii]);
2349:     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]);
2350:   }
2351: #endif

2353:   MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
2354:   MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
2355:   return(0);
2356: }

2358: PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqSBAIJ(MPI_Comm comm,Mat inmat,PetscInt n,MatReuse scall,Mat *outmat)
2359: {
2361:   PetscMPIInt    size;

2364:   MPI_Comm_size(comm,&size);
2365:   if (size == 1 && scall == MAT_REUSE_MATRIX) {
2366:     MatCopy(inmat,*outmat,SAME_NONZERO_PATTERN);
2367:   } else {
2368:     MatCreateMPIMatConcatenateSeqMat_MPISBAIJ(comm,inmat,n,scall,outmat);
2369:   }
2370:   return(0);
2371: }