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: #undef __FUNCT__
36: int MatGetInfo_SeqBDiag(Mat A,MatInfoType flag,MatInfo *info)
37: {
38: Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
41: info->rows_global = (double)A->m;
42: info->columns_global = (double)A->n;
43: info->rows_local = (double)A->m;
44: info->columns_local = (double)A->n;
45: info->block_size = a->bs;
46: info->nz_allocated = (double)a->maxnz;
47: info->nz_used = (double)a->nz;
48: info->nz_unneeded = (double)(a->maxnz - a->nz);
49: info->assemblies = (double)A->num_ass;
50: info->mallocs = (double)a->reallocs;
51: info->memory = A->mem;
52: info->fill_ratio_given = 0; /* supports ILU(0) only */
53: info->fill_ratio_needed = 0;
54: info->factor_mallocs = 0;
55: return(0);
56: }
58: /*
59: Note: this currently will generate different answers if you request
60: all items or a subset. If you request all items it checks if the value is
61: nonzero and only includes it if it is nonzero; if you check a subset of items
62: it returns a list of all active columns in the row (some which may contain
63: a zero)
64: */
65: #undef __FUNCT__
67: int MatGetRow_SeqBDiag(Mat A,int row,int *nz,int **col,PetscScalar **v)
68: {
69: Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
70: int nd = a->nd,bs = a->bs;
71: int nc = A->n,*diag = a->diag,pcol,shift,i,j,k;
74: /* For efficiency, if ((nz) && (col) && (v)) then do all at once */
75: if ((nz) && (col) && (v)) {
76: *col = a->colloc;
77: *v = a->dvalue;
78: k = 0;
79: if (bs == 1) {
80: for (j=0; j<nd; j++) {
81: pcol = row - diag[j];
82: #if defined(PETSC_USE_COMPLEX)
83: if (pcol > -1 && pcol < nc && PetscAbsScalar((a->diagv[j])[row])) {
84: #else
85: if (pcol > -1 && pcol < nc && (a->diagv[j])[row]) {
86: #endif
87: (*v)[k] = (a->diagv[j])[row];
88: (*col)[k] = pcol;
89: k++;
90: }
91: }
92: *nz = k;
93: } else {
94: shift = (row/bs)*bs*bs + row%bs;
95: for (j=0; j<nd; j++) {
96: pcol = bs * (row/bs - diag[j]);
97: if (pcol > -1 && pcol < nc) {
98: for (i=0; i<bs; i++) {
99: #if defined(PETSC_USE_COMPLEX)
100: if (PetscAbsScalar((a->diagv[j])[shift + i*bs])) {
101: #else
102: if ((a->diagv[j])[shift + i*bs]) {
103: #endif
104: (*v)[k] = (a->diagv[j])[shift + i*bs];
105: (*col)[k] = pcol + i;
106: k++;
107: }
108: }
109: }
110: }
111: *nz = k;
112: }
113: } else {
114: if (bs == 1) {
115: if (nz) {
116: k = 0;
117: for (j=0; j<nd; j++) {
118: pcol = row - diag[j];
119: if (pcol > -1 && pcol < nc) k++;
120: }
121: *nz = k;
122: }
123: if (col) {
124: *col = a->colloc;
125: k = 0;
126: for (j=0; j<nd; j++) {
127: pcol = row - diag[j];
128: if (pcol > -1 && pcol < nc) {
129: (*col)[k] = pcol; k++;
130: }
131: }
132: }
133: if (v) {
134: *v = a->dvalue;
135: k = 0;
136: for (j=0; j<nd; j++) {
137: pcol = row - diag[j];
138: if (pcol > -1 && pcol < nc) {
139: (*v)[k] = (a->diagv[j])[row]; k++;
140: }
141: }
142: }
143: } else {
144: if (nz) {
145: k = 0;
146: for (j=0; j<nd; j++) {
147: pcol = bs * (row/bs- diag[j]);
148: if (pcol > -1 && pcol < nc) k += bs;
149: }
150: *nz = k;
151: }
152: if (col) {
153: *col = a->colloc;
154: k = 0;
155: for (j=0; j<nd; j++) {
156: pcol = bs * (row/bs - diag[j]);
157: if (pcol > -1 && pcol < nc) {
158: for (i=0; i<bs; i++) {
159: (*col)[k+i] = pcol + i;
160: }
161: k += bs;
162: }
163: }
164: }
165: if (v) {
166: shift = (row/bs)*bs*bs + row%bs;
167: *v = a->dvalue;
168: k = 0;
169: for (j=0; j<nd; j++) {
170: pcol = bs * (row/bs - diag[j]);
171: if (pcol > -1 && pcol < nc) {
172: for (i=0; i<bs; i++) {
173: (*v)[k+i] = (a->diagv[j])[shift + i*bs];
174: }
175: k += bs;
176: }
177: }
178: }
179: }
180: }
181: return(0);
182: }
184: #undef __FUNCT__
186: int MatRestoreRow_SeqBDiag(Mat A,int row,int *ncols,int **cols,PetscScalar **vals)
187: {
189: /* Work space is allocated during matrix creation and freed
190: when matrix is destroyed */
191: return(0);
192: }
194: /*
195: MatNorm_SeqBDiag_Columns - Computes the column norms of a block diagonal
196: matrix. We code this separately from MatNorm_SeqBDiag() so that the
197: routine can be used for the parallel version as well.
198: */
199: #undef __FUNCT__
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: PetscScalar *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: #undef __FUNCT__
256: int MatNorm_SeqBDiag(Mat A,NormType type,PetscReal *nrm)
257: {
258: Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
259: PetscReal sum = 0.0,*tmp;
260: int ierr,d,i,j,k,nd = a->nd,bs = a->bs,diag,kshift,kloc,len;
261: PetscScalar *dv;
264: if (type == NORM_FROBENIUS) {
265: for (d=0; d<nd; d++) {
266: dv = a->diagv[d];
267: len = a->bdlen[d]*bs*bs;
268: diag = a->diag[d];
269: if (diag > 0) {
270: for (i=0; i<len; i++) {
271: #if defined(PETSC_USE_COMPLEX)
272: sum += PetscRealPart(PetscConj(dv[i+diag])*dv[i+diag]);
273: #else
274: sum += dv[i+diag]*dv[i+diag];
275: #endif
276: }
277: } else {
278: for (i=0; i<len; i++) {
279: #if defined(PETSC_USE_COMPLEX)
280: sum += PetscRealPart(PetscConj(dv[i])*dv[i]);
281: #else
282: sum += dv[i]*dv[i];
283: #endif
284: }
285: }
286: }
287: *nrm = sqrt(sum);
288: } else if (type == NORM_1) { /* max column norm */
289: PetscMalloc((A->n+1)*sizeof(PetscReal),&tmp);
290: MatNorm_SeqBDiag_Columns(A,tmp,A->n);
291: *nrm = 0.0;
292: for (j=0; j<A->n; j++) {
293: if (tmp[j] > *nrm) *nrm = tmp[j];
294: }
295: PetscFree(tmp);
296: } else if (type == NORM_INFINITY) { /* max row norm */
297: PetscMalloc((A->m+1)*sizeof(PetscReal),&tmp);
298: PetscMemzero(tmp,A->m*sizeof(PetscReal));
299: *nrm = 0.0;
300: if (bs == 1) {
301: for (d=0; d<nd; d++) {
302: dv = a->diagv[d];
303: diag = a->diag[d];
304: len = a->bdlen[d];
305: if (diag > 0) { /* lower triangle */
306: for (i=0; i<len; i++) {
307: tmp[i+diag] += PetscAbsScalar(dv[i+diag]);
308: }
309: } else { /* upper triangle */
310: for (i=0; i<len; i++) {
311: tmp[i] += PetscAbsScalar(dv[i]);
312: }
313: }
314: }
315: } else {
316: for (d=0; d<nd; d++) {
317: dv = a->diagv[d];
318: diag = a->diag[d];
319: len = a->bdlen[d];
320: if (diag > 0) {
321: for (k=0; k<len; k++) {
322: kloc = k*bs; kshift = kloc*bs + bs*diag;
323: for (i=0; i<bs; i++) { /* i = local row */
324: for (j=0; j<bs; j++) { /* j = local column */
325: tmp[kloc + i + bs*diag] += PetscAbsScalar(dv[kshift+j*bs+i]);
326: }
327: }
328: }
329: } else {
330: for (k=0; k<len; k++) {
331: kloc = k*bs; kshift = kloc*bs;
332: for (i=0; i<bs; i++) { /* i = local row */
333: for (j=0; j<bs; j++) { /* j = local column */
334: tmp[kloc + i] += PetscAbsScalar(dv[kshift + j*bs + i]);
335: }
336: }
337: }
338: }
339: }
340: }
341: for (j=0; j<A->m; j++) {
342: if (tmp[j] > *nrm) *nrm = tmp[j];
343: }
344: PetscFree(tmp);
345: } else {
346: SETERRQ(PETSC_ERR_SUP,"No support for two norm");
347: }
348: return(0);
349: }
351: #undef __FUNCT__
353: int MatTranspose_SeqBDiag(Mat A,Mat *matout)
354: {
355: Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data,*anew;
356: Mat tmat;
357: int i,j,k,d,ierr,nd = a->nd,*diag = a->diag,*diagnew;
358: int bs = a->bs,kshift,shifto,shiftn;
359: PetscScalar *dwork,*dvnew;
362: PetscMalloc((nd+1)*sizeof(int),&diagnew);
363: for (i=0; i<nd; i++) {
364: diagnew[i] = -diag[nd-i-1]; /* assume sorted in descending order */
365: }
366: MatCreateSeqBDiag(A->comm,A->n,A->m,nd,bs,diagnew,0,&tmat);
367: PetscFree(diagnew);
368: anew = (Mat_SeqBDiag*)tmat->data;
369: for (d=0; d<nd; d++) {
370: dvnew = anew->diagv[d];
371: dwork = a->diagv[nd-d-1];
372: if (anew->bdlen[d] != a->bdlen[nd-d-1]) SETERRQ(PETSC_ERR_ARG_SIZ,"Incompatible diagonal lengths");
373: shifto = a->diag[nd-d-1];
374: shiftn = anew->diag[d];
375: if (shifto > 0) shifto = bs*bs*shifto; else shifto = 0;
376: if (shiftn > 0) shiftn = bs*bs*shiftn; else shiftn = 0;
377: if (bs == 1) {
378: for (k=0; k<anew->bdlen[d]; k++) dvnew[shiftn+k] = dwork[shifto+k];
379: } else {
380: for (k=0; k<anew->bdlen[d]; k++) {
381: kshift = k*bs*bs;
382: for (i=0; i<bs; i++) { /* i = local row */
383: for (j=0; j<bs; j++) { /* j = local column */
384: dvnew[shiftn + kshift + j + i*bs] = dwork[shifto + kshift + j*bs + i];
385: }
386: }
387: }
388: }
389: }
390: MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);
391: MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);
392: if (matout) {
393: *matout = tmat;
394: } else {
395: /* This isn't really an in-place transpose ... but free data
396: structures from a. We should fix this. */
397: if (!a->user_alloc) { /* Free the actual diagonals */
398: for (i=0; i<a->nd; i++) {
399: if (a->diag[i] > 0) {
400: PetscFree(a->diagv[i] + bs*bs*a->diag[i]);
401: } else {
402: PetscFree(a->diagv[i]);
403: }
404: }
405: }
406: if (a->pivot) {PetscFree(a->pivot);}
407: PetscFree(a->diagv);
408: PetscFree(a->diag);
409: PetscFree(a->colloc);
410: PetscFree(a->dvalue);
411: PetscFree(a);
412: PetscMemcpy(A,tmat,sizeof(struct _p_Mat));
413: PetscHeaderDestroy(tmat);
414: }
415: return(0);
416: }
418: /* ----------------------------------------------------------------*/
421: #undef __FUNCT__
423: int MatView_SeqBDiag_Binary(Mat A,PetscViewer viewer)
424: {
425: Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
426: int i,ict,fd,*col_lens,*cval,*col,ierr,nz;
427: PetscScalar *anonz,*val;
430: PetscViewerBinaryGetDescriptor(viewer,&fd);
432: /* For MATSEQBDIAG format,maxnz = nz */
433: ierr = PetscMalloc((4+A->m)*sizeof(int),&col_lens);
434: col_lens[0] = MAT_FILE_COOKIE;
435: col_lens[1] = A->m;
436: col_lens[2] = A->n;
437: col_lens[3] = a->maxnz;
439: /* Should do translation using less memory; this is just a quick initial version */
440: PetscMalloc((a->maxnz)*sizeof(int),&cval);
441: PetscMalloc((a->maxnz)*sizeof(PetscScalar),&anonz);
443: ict = 0;
444: for (i=0; i<A->m; i++) {
445: MatGetRow_SeqBDiag(A,i,&nz,&col,&val);
446: col_lens[4+i] = nz;
447: PetscMemcpy(&cval[ict],col,nz*sizeof(int));
448: PetscMemcpy(&anonz[ict],anonz,nz*sizeof(PetscScalar));
449: MatRestoreRow_SeqBDiag(A,i,&nz,&col,&val);
450: ict += nz;
451: }
452: if (ict != a->maxnz) SETERRQ(PETSC_ERR_PLIB,"Error in nonzero count");
454: /* Store lengths of each row and write (including header) to file */
455: PetscBinaryWrite(fd,col_lens,4+A->m,PETSC_INT,1);
456: PetscFree(col_lens);
458: /* Store column indices (zero start index) */
459: PetscBinaryWrite(fd,cval,a->maxnz,PETSC_INT,0);
461: /* Store nonzero values */
462: PetscBinaryWrite(fd,anonz,a->maxnz,PETSC_SCALAR,0);
463: return(0);
464: }
466: #undef __FUNCT__
468: int MatView_SeqBDiag_ASCII(Mat A,PetscViewer viewer)
469: {
470: Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
471: char *name;
472: int ierr,*col,i,j,len,diag,nr = A->m,bs = a->bs,iprint,nz;
473: PetscScalar *val,*dv,zero = 0.0;
474: PetscViewerFormat format;
477: PetscObjectGetName((PetscObject)A,&name);
478: PetscViewerGetFormat(viewer,&format);
479: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
480: int nline = PetscMin(10,a->nd),k,nk,np;
481: if (a->user_alloc) {
482: PetscViewerASCIIPrintf(viewer,"block size=%d, number of diagonals=%d, user-allocated storagen",bs,a->nd);
483: } else {
484: PetscViewerASCIIPrintf(viewer,"block size=%d, number of diagonals=%d, PETSc-allocated storagen",bs,a->nd);
485: }
486: nk = (a->nd-1)/nline + 1;
487: for (k=0; k<nk; k++) {
488: PetscViewerASCIIPrintf(viewer,"diag numbers:");
489: np = PetscMin(nline,a->nd - nline*k);
490: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
491: for (i=0; i<np; i++) {
492: PetscViewerASCIIPrintf(viewer," %d",a->diag[i+nline*k]);
493: }
494: PetscViewerASCIIPrintf(viewer,"n");
495: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
496: }
497: } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
498: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
499: PetscViewerASCIIPrintf(viewer,"%% Size = %d %d n",nr,A->n);
500: PetscViewerASCIIPrintf(viewer,"%% Nonzeros = %d n",a->nz);
501: PetscViewerASCIIPrintf(viewer,"zzz = zeros(%d,3);n",a->nz);
502: PetscViewerASCIIPrintf(viewer,"zzz = [n");
503: for (i=0; i<A->m; i++) {
504: MatGetRow_SeqBDiag(A,i,&nz,&col,&val);
505: for (j=0; j<nz; j++) {
506: if (val[j] != zero) {
507: #if defined(PETSC_USE_COMPLEX)
508: PetscViewerASCIIPrintf(viewer,"%d %d %18.16e %18.16ei n",
509: i+1,col[j]+1,PetscRealPart(val[j]),PetscImaginaryPart(val[j]));
510: #else
511: PetscViewerASCIIPrintf(viewer,"%d %d %18.16ei n",i+1,col[j]+1,val[j]);
512: #endif
513: }
514: }
515: MatRestoreRow_SeqBDiag(A,i,&nz,&col,&val);
516: }
517: PetscViewerASCIIPrintf(viewer,"];n %s = spconvert(zzz);n",name);
518: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
519: } else if (format == PETSC_VIEWER_ASCII_IMPL) {
520: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
521: if (bs == 1) { /* diagonal format */
522: for (i=0; i<a->nd; i++) {
523: dv = a->diagv[i];
524: diag = a->diag[i];
525: PetscViewerASCIIPrintf(viewer,"n<diagonal %d>n",diag);
526: /* diag[i] is (row-col)/bs */
527: if (diag > 0) { /* lower triangle */
528: len = a->bdlen[i];
529: for (j=0; j<len; j++) {
530: if (dv[diag+j] != zero) {
531: #if defined(PETSC_USE_COMPLEX)
532: if (PetscImaginaryPart(dv[diag+j]) != 0.0) {
533: PetscViewerASCIIPrintf(viewer,"A[ %d , %d ] = %e + %e in",
534: j+diag,j,PetscRealPart(dv[diag+j]),PetscImaginaryPart(dv[diag+j]));
535: } else {
536: PetscViewerASCIIPrintf(viewer,"A[ %d , %d ] = %en",j+diag,j,PetscRealPart(dv[diag+j]));
537: }
538: #else
539: PetscViewerASCIIPrintf(viewer,"A[ %d , %d ] = %en",j+diag,j,dv[diag+j]);
541: #endif
542: }
543: }
544: } else { /* upper triangle, including main diagonal */
545: len = a->bdlen[i];
546: for (j=0; j<len; j++) {
547: if (dv[j] != zero) {
548: #if defined(PETSC_USE_COMPLEX)
549: if (PetscImaginaryPart(dv[j]) != 0.0) {
550: PetscViewerASCIIPrintf(viewer,"A[ %d , %d ] = %g + %g in",
551: j,j-diag,PetscRealPart(dv[j]),PetscImaginaryPart(dv[j]));
552: } else {
553: PetscViewerASCIIPrintf(viewer,"A[ %d , %d ] = %gn",j,j-diag,PetscRealPart(dv[j]));
554: }
555: #else
556: PetscViewerASCIIPrintf(viewer,"A[ %d , %d ] = %gn",j,j-diag,dv[j]);
557: #endif
558: }
559: }
560: }
561: }
562: } else { /* Block diagonals */
563: int d,k,kshift;
564: for (d=0; d< a->nd; d++) {
565: dv = a->diagv[d];
566: diag = a->diag[d];
567: len = a->bdlen[d];
568: PetscViewerASCIIPrintf(viewer,"n<diagonal %d>n", diag);
569: if (diag > 0) { /* lower triangle */
570: for (k=0; k<len; k++) {
571: kshift = (diag+k)*bs*bs;
572: for (i=0; i<bs; i++) {
573: iprint = 0;
574: for (j=0; j<bs; j++) {
575: if (dv[kshift + j*bs + i] != zero) {
576: iprint = 1;
577: #if defined(PETSC_USE_COMPLEX)
578: if (PetscImaginaryPart(dv[kshift + j*bs + i])){
579: PetscViewerASCIIPrintf(viewer,"A[%d,%d]=%5.2e + %5.2e i ",(k+diag)*bs+i,k*bs+j,
580: PetscRealPart(dv[kshift + j*bs + i]),PetscImaginaryPart(dv[kshift + j*bs + i]));
581: } else {
582: PetscViewerASCIIPrintf(viewer,"A[%d,%d]=%5.2e ",(k+diag)*bs+i,k*bs+j,
583: PetscRealPart(dv[kshift + j*bs + i]));
584: }
585: #else
586: PetscViewerASCIIPrintf(viewer,"A[%d,%d]=%5.2e ",(k+diag)*bs+i,k*bs+j,
587: dv[kshift + j*bs + i]);
588: #endif
589: }
590: }
591: if (iprint) {PetscViewerASCIIPrintf(viewer,"n");}
592: }
593: }
594: } else { /* upper triangle, including main diagonal */
595: for (k=0; k<len; k++) {
596: kshift = k*bs*bs;
597: for (i=0; i<bs; i++) {
598: iprint = 0;
599: for (j=0; j<bs; j++) {
600: if (dv[kshift + j*bs + i] != zero) {
601: iprint = 1;
602: #if defined(PETSC_USE_COMPLEX)
603: if (PetscImaginaryPart(dv[kshift + j*bs + i])){
604: PetscViewerASCIIPrintf(viewer,"A[%d,%d]=%5.2e + %5.2e i ",k*bs+i,(k-diag)*bs+j,
605: PetscRealPart(dv[kshift + j*bs + i]),PetscImaginaryPart(dv[kshift + j*bs + i]));
606: } else {
607: PetscViewerASCIIPrintf(viewer,"A[%d,%d]=%5.2e ",k*bs+i,(k-diag)*bs+j,
608: PetscRealPart(dv[kshift + j*bs + i]));
609: }
610: #else
611: PetscViewerASCIIPrintf(viewer,"A[%d,%d]=%5.2e ",k*bs+i,(k-diag)*bs+j,
612: dv[kshift + j*bs + i]);
613: #endif
614: }
615: }
616: if (iprint) {PetscViewerASCIIPrintf(viewer,"n");}
617: }
618: }
619: }
620: }
621: }
622: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
623: } else {
624: PetscViewerASCIIUseTabs(viewer,PETSC_NO);
625: /* the usual row format (PETSC_VIEWER_ASCII_NONZERO_ONLY) */
626: for (i=0; i<A->m; i++) {
627: PetscViewerASCIIPrintf(viewer,"row %d:",i);
628: MatGetRow_SeqBDiag(A,i,&nz,&col,&val);
629: for (j=0; j<nz; j++) {
630: #if defined(PETSC_USE_COMPLEX)
631: if (PetscImaginaryPart(val[j]) != 0.0 && PetscRealPart(val[j]) != 0.0) {
632: PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",col[j],PetscRealPart(val[j]),PetscImaginaryPart(val[j]));
633: } else if (PetscRealPart(val[j]) != 0.0) {
634: PetscViewerASCIIPrintf(viewer," (%d, %g) ",col[j],PetscRealPart(val[j]));
635: }
636: #else
637: if (val[j] != 0.0) {PetscViewerASCIIPrintf(viewer," (%d, %g) ",col[j],val[j]);}
638: #endif
639: }
640: PetscViewerASCIIPrintf(viewer,"n");
641: MatRestoreRow_SeqBDiag(A,i,&nz,&col,&val);
642: }
643: PetscViewerASCIIUseTabs(viewer,PETSC_YES);
644: }
645: PetscViewerFlush(viewer);
646: return(0);
647: }
649: #undef __FUNCT__
651: static int MatView_SeqBDiag_Draw(Mat A,PetscViewer viewer)
652: {
653: PetscDraw draw;
654: PetscReal xl,yl,xr,yr,w,h;
655: int ierr,nz,*col,i,j,nr = A->m;
656: PetscTruth isnull;
659: PetscViewerDrawGetDraw(viewer,0,&draw);
660: PetscDrawIsNull(draw,&isnull); if (isnull) return(0);
662: xr = A->n; yr = A->m; h = yr/10.0; w = xr/10.0;
663: xr += w; yr += h; xl = -w; yl = -h;
664: PetscDrawSetCoordinates(draw,xl,yl,xr,yr);
666: /* loop over matrix elements drawing boxes; we really should do this
667: by diagonals. What do we really want to draw here: nonzeros,
668: allocated space? */
669: for (i=0; i<nr; i++) {
670: yl = nr - i - 1.0; yr = yl + 1.0;
671: MatGetRow_SeqBDiag(A,i,&nz,&col,0);
672: for (j=0; j<nz; j++) {
673: xl = col[j]; xr = xl + 1.0;
674: PetscDrawRectangle(draw,xl,yl,xr,yr,PETSC_DRAW_BLACK,PETSC_DRAW_BLACK,
675: PETSC_DRAW_BLACK,PETSC_DRAW_BLACK);
676: }
677: MatRestoreRow_SeqBDiag(A,i,&nz,&col,0);
678: }
679: PetscDrawFlush(draw);
680: PetscDrawPause(draw);
681: return(0);
682: }
684: #undef __FUNCT__
686: int MatView_SeqBDiag(Mat A,PetscViewer viewer)
687: {
688: int ierr;
689: PetscTruth isascii,isbinary,isdraw;
692: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
693: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
694: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
695: if (isascii) {
696: MatView_SeqBDiag_ASCII(A,viewer);
697: } else if (isbinary) {
698: MatView_SeqBDiag_Binary(A,viewer);
699: } else if (isdraw) {
700: MatView_SeqBDiag_Draw(A,viewer);
701: } else {
702: SETERRQ1(1,"Viewer type %s not supported by BDiag matrices",((PetscObject)viewer)->type_name);
703: }
704: return(0);
705: }