Actual source code: bdiag.c

  1: /*$Id: bdiag.c,v 1.193 2001/03/23 23:22:03 balay Exp $*/

  3: /* Block diagonal matrix format */

  5: #include "petscsys.h"
  6: #include "src/mat/impls/bdiag/seq/bdiag.h"
  7: #include "src/vec/vecimpl.h"
  8: #include "src/inline/ilu.h"

 10: EXTERN int MatSetValues_SeqBDiag_1(Mat,int,int *,int,int *,Scalar *,InsertMode);
 11: EXTERN int MatSetValues_SeqBDiag_N(Mat,int,int *,int,int *,Scalar *,InsertMode);
 12: EXTERN int MatGetValues_SeqBDiag_1(Mat,int,int *,int,int *,Scalar *);
 13: EXTERN int MatGetValues_SeqBDiag_N(Mat,int,int *,int,int *,Scalar *);
 14: EXTERN int MatMult_SeqBDiag_1(Mat,Vec,Vec);
 15: EXTERN int MatMult_SeqBDiag_2(Mat,Vec,Vec);
 16: EXTERN int MatMult_SeqBDiag_3(Mat,Vec,Vec);
 17: EXTERN int MatMult_SeqBDiag_4(Mat,Vec,Vec);
 18: EXTERN int MatMult_SeqBDiag_5(Mat,Vec,Vec);
 19: EXTERN int MatMult_SeqBDiag_N(Mat,Vec,Vec);
 20: EXTERN int MatMultAdd_SeqBDiag_1(Mat,Vec,Vec,Vec);
 21: EXTERN int MatMultAdd_SeqBDiag_2(Mat,Vec,Vec,Vec);
 22: EXTERN int MatMultAdd_SeqBDiag_3(Mat,Vec,Vec,Vec);
 23: EXTERN int MatMultAdd_SeqBDiag_4(Mat,Vec,Vec,Vec);
 24: EXTERN int MatMultAdd_SeqBDiag_5(Mat,Vec,Vec,Vec);
 25: EXTERN int MatMultAdd_SeqBDiag_N(Mat,Vec,Vec,Vec);
 26: EXTERN int MatMultTranspose_SeqBDiag_1(Mat,Vec,Vec);
 27: EXTERN int MatMultTranspose_SeqBDiag_N(Mat,Vec,Vec);
 28: EXTERN int MatMultTransposeAdd_SeqBDiag_1(Mat,Vec,Vec,Vec);
 29: EXTERN int MatMultTransposeAdd_SeqBDiag_N(Mat,Vec,Vec,Vec);
 30: EXTERN int MatRelax_SeqBDiag_N(Mat,Vec,PetscReal,MatSORType,PetscReal,int,Vec);
 31: EXTERN int MatRelax_SeqBDiag_1(Mat,Vec,PetscReal,MatSORType,PetscReal,int,Vec);
 32: EXTERN int MatView_SeqBDiag(Mat,PetscViewer);
 33: EXTERN int MatGetInfo_SeqBDiag(Mat,MatInfoType,MatInfo*);
 34: EXTERN int MatGetRow_SeqBDiag(Mat,int,int *,int **,Scalar **);
 35: EXTERN int MatRestoreRow_SeqBDiag(Mat,int,int *,int **,Scalar **);
 36: EXTERN int MatTranspose_SeqBDiag(Mat,Mat *);
 37: EXTERN int MatNorm_SeqBDiag(Mat,NormType,PetscReal *);
 38: EXTERN int MatGetOwnershipRange_SeqBDiag(Mat,int*,int *);

 40: int MatDestroy_SeqBDiag(Mat A)
 41: {
 42:   Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
 43:   int          i,bs = a->bs,ierr;

 46: #if defined(PETSC_USE_LOG)
 47:   PetscLogObjectState((PetscObject)A,"Rows=%d, Cols=%d, NZ=%d, BSize=%d, NDiag=%d",A->m,A->n,a->nz,a->bs,a->nd);
 48: #endif
 49:   if (!a->user_alloc) { /* Free the actual diagonals */
 50:     for (i=0; i<a->nd; i++) {
 51:       if (a->diag[i] > 0) {
 52:         PetscFree(a->diagv[i] + bs*bs*a->diag[i]);
 53:       } else {
 54:         PetscFree(a->diagv[i]);
 55:       }
 56:     }
 57:   }
 58:   if (a->pivot) {PetscFree(a->pivot);}
 59:   PetscFree(a->diagv);
 60:   PetscFree(a->diag);
 61:   PetscFree(a->colloc);
 62:   PetscFree(a->dvalue);
 63:   PetscFree(a);
 64:   return(0);
 65: }

 67: int MatAssemblyEnd_SeqBDiag(Mat A,MatAssemblyType mode)
 68: {
 69:   Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
 70:   int          i,k,temp,*diag = a->diag,*bdlen = a->bdlen;
 71:   Scalar       *dtemp,**dv = a->diagv;

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

 76:   /* Sort diagonals */
 77:   for (i=0; i<a->nd; i++) {
 78:     for (k=i+1; k<a->nd; k++) {
 79:       if (diag[i] < diag[k]) {
 80:         temp     = diag[i];
 81:         diag[i]  = diag[k];
 82:         diag[k]  = temp;
 83:         temp     = bdlen[i];
 84:         bdlen[i] = bdlen[k];
 85:         bdlen[k] = temp;
 86:         dtemp    = dv[i];
 87:         dv[i]    = dv[k];
 88:         dv[k]    = dtemp;
 89:       }
 90:     }
 91:   }

 93:   /* Set location of main diagonal */
 94:   for (i=0; i<a->nd; i++) {
 95:     if (!a->diag[i]) {a->mainbd = i; break;}
 96:   }
 97:   PetscLogInfo(A,"MatAssemblyEnd_SeqBDiag:Number diagonals %d,memory used %d, block size %dn",a->nd,a->maxnz,a->bs);
 98:   return(0);
 99: }

101: int MatSetOption_SeqBDiag(Mat A,MatOption op)
102: {
103:   Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;

106:   if (op == MAT_NO_NEW_NONZERO_LOCATIONS)       a->nonew       = 1;
107:   else if (op == MAT_YES_NEW_NONZERO_LOCATIONS) a->nonew       = 0;
108:   else if (op == MAT_NO_NEW_DIAGONALS)          a->nonew_diag  = 1;
109:   else if (op == MAT_YES_NEW_DIAGONALS)         a->nonew_diag  = 0;
110:   else if (op == MAT_COLUMN_ORIENTED)           a->roworiented = PETSC_FALSE;
111:   else if (op == MAT_ROW_ORIENTED)              a->roworiented = PETSC_TRUE;
112:   else if (op == MAT_ROWS_SORTED ||
113:            op == MAT_ROWS_UNSORTED ||
114:            op == MAT_COLUMNS_SORTED ||
115:            op == MAT_COLUMNS_UNSORTED ||
116:            op == MAT_IGNORE_OFF_PROC_ENTRIES ||
117:            op == MAT_NEW_NONZERO_LOCATION_ERR ||
118:            op == MAT_NEW_NONZERO_ALLOCATION_ERR ||
119:            op == MAT_USE_HASH_TABLE)
120:     PetscLogInfo(A,"MatSetOption_SeqBDiag:Option ignoredn");
121:   else {
122:     SETERRQ(PETSC_ERR_SUP,"unknown option");
123:   }
124:   return(0);
125: }

127: int MatPrintHelp_SeqBDiag(Mat A)
128: {
129:   static PetscTruth called = PETSC_FALSE;
130:   MPI_Comm          comm = A->comm;
131:   int               ierr;

134:   if (called) {return(0);} else called = PETSC_TRUE;
135:   (*PetscHelpPrintf)(comm," Options for MATSEQBDIAG and MATMPIBDIAG matrix formats:n");
136:   (*PetscHelpPrintf)(comm,"  -mat_block_size <block_size>n");
137:   (*PetscHelpPrintf)(comm,"  -mat_bdiag_diags <d1,d2,d3,...> (diagonal numbers)n");
138:   (*PetscHelpPrintf)(comm,"   (for example) -mat_bdiag_diags -5,-1,0,1,5n");
139:   return(0);
140: }

142: static int MatGetDiagonal_SeqBDiag_N(Mat A,Vec v)
143: {
144:   Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
145:   int          ierr,i,j,n,len,ibase,bs = a->bs,iloc;
146:   Scalar       *x,*dd,zero = 0.0;

149:   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
150:   VecSet(&zero,v);
151:   VecGetLocalSize(v,&n);
152:   if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
153:   if (a->mainbd == -1) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Main diagonal not set");
154:   len = PetscMin(a->mblock,a->nblock);
155:   dd = a->diagv[a->mainbd];
156:   VecGetArray(v,&x);
157:   for (i=0; i<len; i++) {
158:     ibase = i*bs*bs;  iloc = i*bs;
159:     for (j=0; j<bs; j++) x[j + iloc] = dd[ibase + j*(bs+1)];
160:   }
161:   VecRestoreArray(v,&x);
162:   return(0);
163: }

165: static int MatGetDiagonal_SeqBDiag_1(Mat A,Vec v)
166: {
167:   Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
168:   int          ierr,i,n,len;
169:   Scalar       *x,*dd,zero = 0.0;

172:   VecSet(&zero,v);
173:   VecGetLocalSize(v,&n);
174:   if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
175:   if (a->mainbd == -1) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Main diagonal not set");
176:   dd = a->diagv[a->mainbd];
177:   len = PetscMin(A->m,A->n);
178:   VecGetArray(v,&x);
179:   for (i=0; i<len; i++) x[i] = dd[i];
180:   VecRestoreArray(v,&x);
181:   return(0);
182: }

184: int MatZeroEntries_SeqBDiag(Mat A)
185: {
186:   Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
187:   int          d,i,len,bs = a->bs;
188:   Scalar       *dv;

191:   for (d=0; d<a->nd; d++) {
192:     dv  = a->diagv[d];
193:     if (a->diag[d] > 0) {
194:       dv += bs*bs*a->diag[d];
195:     }
196:     len = a->bdlen[d]*bs*bs;
197:     for (i=0; i<len; i++) dv[i] = 0.0;
198:   }
199:   return(0);
200: }

202: int MatGetBlockSize_SeqBDiag(Mat A,int *bs)
203: {
204:   Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;

207:   *bs = a->bs;
208:   return(0);
209: }

211: int MatZeroRows_SeqBDiag(Mat A,IS is,Scalar *diag)
212: {
213:   Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
214:   int          i,ierr,N,*rows,m = A->m - 1,nz,*col;
215:   Scalar       *dd,*val;

218:   ISGetLocalSize(is,&N);
219:   ISGetIndices(is,&rows);
220:   for (i=0; i<N; i++) {
221:     if (rows[i]<0 || rows[i]>m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"row out of range");
222:     MatGetRow(A,rows[i],&nz,&col,&val);
223:     PetscMemzero(val,nz*sizeof(Scalar));
224:     MatSetValues(A,1,&rows[i],nz,col,val,INSERT_VALUES);
225:     MatRestoreRow(A,rows[i],&nz,&col,&val);
226:   }
227:   if (diag) {
228:     if (a->mainbd == -1) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Main diagonal does not exist");
229:     dd = a->diagv[a->mainbd];
230:     for (i=0; i<N; i++) dd[rows[i]] = *diag;
231:   }
232:   ISRestoreIndices(is,&rows);
233:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
234:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
235:   return(0);
236: }

238: int MatGetSubMatrix_SeqBDiag(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *submat)
239: {
240:   Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
241:   int          nznew,*smap,i,j,ierr,oldcols = A->n;
242:   int          *irow,*icol,newr,newc,*cwork,*col,nz,bs;
243:   Scalar       *vwork,*val;
244:   Mat          newmat;

247:   if (scall == MAT_REUSE_MATRIX) { /* no support for reuse so simply destroy all */
248:     MatDestroy(*submat);
249:   }

251:   ISGetIndices(isrow,&irow);
252:   ISGetIndices(iscol,&icol);
253:   ISGetLocalSize(isrow,&newr);
254:   ISGetLocalSize(iscol,&newc);

256:   PetscMalloc((oldcols+1)*sizeof(int),&smap);
257:   PetscMalloc((newc+1)*sizeof(int),&cwork);
258:   PetscMalloc((newc+1)*sizeof(Scalar),&vwork);
259:   ierr  = PetscMemzero((char*)smap,oldcols*sizeof(int));
260:   for (i=0; i<newc; i++) smap[icol[i]] = i+1;

262:   /* Determine diagonals; then create submatrix */
263:   bs = a->bs; /* Default block size remains the same */
264:   MatCreateSeqBDiag(A->comm,newr,newc,0,bs,0,0,&newmat);

266:   /* Fill new matrix */
267:   for (i=0; i<newr; i++) {
268:     MatGetRow(A,irow[i],&nz,&col,&val);
269:     nznew = 0;
270:     for (j=0; j<nz; j++) {
271:       if (smap[col[j]]) {
272:         cwork[nznew]   = smap[col[j]] - 1;
273:         vwork[nznew++] = val[j];
274:       }
275:     }
276:     MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERT_VALUES);
277:     MatRestoreRow(A,i,&nz,&col,&val);
278:   }
279:   MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);
280:   MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);

282:   /* Free work space */
283:   PetscFree(smap);
284:   PetscFree(cwork);
285:   PetscFree(vwork);
286:   ISRestoreIndices(isrow,&irow);
287:   ISRestoreIndices(iscol,&icol);
288:   *submat = newmat;
289:   return(0);
290: }

292: int MatGetSubMatrices_SeqBDiag(Mat A,int n,IS *irow,IS *icol,MatReuse scall,Mat **B)
293: {
294:   int ierr,i;

297:   if (scall == MAT_INITIAL_MATRIX) {
298:     PetscMalloc((n+1)*sizeof(Mat),B);
299:   }

301:   for (i=0; i<n; i++) {
302:     MatGetSubMatrix_SeqBDiag(A,irow[i],icol[i],scall,&(*B)[i]);
303:   }
304:   return(0);
305: }

307: int MatScale_SeqBDiag(Scalar *alpha,Mat inA)
308: {
309:   Mat_SeqBDiag *a = (Mat_SeqBDiag*)inA->data;
310:   int          one = 1,i,len,bs = a->bs;

313:   for (i=0; i<a->nd; i++) {
314:     len = bs*bs*a->bdlen[i];
315:     if (a->diag[i] > 0) {
316:       BLscal_(&len,alpha,a->diagv[i] + bs*bs*a->diag[i],&one);
317:     } else {
318:       BLscal_(&len,alpha,a->diagv[i],&one);
319:     }
320:   }
321:   PetscLogFlops(a->nz);
322:   return(0);
323: }

325: int MatDiagonalScale_SeqBDiag(Mat A,Vec ll,Vec rr)
326: {
327:   Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
328:   Scalar       *l,*r,*dv;
329:   int          d,j,len,ierr;
330:   int          nd = a->nd,bs = a->bs,diag,m,n;

333:   if (ll) {
334:     VecGetSize(ll,&m);
335:     if (m != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length");
336:     if (bs == 1) {
337:       VecGetArray(ll,&l);
338:       for (d=0; d<nd; d++) {
339:         dv   = a->diagv[d];
340:         diag = a->diag[d];
341:         len  = a->bdlen[d];
342:         if (diag > 0) for (j=0; j<len; j++) dv[j+diag] *= l[j+diag];
343:         else          for (j=0; j<len; j++) dv[j]      *= l[j];
344:       }
345:       VecRestoreArray(ll,&l);
346:       PetscLogFlops(a->nz);
347:     } else SETERRQ(PETSC_ERR_SUP,"Not yet done for bs>1");
348:   }
349:   if (rr) {
350:     VecGetSize(rr,&n);
351:     if (n != A->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vector wrong length");
352:     if (bs == 1) {
353:       VecGetArray(rr,&r);
354:       for (d=0; d<nd; d++) {
355:         dv   = a->diagv[d];
356:         diag = a->diag[d];
357:         len  = a->bdlen[d];
358:         if (diag > 0) for (j=0; j<len; j++) dv[j+diag] *= r[j];
359:         else          for (j=0; j<len; j++) dv[j]      *= r[j-diag];
360:       }
361:       VecRestoreArray(rr,&r);
362:       PetscLogFlops(a->nz);
363:     } else SETERRQ(PETSC_ERR_SUP,"Not yet done for bs>1");
364:   }
365:   return(0);
366: }

368: static int MatDuplicate_SeqBDiag(Mat,MatDuplicateOption,Mat *);
369: EXTERN int MatLUFactorSymbolic_SeqBDiag(Mat,IS,IS,MatLUInfo*,Mat*);
370: EXTERN int MatILUFactorSymbolic_SeqBDiag(Mat,IS,IS,MatILUInfo*,Mat*);
371: EXTERN int MatILUFactor_SeqBDiag(Mat,IS,IS,MatILUInfo*);
372: EXTERN int MatLUFactorNumeric_SeqBDiag_N(Mat,Mat*);
373: EXTERN int MatLUFactorNumeric_SeqBDiag_1(Mat,Mat*);
374: EXTERN int MatSolve_SeqBDiag_1(Mat,Vec,Vec);
375: EXTERN int MatSolve_SeqBDiag_2(Mat,Vec,Vec);
376: EXTERN int MatSolve_SeqBDiag_3(Mat,Vec,Vec);
377: EXTERN int MatSolve_SeqBDiag_4(Mat,Vec,Vec);
378: EXTERN int MatSolve_SeqBDiag_5(Mat,Vec,Vec);
379: EXTERN int MatSolve_SeqBDiag_N(Mat,Vec,Vec);

381: int MatSetUpPreallocation_SeqBDiag(Mat A)
382: {
383:   int        ierr;

386:    MatSeqBDiagSetPreallocation(A,PETSC_DEFAULT,PETSC_DEFAULT,0,0);
387:   return(0);
388: }

390: /* -------------------------------------------------------------------*/
391: static struct _MatOps MatOps_Values = {MatSetValues_SeqBDiag_N,
392:        MatGetRow_SeqBDiag,
393:        MatRestoreRow_SeqBDiag,
394:        MatMult_SeqBDiag_N,
395:        MatMultAdd_SeqBDiag_N,
396:        MatMultTranspose_SeqBDiag_N,
397:        MatMultTransposeAdd_SeqBDiag_N,
398:        MatSolve_SeqBDiag_N,
399:        0,
400:        0,
401:        0,
402:        0,
403:        0,
404:        MatRelax_SeqBDiag_N,
405:        MatTranspose_SeqBDiag,
406:        MatGetInfo_SeqBDiag,
407:        0,
408:        MatGetDiagonal_SeqBDiag_N,
409:        MatDiagonalScale_SeqBDiag,
410:        MatNorm_SeqBDiag,
411:        0,
412:        MatAssemblyEnd_SeqBDiag,
413:        0,
414:        MatSetOption_SeqBDiag,
415:        MatZeroEntries_SeqBDiag,
416:        MatZeroRows_SeqBDiag,
417:        0,
418:        MatLUFactorNumeric_SeqBDiag_N,
419:        0,
420:        0,
421:        MatSetUpPreallocation_SeqBDiag,
422:        0,
423:        MatGetOwnershipRange_SeqBDiag,
424:        MatILUFactorSymbolic_SeqBDiag,
425:        0,
426:        0,
427:        0,
428:        MatDuplicate_SeqBDiag,
429:        0,
430:        0,
431:        MatILUFactor_SeqBDiag,
432:        0,
433:        0,
434:        MatGetSubMatrices_SeqBDiag,
435:        0,
436:        MatGetValues_SeqBDiag_N,
437:        0,
438:        MatPrintHelp_SeqBDiag,
439:        MatScale_SeqBDiag,
440:        0,
441:        0,
442:        0,
443:        MatGetBlockSize_SeqBDiag,
444:        0,
445:        0,
446:        0,
447:        0,
448:        0,
449:        0,
450:        0,
451:        0,
452:        0,
453:        0,
454:        MatDestroy_SeqBDiag,
455:        MatView_SeqBDiag,
456:        MatGetMaps_Petsc};

458: /*@C
459:    MatSeqBDiagSetPreallocation - Sets the nonzero structure and (optionally) arrays.

461:    Collective on MPI_Comm

463:    Input Parameters:
464: +  B - the matrix
465: .  nd - number of block diagonals (optional)
466: .  bs - each element of a diagonal is an bs x bs dense matrix
467: .  diag - optional array of block diagonal numbers (length nd).
468:    For a matrix element A[i,j], where i=row and j=column, the
469:    diagonal number is
470: $     diag = i/bs - j/bs  (integer division)
471:    Set diag=PETSC_NULL on input for PETSc to dynamically allocate memory as 
472:    needed (expensive).
473: -  diagv - pointer to actual diagonals (in same order as diag array), 
474:    if allocated by user.  Otherwise, set diagv=PETSC_NULL on input for PETSc
475:    to control memory allocation.

477:    Options Database Keys:
478: .  -mat_block_size <bs> - Sets blocksize
479: .  -mat_bdiag_diags <s1,s2,s3,...> - Sets diagonal numbers

481:    Notes:
482:    See the users manual for further details regarding this storage format.

484:    Fortran Note:
485:    Fortran programmers cannot set diagv; this value is ignored.

487:    Level: intermediate

489: .keywords: matrix, block, diagonal, sparse

491: .seealso: MatCreate(), MatCreateMPIBDiag(), MatSetValues()
492: @*/
493: int MatSeqBDiagSetPreallocation(Mat B,int nd,int bs,int *diag,Scalar **diagv)
494: {
495:   Mat_SeqBDiag *b;
496:   int          i,nda,sizetot,ierr, nd2 = 128,idiag[128];
497:   PetscTruth   flg1;

500:   PetscTypeCompare((PetscObject)B,MATSEQBDIAG,&flg1);
501:   if (!flg1) return(0);

503:   B->preallocated = PETSC_TRUE;
504:   if (bs == PETSC_DEFAULT) bs = 1;
505:   if (bs == 0) SETERRQ(1,"Blocksize cannot be zero");
506:   if (nd == PETSC_DEFAULT) nd = 0;
507:   PetscOptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,PETSC_NULL);
508:   PetscOptionsGetIntArray(PETSC_NULL,"-mat_bdiag_diags",idiag,&nd2,&flg1);
509:   if (flg1) {
510:     diag = idiag;
511:     nd   = nd2;
512:   }

514:   if ((B->n%bs) || (B->m%bs)) SETERRQ(PETSC_ERR_ARG_SIZ,"Invalid block size");
515:   if (!nd) nda = nd + 1;
516:   else     nda = nd;
517:   b            = (Mat_SeqBDiag*)B->data;

519:   PetscOptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg1);
520:   if (!flg1) {
521:     switch (bs) {
522:       case 1:
523:         B->ops->setvalues       = MatSetValues_SeqBDiag_1;
524:         B->ops->getvalues       = MatGetValues_SeqBDiag_1;
525:         B->ops->getdiagonal     = MatGetDiagonal_SeqBDiag_1;
526:         B->ops->mult            = MatMult_SeqBDiag_1;
527:         B->ops->multadd         = MatMultAdd_SeqBDiag_1;
528:         B->ops->multtranspose   = MatMultTranspose_SeqBDiag_1;
529:         B->ops->multtransposeadd= MatMultTransposeAdd_SeqBDiag_1;
530:         B->ops->relax           = MatRelax_SeqBDiag_1;
531:         B->ops->solve           = MatSolve_SeqBDiag_1;
532:         B->ops->lufactornumeric = MatLUFactorNumeric_SeqBDiag_1;
533:         break;
534:       case 2:
535:         B->ops->mult            = MatMult_SeqBDiag_2;
536:         B->ops->multadd         = MatMultAdd_SeqBDiag_2;
537:         B->ops->solve           = MatSolve_SeqBDiag_2;
538:         break;
539:       case 3:
540:         B->ops->mult            = MatMult_SeqBDiag_3;
541:         B->ops->multadd         = MatMultAdd_SeqBDiag_3;
542:         B->ops->solve           = MatSolve_SeqBDiag_3;
543:         break;
544:       case 4:
545:         B->ops->mult            = MatMult_SeqBDiag_4;
546:         B->ops->multadd         = MatMultAdd_SeqBDiag_4;
547:         B->ops->solve           = MatSolve_SeqBDiag_4;
548:         break;
549:       case 5:
550:         B->ops->mult            = MatMult_SeqBDiag_5;
551:         B->ops->multadd         = MatMultAdd_SeqBDiag_5;
552:         B->ops->solve           = MatSolve_SeqBDiag_5;
553:         break;
554:    }
555:   }

557:   b->mblock = B->m/bs;
558:   b->nblock = B->n/bs;
559:   b->nd     = nd;
560:   b->bs     = bs;
561:   b->ndim   = 0;
562:   b->mainbd = -1;
563:   b->pivot  = 0;

565:   ierr      = PetscMalloc(2*nda*sizeof(int),&b->diag);
566:   b->bdlen  = b->diag + nda;
567:   ierr      = PetscMalloc((B->n+1)*sizeof(int),&b->colloc);
568:   ierr      = PetscMalloc(nda*sizeof(Scalar*),&b->diagv);
569:   sizetot   = 0;

571:   if (diagv) { /* user allocated space */
572:     b->user_alloc = PETSC_TRUE;
573:     for (i=0; i<nd; i++) b->diagv[i] = diagv[i];
574:   } else b->user_alloc = PETSC_FALSE;

576:   for (i=0; i<nd; i++) {
577:     b->diag[i] = diag[i];
578:     if (diag[i] > 0) { /* lower triangular */
579:       b->bdlen[i] = PetscMin(b->nblock,b->mblock - diag[i]);
580:     } else {           /* upper triangular */
581:       b->bdlen[i] = PetscMin(b->mblock,b->nblock + diag[i]);
582:     }
583:     sizetot += b->bdlen[i];
584:   }
585:   sizetot   *= bs*bs;
586:   b->maxnz  =  sizetot;
587:   ierr      = PetscMalloc((B->n+1)*sizeof(Scalar),&b->dvalue);
588:   PetscLogObjectMemory(B,(nda*(bs+2))*sizeof(int) + bs*nda*sizeof(Scalar)
589:                     + nda*sizeof(Scalar*) + sizeof(Mat_SeqBDiag)
590:                     + sizeof(struct _p_Mat) + sizetot*sizeof(Scalar));

592:   if (!b->user_alloc) {
593:     for (i=0; i<nd; i++) {
594:       PetscMalloc(bs*bs*b->bdlen[i]*sizeof(Scalar),&b->diagv[i]);
595:       PetscMemzero(b->diagv[i],bs*bs*b->bdlen[i]*sizeof(Scalar));
596:     }
597:     b->nonew = 0; b->nonew_diag = 0;
598:   } else { /* diagonals are set on input; don't allow dynamic allocation */
599:     b->nonew = 1; b->nonew_diag = 1;
600:   }

602:   /* adjust diagv so one may access rows with diagv[diag][row] for all rows */
603:   for (i=0; i<nd; i++) {
604:     if (diag[i] > 0) {
605:       b->diagv[i] -= bs*bs*diag[i];
606:     }
607:   }

609:   b->nz          = b->maxnz; /* Currently not keeping track of exact count */
610:   b->roworiented = PETSC_TRUE;
611:   B->info.nz_unneeded = (double)b->maxnz;
612:   return(0);
613: }

615: static int MatDuplicate_SeqBDiag(Mat A,MatDuplicateOption cpvalues,Mat *matout)
616: {
617:   Mat_SeqBDiag *newmat,*a = (Mat_SeqBDiag*)A->data;
618:   int          i,ierr,len,diag,bs = a->bs;
619:   Mat          mat;

622:   MatCreateSeqBDiag(A->comm,A->m,A->n,a->nd,bs,a->diag,PETSC_NULL,matout);

624:   /* Copy contents of diagonals */
625:   mat = *matout;
626:   newmat = (Mat_SeqBDiag*)mat->data;
627:   if (cpvalues == MAT_COPY_VALUES) {
628:     for (i=0; i<a->nd; i++) {
629:       len = a->bdlen[i] * bs * bs * sizeof(Scalar);
630:       diag = a->diag[i];
631:       if (diag > 0) {
632:         PetscMemcpy(newmat->diagv[i]+bs*bs*diag,a->diagv[i]+bs*bs*diag,len);
633:       } else {
634:         PetscMemcpy(newmat->diagv[i],a->diagv[i],len);
635:       }
636:     }
637:   }
638:   MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
639:   MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);
640:   return(0);
641: }

643: EXTERN_C_BEGIN
644: int MatLoad_SeqBDiag(PetscViewer viewer,MatType type,Mat *A)
645: {
646:   Mat        B;
647:   int        *scols,i,nz,ierr,fd,header[4],size,nd = 128;
648:   int        bs,*rowlengths = 0,M,N,*cols,extra_rows,*diag = 0;
649:   int        idiag[128];
650:   Scalar     *vals,*svals;
651:   MPI_Comm   comm;
652:   PetscTruth flg;
653: 
655:   PetscObjectGetComm((PetscObject)viewer,&comm);
656:   MPI_Comm_size(comm,&size);
657:   if (size > 1) SETERRQ(PETSC_ERR_ARG_SIZ,"view must have one processor");
658:   PetscViewerBinaryGetDescriptor(viewer,&fd);
659:   PetscBinaryRead(fd,header,4,PETSC_INT);
660:   if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object");
661:   M = header[1]; N = header[2]; nz = header[3];
662:   if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only load square matrices");
663:   if (header[3] < 0) {
664:     SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqBDiag");
665:   }

667:   /* 
668:      This code adds extra rows to make sure the number of rows is 
669:     divisible by the blocksize
670:   */
671:   bs = 1;
672:   PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);
673:   extra_rows = bs - M + bs*(M/bs);
674:   if (extra_rows == bs) extra_rows = 0;
675:   if (extra_rows) {
676:     PetscLogInfo(0,"MatLoad_SeqBDiag:Padding loaded matrix to match blocksizen");
677:   }

679:   /* read row lengths */
680:   PetscMalloc((M+extra_rows)*sizeof(int),&rowlengths);
681:   PetscBinaryRead(fd,rowlengths,M,PETSC_INT);
682:   for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;

684:   /* load information about diagonals */
685:   PetscOptionsGetIntArray(PETSC_NULL,"-matload_bdiag_diags",idiag,&nd,&flg);
686:   if (flg) {
687:     diag = idiag;
688:   }

690:   /* create our matrix */
691:   MatCreateSeqBDiag(comm,M+extra_rows,M+extra_rows,nd,bs,diag,
692:                            PETSC_NULL,A);
693:   B = *A;

695:   /* read column indices and nonzeros */
696:   PetscMalloc(nz*sizeof(int),&scols);
697:   cols = scols;
698:   PetscBinaryRead(fd,cols,nz,PETSC_INT);
699:   PetscMalloc(nz*sizeof(Scalar),&svals);
700:   vals = svals;
701:   PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
702:   /* insert into matrix */

704:   for (i=0; i<M; i++) {
705:     MatSetValues(B,1,&i,rowlengths[i],scols,svals,INSERT_VALUES);
706:     scols += rowlengths[i]; svals += rowlengths[i];
707:   }
708:   vals[0] = 1.0;
709:   for (i=M; i<M+extra_rows; i++) {
710:     MatSetValues(B,1,&i,1,&i,vals,INSERT_VALUES);
711:   }

713:   PetscFree(cols);
714:   PetscFree(vals);
715:   PetscFree(rowlengths);

717:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
718:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
719:   return(0);
720: }
721: EXTERN_C_END

723: EXTERN_C_BEGIN
724: int MatCreate_SeqBDiag(Mat B)
725: {
726:   Mat_SeqBDiag *b;
727:   int          ierr,size;

730:   MPI_Comm_size(B->comm,&size);
731:   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");

733:   B->m = B->M = PetscMax(B->m,B->M);
734:   B->n = B->N = PetscMax(B->n,B->N);

736:   ierr            = PetscNew(Mat_SeqBDiag,&b);
737:   B->data         = (void*)b;
738:   ierr            = PetscMemzero(b,sizeof(Mat_SeqBDiag));
739:   ierr            = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
740:   B->factor       = 0;
741:   B->mapping      = 0;

743:   MapCreateMPI(B->comm,B->m,B->m,&B->rmap);
744:   MapCreateMPI(B->comm,B->n,B->n,&B->cmap);

746:   b->ndim   = 0;
747:   b->mainbd = -1;
748:   b->pivot  = 0;

750:   b->roworiented = PETSC_TRUE;
751:   return(0);
752: }
753: EXTERN_C_END

755: /*@C
756:    MatCreateSeqBDiag - Creates a sequential block diagonal matrix.

758:    Collective on MPI_Comm

760:    Input Parameters:
761: +  comm - MPI communicator, set to PETSC_COMM_SELF
762: .  m - number of rows
763: .  n - number of columns
764: .  nd - number of block diagonals (optional)
765: .  bs - each element of a diagonal is an bs x bs dense matrix
766: .  diag - optional array of block diagonal numbers (length nd).
767:    For a matrix element A[i,j], where i=row and j=column, the
768:    diagonal number is
769: $     diag = i/bs - j/bs  (integer division)
770:    Set diag=PETSC_NULL on input for PETSc to dynamically allocate memory as 
771:    needed (expensive).
772: -  diagv - pointer to actual diagonals (in same order as diag array), 
773:    if allocated by user.  Otherwise, set diagv=PETSC_NULL on input for PETSc
774:    to control memory allocation.

776:    Output Parameters:
777: .  A - the matrix

779:    Options Database Keys:
780: .  -mat_block_size <bs> - Sets blocksize
781: .  -mat_bdiag_diags <s1,s2,s3,...> - Sets diagonal numbers

783:    Notes:
784:    See the users manual for further details regarding this storage format.

786:    Fortran Note:
787:    Fortran programmers cannot set diagv; this value is ignored.

789:    Level: intermediate

791: .keywords: matrix, block, diagonal, sparse

793: .seealso: MatCreate(), MatCreateMPIBDiag(), MatSetValues()
794: @*/
795: int MatCreateSeqBDiag(MPI_Comm comm,int m,int n,int nd,int bs,int *diag,Scalar **diagv,Mat *A)
796: {

800:   MatCreate(comm,m,n,m,n,A);
801:   MatSetType(*A,MATSEQBDIAG);
802:   MatSeqBDiagSetPreallocation(*A,nd,bs,diag,diagv);
803:   return(0);
804: }