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/vec/vecimpl.h
 7:  #include src/inline/ilu.h
 8:  #include petscsys.h

 10: EXTERN int MatSetValues_SeqBDiag_1(Mat,int,int *,int,int *,PetscScalar *,InsertMode);
 11: EXTERN int MatSetValues_SeqBDiag_N(Mat,int,int *,int,int *,PetscScalar *,InsertMode);
 12: EXTERN int MatGetValues_SeqBDiag_1(Mat,int,int *,int,int *,PetscScalar *);
 13: EXTERN int MatGetValues_SeqBDiag_N(Mat,int,int *,int,int *,PetscScalar *);
 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);


 34: int MatGetInfo_SeqBDiag(Mat A,MatInfoType flag,MatInfo *info)
 35: {
 36:   Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;

 39:   info->rows_global       = (double)A->m;
 40:   info->columns_global    = (double)A->n;
 41:   info->rows_local        = (double)A->m;
 42:   info->columns_local     = (double)A->n;
 43:   info->block_size        = a->bs;
 44:   info->nz_allocated      = (double)a->maxnz;
 45:   info->nz_used           = (double)a->nz;
 46:   info->nz_unneeded       = (double)(a->maxnz - a->nz);
 47:   info->assemblies        = (double)A->num_ass;
 48:   info->mallocs           = (double)a->reallocs;
 49:   info->memory            = A->mem;
 50:   info->fill_ratio_given  = 0; /* supports ILU(0) only */
 51:   info->fill_ratio_needed = 0;
 52:   info->factor_mallocs    = 0;
 53:   return(0);
 54: }

 56: /*
 57:      Note: this currently will generate different answers if you request
 58:  all items or a subset. If you request all items it checks if the value is
 59:  nonzero and only includes it if it is nonzero; if you check a subset of items
 60:  it returns a list of all active columns in the row (some which may contain
 61:  a zero)
 62: */
 63: int MatGetRow_SeqBDiag(Mat A,int row,int *nz,int **col,PetscScalar **v)
 64: {
 65:   Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
 66:   int          nd = a->nd,bs = a->bs;
 67:   int          nc = A->n,*diag = a->diag,pcol,shift,i,j,k;

 70:   /* For efficiency, if ((nz) && (col) && (v)) then do all at once */
 71:   if ((nz) && (col) && (v)) {
 72:     *col = a->colloc;
 73:     *v   = a->dvalue;
 74:     k    = 0;
 75:     if (bs == 1) {
 76:       for (j=0; j<nd; j++) {
 77:         pcol = row - diag[j];
 78: #if defined(PETSC_USE_COMPLEX)
 79:         if (pcol > -1 && pcol < nc && PetscAbsScalar((a->diagv[j])[row])) {
 80: #else
 81:         if (pcol > -1 && pcol < nc && (a->diagv[j])[row]) {
 82: #endif
 83:           (*v)[k]   = (a->diagv[j])[row];
 84:           (*col)[k] = pcol;
 85:           k++;
 86:         }
 87:       }
 88:       *nz = k;
 89:     } else {
 90:       shift = (row/bs)*bs*bs + row%bs;
 91:       for (j=0; j<nd; j++) {
 92:         pcol = bs * (row/bs - diag[j]);
 93:         if (pcol > -1 && pcol < nc) {
 94:           for (i=0; i<bs; i++) {
 95: #if defined(PETSC_USE_COMPLEX)
 96:             if (PetscAbsScalar((a->diagv[j])[shift + i*bs])) {
 97: #else
 98:             if ((a->diagv[j])[shift + i*bs]) {
 99: #endif
100:               (*v)[k]   = (a->diagv[j])[shift + i*bs];
101:               (*col)[k] = pcol + i;
102:               k++;
103:             }
104:           }
105:         }
106:       }
107:       *nz = k;
108:     }
109:   } else {
110:     if (bs == 1) {
111:       if (nz) {
112:         k = 0;
113:         for (j=0; j<nd; j++) {
114:           pcol = row - diag[j];
115:           if (pcol > -1 && pcol < nc) k++;
116:         }
117:         *nz = k;
118:       }
119:       if (col) {
120:         *col = a->colloc;
121:         k = 0;
122:         for (j=0; j<nd; j++) {
123:           pcol = row - diag[j];
124:           if (pcol > -1 && pcol < nc) {
125:             (*col)[k] = pcol;  k++;
126:           }
127:         }
128:       }
129:       if (v) {
130:         *v = a->dvalue;
131:         k = 0;
132:         for (j=0; j<nd; j++) {
133:           pcol = row - diag[j];
134:           if (pcol > -1 && pcol < nc) {
135:             (*v)[k] = (a->diagv[j])[row]; k++;
136:           }
137:         }
138:       }
139:     } else {
140:       if (nz) {
141:         k = 0;
142:         for (j=0; j<nd; j++) {
143:           pcol = bs * (row/bs- diag[j]);
144:           if (pcol > -1 && pcol < nc) k += bs;
145:         }
146:         *nz = k;
147:       }
148:       if (col) {
149:         *col = a->colloc;
150:         k = 0;
151:         for (j=0; j<nd; j++) {
152:           pcol = bs * (row/bs - diag[j]);
153:           if (pcol > -1 && pcol < nc) {
154:             for (i=0; i<bs; i++) {
155:               (*col)[k+i] = pcol + i;
156:             }
157:             k += bs;
158:           }
159:         }
160:       }
161:       if (v) {
162:         shift = (row/bs)*bs*bs + row%bs;
163:         *v = a->dvalue;
164:         k = 0;
165:         for (j=0; j<nd; j++) {
166:           pcol = bs * (row/bs - diag[j]);
167:           if (pcol > -1 && pcol < nc) {
168:             for (i=0; i<bs; i++) {
169:              (*v)[k+i] = (a->diagv[j])[shift + i*bs];
170:             }
171:             k += bs;
172:           }
173:         }
174:       }
175:     }
176:   }
177:   return(0);
178: }

180: int MatRestoreRow_SeqBDiag(Mat A,int row,int *ncols,int **cols,PetscScalar **vals)
181: {
183:   /* Work space is allocated during matrix creation and freed
184:      when matrix is destroyed */
185:   return(0);
186: }

188: /* 
189:    MatNorm_SeqBDiag_Columns - Computes the column norms of a block diagonal
190:    matrix.  We code this separately from MatNorm_SeqBDiag() so that the
191:    routine can be used for the parallel version as well.
192:  */
193: int MatNorm_SeqBDiag_Columns(Mat A,PetscReal *tmp,int n)
194: {
195:   Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
196:   int          d,i,j,k,nd = a->nd,bs = a->bs,diag,kshift,kloc,len,ierr;
197:   PetscScalar  *dv;

200:   PetscMemzero(tmp,A->n*sizeof(PetscReal));
201:   if (bs == 1) {
202:     for (d=0; d<nd; d++) {
203:       dv   = a->diagv[d];
204:       diag = a->diag[d];
205:       len  = a->bdlen[d];
206:       if (diag > 0) {        /* lower triangle */
207:         for (i=0; i<len; i++) {
208:           tmp[i] += PetscAbsScalar(dv[i+diag]);
209:         }
210:       } else {        /* upper triangle */
211:         for (i=0; i<len; i++) {
212:           tmp[i-diag] += PetscAbsScalar(dv[i]);
213:         }
214:       }
215:     }
216:   } else {
217:     for (d=0; d<nd; d++) {
218:       dv   = a->diagv[d];
219:       diag = a->diag[d];
220:       len  = a->bdlen[d];

222:       if (diag > 0) {        /* lower triangle */
223:         for (k=0; k<len; k++) {
224:           kloc = k*bs; kshift = kloc*bs + diag*bs;
225:           for (i=0; i<bs; i++) {        /* i = local row */
226:             for (j=0; j<bs; j++) {        /* j = local column */
227:               tmp[kloc + j] += PetscAbsScalar(dv[kshift + j*bs + i]);
228:             }
229:           }
230:         }
231:       } else {        /* upper triangle */
232:         for (k=0; k<len; k++) {
233:           kloc = k*bs; kshift = kloc*bs;
234:           for (i=0; i<bs; i++) {        /* i = local row */
235:             for (j=0; j<bs; j++) {        /* j = local column */
236:               tmp[kloc + j - bs*diag] += PetscAbsScalar(dv[kshift + j*bs + i]);
237:             }
238:           }
239:         }
240:       }
241:     }
242:   }
243:   return(0);
244: }

246: int MatNorm_SeqBDiag(Mat A,NormType type,PetscReal *nrm)
247: {
248:   Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
249:   PetscReal    sum = 0.0,*tmp;
250:   int          ierr,d,i,j,k,nd = a->nd,bs = a->bs,diag,kshift,kloc,len;
251:   PetscScalar  *dv;

254:   if (type == NORM_FROBENIUS) {
255:     for (d=0; d<nd; d++) {
256:       dv   = a->diagv[d];
257:       len  = a->bdlen[d]*bs*bs;
258:       diag = a->diag[d];
259:       if (diag > 0) {
260:         for (i=0; i<len; i++) {
261: #if defined(PETSC_USE_COMPLEX)
262:           sum += PetscRealPart(PetscConj(dv[i+diag])*dv[i+diag]);
263: #else
264:           sum += dv[i+diag]*dv[i+diag];
265: #endif
266:         }
267:       } else {
268:         for (i=0; i<len; i++) {
269: #if defined(PETSC_USE_COMPLEX)
270:           sum += PetscRealPart(PetscConj(dv[i])*dv[i]);
271: #else
272:           sum += dv[i]*dv[i];
273: #endif
274:         }
275:       }
276:     }
277:     *nrm = sqrt(sum);
278:   } else if (type == NORM_1) { /* max column norm */
279:     PetscMalloc((A->n+1)*sizeof(PetscReal),&tmp);
280:     MatNorm_SeqBDiag_Columns(A,tmp,A->n);
281:     *nrm = 0.0;
282:     for (j=0; j<A->n; j++) {
283:       if (tmp[j] > *nrm) *nrm = tmp[j];
284:     }
285:     PetscFree(tmp);
286:   } else if (type == NORM_INFINITY) { /* max row norm */
287:     PetscMalloc((A->m+1)*sizeof(PetscReal),&tmp);
288:     PetscMemzero(tmp,A->m*sizeof(PetscReal));
289:     *nrm = 0.0;
290:     if (bs == 1) {
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) {        /* lower triangle */
296:           for (i=0; i<len; i++) {
297:             tmp[i+diag] += PetscAbsScalar(dv[i+diag]);
298:           }
299:         } else {        /* upper triangle */
300:           for (i=0; i<len; i++) {
301:             tmp[i] += PetscAbsScalar(dv[i]);
302:           }
303:         }
304:       }
305:     } else {
306:       for (d=0; d<nd; d++) {
307:         dv   = a->diagv[d];
308:         diag = a->diag[d];
309:         len  = a->bdlen[d];
310:         if (diag > 0) {
311:           for (k=0; k<len; k++) {
312:             kloc = k*bs; kshift = kloc*bs + bs*diag;
313:             for (i=0; i<bs; i++) {        /* i = local row */
314:               for (j=0; j<bs; j++) {        /* j = local column */
315:                 tmp[kloc + i + bs*diag] += PetscAbsScalar(dv[kshift+j*bs+i]);
316:               }
317:             }
318:           }
319:         } else {
320:           for (k=0; k<len; k++) {
321:             kloc = k*bs; kshift = kloc*bs;
322:             for (i=0; i<bs; i++) {        /* i = local row */
323:               for (j=0; j<bs; j++) {        /* j = local column */
324:                 tmp[kloc + i] += PetscAbsScalar(dv[kshift + j*bs + i]);
325:               }
326:             }
327:           }
328:         }
329:       }
330:     }
331:     for (j=0; j<A->m; j++) {
332:       if (tmp[j] > *nrm) *nrm = tmp[j];
333:     }
334:     PetscFree(tmp);
335:   } else {
336:     SETERRQ(PETSC_ERR_SUP,"No support for two norm");
337:   }
338:   return(0);
339: }

341: int MatTranspose_SeqBDiag(Mat A,Mat *matout)
342: {
343:   Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data,*anew;
344:   Mat          tmat;
345:   int          i,j,k,d,ierr,nd = a->nd,*diag = a->diag,*diagnew;
346:   int          bs = a->bs,kshift,shifto,shiftn;
347:   PetscScalar  *dwork,*dvnew;

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

406: /* ----------------------------------------------------------------*/


409: int MatView_SeqBDiag_Binary(Mat A,PetscViewer viewer)
410: {
411:   Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
412:   int          i,ict,fd,*col_lens,*cval,*col,ierr,nz;
413:   PetscScalar  *anonz,*val;

416:   PetscViewerBinaryGetDescriptor(viewer,&fd);

418:   /* For MATSEQBDIAG format,maxnz = nz */
419:   ierr        = PetscMalloc((4+A->m)*sizeof(int),&col_lens);
420:   col_lens[0] = MAT_FILE_COOKIE;
421:   col_lens[1] = A->m;
422:   col_lens[2] = A->n;
423:   col_lens[3] = a->maxnz;

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

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

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

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

447:   /* Store nonzero values */
448:   PetscBinaryWrite(fd,anonz,a->maxnz,PETSC_SCALAR,0);
449:   return(0);
450: }

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

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

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

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

641:   PetscViewerDrawGetDraw(viewer,0,&draw);
642:   PetscDrawIsNull(draw,&isnull); if (isnull) return(0);

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

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

666: int MatView_SeqBDiag(Mat A,PetscViewer viewer)
667: {
668:   int        ierr;
669:   PetscTruth isascii,isbinary,isdraw;

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