Actual source code: sbaij.c
petsc-dev 2014-02-02
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: PetscMalloc1(a->mbs,&a->diag);
58: PetscLogObjectMemory((PetscObject)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,const PetscInt *inia[],const PetscInt *inja[],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;
71: PetscInt **ia = (PetscInt**)inia,**ja = (PetscInt**)inja;
75: *nn = n;
76: if (!ia) return(0);
77: if (!blockcompressed) {
78: /* malloc & create the natural set of indices */
79: PetscMalloc2((n+1)*bs,ia,nz*bs,ja);
80: for (i=0; i<n+1; i++) {
81: for (j=0; j<bs; j++) {
82: *ia[i*bs+j] = a->i[i]*bs+j+oshift;
83: }
84: }
85: for (i=0; i<nz; i++) {
86: for (j=0; j<bs; j++) {
87: *ja[i*bs+j] = a->j[i]*bs+j+oshift;
88: }
89: }
90: } else { /* blockcompressed */
91: if (oshift == 1) {
92: /* temporarily add 1 to i and j indices */
93: for (i=0; i<nz; i++) a->j[i]++;
94: for (i=0; i<n+1; i++) a->i[i]++;
95: }
96: *ia = a->i; *ja = a->j;
97: }
98: return(0);
99: }
103: static PetscErrorCode MatRestoreRowIJ_SeqSBAIJ(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *nn,const PetscInt *ia[],const 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: }
118: return(0);
119: }
123: PetscErrorCode MatDestroy_SeqSBAIJ(Mat A)
124: {
125: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
129: #if defined(PETSC_USE_LOG)
130: PetscLogObjectState((PetscObject)A,"Rows=%D, NZ=%D",A->rmap->N,a->nz);
131: #endif
132: MatSeqXAIJFreeAIJ(A,&a->a,&a->j,&a->i);
133: if (a->free_diag) {PetscFree(a->diag);}
134: ISDestroy(&a->row);
135: ISDestroy(&a->col);
136: ISDestroy(&a->icol);
137: PetscFree(a->idiag);
138: PetscFree(a->inode.size);
139: if (a->free_imax_ilen) {PetscFree2(a->imax,a->ilen);}
140: PetscFree(a->solve_work);
141: PetscFree(a->sor_work);
142: PetscFree(a->solves_work);
143: PetscFree(a->mult_work);
144: PetscFree(a->saved_values);
145: PetscFree(a->xtoy);
146: if (a->free_jshort) {PetscFree(a->jshort);}
147: PetscFree(a->inew);
148: MatDestroy(&a->parent);
149: PetscFree(A->data);
151: PetscObjectChangeTypeName((PetscObject)A,0);
152: PetscObjectComposeFunction((PetscObject)A,"MatStoreValues_C",NULL);
153: PetscObjectComposeFunction((PetscObject)A,"MatRetrieveValues_C",NULL);
154: PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetColumnIndices_C",NULL);
155: PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqaij_C",NULL);
156: PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqbaij_C",NULL);
157: PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetPreallocation_C",NULL);
158: PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetPreallocationCSR_C",NULL);
159: PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqsbstrm_C",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: /* These options are handled directly by MatSetOption() */
205: break;
206: case MAT_SYMMETRIC:
207: case MAT_STRUCTURALLY_SYMMETRIC:
208: case MAT_SYMMETRY_ETERNAL:
209: if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix must be symmetric");
210: PetscInfo1(A,"Option %s not relevent\n",MatOptions[op]);
211: break;
212: case MAT_IGNORE_LOWER_TRIANGULAR:
213: a->ignore_ltriangular = flg;
214: break;
215: case MAT_ERROR_LOWER_TRIANGULAR:
216: a->ignore_ltriangular = flg;
217: break;
218: case MAT_GETROW_UPPERTRIANGULAR:
219: a->getrow_utriangular = flg;
220: break;
221: default:
222: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %d",op);
223: }
224: return(0);
225: }
229: PetscErrorCode MatGetRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **v)
230: {
231: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
233: PetscInt itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*cols_i,bs2;
234: MatScalar *aa,*aa_i;
235: PetscScalar *v_i;
238: 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()");
239: /* Get the upper triangular part of the row */
240: bs = A->rmap->bs;
241: ai = a->i;
242: aj = a->j;
243: aa = a->a;
244: bs2 = a->bs2;
246: if (row < 0 || row >= A->rmap->N) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Row %D out of range", row);
248: bn = row/bs; /* Block number */
249: bp = row % bs; /* Block position */
250: M = ai[bn+1] - ai[bn];
251: *ncols = bs*M;
253: if (v) {
254: *v = 0;
255: if (*ncols) {
256: PetscMalloc1((*ncols+row),v);
257: for (i=0; i<M; i++) { /* for each block in the block row */
258: v_i = *v + i*bs;
259: aa_i = aa + bs2*(ai[bn] + i);
260: for (j=bp,k=0; j<bs2; j+=bs,k++) v_i[k] = aa_i[j];
261: }
262: }
263: }
265: if (cols) {
266: *cols = 0;
267: if (*ncols) {
268: PetscMalloc1((*ncols+row),cols);
269: for (i=0; i<M; i++) { /* for each block in the block row */
270: cols_i = *cols + i*bs;
271: itmp = bs*aj[ai[bn] + i];
272: for (j=0; j<bs; j++) cols_i[j] = itmp++;
273: }
274: }
275: }
277: /*search column A(0:row-1,row) (=A(row,0:row-1)). Could be expensive! */
278: /* this segment is currently removed, so only entries in the upper triangle are obtained */
279: #if defined(column_search)
280: v_i = *v + M*bs;
281: cols_i = *cols + M*bs;
282: for (i=0; i<bn; i++) { /* for each block row */
283: M = ai[i+1] - ai[i];
284: for (j=0; j<M; j++) {
285: itmp = aj[ai[i] + j]; /* block column value */
286: if (itmp == bn) {
287: aa_i = aa + bs2*(ai[i] + j) + bs*bp;
288: for (k=0; k<bs; k++) {
289: *cols_i++ = i*bs+k;
290: *v_i++ = aa_i[k];
291: }
292: *ncols += bs;
293: break;
294: }
295: }
296: }
297: #endif
298: return(0);
299: }
303: PetscErrorCode MatRestoreRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
304: {
308: if (idx) {PetscFree(*idx);}
309: if (v) {PetscFree(*v);}
310: return(0);
311: }
315: PetscErrorCode MatGetRowUpperTriangular_SeqSBAIJ(Mat A)
316: {
317: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
320: a->getrow_utriangular = PETSC_TRUE;
321: return(0);
322: }
325: PetscErrorCode MatRestoreRowUpperTriangular_SeqSBAIJ(Mat A)
326: {
327: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
330: a->getrow_utriangular = PETSC_FALSE;
331: return(0);
332: }
336: PetscErrorCode MatTranspose_SeqSBAIJ(Mat A,MatReuse reuse,Mat *B)
337: {
341: if (reuse == MAT_INITIAL_MATRIX || *B != A) {
342: MatDuplicate(A,MAT_COPY_VALUES,B);
343: }
344: return(0);
345: }
349: static PetscErrorCode MatView_SeqSBAIJ_ASCII(Mat A,PetscViewer viewer)
350: {
351: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
352: PetscErrorCode ierr;
353: PetscInt i,j,bs = A->rmap->bs,k,l,bs2=a->bs2;
354: PetscViewerFormat format;
355: PetscInt *diag;
358: PetscViewerGetFormat(viewer,&format);
359: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
360: PetscViewerASCIIPrintf(viewer," block size is %D\n",bs);
361: } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
362: Mat aij;
363: if (A->factortype && bs>1) {
364: PetscPrintf(PETSC_COMM_SELF,"Warning: matrix is factored with bs>1. MatView() with PETSC_VIEWER_ASCII_MATLAB is not supported and ignored!\n");
365: return(0);
366: }
367: MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&aij);
368: MatView(aij,viewer);
369: MatDestroy(&aij);
370: } else if (format == PETSC_VIEWER_ASCII_COMMON) {
371: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
372: for (i=0; i<a->mbs; i++) {
373: for (j=0; j<bs; j++) {
374: PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);
375: for (k=a->i[i]; k<a->i[i+1]; k++) {
376: for (l=0; l<bs; l++) {
377: #if defined(PETSC_USE_COMPLEX)
378: if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
379: PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",bs*a->j[k]+l,
380: (double)PetscRealPart(a->a[bs2*k + l*bs + j]),(double)PetscImaginaryPart(a->a[bs2*k + l*bs + j]));
381: } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
382: PetscViewerASCIIPrintf(viewer," (%D, %g - %g i) ",bs*a->j[k]+l,
383: (double)PetscRealPart(a->a[bs2*k + l*bs + j]),-(double)PetscImaginaryPart(a->a[bs2*k + l*bs + j]));
384: } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
385: PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,(double)PetscRealPart(a->a[bs2*k + l*bs + j]));
386: }
387: #else
388: if (a->a[bs2*k + l*bs + j] != 0.0) {
389: PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,(double)a->a[bs2*k + l*bs + j]);
390: }
391: #endif
392: }
393: }
394: PetscViewerASCIIPrintf(viewer,"\n");
395: }
396: }
397: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
398: } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
399: return(0);
400: } else {
401: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
402: PetscObjectPrintClassNamePrefixType((PetscObject)A,viewer);
403: if (A->factortype) { /* for factored matrix */
404: if (bs>1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"matrix is factored with bs>1. Not implemented yet");
406: diag=a->diag;
407: for (i=0; i<a->mbs; i++) { /* for row block i */
408: PetscViewerASCIIPrintf(viewer,"row %D:",i);
409: /* diagonal entry */
410: #if defined(PETSC_USE_COMPLEX)
411: if (PetscImaginaryPart(a->a[diag[i]]) > 0.0) {
412: PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",a->j[diag[i]],(double)PetscRealPart(1.0/a->a[diag[i]]),(double)PetscImaginaryPart(1.0/a->a[diag[i]]));
413: } else if (PetscImaginaryPart(a->a[diag[i]]) < 0.0) {
414: PetscViewerASCIIPrintf(viewer," (%D, %g - %g i) ",a->j[diag[i]],(double)PetscRealPart(1.0/a->a[diag[i]]),-(double)PetscImaginaryPart(1.0/a->a[diag[i]]));
415: } else {
416: PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[diag[i]],(double)PetscRealPart(1.0/a->a[diag[i]]));
417: }
418: #else
419: PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[diag[i]],(double)(1.0/a->a[diag[i]]));
420: #endif
421: /* off-diagonal entries */
422: for (k=a->i[i]; k<a->i[i+1]-1; k++) {
423: #if defined(PETSC_USE_COMPLEX)
424: if (PetscImaginaryPart(a->a[k]) > 0.0) {
425: PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",bs*a->j[k],(double)PetscRealPart(a->a[k]),(double)PetscImaginaryPart(a->a[k]));
426: } else if (PetscImaginaryPart(a->a[k]) < 0.0) {
427: PetscViewerASCIIPrintf(viewer," (%D, %g - %g i) ",bs*a->j[k],(double)PetscRealPart(a->a[k]),-(double)PetscImaginaryPart(a->a[k]));
428: } else {
429: PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k],(double)PetscRealPart(a->a[k]));
430: }
431: #else
432: PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[k],(double)a->a[k]);
433: #endif
434: }
435: PetscViewerASCIIPrintf(viewer,"\n");
436: }
438: } else { /* for non-factored matrix */
439: for (i=0; i<a->mbs; i++) { /* for row block i */
440: for (j=0; j<bs; j++) { /* for row bs*i + j */
441: PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);
442: for (k=a->i[i]; k<a->i[i+1]; k++) { /* for column block */
443: for (l=0; l<bs; l++) { /* for column */
444: #if defined(PETSC_USE_COMPLEX)
445: if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) {
446: PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",bs*a->j[k]+l,
447: (double)PetscRealPart(a->a[bs2*k + l*bs + j]),(double)PetscImaginaryPart(a->a[bs2*k + l*bs + j]));
448: } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) {
449: PetscViewerASCIIPrintf(viewer," (%D, %g - %g i) ",bs*a->j[k]+l,
450: (double)PetscRealPart(a->a[bs2*k + l*bs + j]),-(double)PetscImaginaryPart(a->a[bs2*k + l*bs + j]));
451: } else {
452: PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,(double)PetscRealPart(a->a[bs2*k + l*bs + j]));
453: }
454: #else
455: PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,(double)a->a[bs2*k + l*bs + j]);
456: #endif
457: }
458: }
459: PetscViewerASCIIPrintf(viewer,"\n");
460: }
461: }
462: }
463: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
464: }
465: PetscViewerFlush(viewer);
466: return(0);
467: }
469: #include <petscdraw.h>
472: static PetscErrorCode MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw,void *Aa)
473: {
474: Mat A = (Mat) Aa;
475: Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data;
477: PetscInt row,i,j,k,l,mbs=a->mbs,color,bs=A->rmap->bs,bs2=a->bs2;
478: PetscMPIInt rank;
479: PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r;
480: MatScalar *aa;
481: MPI_Comm comm;
482: PetscViewer viewer;
485: /*
486: This is nasty. If this is called from an originally parallel matrix
487: then all processes call this,but only the first has the matrix so the
488: rest should return immediately.
489: */
490: PetscObjectGetComm((PetscObject)draw,&comm);
491: MPI_Comm_rank(comm,&rank);
492: if (rank) return(0);
494: PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);
496: PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);
497: PetscDrawString(draw, .3*(xl+xr), .3*(yl+yr), PETSC_DRAW_BLACK, "symmetric");
499: /* loop over matrix elements drawing boxes */
500: color = PETSC_DRAW_BLUE;
501: for (i=0,row=0; i<mbs; i++,row+=bs) {
502: for (j=a->i[i]; j<a->i[i+1]; j++) {
503: y_l = A->rmap->N - row - 1.0; y_r = y_l + 1.0;
504: x_l = a->j[j]*bs; x_r = x_l + 1.0;
505: aa = a->a + j*bs2;
506: for (k=0; k<bs; k++) {
507: for (l=0; l<bs; l++) {
508: if (PetscRealPart(*aa++) >= 0.) continue;
509: PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
510: }
511: }
512: }
513: }
514: color = PETSC_DRAW_CYAN;
515: for (i=0,row=0; i<mbs; i++,row+=bs) {
516: for (j=a->i[i]; j<a->i[i+1]; j++) {
517: y_l = A->rmap->N - row - 1.0; y_r = y_l + 1.0;
518: x_l = a->j[j]*bs; x_r = x_l + 1.0;
519: aa = a->a + j*bs2;
520: for (k=0; k<bs; k++) {
521: for (l=0; l<bs; l++) {
522: if (PetscRealPart(*aa++) != 0.) continue;
523: PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
524: }
525: }
526: }
527: }
529: color = PETSC_DRAW_RED;
530: for (i=0,row=0; i<mbs; i++,row+=bs) {
531: for (j=a->i[i]; j<a->i[i+1]; j++) {
532: y_l = A->rmap->N - row - 1.0; y_r = y_l + 1.0;
533: x_l = a->j[j]*bs; x_r = x_l + 1.0;
534: aa = a->a + j*bs2;
535: for (k=0; k<bs; k++) {
536: for (l=0; l<bs; l++) {
537: if (PetscRealPart(*aa++) <= 0.) continue;
538: PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
539: }
540: }
541: }
542: }
543: return(0);
544: }
548: static PetscErrorCode MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer)
549: {
551: PetscReal xl,yl,xr,yr,w,h;
552: PetscDraw draw;
553: PetscBool isnull;
556: PetscViewerDrawGetDraw(viewer,0,&draw);
557: PetscDrawIsNull(draw,&isnull); if (isnull) return(0);
559: PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);
560: xr = A->rmap->N; yr = A->rmap->N; h = yr/10.0; w = xr/10.0;
561: xr += w; yr += h; xl = -w; yl = -h;
562: PetscDrawSetCoordinates(draw,xl,yl,xr,yr);
563: PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);
564: PetscObjectCompose((PetscObject)A,"Zoomviewer",NULL);
565: return(0);
566: }
570: PetscErrorCode MatView_SeqSBAIJ(Mat A,PetscViewer viewer)
571: {
573: PetscBool iascii,isdraw;
574: FILE *file = 0;
577: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
578: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
579: if (iascii) {
580: MatView_SeqSBAIJ_ASCII(A,viewer);
581: } else if (isdraw) {
582: MatView_SeqSBAIJ_Draw(A,viewer);
583: } else {
584: Mat B;
585: MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&B);
586: MatView(B,viewer);
587: MatDestroy(&B);
588: PetscViewerBinaryGetInfoPointer(viewer,&file);
589: if (file) {
590: fprintf(file,"-matload_block_size %d\n",(int)A->rmap->bs);
591: }
592: }
593: return(0);
594: }
599: PetscErrorCode MatGetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],PetscScalar v[])
600: {
601: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
602: PetscInt *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
603: PetscInt *ai = a->i,*ailen = a->ilen;
604: PetscInt brow,bcol,ridx,cidx,bs=A->rmap->bs,bs2=a->bs2;
605: MatScalar *ap,*aa = a->a;
608: for (k=0; k<m; k++) { /* loop over rows */
609: row = im[k]; brow = row/bs;
610: if (row < 0) {v += n; continue;} /* SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",row); */
611: 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);
612: rp = aj + ai[brow]; ap = aa + bs2*ai[brow];
613: nrow = ailen[brow];
614: for (l=0; l<n; l++) { /* loop over columns */
615: if (in[l] < 0) {v++; continue;} /* SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %D",in[l]); */
616: 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);
617: col = in[l];
618: bcol = col/bs;
619: cidx = col%bs;
620: ridx = row%bs;
621: high = nrow;
622: low = 0; /* assume unsorted */
623: while (high-low > 5) {
624: t = (low+high)/2;
625: if (rp[t] > bcol) high = t;
626: else low = t;
627: }
628: for (i=low; i<high; i++) {
629: if (rp[i] > bcol) break;
630: if (rp[i] == bcol) {
631: *v++ = ap[bs2*i+bs*cidx+ridx];
632: goto finished;
633: }
634: }
635: *v++ = 0.0;
636: finished:;
637: }
638: }
639: return(0);
640: }
645: PetscErrorCode MatSetValuesBlocked_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
646: {
647: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
648: PetscErrorCode ierr;
649: PetscInt *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,lastcol = -1;
650: PetscInt *imax =a->imax,*ai=a->i,*ailen=a->ilen;
651: PetscInt *aj =a->j,nonew=a->nonew,bs2=a->bs2,bs=A->rmap->bs,stepval;
652: PetscBool roworiented=a->roworiented;
653: const PetscScalar *value = v;
654: MatScalar *ap,*aa = a->a,*bap;
657: if (roworiented) stepval = (n-1)*bs;
658: else stepval = (m-1)*bs;
660: for (k=0; k<m; k++) { /* loop over added rows */
661: row = im[k];
662: if (row < 0) continue;
663: #if defined(PETSC_USE_DEBUG)
664: if (row >= a->mbs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,a->mbs-1);
665: #endif
666: rp = aj + ai[row];
667: ap = aa + bs2*ai[row];
668: rmax = imax[row];
669: nrow = ailen[row];
670: low = 0;
671: high = nrow;
672: for (l=0; l<n; l++) { /* loop over added columns */
673: if (in[l] < 0) continue;
674: col = in[l];
675: #if defined(PETSC_USE_DEBUG)
676: if (col >= a->nbs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",col,a->nbs-1);
677: #endif
678: if (col < row) {
679: if (a->ignore_ltriangular) continue; /* ignore lower triangular block */
680: else 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)");
681: }
682: if (roworiented) value = v + k*(stepval+bs)*bs + l*bs;
683: else value = v + l*(stepval+bs)*bs + k*bs;
685: if (col <= lastcol) low = 0;
686: else high = nrow;
688: lastcol = col;
689: while (high-low > 7) {
690: t = (low+high)/2;
691: if (rp[t] > col) high = t;
692: else low = t;
693: }
694: for (i=low; i<high; i++) {
695: if (rp[i] > col) break;
696: if (rp[i] == col) {
697: bap = ap + bs2*i;
698: if (roworiented) {
699: if (is == ADD_VALUES) {
700: for (ii=0; ii<bs; ii++,value+=stepval) {
701: for (jj=ii; jj<bs2; jj+=bs) {
702: bap[jj] += *value++;
703: }
704: }
705: } else {
706: for (ii=0; ii<bs; ii++,value+=stepval) {
707: for (jj=ii; jj<bs2; jj+=bs) {
708: bap[jj] = *value++;
709: }
710: }
711: }
712: } else {
713: if (is == ADD_VALUES) {
714: for (ii=0; ii<bs; ii++,value+=stepval) {
715: for (jj=0; jj<bs; jj++) {
716: *bap++ += *value++;
717: }
718: }
719: } else {
720: for (ii=0; ii<bs; ii++,value+=stepval) {
721: for (jj=0; jj<bs; jj++) {
722: *bap++ = *value++;
723: }
724: }
725: }
726: }
727: goto noinsert2;
728: }
729: }
730: if (nonew == 1) goto noinsert2;
731: if (nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col);
732: MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar);
733: N = nrow++ - 1; high++;
734: /* shift up all the later entries in this row */
735: for (ii=N; ii>=i; ii--) {
736: rp[ii+1] = rp[ii];
737: PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));
738: }
739: if (N >= i) {
740: PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));
741: }
742: rp[i] = col;
743: bap = ap + bs2*i;
744: if (roworiented) {
745: for (ii=0; ii<bs; ii++,value+=stepval) {
746: for (jj=ii; jj<bs2; jj+=bs) {
747: bap[jj] = *value++;
748: }
749: }
750: } else {
751: for (ii=0; ii<bs; ii++,value+=stepval) {
752: for (jj=0; jj<bs; jj++) {
753: *bap++ = *value++;
754: }
755: }
756: }
757: noinsert2:;
758: low = i;
759: }
760: ailen[row] = nrow;
761: }
762: return(0);
763: }
765: /*
766: This is not yet used
767: */
770: PetscErrorCode MatAssemblyEnd_SeqSBAIJ_SeqAIJ_Inode(Mat A)
771: {
772: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
774: const PetscInt *ai = a->i, *aj = a->j,*cols;
775: PetscInt i = 0,j,blk_size,m = A->rmap->n,node_count = 0,nzx,nzy,*ns,row,nz,cnt,cnt2,*counts;
776: PetscBool flag;
779: PetscMalloc1(m,&ns);
780: while (i < m) {
781: nzx = ai[i+1] - ai[i]; /* Number of nonzeros */
782: /* Limits the number of elements in a node to 'a->inode.limit' */
783: for (j=i+1,blk_size=1; j<m && blk_size <a->inode.limit; ++j,++blk_size) {
784: nzy = ai[j+1] - ai[j];
785: if (nzy != (nzx - j + i)) break;
786: PetscMemcmp(aj + ai[i] + j - i,aj + ai[j],nzy*sizeof(PetscInt),&flag);
787: if (!flag) break;
788: }
789: ns[node_count++] = blk_size;
791: i = j;
792: }
793: if (!a->inode.size && m && node_count > .9*m) {
794: PetscFree(ns);
795: PetscInfo2(A,"Found %D nodes out of %D rows. Not using Inode routines\n",node_count,m);
796: } else {
797: a->inode.node_count = node_count;
799: PetscMalloc1(node_count,&a->inode.size);
800: PetscLogObjectMemory((PetscObject)A,node_count*sizeof(PetscInt));
801: PetscMemcpy(a->inode.size,ns,node_count*sizeof(PetscInt));
802: PetscFree(ns);
803: PetscInfo3(A,"Found %D nodes of %D. Limit used: %D. Using Inode routines\n",node_count,m,a->inode.limit);
805: /* count collections of adjacent columns in each inode */
806: row = 0;
807: cnt = 0;
808: for (i=0; i<node_count; i++) {
809: cols = aj + ai[row] + a->inode.size[i];
810: nz = ai[row+1] - ai[row] - a->inode.size[i];
811: for (j=1; j<nz; j++) {
812: if (cols[j] != cols[j-1]+1) cnt++;
813: }
814: cnt++;
815: row += a->inode.size[i];
816: }
817: PetscMalloc1(2*cnt,&counts);
818: cnt = 0;
819: row = 0;
820: for (i=0; i<node_count; i++) {
821: cols = aj + ai[row] + a->inode.size[i];
822: counts[2*cnt] = cols[0];
823: nz = ai[row+1] - ai[row] - a->inode.size[i];
824: cnt2 = 1;
825: for (j=1; j<nz; j++) {
826: if (cols[j] != cols[j-1]+1) {
827: counts[2*(cnt++)+1] = cnt2;
828: counts[2*cnt] = cols[j];
829: cnt2 = 1;
830: } else cnt2++;
831: }
832: counts[2*(cnt++)+1] = cnt2;
833: row += a->inode.size[i];
834: }
835: PetscIntView(2*cnt,counts,0);
836: }
837: return(0);
838: }
842: PetscErrorCode MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode)
843: {
844: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
846: PetscInt fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
847: PetscInt m = A->rmap->N,*ip,N,*ailen = a->ilen;
848: PetscInt mbs = a->mbs,bs2 = a->bs2,rmax = 0;
849: MatScalar *aa = a->a,*ap;
852: if (mode == MAT_FLUSH_ASSEMBLY) return(0);
854: if (m) rmax = ailen[0];
855: for (i=1; i<mbs; i++) {
856: /* move each row back by the amount of empty slots (fshift) before it*/
857: fshift += imax[i-1] - ailen[i-1];
858: rmax = PetscMax(rmax,ailen[i]);
859: if (fshift) {
860: ip = aj + ai[i]; ap = aa + bs2*ai[i];
861: N = ailen[i];
862: for (j=0; j<N; j++) {
863: ip[j-fshift] = ip[j];
864: PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));
865: }
866: }
867: ai[i] = ai[i-1] + ailen[i-1];
868: }
869: if (mbs) {
870: fshift += imax[mbs-1] - ailen[mbs-1];
871: ai[mbs] = ai[mbs-1] + ailen[mbs-1];
872: }
873: /* reset ilen and imax for each row */
874: for (i=0; i<mbs; i++) {
875: ailen[i] = imax[i] = ai[i+1] - ai[i];
876: }
877: a->nz = ai[mbs];
879: /* diagonals may have moved, reset it */
880: if (a->diag) {
881: PetscMemcpy(a->diag,ai,mbs*sizeof(PetscInt));
882: }
883: if (fshift && a->nounused == -1) 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);
885: 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);
886: PetscInfo1(A,"Number of mallocs during MatSetValues is %D\n",a->reallocs);
887: PetscInfo1(A,"Most nonzeros blocks in any row is %D\n",rmax);
889: A->info.mallocs += a->reallocs;
890: a->reallocs = 0;
891: A->info.nz_unneeded = (PetscReal)fshift*bs2;
892: a->idiagvalid = PETSC_FALSE;
894: if (A->cmap->n < 65536 && A->cmap->bs == 1) {
895: if (a->jshort && a->free_jshort) {
896: /* when matrix data structure is changed, previous jshort must be replaced */
897: PetscFree(a->jshort);
898: }
899: PetscMalloc1(a->i[A->rmap->n],&a->jshort);
900: PetscLogObjectMemory((PetscObject)A,a->i[A->rmap->n]*sizeof(unsigned short));
901: for (i=0; i<a->i[A->rmap->n]; i++) a->jshort[i] = a->j[i];
902: A->ops->mult = MatMult_SeqSBAIJ_1_ushort;
903: A->ops->sor = MatSOR_SeqSBAIJ_ushort;
904: a->free_jshort = PETSC_TRUE;
905: }
906: return(0);
907: }
909: /*
910: This function returns an array of flags which indicate the locations of contiguous
911: blocks that should be zeroed. for eg: if bs = 3 and is = [0,1,2,3,5,6,7,8,9]
912: then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)]
913: Assume: sizes should be long enough to hold all the values.
914: */
917: PetscErrorCode MatZeroRows_SeqSBAIJ_Check_Blocks(PetscInt idx[],PetscInt n,PetscInt bs,PetscInt sizes[], PetscInt *bs_max)
918: {
919: PetscInt i,j,k,row;
920: PetscBool flg;
923: for (i=0,j=0; i<n; j++) {
924: row = idx[i];
925: if (row%bs!=0) { /* Not the begining of a block */
926: sizes[j] = 1;
927: i++;
928: } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */
929: sizes[j] = 1; /* Also makes sure atleast 'bs' values exist for next else */
930: i++;
931: } else { /* Begining of the block, so check if the complete block exists */
932: flg = PETSC_TRUE;
933: for (k=1; k<bs; k++) {
934: if (row+k != idx[i+k]) { /* break in the block */
935: flg = PETSC_FALSE;
936: break;
937: }
938: }
939: if (flg) { /* No break in the bs */
940: sizes[j] = bs;
941: i += bs;
942: } else {
943: sizes[j] = 1;
944: i++;
945: }
946: }
947: }
948: *bs_max = j;
949: return(0);
950: }
953: /* Only add/insert a(i,j) with i<=j (blocks).
954: Any a(i,j) with i>j input by user is ingored.
955: */
959: PetscErrorCode MatSetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
960: {
961: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
963: PetscInt *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,lastcol = -1;
964: PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented;
965: PetscInt *aj =a->j,nonew=a->nonew,bs=A->rmap->bs,brow,bcol;
966: PetscInt ridx,cidx,bs2=a->bs2;
967: MatScalar *ap,value,*aa=a->a,*bap;
971: for (k=0; k<m; k++) { /* loop over added rows */
972: row = im[k]; /* row number */
973: brow = row/bs; /* block row number */
974: if (row < 0) continue;
975: #if defined(PETSC_USE_DEBUG)
976: 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);
977: #endif
978: rp = aj + ai[brow]; /*ptr to beginning of column value of the row block*/
979: ap = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/
980: rmax = imax[brow]; /* maximum space allocated for this row */
981: nrow = ailen[brow]; /* actual length of this row */
982: low = 0;
984: for (l=0; l<n; l++) { /* loop over added columns */
985: if (in[l] < 0) continue;
986: #if defined(PETSC_USE_DEBUG)
987: 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);
988: #endif
989: col = in[l];
990: bcol = col/bs; /* block col number */
992: if (brow > bcol) {
993: if (a->ignore_ltriangular) continue; /* ignore lower triangular values */
994: else 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)");
995: }
997: ridx = row % bs; cidx = col % bs; /*row and col index inside the block */
998: if ((brow==bcol && ridx<=cidx) || (brow<bcol)) {
999: /* element value a(k,l) */
1000: if (roworiented) value = v[l + k*n];
1001: else value = v[k + l*m];
1003: /* move pointer bap to a(k,l) quickly and add/insert value */
1004: if (col <= lastcol) low = 0;
1005: high = nrow;
1006: lastcol = col;
1007: while (high-low > 7) {
1008: t = (low+high)/2;
1009: if (rp[t] > bcol) high = t;
1010: else low = t;
1011: }
1012: for (i=low; i<high; i++) {
1013: if (rp[i] > bcol) break;
1014: if (rp[i] == bcol) {
1015: bap = ap + bs2*i + bs*cidx + ridx;
1016: if (is == ADD_VALUES) *bap += value;
1017: else *bap = value;
1018: /* for diag block, add/insert its symmetric element a(cidx,ridx) */
1019: if (brow == bcol && ridx < cidx) {
1020: bap = ap + bs2*i + bs*ridx + cidx;
1021: if (is == ADD_VALUES) *bap += value;
1022: else *bap = value;
1023: }
1024: goto noinsert1;
1025: }
1026: }
1028: if (nonew == 1) goto noinsert1;
1029: if (nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col);
1030: MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar);
1032: N = nrow++ - 1; high++;
1033: /* shift up all the later entries in this row */
1034: for (ii=N; ii>=i; ii--) {
1035: rp[ii+1] = rp[ii];
1036: PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));
1037: }
1038: if (N>=i) {
1039: PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));
1040: }
1041: rp[i] = bcol;
1042: ap[bs2*i + bs*cidx + ridx] = value;
1043: noinsert1:;
1044: low = i;
1045: }
1046: } /* end of loop over added columns */
1047: ailen[brow] = nrow;
1048: } /* end of loop over added rows */
1049: return(0);
1050: }
1054: PetscErrorCode MatICCFactor_SeqSBAIJ(Mat inA,IS row,const MatFactorInfo *info)
1055: {
1056: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data;
1057: Mat outA;
1059: PetscBool row_identity;
1062: if (info->levels != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only levels=0 is supported for in-place icc");
1063: ISIdentity(row,&row_identity);
1064: if (!row_identity) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix reordering is not supported");
1065: 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()! */
1067: outA = inA;
1068: inA->factortype = MAT_FACTOR_ICC;
1070: MatMarkDiagonal_SeqSBAIJ(inA);
1071: MatSeqSBAIJSetNumericFactorization_inplace(inA,row_identity);
1073: PetscObjectReference((PetscObject)row);
1074: ISDestroy(&a->row);
1075: a->row = row;
1076: PetscObjectReference((PetscObject)row);
1077: ISDestroy(&a->col);
1078: a->col = row;
1080: /* Create the invert permutation so that it can be used in MatCholeskyFactorNumeric() */
1081: if (a->icol) {ISInvertPermutation(row,PETSC_DECIDE, &a->icol);}
1082: PetscLogObjectParent((PetscObject)inA,(PetscObject)a->icol);
1084: if (!a->solve_work) {
1085: PetscMalloc1((inA->rmap->N+inA->rmap->bs),&a->solve_work);
1086: PetscLogObjectMemory((PetscObject)inA,(inA->rmap->N+inA->rmap->bs)*sizeof(PetscScalar));
1087: }
1089: MatCholeskyFactorNumeric(outA,inA,info);
1090: return(0);
1091: }
1095: PetscErrorCode MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,PetscInt *indices)
1096: {
1097: Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ*)mat->data;
1098: PetscInt i,nz,n;
1102: nz = baij->maxnz;
1103: n = mat->cmap->n;
1104: for (i=0; i<nz; i++) baij->j[i] = indices[i];
1106: baij->nz = nz;
1107: for (i=0; i<n; i++) baij->ilen[i] = baij->imax[i];
1109: MatSetOption(mat,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
1110: return(0);
1111: }
1115: /*@
1116: MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows
1117: in the matrix.
1119: Input Parameters:
1120: + mat - the SeqSBAIJ matrix
1121: - indices - the column indices
1123: Level: advanced
1125: Notes:
1126: This can be called if you have precomputed the nonzero structure of the
1127: matrix and want to provide it to the matrix object to improve the performance
1128: of the MatSetValues() operation.
1130: You MUST have set the correct numbers of nonzeros per row in the call to
1131: MatCreateSeqSBAIJ(), and the columns indices MUST be sorted.
1133: MUST be called before any calls to MatSetValues()
1135: .seealso: MatCreateSeqSBAIJ
1136: @*/
1137: PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat mat,PetscInt *indices)
1138: {
1144: PetscUseMethod(mat,"MatSeqSBAIJSetColumnIndices_C",(Mat,PetscInt*),(mat,indices));
1145: return(0);
1146: }
1150: PetscErrorCode MatCopy_SeqSBAIJ(Mat A,Mat B,MatStructure str)
1151: {
1155: /* If the two matrices have the same copy implementation, use fast copy. */
1156: if (str == SAME_NONZERO_PATTERN && (A->ops->copy == B->ops->copy)) {
1157: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1158: Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data;
1160: 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");
1161: PetscMemcpy(b->a,a->a,(a->i[A->rmap->N])*sizeof(PetscScalar));
1162: } else {
1163: MatGetRowUpperTriangular(A);
1164: MatCopy_Basic(A,B,str);
1165: MatRestoreRowUpperTriangular(A);
1166: }
1167: return(0);
1168: }
1172: PetscErrorCode MatSetUp_SeqSBAIJ(Mat A)
1173: {
1177: MatSeqSBAIJSetPreallocation_SeqSBAIJ(A,A->rmap->bs,PETSC_DEFAULT,0);
1178: return(0);
1179: }
1183: PetscErrorCode MatSeqSBAIJGetArray_SeqSBAIJ(Mat A,PetscScalar *array[])
1184: {
1185: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1188: *array = a->a;
1189: return(0);
1190: }
1194: PetscErrorCode MatSeqSBAIJRestoreArray_SeqSBAIJ(Mat A,PetscScalar *array[])
1195: {
1197: return(0);
1198: }
1202: PetscErrorCode MatAXPY_SeqSBAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str)
1203: {
1204: Mat_SeqSBAIJ *x=(Mat_SeqSBAIJ*)X->data, *y=(Mat_SeqSBAIJ*)Y->data;
1206: PetscInt i,bs=Y->rmap->bs,bs2=bs*bs,j;
1207: PetscBLASInt one = 1;
1210: if (str == SAME_NONZERO_PATTERN) {
1211: PetscScalar alpha = a;
1212: PetscBLASInt bnz;
1213: PetscBLASIntCast(x->nz*bs2,&bnz);
1214: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&bnz,&alpha,x->a,&one,y->a,&one));
1215: } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
1216: if (y->xtoy && y->XtoY != X) {
1217: PetscFree(y->xtoy);
1218: MatDestroy(&y->XtoY);
1219: }
1220: if (!y->xtoy) { /* get xtoy */
1221: MatAXPYGetxtoy_Private(x->mbs,x->i,x->j,NULL, y->i,y->j,NULL, &y->xtoy);
1222: y->XtoY = X;
1223: }
1224: for (i=0; i<x->nz; i++) {
1225: j = 0;
1226: while (j < bs2) {
1227: y->a[bs2*y->xtoy[i]+j] += a*(x->a[bs2*i+j]);
1228: j++;
1229: }
1230: }
1231: PetscInfo3(Y,"ratio of nnz_s(X)/nnz_s(Y): %D/%D = %g\n",bs2*x->nz,bs2*y->nz,(double)((PetscReal)(bs2*x->nz)/(PetscReal)(bs2*y->nz)));
1232: } else {
1233: MatGetRowUpperTriangular(X);
1234: MatAXPY_Basic(Y,a,X,str);
1235: MatRestoreRowUpperTriangular(X);
1236: }
1237: return(0);
1238: }
1242: PetscErrorCode MatIsSymmetric_SeqSBAIJ(Mat A,PetscReal tol,PetscBool *flg)
1243: {
1245: *flg = PETSC_TRUE;
1246: return(0);
1247: }
1251: PetscErrorCode MatIsStructurallySymmetric_SeqSBAIJ(Mat A,PetscBool *flg)
1252: {
1254: *flg = PETSC_TRUE;
1255: return(0);
1256: }
1260: PetscErrorCode MatIsHermitian_SeqSBAIJ(Mat A,PetscReal tol,PetscBool *flg)
1261: {
1263: *flg = PETSC_FALSE;
1264: return(0);
1265: }
1269: PetscErrorCode MatRealPart_SeqSBAIJ(Mat A)
1270: {
1271: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1272: PetscInt i,nz = a->bs2*a->i[a->mbs];
1273: MatScalar *aa = a->a;
1276: for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]);
1277: return(0);
1278: }
1282: PetscErrorCode MatImaginaryPart_SeqSBAIJ(Mat A)
1283: {
1284: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1285: PetscInt i,nz = a->bs2*a->i[a->mbs];
1286: MatScalar *aa = a->a;
1289: for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]);
1290: return(0);
1291: }
1295: PetscErrorCode MatZeroRowsColumns_SeqSBAIJ(Mat A,PetscInt is_n,const PetscInt is_idx[],PetscScalar diag,Vec x, Vec b)
1296: {
1297: Mat_SeqSBAIJ *baij=(Mat_SeqSBAIJ*)A->data;
1298: PetscErrorCode ierr;
1299: PetscInt i,j,k,count;
1300: PetscInt bs =A->rmap->bs,bs2=baij->bs2,row,col;
1301: PetscScalar zero = 0.0;
1302: MatScalar *aa;
1303: const PetscScalar *xx;
1304: PetscScalar *bb;
1305: PetscBool *zeroed,vecs = PETSC_FALSE;
1308: /* fix right hand side if needed */
1309: if (x && b) {
1310: VecGetArrayRead(x,&xx);
1311: VecGetArray(b,&bb);
1312: vecs = PETSC_TRUE;
1313: }
1314: A->same_nonzero = PETSC_TRUE;
1316: /* zero the columns */
1317: PetscCalloc1(A->rmap->n,&zeroed);
1318: for (i=0; i<is_n; i++) {
1319: 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]);
1320: zeroed[is_idx[i]] = PETSC_TRUE;
1321: }
1322: if (vecs) {
1323: for (i=0; i<A->rmap->N; i++) {
1324: row = i/bs;
1325: for (j=baij->i[row]; j<baij->i[row+1]; j++) {
1326: for (k=0; k<bs; k++) {
1327: col = bs*baij->j[j] + k;
1328: if (col <= i) continue;
1329: aa = ((MatScalar*)(baij->a)) + j*bs2 + (i%bs) + bs*k;
1330: if (!zeroed[i] && zeroed[col]) bb[i] -= aa[0]*xx[col];
1331: if (zeroed[i] && !zeroed[col]) bb[col] -= aa[0]*xx[i];
1332: }
1333: }
1334: }
1335: for (i=0; i<is_n; i++) bb[is_idx[i]] = diag*xx[is_idx[i]];
1336: }
1338: for (i=0; i<A->rmap->N; i++) {
1339: if (!zeroed[i]) {
1340: row = i/bs;
1341: for (j=baij->i[row]; j<baij->i[row+1]; j++) {
1342: for (k=0; k<bs; k++) {
1343: col = bs*baij->j[j] + k;
1344: if (zeroed[col]) {
1345: aa = ((MatScalar*)(baij->a)) + j*bs2 + (i%bs) + bs*k;
1346: aa[0] = 0.0;
1347: }
1348: }
1349: }
1350: }
1351: }
1352: PetscFree(zeroed);
1353: if (vecs) {
1354: VecRestoreArrayRead(x,&xx);
1355: VecRestoreArray(b,&bb);
1356: }
1358: /* zero the rows */
1359: for (i=0; i<is_n; i++) {
1360: row = is_idx[i];
1361: count = (baij->i[row/bs +1] - baij->i[row/bs])*bs;
1362: aa = ((MatScalar*)(baij->a)) + baij->i[row/bs]*bs2 + (row%bs);
1363: for (k=0; k<count; k++) {
1364: aa[0] = zero;
1365: aa += bs;
1366: }
1367: if (diag != 0.0) {
1368: (*A->ops->setvalues)(A,1,&row,1,&row,&diag,INSERT_VALUES);
1369: }
1370: }
1371: MatAssemblyEnd_SeqSBAIJ(A,MAT_FINAL_ASSEMBLY);
1372: return(0);
1373: }
1375: /* -------------------------------------------------------------------*/
1376: static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ,
1377: MatGetRow_SeqSBAIJ,
1378: MatRestoreRow_SeqSBAIJ,
1379: MatMult_SeqSBAIJ_N,
1380: /* 4*/ MatMultAdd_SeqSBAIJ_N,
1381: MatMult_SeqSBAIJ_N, /* transpose versions are same as non-transpose versions */
1382: MatMultAdd_SeqSBAIJ_N,
1383: 0,
1384: 0,
1385: 0,
1386: /* 10*/ 0,
1387: 0,
1388: MatCholeskyFactor_SeqSBAIJ,
1389: MatSOR_SeqSBAIJ,
1390: MatTranspose_SeqSBAIJ,
1391: /* 15*/ MatGetInfo_SeqSBAIJ,
1392: MatEqual_SeqSBAIJ,
1393: MatGetDiagonal_SeqSBAIJ,
1394: MatDiagonalScale_SeqSBAIJ,
1395: MatNorm_SeqSBAIJ,
1396: /* 20*/ 0,
1397: MatAssemblyEnd_SeqSBAIJ,
1398: MatSetOption_SeqSBAIJ,
1399: MatZeroEntries_SeqSBAIJ,
1400: /* 24*/ 0,
1401: 0,
1402: 0,
1403: 0,
1404: 0,
1405: /* 29*/ MatSetUp_SeqSBAIJ,
1406: 0,
1407: 0,
1408: 0,
1409: 0,
1410: /* 34*/ MatDuplicate_SeqSBAIJ,
1411: 0,
1412: 0,
1413: 0,
1414: MatICCFactor_SeqSBAIJ,
1415: /* 39*/ MatAXPY_SeqSBAIJ,
1416: MatGetSubMatrices_SeqSBAIJ,
1417: MatIncreaseOverlap_SeqSBAIJ,
1418: MatGetValues_SeqSBAIJ,
1419: MatCopy_SeqSBAIJ,
1420: /* 44*/ 0,
1421: MatScale_SeqSBAIJ,
1422: 0,
1423: 0,
1424: MatZeroRowsColumns_SeqSBAIJ,
1425: /* 49*/ 0,
1426: MatGetRowIJ_SeqSBAIJ,
1427: MatRestoreRowIJ_SeqSBAIJ,
1428: 0,
1429: 0,
1430: /* 54*/ 0,
1431: 0,
1432: 0,
1433: 0,
1434: MatSetValuesBlocked_SeqSBAIJ,
1435: /* 59*/ MatGetSubMatrix_SeqSBAIJ,
1436: 0,
1437: 0,
1438: 0,
1439: 0,
1440: /* 64*/ 0,
1441: 0,
1442: 0,
1443: 0,
1444: 0,
1445: /* 69*/ MatGetRowMaxAbs_SeqSBAIJ,
1446: 0,
1447: 0,
1448: 0,
1449: 0,
1450: /* 74*/ 0,
1451: 0,
1452: 0,
1453: 0,
1454: 0,
1455: /* 79*/ 0,
1456: 0,
1457: 0,
1458: MatGetInertia_SeqSBAIJ,
1459: MatLoad_SeqSBAIJ,
1460: /* 84*/ MatIsSymmetric_SeqSBAIJ,
1461: MatIsHermitian_SeqSBAIJ,
1462: MatIsStructurallySymmetric_SeqSBAIJ,
1463: 0,
1464: 0,
1465: /* 89*/ 0,
1466: 0,
1467: 0,
1468: 0,
1469: 0,
1470: /* 94*/ 0,
1471: 0,
1472: 0,
1473: 0,
1474: 0,
1475: /* 99*/ 0,
1476: 0,
1477: 0,
1478: 0,
1479: 0,
1480: /*104*/ 0,
1481: MatRealPart_SeqSBAIJ,
1482: MatImaginaryPart_SeqSBAIJ,
1483: MatGetRowUpperTriangular_SeqSBAIJ,
1484: MatRestoreRowUpperTriangular_SeqSBAIJ,
1485: /*109*/ 0,
1486: 0,
1487: 0,
1488: 0,
1489: MatMissingDiagonal_SeqSBAIJ,
1490: /*114*/ 0,
1491: 0,
1492: 0,
1493: 0,
1494: 0,
1495: /*119*/ 0,
1496: 0,
1497: 0,
1498: 0,
1499: 0,
1500: /*124*/ 0,
1501: 0,
1502: 0,
1503: 0,
1504: 0,
1505: /*129*/ 0,
1506: 0,
1507: 0,
1508: 0,
1509: 0,
1510: /*134*/ 0,
1511: 0,
1512: 0,
1513: 0,
1514: 0,
1515: /*139*/ 0,
1516: 0,
1517: 0
1518: };
1522: PetscErrorCode MatStoreValues_SeqSBAIJ(Mat mat)
1523: {
1524: Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ*)mat->data;
1525: PetscInt nz = aij->i[mat->rmap->N]*mat->rmap->bs*aij->bs2;
1529: if (aij->nonew != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first");
1531: /* allocate space for values if not already there */
1532: if (!aij->saved_values) {
1533: PetscMalloc1((nz+1),&aij->saved_values);
1534: }
1536: /* copy values over */
1537: PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));
1538: return(0);
1539: }
1543: PetscErrorCode MatRetrieveValues_SeqSBAIJ(Mat mat)
1544: {
1545: Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ*)mat->data;
1547: PetscInt nz = aij->i[mat->rmap->N]*mat->rmap->bs*aij->bs2;
1550: if (aij->nonew != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first");
1551: if (!aij->saved_values) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Must call MatStoreValues(A);first");
1553: /* copy values over */
1554: PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));
1555: return(0);
1556: }
1560: PetscErrorCode MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat B,PetscInt bs,PetscInt nz,PetscInt *nnz)
1561: {
1562: Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data;
1564: PetscInt i,mbs,bs2;
1565: PetscBool skipallocation = PETSC_FALSE,flg = PETSC_FALSE,realalloc = PETSC_FALSE;
1568: if (nz >= 0 || nnz) realalloc = PETSC_TRUE;
1569: B->preallocated = PETSC_TRUE;
1571: PetscLayoutSetBlockSize(B->rmap,bs);
1572: PetscLayoutSetBlockSize(B->cmap,bs);
1573: PetscLayoutSetUp(B->rmap);
1574: PetscLayoutSetUp(B->cmap);
1575: PetscLayoutGetBlockSize(B->rmap,&bs);
1577: mbs = B->rmap->N/bs;
1578: bs2 = bs*bs;
1580: if (mbs*bs != B->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize");
1582: if (nz == MAT_SKIP_ALLOCATION) {
1583: skipallocation = PETSC_TRUE;
1584: nz = 0;
1585: }
1587: if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3;
1588: if (nz < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %D",nz);
1589: if (nnz) {
1590: for (i=0; i<mbs; i++) {
1591: 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]);
1592: 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);
1593: }
1594: }
1596: B->ops->mult = MatMult_SeqSBAIJ_N;
1597: B->ops->multadd = MatMultAdd_SeqSBAIJ_N;
1598: B->ops->multtranspose = MatMult_SeqSBAIJ_N;
1599: B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_N;
1601: PetscOptionsGetBool(((PetscObject)B)->prefix,"-mat_no_unroll",&flg,NULL);
1602: if (!flg) {
1603: switch (bs) {
1604: case 1:
1605: B->ops->mult = MatMult_SeqSBAIJ_1;
1606: B->ops->multadd = MatMultAdd_SeqSBAIJ_1;
1607: B->ops->multtranspose = MatMult_SeqSBAIJ_1;
1608: B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_1;
1609: break;
1610: case 2:
1611: B->ops->mult = MatMult_SeqSBAIJ_2;
1612: B->ops->multadd = MatMultAdd_SeqSBAIJ_2;
1613: B->ops->multtranspose = MatMult_SeqSBAIJ_2;
1614: B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_2;
1615: break;
1616: case 3:
1617: B->ops->mult = MatMult_SeqSBAIJ_3;
1618: B->ops->multadd = MatMultAdd_SeqSBAIJ_3;
1619: B->ops->multtranspose = MatMult_SeqSBAIJ_3;
1620: B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_3;
1621: break;
1622: case 4:
1623: B->ops->mult = MatMult_SeqSBAIJ_4;
1624: B->ops->multadd = MatMultAdd_SeqSBAIJ_4;
1625: B->ops->multtranspose = MatMult_SeqSBAIJ_4;
1626: B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_4;
1627: break;
1628: case 5:
1629: B->ops->mult = MatMult_SeqSBAIJ_5;
1630: B->ops->multadd = MatMultAdd_SeqSBAIJ_5;
1631: B->ops->multtranspose = MatMult_SeqSBAIJ_5;
1632: B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_5;
1633: break;
1634: case 6:
1635: B->ops->mult = MatMult_SeqSBAIJ_6;
1636: B->ops->multadd = MatMultAdd_SeqSBAIJ_6;
1637: B->ops->multtranspose = MatMult_SeqSBAIJ_6;
1638: B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_6;
1639: break;
1640: case 7:
1641: B->ops->mult = MatMult_SeqSBAIJ_7;
1642: B->ops->multadd = MatMultAdd_SeqSBAIJ_7;
1643: B->ops->multtranspose = MatMult_SeqSBAIJ_7;
1644: B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_7;
1645: break;
1646: }
1647: }
1649: b->mbs = mbs;
1650: b->nbs = mbs;
1651: if (!skipallocation) {
1652: if (!b->imax) {
1653: PetscMalloc2(mbs,&b->imax,mbs,&b->ilen);
1655: b->free_imax_ilen = PETSC_TRUE;
1657: PetscLogObjectMemory((PetscObject)B,2*mbs*sizeof(PetscInt));
1658: }
1659: if (!nnz) {
1660: if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
1661: else if (nz <= 0) nz = 1;
1662: for (i=0; i<mbs; i++) b->imax[i] = nz;
1663: nz = nz*mbs; /* total nz */
1664: } else {
1665: nz = 0;
1666: for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
1667: }
1668: /* b->ilen will count nonzeros in each block row so far. */
1669: for (i=0; i<mbs; i++) b->ilen[i] = 0;
1670: /* nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */
1672: /* allocate the matrix space */
1673: MatSeqXAIJFreeAIJ(B,&b->a,&b->j,&b->i);
1674: PetscMalloc3(bs2*nz,&b->a,nz,&b->j,B->rmap->N+1,&b->i);
1675: PetscLogObjectMemory((PetscObject)B,(B->rmap->N+1)*sizeof(PetscInt)+nz*(bs2*sizeof(PetscScalar)+sizeof(PetscInt)));
1676: PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));
1677: PetscMemzero(b->j,nz*sizeof(PetscInt));
1679: b->singlemalloc = PETSC_TRUE;
1681: /* pointer to beginning of each row */
1682: b->i[0] = 0;
1683: for (i=1; i<mbs+1; i++) b->i[i] = b->i[i-1] + b->imax[i-1];
1685: b->free_a = PETSC_TRUE;
1686: b->free_ij = PETSC_TRUE;
1687: } else {
1688: b->free_a = PETSC_FALSE;
1689: b->free_ij = PETSC_FALSE;
1690: }
1692: B->rmap->bs = bs;
1693: b->bs2 = bs2;
1694: b->nz = 0;
1695: b->maxnz = nz;
1697: b->inew = 0;
1698: b->jnew = 0;
1699: b->anew = 0;
1700: b->a2anew = 0;
1701: b->permute = PETSC_FALSE;
1702: if (realalloc) {MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);}
1703: return(0);
1704: }
1708: PetscErrorCode MatSeqSBAIJSetPreallocationCSR_SeqSBAIJ(Mat B,PetscInt bs,const PetscInt ii[],const PetscInt jj[], const PetscScalar V[])
1709: {
1710: PetscInt i,j,m,nz,nz_max=0,*nnz;
1711: PetscScalar *values=0;
1712: PetscBool roworiented = ((Mat_SeqSBAIJ*)B->data)->roworiented;
1715: if (bs < 1) SETERRQ1(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_OUTOFRANGE,"Invalid block size specified, must be positive but it is %D",bs);
1716: PetscLayoutSetBlockSize(B->rmap,bs);
1717: PetscLayoutSetBlockSize(B->cmap,bs);
1718: PetscLayoutSetUp(B->rmap);
1719: PetscLayoutSetUp(B->cmap);
1720: PetscLayoutGetBlockSize(B->rmap,&bs);
1721: m = B->rmap->n/bs;
1723: if (ii[0]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"ii[0] must be 0 but it is %D",ii[0]);
1724: PetscMalloc1((m+1),&nnz);
1725: for (i=0; i<m; i++) {
1726: nz = ii[i+1] - ii[i];
1727: if (nz < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row %D has a negative number of columns %D",i,nz);
1728: nz_max = PetscMax(nz_max,nz);
1729: nnz[i] = nz;
1730: }
1731: MatSeqSBAIJSetPreallocation(B,bs,0,nnz);
1732: PetscFree(nnz);
1734: values = (PetscScalar*)V;
1735: if (!values) {
1736: PetscCalloc1(bs*bs*nz_max,&values);
1737: }
1738: for (i=0; i<m; i++) {
1739: PetscInt ncols = ii[i+1] - ii[i];
1740: const PetscInt *icols = jj + ii[i];
1741: if (!roworiented || bs == 1) {
1742: const PetscScalar *svals = values + (V ? (bs*bs*ii[i]) : 0);
1743: MatSetValuesBlocked_SeqSBAIJ(B,1,&i,ncols,icols,svals,INSERT_VALUES);
1744: } else {
1745: for (j=0; j<ncols; j++) {
1746: const PetscScalar *svals = values + (V ? (bs*bs*(ii[i]+j)) : 0);
1747: MatSetValuesBlocked_SeqSBAIJ(B,1,&i,1,&icols[j],svals,INSERT_VALUES);
1748: }
1749: }
1750: }
1751: if (!V) { PetscFree(values); }
1752: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
1753: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
1754: MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
1755: return(0);
1756: }
1758: /*
1759: This is used to set the numeric factorization for both Cholesky and ICC symbolic factorization
1760: */
1763: PetscErrorCode MatSeqSBAIJSetNumericFactorization_inplace(Mat B,PetscBool natural)
1764: {
1766: PetscBool flg = PETSC_FALSE;
1767: PetscInt bs = B->rmap->bs;
1770: PetscOptionsGetBool(((PetscObject)B)->prefix,"-mat_no_unroll",&flg,NULL);
1771: if (flg) bs = 8;
1773: if (!natural) {
1774: switch (bs) {
1775: case 1:
1776: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_inplace;
1777: break;
1778: case 2:
1779: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2;
1780: break;
1781: case 3:
1782: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3;
1783: break;
1784: case 4:
1785: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4;
1786: break;
1787: case 5:
1788: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5;
1789: break;
1790: case 6:
1791: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6;
1792: break;
1793: case 7:
1794: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7;
1795: break;
1796: default:
1797: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N;
1798: break;
1799: }
1800: } else {
1801: switch (bs) {
1802: case 1:
1803: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering_inplace;
1804: break;
1805: case 2:
1806: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering;
1807: break;
1808: case 3:
1809: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering;
1810: break;
1811: case 4:
1812: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering;
1813: break;
1814: case 5:
1815: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering;
1816: break;
1817: case 6:
1818: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering;
1819: break;
1820: case 7:
1821: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering;
1822: break;
1823: default:
1824: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering;
1825: break;
1826: }
1827: }
1828: return(0);
1829: }
1831: PETSC_EXTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqAIJ(Mat, MatType,MatReuse,Mat*);
1832: PETSC_EXTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqBAIJ(Mat, MatType,MatReuse,Mat*);
1836: PETSC_EXTERN PetscErrorCode MatGetFactor_seqsbaij_petsc(Mat A,MatFactorType ftype,Mat *B)
1837: {
1838: PetscInt n = A->rmap->n;
1842: #if defined(PETSC_USE_COMPLEX)
1843: if (A->hermitian) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Hermitian Factor is not supported");
1844: #endif
1845: MatCreate(PetscObjectComm((PetscObject)A),B);
1846: MatSetSizes(*B,n,n,n,n);
1847: if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) {
1848: MatSetType(*B,MATSEQSBAIJ);
1849: MatSeqSBAIJSetPreallocation(*B,A->rmap->bs,MAT_SKIP_ALLOCATION,NULL);
1851: (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqSBAIJ;
1852: (*B)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqSBAIJ;
1853: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
1854: (*B)->factortype = ftype;
1855: return(0);
1856: }
1860: PetscErrorCode MatGetFactorAvailable_seqsbaij_petsc(Mat A,MatFactorType ftype,PetscBool *flg)
1861: {
1863: if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) {
1864: *flg = PETSC_TRUE;
1865: } else {
1866: *flg = PETSC_FALSE;
1867: }
1868: return(0);
1869: }
1871: #if defined(PETSC_HAVE_MUMPS)
1872: PETSC_EXTERN PetscErrorCode MatGetFactor_sbaij_mumps(Mat,MatFactorType,Mat*);
1873: #endif
1874: #if defined(PETSC_HAVE_PASTIX)
1875: PETSC_EXTERN PetscErrorCode MatGetFactor_seqsbaij_pastix(Mat,MatFactorType,Mat*);
1876: #endif
1877: #if defined(PETSC_HAVE_CHOLMOD)
1878: PETSC_EXTERN PetscErrorCode MatGetFactor_seqsbaij_cholmod(Mat,MatFactorType,Mat*);
1879: #endif
1880: PETSC_EXTERN PetscErrorCode MatGetFactor_seqsbaij_sbstrm(Mat,MatFactorType,Mat*);
1882: /*MC
1883: MATSEQSBAIJ - MATSEQSBAIJ = "seqsbaij" - A matrix type to be used for sequential symmetric block sparse matrices,
1884: based on block compressed sparse row format. Only the upper triangular portion of the matrix is stored.
1886: For complex numbers by default this matrix is symmetric, NOT Hermitian symmetric. To make it Hermitian symmetric you
1887: can call MatSetOption(Mat, MAT_HERMITIAN); after MatAssemblyEnd()
1889: Options Database Keys:
1890: . -mat_type seqsbaij - sets the matrix type to "seqsbaij" during a call to MatSetFromOptions()
1892: Notes: By default if you insert values into the lower triangular part of the matrix they are simply ignored (since they are not
1893: stored and it is assumed they symmetric to the upper triangular). If you call MatSetOption(Mat,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_FALSE) or use
1894: 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.
1897: Level: beginner
1899: .seealso: MatCreateSeqSBAIJ
1900: M*/
1902: PETSC_EXTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqSBSTRM(Mat, MatType,MatReuse,Mat*);
1906: PETSC_EXTERN PetscErrorCode MatCreate_SeqSBAIJ(Mat B)
1907: {
1908: Mat_SeqSBAIJ *b;
1910: PetscMPIInt size;
1911: PetscBool no_unroll = PETSC_FALSE,no_inode = PETSC_FALSE;
1914: MPI_Comm_size(PetscObjectComm((PetscObject)B),&size);
1915: if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
1917: PetscNewLog(B,&b);
1918: B->data = (void*)b;
1919: PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
1921: B->ops->destroy = MatDestroy_SeqSBAIJ;
1922: B->ops->view = MatView_SeqSBAIJ;
1923: b->row = 0;
1924: b->icol = 0;
1925: b->reallocs = 0;
1926: b->saved_values = 0;
1927: b->inode.limit = 5;
1928: b->inode.max_limit = 5;
1930: b->roworiented = PETSC_TRUE;
1931: b->nonew = 0;
1932: b->diag = 0;
1933: b->solve_work = 0;
1934: b->mult_work = 0;
1935: B->spptr = 0;
1936: B->info.nz_unneeded = (PetscReal)b->maxnz*b->bs2;
1937: b->keepnonzeropattern = PETSC_FALSE;
1938: b->xtoy = 0;
1939: b->XtoY = 0;
1941: b->inew = 0;
1942: b->jnew = 0;
1943: b->anew = 0;
1944: b->a2anew = 0;
1945: b->permute = PETSC_FALSE;
1947: b->ignore_ltriangular = PETSC_TRUE;
1949: PetscOptionsGetBool(((PetscObject)B)->prefix,"-mat_ignore_lower_triangular",&b->ignore_ltriangular,NULL);
1951: b->getrow_utriangular = PETSC_FALSE;
1953: PetscOptionsGetBool(((PetscObject)B)->prefix,"-mat_getrow_uppertriangular",&b->getrow_utriangular,NULL);
1955: #if defined(PETSC_HAVE_PASTIX)
1956: PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_pastix_C",MatGetFactor_seqsbaij_pastix);
1957: #endif
1958: #if defined(PETSC_HAVE_MUMPS)
1959: PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_mumps_C",MatGetFactor_sbaij_mumps);
1960: #endif
1961: #if defined(PETSC_HAVE_CHOLMOD)
1962: PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_cholmod_C",MatGetFactor_seqsbaij_cholmod);
1963: #endif
1964: PetscObjectComposeFunction((PetscObject)B,"MatGetFactorAvailable_petsc_C",MatGetFactorAvailable_seqsbaij_petsc);
1965: PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_petsc_C",MatGetFactor_seqsbaij_petsc);
1966: PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_sbstrm_C",MatGetFactor_seqsbaij_sbstrm);
1967: PetscObjectComposeFunction((PetscObject)B,"MatStoreValues_C",MatStoreValues_SeqSBAIJ);
1968: PetscObjectComposeFunction((PetscObject)B,"MatRetrieveValues_C",MatRetrieveValues_SeqSBAIJ);
1969: PetscObjectComposeFunction((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C",MatSeqSBAIJSetColumnIndices_SeqSBAIJ);
1970: PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_seqaij_C",MatConvert_SeqSBAIJ_SeqAIJ);
1971: PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_seqbaij_C",MatConvert_SeqSBAIJ_SeqBAIJ);
1972: PetscObjectComposeFunction((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",MatSeqSBAIJSetPreallocation_SeqSBAIJ);
1973: PetscObjectComposeFunction((PetscObject)B,"MatSeqSBAIJSetPreallocationCSR_C",MatSeqSBAIJSetPreallocationCSR_SeqSBAIJ);
1974: PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_seqsbstrm_C",MatConvert_SeqSBAIJ_SeqSBSTRM);
1976: B->symmetric = PETSC_TRUE;
1977: B->structurally_symmetric = PETSC_TRUE;
1978: B->symmetric_set = PETSC_TRUE;
1979: B->structurally_symmetric_set = PETSC_TRUE;
1981: PetscObjectChangeTypeName((PetscObject)B,MATSEQSBAIJ);
1983: PetscOptionsBegin(PetscObjectComm((PetscObject)B),((PetscObject)B)->prefix,"Options for SEQSBAIJ matrix","Mat");
1984: PetscOptionsBool("-mat_no_unroll","Do not optimize for inodes (slower)",NULL,no_unroll,&no_unroll,NULL);
1985: if (no_unroll) {
1986: PetscInfo(B,"Not using Inode routines due to -mat_no_unroll\n");
1987: }
1988: PetscOptionsBool("-mat_no_inode","Do not optimize for inodes (slower)",NULL,no_inode,&no_inode,NULL);
1989: if (no_inode) {
1990: PetscInfo(B,"Not using Inode routines due to -mat_no_inode\n");
1991: }
1992: PetscOptionsInt("-mat_inode_limit","Do not use inodes larger then this value",NULL,b->inode.limit,&b->inode.limit,NULL);
1993: PetscOptionsEnd();
1994: b->inode.use = (PetscBool)(!(no_unroll || no_inode));
1995: if (b->inode.limit > b->inode.max_limit) b->inode.limit = b->inode.max_limit;
1996: return(0);
1997: }
2001: /*@C
2002: MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block
2003: compressed row) format. For good matrix assembly performance the
2004: user should preallocate the matrix storage by setting the parameter nz
2005: (or the array nnz). By setting these parameters accurately, performance
2006: during matrix assembly can be increased by more than a factor of 50.
2008: Collective on Mat
2010: Input Parameters:
2011: + A - the symmetric matrix
2012: . bs - size of block
2013: . nz - number of block nonzeros per block row (same for all rows)
2014: - nnz - array containing the number of block nonzeros in the upper triangular plus
2015: diagonal portion of each block (possibly different for each block row) or NULL
2017: Options Database Keys:
2018: . -mat_no_unroll - uses code that does not unroll the loops in the
2019: block calculations (much slower)
2020: . -mat_block_size - size of the blocks to use (only works if a negative bs is passed in
2022: Level: intermediate
2024: Notes:
2025: Specify the preallocated storage with either nz or nnz (not both).
2026: Set nz=PETSC_DEFAULT and nnz=NULL for PETSc to control dynamic memory
2027: allocation. See the <a href="../../docs/manual.pdf#nameddest=ch_mat">Mat chapter of the users manual</a> for details.
2029: You can call MatGetInfo() to get information on how effective the preallocation was;
2030: for example the fields mallocs,nz_allocated,nz_used,nz_unneeded;
2031: You can also run with the option -info and look for messages with the string
2032: malloc in them to see if additional memory allocation was needed.
2034: If the nnz parameter is given then the nz parameter is ignored
2037: .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateSBAIJ()
2038: @*/
2039: PetscErrorCode MatSeqSBAIJSetPreallocation(Mat B,PetscInt bs,PetscInt nz,const PetscInt nnz[])
2040: {
2047: PetscTryMethod(B,"MatSeqSBAIJSetPreallocation_C",(Mat,PetscInt,PetscInt,const PetscInt[]),(B,bs,nz,nnz));
2048: return(0);
2049: }
2051: #undef __FUNCT__
2053: /*@C
2054: MatSeqSBAIJSetPreallocationCSR - Allocates memory for a sparse sequential matrix in symmetric block AIJ format.
2056: Input Parameters:
2057: + A - the matrix
2058: . i - the indices into j for the start of each local row (starts with zero)
2059: . j - the column indices for each local row (starts with zero) these must be sorted for each row
2060: - v - optional values in the matrix
2062: Level: developer
2064: Notes:
2065: The order of the entries in values is specified by the MatOption MAT_ROW_ORIENTED. For example, C programs
2066: may want to use the default MAT_ROW_ORIENTED=PETSC_TRUE and use an array v[nnz][bs][bs] where the second index is
2067: over rows within a block and the last index is over columns within a block row. Fortran programs will likely set
2068: MAT_ROW_ORIENTED=PETSC_FALSE and use a Fortran array v(bs,bs,nnz) in which the first index is over rows within a
2069: block column and the second index is over columns within a block.
2071: .keywords: matrix, block, aij, compressed row, sparse
2073: .seealso: MatCreate(), MatCreateSeqSBAIJ(), MatSetValuesBlocked(), MatSeqSBAIJSetPreallocation(), MATSEQSBAIJ
2074: @*/
2075: PetscErrorCode MatSeqSBAIJSetPreallocationCSR(Mat B,PetscInt bs,const PetscInt i[],const PetscInt j[], const PetscScalar v[])
2076: {
2083: PetscTryMethod(B,"MatSeqSBAIJSetPreallocationCSR_C",(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]),(B,bs,i,j,v));
2084: return(0);
2085: }
2089: /*@C
2090: MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block
2091: compressed row) format. For good matrix assembly performance the
2092: user should preallocate the matrix storage by setting the parameter nz
2093: (or the array nnz). By setting these parameters accurately, performance
2094: during matrix assembly can be increased by more than a factor of 50.
2096: Collective on MPI_Comm
2098: Input Parameters:
2099: + comm - MPI communicator, set to PETSC_COMM_SELF
2100: . bs - size of block
2101: . m - number of rows, or number of columns
2102: . nz - number of block nonzeros per block row (same for all rows)
2103: - nnz - array containing the number of block nonzeros in the upper triangular plus
2104: diagonal portion of each block (possibly different for each block row) or NULL
2106: Output Parameter:
2107: . A - the symmetric matrix
2109: Options Database Keys:
2110: . -mat_no_unroll - uses code that does not unroll the loops in the
2111: block calculations (much slower)
2112: . -mat_block_size - size of the blocks to use
2114: Level: intermediate
2116: It is recommended that one use the MatCreate(), MatSetType() and/or MatSetFromOptions(),
2117: MatXXXXSetPreallocation() paradgm instead of this routine directly.
2118: [MatXXXXSetPreallocation() is, for example, MatSeqAIJSetPreallocation]
2120: Notes:
2121: The number of rows and columns must be divisible by blocksize.
2122: This matrix type does not support complex Hermitian operation.
2124: Specify the preallocated storage with either nz or nnz (not both).
2125: Set nz=PETSC_DEFAULT and nnz=NULL for PETSc to control dynamic memory
2126: allocation. See the <a href="../../docs/manual.pdf#nameddest=ch_mat">Mat chapter of the users manual</a> for details.
2128: If the nnz parameter is given then the nz parameter is ignored
2130: .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateSBAIJ()
2131: @*/
2132: PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A)
2133: {
2137: MatCreate(comm,A);
2138: MatSetSizes(*A,m,n,m,n);
2139: MatSetType(*A,MATSEQSBAIJ);
2140: MatSeqSBAIJSetPreallocation_SeqSBAIJ(*A,bs,nz,(PetscInt*)nnz);
2141: return(0);
2142: }
2146: PetscErrorCode MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
2147: {
2148: Mat C;
2149: Mat_SeqSBAIJ *c,*a = (Mat_SeqSBAIJ*)A->data;
2151: PetscInt i,mbs = a->mbs,nz = a->nz,bs2 =a->bs2;
2154: if (a->i[mbs] != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Corrupt matrix");
2156: *B = 0;
2157: MatCreate(PetscObjectComm((PetscObject)A),&C);
2158: MatSetSizes(C,A->rmap->N,A->cmap->n,A->rmap->N,A->cmap->n);
2159: MatSetType(C,MATSEQSBAIJ);
2160: PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));
2161: c = (Mat_SeqSBAIJ*)C->data;
2163: C->preallocated = PETSC_TRUE;
2164: C->factortype = A->factortype;
2165: c->row = 0;
2166: c->icol = 0;
2167: c->saved_values = 0;
2168: c->keepnonzeropattern = a->keepnonzeropattern;
2169: C->assembled = PETSC_TRUE;
2171: PetscLayoutReference(A->rmap,&C->rmap);
2172: PetscLayoutReference(A->cmap,&C->cmap);
2173: c->bs2 = a->bs2;
2174: c->mbs = a->mbs;
2175: c->nbs = a->nbs;
2177: if (cpvalues == MAT_SHARE_NONZERO_PATTERN) {
2178: c->imax = a->imax;
2179: c->ilen = a->ilen;
2180: c->free_imax_ilen = PETSC_FALSE;
2181: } else {
2182: PetscMalloc2((mbs+1),&c->imax,(mbs+1),&c->ilen);
2183: PetscLogObjectMemory((PetscObject)C,2*(mbs+1)*sizeof(PetscInt));
2184: for (i=0; i<mbs; i++) {
2185: c->imax[i] = a->imax[i];
2186: c->ilen[i] = a->ilen[i];
2187: }
2188: c->free_imax_ilen = PETSC_TRUE;
2189: }
2191: /* allocate the matrix space */
2192: if (cpvalues == MAT_SHARE_NONZERO_PATTERN) {
2193: PetscMalloc1(bs2*nz,&c->a);
2194: PetscLogObjectMemory((PetscObject)C,nz*bs2*sizeof(MatScalar));
2195: c->i = a->i;
2196: c->j = a->j;
2197: c->singlemalloc = PETSC_FALSE;
2198: c->free_a = PETSC_TRUE;
2199: c->free_ij = PETSC_FALSE;
2200: c->parent = A;
2201: PetscObjectReference((PetscObject)A);
2202: MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
2203: MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
2204: } else {
2205: PetscMalloc3(bs2*nz,&c->a,nz,&c->j,mbs+1,&c->i);
2206: PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(PetscInt));
2207: PetscLogObjectMemory((PetscObject)C,(mbs+1)*sizeof(PetscInt) + nz*(bs2*sizeof(MatScalar) + sizeof(PetscInt)));
2208: c->singlemalloc = PETSC_TRUE;
2209: c->free_a = PETSC_TRUE;
2210: c->free_ij = PETSC_TRUE;
2211: }
2212: if (mbs > 0) {
2213: if (cpvalues != MAT_SHARE_NONZERO_PATTERN) {
2214: PetscMemcpy(c->j,a->j,nz*sizeof(PetscInt));
2215: }
2216: if (cpvalues == MAT_COPY_VALUES) {
2217: PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));
2218: } else {
2219: PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));
2220: }
2221: if (a->jshort) {
2222: /* cannot share jshort, it is reallocated in MatAssemblyEnd_SeqSBAIJ() */
2223: /* if the parent matrix is reassembled, this child matrix will never notice */
2224: PetscMalloc1(nz,&c->jshort);
2225: PetscLogObjectMemory((PetscObject)C,nz*sizeof(unsigned short));
2226: PetscMemcpy(c->jshort,a->jshort,nz*sizeof(unsigned short));
2228: c->free_jshort = PETSC_TRUE;
2229: }
2230: }
2232: c->roworiented = a->roworiented;
2233: c->nonew = a->nonew;
2235: if (a->diag) {
2236: if (cpvalues == MAT_SHARE_NONZERO_PATTERN) {
2237: c->diag = a->diag;
2238: c->free_diag = PETSC_FALSE;
2239: } else {
2240: PetscMalloc1(mbs,&c->diag);
2241: PetscLogObjectMemory((PetscObject)C,mbs*sizeof(PetscInt));
2242: for (i=0; i<mbs; i++) c->diag[i] = a->diag[i];
2243: c->free_diag = PETSC_TRUE;
2244: }
2245: }
2246: c->nz = a->nz;
2247: c->maxnz = a->nz; /* Since we allocate exactly the right amount */
2248: c->solve_work = 0;
2249: c->mult_work = 0;
2251: *B = C;
2252: PetscFunctionListDuplicate(((PetscObject)A)->qlist,&((PetscObject)C)->qlist);
2253: return(0);
2254: }
2258: PetscErrorCode MatLoad_SeqSBAIJ(Mat newmat,PetscViewer viewer)
2259: {
2260: Mat_SeqSBAIJ *a;
2262: int fd;
2263: PetscMPIInt size;
2264: PetscInt i,nz,header[4],*rowlengths=0,M,N,bs=1;
2265: PetscInt *mask,mbs,*jj,j,rowcount,nzcount,k,*s_browlengths,maskcount;
2266: PetscInt kmax,jcount,block,idx,point,nzcountb,extra_rows,rows,cols;
2267: PetscInt *masked,nmask,tmp,bs2,ishift;
2268: PetscScalar *aa;
2269: MPI_Comm comm;
2272: PetscObjectGetComm((PetscObject)viewer,&comm);
2273: PetscOptionsGetInt(((PetscObject)newmat)->prefix,"-matload_block_size",&bs,NULL);
2274: bs2 = bs*bs;
2276: MPI_Comm_size(comm,&size);
2277: if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"view must have one processor");
2278: PetscViewerBinaryGetDescriptor(viewer,&fd);
2279: PetscBinaryRead(fd,header,4,PETSC_INT);
2280: if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"not Mat object");
2281: M = header[1]; N = header[2]; nz = header[3];
2283: if (header[3] < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ");
2285: if (M != N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Can only do square matrices");
2287: /*
2288: This code adds extra rows to make sure the number of rows is
2289: divisible by the blocksize
2290: */
2291: mbs = M/bs;
2292: extra_rows = bs - M + bs*(mbs);
2293: if (extra_rows == bs) extra_rows = 0;
2294: else mbs++;
2295: if (extra_rows) {
2296: PetscInfo(viewer,"Padding loaded matrix to match blocksize\n");
2297: }
2299: /* Set global sizes if not already set */
2300: if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) {
2301: MatSetSizes(newmat,PETSC_DECIDE,PETSC_DECIDE,M+extra_rows,N+extra_rows);
2302: } else { /* Check if the matrix global sizes are correct */
2303: MatGetSize(newmat,&rows,&cols);
2304: 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);
2305: }
2307: /* read in row lengths */
2308: PetscMalloc1((M+extra_rows),&rowlengths);
2309: PetscBinaryRead(fd,rowlengths,M,PETSC_INT);
2310: for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;
2312: /* read in column indices */
2313: PetscMalloc1((nz+extra_rows),&jj);
2314: PetscBinaryRead(fd,jj,nz,PETSC_INT);
2315: for (i=0; i<extra_rows; i++) jj[nz+i] = M+i;
2317: /* loop over row lengths determining block row lengths */
2318: PetscCalloc1(mbs,&s_browlengths);
2319: PetscMalloc2(mbs,&mask,mbs,&masked);
2320: PetscMemzero(mask,mbs*sizeof(PetscInt));
2321: rowcount = 0;
2322: nzcount = 0;
2323: for (i=0; i<mbs; i++) {
2324: nmask = 0;
2325: for (j=0; j<bs; j++) {
2326: kmax = rowlengths[rowcount];
2327: for (k=0; k<kmax; k++) {
2328: tmp = jj[nzcount++]/bs; /* block col. index */
2329: if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;}
2330: }
2331: rowcount++;
2332: }
2333: s_browlengths[i] += nmask;
2335: /* zero out the mask elements we set */
2336: for (j=0; j<nmask; j++) mask[masked[j]] = 0;
2337: }
2339: /* Do preallocation */
2340: MatSeqSBAIJSetPreallocation_SeqSBAIJ(newmat,bs,0,s_browlengths);
2341: a = (Mat_SeqSBAIJ*)newmat->data;
2343: /* set matrix "i" values */
2344: a->i[0] = 0;
2345: for (i=1; i<= mbs; i++) {
2346: a->i[i] = a->i[i-1] + s_browlengths[i-1];
2347: a->ilen[i-1] = s_browlengths[i-1];
2348: }
2349: a->nz = a->i[mbs];
2351: /* read in nonzero values */
2352: PetscMalloc1((nz+extra_rows),&aa);
2353: PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);
2354: for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0;
2356: /* set "a" and "j" values into matrix */
2357: nzcount = 0; jcount = 0;
2358: for (i=0; i<mbs; i++) {
2359: nzcountb = nzcount;
2360: nmask = 0;
2361: for (j=0; j<bs; j++) {
2362: kmax = rowlengths[i*bs+j];
2363: for (k=0; k<kmax; k++) {
2364: tmp = jj[nzcount++]/bs; /* block col. index */
2365: if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;}
2366: }
2367: }
2368: /* sort the masked values */
2369: PetscSortInt(nmask,masked);
2371: /* set "j" values into matrix */
2372: maskcount = 1;
2373: for (j=0; j<nmask; j++) {
2374: a->j[jcount++] = masked[j];
2375: mask[masked[j]] = maskcount++;
2376: }
2378: /* set "a" values into matrix */
2379: ishift = bs2*a->i[i];
2380: for (j=0; j<bs; j++) {
2381: kmax = rowlengths[i*bs+j];
2382: for (k=0; k<kmax; k++) {
2383: tmp = jj[nzcountb]/bs; /* block col. index */
2384: if (tmp >= i) {
2385: block = mask[tmp] - 1;
2386: point = jj[nzcountb] - bs*tmp;
2387: idx = ishift + bs2*block + j + bs*point;
2388: a->a[idx] = aa[nzcountb];
2389: }
2390: nzcountb++;
2391: }
2392: }
2393: /* zero out the mask elements we set */
2394: for (j=0; j<nmask; j++) mask[masked[j]] = 0;
2395: }
2396: if (jcount != a->nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix");
2398: PetscFree(rowlengths);
2399: PetscFree(s_browlengths);
2400: PetscFree(aa);
2401: PetscFree(jj);
2402: PetscFree2(mask,masked);
2404: MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);
2405: MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);
2406: return(0);
2407: }
2411: /*@
2412: MatCreateSeqSBAIJWithArrays - Creates an sequential SBAIJ matrix using matrix elements
2413: (upper triangular entries in CSR format) provided by the user.
2415: Collective on MPI_Comm
2417: Input Parameters:
2418: + comm - must be an MPI communicator of size 1
2419: . bs - size of block
2420: . m - number of rows
2421: . n - number of columns
2422: . i - row indices
2423: . j - column indices
2424: - a - matrix values
2426: Output Parameter:
2427: . mat - the matrix
2429: Level: advanced
2431: Notes:
2432: The i, j, and a arrays are not copied by this routine, the user must free these arrays
2433: once the matrix is destroyed
2435: You cannot set new nonzero locations into this matrix, that will generate an error.
2437: The i and j indices are 0 based
2439: 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
2440: it is the regular CSR format excluding the lower triangular elements.
2442: .seealso: MatCreate(), MatCreateSBAIJ(), MatCreateSeqSBAIJ()
2444: @*/
2445: PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt *i,PetscInt *j,PetscScalar *a,Mat *mat)
2446: {
2448: PetscInt ii;
2449: Mat_SeqSBAIJ *sbaij;
2452: if (bs != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"block size %D > 1 is not supported yet",bs);
2453: if (i[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"i (row indices) must start with 0");
2455: MatCreate(comm,mat);
2456: MatSetSizes(*mat,m,n,m,n);
2457: MatSetType(*mat,MATSEQSBAIJ);
2458: MatSeqSBAIJSetPreallocation_SeqSBAIJ(*mat,bs,MAT_SKIP_ALLOCATION,0);
2459: sbaij = (Mat_SeqSBAIJ*)(*mat)->data;
2460: PetscMalloc2(m,&sbaij->imax,m,&sbaij->ilen);
2461: PetscLogObjectMemory((PetscObject)*mat,2*m*sizeof(PetscInt));
2463: sbaij->i = i;
2464: sbaij->j = j;
2465: sbaij->a = a;
2467: sbaij->singlemalloc = PETSC_FALSE;
2468: sbaij->nonew = -1; /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/
2469: sbaij->free_a = PETSC_FALSE;
2470: sbaij->free_ij = PETSC_FALSE;
2472: for (ii=0; ii<m; ii++) {
2473: sbaij->ilen[ii] = sbaij->imax[ii] = i[ii+1] - i[ii];
2474: #if defined(PETSC_USE_DEBUG)
2475: 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]);
2476: #endif
2477: }
2478: #if defined(PETSC_USE_DEBUG)
2479: for (ii=0; ii<sbaij->i[m]; ii++) {
2480: if (j[ii] < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column index at location = %d index = %d",ii,j[ii]);
2481: 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]);
2482: }
2483: #endif
2485: MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
2486: MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
2487: return(0);
2488: }