Actual source code: mpisbaij.c
petsc-3.7.0 2016-04-25
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: }
637: /* Bmat */
638: v = bmat->a; jj = bmat->j;
639: for (brow=0; brow<mbs; brow++) {
640: grow = bs*(rstart + brow);
641: nz = bmat->i[brow+1] - bmat->i[brow];
642: for (bcol=0; bcol<nz; bcol++) {
643: gcol = bs*garray[*jj]; jj++;
644: for (col=0; col<bs; col++) {
645: for (row=0; row<bs; row++) {
646: vabs = PetscAbsScalar(*v); v++;
647: rsum[gcol+col] += vabs;
648: rsum[grow+row] += vabs;
649: }
650: }
651: }
652: }
653: MPIU_Allreduce(rsum,rsum2,mat->cmap->N,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)mat));
654: *norm = 0.0;
655: for (col=0; col<mat->cmap->N; col++) {
656: if (rsum2[col] > *norm) *norm = rsum2[col];
657: }
658: PetscFree2(rsum,rsum2);
659: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for this norm yet");
660: }
661: return(0);
662: }
666: PetscErrorCode MatAssemblyBegin_MPISBAIJ(Mat mat,MatAssemblyType mode)
667: {
668: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
670: PetscInt nstash,reallocs;
673: if (baij->donotstash || mat->nooffprocentries) return(0);
675: MatStashScatterBegin_Private(mat,&mat->stash,mat->rmap->range);
676: MatStashScatterBegin_Private(mat,&mat->bstash,baij->rangebs);
677: MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);
678: PetscInfo2(mat,"Stash has %D entries,uses %D mallocs.\n",nstash,reallocs);
679: MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);
680: PetscInfo2(mat,"Block-Stash has %D entries, uses %D mallocs.\n",nstash,reallocs);
681: return(0);
682: }
686: PetscErrorCode MatAssemblyEnd_MPISBAIJ(Mat mat,MatAssemblyType mode)
687: {
688: Mat_MPISBAIJ *baij=(Mat_MPISBAIJ*)mat->data;
689: Mat_SeqSBAIJ *a =(Mat_SeqSBAIJ*)baij->A->data;
691: PetscInt i,j,rstart,ncols,flg,bs2=baij->bs2;
692: PetscInt *row,*col;
693: PetscBool other_disassembled;
694: PetscMPIInt n;
695: PetscBool r1,r2,r3;
696: MatScalar *val;
698: /* do not use 'b=(Mat_SeqBAIJ*)baij->B->data' as B can be reset in disassembly */
700: if (!baij->donotstash && !mat->nooffprocentries) {
701: while (1) {
702: MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);
703: if (!flg) break;
705: for (i=0; i<n;) {
706: /* Now identify the consecutive vals belonging to the same row */
707: for (j=i,rstart=row[j]; j<n; j++) {
708: if (row[j] != rstart) break;
709: }
710: if (j < n) ncols = j-i;
711: else ncols = n-i;
712: /* Now assemble all these values with a single function call */
713: MatSetValues_MPISBAIJ(mat,1,row+i,ncols,col+i,val+i,mat->insertmode);
714: i = j;
715: }
716: }
717: MatStashScatterEnd_Private(&mat->stash);
718: /* Now process the block-stash. Since the values are stashed column-oriented,
719: set the roworiented flag to column oriented, and after MatSetValues()
720: restore the original flags */
721: r1 = baij->roworiented;
722: r2 = a->roworiented;
723: r3 = ((Mat_SeqBAIJ*)baij->B->data)->roworiented;
725: baij->roworiented = PETSC_FALSE;
726: a->roworiented = PETSC_FALSE;
728: ((Mat_SeqBAIJ*)baij->B->data)->roworiented = PETSC_FALSE; /* b->roworinted */
729: while (1) {
730: MatStashScatterGetMesg_Private(&mat->bstash,&n,&row,&col,&val,&flg);
731: if (!flg) break;
733: for (i=0; i<n;) {
734: /* Now identify the consecutive vals belonging to the same row */
735: for (j=i,rstart=row[j]; j<n; j++) {
736: if (row[j] != rstart) break;
737: }
738: if (j < n) ncols = j-i;
739: else ncols = n-i;
740: MatSetValuesBlocked_MPISBAIJ(mat,1,row+i,ncols,col+i,val+i*bs2,mat->insertmode);
741: i = j;
742: }
743: }
744: MatStashScatterEnd_Private(&mat->bstash);
746: baij->roworiented = r1;
747: a->roworiented = r2;
749: ((Mat_SeqBAIJ*)baij->B->data)->roworiented = r3; /* b->roworinted */
750: }
752: MatAssemblyBegin(baij->A,mode);
753: MatAssemblyEnd(baij->A,mode);
755: /* determine if any processor has disassembled, if so we must
756: also disassemble ourselfs, in order that we may reassemble. */
757: /*
758: if nonzero structure of submatrix B cannot change then we know that
759: no processor disassembled thus we can skip this stuff
760: */
761: if (!((Mat_SeqBAIJ*)baij->B->data)->nonew) {
762: MPIU_Allreduce(&mat->was_assembled,&other_disassembled,1,MPIU_BOOL,MPI_PROD,PetscObjectComm((PetscObject)mat));
763: if (mat->was_assembled && !other_disassembled) {
764: MatDisAssemble_MPISBAIJ(mat);
765: }
766: }
768: if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
769: MatSetUpMultiply_MPISBAIJ(mat); /* setup Mvctx and sMvctx */
770: }
771: MatAssemblyBegin(baij->B,mode);
772: MatAssemblyEnd(baij->B,mode);
774: PetscFree2(baij->rowvalues,baij->rowindices);
776: baij->rowvalues = 0;
778: /* if no new nonzero locations are allowed in matrix then only set the matrix state the first time through */
779: if ((!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) || !((Mat_SeqBAIJ*)(baij->A->data))->nonew) {
780: PetscObjectState state = baij->A->nonzerostate + baij->B->nonzerostate;
781: MPIU_Allreduce(&state,&mat->nonzerostate,1,MPIU_INT64,MPI_SUM,PetscObjectComm((PetscObject)mat));
782: }
783: return(0);
784: }
786: extern PetscErrorCode MatView_SeqSBAIJ(Mat,PetscViewer);
787: extern PetscErrorCode MatSetValues_MPIBAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
788: #include <petscdraw.h>
791: static PetscErrorCode MatView_MPISBAIJ_ASCIIorDraworSocket(Mat mat,PetscViewer viewer)
792: {
793: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
794: PetscErrorCode ierr;
795: PetscInt bs = mat->rmap->bs;
796: PetscMPIInt rank = baij->rank;
797: PetscBool iascii,isdraw;
798: PetscViewer sviewer;
799: PetscViewerFormat format;
802: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
803: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
804: if (iascii) {
805: PetscViewerGetFormat(viewer,&format);
806: if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
807: MatInfo info;
808: MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&rank);
809: MatGetInfo(mat,MAT_LOCAL,&info);
810: PetscViewerASCIIPushSynchronized(viewer);
811: 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);
812: MatGetInfo(baij->A,MAT_LOCAL,&info);
813: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] on-diagonal part: nz %D \n",rank,(PetscInt)info.nz_used);
814: MatGetInfo(baij->B,MAT_LOCAL,&info);
815: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] off-diagonal part: nz %D \n",rank,(PetscInt)info.nz_used);
816: PetscViewerFlush(viewer);
817: PetscViewerASCIIPopSynchronized(viewer);
818: PetscViewerASCIIPrintf(viewer,"Information on VecScatter used in matrix-vector product: \n");
819: VecScatterView(baij->Mvctx,viewer);
820: return(0);
821: } else if (format == PETSC_VIEWER_ASCII_INFO) {
822: PetscViewerASCIIPrintf(viewer," block size is %D\n",bs);
823: return(0);
824: } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
825: return(0);
826: }
827: }
829: if (isdraw) {
830: PetscDraw draw;
831: PetscBool isnull;
832: PetscViewerDrawGetDraw(viewer,0,&draw);
833: PetscDrawIsNull(draw,&isnull);
834: if (isnull) return(0);
835: }
837: {
838: /* assemble the entire matrix onto first processor. */
839: Mat A;
840: Mat_SeqSBAIJ *Aloc;
841: Mat_SeqBAIJ *Bloc;
842: PetscInt M = mat->rmap->N,N = mat->cmap->N,*ai,*aj,col,i,j,k,*rvals,mbs = baij->mbs;
843: MatScalar *a;
844: const char *matname;
846: /* Should this be the same type as mat? */
847: MatCreate(PetscObjectComm((PetscObject)mat),&A);
848: if (!rank) {
849: MatSetSizes(A,M,N,M,N);
850: } else {
851: MatSetSizes(A,0,0,M,N);
852: }
853: MatSetType(A,MATMPISBAIJ);
854: MatMPISBAIJSetPreallocation(A,mat->rmap->bs,0,NULL,0,NULL);
855: MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);
856: PetscLogObjectParent((PetscObject)mat,(PetscObject)A);
858: /* copy over the A part */
859: Aloc = (Mat_SeqSBAIJ*)baij->A->data;
860: ai = Aloc->i; aj = Aloc->j; a = Aloc->a;
861: PetscMalloc1(bs,&rvals);
863: for (i=0; i<mbs; i++) {
864: rvals[0] = bs*(baij->rstartbs + i);
865: for (j=1; j<bs; j++) rvals[j] = rvals[j-1] + 1;
866: for (j=ai[i]; j<ai[i+1]; j++) {
867: col = (baij->cstartbs+aj[j])*bs;
868: for (k=0; k<bs; k++) {
869: MatSetValues_MPISBAIJ(A,bs,rvals,1,&col,a,INSERT_VALUES);
870: col++;
871: a += bs;
872: }
873: }
874: }
875: /* copy over the B part */
876: Bloc = (Mat_SeqBAIJ*)baij->B->data;
877: ai = Bloc->i; aj = Bloc->j; a = Bloc->a;
878: for (i=0; i<mbs; i++) {
880: rvals[0] = bs*(baij->rstartbs + i);
881: for (j=1; j<bs; j++) rvals[j] = rvals[j-1] + 1;
882: for (j=ai[i]; j<ai[i+1]; j++) {
883: col = baij->garray[aj[j]]*bs;
884: for (k=0; k<bs; k++) {
885: MatSetValues_MPIBAIJ(A,bs,rvals,1,&col,a,INSERT_VALUES);
886: col++;
887: a += bs;
888: }
889: }
890: }
891: PetscFree(rvals);
892: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
893: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
894: /*
895: Everyone has to call to draw the matrix since the graphics waits are
896: synchronized across all processors that share the PetscDraw object
897: */
898: PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);
899: PetscObjectGetName((PetscObject)mat,&matname);
900: if (!rank) {
901: PetscObjectSetName((PetscObject)((Mat_MPISBAIJ*)(A->data))->A,matname);
902: MatView_SeqSBAIJ(((Mat_MPISBAIJ*)(A->data))->A,sviewer);
903: }
904: PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);
905: PetscViewerFlush(viewer);
906: MatDestroy(&A);
907: }
908: return(0);
909: }
913: PetscErrorCode MatView_MPISBAIJ(Mat mat,PetscViewer viewer)
914: {
916: PetscBool iascii,isdraw,issocket,isbinary;
919: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
920: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
921: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSOCKET,&issocket);
922: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
923: if (iascii || isdraw || issocket || isbinary) {
924: MatView_MPISBAIJ_ASCIIorDraworSocket(mat,viewer);
925: }
926: return(0);
927: }
931: PetscErrorCode MatDestroy_MPISBAIJ(Mat mat)
932: {
933: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
937: #if defined(PETSC_USE_LOG)
938: PetscLogObjectState((PetscObject)mat,"Rows=%D,Cols=%D",mat->rmap->N,mat->cmap->N);
939: #endif
940: MatStashDestroy_Private(&mat->stash);
941: MatStashDestroy_Private(&mat->bstash);
942: MatDestroy(&baij->A);
943: MatDestroy(&baij->B);
944: #if defined(PETSC_USE_CTABLE)
945: PetscTableDestroy(&baij->colmap);
946: #else
947: PetscFree(baij->colmap);
948: #endif
949: PetscFree(baij->garray);
950: VecDestroy(&baij->lvec);
951: VecScatterDestroy(&baij->Mvctx);
952: VecDestroy(&baij->slvec0);
953: VecDestroy(&baij->slvec0b);
954: VecDestroy(&baij->slvec1);
955: VecDestroy(&baij->slvec1a);
956: VecDestroy(&baij->slvec1b);
957: VecScatterDestroy(&baij->sMvctx);
958: PetscFree2(baij->rowvalues,baij->rowindices);
959: PetscFree(baij->barray);
960: PetscFree(baij->hd);
961: VecDestroy(&baij->diag);
962: VecDestroy(&baij->bb1);
963: VecDestroy(&baij->xx1);
964: #if defined(PETSC_USE_REAL_MAT_SINGLE)
965: PetscFree(baij->setvaluescopy);
966: #endif
967: PetscFree(baij->in_loc);
968: PetscFree(baij->v_loc);
969: PetscFree(baij->rangebs);
970: PetscFree(mat->data);
972: PetscObjectChangeTypeName((PetscObject)mat,0);
973: PetscObjectComposeFunction((PetscObject)mat,"MatStoreValues_C",NULL);
974: PetscObjectComposeFunction((PetscObject)mat,"MatRetrieveValues_C",NULL);
975: PetscObjectComposeFunction((PetscObject)mat,"MatGetDiagonalBlock_C",NULL);
976: PetscObjectComposeFunction((PetscObject)mat,"MatMPISBAIJSetPreallocation_C",NULL);
977: PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpisbaij_mpisbstrm_C",NULL);
978: #if defined(PETSC_HAVE_ELEMENTAL)
979: PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpisbaij_elemental_C",NULL);
980: #endif
981: return(0);
982: }
986: PetscErrorCode MatMult_MPISBAIJ_Hermitian(Mat A,Vec xx,Vec yy)
987: {
988: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
989: PetscErrorCode ierr;
990: PetscInt nt,mbs=a->mbs,bs=A->rmap->bs;
991: PetscScalar *from;
992: const PetscScalar *x;
995: VecGetLocalSize(xx,&nt);
996: if (nt != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible partition of A and xx");
998: /* diagonal part */
999: (*a->A->ops->mult)(a->A,xx,a->slvec1a);
1000: VecSet(a->slvec1b,0.0);
1002: /* subdiagonal part */
1003: (*a->B->ops->multhermitiantranspose)(a->B,xx,a->slvec0b);
1005: /* copy x into the vec slvec0 */
1006: VecGetArray(a->slvec0,&from);
1007: VecGetArrayRead(xx,&x);
1009: PetscMemcpy(from,x,bs*mbs*sizeof(MatScalar));
1010: VecRestoreArray(a->slvec0,&from);
1011: VecRestoreArrayRead(xx,&x);
1013: VecScatterBegin(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
1014: VecScatterEnd(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
1015: /* supperdiagonal part */
1016: (*a->B->ops->multadd)(a->B,a->slvec1b,a->slvec1a,yy);
1017: return(0);
1018: }
1022: PetscErrorCode MatMult_MPISBAIJ(Mat A,Vec xx,Vec yy)
1023: {
1024: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1025: PetscErrorCode ierr;
1026: PetscInt nt,mbs=a->mbs,bs=A->rmap->bs;
1027: PetscScalar *from;
1028: const PetscScalar *x;
1031: VecGetLocalSize(xx,&nt);
1032: if (nt != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible partition of A and xx");
1034: /* diagonal part */
1035: (*a->A->ops->mult)(a->A,xx,a->slvec1a);
1036: VecSet(a->slvec1b,0.0);
1038: /* subdiagonal part */
1039: (*a->B->ops->multtranspose)(a->B,xx,a->slvec0b);
1041: /* copy x into the vec slvec0 */
1042: VecGetArray(a->slvec0,&from);
1043: VecGetArrayRead(xx,&x);
1045: PetscMemcpy(from,x,bs*mbs*sizeof(MatScalar));
1046: VecRestoreArray(a->slvec0,&from);
1047: VecRestoreArrayRead(xx,&x);
1049: VecScatterBegin(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
1050: VecScatterEnd(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
1051: /* supperdiagonal part */
1052: (*a->B->ops->multadd)(a->B,a->slvec1b,a->slvec1a,yy);
1053: return(0);
1054: }
1058: PetscErrorCode MatMult_MPISBAIJ_2comm(Mat A,Vec xx,Vec yy)
1059: {
1060: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1062: PetscInt nt;
1065: VecGetLocalSize(xx,&nt);
1066: if (nt != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible partition of A and xx");
1068: VecGetLocalSize(yy,&nt);
1069: if (nt != A->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible parition of A and yy");
1071: VecScatterBegin(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);
1072: /* do diagonal part */
1073: (*a->A->ops->mult)(a->A,xx,yy);
1074: /* do supperdiagonal part */
1075: VecScatterEnd(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);
1076: (*a->B->ops->multadd)(a->B,a->lvec,yy,yy);
1077: /* do subdiagonal part */
1078: (*a->B->ops->multtranspose)(a->B,xx,a->lvec);
1079: VecScatterBegin(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);
1080: VecScatterEnd(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);
1081: return(0);
1082: }
1086: PetscErrorCode MatMultAdd_MPISBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
1087: {
1088: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1089: PetscErrorCode ierr;
1090: PetscInt mbs=a->mbs,bs=A->rmap->bs;
1091: PetscScalar *from,zero=0.0;
1092: const PetscScalar *x;
1095: /*
1096: PetscSynchronizedPrintf(PetscObjectComm((PetscObject)A)," MatMultAdd is called ...\n");
1097: PetscSynchronizedFlush(PetscObjectComm((PetscObject)A),PETSC_STDOUT);
1098: */
1099: /* diagonal part */
1100: (*a->A->ops->multadd)(a->A,xx,yy,a->slvec1a);
1101: VecSet(a->slvec1b,zero);
1103: /* subdiagonal part */
1104: (*a->B->ops->multtranspose)(a->B,xx,a->slvec0b);
1106: /* copy x into the vec slvec0 */
1107: VecGetArray(a->slvec0,&from);
1108: VecGetArrayRead(xx,&x);
1109: PetscMemcpy(from,x,bs*mbs*sizeof(MatScalar));
1110: VecRestoreArray(a->slvec0,&from);
1112: VecScatterBegin(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
1113: VecRestoreArrayRead(xx,&x);
1114: VecScatterEnd(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
1116: /* supperdiagonal part */
1117: (*a->B->ops->multadd)(a->B,a->slvec1b,a->slvec1a,zz);
1118: return(0);
1119: }
1123: PetscErrorCode MatMultAdd_MPISBAIJ_2comm(Mat A,Vec xx,Vec yy,Vec zz)
1124: {
1125: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1129: VecScatterBegin(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);
1130: /* do diagonal part */
1131: (*a->A->ops->multadd)(a->A,xx,yy,zz);
1132: /* do supperdiagonal part */
1133: VecScatterEnd(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);
1134: (*a->B->ops->multadd)(a->B,a->lvec,zz,zz);
1136: /* do subdiagonal part */
1137: (*a->B->ops->multtranspose)(a->B,xx,a->lvec);
1138: VecScatterBegin(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);
1139: VecScatterEnd(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);
1140: return(0);
1141: }
1143: /*
1144: This only works correctly for square matrices where the subblock A->A is the
1145: diagonal block
1146: */
1149: PetscErrorCode MatGetDiagonal_MPISBAIJ(Mat A,Vec v)
1150: {
1151: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1155: /* if (a->rmap->N != a->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Supports only square matrix where A->A is diag block"); */
1156: MatGetDiagonal(a->A,v);
1157: return(0);
1158: }
1162: PetscErrorCode MatScale_MPISBAIJ(Mat A,PetscScalar aa)
1163: {
1164: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1168: MatScale(a->A,aa);
1169: MatScale(a->B,aa);
1170: return(0);
1171: }
1175: PetscErrorCode MatGetRow_MPISBAIJ(Mat matin,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
1176: {
1177: Mat_MPISBAIJ *mat = (Mat_MPISBAIJ*)matin->data;
1178: PetscScalar *vworkA,*vworkB,**pvA,**pvB,*v_p;
1180: PetscInt bs = matin->rmap->bs,bs2 = mat->bs2,i,*cworkA,*cworkB,**pcA,**pcB;
1181: PetscInt nztot,nzA,nzB,lrow,brstart = matin->rmap->rstart,brend = matin->rmap->rend;
1182: PetscInt *cmap,*idx_p,cstart = mat->rstartbs;
1185: if (mat->getrowactive) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Already active");
1186: mat->getrowactive = PETSC_TRUE;
1188: if (!mat->rowvalues && (idx || v)) {
1189: /*
1190: allocate enough space to hold information from the longest row.
1191: */
1192: Mat_SeqSBAIJ *Aa = (Mat_SeqSBAIJ*)mat->A->data;
1193: Mat_SeqBAIJ *Ba = (Mat_SeqBAIJ*)mat->B->data;
1194: PetscInt max = 1,mbs = mat->mbs,tmp;
1195: for (i=0; i<mbs; i++) {
1196: tmp = Aa->i[i+1] - Aa->i[i] + Ba->i[i+1] - Ba->i[i]; /* row length */
1197: if (max < tmp) max = tmp;
1198: }
1199: PetscMalloc2(max*bs2,&mat->rowvalues,max*bs2,&mat->rowindices);
1200: }
1202: if (row < brstart || row >= brend) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only local rows");
1203: lrow = row - brstart; /* local row index */
1205: pvA = &vworkA; pcA = &cworkA; pvB = &vworkB; pcB = &cworkB;
1206: if (!v) {pvA = 0; pvB = 0;}
1207: if (!idx) {pcA = 0; if (!v) pcB = 0;}
1208: (*mat->A->ops->getrow)(mat->A,lrow,&nzA,pcA,pvA);
1209: (*mat->B->ops->getrow)(mat->B,lrow,&nzB,pcB,pvB);
1210: nztot = nzA + nzB;
1212: cmap = mat->garray;
1213: if (v || idx) {
1214: if (nztot) {
1215: /* Sort by increasing column numbers, assuming A and B already sorted */
1216: PetscInt imark = -1;
1217: if (v) {
1218: *v = v_p = mat->rowvalues;
1219: for (i=0; i<nzB; i++) {
1220: if (cmap[cworkB[i]/bs] < cstart) v_p[i] = vworkB[i];
1221: else break;
1222: }
1223: imark = i;
1224: for (i=0; i<nzA; i++) v_p[imark+i] = vworkA[i];
1225: for (i=imark; i<nzB; i++) v_p[nzA+i] = vworkB[i];
1226: }
1227: if (idx) {
1228: *idx = idx_p = mat->rowindices;
1229: if (imark > -1) {
1230: for (i=0; i<imark; i++) {
1231: idx_p[i] = cmap[cworkB[i]/bs]*bs + cworkB[i]%bs;
1232: }
1233: } else {
1234: for (i=0; i<nzB; i++) {
1235: if (cmap[cworkB[i]/bs] < cstart) idx_p[i] = cmap[cworkB[i]/bs]*bs + cworkB[i]%bs;
1236: else break;
1237: }
1238: imark = i;
1239: }
1240: for (i=0; i<nzA; i++) idx_p[imark+i] = cstart*bs + cworkA[i];
1241: for (i=imark; i<nzB; i++) idx_p[nzA+i] = cmap[cworkB[i]/bs]*bs + cworkB[i]%bs ;
1242: }
1243: } else {
1244: if (idx) *idx = 0;
1245: if (v) *v = 0;
1246: }
1247: }
1248: *nz = nztot;
1249: (*mat->A->ops->restorerow)(mat->A,lrow,&nzA,pcA,pvA);
1250: (*mat->B->ops->restorerow)(mat->B,lrow,&nzB,pcB,pvB);
1251: return(0);
1252: }
1256: PetscErrorCode MatRestoreRow_MPISBAIJ(Mat mat,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
1257: {
1258: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
1261: if (!baij->getrowactive) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"MatGetRow() must be called first");
1262: baij->getrowactive = PETSC_FALSE;
1263: return(0);
1264: }
1268: PetscErrorCode MatGetRowUpperTriangular_MPISBAIJ(Mat A)
1269: {
1270: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1271: Mat_SeqSBAIJ *aA = (Mat_SeqSBAIJ*)a->A->data;
1274: aA->getrow_utriangular = PETSC_TRUE;
1275: return(0);
1276: }
1279: PetscErrorCode MatRestoreRowUpperTriangular_MPISBAIJ(Mat A)
1280: {
1281: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1282: Mat_SeqSBAIJ *aA = (Mat_SeqSBAIJ*)a->A->data;
1285: aA->getrow_utriangular = PETSC_FALSE;
1286: return(0);
1287: }
1291: PetscErrorCode MatRealPart_MPISBAIJ(Mat A)
1292: {
1293: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1297: MatRealPart(a->A);
1298: MatRealPart(a->B);
1299: return(0);
1300: }
1304: PetscErrorCode MatImaginaryPart_MPISBAIJ(Mat A)
1305: {
1306: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1310: MatImaginaryPart(a->A);
1311: MatImaginaryPart(a->B);
1312: return(0);
1313: }
1315: /* Check if isrow is a subset of iscol_local, called by MatGetSubMatrix_MPISBAIJ()
1316: Input: isrow - distributed(parallel),
1317: iscol_local - locally owned (seq)
1318: */
1321: PetscErrorCode ISEqual_private(IS isrow,IS iscol_local,PetscBool *flg)
1322: {
1324: PetscInt sz1,sz2,*a1,*a2,i,j,k,nmatch;
1325: const PetscInt *ptr1,*ptr2;
1328: ISGetLocalSize(isrow,&sz1);
1329: ISGetLocalSize(iscol_local,&sz2);
1330: if (sz1 > sz2) {
1331: *flg = PETSC_FALSE;
1332: return(0);
1333: }
1335: ISGetIndices(isrow,&ptr1);
1336: ISGetIndices(iscol_local,&ptr2);
1338: PetscMalloc1(sz1,&a1);
1339: PetscMalloc1(sz2,&a2);
1340: PetscMemcpy(a1,ptr1,sz1*sizeof(PetscInt));
1341: PetscMemcpy(a2,ptr2,sz2*sizeof(PetscInt));
1342: PetscSortInt(sz1,a1);
1343: PetscSortInt(sz2,a2);
1345: nmatch=0;
1346: k = 0;
1347: for (i=0; i<sz1; i++){
1348: for (j=k; j<sz2; j++){
1349: if (a1[i] == a2[j]) {
1350: k = j; nmatch++;
1351: break;
1352: }
1353: }
1354: }
1355: ISRestoreIndices(isrow,&ptr1);
1356: ISRestoreIndices(iscol_local,&ptr2);
1357: PetscFree(a1);
1358: PetscFree(a2);
1359: if (nmatch < sz1) {
1360: *flg = PETSC_FALSE;
1361: } else {
1362: *flg = PETSC_TRUE;
1363: }
1364: return(0);
1365: }
1369: PetscErrorCode MatGetSubMatrix_MPISBAIJ(Mat mat,IS isrow,IS iscol,MatReuse call,Mat *newmat)
1370: {
1372: IS iscol_local;
1373: PetscInt csize;
1374: PetscBool isequal;
1377: ISGetLocalSize(iscol,&csize);
1378: if (call == MAT_REUSE_MATRIX) {
1379: PetscObjectQuery((PetscObject)*newmat,"ISAllGather",(PetscObject*)&iscol_local);
1380: if (!iscol_local) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Submatrix passed in was not used before, cannot reuse");
1381: } else {
1382: ISAllGather(iscol,&iscol_local);
1383: ISEqual_private(isrow,iscol_local,&isequal);
1384: if (!isequal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"For symmetric format, iscol must equal isrow");
1385: }
1387: /* now call MatGetSubMatrix_MPIBAIJ() */
1388: MatGetSubMatrix_MPIBAIJ_Private(mat,isrow,iscol_local,csize,call,newmat);
1389: if (call == MAT_INITIAL_MATRIX) {
1390: PetscObjectCompose((PetscObject)*newmat,"ISAllGather",(PetscObject)iscol_local);
1391: ISDestroy(&iscol_local);
1392: }
1393: return(0);
1394: }
1398: PetscErrorCode MatZeroEntries_MPISBAIJ(Mat A)
1399: {
1400: Mat_MPISBAIJ *l = (Mat_MPISBAIJ*)A->data;
1404: MatZeroEntries(l->A);
1405: MatZeroEntries(l->B);
1406: return(0);
1407: }
1411: PetscErrorCode MatGetInfo_MPISBAIJ(Mat matin,MatInfoType flag,MatInfo *info)
1412: {
1413: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)matin->data;
1414: Mat A = a->A,B = a->B;
1416: PetscReal isend[5],irecv[5];
1419: info->block_size = (PetscReal)matin->rmap->bs;
1421: MatGetInfo(A,MAT_LOCAL,info);
1423: isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
1424: isend[3] = info->memory; isend[4] = info->mallocs;
1426: MatGetInfo(B,MAT_LOCAL,info);
1428: isend[0] += info->nz_used; isend[1] += info->nz_allocated; isend[2] += info->nz_unneeded;
1429: isend[3] += info->memory; isend[4] += info->mallocs;
1430: if (flag == MAT_LOCAL) {
1431: info->nz_used = isend[0];
1432: info->nz_allocated = isend[1];
1433: info->nz_unneeded = isend[2];
1434: info->memory = isend[3];
1435: info->mallocs = isend[4];
1436: } else if (flag == MAT_GLOBAL_MAX) {
1437: MPIU_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)matin));
1439: info->nz_used = irecv[0];
1440: info->nz_allocated = irecv[1];
1441: info->nz_unneeded = irecv[2];
1442: info->memory = irecv[3];
1443: info->mallocs = irecv[4];
1444: } else if (flag == MAT_GLOBAL_SUM) {
1445: MPIU_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)matin));
1447: info->nz_used = irecv[0];
1448: info->nz_allocated = irecv[1];
1449: info->nz_unneeded = irecv[2];
1450: info->memory = irecv[3];
1451: info->mallocs = irecv[4];
1452: } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown MatInfoType argument %d",(int)flag);
1453: info->fill_ratio_given = 0; /* no parallel LU/ILU/Cholesky */
1454: info->fill_ratio_needed = 0;
1455: info->factor_mallocs = 0;
1456: return(0);
1457: }
1461: PetscErrorCode MatSetOption_MPISBAIJ(Mat A,MatOption op,PetscBool flg)
1462: {
1463: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1464: Mat_SeqSBAIJ *aA = (Mat_SeqSBAIJ*)a->A->data;
1468: switch (op) {
1469: case MAT_NEW_NONZERO_LOCATIONS:
1470: case MAT_NEW_NONZERO_ALLOCATION_ERR:
1471: case MAT_UNUSED_NONZERO_LOCATION_ERR:
1472: case MAT_KEEP_NONZERO_PATTERN:
1473: case MAT_NEW_NONZERO_LOCATION_ERR:
1474: MatCheckPreallocated(A,1);
1475: MatSetOption(a->A,op,flg);
1476: MatSetOption(a->B,op,flg);
1477: break;
1478: case MAT_ROW_ORIENTED:
1479: MatCheckPreallocated(A,1);
1480: a->roworiented = flg;
1482: MatSetOption(a->A,op,flg);
1483: MatSetOption(a->B,op,flg);
1484: break;
1485: case MAT_NEW_DIAGONALS:
1486: PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
1487: break;
1488: case MAT_IGNORE_OFF_PROC_ENTRIES:
1489: a->donotstash = flg;
1490: break;
1491: case MAT_USE_HASH_TABLE:
1492: a->ht_flag = flg;
1493: break;
1494: case MAT_HERMITIAN:
1495: MatCheckPreallocated(A,1);
1496: if (!A->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatAssemblyEnd() first");
1497: MatSetOption(a->A,op,flg);
1499: A->ops->mult = MatMult_MPISBAIJ_Hermitian;
1500: break;
1501: case MAT_SPD:
1502: A->spd_set = PETSC_TRUE;
1503: A->spd = flg;
1504: if (flg) {
1505: A->symmetric = PETSC_TRUE;
1506: A->structurally_symmetric = PETSC_TRUE;
1507: A->symmetric_set = PETSC_TRUE;
1508: A->structurally_symmetric_set = PETSC_TRUE;
1509: }
1510: break;
1511: case MAT_SYMMETRIC:
1512: MatCheckPreallocated(A,1);
1513: MatSetOption(a->A,op,flg);
1514: break;
1515: case MAT_STRUCTURALLY_SYMMETRIC:
1516: MatCheckPreallocated(A,1);
1517: MatSetOption(a->A,op,flg);
1518: break;
1519: case MAT_SYMMETRY_ETERNAL:
1520: if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix must be symmetric");
1521: PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
1522: break;
1523: case MAT_IGNORE_LOWER_TRIANGULAR:
1524: aA->ignore_ltriangular = flg;
1525: break;
1526: case MAT_ERROR_LOWER_TRIANGULAR:
1527: aA->ignore_ltriangular = flg;
1528: break;
1529: case MAT_GETROW_UPPERTRIANGULAR:
1530: aA->getrow_utriangular = flg;
1531: break;
1532: default:
1533: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %d",op);
1534: }
1535: return(0);
1536: }
1540: PetscErrorCode MatTranspose_MPISBAIJ(Mat A,MatReuse reuse,Mat *B)
1541: {
1545: if (MAT_INITIAL_MATRIX || *B != A) {
1546: MatDuplicate(A,MAT_COPY_VALUES,B);
1547: }
1548: return(0);
1549: }
1553: PetscErrorCode MatDiagonalScale_MPISBAIJ(Mat mat,Vec ll,Vec rr)
1554: {
1555: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
1556: Mat a = baij->A, b=baij->B;
1558: PetscInt nv,m,n;
1559: PetscBool flg;
1562: if (ll != rr) {
1563: VecEqual(ll,rr,&flg);
1564: if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"For symmetric format, left and right scaling vectors must be same\n");
1565: }
1566: if (!ll) return(0);
1568: MatGetLocalSize(mat,&m,&n);
1569: if (m != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"For symmetric format, local size %d %d must be same",m,n);
1571: VecGetLocalSize(rr,&nv);
1572: if (nv!=n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left and right vector non-conforming local size");
1574: VecScatterBegin(baij->Mvctx,rr,baij->lvec,INSERT_VALUES,SCATTER_FORWARD);
1576: /* left diagonalscale the off-diagonal part */
1577: (*b->ops->diagonalscale)(b,ll,NULL);
1579: /* scale the diagonal part */
1580: (*a->ops->diagonalscale)(a,ll,rr);
1582: /* right diagonalscale the off-diagonal part */
1583: VecScatterEnd(baij->Mvctx,rr,baij->lvec,INSERT_VALUES,SCATTER_FORWARD);
1584: (*b->ops->diagonalscale)(b,NULL,baij->lvec);
1585: return(0);
1586: }
1590: PetscErrorCode MatSetUnfactored_MPISBAIJ(Mat A)
1591: {
1592: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1596: MatSetUnfactored(a->A);
1597: return(0);
1598: }
1600: static PetscErrorCode MatDuplicate_MPISBAIJ(Mat,MatDuplicateOption,Mat*);
1604: PetscErrorCode MatEqual_MPISBAIJ(Mat A,Mat B,PetscBool *flag)
1605: {
1606: Mat_MPISBAIJ *matB = (Mat_MPISBAIJ*)B->data,*matA = (Mat_MPISBAIJ*)A->data;
1607: Mat a,b,c,d;
1608: PetscBool flg;
1612: a = matA->A; b = matA->B;
1613: c = matB->A; d = matB->B;
1615: MatEqual(a,c,&flg);
1616: if (flg) {
1617: MatEqual(b,d,&flg);
1618: }
1619: MPIU_Allreduce(&flg,flag,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));
1620: return(0);
1621: }
1625: PetscErrorCode MatCopy_MPISBAIJ(Mat A,Mat B,MatStructure str)
1626: {
1628: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1629: Mat_MPISBAIJ *b = (Mat_MPISBAIJ*)B->data;
1632: /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
1633: if ((str != SAME_NONZERO_PATTERN) || (A->ops->copy != B->ops->copy)) {
1634: MatGetRowUpperTriangular(A);
1635: MatCopy_Basic(A,B,str);
1636: MatRestoreRowUpperTriangular(A);
1637: } else {
1638: MatCopy(a->A,b->A,str);
1639: MatCopy(a->B,b->B,str);
1640: }
1641: return(0);
1642: }
1646: PetscErrorCode MatSetUp_MPISBAIJ(Mat A)
1647: {
1651: MatMPISBAIJSetPreallocation(A,A->rmap->bs,PETSC_DEFAULT,0,PETSC_DEFAULT,0);
1652: return(0);
1653: }
1657: PetscErrorCode MatAXPY_MPISBAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str)
1658: {
1660: Mat_MPISBAIJ *xx=(Mat_MPISBAIJ*)X->data,*yy=(Mat_MPISBAIJ*)Y->data;
1661: PetscBLASInt bnz,one=1;
1662: Mat_SeqSBAIJ *xa,*ya;
1663: Mat_SeqBAIJ *xb,*yb;
1666: if (str == SAME_NONZERO_PATTERN) {
1667: PetscScalar alpha = a;
1668: xa = (Mat_SeqSBAIJ*)xx->A->data;
1669: ya = (Mat_SeqSBAIJ*)yy->A->data;
1670: PetscBLASIntCast(xa->nz,&bnz);
1671: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&bnz,&alpha,xa->a,&one,ya->a,&one));
1672: xb = (Mat_SeqBAIJ*)xx->B->data;
1673: yb = (Mat_SeqBAIJ*)yy->B->data;
1674: PetscBLASIntCast(xb->nz,&bnz);
1675: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&bnz,&alpha,xb->a,&one,yb->a,&one));
1676: PetscObjectStateIncrease((PetscObject)Y);
1677: } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
1678: MatSetOption(X,MAT_GETROW_UPPERTRIANGULAR,PETSC_TRUE);
1679: MatAXPY_Basic(Y,a,X,str);
1680: MatSetOption(X,MAT_GETROW_UPPERTRIANGULAR,PETSC_FALSE);
1681: } else {
1682: Mat B;
1683: PetscInt *nnz_d,*nnz_o,bs=Y->rmap->bs;
1684: if (bs != X->rmap->bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrices must have same block size");
1685: MatGetRowUpperTriangular(X);
1686: MatGetRowUpperTriangular(Y);
1687: PetscMalloc1(yy->A->rmap->N,&nnz_d);
1688: PetscMalloc1(yy->B->rmap->N,&nnz_o);
1689: MatCreate(PetscObjectComm((PetscObject)Y),&B);
1690: PetscObjectSetName((PetscObject)B,((PetscObject)Y)->name);
1691: MatSetSizes(B,Y->rmap->n,Y->cmap->n,Y->rmap->N,Y->cmap->N);
1692: MatSetBlockSizesFromMats(B,Y,Y);
1693: MatSetType(B,MATMPISBAIJ);
1694: MatAXPYGetPreallocation_SeqSBAIJ(yy->A,xx->A,nnz_d);
1695: MatAXPYGetPreallocation_MPIBAIJ(yy->B,yy->garray,xx->B,xx->garray,nnz_o);
1696: MatMPISBAIJSetPreallocation(B,bs,0,nnz_d,0,nnz_o);
1697: MatAXPY_BasicWithPreallocation(B,Y,a,X,str);
1698: MatHeaderReplace(Y,&B);
1699: PetscFree(nnz_d);
1700: PetscFree(nnz_o);
1701: MatRestoreRowUpperTriangular(X);
1702: MatRestoreRowUpperTriangular(Y);
1703: }
1704: return(0);
1705: }
1709: PetscErrorCode MatGetSubMatrices_MPISBAIJ(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1710: {
1712: PetscInt i;
1713: PetscBool flg;
1716: MatGetSubMatrices_MPIBAIJ(A,n,irow,icol,scall,B); /* B[] are sbaij matrices */
1717: for (i=0; i<n; i++) {
1718: ISEqual(irow[i],icol[i],&flg);
1719: if (!flg) {
1720: MatSeqSBAIJZeroOps_Private(*B[i]);
1721: }
1722: }
1723: return(0);
1724: }
1728: PetscErrorCode MatShift_MPISBAIJ(Mat Y,PetscScalar a)
1729: {
1731: Mat_MPISBAIJ *maij = (Mat_MPISBAIJ*)Y->data;
1732: Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ*)maij->A->data;
1735: if (!Y->preallocated) {
1736: MatMPISBAIJSetPreallocation(Y,Y->rmap->bs,1,NULL,0,NULL);
1737: } else if (!aij->nz) {
1738: PetscInt nonew = aij->nonew;
1739: MatSeqSBAIJSetPreallocation(maij->A,Y->rmap->bs,1,NULL);
1740: aij->nonew = nonew;
1741: }
1742: MatShift_Basic(Y,a);
1743: return(0);
1744: }
1748: PetscErrorCode MatMissingDiagonal_MPISBAIJ(Mat A,PetscBool *missing,PetscInt *d)
1749: {
1750: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1754: if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only works for square matrices");
1755: MatMissingDiagonal(a->A,missing,d);
1756: if (d) {
1757: PetscInt rstart;
1758: MatGetOwnershipRange(A,&rstart,NULL);
1759: *d += rstart/A->rmap->bs;
1761: }
1762: return(0);
1763: }
1766: /* -------------------------------------------------------------------*/
1767: static struct _MatOps MatOps_Values = {MatSetValues_MPISBAIJ,
1768: MatGetRow_MPISBAIJ,
1769: MatRestoreRow_MPISBAIJ,
1770: MatMult_MPISBAIJ,
1771: /* 4*/ MatMultAdd_MPISBAIJ,
1772: MatMult_MPISBAIJ, /* transpose versions are same as non-transpose */
1773: MatMultAdd_MPISBAIJ,
1774: 0,
1775: 0,
1776: 0,
1777: /* 10*/ 0,
1778: 0,
1779: 0,
1780: MatSOR_MPISBAIJ,
1781: MatTranspose_MPISBAIJ,
1782: /* 15*/ MatGetInfo_MPISBAIJ,
1783: MatEqual_MPISBAIJ,
1784: MatGetDiagonal_MPISBAIJ,
1785: MatDiagonalScale_MPISBAIJ,
1786: MatNorm_MPISBAIJ,
1787: /* 20*/ MatAssemblyBegin_MPISBAIJ,
1788: MatAssemblyEnd_MPISBAIJ,
1789: MatSetOption_MPISBAIJ,
1790: MatZeroEntries_MPISBAIJ,
1791: /* 24*/ 0,
1792: 0,
1793: 0,
1794: 0,
1795: 0,
1796: /* 29*/ MatSetUp_MPISBAIJ,
1797: 0,
1798: 0,
1799: 0,
1800: 0,
1801: /* 34*/ MatDuplicate_MPISBAIJ,
1802: 0,
1803: 0,
1804: 0,
1805: 0,
1806: /* 39*/ MatAXPY_MPISBAIJ,
1807: MatGetSubMatrices_MPISBAIJ,
1808: MatIncreaseOverlap_MPISBAIJ,
1809: MatGetValues_MPISBAIJ,
1810: MatCopy_MPISBAIJ,
1811: /* 44*/ 0,
1812: MatScale_MPISBAIJ,
1813: MatShift_MPISBAIJ,
1814: 0,
1815: 0,
1816: /* 49*/ 0,
1817: 0,
1818: 0,
1819: 0,
1820: 0,
1821: /* 54*/ 0,
1822: 0,
1823: MatSetUnfactored_MPISBAIJ,
1824: 0,
1825: MatSetValuesBlocked_MPISBAIJ,
1826: /* 59*/ MatGetSubMatrix_MPISBAIJ,
1827: 0,
1828: 0,
1829: 0,
1830: 0,
1831: /* 64*/ 0,
1832: 0,
1833: 0,
1834: 0,
1835: 0,
1836: /* 69*/ MatGetRowMaxAbs_MPISBAIJ,
1837: 0,
1838: 0,
1839: 0,
1840: 0,
1841: /* 74*/ 0,
1842: 0,
1843: 0,
1844: 0,
1845: 0,
1846: /* 79*/ 0,
1847: 0,
1848: 0,
1849: 0,
1850: MatLoad_MPISBAIJ,
1851: /* 84*/ 0,
1852: 0,
1853: 0,
1854: 0,
1855: 0,
1856: /* 89*/ 0,
1857: 0,
1858: 0,
1859: 0,
1860: 0,
1861: /* 94*/ 0,
1862: 0,
1863: 0,
1864: 0,
1865: 0,
1866: /* 99*/ 0,
1867: 0,
1868: 0,
1869: 0,
1870: 0,
1871: /*104*/ 0,
1872: MatRealPart_MPISBAIJ,
1873: MatImaginaryPart_MPISBAIJ,
1874: MatGetRowUpperTriangular_MPISBAIJ,
1875: MatRestoreRowUpperTriangular_MPISBAIJ,
1876: /*109*/ 0,
1877: 0,
1878: 0,
1879: 0,
1880: MatMissingDiagonal_MPISBAIJ,
1881: /*114*/ 0,
1882: 0,
1883: 0,
1884: 0,
1885: 0,
1886: /*119*/ 0,
1887: 0,
1888: 0,
1889: 0,
1890: 0,
1891: /*124*/ 0,
1892: 0,
1893: 0,
1894: 0,
1895: 0,
1896: /*129*/ 0,
1897: 0,
1898: 0,
1899: 0,
1900: 0,
1901: /*134*/ 0,
1902: 0,
1903: 0,
1904: 0,
1905: 0,
1906: /*139*/ 0,
1907: 0,
1908: 0,
1909: 0,
1910: 0,
1911: /*144*/MatCreateMPIMatConcatenateSeqMat_MPISBAIJ
1912: };
1916: PetscErrorCode MatGetDiagonalBlock_MPISBAIJ(Mat A,Mat *a)
1917: {
1919: *a = ((Mat_MPISBAIJ*)A->data)->A;
1920: return(0);
1921: }
1925: PetscErrorCode MatMPISBAIJSetPreallocation_MPISBAIJ(Mat B,PetscInt bs,PetscInt d_nz,const PetscInt *d_nnz,PetscInt o_nz,const PetscInt *o_nnz)
1926: {
1927: Mat_MPISBAIJ *b;
1929: PetscInt i,mbs,Mbs;
1932: MatSetBlockSize(B,PetscAbs(bs));
1933: PetscLayoutSetUp(B->rmap);
1934: PetscLayoutSetUp(B->cmap);
1935: PetscLayoutGetBlockSize(B->rmap,&bs);
1937: b = (Mat_MPISBAIJ*)B->data;
1938: mbs = B->rmap->n/bs;
1939: Mbs = B->rmap->N/bs;
1940: 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);
1942: B->rmap->bs = bs;
1943: b->bs2 = bs*bs;
1944: b->mbs = mbs;
1945: b->Mbs = Mbs;
1946: b->nbs = B->cmap->n/bs;
1947: b->Nbs = B->cmap->N/bs;
1949: for (i=0; i<=b->size; i++) {
1950: b->rangebs[i] = B->rmap->range[i]/bs;
1951: }
1952: b->rstartbs = B->rmap->rstart/bs;
1953: b->rendbs = B->rmap->rend/bs;
1955: b->cstartbs = B->cmap->rstart/bs;
1956: b->cendbs = B->cmap->rend/bs;
1958: if (!B->preallocated) {
1959: MatCreate(PETSC_COMM_SELF,&b->A);
1960: MatSetSizes(b->A,B->rmap->n,B->cmap->n,B->rmap->n,B->cmap->n);
1961: MatSetType(b->A,MATSEQSBAIJ);
1962: PetscLogObjectParent((PetscObject)B,(PetscObject)b->A);
1963: MatCreate(PETSC_COMM_SELF,&b->B);
1964: MatSetSizes(b->B,B->rmap->n,B->cmap->N,B->rmap->n,B->cmap->N);
1965: MatSetType(b->B,MATSEQBAIJ);
1966: PetscLogObjectParent((PetscObject)B,(PetscObject)b->B);
1967: MatStashCreate_Private(PetscObjectComm((PetscObject)B),bs,&B->bstash);
1968: }
1970: MatSeqSBAIJSetPreallocation(b->A,bs,d_nz,d_nnz);
1971: MatSeqBAIJSetPreallocation(b->B,bs,o_nz,o_nnz);
1973: B->preallocated = PETSC_TRUE;
1974: return(0);
1975: }
1979: PetscErrorCode MatMPISBAIJSetPreallocationCSR_MPISBAIJ(Mat B,PetscInt bs,const PetscInt ii[],const PetscInt jj[],const PetscScalar V[])
1980: {
1981: PetscInt m,rstart,cstart,cend;
1982: PetscInt i,j,d,nz,nz_max=0,*d_nnz=0,*o_nnz=0;
1983: const PetscInt *JJ =0;
1984: PetscScalar *values=0;
1988: if (bs < 1) SETERRQ1(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_OUTOFRANGE,"Invalid block size specified, must be positive but it is %D",bs);
1989: PetscLayoutSetBlockSize(B->rmap,bs);
1990: PetscLayoutSetBlockSize(B->cmap,bs);
1991: PetscLayoutSetUp(B->rmap);
1992: PetscLayoutSetUp(B->cmap);
1993: PetscLayoutGetBlockSize(B->rmap,&bs);
1994: m = B->rmap->n/bs;
1995: rstart = B->rmap->rstart/bs;
1996: cstart = B->cmap->rstart/bs;
1997: cend = B->cmap->rend/bs;
1999: if (ii[0]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"ii[0] must be 0 but it is %D",ii[0]);
2000: PetscMalloc2(m,&d_nnz,m,&o_nnz);
2001: for (i=0; i<m; i++) {
2002: nz = ii[i+1] - ii[i];
2003: if (nz < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local row %D has a negative number of columns %D",i,nz);
2004: nz_max = PetscMax(nz_max,nz);
2005: JJ = jj + ii[i];
2006: for (j=0; j<nz; j++) {
2007: if (*JJ >= cstart) break;
2008: JJ++;
2009: }
2010: d = 0;
2011: for (; j<nz; j++) {
2012: if (*JJ++ >= cend) break;
2013: d++;
2014: }
2015: d_nnz[i] = d;
2016: o_nnz[i] = nz - d;
2017: }
2018: MatMPISBAIJSetPreallocation(B,bs,0,d_nnz,0,o_nnz);
2019: PetscFree2(d_nnz,o_nnz);
2021: values = (PetscScalar*)V;
2022: if (!values) {
2023: PetscMalloc1(bs*bs*nz_max,&values);
2024: PetscMemzero(values,bs*bs*nz_max*sizeof(PetscScalar));
2025: }
2026: for (i=0; i<m; i++) {
2027: PetscInt row = i + rstart;
2028: PetscInt ncols = ii[i+1] - ii[i];
2029: const PetscInt *icols = jj + ii[i];
2030: const PetscScalar *svals = values + (V ? (bs*bs*ii[i]) : 0);
2031: MatSetValuesBlocked_MPISBAIJ(B,1,&row,ncols,icols,svals,INSERT_VALUES);
2032: }
2034: if (!V) { PetscFree(values); }
2035: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
2036: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
2037: MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
2038: return(0);
2039: }
2041: /*MC
2042: MATMPISBAIJ - MATMPISBAIJ = "mpisbaij" - A matrix type to be used for distributed symmetric sparse block matrices,
2043: based on block compressed sparse row format. Only the upper triangular portion of the "diagonal" portion of
2044: the matrix is stored.
2046: For complex numbers by default this matrix is symmetric, NOT Hermitian symmetric. To make it Hermitian symmetric you
2047: can call MatSetOption(Mat, MAT_HERMITIAN);
2049: Options Database Keys:
2050: . -mat_type mpisbaij - sets the matrix type to "mpisbaij" during a call to MatSetFromOptions()
2052: Level: beginner
2054: .seealso: MatCreateMPISBAIJ
2055: M*/
2057: PETSC_INTERN PetscErrorCode MatConvert_MPISBAIJ_MPISBSTRM(Mat,MatType,MatReuse,Mat*);
2061: PETSC_EXTERN PetscErrorCode MatCreate_MPISBAIJ(Mat B)
2062: {
2063: Mat_MPISBAIJ *b;
2065: PetscBool flg = PETSC_FALSE;
2068: PetscNewLog(B,&b);
2069: B->data = (void*)b;
2070: PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
2072: B->ops->destroy = MatDestroy_MPISBAIJ;
2073: B->ops->view = MatView_MPISBAIJ;
2074: B->assembled = PETSC_FALSE;
2075: B->insertmode = NOT_SET_VALUES;
2077: MPI_Comm_rank(PetscObjectComm((PetscObject)B),&b->rank);
2078: MPI_Comm_size(PetscObjectComm((PetscObject)B),&b->size);
2080: /* build local table of row and column ownerships */
2081: PetscMalloc1(b->size+2,&b->rangebs);
2083: /* build cache for off array entries formed */
2084: MatStashCreate_Private(PetscObjectComm((PetscObject)B),1,&B->stash);
2086: b->donotstash = PETSC_FALSE;
2087: b->colmap = NULL;
2088: b->garray = NULL;
2089: b->roworiented = PETSC_TRUE;
2091: /* stuff used in block assembly */
2092: b->barray = 0;
2094: /* stuff used for matrix vector multiply */
2095: b->lvec = 0;
2096: b->Mvctx = 0;
2097: b->slvec0 = 0;
2098: b->slvec0b = 0;
2099: b->slvec1 = 0;
2100: b->slvec1a = 0;
2101: b->slvec1b = 0;
2102: b->sMvctx = 0;
2104: /* stuff for MatGetRow() */
2105: b->rowindices = 0;
2106: b->rowvalues = 0;
2107: b->getrowactive = PETSC_FALSE;
2109: /* hash table stuff */
2110: b->ht = 0;
2111: b->hd = 0;
2112: b->ht_size = 0;
2113: b->ht_flag = PETSC_FALSE;
2114: b->ht_fact = 0;
2115: b->ht_total_ct = 0;
2116: b->ht_insert_ct = 0;
2118: /* stuff for MatGetSubMatrices_MPIBAIJ_local() */
2119: b->ijonly = PETSC_FALSE;
2121: b->in_loc = 0;
2122: b->v_loc = 0;
2123: b->n_loc = 0;
2125: PetscObjectComposeFunction((PetscObject)B,"MatStoreValues_C",MatStoreValues_MPISBAIJ);
2126: PetscObjectComposeFunction((PetscObject)B,"MatRetrieveValues_C",MatRetrieveValues_MPISBAIJ);
2127: PetscObjectComposeFunction((PetscObject)B,"MatGetDiagonalBlock_C",MatGetDiagonalBlock_MPISBAIJ);
2128: PetscObjectComposeFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",MatMPISBAIJSetPreallocation_MPISBAIJ);
2129: PetscObjectComposeFunction((PetscObject)B,"MatMPISBAIJSetPreallocationCSR_C",MatMPISBAIJSetPreallocationCSR_MPISBAIJ);
2130: PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpisbaij_mpisbstrm_C",MatConvert_MPISBAIJ_MPISBSTRM);
2131: #if defined(PETSC_HAVE_ELEMENTAL)
2132: PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpisbaij_elemental_C",MatConvert_MPISBAIJ_Elemental);
2133: #endif
2135: B->symmetric = PETSC_TRUE;
2136: B->structurally_symmetric = PETSC_TRUE;
2137: B->symmetric_set = PETSC_TRUE;
2138: B->structurally_symmetric_set = PETSC_TRUE;
2140: PetscObjectChangeTypeName((PetscObject)B,MATMPISBAIJ);
2141: PetscOptionsBegin(PetscObjectComm((PetscObject)B),NULL,"Options for loading MPISBAIJ matrix 1","Mat");
2142: PetscOptionsBool("-mat_use_hash_table","Use hash table to save memory in constructing matrix","MatSetOption",flg,&flg,NULL);
2143: if (flg) {
2144: PetscReal fact = 1.39;
2145: MatSetOption(B,MAT_USE_HASH_TABLE,PETSC_TRUE);
2146: PetscOptionsReal("-mat_use_hash_table","Use hash table factor","MatMPIBAIJSetHashTableFactor",fact,&fact,NULL);
2147: if (fact <= 1.0) fact = 1.39;
2148: MatMPIBAIJSetHashTableFactor(B,fact);
2149: PetscInfo1(B,"Hash table Factor used %5.2f\n",fact);
2150: }
2151: PetscOptionsEnd();
2152: return(0);
2153: }
2155: /*MC
2156: MATSBAIJ - MATSBAIJ = "sbaij" - A matrix type to be used for symmetric block sparse matrices.
2158: This matrix type is identical to MATSEQSBAIJ when constructed with a single process communicator,
2159: and MATMPISBAIJ otherwise.
2161: Options Database Keys:
2162: . -mat_type sbaij - sets the matrix type to "sbaij" during a call to MatSetFromOptions()
2164: Level: beginner
2166: .seealso: MatCreateMPISBAIJ,MATSEQSBAIJ,MATMPISBAIJ
2167: M*/
2171: /*@C
2172: MatMPISBAIJSetPreallocation - For good matrix assembly performance
2173: the user should preallocate the matrix storage by setting the parameters
2174: d_nz (or d_nnz) and o_nz (or o_nnz). By setting these parameters accurately,
2175: performance can be increased by more than a factor of 50.
2177: Collective on Mat
2179: Input Parameters:
2180: + B - the matrix
2181: . bs - size of block, the blocks are ALWAYS square. One can use MatSetBlockSizes() to set a different row and column blocksize but the row
2182: blocksize always defines the size of the blocks. The column blocksize sets the blocksize of the vectors obtained with MatCreateVecs()
2183: . d_nz - number of block nonzeros per block row in diagonal portion of local
2184: submatrix (same for all local rows)
2185: . d_nnz - array containing the number of block nonzeros in the various block rows
2186: in the upper triangular and diagonal part of the in diagonal portion of the local
2187: (possibly different for each block row) or NULL. If you plan to factor the matrix you must leave room
2188: for the diagonal entry and set a value even if it is zero.
2189: . o_nz - number of block nonzeros per block row in the off-diagonal portion of local
2190: submatrix (same for all local rows).
2191: - o_nnz - array containing the number of nonzeros in the various block rows of the
2192: off-diagonal portion of the local submatrix that is right of the diagonal
2193: (possibly different for each block row) or NULL.
2196: Options Database Keys:
2197: . -mat_no_unroll - uses code that does not unroll the loops in the
2198: block calculations (much slower)
2199: . -mat_block_size - size of the blocks to use
2201: Notes:
2203: If PETSC_DECIDE or PETSC_DETERMINE is used for a particular argument on one processor
2204: than it must be used on all processors that share the object for that argument.
2206: If the *_nnz parameter is given then the *_nz parameter is ignored
2208: Storage Information:
2209: For a square global matrix we define each processor's diagonal portion
2210: to be its local rows and the corresponding columns (a square submatrix);
2211: each processor's off-diagonal portion encompasses the remainder of the
2212: local matrix (a rectangular submatrix).
2214: The user can specify preallocated storage for the diagonal part of
2215: the local submatrix with either d_nz or d_nnz (not both). Set
2216: d_nz=PETSC_DEFAULT and d_nnz=NULL for PETSc to control dynamic
2217: memory allocation. Likewise, specify preallocated storage for the
2218: off-diagonal part of the local submatrix with o_nz or o_nnz (not both).
2220: You can call MatGetInfo() to get information on how effective the preallocation was;
2221: for example the fields mallocs,nz_allocated,nz_used,nz_unneeded;
2222: You can also run with the option -info and look for messages with the string
2223: malloc in them to see if additional memory allocation was needed.
2225: Consider a processor that owns rows 3, 4 and 5 of a parallel matrix. In
2226: the figure below we depict these three local rows and all columns (0-11).
2228: .vb
2229: 0 1 2 3 4 5 6 7 8 9 10 11
2230: --------------------------
2231: row 3 |. . . d d d o o o o o o
2232: row 4 |. . . d d d o o o o o o
2233: row 5 |. . . d d d o o o o o o
2234: --------------------------
2235: .ve
2237: Thus, any entries in the d locations are stored in the d (diagonal)
2238: submatrix, and any entries in the o locations are stored in the
2239: o (off-diagonal) submatrix. Note that the d matrix is stored in
2240: MatSeqSBAIJ format and the o submatrix in MATSEQBAIJ format.
2242: Now d_nz should indicate the number of block nonzeros per row in the upper triangular
2243: plus the diagonal part of the d matrix,
2244: and o_nz should indicate the number of block nonzeros per row in the o matrix
2246: In general, for PDE problems in which most nonzeros are near the diagonal,
2247: one expects d_nz >> o_nz. For large problems you MUST preallocate memory
2248: or you will get TERRIBLE performance; see the users' manual chapter on
2249: matrices.
2251: Level: intermediate
2253: .keywords: matrix, block, aij, compressed row, sparse, parallel
2255: .seealso: MatCreate(), MatCreateSeqSBAIJ(), MatSetValues(), MatCreateBAIJ(), PetscSplitOwnership()
2256: @*/
2257: PetscErrorCode MatMPISBAIJSetPreallocation(Mat B,PetscInt bs,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
2258: {
2265: PetscTryMethod(B,"MatMPISBAIJSetPreallocation_C",(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,bs,d_nz,d_nnz,o_nz,o_nnz));
2266: return(0);
2267: }
2271: /*@C
2272: MatCreateSBAIJ - Creates a sparse parallel matrix in symmetric block AIJ format
2273: (block compressed row). For good matrix assembly performance
2274: the user should preallocate the matrix storage by setting the parameters
2275: d_nz (or d_nnz) and o_nz (or o_nnz). By setting these parameters accurately,
2276: performance can be increased by more than a factor of 50.
2278: Collective on MPI_Comm
2280: Input Parameters:
2281: + comm - MPI communicator
2282: . bs - size of block, the blocks are ALWAYS square. One can use MatSetBlockSizes() to set a different row and column blocksize but the row
2283: blocksize always defines the size of the blocks. The column blocksize sets the blocksize of the vectors obtained with MatCreateVecs()
2284: . m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
2285: This value should be the same as the local size used in creating the
2286: y vector for the matrix-vector product y = Ax.
2287: . n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
2288: This value should be the same as the local size used in creating the
2289: x vector for the matrix-vector product y = Ax.
2290: . M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
2291: . N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
2292: . d_nz - number of block nonzeros per block row in diagonal portion of local
2293: submatrix (same for all local rows)
2294: . d_nnz - array containing the number of block nonzeros in the various block rows
2295: in the upper triangular portion of the in diagonal portion of the local
2296: (possibly different for each block block row) or NULL.
2297: If you plan to factor the matrix you must leave room for the diagonal entry and
2298: set its value even if it is zero.
2299: . o_nz - number of block nonzeros per block row in the off-diagonal portion of local
2300: submatrix (same for all local rows).
2301: - o_nnz - array containing the number of nonzeros in the various block rows of the
2302: off-diagonal portion of the local submatrix (possibly different for
2303: each block row) or NULL.
2305: Output Parameter:
2306: . A - the matrix
2308: Options Database Keys:
2309: . -mat_no_unroll - uses code that does not unroll the loops in the
2310: block calculations (much slower)
2311: . -mat_block_size - size of the blocks to use
2312: . -mat_mpi - use the parallel matrix data structures even on one processor
2313: (defaults to using SeqBAIJ format on one processor)
2315: It is recommended that one use the MatCreate(), MatSetType() and/or MatSetFromOptions(),
2316: MatXXXXSetPreallocation() paradgm instead of this routine directly.
2317: [MatXXXXSetPreallocation() is, for example, MatSeqAIJSetPreallocation]
2319: Notes:
2320: The number of rows and columns must be divisible by blocksize.
2321: This matrix type does not support complex Hermitian operation.
2323: The user MUST specify either the local or global matrix dimensions
2324: (possibly both).
2326: If PETSC_DECIDE or PETSC_DETERMINE is used for a particular argument on one processor
2327: than it must be used on all processors that share the object for that argument.
2329: If the *_nnz parameter is given then the *_nz parameter is ignored
2331: Storage Information:
2332: For a square global matrix we define each processor's diagonal portion
2333: to be its local rows and the corresponding columns (a square submatrix);
2334: each processor's off-diagonal portion encompasses the remainder of the
2335: local matrix (a rectangular submatrix).
2337: The user can specify preallocated storage for the diagonal part of
2338: the local submatrix with either d_nz or d_nnz (not both). Set
2339: d_nz=PETSC_DEFAULT and d_nnz=NULL for PETSc to control dynamic
2340: memory allocation. Likewise, specify preallocated storage for the
2341: off-diagonal part of the local submatrix with o_nz or o_nnz (not both).
2343: Consider a processor that owns rows 3, 4 and 5 of a parallel matrix. In
2344: the figure below we depict these three local rows and all columns (0-11).
2346: .vb
2347: 0 1 2 3 4 5 6 7 8 9 10 11
2348: --------------------------
2349: row 3 |. . . d d d o o o o o o
2350: row 4 |. . . d d d o o o o o o
2351: row 5 |. . . d d d o o o o o o
2352: --------------------------
2353: .ve
2355: Thus, any entries in the d locations are stored in the d (diagonal)
2356: submatrix, and any entries in the o locations are stored in the
2357: o (off-diagonal) submatrix. Note that the d matrix is stored in
2358: MatSeqSBAIJ format and the o submatrix in MATSEQBAIJ format.
2360: Now d_nz should indicate the number of block nonzeros per row in the upper triangular
2361: plus the diagonal part of the d matrix,
2362: and o_nz should indicate the number of block nonzeros per row in the o matrix.
2363: In general, for PDE problems in which most nonzeros are near the diagonal,
2364: one expects d_nz >> o_nz. For large problems you MUST preallocate memory
2365: or you will get TERRIBLE performance; see the users' manual chapter on
2366: matrices.
2368: Level: intermediate
2370: .keywords: matrix, block, aij, compressed row, sparse, parallel
2372: .seealso: MatCreate(), MatCreateSeqSBAIJ(), MatSetValues(), MatCreateBAIJ()
2373: @*/
2375: 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)
2376: {
2378: PetscMPIInt size;
2381: MatCreate(comm,A);
2382: MatSetSizes(*A,m,n,M,N);
2383: MPI_Comm_size(comm,&size);
2384: if (size > 1) {
2385: MatSetType(*A,MATMPISBAIJ);
2386: MatMPISBAIJSetPreallocation(*A,bs,d_nz,d_nnz,o_nz,o_nnz);
2387: } else {
2388: MatSetType(*A,MATSEQSBAIJ);
2389: MatSeqSBAIJSetPreallocation(*A,bs,d_nz,d_nnz);
2390: }
2391: return(0);
2392: }
2397: static PetscErrorCode MatDuplicate_MPISBAIJ(Mat matin,MatDuplicateOption cpvalues,Mat *newmat)
2398: {
2399: Mat mat;
2400: Mat_MPISBAIJ *a,*oldmat = (Mat_MPISBAIJ*)matin->data;
2402: PetscInt len=0,nt,bs=matin->rmap->bs,mbs=oldmat->mbs;
2403: PetscScalar *array;
2406: *newmat = 0;
2408: MatCreate(PetscObjectComm((PetscObject)matin),&mat);
2409: MatSetSizes(mat,matin->rmap->n,matin->cmap->n,matin->rmap->N,matin->cmap->N);
2410: MatSetType(mat,((PetscObject)matin)->type_name);
2411: PetscMemcpy(mat->ops,matin->ops,sizeof(struct _MatOps));
2412: PetscLayoutReference(matin->rmap,&mat->rmap);
2413: PetscLayoutReference(matin->cmap,&mat->cmap);
2415: mat->factortype = matin->factortype;
2416: mat->preallocated = PETSC_TRUE;
2417: mat->assembled = PETSC_TRUE;
2418: mat->insertmode = NOT_SET_VALUES;
2420: a = (Mat_MPISBAIJ*)mat->data;
2421: a->bs2 = oldmat->bs2;
2422: a->mbs = oldmat->mbs;
2423: a->nbs = oldmat->nbs;
2424: a->Mbs = oldmat->Mbs;
2425: a->Nbs = oldmat->Nbs;
2428: a->size = oldmat->size;
2429: a->rank = oldmat->rank;
2430: a->donotstash = oldmat->donotstash;
2431: a->roworiented = oldmat->roworiented;
2432: a->rowindices = 0;
2433: a->rowvalues = 0;
2434: a->getrowactive = PETSC_FALSE;
2435: a->barray = 0;
2436: a->rstartbs = oldmat->rstartbs;
2437: a->rendbs = oldmat->rendbs;
2438: a->cstartbs = oldmat->cstartbs;
2439: a->cendbs = oldmat->cendbs;
2441: /* hash table stuff */
2442: a->ht = 0;
2443: a->hd = 0;
2444: a->ht_size = 0;
2445: a->ht_flag = oldmat->ht_flag;
2446: a->ht_fact = oldmat->ht_fact;
2447: a->ht_total_ct = 0;
2448: a->ht_insert_ct = 0;
2450: PetscMemcpy(a->rangebs,oldmat->rangebs,(a->size+2)*sizeof(PetscInt));
2451: if (oldmat->colmap) {
2452: #if defined(PETSC_USE_CTABLE)
2453: PetscTableCreateCopy(oldmat->colmap,&a->colmap);
2454: #else
2455: PetscMalloc1(a->Nbs,&a->colmap);
2456: PetscLogObjectMemory((PetscObject)mat,(a->Nbs)*sizeof(PetscInt));
2457: PetscMemcpy(a->colmap,oldmat->colmap,(a->Nbs)*sizeof(PetscInt));
2458: #endif
2459: } else a->colmap = 0;
2461: if (oldmat->garray && (len = ((Mat_SeqBAIJ*)(oldmat->B->data))->nbs)) {
2462: PetscMalloc1(len,&a->garray);
2463: PetscLogObjectMemory((PetscObject)mat,len*sizeof(PetscInt));
2464: PetscMemcpy(a->garray,oldmat->garray,len*sizeof(PetscInt));
2465: } else a->garray = 0;
2467: MatStashCreate_Private(PetscObjectComm((PetscObject)matin),matin->rmap->bs,&mat->bstash);
2468: VecDuplicate(oldmat->lvec,&a->lvec);
2469: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->lvec);
2470: VecScatterCopy(oldmat->Mvctx,&a->Mvctx);
2471: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->Mvctx);
2473: VecDuplicate(oldmat->slvec0,&a->slvec0);
2474: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec0);
2475: VecDuplicate(oldmat->slvec1,&a->slvec1);
2476: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec1);
2478: VecGetLocalSize(a->slvec1,&nt);
2479: VecGetArray(a->slvec1,&array);
2480: VecCreateSeqWithArray(PETSC_COMM_SELF,1,bs*mbs,array,&a->slvec1a);
2481: VecCreateSeqWithArray(PETSC_COMM_SELF,1,nt-bs*mbs,array+bs*mbs,&a->slvec1b);
2482: VecRestoreArray(a->slvec1,&array);
2483: VecGetArray(a->slvec0,&array);
2484: VecCreateSeqWithArray(PETSC_COMM_SELF,1,nt-bs*mbs,array+bs*mbs,&a->slvec0b);
2485: VecRestoreArray(a->slvec0,&array);
2486: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec0);
2487: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec1);
2488: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec0b);
2489: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec1a);
2490: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec1b);
2492: /* VecScatterCopy(oldmat->sMvctx,&a->sMvctx); - not written yet, replaced by the lazy trick: */
2493: PetscObjectReference((PetscObject)oldmat->sMvctx);
2494: a->sMvctx = oldmat->sMvctx;
2495: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->sMvctx);
2497: MatDuplicate(oldmat->A,cpvalues,&a->A);
2498: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->A);
2499: MatDuplicate(oldmat->B,cpvalues,&a->B);
2500: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->B);
2501: PetscFunctionListDuplicate(((PetscObject)matin)->qlist,&((PetscObject)mat)->qlist);
2502: *newmat = mat;
2503: return(0);
2504: }
2508: PetscErrorCode MatLoad_MPISBAIJ(Mat newmat,PetscViewer viewer)
2509: {
2511: PetscInt i,nz,j,rstart,rend;
2512: PetscScalar *vals,*buf;
2513: MPI_Comm comm;
2514: MPI_Status status;
2515: PetscMPIInt rank,size,tag = ((PetscObject)viewer)->tag,*sndcounts = 0,*browners,maxnz,*rowners,mmbs;
2516: PetscInt header[4],*rowlengths = 0,M,N,m,*cols,*locrowlens;
2517: PetscInt *procsnz = 0,jj,*mycols,*ibuf;
2518: PetscInt bs = newmat->rmap->bs,Mbs,mbs,extra_rows;
2519: PetscInt *dlens,*odlens,*mask,*masked1,*masked2,rowcount,odcount;
2520: PetscInt dcount,kmax,k,nzcount,tmp;
2521: int fd;
2524: /* force binary viewer to load .info file if it has not yet done so */
2525: PetscViewerSetUp(viewer);
2526: PetscObjectGetComm((PetscObject)viewer,&comm);
2527: PetscOptionsBegin(comm,NULL,"Options for loading MPISBAIJ matrix 2","Mat");
2528: PetscOptionsInt("-matload_block_size","Set the blocksize used to store the matrix","MatLoad",bs,&bs,NULL);
2529: PetscOptionsEnd();
2530: if (bs < 0) bs = 1;
2532: MPI_Comm_size(comm,&size);
2533: MPI_Comm_rank(comm,&rank);
2534: PetscViewerBinaryGetDescriptor(viewer,&fd);
2535: if (!rank) {
2536: PetscBinaryRead(fd,(char*)header,4,PETSC_INT);
2537: if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"not matrix object");
2538: if (header[3] < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as MPISBAIJ");
2539: }
2541: MPI_Bcast(header+1,3,MPIU_INT,0,comm);
2542: M = header[1];
2543: N = header[2];
2545: /* If global sizes are set, check if they are consistent with that given in the file */
2546: 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);
2547: 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);
2549: if (M != N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Can only do square matrices");
2551: /*
2552: This code adds extra rows to make sure the number of rows is
2553: divisible by the blocksize
2554: */
2555: Mbs = M/bs;
2556: extra_rows = bs - M + bs*(Mbs);
2557: if (extra_rows == bs) extra_rows = 0;
2558: else Mbs++;
2559: if (extra_rows &&!rank) {
2560: PetscInfo(viewer,"Padding loaded matrix to match blocksize\n");
2561: }
2563: /* determine ownership of all rows */
2564: if (newmat->rmap->n < 0) { /* PETSC_DECIDE */
2565: mbs = Mbs/size + ((Mbs % size) > rank);
2566: m = mbs*bs;
2567: } else { /* User Set */
2568: m = newmat->rmap->n;
2569: mbs = m/bs;
2570: }
2571: PetscMalloc2(size+1,&rowners,size+1,&browners);
2572: PetscMPIIntCast(mbs,&mmbs);
2573: MPI_Allgather(&mmbs,1,MPI_INT,rowners+1,1,MPI_INT,comm);
2574: rowners[0] = 0;
2575: for (i=2; i<=size; i++) rowners[i] += rowners[i-1];
2576: for (i=0; i<=size; i++) browners[i] = rowners[i]*bs;
2577: rstart = rowners[rank];
2578: rend = rowners[rank+1];
2580: /* distribute row lengths to all processors */
2581: PetscMalloc1((rend-rstart)*bs,&locrowlens);
2582: if (!rank) {
2583: PetscMalloc1(M+extra_rows,&rowlengths);
2584: PetscBinaryRead(fd,rowlengths,M,PETSC_INT);
2585: for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;
2586: PetscMalloc1(size,&sndcounts);
2587: for (i=0; i<size; i++) sndcounts[i] = browners[i+1] - browners[i];
2588: MPI_Scatterv(rowlengths,sndcounts,browners,MPIU_INT,locrowlens,(rend-rstart)*bs,MPIU_INT,0,comm);
2589: PetscFree(sndcounts);
2590: } else {
2591: MPI_Scatterv(0,0,0,MPIU_INT,locrowlens,(rend-rstart)*bs,MPIU_INT,0,comm);
2592: }
2594: if (!rank) { /* procs[0] */
2595: /* calculate the number of nonzeros on each processor */
2596: PetscMalloc1(size,&procsnz);
2597: PetscMemzero(procsnz,size*sizeof(PetscInt));
2598: for (i=0; i<size; i++) {
2599: for (j=rowners[i]*bs; j< rowners[i+1]*bs; j++) {
2600: procsnz[i] += rowlengths[j];
2601: }
2602: }
2603: PetscFree(rowlengths);
2605: /* determine max buffer needed and allocate it */
2606: maxnz = 0;
2607: for (i=0; i<size; i++) {
2608: maxnz = PetscMax(maxnz,procsnz[i]);
2609: }
2610: PetscMalloc1(maxnz,&cols);
2612: /* read in my part of the matrix column indices */
2613: nz = procsnz[0];
2614: PetscMalloc1(nz,&ibuf);
2615: mycols = ibuf;
2616: if (size == 1) nz -= extra_rows;
2617: PetscBinaryRead(fd,mycols,nz,PETSC_INT);
2618: if (size == 1) {
2619: for (i=0; i< extra_rows; i++) mycols[nz+i] = M+i;
2620: }
2622: /* read in every ones (except the last) and ship off */
2623: for (i=1; i<size-1; i++) {
2624: nz = procsnz[i];
2625: PetscBinaryRead(fd,cols,nz,PETSC_INT);
2626: MPI_Send(cols,nz,MPIU_INT,i,tag,comm);
2627: }
2628: /* read in the stuff for the last proc */
2629: if (size != 1) {
2630: nz = procsnz[size-1] - extra_rows; /* the extra rows are not on the disk */
2631: PetscBinaryRead(fd,cols,nz,PETSC_INT);
2632: for (i=0; i<extra_rows; i++) cols[nz+i] = M+i;
2633: MPI_Send(cols,nz+extra_rows,MPIU_INT,size-1,tag,comm);
2634: }
2635: PetscFree(cols);
2636: } else { /* procs[i], i>0 */
2637: /* determine buffer space needed for message */
2638: nz = 0;
2639: for (i=0; i<m; i++) nz += locrowlens[i];
2640: PetscMalloc1(nz,&ibuf);
2641: mycols = ibuf;
2642: /* receive message of column indices*/
2643: MPI_Recv(mycols,nz,MPIU_INT,0,tag,comm,&status);
2644: MPI_Get_count(&status,MPIU_INT,&maxnz);
2645: if (maxnz != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
2646: }
2648: /* loop over local rows, determining number of off diagonal entries */
2649: PetscMalloc2(rend-rstart,&dlens,rend-rstart,&odlens);
2650: PetscMalloc3(Mbs,&mask,Mbs,&masked1,Mbs,&masked2);
2651: PetscMemzero(mask,Mbs*sizeof(PetscInt));
2652: PetscMemzero(masked1,Mbs*sizeof(PetscInt));
2653: PetscMemzero(masked2,Mbs*sizeof(PetscInt));
2654: rowcount = 0;
2655: nzcount = 0;
2656: for (i=0; i<mbs; i++) {
2657: dcount = 0;
2658: odcount = 0;
2659: for (j=0; j<bs; j++) {
2660: kmax = locrowlens[rowcount];
2661: for (k=0; k<kmax; k++) {
2662: tmp = mycols[nzcount++]/bs; /* block col. index */
2663: if (!mask[tmp]) {
2664: mask[tmp] = 1;
2665: if (tmp < rstart || tmp >= rend) masked2[odcount++] = tmp; /* entry in off-diag portion */
2666: else masked1[dcount++] = tmp; /* entry in diag portion */
2667: }
2668: }
2669: rowcount++;
2670: }
2672: dlens[i] = dcount; /* d_nzz[i] */
2673: odlens[i] = odcount; /* o_nzz[i] */
2675: /* zero out the mask elements we set */
2676: for (j=0; j<dcount; j++) mask[masked1[j]] = 0;
2677: for (j=0; j<odcount; j++) mask[masked2[j]] = 0;
2678: }
2679: MatSetSizes(newmat,m,m,M+extra_rows,N+extra_rows);
2680: MatMPISBAIJSetPreallocation(newmat,bs,0,dlens,0,odlens);
2681: MatSetOption(newmat,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);
2683: if (!rank) {
2684: PetscMalloc1(maxnz,&buf);
2685: /* read in my part of the matrix numerical values */
2686: nz = procsnz[0];
2687: vals = buf;
2688: mycols = ibuf;
2689: if (size == 1) nz -= extra_rows;
2690: PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
2691: if (size == 1) {
2692: for (i=0; i< extra_rows; i++) vals[nz+i] = 1.0;
2693: }
2695: /* insert into matrix */
2696: jj = rstart*bs;
2697: for (i=0; i<m; i++) {
2698: MatSetValues(newmat,1,&jj,locrowlens[i],mycols,vals,INSERT_VALUES);
2699: mycols += locrowlens[i];
2700: vals += locrowlens[i];
2701: jj++;
2702: }
2704: /* read in other processors (except the last one) and ship out */
2705: for (i=1; i<size-1; i++) {
2706: nz = procsnz[i];
2707: vals = buf;
2708: PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
2709: MPI_Send(vals,nz,MPIU_SCALAR,i,((PetscObject)newmat)->tag,comm);
2710: }
2711: /* the last proc */
2712: if (size != 1) {
2713: nz = procsnz[i] - extra_rows;
2714: vals = buf;
2715: PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
2716: for (i=0; i<extra_rows; i++) vals[nz+i] = 1.0;
2717: MPI_Send(vals,nz+extra_rows,MPIU_SCALAR,size-1,((PetscObject)newmat)->tag,comm);
2718: }
2719: PetscFree(procsnz);
2721: } else {
2722: /* receive numeric values */
2723: PetscMalloc1(nz,&buf);
2725: /* receive message of values*/
2726: vals = buf;
2727: mycols = ibuf;
2728: MPI_Recv(vals,nz,MPIU_SCALAR,0,((PetscObject)newmat)->tag,comm,&status);
2729: MPI_Get_count(&status,MPIU_SCALAR,&maxnz);
2730: if (maxnz != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
2732: /* insert into matrix */
2733: jj = rstart*bs;
2734: for (i=0; i<m; i++) {
2735: MatSetValues_MPISBAIJ(newmat,1,&jj,locrowlens[i],mycols,vals,INSERT_VALUES);
2736: mycols += locrowlens[i];
2737: vals += locrowlens[i];
2738: jj++;
2739: }
2740: }
2742: PetscFree(locrowlens);
2743: PetscFree(buf);
2744: PetscFree(ibuf);
2745: PetscFree2(rowners,browners);
2746: PetscFree2(dlens,odlens);
2747: PetscFree3(mask,masked1,masked2);
2748: MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);
2749: MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);
2750: return(0);
2751: }
2755: /*XXXXX@
2756: MatMPISBAIJSetHashTableFactor - Sets the factor required to compute the size of the HashTable.
2758: Input Parameters:
2759: . mat - the matrix
2760: . fact - factor
2762: Not Collective on Mat, each process can have a different hash factor
2764: Level: advanced
2766: Notes:
2767: This can also be set by the command line option: -mat_use_hash_table fact
2769: .keywords: matrix, hashtable, factor, HT
2771: .seealso: MatSetOption()
2772: @XXXXX*/
2777: PetscErrorCode MatGetRowMaxAbs_MPISBAIJ(Mat A,Vec v,PetscInt idx[])
2778: {
2779: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
2780: Mat_SeqBAIJ *b = (Mat_SeqBAIJ*)(a->B)->data;
2781: PetscReal atmp;
2782: PetscReal *work,*svalues,*rvalues;
2784: PetscInt i,bs,mbs,*bi,*bj,brow,j,ncols,krow,kcol,col,row,Mbs,bcol;
2785: PetscMPIInt rank,size;
2786: PetscInt *rowners_bs,dest,count,source;
2787: PetscScalar *va;
2788: MatScalar *ba;
2789: MPI_Status stat;
2792: if (idx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Send email to petsc-maint@mcs.anl.gov");
2793: MatGetRowMaxAbs(a->A,v,NULL);
2794: VecGetArray(v,&va);
2796: MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);
2797: MPI_Comm_rank(PetscObjectComm((PetscObject)A),&rank);
2799: bs = A->rmap->bs;
2800: mbs = a->mbs;
2801: Mbs = a->Mbs;
2802: ba = b->a;
2803: bi = b->i;
2804: bj = b->j;
2806: /* find ownerships */
2807: rowners_bs = A->rmap->range;
2809: /* each proc creates an array to be distributed */
2810: PetscMalloc1(bs*Mbs,&work);
2811: PetscMemzero(work,bs*Mbs*sizeof(PetscReal));
2813: /* row_max for B */
2814: if (rank != size-1) {
2815: for (i=0; i<mbs; i++) {
2816: ncols = bi[1] - bi[0]; bi++;
2817: brow = bs*i;
2818: for (j=0; j<ncols; j++) {
2819: bcol = bs*(*bj);
2820: for (kcol=0; kcol<bs; kcol++) {
2821: col = bcol + kcol; /* local col index */
2822: col += rowners_bs[rank+1]; /* global col index */
2823: for (krow=0; krow<bs; krow++) {
2824: atmp = PetscAbsScalar(*ba); ba++;
2825: row = brow + krow; /* local row index */
2826: if (PetscRealPart(va[row]) < atmp) va[row] = atmp;
2827: if (work[col] < atmp) work[col] = atmp;
2828: }
2829: }
2830: bj++;
2831: }
2832: }
2834: /* send values to its owners */
2835: for (dest=rank+1; dest<size; dest++) {
2836: svalues = work + rowners_bs[dest];
2837: count = rowners_bs[dest+1]-rowners_bs[dest];
2838: MPI_Send(svalues,count,MPIU_REAL,dest,rank,PetscObjectComm((PetscObject)A));
2839: }
2840: }
2842: /* receive values */
2843: if (rank) {
2844: rvalues = work;
2845: count = rowners_bs[rank+1]-rowners_bs[rank];
2846: for (source=0; source<rank; source++) {
2847: MPI_Recv(rvalues,count,MPIU_REAL,MPI_ANY_SOURCE,MPI_ANY_TAG,PetscObjectComm((PetscObject)A),&stat);
2848: /* process values */
2849: for (i=0; i<count; i++) {
2850: if (PetscRealPart(va[i]) < rvalues[i]) va[i] = rvalues[i];
2851: }
2852: }
2853: }
2855: VecRestoreArray(v,&va);
2856: PetscFree(work);
2857: return(0);
2858: }
2862: PetscErrorCode MatSOR_MPISBAIJ(Mat matin,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
2863: {
2864: Mat_MPISBAIJ *mat = (Mat_MPISBAIJ*)matin->data;
2865: PetscErrorCode ierr;
2866: PetscInt mbs=mat->mbs,bs=matin->rmap->bs;
2867: PetscScalar *x,*ptr,*from;
2868: Vec bb1;
2869: const PetscScalar *b;
2872: 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);
2873: if (bs > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented");
2875: if (flag == SOR_APPLY_UPPER) {
2876: (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx);
2877: return(0);
2878: }
2880: if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP) {
2881: if (flag & SOR_ZERO_INITIAL_GUESS) {
2882: (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,lits,xx);
2883: its--;
2884: }
2886: VecDuplicate(bb,&bb1);
2887: while (its--) {
2889: /* lower triangular part: slvec0b = - B^T*xx */
2890: (*mat->B->ops->multtranspose)(mat->B,xx,mat->slvec0b);
2892: /* copy xx into slvec0a */
2893: VecGetArray(mat->slvec0,&ptr);
2894: VecGetArray(xx,&x);
2895: PetscMemcpy(ptr,x,bs*mbs*sizeof(MatScalar));
2896: VecRestoreArray(mat->slvec0,&ptr);
2898: VecScale(mat->slvec0,-1.0);
2900: /* copy bb into slvec1a */
2901: VecGetArray(mat->slvec1,&ptr);
2902: VecGetArrayRead(bb,&b);
2903: PetscMemcpy(ptr,b,bs*mbs*sizeof(MatScalar));
2904: VecRestoreArray(mat->slvec1,&ptr);
2906: /* set slvec1b = 0 */
2907: VecSet(mat->slvec1b,0.0);
2909: VecScatterBegin(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD);
2910: VecRestoreArray(xx,&x);
2911: VecRestoreArrayRead(bb,&b);
2912: VecScatterEnd(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD);
2914: /* upper triangular part: bb1 = bb1 - B*x */
2915: (*mat->B->ops->multadd)(mat->B,mat->slvec1b,mat->slvec1a,bb1);
2917: /* local diagonal sweep */
2918: (*mat->A->ops->sor)(mat->A,bb1,omega,SOR_SYMMETRIC_SWEEP,fshift,lits,lits,xx);
2919: }
2920: VecDestroy(&bb1);
2921: } else if ((flag & SOR_LOCAL_FORWARD_SWEEP) && (its == 1) && (flag & SOR_ZERO_INITIAL_GUESS)) {
2922: (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx);
2923: } else if ((flag & SOR_LOCAL_BACKWARD_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_EISENSTAT) {
2926: Vec xx1;
2927: PetscBool hasop;
2928: const PetscScalar *diag;
2929: PetscScalar *sl,scale = (omega - 2.0)/omega;
2930: PetscInt i,n;
2932: if (!mat->xx1) {
2933: VecDuplicate(bb,&mat->xx1);
2934: VecDuplicate(bb,&mat->bb1);
2935: }
2936: xx1 = mat->xx1;
2937: bb1 = mat->bb1;
2939: (*mat->A->ops->sor)(mat->A,bb,omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_BACKWARD_SWEEP),fshift,lits,1,xx);
2941: if (!mat->diag) {
2942: /* this is wrong for same matrix with new nonzero values */
2943: MatCreateVecs(matin,&mat->diag,NULL);
2944: MatGetDiagonal(matin,mat->diag);
2945: }
2946: MatHasOperation(matin,MATOP_MULT_DIAGONAL_BLOCK,&hasop);
2948: if (hasop) {
2949: MatMultDiagonalBlock(matin,xx,bb1);
2950: VecAYPX(mat->slvec1a,scale,bb);
2951: } else {
2952: /*
2953: These two lines are replaced by code that may be a bit faster for a good compiler
2954: VecPointwiseMult(mat->slvec1a,mat->diag,xx);
2955: VecAYPX(mat->slvec1a,scale,bb);
2956: */
2957: VecGetArray(mat->slvec1a,&sl);
2958: VecGetArrayRead(mat->diag,&diag);
2959: VecGetArrayRead(bb,&b);
2960: VecGetArray(xx,&x);
2961: VecGetLocalSize(xx,&n);
2962: if (omega == 1.0) {
2963: for (i=0; i<n; i++) sl[i] = b[i] - diag[i]*x[i];
2964: PetscLogFlops(2.0*n);
2965: } else {
2966: for (i=0; i<n; i++) sl[i] = b[i] + scale*diag[i]*x[i];
2967: PetscLogFlops(3.0*n);
2968: }
2969: VecRestoreArray(mat->slvec1a,&sl);
2970: VecRestoreArrayRead(mat->diag,&diag);
2971: VecRestoreArrayRead(bb,&b);
2972: VecRestoreArray(xx,&x);
2973: }
2975: /* multiply off-diagonal portion of matrix */
2976: VecSet(mat->slvec1b,0.0);
2977: (*mat->B->ops->multtranspose)(mat->B,xx,mat->slvec0b);
2978: VecGetArray(mat->slvec0,&from);
2979: VecGetArray(xx,&x);
2980: PetscMemcpy(from,x,bs*mbs*sizeof(MatScalar));
2981: VecRestoreArray(mat->slvec0,&from);
2982: VecRestoreArray(xx,&x);
2983: VecScatterBegin(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD);
2984: VecScatterEnd(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD);
2985: (*mat->B->ops->multadd)(mat->B,mat->slvec1b,mat->slvec1a,mat->slvec1a);
2987: /* local sweep */
2988: (*mat->A->ops->sor)(mat->A,mat->slvec1a,omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_FORWARD_SWEEP),fshift,lits,1,xx1);
2989: VecAXPY(xx,1.0,xx1);
2990: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatSORType is not supported for SBAIJ matrix format");
2991: return(0);
2992: }
2996: PetscErrorCode MatSOR_MPISBAIJ_2comm(Mat matin,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
2997: {
2998: Mat_MPISBAIJ *mat = (Mat_MPISBAIJ*)matin->data;
3000: Vec lvec1,bb1;
3003: 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);
3004: if (matin->rmap->bs > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented");
3006: if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP) {
3007: if (flag & SOR_ZERO_INITIAL_GUESS) {
3008: (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,lits,xx);
3009: its--;
3010: }
3012: VecDuplicate(mat->lvec,&lvec1);
3013: VecDuplicate(bb,&bb1);
3014: while (its--) {
3015: VecScatterBegin(mat->Mvctx,xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD);
3017: /* lower diagonal part: bb1 = bb - B^T*xx */
3018: (*mat->B->ops->multtranspose)(mat->B,xx,lvec1);
3019: VecScale(lvec1,-1.0);
3021: VecScatterEnd(mat->Mvctx,xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD);
3022: VecCopy(bb,bb1);
3023: VecScatterBegin(mat->Mvctx,lvec1,bb1,ADD_VALUES,SCATTER_REVERSE);
3025: /* upper diagonal part: bb1 = bb1 - B*x */
3026: VecScale(mat->lvec,-1.0);
3027: (*mat->B->ops->multadd)(mat->B,mat->lvec,bb1,bb1);
3029: VecScatterEnd(mat->Mvctx,lvec1,bb1,ADD_VALUES,SCATTER_REVERSE);
3031: /* diagonal sweep */
3032: (*mat->A->ops->sor)(mat->A,bb1,omega,SOR_SYMMETRIC_SWEEP,fshift,lits,lits,xx);
3033: }
3034: VecDestroy(&lvec1);
3035: VecDestroy(&bb1);
3036: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatSORType is not supported for SBAIJ matrix format");
3037: return(0);
3038: }
3042: /*@
3043: MatCreateMPISBAIJWithArrays - creates a MPI SBAIJ matrix using arrays that contain in standard
3044: CSR format the local rows.
3046: Collective on MPI_Comm
3048: Input Parameters:
3049: + comm - MPI communicator
3050: . bs - the block size, only a block size of 1 is supported
3051: . m - number of local rows (Cannot be PETSC_DECIDE)
3052: . n - This value should be the same as the local size used in creating the
3053: x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have
3054: calculated if N is given) For square matrices n is almost always m.
3055: . M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
3056: . N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
3057: . i - row indices
3058: . j - column indices
3059: - a - matrix values
3061: Output Parameter:
3062: . mat - the matrix
3064: Level: intermediate
3066: Notes:
3067: The i, j, and a arrays ARE copied by this routine into the internal format used by PETSc;
3068: thus you CANNOT change the matrix entries by changing the values of a[] after you have
3069: called this routine. Use MatCreateMPIAIJWithSplitArrays() to avoid needing to copy the arrays.
3071: The i and j indices are 0 based, and i indices are indices corresponding to the local j array.
3073: .keywords: matrix, aij, compressed row, sparse, parallel
3075: .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatMPIAIJSetPreallocation(), MatMPIAIJSetPreallocationCSR(),
3076: MPIAIJ, MatCreateAIJ(), MatCreateMPIAIJWithSplitArrays()
3077: @*/
3078: 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)
3079: {
3084: if (i[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"i (row indices) must start with 0");
3085: if (m < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"local number of rows (m) cannot be PETSC_DECIDE, or negative");
3086: MatCreate(comm,mat);
3087: MatSetSizes(*mat,m,n,M,N);
3088: MatSetType(*mat,MATMPISBAIJ);
3089: MatMPISBAIJSetPreallocationCSR(*mat,bs,i,j,a);
3090: return(0);
3091: }
3096: /*@C
3097: MatMPISBAIJSetPreallocationCSR - Allocates memory for a sparse parallel matrix in BAIJ format
3098: (the default parallel PETSc format).
3100: Collective on MPI_Comm
3102: Input Parameters:
3103: + B - the matrix
3104: . bs - the block size
3105: . i - the indices into j for the start of each local row (starts with zero)
3106: . j - the column indices for each local row (starts with zero) these must be sorted for each row
3107: - v - optional values in the matrix
3109: Level: developer
3111: .keywords: matrix, aij, compressed row, sparse, parallel
3113: .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatMPIBAIJSetPreallocation(), MatCreateAIJ(), MPIAIJ
3114: @*/
3115: PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat B,PetscInt bs,const PetscInt i[],const PetscInt j[], const PetscScalar v[])
3116: {
3120: PetscTryMethod(B,"MatMPISBAIJSetPreallocationCSR_C",(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]),(B,bs,i,j,v));
3121: return(0);
3122: }
3126: PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPISBAIJ(MPI_Comm comm,Mat inmat,PetscInt n,MatReuse scall,Mat *outmat)
3127: {
3129: PetscInt m,N,i,rstart,nnz,Ii,bs,cbs;
3130: PetscInt *indx;
3131: PetscScalar *values;
3134: MatGetSize(inmat,&m,&N);
3135: if (scall == MAT_INITIAL_MATRIX) { /* symbolic phase */
3136: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inmat->data;
3137: PetscInt *dnz,*onz,sum,bs,cbs,mbs,Nbs;
3138: PetscInt *bindx,rmax=a->rmax,j;
3139:
3140: MatGetBlockSizes(inmat,&bs,&cbs);
3141: mbs = m/bs; Nbs = N/cbs;
3142: if (n == PETSC_DECIDE) {
3143: PetscSplitOwnership(comm,&n,&Nbs);
3144: }
3145: /* Check sum(n) = Nbs */
3146: MPIU_Allreduce(&n,&sum,1,MPIU_INT,MPI_SUM,comm);
3147: if (sum != Nbs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Sum of local columns != global columns %d",Nbs);
3149: MPI_Scan(&mbs, &rstart,1,MPIU_INT,MPI_SUM,comm);
3150: rstart -= mbs;
3152: PetscMalloc1(rmax,&bindx);
3153: MatPreallocateInitialize(comm,mbs,n,dnz,onz);
3154: MatSetOption(inmat,MAT_GETROW_UPPERTRIANGULAR,PETSC_TRUE);
3155: for (i=0; i<mbs; i++) {
3156: MatGetRow_SeqSBAIJ(inmat,i*bs,&nnz,&indx,NULL); /* non-blocked nnz and indx */
3157: nnz = nnz/bs;
3158: for (j=0; j<nnz; j++) bindx[j] = indx[j*bs]/bs;
3159: MatPreallocateSet(i+rstart,nnz,bindx,dnz,onz);
3160: MatRestoreRow_SeqSBAIJ(inmat,i*bs,&nnz,&indx,NULL);
3161: }
3162: MatSetOption(inmat,MAT_GETROW_UPPERTRIANGULAR,PETSC_FALSE);
3163: PetscFree(bindx);
3165: MatCreate(comm,outmat);
3166: MatSetSizes(*outmat,m,n*bs,PETSC_DETERMINE,PETSC_DETERMINE);
3167: MatSetBlockSizes(*outmat,bs,cbs);
3168: MatSetType(*outmat,MATMPISBAIJ);
3169: MatMPISBAIJSetPreallocation(*outmat,bs,0,dnz,0,onz);
3170: MatPreallocateFinalize(dnz,onz);
3171: }
3172:
3173: /* numeric phase */
3174: MatGetBlockSizes(inmat,&bs,&cbs);
3175: MatGetOwnershipRange(*outmat,&rstart,NULL);
3177: MatSetOption(inmat,MAT_GETROW_UPPERTRIANGULAR,PETSC_TRUE);
3178: for (i=0; i<m; i++) {
3179: MatGetRow_SeqSBAIJ(inmat,i,&nnz,&indx,&values);
3180: Ii = i + rstart;
3181: MatSetValues(*outmat,1,&Ii,nnz,indx,values,INSERT_VALUES);
3182: MatRestoreRow_SeqSBAIJ(inmat,i,&nnz,&indx,&values);
3183: }
3184: MatSetOption(inmat,MAT_GETROW_UPPERTRIANGULAR,PETSC_FALSE);
3185: MatAssemblyBegin(*outmat,MAT_FINAL_ASSEMBLY);
3186: MatAssemblyEnd(*outmat,MAT_FINAL_ASSEMBLY);
3187: return(0);
3188: }