Actual source code: bdiag3.c

  2: /* Block diagonal matrix format */

 4:  #include src/mat/impls/bdiag/seq/bdiag.h
 5:  #include src/inline/ilu.h
 6:  #include petscsys.h

 10: PetscErrorCode MatGetInfo_SeqBDiag(Mat A,MatInfoType flag,MatInfo *info)
 11: {
 12:   Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;

 15:   info->rows_global       = (double)A->m;
 16:   info->columns_global    = (double)A->n;
 17:   info->rows_local        = (double)A->m;
 18:   info->columns_local     = (double)A->n;
 19:   info->block_size        = A->bs;
 20:   info->nz_allocated      = (double)a->maxnz;
 21:   info->nz_used           = (double)a->nz;
 22:   info->nz_unneeded       = (double)(a->maxnz - a->nz);
 23:   info->assemblies        = (double)A->num_ass;
 24:   info->mallocs           = (double)a->reallocs;
 25:   info->memory            = A->mem;
 26:   info->fill_ratio_given  = 0; /* supports ILU(0) only */
 27:   info->fill_ratio_needed = 0;
 28:   info->factor_mallocs    = 0;
 29:   return(0);
 30: }

 32: /*
 33:      Note: this currently will generate different answers if you request
 34:  all items or a subset. If you request all items it checks if the value is
 35:  nonzero and only includes it if it is nonzero; if you check a subset of items
 36:  it returns a list of all active columns in the row (some which may contain
 37:  a zero)
 38: */
 41: PetscErrorCode MatGetRow_SeqBDiag(Mat A,PetscInt row,PetscInt *nz,PetscInt **col,PetscScalar **v)
 42: {
 43:   Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
 44:   PetscInt     nd = a->nd,bs = A->bs;
 45:   PetscInt     nc = A->n,*diag = a->diag,pcol,shift,i,j,k;

 48:   /* For efficiency, if ((nz) && (col) && (v)) then do all at once */
 49:   if ((nz) && (col) && (v)) {
 50:     *col = a->colloc;
 51:     *v   = a->dvalue;
 52:     k    = 0;
 53:     if (bs == 1) {
 54:       for (j=0; j<nd; j++) {
 55:         pcol = row - diag[j];
 56: #if defined(PETSC_USE_COMPLEX)
 57:         if (pcol > -1 && pcol < nc && PetscAbsScalar((a->diagv[j])[row])) {
 58: #else
 59:         if (pcol > -1 && pcol < nc && (a->diagv[j])[row]) {
 60: #endif
 61:           (*v)[k]   = (a->diagv[j])[row];
 62:           (*col)[k] = pcol;
 63:           k++;
 64:         }
 65:       }
 66:       *nz = k;
 67:     } else {
 68:       shift = (row/bs)*bs*bs + row%bs;
 69:       for (j=0; j<nd; j++) {
 70:         pcol = bs * (row/bs - diag[j]);
 71:         if (pcol > -1 && pcol < nc) {
 72:           for (i=0; i<bs; i++) {
 73: #if defined(PETSC_USE_COMPLEX)
 74:             if (PetscAbsScalar((a->diagv[j])[shift + i*bs])) {
 75: #else
 76:             if ((a->diagv[j])[shift + i*bs]) {
 77: #endif
 78:               (*v)[k]   = (a->diagv[j])[shift + i*bs];
 79:               (*col)[k] = pcol + i;
 80:               k++;
 81:             }
 82:           }
 83:         }
 84:       }
 85:       *nz = k;
 86:     }
 87:   } else {
 88:     if (bs == 1) {
 89:       if (nz) {
 90:         k = 0;
 91:         for (j=0; j<nd; j++) {
 92:           pcol = row - diag[j];
 93:           if (pcol > -1 && pcol < nc) k++;
 94:         }
 95:         *nz = k;
 96:       }
 97:       if (col) {
 98:         *col = a->colloc;
 99:         k = 0;
100:         for (j=0; j<nd; j++) {
101:           pcol = row - diag[j];
102:           if (pcol > -1 && pcol < nc) {
103:             (*col)[k] = pcol;  k++;
104:           }
105:         }
106:       }
107:       if (v) {
108:         *v = a->dvalue;
109:         k = 0;
110:         for (j=0; j<nd; j++) {
111:           pcol = row - diag[j];
112:           if (pcol > -1 && pcol < nc) {
113:             (*v)[k] = (a->diagv[j])[row]; k++;
114:           }
115:         }
116:       }
117:     } else {
118:       if (nz) {
119:         k = 0;
120:         for (j=0; j<nd; j++) {
121:           pcol = bs * (row/bs- diag[j]);
122:           if (pcol > -1 && pcol < nc) k += bs;
123:         }
124:         *nz = k;
125:       }
126:       if (col) {
127:         *col = a->colloc;
128:         k = 0;
129:         for (j=0; j<nd; j++) {
130:           pcol = bs * (row/bs - diag[j]);
131:           if (pcol > -1 && pcol < nc) {
132:             for (i=0; i<bs; i++) {
133:               (*col)[k+i] = pcol + i;
134:             }
135:             k += bs;
136:           }
137:         }
138:       }
139:       if (v) {
140:         shift = (row/bs)*bs*bs + row%bs;
141:         *v = a->dvalue;
142:         k = 0;
143:         for (j=0; j<nd; j++) {
144:           pcol = bs * (row/bs - diag[j]);
145:           if (pcol > -1 && pcol < nc) {
146:             for (i=0; i<bs; i++) {
147:              (*v)[k+i] = (a->diagv[j])[shift + i*bs];
148:             }
149:             k += bs;
150:           }
151:         }
152:       }
153:     }
154:   }
155:   return(0);
156: }

160: PetscErrorCode MatRestoreRow_SeqBDiag(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals)
161: {
163:   /* Work space is allocated during matrix creation and freed
164:      when matrix is destroyed */
165:   return(0);
166: }

168: /* 
169:    MatNorm_SeqBDiag_Columns - Computes the column norms of a block diagonal
170:    matrix.  We code this separately from MatNorm_SeqBDiag() so that the
171:    routine can be used for the parallel version as well.
172:  */
175: PetscErrorCode MatNorm_SeqBDiag_Columns(Mat A,PetscReal *tmp,PetscInt n)
176: {
177:   Mat_SeqBDiag   *a = (Mat_SeqBDiag*)A->data;
179:   PetscInt       d,i,j,k,nd = a->nd,bs = A->bs,diag,kshift,kloc,len;
180:   PetscScalar    *dv;

183:   PetscMemzero(tmp,A->n*sizeof(PetscReal));
184:   if (bs == 1) {
185:     for (d=0; d<nd; d++) {
186:       dv   = a->diagv[d];
187:       diag = a->diag[d];
188:       len  = a->bdlen[d];
189:       if (diag > 0) {        /* lower triangle */
190:         for (i=0; i<len; i++) {
191:           tmp[i] += PetscAbsScalar(dv[i+diag]);
192:         }
193:       } else {        /* upper triangle */
194:         for (i=0; i<len; i++) {
195:           tmp[i-diag] += PetscAbsScalar(dv[i]);
196:         }
197:       }
198:     }
199:   } else {
200:     for (d=0; d<nd; d++) {
201:       dv   = a->diagv[d];
202:       diag = a->diag[d];
203:       len  = a->bdlen[d];

205:       if (diag > 0) {        /* lower triangle */
206:         for (k=0; k<len; k++) {
207:           kloc = k*bs; kshift = kloc*bs + diag*bs;
208:           for (i=0; i<bs; i++) {        /* i = local row */
209:             for (j=0; j<bs; j++) {        /* j = local column */
210:               tmp[kloc + j] += PetscAbsScalar(dv[kshift + j*bs + i]);
211:             }
212:           }
213:         }
214:       } else {        /* upper triangle */
215:         for (k=0; k<len; k++) {
216:           kloc = k*bs; kshift = kloc*bs;
217:           for (i=0; i<bs; i++) {        /* i = local row */
218:             for (j=0; j<bs; j++) {        /* j = local column */
219:               tmp[kloc + j - bs*diag] += PetscAbsScalar(dv[kshift + j*bs + i]);
220:             }
221:           }
222:         }
223:       }
224:     }
225:   }
226:   return(0);
227: }

231: PetscErrorCode MatNorm_SeqBDiag(Mat A,NormType type,PetscReal *nrm)
232: {
233:   Mat_SeqBDiag   *a = (Mat_SeqBDiag*)A->data;
234:   PetscReal      sum = 0.0,*tmp;
236:   PetscInt       d,i,j,k,nd = a->nd,bs = A->bs,diag,kshift,kloc,len;
237:   PetscScalar    *dv;

240:   if (type == NORM_FROBENIUS) {
241:     for (d=0; d<nd; d++) {
242:       dv   = a->diagv[d];
243:       len  = a->bdlen[d]*bs*bs;
244:       diag = a->diag[d];
245:       if (diag > 0) {
246:         for (i=0; i<len; i++) {
247: #if defined(PETSC_USE_COMPLEX)
248:           sum += PetscRealPart(PetscConj(dv[i+diag])*dv[i+diag]);
249: #else
250:           sum += dv[i+diag]*dv[i+diag];
251: #endif
252:         }
253:       } else {
254:         for (i=0; i<len; i++) {
255: #if defined(PETSC_USE_COMPLEX)
256:           sum += PetscRealPart(PetscConj(dv[i])*dv[i]);
257: #else
258:           sum += dv[i]*dv[i];
259: #endif
260:         }
261:       }
262:     }
263:     *nrm = sqrt(sum);
264:   } else if (type == NORM_1) { /* max column norm */
265:     PetscMalloc((A->n+1)*sizeof(PetscReal),&tmp);
266:     MatNorm_SeqBDiag_Columns(A,tmp,A->n);
267:     *nrm = 0.0;
268:     for (j=0; j<A->n; j++) {
269:       if (tmp[j] > *nrm) *nrm = tmp[j];
270:     }
271:     PetscFree(tmp);
272:   } else if (type == NORM_INFINITY) { /* max row norm */
273:     PetscMalloc((A->m+1)*sizeof(PetscReal),&tmp);
274:     PetscMemzero(tmp,A->m*sizeof(PetscReal));
275:     *nrm = 0.0;
276:     if (bs == 1) {
277:       for (d=0; d<nd; d++) {
278:         dv   = a->diagv[d];
279:         diag = a->diag[d];
280:         len  = a->bdlen[d];
281:         if (diag > 0) {        /* lower triangle */
282:           for (i=0; i<len; i++) {
283:             tmp[i+diag] += PetscAbsScalar(dv[i+diag]);
284:           }
285:         } else {        /* upper triangle */
286:           for (i=0; i<len; i++) {
287:             tmp[i] += PetscAbsScalar(dv[i]);
288:           }
289:         }
290:       }
291:     } else {
292:       for (d=0; d<nd; d++) {
293:         dv   = a->diagv[d];
294:         diag = a->diag[d];
295:         len  = a->bdlen[d];
296:         if (diag > 0) {
297:           for (k=0; k<len; k++) {
298:             kloc = k*bs; kshift = kloc*bs + bs*diag;
299:             for (i=0; i<bs; i++) {        /* i = local row */
300:               for (j=0; j<bs; j++) {        /* j = local column */
301:                 tmp[kloc + i + bs*diag] += PetscAbsScalar(dv[kshift+j*bs+i]);
302:               }
303:             }
304:           }
305:         } else {
306:           for (k=0; k<len; k++) {
307:             kloc = k*bs; kshift = kloc*bs;
308:             for (i=0; i<bs; i++) {        /* i = local row */
309:               for (j=0; j<bs; j++) {        /* j = local column */
310:                 tmp[kloc + i] += PetscAbsScalar(dv[kshift + j*bs + i]);
311:               }
312:             }
313:           }
314:         }
315:       }
316:     }
317:     for (j=0; j<A->m; j++) {
318:       if (tmp[j] > *nrm) *nrm = tmp[j];
319:     }
320:     PetscFree(tmp);
321:   } else {
322:     SETERRQ(PETSC_ERR_SUP,"No support for two norm");
323:   }
324:   return(0);
325: }

329: PetscErrorCode MatTranspose_SeqBDiag(Mat A,Mat *matout)
330: {
331:   Mat_SeqBDiag   *a = (Mat_SeqBDiag*)A->data,*anew;
332:   Mat            tmat;
334:   PetscInt       i,j,k,d,nd = a->nd,*diag = a->diag,*diagnew;
335:   PetscInt       bs = A->bs,kshift,shifto,shiftn;
336:   PetscScalar    *dwork,*dvnew;

339:   PetscMalloc((nd+1)*sizeof(PetscInt),&diagnew);
340:   for (i=0; i<nd; i++) {
341:     diagnew[i] = -diag[nd-i-1]; /* assume sorted in descending order */
342:   }
343:   MatCreate(A->comm,A->n,A->m,A->n,A->m,&tmat);
344:   MatSetType(tmat,A->type_name);
345:   MatSeqBDiagSetPreallocation(tmat,nd,bs,diagnew,PETSC_NULL);
346:   PetscFree(diagnew);
347:   anew = (Mat_SeqBDiag*)tmat->data;
348:   for (d=0; d<nd; d++) {
349:     dvnew = anew->diagv[d];
350:     dwork = a->diagv[nd-d-1];
351:     if (anew->bdlen[d] != a->bdlen[nd-d-1]) SETERRQ(PETSC_ERR_ARG_SIZ,"Incompatible diagonal lengths");
352:     shifto = a->diag[nd-d-1];
353:     shiftn = anew->diag[d];
354:     if (shifto > 0)  shifto = bs*bs*shifto; else shifto = 0;
355:     if (shiftn > 0)  shiftn = bs*bs*shiftn; else shiftn = 0;
356:     if (bs == 1) {
357:       for (k=0; k<anew->bdlen[d]; k++) dvnew[shiftn+k] = dwork[shifto+k];
358:     } else {
359:       for (k=0; k<anew->bdlen[d]; k++) {
360:         kshift = k*bs*bs;
361:         for (i=0; i<bs; i++) {        /* i = local row */
362:           for (j=0; j<bs; j++) {        /* j = local column */
363:             dvnew[shiftn + kshift + j + i*bs] = dwork[shifto + kshift + j*bs + i];
364:           }
365:         }
366:       }
367:     }
368:   }
369:   MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);
370:   MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);
371:   if (matout) {
372:     *matout = tmat;
373:   } else {
374:     /* This isn't really an in-place transpose ... but free data 
375:        structures from a.  We should fix this. */
376:     if (!a->user_alloc) { /* Free the actual diagonals */
377:       for (i=0; i<a->nd; i++) {
378:         if (a->diag[i] > 0) {
379:           PetscScalar *dummy = a->diagv[i] + bs*bs*a->diag[i];
380:           PetscFree(dummy);
381:         } else {
382:           PetscFree(a->diagv[i]);
383:         }
384:       }
385:     }
386:     if (a->pivot) {PetscFree(a->pivot);}
387:     PetscFree(a->diagv);
388:     PetscFree(a->diag);
389:     PetscFree(a->colloc);
390:     PetscFree(a->dvalue);
391:     PetscFree(a);
392:     PetscMemcpy(A,tmat,sizeof(struct _p_Mat));
393:     PetscHeaderDestroy(tmat);
394:   }
395:   return(0);
396: }

398: /* ----------------------------------------------------------------*/


403: PetscErrorCode MatView_SeqBDiag_Binary(Mat A,PetscViewer viewer)
404: {
405:   Mat_SeqBDiag   *a = (Mat_SeqBDiag*)A->data;
407:   PetscInt       i,ict,*col_lens,*cval,*col,nz;
408:   PetscScalar    *anonz,*val;
409:   int            fd;

412:   PetscViewerBinaryGetDescriptor(viewer,&fd);

414:   /* For MATSEQBDIAG format,maxnz = nz */
415:   PetscMalloc((4+A->m)*sizeof(PetscInt),&col_lens);
416:   col_lens[0] = MAT_FILE_COOKIE;
417:   col_lens[1] = A->m;
418:   col_lens[2] = A->n;
419:   col_lens[3] = a->maxnz;

421:   /* Should do translation using less memory; this is just a quick initial version */
422:   PetscMalloc((a->maxnz)*sizeof(PetscInt),&cval);
423:   PetscMalloc((a->maxnz)*sizeof(PetscScalar),&anonz);

425:   ict = 0;
426:   for (i=0; i<A->m; i++) {
427:     MatGetRow_SeqBDiag(A,i,&nz,&col,&val);
428:     col_lens[4+i] = nz;
429:     PetscMemcpy(&cval[ict],col,nz*sizeof(PetscInt));
430:     PetscMemcpy(&anonz[ict],anonz,nz*sizeof(PetscScalar));
431:     MatRestoreRow_SeqBDiag(A,i,&nz,&col,&val);
432:     ict += nz;
433:   }
434:   if (ict != a->maxnz) SETERRQ(PETSC_ERR_PLIB,"Error in nonzero count");

436:   /* Store lengths of each row and write (including header) to file */
437:   PetscBinaryWrite(fd,col_lens,4+A->m,PETSC_INT,PETSC_TRUE);
438:   PetscFree(col_lens);

440:   /* Store column indices (zero start index) */
441:   PetscBinaryWrite(fd,cval,a->maxnz,PETSC_INT,PETSC_FALSE);

443:   /* Store nonzero values */
444:   PetscBinaryWrite(fd,anonz,a->maxnz,PETSC_SCALAR,PETSC_FALSE);
445:   return(0);
446: }

450: PetscErrorCode MatView_SeqBDiag_ASCII(Mat A,PetscViewer viewer)
451: {
452:   Mat_SeqBDiag      *a = (Mat_SeqBDiag*)A->data;
453:   char              *name;
454:   PetscErrorCode    ierr;
455:   PetscInt          *col,i,j,len,diag,nr = A->m,bs = A->bs,iprint,nz;
456:   PetscScalar       *val,*dv,zero = 0.0;
457:   PetscViewerFormat format;

460:   PetscObjectGetName((PetscObject)A,&name);
461:   PetscViewerGetFormat(viewer,&format);
462:   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
463:     PetscInt nline = PetscMin(10,a->nd),k,nk,np;
464:     if (a->user_alloc) {
465:       PetscViewerASCIIPrintf(viewer,"block size=%D, number of diagonals=%D, user-allocated storage\n",bs,a->nd);
466:     } else {
467:       PetscViewerASCIIPrintf(viewer,"block size=%D, number of diagonals=%D, PETSc-allocated storage\n",bs,a->nd);
468:     }
469:     nk = (a->nd-1)/nline + 1;
470:     for (k=0; k<nk; k++) {
471:       PetscViewerASCIIPrintf(viewer,"diag numbers:");
472:       np = PetscMin(nline,a->nd - nline*k);
473:       PetscViewerASCIIUseTabs(viewer,PETSC_NO);
474:       for (i=0; i<np; i++) {
475:         PetscViewerASCIIPrintf(viewer,"  %D",a->diag[i+nline*k]);
476:       }
477:       PetscViewerASCIIPrintf(viewer,"\n");
478:       PetscViewerASCIIUseTabs(viewer,PETSC_YES);
479:     }
480:   } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
481:     PetscViewerASCIIUseTabs(viewer,PETSC_NO);
482:     PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",nr,A->n);
483:     PetscViewerASCIIPrintf(viewer,"%% Nonzeros = %D \n",a->nz);
484:     PetscViewerASCIIPrintf(viewer,"zzz = zeros(%D,3);\n",a->nz);
485:     PetscViewerASCIIPrintf(viewer,"zzz = [\n");
486:     for (i=0; i<A->m; i++) {
487:       MatGetRow_SeqBDiag(A,i,&nz,&col,&val);
488:       for (j=0; j<nz; j++) {
489:         if (val[j] != zero) {
490: #if defined(PETSC_USE_COMPLEX)
491:           PetscViewerASCIIPrintf(viewer,"%D %D  %18.16e  %18.16ei \n",
492:              i+1,col[j]+1,PetscRealPart(val[j]),PetscImaginaryPart(val[j]));
493: #else
494:           PetscViewerASCIIPrintf(viewer,"%D %D  %18.16ei \n",i+1,col[j]+1,val[j]);
495: #endif
496:         }
497:       }
498:       MatRestoreRow_SeqBDiag(A,i,&nz,&col,&val);
499:     }
500:     PetscViewerASCIIPrintf(viewer,"];\n %s = spconvert(zzz);\n",name);
501:     PetscViewerASCIIUseTabs(viewer,PETSC_YES);
502:   } else if (format == PETSC_VIEWER_ASCII_IMPL) {
503:     PetscViewerASCIIUseTabs(viewer,PETSC_NO);
504:     if (bs == 1) { /* diagonal format */
505:       for (i=0; i<a->nd; i++) {
506:         dv   = a->diagv[i];
507:         diag = a->diag[i];
508:         PetscViewerASCIIPrintf(viewer,"\n<diagonal %D>\n",diag);
509:         /* diag[i] is (row-col)/bs */
510:         if (diag > 0) {  /* lower triangle */
511:           len  = a->bdlen[i];
512:           for (j=0; j<len; j++) {
513:             if (dv[diag+j] != zero) {
514: #if defined(PETSC_USE_COMPLEX)
515:               if (PetscImaginaryPart(dv[diag+j]) != 0.0) {
516:                 PetscViewerASCIIPrintf(viewer,"A[ %D , %D ] = %e + %e i\n",
517:                           j+diag,j,PetscRealPart(dv[diag+j]),PetscImaginaryPart(dv[diag+j]));
518:               } else {
519:                 PetscViewerASCIIPrintf(viewer,"A[ %D , %D ] = %e\n",j+diag,j,PetscRealPart(dv[diag+j]));
520:               }
521: #else
522:               PetscViewerASCIIPrintf(viewer,"A[ %D , %D ] = %e\n",j+diag,j,dv[diag+j]);

524: #endif
525:             }
526:           }
527:         } else {         /* upper triangle, including main diagonal */
528:           len  = a->bdlen[i];
529:           for (j=0; j<len; j++) {
530:             if (dv[j] != zero) {
531: #if defined(PETSC_USE_COMPLEX)
532:               if (PetscImaginaryPart(dv[j]) != 0.0) {
533:                 PetscViewerASCIIPrintf(viewer,"A[ %D , %D ] = %g + %g i\n",
534:                                          j,j-diag,PetscRealPart(dv[j]),PetscImaginaryPart(dv[j]));
535:               } else {
536:                 PetscViewerASCIIPrintf(viewer,"A[ %D , %D ] = %g\n",j,j-diag,PetscRealPart(dv[j]));
537:               }
538: #else
539:               PetscViewerASCIIPrintf(viewer,"A[ %D , %D ] = %g\n",j,j-diag,dv[j]);
540: #endif
541:             }
542:           }
543:         }
544:       }
545:     } else {  /* Block diagonals */
546:       PetscInt d,k,kshift;
547:       for (d=0; d< a->nd; d++) {
548:         dv   = a->diagv[d];
549:         diag = a->diag[d];
550:         len  = a->bdlen[d];
551:         PetscViewerASCIIPrintf(viewer,"\n<diagonal %D>\n", diag);
552:         if (diag > 0) {                /* lower triangle */
553:           for (k=0; k<len; k++) {
554:             kshift = (diag+k)*bs*bs;
555:             for (i=0; i<bs; i++) {
556:               iprint = 0;
557:               for (j=0; j<bs; j++) {
558:                 if (dv[kshift + j*bs + i] != zero) {
559:                   iprint = 1;
560: #if defined(PETSC_USE_COMPLEX)
561:                   if (PetscImaginaryPart(dv[kshift + j*bs + i])){
562:                     PetscViewerASCIIPrintf(viewer,"A[%D,%D]=%5.2e + %5.2e i  ",(k+diag)*bs+i,k*bs+j,
563:                       PetscRealPart(dv[kshift + j*bs + i]),PetscImaginaryPart(dv[kshift + j*bs + i]));
564:                   } else {
565:                     PetscViewerASCIIPrintf(viewer,"A[%D,%D]=%5.2e   ",(k+diag)*bs+i,k*bs+j,
566:                       PetscRealPart(dv[kshift + j*bs + i]));
567:                   }
568: #else
569:                   PetscViewerASCIIPrintf(viewer,"A[%D,%D]=%5.2e   ",(k+diag)*bs+i,k*bs+j,
570:                       dv[kshift + j*bs + i]);
571: #endif
572:                 }
573:               }
574:               if (iprint) {PetscViewerASCIIPrintf(viewer,"\n");}
575:             }
576:           }
577:         } else {                /* upper triangle, including main diagonal */
578:           for (k=0; k<len; k++) {
579:             kshift = k*bs*bs;
580:             for (i=0; i<bs; i++) {
581:               iprint = 0;
582:               for (j=0; j<bs; j++) {
583:                 if (dv[kshift + j*bs + i] != zero) {
584:                   iprint = 1;
585: #if defined(PETSC_USE_COMPLEX)
586:                   if (PetscImaginaryPart(dv[kshift + j*bs + i])){
587:                     PetscViewerASCIIPrintf(viewer,"A[%D,%D]=%5.2e + %5.2e i  ",k*bs+i,(k-diag)*bs+j,
588:                        PetscRealPart(dv[kshift + j*bs + i]),PetscImaginaryPart(dv[kshift + j*bs + i]));
589:                   } else {
590:                     PetscViewerASCIIPrintf(viewer,"A[%D,%D]=%5.2e   ",k*bs+i,(k-diag)*bs+j,
591:                        PetscRealPart(dv[kshift + j*bs + i]));
592:                   }
593: #else
594:                   PetscViewerASCIIPrintf(viewer,"A[%D,%D]=%5.2e   ",k*bs+i,(k-diag)*bs+j,
595:                      dv[kshift + j*bs + i]);
596: #endif
597:                 }
598:               }
599:               if (iprint) {PetscViewerASCIIPrintf(viewer,"\n");}
600:             }
601:           }
602:         }
603:       }
604:     }
605:     PetscViewerASCIIUseTabs(viewer,PETSC_YES);
606:   } else {
607:     PetscViewerASCIIUseTabs(viewer,PETSC_NO);
608:     /* the usual row format (PETSC_VIEWER_ASCII_NONZERO_ONLY) */
609:     for (i=0; i<A->m; i++) {
610:       PetscViewerASCIIPrintf(viewer,"row %D:",i);
611:       MatGetRow_SeqBDiag(A,i,&nz,&col,&val);
612:       for (j=0; j<nz; j++) {
613: #if defined(PETSC_USE_COMPLEX)
614:         if (PetscImaginaryPart(val[j]) != 0.0 && PetscRealPart(val[j]) != 0.0) {
615:           PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",col[j],PetscRealPart(val[j]),PetscImaginaryPart(val[j]));
616:         } else if (PetscRealPart(val[j]) != 0.0) {
617:           PetscViewerASCIIPrintf(viewer," (%D, %g) ",col[j],PetscRealPart(val[j]));
618:         }
619: #else
620:         if (val[j] != 0.0) {PetscViewerASCIIPrintf(viewer," (%D, %g) ",col[j],val[j]);}
621: #endif
622:       }
623:       PetscViewerASCIIPrintf(viewer,"\n");
624:       MatRestoreRow_SeqBDiag(A,i,&nz,&col,&val);
625:     }
626:     PetscViewerASCIIUseTabs(viewer,PETSC_YES);
627:   }
628:   PetscViewerFlush(viewer);
629:   return(0);
630: }

634: static PetscErrorCode MatView_SeqBDiag_Draw(Mat A,PetscViewer viewer)
635: {
636:   PetscDraw      draw;
637:   PetscReal      xl,yl,xr,yr,w,h;
639:   PetscInt       nz,*col,i,j,nr = A->m;
640:   PetscTruth     isnull;

643:   PetscViewerDrawGetDraw(viewer,0,&draw);
644:   PetscDrawIsNull(draw,&isnull); if (isnull) return(0);

646:   xr = A->n; yr = A->m; h = yr/10.0; w = xr/10.0;
647:   xr += w; yr += h; xl = -w; yl = -h;
648:   PetscDrawSetCoordinates(draw,xl,yl,xr,yr);

650:   /* loop over matrix elements drawing boxes; we really should do this
651:      by diagonals.  What do we really want to draw here: nonzeros, 
652:      allocated space? */
653:   for (i=0; i<nr; i++) {
654:     yl = nr - i - 1.0; yr = yl + 1.0;
655:     MatGetRow_SeqBDiag(A,i,&nz,&col,0);
656:     for (j=0; j<nz; j++) {
657:       xl = col[j]; xr = xl + 1.0;
658:       PetscDrawRectangle(draw,xl,yl,xr,yr,PETSC_DRAW_BLACK,PETSC_DRAW_BLACK,
659:                            PETSC_DRAW_BLACK,PETSC_DRAW_BLACK);
660:     }
661:     MatRestoreRow_SeqBDiag(A,i,&nz,&col,0);
662:   }
663:   PetscDrawFlush(draw);
664:   PetscDrawPause(draw);
665:   return(0);
666: }

670: PetscErrorCode MatView_SeqBDiag(Mat A,PetscViewer viewer)
671: {
673:   PetscTruth     iascii,isbinary,isdraw;

676:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
677:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
678:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
679:   if (iascii) {
680:     MatView_SeqBDiag_ASCII(A,viewer);
681:   } else if (isbinary) {
682:     MatView_SeqBDiag_Binary(A,viewer);
683:   } else if (isdraw) {
684:     MatView_SeqBDiag_Draw(A,viewer);
685:   } else {
686:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by BDiag matrices",((PetscObject)viewer)->type_name);
687:   }
688:   return(0);
689: }