Actual source code: sbaij.c
petsc-3.3-p7 2013-05-11
2: /*
3: Defines the basic matrix operations for the SBAIJ (compressed row)
4: matrix storage format.
5: */
6: #include <../src/mat/impls/baij/seq/baij.h> /*I "petscmat.h" I*/
7: #include <../src/mat/impls/sbaij/seq/sbaij.h>
8: #include <petscblaslapack.h>
10: #include <../src/mat/impls/sbaij/seq/relax.h>
11: #define USESHORT
12: #include <../src/mat/impls/sbaij/seq/relax.h>
14: extern PetscErrorCode MatSeqSBAIJSetNumericFactorization_inplace(Mat,PetscBool );
16: /*
17: Checks for missing diagonals
18: */
21: PetscErrorCode MatMissingDiagonal_SeqSBAIJ(Mat A,PetscBool *missing,PetscInt *dd)
22: {
23: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
25: PetscInt *diag,*jj = a->j,i;
28: MatMarkDiagonal_SeqSBAIJ(A);
29: *missing = PETSC_FALSE;
30: if(A->rmap->n > 0 && !jj) {
31: *missing = PETSC_TRUE;
32: if (dd) *dd = 0;
33: PetscInfo(A,"Matrix has no entries therefore is missing diagonal");
34: } else {
35: diag = a->diag;
36: for (i=0; i<a->mbs; i++) {
37: if (jj[diag[i]] != i) {
38: *missing = PETSC_TRUE;
39: if (dd) *dd = i;
40: break;
41: }
42: }
43: }
44: return(0);
45: }
49: PetscErrorCode MatMarkDiagonal_SeqSBAIJ(Mat A)
50: {
51: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
53: PetscInt i;
56: if (!a->diag) {
57: PetscMalloc(a->mbs*sizeof(PetscInt),&a->diag);
58: PetscLogObjectMemory(A,a->mbs*sizeof(PetscInt));
59: a->free_diag = PETSC_TRUE;
60: }
61: for (i=0; i<a->mbs; i++) a->diag[i] = a->i[i];
62: return(0);
63: }
67: static PetscErrorCode MatGetRowIJ_SeqSBAIJ(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscBool *done)
68: {
69: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
70: PetscInt i,j,n = a->mbs,nz = a->i[n],bs = A->rmap->bs;
74: *nn = n;
75: if (!ia) return(0);
76: if (!blockcompressed) {
77: /* malloc & create the natural set of indices */
78: PetscMalloc2((n+1)*bs,PetscInt,ia,nz*bs,PetscInt,ja);
79: for (i=0; i<n+1; i++) {
80: for (j=0; j<bs; j++) {
81: *ia[i*bs+j] = a->i[i]*bs+j+oshift;
82: }
83: }
84: for (i=0; i<nz; i++) {
85: for (j=0; j<bs; j++) {
86: *ja[i*bs+j] = a->j[i]*bs+j+oshift;
87: }
88: }
89: } else { /* blockcompressed */
90: if (oshift == 1) {
91: /* temporarily add 1 to i and j indices */
92: for (i=0; i<nz; i++) a->j[i]++;
93: for (i=0; i<n+1; i++) a->i[i]++;
94: }
95: *ia = a->i; *ja = a->j;
96: }
98: return(0);
99: }
103: static PetscErrorCode MatRestoreRowIJ_SeqSBAIJ(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscBool *done)
104: {
105: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
106: PetscInt i,n = a->mbs,nz = a->i[n];
110: if (!ia) return(0);
112: if (!blockcompressed) {
113: PetscFree2(*ia,*ja);
114: } else if (oshift == 1) { /* blockcompressed */
115: for (i=0; i<nz; i++) a->j[i]--;
116: for (i=0; i<n+1; i++) a->i[i]--;
117: }
119: return(0);
120: }
124: PetscErrorCode MatDestroy_SeqSBAIJ(Mat A)
125: {
126: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
130: #if defined(PETSC_USE_LOG)
131: PetscLogObjectState((PetscObject)A,"Rows=%D, NZ=%D",A->rmap->N,a->nz);
132: #endif
133: MatSeqXAIJFreeAIJ(A,&a->a,&a->j,&a->i);
134: if (a->free_diag){PetscFree(a->diag);}
135: ISDestroy(&a->row);
136: ISDestroy(&a->col);
137: ISDestroy(&a->icol);
138: PetscFree(a->idiag);
139: PetscFree(a->inode.size);
140: if (a->free_imax_ilen) {PetscFree2(a->imax,a->ilen);}
141: PetscFree(a->solve_work);
142: PetscFree(a->sor_work);
143: PetscFree(a->solves_work);
144: PetscFree(a->mult_work);
145: PetscFree(a->saved_values);
146: PetscFree(a->xtoy);
147: if (a->free_jshort) {PetscFree(a->jshort);}
148: PetscFree(a->inew);
149: MatDestroy(&a->parent);
150: PetscFree(A->data);
152: PetscObjectChangeTypeName((PetscObject)A,0);
153: PetscObjectComposeFunction((PetscObject)A,"MatStoreValues_C","",PETSC_NULL);
154: PetscObjectComposeFunction((PetscObject)A,"MatRetrieveValues_C","",PETSC_NULL);
155: PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetColumnIndices_C","",PETSC_NULL);
156: PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqaij_C","",PETSC_NULL);
157: PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqbaij_C","",PETSC_NULL);
158: PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetPreallocation_C","",PETSC_NULL);
159: PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqsbstrm_C","",PETSC_NULL);
160: return(0);
161: }
165: PetscErrorCode MatSetOption_SeqSBAIJ(Mat A,MatOption op,PetscBool flg)
166: {
167: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
171: switch (op) {
172: case MAT_ROW_ORIENTED:
173: a->roworiented = flg;
174: break;
175: case MAT_KEEP_NONZERO_PATTERN:
176: a->keepnonzeropattern = flg;
177: break;
178: case MAT_NEW_NONZERO_LOCATIONS:
179: a->nonew = (flg ? 0 : 1);
180: break;
181: case MAT_NEW_NONZERO_LOCATION_ERR:
182: a->nonew = (flg ? -1 : 0);
183: break;
184: case MAT_NEW_NONZERO_ALLOCATION_ERR:
185: a->nonew = (flg ? -2 : 0);
186: break;
187: case MAT_UNUSED_NONZERO_LOCATION_ERR:
188: a->nounused = (flg ? -1 : 0);
189: break;
190: case MAT_NEW_DIAGONALS:
191: case MAT_IGNORE_OFF_PROC_ENTRIES:
192: case MAT_USE_HASH_TABLE:
193: PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
194: break;
195: case MAT_HERMITIAN:
196: if (!A->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatAssemblyEnd() first");
197: if (A->cmap->n < 65536 && A->cmap->bs == 1) {
198: A->ops->mult = MatMult_SeqSBAIJ_1_Hermitian_ushort;
199: } else if (A->cmap->bs == 1) {
200: A->ops->mult = MatMult_SeqSBAIJ_1_Hermitian;
201: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for Hermitian with block size greater than 1");
202: break;
203: case MAT_SPD:
204: A->spd_set = PETSC_TRUE;
205: A->spd = flg;
206: if (flg) {
207: A->symmetric = PETSC_TRUE;
208: A->structurally_symmetric = PETSC_TRUE;
209: A->symmetric_set = PETSC_TRUE;
210: A->structurally_symmetric_set = PETSC_TRUE;
211: }
212: break;
213: case MAT_SYMMETRIC:
214: case MAT_STRUCTURALLY_SYMMETRIC:
215: case MAT_SYMMETRY_ETERNAL:
216: if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix must be symmetric");
217: PetscInfo1(A,"Option %s not relevent\n",MatOptions[op]);
218: break;
219: case MAT_IGNORE_LOWER_TRIANGULAR:
220: a->ignore_ltriangular = flg;
221: break;
222: case MAT_ERROR_LOWER_TRIANGULAR:
223: a->ignore_ltriangular = flg;
224: break;
225: case MAT_GETROW_UPPERTRIANGULAR:
226: a->getrow_utriangular = flg;
227: break;
228: default:
229: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %d",op);
230: }
231: return(0);
232: }
236: PetscErrorCode MatGetRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **v)
237: {
238: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
240: PetscInt itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*cols_i,bs2;
241: MatScalar *aa,*aa_i;
242: PetscScalar *v_i;
245: if (A && !a->getrow_utriangular) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatGetRow is not supported for SBAIJ matrix format. Getting the upper triangular part of row, run with -mat_getrow_uppertriangular, call MatSetOption(mat,MAT_GETROW_UPPERTRIANGULAR,PETSC_TRUE) or MatGetRowUpperTriangular()");
246: /* Get the upper triangular part of the row */
247: bs = A->rmap->bs;
248: ai = a->i;
249: aj = a->j;
250: aa = a->a;
251: bs2 = a->bs2;
252:
253: if (row < 0 || row >= A->rmap->N) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Row %D out of range", row);
254:
255: bn = row/bs; /* Block number */
256: bp = row % bs; /* Block position */
257: M = ai[bn+1] - ai[bn];
258: *ncols = bs*M;
259:
260: if (v) {
261: *v = 0;
262: if (*ncols) {
263: PetscMalloc((*ncols+row)*sizeof(PetscScalar),v);
264: for (i=0; i<M; i++) { /* for each block in the block row */
265: v_i = *v + i*bs;
266: aa_i = aa + bs2*(ai[bn] + i);
267: for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];}
268: }
269: }
270: }
271:
272: if (cols) {
273: *cols = 0;
274: if (*ncols) {
275: PetscMalloc((*ncols+row)*sizeof(PetscInt),cols);
276: for (i=0; i<M; i++) { /* for each block in the block row */
277: cols_i = *cols + i*bs;
278: itmp = bs*aj[ai[bn] + i];
279: for (j=0; j<bs; j++) {cols_i[j] = itmp++;}
280: }
281: }
282: }
283:
284: /*search column A(0:row-1,row) (=A(row,0:row-1)). Could be expensive! */
285: /* this segment is currently removed, so only entries in the upper triangle are obtained */
286: #ifdef column_search
287: v_i = *v + M*bs;
288: cols_i = *cols + M*bs;
289: for (i=0; i<bn; i++){ /* for each block row */
290: M = ai[i+1] - ai[i];
291: for (j=0; j<M; j++){
292: itmp = aj[ai[i] + j]; /* block column value */
293: if (itmp == bn){
294: aa_i = aa + bs2*(ai[i] + j) + bs*bp;
295: for (k=0; k<bs; k++) {
296: *cols_i++ = i*bs+k;
297: *v_i++ = aa_i[k];
298: }
299: *ncols += bs;
300: break;
301: }
302: }
303: }
304: #endif
305: return(0);
306: }
310: PetscErrorCode MatRestoreRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
311: {
313:
315: if (idx) {PetscFree(*idx);}
316: if (v) {PetscFree(*v);}
317: return(0);
318: }
322: PetscErrorCode MatGetRowUpperTriangular_SeqSBAIJ(Mat A)
323: {
324: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
327: a->getrow_utriangular = PETSC_TRUE;
328: return(0);
329: }
332: PetscErrorCode MatRestoreRowUpperTriangular_SeqSBAIJ(Mat A)
333: {
334: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
337: a->getrow_utriangular = PETSC_FALSE;
338: return(0);
339: }
343: PetscErrorCode MatTranspose_SeqSBAIJ(Mat A,MatReuse reuse,Mat *B)
344: {
347: if (reuse == MAT_INITIAL_MATRIX || *B != A) {
348: MatDuplicate(A,MAT_COPY_VALUES,B);
349: }
350: return(0);
351: }
355: static PetscErrorCode MatView_SeqSBAIJ_ASCII(Mat A,PetscViewer viewer)
356: {
357: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
358: PetscErrorCode ierr;
359: PetscInt i,j,bs = A->rmap->bs,k,l,bs2=a->bs2;
360: PetscViewerFormat format;
361: PetscInt *diag;
362:
364: PetscViewerGetFormat(viewer,&format);
365: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
366: PetscViewerASCIIPrintf(viewer," block size is %D\n",bs);
367: } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
368: Mat aij;
369: if (A->factortype && bs>1){
370: PetscPrintf(PETSC_COMM_SELF,"Warning: matrix is factored with bs>1. MatView() with PETSC_VIEWER_ASCII_MATLAB is not supported and ignored!\n");
371: return(0);
372: }
373: MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&aij);
374: MatView(aij,viewer);
375: MatDestroy(&aij);
376: } else if (format == PETSC_VIEWER_ASCII_COMMON) {
377: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
378: for (i=0; i<a->mbs; i++) {
379: for (j=0; j<bs; j++) {
380: PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);
381: for (k=a->i[i]; k<a->i[i+1]; k++) {
382: for (l=0; l<bs; l++) {
383: #if defined(PETSC_USE_COMPLEX)
384: if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
385: PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",bs*a->j[k]+l,
386: PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));
387: } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
388: PetscViewerASCIIPrintf(viewer," (%D, %G - %G i) ",bs*a->j[k]+l,
389: PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));
390: } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
391: PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));
392: }
393: #else
394: if (a->a[bs2*k + l*bs + j] != 0.0) {
395: PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);
396: }
397: #endif
398: }
399: }
400: PetscViewerASCIIPrintf(viewer,"\n");
401: }
402: }
403: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
404: } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
405: return(0);
406: } else {
407: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
408: PetscObjectPrintClassNamePrefixType((PetscObject)A,viewer,"Matrix Object");
409: if (A->factortype){ /* for factored matrix */
410: if (bs>1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"matrix is factored with bs>1. Not implemented yet");
412: diag=a->diag;
413: for (i=0; i<a->mbs; i++) { /* for row block i */
414: PetscViewerASCIIPrintf(viewer,"row %D:",i);
415: /* diagonal entry */
416: #if defined(PETSC_USE_COMPLEX)
417: if (PetscImaginaryPart(a->a[diag[i]]) > 0.0) {
418: PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",a->j[diag[i]],PetscRealPart(1.0/a->a[diag[i]]),PetscImaginaryPart(1.0/a->a[diag[i]]));
419: } else if (PetscImaginaryPart(a->a[diag[i]]) < 0.0) {
420: PetscViewerASCIIPrintf(viewer," (%D, %G - %G i) ",a->j[diag[i]],PetscRealPart(1.0/a->a[diag[i]]),-PetscImaginaryPart(1.0/a->a[diag[i]]));
421: } else {
422: PetscViewerASCIIPrintf(viewer," (%D, %G) ",a->j[diag[i]],PetscRealPart(1.0/a->a[diag[i]]));
423: }
424: #else
425: PetscViewerASCIIPrintf(viewer," (%D, %G) ",a->j[diag[i]],1.0/a->a[diag[i]]);
426: #endif
427: /* off-diagonal entries */
428: for (k=a->i[i]; k<a->i[i+1]-1; k++) {
429: #if defined(PETSC_USE_COMPLEX)
430: if (PetscImaginaryPart(a->a[k]) > 0.0) {
431: PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",bs*a->j[k],PetscRealPart(a->a[k]),PetscImaginaryPart(a->a[k]));
432: } else if (PetscImaginaryPart(a->a[k]) < 0.0) {
433: PetscViewerASCIIPrintf(viewer," (%D, %G - %G i) ",bs*a->j[k],PetscRealPart(a->a[k]),-PetscImaginaryPart(a->a[k]));
434: } else {
435: PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k],PetscRealPart(a->a[k]));
436: }
437: #else
438: PetscViewerASCIIPrintf(viewer," (%D, %G) ",a->j[k],a->a[k]);
439: #endif
440: }
441: PetscViewerASCIIPrintf(viewer,"\n");
442: }
443:
444: } else { /* for non-factored matrix */
445: for (i=0; i<a->mbs; i++) { /* for row block i */
446: for (j=0; j<bs; j++) { /* for row bs*i + j */
447: PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);
448: for (k=a->i[i]; k<a->i[i+1]; k++) { /* for column block */
449: for (l=0; l<bs; l++) { /* for column */
450: #if defined(PETSC_USE_COMPLEX)
451: if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) {
452: PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",bs*a->j[k]+l,
453: PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));
454: } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) {
455: PetscViewerASCIIPrintf(viewer," (%D, %G - %G i) ",bs*a->j[k]+l,
456: PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));
457: } else {
458: PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));
459: }
460: #else
461: PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);
462: #endif
463: }
464: }
465: PetscViewerASCIIPrintf(viewer,"\n");
466: }
467: }
468: }
469: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
470: }
471: PetscViewerFlush(viewer);
472: return(0);
473: }
477: static PetscErrorCode MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw,void *Aa)
478: {
479: Mat A = (Mat) Aa;
480: Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data;
482: PetscInt row,i,j,k,l,mbs=a->mbs,color,bs=A->rmap->bs,bs2=a->bs2;
483: PetscMPIInt rank;
484: PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r;
485: MatScalar *aa;
486: MPI_Comm comm;
487: PetscViewer viewer;
488:
490: /*
491: This is nasty. If this is called from an originally parallel matrix
492: then all processes call this,but only the first has the matrix so the
493: rest should return immediately.
494: */
495: PetscObjectGetComm((PetscObject)draw,&comm);
496: MPI_Comm_rank(comm,&rank);
497: if (rank) return(0);
498:
499: PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);
500:
501: PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);
502: PetscDrawString(draw, .3*(xl+xr), .3*(yl+yr), PETSC_DRAW_BLACK, "symmetric");
503:
504: /* loop over matrix elements drawing boxes */
505: color = PETSC_DRAW_BLUE;
506: for (i=0,row=0; i<mbs; i++,row+=bs) {
507: for (j=a->i[i]; j<a->i[i+1]; j++) {
508: y_l = A->rmap->N - row - 1.0; y_r = y_l + 1.0;
509: x_l = a->j[j]*bs; x_r = x_l + 1.0;
510: aa = a->a + j*bs2;
511: for (k=0; k<bs; k++) {
512: for (l=0; l<bs; l++) {
513: if (PetscRealPart(*aa++) >= 0.) continue;
514: PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
515: }
516: }
517: }
518: }
519: color = PETSC_DRAW_CYAN;
520: for (i=0,row=0; i<mbs; i++,row+=bs) {
521: for (j=a->i[i]; j<a->i[i+1]; j++) {
522: y_l = A->rmap->N - row - 1.0; y_r = y_l + 1.0;
523: x_l = a->j[j]*bs; x_r = x_l + 1.0;
524: aa = a->a + j*bs2;
525: for (k=0; k<bs; k++) {
526: for (l=0; l<bs; l++) {
527: if (PetscRealPart(*aa++) != 0.) continue;
528: PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
529: }
530: }
531: }
532: }
533:
534: color = PETSC_DRAW_RED;
535: for (i=0,row=0; i<mbs; i++,row+=bs) {
536: for (j=a->i[i]; j<a->i[i+1]; j++) {
537: y_l = A->rmap->N - row - 1.0; y_r = y_l + 1.0;
538: x_l = a->j[j]*bs; x_r = x_l + 1.0;
539: aa = a->a + j*bs2;
540: for (k=0; k<bs; k++) {
541: for (l=0; l<bs; l++) {
542: if (PetscRealPart(*aa++) <= 0.) continue;
543: PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
544: }
545: }
546: }
547: }
548: return(0);
549: }
553: static PetscErrorCode MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer)
554: {
556: PetscReal xl,yl,xr,yr,w,h;
557: PetscDraw draw;
558: PetscBool isnull;
559:
561: PetscViewerDrawGetDraw(viewer,0,&draw);
562: PetscDrawIsNull(draw,&isnull); if (isnull) return(0);
563:
564: PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);
565: xr = A->rmap->N; yr = A->rmap->N; h = yr/10.0; w = xr/10.0;
566: xr += w; yr += h; xl = -w; yl = -h;
567: PetscDrawSetCoordinates(draw,xl,yl,xr,yr);
568: PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);
569: PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);
570: return(0);
571: }
575: PetscErrorCode MatView_SeqSBAIJ(Mat A,PetscViewer viewer)
576: {
578: PetscBool iascii,isdraw;
579: FILE *file = 0;
582: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
583: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
584: if (iascii){
585: MatView_SeqSBAIJ_ASCII(A,viewer);
586: } else if (isdraw) {
587: MatView_SeqSBAIJ_Draw(A,viewer);
588: } else {
589: Mat B;
590: MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&B);
591: MatView(B,viewer);
592: MatDestroy(&B);
593: PetscViewerBinaryGetInfoPointer(viewer,&file);
594: if (file) {
595: fprintf(file,"-matload_block_size %d\n",(int)A->rmap->bs);
596: }
597: }
598: return(0);
599: }
604: PetscErrorCode MatGetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],PetscScalar v[])
605: {
606: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
607: PetscInt *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
608: PetscInt *ai = a->i,*ailen = a->ilen;
609: PetscInt brow,bcol,ridx,cidx,bs=A->rmap->bs,bs2=a->bs2;
610: MatScalar *ap,*aa = a->a;
611:
613: for (k=0; k<m; k++) { /* loop over rows */
614: row = im[k]; brow = row/bs;
615: if (row < 0) {v += n; continue;} /* SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",row); */
616: if (row >= A->rmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap->N-1);
617: rp = aj + ai[brow] ; ap = aa + bs2*ai[brow] ;
618: nrow = ailen[brow];
619: for (l=0; l<n; l++) { /* loop over columns */
620: if (in[l] < 0) {v++; continue;} /* SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %D",in[l]); */
621: if (in[l] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->cmap->n-1);
622: col = in[l] ;
623: bcol = col/bs;
624: cidx = col%bs;
625: ridx = row%bs;
626: high = nrow;
627: low = 0; /* assume unsorted */
628: while (high-low > 5) {
629: t = (low+high)/2;
630: if (rp[t] > bcol) high = t;
631: else low = t;
632: }
633: for (i=low; i<high; i++) {
634: if (rp[i] > bcol) break;
635: if (rp[i] == bcol) {
636: *v++ = ap[bs2*i+bs*cidx+ridx];
637: goto finished;
638: }
639: }
640: *v++ = 0.0;
641: finished:;
642: }
643: }
644: return(0);
645: }
650: PetscErrorCode MatSetValuesBlocked_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
651: {
652: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
653: PetscErrorCode ierr;
654: PetscInt *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,lastcol = -1;
655: PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen;
656: PetscInt *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=A->rmap->bs,stepval;
657: PetscBool roworiented=a->roworiented;
658: const PetscScalar *value = v;
659: MatScalar *ap,*aa = a->a,*bap;
660:
662: if (roworiented) {
663: stepval = (n-1)*bs;
664: } else {
665: stepval = (m-1)*bs;
666: }
667: for (k=0; k<m; k++) { /* loop over added rows */
668: row = im[k];
669: if (row < 0) continue;
670: #if defined(PETSC_USE_DEBUG)
671: if (row >= a->mbs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,a->mbs-1);
672: #endif
673: rp = aj + ai[row];
674: ap = aa + bs2*ai[row];
675: rmax = imax[row];
676: nrow = ailen[row];
677: low = 0;
678: high = nrow;
679: for (l=0; l<n; l++) { /* loop over added columns */
680: if (in[l] < 0) continue;
681: col = in[l];
682: #if defined(PETSC_USE_DEBUG)
683: if (col >= a->nbs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",col,a->nbs-1);
684: #endif
685: if (col < row) {
686: if (a->ignore_ltriangular) {
687: continue; /* ignore lower triangular block */
688: } else {
689: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Lower triangular value cannot be set for sbaij format. Ignoring these values, run with -mat_ignore_lower_triangular or call MatSetOption(mat,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE)");
690: }
691: }
692: if (roworiented) {
693: value = v + k*(stepval+bs)*bs + l*bs;
694: } else {
695: value = v + l*(stepval+bs)*bs + k*bs;
696: }
697: if (col <= lastcol) low = 0; else high = nrow;
698: lastcol = col;
699: while (high-low > 7) {
700: t = (low+high)/2;
701: if (rp[t] > col) high = t;
702: else low = t;
703: }
704: for (i=low; i<high; i++) {
705: if (rp[i] > col) break;
706: if (rp[i] == col) {
707: bap = ap + bs2*i;
708: if (roworiented) {
709: if (is == ADD_VALUES) {
710: for (ii=0; ii<bs; ii++,value+=stepval) {
711: for (jj=ii; jj<bs2; jj+=bs) {
712: bap[jj] += *value++;
713: }
714: }
715: } else {
716: for (ii=0; ii<bs; ii++,value+=stepval) {
717: for (jj=ii; jj<bs2; jj+=bs) {
718: bap[jj] = *value++;
719: }
720: }
721: }
722: } else {
723: if (is == ADD_VALUES) {
724: for (ii=0; ii<bs; ii++,value+=stepval) {
725: for (jj=0; jj<bs; jj++) {
726: *bap++ += *value++;
727: }
728: }
729: } else {
730: for (ii=0; ii<bs; ii++,value+=stepval) {
731: for (jj=0; jj<bs; jj++) {
732: *bap++ = *value++;
733: }
734: }
735: }
736: }
737: goto noinsert2;
738: }
739: }
740: if (nonew == 1) goto noinsert2;
741: if (nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col);
742: MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar);
743: N = nrow++ - 1; high++;
744: /* shift up all the later entries in this row */
745: for (ii=N; ii>=i; ii--) {
746: rp[ii+1] = rp[ii];
747: PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));
748: }
749: if (N >= i) {
750: PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));
751: }
752: rp[i] = col;
753: bap = ap + bs2*i;
754: if (roworiented) {
755: for (ii=0; ii<bs; ii++,value+=stepval) {
756: for (jj=ii; jj<bs2; jj+=bs) {
757: bap[jj] = *value++;
758: }
759: }
760: } else {
761: for (ii=0; ii<bs; ii++,value+=stepval) {
762: for (jj=0; jj<bs; jj++) {
763: *bap++ = *value++;
764: }
765: }
766: }
767: noinsert2:;
768: low = i;
769: }
770: ailen[row] = nrow;
771: }
772: return(0);
773: }
775: /*
776: This is not yet used
777: */
780: PetscErrorCode MatAssemblyEnd_SeqSBAIJ_SeqAIJ_Inode(Mat A)
781: {
782: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
783: PetscErrorCode ierr;
784: const PetscInt *ai = a->i, *aj = a->j,*cols;
785: PetscInt i = 0,j,blk_size,m = A->rmap->n,node_count = 0,nzx,nzy,*ns,row,nz,cnt,cnt2,*counts;
786: PetscBool flag;
789: PetscMalloc(m*sizeof(PetscInt),&ns);
790: while (i < m){
791: nzx = ai[i+1] - ai[i]; /* Number of nonzeros */
792: /* Limits the number of elements in a node to 'a->inode.limit' */
793: for (j=i+1,blk_size=1; j<m && blk_size <a->inode.limit; ++j,++blk_size) {
794: nzy = ai[j+1] - ai[j];
795: if (nzy != (nzx - j + i)) break;
796: PetscMemcmp(aj + ai[i] + j - i,aj + ai[j],nzy*sizeof(PetscInt),&flag);
797: if (!flag) break;
798: }
799: ns[node_count++] = blk_size;
800: i = j;
801: }
802: if (!a->inode.size && m && node_count > .9*m) {
803: PetscFree(ns);
804: PetscInfo2(A,"Found %D nodes out of %D rows. Not using Inode routines\n",node_count,m);
805: } else {
806: a->inode.node_count = node_count;
807: PetscMalloc(node_count*sizeof(PetscInt),&a->inode.size);
808: PetscLogObjectMemory(A,node_count*sizeof(PetscInt));
809: PetscMemcpy(a->inode.size,ns,node_count*sizeof(PetscInt));
810: PetscFree(ns);
811: PetscInfo3(A,"Found %D nodes of %D. Limit used: %D. Using Inode routines\n",node_count,m,a->inode.limit);
812:
813: /* count collections of adjacent columns in each inode */
814: row = 0;
815: cnt = 0;
816: for (i=0; i<node_count; i++) {
817: cols = aj + ai[row] + a->inode.size[i];
818: nz = ai[row+1] - ai[row] - a->inode.size[i];
819: for (j=1; j<nz; j++) {
820: if (cols[j] != cols[j-1]+1) {
821: cnt++;
822: }
823: }
824: cnt++;
825: row += a->inode.size[i];
826: }
827: PetscMalloc(2*cnt*sizeof(PetscInt),&counts);
828: cnt = 0;
829: row = 0;
830: for (i=0; i<node_count; i++) {
831: cols = aj + ai[row] + a->inode.size[i];
832: CHKMEMQ;
833: counts[2*cnt] = cols[0];
834: CHKMEMQ;
835: nz = ai[row+1] - ai[row] - a->inode.size[i];
836: cnt2 = 1;
837: for (j=1; j<nz; j++) {
838: if (cols[j] != cols[j-1]+1) {
839: CHKMEMQ;
840: counts[2*(cnt++)+1] = cnt2;
841: counts[2*cnt] = cols[j];
842: CHKMEMQ;
843: cnt2 = 1;
844: } else cnt2++;
845: }
846: CHKMEMQ;
847: counts[2*(cnt++)+1] = cnt2;
848: CHKMEMQ;
849: row += a->inode.size[i];
850: }
851: PetscIntView(2*cnt,counts,0);
852: }
853: return(0);
854: }
858: PetscErrorCode MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode)
859: {
860: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
862: PetscInt fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
863: PetscInt m = A->rmap->N,*ip,N,*ailen = a->ilen;
864: PetscInt mbs = a->mbs,bs2 = a->bs2,rmax = 0;
865: MatScalar *aa = a->a,*ap;
866:
868: if (mode == MAT_FLUSH_ASSEMBLY) return(0);
869:
870: if (m) rmax = ailen[0];
871: for (i=1; i<mbs; i++) {
872: /* move each row back by the amount of empty slots (fshift) before it*/
873: fshift += imax[i-1] - ailen[i-1];
874: rmax = PetscMax(rmax,ailen[i]);
875: if (fshift) {
876: ip = aj + ai[i]; ap = aa + bs2*ai[i];
877: N = ailen[i];
878: for (j=0; j<N; j++) {
879: ip[j-fshift] = ip[j];
880: PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));
881: }
882: }
883: ai[i] = ai[i-1] + ailen[i-1];
884: }
885: if (mbs) {
886: fshift += imax[mbs-1] - ailen[mbs-1];
887: ai[mbs] = ai[mbs-1] + ailen[mbs-1];
888: }
889: /* reset ilen and imax for each row */
890: for (i=0; i<mbs; i++) {
891: ailen[i] = imax[i] = ai[i+1] - ai[i];
892: }
893: a->nz = ai[mbs];
894:
895: /* diagonals may have moved, reset it */
896: if (a->diag) {
897: PetscMemcpy(a->diag,ai,mbs*sizeof(PetscInt));
898: }
899: if (fshift && a->nounused == -1) {
900: SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_PLIB, "Unused space detected in matrix: %D X %D block size %D, %D unneeded", m, A->cmap->n, A->rmap->bs, fshift*bs2);
901: }
902: PetscInfo5(A,"Matrix size: %D X %D, block size %D; storage space: %D unneeded, %D used\n",m,A->rmap->N,A->rmap->bs,fshift*bs2,a->nz*bs2);
903: PetscInfo1(A,"Number of mallocs during MatSetValues is %D\n",a->reallocs);
904: PetscInfo1(A,"Most nonzeros blocks in any row is %D\n",rmax);
905: A->info.mallocs += a->reallocs;
906: a->reallocs = 0;
907: A->info.nz_unneeded = (PetscReal)fshift*bs2;
908: a->idiagvalid = PETSC_FALSE;
910: if (A->cmap->n < 65536 && A->cmap->bs == 1) {
911: if (a->jshort && a->free_jshort){
912: /* when matrix data structure is changed, previous jshort must be replaced */
913: PetscFree(a->jshort);
914: }
915: PetscMalloc(a->i[A->rmap->n]*sizeof(unsigned short),&a->jshort);
916: PetscLogObjectMemory(A,a->i[A->rmap->n]*sizeof(unsigned short));
917: for (i=0; i<a->i[A->rmap->n]; i++) a->jshort[i] = a->j[i];
918: A->ops->mult = MatMult_SeqSBAIJ_1_ushort;
919: A->ops->sor = MatSOR_SeqSBAIJ_ushort;
920: a->free_jshort = PETSC_TRUE;
921: }
922: return(0);
923: }
925: /*
926: This function returns an array of flags which indicate the locations of contiguous
927: blocks that should be zeroed. for eg: if bs = 3 and is = [0,1,2,3,5,6,7,8,9]
928: then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)]
929: Assume: sizes should be long enough to hold all the values.
930: */
933: PetscErrorCode MatZeroRows_SeqSBAIJ_Check_Blocks(PetscInt idx[],PetscInt n,PetscInt bs,PetscInt sizes[], PetscInt *bs_max)
934: {
935: PetscInt i,j,k,row;
936: PetscBool flg;
937:
939: for (i=0,j=0; i<n; j++) {
940: row = idx[i];
941: if (row%bs!=0) { /* Not the begining of a block */
942: sizes[j] = 1;
943: i++;
944: } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */
945: sizes[j] = 1; /* Also makes sure atleast 'bs' values exist for next else */
946: i++;
947: } else { /* Begining of the block, so check if the complete block exists */
948: flg = PETSC_TRUE;
949: for (k=1; k<bs; k++) {
950: if (row+k != idx[i+k]) { /* break in the block */
951: flg = PETSC_FALSE;
952: break;
953: }
954: }
955: if (flg) { /* No break in the bs */
956: sizes[j] = bs;
957: i+= bs;
958: } else {
959: sizes[j] = 1;
960: i++;
961: }
962: }
963: }
964: *bs_max = j;
965: return(0);
966: }
969: /* Only add/insert a(i,j) with i<=j (blocks).
970: Any a(i,j) with i>j input by user is ingored.
971: */
975: PetscErrorCode MatSetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
976: {
977: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
979: PetscInt *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,lastcol = -1;
980: PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented;
981: PetscInt *aj=a->j,nonew=a->nonew,bs=A->rmap->bs,brow,bcol;
982: PetscInt ridx,cidx,bs2=a->bs2;
983: MatScalar *ap,value,*aa=a->a,*bap;
984:
987: for (k=0; k<m; k++) { /* loop over added rows */
988: row = im[k]; /* row number */
989: brow = row/bs; /* block row number */
990: if (row < 0) continue;
991: #if defined(PETSC_USE_DEBUG)
992: if (row >= A->rmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap->N-1);
993: #endif
994: rp = aj + ai[brow]; /*ptr to beginning of column value of the row block*/
995: ap = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/
996: rmax = imax[brow]; /* maximum space allocated for this row */
997: nrow = ailen[brow]; /* actual length of this row */
998: low = 0;
999:
1000: for (l=0; l<n; l++) { /* loop over added columns */
1001: if (in[l] < 0) continue;
1002: #if defined(PETSC_USE_DEBUG)
1003: if (in[l] >= A->rmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->rmap->N-1);
1004: #endif
1005: col = in[l];
1006: bcol = col/bs; /* block col number */
1007:
1008: if (brow > bcol) {
1009: if (a->ignore_ltriangular){
1010: continue; /* ignore lower triangular values */
1011: } else {
1012: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Lower triangular value cannot be set for sbaij format. Ignoring these values, run with -mat_ignore_lower_triangular or call MatSetOption(mat,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE)");
1013: }
1014: }
1015:
1016: ridx = row % bs; cidx = col % bs; /*row and col index inside the block */
1017: if ((brow==bcol && ridx<=cidx) || (brow<bcol)){
1018: /* element value a(k,l) */
1019: if (roworiented) {
1020: value = v[l + k*n];
1021: } else {
1022: value = v[k + l*m];
1023: }
1024:
1025: /* move pointer bap to a(k,l) quickly and add/insert value */
1026: if (col <= lastcol) low = 0; high = nrow;
1027: lastcol = col;
1028: while (high-low > 7) {
1029: t = (low+high)/2;
1030: if (rp[t] > bcol) high = t;
1031: else low = t;
1032: }
1033: for (i=low; i<high; i++) {
1034: if (rp[i] > bcol) break;
1035: if (rp[i] == bcol) {
1036: bap = ap + bs2*i + bs*cidx + ridx;
1037: if (is == ADD_VALUES) *bap += value;
1038: else *bap = value;
1039: /* for diag block, add/insert its symmetric element a(cidx,ridx) */
1040: if (brow == bcol && ridx < cidx){
1041: bap = ap + bs2*i + bs*ridx + cidx;
1042: if (is == ADD_VALUES) *bap += value;
1043: else *bap = value;
1044: }
1045: goto noinsert1;
1046: }
1047: }
1048:
1049: if (nonew == 1) goto noinsert1;
1050: if (nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col);
1051: MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar);
1052:
1053: N = nrow++ - 1; high++;
1054: /* shift up all the later entries in this row */
1055: for (ii=N; ii>=i; ii--) {
1056: rp[ii+1] = rp[ii];
1057: PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));
1058: }
1059: if (N>=i) {
1060: PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));
1061: }
1062: rp[i] = bcol;
1063: ap[bs2*i + bs*cidx + ridx] = value;
1064: noinsert1:;
1065: low = i;
1066: }
1067: } /* end of loop over added columns */
1068: ailen[brow] = nrow;
1069: } /* end of loop over added rows */
1070: return(0);
1071: }
1075: PetscErrorCode MatICCFactor_SeqSBAIJ(Mat inA,IS row,const MatFactorInfo *info)
1076: {
1077: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data;
1078: Mat outA;
1080: PetscBool row_identity;
1081:
1083: if (info->levels != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only levels=0 is supported for in-place icc");
1084: ISIdentity(row,&row_identity);
1085: if (!row_identity) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix reordering is not supported");
1086: if (inA->rmap->bs != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix block size %D is not supported",inA->rmap->bs); /* Need to replace MatCholeskyFactorSymbolic_SeqSBAIJ_MSR()! */
1088: outA = inA;
1089: inA->factortype = MAT_FACTOR_ICC;
1090:
1091: MatMarkDiagonal_SeqSBAIJ(inA);
1092: MatSeqSBAIJSetNumericFactorization_inplace(inA,row_identity);
1094: PetscObjectReference((PetscObject)row);
1095: ISDestroy(&a->row);
1096: a->row = row;
1097: PetscObjectReference((PetscObject)row);
1098: ISDestroy(&a->col);
1099: a->col = row;
1100:
1101: /* Create the invert permutation so that it can be used in MatCholeskyFactorNumeric() */
1102: if (a->icol) {ISInvertPermutation(row,PETSC_DECIDE, &a->icol);}
1103: PetscLogObjectParent(inA,a->icol);
1104:
1105: if (!a->solve_work) {
1106: PetscMalloc((inA->rmap->N+inA->rmap->bs)*sizeof(PetscScalar),&a->solve_work);
1107: PetscLogObjectMemory(inA,(inA->rmap->N+inA->rmap->bs)*sizeof(PetscScalar));
1108: }
1109:
1110: MatCholeskyFactorNumeric(outA,inA,info);
1111: return(0);
1112: }
1114: EXTERN_C_BEGIN
1117: PetscErrorCode MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,PetscInt *indices)
1118: {
1119: Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data;
1120: PetscInt i,nz,n;
1124: nz = baij->maxnz;
1125: n = mat->cmap->n;
1126: for (i=0; i<nz; i++) {
1127: baij->j[i] = indices[i];
1128: }
1129: baij->nz = nz;
1130: for (i=0; i<n; i++) {
1131: baij->ilen[i] = baij->imax[i];
1132: }
1133: MatSetOption(mat,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
1134: return(0);
1135: }
1136: EXTERN_C_END
1140: /*@
1141: MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows
1142: in the matrix.
1143:
1144: Input Parameters:
1145: + mat - the SeqSBAIJ matrix
1146: - indices - the column indices
1147:
1148: Level: advanced
1149:
1150: Notes:
1151: This can be called if you have precomputed the nonzero structure of the
1152: matrix and want to provide it to the matrix object to improve the performance
1153: of the MatSetValues() operation.
1154:
1155: You MUST have set the correct numbers of nonzeros per row in the call to
1156: MatCreateSeqSBAIJ(), and the columns indices MUST be sorted.
1157:
1158: MUST be called before any calls to MatSetValues()
1159:
1160: .seealso: MatCreateSeqSBAIJ
1161: @*/
1162: PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat mat,PetscInt *indices)
1163: {
1165:
1169: PetscUseMethod(mat,"MatSeqSBAIJSetColumnIndices_C",(Mat,PetscInt *),(mat,indices));
1170: return(0);
1171: }
1175: PetscErrorCode MatCopy_SeqSBAIJ(Mat A,Mat B,MatStructure str)
1176: {
1180: /* If the two matrices have the same copy implementation, use fast copy. */
1181: if (str == SAME_NONZERO_PATTERN && (A->ops->copy == B->ops->copy)) {
1182: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1183: Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data;
1185: if (a->i[A->rmap->N] != b->i[B->rmap->N]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of nonzeros in two matrices are different");
1186: PetscMemcpy(b->a,a->a,(a->i[A->rmap->N])*sizeof(PetscScalar));
1187: } else {
1188: MatGetRowUpperTriangular(A);
1189: MatCopy_Basic(A,B,str);
1190: MatRestoreRowUpperTriangular(A);
1191: }
1192: return(0);
1193: }
1197: PetscErrorCode MatSetUp_SeqSBAIJ(Mat A)
1198: {
1202: MatSeqSBAIJSetPreallocation_SeqSBAIJ(A,A->rmap->bs,PETSC_DEFAULT,0);
1203: return(0);
1204: }
1208: PetscErrorCode MatGetArray_SeqSBAIJ(Mat A,PetscScalar *array[])
1209: {
1210: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1212: *array = a->a;
1213: return(0);
1214: }
1218: PetscErrorCode MatRestoreArray_SeqSBAIJ(Mat A,PetscScalar *array[])
1219: {
1221: return(0);
1222: }
1226: PetscErrorCode MatAXPY_SeqSBAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str)
1227: {
1228: Mat_SeqSBAIJ *x=(Mat_SeqSBAIJ *)X->data, *y=(Mat_SeqSBAIJ *)Y->data;
1230: PetscInt i,bs=Y->rmap->bs,bs2=bs*bs,j;
1231: PetscBLASInt one = 1;
1232:
1234: if (str == SAME_NONZERO_PATTERN) {
1235: PetscScalar alpha = a;
1236: PetscBLASInt bnz = PetscBLASIntCast(x->nz*bs2);
1237: BLASaxpy_(&bnz,&alpha,x->a,&one,y->a,&one);
1238: } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
1239: if (y->xtoy && y->XtoY != X) {
1240: PetscFree(y->xtoy);
1241: MatDestroy(&y->XtoY);
1242: }
1243: if (!y->xtoy) { /* get xtoy */
1244: MatAXPYGetxtoy_Private(x->mbs,x->i,x->j,PETSC_NULL, y->i,y->j,PETSC_NULL, &y->xtoy);
1245: y->XtoY = X;
1246: }
1247: for (i=0; i<x->nz; i++) {
1248: j = 0;
1249: while (j < bs2){
1250: y->a[bs2*y->xtoy[i]+j] += a*(x->a[bs2*i+j]);
1251: j++;
1252: }
1253: }
1254: PetscInfo3(Y,"ratio of nnz_s(X)/nnz_s(Y): %D/%D = %G\n",bs2*x->nz,bs2*y->nz,(PetscReal)(bs2*x->nz)/(bs2*y->nz));
1255: } else {
1256: MatGetRowUpperTriangular(X);
1257: MatAXPY_Basic(Y,a,X,str);
1258: MatRestoreRowUpperTriangular(X);
1259: }
1260: return(0);
1261: }
1265: PetscErrorCode MatIsSymmetric_SeqSBAIJ(Mat A,PetscReal tol,PetscBool *flg)
1266: {
1268: *flg = PETSC_TRUE;
1269: return(0);
1270: }
1274: PetscErrorCode MatIsStructurallySymmetric_SeqSBAIJ(Mat A,PetscBool *flg)
1275: {
1277: *flg = PETSC_TRUE;
1278: return(0);
1279: }
1283: PetscErrorCode MatIsHermitian_SeqSBAIJ(Mat A,PetscReal tol,PetscBool *flg)
1284: {
1286: *flg = PETSC_FALSE;
1287: return(0);
1288: }
1292: PetscErrorCode MatRealPart_SeqSBAIJ(Mat A)
1293: {
1294: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1295: PetscInt i,nz = a->bs2*a->i[a->mbs];
1296: MatScalar *aa = a->a;
1299: for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]);
1300: return(0);
1301: }
1305: PetscErrorCode MatImaginaryPart_SeqSBAIJ(Mat A)
1306: {
1307: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1308: PetscInt i,nz = a->bs2*a->i[a->mbs];
1309: MatScalar *aa = a->a;
1312: for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]);
1313: return(0);
1314: }
1318: PetscErrorCode MatZeroRowsColumns_SeqSBAIJ(Mat A,PetscInt is_n,const PetscInt is_idx[],PetscScalar diag,Vec x, Vec b)
1319: {
1320: Mat_SeqSBAIJ *baij=(Mat_SeqSBAIJ*)A->data;
1321: PetscErrorCode ierr;
1322: PetscInt i,j,k,count;
1323: PetscInt bs=A->rmap->bs,bs2=baij->bs2,row,col;
1324: PetscScalar zero = 0.0;
1325: MatScalar *aa;
1326: const PetscScalar *xx;
1327: PetscScalar *bb;
1328: PetscBool *zeroed,vecs = PETSC_FALSE;
1331: /* fix right hand side if needed */
1332: if (x && b) {
1333: VecGetArrayRead(x,&xx);
1334: VecGetArray(b,&bb);
1335: vecs = PETSC_TRUE;
1336: }
1337: A->same_nonzero = PETSC_TRUE;
1339: /* zero the columns */
1340: PetscMalloc(A->rmap->n*sizeof(PetscBool),&zeroed);
1341: PetscMemzero(zeroed,A->rmap->n*sizeof(PetscBool));
1342: for (i=0; i<is_n; i++) {
1343: if (is_idx[i] < 0 || is_idx[i] >= A->rmap->N) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range",is_idx[i]);
1344: zeroed[is_idx[i]] = PETSC_TRUE;
1345: }
1346: if (vecs) {
1347: for (i=0; i<A->rmap->N; i++) {
1348: row = i/bs;
1349: for (j=baij->i[row]; j<baij->i[row+1]; j++) {
1350: for (k=0; k<bs; k++) {
1351: col = bs*baij->j[j] + k;
1352: if (col <= i) continue;
1353: aa = ((MatScalar*)(baij->a)) + j*bs2 + (i%bs) + bs*k;
1354: if (!zeroed[i] && zeroed[col]) {
1355: bb[i] -= aa[0]*xx[col];
1356: }
1357: if (zeroed[i] && !zeroed[col]) {
1358: bb[col] -= aa[0]*xx[i];
1359: }
1360: }
1361: }
1362: }
1363: for (i=0; i<is_n; i++) {
1364: bb[is_idx[i]] = diag*xx[is_idx[i]];
1365: }
1366: }
1368: for (i=0; i<A->rmap->N; i++) {
1369: if (!zeroed[i]) {
1370: row = i/bs;
1371: for (j=baij->i[row]; j<baij->i[row+1]; j++) {
1372: for (k=0; k<bs; k++) {
1373: col = bs*baij->j[j] + k;
1374: if (zeroed[col]) {
1375: aa = ((MatScalar*)(baij->a)) + j*bs2 + (i%bs) + bs*k;
1376: aa[0] = 0.0;
1377: }
1378: }
1379: }
1380: }
1381: }
1382: PetscFree(zeroed);
1383: if (vecs) {
1384: VecRestoreArrayRead(x,&xx);
1385: VecRestoreArray(b,&bb);
1386: }
1388: /* zero the rows */
1389: for (i=0; i<is_n; i++) {
1390: row = is_idx[i];
1391: count = (baij->i[row/bs +1] - baij->i[row/bs])*bs;
1392: aa = ((MatScalar*)(baij->a)) + baij->i[row/bs]*bs2 + (row%bs);
1393: for (k=0; k<count; k++) {
1394: aa[0] = zero;
1395: aa += bs;
1396: }
1397: if (diag != 0.0) {
1398: (*A->ops->setvalues)(A,1,&row,1,&row,&diag,INSERT_VALUES);
1399: }
1400: }
1401: MatAssemblyEnd_SeqSBAIJ(A,MAT_FINAL_ASSEMBLY);
1402: return(0);
1403: }
1405: /* -------------------------------------------------------------------*/
1406: static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ,
1407: MatGetRow_SeqSBAIJ,
1408: MatRestoreRow_SeqSBAIJ,
1409: MatMult_SeqSBAIJ_N,
1410: /* 4*/ MatMultAdd_SeqSBAIJ_N,
1411: MatMult_SeqSBAIJ_N, /* transpose versions are same as non-transpose versions */
1412: MatMultAdd_SeqSBAIJ_N,
1413: 0,
1414: 0,
1415: 0,
1416: /*10*/ 0,
1417: 0,
1418: MatCholeskyFactor_SeqSBAIJ,
1419: MatSOR_SeqSBAIJ,
1420: MatTranspose_SeqSBAIJ,
1421: /*15*/ MatGetInfo_SeqSBAIJ,
1422: MatEqual_SeqSBAIJ,
1423: MatGetDiagonal_SeqSBAIJ,
1424: MatDiagonalScale_SeqSBAIJ,
1425: MatNorm_SeqSBAIJ,
1426: /*20*/ 0,
1427: MatAssemblyEnd_SeqSBAIJ,
1428: MatSetOption_SeqSBAIJ,
1429: MatZeroEntries_SeqSBAIJ,
1430: /*24*/ 0,
1431: 0,
1432: 0,
1433: 0,
1434: 0,
1435: /*29*/ MatSetUp_SeqSBAIJ,
1436: 0,
1437: 0,
1438: MatGetArray_SeqSBAIJ,
1439: MatRestoreArray_SeqSBAIJ,
1440: /*34*/ MatDuplicate_SeqSBAIJ,
1441: 0,
1442: 0,
1443: 0,
1444: MatICCFactor_SeqSBAIJ,
1445: /*39*/ MatAXPY_SeqSBAIJ,
1446: MatGetSubMatrices_SeqSBAIJ,
1447: MatIncreaseOverlap_SeqSBAIJ,
1448: MatGetValues_SeqSBAIJ,
1449: MatCopy_SeqSBAIJ,
1450: /*44*/ 0,
1451: MatScale_SeqSBAIJ,
1452: 0,
1453: 0,
1454: MatZeroRowsColumns_SeqSBAIJ,
1455: /*49*/ 0,
1456: MatGetRowIJ_SeqSBAIJ,
1457: MatRestoreRowIJ_SeqSBAIJ,
1458: 0,
1459: 0,
1460: /*54*/ 0,
1461: 0,
1462: 0,
1463: 0,
1464: MatSetValuesBlocked_SeqSBAIJ,
1465: /*59*/ MatGetSubMatrix_SeqSBAIJ,
1466: 0,
1467: 0,
1468: 0,
1469: 0,
1470: /*64*/ 0,
1471: 0,
1472: 0,
1473: 0,
1474: 0,
1475: /*69*/ MatGetRowMaxAbs_SeqSBAIJ,
1476: 0,
1477: 0,
1478: 0,
1479: 0,
1480: /*74*/ 0,
1481: 0,
1482: 0,
1483: 0,
1484: 0,
1485: /*79*/ 0,
1486: 0,
1487: 0,
1488: MatGetInertia_SeqSBAIJ,
1489: MatLoad_SeqSBAIJ,
1490: /*84*/ MatIsSymmetric_SeqSBAIJ,
1491: MatIsHermitian_SeqSBAIJ,
1492: MatIsStructurallySymmetric_SeqSBAIJ,
1493: 0,
1494: 0,
1495: /*89*/ 0,
1496: 0,
1497: 0,
1498: 0,
1499: 0,
1500: /*94*/ 0,
1501: 0,
1502: 0,
1503: 0,
1504: 0,
1505: /*99*/ 0,
1506: 0,
1507: 0,
1508: 0,
1509: 0,
1510: /*104*/0,
1511: MatRealPart_SeqSBAIJ,
1512: MatImaginaryPart_SeqSBAIJ,
1513: MatGetRowUpperTriangular_SeqSBAIJ,
1514: MatRestoreRowUpperTriangular_SeqSBAIJ,
1515: /*109*/0,
1516: 0,
1517: 0,
1518: 0,
1519: MatMissingDiagonal_SeqSBAIJ,
1520: /*114*/0,
1521: 0,
1522: 0,
1523: 0,
1524: 0,
1525: /*119*/0,
1526: 0,
1527: 0,
1528: 0
1529: };
1531: EXTERN_C_BEGIN
1534: PetscErrorCode MatStoreValues_SeqSBAIJ(Mat mat)
1535: {
1536: Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data;
1537: PetscInt nz = aij->i[mat->rmap->N]*mat->rmap->bs*aij->bs2;
1539:
1541: if (aij->nonew != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first");
1542:
1543: /* allocate space for values if not already there */
1544: if (!aij->saved_values) {
1545: PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);
1546: }
1547:
1548: /* copy values over */
1549: PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));
1550: return(0);
1551: }
1552: EXTERN_C_END
1554: EXTERN_C_BEGIN
1557: PetscErrorCode MatRetrieveValues_SeqSBAIJ(Mat mat)
1558: {
1559: Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data;
1561: PetscInt nz = aij->i[mat->rmap->N]*mat->rmap->bs*aij->bs2;
1562:
1564: if (aij->nonew != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first");
1565: if (!aij->saved_values) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Must call MatStoreValues(A);first");
1566:
1567: /* copy values over */
1568: PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));
1569: return(0);
1570: }
1571: EXTERN_C_END
1573: EXTERN_C_BEGIN
1576: PetscErrorCode MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat B,PetscInt bs,PetscInt nz,PetscInt *nnz)
1577: {
1578: Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data;
1580: PetscInt i,mbs,bs2;
1581: PetscBool skipallocation = PETSC_FALSE,flg = PETSC_FALSE,realalloc = PETSC_FALSE;
1584: if (nz >= 0 || nnz) realalloc = PETSC_TRUE;
1585: B->preallocated = PETSC_TRUE;
1587: PetscLayoutSetBlockSize(B->rmap,bs);
1588: PetscLayoutSetBlockSize(B->cmap,bs);
1589: PetscLayoutSetUp(B->rmap);
1590: PetscLayoutSetUp(B->cmap);
1591: PetscLayoutGetBlockSize(B->rmap,&bs);
1593: mbs = B->rmap->N/bs;
1594: bs2 = bs*bs;
1595:
1596: if (mbs*bs != B->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize");
1597:
1598: if (nz == MAT_SKIP_ALLOCATION) {
1599: skipallocation = PETSC_TRUE;
1600: nz = 0;
1601: }
1603: if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3;
1604: if (nz < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %D",nz);
1605: if (nnz) {
1606: for (i=0; i<mbs; i++) {
1607: if (nnz[i] < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %D value %D",i,nnz[i]);
1608: if (nnz[i] > mbs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be greater than block row length: local row %D value %D rowlength %D",i,nnz[i],mbs);
1609: }
1610: }
1611:
1612: B->ops->mult = MatMult_SeqSBAIJ_N;
1613: B->ops->multadd = MatMultAdd_SeqSBAIJ_N;
1614: B->ops->multtranspose = MatMult_SeqSBAIJ_N;
1615: B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_N;
1616: PetscOptionsGetBool(((PetscObject)B)->prefix,"-mat_no_unroll",&flg,PETSC_NULL);
1617: if (!flg) {
1618: switch (bs) {
1619: case 1:
1620: B->ops->mult = MatMult_SeqSBAIJ_1;
1621: B->ops->multadd = MatMultAdd_SeqSBAIJ_1;
1622: B->ops->multtranspose = MatMult_SeqSBAIJ_1;
1623: B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_1;
1624: break;
1625: case 2:
1626: B->ops->mult = MatMult_SeqSBAIJ_2;
1627: B->ops->multadd = MatMultAdd_SeqSBAIJ_2;
1628: B->ops->multtranspose = MatMult_SeqSBAIJ_2;
1629: B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_2;
1630: break;
1631: case 3:
1632: B->ops->mult = MatMult_SeqSBAIJ_3;
1633: B->ops->multadd = MatMultAdd_SeqSBAIJ_3;
1634: B->ops->multtranspose = MatMult_SeqSBAIJ_3;
1635: B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_3;
1636: break;
1637: case 4:
1638: B->ops->mult = MatMult_SeqSBAIJ_4;
1639: B->ops->multadd = MatMultAdd_SeqSBAIJ_4;
1640: B->ops->multtranspose = MatMult_SeqSBAIJ_4;
1641: B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_4;
1642: break;
1643: case 5:
1644: B->ops->mult = MatMult_SeqSBAIJ_5;
1645: B->ops->multadd = MatMultAdd_SeqSBAIJ_5;
1646: B->ops->multtranspose = MatMult_SeqSBAIJ_5;
1647: B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_5;
1648: break;
1649: case 6:
1650: B->ops->mult = MatMult_SeqSBAIJ_6;
1651: B->ops->multadd = MatMultAdd_SeqSBAIJ_6;
1652: B->ops->multtranspose = MatMult_SeqSBAIJ_6;
1653: B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_6;
1654: break;
1655: case 7:
1656: B->ops->mult = MatMult_SeqSBAIJ_7;
1657: B->ops->multadd = MatMultAdd_SeqSBAIJ_7;
1658: B->ops->multtranspose = MatMult_SeqSBAIJ_7;
1659: B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_7;
1660: break;
1661: }
1662: }
1663:
1664: b->mbs = mbs;
1665: b->nbs = mbs;
1666: if (!skipallocation) {
1667: if (!b->imax) {
1668: PetscMalloc2(mbs,PetscInt,&b->imax,mbs,PetscInt,&b->ilen);
1669: b->free_imax_ilen = PETSC_TRUE;
1670: PetscLogObjectMemory(B,2*mbs*sizeof(PetscInt));
1671: }
1672: if (!nnz) {
1673: if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
1674: else if (nz <= 0) nz = 1;
1675: for (i=0; i<mbs; i++) {
1676: b->imax[i] = nz;
1677: }
1678: nz = nz*mbs; /* total nz */
1679: } else {
1680: nz = 0;
1681: for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
1682: }
1683: /* b->ilen will count nonzeros in each block row so far. */
1684: for (i=0; i<mbs; i++) { b->ilen[i] = 0;}
1685: /* nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */
1686:
1687: /* allocate the matrix space */
1688: MatSeqXAIJFreeAIJ(B,&b->a,&b->j,&b->i);
1689: PetscMalloc3(bs2*nz,PetscScalar,&b->a,nz,PetscInt,&b->j,B->rmap->N+1,PetscInt,&b->i);
1690: PetscLogObjectMemory(B,(B->rmap->N+1)*sizeof(PetscInt)+nz*(bs2*sizeof(PetscScalar)+sizeof(PetscInt)));
1691: PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));
1692: PetscMemzero(b->j,nz*sizeof(PetscInt));
1693: b->singlemalloc = PETSC_TRUE;
1694:
1695: /* pointer to beginning of each row */
1696: b->i[0] = 0;
1697: for (i=1; i<mbs+1; i++) {
1698: b->i[i] = b->i[i-1] + b->imax[i-1];
1699: }
1700: b->free_a = PETSC_TRUE;
1701: b->free_ij = PETSC_TRUE;
1702: } else {
1703: b->free_a = PETSC_FALSE;
1704: b->free_ij = PETSC_FALSE;
1705: }
1706:
1707: B->rmap->bs = bs;
1708: b->bs2 = bs2;
1709: b->nz = 0;
1710: b->maxnz = nz;
1711:
1712: b->inew = 0;
1713: b->jnew = 0;
1714: b->anew = 0;
1715: b->a2anew = 0;
1716: b->permute = PETSC_FALSE;
1717: if (realalloc) {MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);}
1718: return(0);
1719: }
1720: EXTERN_C_END
1722: /*
1723: This is used to set the numeric factorization for both Cholesky and ICC symbolic factorization
1724: */
1727: PetscErrorCode MatSeqSBAIJSetNumericFactorization_inplace(Mat B,PetscBool natural)
1728: {
1730: PetscBool flg = PETSC_FALSE;
1731: PetscInt bs = B->rmap->bs;
1734: PetscOptionsGetBool(((PetscObject)B)->prefix,"-mat_no_unroll",&flg,PETSC_NULL);
1735: if (flg) bs = 8;
1737: if (!natural) {
1738: switch (bs) {
1739: case 1:
1740: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_inplace;
1741: break;
1742: case 2:
1743: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2;
1744: break;
1745: case 3:
1746: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3;
1747: break;
1748: case 4:
1749: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4;
1750: break;
1751: case 5:
1752: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5;
1753: break;
1754: case 6:
1755: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6;
1756: break;
1757: case 7:
1758: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7;
1759: break;
1760: default:
1761: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N;
1762: break;
1763: }
1764: } else {
1765: switch (bs) {
1766: case 1:
1767: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering_inplace;
1768: break;
1769: case 2:
1770: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering;
1771: break;
1772: case 3:
1773: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering;
1774: break;
1775: case 4:
1776: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering;
1777: break;
1778: case 5:
1779: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering;
1780: break;
1781: case 6:
1782: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering;
1783: break;
1784: case 7:
1785: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering;
1786: break;
1787: default:
1788: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering;
1789: break;
1790: }
1791: }
1792: return(0);
1793: }
1795: EXTERN_C_BEGIN
1796: extern PetscErrorCode MatConvert_SeqSBAIJ_SeqAIJ(Mat, MatType,MatReuse,Mat*);
1797: extern PetscErrorCode MatConvert_SeqSBAIJ_SeqBAIJ(Mat, MatType,MatReuse,Mat*);
1798: EXTERN_C_END
1800:
1801: EXTERN_C_BEGIN
1804: PetscErrorCode MatGetFactor_seqsbaij_petsc(Mat A,MatFactorType ftype,Mat *B)
1805: {
1806: PetscInt n = A->rmap->n;
1807: PetscErrorCode ierr;
1810: #if defined(PETSC_USE_COMPLEX)
1811: if (A->hermitian)SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Hermitian Factor is not supported");
1812: #endif
1813: MatCreate(((PetscObject)A)->comm,B);
1814: MatSetSizes(*B,n,n,n,n);
1815: if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) {
1816: MatSetType(*B,MATSEQSBAIJ);
1817: MatSeqSBAIJSetPreallocation(*B,A->rmap->bs,MAT_SKIP_ALLOCATION,PETSC_NULL);
1818: (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqSBAIJ;
1819: (*B)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqSBAIJ;
1820: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
1821: (*B)->factortype = ftype;
1822: return(0);
1823: }
1824: EXTERN_C_END
1826: EXTERN_C_BEGIN
1829: PetscErrorCode MatGetFactorAvailable_seqsbaij_petsc(Mat A,MatFactorType ftype,PetscBool *flg)
1830: {
1832: if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) {
1833: *flg = PETSC_TRUE;
1834: } else {
1835: *flg = PETSC_FALSE;
1836: }
1837: return(0);
1838: }
1839: EXTERN_C_END
1841: EXTERN_C_BEGIN
1842: #if defined(PETSC_HAVE_MUMPS)
1843: extern PetscErrorCode MatGetFactor_sbaij_mumps(Mat,MatFactorType,Mat*);
1844: #endif
1845: #if defined(PETSC_HAVE_SPOOLES)
1846: extern PetscErrorCode MatGetFactor_seqsbaij_spooles(Mat,MatFactorType,Mat*);
1847: #endif
1848: #if defined(PETSC_HAVE_PASTIX)
1849: extern PetscErrorCode MatGetFactor_seqsbaij_pastix(Mat,MatFactorType,Mat*);
1850: #endif
1851: #if defined(PETSC_HAVE_CHOLMOD)
1852: extern PetscErrorCode MatGetFactor_seqsbaij_cholmod(Mat,MatFactorType,Mat*);
1853: #endif
1854: extern PetscErrorCode MatGetFactor_seqsbaij_sbstrm(Mat,MatFactorType,Mat*);
1855: EXTERN_C_END
1857: /*MC
1858: MATSEQSBAIJ - MATSEQSBAIJ = "seqsbaij" - A matrix type to be used for sequential symmetric block sparse matrices,
1859: based on block compressed sparse row format. Only the upper triangular portion of the matrix is stored.
1861: For complex numbers by default this matrix is symmetric, NOT Hermitian symmetric. To make it Hermitian symmetric you
1862: can call MatSetOption(Mat, MAT_HERMITIAN); after MatAssemblyEnd()
1864: Options Database Keys:
1865: . -mat_type seqsbaij - sets the matrix type to "seqsbaij" during a call to MatSetFromOptions()
1866:
1867: Notes: By default if you insert values into the lower triangular part of the matrix they are simply ignored (since they are not
1868: stored and it is assumed they symmetric to the upper triangular). If you call MatSetOption(Mat,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_FALSE) or use
1869: the options database -mat_ignore_lower_triangular false it will generate an error if you try to set a value in the lower triangular portion.
1872: Level: beginner
1873:
1874: .seealso: MatCreateSeqSBAIJ
1875: M*/
1877: EXTERN_C_BEGIN
1878: extern PetscErrorCode MatConvert_SeqSBAIJ_SeqSBSTRM(Mat, MatType,MatReuse,Mat*);
1879: EXTERN_C_END
1882: EXTERN_C_BEGIN
1885: PetscErrorCode MatCreate_SeqSBAIJ(Mat B)
1886: {
1887: Mat_SeqSBAIJ *b;
1889: PetscMPIInt size;
1890: PetscBool no_unroll = PETSC_FALSE,no_inode = PETSC_FALSE;
1893: MPI_Comm_size(((PetscObject)B)->comm,&size);
1894: if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
1895:
1896: PetscNewLog(B,Mat_SeqSBAIJ,&b);
1897: B->data = (void*)b;
1898: PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
1899: B->ops->destroy = MatDestroy_SeqSBAIJ;
1900: B->ops->view = MatView_SeqSBAIJ;
1901: b->row = 0;
1902: b->icol = 0;
1903: b->reallocs = 0;
1904: b->saved_values = 0;
1905: b->inode.limit = 5;
1906: b->inode.max_limit = 5;
1907:
1908: b->roworiented = PETSC_TRUE;
1909: b->nonew = 0;
1910: b->diag = 0;
1911: b->solve_work = 0;
1912: b->mult_work = 0;
1913: B->spptr = 0;
1914: B->info.nz_unneeded = (PetscReal)b->maxnz*b->bs2;
1915: b->keepnonzeropattern = PETSC_FALSE;
1916: b->xtoy = 0;
1917: b->XtoY = 0;
1918:
1919: b->inew = 0;
1920: b->jnew = 0;
1921: b->anew = 0;
1922: b->a2anew = 0;
1923: b->permute = PETSC_FALSE;
1925: b->ignore_ltriangular = PETSC_TRUE;
1926: PetscOptionsGetBool(((PetscObject)B)->prefix,"-mat_ignore_lower_triangular",&b->ignore_ltriangular,PETSC_NULL);
1928: b->getrow_utriangular = PETSC_FALSE;
1929: PetscOptionsGetBool(((PetscObject)B)->prefix,"-mat_getrow_uppertriangular",&b->getrow_utriangular,PETSC_NULL);
1931: #if defined(PETSC_HAVE_PASTIX)
1932: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_pastix_C",
1933: "MatGetFactor_seqsbaij_pastix",
1934: MatGetFactor_seqsbaij_pastix);
1935: #endif
1936: #if defined(PETSC_HAVE_SPOOLES)
1937: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_spooles_C",
1938: "MatGetFactor_seqsbaij_spooles",
1939: MatGetFactor_seqsbaij_spooles);
1940: #endif
1941: #if defined(PETSC_HAVE_MUMPS)
1942: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_mumps_C",
1943: "MatGetFactor_sbaij_mumps",
1944: MatGetFactor_sbaij_mumps);
1945: #endif
1946: #if defined(PETSC_HAVE_CHOLMOD)
1947: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_cholmod_C",
1948: "MatGetFactor_seqsbaij_cholmod",
1949: MatGetFactor_seqsbaij_cholmod);
1950: #endif
1951: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactorAvailable_petsc_C",
1952: "MatGetFactorAvailable_seqsbaij_petsc",
1953: MatGetFactorAvailable_seqsbaij_petsc);
1954: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_petsc_C",
1955: "MatGetFactor_seqsbaij_petsc",
1956: MatGetFactor_seqsbaij_petsc);
1957: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_sbstrm_C",
1958: "MatGetFactor_seqsbaij_sbstrm",
1959: MatGetFactor_seqsbaij_sbstrm);
1960: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
1961: "MatStoreValues_SeqSBAIJ",
1962: MatStoreValues_SeqSBAIJ);
1963: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
1964: "MatRetrieveValues_SeqSBAIJ",
1965: (void*)MatRetrieveValues_SeqSBAIJ);
1966: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C",
1967: "MatSeqSBAIJSetColumnIndices_SeqSBAIJ",
1968: MatSeqSBAIJSetColumnIndices_SeqSBAIJ);
1969: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_seqaij_C",
1970: "MatConvert_SeqSBAIJ_SeqAIJ",
1971: MatConvert_SeqSBAIJ_SeqAIJ);
1972: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_seqbaij_C",
1973: "MatConvert_SeqSBAIJ_SeqBAIJ",
1974: MatConvert_SeqSBAIJ_SeqBAIJ);
1975: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",
1976: "MatSeqSBAIJSetPreallocation_SeqSBAIJ",
1977: MatSeqSBAIJSetPreallocation_SeqSBAIJ);
1978: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_seqsbstrm_C",
1979: "MatConvert_SeqSBAIJ_SeqSBSTRM",
1980: MatConvert_SeqSBAIJ_SeqSBSTRM);
1982: B->symmetric = PETSC_TRUE;
1983: B->structurally_symmetric = PETSC_TRUE;
1984: B->symmetric_set = PETSC_TRUE;
1985: B->structurally_symmetric_set = PETSC_TRUE;
1986: PetscObjectChangeTypeName((PetscObject)B,MATSEQSBAIJ);
1988: PetscOptionsBegin(((PetscObject)B)->comm,((PetscObject)B)->prefix,"Options for SEQSBAIJ matrix","Mat");
1989: PetscOptionsBool("-mat_no_unroll","Do not optimize for inodes (slower)",PETSC_NULL,no_unroll,&no_unroll,PETSC_NULL);
1990: if (no_unroll) {PetscInfo(B,"Not using Inode routines due to -mat_no_unroll\n");}
1991: PetscOptionsBool("-mat_no_inode","Do not optimize for inodes (slower)",PETSC_NULL,no_inode,&no_inode,PETSC_NULL);
1992: if (no_inode) {PetscInfo(B,"Not using Inode routines due to -mat_no_inode\n");}
1993: PetscOptionsInt("-mat_inode_limit","Do not use inodes larger then this value",PETSC_NULL,b->inode.limit,&b->inode.limit,PETSC_NULL);
1994: PetscOptionsEnd();
1995: b->inode.use = (PetscBool)(!(no_unroll || no_inode));
1996: if (b->inode.limit > b->inode.max_limit) b->inode.limit = b->inode.max_limit;
1998: return(0);
1999: }
2000: EXTERN_C_END
2004: /*@C
2005: MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block
2006: compressed row) format. For good matrix assembly performance the
2007: user should preallocate the matrix storage by setting the parameter nz
2008: (or the array nnz). By setting these parameters accurately, performance
2009: during matrix assembly can be increased by more than a factor of 50.
2011: Collective on Mat
2013: Input Parameters:
2014: + A - the symmetric matrix
2015: . bs - size of block
2016: . nz - number of block nonzeros per block row (same for all rows)
2017: - nnz - array containing the number of block nonzeros in the upper triangular plus
2018: diagonal portion of each block (possibly different for each block row) or PETSC_NULL
2020: Options Database Keys:
2021: . -mat_no_unroll - uses code that does not unroll the loops in the
2022: block calculations (much slower)
2023: . -mat_block_size - size of the blocks to use (only works if a negative bs is passed in
2025: Level: intermediate
2027: Notes:
2028: Specify the preallocated storage with either nz or nnz (not both).
2029: Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
2030: allocation. See the <a href="../../docs/manual.pdf#nameddest=ch_mat">Mat chapter of the users manual</a> for details.
2032: You can call MatGetInfo() to get information on how effective the preallocation was;
2033: for example the fields mallocs,nz_allocated,nz_used,nz_unneeded;
2034: You can also run with the option -info and look for messages with the string
2035: malloc in them to see if additional memory allocation was needed.
2037: If the nnz parameter is given then the nz parameter is ignored
2040: .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateSBAIJ()
2041: @*/
2042: PetscErrorCode MatSeqSBAIJSetPreallocation(Mat B,PetscInt bs,PetscInt nz,const PetscInt nnz[])
2043: {
2050: PetscTryMethod(B,"MatSeqSBAIJSetPreallocation_C",(Mat,PetscInt,PetscInt,const PetscInt[]),(B,bs,nz,nnz));
2051: return(0);
2052: }
2056: /*@C
2057: MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block
2058: compressed row) format. For good matrix assembly performance the
2059: user should preallocate the matrix storage by setting the parameter nz
2060: (or the array nnz). By setting these parameters accurately, performance
2061: during matrix assembly can be increased by more than a factor of 50.
2063: Collective on MPI_Comm
2065: Input Parameters:
2066: + comm - MPI communicator, set to PETSC_COMM_SELF
2067: . bs - size of block
2068: . m - number of rows, or number of columns
2069: . nz - number of block nonzeros per block row (same for all rows)
2070: - nnz - array containing the number of block nonzeros in the upper triangular plus
2071: diagonal portion of each block (possibly different for each block row) or PETSC_NULL
2073: Output Parameter:
2074: . A - the symmetric matrix
2076: Options Database Keys:
2077: . -mat_no_unroll - uses code that does not unroll the loops in the
2078: block calculations (much slower)
2079: . -mat_block_size - size of the blocks to use
2081: Level: intermediate
2083: It is recommended that one use the MatCreate(), MatSetType() and/or MatSetFromOptions(),
2084: MatXXXXSetPreallocation() paradgm instead of this routine directly.
2085: [MatXXXXSetPreallocation() is, for example, MatSeqAIJSetPreallocation]
2087: Notes:
2088: The number of rows and columns must be divisible by blocksize.
2089: This matrix type does not support complex Hermitian operation.
2091: Specify the preallocated storage with either nz or nnz (not both).
2092: Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
2093: allocation. See the <a href="../../docs/manual.pdf#nameddest=ch_mat">Mat chapter of the users manual</a> for details.
2095: If the nnz parameter is given then the nz parameter is ignored
2097: .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateSBAIJ()
2098: @*/
2099: PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A)
2100: {
2102:
2104: MatCreate(comm,A);
2105: MatSetSizes(*A,m,n,m,n);
2106: MatSetType(*A,MATSEQSBAIJ);
2107: MatSeqSBAIJSetPreallocation_SeqSBAIJ(*A,bs,nz,(PetscInt*)nnz);
2108: return(0);
2109: }
2113: PetscErrorCode MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
2114: {
2115: Mat C;
2116: Mat_SeqSBAIJ *c,*a = (Mat_SeqSBAIJ*)A->data;
2118: PetscInt i,mbs = a->mbs,nz = a->nz,bs2 =a->bs2;
2121: if (a->i[mbs] != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Corrupt matrix");
2123: *B = 0;
2124: MatCreate(((PetscObject)A)->comm,&C);
2125: MatSetSizes(C,A->rmap->N,A->cmap->n,A->rmap->N,A->cmap->n);
2126: MatSetType(C,MATSEQSBAIJ);
2127: PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));
2128: c = (Mat_SeqSBAIJ*)C->data;
2130: C->preallocated = PETSC_TRUE;
2131: C->factortype = A->factortype;
2132: c->row = 0;
2133: c->icol = 0;
2134: c->saved_values = 0;
2135: c->keepnonzeropattern = a->keepnonzeropattern;
2136: C->assembled = PETSC_TRUE;
2138: PetscLayoutReference(A->rmap,&C->rmap);
2139: PetscLayoutReference(A->cmap,&C->cmap);
2140: c->bs2 = a->bs2;
2141: c->mbs = a->mbs;
2142: c->nbs = a->nbs;
2144: if (cpvalues == MAT_SHARE_NONZERO_PATTERN) {
2145: c->imax = a->imax;
2146: c->ilen = a->ilen;
2147: c->free_imax_ilen = PETSC_FALSE;
2148: } else {
2149: PetscMalloc2((mbs+1),PetscInt,&c->imax,(mbs+1),PetscInt,&c->ilen);
2150: PetscLogObjectMemory(C,2*(mbs+1)*sizeof(PetscInt));
2151: for (i=0; i<mbs; i++) {
2152: c->imax[i] = a->imax[i];
2153: c->ilen[i] = a->ilen[i];
2154: }
2155: c->free_imax_ilen = PETSC_TRUE;
2156: }
2158: /* allocate the matrix space */
2159: if (cpvalues == MAT_SHARE_NONZERO_PATTERN) {
2160: PetscMalloc(bs2*nz*sizeof(MatScalar),&c->a);
2161: PetscLogObjectMemory(C,nz*bs2*sizeof(MatScalar));
2162: c->i = a->i;
2163: c->j = a->j;
2164: c->singlemalloc = PETSC_FALSE;
2165: c->free_a = PETSC_TRUE;
2166: c->free_ij = PETSC_FALSE;
2167: c->parent = A;
2168: PetscObjectReference((PetscObject)A);
2169: MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
2170: MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
2171: } else {
2172: PetscMalloc3(bs2*nz,MatScalar,&c->a,nz,PetscInt,&c->j,mbs+1,PetscInt,&c->i);
2173: PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(PetscInt));
2174: PetscLogObjectMemory(C,(mbs+1)*sizeof(PetscInt) + nz*(bs2*sizeof(MatScalar) + sizeof(PetscInt)));
2175: c->singlemalloc = PETSC_TRUE;
2176: c->free_a = PETSC_TRUE;
2177: c->free_ij = PETSC_TRUE;
2178: }
2179: if (mbs > 0) {
2180: if (cpvalues != MAT_SHARE_NONZERO_PATTERN) {
2181: PetscMemcpy(c->j,a->j,nz*sizeof(PetscInt));
2182: }
2183: if (cpvalues == MAT_COPY_VALUES) {
2184: PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));
2185: } else {
2186: PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));
2187: }
2188: if (a->jshort) {
2189: /* cannot share jshort, it is reallocated in MatAssemblyEnd_SeqSBAIJ() */
2190: /* if the parent matrix is reassembled, this child matrix will never notice */
2191: PetscMalloc(nz*sizeof(unsigned short),&c->jshort);
2192: PetscLogObjectMemory(C,nz*sizeof(unsigned short));
2193: PetscMemcpy(c->jshort,a->jshort,nz*sizeof(unsigned short));
2194: c->free_jshort = PETSC_TRUE;
2195: }
2196: }
2198: c->roworiented = a->roworiented;
2199: c->nonew = a->nonew;
2201: if (a->diag) {
2202: if (cpvalues == MAT_SHARE_NONZERO_PATTERN) {
2203: c->diag = a->diag;
2204: c->free_diag = PETSC_FALSE;
2205: } else {
2206: PetscMalloc(mbs*sizeof(PetscInt),&c->diag);
2207: PetscLogObjectMemory(C,mbs*sizeof(PetscInt));
2208: for (i=0; i<mbs; i++) {
2209: c->diag[i] = a->diag[i];
2210: }
2211: c->free_diag = PETSC_TRUE;
2212: }
2213: }
2214: c->nz = a->nz;
2215: c->maxnz = a->nz; /* Since we allocate exactly the right amount */
2216: c->solve_work = 0;
2217: c->mult_work = 0;
2218: *B = C;
2219: PetscFListDuplicate(((PetscObject)A)->qlist,&((PetscObject)C)->qlist);
2220: return(0);
2221: }
2225: PetscErrorCode MatLoad_SeqSBAIJ(Mat newmat,PetscViewer viewer)
2226: {
2227: Mat_SeqSBAIJ *a;
2229: int fd;
2230: PetscMPIInt size;
2231: PetscInt i,nz,header[4],*rowlengths=0,M,N,bs=1;
2232: PetscInt *mask,mbs,*jj,j,rowcount,nzcount,k,*s_browlengths,maskcount;
2233: PetscInt kmax,jcount,block,idx,point,nzcountb,extra_rows,rows,cols;
2234: PetscInt *masked,nmask,tmp,bs2,ishift;
2235: PetscScalar *aa;
2236: MPI_Comm comm = ((PetscObject)viewer)->comm;
2239: PetscOptionsGetInt(((PetscObject)newmat)->prefix,"-matload_block_size",&bs,PETSC_NULL);
2240: bs2 = bs*bs;
2242: MPI_Comm_size(comm,&size);
2243: if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"view must have one processor");
2244: PetscViewerBinaryGetDescriptor(viewer,&fd);
2245: PetscBinaryRead(fd,header,4,PETSC_INT);
2246: if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"not Mat object");
2247: M = header[1]; N = header[2]; nz = header[3];
2249: if (header[3] < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ");
2251: if (M != N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Can only do square matrices");
2253: /*
2254: This code adds extra rows to make sure the number of rows is
2255: divisible by the blocksize
2256: */
2257: mbs = M/bs;
2258: extra_rows = bs - M + bs*(mbs);
2259: if (extra_rows == bs) extra_rows = 0;
2260: else mbs++;
2261: if (extra_rows) {
2262: PetscInfo(viewer,"Padding loaded matrix to match blocksize\n");
2263: }
2265: /* Set global sizes if not already set */
2266: if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) {
2267: MatSetSizes(newmat,PETSC_DECIDE,PETSC_DECIDE,M+extra_rows,N+extra_rows);
2268: } else { /* Check if the matrix global sizes are correct */
2269: MatGetSize(newmat,&rows,&cols);
2270: if (M != rows || N != cols) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Matrix in file of different length (%d, %d) than the input matrix (%d, %d)",M,N,rows,cols);
2271: }
2273: /* read in row lengths */
2274: PetscMalloc((M+extra_rows)*sizeof(PetscInt),&rowlengths);
2275: PetscBinaryRead(fd,rowlengths,M,PETSC_INT);
2276: for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;
2278: /* read in column indices */
2279: PetscMalloc((nz+extra_rows)*sizeof(PetscInt),&jj);
2280: PetscBinaryRead(fd,jj,nz,PETSC_INT);
2281: for (i=0; i<extra_rows; i++) jj[nz+i] = M+i;
2283: /* loop over row lengths determining block row lengths */
2284: PetscMalloc(mbs*sizeof(PetscInt),&s_browlengths);
2285: PetscMemzero(s_browlengths,mbs*sizeof(PetscInt));
2286: PetscMalloc2(mbs,PetscInt,&mask,mbs,PetscInt,&masked);
2287: PetscMemzero(mask,mbs*sizeof(PetscInt));
2288: rowcount = 0;
2289: nzcount = 0;
2290: for (i=0; i<mbs; i++) {
2291: nmask = 0;
2292: for (j=0; j<bs; j++) {
2293: kmax = rowlengths[rowcount];
2294: for (k=0; k<kmax; k++) {
2295: tmp = jj[nzcount++]/bs; /* block col. index */
2296: if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;}
2297: }
2298: rowcount++;
2299: }
2300: s_browlengths[i] += nmask;
2301:
2302: /* zero out the mask elements we set */
2303: for (j=0; j<nmask; j++) mask[masked[j]] = 0;
2304: }
2306: /* Do preallocation */
2307: MatSeqSBAIJSetPreallocation_SeqSBAIJ(newmat,bs,0,s_browlengths);
2308: a = (Mat_SeqSBAIJ*)newmat->data;
2310: /* set matrix "i" values */
2311: a->i[0] = 0;
2312: for (i=1; i<= mbs; i++) {
2313: a->i[i] = a->i[i-1] + s_browlengths[i-1];
2314: a->ilen[i-1] = s_browlengths[i-1];
2315: }
2316: a->nz = a->i[mbs];
2318: /* read in nonzero values */
2319: PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);
2320: PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);
2321: for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0;
2323: /* set "a" and "j" values into matrix */
2324: nzcount = 0; jcount = 0;
2325: for (i=0; i<mbs; i++) {
2326: nzcountb = nzcount;
2327: nmask = 0;
2328: for (j=0; j<bs; j++) {
2329: kmax = rowlengths[i*bs+j];
2330: for (k=0; k<kmax; k++) {
2331: tmp = jj[nzcount++]/bs; /* block col. index */
2332: if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;}
2333: }
2334: }
2335: /* sort the masked values */
2336: PetscSortInt(nmask,masked);
2338: /* set "j" values into matrix */
2339: maskcount = 1;
2340: for (j=0; j<nmask; j++) {
2341: a->j[jcount++] = masked[j];
2342: mask[masked[j]] = maskcount++;
2343: }
2345: /* set "a" values into matrix */
2346: ishift = bs2*a->i[i];
2347: for (j=0; j<bs; j++) {
2348: kmax = rowlengths[i*bs+j];
2349: for (k=0; k<kmax; k++) {
2350: tmp = jj[nzcountb]/bs ; /* block col. index */
2351: if (tmp >= i){
2352: block = mask[tmp] - 1;
2353: point = jj[nzcountb] - bs*tmp;
2354: idx = ishift + bs2*block + j + bs*point;
2355: a->a[idx] = aa[nzcountb];
2356: }
2357: nzcountb++;
2358: }
2359: }
2360: /* zero out the mask elements we set */
2361: for (j=0; j<nmask; j++) mask[masked[j]] = 0;
2362: }
2363: if (jcount != a->nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix");
2365: PetscFree(rowlengths);
2366: PetscFree(s_browlengths);
2367: PetscFree(aa);
2368: PetscFree(jj);
2369: PetscFree2(mask,masked);
2371: MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);
2372: MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);
2373: MatView_Private(newmat);
2374: return(0);
2375: }
2379: /*@
2380: MatCreateSeqSBAIJWithArrays - Creates an sequential SBAIJ matrix using matrix elements
2381: (upper triangular entries in CSR format) provided by the user.
2383: Collective on MPI_Comm
2385: Input Parameters:
2386: + comm - must be an MPI communicator of size 1
2387: . bs - size of block
2388: . m - number of rows
2389: . n - number of columns
2390: . i - row indices
2391: . j - column indices
2392: - a - matrix values
2394: Output Parameter:
2395: . mat - the matrix
2397: Level: advanced
2399: Notes:
2400: The i, j, and a arrays are not copied by this routine, the user must free these arrays
2401: once the matrix is destroyed
2403: You cannot set new nonzero locations into this matrix, that will generate an error.
2405: The i and j indices are 0 based
2407: When block size is greater than 1 the matrix values must be stored using the SBAIJ storage format (see the SBAIJ code to determine this). For block size of 1
2408: it is the regular CSR format excluding the lower triangular elements.
2410: .seealso: MatCreate(), MatCreateSBAIJ(), MatCreateSeqSBAIJ()
2412: @*/
2413: PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt* i,PetscInt*j,PetscScalar *a,Mat *mat)
2414: {
2416: PetscInt ii;
2417: Mat_SeqSBAIJ *sbaij;
2420: if (bs != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"block size %D > 1 is not supported yet",bs);
2421: if (i[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"i (row indices) must start with 0");
2422:
2423: MatCreate(comm,mat);
2424: MatSetSizes(*mat,m,n,m,n);
2425: MatSetType(*mat,MATSEQSBAIJ);
2426: MatSeqSBAIJSetPreallocation_SeqSBAIJ(*mat,bs,MAT_SKIP_ALLOCATION,0);
2427: sbaij = (Mat_SeqSBAIJ*)(*mat)->data;
2428: PetscMalloc2(m,PetscInt,&sbaij->imax,m,PetscInt,&sbaij->ilen);
2429: PetscLogObjectMemory(*mat,2*m*sizeof(PetscInt));
2431: sbaij->i = i;
2432: sbaij->j = j;
2433: sbaij->a = a;
2434: sbaij->singlemalloc = PETSC_FALSE;
2435: sbaij->nonew = -1; /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/
2436: sbaij->free_a = PETSC_FALSE;
2437: sbaij->free_ij = PETSC_FALSE;
2439: for (ii=0; ii<m; ii++) {
2440: sbaij->ilen[ii] = sbaij->imax[ii] = i[ii+1] - i[ii];
2441: #if defined(PETSC_USE_DEBUG)
2442: if (i[ii+1] - i[ii] < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row length in i (row indices) row = %d length = %d",ii,i[ii+1] - i[ii]);
2443: #endif
2444: }
2445: #if defined(PETSC_USE_DEBUG)
2446: for (ii=0; ii<sbaij->i[m]; ii++) {
2447: if (j[ii] < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column index at location = %d index = %d",ii,j[ii]);
2448: if (j[ii] > n - 1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column index to large at location = %d index = %d",ii,j[ii]);
2449: }
2450: #endif
2452: MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
2453: MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
2454: return(0);
2455: }