Actual source code: bdiag3.c

  1: /*$Id: bdiag3.c,v 1.33 2001/08/07 03:02:53 balay Exp $*/

  3: /* Block diagonal matrix format */

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

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

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

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

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

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

169: /* 
170:    MatNorm_SeqBDiag_Columns - Computes the column norms of a block diagonal
171:    matrix.  We code this separately from MatNorm_SeqBDiag() so that the
172:    routine can be used for the parallel version as well.
173:  */
176: int MatNorm_SeqBDiag_Columns(Mat A,PetscReal *tmp,int n)
177: {
178:   Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
179:   int          d,i,j,k,nd = a->nd,bs = a->bs,diag,kshift,kloc,len,ierr;
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: int MatNorm_SeqBDiag(Mat A,NormType type,PetscReal *nrm)
232: {
233:   Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
234:   PetscReal    sum = 0.0,*tmp;
235:   int          ierr,d,i,j,k,nd = a->nd,bs = a->bs,diag,kshift,kloc,len;
236:   PetscScalar  *dv;

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

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

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

395: /* ----------------------------------------------------------------*/


400: int MatView_SeqBDiag_Binary(Mat A,PetscViewer viewer)
401: {
402:   Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
403:   int          i,ict,fd,*col_lens,*cval,*col,ierr,nz;
404:   PetscScalar  *anonz,*val;

407:   PetscViewerBinaryGetDescriptor(viewer,&fd);

409:   /* For MATSEQBDIAG format,maxnz = nz */
410:   PetscMalloc((4+A->m)*sizeof(int),&col_lens);
411:   col_lens[0] = MAT_FILE_COOKIE;
412:   col_lens[1] = A->m;
413:   col_lens[2] = A->n;
414:   col_lens[3] = a->maxnz;

416:   /* Should do translation using less memory; this is just a quick initial version */
417:   PetscMalloc((a->maxnz)*sizeof(int),&cval);
418:   PetscMalloc((a->maxnz)*sizeof(PetscScalar),&anonz);

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

431:   /* Store lengths of each row and write (including header) to file */
432:   PetscBinaryWrite(fd,col_lens,4+A->m,PETSC_INT,1);
433:   PetscFree(col_lens);

435:   /* Store column indices (zero start index) */
436:   PetscBinaryWrite(fd,cval,a->maxnz,PETSC_INT,0);

438:   /* Store nonzero values */
439:   PetscBinaryWrite(fd,anonz,a->maxnz,PETSC_SCALAR,0);
440:   return(0);
441: }

445: int MatView_SeqBDiag_ASCII(Mat A,PetscViewer viewer)
446: {
447:   Mat_SeqBDiag      *a = (Mat_SeqBDiag*)A->data;
448:   char              *name;
449:   int               ierr,*col,i,j,len,diag,nr = A->m,bs = a->bs,iprint,nz;
450:   PetscScalar       *val,*dv,zero = 0.0;
451:   PetscViewerFormat format;

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

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

628: static int MatView_SeqBDiag_Draw(Mat A,PetscViewer viewer)
629: {
630:   PetscDraw     draw;
631:   PetscReal     xl,yl,xr,yr,w,h;
632:   int           ierr,nz,*col,i,j,nr = A->m;
633:   PetscTruth    isnull;

636:   PetscViewerDrawGetDraw(viewer,0,&draw);
637:   PetscDrawIsNull(draw,&isnull); if (isnull) return(0);

639:   xr = A->n; yr = A->m; h = yr/10.0; w = xr/10.0;
640:   xr += w; yr += h; xl = -w; yl = -h;
641:   PetscDrawSetCoordinates(draw,xl,yl,xr,yr);

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

663: int MatView_SeqBDiag(Mat A,PetscViewer viewer)
664: {
665:   int        ierr;
666:   PetscTruth isascii,isbinary,isdraw;

669:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
670:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
671:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
672:   if (isascii) {
673:     MatView_SeqBDiag_ASCII(A,viewer);
674:   } else if (isbinary) {
675:     MatView_SeqBDiag_Binary(A,viewer);
676:   } else if (isdraw) {
677:     MatView_SeqBDiag_Draw(A,viewer);
678:   } else {
679:     SETERRQ1(1,"Viewer type %s not supported by BDiag matrices",((PetscObject)viewer)->type_name);
680:   }
681:   return(0);
682: }