Actual source code: blockmat.c
petsc-3.5.4 2015-05-23
2: /*
3: This provides a matrix that consists of Mats
4: */
6: #include <petsc-private/matimpl.h> /*I "petscmat.h" I*/
7: #include <../src/mat/impls/baij/seq/baij.h> /* use the common AIJ data-structure */
9: typedef struct {
10: SEQAIJHEADER(Mat);
11: SEQBAIJHEADER;
12: Mat *diags;
14: Vec left,right,middle,workb; /* dummy vectors to perform local parts of product */
15: } Mat_BlockMat;
17: extern PetscErrorCode MatBlockMatSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt*);
21: PetscErrorCode MatSOR_BlockMat_Symmetric(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
22: {
23: Mat_BlockMat *a = (Mat_BlockMat*)A->data;
24: PetscScalar *x;
25: const Mat *v;
26: const PetscScalar *b;
27: PetscErrorCode ierr;
28: PetscInt n = A->cmap->n,i,mbs = n/A->rmap->bs,j,bs = A->rmap->bs;
29: const PetscInt *idx;
30: IS row,col;
31: MatFactorInfo info;
32: Vec left = a->left,right = a->right, middle = a->middle;
33: Mat *diag;
36: its = its*lits;
37: if (its <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
38: if (flag & SOR_EISENSTAT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for Eisenstat");
39: if (omega != 1.0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for omega not equal to 1.0");
40: if (fshift) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for fshift");
41: if ((flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) && !(flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP)) {
42: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot do backward sweep without forward sweep");
43: }
45: if (!a->diags) {
46: PetscMalloc1(mbs,&a->diags);
47: MatFactorInfoInitialize(&info);
48: for (i=0; i<mbs; i++) {
49: MatGetOrdering(a->a[a->diag[i]], MATORDERINGND,&row,&col);
50: MatCholeskyFactorSymbolic(a->diags[i],a->a[a->diag[i]],row,&info);
51: MatCholeskyFactorNumeric(a->diags[i],a->a[a->diag[i]],&info);
52: ISDestroy(&row);
53: ISDestroy(&col);
54: }
55: VecDuplicate(bb,&a->workb);
56: }
57: diag = a->diags;
59: VecSet(xx,0.0);
60: VecGetArray(xx,&x);
61: /* copy right hand side because it must be modified during iteration */
62: VecCopy(bb,a->workb);
63: VecGetArrayRead(a->workb,&b);
65: /* need to add code for when initial guess is zero, see MatSOR_SeqAIJ */
66: while (its--) {
67: if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
69: for (i=0; i<mbs; i++) {
70: n = a->i[i+1] - a->i[i] - 1;
71: idx = a->j + a->i[i] + 1;
72: v = a->a + a->i[i] + 1;
74: VecSet(left,0.0);
75: for (j=0; j<n; j++) {
76: VecPlaceArray(right,x + idx[j]*bs);
77: MatMultAdd(v[j],right,left,left);
78: VecResetArray(right);
79: }
80: VecPlaceArray(right,b + i*bs);
81: VecAYPX(left,-1.0,right);
82: VecResetArray(right);
84: VecPlaceArray(right,x + i*bs);
85: MatSolve(diag[i],left,right);
87: /* now adjust right hand side, see MatSOR_SeqSBAIJ */
88: for (j=0; j<n; j++) {
89: MatMultTranspose(v[j],right,left);
90: VecPlaceArray(middle,b + idx[j]*bs);
91: VecAXPY(middle,-1.0,left);
92: VecResetArray(middle);
93: }
94: VecResetArray(right);
96: }
97: }
98: if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
100: for (i=mbs-1; i>=0; i--) {
101: n = a->i[i+1] - a->i[i] - 1;
102: idx = a->j + a->i[i] + 1;
103: v = a->a + a->i[i] + 1;
105: VecSet(left,0.0);
106: for (j=0; j<n; j++) {
107: VecPlaceArray(right,x + idx[j]*bs);
108: MatMultAdd(v[j],right,left,left);
109: VecResetArray(right);
110: }
111: VecPlaceArray(right,b + i*bs);
112: VecAYPX(left,-1.0,right);
113: VecResetArray(right);
115: VecPlaceArray(right,x + i*bs);
116: MatSolve(diag[i],left,right);
117: VecResetArray(right);
119: }
120: }
121: }
122: VecRestoreArray(xx,&x);
123: VecRestoreArrayRead(a->workb,&b);
124: return(0);
125: }
129: PetscErrorCode MatSOR_BlockMat(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
130: {
131: Mat_BlockMat *a = (Mat_BlockMat*)A->data;
132: PetscScalar *x;
133: const Mat *v;
134: const PetscScalar *b;
135: PetscErrorCode ierr;
136: PetscInt n = A->cmap->n,i,mbs = n/A->rmap->bs,j,bs = A->rmap->bs;
137: const PetscInt *idx;
138: IS row,col;
139: MatFactorInfo info;
140: Vec left = a->left,right = a->right;
141: Mat *diag;
144: its = its*lits;
145: if (its <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
146: if (flag & SOR_EISENSTAT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for Eisenstat");
147: if (omega != 1.0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for omega not equal to 1.0");
148: if (fshift) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for fshift");
150: if (!a->diags) {
151: PetscMalloc1(mbs,&a->diags);
152: MatFactorInfoInitialize(&info);
153: for (i=0; i<mbs; i++) {
154: MatGetOrdering(a->a[a->diag[i]], MATORDERINGND,&row,&col);
155: MatLUFactorSymbolic(a->diags[i],a->a[a->diag[i]],row,col,&info);
156: MatLUFactorNumeric(a->diags[i],a->a[a->diag[i]],&info);
157: ISDestroy(&row);
158: ISDestroy(&col);
159: }
160: }
161: diag = a->diags;
163: VecSet(xx,0.0);
164: VecGetArray(xx,&x);
165: VecGetArrayRead(bb,&b);
167: /* need to add code for when initial guess is zero, see MatSOR_SeqAIJ */
168: while (its--) {
169: if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) {
171: for (i=0; i<mbs; i++) {
172: n = a->i[i+1] - a->i[i];
173: idx = a->j + a->i[i];
174: v = a->a + a->i[i];
176: VecSet(left,0.0);
177: for (j=0; j<n; j++) {
178: if (idx[j] != i) {
179: VecPlaceArray(right,x + idx[j]*bs);
180: MatMultAdd(v[j],right,left,left);
181: VecResetArray(right);
182: }
183: }
184: VecPlaceArray(right,b + i*bs);
185: VecAYPX(left,-1.0,right);
186: VecResetArray(right);
188: VecPlaceArray(right,x + i*bs);
189: MatSolve(diag[i],left,right);
190: VecResetArray(right);
191: }
192: }
193: if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) {
195: for (i=mbs-1; i>=0; i--) {
196: n = a->i[i+1] - a->i[i];
197: idx = a->j + a->i[i];
198: v = a->a + a->i[i];
200: VecSet(left,0.0);
201: for (j=0; j<n; j++) {
202: if (idx[j] != i) {
203: VecPlaceArray(right,x + idx[j]*bs);
204: MatMultAdd(v[j],right,left,left);
205: VecResetArray(right);
206: }
207: }
208: VecPlaceArray(right,b + i*bs);
209: VecAYPX(left,-1.0,right);
210: VecResetArray(right);
212: VecPlaceArray(right,x + i*bs);
213: MatSolve(diag[i],left,right);
214: VecResetArray(right);
216: }
217: }
218: }
219: VecRestoreArray(xx,&x);
220: VecRestoreArrayRead(bb,&b);
221: return(0);
222: }
226: PetscErrorCode MatSetValues_BlockMat(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
227: {
228: Mat_BlockMat *a = (Mat_BlockMat*)A->data;
229: PetscInt *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,lastcol = -1;
230: PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen;
231: PetscInt *aj =a->j,nonew=a->nonew,bs=A->rmap->bs,brow,bcol;
233: PetscInt ridx,cidx;
234: PetscBool roworiented=a->roworiented;
235: MatScalar value;
236: Mat *ap,*aa = a->a;
239: for (k=0; k<m; k++) { /* loop over added rows */
240: row = im[k];
241: brow = row/bs;
242: if (row < 0) continue;
243: #if defined(PETSC_USE_DEBUG)
244: 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);
245: #endif
246: rp = aj + ai[brow];
247: ap = aa + ai[brow];
248: rmax = imax[brow];
249: nrow = ailen[brow];
250: low = 0;
251: high = nrow;
252: for (l=0; l<n; l++) { /* loop over added columns */
253: if (in[l] < 0) continue;
254: #if defined(PETSC_USE_DEBUG)
255: 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);
256: #endif
257: col = in[l]; bcol = col/bs;
258: if (A->symmetric && brow > bcol) continue;
259: ridx = row % bs; cidx = col % bs;
260: if (roworiented) value = v[l + k*n];
261: else value = v[k + l*m];
263: if (col <= lastcol) low = 0;
264: else high = nrow;
265: lastcol = col;
266: while (high-low > 7) {
267: t = (low+high)/2;
268: if (rp[t] > bcol) high = t;
269: else low = t;
270: }
271: for (i=low; i<high; i++) {
272: if (rp[i] > bcol) break;
273: if (rp[i] == bcol) goto noinsert1;
274: }
275: if (nonew == 1) goto noinsert1;
276: if (nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col);
277: MatSeqXAIJReallocateAIJ(A,a->mbs,1,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,imax,nonew,Mat);
278: N = nrow++ - 1; high++;
279: /* shift up all the later entries in this row */
280: for (ii=N; ii>=i; ii--) {
281: rp[ii+1] = rp[ii];
282: ap[ii+1] = ap[ii];
283: }
284: if (N>=i) ap[i] = 0;
285: rp[i] = bcol;
286: a->nz++;
287: A->nonzerostate++;
288: noinsert1:;
289: if (!*(ap+i)) {
290: MatCreateSeqAIJ(PETSC_COMM_SELF,bs,bs,0,0,ap+i);
291: }
292: MatSetValues(ap[i],1,&ridx,1,&cidx,&value,is);
293: low = i;
294: }
295: ailen[brow] = nrow;
296: }
297: return(0);
298: }
302: PetscErrorCode MatLoad_BlockMat(Mat newmat, PetscViewer viewer)
303: {
304: PetscErrorCode ierr;
305: Mat tmpA;
306: PetscInt i,j,m,n,bs = 1,ncols,*lens,currentcol,mbs,**ii,*ilens,nextcol,*llens,cnt = 0;
307: const PetscInt *cols;
308: const PetscScalar *values;
309: PetscBool flg = PETSC_FALSE,notdone;
310: Mat_SeqAIJ *a;
311: Mat_BlockMat *amat;
314: MatCreate(PETSC_COMM_SELF,&tmpA);
315: MatSetType(tmpA,MATSEQAIJ);
316: MatLoad_SeqAIJ(tmpA,viewer);
318: MatGetLocalSize(tmpA,&m,&n);
319: PetscOptionsBegin(PETSC_COMM_SELF,NULL,"Options for loading BlockMat matrix 1","Mat");
320: PetscOptionsInt("-matload_block_size","Set the blocksize used to store the matrix","MatLoad",bs,&bs,NULL);
321: PetscOptionsBool("-matload_symmetric","Store the matrix as symmetric","MatLoad",flg,&flg,NULL);
322: PetscOptionsEnd();
324: /* Determine number of nonzero blocks for each block row */
325: a = (Mat_SeqAIJ*) tmpA->data;
326: mbs = m/bs;
327: PetscMalloc3(mbs,&lens,bs,&ii,bs,&ilens);
328: PetscMemzero(lens,mbs*sizeof(PetscInt));
330: for (i=0; i<mbs; i++) {
331: for (j=0; j<bs; j++) {
332: ii[j] = a->j + a->i[i*bs + j];
333: ilens[j] = a->i[i*bs + j + 1] - a->i[i*bs + j];
334: }
336: currentcol = -1;
337: notdone = PETSC_TRUE;
338: while (PETSC_TRUE) {
339: notdone = PETSC_FALSE;
340: nextcol = 1000000000;
341: for (j=0; j<bs; j++) {
342: while ((ilens[j] > 0 && ii[j][0]/bs <= currentcol)) {
343: ii[j]++;
344: ilens[j]--;
345: }
346: if (ilens[j] > 0) {
347: notdone = PETSC_TRUE;
348: nextcol = PetscMin(nextcol,ii[j][0]/bs);
349: }
350: }
351: if (!notdone) break;
352: if (!flg || (nextcol >= i)) lens[i]++;
353: currentcol = nextcol;
354: }
355: }
357: if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) {
358: MatSetSizes(newmat,m,n,PETSC_DETERMINE,PETSC_DETERMINE);
359: }
360: MatBlockMatSetPreallocation(newmat,bs,0,lens);
361: if (flg) {
362: MatSetOption(newmat,MAT_SYMMETRIC,PETSC_TRUE);
363: }
364: amat = (Mat_BlockMat*)(newmat)->data;
366: /* preallocate the submatrices */
367: PetscMalloc1(bs,&llens);
368: for (i=0; i<mbs; i++) { /* loops for block rows */
369: for (j=0; j<bs; j++) {
370: ii[j] = a->j + a->i[i*bs + j];
371: ilens[j] = a->i[i*bs + j + 1] - a->i[i*bs + j];
372: }
374: currentcol = 1000000000;
375: for (j=0; j<bs; j++) { /* loop over rows in block finding first nonzero block */
376: if (ilens[j] > 0) {
377: currentcol = PetscMin(currentcol,ii[j][0]/bs);
378: }
379: }
381: notdone = PETSC_TRUE;
382: while (PETSC_TRUE) { /* loops over blocks in block row */
384: notdone = PETSC_FALSE;
385: nextcol = 1000000000;
386: PetscMemzero(llens,bs*sizeof(PetscInt));
387: for (j=0; j<bs; j++) { /* loop over rows in block */
388: while ((ilens[j] > 0 && ii[j][0]/bs <= currentcol)) { /* loop over columns in row */
389: ii[j]++;
390: ilens[j]--;
391: llens[j]++;
392: }
393: if (ilens[j] > 0) {
394: notdone = PETSC_TRUE;
395: nextcol = PetscMin(nextcol,ii[j][0]/bs);
396: }
397: }
398: if (cnt >= amat->maxnz) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Number of blocks found greater than expected %D",cnt);
399: if (!flg || currentcol >= i) {
400: amat->j[cnt] = currentcol;
401: MatCreateSeqAIJ(PETSC_COMM_SELF,bs,bs,0,llens,amat->a+cnt++);
402: }
404: if (!notdone) break;
405: currentcol = nextcol;
406: }
407: amat->ilen[i] = lens[i];
408: }
410: PetscFree3(lens,ii,ilens);
411: PetscFree(llens);
413: /* copy over the matrix, one row at a time */
414: for (i=0; i<m; i++) {
415: MatGetRow(tmpA,i,&ncols,&cols,&values);
416: MatSetValues(newmat,1,&i,ncols,cols,values,INSERT_VALUES);
417: MatRestoreRow(tmpA,i,&ncols,&cols,&values);
418: }
419: MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);
420: MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);
421: return(0);
422: }
426: PetscErrorCode MatView_BlockMat(Mat A,PetscViewer viewer)
427: {
428: Mat_BlockMat *a = (Mat_BlockMat*)A->data;
429: PetscErrorCode ierr;
430: const char *name;
431: PetscViewerFormat format;
434: PetscObjectGetName((PetscObject)A,&name);
435: PetscViewerGetFormat(viewer,&format);
436: if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) {
437: PetscViewerASCIIPrintf(viewer,"Nonzero block matrices = %D \n",a->nz);
438: if (A->symmetric) {
439: PetscViewerASCIIPrintf(viewer,"Only upper triangular part of symmetric matrix is stored\n");
440: }
441: }
442: return(0);
443: }
447: PetscErrorCode MatDestroy_BlockMat(Mat mat)
448: {
450: Mat_BlockMat *bmat = (Mat_BlockMat*)mat->data;
451: PetscInt i;
454: VecDestroy(&bmat->right);
455: VecDestroy(&bmat->left);
456: VecDestroy(&bmat->middle);
457: VecDestroy(&bmat->workb);
458: if (bmat->diags) {
459: for (i=0; i<mat->rmap->n/mat->rmap->bs; i++) {
460: MatDestroy(&bmat->diags[i]);
461: }
462: }
463: if (bmat->a) {
464: for (i=0; i<bmat->nz; i++) {
465: MatDestroy(&bmat->a[i]);
466: }
467: }
468: MatSeqXAIJFreeAIJ(mat,(PetscScalar**)&bmat->a,&bmat->j,&bmat->i);
469: PetscFree(mat->data);
470: return(0);
471: }
475: PetscErrorCode MatMult_BlockMat(Mat A,Vec x,Vec y)
476: {
477: Mat_BlockMat *bmat = (Mat_BlockMat*)A->data;
479: PetscScalar *xx,*yy;
480: PetscInt *aj,i,*ii,jrow,m = A->rmap->n/A->rmap->bs,bs = A->rmap->bs,n,j;
481: Mat *aa;
484: /*
485: Standard CSR multiply except each entry is a Mat
486: */
487: VecGetArray(x,&xx);
489: VecSet(y,0.0);
490: VecGetArray(y,&yy);
491: aj = bmat->j;
492: aa = bmat->a;
493: ii = bmat->i;
494: for (i=0; i<m; i++) {
495: jrow = ii[i];
496: VecPlaceArray(bmat->left,yy + bs*i);
497: n = ii[i+1] - jrow;
498: for (j=0; j<n; j++) {
499: VecPlaceArray(bmat->right,xx + bs*aj[jrow]);
500: MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);
501: VecResetArray(bmat->right);
502: jrow++;
503: }
504: VecResetArray(bmat->left);
505: }
506: VecRestoreArray(x,&xx);
507: VecRestoreArray(y,&yy);
508: return(0);
509: }
513: PetscErrorCode MatMult_BlockMat_Symmetric(Mat A,Vec x,Vec y)
514: {
515: Mat_BlockMat *bmat = (Mat_BlockMat*)A->data;
517: PetscScalar *xx,*yy;
518: PetscInt *aj,i,*ii,jrow,m = A->rmap->n/A->rmap->bs,bs = A->rmap->bs,n,j;
519: Mat *aa;
522: /*
523: Standard CSR multiply except each entry is a Mat
524: */
525: VecGetArray(x,&xx);
527: VecSet(y,0.0);
528: VecGetArray(y,&yy);
529: aj = bmat->j;
530: aa = bmat->a;
531: ii = bmat->i;
532: for (i=0; i<m; i++) {
533: jrow = ii[i];
534: n = ii[i+1] - jrow;
535: VecPlaceArray(bmat->left,yy + bs*i);
536: VecPlaceArray(bmat->middle,xx + bs*i);
537: /* if we ALWAYS required a diagonal entry then could remove this if test */
538: if (aj[jrow] == i) {
539: VecPlaceArray(bmat->right,xx + bs*aj[jrow]);
540: MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);
541: VecResetArray(bmat->right);
542: jrow++;
543: n--;
544: }
545: for (j=0; j<n; j++) {
546: VecPlaceArray(bmat->right,xx + bs*aj[jrow]); /* upper triangular part */
547: MatMultAdd(aa[jrow],bmat->right,bmat->left,bmat->left);
548: VecResetArray(bmat->right);
550: VecPlaceArray(bmat->right,yy + bs*aj[jrow]); /* lower triangular part */
551: MatMultTransposeAdd(aa[jrow],bmat->middle,bmat->right,bmat->right);
552: VecResetArray(bmat->right);
553: jrow++;
554: }
555: VecResetArray(bmat->left);
556: VecResetArray(bmat->middle);
557: }
558: VecRestoreArray(x,&xx);
559: VecRestoreArray(y,&yy);
560: return(0);
561: }
565: PetscErrorCode MatMultAdd_BlockMat(Mat A,Vec x,Vec y,Vec z)
566: {
568: return(0);
569: }
573: PetscErrorCode MatMultTranspose_BlockMat(Mat A,Vec x,Vec y)
574: {
576: return(0);
577: }
581: PetscErrorCode MatMultTransposeAdd_BlockMat(Mat A,Vec x,Vec y,Vec z)
582: {
584: return(0);
585: }
587: /*
588: Adds diagonal pointers to sparse matrix structure.
589: */
592: PetscErrorCode MatMarkDiagonal_BlockMat(Mat A)
593: {
594: Mat_BlockMat *a = (Mat_BlockMat*)A->data;
596: PetscInt i,j,mbs = A->rmap->n/A->rmap->bs;
599: if (!a->diag) {
600: PetscMalloc1(mbs,&a->diag);
601: }
602: for (i=0; i<mbs; i++) {
603: a->diag[i] = a->i[i+1];
604: for (j=a->i[i]; j<a->i[i+1]; j++) {
605: if (a->j[j] == i) {
606: a->diag[i] = j;
607: break;
608: }
609: }
610: }
611: return(0);
612: }
616: PetscErrorCode MatGetSubMatrix_BlockMat(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B)
617: {
618: Mat_BlockMat *a = (Mat_BlockMat*)A->data;
619: Mat_SeqAIJ *c;
621: PetscInt i,k,first,step,lensi,nrows,ncols;
622: PetscInt *j_new,*i_new,*aj = a->j,*ailen = a->ilen;
623: PetscScalar *a_new;
624: Mat C,*aa = a->a;
625: PetscBool stride,equal;
628: ISEqual(isrow,iscol,&equal);
629: if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only for idential column and row indices");
630: PetscObjectTypeCompare((PetscObject)iscol,ISSTRIDE,&stride);
631: if (!stride) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only for stride indices");
632: ISStrideGetInfo(iscol,&first,&step);
633: if (step != A->rmap->bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Can only select one entry from each block");
635: ISGetLocalSize(isrow,&nrows);
636: ncols = nrows;
638: /* create submatrix */
639: if (scall == MAT_REUSE_MATRIX) {
640: PetscInt n_cols,n_rows;
641: C = *B;
642: MatGetSize(C,&n_rows,&n_cols);
643: if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size");
644: MatZeroEntries(C);
645: } else {
646: MatCreate(PetscObjectComm((PetscObject)A),&C);
647: MatSetSizes(C,nrows,ncols,PETSC_DETERMINE,PETSC_DETERMINE);
648: if (A->symmetric) {
649: MatSetType(C,MATSEQSBAIJ);
650: } else {
651: MatSetType(C,MATSEQAIJ);
652: }
653: MatSeqAIJSetPreallocation(C,0,ailen);
654: MatSeqSBAIJSetPreallocation(C,1,0,ailen);
655: }
656: c = (Mat_SeqAIJ*)C->data;
658: /* loop over rows inserting into submatrix */
659: a_new = c->a;
660: j_new = c->j;
661: i_new = c->i;
663: for (i=0; i<nrows; i++) {
664: lensi = ailen[i];
665: for (k=0; k<lensi; k++) {
666: *j_new++ = *aj++;
667: MatGetValue(*aa++,first,first,a_new++);
668: }
669: i_new[i+1] = i_new[i] + lensi;
670: c->ilen[i] = lensi;
671: }
673: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
674: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
675: *B = C;
676: return(0);
677: }
681: PetscErrorCode MatAssemblyEnd_BlockMat(Mat A,MatAssemblyType mode)
682: {
683: Mat_BlockMat *a = (Mat_BlockMat*)A->data;
685: PetscInt fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
686: PetscInt m = a->mbs,*ip,N,*ailen = a->ilen,rmax = 0;
687: Mat *aa = a->a,*ap;
690: if (mode == MAT_FLUSH_ASSEMBLY) return(0);
692: if (m) rmax = ailen[0]; /* determine row with most nonzeros */
693: for (i=1; i<m; i++) {
694: /* move each row back by the amount of empty slots (fshift) before it*/
695: fshift += imax[i-1] - ailen[i-1];
696: rmax = PetscMax(rmax,ailen[i]);
697: if (fshift) {
698: ip = aj + ai[i];
699: ap = aa + ai[i];
700: N = ailen[i];
701: for (j=0; j<N; j++) {
702: ip[j-fshift] = ip[j];
703: ap[j-fshift] = ap[j];
704: }
705: }
706: ai[i] = ai[i-1] + ailen[i-1];
707: }
708: if (m) {
709: fshift += imax[m-1] - ailen[m-1];
710: ai[m] = ai[m-1] + ailen[m-1];
711: }
712: /* reset ilen and imax for each row */
713: for (i=0; i<m; i++) {
714: ailen[i] = imax[i] = ai[i+1] - ai[i];
715: }
716: a->nz = ai[m];
717: for (i=0; i<a->nz; i++) {
718: #if defined(PETSC_USE_DEBUG)
719: if (!aa[i]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Null matrix at location %D column %D nz %D",i,aj[i],a->nz);
720: #endif
721: MatAssemblyBegin(aa[i],MAT_FINAL_ASSEMBLY);
722: MatAssemblyEnd(aa[i],MAT_FINAL_ASSEMBLY);
723: }
724: PetscInfo4(A,"Matrix size: %D X %D; storage space: %D unneeded,%D used\n",m,A->cmap->n/A->cmap->bs,fshift,a->nz);
725: PetscInfo1(A,"Number of mallocs during MatSetValues() is %D\n",a->reallocs);
726: PetscInfo1(A,"Maximum nonzeros in any row is %D\n",rmax);
728: A->info.mallocs += a->reallocs;
729: a->reallocs = 0;
730: A->info.nz_unneeded = (double)fshift;
731: a->rmax = rmax;
732: MatMarkDiagonal_BlockMat(A);
733: return(0);
734: }
738: PetscErrorCode MatSetOption_BlockMat(Mat A,MatOption opt,PetscBool flg)
739: {
741: if (opt == MAT_SYMMETRIC && flg) {
742: A->ops->sor = MatSOR_BlockMat_Symmetric;
743: A->ops->mult = MatMult_BlockMat_Symmetric;
744: } else {
745: PetscInfo1(A,"Unused matrix option %s\n",MatOptions[opt]);
746: }
747: return(0);
748: }
751: static struct _MatOps MatOps_Values = {MatSetValues_BlockMat,
752: 0,
753: 0,
754: MatMult_BlockMat,
755: /* 4*/ MatMultAdd_BlockMat,
756: MatMultTranspose_BlockMat,
757: MatMultTransposeAdd_BlockMat,
758: 0,
759: 0,
760: 0,
761: /* 10*/ 0,
762: 0,
763: 0,
764: MatSOR_BlockMat,
765: 0,
766: /* 15*/ 0,
767: 0,
768: 0,
769: 0,
770: 0,
771: /* 20*/ 0,
772: MatAssemblyEnd_BlockMat,
773: MatSetOption_BlockMat,
774: 0,
775: /* 24*/ 0,
776: 0,
777: 0,
778: 0,
779: 0,
780: /* 29*/ 0,
781: 0,
782: 0,
783: 0,
784: 0,
785: /* 34*/ 0,
786: 0,
787: 0,
788: 0,
789: 0,
790: /* 39*/ 0,
791: 0,
792: 0,
793: 0,
794: 0,
795: /* 44*/ 0,
796: 0,
797: 0,
798: 0,
799: 0,
800: /* 49*/ 0,
801: 0,
802: 0,
803: 0,
804: 0,
805: /* 54*/ 0,
806: 0,
807: 0,
808: 0,
809: 0,
810: /* 59*/ MatGetSubMatrix_BlockMat,
811: MatDestroy_BlockMat,
812: MatView_BlockMat,
813: 0,
814: 0,
815: /* 64*/ 0,
816: 0,
817: 0,
818: 0,
819: 0,
820: /* 69*/ 0,
821: 0,
822: 0,
823: 0,
824: 0,
825: /* 74*/ 0,
826: 0,
827: 0,
828: 0,
829: 0,
830: /* 79*/ 0,
831: 0,
832: 0,
833: 0,
834: MatLoad_BlockMat,
835: /* 84*/ 0,
836: 0,
837: 0,
838: 0,
839: 0,
840: /* 89*/ 0,
841: 0,
842: 0,
843: 0,
844: 0,
845: /* 94*/ 0,
846: 0,
847: 0,
848: 0,
849: 0,
850: /* 99*/ 0,
851: 0,
852: 0,
853: 0,
854: 0,
855: /*104*/ 0,
856: 0,
857: 0,
858: 0,
859: 0,
860: /*109*/ 0,
861: 0,
862: 0,
863: 0,
864: 0,
865: /*114*/ 0,
866: 0,
867: 0,
868: 0,
869: 0,
870: /*119*/ 0,
871: 0,
872: 0,
873: 0,
874: 0,
875: /*124*/ 0,
876: 0,
877: 0,
878: 0,
879: 0,
880: /*129*/ 0,
881: 0,
882: 0,
883: 0,
884: 0,
885: /*134*/ 0,
886: 0,
887: 0,
888: 0,
889: 0,
890: /*139*/ 0,
891: 0,
892: 0
893: };
897: /*@C
898: MatBlockMatSetPreallocation - For good matrix assembly performance
899: the user should preallocate the matrix storage by setting the parameter nz
900: (or the array nnz). By setting these parameters accurately, performance
901: during matrix assembly can be increased by more than a factor of 50.
903: Collective on MPI_Comm
905: Input Parameters:
906: + B - The matrix
907: . bs - size of each block in matrix
908: . nz - number of nonzeros per block row (same for all rows)
909: - nnz - array containing the number of nonzeros in the various block rows
910: (possibly different for each row) or NULL
912: Notes:
913: If nnz is given then nz is ignored
915: Specify the preallocated storage with either nz or nnz (not both).
916: Set nz=PETSC_DEFAULT and nnz=NULL for PETSc to control dynamic memory
917: allocation. For large problems you MUST preallocate memory or you
918: will get TERRIBLE performance, see the users' manual chapter on matrices.
920: Level: intermediate
922: .seealso: MatCreate(), MatCreateBlockMat(), MatSetValues()
924: @*/
925: PetscErrorCode MatBlockMatSetPreallocation(Mat B,PetscInt bs,PetscInt nz,const PetscInt nnz[])
926: {
930: PetscTryMethod(B,"MatBlockMatSetPreallocation_C",(Mat,PetscInt,PetscInt,const PetscInt[]),(B,bs,nz,nnz));
931: return(0);
932: }
936: PetscErrorCode MatBlockMatSetPreallocation_BlockMat(Mat A,PetscInt bs,PetscInt nz,PetscInt *nnz)
937: {
938: Mat_BlockMat *bmat = (Mat_BlockMat*)A->data;
940: PetscInt i;
943: PetscLayoutSetBlockSize(A->rmap,bs);
944: PetscLayoutSetBlockSize(A->cmap,bs);
945: PetscLayoutSetUp(A->rmap);
946: PetscLayoutSetUp(A->cmap);
947: PetscLayoutGetBlockSize(A->rmap,&bs);
949: if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
950: if (nz < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz);
951: if (nnz) {
952: for (i=0; i<A->rmap->n/bs; i++) {
953: 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]);
954: if (nnz[i] > A->cmap->n/bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be greater than row length: local row %d value %d rowlength %d",i,nnz[i],A->cmap->n/bs);
955: }
956: }
957: bmat->mbs = A->rmap->n/bs;
959: VecCreateSeqWithArray(PETSC_COMM_SELF,1,bs,NULL,&bmat->right);
960: VecCreateSeqWithArray(PETSC_COMM_SELF,1,bs,NULL,&bmat->middle);
961: VecCreateSeq(PETSC_COMM_SELF,bs,&bmat->left);
963: if (!bmat->imax) {
964: PetscMalloc2(A->rmap->n,&bmat->imax,A->rmap->n,&bmat->ilen);
965: PetscLogObjectMemory((PetscObject)A,2*A->rmap->n*sizeof(PetscInt));
966: }
967: if (nnz) {
968: nz = 0;
969: for (i=0; i<A->rmap->n/A->rmap->bs; i++) {
970: bmat->imax[i] = nnz[i];
971: nz += nnz[i];
972: }
973: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Currently requires block row by row preallocation");
975: /* bmat->ilen will count nonzeros in each row so far. */
976: for (i=0; i<bmat->mbs; i++) bmat->ilen[i] = 0;
978: /* allocate the matrix space */
979: MatSeqXAIJFreeAIJ(A,(PetscScalar**)&bmat->a,&bmat->j,&bmat->i);
980: PetscMalloc3(nz,&bmat->a,nz,&bmat->j,A->rmap->n+1,&bmat->i);
981: PetscLogObjectMemory((PetscObject)A,(A->rmap->n+1)*sizeof(PetscInt)+nz*(sizeof(PetscScalar)+sizeof(PetscInt)));
982: bmat->i[0] = 0;
983: for (i=1; i<bmat->mbs+1; i++) {
984: bmat->i[i] = bmat->i[i-1] + bmat->imax[i-1];
985: }
986: bmat->singlemalloc = PETSC_TRUE;
987: bmat->free_a = PETSC_TRUE;
988: bmat->free_ij = PETSC_TRUE;
990: bmat->nz = 0;
991: bmat->maxnz = nz;
992: A->info.nz_unneeded = (double)bmat->maxnz;
993: MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);
994: return(0);
995: }
997: /*MC
998: MATBLOCKMAT - A matrix that is defined by a set of Mat's that represents a sparse block matrix
999: consisting of (usually) sparse blocks.
1001: Level: advanced
1003: .seealso: MatCreateBlockMat()
1005: M*/
1009: PETSC_EXTERN PetscErrorCode MatCreate_BlockMat(Mat A)
1010: {
1011: Mat_BlockMat *b;
1015: PetscNewLog(A,&b);
1016: A->data = (void*)b;
1017: PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));
1019: A->assembled = PETSC_TRUE;
1020: A->preallocated = PETSC_FALSE;
1021: PetscObjectChangeTypeName((PetscObject)A,MATBLOCKMAT);
1023: PetscObjectComposeFunction((PetscObject)A,"MatBlockMatSetPreallocation_C",MatBlockMatSetPreallocation_BlockMat);
1024: return(0);
1025: }
1029: /*@C
1030: MatCreateBlockMat - Creates a new matrix in which each block contains a uniform-size sequential Mat object
1032: Collective on MPI_Comm
1034: Input Parameters:
1035: + comm - MPI communicator
1036: . m - number of rows
1037: . n - number of columns
1038: . bs - size of each submatrix
1039: . nz - expected maximum number of nonzero blocks in row (use PETSC_DEFAULT if not known)
1040: - nnz - expected number of nonzers per block row if known (use NULL otherwise)
1043: Output Parameter:
1044: . A - the matrix
1046: Level: intermediate
1048: Notes: Matrices of this type are nominally-sparse matrices in which each "entry" is a Mat object. Each Mat must
1049: have the same size and be sequential. The local and global sizes must be compatible with this decomposition.
1051: For matrices containing parallel submatrices and variable block sizes, see MATNEST.
1053: .keywords: matrix, bmat, create
1055: .seealso: MATBLOCKMAT, MatCreateNest()
1056: @*/
1057: PetscErrorCode MatCreateBlockMat(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt bs,PetscInt nz,PetscInt *nnz, Mat *A)
1058: {
1062: MatCreate(comm,A);
1063: MatSetSizes(*A,m,n,PETSC_DETERMINE,PETSC_DETERMINE);
1064: MatSetType(*A,MATBLOCKMAT);
1065: MatBlockMatSetPreallocation(*A,bs,nz,nnz);
1066: return(0);
1067: }