Actual source code: sbaij.c
petsc-3.5.2 2014-09-08
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,*ii = a->i,i;
28: MatMarkDiagonal_SeqSBAIJ(A);
29: *missing = PETSC_FALSE;
30: if (A->rmap->n > 0 && !ii) {
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 (diag[i] >= ii[i+1]) {
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,j;
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++) {
62: a->diag[i] = a->i[i+1];
63: for (j=a->i[i]; j<a->i[i+1]; j++) {
64: if (a->j[j] == i) {
65: a->diag[i] = j;
66: break;
67: }
68: }
69: }
70: return(0);
71: }
75: static PetscErrorCode MatGetRowIJ_SeqSBAIJ(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *nn,const PetscInt *inia[],const PetscInt *inja[],PetscBool *done)
76: {
77: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
78: PetscInt i,j,n = a->mbs,nz = a->i[n],bs = A->rmap->bs;
79: PetscInt **ia = (PetscInt**)inia,**ja = (PetscInt**)inja;
83: *nn = n;
84: if (!ia) return(0);
85: if (!blockcompressed) {
86: /* malloc & create the natural set of indices */
87: PetscMalloc2((n+1)*bs,ia,nz*bs,ja);
88: for (i=0; i<n+1; i++) {
89: for (j=0; j<bs; j++) {
90: *ia[i*bs+j] = a->i[i]*bs+j+oshift;
91: }
92: }
93: for (i=0; i<nz; i++) {
94: for (j=0; j<bs; j++) {
95: *ja[i*bs+j] = a->j[i]*bs+j+oshift;
96: }
97: }
98: } else { /* blockcompressed */
99: if (oshift == 1) {
100: /* temporarily add 1 to i and j indices */
101: for (i=0; i<nz; i++) a->j[i]++;
102: for (i=0; i<n+1; i++) a->i[i]++;
103: }
104: *ia = a->i; *ja = a->j;
105: }
106: return(0);
107: }
111: static PetscErrorCode MatRestoreRowIJ_SeqSBAIJ(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *nn,const PetscInt *ia[],const PetscInt *ja[],PetscBool *done)
112: {
113: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
114: PetscInt i,n = a->mbs,nz = a->i[n];
118: if (!ia) return(0);
120: if (!blockcompressed) {
121: PetscFree2(*ia,*ja);
122: } else if (oshift == 1) { /* blockcompressed */
123: for (i=0; i<nz; i++) a->j[i]--;
124: for (i=0; i<n+1; i++) a->i[i]--;
125: }
126: return(0);
127: }
131: PetscErrorCode MatDestroy_SeqSBAIJ(Mat A)
132: {
133: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
137: #if defined(PETSC_USE_LOG)
138: PetscLogObjectState((PetscObject)A,"Rows=%D, NZ=%D",A->rmap->N,a->nz);
139: #endif
140: MatSeqXAIJFreeAIJ(A,&a->a,&a->j,&a->i);
141: if (a->free_diag) {PetscFree(a->diag);}
142: ISDestroy(&a->row);
143: ISDestroy(&a->col);
144: ISDestroy(&a->icol);
145: PetscFree(a->idiag);
146: PetscFree(a->inode.size);
147: if (a->free_imax_ilen) {PetscFree2(a->imax,a->ilen);}
148: PetscFree(a->solve_work);
149: PetscFree(a->sor_work);
150: PetscFree(a->solves_work);
151: PetscFree(a->mult_work);
152: PetscFree(a->saved_values);
153: PetscFree(a->xtoy);
154: if (a->free_jshort) {PetscFree(a->jshort);}
155: PetscFree(a->inew);
156: MatDestroy(&a->parent);
157: PetscFree(A->data);
159: PetscObjectChangeTypeName((PetscObject)A,0);
160: PetscObjectComposeFunction((PetscObject)A,"MatStoreValues_C",NULL);
161: PetscObjectComposeFunction((PetscObject)A,"MatRetrieveValues_C",NULL);
162: PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetColumnIndices_C",NULL);
163: PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqaij_C",NULL);
164: PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqbaij_C",NULL);
165: PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetPreallocation_C",NULL);
166: PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetPreallocationCSR_C",NULL);
167: PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqsbstrm_C",NULL);
168: return(0);
169: }
173: PetscErrorCode MatSetOption_SeqSBAIJ(Mat A,MatOption op,PetscBool flg)
174: {
175: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
179: switch (op) {
180: case MAT_ROW_ORIENTED:
181: a->roworiented = flg;
182: break;
183: case MAT_KEEP_NONZERO_PATTERN:
184: a->keepnonzeropattern = flg;
185: break;
186: case MAT_NEW_NONZERO_LOCATIONS:
187: a->nonew = (flg ? 0 : 1);
188: break;
189: case MAT_NEW_NONZERO_LOCATION_ERR:
190: a->nonew = (flg ? -1 : 0);
191: break;
192: case MAT_NEW_NONZERO_ALLOCATION_ERR:
193: a->nonew = (flg ? -2 : 0);
194: break;
195: case MAT_UNUSED_NONZERO_LOCATION_ERR:
196: a->nounused = (flg ? -1 : 0);
197: break;
198: case MAT_NEW_DIAGONALS:
199: case MAT_IGNORE_OFF_PROC_ENTRIES:
200: case MAT_USE_HASH_TABLE:
201: PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
202: break;
203: case MAT_HERMITIAN:
204: if (!A->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatAssemblyEnd() first");
205: if (A->cmap->n < 65536 && A->cmap->bs == 1) {
206: A->ops->mult = MatMult_SeqSBAIJ_1_Hermitian_ushort;
207: } else if (A->cmap->bs == 1) {
208: A->ops->mult = MatMult_SeqSBAIJ_1_Hermitian;
209: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for Hermitian with block size greater than 1");
210: break;
211: case MAT_SPD:
212: /* These options are handled directly by MatSetOption() */
213: break;
214: case MAT_SYMMETRIC:
215: case MAT_STRUCTURALLY_SYMMETRIC:
216: case MAT_SYMMETRY_ETERNAL:
217: if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix must be symmetric");
218: PetscInfo1(A,"Option %s not relevent\n",MatOptions[op]);
219: break;
220: case MAT_IGNORE_LOWER_TRIANGULAR:
221: a->ignore_ltriangular = flg;
222: break;
223: case MAT_ERROR_LOWER_TRIANGULAR:
224: a->ignore_ltriangular = flg;
225: break;
226: case MAT_GETROW_UPPERTRIANGULAR:
227: a->getrow_utriangular = flg;
228: break;
229: default:
230: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %d",op);
231: }
232: return(0);
233: }
237: PetscErrorCode MatGetRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
238: {
239: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
243: 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()");
245: /* Get the upper triangular part of the row */
246: MatGetRow_SeqBAIJ_private(A,row,nz,idx,v,a->i,a->j,a->a);
247: return(0);
248: }
252: PetscErrorCode MatRestoreRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
253: {
257: if (idx) {PetscFree(*idx);}
258: if (v) {PetscFree(*v);}
259: return(0);
260: }
264: PetscErrorCode MatGetRowUpperTriangular_SeqSBAIJ(Mat A)
265: {
266: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
269: a->getrow_utriangular = PETSC_TRUE;
270: return(0);
271: }
274: PetscErrorCode MatRestoreRowUpperTriangular_SeqSBAIJ(Mat A)
275: {
276: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
279: a->getrow_utriangular = PETSC_FALSE;
280: return(0);
281: }
285: PetscErrorCode MatTranspose_SeqSBAIJ(Mat A,MatReuse reuse,Mat *B)
286: {
290: if (reuse == MAT_INITIAL_MATRIX || *B != A) {
291: MatDuplicate(A,MAT_COPY_VALUES,B);
292: }
293: return(0);
294: }
298: PetscErrorCode MatView_SeqSBAIJ_ASCII(Mat A,PetscViewer viewer)
299: {
300: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
301: PetscErrorCode ierr;
302: PetscInt i,j,bs = A->rmap->bs,k,l,bs2=a->bs2;
303: PetscViewerFormat format;
304: PetscInt *diag;
307: PetscViewerGetFormat(viewer,&format);
308: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
309: PetscViewerASCIIPrintf(viewer," block size is %D\n",bs);
310: } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
311: Mat aij;
312: if (A->factortype && bs>1) {
313: PetscPrintf(PETSC_COMM_SELF,"Warning: matrix is factored with bs>1. MatView() with PETSC_VIEWER_ASCII_MATLAB is not supported and ignored!\n");
314: return(0);
315: }
316: MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&aij);
317: MatView(aij,viewer);
318: MatDestroy(&aij);
319: } else if (format == PETSC_VIEWER_ASCII_COMMON) {
320: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
321: for (i=0; i<a->mbs; i++) {
322: for (j=0; j<bs; j++) {
323: PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);
324: for (k=a->i[i]; k<a->i[i+1]; k++) {
325: for (l=0; l<bs; l++) {
326: #if defined(PETSC_USE_COMPLEX)
327: if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
328: PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",bs*a->j[k]+l,
329: (double)PetscRealPart(a->a[bs2*k + l*bs + j]),(double)PetscImaginaryPart(a->a[bs2*k + l*bs + j]));
330: } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
331: PetscViewerASCIIPrintf(viewer," (%D, %g - %g i) ",bs*a->j[k]+l,
332: (double)PetscRealPart(a->a[bs2*k + l*bs + j]),-(double)PetscImaginaryPart(a->a[bs2*k + l*bs + j]));
333: } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
334: PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,(double)PetscRealPart(a->a[bs2*k + l*bs + j]));
335: }
336: #else
337: if (a->a[bs2*k + l*bs + j] != 0.0) {
338: PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,(double)a->a[bs2*k + l*bs + j]);
339: }
340: #endif
341: }
342: }
343: PetscViewerASCIIPrintf(viewer,"\n");
344: }
345: }
346: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
347: } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
348: return(0);
349: } else {
350: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
351: if (A->factortype) { /* for factored matrix */
352: if (bs>1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"matrix is factored with bs>1. Not implemented yet");
354: diag=a->diag;
355: for (i=0; i<a->mbs; i++) { /* for row block i */
356: PetscViewerASCIIPrintf(viewer,"row %D:",i);
357: /* diagonal entry */
358: #if defined(PETSC_USE_COMPLEX)
359: if (PetscImaginaryPart(a->a[diag[i]]) > 0.0) {
360: 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]]));
361: } else if (PetscImaginaryPart(a->a[diag[i]]) < 0.0) {
362: 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]]));
363: } else {
364: PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[diag[i]],(double)PetscRealPart(1.0/a->a[diag[i]]));
365: }
366: #else
367: PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[diag[i]],(double)(1.0/a->a[diag[i]]));
368: #endif
369: /* off-diagonal entries */
370: for (k=a->i[i]; k<a->i[i+1]-1; k++) {
371: #if defined(PETSC_USE_COMPLEX)
372: if (PetscImaginaryPart(a->a[k]) > 0.0) {
373: PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",bs*a->j[k],(double)PetscRealPart(a->a[k]),(double)PetscImaginaryPart(a->a[k]));
374: } else if (PetscImaginaryPart(a->a[k]) < 0.0) {
375: PetscViewerASCIIPrintf(viewer," (%D, %g - %g i) ",bs*a->j[k],(double)PetscRealPart(a->a[k]),-(double)PetscImaginaryPart(a->a[k]));
376: } else {
377: PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k],(double)PetscRealPart(a->a[k]));
378: }
379: #else
380: PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[k],(double)a->a[k]);
381: #endif
382: }
383: PetscViewerASCIIPrintf(viewer,"\n");
384: }
386: } else { /* for non-factored matrix */
387: for (i=0; i<a->mbs; i++) { /* for row block i */
388: for (j=0; j<bs; j++) { /* for row bs*i + j */
389: PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);
390: for (k=a->i[i]; k<a->i[i+1]; k++) { /* for column block */
391: for (l=0; l<bs; l++) { /* for column */
392: #if defined(PETSC_USE_COMPLEX)
393: if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) {
394: PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",bs*a->j[k]+l,
395: (double)PetscRealPart(a->a[bs2*k + l*bs + j]),(double)PetscImaginaryPart(a->a[bs2*k + l*bs + j]));
396: } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) {
397: PetscViewerASCIIPrintf(viewer," (%D, %g - %g i) ",bs*a->j[k]+l,
398: (double)PetscRealPart(a->a[bs2*k + l*bs + j]),-(double)PetscImaginaryPart(a->a[bs2*k + l*bs + j]));
399: } else {
400: PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,(double)PetscRealPart(a->a[bs2*k + l*bs + j]));
401: }
402: #else
403: PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,(double)a->a[bs2*k + l*bs + j]);
404: #endif
405: }
406: }
407: PetscViewerASCIIPrintf(viewer,"\n");
408: }
409: }
410: }
411: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
412: }
413: PetscViewerFlush(viewer);
414: return(0);
415: }
417: #include <petscdraw.h>
420: static PetscErrorCode MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw,void *Aa)
421: {
422: Mat A = (Mat) Aa;
423: Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data;
425: PetscInt row,i,j,k,l,mbs=a->mbs,color,bs=A->rmap->bs,bs2=a->bs2;
426: PetscMPIInt rank;
427: PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r;
428: MatScalar *aa;
429: MPI_Comm comm;
430: PetscViewer viewer;
433: /*
434: This is nasty. If this is called from an originally parallel matrix
435: then all processes call this,but only the first has the matrix so the
436: rest should return immediately.
437: */
438: PetscObjectGetComm((PetscObject)draw,&comm);
439: MPI_Comm_rank(comm,&rank);
440: if (rank) return(0);
442: PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);
444: PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);
445: PetscDrawString(draw, .3*(xl+xr), .3*(yl+yr), PETSC_DRAW_BLACK, "symmetric");
447: /* loop over matrix elements drawing boxes */
448: color = PETSC_DRAW_BLUE;
449: for (i=0,row=0; i<mbs; i++,row+=bs) {
450: for (j=a->i[i]; j<a->i[i+1]; j++) {
451: y_l = A->rmap->N - row - 1.0; y_r = y_l + 1.0;
452: x_l = a->j[j]*bs; x_r = x_l + 1.0;
453: aa = a->a + j*bs2;
454: for (k=0; k<bs; k++) {
455: for (l=0; l<bs; l++) {
456: if (PetscRealPart(*aa++) >= 0.) continue;
457: PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
458: }
459: }
460: }
461: }
462: color = PETSC_DRAW_CYAN;
463: for (i=0,row=0; i<mbs; i++,row+=bs) {
464: for (j=a->i[i]; j<a->i[i+1]; j++) {
465: y_l = A->rmap->N - row - 1.0; y_r = y_l + 1.0;
466: x_l = a->j[j]*bs; x_r = x_l + 1.0;
467: aa = a->a + j*bs2;
468: for (k=0; k<bs; k++) {
469: for (l=0; l<bs; l++) {
470: if (PetscRealPart(*aa++) != 0.) continue;
471: PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
472: }
473: }
474: }
475: }
477: color = PETSC_DRAW_RED;
478: for (i=0,row=0; i<mbs; i++,row+=bs) {
479: for (j=a->i[i]; j<a->i[i+1]; j++) {
480: y_l = A->rmap->N - row - 1.0; y_r = y_l + 1.0;
481: x_l = a->j[j]*bs; x_r = x_l + 1.0;
482: aa = a->a + j*bs2;
483: for (k=0; k<bs; k++) {
484: for (l=0; l<bs; l++) {
485: if (PetscRealPart(*aa++) <= 0.) continue;
486: PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
487: }
488: }
489: }
490: }
491: return(0);
492: }
496: static PetscErrorCode MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer)
497: {
499: PetscReal xl,yl,xr,yr,w,h;
500: PetscDraw draw;
501: PetscBool isnull;
504: PetscViewerDrawGetDraw(viewer,0,&draw);
505: PetscDrawIsNull(draw,&isnull); if (isnull) return(0);
507: PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);
508: xr = A->rmap->N; yr = A->rmap->N; h = yr/10.0; w = xr/10.0;
509: xr += w; yr += h; xl = -w; yl = -h;
510: PetscDrawSetCoordinates(draw,xl,yl,xr,yr);
511: PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);
512: PetscObjectCompose((PetscObject)A,"Zoomviewer",NULL);
513: return(0);
514: }
518: PetscErrorCode MatView_SeqSBAIJ(Mat A,PetscViewer viewer)
519: {
521: PetscBool iascii,isdraw;
522: FILE *file = 0;
525: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
526: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
527: if (iascii) {
528: MatView_SeqSBAIJ_ASCII(A,viewer);
529: } else if (isdraw) {
530: MatView_SeqSBAIJ_Draw(A,viewer);
531: } else {
532: Mat B;
533: MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&B);
534: MatView(B,viewer);
535: MatDestroy(&B);
536: PetscViewerBinaryGetInfoPointer(viewer,&file);
537: if (file) {
538: fprintf(file,"-matload_block_size %d\n",(int)A->rmap->bs);
539: }
540: }
541: return(0);
542: }
547: PetscErrorCode MatGetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],PetscScalar v[])
548: {
549: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
550: PetscInt *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
551: PetscInt *ai = a->i,*ailen = a->ilen;
552: PetscInt brow,bcol,ridx,cidx,bs=A->rmap->bs,bs2=a->bs2;
553: MatScalar *ap,*aa = a->a;
556: for (k=0; k<m; k++) { /* loop over rows */
557: row = im[k]; brow = row/bs;
558: if (row < 0) {v += n; continue;} /* SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",row); */
559: 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);
560: rp = aj + ai[brow]; ap = aa + bs2*ai[brow];
561: nrow = ailen[brow];
562: for (l=0; l<n; l++) { /* loop over columns */
563: if (in[l] < 0) {v++; continue;} /* SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %D",in[l]); */
564: 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);
565: col = in[l];
566: bcol = col/bs;
567: cidx = col%bs;
568: ridx = row%bs;
569: high = nrow;
570: low = 0; /* assume unsorted */
571: while (high-low > 5) {
572: t = (low+high)/2;
573: if (rp[t] > bcol) high = t;
574: else low = t;
575: }
576: for (i=low; i<high; i++) {
577: if (rp[i] > bcol) break;
578: if (rp[i] == bcol) {
579: *v++ = ap[bs2*i+bs*cidx+ridx];
580: goto finished;
581: }
582: }
583: *v++ = 0.0;
584: finished:;
585: }
586: }
587: return(0);
588: }
593: PetscErrorCode MatSetValuesBlocked_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
594: {
595: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
596: PetscErrorCode ierr;
597: PetscInt *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,lastcol = -1;
598: PetscInt *imax =a->imax,*ai=a->i,*ailen=a->ilen;
599: PetscInt *aj =a->j,nonew=a->nonew,bs2=a->bs2,bs=A->rmap->bs,stepval;
600: PetscBool roworiented=a->roworiented;
601: const PetscScalar *value = v;
602: MatScalar *ap,*aa = a->a,*bap;
605: if (roworiented) stepval = (n-1)*bs;
606: else stepval = (m-1)*bs;
608: for (k=0; k<m; k++) { /* loop over added rows */
609: row = im[k];
610: if (row < 0) continue;
611: #if defined(PETSC_USE_DEBUG)
612: if (row >= a->mbs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,a->mbs-1);
613: #endif
614: rp = aj + ai[row];
615: ap = aa + bs2*ai[row];
616: rmax = imax[row];
617: nrow = ailen[row];
618: low = 0;
619: high = nrow;
620: for (l=0; l<n; l++) { /* loop over added columns */
621: if (in[l] < 0) continue;
622: col = in[l];
623: #if defined(PETSC_USE_DEBUG)
624: if (col >= a->nbs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",col,a->nbs-1);
625: #endif
626: if (col < row) {
627: if (a->ignore_ltriangular) continue; /* ignore lower triangular block */
628: 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)");
629: }
630: if (roworiented) value = v + k*(stepval+bs)*bs + l*bs;
631: else value = v + l*(stepval+bs)*bs + k*bs;
633: if (col <= lastcol) low = 0;
634: else high = nrow;
636: lastcol = col;
637: while (high-low > 7) {
638: t = (low+high)/2;
639: if (rp[t] > col) high = t;
640: else low = t;
641: }
642: for (i=low; i<high; i++) {
643: if (rp[i] > col) break;
644: if (rp[i] == col) {
645: bap = ap + bs2*i;
646: if (roworiented) {
647: if (is == ADD_VALUES) {
648: for (ii=0; ii<bs; ii++,value+=stepval) {
649: for (jj=ii; jj<bs2; jj+=bs) {
650: bap[jj] += *value++;
651: }
652: }
653: } else {
654: for (ii=0; ii<bs; ii++,value+=stepval) {
655: for (jj=ii; jj<bs2; jj+=bs) {
656: bap[jj] = *value++;
657: }
658: }
659: }
660: } else {
661: if (is == ADD_VALUES) {
662: for (ii=0; ii<bs; ii++,value+=stepval) {
663: for (jj=0; jj<bs; jj++) {
664: *bap++ += *value++;
665: }
666: }
667: } else {
668: for (ii=0; ii<bs; ii++,value+=stepval) {
669: for (jj=0; jj<bs; jj++) {
670: *bap++ = *value++;
671: }
672: }
673: }
674: }
675: goto noinsert2;
676: }
677: }
678: if (nonew == 1) goto noinsert2;
679: if (nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col);
680: MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar);
681: N = nrow++ - 1; high++;
682: /* shift up all the later entries in this row */
683: for (ii=N; ii>=i; ii--) {
684: rp[ii+1] = rp[ii];
685: PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));
686: }
687: if (N >= i) {
688: PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));
689: }
690: rp[i] = col;
691: bap = ap + bs2*i;
692: if (roworiented) {
693: for (ii=0; ii<bs; ii++,value+=stepval) {
694: for (jj=ii; jj<bs2; jj+=bs) {
695: bap[jj] = *value++;
696: }
697: }
698: } else {
699: for (ii=0; ii<bs; ii++,value+=stepval) {
700: for (jj=0; jj<bs; jj++) {
701: *bap++ = *value++;
702: }
703: }
704: }
705: noinsert2:;
706: low = i;
707: }
708: ailen[row] = nrow;
709: }
710: return(0);
711: }
713: /*
714: This is not yet used
715: */
718: PetscErrorCode MatAssemblyEnd_SeqSBAIJ_SeqAIJ_Inode(Mat A)
719: {
720: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
722: const PetscInt *ai = a->i, *aj = a->j,*cols;
723: PetscInt i = 0,j,blk_size,m = A->rmap->n,node_count = 0,nzx,nzy,*ns,row,nz,cnt,cnt2,*counts;
724: PetscBool flag;
727: PetscMalloc1(m,&ns);
728: while (i < m) {
729: nzx = ai[i+1] - ai[i]; /* Number of nonzeros */
730: /* Limits the number of elements in a node to 'a->inode.limit' */
731: for (j=i+1,blk_size=1; j<m && blk_size <a->inode.limit; ++j,++blk_size) {
732: nzy = ai[j+1] - ai[j];
733: if (nzy != (nzx - j + i)) break;
734: PetscMemcmp(aj + ai[i] + j - i,aj + ai[j],nzy*sizeof(PetscInt),&flag);
735: if (!flag) break;
736: }
737: ns[node_count++] = blk_size;
739: i = j;
740: }
741: if (!a->inode.size && m && node_count > .9*m) {
742: PetscFree(ns);
743: PetscInfo2(A,"Found %D nodes out of %D rows. Not using Inode routines\n",node_count,m);
744: } else {
745: a->inode.node_count = node_count;
747: PetscMalloc1(node_count,&a->inode.size);
748: PetscLogObjectMemory((PetscObject)A,node_count*sizeof(PetscInt));
749: PetscMemcpy(a->inode.size,ns,node_count*sizeof(PetscInt));
750: PetscFree(ns);
751: PetscInfo3(A,"Found %D nodes of %D. Limit used: %D. Using Inode routines\n",node_count,m,a->inode.limit);
753: /* count collections of adjacent columns in each inode */
754: row = 0;
755: cnt = 0;
756: for (i=0; i<node_count; i++) {
757: cols = aj + ai[row] + a->inode.size[i];
758: nz = ai[row+1] - ai[row] - a->inode.size[i];
759: for (j=1; j<nz; j++) {
760: if (cols[j] != cols[j-1]+1) cnt++;
761: }
762: cnt++;
763: row += a->inode.size[i];
764: }
765: PetscMalloc1(2*cnt,&counts);
766: cnt = 0;
767: row = 0;
768: for (i=0; i<node_count; i++) {
769: cols = aj + ai[row] + a->inode.size[i];
770: counts[2*cnt] = cols[0];
771: nz = ai[row+1] - ai[row] - a->inode.size[i];
772: cnt2 = 1;
773: for (j=1; j<nz; j++) {
774: if (cols[j] != cols[j-1]+1) {
775: counts[2*(cnt++)+1] = cnt2;
776: counts[2*cnt] = cols[j];
777: cnt2 = 1;
778: } else cnt2++;
779: }
780: counts[2*(cnt++)+1] = cnt2;
781: row += a->inode.size[i];
782: }
783: PetscIntView(2*cnt,counts,0);
784: }
785: return(0);
786: }
790: PetscErrorCode MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode)
791: {
792: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
794: PetscInt fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
795: PetscInt m = A->rmap->N,*ip,N,*ailen = a->ilen;
796: PetscInt mbs = a->mbs,bs2 = a->bs2,rmax = 0;
797: MatScalar *aa = a->a,*ap;
800: if (mode == MAT_FLUSH_ASSEMBLY) return(0);
802: if (m) rmax = ailen[0];
803: for (i=1; i<mbs; i++) {
804: /* move each row back by the amount of empty slots (fshift) before it*/
805: fshift += imax[i-1] - ailen[i-1];
806: rmax = PetscMax(rmax,ailen[i]);
807: if (fshift) {
808: ip = aj + ai[i]; ap = aa + bs2*ai[i];
809: N = ailen[i];
810: for (j=0; j<N; j++) {
811: ip[j-fshift] = ip[j];
812: PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));
813: }
814: }
815: ai[i] = ai[i-1] + ailen[i-1];
816: }
817: if (mbs) {
818: fshift += imax[mbs-1] - ailen[mbs-1];
819: ai[mbs] = ai[mbs-1] + ailen[mbs-1];
820: }
821: /* reset ilen and imax for each row */
822: for (i=0; i<mbs; i++) {
823: ailen[i] = imax[i] = ai[i+1] - ai[i];
824: }
825: a->nz = ai[mbs];
827: /* diagonals may have moved, reset it */
828: if (a->diag) {
829: PetscMemcpy(a->diag,ai,mbs*sizeof(PetscInt));
830: }
831: 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);
833: 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);
834: PetscInfo1(A,"Number of mallocs during MatSetValues is %D\n",a->reallocs);
835: PetscInfo1(A,"Most nonzeros blocks in any row is %D\n",rmax);
837: A->info.mallocs += a->reallocs;
838: a->reallocs = 0;
839: A->info.nz_unneeded = (PetscReal)fshift*bs2;
840: a->idiagvalid = PETSC_FALSE;
842: if (A->cmap->n < 65536 && A->cmap->bs == 1) {
843: if (a->jshort && a->free_jshort) {
844: /* when matrix data structure is changed, previous jshort must be replaced */
845: PetscFree(a->jshort);
846: }
847: PetscMalloc1(a->i[A->rmap->n],&a->jshort);
848: PetscLogObjectMemory((PetscObject)A,a->i[A->rmap->n]*sizeof(unsigned short));
849: for (i=0; i<a->i[A->rmap->n]; i++) a->jshort[i] = a->j[i];
850: A->ops->mult = MatMult_SeqSBAIJ_1_ushort;
851: A->ops->sor = MatSOR_SeqSBAIJ_ushort;
852: a->free_jshort = PETSC_TRUE;
853: }
854: return(0);
855: }
857: /*
858: This function returns an array of flags which indicate the locations of contiguous
859: blocks that should be zeroed. for eg: if bs = 3 and is = [0,1,2,3,5,6,7,8,9]
860: then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)]
861: Assume: sizes should be long enough to hold all the values.
862: */
865: PetscErrorCode MatZeroRows_SeqSBAIJ_Check_Blocks(PetscInt idx[],PetscInt n,PetscInt bs,PetscInt sizes[], PetscInt *bs_max)
866: {
867: PetscInt i,j,k,row;
868: PetscBool flg;
871: for (i=0,j=0; i<n; j++) {
872: row = idx[i];
873: if (row%bs!=0) { /* Not the begining of a block */
874: sizes[j] = 1;
875: i++;
876: } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */
877: sizes[j] = 1; /* Also makes sure atleast 'bs' values exist for next else */
878: i++;
879: } else { /* Begining of the block, so check if the complete block exists */
880: flg = PETSC_TRUE;
881: for (k=1; k<bs; k++) {
882: if (row+k != idx[i+k]) { /* break in the block */
883: flg = PETSC_FALSE;
884: break;
885: }
886: }
887: if (flg) { /* No break in the bs */
888: sizes[j] = bs;
889: i += bs;
890: } else {
891: sizes[j] = 1;
892: i++;
893: }
894: }
895: }
896: *bs_max = j;
897: return(0);
898: }
901: /* Only add/insert a(i,j) with i<=j (blocks).
902: Any a(i,j) with i>j input by user is ingored.
903: */
907: PetscErrorCode MatSetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
908: {
909: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
911: PetscInt *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,lastcol = -1;
912: PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented;
913: PetscInt *aj =a->j,nonew=a->nonew,bs=A->rmap->bs,brow,bcol;
914: PetscInt ridx,cidx,bs2=a->bs2;
915: MatScalar *ap,value,*aa=a->a,*bap;
918: for (k=0; k<m; k++) { /* loop over added rows */
919: row = im[k]; /* row number */
920: brow = row/bs; /* block row number */
921: if (row < 0) continue;
922: #if defined(PETSC_USE_DEBUG)
923: 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);
924: #endif
925: rp = aj + ai[brow]; /*ptr to beginning of column value of the row block*/
926: ap = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/
927: rmax = imax[brow]; /* maximum space allocated for this row */
928: nrow = ailen[brow]; /* actual length of this row */
929: low = 0;
931: for (l=0; l<n; l++) { /* loop over added columns */
932: if (in[l] < 0) continue;
933: #if defined(PETSC_USE_DEBUG)
934: 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);
935: #endif
936: col = in[l];
937: bcol = col/bs; /* block col number */
939: if (brow > bcol) {
940: if (a->ignore_ltriangular) continue; /* ignore lower triangular values */
941: 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)");
942: }
944: ridx = row % bs; cidx = col % bs; /*row and col index inside the block */
945: if ((brow==bcol && ridx<=cidx) || (brow<bcol)) {
946: /* element value a(k,l) */
947: if (roworiented) value = v[l + k*n];
948: else value = v[k + l*m];
950: /* move pointer bap to a(k,l) quickly and add/insert value */
951: if (col <= lastcol) low = 0;
952: high = nrow;
953: lastcol = col;
954: while (high-low > 7) {
955: t = (low+high)/2;
956: if (rp[t] > bcol) high = t;
957: else low = t;
958: }
959: for (i=low; i<high; i++) {
960: if (rp[i] > bcol) break;
961: if (rp[i] == bcol) {
962: bap = ap + bs2*i + bs*cidx + ridx;
963: if (is == ADD_VALUES) *bap += value;
964: else *bap = value;
965: /* for diag block, add/insert its symmetric element a(cidx,ridx) */
966: if (brow == bcol && ridx < cidx) {
967: bap = ap + bs2*i + bs*ridx + cidx;
968: if (is == ADD_VALUES) *bap += value;
969: else *bap = value;
970: }
971: goto noinsert1;
972: }
973: }
975: if (nonew == 1) goto noinsert1;
976: if (nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col);
977: MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar);
979: N = nrow++ - 1; high++;
980: /* shift up all the later entries in this row */
981: for (ii=N; ii>=i; ii--) {
982: rp[ii+1] = rp[ii];
983: PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));
984: }
985: if (N>=i) {
986: PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));
987: }
988: rp[i] = bcol;
989: ap[bs2*i + bs*cidx + ridx] = value;
990: A->nonzerostate++;
991: noinsert1:;
992: low = i;
993: }
994: } /* end of loop over added columns */
995: ailen[brow] = nrow;
996: } /* end of loop over added rows */
997: return(0);
998: }
1002: PetscErrorCode MatICCFactor_SeqSBAIJ(Mat inA,IS row,const MatFactorInfo *info)
1003: {
1004: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data;
1005: Mat outA;
1007: PetscBool row_identity;
1010: if (info->levels != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only levels=0 is supported for in-place icc");
1011: ISIdentity(row,&row_identity);
1012: if (!row_identity) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix reordering is not supported");
1013: 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()! */
1015: outA = inA;
1016: inA->factortype = MAT_FACTOR_ICC;
1018: MatMarkDiagonal_SeqSBAIJ(inA);
1019: MatSeqSBAIJSetNumericFactorization_inplace(inA,row_identity);
1021: PetscObjectReference((PetscObject)row);
1022: ISDestroy(&a->row);
1023: a->row = row;
1024: PetscObjectReference((PetscObject)row);
1025: ISDestroy(&a->col);
1026: a->col = row;
1028: /* Create the invert permutation so that it can be used in MatCholeskyFactorNumeric() */
1029: if (a->icol) {ISInvertPermutation(row,PETSC_DECIDE, &a->icol);}
1030: PetscLogObjectParent((PetscObject)inA,(PetscObject)a->icol);
1032: if (!a->solve_work) {
1033: PetscMalloc1((inA->rmap->N+inA->rmap->bs),&a->solve_work);
1034: PetscLogObjectMemory((PetscObject)inA,(inA->rmap->N+inA->rmap->bs)*sizeof(PetscScalar));
1035: }
1037: MatCholeskyFactorNumeric(outA,inA,info);
1038: return(0);
1039: }
1043: PetscErrorCode MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,PetscInt *indices)
1044: {
1045: Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ*)mat->data;
1046: PetscInt i,nz,n;
1050: nz = baij->maxnz;
1051: n = mat->cmap->n;
1052: for (i=0; i<nz; i++) baij->j[i] = indices[i];
1054: baij->nz = nz;
1055: for (i=0; i<n; i++) baij->ilen[i] = baij->imax[i];
1057: MatSetOption(mat,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
1058: return(0);
1059: }
1063: /*@
1064: MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows
1065: in the matrix.
1067: Input Parameters:
1068: + mat - the SeqSBAIJ matrix
1069: - indices - the column indices
1071: Level: advanced
1073: Notes:
1074: This can be called if you have precomputed the nonzero structure of the
1075: matrix and want to provide it to the matrix object to improve the performance
1076: of the MatSetValues() operation.
1078: You MUST have set the correct numbers of nonzeros per row in the call to
1079: MatCreateSeqSBAIJ(), and the columns indices MUST be sorted.
1081: MUST be called before any calls to MatSetValues()
1083: .seealso: MatCreateSeqSBAIJ
1084: @*/
1085: PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat mat,PetscInt *indices)
1086: {
1092: PetscUseMethod(mat,"MatSeqSBAIJSetColumnIndices_C",(Mat,PetscInt*),(mat,indices));
1093: return(0);
1094: }
1098: PetscErrorCode MatCopy_SeqSBAIJ(Mat A,Mat B,MatStructure str)
1099: {
1103: /* If the two matrices have the same copy implementation, use fast copy. */
1104: if (str == SAME_NONZERO_PATTERN && (A->ops->copy == B->ops->copy)) {
1105: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1106: Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data;
1108: 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");
1109: PetscMemcpy(b->a,a->a,(a->i[A->rmap->N])*sizeof(PetscScalar));
1110: } else {
1111: MatGetRowUpperTriangular(A);
1112: MatCopy_Basic(A,B,str);
1113: MatRestoreRowUpperTriangular(A);
1114: }
1115: return(0);
1116: }
1120: PetscErrorCode MatSetUp_SeqSBAIJ(Mat A)
1121: {
1125: MatSeqSBAIJSetPreallocation_SeqSBAIJ(A,A->rmap->bs,PETSC_DEFAULT,0);
1126: return(0);
1127: }
1131: PetscErrorCode MatSeqSBAIJGetArray_SeqSBAIJ(Mat A,PetscScalar *array[])
1132: {
1133: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1136: *array = a->a;
1137: return(0);
1138: }
1142: PetscErrorCode MatSeqSBAIJRestoreArray_SeqSBAIJ(Mat A,PetscScalar *array[])
1143: {
1145: return(0);
1146: }
1150: PetscErrorCode MatAXPYGetPreallocation_SeqSBAIJ(Mat Y,Mat X,PetscInt *nnz)
1151: {
1152: PetscInt bs = Y->rmap->bs,mbs = Y->rmap->N/bs;
1153: Mat_SeqSBAIJ *x = (Mat_SeqSBAIJ*)X->data;
1154: Mat_SeqSBAIJ *y = (Mat_SeqSBAIJ*)Y->data;
1158: /* Set the number of nonzeros in the new matrix */
1159: MatAXPYGetPreallocation_SeqX_private(mbs,x->i,x->j,y->i,y->j,nnz);
1160: return(0);
1161: }
1165: PetscErrorCode MatAXPY_SeqSBAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str)
1166: {
1167: Mat_SeqSBAIJ *x=(Mat_SeqSBAIJ*)X->data, *y=(Mat_SeqSBAIJ*)Y->data;
1169: PetscInt i,bs=Y->rmap->bs,bs2=bs*bs,j;
1170: PetscBLASInt one = 1;
1173: if (str == SAME_NONZERO_PATTERN) {
1174: PetscScalar alpha = a;
1175: PetscBLASInt bnz;
1176: PetscBLASIntCast(x->nz*bs2,&bnz);
1177: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&bnz,&alpha,x->a,&one,y->a,&one));
1178: PetscObjectStateIncrease((PetscObject)Y);
1179: } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
1180: if (y->xtoy && y->XtoY != X) {
1181: PetscFree(y->xtoy);
1182: MatDestroy(&y->XtoY);
1183: }
1184: if (!y->xtoy) { /* get xtoy */
1185: MatAXPYGetxtoy_Private(x->mbs,x->i,x->j,NULL, y->i,y->j,NULL, &y->xtoy);
1186: y->XtoY = X;
1187: }
1188: for (i=0; i<x->nz; i++) {
1189: j = 0;
1190: while (j < bs2) {
1191: y->a[bs2*y->xtoy[i]+j] += a*(x->a[bs2*i+j]);
1192: j++;
1193: }
1194: }
1195: PetscObjectStateIncrease((PetscObject)Y);
1196: 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)));
1197: } else {
1198: Mat B;
1199: PetscInt *nnz;
1200: if (bs != X->rmap->bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrices must have same block size");
1201: MatGetRowUpperTriangular(X);
1202: MatGetRowUpperTriangular(Y);
1203: PetscMalloc1(Y->rmap->N,&nnz);
1204: MatCreate(PetscObjectComm((PetscObject)Y),&B);
1205: PetscObjectSetName((PetscObject)B,((PetscObject)Y)->name);
1206: MatSetSizes(B,Y->rmap->n,Y->cmap->n,Y->rmap->N,Y->cmap->N);
1207: MatSetBlockSizesFromMats(B,Y,Y);
1208: MatSetType(B,(MatType) ((PetscObject)Y)->type_name);
1209: MatAXPYGetPreallocation_SeqSBAIJ(Y,X,nnz);
1210: MatSeqSBAIJSetPreallocation(B,bs,0,nnz);
1211:
1212: MatAXPY_BasicWithPreallocation(B,Y,a,X,str);
1213:
1214: MatHeaderReplace(Y,B);
1215: PetscFree(nnz);
1216: MatRestoreRowUpperTriangular(X);
1217: MatRestoreRowUpperTriangular(Y);
1218: }
1219: return(0);
1220: }
1224: PetscErrorCode MatIsSymmetric_SeqSBAIJ(Mat A,PetscReal tol,PetscBool *flg)
1225: {
1227: *flg = PETSC_TRUE;
1228: return(0);
1229: }
1233: PetscErrorCode MatIsStructurallySymmetric_SeqSBAIJ(Mat A,PetscBool *flg)
1234: {
1236: *flg = PETSC_TRUE;
1237: return(0);
1238: }
1242: PetscErrorCode MatIsHermitian_SeqSBAIJ(Mat A,PetscReal tol,PetscBool *flg)
1243: {
1245: *flg = PETSC_FALSE;
1246: return(0);
1247: }
1251: PetscErrorCode MatRealPart_SeqSBAIJ(Mat A)
1252: {
1253: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1254: PetscInt i,nz = a->bs2*a->i[a->mbs];
1255: MatScalar *aa = a->a;
1258: for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]);
1259: return(0);
1260: }
1264: PetscErrorCode MatImaginaryPart_SeqSBAIJ(Mat A)
1265: {
1266: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1267: PetscInt i,nz = a->bs2*a->i[a->mbs];
1268: MatScalar *aa = a->a;
1271: for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]);
1272: return(0);
1273: }
1277: PetscErrorCode MatZeroRowsColumns_SeqSBAIJ(Mat A,PetscInt is_n,const PetscInt is_idx[],PetscScalar diag,Vec x, Vec b)
1278: {
1279: Mat_SeqSBAIJ *baij=(Mat_SeqSBAIJ*)A->data;
1280: PetscErrorCode ierr;
1281: PetscInt i,j,k,count;
1282: PetscInt bs =A->rmap->bs,bs2=baij->bs2,row,col;
1283: PetscScalar zero = 0.0;
1284: MatScalar *aa;
1285: const PetscScalar *xx;
1286: PetscScalar *bb;
1287: PetscBool *zeroed,vecs = PETSC_FALSE;
1290: /* fix right hand side if needed */
1291: if (x && b) {
1292: VecGetArrayRead(x,&xx);
1293: VecGetArray(b,&bb);
1294: vecs = PETSC_TRUE;
1295: }
1297: /* zero the columns */
1298: PetscCalloc1(A->rmap->n,&zeroed);
1299: for (i=0; i<is_n; i++) {
1300: 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]);
1301: zeroed[is_idx[i]] = PETSC_TRUE;
1302: }
1303: if (vecs) {
1304: for (i=0; i<A->rmap->N; i++) {
1305: row = i/bs;
1306: for (j=baij->i[row]; j<baij->i[row+1]; j++) {
1307: for (k=0; k<bs; k++) {
1308: col = bs*baij->j[j] + k;
1309: if (col <= i) continue;
1310: aa = ((MatScalar*)(baij->a)) + j*bs2 + (i%bs) + bs*k;
1311: if (!zeroed[i] && zeroed[col]) bb[i] -= aa[0]*xx[col];
1312: if (zeroed[i] && !zeroed[col]) bb[col] -= aa[0]*xx[i];
1313: }
1314: }
1315: }
1316: for (i=0; i<is_n; i++) bb[is_idx[i]] = diag*xx[is_idx[i]];
1317: }
1319: for (i=0; i<A->rmap->N; i++) {
1320: if (!zeroed[i]) {
1321: row = i/bs;
1322: for (j=baij->i[row]; j<baij->i[row+1]; j++) {
1323: for (k=0; k<bs; k++) {
1324: col = bs*baij->j[j] + k;
1325: if (zeroed[col]) {
1326: aa = ((MatScalar*)(baij->a)) + j*bs2 + (i%bs) + bs*k;
1327: aa[0] = 0.0;
1328: }
1329: }
1330: }
1331: }
1332: }
1333: PetscFree(zeroed);
1334: if (vecs) {
1335: VecRestoreArrayRead(x,&xx);
1336: VecRestoreArray(b,&bb);
1337: }
1339: /* zero the rows */
1340: for (i=0; i<is_n; i++) {
1341: row = is_idx[i];
1342: count = (baij->i[row/bs +1] - baij->i[row/bs])*bs;
1343: aa = ((MatScalar*)(baij->a)) + baij->i[row/bs]*bs2 + (row%bs);
1344: for (k=0; k<count; k++) {
1345: aa[0] = zero;
1346: aa += bs;
1347: }
1348: if (diag != 0.0) {
1349: (*A->ops->setvalues)(A,1,&row,1,&row,&diag,INSERT_VALUES);
1350: }
1351: }
1352: MatAssemblyEnd_SeqSBAIJ(A,MAT_FINAL_ASSEMBLY);
1353: return(0);
1354: }
1356: /* -------------------------------------------------------------------*/
1357: static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ,
1358: MatGetRow_SeqSBAIJ,
1359: MatRestoreRow_SeqSBAIJ,
1360: MatMult_SeqSBAIJ_N,
1361: /* 4*/ MatMultAdd_SeqSBAIJ_N,
1362: MatMult_SeqSBAIJ_N, /* transpose versions are same as non-transpose versions */
1363: MatMultAdd_SeqSBAIJ_N,
1364: 0,
1365: 0,
1366: 0,
1367: /* 10*/ 0,
1368: 0,
1369: MatCholeskyFactor_SeqSBAIJ,
1370: MatSOR_SeqSBAIJ,
1371: MatTranspose_SeqSBAIJ,
1372: /* 15*/ MatGetInfo_SeqSBAIJ,
1373: MatEqual_SeqSBAIJ,
1374: MatGetDiagonal_SeqSBAIJ,
1375: MatDiagonalScale_SeqSBAIJ,
1376: MatNorm_SeqSBAIJ,
1377: /* 20*/ 0,
1378: MatAssemblyEnd_SeqSBAIJ,
1379: MatSetOption_SeqSBAIJ,
1380: MatZeroEntries_SeqSBAIJ,
1381: /* 24*/ 0,
1382: 0,
1383: 0,
1384: 0,
1385: 0,
1386: /* 29*/ MatSetUp_SeqSBAIJ,
1387: 0,
1388: 0,
1389: 0,
1390: 0,
1391: /* 34*/ MatDuplicate_SeqSBAIJ,
1392: 0,
1393: 0,
1394: 0,
1395: MatICCFactor_SeqSBAIJ,
1396: /* 39*/ MatAXPY_SeqSBAIJ,
1397: MatGetSubMatrices_SeqSBAIJ,
1398: MatIncreaseOverlap_SeqSBAIJ,
1399: MatGetValues_SeqSBAIJ,
1400: MatCopy_SeqSBAIJ,
1401: /* 44*/ 0,
1402: MatScale_SeqSBAIJ,
1403: 0,
1404: 0,
1405: MatZeroRowsColumns_SeqSBAIJ,
1406: /* 49*/ 0,
1407: MatGetRowIJ_SeqSBAIJ,
1408: MatRestoreRowIJ_SeqSBAIJ,
1409: 0,
1410: 0,
1411: /* 54*/ 0,
1412: 0,
1413: 0,
1414: 0,
1415: MatSetValuesBlocked_SeqSBAIJ,
1416: /* 59*/ MatGetSubMatrix_SeqSBAIJ,
1417: 0,
1418: 0,
1419: 0,
1420: 0,
1421: /* 64*/ 0,
1422: 0,
1423: 0,
1424: 0,
1425: 0,
1426: /* 69*/ MatGetRowMaxAbs_SeqSBAIJ,
1427: 0,
1428: 0,
1429: 0,
1430: 0,
1431: /* 74*/ 0,
1432: 0,
1433: 0,
1434: 0,
1435: 0,
1436: /* 79*/ 0,
1437: 0,
1438: 0,
1439: MatGetInertia_SeqSBAIJ,
1440: MatLoad_SeqSBAIJ,
1441: /* 84*/ MatIsSymmetric_SeqSBAIJ,
1442: MatIsHermitian_SeqSBAIJ,
1443: MatIsStructurallySymmetric_SeqSBAIJ,
1444: 0,
1445: 0,
1446: /* 89*/ 0,
1447: 0,
1448: 0,
1449: 0,
1450: 0,
1451: /* 94*/ 0,
1452: 0,
1453: 0,
1454: 0,
1455: 0,
1456: /* 99*/ 0,
1457: 0,
1458: 0,
1459: 0,
1460: 0,
1461: /*104*/ 0,
1462: MatRealPart_SeqSBAIJ,
1463: MatImaginaryPart_SeqSBAIJ,
1464: MatGetRowUpperTriangular_SeqSBAIJ,
1465: MatRestoreRowUpperTriangular_SeqSBAIJ,
1466: /*109*/ 0,
1467: 0,
1468: 0,
1469: 0,
1470: MatMissingDiagonal_SeqSBAIJ,
1471: /*114*/ 0,
1472: 0,
1473: 0,
1474: 0,
1475: 0,
1476: /*119*/ 0,
1477: 0,
1478: 0,
1479: 0,
1480: 0,
1481: /*124*/ 0,
1482: 0,
1483: 0,
1484: 0,
1485: 0,
1486: /*129*/ 0,
1487: 0,
1488: 0,
1489: 0,
1490: 0,
1491: /*134*/ 0,
1492: 0,
1493: 0,
1494: 0,
1495: 0,
1496: /*139*/ 0,
1497: 0,
1498: 0
1499: };
1503: PetscErrorCode MatStoreValues_SeqSBAIJ(Mat mat)
1504: {
1505: Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ*)mat->data;
1506: PetscInt nz = aij->i[mat->rmap->N]*mat->rmap->bs*aij->bs2;
1510: if (aij->nonew != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first");
1512: /* allocate space for values if not already there */
1513: if (!aij->saved_values) {
1514: PetscMalloc1((nz+1),&aij->saved_values);
1515: }
1517: /* copy values over */
1518: PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));
1519: return(0);
1520: }
1524: PetscErrorCode MatRetrieveValues_SeqSBAIJ(Mat mat)
1525: {
1526: Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ*)mat->data;
1528: PetscInt nz = aij->i[mat->rmap->N]*mat->rmap->bs*aij->bs2;
1531: if (aij->nonew != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first");
1532: if (!aij->saved_values) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Must call MatStoreValues(A);first");
1534: /* copy values over */
1535: PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));
1536: return(0);
1537: }
1541: PetscErrorCode MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat B,PetscInt bs,PetscInt nz,PetscInt *nnz)
1542: {
1543: Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data;
1545: PetscInt i,mbs,bs2;
1546: PetscBool skipallocation = PETSC_FALSE,flg = PETSC_FALSE,realalloc = PETSC_FALSE;
1549: if (nz >= 0 || nnz) realalloc = PETSC_TRUE;
1550: B->preallocated = PETSC_TRUE;
1552: MatSetBlockSize(B,PetscAbs(bs));
1553: PetscLayoutSetUp(B->rmap);
1554: PetscLayoutSetUp(B->cmap);
1555: PetscLayoutGetBlockSize(B->rmap,&bs);
1557: mbs = B->rmap->N/bs;
1558: bs2 = bs*bs;
1560: if (mbs*bs != B->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize");
1562: if (nz == MAT_SKIP_ALLOCATION) {
1563: skipallocation = PETSC_TRUE;
1564: nz = 0;
1565: }
1567: if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3;
1568: if (nz < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %D",nz);
1569: if (nnz) {
1570: for (i=0; i<mbs; i++) {
1571: 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]);
1572: 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);
1573: }
1574: }
1576: B->ops->mult = MatMult_SeqSBAIJ_N;
1577: B->ops->multadd = MatMultAdd_SeqSBAIJ_N;
1578: B->ops->multtranspose = MatMult_SeqSBAIJ_N;
1579: B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_N;
1581: PetscOptionsGetBool(((PetscObject)B)->prefix,"-mat_no_unroll",&flg,NULL);
1582: if (!flg) {
1583: switch (bs) {
1584: case 1:
1585: B->ops->mult = MatMult_SeqSBAIJ_1;
1586: B->ops->multadd = MatMultAdd_SeqSBAIJ_1;
1587: B->ops->multtranspose = MatMult_SeqSBAIJ_1;
1588: B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_1;
1589: break;
1590: case 2:
1591: B->ops->mult = MatMult_SeqSBAIJ_2;
1592: B->ops->multadd = MatMultAdd_SeqSBAIJ_2;
1593: B->ops->multtranspose = MatMult_SeqSBAIJ_2;
1594: B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_2;
1595: break;
1596: case 3:
1597: B->ops->mult = MatMult_SeqSBAIJ_3;
1598: B->ops->multadd = MatMultAdd_SeqSBAIJ_3;
1599: B->ops->multtranspose = MatMult_SeqSBAIJ_3;
1600: B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_3;
1601: break;
1602: case 4:
1603: B->ops->mult = MatMult_SeqSBAIJ_4;
1604: B->ops->multadd = MatMultAdd_SeqSBAIJ_4;
1605: B->ops->multtranspose = MatMult_SeqSBAIJ_4;
1606: B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_4;
1607: break;
1608: case 5:
1609: B->ops->mult = MatMult_SeqSBAIJ_5;
1610: B->ops->multadd = MatMultAdd_SeqSBAIJ_5;
1611: B->ops->multtranspose = MatMult_SeqSBAIJ_5;
1612: B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_5;
1613: break;
1614: case 6:
1615: B->ops->mult = MatMult_SeqSBAIJ_6;
1616: B->ops->multadd = MatMultAdd_SeqSBAIJ_6;
1617: B->ops->multtranspose = MatMult_SeqSBAIJ_6;
1618: B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_6;
1619: break;
1620: case 7:
1621: B->ops->mult = MatMult_SeqSBAIJ_7;
1622: B->ops->multadd = MatMultAdd_SeqSBAIJ_7;
1623: B->ops->multtranspose = MatMult_SeqSBAIJ_7;
1624: B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_7;
1625: break;
1626: }
1627: }
1629: b->mbs = mbs;
1630: b->nbs = mbs;
1631: if (!skipallocation) {
1632: if (!b->imax) {
1633: PetscMalloc2(mbs,&b->imax,mbs,&b->ilen);
1635: b->free_imax_ilen = PETSC_TRUE;
1637: PetscLogObjectMemory((PetscObject)B,2*mbs*sizeof(PetscInt));
1638: }
1639: if (!nnz) {
1640: if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
1641: else if (nz <= 0) nz = 1;
1642: for (i=0; i<mbs; i++) b->imax[i] = nz;
1643: nz = nz*mbs; /* total nz */
1644: } else {
1645: nz = 0;
1646: for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
1647: }
1648: /* b->ilen will count nonzeros in each block row so far. */
1649: for (i=0; i<mbs; i++) b->ilen[i] = 0;
1650: /* nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */
1652: /* allocate the matrix space */
1653: MatSeqXAIJFreeAIJ(B,&b->a,&b->j,&b->i);
1654: PetscMalloc3(bs2*nz,&b->a,nz,&b->j,B->rmap->N+1,&b->i);
1655: PetscLogObjectMemory((PetscObject)B,(B->rmap->N+1)*sizeof(PetscInt)+nz*(bs2*sizeof(PetscScalar)+sizeof(PetscInt)));
1656: PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));
1657: PetscMemzero(b->j,nz*sizeof(PetscInt));
1659: b->singlemalloc = PETSC_TRUE;
1661: /* pointer to beginning of each row */
1662: b->i[0] = 0;
1663: for (i=1; i<mbs+1; i++) b->i[i] = b->i[i-1] + b->imax[i-1];
1665: b->free_a = PETSC_TRUE;
1666: b->free_ij = PETSC_TRUE;
1667: } else {
1668: b->free_a = PETSC_FALSE;
1669: b->free_ij = PETSC_FALSE;
1670: }
1672: B->rmap->bs = bs;
1673: b->bs2 = bs2;
1674: b->nz = 0;
1675: b->maxnz = nz;
1677: b->inew = 0;
1678: b->jnew = 0;
1679: b->anew = 0;
1680: b->a2anew = 0;
1681: b->permute = PETSC_FALSE;
1682: if (realalloc) {MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);}
1683: return(0);
1684: }
1688: PetscErrorCode MatSeqSBAIJSetPreallocationCSR_SeqSBAIJ(Mat B,PetscInt bs,const PetscInt ii[],const PetscInt jj[], const PetscScalar V[])
1689: {
1690: PetscInt i,j,m,nz,nz_max=0,*nnz;
1691: PetscScalar *values=0;
1692: PetscBool roworiented = ((Mat_SeqSBAIJ*)B->data)->roworiented;
1695: if (bs < 1) SETERRQ1(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_OUTOFRANGE,"Invalid block size specified, must be positive but it is %D",bs);
1696: PetscLayoutSetBlockSize(B->rmap,bs);
1697: PetscLayoutSetBlockSize(B->cmap,bs);
1698: PetscLayoutSetUp(B->rmap);
1699: PetscLayoutSetUp(B->cmap);
1700: PetscLayoutGetBlockSize(B->rmap,&bs);
1701: m = B->rmap->n/bs;
1703: if (ii[0]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"ii[0] must be 0 but it is %D",ii[0]);
1704: PetscMalloc1((m+1),&nnz);
1705: for (i=0; i<m; i++) {
1706: nz = ii[i+1] - ii[i];
1707: if (nz < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row %D has a negative number of columns %D",i,nz);
1708: nz_max = PetscMax(nz_max,nz);
1709: nnz[i] = nz;
1710: }
1711: MatSeqSBAIJSetPreallocation(B,bs,0,nnz);
1712: PetscFree(nnz);
1714: values = (PetscScalar*)V;
1715: if (!values) {
1716: PetscCalloc1(bs*bs*nz_max,&values);
1717: }
1718: for (i=0; i<m; i++) {
1719: PetscInt ncols = ii[i+1] - ii[i];
1720: const PetscInt *icols = jj + ii[i];
1721: if (!roworiented || bs == 1) {
1722: const PetscScalar *svals = values + (V ? (bs*bs*ii[i]) : 0);
1723: MatSetValuesBlocked_SeqSBAIJ(B,1,&i,ncols,icols,svals,INSERT_VALUES);
1724: } else {
1725: for (j=0; j<ncols; j++) {
1726: const PetscScalar *svals = values + (V ? (bs*bs*(ii[i]+j)) : 0);
1727: MatSetValuesBlocked_SeqSBAIJ(B,1,&i,1,&icols[j],svals,INSERT_VALUES);
1728: }
1729: }
1730: }
1731: if (!V) { PetscFree(values); }
1732: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
1733: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
1734: MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
1735: return(0);
1736: }
1738: /*
1739: This is used to set the numeric factorization for both Cholesky and ICC symbolic factorization
1740: */
1743: PetscErrorCode MatSeqSBAIJSetNumericFactorization_inplace(Mat B,PetscBool natural)
1744: {
1746: PetscBool flg = PETSC_FALSE;
1747: PetscInt bs = B->rmap->bs;
1750: PetscOptionsGetBool(((PetscObject)B)->prefix,"-mat_no_unroll",&flg,NULL);
1751: if (flg) bs = 8;
1753: if (!natural) {
1754: switch (bs) {
1755: case 1:
1756: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_inplace;
1757: break;
1758: case 2:
1759: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2;
1760: break;
1761: case 3:
1762: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3;
1763: break;
1764: case 4:
1765: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4;
1766: break;
1767: case 5:
1768: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5;
1769: break;
1770: case 6:
1771: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6;
1772: break;
1773: case 7:
1774: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7;
1775: break;
1776: default:
1777: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N;
1778: break;
1779: }
1780: } else {
1781: switch (bs) {
1782: case 1:
1783: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering_inplace;
1784: break;
1785: case 2:
1786: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering;
1787: break;
1788: case 3:
1789: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering;
1790: break;
1791: case 4:
1792: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering;
1793: break;
1794: case 5:
1795: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering;
1796: break;
1797: case 6:
1798: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering;
1799: break;
1800: case 7:
1801: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering;
1802: break;
1803: default:
1804: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering;
1805: break;
1806: }
1807: }
1808: return(0);
1809: }
1811: PETSC_EXTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqAIJ(Mat, MatType,MatReuse,Mat*);
1812: PETSC_EXTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqBAIJ(Mat, MatType,MatReuse,Mat*);
1816: PETSC_EXTERN PetscErrorCode MatGetFactor_seqsbaij_petsc(Mat A,MatFactorType ftype,Mat *B)
1817: {
1818: PetscInt n = A->rmap->n;
1822: #if defined(PETSC_USE_COMPLEX)
1823: if (A->hermitian) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Hermitian Factor is not supported");
1824: #endif
1825: MatCreate(PetscObjectComm((PetscObject)A),B);
1826: MatSetSizes(*B,n,n,n,n);
1827: if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) {
1828: MatSetType(*B,MATSEQSBAIJ);
1829: MatSeqSBAIJSetPreallocation(*B,A->rmap->bs,MAT_SKIP_ALLOCATION,NULL);
1831: (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqSBAIJ;
1832: (*B)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqSBAIJ;
1833: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
1834: (*B)->factortype = ftype;
1835: return(0);
1836: }
1840: PetscErrorCode MatGetFactorAvailable_seqsbaij_petsc(Mat A,MatFactorType ftype,PetscBool *flg)
1841: {
1843: if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) {
1844: *flg = PETSC_TRUE;
1845: } else {
1846: *flg = PETSC_FALSE;
1847: }
1848: return(0);
1849: }
1851: #if defined(PETSC_HAVE_MUMPS)
1852: PETSC_EXTERN PetscErrorCode MatGetFactor_sbaij_mumps(Mat,MatFactorType,Mat*);
1853: #endif
1854: #if defined(PETSC_HAVE_PASTIX)
1855: PETSC_EXTERN PetscErrorCode MatGetFactor_seqsbaij_pastix(Mat,MatFactorType,Mat*);
1856: #endif
1857: #if defined(PETSC_HAVE_SUITESPARSE)
1858: PETSC_EXTERN PetscErrorCode MatGetFactor_seqsbaij_cholmod(Mat,MatFactorType,Mat*);
1859: #endif
1860: PETSC_EXTERN PetscErrorCode MatGetFactor_seqsbaij_sbstrm(Mat,MatFactorType,Mat*);
1862: /*MC
1863: MATSEQSBAIJ - MATSEQSBAIJ = "seqsbaij" - A matrix type to be used for sequential symmetric block sparse matrices,
1864: based on block compressed sparse row format. Only the upper triangular portion of the matrix is stored.
1866: For complex numbers by default this matrix is symmetric, NOT Hermitian symmetric. To make it Hermitian symmetric you
1867: can call MatSetOption(Mat, MAT_HERMITIAN); after MatAssemblyEnd()
1869: Options Database Keys:
1870: . -mat_type seqsbaij - sets the matrix type to "seqsbaij" during a call to MatSetFromOptions()
1872: Notes: By default if you insert values into the lower triangular part of the matrix they are simply ignored (since they are not
1873: stored and it is assumed they symmetric to the upper triangular). If you call MatSetOption(Mat,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_FALSE) or use
1874: 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.
1877: Level: beginner
1879: .seealso: MatCreateSeqSBAIJ
1880: M*/
1882: PETSC_EXTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqSBSTRM(Mat, MatType,MatReuse,Mat*);
1886: PETSC_EXTERN PetscErrorCode MatCreate_SeqSBAIJ(Mat B)
1887: {
1888: Mat_SeqSBAIJ *b;
1890: PetscMPIInt size;
1891: PetscBool no_unroll = PETSC_FALSE,no_inode = PETSC_FALSE;
1894: MPI_Comm_size(PetscObjectComm((PetscObject)B),&size);
1895: if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
1897: PetscNewLog(B,&b);
1898: B->data = (void*)b;
1899: PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
1901: B->ops->destroy = MatDestroy_SeqSBAIJ;
1902: B->ops->view = MatView_SeqSBAIJ;
1903: b->row = 0;
1904: b->icol = 0;
1905: b->reallocs = 0;
1906: b->saved_values = 0;
1907: b->inode.limit = 5;
1908: b->inode.max_limit = 5;
1910: b->roworiented = PETSC_TRUE;
1911: b->nonew = 0;
1912: b->diag = 0;
1913: b->solve_work = 0;
1914: b->mult_work = 0;
1915: B->spptr = 0;
1916: B->info.nz_unneeded = (PetscReal)b->maxnz*b->bs2;
1917: b->keepnonzeropattern = PETSC_FALSE;
1918: b->xtoy = 0;
1919: b->XtoY = 0;
1921: b->inew = 0;
1922: b->jnew = 0;
1923: b->anew = 0;
1924: b->a2anew = 0;
1925: b->permute = PETSC_FALSE;
1927: b->ignore_ltriangular = PETSC_TRUE;
1929: PetscOptionsGetBool(((PetscObject)B)->prefix,"-mat_ignore_lower_triangular",&b->ignore_ltriangular,NULL);
1931: b->getrow_utriangular = PETSC_FALSE;
1933: PetscOptionsGetBool(((PetscObject)B)->prefix,"-mat_getrow_uppertriangular",&b->getrow_utriangular,NULL);
1935: #if defined(PETSC_HAVE_PASTIX)
1936: PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_pastix_C",MatGetFactor_seqsbaij_pastix);
1937: #endif
1938: #if defined(PETSC_HAVE_MUMPS)
1939: PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_mumps_C",MatGetFactor_sbaij_mumps);
1940: #endif
1941: #if defined(PETSC_HAVE_SUITESPARSE)
1942: PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_cholmod_C",MatGetFactor_seqsbaij_cholmod);
1943: #endif
1944: PetscObjectComposeFunction((PetscObject)B,"MatGetFactorAvailable_petsc_C",MatGetFactorAvailable_seqsbaij_petsc);
1945: PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_petsc_C",MatGetFactor_seqsbaij_petsc);
1946: PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_sbstrm_C",MatGetFactor_seqsbaij_sbstrm);
1947: PetscObjectComposeFunction((PetscObject)B,"MatStoreValues_C",MatStoreValues_SeqSBAIJ);
1948: PetscObjectComposeFunction((PetscObject)B,"MatRetrieveValues_C",MatRetrieveValues_SeqSBAIJ);
1949: PetscObjectComposeFunction((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C",MatSeqSBAIJSetColumnIndices_SeqSBAIJ);
1950: PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_seqaij_C",MatConvert_SeqSBAIJ_SeqAIJ);
1951: PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_seqbaij_C",MatConvert_SeqSBAIJ_SeqBAIJ);
1952: PetscObjectComposeFunction((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",MatSeqSBAIJSetPreallocation_SeqSBAIJ);
1953: PetscObjectComposeFunction((PetscObject)B,"MatSeqSBAIJSetPreallocationCSR_C",MatSeqSBAIJSetPreallocationCSR_SeqSBAIJ);
1954: PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_seqsbstrm_C",MatConvert_SeqSBAIJ_SeqSBSTRM);
1956: B->symmetric = PETSC_TRUE;
1957: B->structurally_symmetric = PETSC_TRUE;
1958: B->symmetric_set = PETSC_TRUE;
1959: B->structurally_symmetric_set = PETSC_TRUE;
1961: PetscObjectChangeTypeName((PetscObject)B,MATSEQSBAIJ);
1963: PetscOptionsBegin(PetscObjectComm((PetscObject)B),((PetscObject)B)->prefix,"Options for SEQSBAIJ matrix","Mat");
1964: PetscOptionsBool("-mat_no_unroll","Do not optimize for inodes (slower)",NULL,no_unroll,&no_unroll,NULL);
1965: if (no_unroll) {
1966: PetscInfo(B,"Not using Inode routines due to -mat_no_unroll\n");
1967: }
1968: PetscOptionsBool("-mat_no_inode","Do not optimize for inodes (slower)",NULL,no_inode,&no_inode,NULL);
1969: if (no_inode) {
1970: PetscInfo(B,"Not using Inode routines due to -mat_no_inode\n");
1971: }
1972: PetscOptionsInt("-mat_inode_limit","Do not use inodes larger then this value",NULL,b->inode.limit,&b->inode.limit,NULL);
1973: PetscOptionsEnd();
1974: b->inode.use = (PetscBool)(!(no_unroll || no_inode));
1975: if (b->inode.limit > b->inode.max_limit) b->inode.limit = b->inode.max_limit;
1976: return(0);
1977: }
1981: /*@C
1982: MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block
1983: compressed row) format. For good matrix assembly performance the
1984: user should preallocate the matrix storage by setting the parameter nz
1985: (or the array nnz). By setting these parameters accurately, performance
1986: during matrix assembly can be increased by more than a factor of 50.
1988: Collective on Mat
1990: Input Parameters:
1991: + B - the symmetric matrix
1992: . bs - size of block
1993: . nz - number of block nonzeros per block row (same for all rows)
1994: - nnz - array containing the number of block nonzeros in the upper triangular plus
1995: diagonal portion of each block (possibly different for each block row) or NULL
1997: Options Database Keys:
1998: . -mat_no_unroll - uses code that does not unroll the loops in the
1999: block calculations (much slower)
2000: . -mat_block_size - size of the blocks to use (only works if a negative bs is passed in
2002: Level: intermediate
2004: Notes:
2005: Specify the preallocated storage with either nz or nnz (not both).
2006: Set nz=PETSC_DEFAULT and nnz=NULL for PETSc to control dynamic memory
2007: allocation. See Users-Manual: ch_mat for details.
2009: You can call MatGetInfo() to get information on how effective the preallocation was;
2010: for example the fields mallocs,nz_allocated,nz_used,nz_unneeded;
2011: You can also run with the option -info and look for messages with the string
2012: malloc in them to see if additional memory allocation was needed.
2014: If the nnz parameter is given then the nz parameter is ignored
2017: .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateSBAIJ()
2018: @*/
2019: PetscErrorCode MatSeqSBAIJSetPreallocation(Mat B,PetscInt bs,PetscInt nz,const PetscInt nnz[])
2020: {
2027: PetscTryMethod(B,"MatSeqSBAIJSetPreallocation_C",(Mat,PetscInt,PetscInt,const PetscInt[]),(B,bs,nz,nnz));
2028: return(0);
2029: }
2031: #undef __FUNCT__
2033: /*@C
2034: MatSeqSBAIJSetPreallocationCSR - Allocates memory for a sparse sequential matrix in symmetric block AIJ format.
2036: Input Parameters:
2037: + B - the matrix
2038: . i - the indices into j for the start of each local row (starts with zero)
2039: . j - the column indices for each local row (starts with zero) these must be sorted for each row
2040: - v - optional values in the matrix
2042: Level: developer
2044: Notes:
2045: The order of the entries in values is specified by the MatOption MAT_ROW_ORIENTED. For example, C programs
2046: may want to use the default MAT_ROW_ORIENTED=PETSC_TRUE and use an array v[nnz][bs][bs] where the second index is
2047: over rows within a block and the last index is over columns within a block row. Fortran programs will likely set
2048: MAT_ROW_ORIENTED=PETSC_FALSE and use a Fortran array v(bs,bs,nnz) in which the first index is over rows within a
2049: block column and the second index is over columns within a block.
2051: .keywords: matrix, block, aij, compressed row, sparse
2053: .seealso: MatCreate(), MatCreateSeqSBAIJ(), MatSetValuesBlocked(), MatSeqSBAIJSetPreallocation(), MATSEQSBAIJ
2054: @*/
2055: PetscErrorCode MatSeqSBAIJSetPreallocationCSR(Mat B,PetscInt bs,const PetscInt i[],const PetscInt j[], const PetscScalar v[])
2056: {
2063: PetscTryMethod(B,"MatSeqSBAIJSetPreallocationCSR_C",(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]),(B,bs,i,j,v));
2064: return(0);
2065: }
2069: /*@C
2070: MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block
2071: compressed row) format. For good matrix assembly performance the
2072: user should preallocate the matrix storage by setting the parameter nz
2073: (or the array nnz). By setting these parameters accurately, performance
2074: during matrix assembly can be increased by more than a factor of 50.
2076: Collective on MPI_Comm
2078: Input Parameters:
2079: + comm - MPI communicator, set to PETSC_COMM_SELF
2080: . bs - size of block
2081: . m - number of rows, or number of columns
2082: . nz - number of block nonzeros per block row (same for all rows)
2083: - nnz - array containing the number of block nonzeros in the upper triangular plus
2084: diagonal portion of each block (possibly different for each block row) or NULL
2086: Output Parameter:
2087: . A - the symmetric matrix
2089: Options Database Keys:
2090: . -mat_no_unroll - uses code that does not unroll the loops in the
2091: block calculations (much slower)
2092: . -mat_block_size - size of the blocks to use
2094: Level: intermediate
2096: It is recommended that one use the MatCreate(), MatSetType() and/or MatSetFromOptions(),
2097: MatXXXXSetPreallocation() paradgm instead of this routine directly.
2098: [MatXXXXSetPreallocation() is, for example, MatSeqAIJSetPreallocation]
2100: Notes:
2101: The number of rows and columns must be divisible by blocksize.
2102: This matrix type does not support complex Hermitian operation.
2104: Specify the preallocated storage with either nz or nnz (not both).
2105: Set nz=PETSC_DEFAULT and nnz=NULL for PETSc to control dynamic memory
2106: allocation. See Users-Manual: ch_mat for details.
2108: If the nnz parameter is given then the nz parameter is ignored
2110: .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateSBAIJ()
2111: @*/
2112: PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A)
2113: {
2117: MatCreate(comm,A);
2118: MatSetSizes(*A,m,n,m,n);
2119: MatSetType(*A,MATSEQSBAIJ);
2120: MatSeqSBAIJSetPreallocation_SeqSBAIJ(*A,bs,nz,(PetscInt*)nnz);
2121: return(0);
2122: }
2126: PetscErrorCode MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
2127: {
2128: Mat C;
2129: Mat_SeqSBAIJ *c,*a = (Mat_SeqSBAIJ*)A->data;
2131: PetscInt i,mbs = a->mbs,nz = a->nz,bs2 =a->bs2;
2134: if (a->i[mbs] != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Corrupt matrix");
2136: *B = 0;
2137: MatCreate(PetscObjectComm((PetscObject)A),&C);
2138: MatSetSizes(C,A->rmap->N,A->cmap->n,A->rmap->N,A->cmap->n);
2139: MatSetType(C,MATSEQSBAIJ);
2140: PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));
2141: c = (Mat_SeqSBAIJ*)C->data;
2143: C->preallocated = PETSC_TRUE;
2144: C->factortype = A->factortype;
2145: c->row = 0;
2146: c->icol = 0;
2147: c->saved_values = 0;
2148: c->keepnonzeropattern = a->keepnonzeropattern;
2149: C->assembled = PETSC_TRUE;
2151: PetscLayoutReference(A->rmap,&C->rmap);
2152: PetscLayoutReference(A->cmap,&C->cmap);
2153: c->bs2 = a->bs2;
2154: c->mbs = a->mbs;
2155: c->nbs = a->nbs;
2157: if (cpvalues == MAT_SHARE_NONZERO_PATTERN) {
2158: c->imax = a->imax;
2159: c->ilen = a->ilen;
2160: c->free_imax_ilen = PETSC_FALSE;
2161: } else {
2162: PetscMalloc2((mbs+1),&c->imax,(mbs+1),&c->ilen);
2163: PetscLogObjectMemory((PetscObject)C,2*(mbs+1)*sizeof(PetscInt));
2164: for (i=0; i<mbs; i++) {
2165: c->imax[i] = a->imax[i];
2166: c->ilen[i] = a->ilen[i];
2167: }
2168: c->free_imax_ilen = PETSC_TRUE;
2169: }
2171: /* allocate the matrix space */
2172: if (cpvalues == MAT_SHARE_NONZERO_PATTERN) {
2173: PetscMalloc1(bs2*nz,&c->a);
2174: PetscLogObjectMemory((PetscObject)C,nz*bs2*sizeof(MatScalar));
2175: c->i = a->i;
2176: c->j = a->j;
2177: c->singlemalloc = PETSC_FALSE;
2178: c->free_a = PETSC_TRUE;
2179: c->free_ij = PETSC_FALSE;
2180: c->parent = A;
2181: PetscObjectReference((PetscObject)A);
2182: MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
2183: MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
2184: } else {
2185: PetscMalloc3(bs2*nz,&c->a,nz,&c->j,mbs+1,&c->i);
2186: PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(PetscInt));
2187: PetscLogObjectMemory((PetscObject)C,(mbs+1)*sizeof(PetscInt) + nz*(bs2*sizeof(MatScalar) + sizeof(PetscInt)));
2188: c->singlemalloc = PETSC_TRUE;
2189: c->free_a = PETSC_TRUE;
2190: c->free_ij = PETSC_TRUE;
2191: }
2192: if (mbs > 0) {
2193: if (cpvalues != MAT_SHARE_NONZERO_PATTERN) {
2194: PetscMemcpy(c->j,a->j,nz*sizeof(PetscInt));
2195: }
2196: if (cpvalues == MAT_COPY_VALUES) {
2197: PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));
2198: } else {
2199: PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));
2200: }
2201: if (a->jshort) {
2202: /* cannot share jshort, it is reallocated in MatAssemblyEnd_SeqSBAIJ() */
2203: /* if the parent matrix is reassembled, this child matrix will never notice */
2204: PetscMalloc1(nz,&c->jshort);
2205: PetscLogObjectMemory((PetscObject)C,nz*sizeof(unsigned short));
2206: PetscMemcpy(c->jshort,a->jshort,nz*sizeof(unsigned short));
2208: c->free_jshort = PETSC_TRUE;
2209: }
2210: }
2212: c->roworiented = a->roworiented;
2213: c->nonew = a->nonew;
2215: if (a->diag) {
2216: if (cpvalues == MAT_SHARE_NONZERO_PATTERN) {
2217: c->diag = a->diag;
2218: c->free_diag = PETSC_FALSE;
2219: } else {
2220: PetscMalloc1(mbs,&c->diag);
2221: PetscLogObjectMemory((PetscObject)C,mbs*sizeof(PetscInt));
2222: for (i=0; i<mbs; i++) c->diag[i] = a->diag[i];
2223: c->free_diag = PETSC_TRUE;
2224: }
2225: }
2226: c->nz = a->nz;
2227: c->maxnz = a->nz; /* Since we allocate exactly the right amount */
2228: c->solve_work = 0;
2229: c->mult_work = 0;
2231: *B = C;
2232: PetscFunctionListDuplicate(((PetscObject)A)->qlist,&((PetscObject)C)->qlist);
2233: return(0);
2234: }
2238: PetscErrorCode MatLoad_SeqSBAIJ(Mat newmat,PetscViewer viewer)
2239: {
2240: Mat_SeqSBAIJ *a;
2242: int fd;
2243: PetscMPIInt size;
2244: PetscInt i,nz,header[4],*rowlengths=0,M,N,bs=1;
2245: PetscInt *mask,mbs,*jj,j,rowcount,nzcount,k,*s_browlengths,maskcount;
2246: PetscInt kmax,jcount,block,idx,point,nzcountb,extra_rows,rows,cols;
2247: PetscInt *masked,nmask,tmp,bs2,ishift;
2248: PetscScalar *aa;
2249: MPI_Comm comm;
2252: PetscObjectGetComm((PetscObject)viewer,&comm);
2253: PetscOptionsGetInt(((PetscObject)newmat)->prefix,"-matload_block_size",&bs,NULL);
2254: bs2 = bs*bs;
2256: MPI_Comm_size(comm,&size);
2257: if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"view must have one processor");
2258: PetscViewerBinaryGetDescriptor(viewer,&fd);
2259: PetscBinaryRead(fd,header,4,PETSC_INT);
2260: if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"not Mat object");
2261: M = header[1]; N = header[2]; nz = header[3];
2263: if (header[3] < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ");
2265: if (M != N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Can only do square matrices");
2267: /*
2268: This code adds extra rows to make sure the number of rows is
2269: divisible by the blocksize
2270: */
2271: mbs = M/bs;
2272: extra_rows = bs - M + bs*(mbs);
2273: if (extra_rows == bs) extra_rows = 0;
2274: else mbs++;
2275: if (extra_rows) {
2276: PetscInfo(viewer,"Padding loaded matrix to match blocksize\n");
2277: }
2279: /* Set global sizes if not already set */
2280: if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) {
2281: MatSetSizes(newmat,PETSC_DECIDE,PETSC_DECIDE,M+extra_rows,N+extra_rows);
2282: } else { /* Check if the matrix global sizes are correct */
2283: MatGetSize(newmat,&rows,&cols);
2284: 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);
2285: }
2287: /* read in row lengths */
2288: PetscMalloc1((M+extra_rows),&rowlengths);
2289: PetscBinaryRead(fd,rowlengths,M,PETSC_INT);
2290: for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;
2292: /* read in column indices */
2293: PetscMalloc1((nz+extra_rows),&jj);
2294: PetscBinaryRead(fd,jj,nz,PETSC_INT);
2295: for (i=0; i<extra_rows; i++) jj[nz+i] = M+i;
2297: /* loop over row lengths determining block row lengths */
2298: PetscCalloc1(mbs,&s_browlengths);
2299: PetscMalloc2(mbs,&mask,mbs,&masked);
2300: PetscMemzero(mask,mbs*sizeof(PetscInt));
2301: rowcount = 0;
2302: nzcount = 0;
2303: for (i=0; i<mbs; i++) {
2304: nmask = 0;
2305: for (j=0; j<bs; j++) {
2306: kmax = rowlengths[rowcount];
2307: for (k=0; k<kmax; k++) {
2308: tmp = jj[nzcount++]/bs; /* block col. index */
2309: if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;}
2310: }
2311: rowcount++;
2312: }
2313: s_browlengths[i] += nmask;
2315: /* zero out the mask elements we set */
2316: for (j=0; j<nmask; j++) mask[masked[j]] = 0;
2317: }
2319: /* Do preallocation */
2320: MatSeqSBAIJSetPreallocation_SeqSBAIJ(newmat,bs,0,s_browlengths);
2321: a = (Mat_SeqSBAIJ*)newmat->data;
2323: /* set matrix "i" values */
2324: a->i[0] = 0;
2325: for (i=1; i<= mbs; i++) {
2326: a->i[i] = a->i[i-1] + s_browlengths[i-1];
2327: a->ilen[i-1] = s_browlengths[i-1];
2328: }
2329: a->nz = a->i[mbs];
2331: /* read in nonzero values */
2332: PetscMalloc1((nz+extra_rows),&aa);
2333: PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);
2334: for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0;
2336: /* set "a" and "j" values into matrix */
2337: nzcount = 0; jcount = 0;
2338: for (i=0; i<mbs; i++) {
2339: nzcountb = nzcount;
2340: nmask = 0;
2341: for (j=0; j<bs; j++) {
2342: kmax = rowlengths[i*bs+j];
2343: for (k=0; k<kmax; k++) {
2344: tmp = jj[nzcount++]/bs; /* block col. index */
2345: if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;}
2346: }
2347: }
2348: /* sort the masked values */
2349: PetscSortInt(nmask,masked);
2351: /* set "j" values into matrix */
2352: maskcount = 1;
2353: for (j=0; j<nmask; j++) {
2354: a->j[jcount++] = masked[j];
2355: mask[masked[j]] = maskcount++;
2356: }
2358: /* set "a" values into matrix */
2359: ishift = bs2*a->i[i];
2360: for (j=0; j<bs; j++) {
2361: kmax = rowlengths[i*bs+j];
2362: for (k=0; k<kmax; k++) {
2363: tmp = jj[nzcountb]/bs; /* block col. index */
2364: if (tmp >= i) {
2365: block = mask[tmp] - 1;
2366: point = jj[nzcountb] - bs*tmp;
2367: idx = ishift + bs2*block + j + bs*point;
2368: a->a[idx] = aa[nzcountb];
2369: }
2370: nzcountb++;
2371: }
2372: }
2373: /* zero out the mask elements we set */
2374: for (j=0; j<nmask; j++) mask[masked[j]] = 0;
2375: }
2376: if (jcount != a->nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix");
2378: PetscFree(rowlengths);
2379: PetscFree(s_browlengths);
2380: PetscFree(aa);
2381: PetscFree(jj);
2382: PetscFree2(mask,masked);
2384: MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);
2385: MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);
2386: return(0);
2387: }
2391: /*@
2392: MatCreateSeqSBAIJWithArrays - Creates an sequential SBAIJ matrix using matrix elements
2393: (upper triangular entries in CSR format) provided by the user.
2395: Collective on MPI_Comm
2397: Input Parameters:
2398: + comm - must be an MPI communicator of size 1
2399: . bs - size of block
2400: . m - number of rows
2401: . n - number of columns
2402: . i - row indices
2403: . j - column indices
2404: - a - matrix values
2406: Output Parameter:
2407: . mat - the matrix
2409: Level: advanced
2411: Notes:
2412: The i, j, and a arrays are not copied by this routine, the user must free these arrays
2413: once the matrix is destroyed
2415: You cannot set new nonzero locations into this matrix, that will generate an error.
2417: The i and j indices are 0 based
2419: 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
2420: it is the regular CSR format excluding the lower triangular elements.
2422: .seealso: MatCreate(), MatCreateSBAIJ(), MatCreateSeqSBAIJ()
2424: @*/
2425: PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt *i,PetscInt *j,PetscScalar *a,Mat *mat)
2426: {
2428: PetscInt ii;
2429: Mat_SeqSBAIJ *sbaij;
2432: if (bs != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"block size %D > 1 is not supported yet",bs);
2433: if (i[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"i (row indices) must start with 0");
2435: MatCreate(comm,mat);
2436: MatSetSizes(*mat,m,n,m,n);
2437: MatSetType(*mat,MATSEQSBAIJ);
2438: MatSeqSBAIJSetPreallocation_SeqSBAIJ(*mat,bs,MAT_SKIP_ALLOCATION,0);
2439: sbaij = (Mat_SeqSBAIJ*)(*mat)->data;
2440: PetscMalloc2(m,&sbaij->imax,m,&sbaij->ilen);
2441: PetscLogObjectMemory((PetscObject)*mat,2*m*sizeof(PetscInt));
2443: sbaij->i = i;
2444: sbaij->j = j;
2445: sbaij->a = a;
2447: sbaij->singlemalloc = PETSC_FALSE;
2448: sbaij->nonew = -1; /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/
2449: sbaij->free_a = PETSC_FALSE;
2450: sbaij->free_ij = PETSC_FALSE;
2452: for (ii=0; ii<m; ii++) {
2453: sbaij->ilen[ii] = sbaij->imax[ii] = i[ii+1] - i[ii];
2454: #if defined(PETSC_USE_DEBUG)
2455: 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]);
2456: #endif
2457: }
2458: #if defined(PETSC_USE_DEBUG)
2459: for (ii=0; ii<sbaij->i[m]; ii++) {
2460: if (j[ii] < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column index at location = %d index = %d",ii,j[ii]);
2461: 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]);
2462: }
2463: #endif
2465: MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
2466: MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
2467: return(0);
2468: }