Actual source code: bdiag.c
1: /*$Id: bdiag.c,v 1.198 2001/08/07 03:02:53 balay Exp $*/
3: /* Block diagonal matrix format */
5: #include src/mat/impls/bdiag/seq/bdiag.h
6: #include src/vec/vecimpl.h
7: #include src/inline/ilu.h
9: EXTERN int MatSetValues_SeqBDiag_1(Mat,int,int *,int,int *,PetscScalar *,InsertMode);
10: EXTERN int MatSetValues_SeqBDiag_N(Mat,int,int *,int,int *,PetscScalar *,InsertMode);
11: EXTERN int MatGetValues_SeqBDiag_1(Mat,int,int *,int,int *,PetscScalar *);
12: EXTERN int MatGetValues_SeqBDiag_N(Mat,int,int *,int,int *,PetscScalar *);
13: EXTERN int MatMult_SeqBDiag_1(Mat,Vec,Vec);
14: EXTERN int MatMult_SeqBDiag_2(Mat,Vec,Vec);
15: EXTERN int MatMult_SeqBDiag_3(Mat,Vec,Vec);
16: EXTERN int MatMult_SeqBDiag_4(Mat,Vec,Vec);
17: EXTERN int MatMult_SeqBDiag_5(Mat,Vec,Vec);
18: EXTERN int MatMult_SeqBDiag_N(Mat,Vec,Vec);
19: EXTERN int MatMultAdd_SeqBDiag_1(Mat,Vec,Vec,Vec);
20: EXTERN int MatMultAdd_SeqBDiag_2(Mat,Vec,Vec,Vec);
21: EXTERN int MatMultAdd_SeqBDiag_3(Mat,Vec,Vec,Vec);
22: EXTERN int MatMultAdd_SeqBDiag_4(Mat,Vec,Vec,Vec);
23: EXTERN int MatMultAdd_SeqBDiag_5(Mat,Vec,Vec,Vec);
24: EXTERN int MatMultAdd_SeqBDiag_N(Mat,Vec,Vec,Vec);
25: EXTERN int MatMultTranspose_SeqBDiag_1(Mat,Vec,Vec);
26: EXTERN int MatMultTranspose_SeqBDiag_N(Mat,Vec,Vec);
27: EXTERN int MatMultTransposeAdd_SeqBDiag_1(Mat,Vec,Vec,Vec);
28: EXTERN int MatMultTransposeAdd_SeqBDiag_N(Mat,Vec,Vec,Vec);
29: EXTERN int MatRelax_SeqBDiag_N(Mat,Vec,PetscReal,MatSORType,PetscReal,int,int,Vec);
30: EXTERN int MatRelax_SeqBDiag_1(Mat,Vec,PetscReal,MatSORType,PetscReal,int,int,Vec);
31: EXTERN int MatView_SeqBDiag(Mat,PetscViewer);
32: EXTERN int MatGetInfo_SeqBDiag(Mat,MatInfoType,MatInfo*);
33: EXTERN int MatGetRow_SeqBDiag(Mat,int,int *,int **,PetscScalar **);
34: EXTERN int MatRestoreRow_SeqBDiag(Mat,int,int *,int **,PetscScalar **);
35: EXTERN int MatTranspose_SeqBDiag(Mat,Mat *);
36: EXTERN int MatNorm_SeqBDiag(Mat,NormType,PetscReal *);
38: int MatDestroy_SeqBDiag(Mat A)
39: {
40: Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
41: int i,bs = a->bs,ierr;
44: #if defined(PETSC_USE_LOG)
45: PetscLogObjectState((PetscObject)A,"Rows=%d, Cols=%d, NZ=%d, BSize=%d, NDiag=%d",A->m,A->n,a->nz,a->bs,a->nd);
46: #endif
47: if (!a->user_alloc) { /* Free the actual diagonals */
48: for (i=0; i<a->nd; i++) {
49: if (a->diag[i] > 0) {
50: PetscFree(a->diagv[i] + bs*bs*a->diag[i]);
51: } else {
52: PetscFree(a->diagv[i]);
53: }
54: }
55: }
56: if (a->pivot) {PetscFree(a->pivot);}
57: PetscFree(a->diagv);
58: PetscFree(a->diag);
59: PetscFree(a->colloc);
60: PetscFree(a->dvalue);
61: PetscFree(a);
62: return(0);
63: }
65: int MatAssemblyEnd_SeqBDiag(Mat A,MatAssemblyType mode)
66: {
67: Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
68: int i,k,temp,*diag = a->diag,*bdlen = a->bdlen;
69: PetscScalar *dtemp,**dv = a->diagv;
72: if (mode == MAT_FLUSH_ASSEMBLY) return(0);
74: /* Sort diagonals */
75: for (i=0; i<a->nd; i++) {
76: for (k=i+1; k<a->nd; k++) {
77: if (diag[i] < diag[k]) {
78: temp = diag[i];
79: diag[i] = diag[k];
80: diag[k] = temp;
81: temp = bdlen[i];
82: bdlen[i] = bdlen[k];
83: bdlen[k] = temp;
84: dtemp = dv[i];
85: dv[i] = dv[k];
86: dv[k] = dtemp;
87: }
88: }
89: }
91: /* Set location of main diagonal */
92: for (i=0; i<a->nd; i++) {
93: if (!a->diag[i]) {a->mainbd = i; break;}
94: }
95: PetscLogInfo(A,"MatAssemblyEnd_SeqBDiag:Number diagonals %d,memory used %d, block size %dn",a->nd,a->maxnz,a->bs);
96: return(0);
97: }
99: int MatSetOption_SeqBDiag(Mat A,MatOption op)
100: {
101: Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
104: switch (op) {
105: case MAT_NO_NEW_NONZERO_LOCATIONS:
106: a->nonew = 1;
107: break;
108: case MAT_YES_NEW_NONZERO_LOCATIONS:
109: a->nonew = 0;
110: break;
111: case MAT_NO_NEW_DIAGONALS:
112: a->nonew_diag = 1;
113: break;
114: case MAT_YES_NEW_DIAGONALS:
115: a->nonew_diag = 0;
116: break;
117: case MAT_COLUMN_ORIENTED:
118: a->roworiented = PETSC_FALSE;
119: break;
120: case MAT_ROW_ORIENTED:
121: a->roworiented = PETSC_TRUE;
122: break;
123: case MAT_ROWS_SORTED:
124: case MAT_ROWS_UNSORTED:
125: case MAT_COLUMNS_SORTED:
126: case MAT_COLUMNS_UNSORTED:
127: case MAT_IGNORE_OFF_PROC_ENTRIES:
128: case MAT_NEW_NONZERO_LOCATION_ERR:
129: case MAT_NEW_NONZERO_ALLOCATION_ERR:
130: case MAT_USE_HASH_TABLE:
131: case MAT_USE_SINGLE_PRECISION_SOLVES:
132: PetscLogInfo(A,"MatSetOption_SeqBDiag:Option ignoredn");
133: break;
134: default:
135: SETERRQ(PETSC_ERR_SUP,"unknown option");
136: }
137: return(0);
138: }
140: int MatPrintHelp_SeqBDiag(Mat A)
141: {
142: static PetscTruth called = PETSC_FALSE;
143: MPI_Comm comm = A->comm;
144: int ierr;
147: if (called) {return(0);} else called = PETSC_TRUE;
148: (*PetscHelpPrintf)(comm," Options for MATSEQBDIAG and MATMPIBDIAG matrix formats:n");
149: (*PetscHelpPrintf)(comm," -mat_block_size <block_size>n");
150: (*PetscHelpPrintf)(comm," -mat_bdiag_diags <d1,d2,d3,...> (diagonal numbers)n");
151: (*PetscHelpPrintf)(comm," (for example) -mat_bdiag_diags -5,-1,0,1,5n");
152: return(0);
153: }
155: static int MatGetDiagonal_SeqBDiag_N(Mat A,Vec v)
156: {
157: Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
158: int ierr,i,j,n,len,ibase,bs = a->bs,iloc;
159: PetscScalar *x,*dd,zero = 0.0;
162: if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
163: VecSet(&zero,v);
164: VecGetLocalSize(v,&n);
165: if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
166: if (a->mainbd == -1) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Main diagonal not set");
167: len = PetscMin(a->mblock,a->nblock);
168: dd = a->diagv[a->mainbd];
169: VecGetArray(v,&x);
170: for (i=0; i<len; i++) {
171: ibase = i*bs*bs; iloc = i*bs;
172: for (j=0; j<bs; j++) x[j + iloc] = dd[ibase + j*(bs+1)];
173: }
174: VecRestoreArray(v,&x);
175: return(0);
176: }
178: static int MatGetDiagonal_SeqBDiag_1(Mat A,Vec v)
179: {
180: Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
181: int ierr,i,n,len;
182: PetscScalar *x,*dd,zero = 0.0;
185: VecSet(&zero,v);
186: VecGetLocalSize(v,&n);
187: if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
188: if (a->mainbd == -1) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Main diagonal not set");
189: dd = a->diagv[a->mainbd];
190: len = PetscMin(A->m,A->n);
191: VecGetArray(v,&x);
192: for (i=0; i<len; i++) x[i] = dd[i];
193: VecRestoreArray(v,&x);
194: return(0);
195: }
197: int MatZeroEntries_SeqBDiag(Mat A)
198: {
199: Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
200: int d,i,len,bs = a->bs;
201: PetscScalar *dv;
204: for (d=0; d<a->nd; d++) {
205: dv = a->diagv[d];
206: if (a->diag[d] > 0) {
207: dv += bs*bs*a->diag[d];
208: }
209: len = a->bdlen[d]*bs*bs;
210: for (i=0; i<len; i++) dv[i] = 0.0;
211: }
212: return(0);
213: }
215: int MatGetBlockSize_SeqBDiag(Mat A,int *bs)
216: {
217: Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
220: *bs = a->bs;
221: return(0);
222: }
224: int MatZeroRows_SeqBDiag(Mat A,IS is,PetscScalar *diag)
225: {
226: Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
227: int i,ierr,N,*rows,m = A->m - 1,nz,*col;
228: PetscScalar *dd,*val;
231: ISGetLocalSize(is,&N);
232: ISGetIndices(is,&rows);
233: for (i=0; i<N; i++) {
234: if (rows[i]<0 || rows[i]>m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"row out of range");
235: MatGetRow(A,rows[i],&nz,&col,&val);
236: PetscMemzero(val,nz*sizeof(PetscScalar));
237: MatSetValues(A,1,&rows[i],nz,col,val,INSERT_VALUES);
238: MatRestoreRow(A,rows[i],&nz,&col,&val);
239: }
240: if (diag) {
241: if (a->mainbd == -1) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Main diagonal does not exist");
242: dd = a->diagv[a->mainbd];
243: for (i=0; i<N; i++) dd[rows[i]] = *diag;
244: }
245: ISRestoreIndices(is,&rows);
246: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
247: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
248: return(0);
249: }
251: int MatGetSubMatrix_SeqBDiag(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *submat)
252: {
253: Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
254: int nznew,*smap,i,j,ierr,oldcols = A->n;
255: int *irow,*icol,newr,newc,*cwork,*col,nz,bs;
256: PetscScalar *vwork,*val;
257: Mat newmat;
260: if (scall == MAT_REUSE_MATRIX) { /* no support for reuse so simply destroy all */
261: MatDestroy(*submat);
262: }
264: ISGetIndices(isrow,&irow);
265: ISGetIndices(iscol,&icol);
266: ISGetLocalSize(isrow,&newr);
267: ISGetLocalSize(iscol,&newc);
269: PetscMalloc((oldcols+1)*sizeof(int),&smap);
270: PetscMalloc((newc+1)*sizeof(int),&cwork);
271: PetscMalloc((newc+1)*sizeof(PetscScalar),&vwork);
272: ierr = PetscMemzero((char*)smap,oldcols*sizeof(int));
273: for (i=0; i<newc; i++) smap[icol[i]] = i+1;
275: /* Determine diagonals; then create submatrix */
276: bs = a->bs; /* Default block size remains the same */
277: MatCreateSeqBDiag(A->comm,newr,newc,0,bs,0,0,&newmat);
279: /* Fill new matrix */
280: for (i=0; i<newr; i++) {
281: MatGetRow(A,irow[i],&nz,&col,&val);
282: nznew = 0;
283: for (j=0; j<nz; j++) {
284: if (smap[col[j]]) {
285: cwork[nznew] = smap[col[j]] - 1;
286: vwork[nznew++] = val[j];
287: }
288: }
289: MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERT_VALUES);
290: MatRestoreRow(A,i,&nz,&col,&val);
291: }
292: MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);
293: MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);
295: /* Free work space */
296: PetscFree(smap);
297: PetscFree(cwork);
298: PetscFree(vwork);
299: ISRestoreIndices(isrow,&irow);
300: ISRestoreIndices(iscol,&icol);
301: *submat = newmat;
302: return(0);
303: }
305: int MatGetSubMatrices_SeqBDiag(Mat A,int n,IS *irow,IS *icol,MatReuse scall,Mat **B)
306: {
307: int ierr,i;
310: if (scall == MAT_INITIAL_MATRIX) {
311: PetscMalloc((n+1)*sizeof(Mat),B);
312: }
314: for (i=0; i<n; i++) {
315: MatGetSubMatrix_SeqBDiag(A,irow[i],icol[i],scall,&(*B)[i]);
316: }
317: return(0);
318: }
320: int MatScale_SeqBDiag(PetscScalar *alpha,Mat inA)
321: {
322: Mat_SeqBDiag *a = (Mat_SeqBDiag*)inA->data;
323: int one = 1,i,len,bs = a->bs;
326: for (i=0; i<a->nd; i++) {
327: len = bs*bs*a->bdlen[i];
328: if (a->diag[i] > 0) {
329: BLscal_(&len,alpha,a->diagv[i] + bs*bs*a->diag[i],&one);
330: } else {
331: BLscal_(&len,alpha,a->diagv[i],&one);
332: }
333: }
334: PetscLogFlops(a->nz);
335: return(0);
336: }
338: int MatDiagonalScale_SeqBDiag(Mat A,Vec ll,Vec rr)
339: {
340: Mat_SeqBDiag *a = (Mat_SeqBDiag*)A->data;
341: PetscScalar *l,*r,*dv;
342: int d,j,len,ierr;
343: int nd = a->nd,bs = a->bs,diag,m,n;
346: if (ll) {
347: VecGetSize(ll,&m);
348: if (m != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length");
349: if (bs == 1) {
350: VecGetArray(ll,&l);
351: for (d=0; d<nd; d++) {
352: dv = a->diagv[d];
353: diag = a->diag[d];
354: len = a->bdlen[d];
355: if (diag > 0) for (j=0; j<len; j++) dv[j+diag] *= l[j+diag];
356: else for (j=0; j<len; j++) dv[j] *= l[j];
357: }
358: VecRestoreArray(ll,&l);
359: PetscLogFlops(a->nz);
360: } else SETERRQ(PETSC_ERR_SUP,"Not yet done for bs>1");
361: }
362: if (rr) {
363: VecGetSize(rr,&n);
364: if (n != A->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vector wrong length");
365: if (bs == 1) {
366: VecGetArray(rr,&r);
367: for (d=0; d<nd; d++) {
368: dv = a->diagv[d];
369: diag = a->diag[d];
370: len = a->bdlen[d];
371: if (diag > 0) for (j=0; j<len; j++) dv[j+diag] *= r[j];
372: else for (j=0; j<len; j++) dv[j] *= r[j-diag];
373: }
374: VecRestoreArray(rr,&r);
375: PetscLogFlops(a->nz);
376: } else SETERRQ(PETSC_ERR_SUP,"Not yet done for bs>1");
377: }
378: return(0);
379: }
381: static int MatDuplicate_SeqBDiag(Mat,MatDuplicateOption,Mat *);
382: EXTERN int MatLUFactorSymbolic_SeqBDiag(Mat,IS,IS,MatLUInfo*,Mat*);
383: EXTERN int MatILUFactorSymbolic_SeqBDiag(Mat,IS,IS,MatILUInfo*,Mat*);
384: EXTERN int MatILUFactor_SeqBDiag(Mat,IS,IS,MatILUInfo*);
385: EXTERN int MatLUFactorNumeric_SeqBDiag_N(Mat,Mat*);
386: EXTERN int MatLUFactorNumeric_SeqBDiag_1(Mat,Mat*);
387: EXTERN int MatSolve_SeqBDiag_1(Mat,Vec,Vec);
388: EXTERN int MatSolve_SeqBDiag_2(Mat,Vec,Vec);
389: EXTERN int MatSolve_SeqBDiag_3(Mat,Vec,Vec);
390: EXTERN int MatSolve_SeqBDiag_4(Mat,Vec,Vec);
391: EXTERN int MatSolve_SeqBDiag_5(Mat,Vec,Vec);
392: EXTERN int MatSolve_SeqBDiag_N(Mat,Vec,Vec);
394: int MatSetUpPreallocation_SeqBDiag(Mat A)
395: {
396: int ierr;
399: MatSeqBDiagSetPreallocation(A,PETSC_DEFAULT,PETSC_DEFAULT,0,0);
400: return(0);
401: }
403: /* -------------------------------------------------------------------*/
404: static struct _MatOps MatOps_Values = {MatSetValues_SeqBDiag_N,
405: MatGetRow_SeqBDiag,
406: MatRestoreRow_SeqBDiag,
407: MatMult_SeqBDiag_N,
408: MatMultAdd_SeqBDiag_N,
409: MatMultTranspose_SeqBDiag_N,
410: MatMultTransposeAdd_SeqBDiag_N,
411: MatSolve_SeqBDiag_N,
412: 0,
413: 0,
414: 0,
415: 0,
416: 0,
417: MatRelax_SeqBDiag_N,
418: MatTranspose_SeqBDiag,
419: MatGetInfo_SeqBDiag,
420: 0,
421: MatGetDiagonal_SeqBDiag_N,
422: MatDiagonalScale_SeqBDiag,
423: MatNorm_SeqBDiag,
424: 0,
425: MatAssemblyEnd_SeqBDiag,
426: 0,
427: MatSetOption_SeqBDiag,
428: MatZeroEntries_SeqBDiag,
429: MatZeroRows_SeqBDiag,
430: 0,
431: MatLUFactorNumeric_SeqBDiag_N,
432: 0,
433: 0,
434: MatSetUpPreallocation_SeqBDiag,
435: MatILUFactorSymbolic_SeqBDiag,
436: 0,
437: 0,
438: 0,
439: MatDuplicate_SeqBDiag,
440: 0,
441: 0,
442: MatILUFactor_SeqBDiag,
443: 0,
444: 0,
445: MatGetSubMatrices_SeqBDiag,
446: 0,
447: MatGetValues_SeqBDiag_N,
448: 0,
449: MatPrintHelp_SeqBDiag,
450: MatScale_SeqBDiag,
451: 0,
452: 0,
453: 0,
454: MatGetBlockSize_SeqBDiag,
455: 0,
456: 0,
457: 0,
458: 0,
459: 0,
460: 0,
461: 0,
462: 0,
463: 0,
464: 0,
465: MatDestroy_SeqBDiag,
466: MatView_SeqBDiag,
467: MatGetPetscMaps_Petsc};
469: /*@C
470: MatSeqBDiagSetPreallocation - Sets the nonzero structure and (optionally) arrays.
472: Collective on MPI_Comm
474: Input Parameters:
475: + B - the matrix
476: . nd - number of block diagonals (optional)
477: . bs - each element of a diagonal is an bs x bs dense matrix
478: . diag - optional array of block diagonal numbers (length nd).
479: For a matrix element A[i,j], where i=row and j=column, the
480: diagonal number is
481: $ diag = i/bs - j/bs (integer division)
482: Set diag=PETSC_NULL on input for PETSc to dynamically allocate memory as
483: needed (expensive).
484: - diagv - pointer to actual diagonals (in same order as diag array),
485: if allocated by user. Otherwise, set diagv=PETSC_NULL on input for PETSc
486: to control memory allocation.
488: Options Database Keys:
489: . -mat_block_size <bs> - Sets blocksize
490: . -mat_bdiag_diags <s1,s2,s3,...> - Sets diagonal numbers
492: Notes:
493: See the users manual for further details regarding this storage format.
495: Fortran Note:
496: Fortran programmers cannot set diagv; this value is ignored.
498: Level: intermediate
500: .keywords: matrix, block, diagonal, sparse
502: .seealso: MatCreate(), MatCreateMPIBDiag(), MatSetValues()
503: @*/
504: int MatSeqBDiagSetPreallocation(Mat B,int nd,int bs,int *diag,PetscScalar **diagv)
505: {
506: Mat_SeqBDiag *b;
507: int i,nda,sizetot,ierr, nd2 = 128,idiag[128];
508: PetscTruth flg1;
511: PetscTypeCompare((PetscObject)B,MATSEQBDIAG,&flg1);
512: if (!flg1) return(0);
514: B->preallocated = PETSC_TRUE;
515: if (bs == PETSC_DEFAULT) bs = 1;
516: if (bs == 0) SETERRQ(1,"Blocksize cannot be zero");
517: if (nd == PETSC_DEFAULT) nd = 0;
518: PetscOptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,PETSC_NULL);
519: PetscOptionsGetIntArray(PETSC_NULL,"-mat_bdiag_diags",idiag,&nd2,&flg1);
520: if (flg1) {
521: diag = idiag;
522: nd = nd2;
523: }
525: if ((B->n%bs) || (B->m%bs)) SETERRQ(PETSC_ERR_ARG_SIZ,"Invalid block size");
526: if (!nd) nda = nd + 1;
527: else nda = nd;
528: b = (Mat_SeqBDiag*)B->data;
530: PetscOptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg1);
531: if (!flg1) {
532: switch (bs) {
533: case 1:
534: B->ops->setvalues = MatSetValues_SeqBDiag_1;
535: B->ops->getvalues = MatGetValues_SeqBDiag_1;
536: B->ops->getdiagonal = MatGetDiagonal_SeqBDiag_1;
537: B->ops->mult = MatMult_SeqBDiag_1;
538: B->ops->multadd = MatMultAdd_SeqBDiag_1;
539: B->ops->multtranspose = MatMultTranspose_SeqBDiag_1;
540: B->ops->multtransposeadd= MatMultTransposeAdd_SeqBDiag_1;
541: B->ops->relax = MatRelax_SeqBDiag_1;
542: B->ops->solve = MatSolve_SeqBDiag_1;
543: B->ops->lufactornumeric = MatLUFactorNumeric_SeqBDiag_1;
544: break;
545: case 2:
546: B->ops->mult = MatMult_SeqBDiag_2;
547: B->ops->multadd = MatMultAdd_SeqBDiag_2;
548: B->ops->solve = MatSolve_SeqBDiag_2;
549: break;
550: case 3:
551: B->ops->mult = MatMult_SeqBDiag_3;
552: B->ops->multadd = MatMultAdd_SeqBDiag_3;
553: B->ops->solve = MatSolve_SeqBDiag_3;
554: break;
555: case 4:
556: B->ops->mult = MatMult_SeqBDiag_4;
557: B->ops->multadd = MatMultAdd_SeqBDiag_4;
558: B->ops->solve = MatSolve_SeqBDiag_4;
559: break;
560: case 5:
561: B->ops->mult = MatMult_SeqBDiag_5;
562: B->ops->multadd = MatMultAdd_SeqBDiag_5;
563: B->ops->solve = MatSolve_SeqBDiag_5;
564: break;
565: }
566: }
568: b->mblock = B->m/bs;
569: b->nblock = B->n/bs;
570: b->nd = nd;
571: b->bs = bs;
572: b->ndim = 0;
573: b->mainbd = -1;
574: b->pivot = 0;
576: ierr = PetscMalloc(2*nda*sizeof(int),&b->diag);
577: b->bdlen = b->diag + nda;
578: ierr = PetscMalloc((B->n+1)*sizeof(int),&b->colloc);
579: ierr = PetscMalloc(nda*sizeof(PetscScalar*),&b->diagv);
580: sizetot = 0;
582: if (diagv) { /* user allocated space */
583: b->user_alloc = PETSC_TRUE;
584: for (i=0; i<nd; i++) b->diagv[i] = diagv[i];
585: } else b->user_alloc = PETSC_FALSE;
587: for (i=0; i<nd; i++) {
588: b->diag[i] = diag[i];
589: if (diag[i] > 0) { /* lower triangular */
590: b->bdlen[i] = PetscMin(b->nblock,b->mblock - diag[i]);
591: } else { /* upper triangular */
592: b->bdlen[i] = PetscMin(b->mblock,b->nblock + diag[i]);
593: }
594: sizetot += b->bdlen[i];
595: }
596: sizetot *= bs*bs;
597: b->maxnz = sizetot;
598: ierr = PetscMalloc((B->n+1)*sizeof(PetscScalar),&b->dvalue);
599: PetscLogObjectMemory(B,(nda*(bs+2))*sizeof(int) + bs*nda*sizeof(PetscScalar)
600: + nda*sizeof(PetscScalar*) + sizeof(Mat_SeqBDiag)
601: + sizeof(struct _p_Mat) + sizetot*sizeof(PetscScalar));
603: if (!b->user_alloc) {
604: for (i=0; i<nd; i++) {
605: PetscMalloc(bs*bs*b->bdlen[i]*sizeof(PetscScalar),&b->diagv[i]);
606: PetscMemzero(b->diagv[i],bs*bs*b->bdlen[i]*sizeof(PetscScalar));
607: }
608: b->nonew = 0; b->nonew_diag = 0;
609: } else { /* diagonals are set on input; don't allow dynamic allocation */
610: b->nonew = 1; b->nonew_diag = 1;
611: }
613: /* adjust diagv so one may access rows with diagv[diag][row] for all rows */
614: for (i=0; i<nd; i++) {
615: if (diag[i] > 0) {
616: b->diagv[i] -= bs*bs*diag[i];
617: }
618: }
620: b->nz = b->maxnz; /* Currently not keeping track of exact count */
621: b->roworiented = PETSC_TRUE;
622: B->info.nz_unneeded = (double)b->maxnz;
623: return(0);
624: }
626: static int MatDuplicate_SeqBDiag(Mat A,MatDuplicateOption cpvalues,Mat *matout)
627: {
628: Mat_SeqBDiag *newmat,*a = (Mat_SeqBDiag*)A->data;
629: int i,ierr,len,diag,bs = a->bs;
630: Mat mat;
633: MatCreateSeqBDiag(A->comm,A->m,A->n,a->nd,bs,a->diag,PETSC_NULL,matout);
635: /* Copy contents of diagonals */
636: mat = *matout;
637: newmat = (Mat_SeqBDiag*)mat->data;
638: if (cpvalues == MAT_COPY_VALUES) {
639: for (i=0; i<a->nd; i++) {
640: len = a->bdlen[i] * bs * bs * sizeof(PetscScalar);
641: diag = a->diag[i];
642: if (diag > 0) {
643: PetscMemcpy(newmat->diagv[i]+bs*bs*diag,a->diagv[i]+bs*bs*diag,len);
644: } else {
645: PetscMemcpy(newmat->diagv[i],a->diagv[i],len);
646: }
647: }
648: }
649: MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
650: MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);
651: return(0);
652: }
654: EXTERN_C_BEGIN
655: int MatLoad_SeqBDiag(PetscViewer viewer,MatType type,Mat *A)
656: {
657: Mat B;
658: int *scols,i,nz,ierr,fd,header[4],size,nd = 128;
659: int bs,*rowlengths = 0,M,N,*cols,extra_rows,*diag = 0;
660: int idiag[128];
661: PetscScalar *vals,*svals;
662: MPI_Comm comm;
663: PetscTruth flg;
664:
666: PetscObjectGetComm((PetscObject)viewer,&comm);
667: MPI_Comm_size(comm,&size);
668: if (size > 1) SETERRQ(PETSC_ERR_ARG_SIZ,"view must have one processor");
669: PetscViewerBinaryGetDescriptor(viewer,&fd);
670: PetscBinaryRead(fd,header,4,PETSC_INT);
671: if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object");
672: M = header[1]; N = header[2]; nz = header[3];
673: if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only load square matrices");
674: if (header[3] < 0) {
675: SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqBDiag");
676: }
678: /*
679: This code adds extra rows to make sure the number of rows is
680: divisible by the blocksize
681: */
682: bs = 1;
683: PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);
684: extra_rows = bs - M + bs*(M/bs);
685: if (extra_rows == bs) extra_rows = 0;
686: if (extra_rows) {
687: PetscLogInfo(0,"MatLoad_SeqBDiag:Padding loaded matrix to match blocksizen");
688: }
690: /* read row lengths */
691: PetscMalloc((M+extra_rows)*sizeof(int),&rowlengths);
692: PetscBinaryRead(fd,rowlengths,M,PETSC_INT);
693: for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;
695: /* load information about diagonals */
696: PetscOptionsGetIntArray(PETSC_NULL,"-matload_bdiag_diags",idiag,&nd,&flg);
697: if (flg) {
698: diag = idiag;
699: }
701: /* create our matrix */
702: MatCreateSeqBDiag(comm,M+extra_rows,M+extra_rows,nd,bs,diag,
703: PETSC_NULL,A);
704: B = *A;
706: /* read column indices and nonzeros */
707: PetscMalloc(nz*sizeof(int),&scols);
708: cols = scols;
709: PetscBinaryRead(fd,cols,nz,PETSC_INT);
710: PetscMalloc(nz*sizeof(PetscScalar),&svals);
711: vals = svals;
712: PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
713: /* insert into matrix */
715: for (i=0; i<M; i++) {
716: MatSetValues(B,1,&i,rowlengths[i],scols,svals,INSERT_VALUES);
717: scols += rowlengths[i]; svals += rowlengths[i];
718: }
719: vals[0] = 1.0;
720: for (i=M; i<M+extra_rows; i++) {
721: MatSetValues(B,1,&i,1,&i,vals,INSERT_VALUES);
722: }
724: PetscFree(cols);
725: PetscFree(vals);
726: PetscFree(rowlengths);
728: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
729: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
730: return(0);
731: }
732: EXTERN_C_END
734: EXTERN_C_BEGIN
735: int MatCreate_SeqBDiag(Mat B)
736: {
737: Mat_SeqBDiag *b;
738: int ierr,size;
741: MPI_Comm_size(B->comm,&size);
742: if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
744: B->m = B->M = PetscMax(B->m,B->M);
745: B->n = B->N = PetscMax(B->n,B->N);
747: ierr = PetscNew(Mat_SeqBDiag,&b);
748: B->data = (void*)b;
749: ierr = PetscMemzero(b,sizeof(Mat_SeqBDiag));
750: ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
751: B->factor = 0;
752: B->mapping = 0;
754: PetscMapCreateMPI(B->comm,B->m,B->m,&B->rmap);
755: PetscMapCreateMPI(B->comm,B->n,B->n,&B->cmap);
757: b->ndim = 0;
758: b->mainbd = -1;
759: b->pivot = 0;
761: b->roworiented = PETSC_TRUE;
762: return(0);
763: }
764: EXTERN_C_END
766: /*@C
767: MatCreateSeqBDiag - Creates a sequential block diagonal matrix.
769: Collective on MPI_Comm
771: Input Parameters:
772: + comm - MPI communicator, set to PETSC_COMM_SELF
773: . m - number of rows
774: . n - number of columns
775: . nd - number of block diagonals (optional)
776: . bs - each element of a diagonal is an bs x bs dense matrix
777: . diag - optional array of block diagonal numbers (length nd).
778: For a matrix element A[i,j], where i=row and j=column, the
779: diagonal number is
780: $ diag = i/bs - j/bs (integer division)
781: Set diag=PETSC_NULL on input for PETSc to dynamically allocate memory as
782: needed (expensive).
783: - diagv - pointer to actual diagonals (in same order as diag array),
784: if allocated by user. Otherwise, set diagv=PETSC_NULL on input for PETSc
785: to control memory allocation.
787: Output Parameters:
788: . A - the matrix
790: Options Database Keys:
791: . -mat_block_size <bs> - Sets blocksize
792: . -mat_bdiag_diags <s1,s2,s3,...> - Sets diagonal numbers
794: Notes:
795: See the users manual for further details regarding this storage format.
797: Fortran Note:
798: Fortran programmers cannot set diagv; this value is ignored.
800: Level: intermediate
802: .keywords: matrix, block, diagonal, sparse
804: .seealso: MatCreate(), MatCreateMPIBDiag(), MatSetValues()
805: @*/
806: int MatCreateSeqBDiag(MPI_Comm comm,int m,int n,int nd,int bs,int *diag,PetscScalar **diagv,Mat *A)
807: {
811: MatCreate(comm,m,n,m,n,A);
812: MatSetType(*A,MATSEQBDIAG);
813: MatSeqBDiagSetPreallocation(*A,nd,bs,diag,diagv);
814: return(0);
815: }