Actual source code: mpisbaij.c
petsc-3.7.7 2017-09-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: PetscLogFlops(nz*bs*bs);
637: }
638: /* Bmat */
639: v = bmat->a; jj = bmat->j;
640: for (brow=0; brow<mbs; brow++) {
641: grow = bs*(rstart + brow);
642: nz = bmat->i[brow+1] - bmat->i[brow];
643: for (bcol=0; bcol<nz; bcol++) {
644: gcol = bs*garray[*jj]; jj++;
645: for (col=0; col<bs; col++) {
646: for (row=0; row<bs; row++) {
647: vabs = PetscAbsScalar(*v); v++;
648: rsum[gcol+col] += vabs;
649: rsum[grow+row] += vabs;
650: }
651: }
652: }
653: PetscLogFlops(nz*bs*bs);
654: }
655: MPIU_Allreduce(rsum,rsum2,mat->cmap->N,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)mat));
656: *norm = 0.0;
657: for (col=0; col<mat->cmap->N; col++) {
658: if (rsum2[col] > *norm) *norm = rsum2[col];
659: }
660: PetscFree2(rsum,rsum2);
661: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for this norm yet");
662: }
663: return(0);
664: }
668: PetscErrorCode MatAssemblyBegin_MPISBAIJ(Mat mat,MatAssemblyType mode)
669: {
670: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
672: PetscInt nstash,reallocs;
675: if (baij->donotstash || mat->nooffprocentries) return(0);
677: MatStashScatterBegin_Private(mat,&mat->stash,mat->rmap->range);
678: MatStashScatterBegin_Private(mat,&mat->bstash,baij->rangebs);
679: MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);
680: PetscInfo2(mat,"Stash has %D entries,uses %D mallocs.\n",nstash,reallocs);
681: MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);
682: PetscInfo2(mat,"Block-Stash has %D entries, uses %D mallocs.\n",nstash,reallocs);
683: return(0);
684: }
688: PetscErrorCode MatAssemblyEnd_MPISBAIJ(Mat mat,MatAssemblyType mode)
689: {
690: Mat_MPISBAIJ *baij=(Mat_MPISBAIJ*)mat->data;
691: Mat_SeqSBAIJ *a =(Mat_SeqSBAIJ*)baij->A->data;
693: PetscInt i,j,rstart,ncols,flg,bs2=baij->bs2;
694: PetscInt *row,*col;
695: PetscBool other_disassembled;
696: PetscMPIInt n;
697: PetscBool r1,r2,r3;
698: MatScalar *val;
700: /* do not use 'b=(Mat_SeqBAIJ*)baij->B->data' as B can be reset in disassembly */
702: if (!baij->donotstash && !mat->nooffprocentries) {
703: while (1) {
704: MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);
705: if (!flg) break;
707: for (i=0; i<n;) {
708: /* Now identify the consecutive vals belonging to the same row */
709: for (j=i,rstart=row[j]; j<n; j++) {
710: if (row[j] != rstart) break;
711: }
712: if (j < n) ncols = j-i;
713: else ncols = n-i;
714: /* Now assemble all these values with a single function call */
715: MatSetValues_MPISBAIJ(mat,1,row+i,ncols,col+i,val+i,mat->insertmode);
716: i = j;
717: }
718: }
719: MatStashScatterEnd_Private(&mat->stash);
720: /* Now process the block-stash. Since the values are stashed column-oriented,
721: set the roworiented flag to column oriented, and after MatSetValues()
722: restore the original flags */
723: r1 = baij->roworiented;
724: r2 = a->roworiented;
725: r3 = ((Mat_SeqBAIJ*)baij->B->data)->roworiented;
727: baij->roworiented = PETSC_FALSE;
728: a->roworiented = PETSC_FALSE;
730: ((Mat_SeqBAIJ*)baij->B->data)->roworiented = PETSC_FALSE; /* b->roworinted */
731: while (1) {
732: MatStashScatterGetMesg_Private(&mat->bstash,&n,&row,&col,&val,&flg);
733: if (!flg) break;
735: for (i=0; i<n;) {
736: /* Now identify the consecutive vals belonging to the same row */
737: for (j=i,rstart=row[j]; j<n; j++) {
738: if (row[j] != rstart) break;
739: }
740: if (j < n) ncols = j-i;
741: else ncols = n-i;
742: MatSetValuesBlocked_MPISBAIJ(mat,1,row+i,ncols,col+i,val+i*bs2,mat->insertmode);
743: i = j;
744: }
745: }
746: MatStashScatterEnd_Private(&mat->bstash);
748: baij->roworiented = r1;
749: a->roworiented = r2;
751: ((Mat_SeqBAIJ*)baij->B->data)->roworiented = r3; /* b->roworinted */
752: }
754: MatAssemblyBegin(baij->A,mode);
755: MatAssemblyEnd(baij->A,mode);
757: /* determine if any processor has disassembled, if so we must
758: also disassemble ourselfs, in order that we may reassemble. */
759: /*
760: if nonzero structure of submatrix B cannot change then we know that
761: no processor disassembled thus we can skip this stuff
762: */
763: if (!((Mat_SeqBAIJ*)baij->B->data)->nonew) {
764: MPIU_Allreduce(&mat->was_assembled,&other_disassembled,1,MPIU_BOOL,MPI_PROD,PetscObjectComm((PetscObject)mat));
765: if (mat->was_assembled && !other_disassembled) {
766: MatDisAssemble_MPISBAIJ(mat);
767: }
768: }
770: if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
771: MatSetUpMultiply_MPISBAIJ(mat); /* setup Mvctx and sMvctx */
772: }
773: MatAssemblyBegin(baij->B,mode);
774: MatAssemblyEnd(baij->B,mode);
776: PetscFree2(baij->rowvalues,baij->rowindices);
778: baij->rowvalues = 0;
780: /* if no new nonzero locations are allowed in matrix then only set the matrix state the first time through */
781: if ((!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) || !((Mat_SeqBAIJ*)(baij->A->data))->nonew) {
782: PetscObjectState state = baij->A->nonzerostate + baij->B->nonzerostate;
783: MPIU_Allreduce(&state,&mat->nonzerostate,1,MPIU_INT64,MPI_SUM,PetscObjectComm((PetscObject)mat));
784: }
785: return(0);
786: }
788: extern PetscErrorCode MatView_SeqSBAIJ(Mat,PetscViewer);
789: extern PetscErrorCode MatSetValues_MPIBAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
790: #include <petscdraw.h>
793: static PetscErrorCode MatView_MPISBAIJ_ASCIIorDraworSocket(Mat mat,PetscViewer viewer)
794: {
795: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
796: PetscErrorCode ierr;
797: PetscInt bs = mat->rmap->bs;
798: PetscMPIInt rank = baij->rank;
799: PetscBool iascii,isdraw;
800: PetscViewer sviewer;
801: PetscViewerFormat format;
804: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
805: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
806: if (iascii) {
807: PetscViewerGetFormat(viewer,&format);
808: if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
809: MatInfo info;
810: MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&rank);
811: MatGetInfo(mat,MAT_LOCAL,&info);
812: PetscViewerASCIIPushSynchronized(viewer);
813: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Local rows %D nz %D nz alloced %D bs %D mem %D\n",rank,mat->rmap->n,(PetscInt)info.nz_used,(PetscInt)info.nz_allocated,mat->rmap->bs,(PetscInt)info.memory);
814: MatGetInfo(baij->A,MAT_LOCAL,&info);
815: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] on-diagonal part: nz %D \n",rank,(PetscInt)info.nz_used);
816: MatGetInfo(baij->B,MAT_LOCAL,&info);
817: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] off-diagonal part: nz %D \n",rank,(PetscInt)info.nz_used);
818: PetscViewerFlush(viewer);
819: PetscViewerASCIIPopSynchronized(viewer);
820: PetscViewerASCIIPrintf(viewer,"Information on VecScatter used in matrix-vector product: \n");
821: VecScatterView(baij->Mvctx,viewer);
822: return(0);
823: } else if (format == PETSC_VIEWER_ASCII_INFO) {
824: PetscViewerASCIIPrintf(viewer," block size is %D\n",bs);
825: return(0);
826: } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
827: return(0);
828: }
829: }
831: if (isdraw) {
832: PetscDraw draw;
833: PetscBool isnull;
834: PetscViewerDrawGetDraw(viewer,0,&draw);
835: PetscDrawIsNull(draw,&isnull);
836: if (isnull) return(0);
837: }
839: {
840: /* assemble the entire matrix onto first processor. */
841: Mat A;
842: Mat_SeqSBAIJ *Aloc;
843: Mat_SeqBAIJ *Bloc;
844: PetscInt M = mat->rmap->N,N = mat->cmap->N,*ai,*aj,col,i,j,k,*rvals,mbs = baij->mbs;
845: MatScalar *a;
846: const char *matname;
848: /* Should this be the same type as mat? */
849: MatCreate(PetscObjectComm((PetscObject)mat),&A);
850: if (!rank) {
851: MatSetSizes(A,M,N,M,N);
852: } else {
853: MatSetSizes(A,0,0,M,N);
854: }
855: MatSetType(A,MATMPISBAIJ);
856: MatMPISBAIJSetPreallocation(A,mat->rmap->bs,0,NULL,0,NULL);
857: MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);
858: PetscLogObjectParent((PetscObject)mat,(PetscObject)A);
860: /* copy over the A part */
861: Aloc = (Mat_SeqSBAIJ*)baij->A->data;
862: ai = Aloc->i; aj = Aloc->j; a = Aloc->a;
863: PetscMalloc1(bs,&rvals);
865: for (i=0; i<mbs; i++) {
866: rvals[0] = bs*(baij->rstartbs + i);
867: for (j=1; j<bs; j++) rvals[j] = rvals[j-1] + 1;
868: for (j=ai[i]; j<ai[i+1]; j++) {
869: col = (baij->cstartbs+aj[j])*bs;
870: for (k=0; k<bs; k++) {
871: MatSetValues_MPISBAIJ(A,bs,rvals,1,&col,a,INSERT_VALUES);
872: col++;
873: a += bs;
874: }
875: }
876: }
877: /* copy over the B part */
878: Bloc = (Mat_SeqBAIJ*)baij->B->data;
879: ai = Bloc->i; aj = Bloc->j; a = Bloc->a;
880: for (i=0; i<mbs; i++) {
882: rvals[0] = bs*(baij->rstartbs + i);
883: for (j=1; j<bs; j++) rvals[j] = rvals[j-1] + 1;
884: for (j=ai[i]; j<ai[i+1]; j++) {
885: col = baij->garray[aj[j]]*bs;
886: for (k=0; k<bs; k++) {
887: MatSetValues_MPIBAIJ(A,bs,rvals,1,&col,a,INSERT_VALUES);
888: col++;
889: a += bs;
890: }
891: }
892: }
893: PetscFree(rvals);
894: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
895: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
896: /*
897: Everyone has to call to draw the matrix since the graphics waits are
898: synchronized across all processors that share the PetscDraw object
899: */
900: PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);
901: PetscObjectGetName((PetscObject)mat,&matname);
902: if (!rank) {
903: PetscObjectSetName((PetscObject)((Mat_MPISBAIJ*)(A->data))->A,matname);
904: MatView_SeqSBAIJ(((Mat_MPISBAIJ*)(A->data))->A,sviewer);
905: }
906: PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);
907: PetscViewerFlush(viewer);
908: MatDestroy(&A);
909: }
910: return(0);
911: }
915: static PetscErrorCode MatView_MPISBAIJ_Binary(Mat mat,PetscViewer viewer)
916: {
917: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)mat->data;
918: Mat_SeqSBAIJ *A = (Mat_SeqSBAIJ*)a->A->data;
919: Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)a->B->data;
921: PetscInt i,*row_lens,*crow_lens,bs = mat->rmap->bs,j,k,bs2=a->bs2,header[4],nz,rlen;
922: PetscInt *range=0,nzmax,*column_indices,cnt,col,*garray = a->garray,cstart = mat->cmap->rstart/bs,len,pcnt,l,ll;
923: int fd;
924: PetscScalar *column_values;
925: FILE *file;
926: PetscMPIInt rank,size,tag = ((PetscObject)viewer)->tag;
927: PetscInt message_count,flowcontrolcount;
930: MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&rank);
931: MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);
932: nz = bs2*(A->nz + B->nz);
933: rlen = mat->rmap->n;
934: PetscViewerBinaryGetDescriptor(viewer,&fd);
935: if (!rank) {
936: header[0] = MAT_FILE_CLASSID;
937: header[1] = mat->rmap->N;
938: header[2] = mat->cmap->N;
940: MPI_Reduce(&nz,&header[3],1,MPIU_INT,MPI_SUM,0,PetscObjectComm((PetscObject)mat));
941: PetscBinaryWrite(fd,header,4,PETSC_INT,PETSC_TRUE);
942: /* get largest number of rows any processor has */
943: range = mat->rmap->range;
944: for (i=1; i<size; i++) {
945: rlen = PetscMax(rlen,range[i+1] - range[i]);
946: }
947: } else {
948: MPI_Reduce(&nz,0,1,MPIU_INT,MPI_SUM,0,PetscObjectComm((PetscObject)mat));
949: }
951: PetscMalloc1(rlen/bs,&crow_lens);
952: /* compute lengths of each row */
953: for (i=0; i<a->mbs; i++) {
954: crow_lens[i] = A->i[i+1] - A->i[i] + B->i[i+1] - B->i[i];
955: }
956: /* store the row lengths to the file */
957: PetscViewerFlowControlStart(viewer,&message_count,&flowcontrolcount);
958: if (!rank) {
959: MPI_Status status;
960: PetscMalloc1(rlen,&row_lens);
961: rlen = (range[1] - range[0])/bs;
962: for (i=0; i<rlen; i++) {
963: for (j=0; j<bs; j++) {
964: row_lens[i*bs+j] = bs*crow_lens[i];
965: }
966: }
967: PetscBinaryWrite(fd,row_lens,bs*rlen,PETSC_INT,PETSC_TRUE);
968: for (i=1; i<size; i++) {
969: rlen = (range[i+1] - range[i])/bs;
970: PetscViewerFlowControlStepMaster(viewer,i,&message_count,flowcontrolcount);
971: MPI_Recv(crow_lens,rlen,MPIU_INT,i,tag,PetscObjectComm((PetscObject)mat),&status);
972: for (k=0; k<rlen; k++) {
973: for (j=0; j<bs; j++) {
974: row_lens[k*bs+j] = bs*crow_lens[k];
975: }
976: }
977: PetscBinaryWrite(fd,row_lens,bs*rlen,PETSC_INT,PETSC_TRUE);
978: }
979: PetscViewerFlowControlEndMaster(viewer,&message_count);
980: PetscFree(row_lens);
981: } else {
982: PetscViewerFlowControlStepWorker(viewer,rank,&message_count);
983: MPI_Send(crow_lens,mat->rmap->n/bs,MPIU_INT,0,tag,PetscObjectComm((PetscObject)mat));
984: PetscViewerFlowControlEndWorker(viewer,&message_count);
985: }
986: PetscFree(crow_lens);
988: /* load up the local column indices. Include for all rows not just one for each block row since process 0 does not have the
989: information needed to make it for each row from a block row. This does require more communication but still not more than
990: the communication needed for the nonzero values */
991: nzmax = nz; /* space a largest processor needs */
992: MPI_Reduce(&nz,&nzmax,1,MPIU_INT,MPI_MAX,0,PetscObjectComm((PetscObject)mat));
993: PetscMalloc1(nzmax,&column_indices);
994: cnt = 0;
995: for (i=0; i<a->mbs; i++) {
996: pcnt = cnt;
997: for (j=B->i[i]; j<B->i[i+1]; j++) {
998: if ((col = garray[B->j[j]]) > cstart) break;
999: for (l=0; l<bs; l++) {
1000: column_indices[cnt++] = bs*col+l;
1001: }
1002: }
1003: for (k=A->i[i]; k<A->i[i+1]; k++) {
1004: for (l=0; l<bs; l++) {
1005: column_indices[cnt++] = bs*(A->j[k] + cstart)+l;
1006: }
1007: }
1008: for (; j<B->i[i+1]; j++) {
1009: for (l=0; l<bs; l++) {
1010: column_indices[cnt++] = bs*garray[B->j[j]]+l;
1011: }
1012: }
1013: len = cnt - pcnt;
1014: for (k=1; k<bs; k++) {
1015: PetscMemcpy(&column_indices[cnt],&column_indices[pcnt],len*sizeof(PetscInt));
1016: cnt += len;
1017: }
1018: }
1019: if (cnt != nz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB,"Internal PETSc error: cnt = %D nz = %D",cnt,nz);
1021: /* store the columns to the file */
1022: PetscViewerFlowControlStart(viewer,&message_count,&flowcontrolcount);
1023: if (!rank) {
1024: MPI_Status status;
1025: PetscBinaryWrite(fd,column_indices,nz,PETSC_INT,PETSC_TRUE);
1026: for (i=1; i<size; i++) {
1027: PetscViewerFlowControlStepMaster(viewer,i,&message_count,flowcontrolcount);
1028: MPI_Recv(&cnt,1,MPIU_INT,i,tag,PetscObjectComm((PetscObject)mat),&status);
1029: MPI_Recv(column_indices,cnt,MPIU_INT,i,tag,PetscObjectComm((PetscObject)mat),&status);
1030: PetscBinaryWrite(fd,column_indices,cnt,PETSC_INT,PETSC_TRUE);
1031: }
1032: PetscViewerFlowControlEndMaster(viewer,&message_count);
1033: } else {
1034: PetscViewerFlowControlStepWorker(viewer,rank,&message_count);
1035: MPI_Send(&cnt,1,MPIU_INT,0,tag,PetscObjectComm((PetscObject)mat));
1036: MPI_Send(column_indices,cnt,MPIU_INT,0,tag,PetscObjectComm((PetscObject)mat));
1037: PetscViewerFlowControlEndWorker(viewer,&message_count);
1038: }
1039: PetscFree(column_indices);
1041: /* load up the numerical values */
1042: PetscMalloc1(nzmax,&column_values);
1043: cnt = 0;
1044: for (i=0; i<a->mbs; i++) {
1045: rlen = bs*(B->i[i+1] - B->i[i] + A->i[i+1] - A->i[i]);
1046: for (j=B->i[i]; j<B->i[i+1]; j++) {
1047: if (garray[B->j[j]] > cstart) break;
1048: for (l=0; l<bs; l++) {
1049: for (ll=0; ll<bs; ll++) {
1050: column_values[cnt + l*rlen + ll] = B->a[bs2*j+l+bs*ll];
1051: }
1052: }
1053: cnt += bs;
1054: }
1055: for (k=A->i[i]; k<A->i[i+1]; k++) {
1056: for (l=0; l<bs; l++) {
1057: for (ll=0; ll<bs; ll++) {
1058: column_values[cnt + l*rlen + ll] = A->a[bs2*k+l+bs*ll];
1059: }
1060: }
1061: cnt += bs;
1062: }
1063: for (; j<B->i[i+1]; j++) {
1064: for (l=0; l<bs; l++) {
1065: for (ll=0; ll<bs; ll++) {
1066: column_values[cnt + l*rlen + ll] = B->a[bs2*j+l+bs*ll];
1067: }
1068: }
1069: cnt += bs;
1070: }
1071: cnt += (bs-1)*rlen;
1072: }
1073: if (cnt != nz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Internal PETSc error: cnt = %D nz = %D",cnt,nz);
1075: /* store the column values to the file */
1076: PetscViewerFlowControlStart(viewer,&message_count,&flowcontrolcount);
1077: if (!rank) {
1078: MPI_Status status;
1079: PetscBinaryWrite(fd,column_values,nz,PETSC_SCALAR,PETSC_TRUE);
1080: for (i=1; i<size; i++) {
1081: PetscViewerFlowControlStepMaster(viewer,i,&message_count,flowcontrolcount);
1082: MPI_Recv(&cnt,1,MPIU_INT,i,tag,PetscObjectComm((PetscObject)mat),&status);
1083: MPI_Recv(column_values,cnt,MPIU_SCALAR,i,tag,PetscObjectComm((PetscObject)mat),&status);
1084: PetscBinaryWrite(fd,column_values,cnt,PETSC_SCALAR,PETSC_TRUE);
1085: }
1086: PetscViewerFlowControlEndMaster(viewer,&message_count);
1087: } else {
1088: PetscViewerFlowControlStepWorker(viewer,rank,&message_count);
1089: MPI_Send(&nz,1,MPIU_INT,0,tag,PetscObjectComm((PetscObject)mat));
1090: MPI_Send(column_values,nz,MPIU_SCALAR,0,tag,PetscObjectComm((PetscObject)mat));
1091: PetscViewerFlowControlEndWorker(viewer,&message_count);
1092: }
1093: PetscFree(column_values);
1095: PetscViewerBinaryGetInfoPointer(viewer,&file);
1096: if (file) {
1097: fprintf(file,"-matload_block_size %d\n",(int)mat->rmap->bs);
1098: }
1099: return(0);
1100: }
1104: PetscErrorCode MatView_MPISBAIJ(Mat mat,PetscViewer viewer)
1105: {
1107: PetscBool iascii,isdraw,issocket,isbinary;
1110: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1111: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
1112: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSOCKET,&issocket);
1113: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
1114: if (iascii || isdraw || issocket) {
1115: MatView_MPISBAIJ_ASCIIorDraworSocket(mat,viewer);
1116: } else if (isbinary) {
1117: MatView_MPISBAIJ_Binary(mat,viewer);
1118: }
1119: return(0);
1120: }
1124: PetscErrorCode MatDestroy_MPISBAIJ(Mat mat)
1125: {
1126: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
1130: #if defined(PETSC_USE_LOG)
1131: PetscLogObjectState((PetscObject)mat,"Rows=%D,Cols=%D",mat->rmap->N,mat->cmap->N);
1132: #endif
1133: MatStashDestroy_Private(&mat->stash);
1134: MatStashDestroy_Private(&mat->bstash);
1135: MatDestroy(&baij->A);
1136: MatDestroy(&baij->B);
1137: #if defined(PETSC_USE_CTABLE)
1138: PetscTableDestroy(&baij->colmap);
1139: #else
1140: PetscFree(baij->colmap);
1141: #endif
1142: PetscFree(baij->garray);
1143: VecDestroy(&baij->lvec);
1144: VecScatterDestroy(&baij->Mvctx);
1145: VecDestroy(&baij->slvec0);
1146: VecDestroy(&baij->slvec0b);
1147: VecDestroy(&baij->slvec1);
1148: VecDestroy(&baij->slvec1a);
1149: VecDestroy(&baij->slvec1b);
1150: VecScatterDestroy(&baij->sMvctx);
1151: PetscFree2(baij->rowvalues,baij->rowindices);
1152: PetscFree(baij->barray);
1153: PetscFree(baij->hd);
1154: VecDestroy(&baij->diag);
1155: VecDestroy(&baij->bb1);
1156: VecDestroy(&baij->xx1);
1157: #if defined(PETSC_USE_REAL_MAT_SINGLE)
1158: PetscFree(baij->setvaluescopy);
1159: #endif
1160: PetscFree(baij->in_loc);
1161: PetscFree(baij->v_loc);
1162: PetscFree(baij->rangebs);
1163: PetscFree(mat->data);
1165: PetscObjectChangeTypeName((PetscObject)mat,0);
1166: PetscObjectComposeFunction((PetscObject)mat,"MatStoreValues_C",NULL);
1167: PetscObjectComposeFunction((PetscObject)mat,"MatRetrieveValues_C",NULL);
1168: PetscObjectComposeFunction((PetscObject)mat,"MatGetDiagonalBlock_C",NULL);
1169: PetscObjectComposeFunction((PetscObject)mat,"MatMPISBAIJSetPreallocation_C",NULL);
1170: PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpisbaij_mpisbstrm_C",NULL);
1171: #if defined(PETSC_HAVE_ELEMENTAL)
1172: PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpisbaij_elemental_C",NULL);
1173: #endif
1174: return(0);
1175: }
1179: PetscErrorCode MatMult_MPISBAIJ_Hermitian(Mat A,Vec xx,Vec yy)
1180: {
1181: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1182: PetscErrorCode ierr;
1183: PetscInt nt,mbs=a->mbs,bs=A->rmap->bs;
1184: PetscScalar *from;
1185: const PetscScalar *x;
1188: VecGetLocalSize(xx,&nt);
1189: if (nt != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible partition of A and xx");
1191: /* diagonal part */
1192: (*a->A->ops->mult)(a->A,xx,a->slvec1a);
1193: VecSet(a->slvec1b,0.0);
1195: /* subdiagonal part */
1196: (*a->B->ops->multhermitiantranspose)(a->B,xx,a->slvec0b);
1198: /* copy x into the vec slvec0 */
1199: VecGetArray(a->slvec0,&from);
1200: VecGetArrayRead(xx,&x);
1202: PetscMemcpy(from,x,bs*mbs*sizeof(MatScalar));
1203: VecRestoreArray(a->slvec0,&from);
1204: VecRestoreArrayRead(xx,&x);
1206: VecScatterBegin(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
1207: VecScatterEnd(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
1208: /* supperdiagonal part */
1209: (*a->B->ops->multadd)(a->B,a->slvec1b,a->slvec1a,yy);
1210: return(0);
1211: }
1215: PetscErrorCode MatMult_MPISBAIJ(Mat A,Vec xx,Vec yy)
1216: {
1217: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1218: PetscErrorCode ierr;
1219: PetscInt nt,mbs=a->mbs,bs=A->rmap->bs;
1220: PetscScalar *from;
1221: const PetscScalar *x;
1224: VecGetLocalSize(xx,&nt);
1225: if (nt != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible partition of A and xx");
1227: /* diagonal part */
1228: (*a->A->ops->mult)(a->A,xx,a->slvec1a);
1229: VecSet(a->slvec1b,0.0);
1231: /* subdiagonal part */
1232: (*a->B->ops->multtranspose)(a->B,xx,a->slvec0b);
1234: /* copy x into the vec slvec0 */
1235: VecGetArray(a->slvec0,&from);
1236: VecGetArrayRead(xx,&x);
1238: PetscMemcpy(from,x,bs*mbs*sizeof(MatScalar));
1239: VecRestoreArray(a->slvec0,&from);
1240: VecRestoreArrayRead(xx,&x);
1242: VecScatterBegin(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
1243: VecScatterEnd(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
1244: /* supperdiagonal part */
1245: (*a->B->ops->multadd)(a->B,a->slvec1b,a->slvec1a,yy);
1246: return(0);
1247: }
1251: PetscErrorCode MatMult_MPISBAIJ_2comm(Mat A,Vec xx,Vec yy)
1252: {
1253: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1255: PetscInt nt;
1258: VecGetLocalSize(xx,&nt);
1259: if (nt != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible partition of A and xx");
1261: VecGetLocalSize(yy,&nt);
1262: if (nt != A->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible parition of A and yy");
1264: VecScatterBegin(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);
1265: /* do diagonal part */
1266: (*a->A->ops->mult)(a->A,xx,yy);
1267: /* do supperdiagonal part */
1268: VecScatterEnd(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);
1269: (*a->B->ops->multadd)(a->B,a->lvec,yy,yy);
1270: /* do subdiagonal part */
1271: (*a->B->ops->multtranspose)(a->B,xx,a->lvec);
1272: VecScatterBegin(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);
1273: VecScatterEnd(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);
1274: return(0);
1275: }
1279: PetscErrorCode MatMultAdd_MPISBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
1280: {
1281: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1282: PetscErrorCode ierr;
1283: PetscInt mbs=a->mbs,bs=A->rmap->bs;
1284: PetscScalar *from,zero=0.0;
1285: const PetscScalar *x;
1288: /*
1289: PetscSynchronizedPrintf(PetscObjectComm((PetscObject)A)," MatMultAdd is called ...\n");
1290: PetscSynchronizedFlush(PetscObjectComm((PetscObject)A),PETSC_STDOUT);
1291: */
1292: /* diagonal part */
1293: (*a->A->ops->multadd)(a->A,xx,yy,a->slvec1a);
1294: VecSet(a->slvec1b,zero);
1296: /* subdiagonal part */
1297: (*a->B->ops->multtranspose)(a->B,xx,a->slvec0b);
1299: /* copy x into the vec slvec0 */
1300: VecGetArray(a->slvec0,&from);
1301: VecGetArrayRead(xx,&x);
1302: PetscMemcpy(from,x,bs*mbs*sizeof(MatScalar));
1303: VecRestoreArray(a->slvec0,&from);
1305: VecScatterBegin(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
1306: VecRestoreArrayRead(xx,&x);
1307: VecScatterEnd(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
1309: /* supperdiagonal part */
1310: (*a->B->ops->multadd)(a->B,a->slvec1b,a->slvec1a,zz);
1311: return(0);
1312: }
1316: PetscErrorCode MatMultAdd_MPISBAIJ_2comm(Mat A,Vec xx,Vec yy,Vec zz)
1317: {
1318: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1322: VecScatterBegin(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);
1323: /* do diagonal part */
1324: (*a->A->ops->multadd)(a->A,xx,yy,zz);
1325: /* do supperdiagonal part */
1326: VecScatterEnd(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);
1327: (*a->B->ops->multadd)(a->B,a->lvec,zz,zz);
1329: /* do subdiagonal part */
1330: (*a->B->ops->multtranspose)(a->B,xx,a->lvec);
1331: VecScatterBegin(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);
1332: VecScatterEnd(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);
1333: return(0);
1334: }
1336: /*
1337: This only works correctly for square matrices where the subblock A->A is the
1338: diagonal block
1339: */
1342: PetscErrorCode MatGetDiagonal_MPISBAIJ(Mat A,Vec v)
1343: {
1344: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1348: /* if (a->rmap->N != a->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Supports only square matrix where A->A is diag block"); */
1349: MatGetDiagonal(a->A,v);
1350: return(0);
1351: }
1355: PetscErrorCode MatScale_MPISBAIJ(Mat A,PetscScalar aa)
1356: {
1357: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1361: MatScale(a->A,aa);
1362: MatScale(a->B,aa);
1363: return(0);
1364: }
1368: PetscErrorCode MatGetRow_MPISBAIJ(Mat matin,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
1369: {
1370: Mat_MPISBAIJ *mat = (Mat_MPISBAIJ*)matin->data;
1371: PetscScalar *vworkA,*vworkB,**pvA,**pvB,*v_p;
1373: PetscInt bs = matin->rmap->bs,bs2 = mat->bs2,i,*cworkA,*cworkB,**pcA,**pcB;
1374: PetscInt nztot,nzA,nzB,lrow,brstart = matin->rmap->rstart,brend = matin->rmap->rend;
1375: PetscInt *cmap,*idx_p,cstart = mat->rstartbs;
1378: if (mat->getrowactive) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Already active");
1379: mat->getrowactive = PETSC_TRUE;
1381: if (!mat->rowvalues && (idx || v)) {
1382: /*
1383: allocate enough space to hold information from the longest row.
1384: */
1385: Mat_SeqSBAIJ *Aa = (Mat_SeqSBAIJ*)mat->A->data;
1386: Mat_SeqBAIJ *Ba = (Mat_SeqBAIJ*)mat->B->data;
1387: PetscInt max = 1,mbs = mat->mbs,tmp;
1388: for (i=0; i<mbs; i++) {
1389: tmp = Aa->i[i+1] - Aa->i[i] + Ba->i[i+1] - Ba->i[i]; /* row length */
1390: if (max < tmp) max = tmp;
1391: }
1392: PetscMalloc2(max*bs2,&mat->rowvalues,max*bs2,&mat->rowindices);
1393: }
1395: if (row < brstart || row >= brend) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only local rows");
1396: lrow = row - brstart; /* local row index */
1398: pvA = &vworkA; pcA = &cworkA; pvB = &vworkB; pcB = &cworkB;
1399: if (!v) {pvA = 0; pvB = 0;}
1400: if (!idx) {pcA = 0; if (!v) pcB = 0;}
1401: (*mat->A->ops->getrow)(mat->A,lrow,&nzA,pcA,pvA);
1402: (*mat->B->ops->getrow)(mat->B,lrow,&nzB,pcB,pvB);
1403: nztot = nzA + nzB;
1405: cmap = mat->garray;
1406: if (v || idx) {
1407: if (nztot) {
1408: /* Sort by increasing column numbers, assuming A and B already sorted */
1409: PetscInt imark = -1;
1410: if (v) {
1411: *v = v_p = mat->rowvalues;
1412: for (i=0; i<nzB; i++) {
1413: if (cmap[cworkB[i]/bs] < cstart) v_p[i] = vworkB[i];
1414: else break;
1415: }
1416: imark = i;
1417: for (i=0; i<nzA; i++) v_p[imark+i] = vworkA[i];
1418: for (i=imark; i<nzB; i++) v_p[nzA+i] = vworkB[i];
1419: }
1420: if (idx) {
1421: *idx = idx_p = mat->rowindices;
1422: if (imark > -1) {
1423: for (i=0; i<imark; i++) {
1424: idx_p[i] = cmap[cworkB[i]/bs]*bs + cworkB[i]%bs;
1425: }
1426: } else {
1427: for (i=0; i<nzB; i++) {
1428: if (cmap[cworkB[i]/bs] < cstart) idx_p[i] = cmap[cworkB[i]/bs]*bs + cworkB[i]%bs;
1429: else break;
1430: }
1431: imark = i;
1432: }
1433: for (i=0; i<nzA; i++) idx_p[imark+i] = cstart*bs + cworkA[i];
1434: for (i=imark; i<nzB; i++) idx_p[nzA+i] = cmap[cworkB[i]/bs]*bs + cworkB[i]%bs ;
1435: }
1436: } else {
1437: if (idx) *idx = 0;
1438: if (v) *v = 0;
1439: }
1440: }
1441: *nz = nztot;
1442: (*mat->A->ops->restorerow)(mat->A,lrow,&nzA,pcA,pvA);
1443: (*mat->B->ops->restorerow)(mat->B,lrow,&nzB,pcB,pvB);
1444: return(0);
1445: }
1449: PetscErrorCode MatRestoreRow_MPISBAIJ(Mat mat,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
1450: {
1451: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
1454: if (!baij->getrowactive) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"MatGetRow() must be called first");
1455: baij->getrowactive = PETSC_FALSE;
1456: return(0);
1457: }
1461: PetscErrorCode MatGetRowUpperTriangular_MPISBAIJ(Mat A)
1462: {
1463: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1464: Mat_SeqSBAIJ *aA = (Mat_SeqSBAIJ*)a->A->data;
1467: aA->getrow_utriangular = PETSC_TRUE;
1468: return(0);
1469: }
1472: PetscErrorCode MatRestoreRowUpperTriangular_MPISBAIJ(Mat A)
1473: {
1474: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1475: Mat_SeqSBAIJ *aA = (Mat_SeqSBAIJ*)a->A->data;
1478: aA->getrow_utriangular = PETSC_FALSE;
1479: return(0);
1480: }
1484: PetscErrorCode MatRealPart_MPISBAIJ(Mat A)
1485: {
1486: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1490: MatRealPart(a->A);
1491: MatRealPart(a->B);
1492: return(0);
1493: }
1497: PetscErrorCode MatImaginaryPart_MPISBAIJ(Mat A)
1498: {
1499: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1503: MatImaginaryPart(a->A);
1504: MatImaginaryPart(a->B);
1505: return(0);
1506: }
1508: /* Check if isrow is a subset of iscol_local, called by MatGetSubMatrix_MPISBAIJ()
1509: Input: isrow - distributed(parallel),
1510: iscol_local - locally owned (seq)
1511: */
1514: PetscErrorCode ISEqual_private(IS isrow,IS iscol_local,PetscBool *flg)
1515: {
1517: PetscInt sz1,sz2,*a1,*a2,i,j,k,nmatch;
1518: const PetscInt *ptr1,*ptr2;
1521: ISGetLocalSize(isrow,&sz1);
1522: ISGetLocalSize(iscol_local,&sz2);
1523: if (sz1 > sz2) {
1524: *flg = PETSC_FALSE;
1525: return(0);
1526: }
1528: ISGetIndices(isrow,&ptr1);
1529: ISGetIndices(iscol_local,&ptr2);
1531: PetscMalloc1(sz1,&a1);
1532: PetscMalloc1(sz2,&a2);
1533: PetscMemcpy(a1,ptr1,sz1*sizeof(PetscInt));
1534: PetscMemcpy(a2,ptr2,sz2*sizeof(PetscInt));
1535: PetscSortInt(sz1,a1);
1536: PetscSortInt(sz2,a2);
1538: nmatch=0;
1539: k = 0;
1540: for (i=0; i<sz1; i++){
1541: for (j=k; j<sz2; j++){
1542: if (a1[i] == a2[j]) {
1543: k = j; nmatch++;
1544: break;
1545: }
1546: }
1547: }
1548: ISRestoreIndices(isrow,&ptr1);
1549: ISRestoreIndices(iscol_local,&ptr2);
1550: PetscFree(a1);
1551: PetscFree(a2);
1552: if (nmatch < sz1) {
1553: *flg = PETSC_FALSE;
1554: } else {
1555: *flg = PETSC_TRUE;
1556: }
1557: return(0);
1558: }
1562: PetscErrorCode MatGetSubMatrix_MPISBAIJ(Mat mat,IS isrow,IS iscol,MatReuse call,Mat *newmat)
1563: {
1565: IS iscol_local;
1566: PetscInt csize;
1567: PetscBool isequal;
1570: ISGetLocalSize(iscol,&csize);
1571: if (call == MAT_REUSE_MATRIX) {
1572: PetscObjectQuery((PetscObject)*newmat,"ISAllGather",(PetscObject*)&iscol_local);
1573: if (!iscol_local) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Submatrix passed in was not used before, cannot reuse");
1574: } else {
1575: ISAllGather(iscol,&iscol_local);
1576: ISEqual_private(isrow,iscol_local,&isequal);
1577: if (!isequal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"For symmetric format, iscol must equal isrow");
1578: }
1580: /* now call MatGetSubMatrix_MPIBAIJ() */
1581: MatGetSubMatrix_MPIBAIJ_Private(mat,isrow,iscol_local,csize,call,newmat);
1582: if (call == MAT_INITIAL_MATRIX) {
1583: PetscObjectCompose((PetscObject)*newmat,"ISAllGather",(PetscObject)iscol_local);
1584: ISDestroy(&iscol_local);
1585: }
1586: return(0);
1587: }
1591: PetscErrorCode MatZeroEntries_MPISBAIJ(Mat A)
1592: {
1593: Mat_MPISBAIJ *l = (Mat_MPISBAIJ*)A->data;
1597: MatZeroEntries(l->A);
1598: MatZeroEntries(l->B);
1599: return(0);
1600: }
1604: PetscErrorCode MatGetInfo_MPISBAIJ(Mat matin,MatInfoType flag,MatInfo *info)
1605: {
1606: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)matin->data;
1607: Mat A = a->A,B = a->B;
1609: PetscReal isend[5],irecv[5];
1612: info->block_size = (PetscReal)matin->rmap->bs;
1614: MatGetInfo(A,MAT_LOCAL,info);
1616: isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
1617: isend[3] = info->memory; isend[4] = info->mallocs;
1619: MatGetInfo(B,MAT_LOCAL,info);
1621: isend[0] += info->nz_used; isend[1] += info->nz_allocated; isend[2] += info->nz_unneeded;
1622: isend[3] += info->memory; isend[4] += info->mallocs;
1623: if (flag == MAT_LOCAL) {
1624: info->nz_used = isend[0];
1625: info->nz_allocated = isend[1];
1626: info->nz_unneeded = isend[2];
1627: info->memory = isend[3];
1628: info->mallocs = isend[4];
1629: } else if (flag == MAT_GLOBAL_MAX) {
1630: MPIU_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)matin));
1632: info->nz_used = irecv[0];
1633: info->nz_allocated = irecv[1];
1634: info->nz_unneeded = irecv[2];
1635: info->memory = irecv[3];
1636: info->mallocs = irecv[4];
1637: } else if (flag == MAT_GLOBAL_SUM) {
1638: MPIU_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)matin));
1640: info->nz_used = irecv[0];
1641: info->nz_allocated = irecv[1];
1642: info->nz_unneeded = irecv[2];
1643: info->memory = irecv[3];
1644: info->mallocs = irecv[4];
1645: } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown MatInfoType argument %d",(int)flag);
1646: info->fill_ratio_given = 0; /* no parallel LU/ILU/Cholesky */
1647: info->fill_ratio_needed = 0;
1648: info->factor_mallocs = 0;
1649: return(0);
1650: }
1654: PetscErrorCode MatSetOption_MPISBAIJ(Mat A,MatOption op,PetscBool flg)
1655: {
1656: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1657: Mat_SeqSBAIJ *aA = (Mat_SeqSBAIJ*)a->A->data;
1661: switch (op) {
1662: case MAT_NEW_NONZERO_LOCATIONS:
1663: case MAT_NEW_NONZERO_ALLOCATION_ERR:
1664: case MAT_UNUSED_NONZERO_LOCATION_ERR:
1665: case MAT_KEEP_NONZERO_PATTERN:
1666: case MAT_NEW_NONZERO_LOCATION_ERR:
1667: MatCheckPreallocated(A,1);
1668: MatSetOption(a->A,op,flg);
1669: MatSetOption(a->B,op,flg);
1670: break;
1671: case MAT_ROW_ORIENTED:
1672: MatCheckPreallocated(A,1);
1673: a->roworiented = flg;
1675: MatSetOption(a->A,op,flg);
1676: MatSetOption(a->B,op,flg);
1677: break;
1678: case MAT_NEW_DIAGONALS:
1679: PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
1680: break;
1681: case MAT_IGNORE_OFF_PROC_ENTRIES:
1682: a->donotstash = flg;
1683: break;
1684: case MAT_USE_HASH_TABLE:
1685: a->ht_flag = flg;
1686: break;
1687: case MAT_HERMITIAN:
1688: MatCheckPreallocated(A,1);
1689: if (!A->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatAssemblyEnd() first");
1690: MatSetOption(a->A,op,flg);
1692: A->ops->mult = MatMult_MPISBAIJ_Hermitian;
1693: break;
1694: case MAT_SPD:
1695: A->spd_set = PETSC_TRUE;
1696: A->spd = flg;
1697: if (flg) {
1698: A->symmetric = PETSC_TRUE;
1699: A->structurally_symmetric = PETSC_TRUE;
1700: A->symmetric_set = PETSC_TRUE;
1701: A->structurally_symmetric_set = PETSC_TRUE;
1702: }
1703: break;
1704: case MAT_SYMMETRIC:
1705: MatCheckPreallocated(A,1);
1706: MatSetOption(a->A,op,flg);
1707: break;
1708: case MAT_STRUCTURALLY_SYMMETRIC:
1709: MatCheckPreallocated(A,1);
1710: MatSetOption(a->A,op,flg);
1711: break;
1712: case MAT_SYMMETRY_ETERNAL:
1713: if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix must be symmetric");
1714: PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
1715: break;
1716: case MAT_IGNORE_LOWER_TRIANGULAR:
1717: aA->ignore_ltriangular = flg;
1718: break;
1719: case MAT_ERROR_LOWER_TRIANGULAR:
1720: aA->ignore_ltriangular = flg;
1721: break;
1722: case MAT_GETROW_UPPERTRIANGULAR:
1723: aA->getrow_utriangular = flg;
1724: break;
1725: default:
1726: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %d",op);
1727: }
1728: return(0);
1729: }
1733: PetscErrorCode MatTranspose_MPISBAIJ(Mat A,MatReuse reuse,Mat *B)
1734: {
1738: if (MAT_INITIAL_MATRIX || *B != A) {
1739: MatDuplicate(A,MAT_COPY_VALUES,B);
1740: }
1741: return(0);
1742: }
1746: PetscErrorCode MatDiagonalScale_MPISBAIJ(Mat mat,Vec ll,Vec rr)
1747: {
1748: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
1749: Mat a = baij->A, b=baij->B;
1751: PetscInt nv,m,n;
1752: PetscBool flg;
1755: if (ll != rr) {
1756: VecEqual(ll,rr,&flg);
1757: if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"For symmetric format, left and right scaling vectors must be same\n");
1758: }
1759: if (!ll) return(0);
1761: MatGetLocalSize(mat,&m,&n);
1762: if (m != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"For symmetric format, local size %d %d must be same",m,n);
1764: VecGetLocalSize(rr,&nv);
1765: if (nv!=n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left and right vector non-conforming local size");
1767: VecScatterBegin(baij->Mvctx,rr,baij->lvec,INSERT_VALUES,SCATTER_FORWARD);
1769: /* left diagonalscale the off-diagonal part */
1770: (*b->ops->diagonalscale)(b,ll,NULL);
1772: /* scale the diagonal part */
1773: (*a->ops->diagonalscale)(a,ll,rr);
1775: /* right diagonalscale the off-diagonal part */
1776: VecScatterEnd(baij->Mvctx,rr,baij->lvec,INSERT_VALUES,SCATTER_FORWARD);
1777: (*b->ops->diagonalscale)(b,NULL,baij->lvec);
1778: return(0);
1779: }
1783: PetscErrorCode MatSetUnfactored_MPISBAIJ(Mat A)
1784: {
1785: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1789: MatSetUnfactored(a->A);
1790: return(0);
1791: }
1793: static PetscErrorCode MatDuplicate_MPISBAIJ(Mat,MatDuplicateOption,Mat*);
1797: PetscErrorCode MatEqual_MPISBAIJ(Mat A,Mat B,PetscBool *flag)
1798: {
1799: Mat_MPISBAIJ *matB = (Mat_MPISBAIJ*)B->data,*matA = (Mat_MPISBAIJ*)A->data;
1800: Mat a,b,c,d;
1801: PetscBool flg;
1805: a = matA->A; b = matA->B;
1806: c = matB->A; d = matB->B;
1808: MatEqual(a,c,&flg);
1809: if (flg) {
1810: MatEqual(b,d,&flg);
1811: }
1812: MPIU_Allreduce(&flg,flag,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));
1813: return(0);
1814: }
1818: PetscErrorCode MatCopy_MPISBAIJ(Mat A,Mat B,MatStructure str)
1819: {
1821: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1822: Mat_MPISBAIJ *b = (Mat_MPISBAIJ*)B->data;
1825: /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
1826: if ((str != SAME_NONZERO_PATTERN) || (A->ops->copy != B->ops->copy)) {
1827: MatGetRowUpperTriangular(A);
1828: MatCopy_Basic(A,B,str);
1829: MatRestoreRowUpperTriangular(A);
1830: } else {
1831: MatCopy(a->A,b->A,str);
1832: MatCopy(a->B,b->B,str);
1833: }
1834: return(0);
1835: }
1839: PetscErrorCode MatSetUp_MPISBAIJ(Mat A)
1840: {
1844: MatMPISBAIJSetPreallocation(A,A->rmap->bs,PETSC_DEFAULT,0,PETSC_DEFAULT,0);
1845: return(0);
1846: }
1850: PetscErrorCode MatAXPY_MPISBAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str)
1851: {
1853: Mat_MPISBAIJ *xx=(Mat_MPISBAIJ*)X->data,*yy=(Mat_MPISBAIJ*)Y->data;
1854: PetscBLASInt bnz,one=1;
1855: Mat_SeqSBAIJ *xa,*ya;
1856: Mat_SeqBAIJ *xb,*yb;
1859: if (str == SAME_NONZERO_PATTERN) {
1860: PetscScalar alpha = a;
1861: xa = (Mat_SeqSBAIJ*)xx->A->data;
1862: ya = (Mat_SeqSBAIJ*)yy->A->data;
1863: PetscBLASIntCast(xa->nz,&bnz);
1864: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&bnz,&alpha,xa->a,&one,ya->a,&one));
1865: xb = (Mat_SeqBAIJ*)xx->B->data;
1866: yb = (Mat_SeqBAIJ*)yy->B->data;
1867: PetscBLASIntCast(xb->nz,&bnz);
1868: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&bnz,&alpha,xb->a,&one,yb->a,&one));
1869: PetscObjectStateIncrease((PetscObject)Y);
1870: } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
1871: MatSetOption(X,MAT_GETROW_UPPERTRIANGULAR,PETSC_TRUE);
1872: MatAXPY_Basic(Y,a,X,str);
1873: MatSetOption(X,MAT_GETROW_UPPERTRIANGULAR,PETSC_FALSE);
1874: } else {
1875: Mat B;
1876: PetscInt *nnz_d,*nnz_o,bs=Y->rmap->bs;
1877: if (bs != X->rmap->bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrices must have same block size");
1878: MatGetRowUpperTriangular(X);
1879: MatGetRowUpperTriangular(Y);
1880: PetscMalloc1(yy->A->rmap->N,&nnz_d);
1881: PetscMalloc1(yy->B->rmap->N,&nnz_o);
1882: MatCreate(PetscObjectComm((PetscObject)Y),&B);
1883: PetscObjectSetName((PetscObject)B,((PetscObject)Y)->name);
1884: MatSetSizes(B,Y->rmap->n,Y->cmap->n,Y->rmap->N,Y->cmap->N);
1885: MatSetBlockSizesFromMats(B,Y,Y);
1886: MatSetType(B,MATMPISBAIJ);
1887: MatAXPYGetPreallocation_SeqSBAIJ(yy->A,xx->A,nnz_d);
1888: MatAXPYGetPreallocation_MPIBAIJ(yy->B,yy->garray,xx->B,xx->garray,nnz_o);
1889: MatMPISBAIJSetPreallocation(B,bs,0,nnz_d,0,nnz_o);
1890: MatAXPY_BasicWithPreallocation(B,Y,a,X,str);
1891: MatHeaderReplace(Y,&B);
1892: PetscFree(nnz_d);
1893: PetscFree(nnz_o);
1894: MatRestoreRowUpperTriangular(X);
1895: MatRestoreRowUpperTriangular(Y);
1896: }
1897: return(0);
1898: }
1902: PetscErrorCode MatGetSubMatrices_MPISBAIJ(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1903: {
1905: PetscInt i;
1906: PetscBool flg;
1909: MatGetSubMatrices_MPIBAIJ(A,n,irow,icol,scall,B); /* B[] are sbaij matrices */
1910: for (i=0; i<n; i++) {
1911: ISEqual(irow[i],icol[i],&flg);
1912: if (!flg) {
1913: MatSeqSBAIJZeroOps_Private(*B[i]);
1914: }
1915: }
1916: return(0);
1917: }
1921: PetscErrorCode MatShift_MPISBAIJ(Mat Y,PetscScalar a)
1922: {
1924: Mat_MPISBAIJ *maij = (Mat_MPISBAIJ*)Y->data;
1925: Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ*)maij->A->data;
1928: if (!Y->preallocated) {
1929: MatMPISBAIJSetPreallocation(Y,Y->rmap->bs,1,NULL,0,NULL);
1930: } else if (!aij->nz) {
1931: PetscInt nonew = aij->nonew;
1932: MatSeqSBAIJSetPreallocation(maij->A,Y->rmap->bs,1,NULL);
1933: aij->nonew = nonew;
1934: }
1935: MatShift_Basic(Y,a);
1936: return(0);
1937: }
1941: PetscErrorCode MatMissingDiagonal_MPISBAIJ(Mat A,PetscBool *missing,PetscInt *d)
1942: {
1943: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1947: if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only works for square matrices");
1948: MatMissingDiagonal(a->A,missing,d);
1949: if (d) {
1950: PetscInt rstart;
1951: MatGetOwnershipRange(A,&rstart,NULL);
1952: *d += rstart/A->rmap->bs;
1954: }
1955: return(0);
1956: }
1959: /* -------------------------------------------------------------------*/
1960: static struct _MatOps MatOps_Values = {MatSetValues_MPISBAIJ,
1961: MatGetRow_MPISBAIJ,
1962: MatRestoreRow_MPISBAIJ,
1963: MatMult_MPISBAIJ,
1964: /* 4*/ MatMultAdd_MPISBAIJ,
1965: MatMult_MPISBAIJ, /* transpose versions are same as non-transpose */
1966: MatMultAdd_MPISBAIJ,
1967: 0,
1968: 0,
1969: 0,
1970: /* 10*/ 0,
1971: 0,
1972: 0,
1973: MatSOR_MPISBAIJ,
1974: MatTranspose_MPISBAIJ,
1975: /* 15*/ MatGetInfo_MPISBAIJ,
1976: MatEqual_MPISBAIJ,
1977: MatGetDiagonal_MPISBAIJ,
1978: MatDiagonalScale_MPISBAIJ,
1979: MatNorm_MPISBAIJ,
1980: /* 20*/ MatAssemblyBegin_MPISBAIJ,
1981: MatAssemblyEnd_MPISBAIJ,
1982: MatSetOption_MPISBAIJ,
1983: MatZeroEntries_MPISBAIJ,
1984: /* 24*/ 0,
1985: 0,
1986: 0,
1987: 0,
1988: 0,
1989: /* 29*/ MatSetUp_MPISBAIJ,
1990: 0,
1991: 0,
1992: 0,
1993: 0,
1994: /* 34*/ MatDuplicate_MPISBAIJ,
1995: 0,
1996: 0,
1997: 0,
1998: 0,
1999: /* 39*/ MatAXPY_MPISBAIJ,
2000: MatGetSubMatrices_MPISBAIJ,
2001: MatIncreaseOverlap_MPISBAIJ,
2002: MatGetValues_MPISBAIJ,
2003: MatCopy_MPISBAIJ,
2004: /* 44*/ 0,
2005: MatScale_MPISBAIJ,
2006: MatShift_MPISBAIJ,
2007: 0,
2008: 0,
2009: /* 49*/ 0,
2010: 0,
2011: 0,
2012: 0,
2013: 0,
2014: /* 54*/ 0,
2015: 0,
2016: MatSetUnfactored_MPISBAIJ,
2017: 0,
2018: MatSetValuesBlocked_MPISBAIJ,
2019: /* 59*/ MatGetSubMatrix_MPISBAIJ,
2020: 0,
2021: 0,
2022: 0,
2023: 0,
2024: /* 64*/ 0,
2025: 0,
2026: 0,
2027: 0,
2028: 0,
2029: /* 69*/ MatGetRowMaxAbs_MPISBAIJ,
2030: 0,
2031: 0,
2032: 0,
2033: 0,
2034: /* 74*/ 0,
2035: 0,
2036: 0,
2037: 0,
2038: 0,
2039: /* 79*/ 0,
2040: 0,
2041: 0,
2042: 0,
2043: MatLoad_MPISBAIJ,
2044: /* 84*/ 0,
2045: 0,
2046: 0,
2047: 0,
2048: 0,
2049: /* 89*/ 0,
2050: 0,
2051: 0,
2052: 0,
2053: 0,
2054: /* 94*/ 0,
2055: 0,
2056: 0,
2057: 0,
2058: 0,
2059: /* 99*/ 0,
2060: 0,
2061: 0,
2062: 0,
2063: 0,
2064: /*104*/ 0,
2065: MatRealPart_MPISBAIJ,
2066: MatImaginaryPart_MPISBAIJ,
2067: MatGetRowUpperTriangular_MPISBAIJ,
2068: MatRestoreRowUpperTriangular_MPISBAIJ,
2069: /*109*/ 0,
2070: 0,
2071: 0,
2072: 0,
2073: MatMissingDiagonal_MPISBAIJ,
2074: /*114*/ 0,
2075: 0,
2076: 0,
2077: 0,
2078: 0,
2079: /*119*/ 0,
2080: 0,
2081: 0,
2082: 0,
2083: 0,
2084: /*124*/ 0,
2085: 0,
2086: 0,
2087: 0,
2088: 0,
2089: /*129*/ 0,
2090: 0,
2091: 0,
2092: 0,
2093: 0,
2094: /*134*/ 0,
2095: 0,
2096: 0,
2097: 0,
2098: 0,
2099: /*139*/ 0,
2100: 0,
2101: 0,
2102: 0,
2103: 0,
2104: /*144*/MatCreateMPIMatConcatenateSeqMat_MPISBAIJ
2105: };
2109: PetscErrorCode MatGetDiagonalBlock_MPISBAIJ(Mat A,Mat *a)
2110: {
2112: *a = ((Mat_MPISBAIJ*)A->data)->A;
2113: return(0);
2114: }
2118: PetscErrorCode MatMPISBAIJSetPreallocation_MPISBAIJ(Mat B,PetscInt bs,PetscInt d_nz,const PetscInt *d_nnz,PetscInt o_nz,const PetscInt *o_nnz)
2119: {
2120: Mat_MPISBAIJ *b;
2122: PetscInt i,mbs,Mbs;
2125: MatSetBlockSize(B,PetscAbs(bs));
2126: PetscLayoutSetUp(B->rmap);
2127: PetscLayoutSetUp(B->cmap);
2128: PetscLayoutGetBlockSize(B->rmap,&bs);
2130: b = (Mat_MPISBAIJ*)B->data;
2131: mbs = B->rmap->n/bs;
2132: Mbs = B->rmap->N/bs;
2133: 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);
2135: B->rmap->bs = bs;
2136: b->bs2 = bs*bs;
2137: b->mbs = mbs;
2138: b->Mbs = Mbs;
2139: b->nbs = B->cmap->n/bs;
2140: b->Nbs = B->cmap->N/bs;
2142: for (i=0; i<=b->size; i++) {
2143: b->rangebs[i] = B->rmap->range[i]/bs;
2144: }
2145: b->rstartbs = B->rmap->rstart/bs;
2146: b->rendbs = B->rmap->rend/bs;
2148: b->cstartbs = B->cmap->rstart/bs;
2149: b->cendbs = B->cmap->rend/bs;
2151: if (!B->preallocated) {
2152: MatCreate(PETSC_COMM_SELF,&b->A);
2153: MatSetSizes(b->A,B->rmap->n,B->cmap->n,B->rmap->n,B->cmap->n);
2154: MatSetType(b->A,MATSEQSBAIJ);
2155: PetscLogObjectParent((PetscObject)B,(PetscObject)b->A);
2156: MatCreate(PETSC_COMM_SELF,&b->B);
2157: MatSetSizes(b->B,B->rmap->n,B->cmap->N,B->rmap->n,B->cmap->N);
2158: MatSetType(b->B,MATSEQBAIJ);
2159: PetscLogObjectParent((PetscObject)B,(PetscObject)b->B);
2160: MatStashCreate_Private(PetscObjectComm((PetscObject)B),bs,&B->bstash);
2161: }
2163: MatSeqSBAIJSetPreallocation(b->A,bs,d_nz,d_nnz);
2164: MatSeqBAIJSetPreallocation(b->B,bs,o_nz,o_nnz);
2166: B->preallocated = PETSC_TRUE;
2167: return(0);
2168: }
2172: PetscErrorCode MatMPISBAIJSetPreallocationCSR_MPISBAIJ(Mat B,PetscInt bs,const PetscInt ii[],const PetscInt jj[],const PetscScalar V[])
2173: {
2174: PetscInt m,rstart,cstart,cend;
2175: PetscInt i,j,d,nz,nz_max=0,*d_nnz=0,*o_nnz=0;
2176: const PetscInt *JJ =0;
2177: PetscScalar *values=0;
2181: if (bs < 1) SETERRQ1(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_OUTOFRANGE,"Invalid block size specified, must be positive but it is %D",bs);
2182: PetscLayoutSetBlockSize(B->rmap,bs);
2183: PetscLayoutSetBlockSize(B->cmap,bs);
2184: PetscLayoutSetUp(B->rmap);
2185: PetscLayoutSetUp(B->cmap);
2186: PetscLayoutGetBlockSize(B->rmap,&bs);
2187: m = B->rmap->n/bs;
2188: rstart = B->rmap->rstart/bs;
2189: cstart = B->cmap->rstart/bs;
2190: cend = B->cmap->rend/bs;
2192: if (ii[0]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"ii[0] must be 0 but it is %D",ii[0]);
2193: PetscMalloc2(m,&d_nnz,m,&o_nnz);
2194: for (i=0; i<m; i++) {
2195: nz = ii[i+1] - ii[i];
2196: if (nz < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local row %D has a negative number of columns %D",i,nz);
2197: nz_max = PetscMax(nz_max,nz);
2198: JJ = jj + ii[i];
2199: for (j=0; j<nz; j++) {
2200: if (*JJ >= cstart) break;
2201: JJ++;
2202: }
2203: d = 0;
2204: for (; j<nz; j++) {
2205: if (*JJ++ >= cend) break;
2206: d++;
2207: }
2208: d_nnz[i] = d;
2209: o_nnz[i] = nz - d;
2210: }
2211: MatMPISBAIJSetPreallocation(B,bs,0,d_nnz,0,o_nnz);
2212: PetscFree2(d_nnz,o_nnz);
2214: values = (PetscScalar*)V;
2215: if (!values) {
2216: PetscMalloc1(bs*bs*nz_max,&values);
2217: PetscMemzero(values,bs*bs*nz_max*sizeof(PetscScalar));
2218: }
2219: for (i=0; i<m; i++) {
2220: PetscInt row = i + rstart;
2221: PetscInt ncols = ii[i+1] - ii[i];
2222: const PetscInt *icols = jj + ii[i];
2223: const PetscScalar *svals = values + (V ? (bs*bs*ii[i]) : 0);
2224: MatSetValuesBlocked_MPISBAIJ(B,1,&row,ncols,icols,svals,INSERT_VALUES);
2225: }
2227: if (!V) { PetscFree(values); }
2228: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
2229: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
2230: MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
2231: return(0);
2232: }
2234: /*MC
2235: MATMPISBAIJ - MATMPISBAIJ = "mpisbaij" - A matrix type to be used for distributed symmetric sparse block matrices,
2236: based on block compressed sparse row format. Only the upper triangular portion of the "diagonal" portion of
2237: the matrix is stored.
2239: For complex numbers by default this matrix is symmetric, NOT Hermitian symmetric. To make it Hermitian symmetric you
2240: can call MatSetOption(Mat, MAT_HERMITIAN);
2242: Options Database Keys:
2243: . -mat_type mpisbaij - sets the matrix type to "mpisbaij" during a call to MatSetFromOptions()
2245: Level: beginner
2247: .seealso: MatCreateMPISBAIJ
2248: M*/
2250: PETSC_INTERN PetscErrorCode MatConvert_MPISBAIJ_MPISBSTRM(Mat,MatType,MatReuse,Mat*);
2254: PETSC_EXTERN PetscErrorCode MatCreate_MPISBAIJ(Mat B)
2255: {
2256: Mat_MPISBAIJ *b;
2258: PetscBool flg = PETSC_FALSE;
2261: PetscNewLog(B,&b);
2262: B->data = (void*)b;
2263: PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
2265: B->ops->destroy = MatDestroy_MPISBAIJ;
2266: B->ops->view = MatView_MPISBAIJ;
2267: B->assembled = PETSC_FALSE;
2268: B->insertmode = NOT_SET_VALUES;
2270: MPI_Comm_rank(PetscObjectComm((PetscObject)B),&b->rank);
2271: MPI_Comm_size(PetscObjectComm((PetscObject)B),&b->size);
2273: /* build local table of row and column ownerships */
2274: PetscMalloc1(b->size+2,&b->rangebs);
2276: /* build cache for off array entries formed */
2277: MatStashCreate_Private(PetscObjectComm((PetscObject)B),1,&B->stash);
2279: b->donotstash = PETSC_FALSE;
2280: b->colmap = NULL;
2281: b->garray = NULL;
2282: b->roworiented = PETSC_TRUE;
2284: /* stuff used in block assembly */
2285: b->barray = 0;
2287: /* stuff used for matrix vector multiply */
2288: b->lvec = 0;
2289: b->Mvctx = 0;
2290: b->slvec0 = 0;
2291: b->slvec0b = 0;
2292: b->slvec1 = 0;
2293: b->slvec1a = 0;
2294: b->slvec1b = 0;
2295: b->sMvctx = 0;
2297: /* stuff for MatGetRow() */
2298: b->rowindices = 0;
2299: b->rowvalues = 0;
2300: b->getrowactive = PETSC_FALSE;
2302: /* hash table stuff */
2303: b->ht = 0;
2304: b->hd = 0;
2305: b->ht_size = 0;
2306: b->ht_flag = PETSC_FALSE;
2307: b->ht_fact = 0;
2308: b->ht_total_ct = 0;
2309: b->ht_insert_ct = 0;
2311: /* stuff for MatGetSubMatrices_MPIBAIJ_local() */
2312: b->ijonly = PETSC_FALSE;
2314: b->in_loc = 0;
2315: b->v_loc = 0;
2316: b->n_loc = 0;
2318: PetscObjectComposeFunction((PetscObject)B,"MatStoreValues_C",MatStoreValues_MPISBAIJ);
2319: PetscObjectComposeFunction((PetscObject)B,"MatRetrieveValues_C",MatRetrieveValues_MPISBAIJ);
2320: PetscObjectComposeFunction((PetscObject)B,"MatGetDiagonalBlock_C",MatGetDiagonalBlock_MPISBAIJ);
2321: PetscObjectComposeFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",MatMPISBAIJSetPreallocation_MPISBAIJ);
2322: PetscObjectComposeFunction((PetscObject)B,"MatMPISBAIJSetPreallocationCSR_C",MatMPISBAIJSetPreallocationCSR_MPISBAIJ);
2323: PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpisbaij_mpisbstrm_C",MatConvert_MPISBAIJ_MPISBSTRM);
2324: #if defined(PETSC_HAVE_ELEMENTAL)
2325: PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpisbaij_elemental_C",MatConvert_MPISBAIJ_Elemental);
2326: #endif
2328: B->symmetric = PETSC_TRUE;
2329: B->structurally_symmetric = PETSC_TRUE;
2330: B->symmetric_set = PETSC_TRUE;
2331: B->structurally_symmetric_set = PETSC_TRUE;
2333: PetscObjectChangeTypeName((PetscObject)B,MATMPISBAIJ);
2334: PetscOptionsBegin(PetscObjectComm((PetscObject)B),NULL,"Options for loading MPISBAIJ matrix 1","Mat");
2335: PetscOptionsBool("-mat_use_hash_table","Use hash table to save memory in constructing matrix","MatSetOption",flg,&flg,NULL);
2336: if (flg) {
2337: PetscReal fact = 1.39;
2338: MatSetOption(B,MAT_USE_HASH_TABLE,PETSC_TRUE);
2339: PetscOptionsReal("-mat_use_hash_table","Use hash table factor","MatMPIBAIJSetHashTableFactor",fact,&fact,NULL);
2340: if (fact <= 1.0) fact = 1.39;
2341: MatMPIBAIJSetHashTableFactor(B,fact);
2342: PetscInfo1(B,"Hash table Factor used %5.2f\n",fact);
2343: }
2344: PetscOptionsEnd();
2345: return(0);
2346: }
2348: /*MC
2349: MATSBAIJ - MATSBAIJ = "sbaij" - A matrix type to be used for symmetric block sparse matrices.
2351: This matrix type is identical to MATSEQSBAIJ when constructed with a single process communicator,
2352: and MATMPISBAIJ otherwise.
2354: Options Database Keys:
2355: . -mat_type sbaij - sets the matrix type to "sbaij" during a call to MatSetFromOptions()
2357: Level: beginner
2359: .seealso: MatCreateMPISBAIJ,MATSEQSBAIJ,MATMPISBAIJ
2360: M*/
2364: /*@C
2365: MatMPISBAIJSetPreallocation - For good matrix assembly performance
2366: the user should preallocate the matrix storage by setting the parameters
2367: d_nz (or d_nnz) and o_nz (or o_nnz). By setting these parameters accurately,
2368: performance can be increased by more than a factor of 50.
2370: Collective on Mat
2372: Input Parameters:
2373: + B - the matrix
2374: . bs - size of block, the blocks are ALWAYS square. One can use MatSetBlockSizes() to set a different row and column blocksize but the row
2375: blocksize always defines the size of the blocks. The column blocksize sets the blocksize of the vectors obtained with MatCreateVecs()
2376: . d_nz - number of block nonzeros per block row in diagonal portion of local
2377: submatrix (same for all local rows)
2378: . d_nnz - array containing the number of block nonzeros in the various block rows
2379: in the upper triangular and diagonal part of the in diagonal portion of the local
2380: (possibly different for each block row) or NULL. If you plan to factor the matrix you must leave room
2381: for the diagonal entry and set a value even if it is zero.
2382: . o_nz - number of block nonzeros per block row in the off-diagonal portion of local
2383: submatrix (same for all local rows).
2384: - o_nnz - array containing the number of nonzeros in the various block rows of the
2385: off-diagonal portion of the local submatrix that is right of the diagonal
2386: (possibly different for each block row) or NULL.
2389: Options Database Keys:
2390: . -mat_no_unroll - uses code that does not unroll the loops in the
2391: block calculations (much slower)
2392: . -mat_block_size - size of the blocks to use
2394: Notes:
2396: If PETSC_DECIDE or PETSC_DETERMINE is used for a particular argument on one processor
2397: than it must be used on all processors that share the object for that argument.
2399: If the *_nnz parameter is given then the *_nz parameter is ignored
2401: Storage Information:
2402: For a square global matrix we define each processor's diagonal portion
2403: to be its local rows and the corresponding columns (a square submatrix);
2404: each processor's off-diagonal portion encompasses the remainder of the
2405: local matrix (a rectangular submatrix).
2407: The user can specify preallocated storage for the diagonal part of
2408: the local submatrix with either d_nz or d_nnz (not both). Set
2409: d_nz=PETSC_DEFAULT and d_nnz=NULL for PETSc to control dynamic
2410: memory allocation. Likewise, specify preallocated storage for the
2411: off-diagonal part of the local submatrix with o_nz or o_nnz (not both).
2413: You can call MatGetInfo() to get information on how effective the preallocation was;
2414: for example the fields mallocs,nz_allocated,nz_used,nz_unneeded;
2415: You can also run with the option -info and look for messages with the string
2416: malloc in them to see if additional memory allocation was needed.
2418: Consider a processor that owns rows 3, 4 and 5 of a parallel matrix. In
2419: the figure below we depict these three local rows and all columns (0-11).
2421: .vb
2422: 0 1 2 3 4 5 6 7 8 9 10 11
2423: --------------------------
2424: row 3 |. . . d d d o o o o o o
2425: row 4 |. . . d d d o o o o o o
2426: row 5 |. . . d d d o o o o o o
2427: --------------------------
2428: .ve
2430: Thus, any entries in the d locations are stored in the d (diagonal)
2431: submatrix, and any entries in the o locations are stored in the
2432: o (off-diagonal) submatrix. Note that the d matrix is stored in
2433: MatSeqSBAIJ format and the o submatrix in MATSEQBAIJ format.
2435: Now d_nz should indicate the number of block nonzeros per row in the upper triangular
2436: plus the diagonal part of the d matrix,
2437: and o_nz should indicate the number of block nonzeros per row in the o matrix
2439: In general, for PDE problems in which most nonzeros are near the diagonal,
2440: one expects d_nz >> o_nz. For large problems you MUST preallocate memory
2441: or you will get TERRIBLE performance; see the users' manual chapter on
2442: matrices.
2444: Level: intermediate
2446: .keywords: matrix, block, aij, compressed row, sparse, parallel
2448: .seealso: MatCreate(), MatCreateSeqSBAIJ(), MatSetValues(), MatCreateBAIJ(), PetscSplitOwnership()
2449: @*/
2450: PetscErrorCode MatMPISBAIJSetPreallocation(Mat B,PetscInt bs,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
2451: {
2458: PetscTryMethod(B,"MatMPISBAIJSetPreallocation_C",(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,bs,d_nz,d_nnz,o_nz,o_nnz));
2459: return(0);
2460: }
2464: /*@C
2465: MatCreateSBAIJ - Creates a sparse parallel matrix in symmetric block AIJ format
2466: (block compressed row). For good matrix assembly performance
2467: the user should preallocate the matrix storage by setting the parameters
2468: d_nz (or d_nnz) and o_nz (or o_nnz). By setting these parameters accurately,
2469: performance can be increased by more than a factor of 50.
2471: Collective on MPI_Comm
2473: Input Parameters:
2474: + comm - MPI communicator
2475: . bs - size of block, the blocks are ALWAYS square. One can use MatSetBlockSizes() to set a different row and column blocksize but the row
2476: blocksize always defines the size of the blocks. The column blocksize sets the blocksize of the vectors obtained with MatCreateVecs()
2477: . m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
2478: This value should be the same as the local size used in creating the
2479: y vector for the matrix-vector product y = Ax.
2480: . n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
2481: This value should be the same as the local size used in creating the
2482: x vector for the matrix-vector product y = Ax.
2483: . M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
2484: . N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
2485: . d_nz - number of block nonzeros per block row in diagonal portion of local
2486: submatrix (same for all local rows)
2487: . d_nnz - array containing the number of block nonzeros in the various block rows
2488: in the upper triangular portion of the in diagonal portion of the local
2489: (possibly different for each block block row) or NULL.
2490: If you plan to factor the matrix you must leave room for the diagonal entry and
2491: set its value even if it is zero.
2492: . o_nz - number of block nonzeros per block row in the off-diagonal portion of local
2493: submatrix (same for all local rows).
2494: - o_nnz - array containing the number of nonzeros in the various block rows of the
2495: off-diagonal portion of the local submatrix (possibly different for
2496: each block row) or NULL.
2498: Output Parameter:
2499: . A - the matrix
2501: Options Database Keys:
2502: . -mat_no_unroll - uses code that does not unroll the loops in the
2503: block calculations (much slower)
2504: . -mat_block_size - size of the blocks to use
2505: . -mat_mpi - use the parallel matrix data structures even on one processor
2506: (defaults to using SeqBAIJ format on one processor)
2508: It is recommended that one use the MatCreate(), MatSetType() and/or MatSetFromOptions(),
2509: MatXXXXSetPreallocation() paradgm instead of this routine directly.
2510: [MatXXXXSetPreallocation() is, for example, MatSeqAIJSetPreallocation]
2512: Notes:
2513: The number of rows and columns must be divisible by blocksize.
2514: This matrix type does not support complex Hermitian operation.
2516: The user MUST specify either the local or global matrix dimensions
2517: (possibly both).
2519: If PETSC_DECIDE or PETSC_DETERMINE is used for a particular argument on one processor
2520: than it must be used on all processors that share the object for that argument.
2522: If the *_nnz parameter is given then the *_nz parameter is ignored
2524: Storage Information:
2525: For a square global matrix we define each processor's diagonal portion
2526: to be its local rows and the corresponding columns (a square submatrix);
2527: each processor's off-diagonal portion encompasses the remainder of the
2528: local matrix (a rectangular submatrix).
2530: The user can specify preallocated storage for the diagonal part of
2531: the local submatrix with either d_nz or d_nnz (not both). Set
2532: d_nz=PETSC_DEFAULT and d_nnz=NULL for PETSc to control dynamic
2533: memory allocation. Likewise, specify preallocated storage for the
2534: off-diagonal part of the local submatrix with o_nz or o_nnz (not both).
2536: Consider a processor that owns rows 3, 4 and 5 of a parallel matrix. In
2537: the figure below we depict these three local rows and all columns (0-11).
2539: .vb
2540: 0 1 2 3 4 5 6 7 8 9 10 11
2541: --------------------------
2542: row 3 |. . . d d d o o o o o o
2543: row 4 |. . . d d d o o o o o o
2544: row 5 |. . . d d d o o o o o o
2545: --------------------------
2546: .ve
2548: Thus, any entries in the d locations are stored in the d (diagonal)
2549: submatrix, and any entries in the o locations are stored in the
2550: o (off-diagonal) submatrix. Note that the d matrix is stored in
2551: MatSeqSBAIJ format and the o submatrix in MATSEQBAIJ format.
2553: Now d_nz should indicate the number of block nonzeros per row in the upper triangular
2554: plus the diagonal part of the d matrix,
2555: and o_nz should indicate the number of block nonzeros per row in the o matrix.
2556: In general, for PDE problems in which most nonzeros are near the diagonal,
2557: one expects d_nz >> o_nz. For large problems you MUST preallocate memory
2558: or you will get TERRIBLE performance; see the users' manual chapter on
2559: matrices.
2561: Level: intermediate
2563: .keywords: matrix, block, aij, compressed row, sparse, parallel
2565: .seealso: MatCreate(), MatCreateSeqSBAIJ(), MatSetValues(), MatCreateBAIJ()
2566: @*/
2568: 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)
2569: {
2571: PetscMPIInt size;
2574: MatCreate(comm,A);
2575: MatSetSizes(*A,m,n,M,N);
2576: MPI_Comm_size(comm,&size);
2577: if (size > 1) {
2578: MatSetType(*A,MATMPISBAIJ);
2579: MatMPISBAIJSetPreallocation(*A,bs,d_nz,d_nnz,o_nz,o_nnz);
2580: } else {
2581: MatSetType(*A,MATSEQSBAIJ);
2582: MatSeqSBAIJSetPreallocation(*A,bs,d_nz,d_nnz);
2583: }
2584: return(0);
2585: }
2590: static PetscErrorCode MatDuplicate_MPISBAIJ(Mat matin,MatDuplicateOption cpvalues,Mat *newmat)
2591: {
2592: Mat mat;
2593: Mat_MPISBAIJ *a,*oldmat = (Mat_MPISBAIJ*)matin->data;
2595: PetscInt len=0,nt,bs=matin->rmap->bs,mbs=oldmat->mbs;
2596: PetscScalar *array;
2599: *newmat = 0;
2601: MatCreate(PetscObjectComm((PetscObject)matin),&mat);
2602: MatSetSizes(mat,matin->rmap->n,matin->cmap->n,matin->rmap->N,matin->cmap->N);
2603: MatSetType(mat,((PetscObject)matin)->type_name);
2604: PetscMemcpy(mat->ops,matin->ops,sizeof(struct _MatOps));
2605: PetscLayoutReference(matin->rmap,&mat->rmap);
2606: PetscLayoutReference(matin->cmap,&mat->cmap);
2608: mat->factortype = matin->factortype;
2609: mat->preallocated = PETSC_TRUE;
2610: mat->assembled = PETSC_TRUE;
2611: mat->insertmode = NOT_SET_VALUES;
2613: a = (Mat_MPISBAIJ*)mat->data;
2614: a->bs2 = oldmat->bs2;
2615: a->mbs = oldmat->mbs;
2616: a->nbs = oldmat->nbs;
2617: a->Mbs = oldmat->Mbs;
2618: a->Nbs = oldmat->Nbs;
2621: a->size = oldmat->size;
2622: a->rank = oldmat->rank;
2623: a->donotstash = oldmat->donotstash;
2624: a->roworiented = oldmat->roworiented;
2625: a->rowindices = 0;
2626: a->rowvalues = 0;
2627: a->getrowactive = PETSC_FALSE;
2628: a->barray = 0;
2629: a->rstartbs = oldmat->rstartbs;
2630: a->rendbs = oldmat->rendbs;
2631: a->cstartbs = oldmat->cstartbs;
2632: a->cendbs = oldmat->cendbs;
2634: /* hash table stuff */
2635: a->ht = 0;
2636: a->hd = 0;
2637: a->ht_size = 0;
2638: a->ht_flag = oldmat->ht_flag;
2639: a->ht_fact = oldmat->ht_fact;
2640: a->ht_total_ct = 0;
2641: a->ht_insert_ct = 0;
2643: PetscMemcpy(a->rangebs,oldmat->rangebs,(a->size+2)*sizeof(PetscInt));
2644: if (oldmat->colmap) {
2645: #if defined(PETSC_USE_CTABLE)
2646: PetscTableCreateCopy(oldmat->colmap,&a->colmap);
2647: #else
2648: PetscMalloc1(a->Nbs,&a->colmap);
2649: PetscLogObjectMemory((PetscObject)mat,(a->Nbs)*sizeof(PetscInt));
2650: PetscMemcpy(a->colmap,oldmat->colmap,(a->Nbs)*sizeof(PetscInt));
2651: #endif
2652: } else a->colmap = 0;
2654: if (oldmat->garray && (len = ((Mat_SeqBAIJ*)(oldmat->B->data))->nbs)) {
2655: PetscMalloc1(len,&a->garray);
2656: PetscLogObjectMemory((PetscObject)mat,len*sizeof(PetscInt));
2657: PetscMemcpy(a->garray,oldmat->garray,len*sizeof(PetscInt));
2658: } else a->garray = 0;
2660: MatStashCreate_Private(PetscObjectComm((PetscObject)matin),matin->rmap->bs,&mat->bstash);
2661: VecDuplicate(oldmat->lvec,&a->lvec);
2662: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->lvec);
2663: VecScatterCopy(oldmat->Mvctx,&a->Mvctx);
2664: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->Mvctx);
2666: VecDuplicate(oldmat->slvec0,&a->slvec0);
2667: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec0);
2668: VecDuplicate(oldmat->slvec1,&a->slvec1);
2669: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec1);
2671: VecGetLocalSize(a->slvec1,&nt);
2672: VecGetArray(a->slvec1,&array);
2673: VecCreateSeqWithArray(PETSC_COMM_SELF,1,bs*mbs,array,&a->slvec1a);
2674: VecCreateSeqWithArray(PETSC_COMM_SELF,1,nt-bs*mbs,array+bs*mbs,&a->slvec1b);
2675: VecRestoreArray(a->slvec1,&array);
2676: VecGetArray(a->slvec0,&array);
2677: VecCreateSeqWithArray(PETSC_COMM_SELF,1,nt-bs*mbs,array+bs*mbs,&a->slvec0b);
2678: VecRestoreArray(a->slvec0,&array);
2679: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec0);
2680: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec1);
2681: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec0b);
2682: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec1a);
2683: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec1b);
2685: /* VecScatterCopy(oldmat->sMvctx,&a->sMvctx); - not written yet, replaced by the lazy trick: */
2686: PetscObjectReference((PetscObject)oldmat->sMvctx);
2687: a->sMvctx = oldmat->sMvctx;
2688: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->sMvctx);
2690: MatDuplicate(oldmat->A,cpvalues,&a->A);
2691: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->A);
2692: MatDuplicate(oldmat->B,cpvalues,&a->B);
2693: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->B);
2694: PetscFunctionListDuplicate(((PetscObject)matin)->qlist,&((PetscObject)mat)->qlist);
2695: *newmat = mat;
2696: return(0);
2697: }
2701: PetscErrorCode MatLoad_MPISBAIJ(Mat newmat,PetscViewer viewer)
2702: {
2704: PetscInt i,nz,j,rstart,rend;
2705: PetscScalar *vals,*buf;
2706: MPI_Comm comm;
2707: MPI_Status status;
2708: PetscMPIInt rank,size,tag = ((PetscObject)viewer)->tag,*sndcounts = 0,*browners,maxnz,*rowners,mmbs;
2709: PetscInt header[4],*rowlengths = 0,M,N,m,*cols,*locrowlens;
2710: PetscInt *procsnz = 0,jj,*mycols,*ibuf;
2711: PetscInt bs = newmat->rmap->bs,Mbs,mbs,extra_rows;
2712: PetscInt *dlens,*odlens,*mask,*masked1,*masked2,rowcount,odcount;
2713: PetscInt dcount,kmax,k,nzcount,tmp;
2714: int fd;
2717: /* force binary viewer to load .info file if it has not yet done so */
2718: PetscViewerSetUp(viewer);
2719: PetscObjectGetComm((PetscObject)viewer,&comm);
2720: PetscOptionsBegin(comm,NULL,"Options for loading MPISBAIJ matrix 2","Mat");
2721: PetscOptionsInt("-matload_block_size","Set the blocksize used to store the matrix","MatLoad",bs,&bs,NULL);
2722: PetscOptionsEnd();
2723: if (bs < 0) bs = 1;
2725: MPI_Comm_size(comm,&size);
2726: MPI_Comm_rank(comm,&rank);
2727: PetscViewerBinaryGetDescriptor(viewer,&fd);
2728: if (!rank) {
2729: PetscBinaryRead(fd,(char*)header,4,PETSC_INT);
2730: if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"not matrix object");
2731: if (header[3] < 0) SETERRQ(PetscObjectComm((PetscObject)newmat),PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as MPISBAIJ");
2732: }
2734: MPI_Bcast(header+1,3,MPIU_INT,0,comm);
2735: M = header[1];
2736: N = header[2];
2738: /* If global sizes are set, check if they are consistent with that given in the file */
2739: 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);
2740: 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);
2742: if (M != N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Can only do square matrices");
2744: /*
2745: This code adds extra rows to make sure the number of rows is
2746: divisible by the blocksize
2747: */
2748: Mbs = M/bs;
2749: extra_rows = bs - M + bs*(Mbs);
2750: if (extra_rows == bs) extra_rows = 0;
2751: else Mbs++;
2752: if (extra_rows &&!rank) {
2753: PetscInfo(viewer,"Padding loaded matrix to match blocksize\n");
2754: }
2756: /* determine ownership of all rows */
2757: if (newmat->rmap->n < 0) { /* PETSC_DECIDE */
2758: mbs = Mbs/size + ((Mbs % size) > rank);
2759: m = mbs*bs;
2760: } else { /* User Set */
2761: m = newmat->rmap->n;
2762: mbs = m/bs;
2763: }
2764: PetscMalloc2(size+1,&rowners,size+1,&browners);
2765: PetscMPIIntCast(mbs,&mmbs);
2766: MPI_Allgather(&mmbs,1,MPI_INT,rowners+1,1,MPI_INT,comm);
2767: rowners[0] = 0;
2768: for (i=2; i<=size; i++) rowners[i] += rowners[i-1];
2769: for (i=0; i<=size; i++) browners[i] = rowners[i]*bs;
2770: rstart = rowners[rank];
2771: rend = rowners[rank+1];
2773: /* distribute row lengths to all processors */
2774: PetscMalloc1((rend-rstart)*bs,&locrowlens);
2775: if (!rank) {
2776: PetscMalloc1(M+extra_rows,&rowlengths);
2777: PetscBinaryRead(fd,rowlengths,M,PETSC_INT);
2778: for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;
2779: PetscMalloc1(size,&sndcounts);
2780: for (i=0; i<size; i++) sndcounts[i] = browners[i+1] - browners[i];
2781: MPI_Scatterv(rowlengths,sndcounts,browners,MPIU_INT,locrowlens,(rend-rstart)*bs,MPIU_INT,0,comm);
2782: PetscFree(sndcounts);
2783: } else {
2784: MPI_Scatterv(0,0,0,MPIU_INT,locrowlens,(rend-rstart)*bs,MPIU_INT,0,comm);
2785: }
2787: if (!rank) { /* procs[0] */
2788: /* calculate the number of nonzeros on each processor */
2789: PetscMalloc1(size,&procsnz);
2790: PetscMemzero(procsnz,size*sizeof(PetscInt));
2791: for (i=0; i<size; i++) {
2792: for (j=rowners[i]*bs; j< rowners[i+1]*bs; j++) {
2793: procsnz[i] += rowlengths[j];
2794: }
2795: }
2796: PetscFree(rowlengths);
2798: /* determine max buffer needed and allocate it */
2799: maxnz = 0;
2800: for (i=0; i<size; i++) {
2801: maxnz = PetscMax(maxnz,procsnz[i]);
2802: }
2803: PetscMalloc1(maxnz,&cols);
2805: /* read in my part of the matrix column indices */
2806: nz = procsnz[0];
2807: PetscMalloc1(nz,&ibuf);
2808: mycols = ibuf;
2809: if (size == 1) nz -= extra_rows;
2810: PetscBinaryRead(fd,mycols,nz,PETSC_INT);
2811: if (size == 1) {
2812: for (i=0; i< extra_rows; i++) mycols[nz+i] = M+i;
2813: }
2815: /* read in every ones (except the last) and ship off */
2816: for (i=1; i<size-1; i++) {
2817: nz = procsnz[i];
2818: PetscBinaryRead(fd,cols,nz,PETSC_INT);
2819: MPI_Send(cols,nz,MPIU_INT,i,tag,comm);
2820: }
2821: /* read in the stuff for the last proc */
2822: if (size != 1) {
2823: nz = procsnz[size-1] - extra_rows; /* the extra rows are not on the disk */
2824: PetscBinaryRead(fd,cols,nz,PETSC_INT);
2825: for (i=0; i<extra_rows; i++) cols[nz+i] = M+i;
2826: MPI_Send(cols,nz+extra_rows,MPIU_INT,size-1,tag,comm);
2827: }
2828: PetscFree(cols);
2829: } else { /* procs[i], i>0 */
2830: /* determine buffer space needed for message */
2831: nz = 0;
2832: for (i=0; i<m; i++) nz += locrowlens[i];
2833: PetscMalloc1(nz,&ibuf);
2834: mycols = ibuf;
2835: /* receive message of column indices*/
2836: MPI_Recv(mycols,nz,MPIU_INT,0,tag,comm,&status);
2837: MPI_Get_count(&status,MPIU_INT,&maxnz);
2838: if (maxnz != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
2839: }
2841: /* loop over local rows, determining number of off diagonal entries */
2842: PetscMalloc2(rend-rstart,&dlens,rend-rstart,&odlens);
2843: PetscMalloc3(Mbs,&mask,Mbs,&masked1,Mbs,&masked2);
2844: PetscMemzero(mask,Mbs*sizeof(PetscInt));
2845: PetscMemzero(masked1,Mbs*sizeof(PetscInt));
2846: PetscMemzero(masked2,Mbs*sizeof(PetscInt));
2847: rowcount = 0;
2848: nzcount = 0;
2849: for (i=0; i<mbs; i++) {
2850: dcount = 0;
2851: odcount = 0;
2852: for (j=0; j<bs; j++) {
2853: kmax = locrowlens[rowcount];
2854: for (k=0; k<kmax; k++) {
2855: tmp = mycols[nzcount++]/bs; /* block col. index */
2856: if (!mask[tmp]) {
2857: mask[tmp] = 1;
2858: if (tmp < rstart || tmp >= rend) masked2[odcount++] = tmp; /* entry in off-diag portion */
2859: else masked1[dcount++] = tmp; /* entry in diag portion */
2860: }
2861: }
2862: rowcount++;
2863: }
2865: dlens[i] = dcount; /* d_nzz[i] */
2866: odlens[i] = odcount; /* o_nzz[i] */
2868: /* zero out the mask elements we set */
2869: for (j=0; j<dcount; j++) mask[masked1[j]] = 0;
2870: for (j=0; j<odcount; j++) mask[masked2[j]] = 0;
2871: }
2872: MatSetSizes(newmat,m,m,M+extra_rows,N+extra_rows);
2873: MatMPISBAIJSetPreallocation(newmat,bs,0,dlens,0,odlens);
2874: MatSetOption(newmat,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);
2876: if (!rank) {
2877: PetscMalloc1(maxnz,&buf);
2878: /* read in my part of the matrix numerical values */
2879: nz = procsnz[0];
2880: vals = buf;
2881: mycols = ibuf;
2882: if (size == 1) nz -= extra_rows;
2883: PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
2884: if (size == 1) {
2885: for (i=0; i< extra_rows; i++) vals[nz+i] = 1.0;
2886: }
2888: /* insert into matrix */
2889: jj = rstart*bs;
2890: for (i=0; i<m; i++) {
2891: MatSetValues(newmat,1,&jj,locrowlens[i],mycols,vals,INSERT_VALUES);
2892: mycols += locrowlens[i];
2893: vals += locrowlens[i];
2894: jj++;
2895: }
2897: /* read in other processors (except the last one) and ship out */
2898: for (i=1; i<size-1; i++) {
2899: nz = procsnz[i];
2900: vals = buf;
2901: PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
2902: MPI_Send(vals,nz,MPIU_SCALAR,i,((PetscObject)newmat)->tag,comm);
2903: }
2904: /* the last proc */
2905: if (size != 1) {
2906: nz = procsnz[i] - extra_rows;
2907: vals = buf;
2908: PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
2909: for (i=0; i<extra_rows; i++) vals[nz+i] = 1.0;
2910: MPI_Send(vals,nz+extra_rows,MPIU_SCALAR,size-1,((PetscObject)newmat)->tag,comm);
2911: }
2912: PetscFree(procsnz);
2914: } else {
2915: /* receive numeric values */
2916: PetscMalloc1(nz,&buf);
2918: /* receive message of values*/
2919: vals = buf;
2920: mycols = ibuf;
2921: MPI_Recv(vals,nz,MPIU_SCALAR,0,((PetscObject)newmat)->tag,comm,&status);
2922: MPI_Get_count(&status,MPIU_SCALAR,&maxnz);
2923: if (maxnz != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
2925: /* insert into matrix */
2926: jj = rstart*bs;
2927: for (i=0; i<m; i++) {
2928: MatSetValues_MPISBAIJ(newmat,1,&jj,locrowlens[i],mycols,vals,INSERT_VALUES);
2929: mycols += locrowlens[i];
2930: vals += locrowlens[i];
2931: jj++;
2932: }
2933: }
2935: PetscFree(locrowlens);
2936: PetscFree(buf);
2937: PetscFree(ibuf);
2938: PetscFree2(rowners,browners);
2939: PetscFree2(dlens,odlens);
2940: PetscFree3(mask,masked1,masked2);
2941: MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);
2942: MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);
2943: return(0);
2944: }
2948: /*XXXXX@
2949: MatMPISBAIJSetHashTableFactor - Sets the factor required to compute the size of the HashTable.
2951: Input Parameters:
2952: . mat - the matrix
2953: . fact - factor
2955: Not Collective on Mat, each process can have a different hash factor
2957: Level: advanced
2959: Notes:
2960: This can also be set by the command line option: -mat_use_hash_table fact
2962: .keywords: matrix, hashtable, factor, HT
2964: .seealso: MatSetOption()
2965: @XXXXX*/
2970: PetscErrorCode MatGetRowMaxAbs_MPISBAIJ(Mat A,Vec v,PetscInt idx[])
2971: {
2972: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
2973: Mat_SeqBAIJ *b = (Mat_SeqBAIJ*)(a->B)->data;
2974: PetscReal atmp;
2975: PetscReal *work,*svalues,*rvalues;
2977: PetscInt i,bs,mbs,*bi,*bj,brow,j,ncols,krow,kcol,col,row,Mbs,bcol;
2978: PetscMPIInt rank,size;
2979: PetscInt *rowners_bs,dest,count,source;
2980: PetscScalar *va;
2981: MatScalar *ba;
2982: MPI_Status stat;
2985: if (idx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Send email to petsc-maint@mcs.anl.gov");
2986: MatGetRowMaxAbs(a->A,v,NULL);
2987: VecGetArray(v,&va);
2989: MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);
2990: MPI_Comm_rank(PetscObjectComm((PetscObject)A),&rank);
2992: bs = A->rmap->bs;
2993: mbs = a->mbs;
2994: Mbs = a->Mbs;
2995: ba = b->a;
2996: bi = b->i;
2997: bj = b->j;
2999: /* find ownerships */
3000: rowners_bs = A->rmap->range;
3002: /* each proc creates an array to be distributed */
3003: PetscMalloc1(bs*Mbs,&work);
3004: PetscMemzero(work,bs*Mbs*sizeof(PetscReal));
3006: /* row_max for B */
3007: if (rank != size-1) {
3008: for (i=0; i<mbs; i++) {
3009: ncols = bi[1] - bi[0]; bi++;
3010: brow = bs*i;
3011: for (j=0; j<ncols; j++) {
3012: bcol = bs*(*bj);
3013: for (kcol=0; kcol<bs; kcol++) {
3014: col = bcol + kcol; /* local col index */
3015: col += rowners_bs[rank+1]; /* global col index */
3016: for (krow=0; krow<bs; krow++) {
3017: atmp = PetscAbsScalar(*ba); ba++;
3018: row = brow + krow; /* local row index */
3019: if (PetscRealPart(va[row]) < atmp) va[row] = atmp;
3020: if (work[col] < atmp) work[col] = atmp;
3021: }
3022: }
3023: bj++;
3024: }
3025: }
3027: /* send values to its owners */
3028: for (dest=rank+1; dest<size; dest++) {
3029: svalues = work + rowners_bs[dest];
3030: count = rowners_bs[dest+1]-rowners_bs[dest];
3031: MPI_Send(svalues,count,MPIU_REAL,dest,rank,PetscObjectComm((PetscObject)A));
3032: }
3033: }
3035: /* receive values */
3036: if (rank) {
3037: rvalues = work;
3038: count = rowners_bs[rank+1]-rowners_bs[rank];
3039: for (source=0; source<rank; source++) {
3040: MPI_Recv(rvalues,count,MPIU_REAL,MPI_ANY_SOURCE,MPI_ANY_TAG,PetscObjectComm((PetscObject)A),&stat);
3041: /* process values */
3042: for (i=0; i<count; i++) {
3043: if (PetscRealPart(va[i]) < rvalues[i]) va[i] = rvalues[i];
3044: }
3045: }
3046: }
3048: VecRestoreArray(v,&va);
3049: PetscFree(work);
3050: return(0);
3051: }
3055: PetscErrorCode MatSOR_MPISBAIJ(Mat matin,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
3056: {
3057: Mat_MPISBAIJ *mat = (Mat_MPISBAIJ*)matin->data;
3058: PetscErrorCode ierr;
3059: PetscInt mbs=mat->mbs,bs=matin->rmap->bs;
3060: PetscScalar *x,*ptr,*from;
3061: Vec bb1;
3062: const PetscScalar *b;
3065: 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);
3066: if (bs > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented");
3068: if (flag == SOR_APPLY_UPPER) {
3069: (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx);
3070: return(0);
3071: }
3073: if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP) {
3074: if (flag & SOR_ZERO_INITIAL_GUESS) {
3075: (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,lits,xx);
3076: its--;
3077: }
3079: VecDuplicate(bb,&bb1);
3080: while (its--) {
3082: /* lower triangular part: slvec0b = - B^T*xx */
3083: (*mat->B->ops->multtranspose)(mat->B,xx,mat->slvec0b);
3085: /* copy xx into slvec0a */
3086: VecGetArray(mat->slvec0,&ptr);
3087: VecGetArray(xx,&x);
3088: PetscMemcpy(ptr,x,bs*mbs*sizeof(MatScalar));
3089: VecRestoreArray(mat->slvec0,&ptr);
3091: VecScale(mat->slvec0,-1.0);
3093: /* copy bb into slvec1a */
3094: VecGetArray(mat->slvec1,&ptr);
3095: VecGetArrayRead(bb,&b);
3096: PetscMemcpy(ptr,b,bs*mbs*sizeof(MatScalar));
3097: VecRestoreArray(mat->slvec1,&ptr);
3099: /* set slvec1b = 0 */
3100: VecSet(mat->slvec1b,0.0);
3102: VecScatterBegin(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD);
3103: VecRestoreArray(xx,&x);
3104: VecRestoreArrayRead(bb,&b);
3105: VecScatterEnd(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD);
3107: /* upper triangular part: bb1 = bb1 - B*x */
3108: (*mat->B->ops->multadd)(mat->B,mat->slvec1b,mat->slvec1a,bb1);
3110: /* local diagonal sweep */
3111: (*mat->A->ops->sor)(mat->A,bb1,omega,SOR_SYMMETRIC_SWEEP,fshift,lits,lits,xx);
3112: }
3113: VecDestroy(&bb1);
3114: } else if ((flag & SOR_LOCAL_FORWARD_SWEEP) && (its == 1) && (flag & SOR_ZERO_INITIAL_GUESS)) {
3115: (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx);
3116: } else if ((flag & SOR_LOCAL_BACKWARD_SWEEP) && (its == 1) && (flag & SOR_ZERO_INITIAL_GUESS)) {
3117: (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx);
3118: } else if (flag & SOR_EISENSTAT) {
3119: Vec xx1;
3120: PetscBool hasop;
3121: const PetscScalar *diag;
3122: PetscScalar *sl,scale = (omega - 2.0)/omega;
3123: PetscInt i,n;
3125: if (!mat->xx1) {
3126: VecDuplicate(bb,&mat->xx1);
3127: VecDuplicate(bb,&mat->bb1);
3128: }
3129: xx1 = mat->xx1;
3130: bb1 = mat->bb1;
3132: (*mat->A->ops->sor)(mat->A,bb,omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_BACKWARD_SWEEP),fshift,lits,1,xx);
3134: if (!mat->diag) {
3135: /* this is wrong for same matrix with new nonzero values */
3136: MatCreateVecs(matin,&mat->diag,NULL);
3137: MatGetDiagonal(matin,mat->diag);
3138: }
3139: MatHasOperation(matin,MATOP_MULT_DIAGONAL_BLOCK,&hasop);
3141: if (hasop) {
3142: MatMultDiagonalBlock(matin,xx,bb1);
3143: VecAYPX(mat->slvec1a,scale,bb);
3144: } else {
3145: /*
3146: These two lines are replaced by code that may be a bit faster for a good compiler
3147: VecPointwiseMult(mat->slvec1a,mat->diag,xx);
3148: VecAYPX(mat->slvec1a,scale,bb);
3149: */
3150: VecGetArray(mat->slvec1a,&sl);
3151: VecGetArrayRead(mat->diag,&diag);
3152: VecGetArrayRead(bb,&b);
3153: VecGetArray(xx,&x);
3154: VecGetLocalSize(xx,&n);
3155: if (omega == 1.0) {
3156: for (i=0; i<n; i++) sl[i] = b[i] - diag[i]*x[i];
3157: PetscLogFlops(2.0*n);
3158: } else {
3159: for (i=0; i<n; i++) sl[i] = b[i] + scale*diag[i]*x[i];
3160: PetscLogFlops(3.0*n);
3161: }
3162: VecRestoreArray(mat->slvec1a,&sl);
3163: VecRestoreArrayRead(mat->diag,&diag);
3164: VecRestoreArrayRead(bb,&b);
3165: VecRestoreArray(xx,&x);
3166: }
3168: /* multiply off-diagonal portion of matrix */
3169: VecSet(mat->slvec1b,0.0);
3170: (*mat->B->ops->multtranspose)(mat->B,xx,mat->slvec0b);
3171: VecGetArray(mat->slvec0,&from);
3172: VecGetArray(xx,&x);
3173: PetscMemcpy(from,x,bs*mbs*sizeof(MatScalar));
3174: VecRestoreArray(mat->slvec0,&from);
3175: VecRestoreArray(xx,&x);
3176: VecScatterBegin(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD);
3177: VecScatterEnd(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD);
3178: (*mat->B->ops->multadd)(mat->B,mat->slvec1b,mat->slvec1a,mat->slvec1a);
3180: /* local sweep */
3181: (*mat->A->ops->sor)(mat->A,mat->slvec1a,omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_FORWARD_SWEEP),fshift,lits,1,xx1);
3182: VecAXPY(xx,1.0,xx1);
3183: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatSORType is not supported for SBAIJ matrix format");
3184: return(0);
3185: }
3189: PetscErrorCode MatSOR_MPISBAIJ_2comm(Mat matin,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
3190: {
3191: Mat_MPISBAIJ *mat = (Mat_MPISBAIJ*)matin->data;
3193: Vec lvec1,bb1;
3196: 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);
3197: if (matin->rmap->bs > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented");
3199: if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP) {
3200: if (flag & SOR_ZERO_INITIAL_GUESS) {
3201: (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,lits,xx);
3202: its--;
3203: }
3205: VecDuplicate(mat->lvec,&lvec1);
3206: VecDuplicate(bb,&bb1);
3207: while (its--) {
3208: VecScatterBegin(mat->Mvctx,xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD);
3210: /* lower diagonal part: bb1 = bb - B^T*xx */
3211: (*mat->B->ops->multtranspose)(mat->B,xx,lvec1);
3212: VecScale(lvec1,-1.0);
3214: VecScatterEnd(mat->Mvctx,xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD);
3215: VecCopy(bb,bb1);
3216: VecScatterBegin(mat->Mvctx,lvec1,bb1,ADD_VALUES,SCATTER_REVERSE);
3218: /* upper diagonal part: bb1 = bb1 - B*x */
3219: VecScale(mat->lvec,-1.0);
3220: (*mat->B->ops->multadd)(mat->B,mat->lvec,bb1,bb1);
3222: VecScatterEnd(mat->Mvctx,lvec1,bb1,ADD_VALUES,SCATTER_REVERSE);
3224: /* diagonal sweep */
3225: (*mat->A->ops->sor)(mat->A,bb1,omega,SOR_SYMMETRIC_SWEEP,fshift,lits,lits,xx);
3226: }
3227: VecDestroy(&lvec1);
3228: VecDestroy(&bb1);
3229: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatSORType is not supported for SBAIJ matrix format");
3230: return(0);
3231: }
3235: /*@
3236: MatCreateMPISBAIJWithArrays - creates a MPI SBAIJ matrix using arrays that contain in standard
3237: CSR format the local rows.
3239: Collective on MPI_Comm
3241: Input Parameters:
3242: + comm - MPI communicator
3243: . bs - the block size, only a block size of 1 is supported
3244: . m - number of local rows (Cannot be PETSC_DECIDE)
3245: . n - This value should be the same as the local size used in creating the
3246: x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have
3247: calculated if N is given) For square matrices n is almost always m.
3248: . M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
3249: . N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
3250: . i - row indices
3251: . j - column indices
3252: - a - matrix values
3254: Output Parameter:
3255: . mat - the matrix
3257: Level: intermediate
3259: Notes:
3260: The i, j, and a arrays ARE copied by this routine into the internal format used by PETSc;
3261: thus you CANNOT change the matrix entries by changing the values of a[] after you have
3262: called this routine. Use MatCreateMPIAIJWithSplitArrays() to avoid needing to copy the arrays.
3264: The i and j indices are 0 based, and i indices are indices corresponding to the local j array.
3266: .keywords: matrix, aij, compressed row, sparse, parallel
3268: .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatMPIAIJSetPreallocation(), MatMPIAIJSetPreallocationCSR(),
3269: MPIAIJ, MatCreateAIJ(), MatCreateMPIAIJWithSplitArrays()
3270: @*/
3271: 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)
3272: {
3277: if (i[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"i (row indices) must start with 0");
3278: if (m < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"local number of rows (m) cannot be PETSC_DECIDE, or negative");
3279: MatCreate(comm,mat);
3280: MatSetSizes(*mat,m,n,M,N);
3281: MatSetType(*mat,MATMPISBAIJ);
3282: MatMPISBAIJSetPreallocationCSR(*mat,bs,i,j,a);
3283: return(0);
3284: }
3289: /*@C
3290: MatMPISBAIJSetPreallocationCSR - Allocates memory for a sparse parallel matrix in BAIJ format
3291: (the default parallel PETSc format).
3293: Collective on MPI_Comm
3295: Input Parameters:
3296: + B - the matrix
3297: . bs - the block size
3298: . i - the indices into j for the start of each local row (starts with zero)
3299: . j - the column indices for each local row (starts with zero) these must be sorted for each row
3300: - v - optional values in the matrix
3302: Level: developer
3304: .keywords: matrix, aij, compressed row, sparse, parallel
3306: .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatMPIBAIJSetPreallocation(), MatCreateAIJ(), MPIAIJ
3307: @*/
3308: PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat B,PetscInt bs,const PetscInt i[],const PetscInt j[], const PetscScalar v[])
3309: {
3313: PetscTryMethod(B,"MatMPISBAIJSetPreallocationCSR_C",(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]),(B,bs,i,j,v));
3314: return(0);
3315: }
3319: PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPISBAIJ(MPI_Comm comm,Mat inmat,PetscInt n,MatReuse scall,Mat *outmat)
3320: {
3322: PetscInt m,N,i,rstart,nnz,Ii,bs,cbs;
3323: PetscInt *indx;
3324: PetscScalar *values;
3327: MatGetSize(inmat,&m,&N);
3328: if (scall == MAT_INITIAL_MATRIX) { /* symbolic phase */
3329: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inmat->data;
3330: PetscInt *dnz,*onz,sum,bs,cbs,mbs,Nbs;
3331: PetscInt *bindx,rmax=a->rmax,j;
3332:
3333: MatGetBlockSizes(inmat,&bs,&cbs);
3334: mbs = m/bs; Nbs = N/cbs;
3335: if (n == PETSC_DECIDE) {
3336: PetscSplitOwnership(comm,&n,&Nbs);
3337: }
3338: /* Check sum(n) = Nbs */
3339: MPIU_Allreduce(&n,&sum,1,MPIU_INT,MPI_SUM,comm);
3340: if (sum != Nbs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Sum of local columns != global columns %d",Nbs);
3342: MPI_Scan(&mbs, &rstart,1,MPIU_INT,MPI_SUM,comm);
3343: rstart -= mbs;
3345: PetscMalloc1(rmax,&bindx);
3346: MatPreallocateInitialize(comm,mbs,n,dnz,onz);
3347: MatSetOption(inmat,MAT_GETROW_UPPERTRIANGULAR,PETSC_TRUE);
3348: for (i=0; i<mbs; i++) {
3349: MatGetRow_SeqSBAIJ(inmat,i*bs,&nnz,&indx,NULL); /* non-blocked nnz and indx */
3350: nnz = nnz/bs;
3351: for (j=0; j<nnz; j++) bindx[j] = indx[j*bs]/bs;
3352: MatPreallocateSet(i+rstart,nnz,bindx,dnz,onz);
3353: MatRestoreRow_SeqSBAIJ(inmat,i*bs,&nnz,&indx,NULL);
3354: }
3355: MatSetOption(inmat,MAT_GETROW_UPPERTRIANGULAR,PETSC_FALSE);
3356: PetscFree(bindx);
3358: MatCreate(comm,outmat);
3359: MatSetSizes(*outmat,m,n*bs,PETSC_DETERMINE,PETSC_DETERMINE);
3360: MatSetBlockSizes(*outmat,bs,cbs);
3361: MatSetType(*outmat,MATMPISBAIJ);
3362: MatMPISBAIJSetPreallocation(*outmat,bs,0,dnz,0,onz);
3363: MatPreallocateFinalize(dnz,onz);
3364: }
3365:
3366: /* numeric phase */
3367: MatGetBlockSizes(inmat,&bs,&cbs);
3368: MatGetOwnershipRange(*outmat,&rstart,NULL);
3370: MatSetOption(inmat,MAT_GETROW_UPPERTRIANGULAR,PETSC_TRUE);
3371: for (i=0; i<m; i++) {
3372: MatGetRow_SeqSBAIJ(inmat,i,&nnz,&indx,&values);
3373: Ii = i + rstart;
3374: MatSetValues(*outmat,1,&Ii,nnz,indx,values,INSERT_VALUES);
3375: MatRestoreRow_SeqSBAIJ(inmat,i,&nnz,&indx,&values);
3376: }
3377: MatSetOption(inmat,MAT_GETROW_UPPERTRIANGULAR,PETSC_FALSE);
3378: MatAssemblyBegin(*outmat,MAT_FINAL_ASSEMBLY);
3379: MatAssemblyEnd(*outmat,MAT_FINAL_ASSEMBLY);
3380: return(0);
3381: }