Actual source code: bdiag3.c

  1: /*$Id: bdiag3.c,v 1.31 2001/04/10 19:35:32 bsmith 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);


 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: int MatGetOwnershipRange_SeqBDiag(Mat A,int *m,int *n)
 57: {
 59:   if (m) *m = 0;
 60:   if (n) *n = A->m;
 61:   return(0);
 62: }

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

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

188: int MatRestoreRow_SeqBDiag(Mat A,int row,int *ncols,int **cols,Scalar **vals)
189: {
191:   /* Work space is allocated during matrix creation and freed
192:      when matrix is destroyed */
193:   return(0);
194: }

196: /* 
197:    MatNorm_SeqBDiag_Columns - Computes the column norms of a block diagonal
198:    matrix.  We code this separately from MatNorm_SeqBDiag() so that the
199:    routine can be used for the parallel version as well.
200:  */
201: int MatNorm_SeqBDiag_Columns(Mat A,PetscReal *tmp,int n)
202: {
203:   Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
204:   int          d,i,j,k,nd = a->nd,bs = a->bs,diag,kshift,kloc,len,ierr;
205:   Scalar       *dv;

208:   PetscMemzero(tmp,A->n*sizeof(PetscReal));
209:   if (bs == 1) {
210:     for (d=0; d<nd; d++) {
211:       dv   = a->diagv[d];
212:       diag = a->diag[d];
213:       len  = a->bdlen[d];
214:       if (diag > 0) {        /* lower triangle */
215:         for (i=0; i<len; i++) {
216:           tmp[i] += PetscAbsScalar(dv[i+diag]);
217:         }
218:       } else {        /* upper triangle */
219:         for (i=0; i<len; i++) {
220:           tmp[i-diag] += PetscAbsScalar(dv[i]);
221:         }
222:       }
223:     }
224:   } else {
225:     for (d=0; d<nd; d++) {
226:       dv   = a->diagv[d];
227:       diag = a->diag[d];
228:       len  = a->bdlen[d];

230:       if (diag > 0) {        /* lower triangle */
231:         for (k=0; k<len; k++) {
232:           kloc = k*bs; kshift = kloc*bs + diag*bs;
233:           for (i=0; i<bs; i++) {        /* i = local row */
234:             for (j=0; j<bs; j++) {        /* j = local column */
235:               tmp[kloc + j] += PetscAbsScalar(dv[kshift + j*bs + i]);
236:             }
237:           }
238:         }
239:       } else {        /* upper triangle */
240:         for (k=0; k<len; k++) {
241:           kloc = k*bs; kshift = kloc*bs;
242:           for (i=0; i<bs; i++) {        /* i = local row */
243:             for (j=0; j<bs; j++) {        /* j = local column */
244:               tmp[kloc + j - bs*diag] += PetscAbsScalar(dv[kshift + j*bs + i]);
245:             }
246:           }
247:         }
248:       }
249:     }
250:   }
251:   return(0);
252: }

254: int MatNorm_SeqBDiag(Mat A,NormType type,PetscReal *norm)
255: {
256:   Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
257:   PetscReal    sum = 0.0,*tmp;
258:   int          ierr,d,i,j,k,nd = a->nd,bs = a->bs,diag,kshift,kloc,len;
259:   Scalar       *dv;

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

349: int MatTranspose_SeqBDiag(Mat A,Mat *matout)
350: {
351:   Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data,*anew;
352:   Mat          tmat;
353:   int          i,j,k,d,ierr,nd = a->nd,*diag = a->diag,*diagnew;
354:   int          bs = a->bs,kshift,shifto,shiftn;
355:   Scalar       *dwork,*dvnew;

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

414: /* ----------------------------------------------------------------*/


417: int MatView_SeqBDiag_Binary(Mat A,PetscViewer viewer)
418: {
419:   Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
420:   int          i,ict,fd,*col_lens,*cval,*col,ierr,nz;
421:   Scalar       *anonz,*val;

424:   PetscViewerBinaryGetDescriptor(viewer,&fd);

426:   /* For MATSEQBDIAG format,maxnz = nz */
427:   ierr        = PetscMalloc((4+A->m)*sizeof(int),&col_lens);
428:   col_lens[0] = MAT_COOKIE;
429:   col_lens[1] = A->m;
430:   col_lens[2] = A->n;
431:   col_lens[3] = a->maxnz;

433:   /* Should do translation using less memory; this is just a quick initial version */
434:   PetscMalloc((a->maxnz)*sizeof(int),&cval);
435:   PetscMalloc((a->maxnz)*sizeof(Scalar),&anonz);

437:   ict = 0;
438:   for (i=0; i<A->m; i++) {
439:     MatGetRow_SeqBDiag(A,i,&nz,&col,&val);
440:     col_lens[4+i] = nz;
441:     PetscMemcpy(&cval[ict],col,nz*sizeof(int));
442:     PetscMemcpy(&anonz[ict],anonz,nz*sizeof(Scalar));
443:     MatRestoreRow_SeqBDiag(A,i,&nz,&col,&val);
444:     ict += nz;
445:   }
446:   if (ict != a->maxnz) SETERRQ(PETSC_ERR_PLIB,"Error in nonzero count");

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

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

455:   /* Store nonzero values */
456:   PetscBinaryWrite(fd,anonz,a->maxnz,PETSC_SCALAR,0);
457:   return(0);
458: }

460: int MatView_SeqBDiag_ASCII(Mat A,PetscViewer viewer)
461: {
462:   Mat_SeqBDiag      *a = (Mat_SeqBDiag*)A->data;
463:   char              *name;
464:   int               ierr,*col,i,j,len,diag,nr = A->m,bs = a->bs,iprint,nz;
465:   Scalar            *val,*dv,zero = 0.0;
466:   PetscViewerFormat format;

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

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

641: static int MatView_SeqBDiag_Draw(Mat A,PetscViewer viewer)
642: {
643:   PetscDraw     draw;
644:   PetscReal     xl,yl,xr,yr,w,h;
645:   int           ierr,nz,*col,i,j,nr = A->m;
646:   PetscTruth    isnull;

649:   PetscViewerDrawGetDraw(viewer,0,&draw);
650:   PetscDrawIsNull(draw,&isnull); if (isnull) return(0);

652:   xr = A->n; yr = A->m; h = yr/10.0; w = xr/10.0;
653:   xr += w; yr += h; xl = -w; yl = -h;
654:   PetscDrawSetCoordinates(draw,xl,yl,xr,yr);

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

674: int MatView_SeqBDiag(Mat A,PetscViewer viewer)
675: {
676:   int        ierr;
677:   PetscTruth isascii,isbinary,isdraw;

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