Actual source code: mpisbaij.c
petsc-3.7.3 2016-08-01
2: #include <../src/mat/impls/baij/mpi/mpibaij.h> /*I "petscmat.h" I*/
3: #include <../src/mat/impls/sbaij/mpi/mpisbaij.h>
4: #include <../src/mat/impls/sbaij/seq/sbaij.h>
5: #include <petscblaslapack.h>
7: #if defined(PETSC_HAVE_ELEMENTAL)
8: PETSC_INTERN PetscErrorCode MatConvert_MPISBAIJ_Elemental(Mat,MatType,MatReuse,Mat*);
9: #endif
12: PetscErrorCode MatStoreValues_MPISBAIJ(Mat mat)
13: {
14: Mat_MPISBAIJ *aij = (Mat_MPISBAIJ*)mat->data;
18: MatStoreValues(aij->A);
19: MatStoreValues(aij->B);
20: return(0);
21: }
25: PetscErrorCode MatRetrieveValues_MPISBAIJ(Mat mat)
26: {
27: Mat_MPISBAIJ *aij = (Mat_MPISBAIJ*)mat->data;
31: MatRetrieveValues(aij->A);
32: MatRetrieveValues(aij->B);
33: return(0);
34: }
36: #define MatSetValues_SeqSBAIJ_A_Private(row,col,value,addv,orow,ocol) \
37: { \
38: \
39: brow = row/bs; \
40: rp = aj + ai[brow]; ap = aa + bs2*ai[brow]; \
41: rmax = aimax[brow]; nrow = ailen[brow]; \
42: bcol = col/bs; \
43: ridx = row % bs; cidx = col % bs; \
44: low = 0; high = nrow; \
45: while (high-low > 3) { \
46: t = (low+high)/2; \
47: if (rp[t] > bcol) high = t; \
48: else low = t; \
49: } \
50: for (_i=low; _i<high; _i++) { \
51: if (rp[_i] > bcol) break; \
52: if (rp[_i] == bcol) { \
53: bap = ap + bs2*_i + bs*cidx + ridx; \
54: if (addv == ADD_VALUES) *bap += value; \
55: else *bap = value; \
56: goto a_noinsert; \
57: } \
58: } \
59: if (a->nonew == 1) goto a_noinsert; \
60: if (a->nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero at global row/column (%D, %D) into matrix", orow, ocol); \
61: MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,aimax,a->nonew,MatScalar); \
62: N = nrow++ - 1; \
63: /* shift up all the later entries in this row */ \
64: for (ii=N; ii>=_i; ii--) { \
65: rp[ii+1] = rp[ii]; \
66: PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar)); \
67: } \
68: if (N>=_i) { PetscMemzero(ap+bs2*_i,bs2*sizeof(MatScalar)); } \
69: rp[_i] = bcol; \
70: ap[bs2*_i + bs*cidx + ridx] = value; \
71: A->nonzerostate++;\
72: a_noinsert:; \
73: ailen[brow] = nrow; \
74: }
76: #define MatSetValues_SeqSBAIJ_B_Private(row,col,value,addv,orow,ocol) \
77: { \
78: brow = row/bs; \
79: rp = bj + bi[brow]; ap = ba + bs2*bi[brow]; \
80: rmax = bimax[brow]; nrow = bilen[brow]; \
81: bcol = col/bs; \
82: ridx = row % bs; cidx = col % bs; \
83: low = 0; high = nrow; \
84: while (high-low > 3) { \
85: t = (low+high)/2; \
86: if (rp[t] > bcol) high = t; \
87: else low = t; \
88: } \
89: for (_i=low; _i<high; _i++) { \
90: if (rp[_i] > bcol) break; \
91: if (rp[_i] == bcol) { \
92: bap = ap + bs2*_i + bs*cidx + ridx; \
93: if (addv == ADD_VALUES) *bap += value; \
94: else *bap = value; \
95: goto b_noinsert; \
96: } \
97: } \
98: if (b->nonew == 1) goto b_noinsert; \
99: if (b->nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero at global row/column (%D, %D) into matrix", orow, ocol); \
100: MatSeqXAIJReallocateAIJ(B,b->mbs,bs2,nrow,brow,bcol,rmax,ba,bi,bj,rp,ap,bimax,b->nonew,MatScalar); \
101: N = nrow++ - 1; \
102: /* shift up all the later entries in this row */ \
103: for (ii=N; ii>=_i; ii--) { \
104: rp[ii+1] = rp[ii]; \
105: PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar)); \
106: } \
107: if (N>=_i) { PetscMemzero(ap+bs2*_i,bs2*sizeof(MatScalar));} \
108: rp[_i] = bcol; \
109: ap[bs2*_i + bs*cidx + ridx] = value; \
110: B->nonzerostate++;\
111: b_noinsert:; \
112: bilen[brow] = nrow; \
113: }
115: /* Only add/insert a(i,j) with i<=j (blocks).
116: Any a(i,j) with i>j input by user is ingored.
117: */
120: PetscErrorCode MatSetValues_MPISBAIJ(Mat mat,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode addv)
121: {
122: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
123: MatScalar value;
124: PetscBool roworiented = baij->roworiented;
126: PetscInt i,j,row,col;
127: PetscInt rstart_orig=mat->rmap->rstart;
128: PetscInt rend_orig =mat->rmap->rend,cstart_orig=mat->cmap->rstart;
129: PetscInt cend_orig =mat->cmap->rend,bs=mat->rmap->bs;
131: /* Some Variables required in the macro */
132: Mat A = baij->A;
133: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)(A)->data;
134: PetscInt *aimax=a->imax,*ai=a->i,*ailen=a->ilen,*aj=a->j;
135: MatScalar *aa =a->a;
137: Mat B = baij->B;
138: Mat_SeqBAIJ *b = (Mat_SeqBAIJ*)(B)->data;
139: PetscInt *bimax=b->imax,*bi=b->i,*bilen=b->ilen,*bj=b->j;
140: MatScalar *ba =b->a;
142: PetscInt *rp,ii,nrow,_i,rmax,N,brow,bcol;
143: PetscInt low,high,t,ridx,cidx,bs2=a->bs2;
144: MatScalar *ap,*bap;
146: /* for stash */
147: PetscInt n_loc, *in_loc = NULL;
148: MatScalar *v_loc = NULL;
151: if (!baij->donotstash) {
152: if (n > baij->n_loc) {
153: PetscFree(baij->in_loc);
154: PetscFree(baij->v_loc);
155: PetscMalloc1(n,&baij->in_loc);
156: PetscMalloc1(n,&baij->v_loc);
158: baij->n_loc = n;
159: }
160: in_loc = baij->in_loc;
161: v_loc = baij->v_loc;
162: }
164: for (i=0; i<m; i++) {
165: if (im[i] < 0) continue;
166: #if defined(PETSC_USE_DEBUG)
167: if (im[i] >= mat->rmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",im[i],mat->rmap->N-1);
168: #endif
169: if (im[i] >= rstart_orig && im[i] < rend_orig) { /* this processor entry */
170: row = im[i] - rstart_orig; /* local row index */
171: for (j=0; j<n; j++) {
172: if (im[i]/bs > in[j]/bs) {
173: if (a->ignore_ltriangular) {
174: continue; /* ignore lower triangular blocks */
175: } 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)");
176: }
177: if (in[j] >= cstart_orig && in[j] < cend_orig) { /* diag entry (A) */
178: col = in[j] - cstart_orig; /* local col index */
179: brow = row/bs; bcol = col/bs;
180: if (brow > bcol) continue; /* ignore lower triangular blocks of A */
181: if (roworiented) value = v[i*n+j];
182: else value = v[i+j*m];
183: MatSetValues_SeqSBAIJ_A_Private(row,col,value,addv,im[i],in[j]);
184: /* MatSetValues_SeqBAIJ(baij->A,1,&row,1,&col,&value,addv); */
185: } else if (in[j] < 0) continue;
186: #if defined(PETSC_USE_DEBUG)
187: else if (in[j] >= mat->cmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[j],mat->cmap->N-1);
188: #endif
189: else { /* off-diag entry (B) */
190: if (mat->was_assembled) {
191: if (!baij->colmap) {
192: MatCreateColmap_MPIBAIJ_Private(mat);
193: }
194: #if defined(PETSC_USE_CTABLE)
195: PetscTableFind(baij->colmap,in[j]/bs + 1,&col);
196: col = col - 1;
197: #else
198: col = baij->colmap[in[j]/bs] - 1;
199: #endif
200: if (col < 0 && !((Mat_SeqSBAIJ*)(baij->A->data))->nonew) {
201: MatDisAssemble_MPISBAIJ(mat);
202: col = in[j];
203: /* Reinitialize the variables required by MatSetValues_SeqBAIJ_B_Private() */
204: B = baij->B;
205: b = (Mat_SeqBAIJ*)(B)->data;
206: bimax= b->imax;bi=b->i;bilen=b->ilen;bj=b->j;
207: ba = b->a;
208: } else col += in[j]%bs;
209: } else col = in[j];
210: if (roworiented) value = v[i*n+j];
211: else value = v[i+j*m];
212: MatSetValues_SeqSBAIJ_B_Private(row,col,value,addv,im[i],in[j]);
213: /* MatSetValues_SeqBAIJ(baij->B,1,&row,1,&col,&value,addv); */
214: }
215: }
216: } else { /* off processor entry */
217: if (mat->nooffprocentries) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Setting off process row %D even though MatSetOption(,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE) was set",im[i]);
218: if (!baij->donotstash) {
219: mat->assembled = PETSC_FALSE;
220: n_loc = 0;
221: for (j=0; j<n; j++) {
222: if (im[i]/bs > in[j]/bs) continue; /* ignore lower triangular blocks */
223: in_loc[n_loc] = in[j];
224: if (roworiented) {
225: v_loc[n_loc] = v[i*n+j];
226: } else {
227: v_loc[n_loc] = v[j*m+i];
228: }
229: n_loc++;
230: }
231: MatStashValuesRow_Private(&mat->stash,im[i],n_loc,in_loc,v_loc,PETSC_FALSE);
232: }
233: }
234: }
235: return(0);
236: }
240: PETSC_STATIC_INLINE PetscErrorCode MatSetValuesBlocked_SeqSBAIJ_Inlined(Mat A,PetscInt row,PetscInt col,const PetscScalar v[],InsertMode is,PetscInt orow,PetscInt ocol)
241: {
242: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
243: PetscErrorCode ierr;
244: PetscInt *rp,low,high,t,ii,jj,nrow,i,rmax,N;
245: PetscInt *imax =a->imax,*ai=a->i,*ailen=a->ilen;
246: PetscInt *aj =a->j,nonew=a->nonew,bs2=a->bs2,bs=A->rmap->bs;
247: PetscBool roworiented=a->roworiented;
248: const PetscScalar *value = v;
249: MatScalar *ap,*aa = a->a,*bap;
252: if (col < row) {
253: if (a->ignore_ltriangular) return(0); /* ignore lower triangular block */
254: 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)");
255: }
256: rp = aj + ai[row];
257: ap = aa + bs2*ai[row];
258: rmax = imax[row];
259: nrow = ailen[row];
260: value = v;
261: low = 0;
262: high = nrow;
264: while (high-low > 7) {
265: t = (low+high)/2;
266: if (rp[t] > col) high = t;
267: else low = t;
268: }
269: for (i=low; i<high; i++) {
270: if (rp[i] > col) break;
271: if (rp[i] == col) {
272: bap = ap + bs2*i;
273: if (roworiented) {
274: if (is == ADD_VALUES) {
275: for (ii=0; ii<bs; ii++) {
276: for (jj=ii; jj<bs2; jj+=bs) {
277: bap[jj] += *value++;
278: }
279: }
280: } else {
281: for (ii=0; ii<bs; ii++) {
282: for (jj=ii; jj<bs2; jj+=bs) {
283: bap[jj] = *value++;
284: }
285: }
286: }
287: } else {
288: if (is == ADD_VALUES) {
289: for (ii=0; ii<bs; ii++) {
290: for (jj=0; jj<bs; jj++) {
291: *bap++ += *value++;
292: }
293: }
294: } else {
295: for (ii=0; ii<bs; ii++) {
296: for (jj=0; jj<bs; jj++) {
297: *bap++ = *value++;
298: }
299: }
300: }
301: }
302: goto noinsert2;
303: }
304: }
305: if (nonew == 1) goto noinsert2;
306: if (nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new block index nonzero block (%D, %D) in the matrix", orow, ocol);
307: MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar);
308: N = nrow++ - 1; high++;
309: /* shift up all the later entries in this row */
310: for (ii=N; ii>=i; ii--) {
311: rp[ii+1] = rp[ii];
312: PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));
313: }
314: if (N >= i) {
315: PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));
316: }
317: rp[i] = col;
318: bap = ap + bs2*i;
319: if (roworiented) {
320: for (ii=0; ii<bs; ii++) {
321: for (jj=ii; jj<bs2; jj+=bs) {
322: bap[jj] = *value++;
323: }
324: }
325: } else {
326: for (ii=0; ii<bs; ii++) {
327: for (jj=0; jj<bs; jj++) {
328: *bap++ = *value++;
329: }
330: }
331: }
332: noinsert2:;
333: ailen[row] = nrow;
334: return(0);
335: }
339: /*
340: This routine is exactly duplicated in mpibaij.c
341: */
342: PETSC_STATIC_INLINE PetscErrorCode MatSetValuesBlocked_SeqBAIJ_Inlined(Mat A,PetscInt row,PetscInt col,const PetscScalar v[],InsertMode is,PetscInt orow,PetscInt ocol)
343: {
344: Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
345: PetscInt *rp,low,high,t,ii,jj,nrow,i,rmax,N;
346: PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen;
347: PetscErrorCode ierr;
348: PetscInt *aj =a->j,nonew=a->nonew,bs2=a->bs2,bs=A->rmap->bs;
349: PetscBool roworiented=a->roworiented;
350: const PetscScalar *value = v;
351: MatScalar *ap,*aa = a->a,*bap;
354: rp = aj + ai[row];
355: ap = aa + bs2*ai[row];
356: rmax = imax[row];
357: nrow = ailen[row];
358: low = 0;
359: high = nrow;
360: value = v;
361: while (high-low > 7) {
362: t = (low+high)/2;
363: if (rp[t] > col) high = t;
364: else low = t;
365: }
366: for (i=low; i<high; i++) {
367: if (rp[i] > col) break;
368: if (rp[i] == col) {
369: bap = ap + bs2*i;
370: if (roworiented) {
371: if (is == ADD_VALUES) {
372: for (ii=0; ii<bs; ii++) {
373: for (jj=ii; jj<bs2; jj+=bs) {
374: bap[jj] += *value++;
375: }
376: }
377: } else {
378: for (ii=0; ii<bs; ii++) {
379: for (jj=ii; jj<bs2; jj+=bs) {
380: bap[jj] = *value++;
381: }
382: }
383: }
384: } else {
385: if (is == ADD_VALUES) {
386: for (ii=0; ii<bs; ii++,value+=bs) {
387: for (jj=0; jj<bs; jj++) {
388: bap[jj] += value[jj];
389: }
390: bap += bs;
391: }
392: } else {
393: for (ii=0; ii<bs; ii++,value+=bs) {
394: for (jj=0; jj<bs; jj++) {
395: bap[jj] = value[jj];
396: }
397: bap += bs;
398: }
399: }
400: }
401: goto noinsert2;
402: }
403: }
404: if (nonew == 1) goto noinsert2;
405: if (nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new global block indexed nonzero block (%D, %D) in the matrix", orow, ocol);
406: MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar);
407: N = nrow++ - 1; high++;
408: /* shift up all the later entries in this row */
409: for (ii=N; ii>=i; ii--) {
410: rp[ii+1] = rp[ii];
411: PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));
412: }
413: if (N >= i) {
414: PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));
415: }
416: rp[i] = col;
417: bap = ap + bs2*i;
418: if (roworiented) {
419: for (ii=0; ii<bs; ii++) {
420: for (jj=ii; jj<bs2; jj+=bs) {
421: bap[jj] = *value++;
422: }
423: }
424: } else {
425: for (ii=0; ii<bs; ii++) {
426: for (jj=0; jj<bs; jj++) {
427: *bap++ = *value++;
428: }
429: }
430: }
431: noinsert2:;
432: ailen[row] = nrow;
433: return(0);
434: }
438: /*
439: This routine could be optimized by removing the need for the block copy below and passing stride information
440: to the above inline routines; similarly in MatSetValuesBlocked_MPIBAIJ()
441: */
442: PetscErrorCode MatSetValuesBlocked_MPISBAIJ(Mat mat,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const MatScalar v[],InsertMode addv)
443: {
444: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
445: const MatScalar *value;
446: MatScalar *barray =baij->barray;
447: PetscBool roworiented = baij->roworiented,ignore_ltriangular = ((Mat_SeqSBAIJ*)baij->A->data)->ignore_ltriangular;
448: PetscErrorCode ierr;
449: PetscInt i,j,ii,jj,row,col,rstart=baij->rstartbs;
450: PetscInt rend=baij->rendbs,cstart=baij->rstartbs,stepval;
451: PetscInt cend=baij->rendbs,bs=mat->rmap->bs,bs2=baij->bs2;
454: if (!barray) {
455: PetscMalloc1(bs2,&barray);
456: baij->barray = barray;
457: }
459: if (roworiented) {
460: stepval = (n-1)*bs;
461: } else {
462: stepval = (m-1)*bs;
463: }
464: for (i=0; i<m; i++) {
465: if (im[i] < 0) continue;
466: #if defined(PETSC_USE_DEBUG)
467: if (im[i] >= baij->Mbs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Block indexed row too large %D max %D",im[i],baij->Mbs-1);
468: #endif
469: if (im[i] >= rstart && im[i] < rend) {
470: row = im[i] - rstart;
471: for (j=0; j<n; j++) {
472: if (im[i] > in[j]) {
473: if (ignore_ltriangular) continue; /* ignore lower triangular blocks */
474: 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)");
475: }
476: /* If NumCol = 1 then a copy is not required */
477: if ((roworiented) && (n == 1)) {
478: barray = (MatScalar*) v + i*bs2;
479: } else if ((!roworiented) && (m == 1)) {
480: barray = (MatScalar*) v + j*bs2;
481: } else { /* Here a copy is required */
482: if (roworiented) {
483: value = v + i*(stepval+bs)*bs + j*bs;
484: } else {
485: value = v + j*(stepval+bs)*bs + i*bs;
486: }
487: for (ii=0; ii<bs; ii++,value+=stepval) {
488: for (jj=0; jj<bs; jj++) {
489: *barray++ = *value++;
490: }
491: }
492: barray -=bs2;
493: }
495: if (in[j] >= cstart && in[j] < cend) {
496: col = in[j] - cstart;
497: MatSetValuesBlocked_SeqSBAIJ_Inlined(baij->A,row,col,barray,addv,im[i],in[j]);
498: } else if (in[j] < 0) continue;
499: #if defined(PETSC_USE_DEBUG)
500: else if (in[j] >= baij->Nbs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Block indexed column too large %D max %D",in[j],baij->Nbs-1);
501: #endif
502: else {
503: if (mat->was_assembled) {
504: if (!baij->colmap) {
505: MatCreateColmap_MPIBAIJ_Private(mat);
506: }
508: #if defined(PETSC_USE_DEBUG)
509: #if defined(PETSC_USE_CTABLE)
510: { PetscInt data;
511: PetscTableFind(baij->colmap,in[j]+1,&data);
512: if ((data - 1) % bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect colmap");
513: }
514: #else
515: if ((baij->colmap[in[j]] - 1) % bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect colmap");
516: #endif
517: #endif
518: #if defined(PETSC_USE_CTABLE)
519: PetscTableFind(baij->colmap,in[j]+1,&col);
520: col = (col - 1)/bs;
521: #else
522: col = (baij->colmap[in[j]] - 1)/bs;
523: #endif
524: if (col < 0 && !((Mat_SeqBAIJ*)(baij->A->data))->nonew) {
525: MatDisAssemble_MPISBAIJ(mat);
526: col = in[j];
527: }
528: } else col = in[j];
529: MatSetValuesBlocked_SeqBAIJ_Inlined(baij->B,row,col,barray,addv,im[i],in[j]);
530: }
531: }
532: } else {
533: if (mat->nooffprocentries) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Setting off process block indexed row %D even though MatSetOption(,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE) was set",im[i]);
534: if (!baij->donotstash) {
535: if (roworiented) {
536: MatStashValuesRowBlocked_Private(&mat->bstash,im[i],n,in,v,m,n,i);
537: } else {
538: MatStashValuesColBlocked_Private(&mat->bstash,im[i],n,in,v,m,n,i);
539: }
540: }
541: }
542: }
543: return(0);
544: }
548: PetscErrorCode MatGetValues_MPISBAIJ(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[])
549: {
550: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
552: PetscInt bs = mat->rmap->bs,i,j,bsrstart = mat->rmap->rstart,bsrend = mat->rmap->rend;
553: PetscInt bscstart = mat->cmap->rstart,bscend = mat->cmap->rend,row,col,data;
556: for (i=0; i<m; i++) {
557: if (idxm[i] < 0) continue; /* SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",idxm[i]); */
558: if (idxm[i] >= mat->rmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm[i],mat->rmap->N-1);
559: if (idxm[i] >= bsrstart && idxm[i] < bsrend) {
560: row = idxm[i] - bsrstart;
561: for (j=0; j<n; j++) {
562: if (idxn[j] < 0) continue; /* SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column %D",idxn[j]); */
563: if (idxn[j] >= mat->cmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",idxn[j],mat->cmap->N-1);
564: if (idxn[j] >= bscstart && idxn[j] < bscend) {
565: col = idxn[j] - bscstart;
566: MatGetValues_SeqSBAIJ(baij->A,1,&row,1,&col,v+i*n+j);
567: } else {
568: if (!baij->colmap) {
569: MatCreateColmap_MPIBAIJ_Private(mat);
570: }
571: #if defined(PETSC_USE_CTABLE)
572: PetscTableFind(baij->colmap,idxn[j]/bs+1,&data);
573: data--;
574: #else
575: data = baij->colmap[idxn[j]/bs]-1;
576: #endif
577: if ((data < 0) || (baij->garray[data/bs] != idxn[j]/bs)) *(v+i*n+j) = 0.0;
578: else {
579: col = data + idxn[j]%bs;
580: MatGetValues_SeqBAIJ(baij->B,1,&row,1,&col,v+i*n+j);
581: }
582: }
583: }
584: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only local values currently supported");
585: }
586: return(0);
587: }
591: PetscErrorCode MatNorm_MPISBAIJ(Mat mat,NormType type,PetscReal *norm)
592: {
593: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
595: PetscReal sum[2],*lnorm2;
598: if (baij->size == 1) {
599: MatNorm(baij->A,type,norm);
600: } else {
601: if (type == NORM_FROBENIUS) {
602: PetscMalloc1(2,&lnorm2);
603: MatNorm(baij->A,type,lnorm2);
604: *lnorm2 = (*lnorm2)*(*lnorm2); lnorm2++; /* squar power of norm(A) */
605: MatNorm(baij->B,type,lnorm2);
606: *lnorm2 = (*lnorm2)*(*lnorm2); lnorm2--; /* squar power of norm(B) */
607: MPIU_Allreduce(lnorm2,sum,2,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)mat));
608: *norm = PetscSqrtReal(sum[0] + 2*sum[1]);
609: PetscFree(lnorm2);
610: } else if (type == NORM_INFINITY || type == NORM_1) { /* max row/column sum */
611: Mat_SeqSBAIJ *amat=(Mat_SeqSBAIJ*)baij->A->data;
612: Mat_SeqBAIJ *bmat=(Mat_SeqBAIJ*)baij->B->data;
613: PetscReal *rsum,*rsum2,vabs;
614: PetscInt *jj,*garray=baij->garray,rstart=baij->rstartbs,nz;
615: PetscInt brow,bcol,col,bs=baij->A->rmap->bs,row,grow,gcol,mbs=amat->mbs;
616: MatScalar *v;
618: PetscMalloc2(mat->cmap->N,&rsum,mat->cmap->N,&rsum2);
619: PetscMemzero(rsum,mat->cmap->N*sizeof(PetscReal));
620: /* Amat */
621: v = amat->a; jj = amat->j;
622: for (brow=0; brow<mbs; brow++) {
623: grow = bs*(rstart + brow);
624: nz = amat->i[brow+1] - amat->i[brow];
625: for (bcol=0; bcol<nz; bcol++) {
626: gcol = bs*(rstart + *jj); jj++;
627: for (col=0; col<bs; col++) {
628: for (row=0; row<bs; row++) {
629: vabs = PetscAbsScalar(*v); v++;
630: rsum[gcol+col] += vabs;
631: /* non-diagonal block */
632: if (bcol > 0 && vabs > 0.0) rsum[grow+row] += vabs;
633: }
634: }
635: }
636: PetscLogFlops(nz*bs*bs);
637: }
638: /* Bmat */
639: v = bmat->a; jj = bmat->j;
640: for (brow=0; brow<mbs; brow++) {
641: grow = bs*(rstart + brow);
642: nz = bmat->i[brow+1] - bmat->i[brow];
643: for (bcol=0; bcol<nz; bcol++) {
644: gcol = bs*garray[*jj]; jj++;
645: for (col=0; col<bs; col++) {
646: for (row=0; row<bs; row++) {
647: vabs = PetscAbsScalar(*v); v++;
648: rsum[gcol+col] += vabs;
649: rsum[grow+row] += vabs;
650: }
651: }
652: }
653: PetscLogFlops(nz*bs*bs);
654: }
655: MPIU_Allreduce(rsum,rsum2,mat->cmap->N,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)mat));
656: *norm = 0.0;
657: for (col=0; col<mat->cmap->N; col++) {
658: if (rsum2[col] > *norm) *norm = rsum2[col];
659: }
660: PetscFree2(rsum,rsum2);
661: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for this norm yet");
662: }
663: return(0);
664: }
668: PetscErrorCode MatAssemblyBegin_MPISBAIJ(Mat mat,MatAssemblyType mode)
669: {
670: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
672: PetscInt nstash,reallocs;
675: if (baij->donotstash || mat->nooffprocentries) return(0);
677: MatStashScatterBegin_Private(mat,&mat->stash,mat->rmap->range);
678: MatStashScatterBegin_Private(mat,&mat->bstash,baij->rangebs);
679: MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);
680: PetscInfo2(mat,"Stash has %D entries,uses %D mallocs.\n",nstash,reallocs);
681: MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);
682: PetscInfo2(mat,"Block-Stash has %D entries, uses %D mallocs.\n",nstash,reallocs);
683: return(0);
684: }
688: PetscErrorCode MatAssemblyEnd_MPISBAIJ(Mat mat,MatAssemblyType mode)
689: {
690: Mat_MPISBAIJ *baij=(Mat_MPISBAIJ*)mat->data;
691: Mat_SeqSBAIJ *a =(Mat_SeqSBAIJ*)baij->A->data;
693: PetscInt i,j,rstart,ncols,flg,bs2=baij->bs2;
694: PetscInt *row,*col;
695: PetscBool other_disassembled;
696: PetscMPIInt n;
697: PetscBool r1,r2,r3;
698: MatScalar *val;
700: /* do not use 'b=(Mat_SeqBAIJ*)baij->B->data' as B can be reset in disassembly */
702: if (!baij->donotstash && !mat->nooffprocentries) {
703: while (1) {
704: MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);
705: if (!flg) break;
707: for (i=0; i<n;) {
708: /* Now identify the consecutive vals belonging to the same row */
709: for (j=i,rstart=row[j]; j<n; j++) {
710: if (row[j] != rstart) break;
711: }
712: if (j < n) ncols = j-i;
713: else ncols = n-i;
714: /* Now assemble all these values with a single function call */
715: MatSetValues_MPISBAIJ(mat,1,row+i,ncols,col+i,val+i,mat->insertmode);
716: i = j;
717: }
718: }
719: MatStashScatterEnd_Private(&mat->stash);
720: /* Now process the block-stash. Since the values are stashed column-oriented,
721: set the roworiented flag to column oriented, and after MatSetValues()
722: restore the original flags */
723: r1 = baij->roworiented;
724: r2 = a->roworiented;
725: r3 = ((Mat_SeqBAIJ*)baij->B->data)->roworiented;
727: baij->roworiented = PETSC_FALSE;
728: a->roworiented = PETSC_FALSE;
730: ((Mat_SeqBAIJ*)baij->B->data)->roworiented = PETSC_FALSE; /* b->roworinted */
731: while (1) {
732: MatStashScatterGetMesg_Private(&mat->bstash,&n,&row,&col,&val,&flg);
733: if (!flg) break;
735: for (i=0; i<n;) {
736: /* Now identify the consecutive vals belonging to the same row */
737: for (j=i,rstart=row[j]; j<n; j++) {
738: if (row[j] != rstart) break;
739: }
740: if (j < n) ncols = j-i;
741: else ncols = n-i;
742: MatSetValuesBlocked_MPISBAIJ(mat,1,row+i,ncols,col+i,val+i*bs2,mat->insertmode);
743: i = j;
744: }
745: }
746: MatStashScatterEnd_Private(&mat->bstash);
748: baij->roworiented = r1;
749: a->roworiented = r2;
751: ((Mat_SeqBAIJ*)baij->B->data)->roworiented = r3; /* b->roworinted */
752: }
754: MatAssemblyBegin(baij->A,mode);
755: MatAssemblyEnd(baij->A,mode);
757: /* determine if any processor has disassembled, if so we must
758: also disassemble ourselfs, in order that we may reassemble. */
759: /*
760: if nonzero structure of submatrix B cannot change then we know that
761: no processor disassembled thus we can skip this stuff
762: */
763: if (!((Mat_SeqBAIJ*)baij->B->data)->nonew) {
764: MPIU_Allreduce(&mat->was_assembled,&other_disassembled,1,MPIU_BOOL,MPI_PROD,PetscObjectComm((PetscObject)mat));
765: if (mat->was_assembled && !other_disassembled) {
766: MatDisAssemble_MPISBAIJ(mat);
767: }
768: }
770: if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
771: MatSetUpMultiply_MPISBAIJ(mat); /* setup Mvctx and sMvctx */
772: }
773: MatAssemblyBegin(baij->B,mode);
774: MatAssemblyEnd(baij->B,mode);
776: PetscFree2(baij->rowvalues,baij->rowindices);
778: baij->rowvalues = 0;
780: /* if no new nonzero locations are allowed in matrix then only set the matrix state the first time through */
781: if ((!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) || !((Mat_SeqBAIJ*)(baij->A->data))->nonew) {
782: PetscObjectState state = baij->A->nonzerostate + baij->B->nonzerostate;
783: MPIU_Allreduce(&state,&mat->nonzerostate,1,MPIU_INT64,MPI_SUM,PetscObjectComm((PetscObject)mat));
784: }
785: return(0);
786: }
788: extern PetscErrorCode MatView_SeqSBAIJ(Mat,PetscViewer);
789: extern PetscErrorCode MatSetValues_MPIBAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
790: #include <petscdraw.h>
793: static PetscErrorCode MatView_MPISBAIJ_ASCIIorDraworSocket(Mat mat,PetscViewer viewer)
794: {
795: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
796: PetscErrorCode ierr;
797: PetscInt bs = mat->rmap->bs;
798: PetscMPIInt rank = baij->rank;
799: PetscBool iascii,isdraw;
800: PetscViewer sviewer;
801: PetscViewerFormat format;
804: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
805: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
806: if (iascii) {
807: PetscViewerGetFormat(viewer,&format);
808: if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
809: MatInfo info;
810: MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&rank);
811: MatGetInfo(mat,MAT_LOCAL,&info);
812: PetscViewerASCIIPushSynchronized(viewer);
813: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Local rows %D nz %D nz alloced %D bs %D mem %D\n",rank,mat->rmap->n,(PetscInt)info.nz_used,(PetscInt)info.nz_allocated,mat->rmap->bs,(PetscInt)info.memory);
814: MatGetInfo(baij->A,MAT_LOCAL,&info);
815: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] on-diagonal part: nz %D \n",rank,(PetscInt)info.nz_used);
816: MatGetInfo(baij->B,MAT_LOCAL,&info);
817: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] off-diagonal part: nz %D \n",rank,(PetscInt)info.nz_used);
818: PetscViewerFlush(viewer);
819: PetscViewerASCIIPopSynchronized(viewer);
820: PetscViewerASCIIPrintf(viewer,"Information on VecScatter used in matrix-vector product: \n");
821: VecScatterView(baij->Mvctx,viewer);
822: return(0);
823: } else if (format == PETSC_VIEWER_ASCII_INFO) {
824: PetscViewerASCIIPrintf(viewer," block size is %D\n",bs);
825: return(0);
826: } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
827: return(0);
828: }
829: }
831: if (isdraw) {
832: PetscDraw draw;
833: PetscBool isnull;
834: PetscViewerDrawGetDraw(viewer,0,&draw);
835: PetscDrawIsNull(draw,&isnull);
836: if (isnull) return(0);
837: }
839: {
840: /* assemble the entire matrix onto first processor. */
841: Mat A;
842: Mat_SeqSBAIJ *Aloc;
843: Mat_SeqBAIJ *Bloc;
844: PetscInt M = mat->rmap->N,N = mat->cmap->N,*ai,*aj,col,i,j,k,*rvals,mbs = baij->mbs;
845: MatScalar *a;
846: const char *matname;
848: /* Should this be the same type as mat? */
849: MatCreate(PetscObjectComm((PetscObject)mat),&A);
850: if (!rank) {
851: MatSetSizes(A,M,N,M,N);
852: } else {
853: MatSetSizes(A,0,0,M,N);
854: }
855: MatSetType(A,MATMPISBAIJ);
856: MatMPISBAIJSetPreallocation(A,mat->rmap->bs,0,NULL,0,NULL);
857: MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);
858: PetscLogObjectParent((PetscObject)mat,(PetscObject)A);
860: /* copy over the A part */
861: Aloc = (Mat_SeqSBAIJ*)baij->A->data;
862: ai = Aloc->i; aj = Aloc->j; a = Aloc->a;
863: PetscMalloc1(bs,&rvals);
865: for (i=0; i<mbs; i++) {
866: rvals[0] = bs*(baij->rstartbs + i);
867: for (j=1; j<bs; j++) rvals[j] = rvals[j-1] + 1;
868: for (j=ai[i]; j<ai[i+1]; j++) {
869: col = (baij->cstartbs+aj[j])*bs;
870: for (k=0; k<bs; k++) {
871: MatSetValues_MPISBAIJ(A,bs,rvals,1,&col,a,INSERT_VALUES);
872: col++;
873: a += bs;
874: }
875: }
876: }
877: /* copy over the B part */
878: Bloc = (Mat_SeqBAIJ*)baij->B->data;
879: ai = Bloc->i; aj = Bloc->j; a = Bloc->a;
880: for (i=0; i<mbs; i++) {
882: rvals[0] = bs*(baij->rstartbs + i);
883: for (j=1; j<bs; j++) rvals[j] = rvals[j-1] + 1;
884: for (j=ai[i]; j<ai[i+1]; j++) {
885: col = baij->garray[aj[j]]*bs;
886: for (k=0; k<bs; k++) {
887: MatSetValues_MPIBAIJ(A,bs,rvals,1,&col,a,INSERT_VALUES);
888: col++;
889: a += bs;
890: }
891: }
892: }
893: PetscFree(rvals);
894: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
895: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
896: /*
897: Everyone has to call to draw the matrix since the graphics waits are
898: synchronized across all processors that share the PetscDraw object
899: */
900: PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);
901: PetscObjectGetName((PetscObject)mat,&matname);
902: if (!rank) {
903: PetscObjectSetName((PetscObject)((Mat_MPISBAIJ*)(A->data))->A,matname);
904: MatView_SeqSBAIJ(((Mat_MPISBAIJ*)(A->data))->A,sviewer);
905: }
906: PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);
907: PetscViewerFlush(viewer);
908: MatDestroy(&A);
909: }
910: return(0);
911: }
915: PetscErrorCode MatView_MPISBAIJ(Mat mat,PetscViewer viewer)
916: {
918: PetscBool iascii,isdraw,issocket,isbinary;
921: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
922: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
923: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSOCKET,&issocket);
924: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
925: if (iascii || isdraw || issocket || isbinary) {
926: MatView_MPISBAIJ_ASCIIorDraworSocket(mat,viewer);
927: }
928: return(0);
929: }
933: PetscErrorCode MatDestroy_MPISBAIJ(Mat mat)
934: {
935: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
939: #if defined(PETSC_USE_LOG)
940: PetscLogObjectState((PetscObject)mat,"Rows=%D,Cols=%D",mat->rmap->N,mat->cmap->N);
941: #endif
942: MatStashDestroy_Private(&mat->stash);
943: MatStashDestroy_Private(&mat->bstash);
944: MatDestroy(&baij->A);
945: MatDestroy(&baij->B);
946: #if defined(PETSC_USE_CTABLE)
947: PetscTableDestroy(&baij->colmap);
948: #else
949: PetscFree(baij->colmap);
950: #endif
951: PetscFree(baij->garray);
952: VecDestroy(&baij->lvec);
953: VecScatterDestroy(&baij->Mvctx);
954: VecDestroy(&baij->slvec0);
955: VecDestroy(&baij->slvec0b);
956: VecDestroy(&baij->slvec1);
957: VecDestroy(&baij->slvec1a);
958: VecDestroy(&baij->slvec1b);
959: VecScatterDestroy(&baij->sMvctx);
960: PetscFree2(baij->rowvalues,baij->rowindices);
961: PetscFree(baij->barray);
962: PetscFree(baij->hd);
963: VecDestroy(&baij->diag);
964: VecDestroy(&baij->bb1);
965: VecDestroy(&baij->xx1);
966: #if defined(PETSC_USE_REAL_MAT_SINGLE)
967: PetscFree(baij->setvaluescopy);
968: #endif
969: PetscFree(baij->in_loc);
970: PetscFree(baij->v_loc);
971: PetscFree(baij->rangebs);
972: PetscFree(mat->data);
974: PetscObjectChangeTypeName((PetscObject)mat,0);
975: PetscObjectComposeFunction((PetscObject)mat,"MatStoreValues_C",NULL);
976: PetscObjectComposeFunction((PetscObject)mat,"MatRetrieveValues_C",NULL);
977: PetscObjectComposeFunction((PetscObject)mat,"MatGetDiagonalBlock_C",NULL);
978: PetscObjectComposeFunction((PetscObject)mat,"MatMPISBAIJSetPreallocation_C",NULL);
979: PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpisbaij_mpisbstrm_C",NULL);
980: #if defined(PETSC_HAVE_ELEMENTAL)
981: PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpisbaij_elemental_C",NULL);
982: #endif
983: return(0);
984: }
988: PetscErrorCode MatMult_MPISBAIJ_Hermitian(Mat A,Vec xx,Vec yy)
989: {
990: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
991: PetscErrorCode ierr;
992: PetscInt nt,mbs=a->mbs,bs=A->rmap->bs;
993: PetscScalar *from;
994: const PetscScalar *x;
997: VecGetLocalSize(xx,&nt);
998: if (nt != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible partition of A and xx");
1000: /* diagonal part */
1001: (*a->A->ops->mult)(a->A,xx,a->slvec1a);
1002: VecSet(a->slvec1b,0.0);
1004: /* subdiagonal part */
1005: (*a->B->ops->multhermitiantranspose)(a->B,xx,a->slvec0b);
1007: /* copy x into the vec slvec0 */
1008: VecGetArray(a->slvec0,&from);
1009: VecGetArrayRead(xx,&x);
1011: PetscMemcpy(from,x,bs*mbs*sizeof(MatScalar));
1012: VecRestoreArray(a->slvec0,&from);
1013: VecRestoreArrayRead(xx,&x);
1015: VecScatterBegin(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
1016: VecScatterEnd(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
1017: /* supperdiagonal part */
1018: (*a->B->ops->multadd)(a->B,a->slvec1b,a->slvec1a,yy);
1019: return(0);
1020: }
1024: PetscErrorCode MatMult_MPISBAIJ(Mat A,Vec xx,Vec yy)
1025: {
1026: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1027: PetscErrorCode ierr;
1028: PetscInt nt,mbs=a->mbs,bs=A->rmap->bs;
1029: PetscScalar *from;
1030: const PetscScalar *x;
1033: VecGetLocalSize(xx,&nt);
1034: if (nt != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible partition of A and xx");
1036: /* diagonal part */
1037: (*a->A->ops->mult)(a->A,xx,a->slvec1a);
1038: VecSet(a->slvec1b,0.0);
1040: /* subdiagonal part */
1041: (*a->B->ops->multtranspose)(a->B,xx,a->slvec0b);
1043: /* copy x into the vec slvec0 */
1044: VecGetArray(a->slvec0,&from);
1045: VecGetArrayRead(xx,&x);
1047: PetscMemcpy(from,x,bs*mbs*sizeof(MatScalar));
1048: VecRestoreArray(a->slvec0,&from);
1049: VecRestoreArrayRead(xx,&x);
1051: VecScatterBegin(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
1052: VecScatterEnd(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
1053: /* supperdiagonal part */
1054: (*a->B->ops->multadd)(a->B,a->slvec1b,a->slvec1a,yy);
1055: return(0);
1056: }
1060: PetscErrorCode MatMult_MPISBAIJ_2comm(Mat A,Vec xx,Vec yy)
1061: {
1062: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1064: PetscInt nt;
1067: VecGetLocalSize(xx,&nt);
1068: if (nt != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible partition of A and xx");
1070: VecGetLocalSize(yy,&nt);
1071: if (nt != A->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible parition of A and yy");
1073: VecScatterBegin(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);
1074: /* do diagonal part */
1075: (*a->A->ops->mult)(a->A,xx,yy);
1076: /* do supperdiagonal part */
1077: VecScatterEnd(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);
1078: (*a->B->ops->multadd)(a->B,a->lvec,yy,yy);
1079: /* do subdiagonal part */
1080: (*a->B->ops->multtranspose)(a->B,xx,a->lvec);
1081: VecScatterBegin(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);
1082: VecScatterEnd(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);
1083: return(0);
1084: }
1088: PetscErrorCode MatMultAdd_MPISBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
1089: {
1090: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1091: PetscErrorCode ierr;
1092: PetscInt mbs=a->mbs,bs=A->rmap->bs;
1093: PetscScalar *from,zero=0.0;
1094: const PetscScalar *x;
1097: /*
1098: PetscSynchronizedPrintf(PetscObjectComm((PetscObject)A)," MatMultAdd is called ...\n");
1099: PetscSynchronizedFlush(PetscObjectComm((PetscObject)A),PETSC_STDOUT);
1100: */
1101: /* diagonal part */
1102: (*a->A->ops->multadd)(a->A,xx,yy,a->slvec1a);
1103: VecSet(a->slvec1b,zero);
1105: /* subdiagonal part */
1106: (*a->B->ops->multtranspose)(a->B,xx,a->slvec0b);
1108: /* copy x into the vec slvec0 */
1109: VecGetArray(a->slvec0,&from);
1110: VecGetArrayRead(xx,&x);
1111: PetscMemcpy(from,x,bs*mbs*sizeof(MatScalar));
1112: VecRestoreArray(a->slvec0,&from);
1114: VecScatterBegin(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
1115: VecRestoreArrayRead(xx,&x);
1116: VecScatterEnd(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
1118: /* supperdiagonal part */
1119: (*a->B->ops->multadd)(a->B,a->slvec1b,a->slvec1a,zz);
1120: return(0);
1121: }
1125: PetscErrorCode MatMultAdd_MPISBAIJ_2comm(Mat A,Vec xx,Vec yy,Vec zz)
1126: {
1127: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1131: VecScatterBegin(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);
1132: /* do diagonal part */
1133: (*a->A->ops->multadd)(a->A,xx,yy,zz);
1134: /* do supperdiagonal part */
1135: VecScatterEnd(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);
1136: (*a->B->ops->multadd)(a->B,a->lvec,zz,zz);
1138: /* do subdiagonal part */
1139: (*a->B->ops->multtranspose)(a->B,xx,a->lvec);
1140: VecScatterBegin(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);
1141: VecScatterEnd(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);
1142: return(0);
1143: }
1145: /*
1146: This only works correctly for square matrices where the subblock A->A is the
1147: diagonal block
1148: */
1151: PetscErrorCode MatGetDiagonal_MPISBAIJ(Mat A,Vec v)
1152: {
1153: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1157: /* if (a->rmap->N != a->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Supports only square matrix where A->A is diag block"); */
1158: MatGetDiagonal(a->A,v);
1159: return(0);
1160: }
1164: PetscErrorCode MatScale_MPISBAIJ(Mat A,PetscScalar aa)
1165: {
1166: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1170: MatScale(a->A,aa);
1171: MatScale(a->B,aa);
1172: return(0);
1173: }
1177: PetscErrorCode MatGetRow_MPISBAIJ(Mat matin,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
1178: {
1179: Mat_MPISBAIJ *mat = (Mat_MPISBAIJ*)matin->data;
1180: PetscScalar *vworkA,*vworkB,**pvA,**pvB,*v_p;
1182: PetscInt bs = matin->rmap->bs,bs2 = mat->bs2,i,*cworkA,*cworkB,**pcA,**pcB;
1183: PetscInt nztot,nzA,nzB,lrow,brstart = matin->rmap->rstart,brend = matin->rmap->rend;
1184: PetscInt *cmap,*idx_p,cstart = mat->rstartbs;
1187: if (mat->getrowactive) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Already active");
1188: mat->getrowactive = PETSC_TRUE;
1190: if (!mat->rowvalues && (idx || v)) {
1191: /*
1192: allocate enough space to hold information from the longest row.
1193: */
1194: Mat_SeqSBAIJ *Aa = (Mat_SeqSBAIJ*)mat->A->data;
1195: Mat_SeqBAIJ *Ba = (Mat_SeqBAIJ*)mat->B->data;
1196: PetscInt max = 1,mbs = mat->mbs,tmp;
1197: for (i=0; i<mbs; i++) {
1198: tmp = Aa->i[i+1] - Aa->i[i] + Ba->i[i+1] - Ba->i[i]; /* row length */
1199: if (max < tmp) max = tmp;
1200: }
1201: PetscMalloc2(max*bs2,&mat->rowvalues,max*bs2,&mat->rowindices);
1202: }
1204: if (row < brstart || row >= brend) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only local rows");
1205: lrow = row - brstart; /* local row index */
1207: pvA = &vworkA; pcA = &cworkA; pvB = &vworkB; pcB = &cworkB;
1208: if (!v) {pvA = 0; pvB = 0;}
1209: if (!idx) {pcA = 0; if (!v) pcB = 0;}
1210: (*mat->A->ops->getrow)(mat->A,lrow,&nzA,pcA,pvA);
1211: (*mat->B->ops->getrow)(mat->B,lrow,&nzB,pcB,pvB);
1212: nztot = nzA + nzB;
1214: cmap = mat->garray;
1215: if (v || idx) {
1216: if (nztot) {
1217: /* Sort by increasing column numbers, assuming A and B already sorted */
1218: PetscInt imark = -1;
1219: if (v) {
1220: *v = v_p = mat->rowvalues;
1221: for (i=0; i<nzB; i++) {
1222: if (cmap[cworkB[i]/bs] < cstart) v_p[i] = vworkB[i];
1223: else break;
1224: }
1225: imark = i;
1226: for (i=0; i<nzA; i++) v_p[imark+i] = vworkA[i];
1227: for (i=imark; i<nzB; i++) v_p[nzA+i] = vworkB[i];
1228: }
1229: if (idx) {
1230: *idx = idx_p = mat->rowindices;
1231: if (imark > -1) {
1232: for (i=0; i<imark; i++) {
1233: idx_p[i] = cmap[cworkB[i]/bs]*bs + cworkB[i]%bs;
1234: }
1235: } else {
1236: for (i=0; i<nzB; i++) {
1237: if (cmap[cworkB[i]/bs] < cstart) idx_p[i] = cmap[cworkB[i]/bs]*bs + cworkB[i]%bs;
1238: else break;
1239: }
1240: imark = i;
1241: }
1242: for (i=0; i<nzA; i++) idx_p[imark+i] = cstart*bs + cworkA[i];
1243: for (i=imark; i<nzB; i++) idx_p[nzA+i] = cmap[cworkB[i]/bs]*bs + cworkB[i]%bs ;
1244: }
1245: } else {
1246: if (idx) *idx = 0;
1247: if (v) *v = 0;
1248: }
1249: }
1250: *nz = nztot;
1251: (*mat->A->ops->restorerow)(mat->A,lrow,&nzA,pcA,pvA);
1252: (*mat->B->ops->restorerow)(mat->B,lrow,&nzB,pcB,pvB);
1253: return(0);
1254: }
1258: PetscErrorCode MatRestoreRow_MPISBAIJ(Mat mat,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
1259: {
1260: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
1263: if (!baij->getrowactive) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"MatGetRow() must be called first");
1264: baij->getrowactive = PETSC_FALSE;
1265: return(0);
1266: }
1270: PetscErrorCode MatGetRowUpperTriangular_MPISBAIJ(Mat A)
1271: {
1272: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1273: Mat_SeqSBAIJ *aA = (Mat_SeqSBAIJ*)a->A->data;
1276: aA->getrow_utriangular = PETSC_TRUE;
1277: return(0);
1278: }
1281: PetscErrorCode MatRestoreRowUpperTriangular_MPISBAIJ(Mat A)
1282: {
1283: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1284: Mat_SeqSBAIJ *aA = (Mat_SeqSBAIJ*)a->A->data;
1287: aA->getrow_utriangular = PETSC_FALSE;
1288: return(0);
1289: }
1293: PetscErrorCode MatRealPart_MPISBAIJ(Mat A)
1294: {
1295: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1299: MatRealPart(a->A);
1300: MatRealPart(a->B);
1301: return(0);
1302: }
1306: PetscErrorCode MatImaginaryPart_MPISBAIJ(Mat A)
1307: {
1308: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1312: MatImaginaryPart(a->A);
1313: MatImaginaryPart(a->B);
1314: return(0);
1315: }
1317: /* Check if isrow is a subset of iscol_local, called by MatGetSubMatrix_MPISBAIJ()
1318: Input: isrow - distributed(parallel),
1319: iscol_local - locally owned (seq)
1320: */
1323: PetscErrorCode ISEqual_private(IS isrow,IS iscol_local,PetscBool *flg)
1324: {
1326: PetscInt sz1,sz2,*a1,*a2,i,j,k,nmatch;
1327: const PetscInt *ptr1,*ptr2;
1330: ISGetLocalSize(isrow,&sz1);
1331: ISGetLocalSize(iscol_local,&sz2);
1332: if (sz1 > sz2) {
1333: *flg = PETSC_FALSE;
1334: return(0);
1335: }
1337: ISGetIndices(isrow,&ptr1);
1338: ISGetIndices(iscol_local,&ptr2);
1340: PetscMalloc1(sz1,&a1);
1341: PetscMalloc1(sz2,&a2);
1342: PetscMemcpy(a1,ptr1,sz1*sizeof(PetscInt));
1343: PetscMemcpy(a2,ptr2,sz2*sizeof(PetscInt));
1344: PetscSortInt(sz1,a1);
1345: PetscSortInt(sz2,a2);
1347: nmatch=0;
1348: k = 0;
1349: for (i=0; i<sz1; i++){
1350: for (j=k; j<sz2; j++){
1351: if (a1[i] == a2[j]) {
1352: k = j; nmatch++;
1353: break;
1354: }
1355: }
1356: }
1357: ISRestoreIndices(isrow,&ptr1);
1358: ISRestoreIndices(iscol_local,&ptr2);
1359: PetscFree(a1);
1360: PetscFree(a2);
1361: if (nmatch < sz1) {
1362: *flg = PETSC_FALSE;
1363: } else {
1364: *flg = PETSC_TRUE;
1365: }
1366: return(0);
1367: }
1371: PetscErrorCode MatGetSubMatrix_MPISBAIJ(Mat mat,IS isrow,IS iscol,MatReuse call,Mat *newmat)
1372: {
1374: IS iscol_local;
1375: PetscInt csize;
1376: PetscBool isequal;
1379: ISGetLocalSize(iscol,&csize);
1380: if (call == MAT_REUSE_MATRIX) {
1381: PetscObjectQuery((PetscObject)*newmat,"ISAllGather",(PetscObject*)&iscol_local);
1382: if (!iscol_local) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Submatrix passed in was not used before, cannot reuse");
1383: } else {
1384: ISAllGather(iscol,&iscol_local);
1385: ISEqual_private(isrow,iscol_local,&isequal);
1386: if (!isequal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"For symmetric format, iscol must equal isrow");
1387: }
1389: /* now call MatGetSubMatrix_MPIBAIJ() */
1390: MatGetSubMatrix_MPIBAIJ_Private(mat,isrow,iscol_local,csize,call,newmat);
1391: if (call == MAT_INITIAL_MATRIX) {
1392: PetscObjectCompose((PetscObject)*newmat,"ISAllGather",(PetscObject)iscol_local);
1393: ISDestroy(&iscol_local);
1394: }
1395: return(0);
1396: }
1400: PetscErrorCode MatZeroEntries_MPISBAIJ(Mat A)
1401: {
1402: Mat_MPISBAIJ *l = (Mat_MPISBAIJ*)A->data;
1406: MatZeroEntries(l->A);
1407: MatZeroEntries(l->B);
1408: return(0);
1409: }
1413: PetscErrorCode MatGetInfo_MPISBAIJ(Mat matin,MatInfoType flag,MatInfo *info)
1414: {
1415: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)matin->data;
1416: Mat A = a->A,B = a->B;
1418: PetscReal isend[5],irecv[5];
1421: info->block_size = (PetscReal)matin->rmap->bs;
1423: MatGetInfo(A,MAT_LOCAL,info);
1425: isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
1426: isend[3] = info->memory; isend[4] = info->mallocs;
1428: MatGetInfo(B,MAT_LOCAL,info);
1430: isend[0] += info->nz_used; isend[1] += info->nz_allocated; isend[2] += info->nz_unneeded;
1431: isend[3] += info->memory; isend[4] += info->mallocs;
1432: if (flag == MAT_LOCAL) {
1433: info->nz_used = isend[0];
1434: info->nz_allocated = isend[1];
1435: info->nz_unneeded = isend[2];
1436: info->memory = isend[3];
1437: info->mallocs = isend[4];
1438: } else if (flag == MAT_GLOBAL_MAX) {
1439: MPIU_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)matin));
1441: info->nz_used = irecv[0];
1442: info->nz_allocated = irecv[1];
1443: info->nz_unneeded = irecv[2];
1444: info->memory = irecv[3];
1445: info->mallocs = irecv[4];
1446: } else if (flag == MAT_GLOBAL_SUM) {
1447: MPIU_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)matin));
1449: info->nz_used = irecv[0];
1450: info->nz_allocated = irecv[1];
1451: info->nz_unneeded = irecv[2];
1452: info->memory = irecv[3];
1453: info->mallocs = irecv[4];
1454: } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown MatInfoType argument %d",(int)flag);
1455: info->fill_ratio_given = 0; /* no parallel LU/ILU/Cholesky */
1456: info->fill_ratio_needed = 0;
1457: info->factor_mallocs = 0;
1458: return(0);
1459: }
1463: PetscErrorCode MatSetOption_MPISBAIJ(Mat A,MatOption op,PetscBool flg)
1464: {
1465: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1466: Mat_SeqSBAIJ *aA = (Mat_SeqSBAIJ*)a->A->data;
1470: switch (op) {
1471: case MAT_NEW_NONZERO_LOCATIONS:
1472: case MAT_NEW_NONZERO_ALLOCATION_ERR:
1473: case MAT_UNUSED_NONZERO_LOCATION_ERR:
1474: case MAT_KEEP_NONZERO_PATTERN:
1475: case MAT_NEW_NONZERO_LOCATION_ERR:
1476: MatCheckPreallocated(A,1);
1477: MatSetOption(a->A,op,flg);
1478: MatSetOption(a->B,op,flg);
1479: break;
1480: case MAT_ROW_ORIENTED:
1481: MatCheckPreallocated(A,1);
1482: a->roworiented = flg;
1484: MatSetOption(a->A,op,flg);
1485: MatSetOption(a->B,op,flg);
1486: break;
1487: case MAT_NEW_DIAGONALS:
1488: PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
1489: break;
1490: case MAT_IGNORE_OFF_PROC_ENTRIES:
1491: a->donotstash = flg;
1492: break;
1493: case MAT_USE_HASH_TABLE:
1494: a->ht_flag = flg;
1495: break;
1496: case MAT_HERMITIAN:
1497: MatCheckPreallocated(A,1);
1498: if (!A->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatAssemblyEnd() first");
1499: MatSetOption(a->A,op,flg);
1501: A->ops->mult = MatMult_MPISBAIJ_Hermitian;
1502: break;
1503: case MAT_SPD:
1504: A->spd_set = PETSC_TRUE;
1505: A->spd = flg;
1506: if (flg) {
1507: A->symmetric = PETSC_TRUE;
1508: A->structurally_symmetric = PETSC_TRUE;
1509: A->symmetric_set = PETSC_TRUE;
1510: A->structurally_symmetric_set = PETSC_TRUE;
1511: }
1512: break;
1513: case MAT_SYMMETRIC:
1514: MatCheckPreallocated(A,1);
1515: MatSetOption(a->A,op,flg);
1516: break;
1517: case MAT_STRUCTURALLY_SYMMETRIC:
1518: MatCheckPreallocated(A,1);
1519: MatSetOption(a->A,op,flg);
1520: break;
1521: case MAT_SYMMETRY_ETERNAL:
1522: if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix must be symmetric");
1523: PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
1524: break;
1525: case MAT_IGNORE_LOWER_TRIANGULAR:
1526: aA->ignore_ltriangular = flg;
1527: break;
1528: case MAT_ERROR_LOWER_TRIANGULAR:
1529: aA->ignore_ltriangular = flg;
1530: break;
1531: case MAT_GETROW_UPPERTRIANGULAR:
1532: aA->getrow_utriangular = flg;
1533: break;
1534: default:
1535: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %d",op);
1536: }
1537: return(0);
1538: }
1542: PetscErrorCode MatTranspose_MPISBAIJ(Mat A,MatReuse reuse,Mat *B)
1543: {
1547: if (MAT_INITIAL_MATRIX || *B != A) {
1548: MatDuplicate(A,MAT_COPY_VALUES,B);
1549: }
1550: return(0);
1551: }
1555: PetscErrorCode MatDiagonalScale_MPISBAIJ(Mat mat,Vec ll,Vec rr)
1556: {
1557: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
1558: Mat a = baij->A, b=baij->B;
1560: PetscInt nv,m,n;
1561: PetscBool flg;
1564: if (ll != rr) {
1565: VecEqual(ll,rr,&flg);
1566: if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"For symmetric format, left and right scaling vectors must be same\n");
1567: }
1568: if (!ll) return(0);
1570: MatGetLocalSize(mat,&m,&n);
1571: if (m != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"For symmetric format, local size %d %d must be same",m,n);
1573: VecGetLocalSize(rr,&nv);
1574: if (nv!=n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left and right vector non-conforming local size");
1576: VecScatterBegin(baij->Mvctx,rr,baij->lvec,INSERT_VALUES,SCATTER_FORWARD);
1578: /* left diagonalscale the off-diagonal part */
1579: (*b->ops->diagonalscale)(b,ll,NULL);
1581: /* scale the diagonal part */
1582: (*a->ops->diagonalscale)(a,ll,rr);
1584: /* right diagonalscale the off-diagonal part */
1585: VecScatterEnd(baij->Mvctx,rr,baij->lvec,INSERT_VALUES,SCATTER_FORWARD);
1586: (*b->ops->diagonalscale)(b,NULL,baij->lvec);
1587: return(0);
1588: }
1592: PetscErrorCode MatSetUnfactored_MPISBAIJ(Mat A)
1593: {
1594: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1598: MatSetUnfactored(a->A);
1599: return(0);
1600: }
1602: static PetscErrorCode MatDuplicate_MPISBAIJ(Mat,MatDuplicateOption,Mat*);
1606: PetscErrorCode MatEqual_MPISBAIJ(Mat A,Mat B,PetscBool *flag)
1607: {
1608: Mat_MPISBAIJ *matB = (Mat_MPISBAIJ*)B->data,*matA = (Mat_MPISBAIJ*)A->data;
1609: Mat a,b,c,d;
1610: PetscBool flg;
1614: a = matA->A; b = matA->B;
1615: c = matB->A; d = matB->B;
1617: MatEqual(a,c,&flg);
1618: if (flg) {
1619: MatEqual(b,d,&flg);
1620: }
1621: MPIU_Allreduce(&flg,flag,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));
1622: return(0);
1623: }
1627: PetscErrorCode MatCopy_MPISBAIJ(Mat A,Mat B,MatStructure str)
1628: {
1630: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1631: Mat_MPISBAIJ *b = (Mat_MPISBAIJ*)B->data;
1634: /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
1635: if ((str != SAME_NONZERO_PATTERN) || (A->ops->copy != B->ops->copy)) {
1636: MatGetRowUpperTriangular(A);
1637: MatCopy_Basic(A,B,str);
1638: MatRestoreRowUpperTriangular(A);
1639: } else {
1640: MatCopy(a->A,b->A,str);
1641: MatCopy(a->B,b->B,str);
1642: }
1643: return(0);
1644: }
1648: PetscErrorCode MatSetUp_MPISBAIJ(Mat A)
1649: {
1653: MatMPISBAIJSetPreallocation(A,A->rmap->bs,PETSC_DEFAULT,0,PETSC_DEFAULT,0);
1654: return(0);
1655: }
1659: PetscErrorCode MatAXPY_MPISBAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str)
1660: {
1662: Mat_MPISBAIJ *xx=(Mat_MPISBAIJ*)X->data,*yy=(Mat_MPISBAIJ*)Y->data;
1663: PetscBLASInt bnz,one=1;
1664: Mat_SeqSBAIJ *xa,*ya;
1665: Mat_SeqBAIJ *xb,*yb;
1668: if (str == SAME_NONZERO_PATTERN) {
1669: PetscScalar alpha = a;
1670: xa = (Mat_SeqSBAIJ*)xx->A->data;
1671: ya = (Mat_SeqSBAIJ*)yy->A->data;
1672: PetscBLASIntCast(xa->nz,&bnz);
1673: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&bnz,&alpha,xa->a,&one,ya->a,&one));
1674: xb = (Mat_SeqBAIJ*)xx->B->data;
1675: yb = (Mat_SeqBAIJ*)yy->B->data;
1676: PetscBLASIntCast(xb->nz,&bnz);
1677: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&bnz,&alpha,xb->a,&one,yb->a,&one));
1678: PetscObjectStateIncrease((PetscObject)Y);
1679: } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
1680: MatSetOption(X,MAT_GETROW_UPPERTRIANGULAR,PETSC_TRUE);
1681: MatAXPY_Basic(Y,a,X,str);
1682: MatSetOption(X,MAT_GETROW_UPPERTRIANGULAR,PETSC_FALSE);
1683: } else {
1684: Mat B;
1685: PetscInt *nnz_d,*nnz_o,bs=Y->rmap->bs;
1686: if (bs != X->rmap->bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrices must have same block size");
1687: MatGetRowUpperTriangular(X);
1688: MatGetRowUpperTriangular(Y);
1689: PetscMalloc1(yy->A->rmap->N,&nnz_d);
1690: PetscMalloc1(yy->B->rmap->N,&nnz_o);
1691: MatCreate(PetscObjectComm((PetscObject)Y),&B);
1692: PetscObjectSetName((PetscObject)B,((PetscObject)Y)->name);
1693: MatSetSizes(B,Y->rmap->n,Y->cmap->n,Y->rmap->N,Y->cmap->N);
1694: MatSetBlockSizesFromMats(B,Y,Y);
1695: MatSetType(B,MATMPISBAIJ);
1696: MatAXPYGetPreallocation_SeqSBAIJ(yy->A,xx->A,nnz_d);
1697: MatAXPYGetPreallocation_MPIBAIJ(yy->B,yy->garray,xx->B,xx->garray,nnz_o);
1698: MatMPISBAIJSetPreallocation(B,bs,0,nnz_d,0,nnz_o);
1699: MatAXPY_BasicWithPreallocation(B,Y,a,X,str);
1700: MatHeaderReplace(Y,&B);
1701: PetscFree(nnz_d);
1702: PetscFree(nnz_o);
1703: MatRestoreRowUpperTriangular(X);
1704: MatRestoreRowUpperTriangular(Y);
1705: }
1706: return(0);
1707: }
1711: PetscErrorCode MatGetSubMatrices_MPISBAIJ(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1712: {
1714: PetscInt i;
1715: PetscBool flg;
1718: MatGetSubMatrices_MPIBAIJ(A,n,irow,icol,scall,B); /* B[] are sbaij matrices */
1719: for (i=0; i<n; i++) {
1720: ISEqual(irow[i],icol[i],&flg);
1721: if (!flg) {
1722: MatSeqSBAIJZeroOps_Private(*B[i]);
1723: }
1724: }
1725: return(0);
1726: }
1730: PetscErrorCode MatShift_MPISBAIJ(Mat Y,PetscScalar a)
1731: {
1733: Mat_MPISBAIJ *maij = (Mat_MPISBAIJ*)Y->data;
1734: Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ*)maij->A->data;
1737: if (!Y->preallocated) {
1738: MatMPISBAIJSetPreallocation(Y,Y->rmap->bs,1,NULL,0,NULL);
1739: } else if (!aij->nz) {
1740: PetscInt nonew = aij->nonew;
1741: MatSeqSBAIJSetPreallocation(maij->A,Y->rmap->bs,1,NULL);
1742: aij->nonew = nonew;
1743: }
1744: MatShift_Basic(Y,a);
1745: return(0);
1746: }
1750: PetscErrorCode MatMissingDiagonal_MPISBAIJ(Mat A,PetscBool *missing,PetscInt *d)
1751: {
1752: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1756: if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only works for square matrices");
1757: MatMissingDiagonal(a->A,missing,d);
1758: if (d) {
1759: PetscInt rstart;
1760: MatGetOwnershipRange(A,&rstart,NULL);
1761: *d += rstart/A->rmap->bs;
1763: }
1764: return(0);
1765: }
1768: /* -------------------------------------------------------------------*/
1769: static struct _MatOps MatOps_Values = {MatSetValues_MPISBAIJ,
1770: MatGetRow_MPISBAIJ,
1771: MatRestoreRow_MPISBAIJ,
1772: MatMult_MPISBAIJ,
1773: /* 4*/ MatMultAdd_MPISBAIJ,
1774: MatMult_MPISBAIJ, /* transpose versions are same as non-transpose */
1775: MatMultAdd_MPISBAIJ,
1776: 0,
1777: 0,
1778: 0,
1779: /* 10*/ 0,
1780: 0,
1781: 0,
1782: MatSOR_MPISBAIJ,
1783: MatTranspose_MPISBAIJ,
1784: /* 15*/ MatGetInfo_MPISBAIJ,
1785: MatEqual_MPISBAIJ,
1786: MatGetDiagonal_MPISBAIJ,
1787: MatDiagonalScale_MPISBAIJ,
1788: MatNorm_MPISBAIJ,
1789: /* 20*/ MatAssemblyBegin_MPISBAIJ,
1790: MatAssemblyEnd_MPISBAIJ,
1791: MatSetOption_MPISBAIJ,
1792: MatZeroEntries_MPISBAIJ,
1793: /* 24*/ 0,
1794: 0,
1795: 0,
1796: 0,
1797: 0,
1798: /* 29*/ MatSetUp_MPISBAIJ,
1799: 0,
1800: 0,
1801: 0,
1802: 0,
1803: /* 34*/ MatDuplicate_MPISBAIJ,
1804: 0,
1805: 0,
1806: 0,
1807: 0,
1808: /* 39*/ MatAXPY_MPISBAIJ,
1809: MatGetSubMatrices_MPISBAIJ,
1810: MatIncreaseOverlap_MPISBAIJ,
1811: MatGetValues_MPISBAIJ,
1812: MatCopy_MPISBAIJ,
1813: /* 44*/ 0,
1814: MatScale_MPISBAIJ,
1815: MatShift_MPISBAIJ,
1816: 0,
1817: 0,
1818: /* 49*/ 0,
1819: 0,
1820: 0,
1821: 0,
1822: 0,
1823: /* 54*/ 0,
1824: 0,
1825: MatSetUnfactored_MPISBAIJ,
1826: 0,
1827: MatSetValuesBlocked_MPISBAIJ,
1828: /* 59*/ MatGetSubMatrix_MPISBAIJ,
1829: 0,
1830: 0,
1831: 0,
1832: 0,
1833: /* 64*/ 0,
1834: 0,
1835: 0,
1836: 0,
1837: 0,
1838: /* 69*/ MatGetRowMaxAbs_MPISBAIJ,
1839: 0,
1840: 0,
1841: 0,
1842: 0,
1843: /* 74*/ 0,
1844: 0,
1845: 0,
1846: 0,
1847: 0,
1848: /* 79*/ 0,
1849: 0,
1850: 0,
1851: 0,
1852: MatLoad_MPISBAIJ,
1853: /* 84*/ 0,
1854: 0,
1855: 0,
1856: 0,
1857: 0,
1858: /* 89*/ 0,
1859: 0,
1860: 0,
1861: 0,
1862: 0,
1863: /* 94*/ 0,
1864: 0,
1865: 0,
1866: 0,
1867: 0,
1868: /* 99*/ 0,
1869: 0,
1870: 0,
1871: 0,
1872: 0,
1873: /*104*/ 0,
1874: MatRealPart_MPISBAIJ,
1875: MatImaginaryPart_MPISBAIJ,
1876: MatGetRowUpperTriangular_MPISBAIJ,
1877: MatRestoreRowUpperTriangular_MPISBAIJ,
1878: /*109*/ 0,
1879: 0,
1880: 0,
1881: 0,
1882: MatMissingDiagonal_MPISBAIJ,
1883: /*114*/ 0,
1884: 0,
1885: 0,
1886: 0,
1887: 0,
1888: /*119*/ 0,
1889: 0,
1890: 0,
1891: 0,
1892: 0,
1893: /*124*/ 0,
1894: 0,
1895: 0,
1896: 0,
1897: 0,
1898: /*129*/ 0,
1899: 0,
1900: 0,
1901: 0,
1902: 0,
1903: /*134*/ 0,
1904: 0,
1905: 0,
1906: 0,
1907: 0,
1908: /*139*/ 0,
1909: 0,
1910: 0,
1911: 0,
1912: 0,
1913: /*144*/MatCreateMPIMatConcatenateSeqMat_MPISBAIJ
1914: };
1918: PetscErrorCode MatGetDiagonalBlock_MPISBAIJ(Mat A,Mat *a)
1919: {
1921: *a = ((Mat_MPISBAIJ*)A->data)->A;
1922: return(0);
1923: }
1927: PetscErrorCode MatMPISBAIJSetPreallocation_MPISBAIJ(Mat B,PetscInt bs,PetscInt d_nz,const PetscInt *d_nnz,PetscInt o_nz,const PetscInt *o_nnz)
1928: {
1929: Mat_MPISBAIJ *b;
1931: PetscInt i,mbs,Mbs;
1934: MatSetBlockSize(B,PetscAbs(bs));
1935: PetscLayoutSetUp(B->rmap);
1936: PetscLayoutSetUp(B->cmap);
1937: PetscLayoutGetBlockSize(B->rmap,&bs);
1939: b = (Mat_MPISBAIJ*)B->data;
1940: mbs = B->rmap->n/bs;
1941: Mbs = B->rmap->N/bs;
1942: if (mbs*bs != B->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"No of local rows %D must be divisible by blocksize %D",B->rmap->N,bs);
1944: B->rmap->bs = bs;
1945: b->bs2 = bs*bs;
1946: b->mbs = mbs;
1947: b->Mbs = Mbs;
1948: b->nbs = B->cmap->n/bs;
1949: b->Nbs = B->cmap->N/bs;
1951: for (i=0; i<=b->size; i++) {
1952: b->rangebs[i] = B->rmap->range[i]/bs;
1953: }
1954: b->rstartbs = B->rmap->rstart/bs;
1955: b->rendbs = B->rmap->rend/bs;
1957: b->cstartbs = B->cmap->rstart/bs;
1958: b->cendbs = B->cmap->rend/bs;
1960: if (!B->preallocated) {
1961: MatCreate(PETSC_COMM_SELF,&b->A);
1962: MatSetSizes(b->A,B->rmap->n,B->cmap->n,B->rmap->n,B->cmap->n);
1963: MatSetType(b->A,MATSEQSBAIJ);
1964: PetscLogObjectParent((PetscObject)B,(PetscObject)b->A);
1965: MatCreate(PETSC_COMM_SELF,&b->B);
1966: MatSetSizes(b->B,B->rmap->n,B->cmap->N,B->rmap->n,B->cmap->N);
1967: MatSetType(b->B,MATSEQBAIJ);
1968: PetscLogObjectParent((PetscObject)B,(PetscObject)b->B);
1969: MatStashCreate_Private(PetscObjectComm((PetscObject)B),bs,&B->bstash);
1970: }
1972: MatSeqSBAIJSetPreallocation(b->A,bs,d_nz,d_nnz);
1973: MatSeqBAIJSetPreallocation(b->B,bs,o_nz,o_nnz);
1975: B->preallocated = PETSC_TRUE;
1976: return(0);
1977: }
1981: PetscErrorCode MatMPISBAIJSetPreallocationCSR_MPISBAIJ(Mat B,PetscInt bs,const PetscInt ii[],const PetscInt jj[],const PetscScalar V[])
1982: {
1983: PetscInt m,rstart,cstart,cend;
1984: PetscInt i,j,d,nz,nz_max=0,*d_nnz=0,*o_nnz=0;
1985: const PetscInt *JJ =0;
1986: PetscScalar *values=0;
1990: if (bs < 1) SETERRQ1(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_OUTOFRANGE,"Invalid block size specified, must be positive but it is %D",bs);
1991: PetscLayoutSetBlockSize(B->rmap,bs);
1992: PetscLayoutSetBlockSize(B->cmap,bs);
1993: PetscLayoutSetUp(B->rmap);
1994: PetscLayoutSetUp(B->cmap);
1995: PetscLayoutGetBlockSize(B->rmap,&bs);
1996: m = B->rmap->n/bs;
1997: rstart = B->rmap->rstart/bs;
1998: cstart = B->cmap->rstart/bs;
1999: cend = B->cmap->rend/bs;
2001: if (ii[0]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"ii[0] must be 0 but it is %D",ii[0]);
2002: PetscMalloc2(m,&d_nnz,m,&o_nnz);
2003: for (i=0; i<m; i++) {
2004: nz = ii[i+1] - ii[i];
2005: if (nz < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local row %D has a negative number of columns %D",i,nz);
2006: nz_max = PetscMax(nz_max,nz);
2007: JJ = jj + ii[i];
2008: for (j=0; j<nz; j++) {
2009: if (*JJ >= cstart) break;
2010: JJ++;
2011: }
2012: d = 0;
2013: for (; j<nz; j++) {
2014: if (*JJ++ >= cend) break;
2015: d++;
2016: }
2017: d_nnz[i] = d;
2018: o_nnz[i] = nz - d;
2019: }
2020: MatMPISBAIJSetPreallocation(B,bs,0,d_nnz,0,o_nnz);
2021: PetscFree2(d_nnz,o_nnz);
2023: values = (PetscScalar*)V;
2024: if (!values) {
2025: PetscMalloc1(bs*bs*nz_max,&values);
2026: PetscMemzero(values,bs*bs*nz_max*sizeof(PetscScalar));
2027: }
2028: for (i=0; i<m; i++) {
2029: PetscInt row = i + rstart;
2030: PetscInt ncols = ii[i+1] - ii[i];
2031: const PetscInt *icols = jj + ii[i];
2032: const PetscScalar *svals = values + (V ? (bs*bs*ii[i]) : 0);
2033: MatSetValuesBlocked_MPISBAIJ(B,1,&row,ncols,icols,svals,INSERT_VALUES);
2034: }
2036: if (!V) { PetscFree(values); }
2037: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
2038: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
2039: MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
2040: return(0);
2041: }
2043: /*MC
2044: MATMPISBAIJ - MATMPISBAIJ = "mpisbaij" - A matrix type to be used for distributed symmetric sparse block matrices,
2045: based on block compressed sparse row format. Only the upper triangular portion of the "diagonal" portion of
2046: the matrix is stored.
2048: For complex numbers by default this matrix is symmetric, NOT Hermitian symmetric. To make it Hermitian symmetric you
2049: can call MatSetOption(Mat, MAT_HERMITIAN);
2051: Options Database Keys:
2052: . -mat_type mpisbaij - sets the matrix type to "mpisbaij" during a call to MatSetFromOptions()
2054: Level: beginner
2056: .seealso: MatCreateMPISBAIJ
2057: M*/
2059: PETSC_INTERN PetscErrorCode MatConvert_MPISBAIJ_MPISBSTRM(Mat,MatType,MatReuse,Mat*);
2063: PETSC_EXTERN PetscErrorCode MatCreate_MPISBAIJ(Mat B)
2064: {
2065: Mat_MPISBAIJ *b;
2067: PetscBool flg = PETSC_FALSE;
2070: PetscNewLog(B,&b);
2071: B->data = (void*)b;
2072: PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
2074: B->ops->destroy = MatDestroy_MPISBAIJ;
2075: B->ops->view = MatView_MPISBAIJ;
2076: B->assembled = PETSC_FALSE;
2077: B->insertmode = NOT_SET_VALUES;
2079: MPI_Comm_rank(PetscObjectComm((PetscObject)B),&b->rank);
2080: MPI_Comm_size(PetscObjectComm((PetscObject)B),&b->size);
2082: /* build local table of row and column ownerships */
2083: PetscMalloc1(b->size+2,&b->rangebs);
2085: /* build cache for off array entries formed */
2086: MatStashCreate_Private(PetscObjectComm((PetscObject)B),1,&B->stash);
2088: b->donotstash = PETSC_FALSE;
2089: b->colmap = NULL;
2090: b->garray = NULL;
2091: b->roworiented = PETSC_TRUE;
2093: /* stuff used in block assembly */
2094: b->barray = 0;
2096: /* stuff used for matrix vector multiply */
2097: b->lvec = 0;
2098: b->Mvctx = 0;
2099: b->slvec0 = 0;
2100: b->slvec0b = 0;
2101: b->slvec1 = 0;
2102: b->slvec1a = 0;
2103: b->slvec1b = 0;
2104: b->sMvctx = 0;
2106: /* stuff for MatGetRow() */
2107: b->rowindices = 0;
2108: b->rowvalues = 0;
2109: b->getrowactive = PETSC_FALSE;
2111: /* hash table stuff */
2112: b->ht = 0;
2113: b->hd = 0;
2114: b->ht_size = 0;
2115: b->ht_flag = PETSC_FALSE;
2116: b->ht_fact = 0;
2117: b->ht_total_ct = 0;
2118: b->ht_insert_ct = 0;
2120: /* stuff for MatGetSubMatrices_MPIBAIJ_local() */
2121: b->ijonly = PETSC_FALSE;
2123: b->in_loc = 0;
2124: b->v_loc = 0;
2125: b->n_loc = 0;
2127: PetscObjectComposeFunction((PetscObject)B,"MatStoreValues_C",MatStoreValues_MPISBAIJ);
2128: PetscObjectComposeFunction((PetscObject)B,"MatRetrieveValues_C",MatRetrieveValues_MPISBAIJ);
2129: PetscObjectComposeFunction((PetscObject)B,"MatGetDiagonalBlock_C",MatGetDiagonalBlock_MPISBAIJ);
2130: PetscObjectComposeFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",MatMPISBAIJSetPreallocation_MPISBAIJ);
2131: PetscObjectComposeFunction((PetscObject)B,"MatMPISBAIJSetPreallocationCSR_C",MatMPISBAIJSetPreallocationCSR_MPISBAIJ);
2132: PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpisbaij_mpisbstrm_C",MatConvert_MPISBAIJ_MPISBSTRM);
2133: #if defined(PETSC_HAVE_ELEMENTAL)
2134: PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpisbaij_elemental_C",MatConvert_MPISBAIJ_Elemental);
2135: #endif
2137: B->symmetric = PETSC_TRUE;
2138: B->structurally_symmetric = PETSC_TRUE;
2139: B->symmetric_set = PETSC_TRUE;
2140: B->structurally_symmetric_set = PETSC_TRUE;
2142: PetscObjectChangeTypeName((PetscObject)B,MATMPISBAIJ);
2143: PetscOptionsBegin(PetscObjectComm((PetscObject)B),NULL,"Options for loading MPISBAIJ matrix 1","Mat");
2144: PetscOptionsBool("-mat_use_hash_table","Use hash table to save memory in constructing matrix","MatSetOption",flg,&flg,NULL);
2145: if (flg) {
2146: PetscReal fact = 1.39;
2147: MatSetOption(B,MAT_USE_HASH_TABLE,PETSC_TRUE);
2148: PetscOptionsReal("-mat_use_hash_table","Use hash table factor","MatMPIBAIJSetHashTableFactor",fact,&fact,NULL);
2149: if (fact <= 1.0) fact = 1.39;
2150: MatMPIBAIJSetHashTableFactor(B,fact);
2151: PetscInfo1(B,"Hash table Factor used %5.2f\n",fact);
2152: }
2153: PetscOptionsEnd();
2154: return(0);
2155: }
2157: /*MC
2158: MATSBAIJ - MATSBAIJ = "sbaij" - A matrix type to be used for symmetric block sparse matrices.
2160: This matrix type is identical to MATSEQSBAIJ when constructed with a single process communicator,
2161: and MATMPISBAIJ otherwise.
2163: Options Database Keys:
2164: . -mat_type sbaij - sets the matrix type to "sbaij" during a call to MatSetFromOptions()
2166: Level: beginner
2168: .seealso: MatCreateMPISBAIJ,MATSEQSBAIJ,MATMPISBAIJ
2169: M*/
2173: /*@C
2174: MatMPISBAIJSetPreallocation - For good matrix assembly performance
2175: the user should preallocate the matrix storage by setting the parameters
2176: d_nz (or d_nnz) and o_nz (or o_nnz). By setting these parameters accurately,
2177: performance can be increased by more than a factor of 50.
2179: Collective on Mat
2181: Input Parameters:
2182: + B - the matrix
2183: . bs - size of block, the blocks are ALWAYS square. One can use MatSetBlockSizes() to set a different row and column blocksize but the row
2184: blocksize always defines the size of the blocks. The column blocksize sets the blocksize of the vectors obtained with MatCreateVecs()
2185: . d_nz - number of block nonzeros per block row in diagonal portion of local
2186: submatrix (same for all local rows)
2187: . d_nnz - array containing the number of block nonzeros in the various block rows
2188: in the upper triangular and diagonal part of the in diagonal portion of the local
2189: (possibly different for each block row) or NULL. If you plan to factor the matrix you must leave room
2190: for the diagonal entry and set a value even if it is zero.
2191: . o_nz - number of block nonzeros per block row in the off-diagonal portion of local
2192: submatrix (same for all local rows).
2193: - o_nnz - array containing the number of nonzeros in the various block rows of the
2194: off-diagonal portion of the local submatrix that is right of the diagonal
2195: (possibly different for each block row) or NULL.
2198: Options Database Keys:
2199: . -mat_no_unroll - uses code that does not unroll the loops in the
2200: block calculations (much slower)
2201: . -mat_block_size - size of the blocks to use
2203: Notes:
2205: If PETSC_DECIDE or PETSC_DETERMINE is used for a particular argument on one processor
2206: than it must be used on all processors that share the object for that argument.
2208: If the *_nnz parameter is given then the *_nz parameter is ignored
2210: Storage Information:
2211: For a square global matrix we define each processor's diagonal portion
2212: to be its local rows and the corresponding columns (a square submatrix);
2213: each processor's off-diagonal portion encompasses the remainder of the
2214: local matrix (a rectangular submatrix).
2216: The user can specify preallocated storage for the diagonal part of
2217: the local submatrix with either d_nz or d_nnz (not both). Set
2218: d_nz=PETSC_DEFAULT and d_nnz=NULL for PETSc to control dynamic
2219: memory allocation. Likewise, specify preallocated storage for the
2220: off-diagonal part of the local submatrix with o_nz or o_nnz (not both).
2222: You can call MatGetInfo() to get information on how effective the preallocation was;
2223: for example the fields mallocs,nz_allocated,nz_used,nz_unneeded;
2224: You can also run with the option -info and look for messages with the string
2225: malloc in them to see if additional memory allocation was needed.
2227: Consider a processor that owns rows 3, 4 and 5 of a parallel matrix. In
2228: the figure below we depict these three local rows and all columns (0-11).
2230: .vb
2231: 0 1 2 3 4 5 6 7 8 9 10 11
2232: --------------------------
2233: row 3 |. . . d d d o o o o o o
2234: row 4 |. . . d d d o o o o o o
2235: row 5 |. . . d d d o o o o o o
2236: --------------------------
2237: .ve
2239: Thus, any entries in the d locations are stored in the d (diagonal)
2240: submatrix, and any entries in the o locations are stored in the
2241: o (off-diagonal) submatrix. Note that the d matrix is stored in
2242: MatSeqSBAIJ format and the o submatrix in MATSEQBAIJ format.
2244: Now d_nz should indicate the number of block nonzeros per row in the upper triangular
2245: plus the diagonal part of the d matrix,
2246: and o_nz should indicate the number of block nonzeros per row in the o matrix
2248: In general, for PDE problems in which most nonzeros are near the diagonal,
2249: one expects d_nz >> o_nz. For large problems you MUST preallocate memory
2250: or you will get TERRIBLE performance; see the users' manual chapter on
2251: matrices.
2253: Level: intermediate
2255: .keywords: matrix, block, aij, compressed row, sparse, parallel
2257: .seealso: MatCreate(), MatCreateSeqSBAIJ(), MatSetValues(), MatCreateBAIJ(), PetscSplitOwnership()
2258: @*/
2259: PetscErrorCode MatMPISBAIJSetPreallocation(Mat B,PetscInt bs,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
2260: {
2267: PetscTryMethod(B,"MatMPISBAIJSetPreallocation_C",(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,bs,d_nz,d_nnz,o_nz,o_nnz));
2268: return(0);
2269: }
2273: /*@C
2274: MatCreateSBAIJ - Creates a sparse parallel matrix in symmetric block AIJ format
2275: (block compressed row). For good matrix assembly performance
2276: the user should preallocate the matrix storage by setting the parameters
2277: d_nz (or d_nnz) and o_nz (or o_nnz). By setting these parameters accurately,
2278: performance can be increased by more than a factor of 50.
2280: Collective on MPI_Comm
2282: Input Parameters:
2283: + comm - MPI communicator
2284: . bs - size of block, the blocks are ALWAYS square. One can use MatSetBlockSizes() to set a different row and column blocksize but the row
2285: blocksize always defines the size of the blocks. The column blocksize sets the blocksize of the vectors obtained with MatCreateVecs()
2286: . m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
2287: This value should be the same as the local size used in creating the
2288: y vector for the matrix-vector product y = Ax.
2289: . n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
2290: This value should be the same as the local size used in creating the
2291: x vector for the matrix-vector product y = Ax.
2292: . M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
2293: . N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
2294: . d_nz - number of block nonzeros per block row in diagonal portion of local
2295: submatrix (same for all local rows)
2296: . d_nnz - array containing the number of block nonzeros in the various block rows
2297: in the upper triangular portion of the in diagonal portion of the local
2298: (possibly different for each block block row) or NULL.
2299: If you plan to factor the matrix you must leave room for the diagonal entry and
2300: set its value even if it is zero.
2301: . o_nz - number of block nonzeros per block row in the off-diagonal portion of local
2302: submatrix (same for all local rows).
2303: - o_nnz - array containing the number of nonzeros in the various block rows of the
2304: off-diagonal portion of the local submatrix (possibly different for
2305: each block row) or NULL.
2307: Output Parameter:
2308: . A - the matrix
2310: Options Database Keys:
2311: . -mat_no_unroll - uses code that does not unroll the loops in the
2312: block calculations (much slower)
2313: . -mat_block_size - size of the blocks to use
2314: . -mat_mpi - use the parallel matrix data structures even on one processor
2315: (defaults to using SeqBAIJ format on one processor)
2317: It is recommended that one use the MatCreate(), MatSetType() and/or MatSetFromOptions(),
2318: MatXXXXSetPreallocation() paradgm instead of this routine directly.
2319: [MatXXXXSetPreallocation() is, for example, MatSeqAIJSetPreallocation]
2321: Notes:
2322: The number of rows and columns must be divisible by blocksize.
2323: This matrix type does not support complex Hermitian operation.
2325: The user MUST specify either the local or global matrix dimensions
2326: (possibly both).
2328: If PETSC_DECIDE or PETSC_DETERMINE is used for a particular argument on one processor
2329: than it must be used on all processors that share the object for that argument.
2331: If the *_nnz parameter is given then the *_nz parameter is ignored
2333: Storage Information:
2334: For a square global matrix we define each processor's diagonal portion
2335: to be its local rows and the corresponding columns (a square submatrix);
2336: each processor's off-diagonal portion encompasses the remainder of the
2337: local matrix (a rectangular submatrix).
2339: The user can specify preallocated storage for the diagonal part of
2340: the local submatrix with either d_nz or d_nnz (not both). Set
2341: d_nz=PETSC_DEFAULT and d_nnz=NULL for PETSc to control dynamic
2342: memory allocation. Likewise, specify preallocated storage for the
2343: off-diagonal part of the local submatrix with o_nz or o_nnz (not both).
2345: Consider a processor that owns rows 3, 4 and 5 of a parallel matrix. In
2346: the figure below we depict these three local rows and all columns (0-11).
2348: .vb
2349: 0 1 2 3 4 5 6 7 8 9 10 11
2350: --------------------------
2351: row 3 |. . . d d d o o o o o o
2352: row 4 |. . . d d d o o o o o o
2353: row 5 |. . . d d d o o o o o o
2354: --------------------------
2355: .ve
2357: Thus, any entries in the d locations are stored in the d (diagonal)
2358: submatrix, and any entries in the o locations are stored in the
2359: o (off-diagonal) submatrix. Note that the d matrix is stored in
2360: MatSeqSBAIJ format and the o submatrix in MATSEQBAIJ format.
2362: Now d_nz should indicate the number of block nonzeros per row in the upper triangular
2363: plus the diagonal part of the d matrix,
2364: and o_nz should indicate the number of block nonzeros per row in the o matrix.
2365: In general, for PDE problems in which most nonzeros are near the diagonal,
2366: one expects d_nz >> o_nz. For large problems you MUST preallocate memory
2367: or you will get TERRIBLE performance; see the users' manual chapter on
2368: matrices.
2370: Level: intermediate
2372: .keywords: matrix, block, aij, compressed row, sparse, parallel
2374: .seealso: MatCreate(), MatCreateSeqSBAIJ(), MatSetValues(), MatCreateBAIJ()
2375: @*/
2377: PetscErrorCode MatCreateSBAIJ(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[],Mat *A)
2378: {
2380: PetscMPIInt size;
2383: MatCreate(comm,A);
2384: MatSetSizes(*A,m,n,M,N);
2385: MPI_Comm_size(comm,&size);
2386: if (size > 1) {
2387: MatSetType(*A,MATMPISBAIJ);
2388: MatMPISBAIJSetPreallocation(*A,bs,d_nz,d_nnz,o_nz,o_nnz);
2389: } else {
2390: MatSetType(*A,MATSEQSBAIJ);
2391: MatSeqSBAIJSetPreallocation(*A,bs,d_nz,d_nnz);
2392: }
2393: return(0);
2394: }
2399: static PetscErrorCode MatDuplicate_MPISBAIJ(Mat matin,MatDuplicateOption cpvalues,Mat *newmat)
2400: {
2401: Mat mat;
2402: Mat_MPISBAIJ *a,*oldmat = (Mat_MPISBAIJ*)matin->data;
2404: PetscInt len=0,nt,bs=matin->rmap->bs,mbs=oldmat->mbs;
2405: PetscScalar *array;
2408: *newmat = 0;
2410: MatCreate(PetscObjectComm((PetscObject)matin),&mat);
2411: MatSetSizes(mat,matin->rmap->n,matin->cmap->n,matin->rmap->N,matin->cmap->N);
2412: MatSetType(mat,((PetscObject)matin)->type_name);
2413: PetscMemcpy(mat->ops,matin->ops,sizeof(struct _MatOps));
2414: PetscLayoutReference(matin->rmap,&mat->rmap);
2415: PetscLayoutReference(matin->cmap,&mat->cmap);
2417: mat->factortype = matin->factortype;
2418: mat->preallocated = PETSC_TRUE;
2419: mat->assembled = PETSC_TRUE;
2420: mat->insertmode = NOT_SET_VALUES;
2422: a = (Mat_MPISBAIJ*)mat->data;
2423: a->bs2 = oldmat->bs2;
2424: a->mbs = oldmat->mbs;
2425: a->nbs = oldmat->nbs;
2426: a->Mbs = oldmat->Mbs;
2427: a->Nbs = oldmat->Nbs;
2430: a->size = oldmat->size;
2431: a->rank = oldmat->rank;
2432: a->donotstash = oldmat->donotstash;
2433: a->roworiented = oldmat->roworiented;
2434: a->rowindices = 0;
2435: a->rowvalues = 0;
2436: a->getrowactive = PETSC_FALSE;
2437: a->barray = 0;
2438: a->rstartbs = oldmat->rstartbs;
2439: a->rendbs = oldmat->rendbs;
2440: a->cstartbs = oldmat->cstartbs;
2441: a->cendbs = oldmat->cendbs;
2443: /* hash table stuff */
2444: a->ht = 0;
2445: a->hd = 0;
2446: a->ht_size = 0;
2447: a->ht_flag = oldmat->ht_flag;
2448: a->ht_fact = oldmat->ht_fact;
2449: a->ht_total_ct = 0;
2450: a->ht_insert_ct = 0;
2452: PetscMemcpy(a->rangebs,oldmat->rangebs,(a->size+2)*sizeof(PetscInt));
2453: if (oldmat->colmap) {
2454: #if defined(PETSC_USE_CTABLE)
2455: PetscTableCreateCopy(oldmat->colmap,&a->colmap);
2456: #else
2457: PetscMalloc1(a->Nbs,&a->colmap);
2458: PetscLogObjectMemory((PetscObject)mat,(a->Nbs)*sizeof(PetscInt));
2459: PetscMemcpy(a->colmap,oldmat->colmap,(a->Nbs)*sizeof(PetscInt));
2460: #endif
2461: } else a->colmap = 0;
2463: if (oldmat->garray && (len = ((Mat_SeqBAIJ*)(oldmat->B->data))->nbs)) {
2464: PetscMalloc1(len,&a->garray);
2465: PetscLogObjectMemory((PetscObject)mat,len*sizeof(PetscInt));
2466: PetscMemcpy(a->garray,oldmat->garray,len*sizeof(PetscInt));
2467: } else a->garray = 0;
2469: MatStashCreate_Private(PetscObjectComm((PetscObject)matin),matin->rmap->bs,&mat->bstash);
2470: VecDuplicate(oldmat->lvec,&a->lvec);
2471: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->lvec);
2472: VecScatterCopy(oldmat->Mvctx,&a->Mvctx);
2473: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->Mvctx);
2475: VecDuplicate(oldmat->slvec0,&a->slvec0);
2476: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec0);
2477: VecDuplicate(oldmat->slvec1,&a->slvec1);
2478: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec1);
2480: VecGetLocalSize(a->slvec1,&nt);
2481: VecGetArray(a->slvec1,&array);
2482: VecCreateSeqWithArray(PETSC_COMM_SELF,1,bs*mbs,array,&a->slvec1a);
2483: VecCreateSeqWithArray(PETSC_COMM_SELF,1,nt-bs*mbs,array+bs*mbs,&a->slvec1b);
2484: VecRestoreArray(a->slvec1,&array);
2485: VecGetArray(a->slvec0,&array);
2486: VecCreateSeqWithArray(PETSC_COMM_SELF,1,nt-bs*mbs,array+bs*mbs,&a->slvec0b);
2487: VecRestoreArray(a->slvec0,&array);
2488: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec0);
2489: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec1);
2490: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec0b);
2491: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec1a);
2492: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec1b);
2494: /* VecScatterCopy(oldmat->sMvctx,&a->sMvctx); - not written yet, replaced by the lazy trick: */
2495: PetscObjectReference((PetscObject)oldmat->sMvctx);
2496: a->sMvctx = oldmat->sMvctx;
2497: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->sMvctx);
2499: MatDuplicate(oldmat->A,cpvalues,&a->A);
2500: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->A);
2501: MatDuplicate(oldmat->B,cpvalues,&a->B);
2502: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->B);
2503: PetscFunctionListDuplicate(((PetscObject)matin)->qlist,&((PetscObject)mat)->qlist);
2504: *newmat = mat;
2505: return(0);
2506: }
2510: PetscErrorCode MatLoad_MPISBAIJ(Mat newmat,PetscViewer viewer)
2511: {
2513: PetscInt i,nz,j,rstart,rend;
2514: PetscScalar *vals,*buf;
2515: MPI_Comm comm;
2516: MPI_Status status;
2517: PetscMPIInt rank,size,tag = ((PetscObject)viewer)->tag,*sndcounts = 0,*browners,maxnz,*rowners,mmbs;
2518: PetscInt header[4],*rowlengths = 0,M,N,m,*cols,*locrowlens;
2519: PetscInt *procsnz = 0,jj,*mycols,*ibuf;
2520: PetscInt bs = newmat->rmap->bs,Mbs,mbs,extra_rows;
2521: PetscInt *dlens,*odlens,*mask,*masked1,*masked2,rowcount,odcount;
2522: PetscInt dcount,kmax,k,nzcount,tmp;
2523: int fd;
2526: /* force binary viewer to load .info file if it has not yet done so */
2527: PetscViewerSetUp(viewer);
2528: PetscObjectGetComm((PetscObject)viewer,&comm);
2529: PetscOptionsBegin(comm,NULL,"Options for loading MPISBAIJ matrix 2","Mat");
2530: PetscOptionsInt("-matload_block_size","Set the blocksize used to store the matrix","MatLoad",bs,&bs,NULL);
2531: PetscOptionsEnd();
2532: if (bs < 0) bs = 1;
2534: MPI_Comm_size(comm,&size);
2535: MPI_Comm_rank(comm,&rank);
2536: PetscViewerBinaryGetDescriptor(viewer,&fd);
2537: if (!rank) {
2538: PetscBinaryRead(fd,(char*)header,4,PETSC_INT);
2539: if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"not matrix object");
2540: if (header[3] < 0) SETERRQ(PetscObjectComm((PetscObject)newmat),PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as MPISBAIJ");
2541: }
2543: MPI_Bcast(header+1,3,MPIU_INT,0,comm);
2544: M = header[1];
2545: N = header[2];
2547: /* If global sizes are set, check if they are consistent with that given in the file */
2548: if (newmat->rmap->N >= 0 && newmat->rmap->N != M) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Inconsistent # of rows:Matrix in file has (%D) and input matrix has (%D)",newmat->rmap->N,M);
2549: if (newmat->cmap->N >= 0 && newmat->cmap->N != N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Inconsistent # of cols:Matrix in file has (%D) and input matrix has (%D)",newmat->cmap->N,N);
2551: if (M != N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Can only do square matrices");
2553: /*
2554: This code adds extra rows to make sure the number of rows is
2555: divisible by the blocksize
2556: */
2557: Mbs = M/bs;
2558: extra_rows = bs - M + bs*(Mbs);
2559: if (extra_rows == bs) extra_rows = 0;
2560: else Mbs++;
2561: if (extra_rows &&!rank) {
2562: PetscInfo(viewer,"Padding loaded matrix to match blocksize\n");
2563: }
2565: /* determine ownership of all rows */
2566: if (newmat->rmap->n < 0) { /* PETSC_DECIDE */
2567: mbs = Mbs/size + ((Mbs % size) > rank);
2568: m = mbs*bs;
2569: } else { /* User Set */
2570: m = newmat->rmap->n;
2571: mbs = m/bs;
2572: }
2573: PetscMalloc2(size+1,&rowners,size+1,&browners);
2574: PetscMPIIntCast(mbs,&mmbs);
2575: MPI_Allgather(&mmbs,1,MPI_INT,rowners+1,1,MPI_INT,comm);
2576: rowners[0] = 0;
2577: for (i=2; i<=size; i++) rowners[i] += rowners[i-1];
2578: for (i=0; i<=size; i++) browners[i] = rowners[i]*bs;
2579: rstart = rowners[rank];
2580: rend = rowners[rank+1];
2582: /* distribute row lengths to all processors */
2583: PetscMalloc1((rend-rstart)*bs,&locrowlens);
2584: if (!rank) {
2585: PetscMalloc1(M+extra_rows,&rowlengths);
2586: PetscBinaryRead(fd,rowlengths,M,PETSC_INT);
2587: for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;
2588: PetscMalloc1(size,&sndcounts);
2589: for (i=0; i<size; i++) sndcounts[i] = browners[i+1] - browners[i];
2590: MPI_Scatterv(rowlengths,sndcounts,browners,MPIU_INT,locrowlens,(rend-rstart)*bs,MPIU_INT,0,comm);
2591: PetscFree(sndcounts);
2592: } else {
2593: MPI_Scatterv(0,0,0,MPIU_INT,locrowlens,(rend-rstart)*bs,MPIU_INT,0,comm);
2594: }
2596: if (!rank) { /* procs[0] */
2597: /* calculate the number of nonzeros on each processor */
2598: PetscMalloc1(size,&procsnz);
2599: PetscMemzero(procsnz,size*sizeof(PetscInt));
2600: for (i=0; i<size; i++) {
2601: for (j=rowners[i]*bs; j< rowners[i+1]*bs; j++) {
2602: procsnz[i] += rowlengths[j];
2603: }
2604: }
2605: PetscFree(rowlengths);
2607: /* determine max buffer needed and allocate it */
2608: maxnz = 0;
2609: for (i=0; i<size; i++) {
2610: maxnz = PetscMax(maxnz,procsnz[i]);
2611: }
2612: PetscMalloc1(maxnz,&cols);
2614: /* read in my part of the matrix column indices */
2615: nz = procsnz[0];
2616: PetscMalloc1(nz,&ibuf);
2617: mycols = ibuf;
2618: if (size == 1) nz -= extra_rows;
2619: PetscBinaryRead(fd,mycols,nz,PETSC_INT);
2620: if (size == 1) {
2621: for (i=0; i< extra_rows; i++) mycols[nz+i] = M+i;
2622: }
2624: /* read in every ones (except the last) and ship off */
2625: for (i=1; i<size-1; i++) {
2626: nz = procsnz[i];
2627: PetscBinaryRead(fd,cols,nz,PETSC_INT);
2628: MPI_Send(cols,nz,MPIU_INT,i,tag,comm);
2629: }
2630: /* read in the stuff for the last proc */
2631: if (size != 1) {
2632: nz = procsnz[size-1] - extra_rows; /* the extra rows are not on the disk */
2633: PetscBinaryRead(fd,cols,nz,PETSC_INT);
2634: for (i=0; i<extra_rows; i++) cols[nz+i] = M+i;
2635: MPI_Send(cols,nz+extra_rows,MPIU_INT,size-1,tag,comm);
2636: }
2637: PetscFree(cols);
2638: } else { /* procs[i], i>0 */
2639: /* determine buffer space needed for message */
2640: nz = 0;
2641: for (i=0; i<m; i++) nz += locrowlens[i];
2642: PetscMalloc1(nz,&ibuf);
2643: mycols = ibuf;
2644: /* receive message of column indices*/
2645: MPI_Recv(mycols,nz,MPIU_INT,0,tag,comm,&status);
2646: MPI_Get_count(&status,MPIU_INT,&maxnz);
2647: if (maxnz != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
2648: }
2650: /* loop over local rows, determining number of off diagonal entries */
2651: PetscMalloc2(rend-rstart,&dlens,rend-rstart,&odlens);
2652: PetscMalloc3(Mbs,&mask,Mbs,&masked1,Mbs,&masked2);
2653: PetscMemzero(mask,Mbs*sizeof(PetscInt));
2654: PetscMemzero(masked1,Mbs*sizeof(PetscInt));
2655: PetscMemzero(masked2,Mbs*sizeof(PetscInt));
2656: rowcount = 0;
2657: nzcount = 0;
2658: for (i=0; i<mbs; i++) {
2659: dcount = 0;
2660: odcount = 0;
2661: for (j=0; j<bs; j++) {
2662: kmax = locrowlens[rowcount];
2663: for (k=0; k<kmax; k++) {
2664: tmp = mycols[nzcount++]/bs; /* block col. index */
2665: if (!mask[tmp]) {
2666: mask[tmp] = 1;
2667: if (tmp < rstart || tmp >= rend) masked2[odcount++] = tmp; /* entry in off-diag portion */
2668: else masked1[dcount++] = tmp; /* entry in diag portion */
2669: }
2670: }
2671: rowcount++;
2672: }
2674: dlens[i] = dcount; /* d_nzz[i] */
2675: odlens[i] = odcount; /* o_nzz[i] */
2677: /* zero out the mask elements we set */
2678: for (j=0; j<dcount; j++) mask[masked1[j]] = 0;
2679: for (j=0; j<odcount; j++) mask[masked2[j]] = 0;
2680: }
2681: MatSetSizes(newmat,m,m,M+extra_rows,N+extra_rows);
2682: MatMPISBAIJSetPreallocation(newmat,bs,0,dlens,0,odlens);
2683: MatSetOption(newmat,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);
2685: if (!rank) {
2686: PetscMalloc1(maxnz,&buf);
2687: /* read in my part of the matrix numerical values */
2688: nz = procsnz[0];
2689: vals = buf;
2690: mycols = ibuf;
2691: if (size == 1) nz -= extra_rows;
2692: PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
2693: if (size == 1) {
2694: for (i=0; i< extra_rows; i++) vals[nz+i] = 1.0;
2695: }
2697: /* insert into matrix */
2698: jj = rstart*bs;
2699: for (i=0; i<m; i++) {
2700: MatSetValues(newmat,1,&jj,locrowlens[i],mycols,vals,INSERT_VALUES);
2701: mycols += locrowlens[i];
2702: vals += locrowlens[i];
2703: jj++;
2704: }
2706: /* read in other processors (except the last one) and ship out */
2707: for (i=1; i<size-1; i++) {
2708: nz = procsnz[i];
2709: vals = buf;
2710: PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
2711: MPI_Send(vals,nz,MPIU_SCALAR,i,((PetscObject)newmat)->tag,comm);
2712: }
2713: /* the last proc */
2714: if (size != 1) {
2715: nz = procsnz[i] - extra_rows;
2716: vals = buf;
2717: PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
2718: for (i=0; i<extra_rows; i++) vals[nz+i] = 1.0;
2719: MPI_Send(vals,nz+extra_rows,MPIU_SCALAR,size-1,((PetscObject)newmat)->tag,comm);
2720: }
2721: PetscFree(procsnz);
2723: } else {
2724: /* receive numeric values */
2725: PetscMalloc1(nz,&buf);
2727: /* receive message of values*/
2728: vals = buf;
2729: mycols = ibuf;
2730: MPI_Recv(vals,nz,MPIU_SCALAR,0,((PetscObject)newmat)->tag,comm,&status);
2731: MPI_Get_count(&status,MPIU_SCALAR,&maxnz);
2732: if (maxnz != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
2734: /* insert into matrix */
2735: jj = rstart*bs;
2736: for (i=0; i<m; i++) {
2737: MatSetValues_MPISBAIJ(newmat,1,&jj,locrowlens[i],mycols,vals,INSERT_VALUES);
2738: mycols += locrowlens[i];
2739: vals += locrowlens[i];
2740: jj++;
2741: }
2742: }
2744: PetscFree(locrowlens);
2745: PetscFree(buf);
2746: PetscFree(ibuf);
2747: PetscFree2(rowners,browners);
2748: PetscFree2(dlens,odlens);
2749: PetscFree3(mask,masked1,masked2);
2750: MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);
2751: MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);
2752: return(0);
2753: }
2757: /*XXXXX@
2758: MatMPISBAIJSetHashTableFactor - Sets the factor required to compute the size of the HashTable.
2760: Input Parameters:
2761: . mat - the matrix
2762: . fact - factor
2764: Not Collective on Mat, each process can have a different hash factor
2766: Level: advanced
2768: Notes:
2769: This can also be set by the command line option: -mat_use_hash_table fact
2771: .keywords: matrix, hashtable, factor, HT
2773: .seealso: MatSetOption()
2774: @XXXXX*/
2779: PetscErrorCode MatGetRowMaxAbs_MPISBAIJ(Mat A,Vec v,PetscInt idx[])
2780: {
2781: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
2782: Mat_SeqBAIJ *b = (Mat_SeqBAIJ*)(a->B)->data;
2783: PetscReal atmp;
2784: PetscReal *work,*svalues,*rvalues;
2786: PetscInt i,bs,mbs,*bi,*bj,brow,j,ncols,krow,kcol,col,row,Mbs,bcol;
2787: PetscMPIInt rank,size;
2788: PetscInt *rowners_bs,dest,count,source;
2789: PetscScalar *va;
2790: MatScalar *ba;
2791: MPI_Status stat;
2794: if (idx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Send email to petsc-maint@mcs.anl.gov");
2795: MatGetRowMaxAbs(a->A,v,NULL);
2796: VecGetArray(v,&va);
2798: MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);
2799: MPI_Comm_rank(PetscObjectComm((PetscObject)A),&rank);
2801: bs = A->rmap->bs;
2802: mbs = a->mbs;
2803: Mbs = a->Mbs;
2804: ba = b->a;
2805: bi = b->i;
2806: bj = b->j;
2808: /* find ownerships */
2809: rowners_bs = A->rmap->range;
2811: /* each proc creates an array to be distributed */
2812: PetscMalloc1(bs*Mbs,&work);
2813: PetscMemzero(work,bs*Mbs*sizeof(PetscReal));
2815: /* row_max for B */
2816: if (rank != size-1) {
2817: for (i=0; i<mbs; i++) {
2818: ncols = bi[1] - bi[0]; bi++;
2819: brow = bs*i;
2820: for (j=0; j<ncols; j++) {
2821: bcol = bs*(*bj);
2822: for (kcol=0; kcol<bs; kcol++) {
2823: col = bcol + kcol; /* local col index */
2824: col += rowners_bs[rank+1]; /* global col index */
2825: for (krow=0; krow<bs; krow++) {
2826: atmp = PetscAbsScalar(*ba); ba++;
2827: row = brow + krow; /* local row index */
2828: if (PetscRealPart(va[row]) < atmp) va[row] = atmp;
2829: if (work[col] < atmp) work[col] = atmp;
2830: }
2831: }
2832: bj++;
2833: }
2834: }
2836: /* send values to its owners */
2837: for (dest=rank+1; dest<size; dest++) {
2838: svalues = work + rowners_bs[dest];
2839: count = rowners_bs[dest+1]-rowners_bs[dest];
2840: MPI_Send(svalues,count,MPIU_REAL,dest,rank,PetscObjectComm((PetscObject)A));
2841: }
2842: }
2844: /* receive values */
2845: if (rank) {
2846: rvalues = work;
2847: count = rowners_bs[rank+1]-rowners_bs[rank];
2848: for (source=0; source<rank; source++) {
2849: MPI_Recv(rvalues,count,MPIU_REAL,MPI_ANY_SOURCE,MPI_ANY_TAG,PetscObjectComm((PetscObject)A),&stat);
2850: /* process values */
2851: for (i=0; i<count; i++) {
2852: if (PetscRealPart(va[i]) < rvalues[i]) va[i] = rvalues[i];
2853: }
2854: }
2855: }
2857: VecRestoreArray(v,&va);
2858: PetscFree(work);
2859: return(0);
2860: }
2864: PetscErrorCode MatSOR_MPISBAIJ(Mat matin,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
2865: {
2866: Mat_MPISBAIJ *mat = (Mat_MPISBAIJ*)matin->data;
2867: PetscErrorCode ierr;
2868: PetscInt mbs=mat->mbs,bs=matin->rmap->bs;
2869: PetscScalar *x,*ptr,*from;
2870: Vec bb1;
2871: const PetscScalar *b;
2874: if (its <= 0 || lits <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
2875: if (bs > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented");
2877: if (flag == SOR_APPLY_UPPER) {
2878: (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx);
2879: return(0);
2880: }
2882: if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP) {
2883: if (flag & SOR_ZERO_INITIAL_GUESS) {
2884: (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,lits,xx);
2885: its--;
2886: }
2888: VecDuplicate(bb,&bb1);
2889: while (its--) {
2891: /* lower triangular part: slvec0b = - B^T*xx */
2892: (*mat->B->ops->multtranspose)(mat->B,xx,mat->slvec0b);
2894: /* copy xx into slvec0a */
2895: VecGetArray(mat->slvec0,&ptr);
2896: VecGetArray(xx,&x);
2897: PetscMemcpy(ptr,x,bs*mbs*sizeof(MatScalar));
2898: VecRestoreArray(mat->slvec0,&ptr);
2900: VecScale(mat->slvec0,-1.0);
2902: /* copy bb into slvec1a */
2903: VecGetArray(mat->slvec1,&ptr);
2904: VecGetArrayRead(bb,&b);
2905: PetscMemcpy(ptr,b,bs*mbs*sizeof(MatScalar));
2906: VecRestoreArray(mat->slvec1,&ptr);
2908: /* set slvec1b = 0 */
2909: VecSet(mat->slvec1b,0.0);
2911: VecScatterBegin(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD);
2912: VecRestoreArray(xx,&x);
2913: VecRestoreArrayRead(bb,&b);
2914: VecScatterEnd(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD);
2916: /* upper triangular part: bb1 = bb1 - B*x */
2917: (*mat->B->ops->multadd)(mat->B,mat->slvec1b,mat->slvec1a,bb1);
2919: /* local diagonal sweep */
2920: (*mat->A->ops->sor)(mat->A,bb1,omega,SOR_SYMMETRIC_SWEEP,fshift,lits,lits,xx);
2921: }
2922: VecDestroy(&bb1);
2923: } else if ((flag & SOR_LOCAL_FORWARD_SWEEP) && (its == 1) && (flag & SOR_ZERO_INITIAL_GUESS)) {
2924: (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx);
2925: } else if ((flag & SOR_LOCAL_BACKWARD_SWEEP) && (its == 1) && (flag & SOR_ZERO_INITIAL_GUESS)) {
2926: (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx);
2927: } else if (flag & SOR_EISENSTAT) {
2928: Vec xx1;
2929: PetscBool hasop;
2930: const PetscScalar *diag;
2931: PetscScalar *sl,scale = (omega - 2.0)/omega;
2932: PetscInt i,n;
2934: if (!mat->xx1) {
2935: VecDuplicate(bb,&mat->xx1);
2936: VecDuplicate(bb,&mat->bb1);
2937: }
2938: xx1 = mat->xx1;
2939: bb1 = mat->bb1;
2941: (*mat->A->ops->sor)(mat->A,bb,omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_BACKWARD_SWEEP),fshift,lits,1,xx);
2943: if (!mat->diag) {
2944: /* this is wrong for same matrix with new nonzero values */
2945: MatCreateVecs(matin,&mat->diag,NULL);
2946: MatGetDiagonal(matin,mat->diag);
2947: }
2948: MatHasOperation(matin,MATOP_MULT_DIAGONAL_BLOCK,&hasop);
2950: if (hasop) {
2951: MatMultDiagonalBlock(matin,xx,bb1);
2952: VecAYPX(mat->slvec1a,scale,bb);
2953: } else {
2954: /*
2955: These two lines are replaced by code that may be a bit faster for a good compiler
2956: VecPointwiseMult(mat->slvec1a,mat->diag,xx);
2957: VecAYPX(mat->slvec1a,scale,bb);
2958: */
2959: VecGetArray(mat->slvec1a,&sl);
2960: VecGetArrayRead(mat->diag,&diag);
2961: VecGetArrayRead(bb,&b);
2962: VecGetArray(xx,&x);
2963: VecGetLocalSize(xx,&n);
2964: if (omega == 1.0) {
2965: for (i=0; i<n; i++) sl[i] = b[i] - diag[i]*x[i];
2966: PetscLogFlops(2.0*n);
2967: } else {
2968: for (i=0; i<n; i++) sl[i] = b[i] + scale*diag[i]*x[i];
2969: PetscLogFlops(3.0*n);
2970: }
2971: VecRestoreArray(mat->slvec1a,&sl);
2972: VecRestoreArrayRead(mat->diag,&diag);
2973: VecRestoreArrayRead(bb,&b);
2974: VecRestoreArray(xx,&x);
2975: }
2977: /* multiply off-diagonal portion of matrix */
2978: VecSet(mat->slvec1b,0.0);
2979: (*mat->B->ops->multtranspose)(mat->B,xx,mat->slvec0b);
2980: VecGetArray(mat->slvec0,&from);
2981: VecGetArray(xx,&x);
2982: PetscMemcpy(from,x,bs*mbs*sizeof(MatScalar));
2983: VecRestoreArray(mat->slvec0,&from);
2984: VecRestoreArray(xx,&x);
2985: VecScatterBegin(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD);
2986: VecScatterEnd(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD);
2987: (*mat->B->ops->multadd)(mat->B,mat->slvec1b,mat->slvec1a,mat->slvec1a);
2989: /* local sweep */
2990: (*mat->A->ops->sor)(mat->A,mat->slvec1a,omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_FORWARD_SWEEP),fshift,lits,1,xx1);
2991: VecAXPY(xx,1.0,xx1);
2992: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatSORType is not supported for SBAIJ matrix format");
2993: return(0);
2994: }
2998: PetscErrorCode MatSOR_MPISBAIJ_2comm(Mat matin,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
2999: {
3000: Mat_MPISBAIJ *mat = (Mat_MPISBAIJ*)matin->data;
3002: Vec lvec1,bb1;
3005: if (its <= 0 || lits <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
3006: if (matin->rmap->bs > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented");
3008: if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP) {
3009: if (flag & SOR_ZERO_INITIAL_GUESS) {
3010: (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,lits,xx);
3011: its--;
3012: }
3014: VecDuplicate(mat->lvec,&lvec1);
3015: VecDuplicate(bb,&bb1);
3016: while (its--) {
3017: VecScatterBegin(mat->Mvctx,xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD);
3019: /* lower diagonal part: bb1 = bb - B^T*xx */
3020: (*mat->B->ops->multtranspose)(mat->B,xx,lvec1);
3021: VecScale(lvec1,-1.0);
3023: VecScatterEnd(mat->Mvctx,xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD);
3024: VecCopy(bb,bb1);
3025: VecScatterBegin(mat->Mvctx,lvec1,bb1,ADD_VALUES,SCATTER_REVERSE);
3027: /* upper diagonal part: bb1 = bb1 - B*x */
3028: VecScale(mat->lvec,-1.0);
3029: (*mat->B->ops->multadd)(mat->B,mat->lvec,bb1,bb1);
3031: VecScatterEnd(mat->Mvctx,lvec1,bb1,ADD_VALUES,SCATTER_REVERSE);
3033: /* diagonal sweep */
3034: (*mat->A->ops->sor)(mat->A,bb1,omega,SOR_SYMMETRIC_SWEEP,fshift,lits,lits,xx);
3035: }
3036: VecDestroy(&lvec1);
3037: VecDestroy(&bb1);
3038: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatSORType is not supported for SBAIJ matrix format");
3039: return(0);
3040: }
3044: /*@
3045: MatCreateMPISBAIJWithArrays - creates a MPI SBAIJ matrix using arrays that contain in standard
3046: CSR format the local rows.
3048: Collective on MPI_Comm
3050: Input Parameters:
3051: + comm - MPI communicator
3052: . bs - the block size, only a block size of 1 is supported
3053: . m - number of local rows (Cannot be PETSC_DECIDE)
3054: . n - This value should be the same as the local size used in creating the
3055: x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have
3056: calculated if N is given) For square matrices n is almost always m.
3057: . M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
3058: . N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
3059: . i - row indices
3060: . j - column indices
3061: - a - matrix values
3063: Output Parameter:
3064: . mat - the matrix
3066: Level: intermediate
3068: Notes:
3069: The i, j, and a arrays ARE copied by this routine into the internal format used by PETSc;
3070: thus you CANNOT change the matrix entries by changing the values of a[] after you have
3071: called this routine. Use MatCreateMPIAIJWithSplitArrays() to avoid needing to copy the arrays.
3073: The i and j indices are 0 based, and i indices are indices corresponding to the local j array.
3075: .keywords: matrix, aij, compressed row, sparse, parallel
3077: .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatMPIAIJSetPreallocation(), MatMPIAIJSetPreallocationCSR(),
3078: MPIAIJ, MatCreateAIJ(), MatCreateMPIAIJWithSplitArrays()
3079: @*/
3080: PetscErrorCode MatCreateMPISBAIJWithArrays(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt i[],const PetscInt j[],const PetscScalar a[],Mat *mat)
3081: {
3086: if (i[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"i (row indices) must start with 0");
3087: if (m < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"local number of rows (m) cannot be PETSC_DECIDE, or negative");
3088: MatCreate(comm,mat);
3089: MatSetSizes(*mat,m,n,M,N);
3090: MatSetType(*mat,MATMPISBAIJ);
3091: MatMPISBAIJSetPreallocationCSR(*mat,bs,i,j,a);
3092: return(0);
3093: }
3098: /*@C
3099: MatMPISBAIJSetPreallocationCSR - Allocates memory for a sparse parallel matrix in BAIJ format
3100: (the default parallel PETSc format).
3102: Collective on MPI_Comm
3104: Input Parameters:
3105: + B - the matrix
3106: . bs - the block size
3107: . i - the indices into j for the start of each local row (starts with zero)
3108: . j - the column indices for each local row (starts with zero) these must be sorted for each row
3109: - v - optional values in the matrix
3111: Level: developer
3113: .keywords: matrix, aij, compressed row, sparse, parallel
3115: .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatMPIBAIJSetPreallocation(), MatCreateAIJ(), MPIAIJ
3116: @*/
3117: PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat B,PetscInt bs,const PetscInt i[],const PetscInt j[], const PetscScalar v[])
3118: {
3122: PetscTryMethod(B,"MatMPISBAIJSetPreallocationCSR_C",(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]),(B,bs,i,j,v));
3123: return(0);
3124: }
3128: PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPISBAIJ(MPI_Comm comm,Mat inmat,PetscInt n,MatReuse scall,Mat *outmat)
3129: {
3131: PetscInt m,N,i,rstart,nnz,Ii,bs,cbs;
3132: PetscInt *indx;
3133: PetscScalar *values;
3136: MatGetSize(inmat,&m,&N);
3137: if (scall == MAT_INITIAL_MATRIX) { /* symbolic phase */
3138: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inmat->data;
3139: PetscInt *dnz,*onz,sum,bs,cbs,mbs,Nbs;
3140: PetscInt *bindx,rmax=a->rmax,j;
3141:
3142: MatGetBlockSizes(inmat,&bs,&cbs);
3143: mbs = m/bs; Nbs = N/cbs;
3144: if (n == PETSC_DECIDE) {
3145: PetscSplitOwnership(comm,&n,&Nbs);
3146: }
3147: /* Check sum(n) = Nbs */
3148: MPIU_Allreduce(&n,&sum,1,MPIU_INT,MPI_SUM,comm);
3149: if (sum != Nbs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Sum of local columns != global columns %d",Nbs);
3151: MPI_Scan(&mbs, &rstart,1,MPIU_INT,MPI_SUM,comm);
3152: rstart -= mbs;
3154: PetscMalloc1(rmax,&bindx);
3155: MatPreallocateInitialize(comm,mbs,n,dnz,onz);
3156: MatSetOption(inmat,MAT_GETROW_UPPERTRIANGULAR,PETSC_TRUE);
3157: for (i=0; i<mbs; i++) {
3158: MatGetRow_SeqSBAIJ(inmat,i*bs,&nnz,&indx,NULL); /* non-blocked nnz and indx */
3159: nnz = nnz/bs;
3160: for (j=0; j<nnz; j++) bindx[j] = indx[j*bs]/bs;
3161: MatPreallocateSet(i+rstart,nnz,bindx,dnz,onz);
3162: MatRestoreRow_SeqSBAIJ(inmat,i*bs,&nnz,&indx,NULL);
3163: }
3164: MatSetOption(inmat,MAT_GETROW_UPPERTRIANGULAR,PETSC_FALSE);
3165: PetscFree(bindx);
3167: MatCreate(comm,outmat);
3168: MatSetSizes(*outmat,m,n*bs,PETSC_DETERMINE,PETSC_DETERMINE);
3169: MatSetBlockSizes(*outmat,bs,cbs);
3170: MatSetType(*outmat,MATMPISBAIJ);
3171: MatMPISBAIJSetPreallocation(*outmat,bs,0,dnz,0,onz);
3172: MatPreallocateFinalize(dnz,onz);
3173: }
3174:
3175: /* numeric phase */
3176: MatGetBlockSizes(inmat,&bs,&cbs);
3177: MatGetOwnershipRange(*outmat,&rstart,NULL);
3179: MatSetOption(inmat,MAT_GETROW_UPPERTRIANGULAR,PETSC_TRUE);
3180: for (i=0; i<m; i++) {
3181: MatGetRow_SeqSBAIJ(inmat,i,&nnz,&indx,&values);
3182: Ii = i + rstart;
3183: MatSetValues(*outmat,1,&Ii,nnz,indx,values,INSERT_VALUES);
3184: MatRestoreRow_SeqSBAIJ(inmat,i,&nnz,&indx,&values);
3185: }
3186: MatSetOption(inmat,MAT_GETROW_UPPERTRIANGULAR,PETSC_FALSE);
3187: MatAssemblyBegin(*outmat,MAT_FINAL_ASSEMBLY);
3188: MatAssemblyEnd(*outmat,MAT_FINAL_ASSEMBLY);
3189: return(0);
3190: }