Actual source code: sbaij.c

  1: /*$Id: sbaij.c,v 1.56 2001/04/10 22:35:39 balay Exp $*/

  3: /*
  4:     Defines the basic matrix operations for the BAIJ (compressed row)
  5:   matrix storage format.
  6: */
  7: #include "petscsys.h"                            /*I "petscmat.h" I*/
  8: #include "src/mat/impls/baij/seq/baij.h"
  9: #include "src/vec/vecimpl.h"
 10: #include "src/inline/spops.h"
 11: #include "src/mat/impls/sbaij/seq/sbaij.h"

 13: #define CHUNKSIZE  10

 15: /*
 16:      Checks for missing diagonals
 17: */
 18: int MatMissingDiagonal_SeqSBAIJ(Mat A)
 19: {
 20:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
 21:   int         *diag,*jj = a->j,i,ierr;

 24:   MatMarkDiagonal_SeqSBAIJ(A);
 25:   diag = a->diag;
 26:   for (i=0; i<a->mbs; i++) {
 27:     if (jj[diag[i]] != i) {
 28:       SETERRQ1(1,"Matrix is missing diagonal number %d",i);
 29:     }
 30:   }
 31:   return(0);
 32: }

 34: int MatMarkDiagonal_SeqSBAIJ(Mat A)
 35: {
 36:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
 37:   int          i,j,*diag,m = a->mbs,ierr;

 40:   if (a->diag) return(0);

 42:   PetscMalloc((m+1)*sizeof(int),&diag);
 43:   PetscLogObjectMemory(A,(m+1)*sizeof(int));
 44:   for (i=0; i<m; i++) {
 45:     diag[i] = a->i[i+1];
 46:     for (j=a->i[i]; j<a->i[i+1]; j++) {
 47:       if (a->j[j] == i) {
 48:         diag[i] = j;
 49:         break;
 50:       }
 51:     }
 52:   }
 53:   a->diag = diag;
 54:   return(0);
 55: }


 58: extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**);

 60: static int MatGetRowIJ_SeqSBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done)
 61: {
 62:   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;

 65:   if (ia) {
 66:     SETERRQ(1,"Function not yet written for SBAIJ format, only supports natural ordering");
 67:   }
 68:   *nn = a->mbs;
 69:   return(0);
 70: }

 72: static int MatRestoreRowIJ_SeqSBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done)
 73: {
 75:   if (!ia) return(0);
 76:   SETERRQ(1,"Function not yet written for SBAIJ format, only supports natural ordering");
 77:   /* return(0); */
 78: }

 80: int MatGetBlockSize_SeqSBAIJ(Mat mat,int *bs)
 81: {
 82:   Mat_SeqSBAIJ *sbaij = (Mat_SeqSBAIJ*)mat->data;

 85:   *bs = sbaij->bs;
 86:   return(0);
 87: }

 89: int MatDestroy_SeqSBAIJ(Mat A)
 90: {
 91:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
 92:   int         ierr;

 95: #if defined(PETSC_USE_LOG)
 96:   PetscLogObjectState((PetscObject)A,"Rows=%d, s_NZ=%d",A->m,a->s_nz);
 97: #endif
 98:   PetscFree(a->a);
 99:   if (!a->singlemalloc) {
100:     PetscFree(a->i);
101:     PetscFree(a->j);
102:   }
103:   if (a->row) {
104:     ISDestroy(a->row);
105:   }
106:   if (a->diag) {PetscFree(a->diag);}
107:   if (a->ilen) {PetscFree(a->ilen);}
108:   if (a->imax) {PetscFree(a->imax);}
109:   if (a->solve_work) {PetscFree(a->solve_work);}
110:   if (a->mult_work) {PetscFree(a->mult_work);}
111:   if (a->icol) {ISDestroy(a->icol);}
112:   if (a->saved_values) {PetscFree(a->saved_values);}

114:   if (a->inew){
115:     PetscFree(a->inew);
116:     a->inew = 0;
117:   }
118:   PetscFree(a);
119:   return(0);
120: }

122: int MatSetOption_SeqSBAIJ(Mat A,MatOption op)
123: {
124:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;

127:   if      (op == MAT_ROW_ORIENTED)                 a->roworiented    = PETSC_TRUE;
128:   else if (op == MAT_COLUMN_ORIENTED)              a->roworiented    = PETSC_FALSE;
129:   else if (op == MAT_COLUMNS_SORTED)               a->sorted         = PETSC_TRUE;
130:   else if (op == MAT_COLUMNS_UNSORTED)             a->sorted         = PETSC_FALSE;
131:   else if (op == MAT_KEEP_ZEROED_ROWS)             a->keepzeroedrows = PETSC_TRUE;
132:   else if (op == MAT_NO_NEW_NONZERO_LOCATIONS)     a->nonew          = 1;
133:   else if (op == MAT_NEW_NONZERO_LOCATION_ERR)     a->nonew          = -1;
134:   else if (op == MAT_NEW_NONZERO_ALLOCATION_ERR)   a->nonew          = -2;
135:   else if (op == MAT_YES_NEW_NONZERO_LOCATIONS)    a->nonew          = 0;
136:   else if (op == MAT_ROWS_SORTED ||
137:            op == MAT_ROWS_UNSORTED ||
138:            op == MAT_YES_NEW_DIAGONALS ||
139:            op == MAT_IGNORE_OFF_PROC_ENTRIES ||
140:            op == MAT_USE_HASH_TABLE) {
141:     PetscLogInfo(A,"MatSetOption_SeqSBAIJ:Option ignoredn");
142:   } else if (op == MAT_NO_NEW_DIAGONALS) {
143:     SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
144:   } else {
145:     SETERRQ(PETSC_ERR_SUP,"unknown option");
146:   }
147:   return(0);
148: }

150: int MatGetOwnershipRange_SeqSBAIJ(Mat A,int *m,int *n)
151: {
153:   *m = 0;
154:   *n = A->m;
155:   return(0);
156: }

158: int MatGetRow_SeqSBAIJ(Mat A,int row,int *ncols,int **cols,Scalar **v)
159: {
160:   Mat_SeqSBAIJ  *a = (Mat_SeqSBAIJ*)A->data;
161:   int          itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*cols_i,bs2,ierr;
162:   MatScalar    *aa,*aa_i;
163:   Scalar       *v_i;

166:   bs  = a->bs;
167:   ai  = a->i;
168:   aj  = a->j;
169:   aa  = a->a;
170:   bs2 = a->bs2;
171: 
172:   if (row < 0 || row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row out of range");
173: 
174:   bn  = row/bs;   /* Block number */
175:   bp  = row % bs; /* Block position */
176:   M   = ai[bn+1] - ai[bn];
177:   *ncols = bs*M;
178: 
179:   if (v) {
180:     *v = 0;
181:     if (*ncols) {
182:       PetscMalloc((*ncols+row)*sizeof(Scalar),v);
183:       for (i=0; i<M; i++) { /* for each block in the block row */
184:         v_i  = *v + i*bs;
185:         aa_i = aa + bs2*(ai[bn] + i);
186:         for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];}
187:       }
188:     }
189:   }

191:   if (cols) {
192:     *cols = 0;
193:     if (*ncols) {
194:       PetscMalloc((*ncols+row)*sizeof(int),cols);
195:       for (i=0; i<M; i++) { /* for each block in the block row */
196:         cols_i = *cols + i*bs;
197:         itmp  = bs*aj[ai[bn] + i];
198:         for (j=0; j<bs; j++) {cols_i[j] = itmp++;}
199:       }
200:     }
201:   }

203:   /*search column A(0:row-1,row) (=A(row,0:row-1)). Could be expensive! */
204:   /* this segment is currently removed, so only entries in the upper triangle are obtained */
205: #ifdef column_search
206:   v_i    = *v    + M*bs;
207:   cols_i = *cols + M*bs;
208:   for (i=0; i<bn; i++){ /* for each block row */
209:     M = ai[i+1] - ai[i];
210:     for (j=0; j<M; j++){
211:       itmp = aj[ai[i] + j];    /* block column value */
212:       if (itmp == bn){
213:         aa_i   = aa    + bs2*(ai[i] + j) + bs*bp;
214:         for (k=0; k<bs; k++) {
215:           *cols_i++ = i*bs+k;
216:           *v_i++    = aa_i[k];
217:         }
218:         *ncols += bs;
219:         break;
220:       }
221:     }
222:   }
223: #endif

225:   return(0);
226: }

228: int MatRestoreRow_SeqSBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
229: {

233:   if (idx) {if (*idx) {PetscFree(*idx);}}
234:   if (v)   {if (*v)   {PetscFree(*v);}}
235:   return(0);
236: }

238: int MatTranspose_SeqSBAIJ(Mat A,Mat *B)
239: {
241:   SETERRQ(1,"Matrix is symmetric. MatTranspose() should not be called");
242:   /* return(0); */
243: }

245: static int MatView_SeqSBAIJ_Binary(Mat A,PetscViewer viewer)
246: {
247:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
248:   int          i,fd,*col_lens,ierr,bs = a->bs,count,*jj,j,k,l,bs2=a->bs2;
249:   Scalar       *aa;
250:   FILE         *file;

253:   ierr        = PetscViewerBinaryGetDescriptor(viewer,&fd);
254:   ierr        = PetscMalloc((4+A->m)*sizeof(int),&col_lens);
255:   col_lens[0] = MAT_COOKIE;

257:   col_lens[1] = A->m;
258:   col_lens[2] = A->m;
259:   col_lens[3] = a->s_nz*bs2;

261:   /* store lengths of each row and write (including header) to file */
262:   count = 0;
263:   for (i=0; i<a->mbs; i++) {
264:     for (j=0; j<bs; j++) {
265:       col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]);
266:     }
267:   }
268:   PetscBinaryWrite(fd,col_lens,4+A->m,PETSC_INT,1);
269:   PetscFree(col_lens);

271:   /* store column indices (zero start index) */
272:   ierr  = PetscMalloc((a->s_nz+1)*bs2*sizeof(int),&jj);
273:   count = 0;
274:   for (i=0; i<a->mbs; i++) {
275:     for (j=0; j<bs; j++) {
276:       for (k=a->i[i]; k<a->i[i+1]; k++) {
277:         for (l=0; l<bs; l++) {
278:           jj[count++] = bs*a->j[k] + l;
279:         }
280:       }
281:     }
282:   }
283:   PetscBinaryWrite(fd,jj,bs2*a->s_nz,PETSC_INT,0);
284:   PetscFree(jj);

286:   /* store nonzero values */
287:   ierr  = PetscMalloc((a->s_nz+1)*bs2*sizeof(Scalar),&aa);
288:   count = 0;
289:   for (i=0; i<a->mbs; i++) {
290:     for (j=0; j<bs; j++) {
291:       for (k=a->i[i]; k<a->i[i+1]; k++) {
292:         for (l=0; l<bs; l++) {
293:           aa[count++] = a->a[bs2*k + l*bs + j];
294:         }
295:       }
296:     }
297:   }
298:   PetscBinaryWrite(fd,aa,bs2*a->s_nz,PETSC_SCALAR,0);
299:   PetscFree(aa);

301:   PetscViewerBinaryGetInfoPointer(viewer,&file);
302:   if (file) {
303:     fprintf(file,"-matload_block_size %dn",a->bs);
304:   }
305:   return(0);
306: }

308: static int MatView_SeqSBAIJ_ASCII(Mat A,PetscViewer viewer)
309: {
310:   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
311:   int               ierr,i,j,bs = a->bs,k,l,bs2=a->bs2;
312:   char              *name;
313:   PetscViewerFormat format;

316:   PetscObjectGetName((PetscObject)viewer,&name);
317:   PetscViewerGetFormat(viewer,&format);
318:   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_LONG) {
319:     PetscViewerASCIIPrintf(viewer,"  block size is %dn",bs);
320:   } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
321:     SETERRQ(PETSC_ERR_SUP,"Matlab format not supported");
322:   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
323:     PetscViewerASCIIUseTabs(viewer,PETSC_NO);
324:     for (i=0; i<a->mbs; i++) {
325:       for (j=0; j<bs; j++) {
326:         PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);
327:         for (k=a->i[i]; k<a->i[i+1]; k++) {
328:           for (l=0; l<bs; l++) {
329: #if defined(PETSC_USE_COMPLEX)
330:             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
331:               PetscViewerASCIIPrintf(viewer," %d %g + %g i",bs*a->j[k]+l,
332:                       PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));
333:             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
334:               PetscViewerASCIIPrintf(viewer," %d %g - %g i",bs*a->j[k]+l,
335:                       PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));
336:             } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
337:               PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));
338:             }
339: #else
340:             if (a->a[bs2*k + l*bs + j] != 0.0) {
341:               PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);
342:             }
343: #endif
344:           }
345:         }
346:         PetscViewerASCIIPrintf(viewer,"n");
347:       }
348:     }
349:     PetscViewerASCIIUseTabs(viewer,PETSC_YES);
350:   } else {
351:     PetscViewerASCIIUseTabs(viewer,PETSC_NO);
352:     for (i=0; i<a->mbs; i++) {
353:       for (j=0; j<bs; j++) {
354:         PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);
355:         for (k=a->i[i]; k<a->i[i+1]; k++) {
356:           for (l=0; l<bs; l++) {
357: #if defined(PETSC_USE_COMPLEX)
358:             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) {
359:               PetscViewerASCIIPrintf(viewer," %d %g + %g i",bs*a->j[k]+l,
360:                 PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));
361:             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) {
362:               PetscViewerASCIIPrintf(viewer," %d %g - %g i",bs*a->j[k]+l,
363:                 PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));
364:             } else {
365:               PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));
366:             }
367: #else
368:             PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);
369: #endif
370:           }
371:         }
372:         PetscViewerASCIIPrintf(viewer,"n");
373:       }
374:     }
375:     PetscViewerASCIIUseTabs(viewer,PETSC_YES);
376:   }
377:   PetscViewerFlush(viewer);
378:   return(0);
379: }

381: static int MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw,void *Aa)
382: {
383:   Mat          A = (Mat) Aa;
384:   Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data;
385:   int          row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2,rank;
386:   PetscReal    xl,yl,xr,yr,x_l,x_r,y_l,y_r;
387:   MatScalar    *aa;
388:   MPI_Comm     comm;
389:   PetscViewer  viewer;

392:   /*
393:       This is nasty. If this is called from an originally parallel matrix
394:    then all processes call this,but only the first has the matrix so the
395:    rest should return immediately.
396:   */
397:   PetscObjectGetComm((PetscObject)draw,&comm);
398:   MPI_Comm_rank(comm,&rank);
399:   if (rank) return(0);

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

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

406:   /* loop over matrix elements drawing boxes */
407:   color = PETSC_DRAW_BLUE;
408:   for (i=0,row=0; i<mbs; i++,row+=bs) {
409:     for (j=a->i[i]; j<a->i[i+1]; j++) {
410:       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
411:       x_l = a->j[j]*bs; x_r = x_l + 1.0;
412:       aa = a->a + j*bs2;
413:       for (k=0; k<bs; k++) {
414:         for (l=0; l<bs; l++) {
415:           if (PetscRealPart(*aa++) >=  0.) continue;
416:           PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
417:         }
418:       }
419:     }
420:   }
421:   color = PETSC_DRAW_CYAN;
422:   for (i=0,row=0; i<mbs; i++,row+=bs) {
423:     for (j=a->i[i]; j<a->i[i+1]; j++) {
424:       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
425:       x_l = a->j[j]*bs; x_r = x_l + 1.0;
426:       aa = a->a + j*bs2;
427:       for (k=0; k<bs; k++) {
428:         for (l=0; l<bs; l++) {
429:           if (PetscRealPart(*aa++) != 0.) continue;
430:           PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
431:         }
432:       }
433:     }
434:   }

436:   color = PETSC_DRAW_RED;
437:   for (i=0,row=0; i<mbs; i++,row+=bs) {
438:     for (j=a->i[i]; j<a->i[i+1]; j++) {
439:       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
440:       x_l = a->j[j]*bs; x_r = x_l + 1.0;
441:       aa = a->a + j*bs2;
442:       for (k=0; k<bs; k++) {
443:         for (l=0; l<bs; l++) {
444:           if (PetscRealPart(*aa++) <= 0.) continue;
445:           PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
446:         }
447:       }
448:     }
449:   }
450:   return(0);
451: }

453: static int MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer)
454: {
455:   int          ierr;
456:   PetscReal    xl,yl,xr,yr,w,h;
457:   PetscDraw    draw;
458:   PetscTruth   isnull;


462:   PetscViewerDrawGetDraw(viewer,0,&draw);
463:   PetscDrawIsNull(draw,&isnull); if (isnull) return(0);

465:   PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);
466:   xr  = A->m; yr = A->m; h = yr/10.0; w = xr/10.0;
467:   xr += w;    yr += h;  xl = -w;     yl = -h;
468:   PetscDrawSetCoordinates(draw,xl,yl,xr,yr);
469:   PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);
470:   PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);
471:   return(0);
472: }

474: int MatView_SeqSBAIJ(Mat A,PetscViewer viewer)
475: {
476:   int        ierr;
477:   PetscTruth issocket,isascii,isbinary,isdraw;

480:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);
481:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
482:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
483:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
484:   if (issocket) {
485:     SETERRQ(PETSC_ERR_SUP,"Socket viewer not supported");
486:   } else if (isascii){
487:     MatView_SeqSBAIJ_ASCII(A,viewer);
488:   } else if (isbinary) {
489:     MatView_SeqSBAIJ_Binary(A,viewer);
490:   } else if (isdraw) {
491:     MatView_SeqSBAIJ_Draw(A,viewer);
492:   } else {
493:     SETERRQ1(1,"Viewer type %s not supported by SeqBAIJ matrices",((PetscObject)viewer)->type_name);
494:   }
495:   return(0);
496: }


499: int MatGetValues_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v)
500: {
501:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
502:   int        *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
503:   int        *ai = a->i,*ailen = a->ilen;
504:   int        brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2;
505:   MatScalar  *ap,*aa = a->a,zero = 0.0;

508:   for (k=0; k<m; k++) { /* loop over rows */
509:     row  = im[k]; brow = row/bs;
510:     if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row");
511:     if (row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
512:     rp   = aj + ai[brow] ; ap = aa + bs2*ai[brow] ;
513:     nrow = ailen[brow];
514:     for (l=0; l<n; l++) { /* loop over columns */
515:       if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column");
516:       if (in[l] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
517:       col  = in[l] ;
518:       bcol = col/bs;
519:       cidx = col%bs;
520:       ridx = row%bs;
521:       high = nrow;
522:       low  = 0; /* assume unsorted */
523:       while (high-low > 5) {
524:         t = (low+high)/2;
525:         if (rp[t] > bcol) high = t;
526:         else             low  = t;
527:       }
528:       for (i=low; i<high; i++) {
529:         if (rp[i] > bcol) break;
530:         if (rp[i] == bcol) {
531:           *v++ = ap[bs2*i+bs*cidx+ridx];
532:           goto finished;
533:         }
534:       }
535:       *v++ = zero;
536:       finished:;
537:     }
538:   }
539:   return(0);
540: }


543: int MatSetValuesBlocked_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
544: {
546:   SETERRQ(1,"Function not yet written for SBAIJ format");
547:   /* return(0); */
548: }

550: int MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode)
551: {
552:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
553:   int        fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
554:   int        m = A->m,*ip,N,*ailen = a->ilen;
555:   int        mbs = a->mbs,bs2 = a->bs2,rmax = 0,ierr;
556:   MatScalar  *aa = a->a,*ap;

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

561:   if (m) rmax = ailen[0];
562:   for (i=1; i<mbs; i++) {
563:     /* move each row back by the amount of empty slots (fshift) before it*/
564:     fshift += imax[i-1] - ailen[i-1];
565:     rmax   = PetscMax(rmax,ailen[i]);
566:     if (fshift) {
567:       ip = aj + ai[i]; ap = aa + bs2*ai[i];
568:       N = ailen[i];
569:       for (j=0; j<N; j++) {
570:         ip[j-fshift] = ip[j];
571:         PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));
572:       }
573:     }
574:     ai[i] = ai[i-1] + ailen[i-1];
575:   }
576:   if (mbs) {
577:     fshift += imax[mbs-1] - ailen[mbs-1];
578:     ai[mbs] = ai[mbs-1] + ailen[mbs-1];
579:   }
580:   /* reset ilen and imax for each row */
581:   for (i=0; i<mbs; i++) {
582:     ailen[i] = imax[i] = ai[i+1] - ai[i];
583:   }
584:   a->s_nz = ai[mbs];

586:   /* diagonals may have moved, so kill the diagonal pointers */
587:   if (fshift && a->diag) {
588:     PetscFree(a->diag);
589:     PetscLogObjectMemory(A,-(m+1)*sizeof(int));
590:     a->diag = 0;
591:   }
592:   PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d usedn",
593:            m,A->m,a->bs,fshift*bs2,a->s_nz*bs2);
594:   PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Number of mallocs during MatSetValues is %dn",
595:            a->reallocs);
596:   PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Most nonzeros blocks in any row is %dn",rmax);
597:   a->reallocs          = 0;
598:   A->info.nz_unneeded  = (PetscReal)fshift*bs2;

600:   return(0);
601: }

603: /* 
604:    This function returns an array of flags which indicate the locations of contiguous
605:    blocks that should be zeroed. for eg: if bs = 3  and is = [0,1,2,3,5,6,7,8,9]
606:    then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)]
607:    Assume: sizes should be long enough to hold all the values.
608: */
609: static int MatZeroRows_SeqSBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max)
610: {
611:   int        i,j,k,row;
612:   PetscTruth flg;

615:   for (i=0,j=0; i<n; j++) {
616:     row = idx[i];
617:     if (row%bs!=0) { /* Not the begining of a block */
618:       sizes[j] = 1;
619:       i++;
620:     } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */
621:       sizes[j] = 1;         /* Also makes sure atleast 'bs' values exist for next else */
622:       i++;
623:     } else { /* Begining of the block, so check if the complete block exists */
624:       flg = PETSC_TRUE;
625:       for (k=1; k<bs; k++) {
626:         if (row+k != idx[i+k]) { /* break in the block */
627:           flg = PETSC_FALSE;
628:           break;
629:         }
630:       }
631:       if (flg == PETSC_TRUE) { /* No break in the bs */
632:         sizes[j] = bs;
633:         i+= bs;
634:       } else {
635:         sizes[j] = 1;
636:         i++;
637:       }
638:     }
639:   }
640:   *bs_max = j;
641:   return(0);
642: }
643: 
644: int MatZeroRows_SeqSBAIJ(Mat A,IS is,Scalar *diag)
645: {
646:   Mat_SeqSBAIJ *sbaij=(Mat_SeqSBAIJ*)A->data;
647:   int         ierr,i,j,k,count,is_n,*is_idx,*rows;
648:   int         bs=sbaij->bs,bs2=sbaij->bs2,*sizes,row,bs_max;
649:   Scalar      zero = 0.0;
650:   MatScalar   *aa;

653:   /* Make a copy of the IS and  sort it */
654:   ISGetSize(is,&is_n);
655:   ISGetIndices(is,&is_idx);

657:   /* allocate memory for rows,sizes */
658:   PetscMalloc((3*is_n+1)*sizeof(int),&rows);
659:   sizes = rows + is_n;

661:   /* initialize copy IS values to rows, and sort them */
662:   for (i=0; i<is_n; i++) { rows[i] = is_idx[i]; }
663:   PetscSortInt(is_n,rows);
664:   if (sbaij->keepzeroedrows) { /* do not change nonzero structure */
665:     for (i=0; i<is_n; i++) { sizes[i] = 1; } /* sizes: size of blocks, = 1 or bs */
666:     bs_max = is_n;            /* bs_max: num. of contiguous block row in the row */
667:   } else {
668:     MatZeroRows_SeqSBAIJ_Check_Blocks(rows,is_n,bs,sizes,&bs_max);
669:   }
670:   ISRestoreIndices(is,&is_idx);

672:   for (i=0,j=0; i<bs_max; j+=sizes[i],i++) {
673:     row   = rows[j];                  /* row to be zeroed */
674:     if (row < 0 || row > A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %d out of range",row);
675:     count = (sbaij->i[row/bs +1] - sbaij->i[row/bs])*bs; /* num. of elements in the row */
676:     aa    = sbaij->a + sbaij->i[row/bs]*bs2 + (row%bs);
677:     if (sizes[i] == bs && !sbaij->keepzeroedrows) {
678:       if (diag) {
679:         if (sbaij->ilen[row/bs] > 0) {
680:           sbaij->ilen[row/bs] = 1;
681:           sbaij->j[sbaij->i[row/bs]] = row/bs;
682:           PetscMemzero(aa,count*bs*sizeof(MatScalar));
683:         }
684:         /* Now insert all the diagoanl values for this bs */
685:         for (k=0; k<bs; k++) {
686:           (*A->ops->setvalues)(A,1,rows+j+k,1,rows+j+k,diag,INSERT_VALUES);
687:         }
688:       } else { /* (!diag) */
689:         sbaij->ilen[row/bs] = 0;
690:       } /* end (!diag) */
691:     } else { /* (sizes[i] != bs), broken block */
692: #if defined (PETSC_USE_DEBUG)
693:       if (sizes[i] != 1) SETERRQ(1,"Internal Error. Value should be 1");
694: #endif
695:       for (k=0; k<count; k++) {
696:         aa[0] = zero;
697:         aa+=bs;
698:       }
699:       if (diag) {
700:         (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);
701:       }
702:     }
703:   }

705:   PetscFree(rows);
706:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
707:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
708:   return(0);
709: }

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

715: int MatSetValues_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
716: {
717:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
718:   int         *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
719:   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented;
720:   int         *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol;
721:   int         ridx,cidx,bs2=a->bs2,ierr;
722:   MatScalar   *ap,value,*aa=a->a,*bap;


726:   for (k=0; k<m; k++) { /* loop over added rows */
727:     row  = im[k];       /* row number */
728:     brow = row/bs;      /* block row number */
729:     if (row < 0) continue;
730: #if defined(PETSC_USE_BOPT_g)  
731:     if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,A->m);
732: #endif
733:     rp   = aj + ai[brow]; /*ptr to beginning of column value of the row block*/
734:     ap   = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/
735:     rmax = imax[brow];  /* maximum space allocated for this row */
736:     nrow = ailen[brow]; /* actual length of this row */
737:     low  = 0;

739:     for (l=0; l<n; l++) { /* loop over added columns */
740:       if (in[l] < 0) continue;
741: #if defined(PETSC_USE_BOPT_g)  
742:       if (in[l] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],A->m);
743: #endif
744:       col = in[l];
745:       bcol = col/bs;              /* block col number */

747:       if (brow <= bcol){
748:         ridx = row % bs; cidx = col % bs; /*row and col index inside the block */
749:         if ((brow==bcol && ridx<=cidx) || (brow<bcol)){
750:           /* element value a(k,l) */
751:           if (roworiented) {
752:             value = v[l + k*n];
753:           } else {
754:             value = v[k + l*m];
755:           }

757:           /* move pointer bap to a(k,l) quickly and add/insert value */
758:           if (!sorted) low = 0; high = nrow;
759:           while (high-low > 7) {
760:             t = (low+high)/2;
761:             if (rp[t] > bcol) high = t;
762:             else              low  = t;
763:           }
764:           for (i=low; i<high; i++) {
765:             /* printf("The loop of i=low.., rp[%d]=%dn",i,rp[i]); */
766:             if (rp[i] > bcol) break;
767:             if (rp[i] == bcol) {
768:               bap  = ap +  bs2*i + bs*cidx + ridx;
769:               if (is == ADD_VALUES) *bap += value;
770:               else                  *bap  = value;
771:               /* for diag block, add/insert its symmetric element a(cidx,ridx) */
772:               if (brow == bcol && ridx < cidx){
773:                 bap  = ap +  bs2*i + bs*ridx + cidx;
774:                 if (is == ADD_VALUES) *bap += value;
775:                 else                  *bap  = value;
776:               }
777:               goto noinsert1;
778:             }
779:           }
780: 
781:       if (nonew == 1) goto noinsert1;
782:       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
783:       if (nrow >= rmax) {
784:         /* there is no extra room in row, therefore enlarge */
785:         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
786:         MatScalar *new_a;

788:         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");

790:         /* Malloc new storage space */
791:         len   = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
792:         ierr  = PetscMalloc(len,&new_a);
793:         new_j = (int*)(new_a + bs2*new_nz);
794:         new_i = new_j + new_nz;

796:         /* copy over old data into new slots */
797:         for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];}
798:         for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;}
799:         PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));
800:         len  = (new_nz - CHUNKSIZE - ai[brow] - nrow);
801:         PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));
802:         PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));
803:         PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));
804:         PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));
805:         /* free up old matrix storage */
806:         PetscFree(a->a);
807:         if (!a->singlemalloc) {
808:           PetscFree(a->i);
809:           PetscFree(a->j);
810:         }
811:         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
812:         a->singlemalloc = PETSC_TRUE;

814:         rp   = aj + ai[brow]; ap = aa + bs2*ai[brow];
815:         rmax = imax[brow] = imax[brow] + CHUNKSIZE;
816:         PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
817:         a->s_maxnz += bs2*CHUNKSIZE;
818:         a->reallocs++;
819:         a->s_nz++;
820:       }

822:       N = nrow++ - 1;
823:       /* shift up all the later entries in this row */
824:       for (ii=N; ii>=i; ii--) {
825:         rp[ii+1] = rp[ii];
826:         ierr     = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));
827:       }
828:       if (N>=i) {
829:         PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));
830:       }
831:       rp[i]                      = bcol;
832:       ap[bs2*i + bs*cidx + ridx] = value;
833:       noinsert1:;
834:       low = i;
835:       /* } */
836:         }
837:       } /* end of if .. if..  */
838:     }                     /* end of loop over added columns */
839:     ailen[brow] = nrow;
840:   }                       /* end of loop over added rows */

842:   return(0);
843: }

845: extern int MatCholeskyFactorSymbolic_SeqSBAIJ(Mat,IS,PetscReal,Mat*);
846: extern int MatCholeskyFactor_SeqSBAIJ(Mat,IS,PetscReal);
847: extern int MatIncreaseOverlap_SeqSBAIJ(Mat,int,IS*,int);
848: extern int MatGetSubMatrix_SeqSBAIJ(Mat,IS,IS,int,MatReuse,Mat*);
849: extern int MatGetSubMatrices_SeqSBAIJ(Mat,int,IS*,IS*,MatReuse,Mat**);
850: extern int MatMultTranspose_SeqSBAIJ(Mat,Vec,Vec);
851: extern int MatMultTransposeAdd_SeqSBAIJ(Mat,Vec,Vec,Vec);
852: extern int MatScale_SeqSBAIJ(Scalar*,Mat);
853: extern int MatNorm_SeqSBAIJ(Mat,NormType,PetscReal *);
854: extern int MatEqual_SeqSBAIJ(Mat,Mat,PetscTruth*);
855: extern int MatGetDiagonal_SeqSBAIJ(Mat,Vec);
856: extern int MatDiagonalScale_SeqSBAIJ(Mat,Vec,Vec);
857: extern int MatGetInfo_SeqSBAIJ(Mat,MatInfoType,MatInfo *);
858: extern int MatZeroEntries_SeqSBAIJ(Mat);
859: extern int MatGetRowMax_SeqSBAIJ(Mat,Vec);

861: extern int MatSolve_SeqSBAIJ_N(Mat,Vec,Vec);
862: extern int MatSolve_SeqSBAIJ_1(Mat,Vec,Vec);
863: extern int MatSolve_SeqSBAIJ_2(Mat,Vec,Vec);
864: extern int MatSolve_SeqSBAIJ_3(Mat,Vec,Vec);
865: extern int MatSolve_SeqSBAIJ_4(Mat,Vec,Vec);
866: extern int MatSolve_SeqSBAIJ_5(Mat,Vec,Vec);
867: extern int MatSolve_SeqSBAIJ_6(Mat,Vec,Vec);
868: extern int MatSolve_SeqSBAIJ_7(Mat,Vec,Vec);
869: extern int MatSolveTranspose_SeqSBAIJ_7(Mat,Vec,Vec);
870: extern int MatSolveTranspose_SeqSBAIJ_6(Mat,Vec,Vec);
871: extern int MatSolveTranspose_SeqSBAIJ_5(Mat,Vec,Vec);
872: extern int MatSolveTranspose_SeqSBAIJ_4(Mat,Vec,Vec);
873: extern int MatSolveTranspose_SeqSBAIJ_3(Mat,Vec,Vec);
874: extern int MatSolveTranspose_SeqSBAIJ_2(Mat,Vec,Vec);
875: extern int MatSolveTranspose_SeqSBAIJ_1(Mat,Vec,Vec);

877: extern int MatSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec);
878: extern int MatSolve_SeqSBAIJ_2_NaturalOrdering(Mat,Vec,Vec);
879: extern int MatSolve_SeqSBAIJ_3_NaturalOrdering(Mat,Vec,Vec);
880: extern int MatSolve_SeqSBAIJ_4_NaturalOrdering(Mat,Vec,Vec);
881: extern int MatSolve_SeqSBAIJ_5_NaturalOrdering(Mat,Vec,Vec);
882: extern int MatSolve_SeqSBAIJ_6_NaturalOrdering(Mat,Vec,Vec);
883: extern int MatSolve_SeqSBAIJ_7_NaturalOrdering(Mat,Vec,Vec);
884: extern int MatSolve_SeqSBAIJ_N_NaturalOrdering(Mat,Vec,Vec);

886: extern int MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat,Mat*);
887: extern int MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat,Mat*);
888: extern int MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat,Mat*);
889: extern int MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat,Mat*);
890: extern int MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat,Mat*);
891: extern int MatCholeskyFactorNumeric_SeqSBAIJ_5(Mat,Mat*);
892: extern int MatCholeskyFactorNumeric_SeqSBAIJ_6(Mat,Mat*);
893: extern int MatCholeskyFactorNumeric_SeqSBAIJ_7(Mat,Mat*);

895: extern int MatMult_SeqSBAIJ_1(Mat,Vec,Vec);
896: extern int MatMult_SeqSBAIJ_2(Mat,Vec,Vec);
897: extern int MatMult_SeqSBAIJ_3(Mat,Vec,Vec);
898: extern int MatMult_SeqSBAIJ_4(Mat,Vec,Vec);
899: extern int MatMult_SeqSBAIJ_5(Mat,Vec,Vec);
900: extern int MatMult_SeqSBAIJ_6(Mat,Vec,Vec);
901: extern int MatMult_SeqSBAIJ_7(Mat,Vec,Vec);
902: extern int MatMult_SeqSBAIJ_N(Mat,Vec,Vec);

904: extern int MatMultAdd_SeqSBAIJ_1(Mat,Vec,Vec,Vec);
905: extern int MatMultAdd_SeqSBAIJ_2(Mat,Vec,Vec,Vec);
906: extern int MatMultAdd_SeqSBAIJ_3(Mat,Vec,Vec,Vec);
907: extern int MatMultAdd_SeqSBAIJ_4(Mat,Vec,Vec,Vec);
908: extern int MatMultAdd_SeqSBAIJ_5(Mat,Vec,Vec,Vec);
909: extern int MatMultAdd_SeqSBAIJ_6(Mat,Vec,Vec,Vec);
910: extern int MatMultAdd_SeqSBAIJ_7(Mat,Vec,Vec,Vec);
911: extern int MatMultAdd_SeqSBAIJ_N(Mat,Vec,Vec,Vec);

913: #ifdef HAVE_MatICCFactor
914: /* modefied from MatILUFactor_SeqSBAIJ, needs further work!  */
915: int MatICCFactor_SeqSBAIJ(Mat inA,IS row,PetscReal fill,int level)
916: {
917:   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data;
918:   Mat         outA;
919:   int         ierr;
920:   PetscTruth  row_identity,col_identity;

923:   /*
924:   if (level != 0) SETERRQ(PETSC_ERR_SUP,"Only levels = 0 supported for in-place ILU"); 
925:   ISIdentity(row,&row_identity);
926:   ISIdentity(col,&col_identity);
927:   if (!row_identity || !col_identity) {
928:     SETERRQ(1,"Row and column permutations must be identity for in-place ICC");
929:   }
930:   */

932:   outA          = inA;
933:   inA->factor   = FACTOR_CHOLESKY;

935:   if (!a->diag) {
936:     MatMarkDiagonal_SeqSBAIJ(inA);
937:   }
938:   /*
939:       Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver 
940:       for ILU(0) factorization with natural ordering
941:   */
942:   switch (a->bs) {
943:   case 1:
944:     inA->ops->solve            = MatSolve_SeqSBAIJ_1_NaturalOrdering;
945:     inA->ops->solvetranspose   = MatSolve_SeqSBAIJ_1_NaturalOrdering;
946:     PetscLoginfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering solvetrans BS=1n");
947:   case 2:
948:     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering;
949:     inA->ops->solve           = MatSolve_SeqSBAIJ_2_NaturalOrdering;
950:     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_2_NaturalOrdering;
951:     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=2n");
952:     break;
953:   case 3:
954:     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering;
955:     inA->ops->solve           = MatSolve_SeqSBAIJ_3_NaturalOrdering;
956:     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_3_NaturalOrdering;
957:     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=3n");
958:     break;
959:   case 4:
960:     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering;
961:     inA->ops->solve           = MatSolve_SeqSBAIJ_4_NaturalOrdering;
962:     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_4_NaturalOrdering;
963:     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=4n");
964:     break;
965:   case 5:
966:     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering;
967:     inA->ops->solve           = MatSolve_SeqSBAIJ_5_NaturalOrdering;
968:     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_5_NaturalOrdering;
969:     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=5n");
970:     break;
971:   case 6:
972:     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering;
973:     inA->ops->solve           = MatSolve_SeqSBAIJ_6_NaturalOrdering;
974:     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_6_NaturalOrdering;
975:     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=6n");
976:     break;
977:   case 7:
978:     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering;
979:     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_7_NaturalOrdering;
980:     inA->ops->solve           = MatSolve_SeqSBAIJ_7_NaturalOrdering;
981:     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=7n");
982:     break;
983:   default:
984:     a->row        = row;
985:     a->icol       = col;
986:     ierr          = PetscObjectReference((PetscObject)row);
987:     ierr          = PetscObjectReference((PetscObject)col);

989:     /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */
990:     ISInvertPermutation(col,PETSC_DECIDE, &(a->icol));
991:     PetscLogObjectParent(inA,a->icol);

993:     if (!a->solve_work) {
994:       PetscMalloc((A->m+a->bs)*sizeof(Scalar),&a->solve_work);
995:       PetscLogObjectMemory(inA,(A->m+a->bs)*sizeof(Scalar));
996:     }
997:   }

999:   MatCholeskyFactorNumeric(inA,&outA);

1001:   return(0);
1002: }
1003: #endif

1005: int MatPrintHelp_SeqSBAIJ(Mat A)
1006: {
1007:   static PetscTruth called = PETSC_FALSE;
1008:   MPI_Comm          comm = A->comm;
1009:   int               ierr;

1012:   if (called) {return(0);} else called = PETSC_TRUE;
1013:   (*PetscHelpPrintf)(comm," Options for MATSEQSBAIJ and MATMPISBAIJ matrix formats (the defaults):n");
1014:   (*PetscHelpPrintf)(comm,"  -mat_block_size <block_size>n");
1015:   return(0);
1016: }

1018: EXTERN_C_BEGIN
1019: int MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,int *indices)
1020: {
1021:   Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data;
1022:   int         i,nz,n;

1025:   nz = baij->s_maxnz;
1026:   n  = mat->n;
1027:   for (i=0; i<nz; i++) {
1028:     baij->j[i] = indices[i];
1029:   }
1030:   baij->s_nz = nz;
1031:   for (i=0; i<n; i++) {
1032:     baij->ilen[i] = baij->imax[i];
1033:   }

1035:   return(0);
1036: }
1037: EXTERN_C_END

1039: /*@
1040:     MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows
1041:        in the matrix.

1043:   Input Parameters:
1044: +  mat     - the SeqSBAIJ matrix
1045: -  indices - the column indices

1047:   Level: advanced

1049:   Notes:
1050:     This can be called if you have precomputed the nonzero structure of the 
1051:   matrix and want to provide it to the matrix object to improve the performance
1052:   of the MatSetValues() operation.

1054:     You MUST have set the correct numbers of nonzeros per row in the call to 
1055:   MatCreateSeqSBAIJ().

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

1059: .seealso: MatCreateSeqSBAIJ
1060: @*/
1061: int MatSeqSBAIJSetColumnIndices(Mat mat,int *indices)
1062: {
1063:   int ierr,(*f)(Mat,int *);

1067:   PetscObjectQueryFunction((PetscObject)mat,"MatSeqSBAIJSetColumnIndices_C",(void (**)())&f);
1068:   if (f) {
1069:     (*f)(mat,indices);
1070:   } else {
1071:     SETERRQ(1,"Wrong type of matrix to set column indices");
1072:   }
1073:   return(0);
1074: }

1076: int MatSetUpPreallocation_SeqSBAIJ(Mat A)
1077: {
1078:   int        ierr;

1081:    MatSeqSBAIJSetPreallocation(A,1,PETSC_DEFAULT,0);
1082:   return(0);
1083: }

1085: /* -------------------------------------------------------------------*/
1086: static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ,
1087:        MatGetRow_SeqSBAIJ,
1088:        MatRestoreRow_SeqSBAIJ,
1089:        MatMult_SeqSBAIJ_N,
1090:        MatMultAdd_SeqSBAIJ_N,
1091:        MatMultTranspose_SeqSBAIJ,
1092:        MatMultTransposeAdd_SeqSBAIJ,
1093:        MatSolve_SeqSBAIJ_N,
1094:        0,
1095:        0,
1096:        0,
1097:        0,
1098:        MatCholeskyFactor_SeqSBAIJ,
1099:        0,
1100:        MatTranspose_SeqSBAIJ,
1101:        MatGetInfo_SeqSBAIJ,
1102:        MatEqual_SeqSBAIJ,
1103:        MatGetDiagonal_SeqSBAIJ,
1104:        MatDiagonalScale_SeqSBAIJ,
1105:        MatNorm_SeqSBAIJ,
1106:        0,
1107:        MatAssemblyEnd_SeqSBAIJ,
1108:        0,
1109:        MatSetOption_SeqSBAIJ,
1110:        MatZeroEntries_SeqSBAIJ,
1111:        MatZeroRows_SeqSBAIJ,
1112:        0,
1113:        0,
1114:        MatCholeskyFactorSymbolic_SeqSBAIJ,
1115:        MatCholeskyFactorNumeric_SeqSBAIJ_N,
1116:        MatSetUpPreallocation_SeqSBAIJ,
1117:        0,
1118:        MatGetOwnershipRange_SeqSBAIJ,
1119:        0,
1120:        MatICCFactorSymbolic_SeqSBAIJ,
1121:        0,
1122:        0,
1123:        MatDuplicate_SeqSBAIJ,
1124:        0,
1125:        0,
1126:        0,
1127:        0,
1128:        0,
1129:        MatGetSubMatrices_SeqSBAIJ,
1130:        MatIncreaseOverlap_SeqSBAIJ,
1131:        MatGetValues_SeqSBAIJ,
1132:        0,
1133:        MatPrintHelp_SeqSBAIJ,
1134:        MatScale_SeqSBAIJ,
1135:        0,
1136:        0,
1137:        0,
1138:        MatGetBlockSize_SeqSBAIJ,
1139:        MatGetRowIJ_SeqSBAIJ,
1140:        MatRestoreRowIJ_SeqSBAIJ,
1141:        0,
1142:        0,
1143:        0,
1144:        0,
1145:        0,
1146:        0,
1147:        MatSetValuesBlocked_SeqSBAIJ,
1148:        MatGetSubMatrix_SeqSBAIJ,
1149:        0,
1150:        0,
1151:        MatGetMaps_Petsc,
1152:        0,
1153:        0,
1154:        0,
1155:        0,
1156:        0,
1157:        0,
1158:        MatGetRowMax_SeqSBAIJ};

1160: EXTERN_C_BEGIN
1161: int MatStoreValues_SeqSBAIJ(Mat mat)
1162: {
1163:   Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data;
1164:   int         nz = aij->i[mat->m]*aij->bs*aij->bs2;
1165:   int         ierr;

1168:   if (aij->nonew != 1) {
1169:     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
1170:   }

1172:   /* allocate space for values if not already there */
1173:   if (!aij->saved_values) {
1174:     PetscMalloc((nz+1)*sizeof(Scalar),&aij->saved_values);
1175:   }

1177:   /* copy values over */
1178:   PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(Scalar));
1179:   return(0);
1180: }
1181: EXTERN_C_END

1183: EXTERN_C_BEGIN
1184: int MatRetrieveValues_SeqSBAIJ(Mat mat)
1185: {
1186:   Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data;
1187:   int         nz = aij->i[mat->m]*aij->bs*aij->bs2,ierr;

1190:   if (aij->nonew != 1) {
1191:     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
1192:   }
1193:   if (!aij->saved_values) {
1194:     SETERRQ(1,"Must call MatStoreValues(A);first");
1195:   }

1197:   /* copy values over */
1198:   PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(Scalar));
1199:   return(0);
1200: }
1201: EXTERN_C_END

1203: EXTERN_C_BEGIN
1204: int MatCreate_SeqSBAIJ(Mat B)
1205: {
1206:   Mat_SeqSBAIJ *b;
1207:   int          ierr,size;

1210:   MPI_Comm_size(B->comm,&size);
1211:   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
1212:   B->m = B->M = PetscMax(B->m,B->M);
1213:   B->n = B->N = PetscMax(B->n,B->N);

1215:   ierr    = PetscNew(Mat_SeqSBAIJ,&b);
1216:   B->data = (void*)b;
1217:   ierr    = PetscMemzero(b,sizeof(Mat_SeqSBAIJ));
1218:   ierr    = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
1219:   B->ops->destroy     = MatDestroy_SeqSBAIJ;
1220:   B->ops->view        = MatView_SeqSBAIJ;
1221:   B->factor           = 0;
1222:   B->lupivotthreshold = 1.0;
1223:   B->mapping          = 0;
1224:   b->row              = 0;
1225:   b->icol             = 0;
1226:   b->reallocs         = 0;
1227:   b->saved_values     = 0;
1228: 
1229:   MapCreateMPI(B->comm,B->m,B->M,&B->rmap);
1230:   MapCreateMPI(B->comm,B->n,B->N,&B->cmap);

1232:   b->sorted           = PETSC_FALSE;
1233:   b->roworiented      = PETSC_TRUE;
1234:   b->nonew            = 0;
1235:   b->diag             = 0;
1236:   b->solve_work       = 0;
1237:   b->mult_work        = 0;
1238:   b->spptr            = 0;
1239:   b->keepzeroedrows   = PETSC_FALSE;
1240: 
1241:   b->inew             = 0;
1242:   b->jnew             = 0;
1243:   b->anew             = 0;
1244:   b->a2anew           = 0;
1245:   b->permute          = PETSC_FALSE;

1247:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
1248:                                      "MatStoreValues_SeqSBAIJ",
1249:                                      (void*)MatStoreValues_SeqSBAIJ);
1250:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
1251:                                      "MatRetrieveValues_SeqSBAIJ",
1252:                                      (void*)MatRetrieveValues_SeqSBAIJ);
1253:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C",
1254:                                      "MatSeqSBAIJSetColumnIndices_SeqSBAIJ",
1255:                                      (void*)MatSeqSBAIJSetColumnIndices_SeqSBAIJ);
1256:   return(0);
1257: }
1258: EXTERN_C_END

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

1267:    Collective on Mat

1269:    Input Parameters:
1270: +  A - the symmetric matrix 
1271: .  bs - size of block
1272: .  nz - number of block nonzeros per block row (same for all rows)
1273: -  nnz - array containing the number of block nonzeros in the various block rows 
1274:          (possibly different for each block row) or PETSC_NULL

1276:    Options Database Keys:
1277: .   -mat_no_unroll - uses code that does not unroll the loops in the 
1278:                      block calculations (much slower)
1279: .    -mat_block_size - size of the blocks to use

1281:    Level: intermediate

1283:    Notes:
1284:    The block AIJ format is compatible with standard Fortran 77
1285:    storage.  That is, the stored row and column indices can begin at
1286:    either one (as in Fortran) or zero.  See the users' manual for details.

1288:    Specify the preallocated storage with either nz or nnz (not both).
1289:    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 
1290:    allocation.  For additional details, see the users manual chapter on
1291:    matrices.

1293: .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ()
1294: @*/
1295: int MatSeqSBAIJSetPreallocation(Mat B,int bs,int nz,int *nnz)
1296: {
1297:   Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data;
1298:   int          i,len,ierr,mbs,bs2;
1299:   PetscTruth   flg;
1300:   int          s_nz;

1303:   B->preallocated = PETSC_TRUE;
1304:   PetscOptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,PETSC_NULL);
1305:   mbs  = B->m/bs;
1306:   bs2  = bs*bs;

1308:   if (mbs*bs != B->m) {
1309:     SETERRQ(PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize");
1310:   }

1312:   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3;
1313:   if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz);
1314:   if (nnz) {
1315:     for (i=0; i<mbs; i++) {
1316:       if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]);
1317:     }
1318:   }

1320:   ierr    = PetscOptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg);
1321:   if (!flg) {
1322:     switch (bs) {
1323:     case 1:
1324:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1;
1325:       B->ops->solve           = MatSolve_SeqSBAIJ_1;
1326:       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_1;
1327:       B->ops->mult            = MatMult_SeqSBAIJ_1;
1328:       B->ops->multadd         = MatMultAdd_SeqSBAIJ_1;
1329:       break;
1330:     case 2:
1331:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2;
1332:       B->ops->solve           = MatSolve_SeqSBAIJ_2;
1333:       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_2;
1334:       B->ops->mult            = MatMult_SeqSBAIJ_2;
1335:       B->ops->multadd         = MatMultAdd_SeqSBAIJ_2;
1336:       break;
1337:     case 3:
1338:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3;
1339:       B->ops->solve           = MatSolve_SeqSBAIJ_3;
1340:       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_3;
1341:       B->ops->mult            = MatMult_SeqSBAIJ_3;
1342:       B->ops->multadd         = MatMultAdd_SeqSBAIJ_3;
1343:       break;
1344:     case 4:
1345:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4;
1346:       B->ops->solve           = MatSolve_SeqSBAIJ_4;
1347:       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_4;
1348:       B->ops->mult            = MatMult_SeqSBAIJ_4;
1349:       B->ops->multadd         = MatMultAdd_SeqSBAIJ_4;
1350:       break;
1351:     case 5:
1352:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5;
1353:       B->ops->solve           = MatSolve_SeqSBAIJ_5;
1354:       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_5;
1355:       B->ops->mult            = MatMult_SeqSBAIJ_5;
1356:       B->ops->multadd         = MatMultAdd_SeqSBAIJ_5;
1357:       break;
1358:     case 6:
1359:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6;
1360:       B->ops->solve           = MatSolve_SeqSBAIJ_6;
1361:       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_6;
1362:       B->ops->mult            = MatMult_SeqSBAIJ_6;
1363:       B->ops->multadd         = MatMultAdd_SeqSBAIJ_6;
1364:       break;
1365:     case 7:
1366:       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7;
1367:       B->ops->solve           = MatSolve_SeqSBAIJ_7;
1368:       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_7;
1369:       B->ops->mult            = MatMult_SeqSBAIJ_7;
1370:       B->ops->multadd         = MatMultAdd_SeqSBAIJ_7;
1371:       break;
1372:     }
1373:   }

1375:   b->mbs = mbs;
1376:   b->nbs = mbs;
1377:   ierr   = PetscMalloc((mbs+1)*sizeof(int),&b->imax);
1378:   if (!nnz) {
1379:     if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
1380:     else if (nz <= 0)        nz = 1;
1381:     for (i=0; i<mbs; i++) {
1382:       b->imax[i] = nz;
1383:     }
1384:     nz = nz*mbs; /* total nz */
1385:   } else {
1386:     nz = 0;
1387:     for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
1388:   }
1389:   /* s_nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */
1390:   s_nz = nz;

1392:   /* allocate the matrix space */
1393:   len  = s_nz*sizeof(int) + s_nz*bs2*sizeof(MatScalar) + (B->m+1)*sizeof(int);
1394:   PetscMalloc(len,&b->a);
1395:   PetscMemzero(b->a,s_nz*bs2*sizeof(MatScalar));
1396:   b->j = (int*)(b->a + s_nz*bs2);
1397:   PetscMemzero(b->j,s_nz*sizeof(int));
1398:   b->i = b->j + s_nz;
1399:   b->singlemalloc = PETSC_TRUE;

1401:   /* pointer to beginning of each row */
1402:   b->i[0] = 0;
1403:   for (i=1; i<mbs+1; i++) {
1404:     b->i[i] = b->i[i-1] + b->imax[i-1];
1405:   }

1407:   /* b->ilen will count nonzeros in each block row so far. */
1408:   PetscMalloc((mbs+1)*sizeof(int),&b->ilen);
1409:   PetscLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ));
1410:   for (i=0; i<mbs; i++) { b->ilen[i] = 0;}
1411: 
1412:   b->bs               = bs;
1413:   b->bs2              = bs2;
1414:   b->s_nz             = 0;
1415:   b->s_maxnz          = s_nz*bs2;
1416: 
1417:   b->inew             = 0;
1418:   b->jnew             = 0;
1419:   b->anew             = 0;
1420:   b->a2anew           = 0;
1421:   b->permute          = PETSC_FALSE;
1422:   return(0);
1423: }


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

1433:    Collective on MPI_Comm

1435:    Input Parameters:
1436: +  comm - MPI communicator, set to PETSC_COMM_SELF
1437: .  bs - size of block
1438: .  m - number of rows, or number of columns
1439: .  nz - number of block nonzeros per block row (same for all rows)
1440: -  nnz - array containing the number of block nonzeros in the various block rows 
1441:          (possibly different for each block row) or PETSC_NULL

1443:    Output Parameter:
1444: .  A - the symmetric matrix 

1446:    Options Database Keys:
1447: .   -mat_no_unroll - uses code that does not unroll the loops in the 
1448:                      block calculations (much slower)
1449: .    -mat_block_size - size of the blocks to use

1451:    Level: intermediate

1453:    Notes:
1454:    The block AIJ format is compatible with standard Fortran 77
1455:    storage.  That is, the stored row and column indices can begin at
1456:    either one (as in Fortran) or zero.  See the users' manual for details.

1458:    Specify the preallocated storage with either nz or nnz (not both).
1459:    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 
1460:    allocation.  For additional details, see the users manual chapter on
1461:    matrices.

1463: .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ()
1464: @*/
1465: int MatCreateSeqSBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz,Mat *A)
1466: {
1468: 
1470:   MatCreate(comm,m,n,m,n,A);
1471:   MatSetType(*A,MATSEQSBAIJ);
1472:   MatSeqSBAIJSetPreallocation(*A,bs,nz,nnz);
1473:   return(0);
1474: }

1476: int MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
1477: {
1478:   Mat         C;
1479:   Mat_SeqSBAIJ *c,*a = (Mat_SeqSBAIJ*)A->data;
1480:   int         i,len,mbs = a->mbs,nz = a->s_nz,bs2 =a->bs2,ierr;

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

1485:   *B = 0;
1486:   MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);
1487:   MatSetType(C,MATSEQSBAIJ);
1488:   c    = (Mat_SeqSBAIJ*)C->data;

1490:   ierr              = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));
1491:   C->preallocated   = PETSC_TRUE;
1492:   C->factor         = A->factor;
1493:   c->row            = 0;
1494:   c->icol           = 0;
1495:   c->saved_values   = 0;
1496:   c->keepzeroedrows = a->keepzeroedrows;
1497:   C->assembled      = PETSC_TRUE;

1499:   c->bs         = a->bs;
1500:   c->bs2        = a->bs2;
1501:   c->mbs        = a->mbs;
1502:   c->nbs        = a->nbs;

1504:   PetscMalloc((mbs+1)*sizeof(int),&c->imax);
1505:   PetscMalloc((mbs+1)*sizeof(int),&c->ilen);
1506:   for (i=0; i<mbs; i++) {
1507:     c->imax[i] = a->imax[i];
1508:     c->ilen[i] = a->ilen[i];
1509:   }

1511:   /* allocate the matrix space */
1512:   c->singlemalloc = PETSC_TRUE;
1513:   len  = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int));
1514:   PetscMalloc(len,&c->a);
1515:   c->j = (int*)(c->a + nz*bs2);
1516:   c->i = c->j + nz;
1517:   PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));
1518:   if (mbs > 0) {
1519:     PetscMemcpy(c->j,a->j,nz*sizeof(int));
1520:     if (cpvalues == MAT_COPY_VALUES) {
1521:       PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));
1522:     } else {
1523:       PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));
1524:     }
1525:   }

1527:   PetscLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ));
1528:   c->sorted      = a->sorted;
1529:   c->roworiented = a->roworiented;
1530:   c->nonew       = a->nonew;

1532:   if (a->diag) {
1533:     PetscMalloc((mbs+1)*sizeof(int),&c->diag);
1534:     PetscLogObjectMemory(C,(mbs+1)*sizeof(int));
1535:     for (i=0; i<mbs; i++) {
1536:       c->diag[i] = a->diag[i];
1537:     }
1538:   } else c->diag        = 0;
1539:   c->s_nz               = a->s_nz;
1540:   c->s_maxnz            = a->s_maxnz;
1541:   c->solve_work         = 0;
1542:   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
1543:   c->mult_work          = 0;
1544:   *B = C;
1545:   PetscFListDuplicate(A->qlist,&C->qlist);
1546:   return(0);
1547: }

1549: EXTERN_C_BEGIN
1550: int MatLoad_SeqSBAIJ(PetscViewer viewer,MatType type,Mat *A)
1551: {
1552:   Mat_SeqSBAIJ *a;
1553:   Mat          B;
1554:   int          i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1;
1555:   int          *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,*s_browlengths,maskcount;
1556:   int          kmax,jcount,block,idx,point,nzcountb,extra_rows;
1557:   int          *masked,nmask,tmp,bs2,ishift;
1558:   Scalar       *aa;
1559:   MPI_Comm     comm = ((PetscObject)viewer)->comm;

1562:   PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);
1563:   bs2  = bs*bs;

1565:   MPI_Comm_size(comm,&size);
1566:   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
1567:   PetscViewerBinaryGetDescriptor(viewer,&fd);
1568:   PetscBinaryRead(fd,header,4,PETSC_INT);
1569:   if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object");
1570:   M = header[1]; N = header[2]; nz = header[3];

1572:   if (header[3] < 0) {
1573:     SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ");
1574:   }

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

1578:   /* 
1579:      This code adds extra rows to make sure the number of rows is 
1580:     divisible by the blocksize
1581:   */
1582:   mbs        = M/bs;
1583:   extra_rows = bs - M + bs*(mbs);
1584:   if (extra_rows == bs) extra_rows = 0;
1585:   else                  mbs++;
1586:   if (extra_rows) {
1587:     PetscLogInfo(0,"MatLoad_SeqSBAIJ:Padding loaded matrix to match blocksizen");
1588:   }

1590:   /* read in row lengths */
1591:   PetscMalloc((M+extra_rows)*sizeof(int),&rowlengths);
1592:   PetscBinaryRead(fd,rowlengths,M,PETSC_INT);
1593:   for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;

1595:   /* read in column indices */
1596:   PetscMalloc((nz+extra_rows)*sizeof(int),&jj);
1597:   PetscBinaryRead(fd,jj,nz,PETSC_INT);
1598:   for (i=0; i<extra_rows; i++) jj[nz+i] = M+i;

1600:   /* loop over row lengths determining block row lengths */
1601:   ierr     = PetscMalloc(mbs*sizeof(int),&browlengths);
1602:   ierr     = PetscMalloc(mbs*sizeof(int),&s_browlengths);
1603:   ierr     = PetscMemzero(s_browlengths,mbs*sizeof(int));
1604:   ierr     = PetscMalloc(2*mbs*sizeof(int),&mask);
1605:   ierr     = PetscMemzero(mask,mbs*sizeof(int));
1606:   masked   = mask + mbs;
1607:   rowcount = 0; nzcount = 0;
1608:   for (i=0; i<mbs; i++) {
1609:     nmask = 0;
1610:     for (j=0; j<bs; j++) {
1611:       kmax = rowlengths[rowcount];
1612:       for (k=0; k<kmax; k++) {
1613:         tmp = jj[nzcount++]/bs;   /* block col. index */
1614:         if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;}
1615:       }
1616:       rowcount++;
1617:     }
1618:     s_browlengths[i] += nmask;
1619:     browlengths[i]    = 2*s_browlengths[i];

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

1625:   /* create our matrix */
1626:   MatCreateSeqSBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);
1627:   B = *A;
1628:   a = (Mat_SeqSBAIJ*)B->data;

1630:   /* set matrix "i" values */
1631:   a->i[0] = 0;
1632:   for (i=1; i<= mbs; i++) {
1633:     a->i[i]      = a->i[i-1] + s_browlengths[i-1];
1634:     a->ilen[i-1] = s_browlengths[i-1];
1635:   }
1636:   a->s_nz         = 0;
1637:   for (i=0; i<mbs; i++) a->s_nz += s_browlengths[i];

1639:   /* read in nonzero values */
1640:   PetscMalloc((nz+extra_rows)*sizeof(Scalar),&aa);
1641:   PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);
1642:   for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0;

1644:   /* set "a" and "j" values into matrix */
1645:   nzcount = 0; jcount = 0;
1646:   for (i=0; i<mbs; i++) {
1647:     nzcountb = nzcount;
1648:     nmask    = 0;
1649:     for (j=0; j<bs; j++) {
1650:       kmax = rowlengths[i*bs+j];
1651:       for (k=0; k<kmax; k++) {
1652:         tmp = jj[nzcount++]/bs; /* block col. index */
1653:         if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;}
1654:       }
1655:     }
1656:     /* sort the masked values */
1657:     PetscSortInt(nmask,masked);

1659:     /* set "j" values into matrix */
1660:     maskcount = 1;
1661:     for (j=0; j<nmask; j++) {
1662:       a->j[jcount++]  = masked[j];
1663:       mask[masked[j]] = maskcount++;
1664:     }

1666:     /* set "a" values into matrix */
1667:     ishift = bs2*a->i[i];
1668:     for (j=0; j<bs; j++) {
1669:       kmax = rowlengths[i*bs+j];
1670:       for (k=0; k<kmax; k++) {
1671:         tmp       = jj[nzcountb]/bs ; /* block col. index */
1672:         if (tmp >= i){
1673:           block     = mask[tmp] - 1;
1674:           point     = jj[nzcountb] - bs*tmp;
1675:           idx       = ishift + bs2*block + j + bs*point;
1676:           a->a[idx] = aa[nzcountb];
1677:         }
1678:         nzcountb++;
1679:       }
1680:     }
1681:     /* zero out the mask elements we set */
1682:     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
1683:   }
1684:   if (jcount != a->s_nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix");

1686:   PetscFree(rowlengths);
1687:   PetscFree(browlengths);
1688:   PetscFree(s_browlengths);
1689:   PetscFree(aa);
1690:   PetscFree(jj);
1691:   PetscFree(mask);

1693:   B->assembled = PETSC_TRUE;
1694:   MatView_Private(B);
1695:   return(0);
1696: }
1697: EXTERN_C_END