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
12: int MatGetInfo_SeqBDiag(Mat A,MatInfoType flag,MatInfo *info)
13: {
14: Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
17: info->rows_global = (double)A->m;
18: info->columns_global = (double)A->n;
19: info->rows_local = (double)A->m;
20: info->columns_local = (double)A->n;
21: info->block_size = a->bs;
22: info->nz_allocated = (double)a->maxnz;
23: info->nz_used = (double)a->nz;
24: info->nz_unneeded = (double)(a->maxnz - a->nz);
25: info->assemblies = (double)A->num_ass;
26: info->mallocs = (double)a->reallocs;
27: info->memory = A->mem;
28: info->fill_ratio_given = 0; /* supports ILU(0) only */
29: info->fill_ratio_needed = 0;
30: info->factor_mallocs = 0;
31: return(0);
32: }
34: /*
35: Note: this currently will generate different answers if you request
36: all items or a subset. If you request all items it checks if the value is
37: nonzero and only includes it if it is nonzero; if you check a subset of items
38: it returns a list of all active columns in the row (some which may contain
39: a zero)
40: */
43: int MatGetRow_SeqBDiag(Mat A,int row,int *nz,int **col,PetscScalar **v)
44: {
45: Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
46: int nd = a->nd,bs = a->bs;
47: int nc = A->n,*diag = a->diag,pcol,shift,i,j,k;
50: /* For efficiency, if ((nz) && (col) && (v)) then do all at once */
51: if ((nz) && (col) && (v)) {
52: *col = a->colloc;
53: *v = a->dvalue;
54: k = 0;
55: if (bs == 1) {
56: for (j=0; j<nd; j++) {
57: pcol = row - diag[j];
58: #if defined(PETSC_USE_COMPLEX)
59: if (pcol > -1 && pcol < nc && PetscAbsScalar((a->diagv[j])[row])) {
60: #else
61: if (pcol > -1 && pcol < nc && (a->diagv[j])[row]) {
62: #endif
63: (*v)[k] = (a->diagv[j])[row];
64: (*col)[k] = pcol;
65: k++;
66: }
67: }
68: *nz = k;
69: } else {
70: shift = (row/bs)*bs*bs + row%bs;
71: for (j=0; j<nd; j++) {
72: pcol = bs * (row/bs - diag[j]);
73: if (pcol > -1 && pcol < nc) {
74: for (i=0; i<bs; i++) {
75: #if defined(PETSC_USE_COMPLEX)
76: if (PetscAbsScalar((a->diagv[j])[shift + i*bs])) {
77: #else
78: if ((a->diagv[j])[shift + i*bs]) {
79: #endif
80: (*v)[k] = (a->diagv[j])[shift + i*bs];
81: (*col)[k] = pcol + i;
82: k++;
83: }
84: }
85: }
86: }
87: *nz = k;
88: }
89: } else {
90: if (bs == 1) {
91: if (nz) {
92: k = 0;
93: for (j=0; j<nd; j++) {
94: pcol = row - diag[j];
95: if (pcol > -1 && pcol < nc) k++;
96: }
97: *nz = k;
98: }
99: if (col) {
100: *col = a->colloc;
101: k = 0;
102: for (j=0; j<nd; j++) {
103: pcol = row - diag[j];
104: if (pcol > -1 && pcol < nc) {
105: (*col)[k] = pcol; k++;
106: }
107: }
108: }
109: if (v) {
110: *v = a->dvalue;
111: k = 0;
112: for (j=0; j<nd; j++) {
113: pcol = row - diag[j];
114: if (pcol > -1 && pcol < nc) {
115: (*v)[k] = (a->diagv[j])[row]; k++;
116: }
117: }
118: }
119: } else {
120: if (nz) {
121: k = 0;
122: for (j=0; j<nd; j++) {
123: pcol = bs * (row/bs- diag[j]);
124: if (pcol > -1 && pcol < nc) k += bs;
125: }
126: *nz = k;
127: }
128: if (col) {
129: *col = a->colloc;
130: k = 0;
131: for (j=0; j<nd; j++) {
132: pcol = bs * (row/bs - diag[j]);
133: if (pcol > -1 && pcol < nc) {
134: for (i=0; i<bs; i++) {
135: (*col)[k+i] = pcol + i;
136: }
137: k += bs;
138: }
139: }
140: }
141: if (v) {
142: shift = (row/bs)*bs*bs + row%bs;
143: *v = a->dvalue;
144: k = 0;
145: for (j=0; j<nd; j++) {
146: pcol = bs * (row/bs - diag[j]);
147: if (pcol > -1 && pcol < nc) {
148: for (i=0; i<bs; i++) {
149: (*v)[k+i] = (a->diagv[j])[shift + i*bs];
150: }
151: k += bs;
152: }
153: }
154: }
155: }
156: }
157: return(0);
158: }
162: int MatRestoreRow_SeqBDiag(Mat A,int row,int *ncols,int **cols,PetscScalar **vals)
163: {
165: /* Work space is allocated during matrix creation and freed
166: when matrix is destroyed */
167: return(0);
168: }
170: /*
171: MatNorm_SeqBDiag_Columns - Computes the column norms of a block diagonal
172: matrix. We code this separately from MatNorm_SeqBDiag() so that the
173: routine can be used for the parallel version as well.
174: */
177: int MatNorm_SeqBDiag_Columns(Mat A,PetscReal *tmp,int n)
178: {
179: Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
180: int d,i,j,k,nd = a->nd,bs = a->bs,diag,kshift,kloc,len,ierr;
181: PetscScalar *dv;
184: PetscMemzero(tmp,A->n*sizeof(PetscReal));
185: if (bs == 1) {
186: for (d=0; d<nd; d++) {
187: dv = a->diagv[d];
188: diag = a->diag[d];
189: len = a->bdlen[d];
190: if (diag > 0) { /* lower triangle */
191: for (i=0; i<len; i++) {
192: tmp[i] += PetscAbsScalar(dv[i+diag]);
193: }
194: } else { /* upper triangle */
195: for (i=0; i<len; i++) {
196: tmp[i-diag] += PetscAbsScalar(dv[i]);
197: }
198: }
199: }
200: } else {
201: for (d=0; d<nd; d++) {
202: dv = a->diagv[d];
203: diag = a->diag[d];
204: len = a->bdlen[d];
206: if (diag > 0) { /* lower triangle */
207: for (k=0; k<len; k++) {
208: kloc = k*bs; kshift = kloc*bs + diag*bs;
209: for (i=0; i<bs; i++) { /* i = local row */
210: for (j=0; j<bs; j++) { /* j = local column */
211: tmp[kloc + j] += PetscAbsScalar(dv[kshift + j*bs + i]);
212: }
213: }
214: }
215: } else { /* upper triangle */
216: for (k=0; k<len; k++) {
217: kloc = k*bs; kshift = kloc*bs;
218: for (i=0; i<bs; i++) { /* i = local row */
219: for (j=0; j<bs; j++) { /* j = local column */
220: tmp[kloc + j - bs*diag] += PetscAbsScalar(dv[kshift + j*bs + i]);
221: }
222: }
223: }
224: }
225: }
226: }
227: return(0);
228: }
232: int MatNorm_SeqBDiag(Mat A,NormType type,PetscReal *nrm)
233: {
234: Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
235: PetscReal sum = 0.0,*tmp;
236: int ierr,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: int MatTranspose_SeqBDiag(Mat A,Mat *matout)
330: {
331: Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data,*anew;
332: Mat tmat;
333: int i,j,k,d,ierr,nd = a->nd,*diag = a->diag,*diagnew;
334: int bs = a->bs,kshift,shifto,shiftn;
335: PetscScalar *dwork,*dvnew;
338: PetscMalloc((nd+1)*sizeof(int),&diagnew);
339: for (i=0; i<nd; i++) {
340: diagnew[i] = -diag[nd-i-1]; /* assume sorted in descending order */
341: }
342: MatCreateSeqBDiag(A->comm,A->n,A->m,nd,bs,diagnew,0,&tmat);
343: PetscFree(diagnew);
344: anew = (Mat_SeqBDiag*)tmat->data;
345: for (d=0; d<nd; d++) {
346: dvnew = anew->diagv[d];
347: dwork = a->diagv[nd-d-1];
348: if (anew->bdlen[d] != a->bdlen[nd-d-1]) SETERRQ(PETSC_ERR_ARG_SIZ,"Incompatible diagonal lengths");
349: shifto = a->diag[nd-d-1];
350: shiftn = anew->diag[d];
351: if (shifto > 0) shifto = bs*bs*shifto; else shifto = 0;
352: if (shiftn > 0) shiftn = bs*bs*shiftn; else shiftn = 0;
353: if (bs == 1) {
354: for (k=0; k<anew->bdlen[d]; k++) dvnew[shiftn+k] = dwork[shifto+k];
355: } else {
356: for (k=0; k<anew->bdlen[d]; k++) {
357: kshift = k*bs*bs;
358: for (i=0; i<bs; i++) { /* i = local row */
359: for (j=0; j<bs; j++) { /* j = local column */
360: dvnew[shiftn + kshift + j + i*bs] = dwork[shifto + kshift + j*bs + i];
361: }
362: }
363: }
364: }
365: }
366: MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);
367: MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);
368: if (matout) {
369: *matout = tmat;
370: } else {
371: /* This isn't really an in-place transpose ... but free data
372: structures from a. We should fix this. */
373: if (!a->user_alloc) { /* Free the actual diagonals */
374: for (i=0; i<a->nd; i++) {
375: if (a->diag[i] > 0) {
376: PetscFree(a->diagv[i] + bs*bs*a->diag[i]);
377: } else {
378: PetscFree(a->diagv[i]);
379: }
380: }
381: }
382: if (a->pivot) {PetscFree(a->pivot);}
383: PetscFree(a->diagv);
384: PetscFree(a->diag);
385: PetscFree(a->colloc);
386: PetscFree(a->dvalue);
387: PetscFree(a);
388: PetscMemcpy(A,tmat,sizeof(struct _p_Mat));
389: PetscHeaderDestroy(tmat);
390: }
391: return(0);
392: }
394: /* ----------------------------------------------------------------*/
399: int MatView_SeqBDiag_Binary(Mat A,PetscViewer viewer)
400: {
401: Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
402: int i,ict,fd,*col_lens,*cval,*col,ierr,nz;
403: PetscScalar *anonz,*val;
406: PetscViewerBinaryGetDescriptor(viewer,&fd);
408: /* For MATSEQBDIAG format,maxnz = nz */
409: PetscMalloc((4+A->m)*sizeof(int),&col_lens);
410: col_lens[0] = MAT_FILE_COOKIE;
411: col_lens[1] = A->m;
412: col_lens[2] = A->n;
413: col_lens[3] = a->maxnz;
415: /* Should do translation using less memory; this is just a quick initial version */
416: PetscMalloc((a->maxnz)*sizeof(int),&cval);
417: PetscMalloc((a->maxnz)*sizeof(PetscScalar),&anonz);
419: ict = 0;
420: for (i=0; i<A->m; i++) {
421: MatGetRow_SeqBDiag(A,i,&nz,&col,&val);
422: col_lens[4+i] = nz;
423: PetscMemcpy(&cval[ict],col,nz*sizeof(int));
424: PetscMemcpy(&anonz[ict],anonz,nz*sizeof(PetscScalar));
425: MatRestoreRow_SeqBDiag(A,i,&nz,&col,&val);
426: ict += nz;
427: }
428: if (ict != a->maxnz) SETERRQ(PETSC_ERR_PLIB,"Error in nonzero count");
430: /* Store lengths of each row and write (including header) to file */
431: PetscBinaryWrite(fd,col_lens,4+A->m,PETSC_INT,1);
432: PetscFree(col_lens);
434: /* Store column indices (zero start index) */
435: PetscBinaryWrite(fd,cval,a->maxnz,PETSC_INT,0);
437: /* Store nonzero values */
438: PetscBinaryWrite(fd,anonz,a->maxnz,PETSC_SCALAR,0);
439: return(0);
440: }
444: int MatView_SeqBDiag_ASCII(Mat A,PetscViewer viewer)
445: {
446: Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
447: char *name;
448: int ierr,*col,i,j,len,diag,nr = A->m,bs = a->bs,iprint,nz;
449: PetscScalar *val,*dv,zero = 0.0;
450: PetscViewerFormat format;
453: PetscObjectGetName((PetscObject)A,&name);
454: PetscViewerGetFormat(viewer,&format);
455: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
456: int nline = PetscMin(10,a->nd),k,nk,np;
457: if (a->user_alloc) {
458: PetscViewerASCIIPrintf(viewer,"block size=%d, number of diagonals=%d, user-allocated storage\n",bs,a->nd);
459: } else {
460: PetscViewerASCIIPrintf(viewer,"block size=%d, number of diagonals=%d, PETSc-allocated storage\n",bs,a->nd);
461: }
462: nk = (a->nd-1)/nline + 1;
463: for (k=0; k<nk; k++) {
464: PetscViewerASCIIPrintf(viewer,"diag numbers:");
465: np = PetscMin(nline,a->nd - nline*k);
466: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
467: for (i=0; i<np; i++) {
468: PetscViewerASCIIPrintf(viewer," %d",a->diag[i+nline*k]);
469: }
470: PetscViewerASCIIPrintf(viewer,"\n");
471: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
472: }
473: } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
474: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
475: PetscViewerASCIIPrintf(viewer,"%% Size = %d %d \n",nr,A->n);
476: PetscViewerASCIIPrintf(viewer,"%% Nonzeros = %d \n",a->nz);
477: PetscViewerASCIIPrintf(viewer,"zzz = zeros(%d,3);\n",a->nz);
478: PetscViewerASCIIPrintf(viewer,"zzz = [\n");
479: for (i=0; i<A->m; i++) {
480: MatGetRow_SeqBDiag(A,i,&nz,&col,&val);
481: for (j=0; j<nz; j++) {
482: if (val[j] != zero) {
483: #if defined(PETSC_USE_COMPLEX)
484: PetscViewerASCIIPrintf(viewer,"%d %d %18.16e %18.16ei \n",
485: i+1,col[j]+1,PetscRealPart(val[j]),PetscImaginaryPart(val[j]));
486: #else
487: PetscViewerASCIIPrintf(viewer,"%d %d %18.16ei \n",i+1,col[j]+1,val[j]);
488: #endif
489: }
490: }
491: MatRestoreRow_SeqBDiag(A,i,&nz,&col,&val);
492: }
493: PetscViewerASCIIPrintf(viewer,"];\n %s = spconvert(zzz);\n",name);
494: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
495: } else if (format == PETSC_VIEWER_ASCII_IMPL) {
496: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
497: if (bs == 1) { /* diagonal format */
498: for (i=0; i<a->nd; i++) {
499: dv = a->diagv[i];
500: diag = a->diag[i];
501: PetscViewerASCIIPrintf(viewer,"\n<diagonal %d>\n",diag);
502: /* diag[i] is (row-col)/bs */
503: if (diag > 0) { /* lower triangle */
504: len = a->bdlen[i];
505: for (j=0; j<len; j++) {
506: if (dv[diag+j] != zero) {
507: #if defined(PETSC_USE_COMPLEX)
508: if (PetscImaginaryPart(dv[diag+j]) != 0.0) {
509: PetscViewerASCIIPrintf(viewer,"A[ %d , %d ] = %e + %e i\n",
510: j+diag,j,PetscRealPart(dv[diag+j]),PetscImaginaryPart(dv[diag+j]));
511: } else {
512: PetscViewerASCIIPrintf(viewer,"A[ %d , %d ] = %e\n",j+diag,j,PetscRealPart(dv[diag+j]));
513: }
514: #else
515: PetscViewerASCIIPrintf(viewer,"A[ %d , %d ] = %e\n",j+diag,j,dv[diag+j]);
517: #endif
518: }
519: }
520: } else { /* upper triangle, including main diagonal */
521: len = a->bdlen[i];
522: for (j=0; j<len; j++) {
523: if (dv[j] != zero) {
524: #if defined(PETSC_USE_COMPLEX)
525: if (PetscImaginaryPart(dv[j]) != 0.0) {
526: PetscViewerASCIIPrintf(viewer,"A[ %d , %d ] = %g + %g i\n",
527: j,j-diag,PetscRealPart(dv[j]),PetscImaginaryPart(dv[j]));
528: } else {
529: PetscViewerASCIIPrintf(viewer,"A[ %d , %d ] = %g\n",j,j-diag,PetscRealPart(dv[j]));
530: }
531: #else
532: PetscViewerASCIIPrintf(viewer,"A[ %d , %d ] = %g\n",j,j-diag,dv[j]);
533: #endif
534: }
535: }
536: }
537: }
538: } else { /* Block diagonals */
539: int d,k,kshift;
540: for (d=0; d< a->nd; d++) {
541: dv = a->diagv[d];
542: diag = a->diag[d];
543: len = a->bdlen[d];
544: PetscViewerASCIIPrintf(viewer,"\n<diagonal %d>\n", diag);
545: if (diag > 0) { /* lower triangle */
546: for (k=0; k<len; k++) {
547: kshift = (diag+k)*bs*bs;
548: for (i=0; i<bs; i++) {
549: iprint = 0;
550: for (j=0; j<bs; j++) {
551: if (dv[kshift + j*bs + i] != zero) {
552: iprint = 1;
553: #if defined(PETSC_USE_COMPLEX)
554: if (PetscImaginaryPart(dv[kshift + j*bs + i])){
555: PetscViewerASCIIPrintf(viewer,"A[%d,%d]=%5.2e + %5.2e i ",(k+diag)*bs+i,k*bs+j,
556: PetscRealPart(dv[kshift + j*bs + i]),PetscImaginaryPart(dv[kshift + j*bs + i]));
557: } else {
558: PetscViewerASCIIPrintf(viewer,"A[%d,%d]=%5.2e ",(k+diag)*bs+i,k*bs+j,
559: PetscRealPart(dv[kshift + j*bs + i]));
560: }
561: #else
562: PetscViewerASCIIPrintf(viewer,"A[%d,%d]=%5.2e ",(k+diag)*bs+i,k*bs+j,
563: dv[kshift + j*bs + i]);
564: #endif
565: }
566: }
567: if (iprint) {PetscViewerASCIIPrintf(viewer,"\n");}
568: }
569: }
570: } else { /* upper triangle, including main diagonal */
571: for (k=0; k<len; k++) {
572: kshift = k*bs*bs;
573: for (i=0; i<bs; i++) {
574: iprint = 0;
575: for (j=0; j<bs; j++) {
576: if (dv[kshift + j*bs + i] != zero) {
577: iprint = 1;
578: #if defined(PETSC_USE_COMPLEX)
579: if (PetscImaginaryPart(dv[kshift + j*bs + i])){
580: PetscViewerASCIIPrintf(viewer,"A[%d,%d]=%5.2e + %5.2e i ",k*bs+i,(k-diag)*bs+j,
581: PetscRealPart(dv[kshift + j*bs + i]),PetscImaginaryPart(dv[kshift + j*bs + i]));
582: } else {
583: PetscViewerASCIIPrintf(viewer,"A[%d,%d]=%5.2e ",k*bs+i,(k-diag)*bs+j,
584: PetscRealPart(dv[kshift + j*bs + i]));
585: }
586: #else
587: PetscViewerASCIIPrintf(viewer,"A[%d,%d]=%5.2e ",k*bs+i,(k-diag)*bs+j,
588: dv[kshift + j*bs + i]);
589: #endif
590: }
591: }
592: if (iprint) {PetscViewerASCIIPrintf(viewer,"\n");}
593: }
594: }
595: }
596: }
597: }
598: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
599: } else {
600: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
601: /* the usual row format (PETSC_VIEWER_ASCII_NONZERO_ONLY) */
602: for (i=0; i<A->m; i++) {
603: PetscViewerASCIIPrintf(viewer,"row %d:",i);
604: MatGetRow_SeqBDiag(A,i,&nz,&col,&val);
605: for (j=0; j<nz; j++) {
606: #if defined(PETSC_USE_COMPLEX)
607: if (PetscImaginaryPart(val[j]) != 0.0 && PetscRealPart(val[j]) != 0.0) {
608: PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",col[j],PetscRealPart(val[j]),PetscImaginaryPart(val[j]));
609: } else if (PetscRealPart(val[j]) != 0.0) {
610: PetscViewerASCIIPrintf(viewer," (%d, %g) ",col[j],PetscRealPart(val[j]));
611: }
612: #else
613: if (val[j] != 0.0) {PetscViewerASCIIPrintf(viewer," (%d, %g) ",col[j],val[j]);}
614: #endif
615: }
616: PetscViewerASCIIPrintf(viewer,"\n");
617: MatRestoreRow_SeqBDiag(A,i,&nz,&col,&val);
618: }
619: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
620: }
621: PetscViewerFlush(viewer);
622: return(0);
623: }
627: static int MatView_SeqBDiag_Draw(Mat A,PetscViewer viewer)
628: {
629: PetscDraw draw;
630: PetscReal xl,yl,xr,yr,w,h;
631: int ierr,nz,*col,i,j,nr = A->m;
632: PetscTruth isnull;
635: PetscViewerDrawGetDraw(viewer,0,&draw);
636: PetscDrawIsNull(draw,&isnull); if (isnull) return(0);
638: xr = A->n; yr = A->m; h = yr/10.0; w = xr/10.0;
639: xr += w; yr += h; xl = -w; yl = -h;
640: PetscDrawSetCoordinates(draw,xl,yl,xr,yr);
642: /* loop over matrix elements drawing boxes; we really should do this
643: by diagonals. What do we really want to draw here: nonzeros,
644: allocated space? */
645: for (i=0; i<nr; i++) {
646: yl = nr - i - 1.0; yr = yl + 1.0;
647: MatGetRow_SeqBDiag(A,i,&nz,&col,0);
648: for (j=0; j<nz; j++) {
649: xl = col[j]; xr = xl + 1.0;
650: PetscDrawRectangle(draw,xl,yl,xr,yr,PETSC_DRAW_BLACK,PETSC_DRAW_BLACK,
651: PETSC_DRAW_BLACK,PETSC_DRAW_BLACK);
652: }
653: MatRestoreRow_SeqBDiag(A,i,&nz,&col,0);
654: }
655: PetscDrawFlush(draw);
656: PetscDrawPause(draw);
657: return(0);
658: }
662: int MatView_SeqBDiag(Mat A,PetscViewer viewer)
663: {
664: int ierr;
665: PetscTruth isascii,isbinary,isdraw;
668: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
669: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
670: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
671: if (isascii) {
672: MatView_SeqBDiag_ASCII(A,viewer);
673: } else if (isbinary) {
674: MatView_SeqBDiag_Binary(A,viewer);
675: } else if (isdraw) {
676: MatView_SeqBDiag_Draw(A,viewer);
677: } else {
678: SETERRQ1(1,"Viewer type %s not supported by BDiag matrices",((PetscObject)viewer)->type_name);
679: }
680: return(0);
681: }