Actual source code: sbaij.c

petsc-master 2014-12-27
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>         /*I "petscmat.h" I*/
  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: extern PetscErrorCode MatSeqSBAIJSetNumericFactorization_inplace(Mat,PetscBool);

 16: /*
 17:      Checks for missing diagonals
 18: */
 21: PetscErrorCode MatMissingDiagonal_SeqSBAIJ(Mat A,PetscBool  *missing,PetscInt *dd)
 22: {
 23:   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
 25:   PetscInt       *diag,*ii = a->i,i;

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

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

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

 75: static PetscErrorCode MatGetRowIJ_SeqSBAIJ(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *nn,const PetscInt *inia[],const PetscInt *inja[],PetscBool  *done)
 76: {
 77:   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
 78:   PetscInt       i,j,n = a->mbs,nz = a->i[n],bs = A->rmap->bs;
 79:   PetscInt       **ia = (PetscInt**)inia,**ja = (PetscInt**)inja;

 83:   *nn = n;
 84:   if (!ia) return(0);
 85:   if (!blockcompressed) {
 86:     /* malloc & create the natural set of indices */
 87:     PetscMalloc2((n+1)*bs,ia,nz*bs,ja);
 88:     for (i=0; i<n+1; i++) {
 89:       for (j=0; j<bs; j++) {
 90:         *ia[i*bs+j] = a->i[i]*bs+j+oshift;
 91:       }
 92:     }
 93:     for (i=0; i<nz; i++) {
 94:       for (j=0; j<bs; j++) {
 95:         *ja[i*bs+j] = a->j[i]*bs+j+oshift;
 96:       }
 97:     }
 98:   } else { /* blockcompressed */
 99:     if (oshift == 1) {
100:       /* temporarily add 1 to i and j indices */
101:       for (i=0; i<nz; i++) a->j[i]++;
102:       for (i=0; i<n+1; i++) a->i[i]++;
103:     }
104:     *ia = a->i; *ja = a->j;
105:   }
106:   return(0);
107: }

111: static PetscErrorCode MatRestoreRowIJ_SeqSBAIJ(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *nn,const PetscInt *ia[],const PetscInt *ja[],PetscBool  *done)
112: {
113:   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
114:   PetscInt       i,n = a->mbs,nz = a->i[n];

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

120:   if (!blockcompressed) {
121:     PetscFree2(*ia,*ja);
122:   } else if (oshift == 1) { /* blockcompressed */
123:     for (i=0; i<nz; i++) a->j[i]--;
124:     for (i=0; i<n+1; i++) a->i[i]--;
125:   }
126:   return(0);
127: }

131: PetscErrorCode MatDestroy_SeqSBAIJ(Mat A)
132: {
133:   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;

137: #if defined(PETSC_USE_LOG)
138:   PetscLogObjectState((PetscObject)A,"Rows=%D, NZ=%D",A->rmap->N,a->nz);
139: #endif
140:   MatSeqXAIJFreeAIJ(A,&a->a,&a->j,&a->i);
141:   if (a->free_diag) {PetscFree(a->diag);}
142:   ISDestroy(&a->row);
143:   ISDestroy(&a->col);
144:   ISDestroy(&a->icol);
145:   PetscFree(a->idiag);
146:   PetscFree(a->inode.size);
147:   if (a->free_imax_ilen) {PetscFree2(a->imax,a->ilen);}
148:   PetscFree(a->solve_work);
149:   PetscFree(a->sor_work);
150:   PetscFree(a->solves_work);
151:   PetscFree(a->mult_work);
152:   PetscFree(a->saved_values);
153:   if (a->free_jshort) {PetscFree(a->jshort);}
154:   PetscFree(a->inew);
155:   MatDestroy(&a->parent);
156:   PetscFree(A->data);

158:   PetscObjectChangeTypeName((PetscObject)A,0);
159:   PetscObjectComposeFunction((PetscObject)A,"MatStoreValues_C",NULL);
160:   PetscObjectComposeFunction((PetscObject)A,"MatRetrieveValues_C",NULL);
161:   PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetColumnIndices_C",NULL);
162:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqaij_C",NULL);
163:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqbaij_C",NULL);
164:   PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetPreallocation_C",NULL);
165:   PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetPreallocationCSR_C",NULL);
166:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqsbstrm_C",NULL);
167:   return(0);
168: }

172: PetscErrorCode MatSetOption_SeqSBAIJ(Mat A,MatOption op,PetscBool flg)
173: {
174:   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;

178:   switch (op) {
179:   case MAT_ROW_ORIENTED:
180:     a->roworiented = flg;
181:     break;
182:   case MAT_KEEP_NONZERO_PATTERN:
183:     a->keepnonzeropattern = flg;
184:     break;
185:   case MAT_NEW_NONZERO_LOCATIONS:
186:     a->nonew = (flg ? 0 : 1);
187:     break;
188:   case MAT_NEW_NONZERO_LOCATION_ERR:
189:     a->nonew = (flg ? -1 : 0);
190:     break;
191:   case MAT_NEW_NONZERO_ALLOCATION_ERR:
192:     a->nonew = (flg ? -2 : 0);
193:     break;
194:   case MAT_UNUSED_NONZERO_LOCATION_ERR:
195:     a->nounused = (flg ? -1 : 0);
196:     break;
197:   case MAT_NEW_DIAGONALS:
198:   case MAT_IGNORE_OFF_PROC_ENTRIES:
199:   case MAT_USE_HASH_TABLE:
200:     PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
201:     break;
202:   case MAT_HERMITIAN:
203:     if (!A->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatAssemblyEnd() first");
204:     if (A->cmap->n < 65536 && A->cmap->bs == 1) {
205:       A->ops->mult = MatMult_SeqSBAIJ_1_Hermitian_ushort;
206:     } else if (A->cmap->bs == 1) {
207:       A->ops->mult = MatMult_SeqSBAIJ_1_Hermitian;
208:     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for Hermitian with block size greater than 1");
209:     break;
210:   case MAT_SPD:
211:     /* These options are handled directly by MatSetOption() */
212:     break;
213:   case MAT_SYMMETRIC:
214:   case MAT_STRUCTURALLY_SYMMETRIC:
215:   case MAT_SYMMETRY_ETERNAL:
216:     /* These options are handled directly by MatSetOption() */
217:     break;
218:   case MAT_IGNORE_LOWER_TRIANGULAR:
219:     a->ignore_ltriangular = flg;
220:     break;
221:   case MAT_ERROR_LOWER_TRIANGULAR:
222:     a->ignore_ltriangular = flg;
223:     break;
224:   case MAT_GETROW_UPPERTRIANGULAR:
225:     a->getrow_utriangular = flg;
226:     break;
227:   default:
228:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %d",op);
229:   }
230:   return(0);
231: }

235: PetscErrorCode MatGetRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
236: {
237:   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;

241:   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()");

243:   /* Get the upper triangular part of the row */
244:   MatGetRow_SeqBAIJ_private(A,row,nz,idx,v,a->i,a->j,a->a);
245:   return(0);
246: }

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

255:   if (idx) {PetscFree(*idx);}
256:   if (v)   {PetscFree(*v);}
257:   return(0);
258: }

262: PetscErrorCode MatGetRowUpperTriangular_SeqSBAIJ(Mat A)
263: {
264:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;

267:   a->getrow_utriangular = PETSC_TRUE;
268:   return(0);
269: }
272: PetscErrorCode MatRestoreRowUpperTriangular_SeqSBAIJ(Mat A)
273: {
274:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;

277:   a->getrow_utriangular = PETSC_FALSE;
278:   return(0);
279: }

283: PetscErrorCode MatTranspose_SeqSBAIJ(Mat A,MatReuse reuse,Mat *B)
284: {

288:   if (reuse == MAT_INITIAL_MATRIX || *B != A) {
289:     MatDuplicate(A,MAT_COPY_VALUES,B);
290:   }
291:   return(0);
292: }

296: PetscErrorCode MatView_SeqSBAIJ_ASCII(Mat A,PetscViewer viewer)
297: {
298:   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
299:   PetscErrorCode    ierr;
300:   PetscInt          i,j,bs = A->rmap->bs,k,l,bs2=a->bs2;
301:   PetscViewerFormat format;
302:   PetscInt          *diag;

305:   PetscViewerGetFormat(viewer,&format);
306:   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
307:     PetscViewerASCIIPrintf(viewer,"  block size is %D\n",bs);
308:   } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
309:     Mat aij;
310:     if (A->factortype && bs>1) {
311:       PetscPrintf(PETSC_COMM_SELF,"Warning: matrix is factored with bs>1. MatView() with PETSC_VIEWER_ASCII_MATLAB is not supported and ignored!\n");
312:       return(0);
313:     }
314:     MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&aij);
315:     MatView(aij,viewer);
316:     MatDestroy(&aij);
317:   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
318:     PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
319:     for (i=0; i<a->mbs; i++) {
320:       for (j=0; j<bs; j++) {
321:         PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);
322:         for (k=a->i[i]; k<a->i[i+1]; k++) {
323:           for (l=0; l<bs; l++) {
324: #if defined(PETSC_USE_COMPLEX)
325:             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
326:               PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",bs*a->j[k]+l,
327:                                             (double)PetscRealPart(a->a[bs2*k + l*bs + j]),(double)PetscImaginaryPart(a->a[bs2*k + l*bs + j]));
328:             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
329:               PetscViewerASCIIPrintf(viewer," (%D, %g - %g i) ",bs*a->j[k]+l,
330:                                             (double)PetscRealPart(a->a[bs2*k + l*bs + j]),-(double)PetscImaginaryPart(a->a[bs2*k + l*bs + j]));
331:             } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
332:               PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,(double)PetscRealPart(a->a[bs2*k + l*bs + j]));
333:             }
334: #else
335:             if (a->a[bs2*k + l*bs + j] != 0.0) {
336:               PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,(double)a->a[bs2*k + l*bs + j]);
337:             }
338: #endif
339:           }
340:         }
341:         PetscViewerASCIIPrintf(viewer,"\n");
342:       }
343:     }
344:     PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
345:   } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
346:     return(0);
347:   } else {
348:     PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
349:     if (A->factortype) { /* for factored matrix */
350:       if (bs>1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"matrix is factored with bs>1. Not implemented yet");

352:       diag=a->diag;
353:       for (i=0; i<a->mbs; i++) { /* for row block i */
354:         PetscViewerASCIIPrintf(viewer,"row %D:",i);
355:         /* diagonal entry */
356: #if defined(PETSC_USE_COMPLEX)
357:         if (PetscImaginaryPart(a->a[diag[i]]) > 0.0) {
358:           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]]));
359:         } else if (PetscImaginaryPart(a->a[diag[i]]) < 0.0) {
360:           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]]));
361:         } else {
362:           PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[diag[i]],(double)PetscRealPart(1.0/a->a[diag[i]]));
363:         }
364: #else
365:         PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[diag[i]],(double)(1.0/a->a[diag[i]]));
366: #endif
367:         /* off-diagonal entries */
368:         for (k=a->i[i]; k<a->i[i+1]-1; k++) {
369: #if defined(PETSC_USE_COMPLEX)
370:           if (PetscImaginaryPart(a->a[k]) > 0.0) {
371:             PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",bs*a->j[k],(double)PetscRealPart(a->a[k]),(double)PetscImaginaryPart(a->a[k]));
372:           } else if (PetscImaginaryPart(a->a[k]) < 0.0) {
373:             PetscViewerASCIIPrintf(viewer," (%D, %g - %g i) ",bs*a->j[k],(double)PetscRealPart(a->a[k]),-(double)PetscImaginaryPart(a->a[k]));
374:           } else {
375:             PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k],(double)PetscRealPart(a->a[k]));
376:           }
377: #else
378:           PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[k],(double)a->a[k]);
379: #endif
380:         }
381:         PetscViewerASCIIPrintf(viewer,"\n");
382:       }

384:     } else { /* for non-factored matrix */
385:       for (i=0; i<a->mbs; i++) { /* for row block i */
386:         for (j=0; j<bs; j++) {   /* for row bs*i + j */
387:           PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);
388:           for (k=a->i[i]; k<a->i[i+1]; k++) { /* for column block */
389:             for (l=0; l<bs; l++) {            /* for column */
390: #if defined(PETSC_USE_COMPLEX)
391:               if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) {
392:                 PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",bs*a->j[k]+l,
393:                                               (double)PetscRealPart(a->a[bs2*k + l*bs + j]),(double)PetscImaginaryPart(a->a[bs2*k + l*bs + j]));
394:               } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) {
395:                 PetscViewerASCIIPrintf(viewer," (%D, %g - %g i) ",bs*a->j[k]+l,
396:                                               (double)PetscRealPart(a->a[bs2*k + l*bs + j]),-(double)PetscImaginaryPart(a->a[bs2*k + l*bs + j]));
397:               } else {
398:                 PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,(double)PetscRealPart(a->a[bs2*k + l*bs + j]));
399:               }
400: #else
401:               PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,(double)a->a[bs2*k + l*bs + j]);
402: #endif
403:             }
404:           }
405:           PetscViewerASCIIPrintf(viewer,"\n");
406:         }
407:       }
408:     }
409:     PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
410:   }
411:   PetscViewerFlush(viewer);
412:   return(0);
413: }

415: #include <petscdraw.h>
418: static PetscErrorCode MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw,void *Aa)
419: {
420:   Mat            A = (Mat) Aa;
421:   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data;
423:   PetscInt       row,i,j,k,l,mbs=a->mbs,color,bs=A->rmap->bs,bs2=a->bs2;
424:   PetscMPIInt    rank;
425:   PetscReal      xl,yl,xr,yr,x_l,x_r,y_l,y_r;
426:   MatScalar      *aa;
427:   MPI_Comm       comm;
428:   PetscViewer    viewer;

431:   /*
432:     This is nasty. If this is called from an originally parallel matrix
433:     then all processes call this,but only the first has the matrix so the
434:     rest should return immediately.
435:   */
436:   PetscObjectGetComm((PetscObject)draw,&comm);
437:   MPI_Comm_rank(comm,&rank);
438:   if (rank) return(0);

440:   PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);

442:   PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);
443:   PetscDrawString(draw, .3*(xl+xr), .3*(yl+yr), PETSC_DRAW_BLACK, "symmetric");

445:   /* loop over matrix elements drawing boxes */
446:   color = PETSC_DRAW_BLUE;
447:   for (i=0,row=0; i<mbs; i++,row+=bs) {
448:     for (j=a->i[i]; j<a->i[i+1]; j++) {
449:       y_l = A->rmap->N - row - 1.0; y_r = y_l + 1.0;
450:       x_l = a->j[j]*bs; x_r = x_l + 1.0;
451:       aa  = a->a + j*bs2;
452:       for (k=0; k<bs; k++) {
453:         for (l=0; l<bs; l++) {
454:           if (PetscRealPart(*aa++) >=  0.) continue;
455:           PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
456:         }
457:       }
458:     }
459:   }
460:   color = PETSC_DRAW_CYAN;
461:   for (i=0,row=0; i<mbs; i++,row+=bs) {
462:     for (j=a->i[i]; j<a->i[i+1]; j++) {
463:       y_l = A->rmap->N - row - 1.0; y_r = y_l + 1.0;
464:       x_l = a->j[j]*bs; x_r = x_l + 1.0;
465:       aa = a->a + j*bs2;
466:       for (k=0; k<bs; k++) {
467:         for (l=0; l<bs; l++) {
468:           if (PetscRealPart(*aa++) != 0.) continue;
469:           PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
470:         }
471:       }
472:     }
473:   }

475:   color = PETSC_DRAW_RED;
476:   for (i=0,row=0; i<mbs; i++,row+=bs) {
477:     for (j=a->i[i]; j<a->i[i+1]; j++) {
478:       y_l = A->rmap->N - row - 1.0; y_r = y_l + 1.0;
479:       x_l = a->j[j]*bs; x_r = x_l + 1.0;
480:       aa = a->a + j*bs2;
481:       for (k=0; k<bs; k++) {
482:         for (l=0; l<bs; l++) {
483:           if (PetscRealPart(*aa++) <= 0.) continue;
484:           PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
485:         }
486:       }
487:     }
488:   }
489:   return(0);
490: }

494: static PetscErrorCode MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer)
495: {
497:   PetscReal      xl,yl,xr,yr,w,h;
498:   PetscDraw      draw;
499:   PetscBool      isnull;

502:   PetscViewerDrawGetDraw(viewer,0,&draw);
503:   PetscDrawIsNull(draw,&isnull); if (isnull) return(0);

505:   PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);
506:   xr   = A->rmap->N; yr = A->rmap->N; h = yr/10.0; w = xr/10.0;
507:   xr  += w;    yr += h;  xl = -w;     yl = -h;
508:   PetscDrawSetCoordinates(draw,xl,yl,xr,yr);
509:   PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);
510:   PetscObjectCompose((PetscObject)A,"Zoomviewer",NULL);
511:   return(0);
512: }

516: PetscErrorCode MatView_SeqSBAIJ(Mat A,PetscViewer viewer)
517: {
519:   PetscBool      iascii,isdraw;
520:   FILE           *file = 0;

523:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
524:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
525:   if (iascii) {
526:     MatView_SeqSBAIJ_ASCII(A,viewer);
527:   } else if (isdraw) {
528:     MatView_SeqSBAIJ_Draw(A,viewer);
529:   } else {
530:     Mat B;
531:     MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&B);
532:     MatView(B,viewer);
533:     MatDestroy(&B);
534:     PetscViewerBinaryGetInfoPointer(viewer,&file);
535:     if (file) {
536:       fprintf(file,"-matload_block_size %d\n",(int)A->rmap->bs);
537:     }
538:   }
539:   return(0);
540: }


545: PetscErrorCode MatGetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],PetscScalar v[])
546: {
547:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
548:   PetscInt     *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
549:   PetscInt     *ai = a->i,*ailen = a->ilen;
550:   PetscInt     brow,bcol,ridx,cidx,bs=A->rmap->bs,bs2=a->bs2;
551:   MatScalar    *ap,*aa = a->a;

554:   for (k=0; k<m; k++) { /* loop over rows */
555:     row = im[k]; brow = row/bs;
556:     if (row < 0) {v += n; continue;} /* SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",row); */
557:     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);
558:     rp   = aj + ai[brow]; ap = aa + bs2*ai[brow];
559:     nrow = ailen[brow];
560:     for (l=0; l<n; l++) { /* loop over columns */
561:       if (in[l] < 0) {v++; continue;} /* SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %D",in[l]); */
562:       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);
563:       col  = in[l];
564:       bcol = col/bs;
565:       cidx = col%bs;
566:       ridx = row%bs;
567:       high = nrow;
568:       low  = 0; /* assume unsorted */
569:       while (high-low > 5) {
570:         t = (low+high)/2;
571:         if (rp[t] > bcol) high = t;
572:         else              low  = t;
573:       }
574:       for (i=low; i<high; i++) {
575:         if (rp[i] > bcol) break;
576:         if (rp[i] == bcol) {
577:           *v++ = ap[bs2*i+bs*cidx+ridx];
578:           goto finished;
579:         }
580:       }
581:       *v++ = 0.0;
582: finished:;
583:     }
584:   }
585:   return(0);
586: }


591: PetscErrorCode MatSetValuesBlocked_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
592: {
593:   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
594:   PetscErrorCode    ierr;
595:   PetscInt          *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,lastcol = -1;
596:   PetscInt          *imax      =a->imax,*ai=a->i,*ailen=a->ilen;
597:   PetscInt          *aj        =a->j,nonew=a->nonew,bs2=a->bs2,bs=A->rmap->bs,stepval;
598:   PetscBool         roworiented=a->roworiented;
599:   const PetscScalar *value     = v;
600:   MatScalar         *ap,*aa = a->a,*bap;

603:   if (roworiented) stepval = (n-1)*bs;
604:   else stepval = (m-1)*bs;

606:   for (k=0; k<m; k++) { /* loop over added rows */
607:     row = im[k];
608:     if (row < 0) continue;
609: #if defined(PETSC_USE_DEBUG)
610:     if (row >= a->mbs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,a->mbs-1);
611: #endif
612:     rp   = aj + ai[row];
613:     ap   = aa + bs2*ai[row];
614:     rmax = imax[row];
615:     nrow = ailen[row];
616:     low  = 0;
617:     high = nrow;
618:     for (l=0; l<n; l++) { /* loop over added columns */
619:       if (in[l] < 0) continue;
620:       col = in[l];
621: #if defined(PETSC_USE_DEBUG)
622:       if (col >= a->nbs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",col,a->nbs-1);
623: #endif
624:       if (col < row) {
625:         if (a->ignore_ltriangular) continue; /* ignore lower triangular block */
626:         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)");
627:       }
628:       if (roworiented) value = v + k*(stepval+bs)*bs + l*bs;
629:       else value = v + l*(stepval+bs)*bs + k*bs;

631:       if (col <= lastcol) low = 0;
632:       else high = nrow;

634:       lastcol = col;
635:       while (high-low > 7) {
636:         t = (low+high)/2;
637:         if (rp[t] > col) high = t;
638:         else             low  = t;
639:       }
640:       for (i=low; i<high; i++) {
641:         if (rp[i] > col) break;
642:         if (rp[i] == col) {
643:           bap = ap +  bs2*i;
644:           if (roworiented) {
645:             if (is == ADD_VALUES) {
646:               for (ii=0; ii<bs; ii++,value+=stepval) {
647:                 for (jj=ii; jj<bs2; jj+=bs) {
648:                   bap[jj] += *value++;
649:                 }
650:               }
651:             } else {
652:               for (ii=0; ii<bs; ii++,value+=stepval) {
653:                 for (jj=ii; jj<bs2; jj+=bs) {
654:                   bap[jj] = *value++;
655:                 }
656:                }
657:             }
658:           } else {
659:             if (is == ADD_VALUES) {
660:               for (ii=0; ii<bs; ii++,value+=stepval) {
661:                 for (jj=0; jj<bs; jj++) {
662:                   *bap++ += *value++;
663:                 }
664:               }
665:             } else {
666:               for (ii=0; ii<bs; ii++,value+=stepval) {
667:                 for (jj=0; jj<bs; jj++) {
668:                   *bap++  = *value++;
669:                 }
670:               }
671:             }
672:           }
673:           goto noinsert2;
674:         }
675:       }
676:       if (nonew == 1) goto noinsert2;
677:       if (nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col);
678:       MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar);
679:       N = nrow++ - 1; high++;
680:       /* shift up all the later entries in this row */
681:       for (ii=N; ii>=i; ii--) {
682:         rp[ii+1] = rp[ii];
683:         PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));
684:       }
685:       if (N >= i) {
686:         PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));
687:       }
688:       rp[i] = col;
689:       bap   = ap +  bs2*i;
690:       if (roworiented) {
691:         for (ii=0; ii<bs; ii++,value+=stepval) {
692:           for (jj=ii; jj<bs2; jj+=bs) {
693:             bap[jj] = *value++;
694:           }
695:         }
696:       } else {
697:         for (ii=0; ii<bs; ii++,value+=stepval) {
698:           for (jj=0; jj<bs; jj++) {
699:             *bap++ = *value++;
700:           }
701:         }
702:        }
703:     noinsert2:;
704:       low = i;
705:     }
706:     ailen[row] = nrow;
707:   }
708:   return(0);
709: }

711: /*
712:     This is not yet used
713: */
716: PetscErrorCode MatAssemblyEnd_SeqSBAIJ_SeqAIJ_Inode(Mat A)
717: {
718:   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
720:   const PetscInt *ai = a->i, *aj = a->j,*cols;
721:   PetscInt       i   = 0,j,blk_size,m = A->rmap->n,node_count = 0,nzx,nzy,*ns,row,nz,cnt,cnt2,*counts;
722:   PetscBool      flag;

725:   PetscMalloc1(m,&ns);
726:   while (i < m) {
727:     nzx = ai[i+1] - ai[i];       /* Number of nonzeros */
728:     /* Limits the number of elements in a node to 'a->inode.limit' */
729:     for (j=i+1,blk_size=1; j<m && blk_size <a->inode.limit; ++j,++blk_size) {
730:       nzy = ai[j+1] - ai[j];
731:       if (nzy != (nzx - j + i)) break;
732:       PetscMemcmp(aj + ai[i] + j - i,aj + ai[j],nzy*sizeof(PetscInt),&flag);
733:       if (!flag) break;
734:     }
735:     ns[node_count++] = blk_size;

737:     i = j;
738:   }
739:   if (!a->inode.size && m && node_count > .9*m) {
740:     PetscFree(ns);
741:     PetscInfo2(A,"Found %D nodes out of %D rows. Not using Inode routines\n",node_count,m);
742:   } else {
743:     a->inode.node_count = node_count;

745:     PetscMalloc1(node_count,&a->inode.size);
746:     PetscLogObjectMemory((PetscObject)A,node_count*sizeof(PetscInt));
747:     PetscMemcpy(a->inode.size,ns,node_count*sizeof(PetscInt));
748:     PetscFree(ns);
749:     PetscInfo3(A,"Found %D nodes of %D. Limit used: %D. Using Inode routines\n",node_count,m,a->inode.limit);

751:     /* count collections of adjacent columns in each inode */
752:     row = 0;
753:     cnt = 0;
754:     for (i=0; i<node_count; i++) {
755:       cols = aj + ai[row] + a->inode.size[i];
756:       nz   = ai[row+1] - ai[row] - a->inode.size[i];
757:       for (j=1; j<nz; j++) {
758:         if (cols[j] != cols[j-1]+1) cnt++;
759:       }
760:       cnt++;
761:       row += a->inode.size[i];
762:     }
763:     PetscMalloc1(2*cnt,&counts);
764:     cnt  = 0;
765:     row  = 0;
766:     for (i=0; i<node_count; i++) {
767:       cols = aj + ai[row] + a->inode.size[i];
768:       counts[2*cnt] = cols[0];
769:       nz   = ai[row+1] - ai[row] - a->inode.size[i];
770:       cnt2 = 1;
771:       for (j=1; j<nz; j++) {
772:         if (cols[j] != cols[j-1]+1) {
773:           counts[2*(cnt++)+1] = cnt2;
774:           counts[2*cnt]       = cols[j];
775:           cnt2 = 1;
776:         } else cnt2++;
777:       }
778:       counts[2*(cnt++)+1] = cnt2;
779:       row += a->inode.size[i];
780:     }
781:     PetscIntView(2*cnt,counts,0);
782:   }
783:   return(0);
784: }

788: PetscErrorCode MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode)
789: {
790:   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
792:   PetscInt       fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
793:   PetscInt       m      = A->rmap->N,*ip,N,*ailen = a->ilen;
794:   PetscInt       mbs    = a->mbs,bs2 = a->bs2,rmax = 0;
795:   MatScalar      *aa    = a->a,*ap;

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

800:   if (m) rmax = ailen[0];
801:   for (i=1; i<mbs; i++) {
802:     /* move each row back by the amount of empty slots (fshift) before it*/
803:     fshift += imax[i-1] - ailen[i-1];
804:     rmax    = PetscMax(rmax,ailen[i]);
805:     if (fshift) {
806:       ip = aj + ai[i]; ap = aa + bs2*ai[i];
807:       N  = ailen[i];
808:       for (j=0; j<N; j++) {
809:         ip[j-fshift] = ip[j];
810:         PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));
811:       }
812:     }
813:     ai[i] = ai[i-1] + ailen[i-1];
814:   }
815:   if (mbs) {
816:     fshift += imax[mbs-1] - ailen[mbs-1];
817:     ai[mbs] = ai[mbs-1] + ailen[mbs-1];
818:   }
819:   /* reset ilen and imax for each row */
820:   for (i=0; i<mbs; i++) {
821:     ailen[i] = imax[i] = ai[i+1] - ai[i];
822:   }
823:   a->nz = ai[mbs];

825:   /* diagonals may have moved, reset it */
826:   if (a->diag) {
827:     PetscMemcpy(a->diag,ai,mbs*sizeof(PetscInt));
828:   }
829:   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);

831:   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);
832:   PetscInfo1(A,"Number of mallocs during MatSetValues is %D\n",a->reallocs);
833:   PetscInfo1(A,"Most nonzeros blocks in any row is %D\n",rmax);

835:   A->info.mallocs    += a->reallocs;
836:   a->reallocs         = 0;
837:   A->info.nz_unneeded = (PetscReal)fshift*bs2;
838:   a->idiagvalid       = PETSC_FALSE;
839:   a->rmax             = rmax;

841:   if (A->cmap->n < 65536 && A->cmap->bs == 1) {
842:     if (a->jshort && a->free_jshort) {
843:       /* when matrix data structure is changed, previous jshort must be replaced */
844:       PetscFree(a->jshort);
845:     }
846:     PetscMalloc1(a->i[A->rmap->n],&a->jshort);
847:     PetscLogObjectMemory((PetscObject)A,a->i[A->rmap->n]*sizeof(unsigned short));
848:     for (i=0; i<a->i[A->rmap->n]; i++) a->jshort[i] = a->j[i];
849:     A->ops->mult   = MatMult_SeqSBAIJ_1_ushort;
850:     A->ops->sor    = MatSOR_SeqSBAIJ_ushort;
851:     a->free_jshort = PETSC_TRUE;
852:   }
853:   return(0);
854: }

856: /*
857:    This function returns an array of flags which indicate the locations of contiguous
858:    blocks that should be zeroed. for eg: if bs = 3  and is = [0,1,2,3,5,6,7,8,9]
859:    then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)]
860:    Assume: sizes should be long enough to hold all the values.
861: */
864: PetscErrorCode MatZeroRows_SeqSBAIJ_Check_Blocks(PetscInt idx[],PetscInt n,PetscInt bs,PetscInt sizes[], PetscInt *bs_max)
865: {
866:   PetscInt  i,j,k,row;
867:   PetscBool flg;

870:   for (i=0,j=0; i<n; j++) {
871:     row = idx[i];
872:     if (row%bs!=0) { /* Not the begining of a block */
873:       sizes[j] = 1;
874:       i++;
875:     } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */
876:       sizes[j] = 1;         /* Also makes sure atleast 'bs' values exist for next else */
877:       i++;
878:     } else { /* Begining of the block, so check if the complete block exists */
879:       flg = PETSC_TRUE;
880:       for (k=1; k<bs; k++) {
881:         if (row+k != idx[i+k]) { /* break in the block */
882:           flg = PETSC_FALSE;
883:           break;
884:         }
885:       }
886:       if (flg) { /* No break in the bs */
887:         sizes[j] = bs;
888:         i       += bs;
889:       } else {
890:         sizes[j] = 1;
891:         i++;
892:       }
893:     }
894:   }
895:   *bs_max = j;
896:   return(0);
897: }


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

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

917:   for (k=0; k<m; k++) { /* loop over added rows */
918:     row  = im[k];       /* row number */
919:     brow = row/bs;      /* block row number */
920:     if (row < 0) continue;
921: #if defined(PETSC_USE_DEBUG)
922:     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);
923: #endif
924:     rp   = aj + ai[brow]; /*ptr to beginning of column value of the row block*/
925:     ap   = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/
926:     rmax = imax[brow];  /* maximum space allocated for this row */
927:     nrow = ailen[brow]; /* actual length of this row */
928:     low  = 0;

930:     for (l=0; l<n; l++) { /* loop over added columns */
931:       if (in[l] < 0) continue;
932: #if defined(PETSC_USE_DEBUG)
933:       if (in[l] >= A->rmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->rmap->N-1);
934: #endif
935:       col  = in[l];
936:       bcol = col/bs;              /* block col number */

938:       if (brow > bcol) {
939:         if (a->ignore_ltriangular) continue; /* ignore lower triangular values */
940:         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)");
941:       }

943:       ridx = row % bs; cidx = col % bs; /*row and col index inside the block */
944:       if ((brow==bcol && ridx<=cidx) || (brow<bcol)) {
945:         /* element value a(k,l) */
946:         if (roworiented) value = v[l + k*n];
947:         else value = v[k + l*m];

949:         /* move pointer bap to a(k,l) quickly and add/insert value */
950:         if (col <= lastcol) low = 0;
951:         high = nrow;
952:         lastcol = col;
953:         while (high-low > 7) {
954:           t = (low+high)/2;
955:           if (rp[t] > bcol) high = t;
956:           else              low  = t;
957:         }
958:         for (i=low; i<high; i++) {
959:           if (rp[i] > bcol) break;
960:           if (rp[i] == bcol) {
961:             bap = ap +  bs2*i + bs*cidx + ridx;
962:             if (is == ADD_VALUES) *bap += value;
963:             else                  *bap  = value;
964:             /* for diag block, add/insert its symmetric element a(cidx,ridx) */
965:             if (brow == bcol && ridx < cidx) {
966:               bap = ap +  bs2*i + bs*ridx + cidx;
967:               if (is == ADD_VALUES) *bap += value;
968:               else                  *bap  = value;
969:             }
970:             goto noinsert1;
971:           }
972:         }

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

978:         N = nrow++ - 1; high++;
979:         /* shift up all the later entries in this row */
980:         for (ii=N; ii>=i; ii--) {
981:           rp[ii+1] = rp[ii];
982:           PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));
983:         }
984:         if (N>=i) {
985:           PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));
986:         }
987:         rp[i]                      = bcol;
988:         ap[bs2*i + bs*cidx + ridx] = value;
989:         A->nonzerostate++;
990: noinsert1:;
991:         low = i;
992:       }
993:     }   /* end of loop over added columns */
994:     ailen[brow] = nrow;
995:   }   /* end of loop over added rows */
996:   return(0);
997: }

1001: PetscErrorCode MatICCFactor_SeqSBAIJ(Mat inA,IS row,const MatFactorInfo *info)
1002: {
1003:   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)inA->data;
1004:   Mat            outA;
1006:   PetscBool      row_identity;

1009:   if (info->levels != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only levels=0 is supported for in-place icc");
1010:   ISIdentity(row,&row_identity);
1011:   if (!row_identity) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix reordering is not supported");
1012:   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()! */

1014:   outA            = inA;
1015:   inA->factortype = MAT_FACTOR_ICC;

1017:   MatMarkDiagonal_SeqSBAIJ(inA);
1018:   MatSeqSBAIJSetNumericFactorization_inplace(inA,row_identity);

1020:   PetscObjectReference((PetscObject)row);
1021:   ISDestroy(&a->row);
1022:   a->row = row;
1023:   PetscObjectReference((PetscObject)row);
1024:   ISDestroy(&a->col);
1025:   a->col = row;

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

1031:   if (!a->solve_work) {
1032:     PetscMalloc1(inA->rmap->N+inA->rmap->bs,&a->solve_work);
1033:     PetscLogObjectMemory((PetscObject)inA,(inA->rmap->N+inA->rmap->bs)*sizeof(PetscScalar));
1034:   }

1036:   MatCholeskyFactorNumeric(outA,inA,info);
1037:   return(0);
1038: }

1042: PetscErrorCode  MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,PetscInt *indices)
1043: {
1044:   Mat_SeqSBAIJ   *baij = (Mat_SeqSBAIJ*)mat->data;
1045:   PetscInt       i,nz,n;

1049:   nz = baij->maxnz;
1050:   n  = mat->cmap->n;
1051:   for (i=0; i<nz; i++) baij->j[i] = indices[i];

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

1056:   MatSetOption(mat,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
1057:   return(0);
1058: }

1062: /*@
1063:   MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows
1064:   in the matrix.

1066:   Input Parameters:
1067:   +  mat     - the SeqSBAIJ matrix
1068:   -  indices - the column indices

1070:   Level: advanced

1072:   Notes:
1073:   This can be called if you have precomputed the nonzero structure of the
1074:   matrix and want to provide it to the matrix object to improve the performance
1075:   of the MatSetValues() operation.

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

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

1082:   .seealso: MatCreateSeqSBAIJ
1083: @*/
1084: PetscErrorCode  MatSeqSBAIJSetColumnIndices(Mat mat,PetscInt *indices)
1085: {

1091:   PetscUseMethod(mat,"MatSeqSBAIJSetColumnIndices_C",(Mat,PetscInt*),(mat,indices));
1092:   return(0);
1093: }

1097: PetscErrorCode MatCopy_SeqSBAIJ(Mat A,Mat B,MatStructure str)
1098: {

1102:   /* If the two matrices have the same copy implementation, use fast copy. */
1103:   if (str == SAME_NONZERO_PATTERN && (A->ops->copy == B->ops->copy)) {
1104:     Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1105:     Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data;

1107:     if (a->i[A->rmap->N] != b->i[B->rmap->N]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of nonzeros in two matrices are different");
1108:     PetscMemcpy(b->a,a->a,(a->i[A->rmap->N])*sizeof(PetscScalar));
1109:   } else {
1110:     MatGetRowUpperTriangular(A);
1111:     MatCopy_Basic(A,B,str);
1112:     MatRestoreRowUpperTriangular(A);
1113:   }
1114:   return(0);
1115: }

1119: PetscErrorCode MatSetUp_SeqSBAIJ(Mat A)
1120: {

1124:   MatSeqSBAIJSetPreallocation_SeqSBAIJ(A,A->rmap->bs,PETSC_DEFAULT,0);
1125:   return(0);
1126: }

1130: PetscErrorCode MatSeqSBAIJGetArray_SeqSBAIJ(Mat A,PetscScalar *array[])
1131: {
1132:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;

1135:   *array = a->a;
1136:   return(0);
1137: }

1141: PetscErrorCode MatSeqSBAIJRestoreArray_SeqSBAIJ(Mat A,PetscScalar *array[])
1142: {
1144:   return(0);
1145: }

1149: PetscErrorCode MatAXPYGetPreallocation_SeqSBAIJ(Mat Y,Mat X,PetscInt *nnz)
1150: {
1151:   PetscInt       bs = Y->rmap->bs,mbs = Y->rmap->N/bs;
1152:   Mat_SeqSBAIJ   *x = (Mat_SeqSBAIJ*)X->data;
1153:   Mat_SeqSBAIJ   *y = (Mat_SeqSBAIJ*)Y->data;

1157:   /* Set the number of nonzeros in the new matrix */
1158:   MatAXPYGetPreallocation_SeqX_private(mbs,x->i,x->j,y->i,y->j,nnz);
1159:   return(0);
1160: }

1164: PetscErrorCode MatAXPY_SeqSBAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str)
1165: {
1166:   Mat_SeqSBAIJ   *x=(Mat_SeqSBAIJ*)X->data, *y=(Mat_SeqSBAIJ*)Y->data;
1168:   PetscInt       bs=Y->rmap->bs,bs2=bs*bs;
1169:   PetscBLASInt   one = 1;

1172:   if (str == SAME_NONZERO_PATTERN) {
1173:     PetscScalar  alpha = a;
1174:     PetscBLASInt bnz;
1175:     PetscBLASIntCast(x->nz*bs2,&bnz);
1176:     PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&bnz,&alpha,x->a,&one,y->a,&one));
1177:     PetscObjectStateIncrease((PetscObject)Y);
1178:   } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
1179:     MatSetOption(X,MAT_GETROW_UPPERTRIANGULAR,PETSC_TRUE);
1180:     MatAXPY_Basic(Y,a,X,str);
1181:     MatSetOption(X,MAT_GETROW_UPPERTRIANGULAR,PETSC_FALSE);
1182:   } else {
1183:     Mat      B;
1184:     PetscInt *nnz;
1185:     if (bs != X->rmap->bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrices must have same block size");
1186:     MatGetRowUpperTriangular(X);
1187:     MatGetRowUpperTriangular(Y);
1188:     PetscMalloc1(Y->rmap->N,&nnz);
1189:     MatCreate(PetscObjectComm((PetscObject)Y),&B);
1190:     PetscObjectSetName((PetscObject)B,((PetscObject)Y)->name);
1191:     MatSetSizes(B,Y->rmap->n,Y->cmap->n,Y->rmap->N,Y->cmap->N);
1192:     MatSetBlockSizesFromMats(B,Y,Y);
1193:     MatSetType(B,(MatType) ((PetscObject)Y)->type_name);
1194:     MatAXPYGetPreallocation_SeqSBAIJ(Y,X,nnz);
1195:     MatSeqSBAIJSetPreallocation(B,bs,0,nnz);
1196: 
1197:     MatAXPY_BasicWithPreallocation(B,Y,a,X,str);
1198: 
1199:     MatHeaderReplace(Y,B);
1200:     PetscFree(nnz);
1201:     MatRestoreRowUpperTriangular(X);
1202:     MatRestoreRowUpperTriangular(Y);
1203:   }
1204:   return(0);
1205: }

1209: PetscErrorCode MatIsSymmetric_SeqSBAIJ(Mat A,PetscReal tol,PetscBool  *flg)
1210: {
1212:   *flg = PETSC_TRUE;
1213:   return(0);
1214: }

1218: PetscErrorCode MatIsStructurallySymmetric_SeqSBAIJ(Mat A,PetscBool  *flg)
1219: {
1221:   *flg = PETSC_TRUE;
1222:   return(0);
1223: }

1227: PetscErrorCode MatIsHermitian_SeqSBAIJ(Mat A,PetscReal tol,PetscBool  *flg)
1228: {
1230:   *flg = PETSC_FALSE;
1231:   return(0);
1232: }

1236: PetscErrorCode MatRealPart_SeqSBAIJ(Mat A)
1237: {
1238:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1239:   PetscInt     i,nz = a->bs2*a->i[a->mbs];
1240:   MatScalar    *aa = a->a;

1243:   for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]);
1244:   return(0);
1245: }

1249: PetscErrorCode MatImaginaryPart_SeqSBAIJ(Mat A)
1250: {
1251:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1252:   PetscInt     i,nz = a->bs2*a->i[a->mbs];
1253:   MatScalar    *aa = a->a;

1256:   for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]);
1257:   return(0);
1258: }

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

1275:   /* fix right hand side if needed */
1276:   if (x && b) {
1277:     VecGetArrayRead(x,&xx);
1278:     VecGetArray(b,&bb);
1279:     vecs = PETSC_TRUE;
1280:   }

1282:   /* zero the columns */
1283:   PetscCalloc1(A->rmap->n,&zeroed);
1284:   for (i=0; i<is_n; i++) {
1285:     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]);
1286:     zeroed[is_idx[i]] = PETSC_TRUE;
1287:   }
1288:   if (vecs) {
1289:     for (i=0; i<A->rmap->N; i++) {
1290:       row = i/bs;
1291:       for (j=baij->i[row]; j<baij->i[row+1]; j++) {
1292:         for (k=0; k<bs; k++) {
1293:           col = bs*baij->j[j] + k;
1294:           if (col <= i) continue;
1295:           aa = ((MatScalar*)(baij->a)) + j*bs2 + (i%bs) + bs*k;
1296:           if (!zeroed[i] && zeroed[col]) bb[i]   -= aa[0]*xx[col];
1297:           if (zeroed[i] && !zeroed[col]) bb[col] -= aa[0]*xx[i];
1298:         }
1299:       }
1300:     }
1301:     for (i=0; i<is_n; i++) bb[is_idx[i]] = diag*xx[is_idx[i]];
1302:   }

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

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

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

1491: PetscErrorCode  MatStoreValues_SeqSBAIJ(Mat mat)
1492: {
1493:   Mat_SeqSBAIJ   *aij = (Mat_SeqSBAIJ*)mat->data;
1494:   PetscInt       nz   = aij->i[mat->rmap->N]*mat->rmap->bs*aij->bs2;

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

1500:   /* allocate space for values if not already there */
1501:   if (!aij->saved_values) {
1502:     PetscMalloc1(nz+1,&aij->saved_values);
1503:   }

1505:   /* copy values over */
1506:   PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));
1507:   return(0);
1508: }

1512: PetscErrorCode  MatRetrieveValues_SeqSBAIJ(Mat mat)
1513: {
1514:   Mat_SeqSBAIJ   *aij = (Mat_SeqSBAIJ*)mat->data;
1516:   PetscInt       nz = aij->i[mat->rmap->N]*mat->rmap->bs*aij->bs2;

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

1522:   /* copy values over */
1523:   PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));
1524:   return(0);
1525: }

1529: PetscErrorCode  MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat B,PetscInt bs,PetscInt nz,PetscInt *nnz)
1530: {
1531:   Mat_SeqSBAIJ   *b = (Mat_SeqSBAIJ*)B->data;
1533:   PetscInt       i,mbs,nbs,bs2;
1534:   PetscBool      skipallocation = PETSC_FALSE,flg = PETSC_FALSE,realalloc = PETSC_FALSE;

1537:   if (nz >= 0 || nnz) realalloc = PETSC_TRUE;
1538:   B->preallocated = PETSC_TRUE;

1540:   MatSetBlockSize(B,PetscAbs(bs));
1541:   PetscLayoutSetUp(B->rmap);
1542:   PetscLayoutSetUp(B->cmap);
1543:   PetscLayoutGetBlockSize(B->rmap,&bs);

1545:   mbs = B->rmap->N/bs;
1546:   nbs = B->cmap->n/bs;
1547:   bs2 = bs*bs;

1549:   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");

1551:   if (nz == MAT_SKIP_ALLOCATION) {
1552:     skipallocation = PETSC_TRUE;
1553:     nz             = 0;
1554:   }

1556:   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3;
1557:   if (nz < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %D",nz);
1558:   if (nnz) {
1559:     for (i=0; i<mbs; i++) {
1560:       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]);
1561:       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);
1562:     }
1563:   }

1565:   B->ops->mult             = MatMult_SeqSBAIJ_N;
1566:   B->ops->multadd          = MatMultAdd_SeqSBAIJ_N;
1567:   B->ops->multtranspose    = MatMult_SeqSBAIJ_N;
1568:   B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_N;

1570:   PetscOptionsGetBool(((PetscObject)B)->prefix,"-mat_no_unroll",&flg,NULL);
1571:   if (!flg) {
1572:     switch (bs) {
1573:     case 1:
1574:       B->ops->mult             = MatMult_SeqSBAIJ_1;
1575:       B->ops->multadd          = MatMultAdd_SeqSBAIJ_1;
1576:       B->ops->multtranspose    = MatMult_SeqSBAIJ_1;
1577:       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_1;
1578:       break;
1579:     case 2:
1580:       B->ops->mult             = MatMult_SeqSBAIJ_2;
1581:       B->ops->multadd          = MatMultAdd_SeqSBAIJ_2;
1582:       B->ops->multtranspose    = MatMult_SeqSBAIJ_2;
1583:       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_2;
1584:       break;
1585:     case 3:
1586:       B->ops->mult             = MatMult_SeqSBAIJ_3;
1587:       B->ops->multadd          = MatMultAdd_SeqSBAIJ_3;
1588:       B->ops->multtranspose    = MatMult_SeqSBAIJ_3;
1589:       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_3;
1590:       break;
1591:     case 4:
1592:       B->ops->mult             = MatMult_SeqSBAIJ_4;
1593:       B->ops->multadd          = MatMultAdd_SeqSBAIJ_4;
1594:       B->ops->multtranspose    = MatMult_SeqSBAIJ_4;
1595:       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_4;
1596:       break;
1597:     case 5:
1598:       B->ops->mult             = MatMult_SeqSBAIJ_5;
1599:       B->ops->multadd          = MatMultAdd_SeqSBAIJ_5;
1600:       B->ops->multtranspose    = MatMult_SeqSBAIJ_5;
1601:       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_5;
1602:       break;
1603:     case 6:
1604:       B->ops->mult             = MatMult_SeqSBAIJ_6;
1605:       B->ops->multadd          = MatMultAdd_SeqSBAIJ_6;
1606:       B->ops->multtranspose    = MatMult_SeqSBAIJ_6;
1607:       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_6;
1608:       break;
1609:     case 7:
1610:       B->ops->mult             = MatMult_SeqSBAIJ_7;
1611:       B->ops->multadd          = MatMultAdd_SeqSBAIJ_7;
1612:       B->ops->multtranspose    = MatMult_SeqSBAIJ_7;
1613:       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_7;
1614:       break;
1615:     }
1616:   }

1618:   b->mbs = mbs;
1619:   b->nbs = nbs;
1620:   if (!skipallocation) {
1621:     if (!b->imax) {
1622:       PetscMalloc2(mbs,&b->imax,mbs,&b->ilen);

1624:       b->free_imax_ilen = PETSC_TRUE;

1626:       PetscLogObjectMemory((PetscObject)B,2*mbs*sizeof(PetscInt));
1627:     }
1628:     if (!nnz) {
1629:       if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
1630:       else if (nz <= 0) nz = 1;
1631:       for (i=0; i<mbs; i++) b->imax[i] = nz;
1632:       nz = nz*mbs; /* total nz */
1633:     } else {
1634:       nz = 0;
1635:       for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
1636:     }
1637:     /* b->ilen will count nonzeros in each block row so far. */
1638:     for (i=0; i<mbs; i++) b->ilen[i] = 0;
1639:     /* nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */

1641:     /* allocate the matrix space */
1642:     MatSeqXAIJFreeAIJ(B,&b->a,&b->j,&b->i);
1643:     PetscMalloc3(bs2*nz,&b->a,nz,&b->j,B->rmap->N+1,&b->i);
1644:     PetscLogObjectMemory((PetscObject)B,(B->rmap->N+1)*sizeof(PetscInt)+nz*(bs2*sizeof(PetscScalar)+sizeof(PetscInt)));
1645:     PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));
1646:     PetscMemzero(b->j,nz*sizeof(PetscInt));

1648:     b->singlemalloc = PETSC_TRUE;

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

1654:     b->free_a  = PETSC_TRUE;
1655:     b->free_ij = PETSC_TRUE;
1656:   } else {
1657:     b->free_a  = PETSC_FALSE;
1658:     b->free_ij = PETSC_FALSE;
1659:   }

1661:   B->rmap->bs = bs;
1662:   b->bs2      = bs2;
1663:   b->nz       = 0;
1664:   b->maxnz    = nz;

1666:   b->inew    = 0;
1667:   b->jnew    = 0;
1668:   b->anew    = 0;
1669:   b->a2anew  = 0;
1670:   b->permute = PETSC_FALSE;
1671:   if (realalloc) {MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);}
1672:   return(0);
1673: }

1677: PetscErrorCode MatSeqSBAIJSetPreallocationCSR_SeqSBAIJ(Mat B,PetscInt bs,const PetscInt ii[],const PetscInt jj[], const PetscScalar V[])
1678: {
1679:   PetscInt       i,j,m,nz,nz_max=0,*nnz;
1680:   PetscScalar    *values=0;
1681:   PetscBool      roworiented = ((Mat_SeqSBAIJ*)B->data)->roworiented;
1684:   if (bs < 1) SETERRQ1(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_OUTOFRANGE,"Invalid block size specified, must be positive but it is %D",bs);
1685:   PetscLayoutSetBlockSize(B->rmap,bs);
1686:   PetscLayoutSetBlockSize(B->cmap,bs);
1687:   PetscLayoutSetUp(B->rmap);
1688:   PetscLayoutSetUp(B->cmap);
1689:   PetscLayoutGetBlockSize(B->rmap,&bs);
1690:   m      = B->rmap->n/bs;

1692:   if (ii[0]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"ii[0] must be 0 but it is %D",ii[0]);
1693:   PetscMalloc1(m+1,&nnz);
1694:   for (i=0; i<m; i++) {
1695:     nz = ii[i+1] - ii[i];
1696:     if (nz < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row %D has a negative number of columns %D",i,nz);
1697:     nz_max = PetscMax(nz_max,nz);
1698:     nnz[i] = nz;
1699:   }
1700:   MatSeqSBAIJSetPreallocation(B,bs,0,nnz);
1701:   PetscFree(nnz);

1703:   values = (PetscScalar*)V;
1704:   if (!values) {
1705:     PetscCalloc1(bs*bs*nz_max,&values);
1706:   }
1707:   for (i=0; i<m; i++) {
1708:     PetscInt          ncols  = ii[i+1] - ii[i];
1709:     const PetscInt    *icols = jj + ii[i];
1710:     if (!roworiented || bs == 1) {
1711:       const PetscScalar *svals = values + (V ? (bs*bs*ii[i]) : 0);
1712:       MatSetValuesBlocked_SeqSBAIJ(B,1,&i,ncols,icols,svals,INSERT_VALUES);
1713:     } else {
1714:       for (j=0; j<ncols; j++) {
1715:         const PetscScalar *svals = values + (V ? (bs*bs*(ii[i]+j)) : 0);
1716:         MatSetValuesBlocked_SeqSBAIJ(B,1,&i,1,&icols[j],svals,INSERT_VALUES);
1717:       }
1718:     }
1719:   }
1720:   if (!V) { PetscFree(values); }
1721:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
1722:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
1723:   MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
1724:   return(0);
1725: }

1727: /*
1728:    This is used to set the numeric factorization for both Cholesky and ICC symbolic factorization
1729: */
1732: PetscErrorCode MatSeqSBAIJSetNumericFactorization_inplace(Mat B,PetscBool natural)
1733: {
1735:   PetscBool      flg = PETSC_FALSE;
1736:   PetscInt       bs  = B->rmap->bs;

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

1742:   if (!natural) {
1743:     switch (bs) {
1744:     case 1:
1745:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_inplace;
1746:       break;
1747:     case 2:
1748:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2;
1749:       break;
1750:     case 3:
1751:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3;
1752:       break;
1753:     case 4:
1754:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4;
1755:       break;
1756:     case 5:
1757:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5;
1758:       break;
1759:     case 6:
1760:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6;
1761:       break;
1762:     case 7:
1763:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7;
1764:       break;
1765:     default:
1766:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N;
1767:       break;
1768:     }
1769:   } else {
1770:     switch (bs) {
1771:     case 1:
1772:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering_inplace;
1773:       break;
1774:     case 2:
1775:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering;
1776:       break;
1777:     case 3:
1778:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering;
1779:       break;
1780:     case 4:
1781:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering;
1782:       break;
1783:     case 5:
1784:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering;
1785:       break;
1786:     case 6:
1787:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering;
1788:       break;
1789:     case 7:
1790:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering;
1791:       break;
1792:     default:
1793:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering;
1794:       break;
1795:     }
1796:   }
1797:   return(0);
1798: }

1800: PETSC_EXTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqAIJ(Mat, MatType,MatReuse,Mat*);
1801: PETSC_EXTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqBAIJ(Mat, MatType,MatReuse,Mat*);

1805: PETSC_EXTERN PetscErrorCode MatGetFactor_seqsbaij_petsc(Mat A,MatFactorType ftype,Mat *B)
1806: {
1807:   PetscInt       n = A->rmap->n;

1811: #if defined(PETSC_USE_COMPLEX)
1812:   if (A->hermitian) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Hermitian Factor is not supported");
1813: #endif
1814:   MatCreate(PetscObjectComm((PetscObject)A),B);
1815:   MatSetSizes(*B,n,n,n,n);
1816:   if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) {
1817:     MatSetType(*B,MATSEQSBAIJ);
1818:     MatSeqSBAIJSetPreallocation(*B,A->rmap->bs,MAT_SKIP_ALLOCATION,NULL);

1820:     (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqSBAIJ;
1821:     (*B)->ops->iccfactorsymbolic      = MatICCFactorSymbolic_SeqSBAIJ;
1822:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
1823:   (*B)->factortype = ftype;
1824:   return(0);
1825: }

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

1831:   For complex numbers by default this matrix is symmetric, NOT Hermitian symmetric. To make it Hermitian symmetric you
1832:   can call MatSetOption(Mat, MAT_HERMITIAN); after MatAssemblyEnd()

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

1837:   Notes: By default if you insert values into the lower triangular part of the matrix they are simply ignored (since they are not
1838:      stored and it is assumed they symmetric to the upper triangular). If you call MatSetOption(Mat,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_FALSE) or use
1839:      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.


1842:   Level: beginner

1844:   .seealso: MatCreateSeqSBAIJ
1845: M*/

1847: PETSC_EXTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqSBSTRM(Mat, MatType,MatReuse,Mat*);

1851: PETSC_EXTERN PetscErrorCode MatCreate_SeqSBAIJ(Mat B)
1852: {
1853:   Mat_SeqSBAIJ   *b;
1855:   PetscMPIInt    size;
1856:   PetscBool      no_unroll = PETSC_FALSE,no_inode = PETSC_FALSE;

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

1862:   PetscNewLog(B,&b);
1863:   B->data = (void*)b;
1864:   PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));

1866:   B->ops->destroy    = MatDestroy_SeqSBAIJ;
1867:   B->ops->view       = MatView_SeqSBAIJ;
1868:   b->row             = 0;
1869:   b->icol            = 0;
1870:   b->reallocs        = 0;
1871:   b->saved_values    = 0;
1872:   b->inode.limit     = 5;
1873:   b->inode.max_limit = 5;

1875:   b->roworiented        = PETSC_TRUE;
1876:   b->nonew              = 0;
1877:   b->diag               = 0;
1878:   b->solve_work         = 0;
1879:   b->mult_work          = 0;
1880:   B->spptr              = 0;
1881:   B->info.nz_unneeded   = (PetscReal)b->maxnz*b->bs2;
1882:   b->keepnonzeropattern = PETSC_FALSE;

1884:   b->inew    = 0;
1885:   b->jnew    = 0;
1886:   b->anew    = 0;
1887:   b->a2anew  = 0;
1888:   b->permute = PETSC_FALSE;

1890:   b->ignore_ltriangular = PETSC_TRUE;

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

1894:   b->getrow_utriangular = PETSC_FALSE;

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

1898:   PetscObjectComposeFunction((PetscObject)B,"MatStoreValues_C",MatStoreValues_SeqSBAIJ);
1899:   PetscObjectComposeFunction((PetscObject)B,"MatRetrieveValues_C",MatRetrieveValues_SeqSBAIJ);
1900:   PetscObjectComposeFunction((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C",MatSeqSBAIJSetColumnIndices_SeqSBAIJ);
1901:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_seqaij_C",MatConvert_SeqSBAIJ_SeqAIJ);
1902:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_seqbaij_C",MatConvert_SeqSBAIJ_SeqBAIJ);
1903:   PetscObjectComposeFunction((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",MatSeqSBAIJSetPreallocation_SeqSBAIJ);
1904:   PetscObjectComposeFunction((PetscObject)B,"MatSeqSBAIJSetPreallocationCSR_C",MatSeqSBAIJSetPreallocationCSR_SeqSBAIJ);
1905:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_seqsbstrm_C",MatConvert_SeqSBAIJ_SeqSBSTRM);

1907:   B->symmetric                  = PETSC_TRUE;
1908:   B->structurally_symmetric     = PETSC_TRUE;
1909:   B->symmetric_set              = PETSC_TRUE;
1910:   B->structurally_symmetric_set = PETSC_TRUE;

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

1914:   PetscOptionsBegin(PetscObjectComm((PetscObject)B),((PetscObject)B)->prefix,"Options for SEQSBAIJ matrix","Mat");
1915:   PetscOptionsBool("-mat_no_unroll","Do not optimize for inodes (slower)",NULL,no_unroll,&no_unroll,NULL);
1916:   if (no_unroll) {
1917:     PetscInfo(B,"Not using Inode routines due to -mat_no_unroll\n");
1918:   }
1919:   PetscOptionsBool("-mat_no_inode","Do not optimize for inodes (slower)",NULL,no_inode,&no_inode,NULL);
1920:   if (no_inode) {
1921:     PetscInfo(B,"Not using Inode routines due to -mat_no_inode\n");
1922:   }
1923:   PetscOptionsInt("-mat_inode_limit","Do not use inodes larger then this value",NULL,b->inode.limit,&b->inode.limit,NULL);
1924:   PetscOptionsEnd();
1925:   b->inode.use = (PetscBool)(!(no_unroll || no_inode));
1926:   if (b->inode.limit > b->inode.max_limit) b->inode.limit = b->inode.max_limit;
1927:   return(0);
1928: }

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

1939:    Collective on Mat

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

1949:    Options Database Keys:
1950: .   -mat_no_unroll - uses code that does not unroll the loops in the
1951:                      block calculations (much slower)
1952: .    -mat_block_size - size of the blocks to use (only works if a negative bs is passed in

1954:    Level: intermediate

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

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

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


1969: .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateSBAIJ()
1970: @*/
1971: PetscErrorCode  MatSeqSBAIJSetPreallocation(Mat B,PetscInt bs,PetscInt nz,const PetscInt nnz[])
1972: {

1979:   PetscTryMethod(B,"MatSeqSBAIJSetPreallocation_C",(Mat,PetscInt,PetscInt,const PetscInt[]),(B,bs,nz,nnz));
1980:   return(0);
1981: }

1983: #undef  __FUNCT__
1985: /*@C
1986:    MatSeqSBAIJSetPreallocationCSR - Allocates memory for a sparse sequential matrix in symmetric block AIJ format.

1988:    Input Parameters:
1989: +  B - the matrix
1990: .  i - the indices into j for the start of each local row (starts with zero)
1991: .  j - the column indices for each local row (starts with zero) these must be sorted for each row
1992: -  v - optional values in the matrix

1994:    Level: developer

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

2003: .keywords: matrix, block, aij, compressed row, sparse

2005: .seealso: MatCreate(), MatCreateSeqSBAIJ(), MatSetValuesBlocked(), MatSeqSBAIJSetPreallocation(), MATSEQSBAIJ
2006: @*/
2007: PetscErrorCode MatSeqSBAIJSetPreallocationCSR(Mat B,PetscInt bs,const PetscInt i[],const PetscInt j[], const PetscScalar v[])
2008: {

2015:   PetscTryMethod(B,"MatSeqSBAIJSetPreallocationCSR_C",(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]),(B,bs,i,j,v));
2016:   return(0);
2017: }

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

2028:    Collective on MPI_Comm

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

2039:    Output Parameter:
2040: .  A - the symmetric matrix

2042:    Options Database Keys:
2043: .   -mat_no_unroll - uses code that does not unroll the loops in the
2044:                      block calculations (much slower)
2045: .    -mat_block_size - size of the blocks to use

2047:    Level: intermediate

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

2053:    Notes:
2054:    The number of rows and columns must be divisible by blocksize.
2055:    This matrix type does not support complex Hermitian operation.

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

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

2063: .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateSBAIJ()
2064: @*/
2065: PetscErrorCode  MatCreateSeqSBAIJ(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A)
2066: {

2070:   MatCreate(comm,A);
2071:   MatSetSizes(*A,m,n,m,n);
2072:   MatSetType(*A,MATSEQSBAIJ);
2073:   MatSeqSBAIJSetPreallocation_SeqSBAIJ(*A,bs,nz,(PetscInt*)nnz);
2074:   return(0);
2075: }

2079: PetscErrorCode MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
2080: {
2081:   Mat            C;
2082:   Mat_SeqSBAIJ   *c,*a = (Mat_SeqSBAIJ*)A->data;
2084:   PetscInt       i,mbs = a->mbs,nz = a->nz,bs2 =a->bs2;

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

2089:   *B   = 0;
2090:   MatCreate(PetscObjectComm((PetscObject)A),&C);
2091:   MatSetSizes(C,A->rmap->N,A->cmap->n,A->rmap->N,A->cmap->n);
2092:   MatSetType(C,MATSEQSBAIJ);
2093:   PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));
2094:   c    = (Mat_SeqSBAIJ*)C->data;

2096:   C->preallocated       = PETSC_TRUE;
2097:   C->factortype         = A->factortype;
2098:   c->row                = 0;
2099:   c->icol               = 0;
2100:   c->saved_values       = 0;
2101:   c->keepnonzeropattern = a->keepnonzeropattern;
2102:   C->assembled          = PETSC_TRUE;

2104:   PetscLayoutReference(A->rmap,&C->rmap);
2105:   PetscLayoutReference(A->cmap,&C->cmap);
2106:   c->bs2 = a->bs2;
2107:   c->mbs = a->mbs;
2108:   c->nbs = a->nbs;

2110:   if  (cpvalues == MAT_SHARE_NONZERO_PATTERN) {
2111:     c->imax           = a->imax;
2112:     c->ilen           = a->ilen;
2113:     c->free_imax_ilen = PETSC_FALSE;
2114:   } else {
2115:     PetscMalloc2((mbs+1),&c->imax,(mbs+1),&c->ilen);
2116:     PetscLogObjectMemory((PetscObject)C,2*(mbs+1)*sizeof(PetscInt));
2117:     for (i=0; i<mbs; i++) {
2118:       c->imax[i] = a->imax[i];
2119:       c->ilen[i] = a->ilen[i];
2120:     }
2121:     c->free_imax_ilen = PETSC_TRUE;
2122:   }

2124:   /* allocate the matrix space */
2125:   if (cpvalues == MAT_SHARE_NONZERO_PATTERN) {
2126:     PetscMalloc1(bs2*nz,&c->a);
2127:     PetscLogObjectMemory((PetscObject)C,nz*bs2*sizeof(MatScalar));
2128:     c->i            = a->i;
2129:     c->j            = a->j;
2130:     c->singlemalloc = PETSC_FALSE;
2131:     c->free_a       = PETSC_TRUE;
2132:     c->free_ij      = PETSC_FALSE;
2133:     c->parent       = A;
2134:     PetscObjectReference((PetscObject)A);
2135:     MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
2136:     MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
2137:   } else {
2138:     PetscMalloc3(bs2*nz,&c->a,nz,&c->j,mbs+1,&c->i);
2139:     PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(PetscInt));
2140:     PetscLogObjectMemory((PetscObject)C,(mbs+1)*sizeof(PetscInt) + nz*(bs2*sizeof(MatScalar) + sizeof(PetscInt)));
2141:     c->singlemalloc = PETSC_TRUE;
2142:     c->free_a       = PETSC_TRUE;
2143:     c->free_ij      = PETSC_TRUE;
2144:   }
2145:   if (mbs > 0) {
2146:     if (cpvalues != MAT_SHARE_NONZERO_PATTERN) {
2147:       PetscMemcpy(c->j,a->j,nz*sizeof(PetscInt));
2148:     }
2149:     if (cpvalues == MAT_COPY_VALUES) {
2150:       PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));
2151:     } else {
2152:       PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));
2153:     }
2154:     if (a->jshort) {
2155:       /* cannot share jshort, it is reallocated in MatAssemblyEnd_SeqSBAIJ() */
2156:       /* if the parent matrix is reassembled, this child matrix will never notice */
2157:       PetscMalloc1(nz,&c->jshort);
2158:       PetscLogObjectMemory((PetscObject)C,nz*sizeof(unsigned short));
2159:       PetscMemcpy(c->jshort,a->jshort,nz*sizeof(unsigned short));

2161:       c->free_jshort = PETSC_TRUE;
2162:     }
2163:   }

2165:   c->roworiented = a->roworiented;
2166:   c->nonew       = a->nonew;

2168:   if (a->diag) {
2169:     if (cpvalues == MAT_SHARE_NONZERO_PATTERN) {
2170:       c->diag      = a->diag;
2171:       c->free_diag = PETSC_FALSE;
2172:     } else {
2173:       PetscMalloc1(mbs,&c->diag);
2174:       PetscLogObjectMemory((PetscObject)C,mbs*sizeof(PetscInt));
2175:       for (i=0; i<mbs; i++) c->diag[i] = a->diag[i];
2176:       c->free_diag = PETSC_TRUE;
2177:     }
2178:   }
2179:   c->nz         = a->nz;
2180:   c->maxnz      = a->nz; /* Since we allocate exactly the right amount */
2181:   c->solve_work = 0;
2182:   c->mult_work  = 0;

2184:   *B   = C;
2185:   PetscFunctionListDuplicate(((PetscObject)A)->qlist,&((PetscObject)C)->qlist);
2186:   return(0);
2187: }

2191: PetscErrorCode MatLoad_SeqSBAIJ(Mat newmat,PetscViewer viewer)
2192: {
2193:   Mat_SeqSBAIJ   *a;
2195:   int            fd;
2196:   PetscMPIInt    size;
2197:   PetscInt       i,nz,header[4],*rowlengths=0,M,N,bs = newmat->rmap->bs;
2198:   PetscInt       *mask,mbs,*jj,j,rowcount,nzcount,k,*s_browlengths,maskcount;
2199:   PetscInt       kmax,jcount,block,idx,point,nzcountb,extra_rows,rows,cols;
2200:   PetscInt       *masked,nmask,tmp,bs2,ishift;
2201:   PetscScalar    *aa;
2202:   MPI_Comm       comm;

2205:   PetscObjectGetComm((PetscObject)viewer,&comm);
2206:   PetscOptionsGetInt(((PetscObject)newmat)->prefix,"-matload_block_size",&bs,NULL);
2207:   if (bs < 0) bs = 1;
2208:   bs2  = bs*bs;

2210:   MPI_Comm_size(comm,&size);
2211:   if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"view must have one processor");
2212:   PetscViewerBinaryGetDescriptor(viewer,&fd);
2213:   PetscBinaryRead(fd,header,4,PETSC_INT);
2214:   if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"not Mat object");
2215:   M = header[1]; N = header[2]; nz = header[3];

2217:   if (header[3] < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ");

2219:   if (M != N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Can only do square matrices");

2221:   /*
2222:      This code adds extra rows to make sure the number of rows is
2223:     divisible by the blocksize
2224:   */
2225:   mbs        = M/bs;
2226:   extra_rows = bs - M + bs*(mbs);
2227:   if (extra_rows == bs) extra_rows = 0;
2228:   else                  mbs++;
2229:   if (extra_rows) {
2230:     PetscInfo(viewer,"Padding loaded matrix to match blocksize\n");
2231:   }

2233:   /* Set global sizes if not already set */
2234:   if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) {
2235:     MatSetSizes(newmat,PETSC_DECIDE,PETSC_DECIDE,M+extra_rows,N+extra_rows);
2236:   } else { /* Check if the matrix global sizes are correct */
2237:     MatGetSize(newmat,&rows,&cols);
2238:     if (M != rows ||  N != cols) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Matrix in file of different length (%d, %d) than the input matrix (%d, %d)",M,N,rows,cols);
2239:   }

2241:   /* read in row lengths */
2242:   PetscMalloc1(M+extra_rows,&rowlengths);
2243:   PetscBinaryRead(fd,rowlengths,M,PETSC_INT);
2244:   for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;

2246:   /* read in column indices */
2247:   PetscMalloc1(nz+extra_rows,&jj);
2248:   PetscBinaryRead(fd,jj,nz,PETSC_INT);
2249:   for (i=0; i<extra_rows; i++) jj[nz+i] = M+i;

2251:   /* loop over row lengths determining block row lengths */
2252:   PetscCalloc1(mbs,&s_browlengths);
2253:   PetscMalloc2(mbs,&mask,mbs,&masked);
2254:   PetscMemzero(mask,mbs*sizeof(PetscInt));
2255:   rowcount = 0;
2256:   nzcount  = 0;
2257:   for (i=0; i<mbs; i++) {
2258:     nmask = 0;
2259:     for (j=0; j<bs; j++) {
2260:       kmax = rowlengths[rowcount];
2261:       for (k=0; k<kmax; k++) {
2262:         tmp = jj[nzcount++]/bs;   /* block col. index */
2263:         if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;}
2264:       }
2265:       rowcount++;
2266:     }
2267:     s_browlengths[i] += nmask;

2269:     /* zero out the mask elements we set */
2270:     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
2271:   }

2273:   /* Do preallocation */
2274:   MatSeqSBAIJSetPreallocation_SeqSBAIJ(newmat,bs,0,s_browlengths);
2275:   a    = (Mat_SeqSBAIJ*)newmat->data;

2277:   /* set matrix "i" values */
2278:   a->i[0] = 0;
2279:   for (i=1; i<= mbs; i++) {
2280:     a->i[i]      = a->i[i-1] + s_browlengths[i-1];
2281:     a->ilen[i-1] = s_browlengths[i-1];
2282:   }
2283:   a->nz = a->i[mbs];

2285:   /* read in nonzero values */
2286:   PetscMalloc1(nz+extra_rows,&aa);
2287:   PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);
2288:   for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0;

2290:   /* set "a" and "j" values into matrix */
2291:   nzcount = 0; jcount = 0;
2292:   for (i=0; i<mbs; i++) {
2293:     nzcountb = nzcount;
2294:     nmask    = 0;
2295:     for (j=0; j<bs; j++) {
2296:       kmax = rowlengths[i*bs+j];
2297:       for (k=0; k<kmax; k++) {
2298:         tmp = jj[nzcount++]/bs; /* block col. index */
2299:         if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;}
2300:       }
2301:     }
2302:     /* sort the masked values */
2303:     PetscSortInt(nmask,masked);

2305:     /* set "j" values into matrix */
2306:     maskcount = 1;
2307:     for (j=0; j<nmask; j++) {
2308:       a->j[jcount++]  = masked[j];
2309:       mask[masked[j]] = maskcount++;
2310:     }

2312:     /* set "a" values into matrix */
2313:     ishift = bs2*a->i[i];
2314:     for (j=0; j<bs; j++) {
2315:       kmax = rowlengths[i*bs+j];
2316:       for (k=0; k<kmax; k++) {
2317:         tmp = jj[nzcountb]/bs;        /* block col. index */
2318:         if (tmp >= i) {
2319:           block     = mask[tmp] - 1;
2320:           point     = jj[nzcountb] - bs*tmp;
2321:           idx       = ishift + bs2*block + j + bs*point;
2322:           a->a[idx] = aa[nzcountb];
2323:         }
2324:         nzcountb++;
2325:       }
2326:     }
2327:     /* zero out the mask elements we set */
2328:     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
2329:   }
2330:   if (jcount != a->nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix");

2332:   PetscFree(rowlengths);
2333:   PetscFree(s_browlengths);
2334:   PetscFree(aa);
2335:   PetscFree(jj);
2336:   PetscFree2(mask,masked);

2338:   MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);
2339:   MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);
2340:   return(0);
2341: }

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

2349:      Collective on MPI_Comm

2351:    Input Parameters:
2352: +  comm - must be an MPI communicator of size 1
2353: .  bs - size of block
2354: .  m - number of rows
2355: .  n - number of columns
2356: .  i - row indices
2357: .  j - column indices
2358: -  a - matrix values

2360:    Output Parameter:
2361: .  mat - the matrix

2363:    Level: advanced

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

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

2371:        The i and j indices are 0 based

2373:        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
2374:        it is the regular CSR format excluding the lower triangular elements.

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

2378: @*/
2379: PetscErrorCode  MatCreateSeqSBAIJWithArrays(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt *i,PetscInt *j,PetscScalar *a,Mat *mat)
2380: {
2382:   PetscInt       ii;
2383:   Mat_SeqSBAIJ   *sbaij;

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

2389:   MatCreate(comm,mat);
2390:   MatSetSizes(*mat,m,n,m,n);
2391:   MatSetType(*mat,MATSEQSBAIJ);
2392:   MatSeqSBAIJSetPreallocation_SeqSBAIJ(*mat,bs,MAT_SKIP_ALLOCATION,0);
2393:   sbaij = (Mat_SeqSBAIJ*)(*mat)->data;
2394:   PetscMalloc2(m,&sbaij->imax,m,&sbaij->ilen);
2395:   PetscLogObjectMemory((PetscObject)*mat,2*m*sizeof(PetscInt));

2397:   sbaij->i = i;
2398:   sbaij->j = j;
2399:   sbaij->a = a;

2401:   sbaij->singlemalloc = PETSC_FALSE;
2402:   sbaij->nonew        = -1;             /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/
2403:   sbaij->free_a       = PETSC_FALSE;
2404:   sbaij->free_ij      = PETSC_FALSE;

2406:   for (ii=0; ii<m; ii++) {
2407:     sbaij->ilen[ii] = sbaij->imax[ii] = i[ii+1] - i[ii];
2408: #if defined(PETSC_USE_DEBUG)
2409:     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]);
2410: #endif
2411:   }
2412: #if defined(PETSC_USE_DEBUG)
2413:   for (ii=0; ii<sbaij->i[m]; ii++) {
2414:     if (j[ii] < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column index at location = %d index = %d",ii,j[ii]);
2415:     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]);
2416:   }
2417: #endif

2419:   MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
2420:   MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
2421:   return(0);
2422: }

2426: PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqSBAIJ(MPI_Comm comm,Mat inmat,PetscInt n,MatReuse scall,Mat *outmat)
2427: {

2431:   MatCreateMPIMatConcatenateSeqMat_MPISBAIJ(comm,inmat,n,scall,outmat);
2432:   return(0);
2433: }