Actual source code: sbaij.c

petsc-3.7.3 2016-08-01
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);
 15: #if defined(PETSC_HAVE_ELEMENTAL)
 16: PETSC_INTERN PetscErrorCode MatConvert_SeqSBAIJ_Elemental(Mat,MatType,MatReuse,Mat*);
 17: #endif

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

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

 52: PetscErrorCode MatMarkDiagonal_SeqSBAIJ(Mat A)
 53: {
 54:   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
 56:   PetscInt       i,j;

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

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

 86:   *nn = n;
 87:   if (!ia) return(0);
 88:   if (symmetric) {
 89:     MatToSymmetricIJ_SeqAIJ(n,a->i,a->j,PETSC_FALSE,0,0,&tia,&tja);
 90:     nz   = tia[n];
 91:   } else {
 92:     tia = a->i; tja = a->j;
 93:   }

 95:   if (!blockcompressed && bs > 1) {
 96:     (*nn) *= bs;
 97:     /* malloc & create the natural set of indices */
 98:     PetscMalloc1((n+1)*bs,ia);
 99:     if (n) {
100:       (*ia)[0] = oshift;
101:       for (j=1; j<bs; j++) {
102:         (*ia)[j] = (tia[1]-tia[0])*bs+(*ia)[j-1];
103:       }
104:     }

106:     for (i=1; i<n; i++) {
107:       (*ia)[i*bs] = (tia[i]-tia[i-1])*bs + (*ia)[i*bs-1];
108:       for (j=1; j<bs; j++) {
109:         (*ia)[i*bs+j] = (tia[i+1]-tia[i])*bs + (*ia)[i*bs+j-1];
110:       }
111:     }
112:     if (n) {
113:       (*ia)[n*bs] = (tia[n]-tia[n-1])*bs + (*ia)[n*bs-1];
114:     }

116:     if (inja) {
117:       PetscMalloc1(nz*bs*bs,ja);
118:       cnt = 0;
119:       for (i=0; i<n; i++) {
120:         for (j=0; j<bs; j++) {
121:           for (k=tia[i]; k<tia[i+1]; k++) {
122:             for (l=0; l<bs; l++) {
123:               (*ja)[cnt++] = bs*tja[k] + l;
124:             }
125:           }
126:         }
127:       }
128:     }

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

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

168:   if (!ia) return(0);
169:   if ((!blockcompressed && A->rmap->bs > 1) || (symmetric || oshift == 1)) {
170:     PetscFree(*ia);
171:     if (ja) {PetscFree(*ja);}
172:   }
173:   return(0);
174: }

178: PetscErrorCode MatDestroy_SeqSBAIJ(Mat A)
179: {
180:   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;

184: #if defined(PETSC_USE_LOG)
185:   PetscLogObjectState((PetscObject)A,"Rows=%D, NZ=%D",A->rmap->N,a->nz);
186: #endif
187:   MatSeqXAIJFreeAIJ(A,&a->a,&a->j,&a->i);
188:   if (a->free_diag) {PetscFree(a->diag);}
189:   ISDestroy(&a->row);
190:   ISDestroy(&a->col);
191:   ISDestroy(&a->icol);
192:   PetscFree(a->idiag);
193:   PetscFree(a->inode.size);
194:   if (a->free_imax_ilen) {PetscFree2(a->imax,a->ilen);}
195:   PetscFree(a->solve_work);
196:   PetscFree(a->sor_work);
197:   PetscFree(a->solves_work);
198:   PetscFree(a->mult_work);
199:   PetscFree(a->saved_values);
200:   if (a->free_jshort) {PetscFree(a->jshort);}
201:   PetscFree(a->inew);
202:   MatDestroy(&a->parent);
203:   PetscFree(A->data);

205:   PetscObjectChangeTypeName((PetscObject)A,0);
206:   PetscObjectComposeFunction((PetscObject)A,"MatStoreValues_C",NULL);
207:   PetscObjectComposeFunction((PetscObject)A,"MatRetrieveValues_C",NULL);
208:   PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetColumnIndices_C",NULL);
209:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqaij_C",NULL);
210:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqbaij_C",NULL);
211:   PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetPreallocation_C",NULL);
212:   PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetPreallocationCSR_C",NULL);
213:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqsbstrm_C",NULL);
214: #if defined(PETSC_HAVE_ELEMENTAL)
215:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_elemental_C",NULL);
216: #endif
217:   return(0);
218: }

222: PetscErrorCode MatSetOption_SeqSBAIJ(Mat A,MatOption op,PetscBool flg)
223: {
224:   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;

228:   switch (op) {
229:   case MAT_ROW_ORIENTED:
230:     a->roworiented = flg;
231:     break;
232:   case MAT_KEEP_NONZERO_PATTERN:
233:     a->keepnonzeropattern = flg;
234:     break;
235:   case MAT_NEW_NONZERO_LOCATIONS:
236:     a->nonew = (flg ? 0 : 1);
237:     break;
238:   case MAT_NEW_NONZERO_LOCATION_ERR:
239:     a->nonew = (flg ? -1 : 0);
240:     break;
241:   case MAT_NEW_NONZERO_ALLOCATION_ERR:
242:     a->nonew = (flg ? -2 : 0);
243:     break;
244:   case MAT_UNUSED_NONZERO_LOCATION_ERR:
245:     a->nounused = (flg ? -1 : 0);
246:     break;
247:   case MAT_NEW_DIAGONALS:
248:   case MAT_IGNORE_OFF_PROC_ENTRIES:
249:   case MAT_USE_HASH_TABLE:
250:     PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
251:     break;
252:   case MAT_HERMITIAN:
253:     if (!A->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatAssemblyEnd() first");
254:     if (A->cmap->n < 65536 && A->cmap->bs == 1) {
255:       A->ops->mult = MatMult_SeqSBAIJ_1_Hermitian_ushort;
256:     } else if (A->cmap->bs == 1) {
257:       A->ops->mult = MatMult_SeqSBAIJ_1_Hermitian;
258:     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for Hermitian with block size greater than 1");
259:     break;
260:   case MAT_SPD:
261:     /* These options are handled directly by MatSetOption() */
262:     break;
263:   case MAT_SYMMETRIC:
264:   case MAT_STRUCTURALLY_SYMMETRIC:
265:   case MAT_SYMMETRY_ETERNAL:
266:     /* These options are handled directly by MatSetOption() */
267:     break;
268:   case MAT_IGNORE_LOWER_TRIANGULAR:
269:     a->ignore_ltriangular = flg;
270:     break;
271:   case MAT_ERROR_LOWER_TRIANGULAR:
272:     a->ignore_ltriangular = flg;
273:     break;
274:   case MAT_GETROW_UPPERTRIANGULAR:
275:     a->getrow_utriangular = flg;
276:     break;
277:   default:
278:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %d",op);
279:   }
280:   return(0);
281: }

285: PetscErrorCode MatGetRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
286: {
287:   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;

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

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

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

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

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

317:   a->getrow_utriangular = PETSC_TRUE;
318:   return(0);
319: }
322: PetscErrorCode MatRestoreRowUpperTriangular_SeqSBAIJ(Mat A)
323: {
324:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;

327:   a->getrow_utriangular = PETSC_FALSE;
328:   return(0);
329: }

333: PetscErrorCode MatTranspose_SeqSBAIJ(Mat A,MatReuse reuse,Mat *B)
334: {

338:   if (reuse == MAT_INITIAL_MATRIX || *B != A) {
339:     MatDuplicate(A,MAT_COPY_VALUES,B);
340:   }
341:   return(0);
342: }

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

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

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

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

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

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

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

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

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

539: static PetscErrorCode MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer)
540: {
542:   PetscReal      xl,yl,xr,yr,w,h;
543:   PetscDraw      draw;
544:   PetscBool      isnull;

547:   PetscViewerDrawGetDraw(viewer,0,&draw);
548:   PetscDrawIsNull(draw,&isnull);
549:   if (isnull) return(0);

551:   xr   = A->rmap->N; yr = A->rmap->N; h = yr/10.0; w = xr/10.0;
552:   xr  += w;          yr += h;        xl = -w;     yl = -h;
553:   PetscDrawSetCoordinates(draw,xl,yl,xr,yr);
554:   PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);
555:   PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);
556:   PetscObjectCompose((PetscObject)A,"Zoomviewer",NULL);
557:   PetscDrawSave(draw);
558:   return(0);
559: }

563: PetscErrorCode MatView_SeqSBAIJ(Mat A,PetscViewer viewer)
564: {
566:   PetscBool      iascii,isdraw;
567:   FILE           *file = 0;

570:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
571:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
572:   if (iascii) {
573:     MatView_SeqSBAIJ_ASCII(A,viewer);
574:   } else if (isdraw) {
575:     MatView_SeqSBAIJ_Draw(A,viewer);
576:   } else {
577:     Mat        B;
578:     const char *matname;
579:     MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&B);
580:     PetscObjectGetName((PetscObject)A,&matname);
581:     PetscObjectSetName((PetscObject)B,matname);
582:     MatView(B,viewer);
583:     MatDestroy(&B);
584:     PetscViewerBinaryGetInfoPointer(viewer,&file);
585:     if (file) {
586:       fprintf(file,"-matload_block_size %d\n",(int)A->rmap->bs);
587:     }
588:   }
589:   return(0);
590: }


595: PetscErrorCode MatGetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],PetscScalar v[])
596: {
597:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
598:   PetscInt     *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
599:   PetscInt     *ai = a->i,*ailen = a->ilen;
600:   PetscInt     brow,bcol,ridx,cidx,bs=A->rmap->bs,bs2=a->bs2;
601:   MatScalar    *ap,*aa = a->a;

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


641: PetscErrorCode MatSetValuesBlocked_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
642: {
643:   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
644:   PetscErrorCode    ierr;
645:   PetscInt          *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,lastcol = -1;
646:   PetscInt          *imax      =a->imax,*ai=a->i,*ailen=a->ilen;
647:   PetscInt          *aj        =a->j,nonew=a->nonew,bs2=a->bs2,bs=A->rmap->bs,stepval;
648:   PetscBool         roworiented=a->roworiented;
649:   const PetscScalar *value     = v;
650:   MatScalar         *ap,*aa = a->a,*bap;

653:   if (roworiented) stepval = (n-1)*bs;
654:   else stepval = (m-1)*bs;

656:   for (k=0; k<m; k++) { /* loop over added rows */
657:     row = im[k];
658:     if (row < 0) continue;
659: #if defined(PETSC_USE_DEBUG)
660:     if (row >= a->mbs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Block index row too large %D max %D",row,a->mbs-1);
661: #endif
662:     rp   = aj + ai[row];
663:     ap   = aa + bs2*ai[row];
664:     rmax = imax[row];
665:     nrow = ailen[row];
666:     low  = 0;
667:     high = nrow;
668:     for (l=0; l<n; l++) { /* loop over added columns */
669:       if (in[l] < 0) continue;
670:       col = in[l];
671: #if defined(PETSC_USE_DEBUG)
672:       if (col >= a->nbs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Block index column too large %D max %D",col,a->nbs-1);
673: #endif
674:       if (col < row) {
675:         if (a->ignore_ltriangular) continue; /* ignore lower triangular block */
676:         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)");
677:       }
678:       if (roworiented) value = v + k*(stepval+bs)*bs + l*bs;
679:       else value = v + l*(stepval+bs)*bs + k*bs;

681:       if (col <= lastcol) low = 0;
682:       else high = nrow;

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

761: /*
762:     This is not yet used
763: */
766: PetscErrorCode MatAssemblyEnd_SeqSBAIJ_SeqAIJ_Inode(Mat A)
767: {
768:   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
770:   const PetscInt *ai = a->i, *aj = a->j,*cols;
771:   PetscInt       i   = 0,j,blk_size,m = A->rmap->n,node_count = 0,nzx,nzy,*ns,row,nz,cnt,cnt2,*counts;
772:   PetscBool      flag;

775:   PetscMalloc1(m,&ns);
776:   while (i < m) {
777:     nzx = ai[i+1] - ai[i];       /* Number of nonzeros */
778:     /* Limits the number of elements in a node to 'a->inode.limit' */
779:     for (j=i+1,blk_size=1; j<m && blk_size <a->inode.limit; ++j,++blk_size) {
780:       nzy = ai[j+1] - ai[j];
781:       if (nzy != (nzx - j + i)) break;
782:       PetscMemcmp(aj + ai[i] + j - i,aj + ai[j],nzy*sizeof(PetscInt),&flag);
783:       if (!flag) break;
784:     }
785:     ns[node_count++] = blk_size;

787:     i = j;
788:   }
789:   if (!a->inode.size && m && node_count > .9*m) {
790:     PetscFree(ns);
791:     PetscInfo2(A,"Found %D nodes out of %D rows. Not using Inode routines\n",node_count,m);
792:   } else {
793:     a->inode.node_count = node_count;

795:     PetscMalloc1(node_count,&a->inode.size);
796:     PetscLogObjectMemory((PetscObject)A,node_count*sizeof(PetscInt));
797:     PetscMemcpy(a->inode.size,ns,node_count*sizeof(PetscInt));
798:     PetscFree(ns);
799:     PetscInfo3(A,"Found %D nodes of %D. Limit used: %D. Using Inode routines\n",node_count,m,a->inode.limit);

801:     /* count collections of adjacent columns in each inode */
802:     row = 0;
803:     cnt = 0;
804:     for (i=0; i<node_count; i++) {
805:       cols = aj + ai[row] + a->inode.size[i];
806:       nz   = ai[row+1] - ai[row] - a->inode.size[i];
807:       for (j=1; j<nz; j++) {
808:         if (cols[j] != cols[j-1]+1) cnt++;
809:       }
810:       cnt++;
811:       row += a->inode.size[i];
812:     }
813:     PetscMalloc1(2*cnt,&counts);
814:     cnt  = 0;
815:     row  = 0;
816:     for (i=0; i<node_count; i++) {
817:       cols = aj + ai[row] + a->inode.size[i];
818:       counts[2*cnt] = cols[0];
819:       nz   = ai[row+1] - ai[row] - a->inode.size[i];
820:       cnt2 = 1;
821:       for (j=1; j<nz; j++) {
822:         if (cols[j] != cols[j-1]+1) {
823:           counts[2*(cnt++)+1] = cnt2;
824:           counts[2*cnt]       = cols[j];
825:           cnt2 = 1;
826:         } else cnt2++;
827:       }
828:       counts[2*(cnt++)+1] = cnt2;
829:       row += a->inode.size[i];
830:     }
831:     PetscIntView(2*cnt,counts,0);
832:   }
833:   return(0);
834: }

838: PetscErrorCode MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode)
839: {
840:   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
842:   PetscInt       fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
843:   PetscInt       m      = A->rmap->N,*ip,N,*ailen = a->ilen;
844:   PetscInt       mbs    = a->mbs,bs2 = a->bs2,rmax = 0;
845:   MatScalar      *aa    = a->a,*ap;

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

850:   if (m) rmax = ailen[0];
851:   for (i=1; i<mbs; i++) {
852:     /* move each row back by the amount of empty slots (fshift) before it*/
853:     fshift += imax[i-1] - ailen[i-1];
854:     rmax    = PetscMax(rmax,ailen[i]);
855:     if (fshift) {
856:       ip = aj + ai[i]; ap = aa + bs2*ai[i];
857:       N  = ailen[i];
858:       for (j=0; j<N; j++) {
859:         ip[j-fshift] = ip[j];
860:         PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));
861:       }
862:     }
863:     ai[i] = ai[i-1] + ailen[i-1];
864:   }
865:   if (mbs) {
866:     fshift += imax[mbs-1] - ailen[mbs-1];
867:     ai[mbs] = ai[mbs-1] + ailen[mbs-1];
868:   }
869:   /* reset ilen and imax for each row */
870:   for (i=0; i<mbs; i++) {
871:     ailen[i] = imax[i] = ai[i+1] - ai[i];
872:   }
873:   a->nz = ai[mbs];

875:   /* diagonals may have moved, reset it */
876:   if (a->diag) {
877:     PetscMemcpy(a->diag,ai,mbs*sizeof(PetscInt));
878:   }
879:   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);

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

885:   A->info.mallocs    += a->reallocs;
886:   a->reallocs         = 0;
887:   A->info.nz_unneeded = (PetscReal)fshift*bs2;
888:   a->idiagvalid       = PETSC_FALSE;
889:   a->rmax             = rmax;

891:   if (A->cmap->n < 65536 && A->cmap->bs == 1) {
892:     if (a->jshort && a->free_jshort) {
893:       /* when matrix data structure is changed, previous jshort must be replaced */
894:       PetscFree(a->jshort);
895:     }
896:     PetscMalloc1(a->i[A->rmap->n],&a->jshort);
897:     PetscLogObjectMemory((PetscObject)A,a->i[A->rmap->n]*sizeof(unsigned short));
898:     for (i=0; i<a->i[A->rmap->n]; i++) a->jshort[i] = a->j[i];
899:     A->ops->mult   = MatMult_SeqSBAIJ_1_ushort;
900:     A->ops->sor    = MatSOR_SeqSBAIJ_ushort;
901:     a->free_jshort = PETSC_TRUE;
902:   }
903:   return(0);
904: }

906: /*
907:    This function returns an array of flags which indicate the locations of contiguous
908:    blocks that should be zeroed. for eg: if bs = 3  and is = [0,1,2,3,5,6,7,8,9]
909:    then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)]
910:    Assume: sizes should be long enough to hold all the values.
911: */
914: PetscErrorCode MatZeroRows_SeqSBAIJ_Check_Blocks(PetscInt idx[],PetscInt n,PetscInt bs,PetscInt sizes[], PetscInt *bs_max)
915: {
916:   PetscInt  i,j,k,row;
917:   PetscBool flg;

920:   for (i=0,j=0; i<n; j++) {
921:     row = idx[i];
922:     if (row%bs!=0) { /* Not the begining of a block */
923:       sizes[j] = 1;
924:       i++;
925:     } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */
926:       sizes[j] = 1;         /* Also makes sure atleast 'bs' values exist for next else */
927:       i++;
928:     } else { /* Begining of the block, so check if the complete block exists */
929:       flg = PETSC_TRUE;
930:       for (k=1; k<bs; k++) {
931:         if (row+k != idx[i+k]) { /* break in the block */
932:           flg = PETSC_FALSE;
933:           break;
934:         }
935:       }
936:       if (flg) { /* No break in the bs */
937:         sizes[j] = bs;
938:         i       += bs;
939:       } else {
940:         sizes[j] = 1;
941:         i++;
942:       }
943:     }
944:   }
945:   *bs_max = j;
946:   return(0);
947: }


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

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

967:   for (k=0; k<m; k++) { /* loop over added rows */
968:     row  = im[k];       /* row number */
969:     brow = row/bs;      /* block row number */
970:     if (row < 0) continue;
971: #if defined(PETSC_USE_DEBUG)
972:     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);
973: #endif
974:     rp   = aj + ai[brow]; /*ptr to beginning of column value of the row block*/
975:     ap   = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/
976:     rmax = imax[brow];  /* maximum space allocated for this row */
977:     nrow = ailen[brow]; /* actual length of this row */
978:     low  = 0;

980:     for (l=0; l<n; l++) { /* loop over added columns */
981:       if (in[l] < 0) continue;
982: #if defined(PETSC_USE_DEBUG)
983:       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);
984: #endif
985:       col  = in[l];
986:       bcol = col/bs;              /* block col number */

988:       if (brow > bcol) {
989:         if (a->ignore_ltriangular) continue; /* ignore lower triangular values */
990:         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)");
991:       }

993:       ridx = row % bs; cidx = col % bs; /*row and col index inside the block */
994:       if ((brow==bcol && ridx<=cidx) || (brow<bcol)) {
995:         /* element value a(k,l) */
996:         if (roworiented) value = v[l + k*n];
997:         else value = v[k + l*m];

999:         /* move pointer bap to a(k,l) quickly and add/insert value */
1000:         if (col <= lastcol) low = 0;
1001:         high = nrow;
1002:         lastcol = col;
1003:         while (high-low > 7) {
1004:           t = (low+high)/2;
1005:           if (rp[t] > bcol) high = t;
1006:           else              low  = t;
1007:         }
1008:         for (i=low; i<high; i++) {
1009:           if (rp[i] > bcol) break;
1010:           if (rp[i] == bcol) {
1011:             bap = ap +  bs2*i + bs*cidx + ridx;
1012:             if (is == ADD_VALUES) *bap += value;
1013:             else                  *bap  = value;
1014:             /* for diag block, add/insert its symmetric element a(cidx,ridx) */
1015:             if (brow == bcol && ridx < cidx) {
1016:               bap = ap +  bs2*i + bs*ridx + cidx;
1017:               if (is == ADD_VALUES) *bap += value;
1018:               else                  *bap  = value;
1019:             }
1020:             goto noinsert1;
1021:           }
1022:         }

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

1028:         N = nrow++ - 1; high++;
1029:         /* shift up all the later entries in this row */
1030:         for (ii=N; ii>=i; ii--) {
1031:           rp[ii+1] = rp[ii];
1032:           PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));
1033:         }
1034:         if (N>=i) {
1035:           PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));
1036:         }
1037:         rp[i]                      = bcol;
1038:         ap[bs2*i + bs*cidx + ridx] = value;
1039:         A->nonzerostate++;
1040: noinsert1:;
1041:         low = i;
1042:       }
1043:     }   /* end of loop over added columns */
1044:     ailen[brow] = nrow;
1045:   }   /* end of loop over added rows */
1046:   return(0);
1047: }

1051: PetscErrorCode MatICCFactor_SeqSBAIJ(Mat inA,IS row,const MatFactorInfo *info)
1052: {
1053:   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)inA->data;
1054:   Mat            outA;
1056:   PetscBool      row_identity;

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

1064:   outA            = inA;
1065:   inA->factortype = MAT_FACTOR_ICC;
1066:   PetscFree(inA->solvertype);
1067:   PetscStrallocpy(MATSOLVERPETSC,&inA->solvertype);

1069:   MatMarkDiagonal_SeqSBAIJ(inA);
1070:   MatSeqSBAIJSetNumericFactorization_inplace(inA,row_identity);

1072:   PetscObjectReference((PetscObject)row);
1073:   ISDestroy(&a->row);
1074:   a->row = row;
1075:   PetscObjectReference((PetscObject)row);
1076:   ISDestroy(&a->col);
1077:   a->col = row;

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

1083:   if (!a->solve_work) {
1084:     PetscMalloc1(inA->rmap->N+inA->rmap->bs,&a->solve_work);
1085:     PetscLogObjectMemory((PetscObject)inA,(inA->rmap->N+inA->rmap->bs)*sizeof(PetscScalar));
1086:   }

1088:   MatCholeskyFactorNumeric(outA,inA,info);
1089:   return(0);
1090: }

1094: PetscErrorCode  MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,PetscInt *indices)
1095: {
1096:   Mat_SeqSBAIJ   *baij = (Mat_SeqSBAIJ*)mat->data;
1097:   PetscInt       i,nz,n;

1101:   nz = baij->maxnz;
1102:   n  = mat->cmap->n;
1103:   for (i=0; i<nz; i++) baij->j[i] = indices[i];

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

1108:   MatSetOption(mat,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
1109:   return(0);
1110: }

1114: /*@
1115:   MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows
1116:   in the matrix.

1118:   Input Parameters:
1119:   +  mat     - the SeqSBAIJ matrix
1120:   -  indices - the column indices

1122:   Level: advanced

1124:   Notes:
1125:   This can be called if you have precomputed the nonzero structure of the
1126:   matrix and want to provide it to the matrix object to improve the performance
1127:   of the MatSetValues() operation.

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

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

1134:   .seealso: MatCreateSeqSBAIJ
1135: @*/
1136: PetscErrorCode  MatSeqSBAIJSetColumnIndices(Mat mat,PetscInt *indices)
1137: {

1143:   PetscUseMethod(mat,"MatSeqSBAIJSetColumnIndices_C",(Mat,PetscInt*),(mat,indices));
1144:   return(0);
1145: }

1149: PetscErrorCode MatCopy_SeqSBAIJ(Mat A,Mat B,MatStructure str)
1150: {

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

1159:     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");
1160:     PetscMemcpy(b->a,a->a,(a->i[A->rmap->N])*sizeof(PetscScalar));
1161:   } else {
1162:     MatGetRowUpperTriangular(A);
1163:     MatCopy_Basic(A,B,str);
1164:     MatRestoreRowUpperTriangular(A);
1165:   }
1166:   return(0);
1167: }

1171: PetscErrorCode MatSetUp_SeqSBAIJ(Mat A)
1172: {

1176:   MatSeqSBAIJSetPreallocation_SeqSBAIJ(A,A->rmap->bs,PETSC_DEFAULT,0);
1177:   return(0);
1178: }

1182: PetscErrorCode MatSeqSBAIJGetArray_SeqSBAIJ(Mat A,PetscScalar *array[])
1183: {
1184:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;

1187:   *array = a->a;
1188:   return(0);
1189: }

1193: PetscErrorCode MatSeqSBAIJRestoreArray_SeqSBAIJ(Mat A,PetscScalar *array[])
1194: {
1196:   return(0);
1197: }

1201: PetscErrorCode MatAXPYGetPreallocation_SeqSBAIJ(Mat Y,Mat X,PetscInt *nnz)
1202: {
1203:   PetscInt       bs = Y->rmap->bs,mbs = Y->rmap->N/bs;
1204:   Mat_SeqSBAIJ   *x = (Mat_SeqSBAIJ*)X->data;
1205:   Mat_SeqSBAIJ   *y = (Mat_SeqSBAIJ*)Y->data;

1209:   /* Set the number of nonzeros in the new matrix */
1210:   MatAXPYGetPreallocation_SeqX_private(mbs,x->i,x->j,y->i,y->j,nnz);
1211:   return(0);
1212: }

1216: PetscErrorCode MatAXPY_SeqSBAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str)
1217: {
1218:   Mat_SeqSBAIJ   *x=(Mat_SeqSBAIJ*)X->data, *y=(Mat_SeqSBAIJ*)Y->data;
1220:   PetscInt       bs=Y->rmap->bs,bs2=bs*bs;
1221:   PetscBLASInt   one = 1;

1224:   if (str == SAME_NONZERO_PATTERN) {
1225:     PetscScalar  alpha = a;
1226:     PetscBLASInt bnz;
1227:     PetscBLASIntCast(x->nz*bs2,&bnz);
1228:     PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&bnz,&alpha,x->a,&one,y->a,&one));
1229:     PetscObjectStateIncrease((PetscObject)Y);
1230:   } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
1231:     MatSetOption(X,MAT_GETROW_UPPERTRIANGULAR,PETSC_TRUE);
1232:     MatAXPY_Basic(Y,a,X,str);
1233:     MatSetOption(X,MAT_GETROW_UPPERTRIANGULAR,PETSC_FALSE);
1234:   } else {
1235:     Mat      B;
1236:     PetscInt *nnz;
1237:     if (bs != X->rmap->bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrices must have same block size");
1238:     MatGetRowUpperTriangular(X);
1239:     MatGetRowUpperTriangular(Y);
1240:     PetscMalloc1(Y->rmap->N,&nnz);
1241:     MatCreate(PetscObjectComm((PetscObject)Y),&B);
1242:     PetscObjectSetName((PetscObject)B,((PetscObject)Y)->name);
1243:     MatSetSizes(B,Y->rmap->n,Y->cmap->n,Y->rmap->N,Y->cmap->N);
1244:     MatSetBlockSizesFromMats(B,Y,Y);
1245:     MatSetType(B,(MatType) ((PetscObject)Y)->type_name);
1246:     MatAXPYGetPreallocation_SeqSBAIJ(Y,X,nnz);
1247:     MatSeqSBAIJSetPreallocation(B,bs,0,nnz);

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

1251:     MatHeaderReplace(Y,&B);
1252:     PetscFree(nnz);
1253:     MatRestoreRowUpperTriangular(X);
1254:     MatRestoreRowUpperTriangular(Y);
1255:   }
1256:   return(0);
1257: }

1261: PetscErrorCode MatIsSymmetric_SeqSBAIJ(Mat A,PetscReal tol,PetscBool  *flg)
1262: {
1264:   *flg = PETSC_TRUE;
1265:   return(0);
1266: }

1270: PetscErrorCode MatIsStructurallySymmetric_SeqSBAIJ(Mat A,PetscBool  *flg)
1271: {
1273:   *flg = PETSC_TRUE;
1274:   return(0);
1275: }

1279: PetscErrorCode MatIsHermitian_SeqSBAIJ(Mat A,PetscReal tol,PetscBool  *flg)
1280: {
1282:   *flg = PETSC_FALSE;
1283:   return(0);
1284: }

1288: PetscErrorCode MatRealPart_SeqSBAIJ(Mat A)
1289: {
1290:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1291:   PetscInt     i,nz = a->bs2*a->i[a->mbs];
1292:   MatScalar    *aa = a->a;

1295:   for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]);
1296:   return(0);
1297: }

1301: PetscErrorCode MatImaginaryPart_SeqSBAIJ(Mat A)
1302: {
1303:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1304:   PetscInt     i,nz = a->bs2*a->i[a->mbs];
1305:   MatScalar    *aa = a->a;

1308:   for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]);
1309:   return(0);
1310: }

1314: PetscErrorCode MatZeroRowsColumns_SeqSBAIJ(Mat A,PetscInt is_n,const PetscInt is_idx[],PetscScalar diag,Vec x, Vec b)
1315: {
1316:   Mat_SeqSBAIJ      *baij=(Mat_SeqSBAIJ*)A->data;
1317:   PetscErrorCode    ierr;
1318:   PetscInt          i,j,k,count;
1319:   PetscInt          bs   =A->rmap->bs,bs2=baij->bs2,row,col;
1320:   PetscScalar       zero = 0.0;
1321:   MatScalar         *aa;
1322:   const PetscScalar *xx;
1323:   PetscScalar       *bb;
1324:   PetscBool         *zeroed,vecs = PETSC_FALSE;

1327:   /* fix right hand side if needed */
1328:   if (x && b) {
1329:     VecGetArrayRead(x,&xx);
1330:     VecGetArray(b,&bb);
1331:     vecs = PETSC_TRUE;
1332:   }

1334:   /* zero the columns */
1335:   PetscCalloc1(A->rmap->n,&zeroed);
1336:   for (i=0; i<is_n; i++) {
1337:     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]);
1338:     zeroed[is_idx[i]] = PETSC_TRUE;
1339:   }
1340:   if (vecs) {
1341:     for (i=0; i<A->rmap->N; i++) {
1342:       row = i/bs;
1343:       for (j=baij->i[row]; j<baij->i[row+1]; j++) {
1344:         for (k=0; k<bs; k++) {
1345:           col = bs*baij->j[j] + k;
1346:           if (col <= i) continue;
1347:           aa = ((MatScalar*)(baij->a)) + j*bs2 + (i%bs) + bs*k;
1348:           if (!zeroed[i] && zeroed[col]) bb[i]   -= aa[0]*xx[col];
1349:           if (zeroed[i] && !zeroed[col]) bb[col] -= aa[0]*xx[i];
1350:         }
1351:       }
1352:     }
1353:     for (i=0; i<is_n; i++) bb[is_idx[i]] = diag*xx[is_idx[i]];
1354:   }

1356:   for (i=0; i<A->rmap->N; i++) {
1357:     if (!zeroed[i]) {
1358:       row = i/bs;
1359:       for (j=baij->i[row]; j<baij->i[row+1]; j++) {
1360:         for (k=0; k<bs; k++) {
1361:           col = bs*baij->j[j] + k;
1362:           if (zeroed[col]) {
1363:             aa = ((MatScalar*)(baij->a)) + j*bs2 + (i%bs) + bs*k;
1364:             aa[0] = 0.0;
1365:           }
1366:         }
1367:       }
1368:     }
1369:   }
1370:   PetscFree(zeroed);
1371:   if (vecs) {
1372:     VecRestoreArrayRead(x,&xx);
1373:     VecRestoreArray(b,&bb);
1374:   }

1376:   /* zero the rows */
1377:   for (i=0; i<is_n; i++) {
1378:     row   = is_idx[i];
1379:     count = (baij->i[row/bs +1] - baij->i[row/bs])*bs;
1380:     aa    = ((MatScalar*)(baij->a)) + baij->i[row/bs]*bs2 + (row%bs);
1381:     for (k=0; k<count; k++) {
1382:       aa[0] =  zero;
1383:       aa   += bs;
1384:     }
1385:     if (diag != 0.0) {
1386:       (*A->ops->setvalues)(A,1,&row,1,&row,&diag,INSERT_VALUES);
1387:     }
1388:   }
1389:   MatAssemblyEnd_SeqSBAIJ(A,MAT_FINAL_ASSEMBLY);
1390:   return(0);
1391: }

1395: PetscErrorCode MatShift_SeqSBAIJ(Mat Y,PetscScalar a)
1396: {
1398:   Mat_SeqSBAIJ    *aij = (Mat_SeqSBAIJ*)Y->data;

1401:   if (!Y->preallocated || !aij->nz) {
1402:     MatSeqSBAIJSetPreallocation(Y,Y->rmap->bs,1,NULL);
1403:   }
1404:   MatShift_Basic(Y,a);
1405:   return(0);
1406: }

1408: /* -------------------------------------------------------------------*/
1409: static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ,
1410:                                        MatGetRow_SeqSBAIJ,
1411:                                        MatRestoreRow_SeqSBAIJ,
1412:                                        MatMult_SeqSBAIJ_N,
1413:                                /*  4*/ MatMultAdd_SeqSBAIJ_N,
1414:                                        MatMult_SeqSBAIJ_N,       /* transpose versions are same as non-transpose versions */
1415:                                        MatMultAdd_SeqSBAIJ_N,
1416:                                        0,
1417:                                        0,
1418:                                        0,
1419:                                /* 10*/ 0,
1420:                                        0,
1421:                                        MatCholeskyFactor_SeqSBAIJ,
1422:                                        MatSOR_SeqSBAIJ,
1423:                                        MatTranspose_SeqSBAIJ,
1424:                                /* 15*/ MatGetInfo_SeqSBAIJ,
1425:                                        MatEqual_SeqSBAIJ,
1426:                                        MatGetDiagonal_SeqSBAIJ,
1427:                                        MatDiagonalScale_SeqSBAIJ,
1428:                                        MatNorm_SeqSBAIJ,
1429:                                /* 20*/ 0,
1430:                                        MatAssemblyEnd_SeqSBAIJ,
1431:                                        MatSetOption_SeqSBAIJ,
1432:                                        MatZeroEntries_SeqSBAIJ,
1433:                                /* 24*/ 0,
1434:                                        0,
1435:                                        0,
1436:                                        0,
1437:                                        0,
1438:                                /* 29*/ MatSetUp_SeqSBAIJ,
1439:                                        0,
1440:                                        0,
1441:                                        0,
1442:                                        0,
1443:                                /* 34*/ MatDuplicate_SeqSBAIJ,
1444:                                        0,
1445:                                        0,
1446:                                        0,
1447:                                        MatICCFactor_SeqSBAIJ,
1448:                                /* 39*/ MatAXPY_SeqSBAIJ,
1449:                                        MatGetSubMatrices_SeqSBAIJ,
1450:                                        MatIncreaseOverlap_SeqSBAIJ,
1451:                                        MatGetValues_SeqSBAIJ,
1452:                                        MatCopy_SeqSBAIJ,
1453:                                /* 44*/ 0,
1454:                                        MatScale_SeqSBAIJ,
1455:                                        MatShift_SeqSBAIJ,
1456:                                        0,
1457:                                        MatZeroRowsColumns_SeqSBAIJ,
1458:                                /* 49*/ 0,
1459:                                        MatGetRowIJ_SeqSBAIJ,
1460:                                        MatRestoreRowIJ_SeqSBAIJ,
1461:                                        0,
1462:                                        0,
1463:                                /* 54*/ 0,
1464:                                        0,
1465:                                        0,
1466:                                        0,
1467:                                        MatSetValuesBlocked_SeqSBAIJ,
1468:                                /* 59*/ MatGetSubMatrix_SeqSBAIJ,
1469:                                        0,
1470:                                        0,
1471:                                        0,
1472:                                        0,
1473:                                /* 64*/ 0,
1474:                                        0,
1475:                                        0,
1476:                                        0,
1477:                                        0,
1478:                                /* 69*/ MatGetRowMaxAbs_SeqSBAIJ,
1479:                                        0,
1480:                                        0,
1481:                                        0,
1482:                                        0,
1483:                                /* 74*/ 0,
1484:                                        0,
1485:                                        0,
1486:                                        0,
1487:                                        0,
1488:                                /* 79*/ 0,
1489:                                        0,
1490:                                        0,
1491:                                        MatGetInertia_SeqSBAIJ,
1492:                                        MatLoad_SeqSBAIJ,
1493:                                /* 84*/ MatIsSymmetric_SeqSBAIJ,
1494:                                        MatIsHermitian_SeqSBAIJ,
1495:                                        MatIsStructurallySymmetric_SeqSBAIJ,
1496:                                        0,
1497:                                        0,
1498:                                /* 89*/ 0,
1499:                                        0,
1500:                                        0,
1501:                                        0,
1502:                                        0,
1503:                                /* 94*/ 0,
1504:                                        0,
1505:                                        0,
1506:                                        0,
1507:                                        0,
1508:                                /* 99*/ 0,
1509:                                        0,
1510:                                        0,
1511:                                        0,
1512:                                        0,
1513:                                /*104*/ 0,
1514:                                        MatRealPart_SeqSBAIJ,
1515:                                        MatImaginaryPart_SeqSBAIJ,
1516:                                        MatGetRowUpperTriangular_SeqSBAIJ,
1517:                                        MatRestoreRowUpperTriangular_SeqSBAIJ,
1518:                                /*109*/ 0,
1519:                                        0,
1520:                                        0,
1521:                                        0,
1522:                                        MatMissingDiagonal_SeqSBAIJ,
1523:                                /*114*/ 0,
1524:                                        0,
1525:                                        0,
1526:                                        0,
1527:                                        0,
1528:                                /*119*/ 0,
1529:                                        0,
1530:                                        0,
1531:                                        0,
1532:                                        0,
1533:                                /*124*/ 0,
1534:                                        0,
1535:                                        0,
1536:                                        0,
1537:                                        0,
1538:                                /*129*/ 0,
1539:                                        0,
1540:                                        0,
1541:                                        0,
1542:                                        0,
1543:                                /*134*/ 0,
1544:                                        0,
1545:                                        0,
1546:                                        0,
1547:                                        0,
1548:                                /*139*/ 0,
1549:                                        0,
1550:                                        0,
1551:                                        0,
1552:                                        0,
1553:                                 /*144*/MatCreateMPIMatConcatenateSeqMat_SeqSBAIJ
1554: };

1558: PetscErrorCode  MatStoreValues_SeqSBAIJ(Mat mat)
1559: {
1560:   Mat_SeqSBAIJ   *aij = (Mat_SeqSBAIJ*)mat->data;
1561:   PetscInt       nz   = aij->i[mat->rmap->N]*mat->rmap->bs*aij->bs2;

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

1567:   /* allocate space for values if not already there */
1568:   if (!aij->saved_values) {
1569:     PetscMalloc1(nz+1,&aij->saved_values);
1570:   }

1572:   /* copy values over */
1573:   PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));
1574:   return(0);
1575: }

1579: PetscErrorCode  MatRetrieveValues_SeqSBAIJ(Mat mat)
1580: {
1581:   Mat_SeqSBAIJ   *aij = (Mat_SeqSBAIJ*)mat->data;
1583:   PetscInt       nz = aij->i[mat->rmap->N]*mat->rmap->bs*aij->bs2;

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

1589:   /* copy values over */
1590:   PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));
1591:   return(0);
1592: }

1596: PetscErrorCode  MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat B,PetscInt bs,PetscInt nz,PetscInt *nnz)
1597: {
1598:   Mat_SeqSBAIJ   *b = (Mat_SeqSBAIJ*)B->data;
1600:   PetscInt       i,mbs,nbs,bs2;
1601:   PetscBool      skipallocation = PETSC_FALSE,flg = PETSC_FALSE,realalloc = PETSC_FALSE;

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

1607:   MatSetBlockSize(B,PetscAbs(bs));
1608:   PetscLayoutSetUp(B->rmap);
1609:   PetscLayoutSetUp(B->cmap);
1610:   PetscLayoutGetBlockSize(B->rmap,&bs);

1612:   mbs = B->rmap->N/bs;
1613:   nbs = B->cmap->n/bs;
1614:   bs2 = bs*bs;

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

1618:   if (nz == MAT_SKIP_ALLOCATION) {
1619:     skipallocation = PETSC_TRUE;
1620:     nz             = 0;
1621:   }

1623:   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3;
1624:   if (nz < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %D",nz);
1625:   if (nnz) {
1626:     for (i=0; i<mbs; i++) {
1627:       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]);
1628:       if (nnz[i] > nbs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be greater than block row length: local row %D value %D block rowlength %D",i,nnz[i],nbs);
1629:     }
1630:   }

1632:   B->ops->mult             = MatMult_SeqSBAIJ_N;
1633:   B->ops->multadd          = MatMultAdd_SeqSBAIJ_N;
1634:   B->ops->multtranspose    = MatMult_SeqSBAIJ_N;
1635:   B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_N;

1637:   PetscOptionsGetBool(((PetscObject)B)->options,((PetscObject)B)->prefix,"-mat_no_unroll",&flg,NULL);
1638:   if (!flg) {
1639:     switch (bs) {
1640:     case 1:
1641:       B->ops->mult             = MatMult_SeqSBAIJ_1;
1642:       B->ops->multadd          = MatMultAdd_SeqSBAIJ_1;
1643:       B->ops->multtranspose    = MatMult_SeqSBAIJ_1;
1644:       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_1;
1645:       break;
1646:     case 2:
1647:       B->ops->mult             = MatMult_SeqSBAIJ_2;
1648:       B->ops->multadd          = MatMultAdd_SeqSBAIJ_2;
1649:       B->ops->multtranspose    = MatMult_SeqSBAIJ_2;
1650:       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_2;
1651:       break;
1652:     case 3:
1653:       B->ops->mult             = MatMult_SeqSBAIJ_3;
1654:       B->ops->multadd          = MatMultAdd_SeqSBAIJ_3;
1655:       B->ops->multtranspose    = MatMult_SeqSBAIJ_3;
1656:       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_3;
1657:       break;
1658:     case 4:
1659:       B->ops->mult             = MatMult_SeqSBAIJ_4;
1660:       B->ops->multadd          = MatMultAdd_SeqSBAIJ_4;
1661:       B->ops->multtranspose    = MatMult_SeqSBAIJ_4;
1662:       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_4;
1663:       break;
1664:     case 5:
1665:       B->ops->mult             = MatMult_SeqSBAIJ_5;
1666:       B->ops->multadd          = MatMultAdd_SeqSBAIJ_5;
1667:       B->ops->multtranspose    = MatMult_SeqSBAIJ_5;
1668:       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_5;
1669:       break;
1670:     case 6:
1671:       B->ops->mult             = MatMult_SeqSBAIJ_6;
1672:       B->ops->multadd          = MatMultAdd_SeqSBAIJ_6;
1673:       B->ops->multtranspose    = MatMult_SeqSBAIJ_6;
1674:       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_6;
1675:       break;
1676:     case 7:
1677:       B->ops->mult             = MatMult_SeqSBAIJ_7;
1678:       B->ops->multadd          = MatMultAdd_SeqSBAIJ_7;
1679:       B->ops->multtranspose    = MatMult_SeqSBAIJ_7;
1680:       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_7;
1681:       break;
1682:     }
1683:   }

1685:   b->mbs = mbs;
1686:   b->nbs = nbs;
1687:   if (!skipallocation) {
1688:     if (!b->imax) {
1689:       PetscMalloc2(mbs,&b->imax,mbs,&b->ilen);

1691:       b->free_imax_ilen = PETSC_TRUE;

1693:       PetscLogObjectMemory((PetscObject)B,2*mbs*sizeof(PetscInt));
1694:     }
1695:     if (!nnz) {
1696:       if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
1697:       else if (nz <= 0) nz = 1;
1698:       for (i=0; i<mbs; i++) b->imax[i] = nz;
1699:       nz = nz*mbs; /* total nz */
1700:     } else {
1701:       nz = 0;
1702:       for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
1703:     }
1704:     /* b->ilen will count nonzeros in each block row so far. */
1705:     for (i=0; i<mbs; i++) b->ilen[i] = 0;
1706:     /* nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */

1708:     /* allocate the matrix space */
1709:     MatSeqXAIJFreeAIJ(B,&b->a,&b->j,&b->i);
1710:     PetscMalloc3(bs2*nz,&b->a,nz,&b->j,B->rmap->N+1,&b->i);
1711:     PetscLogObjectMemory((PetscObject)B,(B->rmap->N+1)*sizeof(PetscInt)+nz*(bs2*sizeof(PetscScalar)+sizeof(PetscInt)));
1712:     PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));
1713:     PetscMemzero(b->j,nz*sizeof(PetscInt));

1715:     b->singlemalloc = PETSC_TRUE;

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

1721:     b->free_a  = PETSC_TRUE;
1722:     b->free_ij = PETSC_TRUE;
1723:   } else {
1724:     b->free_a  = PETSC_FALSE;
1725:     b->free_ij = PETSC_FALSE;
1726:   }

1728:   B->rmap->bs = bs;
1729:   b->bs2      = bs2;
1730:   b->nz       = 0;
1731:   b->maxnz    = nz;

1733:   b->inew    = 0;
1734:   b->jnew    = 0;
1735:   b->anew    = 0;
1736:   b->a2anew  = 0;
1737:   b->permute = PETSC_FALSE;
1738:   if (realalloc) {MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);}
1739:   return(0);
1740: }

1744: PetscErrorCode MatSeqSBAIJSetPreallocationCSR_SeqSBAIJ(Mat B,PetscInt bs,const PetscInt ii[],const PetscInt jj[], const PetscScalar V[])
1745: {
1746:   PetscInt       i,j,m,nz,nz_max=0,*nnz;
1747:   PetscScalar    *values=0;
1748:   PetscBool      roworiented = ((Mat_SeqSBAIJ*)B->data)->roworiented;
1751:   if (bs < 1) SETERRQ1(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_OUTOFRANGE,"Invalid block size specified, must be positive but it is %D",bs);
1752:   PetscLayoutSetBlockSize(B->rmap,bs);
1753:   PetscLayoutSetBlockSize(B->cmap,bs);
1754:   PetscLayoutSetUp(B->rmap);
1755:   PetscLayoutSetUp(B->cmap);
1756:   PetscLayoutGetBlockSize(B->rmap,&bs);
1757:   m      = B->rmap->n/bs;

1759:   if (ii[0]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"ii[0] must be 0 but it is %D",ii[0]);
1760:   PetscMalloc1(m+1,&nnz);
1761:   for (i=0; i<m; i++) {
1762:     nz = ii[i+1] - ii[i];
1763:     if (nz < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row %D has a negative number of columns %D",i,nz);
1764:     nz_max = PetscMax(nz_max,nz);
1765:     nnz[i] = nz;
1766:   }
1767:   MatSeqSBAIJSetPreallocation(B,bs,0,nnz);
1768:   PetscFree(nnz);

1770:   values = (PetscScalar*)V;
1771:   if (!values) {
1772:     PetscCalloc1(bs*bs*nz_max,&values);
1773:   }
1774:   for (i=0; i<m; i++) {
1775:     PetscInt          ncols  = ii[i+1] - ii[i];
1776:     const PetscInt    *icols = jj + ii[i];
1777:     if (!roworiented || bs == 1) {
1778:       const PetscScalar *svals = values + (V ? (bs*bs*ii[i]) : 0);
1779:       MatSetValuesBlocked_SeqSBAIJ(B,1,&i,ncols,icols,svals,INSERT_VALUES);
1780:     } else {
1781:       for (j=0; j<ncols; j++) {
1782:         const PetscScalar *svals = values + (V ? (bs*bs*(ii[i]+j)) : 0);
1783:         MatSetValuesBlocked_SeqSBAIJ(B,1,&i,1,&icols[j],svals,INSERT_VALUES);
1784:       }
1785:     }
1786:   }
1787:   if (!V) { PetscFree(values); }
1788:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
1789:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
1790:   MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
1791:   return(0);
1792: }

1794: /*
1795:    This is used to set the numeric factorization for both Cholesky and ICC symbolic factorization
1796: */
1799: PetscErrorCode MatSeqSBAIJSetNumericFactorization_inplace(Mat B,PetscBool natural)
1800: {
1802:   PetscBool      flg = PETSC_FALSE;
1803:   PetscInt       bs  = B->rmap->bs;

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

1809:   if (!natural) {
1810:     switch (bs) {
1811:     case 1:
1812:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_inplace;
1813:       break;
1814:     case 2:
1815:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2;
1816:       break;
1817:     case 3:
1818:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3;
1819:       break;
1820:     case 4:
1821:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4;
1822:       break;
1823:     case 5:
1824:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5;
1825:       break;
1826:     case 6:
1827:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6;
1828:       break;
1829:     case 7:
1830:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7;
1831:       break;
1832:     default:
1833:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N;
1834:       break;
1835:     }
1836:   } else {
1837:     switch (bs) {
1838:     case 1:
1839:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering_inplace;
1840:       break;
1841:     case 2:
1842:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering;
1843:       break;
1844:     case 3:
1845:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering;
1846:       break;
1847:     case 4:
1848:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering;
1849:       break;
1850:     case 5:
1851:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering;
1852:       break;
1853:     case 6:
1854:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering;
1855:       break;
1856:     case 7:
1857:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering;
1858:       break;
1859:     default:
1860:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering;
1861:       break;
1862:     }
1863:   }
1864:   return(0);
1865: }

1867: PETSC_INTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqAIJ(Mat, MatType,MatReuse,Mat*);
1868: PETSC_INTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqBAIJ(Mat, MatType,MatReuse,Mat*);

1872: PETSC_INTERN PetscErrorCode MatGetFactor_seqsbaij_petsc(Mat A,MatFactorType ftype,Mat *B)
1873: {
1874:   PetscInt       n = A->rmap->n;

1878: #if defined(PETSC_USE_COMPLEX)
1879:   if (A->hermitian) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Hermitian Factor is not supported");
1880: #endif
1881:   MatCreate(PetscObjectComm((PetscObject)A),B);
1882:   MatSetSizes(*B,n,n,n,n);
1883:   if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) {
1884:     MatSetType(*B,MATSEQSBAIJ);
1885:     MatSeqSBAIJSetPreallocation(*B,A->rmap->bs,MAT_SKIP_ALLOCATION,NULL);

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

1891:   (*B)->factortype = ftype;
1892:   PetscFree((*B)->solvertype);
1893:   PetscStrallocpy(MATSOLVERPETSC,&(*B)->solvertype);
1894:   return(0);
1895: }

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

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

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

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


1912:   Level: beginner

1914:   .seealso: MatCreateSeqSBAIJ
1915: M*/

1917: PETSC_INTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqSBSTRM(Mat, MatType,MatReuse,Mat*);

1921: PETSC_EXTERN PetscErrorCode MatCreate_SeqSBAIJ(Mat B)
1922: {
1923:   Mat_SeqSBAIJ   *b;
1925:   PetscMPIInt    size;
1926:   PetscBool      no_unroll = PETSC_FALSE,no_inode = PETSC_FALSE;

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

1932:   PetscNewLog(B,&b);
1933:   B->data = (void*)b;
1934:   PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));

1936:   B->ops->destroy    = MatDestroy_SeqSBAIJ;
1937:   B->ops->view       = MatView_SeqSBAIJ;
1938:   b->row             = 0;
1939:   b->icol            = 0;
1940:   b->reallocs        = 0;
1941:   b->saved_values    = 0;
1942:   b->inode.limit     = 5;
1943:   b->inode.max_limit = 5;

1945:   b->roworiented        = PETSC_TRUE;
1946:   b->nonew              = 0;
1947:   b->diag               = 0;
1948:   b->solve_work         = 0;
1949:   b->mult_work          = 0;
1950:   B->spptr              = 0;
1951:   B->info.nz_unneeded   = (PetscReal)b->maxnz*b->bs2;
1952:   b->keepnonzeropattern = PETSC_FALSE;

1954:   b->inew    = 0;
1955:   b->jnew    = 0;
1956:   b->anew    = 0;
1957:   b->a2anew  = 0;
1958:   b->permute = PETSC_FALSE;

1960:   b->ignore_ltriangular = PETSC_TRUE;

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

1964:   b->getrow_utriangular = PETSC_FALSE;

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

1968:   PetscObjectComposeFunction((PetscObject)B,"MatStoreValues_C",MatStoreValues_SeqSBAIJ);
1969:   PetscObjectComposeFunction((PetscObject)B,"MatRetrieveValues_C",MatRetrieveValues_SeqSBAIJ);
1970:   PetscObjectComposeFunction((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C",MatSeqSBAIJSetColumnIndices_SeqSBAIJ);
1971:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_seqaij_C",MatConvert_SeqSBAIJ_SeqAIJ);
1972:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_seqbaij_C",MatConvert_SeqSBAIJ_SeqBAIJ);
1973:   PetscObjectComposeFunction((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",MatSeqSBAIJSetPreallocation_SeqSBAIJ);
1974:   PetscObjectComposeFunction((PetscObject)B,"MatSeqSBAIJSetPreallocationCSR_C",MatSeqSBAIJSetPreallocationCSR_SeqSBAIJ);
1975:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_seqsbstrm_C",MatConvert_SeqSBAIJ_SeqSBSTRM);
1976: #if defined(PETSC_HAVE_ELEMENTAL)
1977:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_elemental_C",MatConvert_SeqSBAIJ_Elemental);
1978: #endif

1980:   B->symmetric                  = PETSC_TRUE;
1981:   B->structurally_symmetric     = PETSC_TRUE;
1982:   B->symmetric_set              = PETSC_TRUE;
1983:   B->structurally_symmetric_set = PETSC_TRUE;

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

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

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

2012:    Collective on Mat

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

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

2027:    Level: intermediate

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

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

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


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

2052:   PetscTryMethod(B,"MatSeqSBAIJSetPreallocation_C",(Mat,PetscInt,PetscInt,const PetscInt[]),(B,bs,nz,nnz));
2053:   return(0);
2054: }

2056: #undef  __FUNCT__
2058: /*@C
2059:    MatSeqSBAIJSetPreallocationCSR - Allocates memory for a sparse sequential matrix in symmetric block AIJ format.

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

2068:    Level: developer

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

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

2079: .seealso: MatCreate(), MatCreateSeqSBAIJ(), MatSetValuesBlocked(), MatSeqSBAIJSetPreallocation(), MATSEQSBAIJ
2080: @*/
2081: PetscErrorCode MatSeqSBAIJSetPreallocationCSR(Mat B,PetscInt bs,const PetscInt i[],const PetscInt j[], const PetscScalar v[])
2082: {

2089:   PetscTryMethod(B,"MatSeqSBAIJSetPreallocationCSR_C",(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]),(B,bs,i,j,v));
2090:   return(0);
2091: }

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

2102:    Collective on MPI_Comm

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

2113:    Output Parameter:
2114: .  A - the symmetric matrix

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

2121:    Level: intermediate

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

2265: PetscErrorCode MatLoad_SeqSBAIJ(Mat newmat,PetscViewer viewer)
2266: {
2267:   Mat_SeqSBAIJ   *a;
2269:   int            fd;
2270:   PetscMPIInt    size;
2271:   PetscInt       i,nz,header[4],*rowlengths=0,M,N,bs = newmat->rmap->bs;
2272:   PetscInt       *mask,mbs,*jj,j,rowcount,nzcount,k,*s_browlengths,maskcount;
2273:   PetscInt       kmax,jcount,block,idx,point,nzcountb,extra_rows,rows,cols;
2274:   PetscInt       *masked,nmask,tmp,bs2,ishift;
2275:   PetscScalar    *aa;
2276:   MPI_Comm       comm;

2279:   /* force binary viewer to load .info file if it has not yet done so */
2280:   PetscViewerSetUp(viewer);
2281:   PetscObjectGetComm((PetscObject)viewer,&comm);
2282:   PetscOptionsGetInt(((PetscObject)newmat)->options,((PetscObject)newmat)->prefix,"-matload_block_size",&bs,NULL);
2283:   if (bs < 0) bs = 1;
2284:   bs2  = bs*bs;

2286:   MPI_Comm_size(comm,&size);
2287:   if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"view must have one processor");
2288:   PetscViewerBinaryGetDescriptor(viewer,&fd);
2289:   PetscBinaryRead(fd,header,4,PETSC_INT);
2290:   if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"not Mat object");
2291:   M = header[1]; N = header[2]; nz = header[3];

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

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

2297:   /*
2298:      This code adds extra rows to make sure the number of rows is
2299:     divisible by the blocksize
2300:   */
2301:   mbs        = M/bs;
2302:   extra_rows = bs - M + bs*(mbs);
2303:   if (extra_rows == bs) extra_rows = 0;
2304:   else                  mbs++;
2305:   if (extra_rows) {
2306:     PetscInfo(viewer,"Padding loaded matrix to match blocksize\n");
2307:   }

2309:   /* Set global sizes if not already set */
2310:   if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) {
2311:     MatSetSizes(newmat,PETSC_DECIDE,PETSC_DECIDE,M+extra_rows,N+extra_rows);
2312:   } else { /* Check if the matrix global sizes are correct */
2313:     MatGetSize(newmat,&rows,&cols);
2314:     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);
2315:   }

2317:   /* read in row lengths */
2318:   PetscMalloc1(M+extra_rows,&rowlengths);
2319:   PetscBinaryRead(fd,rowlengths,M,PETSC_INT);
2320:   for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;

2322:   /* read in column indices */
2323:   PetscMalloc1(nz+extra_rows,&jj);
2324:   PetscBinaryRead(fd,jj,nz,PETSC_INT);
2325:   for (i=0; i<extra_rows; i++) jj[nz+i] = M+i;

2327:   /* loop over row lengths determining block row lengths */
2328:   PetscCalloc1(mbs,&s_browlengths);
2329:   PetscMalloc2(mbs,&mask,mbs,&masked);
2330:   PetscMemzero(mask,mbs*sizeof(PetscInt));
2331:   rowcount = 0;
2332:   nzcount  = 0;
2333:   for (i=0; i<mbs; i++) {
2334:     nmask = 0;
2335:     for (j=0; j<bs; j++) {
2336:       kmax = rowlengths[rowcount];
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:       rowcount++;
2342:     }
2343:     s_browlengths[i] += nmask;

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

2349:   /* Do preallocation */
2350:   MatSeqSBAIJSetPreallocation_SeqSBAIJ(newmat,bs,0,s_browlengths);
2351:   a    = (Mat_SeqSBAIJ*)newmat->data;

2353:   /* set matrix "i" values */
2354:   a->i[0] = 0;
2355:   for (i=1; i<= mbs; i++) {
2356:     a->i[i]      = a->i[i-1] + s_browlengths[i-1];
2357:     a->ilen[i-1] = s_browlengths[i-1];
2358:   }
2359:   a->nz = a->i[mbs];

2361:   /* read in nonzero values */
2362:   PetscMalloc1(nz+extra_rows,&aa);
2363:   PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);
2364:   for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0;

2366:   /* set "a" and "j" values into matrix */
2367:   nzcount = 0; jcount = 0;
2368:   for (i=0; i<mbs; i++) {
2369:     nzcountb = nzcount;
2370:     nmask    = 0;
2371:     for (j=0; j<bs; j++) {
2372:       kmax = rowlengths[i*bs+j];
2373:       for (k=0; k<kmax; k++) {
2374:         tmp = jj[nzcount++]/bs; /* block col. index */
2375:         if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;}
2376:       }
2377:     }
2378:     /* sort the masked values */
2379:     PetscSortInt(nmask,masked);

2381:     /* set "j" values into matrix */
2382:     maskcount = 1;
2383:     for (j=0; j<nmask; j++) {
2384:       a->j[jcount++]  = masked[j];
2385:       mask[masked[j]] = maskcount++;
2386:     }

2388:     /* set "a" values into matrix */
2389:     ishift = bs2*a->i[i];
2390:     for (j=0; j<bs; j++) {
2391:       kmax = rowlengths[i*bs+j];
2392:       for (k=0; k<kmax; k++) {
2393:         tmp = jj[nzcountb]/bs;        /* block col. index */
2394:         if (tmp >= i) {
2395:           block     = mask[tmp] - 1;
2396:           point     = jj[nzcountb] - bs*tmp;
2397:           idx       = ishift + bs2*block + j + bs*point;
2398:           a->a[idx] = aa[nzcountb];
2399:         }
2400:         nzcountb++;
2401:       }
2402:     }
2403:     /* zero out the mask elements we set */
2404:     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
2405:   }
2406:   if (jcount != a->nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix");

2408:   PetscFree(rowlengths);
2409:   PetscFree(s_browlengths);
2410:   PetscFree(aa);
2411:   PetscFree(jj);
2412:   PetscFree2(mask,masked);

2414:   MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);
2415:   MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);
2416:   return(0);
2417: }

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

2425:      Collective on MPI_Comm

2427:    Input Parameters:
2428: +  comm - must be an MPI communicator of size 1
2429: .  bs - size of block
2430: .  m - number of rows
2431: .  n - number of columns
2432: .  i - row indices
2433: .  j - column indices
2434: -  a - matrix values

2436:    Output Parameter:
2437: .  mat - the matrix

2439:    Level: advanced

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

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

2447:        The i and j indices are 0 based

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

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

2454: @*/
2455: PetscErrorCode  MatCreateSeqSBAIJWithArrays(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt *i,PetscInt *j,PetscScalar *a,Mat *mat)
2456: {
2458:   PetscInt       ii;
2459:   Mat_SeqSBAIJ   *sbaij;

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

2465:   MatCreate(comm,mat);
2466:   MatSetSizes(*mat,m,n,m,n);
2467:   MatSetType(*mat,MATSEQSBAIJ);
2468:   MatSeqSBAIJSetPreallocation_SeqSBAIJ(*mat,bs,MAT_SKIP_ALLOCATION,0);
2469:   sbaij = (Mat_SeqSBAIJ*)(*mat)->data;
2470:   PetscMalloc2(m,&sbaij->imax,m,&sbaij->ilen);
2471:   PetscLogObjectMemory((PetscObject)*mat,2*m*sizeof(PetscInt));

2473:   sbaij->i = i;
2474:   sbaij->j = j;
2475:   sbaij->a = a;

2477:   sbaij->singlemalloc = PETSC_FALSE;
2478:   sbaij->nonew        = -1;             /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/
2479:   sbaij->free_a       = PETSC_FALSE;
2480:   sbaij->free_ij      = PETSC_FALSE;

2482:   for (ii=0; ii<m; ii++) {
2483:     sbaij->ilen[ii] = sbaij->imax[ii] = i[ii+1] - i[ii];
2484: #if defined(PETSC_USE_DEBUG)
2485:     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]);
2486: #endif
2487:   }
2488: #if defined(PETSC_USE_DEBUG)
2489:   for (ii=0; ii<sbaij->i[m]; ii++) {
2490:     if (j[ii] < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column index at location = %d index = %d",ii,j[ii]);
2491:     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]);
2492:   }
2493: #endif

2495:   MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
2496:   MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
2497:   return(0);
2498: }

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

2507:   MatCreateMPIMatConcatenateSeqMat_MPISBAIJ(comm,inmat,n,scall,outmat);
2508:   return(0);
2509: }