Actual source code: sbaij.c

petsc-3.5.4 2015-05-23
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");
 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:   PetscFree(a->xtoy);
154:   if (a->free_jshort) {PetscFree(a->jshort);}
155:   PetscFree(a->inew);
156:   MatDestroy(&a->parent);
157:   PetscFree(A->data);

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

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

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

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

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

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

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

257:   if (idx) {PetscFree(*idx);}
258:   if (v)   {PetscFree(*v);}
259:   return(0);
260: }

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

269:   a->getrow_utriangular = PETSC_TRUE;
270:   return(0);
271: }
274: PetscErrorCode MatRestoreRowUpperTriangular_SeqSBAIJ(Mat A)
275: {
276:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;

279:   a->getrow_utriangular = PETSC_FALSE;
280:   return(0);
281: }

285: PetscErrorCode MatTranspose_SeqSBAIJ(Mat A,MatReuse reuse,Mat *B)
286: {

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

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

307:   PetscViewerGetFormat(viewer,&format);
308:   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
309:     PetscViewerASCIIPrintf(viewer,"  block size is %D\n",bs);
310:   } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
311:     Mat        aij;
312:     const char *matname;

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

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

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

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

437:   /*
438:     This is nasty. If this is called from an originally parallel matrix
439:     then all processes call this,but only the first has the matrix so the
440:     rest should return immediately.
441:   */
442:   PetscObjectGetComm((PetscObject)draw,&comm);
443:   MPI_Comm_rank(comm,&rank);
444:   if (rank) return(0);

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

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

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

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

500: static PetscErrorCode MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer)
501: {
503:   PetscReal      xl,yl,xr,yr,w,h;
504:   PetscDraw      draw;
505:   PetscBool      isnull;

508:   PetscViewerDrawGetDraw(viewer,0,&draw);
509:   PetscDrawIsNull(draw,&isnull); if (isnull) return(0);

511:   PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);
512:   xr   = A->rmap->N; yr = A->rmap->N; h = yr/10.0; w = xr/10.0;
513:   xr  += w;    yr += h;  xl = -w;     yl = -h;
514:   PetscDrawSetCoordinates(draw,xl,yl,xr,yr);
515:   PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);
516:   PetscObjectCompose((PetscObject)A,"Zoomviewer",NULL);
517:   return(0);
518: }

522: PetscErrorCode MatView_SeqSBAIJ(Mat A,PetscViewer viewer)
523: {
525:   PetscBool      iascii,isdraw;
526:   FILE           *file = 0;

529:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
530:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
531:   if (iascii) {
532:     MatView_SeqSBAIJ_ASCII(A,viewer);
533:   } else if (isdraw) {
534:     MatView_SeqSBAIJ_Draw(A,viewer);
535:   } else {
536:     Mat        B;
537:     const char *matname;
538:     MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&B);
539:     PetscObjectGetName((PetscObject)A,&matname);
540:     PetscObjectSetName((PetscObject)B,matname);
541:     MatView(B,viewer);
542:     MatDestroy(&B);
543:     PetscViewerBinaryGetInfoPointer(viewer,&file);
544:     if (file) {
545:       fprintf(file,"-matload_block_size %d\n",(int)A->rmap->bs);
546:     }
547:   }
548:   return(0);
549: }


554: PetscErrorCode MatGetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],PetscScalar v[])
555: {
556:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
557:   PetscInt     *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
558:   PetscInt     *ai = a->i,*ailen = a->ilen;
559:   PetscInt     brow,bcol,ridx,cidx,bs=A->rmap->bs,bs2=a->bs2;
560:   MatScalar    *ap,*aa = a->a;

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


600: PetscErrorCode MatSetValuesBlocked_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
601: {
602:   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
603:   PetscErrorCode    ierr;
604:   PetscInt          *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,lastcol = -1;
605:   PetscInt          *imax      =a->imax,*ai=a->i,*ailen=a->ilen;
606:   PetscInt          *aj        =a->j,nonew=a->nonew,bs2=a->bs2,bs=A->rmap->bs,stepval;
607:   PetscBool         roworiented=a->roworiented;
608:   const PetscScalar *value     = v;
609:   MatScalar         *ap,*aa = a->a,*bap;

612:   if (roworiented) stepval = (n-1)*bs;
613:   else stepval = (m-1)*bs;

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

640:       if (col <= lastcol) low = 0;
641:       else high = nrow;

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

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

734:   PetscMalloc1(m,&ns);
735:   while (i < m) {
736:     nzx = ai[i+1] - ai[i];       /* Number of nonzeros */
737:     /* Limits the number of elements in a node to 'a->inode.limit' */
738:     for (j=i+1,blk_size=1; j<m && blk_size <a->inode.limit; ++j,++blk_size) {
739:       nzy = ai[j+1] - ai[j];
740:       if (nzy != (nzx - j + i)) break;
741:       PetscMemcmp(aj + ai[i] + j - i,aj + ai[j],nzy*sizeof(PetscInt),&flag);
742:       if (!flag) break;
743:     }
744:     ns[node_count++] = blk_size;

746:     i = j;
747:   }
748:   if (!a->inode.size && m && node_count > .9*m) {
749:     PetscFree(ns);
750:     PetscInfo2(A,"Found %D nodes out of %D rows. Not using Inode routines\n",node_count,m);
751:   } else {
752:     a->inode.node_count = node_count;

754:     PetscMalloc1(node_count,&a->inode.size);
755:     PetscLogObjectMemory((PetscObject)A,node_count*sizeof(PetscInt));
756:     PetscMemcpy(a->inode.size,ns,node_count*sizeof(PetscInt));
757:     PetscFree(ns);
758:     PetscInfo3(A,"Found %D nodes of %D. Limit used: %D. Using Inode routines\n",node_count,m,a->inode.limit);

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

797: PetscErrorCode MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode)
798: {
799:   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
801:   PetscInt       fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
802:   PetscInt       m      = A->rmap->N,*ip,N,*ailen = a->ilen;
803:   PetscInt       mbs    = a->mbs,bs2 = a->bs2,rmax = 0;
804:   MatScalar      *aa    = a->a,*ap;

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

809:   if (m) rmax = ailen[0];
810:   for (i=1; i<mbs; i++) {
811:     /* move each row back by the amount of empty slots (fshift) before it*/
812:     fshift += imax[i-1] - ailen[i-1];
813:     rmax    = PetscMax(rmax,ailen[i]);
814:     if (fshift) {
815:       ip = aj + ai[i]; ap = aa + bs2*ai[i];
816:       N  = ailen[i];
817:       for (j=0; j<N; j++) {
818:         ip[j-fshift] = ip[j];
819:         PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));
820:       }
821:     }
822:     ai[i] = ai[i-1] + ailen[i-1];
823:   }
824:   if (mbs) {
825:     fshift += imax[mbs-1] - ailen[mbs-1];
826:     ai[mbs] = ai[mbs-1] + ailen[mbs-1];
827:   }
828:   /* reset ilen and imax for each row */
829:   for (i=0; i<mbs; i++) {
830:     ailen[i] = imax[i] = ai[i+1] - ai[i];
831:   }
832:   a->nz = ai[mbs];

834:   /* diagonals may have moved, reset it */
835:   if (a->diag) {
836:     PetscMemcpy(a->diag,ai,mbs*sizeof(PetscInt));
837:   }
838:   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);

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

844:   A->info.mallocs    += a->reallocs;
845:   a->reallocs         = 0;
846:   A->info.nz_unneeded = (PetscReal)fshift*bs2;
847:   a->idiagvalid       = PETSC_FALSE;

849:   if (A->cmap->n < 65536 && A->cmap->bs == 1) {
850:     if (a->jshort && a->free_jshort) {
851:       /* when matrix data structure is changed, previous jshort must be replaced */
852:       PetscFree(a->jshort);
853:     }
854:     PetscMalloc1(a->i[A->rmap->n],&a->jshort);
855:     PetscLogObjectMemory((PetscObject)A,a->i[A->rmap->n]*sizeof(unsigned short));
856:     for (i=0; i<a->i[A->rmap->n]; i++) a->jshort[i] = a->j[i];
857:     A->ops->mult   = MatMult_SeqSBAIJ_1_ushort;
858:     A->ops->sor    = MatSOR_SeqSBAIJ_ushort;
859:     a->free_jshort = PETSC_TRUE;
860:   }
861:   return(0);
862: }

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

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


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

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

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

938:     for (l=0; l<n; l++) { /* loop over added columns */
939:       if (in[l] < 0) continue;
940: #if defined(PETSC_USE_DEBUG)
941:       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);
942: #endif
943:       col  = in[l];
944:       bcol = col/bs;              /* block col number */

946:       if (brow > bcol) {
947:         if (a->ignore_ltriangular) continue; /* ignore lower triangular values */
948:         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)");
949:       }

951:       ridx = row % bs; cidx = col % bs; /*row and col index inside the block */
952:       if ((brow==bcol && ridx<=cidx) || (brow<bcol)) {
953:         /* element value a(k,l) */
954:         if (roworiented) value = v[l + k*n];
955:         else value = v[k + l*m];

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

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

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

1009: PetscErrorCode MatICCFactor_SeqSBAIJ(Mat inA,IS row,const MatFactorInfo *info)
1010: {
1011:   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)inA->data;
1012:   Mat            outA;
1014:   PetscBool      row_identity;

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

1022:   outA            = inA;
1023:   inA->factortype = MAT_FACTOR_ICC;

1025:   MatMarkDiagonal_SeqSBAIJ(inA);
1026:   MatSeqSBAIJSetNumericFactorization_inplace(inA,row_identity);

1028:   PetscObjectReference((PetscObject)row);
1029:   ISDestroy(&a->row);
1030:   a->row = row;
1031:   PetscObjectReference((PetscObject)row);
1032:   ISDestroy(&a->col);
1033:   a->col = row;

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

1039:   if (!a->solve_work) {
1040:     PetscMalloc1((inA->rmap->N+inA->rmap->bs),&a->solve_work);
1041:     PetscLogObjectMemory((PetscObject)inA,(inA->rmap->N+inA->rmap->bs)*sizeof(PetscScalar));
1042:   }

1044:   MatCholeskyFactorNumeric(outA,inA,info);
1045:   return(0);
1046: }

1050: PetscErrorCode  MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,PetscInt *indices)
1051: {
1052:   Mat_SeqSBAIJ   *baij = (Mat_SeqSBAIJ*)mat->data;
1053:   PetscInt       i,nz,n;

1057:   nz = baij->maxnz;
1058:   n  = mat->cmap->n;
1059:   for (i=0; i<nz; i++) baij->j[i] = indices[i];

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

1064:   MatSetOption(mat,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
1065:   return(0);
1066: }

1070: /*@
1071:   MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows
1072:   in the matrix.

1074:   Input Parameters:
1075:   +  mat     - the SeqSBAIJ matrix
1076:   -  indices - the column indices

1078:   Level: advanced

1080:   Notes:
1081:   This can be called if you have precomputed the nonzero structure of the
1082:   matrix and want to provide it to the matrix object to improve the performance
1083:   of the MatSetValues() operation.

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

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

1090:   .seealso: MatCreateSeqSBAIJ
1091: @*/
1092: PetscErrorCode  MatSeqSBAIJSetColumnIndices(Mat mat,PetscInt *indices)
1093: {

1099:   PetscUseMethod(mat,"MatSeqSBAIJSetColumnIndices_C",(Mat,PetscInt*),(mat,indices));
1100:   return(0);
1101: }

1105: PetscErrorCode MatCopy_SeqSBAIJ(Mat A,Mat B,MatStructure str)
1106: {

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

1115:     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");
1116:     PetscMemcpy(b->a,a->a,(a->i[A->rmap->N])*sizeof(PetscScalar));
1117:   } else {
1118:     MatGetRowUpperTriangular(A);
1119:     MatCopy_Basic(A,B,str);
1120:     MatRestoreRowUpperTriangular(A);
1121:   }
1122:   return(0);
1123: }

1127: PetscErrorCode MatSetUp_SeqSBAIJ(Mat A)
1128: {

1132:   MatSeqSBAIJSetPreallocation_SeqSBAIJ(A,A->rmap->bs,PETSC_DEFAULT,0);
1133:   return(0);
1134: }

1138: PetscErrorCode MatSeqSBAIJGetArray_SeqSBAIJ(Mat A,PetscScalar *array[])
1139: {
1140:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;

1143:   *array = a->a;
1144:   return(0);
1145: }

1149: PetscErrorCode MatSeqSBAIJRestoreArray_SeqSBAIJ(Mat A,PetscScalar *array[])
1150: {
1152:   return(0);
1153: }

1157: PetscErrorCode MatAXPYGetPreallocation_SeqSBAIJ(Mat Y,Mat X,PetscInt *nnz)
1158: {
1159:   PetscInt       bs = Y->rmap->bs,mbs = Y->rmap->N/bs;
1160:   Mat_SeqSBAIJ   *x = (Mat_SeqSBAIJ*)X->data;
1161:   Mat_SeqSBAIJ   *y = (Mat_SeqSBAIJ*)Y->data;

1165:   /* Set the number of nonzeros in the new matrix */
1166:   MatAXPYGetPreallocation_SeqX_private(mbs,x->i,x->j,y->i,y->j,nnz);
1167:   return(0);
1168: }

1172: PetscErrorCode MatAXPY_SeqSBAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str)
1173: {
1174:   Mat_SeqSBAIJ   *x=(Mat_SeqSBAIJ*)X->data, *y=(Mat_SeqSBAIJ*)Y->data;
1176:   PetscInt       bs=Y->rmap->bs,bs2=bs*bs;
1177:   PetscBLASInt   one = 1;

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

1217: PetscErrorCode MatIsSymmetric_SeqSBAIJ(Mat A,PetscReal tol,PetscBool  *flg)
1218: {
1220:   *flg = PETSC_TRUE;
1221:   return(0);
1222: }

1226: PetscErrorCode MatIsStructurallySymmetric_SeqSBAIJ(Mat A,PetscBool  *flg)
1227: {
1229:   *flg = PETSC_TRUE;
1230:   return(0);
1231: }

1235: PetscErrorCode MatIsHermitian_SeqSBAIJ(Mat A,PetscReal tol,PetscBool  *flg)
1236: {
1238:   *flg = PETSC_FALSE;
1239:   return(0);
1240: }

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

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

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

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

1283:   /* fix right hand side if needed */
1284:   if (x && b) {
1285:     VecGetArrayRead(x,&xx);
1286:     VecGetArray(b,&bb);
1287:     vecs = PETSC_TRUE;
1288:   }

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

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

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

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

1496: PetscErrorCode  MatStoreValues_SeqSBAIJ(Mat mat)
1497: {
1498:   Mat_SeqSBAIJ   *aij = (Mat_SeqSBAIJ*)mat->data;
1499:   PetscInt       nz   = aij->i[mat->rmap->N]*mat->rmap->bs*aij->bs2;

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

1505:   /* allocate space for values if not already there */
1506:   if (!aij->saved_values) {
1507:     PetscMalloc1((nz+1),&aij->saved_values);
1508:   }

1510:   /* copy values over */
1511:   PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));
1512:   return(0);
1513: }

1517: PetscErrorCode  MatRetrieveValues_SeqSBAIJ(Mat mat)
1518: {
1519:   Mat_SeqSBAIJ   *aij = (Mat_SeqSBAIJ*)mat->data;
1521:   PetscInt       nz = aij->i[mat->rmap->N]*mat->rmap->bs*aij->bs2;

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

1527:   /* copy values over */
1528:   PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));
1529:   return(0);
1530: }

1534: PetscErrorCode  MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat B,PetscInt bs,PetscInt nz,PetscInt *nnz)
1535: {
1536:   Mat_SeqSBAIJ   *b = (Mat_SeqSBAIJ*)B->data;
1538:   PetscInt       i,mbs,bs2;
1539:   PetscBool      skipallocation = PETSC_FALSE,flg = PETSC_FALSE,realalloc = PETSC_FALSE;

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

1545:   MatSetBlockSize(B,PetscAbs(bs));
1546:   PetscLayoutSetUp(B->rmap);
1547:   PetscLayoutSetUp(B->cmap);
1548:   PetscLayoutGetBlockSize(B->rmap,&bs);

1550:   mbs = B->rmap->N/bs;
1551:   bs2 = bs*bs;

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

1555:   if (nz == MAT_SKIP_ALLOCATION) {
1556:     skipallocation = PETSC_TRUE;
1557:     nz             = 0;
1558:   }

1560:   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3;
1561:   if (nz < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %D",nz);
1562:   if (nnz) {
1563:     for (i=0; i<mbs; i++) {
1564:       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]);
1565:       if (nnz[i] > mbs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be greater than block row length: local row %D value %D rowlength %D",i,nnz[i],mbs);
1566:     }
1567:   }

1569:   B->ops->mult             = MatMult_SeqSBAIJ_N;
1570:   B->ops->multadd          = MatMultAdd_SeqSBAIJ_N;
1571:   B->ops->multtranspose    = MatMult_SeqSBAIJ_N;
1572:   B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_N;

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

1622:   b->mbs = mbs;
1623:   b->nbs = mbs;
1624:   if (!skipallocation) {
1625:     if (!b->imax) {
1626:       PetscMalloc2(mbs,&b->imax,mbs,&b->ilen);

1628:       b->free_imax_ilen = PETSC_TRUE;

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

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

1652:     b->singlemalloc = PETSC_TRUE;

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

1658:     b->free_a  = PETSC_TRUE;
1659:     b->free_ij = PETSC_TRUE;
1660:   } else {
1661:     b->free_a  = PETSC_FALSE;
1662:     b->free_ij = PETSC_FALSE;
1663:   }

1665:   B->rmap->bs = bs;
1666:   b->bs2      = bs2;
1667:   b->nz       = 0;
1668:   b->maxnz    = nz;

1670:   b->inew    = 0;
1671:   b->jnew    = 0;
1672:   b->anew    = 0;
1673:   b->a2anew  = 0;
1674:   b->permute = PETSC_FALSE;
1675:   if (realalloc) {MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);}
1676:   return(0);
1677: }

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

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

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

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

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

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

1804: PETSC_EXTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqAIJ(Mat, MatType,MatReuse,Mat*);
1805: PETSC_EXTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqBAIJ(Mat, MatType,MatReuse,Mat*);

1809: PETSC_EXTERN PetscErrorCode MatGetFactor_seqsbaij_petsc(Mat A,MatFactorType ftype,Mat *B)
1810: {
1811:   PetscInt       n = A->rmap->n;

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

1824:     (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqSBAIJ;
1825:     (*B)->ops->iccfactorsymbolic      = MatICCFactorSymbolic_SeqSBAIJ;
1826:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
1827:   (*B)->factortype = ftype;
1828:   return(0);
1829: }

1833: PetscErrorCode MatGetFactorAvailable_seqsbaij_petsc(Mat A,MatFactorType ftype,PetscBool  *flg)
1834: {
1836:   if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) {
1837:     *flg = PETSC_TRUE;
1838:   } else {
1839:     *flg = PETSC_FALSE;
1840:   }
1841:   return(0);
1842: }

1844: #if defined(PETSC_HAVE_MUMPS)
1845: PETSC_EXTERN PetscErrorCode MatGetFactor_sbaij_mumps(Mat,MatFactorType,Mat*);
1846: #endif
1847: #if defined(PETSC_HAVE_PASTIX)
1848: PETSC_EXTERN PetscErrorCode MatGetFactor_seqsbaij_pastix(Mat,MatFactorType,Mat*);
1849: #endif
1850: #if defined(PETSC_HAVE_SUITESPARSE)
1851: PETSC_EXTERN PetscErrorCode MatGetFactor_seqsbaij_cholmod(Mat,MatFactorType,Mat*);
1852: #endif
1853: PETSC_EXTERN PetscErrorCode MatGetFactor_seqsbaij_sbstrm(Mat,MatFactorType,Mat*);

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

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

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

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


1870:   Level: beginner

1872:   .seealso: MatCreateSeqSBAIJ
1873: M*/

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

1879: PETSC_EXTERN PetscErrorCode MatCreate_SeqSBAIJ(Mat B)
1880: {
1881:   Mat_SeqSBAIJ   *b;
1883:   PetscMPIInt    size;
1884:   PetscBool      no_unroll = PETSC_FALSE,no_inode = PETSC_FALSE;

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

1890:   PetscNewLog(B,&b);
1891:   B->data = (void*)b;
1892:   PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));

1894:   B->ops->destroy    = MatDestroy_SeqSBAIJ;
1895:   B->ops->view       = MatView_SeqSBAIJ;
1896:   b->row             = 0;
1897:   b->icol            = 0;
1898:   b->reallocs        = 0;
1899:   b->saved_values    = 0;
1900:   b->inode.limit     = 5;
1901:   b->inode.max_limit = 5;

1903:   b->roworiented        = PETSC_TRUE;
1904:   b->nonew              = 0;
1905:   b->diag               = 0;
1906:   b->solve_work         = 0;
1907:   b->mult_work          = 0;
1908:   B->spptr              = 0;
1909:   B->info.nz_unneeded   = (PetscReal)b->maxnz*b->bs2;
1910:   b->keepnonzeropattern = PETSC_FALSE;
1911:   b->xtoy               = 0;
1912:   b->XtoY               = 0;

1914:   b->inew    = 0;
1915:   b->jnew    = 0;
1916:   b->anew    = 0;
1917:   b->a2anew  = 0;
1918:   b->permute = PETSC_FALSE;

1920:   b->ignore_ltriangular = PETSC_TRUE;

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

1924:   b->getrow_utriangular = PETSC_FALSE;

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

1928: #if defined(PETSC_HAVE_PASTIX)
1929:   PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_pastix_C",MatGetFactor_seqsbaij_pastix);
1930: #endif
1931: #if defined(PETSC_HAVE_MUMPS)
1932:   PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_mumps_C",MatGetFactor_sbaij_mumps);
1933: #endif
1934: #if defined(PETSC_HAVE_SUITESPARSE)
1935:   PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_cholmod_C",MatGetFactor_seqsbaij_cholmod);
1936: #endif
1937:   PetscObjectComposeFunction((PetscObject)B,"MatGetFactorAvailable_petsc_C",MatGetFactorAvailable_seqsbaij_petsc);
1938:   PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_petsc_C",MatGetFactor_seqsbaij_petsc);
1939:   PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_sbstrm_C",MatGetFactor_seqsbaij_sbstrm);
1940:   PetscObjectComposeFunction((PetscObject)B,"MatStoreValues_C",MatStoreValues_SeqSBAIJ);
1941:   PetscObjectComposeFunction((PetscObject)B,"MatRetrieveValues_C",MatRetrieveValues_SeqSBAIJ);
1942:   PetscObjectComposeFunction((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C",MatSeqSBAIJSetColumnIndices_SeqSBAIJ);
1943:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_seqaij_C",MatConvert_SeqSBAIJ_SeqAIJ);
1944:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_seqbaij_C",MatConvert_SeqSBAIJ_SeqBAIJ);
1945:   PetscObjectComposeFunction((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",MatSeqSBAIJSetPreallocation_SeqSBAIJ);
1946:   PetscObjectComposeFunction((PetscObject)B,"MatSeqSBAIJSetPreallocationCSR_C",MatSeqSBAIJSetPreallocationCSR_SeqSBAIJ);
1947:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_seqsbstrm_C",MatConvert_SeqSBAIJ_SeqSBSTRM);

1949:   B->symmetric                  = PETSC_TRUE;
1950:   B->structurally_symmetric     = PETSC_TRUE;
1951:   B->symmetric_set              = PETSC_TRUE;
1952:   B->structurally_symmetric_set = PETSC_TRUE;

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

1956:   PetscOptionsBegin(PetscObjectComm((PetscObject)B),((PetscObject)B)->prefix,"Options for SEQSBAIJ matrix","Mat");
1957:   PetscOptionsBool("-mat_no_unroll","Do not optimize for inodes (slower)",NULL,no_unroll,&no_unroll,NULL);
1958:   if (no_unroll) {
1959:     PetscInfo(B,"Not using Inode routines due to -mat_no_unroll\n");
1960:   }
1961:   PetscOptionsBool("-mat_no_inode","Do not optimize for inodes (slower)",NULL,no_inode,&no_inode,NULL);
1962:   if (no_inode) {
1963:     PetscInfo(B,"Not using Inode routines due to -mat_no_inode\n");
1964:   }
1965:   PetscOptionsInt("-mat_inode_limit","Do not use inodes larger then this value",NULL,b->inode.limit,&b->inode.limit,NULL);
1966:   PetscOptionsEnd();
1967:   b->inode.use = (PetscBool)(!(no_unroll || no_inode));
1968:   if (b->inode.limit > b->inode.max_limit) b->inode.limit = b->inode.max_limit;
1969:   return(0);
1970: }

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

1981:    Collective on Mat

1983:    Input Parameters:
1984: +  B - the symmetric matrix
1985: .  bs - size of block
1986: .  nz - number of block nonzeros per block row (same for all rows)
1987: -  nnz - array containing the number of block nonzeros in the upper triangular plus
1988:          diagonal portion of each block (possibly different for each block row) or NULL

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

1995:    Level: intermediate

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

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

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


2010: .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateSBAIJ()
2011: @*/
2012: PetscErrorCode  MatSeqSBAIJSetPreallocation(Mat B,PetscInt bs,PetscInt nz,const PetscInt nnz[])
2013: {

2020:   PetscTryMethod(B,"MatSeqSBAIJSetPreallocation_C",(Mat,PetscInt,PetscInt,const PetscInt[]),(B,bs,nz,nnz));
2021:   return(0);
2022: }

2024: #undef  __FUNCT__
2026: /*@C
2027:    MatSeqSBAIJSetPreallocationCSR - Allocates memory for a sparse sequential matrix in symmetric block AIJ format.

2029:    Input Parameters:
2030: +  B - the matrix
2031: .  i - the indices into j for the start of each local row (starts with zero)
2032: .  j - the column indices for each local row (starts with zero) these must be sorted for each row
2033: -  v - optional values in the matrix

2035:    Level: developer

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

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

2046: .seealso: MatCreate(), MatCreateSeqSBAIJ(), MatSetValuesBlocked(), MatSeqSBAIJSetPreallocation(), MATSEQSBAIJ
2047: @*/
2048: PetscErrorCode MatSeqSBAIJSetPreallocationCSR(Mat B,PetscInt bs,const PetscInt i[],const PetscInt j[], const PetscScalar v[])
2049: {

2056:   PetscTryMethod(B,"MatSeqSBAIJSetPreallocationCSR_C",(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]),(B,bs,i,j,v));
2057:   return(0);
2058: }

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

2069:    Collective on MPI_Comm

2071:    Input Parameters:
2072: +  comm - MPI communicator, set to PETSC_COMM_SELF
2073: .  bs - size of block
2074: .  m - number of rows, or number of columns
2075: .  nz - number of block nonzeros per block row (same for all rows)
2076: -  nnz - array containing the number of block nonzeros in the upper triangular plus
2077:          diagonal portion of each block (possibly different for each block row) or NULL

2079:    Output Parameter:
2080: .  A - the symmetric matrix

2082:    Options Database Keys:
2083: .   -mat_no_unroll - uses code that does not unroll the loops in the
2084:                      block calculations (much slower)
2085: .    -mat_block_size - size of the blocks to use

2087:    Level: intermediate

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

2093:    Notes:
2094:    The number of rows and columns must be divisible by blocksize.
2095:    This matrix type does not support complex Hermitian operation.

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

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

2103: .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateSBAIJ()
2104: @*/
2105: PetscErrorCode  MatCreateSeqSBAIJ(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A)
2106: {

2110:   MatCreate(comm,A);
2111:   MatSetSizes(*A,m,n,m,n);
2112:   MatSetType(*A,MATSEQSBAIJ);
2113:   MatSeqSBAIJSetPreallocation_SeqSBAIJ(*A,bs,nz,(PetscInt*)nnz);
2114:   return(0);
2115: }

2119: PetscErrorCode MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
2120: {
2121:   Mat            C;
2122:   Mat_SeqSBAIJ   *c,*a = (Mat_SeqSBAIJ*)A->data;
2124:   PetscInt       i,mbs = a->mbs,nz = a->nz,bs2 =a->bs2;

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

2129:   *B   = 0;
2130:   MatCreate(PetscObjectComm((PetscObject)A),&C);
2131:   MatSetSizes(C,A->rmap->N,A->cmap->n,A->rmap->N,A->cmap->n);
2132:   MatSetType(C,MATSEQSBAIJ);
2133:   PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));
2134:   c    = (Mat_SeqSBAIJ*)C->data;

2136:   C->preallocated       = PETSC_TRUE;
2137:   C->factortype         = A->factortype;
2138:   c->row                = 0;
2139:   c->icol               = 0;
2140:   c->saved_values       = 0;
2141:   c->keepnonzeropattern = a->keepnonzeropattern;
2142:   C->assembled          = PETSC_TRUE;

2144:   PetscLayoutReference(A->rmap,&C->rmap);
2145:   PetscLayoutReference(A->cmap,&C->cmap);
2146:   c->bs2 = a->bs2;
2147:   c->mbs = a->mbs;
2148:   c->nbs = a->nbs;

2150:   if  (cpvalues == MAT_SHARE_NONZERO_PATTERN) {
2151:     c->imax           = a->imax;
2152:     c->ilen           = a->ilen;
2153:     c->free_imax_ilen = PETSC_FALSE;
2154:   } else {
2155:     PetscMalloc2((mbs+1),&c->imax,(mbs+1),&c->ilen);
2156:     PetscLogObjectMemory((PetscObject)C,2*(mbs+1)*sizeof(PetscInt));
2157:     for (i=0; i<mbs; i++) {
2158:       c->imax[i] = a->imax[i];
2159:       c->ilen[i] = a->ilen[i];
2160:     }
2161:     c->free_imax_ilen = PETSC_TRUE;
2162:   }

2164:   /* allocate the matrix space */
2165:   if (cpvalues == MAT_SHARE_NONZERO_PATTERN) {
2166:     PetscMalloc1(bs2*nz,&c->a);
2167:     PetscLogObjectMemory((PetscObject)C,nz*bs2*sizeof(MatScalar));
2168:     c->i            = a->i;
2169:     c->j            = a->j;
2170:     c->singlemalloc = PETSC_FALSE;
2171:     c->free_a       = PETSC_TRUE;
2172:     c->free_ij      = PETSC_FALSE;
2173:     c->parent       = A;
2174:     PetscObjectReference((PetscObject)A);
2175:     MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
2176:     MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
2177:   } else {
2178:     PetscMalloc3(bs2*nz,&c->a,nz,&c->j,mbs+1,&c->i);
2179:     PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(PetscInt));
2180:     PetscLogObjectMemory((PetscObject)C,(mbs+1)*sizeof(PetscInt) + nz*(bs2*sizeof(MatScalar) + sizeof(PetscInt)));
2181:     c->singlemalloc = PETSC_TRUE;
2182:     c->free_a       = PETSC_TRUE;
2183:     c->free_ij      = PETSC_TRUE;
2184:   }
2185:   if (mbs > 0) {
2186:     if (cpvalues != MAT_SHARE_NONZERO_PATTERN) {
2187:       PetscMemcpy(c->j,a->j,nz*sizeof(PetscInt));
2188:     }
2189:     if (cpvalues == MAT_COPY_VALUES) {
2190:       PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));
2191:     } else {
2192:       PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));
2193:     }
2194:     if (a->jshort) {
2195:       /* cannot share jshort, it is reallocated in MatAssemblyEnd_SeqSBAIJ() */
2196:       /* if the parent matrix is reassembled, this child matrix will never notice */
2197:       PetscMalloc1(nz,&c->jshort);
2198:       PetscLogObjectMemory((PetscObject)C,nz*sizeof(unsigned short));
2199:       PetscMemcpy(c->jshort,a->jshort,nz*sizeof(unsigned short));

2201:       c->free_jshort = PETSC_TRUE;
2202:     }
2203:   }

2205:   c->roworiented = a->roworiented;
2206:   c->nonew       = a->nonew;

2208:   if (a->diag) {
2209:     if (cpvalues == MAT_SHARE_NONZERO_PATTERN) {
2210:       c->diag      = a->diag;
2211:       c->free_diag = PETSC_FALSE;
2212:     } else {
2213:       PetscMalloc1(mbs,&c->diag);
2214:       PetscLogObjectMemory((PetscObject)C,mbs*sizeof(PetscInt));
2215:       for (i=0; i<mbs; i++) c->diag[i] = a->diag[i];
2216:       c->free_diag = PETSC_TRUE;
2217:     }
2218:   }
2219:   c->nz         = a->nz;
2220:   c->maxnz      = a->nz; /* Since we allocate exactly the right amount */
2221:   c->solve_work = 0;
2222:   c->mult_work  = 0;

2224:   *B   = C;
2225:   PetscFunctionListDuplicate(((PetscObject)A)->qlist,&((PetscObject)C)->qlist);
2226:   return(0);
2227: }

2231: PetscErrorCode MatLoad_SeqSBAIJ(Mat newmat,PetscViewer viewer)
2232: {
2233:   Mat_SeqSBAIJ   *a;
2235:   int            fd;
2236:   PetscMPIInt    size;
2237:   PetscInt       i,nz,header[4],*rowlengths=0,M,N,bs = newmat->rmap->bs;
2238:   PetscInt       *mask,mbs,*jj,j,rowcount,nzcount,k,*s_browlengths,maskcount;
2239:   PetscInt       kmax,jcount,block,idx,point,nzcountb,extra_rows,rows,cols;
2240:   PetscInt       *masked,nmask,tmp,bs2,ishift;
2241:   PetscScalar    *aa;
2242:   MPI_Comm       comm;

2245:   PetscObjectGetComm((PetscObject)viewer,&comm);
2246:   PetscOptionsGetInt(((PetscObject)newmat)->prefix,"-matload_block_size",&bs,NULL);
2247:   if (bs < 0) bs = 1;
2248:   bs2  = bs*bs;

2250:   MPI_Comm_size(comm,&size);
2251:   if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"view must have one processor");
2252:   PetscViewerBinaryGetDescriptor(viewer,&fd);
2253:   PetscBinaryRead(fd,header,4,PETSC_INT);
2254:   if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"not Mat object");
2255:   M = header[1]; N = header[2]; nz = header[3];

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

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

2261:   /*
2262:      This code adds extra rows to make sure the number of rows is
2263:     divisible by the blocksize
2264:   */
2265:   mbs        = M/bs;
2266:   extra_rows = bs - M + bs*(mbs);
2267:   if (extra_rows == bs) extra_rows = 0;
2268:   else                  mbs++;
2269:   if (extra_rows) {
2270:     PetscInfo(viewer,"Padding loaded matrix to match blocksize\n");
2271:   }

2273:   /* Set global sizes if not already set */
2274:   if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) {
2275:     MatSetSizes(newmat,PETSC_DECIDE,PETSC_DECIDE,M+extra_rows,N+extra_rows);
2276:   } else { /* Check if the matrix global sizes are correct */
2277:     MatGetSize(newmat,&rows,&cols);
2278:     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);
2279:   }

2281:   /* read in row lengths */
2282:   PetscMalloc1((M+extra_rows),&rowlengths);
2283:   PetscBinaryRead(fd,rowlengths,M,PETSC_INT);
2284:   for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;

2286:   /* read in column indices */
2287:   PetscMalloc1((nz+extra_rows),&jj);
2288:   PetscBinaryRead(fd,jj,nz,PETSC_INT);
2289:   for (i=0; i<extra_rows; i++) jj[nz+i] = M+i;

2291:   /* loop over row lengths determining block row lengths */
2292:   PetscCalloc1(mbs,&s_browlengths);
2293:   PetscMalloc2(mbs,&mask,mbs,&masked);
2294:   PetscMemzero(mask,mbs*sizeof(PetscInt));
2295:   rowcount = 0;
2296:   nzcount  = 0;
2297:   for (i=0; i<mbs; i++) {
2298:     nmask = 0;
2299:     for (j=0; j<bs; j++) {
2300:       kmax = rowlengths[rowcount];
2301:       for (k=0; k<kmax; k++) {
2302:         tmp = jj[nzcount++]/bs;   /* block col. index */
2303:         if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;}
2304:       }
2305:       rowcount++;
2306:     }
2307:     s_browlengths[i] += nmask;

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

2313:   /* Do preallocation */
2314:   MatSeqSBAIJSetPreallocation_SeqSBAIJ(newmat,bs,0,s_browlengths);
2315:   a    = (Mat_SeqSBAIJ*)newmat->data;

2317:   /* set matrix "i" values */
2318:   a->i[0] = 0;
2319:   for (i=1; i<= mbs; i++) {
2320:     a->i[i]      = a->i[i-1] + s_browlengths[i-1];
2321:     a->ilen[i-1] = s_browlengths[i-1];
2322:   }
2323:   a->nz = a->i[mbs];

2325:   /* read in nonzero values */
2326:   PetscMalloc1((nz+extra_rows),&aa);
2327:   PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);
2328:   for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0;

2330:   /* set "a" and "j" values into matrix */
2331:   nzcount = 0; jcount = 0;
2332:   for (i=0; i<mbs; i++) {
2333:     nzcountb = nzcount;
2334:     nmask    = 0;
2335:     for (j=0; j<bs; j++) {
2336:       kmax = rowlengths[i*bs+j];
2337:       for (k=0; k<kmax; k++) {
2338:         tmp = jj[nzcount++]/bs; /* block col. index */
2339:         if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;}
2340:       }
2341:     }
2342:     /* sort the masked values */
2343:     PetscSortInt(nmask,masked);

2345:     /* set "j" values into matrix */
2346:     maskcount = 1;
2347:     for (j=0; j<nmask; j++) {
2348:       a->j[jcount++]  = masked[j];
2349:       mask[masked[j]] = maskcount++;
2350:     }

2352:     /* set "a" values into matrix */
2353:     ishift = bs2*a->i[i];
2354:     for (j=0; j<bs; j++) {
2355:       kmax = rowlengths[i*bs+j];
2356:       for (k=0; k<kmax; k++) {
2357:         tmp = jj[nzcountb]/bs;        /* block col. index */
2358:         if (tmp >= i) {
2359:           block     = mask[tmp] - 1;
2360:           point     = jj[nzcountb] - bs*tmp;
2361:           idx       = ishift + bs2*block + j + bs*point;
2362:           a->a[idx] = aa[nzcountb];
2363:         }
2364:         nzcountb++;
2365:       }
2366:     }
2367:     /* zero out the mask elements we set */
2368:     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
2369:   }
2370:   if (jcount != a->nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix");

2372:   PetscFree(rowlengths);
2373:   PetscFree(s_browlengths);
2374:   PetscFree(aa);
2375:   PetscFree(jj);
2376:   PetscFree2(mask,masked);

2378:   MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);
2379:   MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);
2380:   return(0);
2381: }

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

2389:      Collective on MPI_Comm

2391:    Input Parameters:
2392: +  comm - must be an MPI communicator of size 1
2393: .  bs - size of block
2394: .  m - number of rows
2395: .  n - number of columns
2396: .  i - row indices
2397: .  j - column indices
2398: -  a - matrix values

2400:    Output Parameter:
2401: .  mat - the matrix

2403:    Level: advanced

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

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

2411:        The i and j indices are 0 based

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

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

2418: @*/
2419: PetscErrorCode  MatCreateSeqSBAIJWithArrays(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt *i,PetscInt *j,PetscScalar *a,Mat *mat)
2420: {
2422:   PetscInt       ii;
2423:   Mat_SeqSBAIJ   *sbaij;

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

2429:   MatCreate(comm,mat);
2430:   MatSetSizes(*mat,m,n,m,n);
2431:   MatSetType(*mat,MATSEQSBAIJ);
2432:   MatSeqSBAIJSetPreallocation_SeqSBAIJ(*mat,bs,MAT_SKIP_ALLOCATION,0);
2433:   sbaij = (Mat_SeqSBAIJ*)(*mat)->data;
2434:   PetscMalloc2(m,&sbaij->imax,m,&sbaij->ilen);
2435:   PetscLogObjectMemory((PetscObject)*mat,2*m*sizeof(PetscInt));

2437:   sbaij->i = i;
2438:   sbaij->j = j;
2439:   sbaij->a = a;

2441:   sbaij->singlemalloc = PETSC_FALSE;
2442:   sbaij->nonew        = -1;             /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/
2443:   sbaij->free_a       = PETSC_FALSE;
2444:   sbaij->free_ij      = PETSC_FALSE;

2446:   for (ii=0; ii<m; ii++) {
2447:     sbaij->ilen[ii] = sbaij->imax[ii] = i[ii+1] - i[ii];
2448: #if defined(PETSC_USE_DEBUG)
2449:     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]);
2450: #endif
2451:   }
2452: #if defined(PETSC_USE_DEBUG)
2453:   for (ii=0; ii<sbaij->i[m]; ii++) {
2454:     if (j[ii] < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column index at location = %d index = %d",ii,j[ii]);
2455:     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]);
2456:   }
2457: #endif

2459:   MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
2460:   MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
2461:   return(0);
2462: }