Actual source code: mpisbaij.c
2: #include <../src/mat/impls/baij/mpi/mpibaij.h>
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
10: #if defined(PETSC_HAVE_SCALAPACK)
11: PETSC_INTERN PetscErrorCode MatConvert_SBAIJ_ScaLAPACK(Mat,MatType,MatReuse,Mat*);
12: #endif
14: /* This could be moved to matimpl.h */
15: static PetscErrorCode MatPreallocateWithMats_Private(Mat B, PetscInt nm, Mat X[], PetscBool symm[], PetscBool fill)
16: {
17: Mat preallocator;
18: PetscInt r,rstart,rend;
19: PetscInt bs,i,m,n,M,N;
20: PetscBool cong = PETSC_TRUE;
26: for (i = 0; i < nm; i++) {
28: PetscLayoutCompare(B->rmap,X[i]->rmap,&cong);
29: if (!cong) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"Not for different layouts");
30: }
32: MatGetBlockSize(B,&bs);
33: MatGetSize(B,&M,&N);
34: MatGetLocalSize(B,&m,&n);
35: MatCreate(PetscObjectComm((PetscObject)B),&preallocator);
36: MatSetType(preallocator,MATPREALLOCATOR);
37: MatSetBlockSize(preallocator,bs);
38: MatSetSizes(preallocator,m,n,M,N);
39: MatSetUp(preallocator);
40: MatGetOwnershipRange(preallocator,&rstart,&rend);
41: for (r = rstart; r < rend; ++r) {
42: PetscInt ncols;
43: const PetscInt *row;
44: const PetscScalar *vals;
46: for (i = 0; i < nm; i++) {
47: MatGetRow(X[i],r,&ncols,&row,&vals);
48: MatSetValues(preallocator,1,&r,ncols,row,vals,INSERT_VALUES);
49: if (symm && symm[i]) {
50: MatSetValues(preallocator,ncols,row,1,&r,vals,INSERT_VALUES);
51: }
52: MatRestoreRow(X[i],r,&ncols,&row,&vals);
53: }
54: }
55: MatAssemblyBegin(preallocator,MAT_FINAL_ASSEMBLY);
56: MatAssemblyEnd(preallocator,MAT_FINAL_ASSEMBLY);
57: MatPreallocatorPreallocate(preallocator,fill,B);
58: MatDestroy(&preallocator);
59: return(0);
60: }
62: PETSC_INTERN PetscErrorCode MatConvert_MPISBAIJ_Basic(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
63: {
64: Mat B;
66: PetscInt r;
69: if (reuse != MAT_REUSE_MATRIX) {
70: PetscBool symm = PETSC_TRUE,isdense;
71: PetscInt bs;
73: MatCreate(PetscObjectComm((PetscObject)A),&B);
74: MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
75: MatSetType(B,newtype);
76: MatGetBlockSize(A,&bs);
77: MatSetBlockSize(B,bs);
78: PetscLayoutSetUp(B->rmap);
79: PetscLayoutSetUp(B->cmap);
80: PetscObjectTypeCompareAny((PetscObject)B,&isdense,MATSEQDENSE,MATMPIDENSE,MATSEQDENSECUDA,"");
81: if (!isdense) {
82: MatGetRowUpperTriangular(A);
83: MatPreallocateWithMats_Private(B,1,&A,&symm,PETSC_TRUE);
84: MatRestoreRowUpperTriangular(A);
85: } else {
86: MatSetUp(B);
87: }
88: } else {
89: B = *newmat;
90: MatZeroEntries(B);
91: }
93: MatGetRowUpperTriangular(A);
94: for (r = A->rmap->rstart; r < A->rmap->rend; r++) {
95: PetscInt ncols;
96: const PetscInt *row;
97: const PetscScalar *vals;
99: MatGetRow(A,r,&ncols,&row,&vals);
100: MatSetValues(B,1,&r,ncols,row,vals,INSERT_VALUES);
101: #if defined(PETSC_USE_COMPLEX)
102: if (A->hermitian) {
103: PetscInt i;
104: for (i = 0; i < ncols; i++) {
105: MatSetValue(B,row[i],r,PetscConj(vals[i]),INSERT_VALUES);
106: }
107: } else {
108: MatSetValues(B,ncols,row,1,&r,vals,INSERT_VALUES);
109: }
110: #else
111: MatSetValues(B,ncols,row,1,&r,vals,INSERT_VALUES);
112: #endif
113: MatRestoreRow(A,r,&ncols,&row,&vals);
114: }
115: MatRestoreRowUpperTriangular(A);
116: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
117: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
119: if (reuse == MAT_INPLACE_MATRIX) {
120: MatHeaderReplace(A,&B);
121: } else {
122: *newmat = B;
123: }
124: return(0);
125: }
127: PetscErrorCode MatStoreValues_MPISBAIJ(Mat mat)
128: {
129: Mat_MPISBAIJ *aij = (Mat_MPISBAIJ*)mat->data;
133: MatStoreValues(aij->A);
134: MatStoreValues(aij->B);
135: return(0);
136: }
138: PetscErrorCode MatRetrieveValues_MPISBAIJ(Mat mat)
139: {
140: Mat_MPISBAIJ *aij = (Mat_MPISBAIJ*)mat->data;
144: MatRetrieveValues(aij->A);
145: MatRetrieveValues(aij->B);
146: return(0);
147: }
149: #define MatSetValues_SeqSBAIJ_A_Private(row,col,value,addv,orow,ocol) \
150: { \
151: brow = row/bs; \
152: rp = aj + ai[brow]; ap = aa + bs2*ai[brow]; \
153: rmax = aimax[brow]; nrow = ailen[brow]; \
154: bcol = col/bs; \
155: ridx = row % bs; cidx = col % bs; \
156: low = 0; high = nrow; \
157: while (high-low > 3) { \
158: t = (low+high)/2; \
159: if (rp[t] > bcol) high = t; \
160: else low = t; \
161: } \
162: for (_i=low; _i<high; _i++) { \
163: if (rp[_i] > bcol) break; \
164: if (rp[_i] == bcol) { \
165: bap = ap + bs2*_i + bs*cidx + ridx; \
166: if (addv == ADD_VALUES) *bap += value; \
167: else *bap = value; \
168: goto a_noinsert; \
169: } \
170: } \
171: if (a->nonew == 1) goto a_noinsert; \
172: 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); \
173: MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,aimax,a->nonew,MatScalar); \
174: N = nrow++ - 1; \
175: /* shift up all the later entries in this row */ \
176: PetscArraymove(rp+_i+1,rp+_i,N-_i+1); \
177: PetscArraymove(ap+bs2*(_i+1),ap+bs2*_i,bs2*(N-_i+1)); \
178: PetscArrayzero(ap+bs2*_i,bs2); \
179: rp[_i] = bcol; \
180: ap[bs2*_i + bs*cidx + ridx] = value; \
181: A->nonzerostate++;\
182: a_noinsert:; \
183: ailen[brow] = nrow; \
184: }
186: #define MatSetValues_SeqSBAIJ_B_Private(row,col,value,addv,orow,ocol) \
187: { \
188: brow = row/bs; \
189: rp = bj + bi[brow]; ap = ba + bs2*bi[brow]; \
190: rmax = bimax[brow]; nrow = bilen[brow]; \
191: bcol = col/bs; \
192: ridx = row % bs; cidx = col % bs; \
193: low = 0; high = nrow; \
194: while (high-low > 3) { \
195: t = (low+high)/2; \
196: if (rp[t] > bcol) high = t; \
197: else low = t; \
198: } \
199: for (_i=low; _i<high; _i++) { \
200: if (rp[_i] > bcol) break; \
201: if (rp[_i] == bcol) { \
202: bap = ap + bs2*_i + bs*cidx + ridx; \
203: if (addv == ADD_VALUES) *bap += value; \
204: else *bap = value; \
205: goto b_noinsert; \
206: } \
207: } \
208: if (b->nonew == 1) goto b_noinsert; \
209: 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); \
210: MatSeqXAIJReallocateAIJ(B,b->mbs,bs2,nrow,brow,bcol,rmax,ba,bi,bj,rp,ap,bimax,b->nonew,MatScalar); \
211: N = nrow++ - 1; \
212: /* shift up all the later entries in this row */ \
213: PetscArraymove(rp+_i+1,rp+_i,N-_i+1); \
214: PetscArraymove(ap+bs2*(_i+1),ap+bs2*_i,bs2*(N-_i+1)); \
215: PetscArrayzero(ap+bs2*_i,bs2); \
216: rp[_i] = bcol; \
217: ap[bs2*_i + bs*cidx + ridx] = value; \
218: B->nonzerostate++;\
219: b_noinsert:; \
220: bilen[brow] = nrow; \
221: }
223: /* Only add/insert a(i,j) with i<=j (blocks).
224: Any a(i,j) with i>j input by user is ingored or generates an error
225: */
226: PetscErrorCode MatSetValues_MPISBAIJ(Mat mat,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode addv)
227: {
228: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
229: MatScalar value;
230: PetscBool roworiented = baij->roworiented;
232: PetscInt i,j,row,col;
233: PetscInt rstart_orig=mat->rmap->rstart;
234: PetscInt rend_orig =mat->rmap->rend,cstart_orig=mat->cmap->rstart;
235: PetscInt cend_orig =mat->cmap->rend,bs=mat->rmap->bs;
237: /* Some Variables required in the macro */
238: Mat A = baij->A;
239: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)(A)->data;
240: PetscInt *aimax=a->imax,*ai=a->i,*ailen=a->ilen,*aj=a->j;
241: MatScalar *aa =a->a;
243: Mat B = baij->B;
244: Mat_SeqBAIJ *b = (Mat_SeqBAIJ*)(B)->data;
245: PetscInt *bimax=b->imax,*bi=b->i,*bilen=b->ilen,*bj=b->j;
246: MatScalar *ba =b->a;
248: PetscInt *rp,ii,nrow,_i,rmax,N,brow,bcol;
249: PetscInt low,high,t,ridx,cidx,bs2=a->bs2;
250: MatScalar *ap,*bap;
252: /* for stash */
253: PetscInt n_loc, *in_loc = NULL;
254: MatScalar *v_loc = NULL;
257: if (!baij->donotstash) {
258: if (n > baij->n_loc) {
259: PetscFree(baij->in_loc);
260: PetscFree(baij->v_loc);
261: PetscMalloc1(n,&baij->in_loc);
262: PetscMalloc1(n,&baij->v_loc);
264: baij->n_loc = n;
265: }
266: in_loc = baij->in_loc;
267: v_loc = baij->v_loc;
268: }
270: for (i=0; i<m; i++) {
271: if (im[i] < 0) continue;
272: if (PetscUnlikely(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);
273: if (im[i] >= rstart_orig && im[i] < rend_orig) { /* this processor entry */
274: row = im[i] - rstart_orig; /* local row index */
275: for (j=0; j<n; j++) {
276: if (im[i]/bs > in[j]/bs) {
277: if (a->ignore_ltriangular) {
278: continue; /* ignore lower triangular blocks */
279: } 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)");
280: }
281: if (in[j] >= cstart_orig && in[j] < cend_orig) { /* diag entry (A) */
282: col = in[j] - cstart_orig; /* local col index */
283: brow = row/bs; bcol = col/bs;
284: if (brow > bcol) continue; /* ignore lower triangular blocks of A */
285: if (roworiented) value = v[i*n+j];
286: else value = v[i+j*m];
287: MatSetValues_SeqSBAIJ_A_Private(row,col,value,addv,im[i],in[j]);
288: /* MatSetValues_SeqBAIJ(baij->A,1,&row,1,&col,&value,addv); */
289: } else if (in[j] < 0) continue;
290: 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);
291: else { /* off-diag entry (B) */
292: if (mat->was_assembled) {
293: if (!baij->colmap) {
294: MatCreateColmap_MPIBAIJ_Private(mat);
295: }
296: #if defined(PETSC_USE_CTABLE)
297: PetscTableFind(baij->colmap,in[j]/bs + 1,&col);
298: col = col - 1;
299: #else
300: col = baij->colmap[in[j]/bs] - 1;
301: #endif
302: if (col < 0 && !((Mat_SeqSBAIJ*)(baij->A->data))->nonew) {
303: MatDisAssemble_MPISBAIJ(mat);
304: col = in[j];
305: /* Reinitialize the variables required by MatSetValues_SeqBAIJ_B_Private() */
306: B = baij->B;
307: b = (Mat_SeqBAIJ*)(B)->data;
308: bimax= b->imax;bi=b->i;bilen=b->ilen;bj=b->j;
309: ba = b->a;
310: } else col += in[j]%bs;
311: } else col = in[j];
312: if (roworiented) value = v[i*n+j];
313: else value = v[i+j*m];
314: MatSetValues_SeqSBAIJ_B_Private(row,col,value,addv,im[i],in[j]);
315: /* MatSetValues_SeqBAIJ(baij->B,1,&row,1,&col,&value,addv); */
316: }
317: }
318: } else { /* off processor entry */
319: 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]);
320: if (!baij->donotstash) {
321: mat->assembled = PETSC_FALSE;
322: n_loc = 0;
323: for (j=0; j<n; j++) {
324: if (im[i]/bs > in[j]/bs) continue; /* ignore lower triangular blocks */
325: in_loc[n_loc] = in[j];
326: if (roworiented) {
327: v_loc[n_loc] = v[i*n+j];
328: } else {
329: v_loc[n_loc] = v[j*m+i];
330: }
331: n_loc++;
332: }
333: MatStashValuesRow_Private(&mat->stash,im[i],n_loc,in_loc,v_loc,PETSC_FALSE);
334: }
335: }
336: }
337: return(0);
338: }
340: PETSC_STATIC_INLINE PetscErrorCode MatSetValuesBlocked_SeqSBAIJ_Inlined(Mat A,PetscInt row,PetscInt col,const PetscScalar v[],InsertMode is,PetscInt orow,PetscInt ocol)
341: {
342: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
343: PetscErrorCode ierr;
344: PetscInt *rp,low,high,t,ii,jj,nrow,i,rmax,N;
345: PetscInt *imax =a->imax,*ai=a->i,*ailen=a->ilen;
346: PetscInt *aj =a->j,nonew=a->nonew,bs2=a->bs2,bs=A->rmap->bs;
347: PetscBool roworiented=a->roworiented;
348: const PetscScalar *value = v;
349: MatScalar *ap,*aa = a->a,*bap;
352: if (col < row) {
353: if (a->ignore_ltriangular) return(0); /* ignore lower triangular block */
354: 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)");
355: }
356: rp = aj + ai[row];
357: ap = aa + bs2*ai[row];
358: rmax = imax[row];
359: nrow = ailen[row];
360: value = v;
361: low = 0;
362: high = nrow;
364: while (high-low > 7) {
365: t = (low+high)/2;
366: if (rp[t] > col) high = t;
367: else low = t;
368: }
369: for (i=low; i<high; i++) {
370: if (rp[i] > col) break;
371: if (rp[i] == col) {
372: bap = ap + bs2*i;
373: if (roworiented) {
374: if (is == ADD_VALUES) {
375: for (ii=0; ii<bs; ii++) {
376: for (jj=ii; jj<bs2; jj+=bs) {
377: bap[jj] += *value++;
378: }
379: }
380: } else {
381: for (ii=0; ii<bs; ii++) {
382: for (jj=ii; jj<bs2; jj+=bs) {
383: bap[jj] = *value++;
384: }
385: }
386: }
387: } else {
388: if (is == ADD_VALUES) {
389: for (ii=0; ii<bs; ii++) {
390: for (jj=0; jj<bs; jj++) {
391: *bap++ += *value++;
392: }
393: }
394: } else {
395: for (ii=0; ii<bs; ii++) {
396: for (jj=0; jj<bs; jj++) {
397: *bap++ = *value++;
398: }
399: }
400: }
401: }
402: goto noinsert2;
403: }
404: }
405: if (nonew == 1) goto noinsert2;
406: 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);
407: MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar);
408: N = nrow++ - 1; high++;
409: /* shift up all the later entries in this row */
410: PetscArraymove(rp+i+1,rp+i,N-i+1);
411: PetscArraymove(ap+bs2*(i+1),ap+bs2*i,bs2*(N-i+1));
412: rp[i] = col;
413: bap = ap + bs2*i;
414: if (roworiented) {
415: for (ii=0; ii<bs; ii++) {
416: for (jj=ii; jj<bs2; jj+=bs) {
417: bap[jj] = *value++;
418: }
419: }
420: } else {
421: for (ii=0; ii<bs; ii++) {
422: for (jj=0; jj<bs; jj++) {
423: *bap++ = *value++;
424: }
425: }
426: }
427: noinsert2:;
428: ailen[row] = nrow;
429: return(0);
430: }
432: /*
433: This routine is exactly duplicated in mpibaij.c
434: */
435: PETSC_STATIC_INLINE PetscErrorCode MatSetValuesBlocked_SeqBAIJ_Inlined(Mat A,PetscInt row,PetscInt col,const PetscScalar v[],InsertMode is,PetscInt orow,PetscInt ocol)
436: {
437: Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
438: PetscInt *rp,low,high,t,ii,jj,nrow,i,rmax,N;
439: PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen;
440: PetscErrorCode ierr;
441: PetscInt *aj =a->j,nonew=a->nonew,bs2=a->bs2,bs=A->rmap->bs;
442: PetscBool roworiented=a->roworiented;
443: const PetscScalar *value = v;
444: MatScalar *ap,*aa = a->a,*bap;
447: rp = aj + ai[row];
448: ap = aa + bs2*ai[row];
449: rmax = imax[row];
450: nrow = ailen[row];
451: low = 0;
452: high = nrow;
453: value = v;
454: while (high-low > 7) {
455: t = (low+high)/2;
456: if (rp[t] > col) high = t;
457: else low = t;
458: }
459: for (i=low; i<high; i++) {
460: if (rp[i] > col) break;
461: if (rp[i] == col) {
462: bap = ap + bs2*i;
463: if (roworiented) {
464: if (is == ADD_VALUES) {
465: for (ii=0; ii<bs; ii++) {
466: for (jj=ii; jj<bs2; jj+=bs) {
467: bap[jj] += *value++;
468: }
469: }
470: } else {
471: for (ii=0; ii<bs; ii++) {
472: for (jj=ii; jj<bs2; jj+=bs) {
473: bap[jj] = *value++;
474: }
475: }
476: }
477: } else {
478: if (is == ADD_VALUES) {
479: for (ii=0; ii<bs; ii++,value+=bs) {
480: for (jj=0; jj<bs; jj++) {
481: bap[jj] += value[jj];
482: }
483: bap += bs;
484: }
485: } else {
486: for (ii=0; ii<bs; ii++,value+=bs) {
487: for (jj=0; jj<bs; jj++) {
488: bap[jj] = value[jj];
489: }
490: bap += bs;
491: }
492: }
493: }
494: goto noinsert2;
495: }
496: }
497: if (nonew == 1) goto noinsert2;
498: 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);
499: MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar);
500: N = nrow++ - 1; high++;
501: /* shift up all the later entries in this row */
502: PetscArraymove(rp+i+1,rp+i,N-i+1);
503: PetscArraymove(ap+bs2*(i+1),ap+bs2*i,bs2*(N-i+1));
504: rp[i] = col;
505: bap = ap + bs2*i;
506: if (roworiented) {
507: for (ii=0; ii<bs; ii++) {
508: for (jj=ii; jj<bs2; jj+=bs) {
509: bap[jj] = *value++;
510: }
511: }
512: } else {
513: for (ii=0; ii<bs; ii++) {
514: for (jj=0; jj<bs; jj++) {
515: *bap++ = *value++;
516: }
517: }
518: }
519: noinsert2:;
520: ailen[row] = nrow;
521: return(0);
522: }
524: /*
525: This routine could be optimized by removing the need for the block copy below and passing stride information
526: to the above inline routines; similarly in MatSetValuesBlocked_MPIBAIJ()
527: */
528: PetscErrorCode MatSetValuesBlocked_MPISBAIJ(Mat mat,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const MatScalar v[],InsertMode addv)
529: {
530: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
531: const MatScalar *value;
532: MatScalar *barray =baij->barray;
533: PetscBool roworiented = baij->roworiented,ignore_ltriangular = ((Mat_SeqSBAIJ*)baij->A->data)->ignore_ltriangular;
534: PetscErrorCode ierr;
535: PetscInt i,j,ii,jj,row,col,rstart=baij->rstartbs;
536: PetscInt rend=baij->rendbs,cstart=baij->cstartbs,stepval;
537: PetscInt cend=baij->cendbs,bs=mat->rmap->bs,bs2=baij->bs2;
540: if (!barray) {
541: PetscMalloc1(bs2,&barray);
542: baij->barray = barray;
543: }
545: if (roworiented) {
546: stepval = (n-1)*bs;
547: } else {
548: stepval = (m-1)*bs;
549: }
550: for (i=0; i<m; i++) {
551: if (im[i] < 0) continue;
552: if (PetscUnlikelyDebug(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);
553: if (im[i] >= rstart && im[i] < rend) {
554: row = im[i] - rstart;
555: for (j=0; j<n; j++) {
556: if (im[i] > in[j]) {
557: if (ignore_ltriangular) continue; /* ignore lower triangular blocks */
558: 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)");
559: }
560: /* If NumCol = 1 then a copy is not required */
561: if ((roworiented) && (n == 1)) {
562: barray = (MatScalar*) v + i*bs2;
563: } else if ((!roworiented) && (m == 1)) {
564: barray = (MatScalar*) v + j*bs2;
565: } else { /* Here a copy is required */
566: if (roworiented) {
567: value = v + i*(stepval+bs)*bs + j*bs;
568: } else {
569: value = v + j*(stepval+bs)*bs + i*bs;
570: }
571: for (ii=0; ii<bs; ii++,value+=stepval) {
572: for (jj=0; jj<bs; jj++) {
573: *barray++ = *value++;
574: }
575: }
576: barray -=bs2;
577: }
579: if (in[j] >= cstart && in[j] < cend) {
580: col = in[j] - cstart;
581: MatSetValuesBlocked_SeqSBAIJ_Inlined(baij->A,row,col,barray,addv,im[i],in[j]);
582: } else if (in[j] < 0) continue;
583: else if (PetscUnlikelyDebug(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);
584: else {
585: if (mat->was_assembled) {
586: if (!baij->colmap) {
587: MatCreateColmap_MPIBAIJ_Private(mat);
588: }
590: #if defined(PETSC_USE_DEBUG)
591: #if defined(PETSC_USE_CTABLE)
592: { PetscInt data;
593: PetscTableFind(baij->colmap,in[j]+1,&data);
594: if ((data - 1) % bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect colmap");
595: }
596: #else
597: if ((baij->colmap[in[j]] - 1) % bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect colmap");
598: #endif
599: #endif
600: #if defined(PETSC_USE_CTABLE)
601: PetscTableFind(baij->colmap,in[j]+1,&col);
602: col = (col - 1)/bs;
603: #else
604: col = (baij->colmap[in[j]] - 1)/bs;
605: #endif
606: if (col < 0 && !((Mat_SeqBAIJ*)(baij->A->data))->nonew) {
607: MatDisAssemble_MPISBAIJ(mat);
608: col = in[j];
609: }
610: } else col = in[j];
611: MatSetValuesBlocked_SeqBAIJ_Inlined(baij->B,row,col,barray,addv,im[i],in[j]);
612: }
613: }
614: } else {
615: 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]);
616: if (!baij->donotstash) {
617: if (roworiented) {
618: MatStashValuesRowBlocked_Private(&mat->bstash,im[i],n,in,v,m,n,i);
619: } else {
620: MatStashValuesColBlocked_Private(&mat->bstash,im[i],n,in,v,m,n,i);
621: }
622: }
623: }
624: }
625: return(0);
626: }
628: PetscErrorCode MatGetValues_MPISBAIJ(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[])
629: {
630: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
632: PetscInt bs = mat->rmap->bs,i,j,bsrstart = mat->rmap->rstart,bsrend = mat->rmap->rend;
633: PetscInt bscstart = mat->cmap->rstart,bscend = mat->cmap->rend,row,col,data;
636: for (i=0; i<m; i++) {
637: if (idxm[i] < 0) continue; /* SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",idxm[i]); */
638: 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);
639: if (idxm[i] >= bsrstart && idxm[i] < bsrend) {
640: row = idxm[i] - bsrstart;
641: for (j=0; j<n; j++) {
642: if (idxn[j] < 0) continue; /* SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column %D",idxn[j]); */
643: 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);
644: if (idxn[j] >= bscstart && idxn[j] < bscend) {
645: col = idxn[j] - bscstart;
646: MatGetValues_SeqSBAIJ(baij->A,1,&row,1,&col,v+i*n+j);
647: } else {
648: if (!baij->colmap) {
649: MatCreateColmap_MPIBAIJ_Private(mat);
650: }
651: #if defined(PETSC_USE_CTABLE)
652: PetscTableFind(baij->colmap,idxn[j]/bs+1,&data);
653: data--;
654: #else
655: data = baij->colmap[idxn[j]/bs]-1;
656: #endif
657: if ((data < 0) || (baij->garray[data/bs] != idxn[j]/bs)) *(v+i*n+j) = 0.0;
658: else {
659: col = data + idxn[j]%bs;
660: MatGetValues_SeqBAIJ(baij->B,1,&row,1,&col,v+i*n+j);
661: }
662: }
663: }
664: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only local values currently supported");
665: }
666: return(0);
667: }
669: PetscErrorCode MatNorm_MPISBAIJ(Mat mat,NormType type,PetscReal *norm)
670: {
671: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
673: PetscReal sum[2],*lnorm2;
676: if (baij->size == 1) {
677: MatNorm(baij->A,type,norm);
678: } else {
679: if (type == NORM_FROBENIUS) {
680: PetscMalloc1(2,&lnorm2);
681: MatNorm(baij->A,type,lnorm2);
682: *lnorm2 = (*lnorm2)*(*lnorm2); lnorm2++; /* squar power of norm(A) */
683: MatNorm(baij->B,type,lnorm2);
684: *lnorm2 = (*lnorm2)*(*lnorm2); lnorm2--; /* squar power of norm(B) */
685: MPIU_Allreduce(lnorm2,sum,2,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)mat));
686: *norm = PetscSqrtReal(sum[0] + 2*sum[1]);
687: PetscFree(lnorm2);
688: } else if (type == NORM_INFINITY || type == NORM_1) { /* max row/column sum */
689: Mat_SeqSBAIJ *amat=(Mat_SeqSBAIJ*)baij->A->data;
690: Mat_SeqBAIJ *bmat=(Mat_SeqBAIJ*)baij->B->data;
691: PetscReal *rsum,*rsum2,vabs;
692: PetscInt *jj,*garray=baij->garray,rstart=baij->rstartbs,nz;
693: PetscInt brow,bcol,col,bs=baij->A->rmap->bs,row,grow,gcol,mbs=amat->mbs;
694: MatScalar *v;
696: PetscMalloc2(mat->cmap->N,&rsum,mat->cmap->N,&rsum2);
697: PetscArrayzero(rsum,mat->cmap->N);
698: /* Amat */
699: v = amat->a; jj = amat->j;
700: for (brow=0; brow<mbs; brow++) {
701: grow = bs*(rstart + brow);
702: nz = amat->i[brow+1] - amat->i[brow];
703: for (bcol=0; bcol<nz; bcol++) {
704: gcol = bs*(rstart + *jj); jj++;
705: for (col=0; col<bs; col++) {
706: for (row=0; row<bs; row++) {
707: vabs = PetscAbsScalar(*v); v++;
708: rsum[gcol+col] += vabs;
709: /* non-diagonal block */
710: if (bcol > 0 && vabs > 0.0) rsum[grow+row] += vabs;
711: }
712: }
713: }
714: PetscLogFlops(nz*bs*bs);
715: }
716: /* Bmat */
717: v = bmat->a; jj = bmat->j;
718: for (brow=0; brow<mbs; brow++) {
719: grow = bs*(rstart + brow);
720: nz = bmat->i[brow+1] - bmat->i[brow];
721: for (bcol=0; bcol<nz; bcol++) {
722: gcol = bs*garray[*jj]; jj++;
723: for (col=0; col<bs; col++) {
724: for (row=0; row<bs; row++) {
725: vabs = PetscAbsScalar(*v); v++;
726: rsum[gcol+col] += vabs;
727: rsum[grow+row] += vabs;
728: }
729: }
730: }
731: PetscLogFlops(nz*bs*bs);
732: }
733: MPIU_Allreduce(rsum,rsum2,mat->cmap->N,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)mat));
734: *norm = 0.0;
735: for (col=0; col<mat->cmap->N; col++) {
736: if (rsum2[col] > *norm) *norm = rsum2[col];
737: }
738: PetscFree2(rsum,rsum2);
739: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for this norm yet");
740: }
741: return(0);
742: }
744: PetscErrorCode MatAssemblyBegin_MPISBAIJ(Mat mat,MatAssemblyType mode)
745: {
746: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
748: PetscInt nstash,reallocs;
751: if (baij->donotstash || mat->nooffprocentries) return(0);
753: MatStashScatterBegin_Private(mat,&mat->stash,mat->rmap->range);
754: MatStashScatterBegin_Private(mat,&mat->bstash,baij->rangebs);
755: MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);
756: PetscInfo2(mat,"Stash has %D entries,uses %D mallocs.\n",nstash,reallocs);
757: MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);
758: PetscInfo2(mat,"Block-Stash has %D entries, uses %D mallocs.\n",nstash,reallocs);
759: return(0);
760: }
762: PetscErrorCode MatAssemblyEnd_MPISBAIJ(Mat mat,MatAssemblyType mode)
763: {
764: Mat_MPISBAIJ *baij=(Mat_MPISBAIJ*)mat->data;
765: Mat_SeqSBAIJ *a =(Mat_SeqSBAIJ*)baij->A->data;
767: PetscInt i,j,rstart,ncols,flg,bs2=baij->bs2;
768: PetscInt *row,*col;
769: PetscBool other_disassembled;
770: PetscMPIInt n;
771: PetscBool r1,r2,r3;
772: MatScalar *val;
774: /* do not use 'b=(Mat_SeqBAIJ*)baij->B->data' as B can be reset in disassembly */
776: if (!baij->donotstash && !mat->nooffprocentries) {
777: while (1) {
778: MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);
779: if (!flg) break;
781: for (i=0; i<n;) {
782: /* Now identify the consecutive vals belonging to the same row */
783: for (j=i,rstart=row[j]; j<n; j++) {
784: if (row[j] != rstart) break;
785: }
786: if (j < n) ncols = j-i;
787: else ncols = n-i;
788: /* Now assemble all these values with a single function call */
789: MatSetValues_MPISBAIJ(mat,1,row+i,ncols,col+i,val+i,mat->insertmode);
790: i = j;
791: }
792: }
793: MatStashScatterEnd_Private(&mat->stash);
794: /* Now process the block-stash. Since the values are stashed column-oriented,
795: set the roworiented flag to column oriented, and after MatSetValues()
796: restore the original flags */
797: r1 = baij->roworiented;
798: r2 = a->roworiented;
799: r3 = ((Mat_SeqBAIJ*)baij->B->data)->roworiented;
801: baij->roworiented = PETSC_FALSE;
802: a->roworiented = PETSC_FALSE;
804: ((Mat_SeqBAIJ*)baij->B->data)->roworiented = PETSC_FALSE; /* b->roworinted */
805: while (1) {
806: MatStashScatterGetMesg_Private(&mat->bstash,&n,&row,&col,&val,&flg);
807: if (!flg) break;
809: for (i=0; i<n;) {
810: /* Now identify the consecutive vals belonging to the same row */
811: for (j=i,rstart=row[j]; j<n; j++) {
812: if (row[j] != rstart) break;
813: }
814: if (j < n) ncols = j-i;
815: else ncols = n-i;
816: MatSetValuesBlocked_MPISBAIJ(mat,1,row+i,ncols,col+i,val+i*bs2,mat->insertmode);
817: i = j;
818: }
819: }
820: MatStashScatterEnd_Private(&mat->bstash);
822: baij->roworiented = r1;
823: a->roworiented = r2;
825: ((Mat_SeqBAIJ*)baij->B->data)->roworiented = r3; /* b->roworinted */
826: }
828: MatAssemblyBegin(baij->A,mode);
829: MatAssemblyEnd(baij->A,mode);
831: /* determine if any processor has disassembled, if so we must
832: also disassemble ourselfs, in order that we may reassemble. */
833: /*
834: if nonzero structure of submatrix B cannot change then we know that
835: no processor disassembled thus we can skip this stuff
836: */
837: if (!((Mat_SeqBAIJ*)baij->B->data)->nonew) {
838: MPIU_Allreduce(&mat->was_assembled,&other_disassembled,1,MPIU_BOOL,MPI_PROD,PetscObjectComm((PetscObject)mat));
839: if (mat->was_assembled && !other_disassembled) {
840: MatDisAssemble_MPISBAIJ(mat);
841: }
842: }
844: if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
845: MatSetUpMultiply_MPISBAIJ(mat); /* setup Mvctx and sMvctx */
846: }
847: MatAssemblyBegin(baij->B,mode);
848: MatAssemblyEnd(baij->B,mode);
850: PetscFree2(baij->rowvalues,baij->rowindices);
852: baij->rowvalues = NULL;
854: /* if no new nonzero locations are allowed in matrix then only set the matrix state the first time through */
855: if ((!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) || !((Mat_SeqBAIJ*)(baij->A->data))->nonew) {
856: PetscObjectState state = baij->A->nonzerostate + baij->B->nonzerostate;
857: MPIU_Allreduce(&state,&mat->nonzerostate,1,MPIU_INT64,MPI_SUM,PetscObjectComm((PetscObject)mat));
858: }
859: return(0);
860: }
862: extern PetscErrorCode MatSetValues_MPIBAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
863: #include <petscdraw.h>
864: static PetscErrorCode MatView_MPISBAIJ_ASCIIorDraworSocket(Mat mat,PetscViewer viewer)
865: {
866: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
867: PetscErrorCode ierr;
868: PetscInt bs = mat->rmap->bs;
869: PetscMPIInt rank = baij->rank;
870: PetscBool iascii,isdraw;
871: PetscViewer sviewer;
872: PetscViewerFormat format;
875: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
876: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
877: if (iascii) {
878: PetscViewerGetFormat(viewer,&format);
879: if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
880: MatInfo info;
881: MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&rank);
882: MatGetInfo(mat,MAT_LOCAL,&info);
883: PetscViewerASCIIPushSynchronized(viewer);
884: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Local rows %D nz %D nz alloced %D bs %D mem %g\n",rank,mat->rmap->n,(PetscInt)info.nz_used,(PetscInt)info.nz_allocated,mat->rmap->bs,(double)info.memory);
885: MatGetInfo(baij->A,MAT_LOCAL,&info);
886: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] on-diagonal part: nz %D \n",rank,(PetscInt)info.nz_used);
887: MatGetInfo(baij->B,MAT_LOCAL,&info);
888: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] off-diagonal part: nz %D \n",rank,(PetscInt)info.nz_used);
889: PetscViewerFlush(viewer);
890: PetscViewerASCIIPopSynchronized(viewer);
891: PetscViewerASCIIPrintf(viewer,"Information on VecScatter used in matrix-vector product: \n");
892: VecScatterView(baij->Mvctx,viewer);
893: return(0);
894: } else if (format == PETSC_VIEWER_ASCII_INFO) {
895: PetscViewerASCIIPrintf(viewer," block size is %D\n",bs);
896: return(0);
897: } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
898: return(0);
899: }
900: }
902: if (isdraw) {
903: PetscDraw draw;
904: PetscBool isnull;
905: PetscViewerDrawGetDraw(viewer,0,&draw);
906: PetscDrawIsNull(draw,&isnull);
907: if (isnull) return(0);
908: }
910: {
911: /* assemble the entire matrix onto first processor. */
912: Mat A;
913: Mat_SeqSBAIJ *Aloc;
914: Mat_SeqBAIJ *Bloc;
915: PetscInt M = mat->rmap->N,N = mat->cmap->N,*ai,*aj,col,i,j,k,*rvals,mbs = baij->mbs;
916: MatScalar *a;
917: const char *matname;
919: /* Should this be the same type as mat? */
920: MatCreate(PetscObjectComm((PetscObject)mat),&A);
921: if (rank == 0) {
922: MatSetSizes(A,M,N,M,N);
923: } else {
924: MatSetSizes(A,0,0,M,N);
925: }
926: MatSetType(A,MATMPISBAIJ);
927: MatMPISBAIJSetPreallocation(A,mat->rmap->bs,0,NULL,0,NULL);
928: MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);
929: PetscLogObjectParent((PetscObject)mat,(PetscObject)A);
931: /* copy over the A part */
932: Aloc = (Mat_SeqSBAIJ*)baij->A->data;
933: ai = Aloc->i; aj = Aloc->j; a = Aloc->a;
934: PetscMalloc1(bs,&rvals);
936: for (i=0; i<mbs; i++) {
937: rvals[0] = bs*(baij->rstartbs + i);
938: for (j=1; j<bs; j++) rvals[j] = rvals[j-1] + 1;
939: for (j=ai[i]; j<ai[i+1]; j++) {
940: col = (baij->cstartbs+aj[j])*bs;
941: for (k=0; k<bs; k++) {
942: MatSetValues_MPISBAIJ(A,bs,rvals,1,&col,a,INSERT_VALUES);
943: col++;
944: a += bs;
945: }
946: }
947: }
948: /* copy over the B part */
949: Bloc = (Mat_SeqBAIJ*)baij->B->data;
950: ai = Bloc->i; aj = Bloc->j; a = Bloc->a;
951: for (i=0; i<mbs; i++) {
953: rvals[0] = bs*(baij->rstartbs + i);
954: for (j=1; j<bs; j++) rvals[j] = rvals[j-1] + 1;
955: for (j=ai[i]; j<ai[i+1]; j++) {
956: col = baij->garray[aj[j]]*bs;
957: for (k=0; k<bs; k++) {
958: MatSetValues_MPIBAIJ(A,bs,rvals,1,&col,a,INSERT_VALUES);
959: col++;
960: a += bs;
961: }
962: }
963: }
964: PetscFree(rvals);
965: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
966: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
967: /*
968: Everyone has to call to draw the matrix since the graphics waits are
969: synchronized across all processors that share the PetscDraw object
970: */
971: PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);
972: PetscObjectGetName((PetscObject)mat,&matname);
973: if (rank == 0) {
974: PetscObjectSetName((PetscObject)((Mat_MPISBAIJ*)(A->data))->A,matname);
975: MatView_SeqSBAIJ(((Mat_MPISBAIJ*)(A->data))->A,sviewer);
976: }
977: PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);
978: PetscViewerFlush(viewer);
979: MatDestroy(&A);
980: }
981: return(0);
982: }
984: /* Used for both MPIBAIJ and MPISBAIJ matrices */
985: #define MatView_MPISBAIJ_Binary MatView_MPIBAIJ_Binary
987: PetscErrorCode MatView_MPISBAIJ(Mat mat,PetscViewer viewer)
988: {
990: PetscBool iascii,isdraw,issocket,isbinary;
993: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
994: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
995: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSOCKET,&issocket);
996: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
997: if (iascii || isdraw || issocket) {
998: MatView_MPISBAIJ_ASCIIorDraworSocket(mat,viewer);
999: } else if (isbinary) {
1000: MatView_MPISBAIJ_Binary(mat,viewer);
1001: }
1002: return(0);
1003: }
1005: PetscErrorCode MatDestroy_MPISBAIJ(Mat mat)
1006: {
1007: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
1011: #if defined(PETSC_USE_LOG)
1012: PetscLogObjectState((PetscObject)mat,"Rows=%D,Cols=%D",mat->rmap->N,mat->cmap->N);
1013: #endif
1014: MatStashDestroy_Private(&mat->stash);
1015: MatStashDestroy_Private(&mat->bstash);
1016: MatDestroy(&baij->A);
1017: MatDestroy(&baij->B);
1018: #if defined(PETSC_USE_CTABLE)
1019: PetscTableDestroy(&baij->colmap);
1020: #else
1021: PetscFree(baij->colmap);
1022: #endif
1023: PetscFree(baij->garray);
1024: VecDestroy(&baij->lvec);
1025: VecScatterDestroy(&baij->Mvctx);
1026: VecDestroy(&baij->slvec0);
1027: VecDestroy(&baij->slvec0b);
1028: VecDestroy(&baij->slvec1);
1029: VecDestroy(&baij->slvec1a);
1030: VecDestroy(&baij->slvec1b);
1031: VecScatterDestroy(&baij->sMvctx);
1032: PetscFree2(baij->rowvalues,baij->rowindices);
1033: PetscFree(baij->barray);
1034: PetscFree(baij->hd);
1035: VecDestroy(&baij->diag);
1036: VecDestroy(&baij->bb1);
1037: VecDestroy(&baij->xx1);
1038: #if defined(PETSC_USE_REAL_MAT_SINGLE)
1039: PetscFree(baij->setvaluescopy);
1040: #endif
1041: PetscFree(baij->in_loc);
1042: PetscFree(baij->v_loc);
1043: PetscFree(baij->rangebs);
1044: PetscFree(mat->data);
1046: PetscObjectChangeTypeName((PetscObject)mat,NULL);
1047: PetscObjectComposeFunction((PetscObject)mat,"MatStoreValues_C",NULL);
1048: PetscObjectComposeFunction((PetscObject)mat,"MatRetrieveValues_C",NULL);
1049: PetscObjectComposeFunction((PetscObject)mat,"MatMPISBAIJSetPreallocation_C",NULL);
1050: PetscObjectComposeFunction((PetscObject)mat,"MatMPISBAIJSetPreallocationCSR_C",NULL);
1051: #if defined(PETSC_HAVE_ELEMENTAL)
1052: PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpisbaij_elemental_C",NULL);
1053: #endif
1054: #if defined(PETSC_HAVE_SCALAPACK)
1055: PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpisbaij_scalapack_C",NULL);
1056: #endif
1057: PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpisbaij_mpiaij_C",NULL);
1058: PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpisbaij_mpibaij_C",NULL);
1059: return(0);
1060: }
1062: PetscErrorCode MatMult_MPISBAIJ_Hermitian(Mat A,Vec xx,Vec yy)
1063: {
1064: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1065: PetscErrorCode ierr;
1066: PetscInt mbs=a->mbs,bs=A->rmap->bs;
1067: PetscScalar *from;
1068: const PetscScalar *x;
1071: /* diagonal part */
1072: (*a->A->ops->mult)(a->A,xx,a->slvec1a);
1073: VecSet(a->slvec1b,0.0);
1075: /* subdiagonal part */
1076: if (!a->B->ops->multhermitiantranspose) SETERRQ1(PetscObjectComm((PetscObject)a->B),PETSC_ERR_SUP,"Not for type %s\n",((PetscObject)a->B)->type_name);
1077: (*a->B->ops->multhermitiantranspose)(a->B,xx,a->slvec0b);
1079: /* copy x into the vec slvec0 */
1080: VecGetArray(a->slvec0,&from);
1081: VecGetArrayRead(xx,&x);
1083: PetscArraycpy(from,x,bs*mbs);
1084: VecRestoreArray(a->slvec0,&from);
1085: VecRestoreArrayRead(xx,&x);
1087: VecScatterBegin(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
1088: VecScatterEnd(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
1089: /* supperdiagonal part */
1090: (*a->B->ops->multadd)(a->B,a->slvec1b,a->slvec1a,yy);
1091: return(0);
1092: }
1094: PetscErrorCode MatMult_MPISBAIJ(Mat A,Vec xx,Vec yy)
1095: {
1096: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1097: PetscErrorCode ierr;
1098: PetscInt mbs=a->mbs,bs=A->rmap->bs;
1099: PetscScalar *from;
1100: const PetscScalar *x;
1103: /* diagonal part */
1104: (*a->A->ops->mult)(a->A,xx,a->slvec1a);
1105: VecSet(a->slvec1b,0.0);
1107: /* subdiagonal part */
1108: (*a->B->ops->multtranspose)(a->B,xx,a->slvec0b);
1110: /* copy x into the vec slvec0 */
1111: VecGetArray(a->slvec0,&from);
1112: VecGetArrayRead(xx,&x);
1114: PetscArraycpy(from,x,bs*mbs);
1115: VecRestoreArray(a->slvec0,&from);
1116: VecRestoreArrayRead(xx,&x);
1118: VecScatterBegin(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
1119: VecScatterEnd(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
1120: /* supperdiagonal part */
1121: (*a->B->ops->multadd)(a->B,a->slvec1b,a->slvec1a,yy);
1122: return(0);
1123: }
1125: PetscErrorCode MatMultAdd_MPISBAIJ_Hermitian(Mat A,Vec xx,Vec yy,Vec zz)
1126: {
1127: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1128: PetscErrorCode ierr;
1129: PetscInt mbs=a->mbs,bs=A->rmap->bs;
1130: PetscScalar *from,zero=0.0;
1131: const PetscScalar *x;
1134: /* diagonal part */
1135: (*a->A->ops->multadd)(a->A,xx,yy,a->slvec1a);
1136: VecSet(a->slvec1b,zero);
1138: /* subdiagonal part */
1139: if (!a->B->ops->multhermitiantranspose) SETERRQ1(PetscObjectComm((PetscObject)a->B),PETSC_ERR_SUP,"Not for type %s\n",((PetscObject)a->B)->type_name);
1140: (*a->B->ops->multhermitiantranspose)(a->B,xx,a->slvec0b);
1142: /* copy x into the vec slvec0 */
1143: VecGetArray(a->slvec0,&from);
1144: VecGetArrayRead(xx,&x);
1145: PetscArraycpy(from,x,bs*mbs);
1146: VecRestoreArray(a->slvec0,&from);
1148: VecScatterBegin(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
1149: VecRestoreArrayRead(xx,&x);
1150: VecScatterEnd(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
1152: /* supperdiagonal part */
1153: (*a->B->ops->multadd)(a->B,a->slvec1b,a->slvec1a,zz);
1154: return(0);
1155: }
1157: PetscErrorCode MatMultAdd_MPISBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
1158: {
1159: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1160: PetscErrorCode ierr;
1161: PetscInt mbs=a->mbs,bs=A->rmap->bs;
1162: PetscScalar *from,zero=0.0;
1163: const PetscScalar *x;
1166: /* diagonal part */
1167: (*a->A->ops->multadd)(a->A,xx,yy,a->slvec1a);
1168: VecSet(a->slvec1b,zero);
1170: /* subdiagonal part */
1171: (*a->B->ops->multtranspose)(a->B,xx,a->slvec0b);
1173: /* copy x into the vec slvec0 */
1174: VecGetArray(a->slvec0,&from);
1175: VecGetArrayRead(xx,&x);
1176: PetscArraycpy(from,x,bs*mbs);
1177: VecRestoreArray(a->slvec0,&from);
1179: VecScatterBegin(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
1180: VecRestoreArrayRead(xx,&x);
1181: VecScatterEnd(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
1183: /* supperdiagonal part */
1184: (*a->B->ops->multadd)(a->B,a->slvec1b,a->slvec1a,zz);
1185: return(0);
1186: }
1188: /*
1189: This only works correctly for square matrices where the subblock A->A is the
1190: diagonal block
1191: */
1192: PetscErrorCode MatGetDiagonal_MPISBAIJ(Mat A,Vec v)
1193: {
1194: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1198: /* if (a->rmap->N != a->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Supports only square matrix where A->A is diag block"); */
1199: MatGetDiagonal(a->A,v);
1200: return(0);
1201: }
1203: PetscErrorCode MatScale_MPISBAIJ(Mat A,PetscScalar aa)
1204: {
1205: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1209: MatScale(a->A,aa);
1210: MatScale(a->B,aa);
1211: return(0);
1212: }
1214: PetscErrorCode MatGetRow_MPISBAIJ(Mat matin,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
1215: {
1216: Mat_MPISBAIJ *mat = (Mat_MPISBAIJ*)matin->data;
1217: PetscScalar *vworkA,*vworkB,**pvA,**pvB,*v_p;
1219: PetscInt bs = matin->rmap->bs,bs2 = mat->bs2,i,*cworkA,*cworkB,**pcA,**pcB;
1220: PetscInt nztot,nzA,nzB,lrow,brstart = matin->rmap->rstart,brend = matin->rmap->rend;
1221: PetscInt *cmap,*idx_p,cstart = mat->rstartbs;
1224: if (mat->getrowactive) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Already active");
1225: mat->getrowactive = PETSC_TRUE;
1227: if (!mat->rowvalues && (idx || v)) {
1228: /*
1229: allocate enough space to hold information from the longest row.
1230: */
1231: Mat_SeqSBAIJ *Aa = (Mat_SeqSBAIJ*)mat->A->data;
1232: Mat_SeqBAIJ *Ba = (Mat_SeqBAIJ*)mat->B->data;
1233: PetscInt max = 1,mbs = mat->mbs,tmp;
1234: for (i=0; i<mbs; i++) {
1235: tmp = Aa->i[i+1] - Aa->i[i] + Ba->i[i+1] - Ba->i[i]; /* row length */
1236: if (max < tmp) max = tmp;
1237: }
1238: PetscMalloc2(max*bs2,&mat->rowvalues,max*bs2,&mat->rowindices);
1239: }
1241: if (row < brstart || row >= brend) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only local rows");
1242: lrow = row - brstart; /* local row index */
1244: pvA = &vworkA; pcA = &cworkA; pvB = &vworkB; pcB = &cworkB;
1245: if (!v) {pvA = NULL; pvB = NULL;}
1246: if (!idx) {pcA = NULL; if (!v) pcB = NULL;}
1247: (*mat->A->ops->getrow)(mat->A,lrow,&nzA,pcA,pvA);
1248: (*mat->B->ops->getrow)(mat->B,lrow,&nzB,pcB,pvB);
1249: nztot = nzA + nzB;
1251: cmap = mat->garray;
1252: if (v || idx) {
1253: if (nztot) {
1254: /* Sort by increasing column numbers, assuming A and B already sorted */
1255: PetscInt imark = -1;
1256: if (v) {
1257: *v = v_p = mat->rowvalues;
1258: for (i=0; i<nzB; i++) {
1259: if (cmap[cworkB[i]/bs] < cstart) v_p[i] = vworkB[i];
1260: else break;
1261: }
1262: imark = i;
1263: for (i=0; i<nzA; i++) v_p[imark+i] = vworkA[i];
1264: for (i=imark; i<nzB; i++) v_p[nzA+i] = vworkB[i];
1265: }
1266: if (idx) {
1267: *idx = idx_p = mat->rowindices;
1268: if (imark > -1) {
1269: for (i=0; i<imark; i++) {
1270: idx_p[i] = cmap[cworkB[i]/bs]*bs + cworkB[i]%bs;
1271: }
1272: } else {
1273: for (i=0; i<nzB; i++) {
1274: if (cmap[cworkB[i]/bs] < cstart) idx_p[i] = cmap[cworkB[i]/bs]*bs + cworkB[i]%bs;
1275: else break;
1276: }
1277: imark = i;
1278: }
1279: for (i=0; i<nzA; i++) idx_p[imark+i] = cstart*bs + cworkA[i];
1280: for (i=imark; i<nzB; i++) idx_p[nzA+i] = cmap[cworkB[i]/bs]*bs + cworkB[i]%bs ;
1281: }
1282: } else {
1283: if (idx) *idx = NULL;
1284: if (v) *v = NULL;
1285: }
1286: }
1287: *nz = nztot;
1288: (*mat->A->ops->restorerow)(mat->A,lrow,&nzA,pcA,pvA);
1289: (*mat->B->ops->restorerow)(mat->B,lrow,&nzB,pcB,pvB);
1290: return(0);
1291: }
1293: PetscErrorCode MatRestoreRow_MPISBAIJ(Mat mat,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
1294: {
1295: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
1298: if (!baij->getrowactive) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"MatGetRow() must be called first");
1299: baij->getrowactive = PETSC_FALSE;
1300: return(0);
1301: }
1303: PetscErrorCode MatGetRowUpperTriangular_MPISBAIJ(Mat A)
1304: {
1305: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1306: Mat_SeqSBAIJ *aA = (Mat_SeqSBAIJ*)a->A->data;
1309: aA->getrow_utriangular = PETSC_TRUE;
1310: return(0);
1311: }
1312: PetscErrorCode MatRestoreRowUpperTriangular_MPISBAIJ(Mat A)
1313: {
1314: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1315: Mat_SeqSBAIJ *aA = (Mat_SeqSBAIJ*)a->A->data;
1318: aA->getrow_utriangular = PETSC_FALSE;
1319: return(0);
1320: }
1322: PetscErrorCode MatConjugate_MPISBAIJ(Mat mat)
1323: {
1324: #if defined(PETSC_USE_COMPLEX)
1326: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)mat->data;
1329: MatConjugate(a->A);
1330: MatConjugate(a->B);
1331: #else
1333: #endif
1334: return(0);
1335: }
1337: PetscErrorCode MatRealPart_MPISBAIJ(Mat A)
1338: {
1339: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1343: MatRealPart(a->A);
1344: MatRealPart(a->B);
1345: return(0);
1346: }
1348: PetscErrorCode MatImaginaryPart_MPISBAIJ(Mat A)
1349: {
1350: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1354: MatImaginaryPart(a->A);
1355: MatImaginaryPart(a->B);
1356: return(0);
1357: }
1359: /* Check if isrow is a subset of iscol_local, called by MatCreateSubMatrix_MPISBAIJ()
1360: Input: isrow - distributed(parallel),
1361: iscol_local - locally owned (seq)
1362: */
1363: PetscErrorCode ISEqual_private(IS isrow,IS iscol_local,PetscBool *flg)
1364: {
1366: PetscInt sz1,sz2,*a1,*a2,i,j,k,nmatch;
1367: const PetscInt *ptr1,*ptr2;
1370: ISGetLocalSize(isrow,&sz1);
1371: ISGetLocalSize(iscol_local,&sz2);
1372: if (sz1 > sz2) {
1373: *flg = PETSC_FALSE;
1374: return(0);
1375: }
1377: ISGetIndices(isrow,&ptr1);
1378: ISGetIndices(iscol_local,&ptr2);
1380: PetscMalloc1(sz1,&a1);
1381: PetscMalloc1(sz2,&a2);
1382: PetscArraycpy(a1,ptr1,sz1);
1383: PetscArraycpy(a2,ptr2,sz2);
1384: PetscSortInt(sz1,a1);
1385: PetscSortInt(sz2,a2);
1387: nmatch=0;
1388: k = 0;
1389: for (i=0; i<sz1; i++) {
1390: for (j=k; j<sz2; j++) {
1391: if (a1[i] == a2[j]) {
1392: k = j; nmatch++;
1393: break;
1394: }
1395: }
1396: }
1397: ISRestoreIndices(isrow,&ptr1);
1398: ISRestoreIndices(iscol_local,&ptr2);
1399: PetscFree(a1);
1400: PetscFree(a2);
1401: if (nmatch < sz1) {
1402: *flg = PETSC_FALSE;
1403: } else {
1404: *flg = PETSC_TRUE;
1405: }
1406: return(0);
1407: }
1409: PetscErrorCode MatCreateSubMatrix_MPISBAIJ(Mat mat,IS isrow,IS iscol,MatReuse call,Mat *newmat)
1410: {
1412: IS iscol_local;
1413: PetscInt csize;
1414: PetscBool isequal;
1417: ISGetLocalSize(iscol,&csize);
1418: if (call == MAT_REUSE_MATRIX) {
1419: PetscObjectQuery((PetscObject)*newmat,"ISAllGather",(PetscObject*)&iscol_local);
1420: if (!iscol_local) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Submatrix passed in was not used before, cannot reuse");
1421: } else {
1422: PetscBool issorted;
1424: ISAllGather(iscol,&iscol_local);
1425: ISEqual_private(isrow,iscol_local,&isequal);
1426: ISSorted(iscol_local, &issorted);
1427: if (!isequal || !issorted) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"For symmetric format, iscol must equal isrow and be sorted");
1428: }
1430: /* now call MatCreateSubMatrix_MPIBAIJ() */
1431: MatCreateSubMatrix_MPIBAIJ_Private(mat,isrow,iscol_local,csize,call,newmat);
1432: if (call == MAT_INITIAL_MATRIX) {
1433: PetscObjectCompose((PetscObject)*newmat,"ISAllGather",(PetscObject)iscol_local);
1434: ISDestroy(&iscol_local);
1435: }
1436: return(0);
1437: }
1439: PetscErrorCode MatZeroEntries_MPISBAIJ(Mat A)
1440: {
1441: Mat_MPISBAIJ *l = (Mat_MPISBAIJ*)A->data;
1445: MatZeroEntries(l->A);
1446: MatZeroEntries(l->B);
1447: return(0);
1448: }
1450: PetscErrorCode MatGetInfo_MPISBAIJ(Mat matin,MatInfoType flag,MatInfo *info)
1451: {
1452: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)matin->data;
1453: Mat A = a->A,B = a->B;
1455: PetscLogDouble isend[5],irecv[5];
1458: info->block_size = (PetscReal)matin->rmap->bs;
1460: MatGetInfo(A,MAT_LOCAL,info);
1462: isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
1463: isend[3] = info->memory; isend[4] = info->mallocs;
1465: MatGetInfo(B,MAT_LOCAL,info);
1467: isend[0] += info->nz_used; isend[1] += info->nz_allocated; isend[2] += info->nz_unneeded;
1468: isend[3] += info->memory; isend[4] += info->mallocs;
1469: if (flag == MAT_LOCAL) {
1470: info->nz_used = isend[0];
1471: info->nz_allocated = isend[1];
1472: info->nz_unneeded = isend[2];
1473: info->memory = isend[3];
1474: info->mallocs = isend[4];
1475: } else if (flag == MAT_GLOBAL_MAX) {
1476: MPIU_Allreduce(isend,irecv,5,MPIU_PETSCLOGDOUBLE,MPI_MAX,PetscObjectComm((PetscObject)matin));
1478: info->nz_used = irecv[0];
1479: info->nz_allocated = irecv[1];
1480: info->nz_unneeded = irecv[2];
1481: info->memory = irecv[3];
1482: info->mallocs = irecv[4];
1483: } else if (flag == MAT_GLOBAL_SUM) {
1484: MPIU_Allreduce(isend,irecv,5,MPIU_PETSCLOGDOUBLE,MPI_SUM,PetscObjectComm((PetscObject)matin));
1486: info->nz_used = irecv[0];
1487: info->nz_allocated = irecv[1];
1488: info->nz_unneeded = irecv[2];
1489: info->memory = irecv[3];
1490: info->mallocs = irecv[4];
1491: } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown MatInfoType argument %d",(int)flag);
1492: info->fill_ratio_given = 0; /* no parallel LU/ILU/Cholesky */
1493: info->fill_ratio_needed = 0;
1494: info->factor_mallocs = 0;
1495: return(0);
1496: }
1498: PetscErrorCode MatSetOption_MPISBAIJ(Mat A,MatOption op,PetscBool flg)
1499: {
1500: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1501: Mat_SeqSBAIJ *aA = (Mat_SeqSBAIJ*)a->A->data;
1505: switch (op) {
1506: case MAT_NEW_NONZERO_LOCATIONS:
1507: case MAT_NEW_NONZERO_ALLOCATION_ERR:
1508: case MAT_UNUSED_NONZERO_LOCATION_ERR:
1509: case MAT_KEEP_NONZERO_PATTERN:
1510: case MAT_SUBMAT_SINGLEIS:
1511: case MAT_NEW_NONZERO_LOCATION_ERR:
1512: MatCheckPreallocated(A,1);
1513: MatSetOption(a->A,op,flg);
1514: MatSetOption(a->B,op,flg);
1515: break;
1516: case MAT_ROW_ORIENTED:
1517: MatCheckPreallocated(A,1);
1518: a->roworiented = flg;
1520: MatSetOption(a->A,op,flg);
1521: MatSetOption(a->B,op,flg);
1522: break;
1523: case MAT_FORCE_DIAGONAL_ENTRIES:
1524: case MAT_SORTED_FULL:
1525: PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
1526: break;
1527: case MAT_IGNORE_OFF_PROC_ENTRIES:
1528: a->donotstash = flg;
1529: break;
1530: case MAT_USE_HASH_TABLE:
1531: a->ht_flag = flg;
1532: break;
1533: case MAT_HERMITIAN:
1534: MatCheckPreallocated(A,1);
1535: MatSetOption(a->A,op,flg);
1536: #if defined(PETSC_USE_COMPLEX)
1537: if (flg) { /* need different mat-vec ops */
1538: A->ops->mult = MatMult_MPISBAIJ_Hermitian;
1539: A->ops->multadd = MatMultAdd_MPISBAIJ_Hermitian;
1540: A->ops->multtranspose = NULL;
1541: A->ops->multtransposeadd = NULL;
1542: A->symmetric = PETSC_FALSE;
1543: }
1544: #endif
1545: break;
1546: case MAT_SPD:
1547: case MAT_SYMMETRIC:
1548: MatCheckPreallocated(A,1);
1549: MatSetOption(a->A,op,flg);
1550: #if defined(PETSC_USE_COMPLEX)
1551: if (flg) { /* restore to use default mat-vec ops */
1552: A->ops->mult = MatMult_MPISBAIJ;
1553: A->ops->multadd = MatMultAdd_MPISBAIJ;
1554: A->ops->multtranspose = MatMult_MPISBAIJ;
1555: A->ops->multtransposeadd = MatMultAdd_MPISBAIJ;
1556: }
1557: #endif
1558: break;
1559: case MAT_STRUCTURALLY_SYMMETRIC:
1560: MatCheckPreallocated(A,1);
1561: MatSetOption(a->A,op,flg);
1562: break;
1563: case MAT_SYMMETRY_ETERNAL:
1564: if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix must be symmetric");
1565: PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
1566: break;
1567: case MAT_IGNORE_LOWER_TRIANGULAR:
1568: aA->ignore_ltriangular = flg;
1569: break;
1570: case MAT_ERROR_LOWER_TRIANGULAR:
1571: aA->ignore_ltriangular = flg;
1572: break;
1573: case MAT_GETROW_UPPERTRIANGULAR:
1574: aA->getrow_utriangular = flg;
1575: break;
1576: default:
1577: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %d",op);
1578: }
1579: return(0);
1580: }
1582: PetscErrorCode MatTranspose_MPISBAIJ(Mat A,MatReuse reuse,Mat *B)
1583: {
1587: if (reuse == MAT_INITIAL_MATRIX) {
1588: MatDuplicate(A,MAT_COPY_VALUES,B);
1589: } else if (reuse == MAT_REUSE_MATRIX) {
1590: MatCopy(A,*B,SAME_NONZERO_PATTERN);
1591: }
1592: return(0);
1593: }
1595: PetscErrorCode MatDiagonalScale_MPISBAIJ(Mat mat,Vec ll,Vec rr)
1596: {
1597: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
1598: Mat a = baij->A, b=baij->B;
1600: PetscInt nv,m,n;
1601: PetscBool flg;
1604: if (ll != rr) {
1605: VecEqual(ll,rr,&flg);
1606: if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"For symmetric format, left and right scaling vectors must be same\n");
1607: }
1608: if (!ll) return(0);
1610: MatGetLocalSize(mat,&m,&n);
1611: if (m != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"For symmetric format, local size %d %d must be same",m,n);
1613: VecGetLocalSize(rr,&nv);
1614: if (nv!=n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left and right vector non-conforming local size");
1616: VecScatterBegin(baij->Mvctx,rr,baij->lvec,INSERT_VALUES,SCATTER_FORWARD);
1618: /* left diagonalscale the off-diagonal part */
1619: (*b->ops->diagonalscale)(b,ll,NULL);
1621: /* scale the diagonal part */
1622: (*a->ops->diagonalscale)(a,ll,rr);
1624: /* right diagonalscale the off-diagonal part */
1625: VecScatterEnd(baij->Mvctx,rr,baij->lvec,INSERT_VALUES,SCATTER_FORWARD);
1626: (*b->ops->diagonalscale)(b,NULL,baij->lvec);
1627: return(0);
1628: }
1630: PetscErrorCode MatSetUnfactored_MPISBAIJ(Mat A)
1631: {
1632: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1636: MatSetUnfactored(a->A);
1637: return(0);
1638: }
1640: static PetscErrorCode MatDuplicate_MPISBAIJ(Mat,MatDuplicateOption,Mat*);
1642: PetscErrorCode MatEqual_MPISBAIJ(Mat A,Mat B,PetscBool *flag)
1643: {
1644: Mat_MPISBAIJ *matB = (Mat_MPISBAIJ*)B->data,*matA = (Mat_MPISBAIJ*)A->data;
1645: Mat a,b,c,d;
1646: PetscBool flg;
1650: a = matA->A; b = matA->B;
1651: c = matB->A; d = matB->B;
1653: MatEqual(a,c,&flg);
1654: if (flg) {
1655: MatEqual(b,d,&flg);
1656: }
1657: MPIU_Allreduce(&flg,flag,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));
1658: return(0);
1659: }
1661: PetscErrorCode MatCopy_MPISBAIJ(Mat A,Mat B,MatStructure str)
1662: {
1664: PetscBool isbaij;
1667: PetscObjectTypeCompareAny((PetscObject)B,&isbaij,MATSEQSBAIJ,MATMPISBAIJ,"");
1668: if (!isbaij) SETERRQ1(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"Not for matrix type %s",((PetscObject)B)->type_name);
1669: /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
1670: if ((str != SAME_NONZERO_PATTERN) || (A->ops->copy != B->ops->copy)) {
1671: MatGetRowUpperTriangular(A);
1672: MatCopy_Basic(A,B,str);
1673: MatRestoreRowUpperTriangular(A);
1674: } else {
1675: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1676: Mat_MPISBAIJ *b = (Mat_MPISBAIJ*)B->data;
1678: MatCopy(a->A,b->A,str);
1679: MatCopy(a->B,b->B,str);
1680: }
1681: PetscObjectStateIncrease((PetscObject)B);
1682: return(0);
1683: }
1685: PetscErrorCode MatSetUp_MPISBAIJ(Mat A)
1686: {
1690: MatMPISBAIJSetPreallocation(A,A->rmap->bs,PETSC_DEFAULT,NULL,PETSC_DEFAULT,NULL);
1691: return(0);
1692: }
1694: PetscErrorCode MatAXPY_MPISBAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str)
1695: {
1697: Mat_MPISBAIJ *xx=(Mat_MPISBAIJ*)X->data,*yy=(Mat_MPISBAIJ*)Y->data;
1698: PetscBLASInt bnz,one=1;
1699: Mat_SeqSBAIJ *xa,*ya;
1700: Mat_SeqBAIJ *xb,*yb;
1703: if (str == SAME_NONZERO_PATTERN) {
1704: PetscScalar alpha = a;
1705: xa = (Mat_SeqSBAIJ*)xx->A->data;
1706: ya = (Mat_SeqSBAIJ*)yy->A->data;
1707: PetscBLASIntCast(xa->nz,&bnz);
1708: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&bnz,&alpha,xa->a,&one,ya->a,&one));
1709: xb = (Mat_SeqBAIJ*)xx->B->data;
1710: yb = (Mat_SeqBAIJ*)yy->B->data;
1711: PetscBLASIntCast(xb->nz,&bnz);
1712: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&bnz,&alpha,xb->a,&one,yb->a,&one));
1713: PetscObjectStateIncrease((PetscObject)Y);
1714: } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
1715: MatSetOption(X,MAT_GETROW_UPPERTRIANGULAR,PETSC_TRUE);
1716: MatAXPY_Basic(Y,a,X,str);
1717: MatSetOption(X,MAT_GETROW_UPPERTRIANGULAR,PETSC_FALSE);
1718: } else {
1719: Mat B;
1720: PetscInt *nnz_d,*nnz_o,bs=Y->rmap->bs;
1721: if (bs != X->rmap->bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrices must have same block size");
1722: MatGetRowUpperTriangular(X);
1723: MatGetRowUpperTriangular(Y);
1724: PetscMalloc1(yy->A->rmap->N,&nnz_d);
1725: PetscMalloc1(yy->B->rmap->N,&nnz_o);
1726: MatCreate(PetscObjectComm((PetscObject)Y),&B);
1727: PetscObjectSetName((PetscObject)B,((PetscObject)Y)->name);
1728: MatSetSizes(B,Y->rmap->n,Y->cmap->n,Y->rmap->N,Y->cmap->N);
1729: MatSetBlockSizesFromMats(B,Y,Y);
1730: MatSetType(B,MATMPISBAIJ);
1731: MatAXPYGetPreallocation_SeqSBAIJ(yy->A,xx->A,nnz_d);
1732: MatAXPYGetPreallocation_MPIBAIJ(yy->B,yy->garray,xx->B,xx->garray,nnz_o);
1733: MatMPISBAIJSetPreallocation(B,bs,0,nnz_d,0,nnz_o);
1734: MatAXPY_BasicWithPreallocation(B,Y,a,X,str);
1735: MatHeaderReplace(Y,&B);
1736: PetscFree(nnz_d);
1737: PetscFree(nnz_o);
1738: MatRestoreRowUpperTriangular(X);
1739: MatRestoreRowUpperTriangular(Y);
1740: }
1741: return(0);
1742: }
1744: PetscErrorCode MatCreateSubMatrices_MPISBAIJ(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1745: {
1747: PetscInt i;
1748: PetscBool flg;
1751: MatCreateSubMatrices_MPIBAIJ(A,n,irow,icol,scall,B); /* B[] are sbaij matrices */
1752: for (i=0; i<n; i++) {
1753: ISEqual(irow[i],icol[i],&flg);
1754: if (!flg) {
1755: MatSeqSBAIJZeroOps_Private(*B[i]);
1756: }
1757: }
1758: return(0);
1759: }
1761: PetscErrorCode MatShift_MPISBAIJ(Mat Y,PetscScalar a)
1762: {
1764: Mat_MPISBAIJ *maij = (Mat_MPISBAIJ*)Y->data;
1765: Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ*)maij->A->data;
1768: if (!Y->preallocated) {
1769: MatMPISBAIJSetPreallocation(Y,Y->rmap->bs,1,NULL,0,NULL);
1770: } else if (!aij->nz) {
1771: PetscInt nonew = aij->nonew;
1772: MatSeqSBAIJSetPreallocation(maij->A,Y->rmap->bs,1,NULL);
1773: aij->nonew = nonew;
1774: }
1775: MatShift_Basic(Y,a);
1776: return(0);
1777: }
1779: PetscErrorCode MatMissingDiagonal_MPISBAIJ(Mat A,PetscBool *missing,PetscInt *d)
1780: {
1781: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1785: if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only works for square matrices");
1786: MatMissingDiagonal(a->A,missing,d);
1787: if (d) {
1788: PetscInt rstart;
1789: MatGetOwnershipRange(A,&rstart,NULL);
1790: *d += rstart/A->rmap->bs;
1792: }
1793: return(0);
1794: }
1796: PetscErrorCode MatGetDiagonalBlock_MPISBAIJ(Mat A,Mat *a)
1797: {
1799: *a = ((Mat_MPISBAIJ*)A->data)->A;
1800: return(0);
1801: }
1803: /* -------------------------------------------------------------------*/
1804: static struct _MatOps MatOps_Values = {MatSetValues_MPISBAIJ,
1805: MatGetRow_MPISBAIJ,
1806: MatRestoreRow_MPISBAIJ,
1807: MatMult_MPISBAIJ,
1808: /* 4*/ MatMultAdd_MPISBAIJ,
1809: MatMult_MPISBAIJ, /* transpose versions are same as non-transpose */
1810: MatMultAdd_MPISBAIJ,
1811: NULL,
1812: NULL,
1813: NULL,
1814: /* 10*/ NULL,
1815: NULL,
1816: NULL,
1817: MatSOR_MPISBAIJ,
1818: MatTranspose_MPISBAIJ,
1819: /* 15*/ MatGetInfo_MPISBAIJ,
1820: MatEqual_MPISBAIJ,
1821: MatGetDiagonal_MPISBAIJ,
1822: MatDiagonalScale_MPISBAIJ,
1823: MatNorm_MPISBAIJ,
1824: /* 20*/ MatAssemblyBegin_MPISBAIJ,
1825: MatAssemblyEnd_MPISBAIJ,
1826: MatSetOption_MPISBAIJ,
1827: MatZeroEntries_MPISBAIJ,
1828: /* 24*/ NULL,
1829: NULL,
1830: NULL,
1831: NULL,
1832: NULL,
1833: /* 29*/ MatSetUp_MPISBAIJ,
1834: NULL,
1835: NULL,
1836: MatGetDiagonalBlock_MPISBAIJ,
1837: NULL,
1838: /* 34*/ MatDuplicate_MPISBAIJ,
1839: NULL,
1840: NULL,
1841: NULL,
1842: NULL,
1843: /* 39*/ MatAXPY_MPISBAIJ,
1844: MatCreateSubMatrices_MPISBAIJ,
1845: MatIncreaseOverlap_MPISBAIJ,
1846: MatGetValues_MPISBAIJ,
1847: MatCopy_MPISBAIJ,
1848: /* 44*/ NULL,
1849: MatScale_MPISBAIJ,
1850: MatShift_MPISBAIJ,
1851: NULL,
1852: NULL,
1853: /* 49*/ NULL,
1854: NULL,
1855: NULL,
1856: NULL,
1857: NULL,
1858: /* 54*/ NULL,
1859: NULL,
1860: MatSetUnfactored_MPISBAIJ,
1861: NULL,
1862: MatSetValuesBlocked_MPISBAIJ,
1863: /* 59*/ MatCreateSubMatrix_MPISBAIJ,
1864: NULL,
1865: NULL,
1866: NULL,
1867: NULL,
1868: /* 64*/ NULL,
1869: NULL,
1870: NULL,
1871: NULL,
1872: NULL,
1873: /* 69*/ MatGetRowMaxAbs_MPISBAIJ,
1874: NULL,
1875: MatConvert_MPISBAIJ_Basic,
1876: NULL,
1877: NULL,
1878: /* 74*/ NULL,
1879: NULL,
1880: NULL,
1881: NULL,
1882: NULL,
1883: /* 79*/ NULL,
1884: NULL,
1885: NULL,
1886: NULL,
1887: MatLoad_MPISBAIJ,
1888: /* 84*/ NULL,
1889: NULL,
1890: NULL,
1891: NULL,
1892: NULL,
1893: /* 89*/ NULL,
1894: NULL,
1895: NULL,
1896: NULL,
1897: NULL,
1898: /* 94*/ NULL,
1899: NULL,
1900: NULL,
1901: NULL,
1902: NULL,
1903: /* 99*/ NULL,
1904: NULL,
1905: NULL,
1906: MatConjugate_MPISBAIJ,
1907: NULL,
1908: /*104*/ NULL,
1909: MatRealPart_MPISBAIJ,
1910: MatImaginaryPart_MPISBAIJ,
1911: MatGetRowUpperTriangular_MPISBAIJ,
1912: MatRestoreRowUpperTriangular_MPISBAIJ,
1913: /*109*/ NULL,
1914: NULL,
1915: NULL,
1916: NULL,
1917: MatMissingDiagonal_MPISBAIJ,
1918: /*114*/ NULL,
1919: NULL,
1920: NULL,
1921: NULL,
1922: NULL,
1923: /*119*/ NULL,
1924: NULL,
1925: NULL,
1926: NULL,
1927: NULL,
1928: /*124*/ NULL,
1929: NULL,
1930: NULL,
1931: NULL,
1932: NULL,
1933: /*129*/ NULL,
1934: NULL,
1935: NULL,
1936: NULL,
1937: NULL,
1938: /*134*/ NULL,
1939: NULL,
1940: NULL,
1941: NULL,
1942: NULL,
1943: /*139*/ MatSetBlockSizes_Default,
1944: NULL,
1945: NULL,
1946: NULL,
1947: NULL,
1948: /*144*/MatCreateMPIMatConcatenateSeqMat_MPISBAIJ
1949: };
1951: PetscErrorCode MatMPISBAIJSetPreallocation_MPISBAIJ(Mat B,PetscInt bs,PetscInt d_nz,const PetscInt *d_nnz,PetscInt o_nz,const PetscInt *o_nnz)
1952: {
1953: Mat_MPISBAIJ *b = (Mat_MPISBAIJ*)B->data;
1955: PetscInt i,mbs,Mbs;
1956: PetscMPIInt size;
1959: MatSetBlockSize(B,PetscAbs(bs));
1960: PetscLayoutSetUp(B->rmap);
1961: PetscLayoutSetUp(B->cmap);
1962: PetscLayoutGetBlockSize(B->rmap,&bs);
1963: if (B->rmap->N > B->cmap->N) SETERRQ2(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"MPISBAIJ matrix cannot have more rows %D than columns %D",B->rmap->N,B->cmap->N);
1964: if (B->rmap->n > B->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"MPISBAIJ matrix cannot have more local rows %D than columns %D",B->rmap->n,B->cmap->n);
1966: mbs = B->rmap->n/bs;
1967: Mbs = B->rmap->N/bs;
1968: 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);
1970: B->rmap->bs = bs;
1971: b->bs2 = bs*bs;
1972: b->mbs = mbs;
1973: b->Mbs = Mbs;
1974: b->nbs = B->cmap->n/bs;
1975: b->Nbs = B->cmap->N/bs;
1977: for (i=0; i<=b->size; i++) {
1978: b->rangebs[i] = B->rmap->range[i]/bs;
1979: }
1980: b->rstartbs = B->rmap->rstart/bs;
1981: b->rendbs = B->rmap->rend/bs;
1983: b->cstartbs = B->cmap->rstart/bs;
1984: b->cendbs = B->cmap->rend/bs;
1986: #if defined(PETSC_USE_CTABLE)
1987: PetscTableDestroy(&b->colmap);
1988: #else
1989: PetscFree(b->colmap);
1990: #endif
1991: PetscFree(b->garray);
1992: VecDestroy(&b->lvec);
1993: VecScatterDestroy(&b->Mvctx);
1994: VecDestroy(&b->slvec0);
1995: VecDestroy(&b->slvec0b);
1996: VecDestroy(&b->slvec1);
1997: VecDestroy(&b->slvec1a);
1998: VecDestroy(&b->slvec1b);
1999: VecScatterDestroy(&b->sMvctx);
2001: /* Because the B will have been resized we simply destroy it and create a new one each time */
2002: MPI_Comm_size(PetscObjectComm((PetscObject)B),&size);
2003: MatDestroy(&b->B);
2004: MatCreate(PETSC_COMM_SELF,&b->B);
2005: MatSetSizes(b->B,B->rmap->n,size > 1 ? B->cmap->N : 0,B->rmap->n,size > 1 ? B->cmap->N : 0);
2006: MatSetType(b->B,MATSEQBAIJ);
2007: PetscLogObjectParent((PetscObject)B,(PetscObject)b->B);
2009: if (!B->preallocated) {
2010: MatCreate(PETSC_COMM_SELF,&b->A);
2011: MatSetSizes(b->A,B->rmap->n,B->cmap->n,B->rmap->n,B->cmap->n);
2012: MatSetType(b->A,MATSEQSBAIJ);
2013: PetscLogObjectParent((PetscObject)B,(PetscObject)b->A);
2014: MatStashCreate_Private(PetscObjectComm((PetscObject)B),bs,&B->bstash);
2015: }
2017: MatSeqSBAIJSetPreallocation(b->A,bs,d_nz,d_nnz);
2018: MatSeqBAIJSetPreallocation(b->B,bs,o_nz,o_nnz);
2020: B->preallocated = PETSC_TRUE;
2021: B->was_assembled = PETSC_FALSE;
2022: B->assembled = PETSC_FALSE;
2023: return(0);
2024: }
2026: PetscErrorCode MatMPISBAIJSetPreallocationCSR_MPISBAIJ(Mat B,PetscInt bs,const PetscInt ii[],const PetscInt jj[],const PetscScalar V[])
2027: {
2028: PetscInt m,rstart,cend;
2029: PetscInt i,j,d,nz,bd, nz_max=0,*d_nnz=NULL,*o_nnz=NULL;
2030: const PetscInt *JJ =NULL;
2031: PetscScalar *values=NULL;
2032: PetscBool roworiented = ((Mat_MPISBAIJ*)B->data)->roworiented;
2034: PetscBool nooffprocentries;
2037: if (bs < 1) SETERRQ1(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_OUTOFRANGE,"Invalid block size specified, must be positive but it is %D",bs);
2038: PetscLayoutSetBlockSize(B->rmap,bs);
2039: PetscLayoutSetBlockSize(B->cmap,bs);
2040: PetscLayoutSetUp(B->rmap);
2041: PetscLayoutSetUp(B->cmap);
2042: PetscLayoutGetBlockSize(B->rmap,&bs);
2043: m = B->rmap->n/bs;
2044: rstart = B->rmap->rstart/bs;
2045: cend = B->cmap->rend/bs;
2047: if (ii[0]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"ii[0] must be 0 but it is %D",ii[0]);
2048: PetscMalloc2(m,&d_nnz,m,&o_nnz);
2049: for (i=0; i<m; i++) {
2050: nz = ii[i+1] - ii[i];
2051: if (nz < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local row %D has a negative number of columns %D",i,nz);
2052: /* count the ones on the diagonal and above, split into diagonal and off diagonal portions. */
2053: JJ = jj + ii[i];
2054: bd = 0;
2055: for (j=0; j<nz; j++) {
2056: if (*JJ >= i + rstart) break;
2057: JJ++;
2058: bd++;
2059: }
2060: d = 0;
2061: for (; j<nz; j++) {
2062: if (*JJ++ >= cend) break;
2063: d++;
2064: }
2065: d_nnz[i] = d;
2066: o_nnz[i] = nz - d - bd;
2067: nz = nz - bd;
2068: nz_max = PetscMax(nz_max,nz);
2069: }
2070: MatMPISBAIJSetPreallocation(B,bs,0,d_nnz,0,o_nnz);
2071: MatSetOption(B,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);
2072: PetscFree2(d_nnz,o_nnz);
2074: values = (PetscScalar*)V;
2075: if (!values) {
2076: PetscCalloc1(bs*bs*nz_max,&values);
2077: }
2078: for (i=0; i<m; i++) {
2079: PetscInt row = i + rstart;
2080: PetscInt ncols = ii[i+1] - ii[i];
2081: const PetscInt *icols = jj + ii[i];
2082: if (bs == 1 || !roworiented) { /* block ordering matches the non-nested layout of MatSetValues so we can insert entire rows */
2083: const PetscScalar *svals = values + (V ? (bs*bs*ii[i]) : 0);
2084: MatSetValuesBlocked_MPISBAIJ(B,1,&row,ncols,icols,svals,INSERT_VALUES);
2085: } else { /* block ordering does not match so we can only insert one block at a time. */
2086: PetscInt j;
2087: for (j=0; j<ncols; j++) {
2088: const PetscScalar *svals = values + (V ? (bs*bs*(ii[i]+j)) : 0);
2089: MatSetValuesBlocked_MPISBAIJ(B,1,&row,1,&icols[j],svals,INSERT_VALUES);
2090: }
2091: }
2092: }
2094: if (!V) { PetscFree(values); }
2095: nooffprocentries = B->nooffprocentries;
2096: B->nooffprocentries = PETSC_TRUE;
2097: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
2098: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
2099: B->nooffprocentries = nooffprocentries;
2101: MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
2102: return(0);
2103: }
2105: /*MC
2106: MATMPISBAIJ - MATMPISBAIJ = "mpisbaij" - A matrix type to be used for distributed symmetric sparse block matrices,
2107: based on block compressed sparse row format. Only the upper triangular portion of the "diagonal" portion of
2108: the matrix is stored.
2110: For complex numbers by default this matrix is symmetric, NOT Hermitian symmetric. To make it Hermitian symmetric you
2111: can call MatSetOption(Mat, MAT_HERMITIAN);
2113: Options Database Keys:
2114: . -mat_type mpisbaij - sets the matrix type to "mpisbaij" during a call to MatSetFromOptions()
2116: Notes:
2117: The number of rows in the matrix must be less than or equal to the number of columns. Similarly the number of rows in the
2118: diagonal portion of the matrix of each process has to less than or equal the number of columns.
2120: Level: beginner
2122: .seealso: MatCreateBAIJ(), MATSEQSBAIJ, MatType
2123: M*/
2125: PETSC_EXTERN PetscErrorCode MatCreate_MPISBAIJ(Mat B)
2126: {
2127: Mat_MPISBAIJ *b;
2129: PetscBool flg = PETSC_FALSE;
2132: PetscNewLog(B,&b);
2133: B->data = (void*)b;
2134: PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
2136: B->ops->destroy = MatDestroy_MPISBAIJ;
2137: B->ops->view = MatView_MPISBAIJ;
2138: B->assembled = PETSC_FALSE;
2139: B->insertmode = NOT_SET_VALUES;
2141: MPI_Comm_rank(PetscObjectComm((PetscObject)B),&b->rank);
2142: MPI_Comm_size(PetscObjectComm((PetscObject)B),&b->size);
2144: /* build local table of row and column ownerships */
2145: PetscMalloc1(b->size+2,&b->rangebs);
2147: /* build cache for off array entries formed */
2148: MatStashCreate_Private(PetscObjectComm((PetscObject)B),1,&B->stash);
2150: b->donotstash = PETSC_FALSE;
2151: b->colmap = NULL;
2152: b->garray = NULL;
2153: b->roworiented = PETSC_TRUE;
2155: /* stuff used in block assembly */
2156: b->barray = NULL;
2158: /* stuff used for matrix vector multiply */
2159: b->lvec = NULL;
2160: b->Mvctx = NULL;
2161: b->slvec0 = NULL;
2162: b->slvec0b = NULL;
2163: b->slvec1 = NULL;
2164: b->slvec1a = NULL;
2165: b->slvec1b = NULL;
2166: b->sMvctx = NULL;
2168: /* stuff for MatGetRow() */
2169: b->rowindices = NULL;
2170: b->rowvalues = NULL;
2171: b->getrowactive = PETSC_FALSE;
2173: /* hash table stuff */
2174: b->ht = NULL;
2175: b->hd = NULL;
2176: b->ht_size = 0;
2177: b->ht_flag = PETSC_FALSE;
2178: b->ht_fact = 0;
2179: b->ht_total_ct = 0;
2180: b->ht_insert_ct = 0;
2182: /* stuff for MatCreateSubMatrices_MPIBAIJ_local() */
2183: b->ijonly = PETSC_FALSE;
2185: b->in_loc = NULL;
2186: b->v_loc = NULL;
2187: b->n_loc = 0;
2189: PetscObjectComposeFunction((PetscObject)B,"MatStoreValues_C",MatStoreValues_MPISBAIJ);
2190: PetscObjectComposeFunction((PetscObject)B,"MatRetrieveValues_C",MatRetrieveValues_MPISBAIJ);
2191: PetscObjectComposeFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",MatMPISBAIJSetPreallocation_MPISBAIJ);
2192: PetscObjectComposeFunction((PetscObject)B,"MatMPISBAIJSetPreallocationCSR_C",MatMPISBAIJSetPreallocationCSR_MPISBAIJ);
2193: #if defined(PETSC_HAVE_ELEMENTAL)
2194: PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpisbaij_elemental_C",MatConvert_MPISBAIJ_Elemental);
2195: #endif
2196: #if defined(PETSC_HAVE_SCALAPACK)
2197: PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpisbaij_scalapack_C",MatConvert_SBAIJ_ScaLAPACK);
2198: #endif
2199: PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpisbaij_mpiaij_C",MatConvert_MPISBAIJ_Basic);
2200: PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpisbaij_mpibaij_C",MatConvert_MPISBAIJ_Basic);
2202: B->symmetric = PETSC_TRUE;
2203: B->structurally_symmetric = PETSC_TRUE;
2204: B->symmetric_set = PETSC_TRUE;
2205: B->structurally_symmetric_set = PETSC_TRUE;
2206: B->symmetric_eternal = PETSC_TRUE;
2207: #if defined(PETSC_USE_COMPLEX)
2208: B->hermitian = PETSC_FALSE;
2209: B->hermitian_set = PETSC_FALSE;
2210: #else
2211: B->hermitian = PETSC_TRUE;
2212: B->hermitian_set = PETSC_TRUE;
2213: #endif
2215: PetscObjectChangeTypeName((PetscObject)B,MATMPISBAIJ);
2216: PetscOptionsBegin(PetscObjectComm((PetscObject)B),NULL,"Options for loading MPISBAIJ matrix 1","Mat");
2217: PetscOptionsBool("-mat_use_hash_table","Use hash table to save memory in constructing matrix","MatSetOption",flg,&flg,NULL);
2218: if (flg) {
2219: PetscReal fact = 1.39;
2220: MatSetOption(B,MAT_USE_HASH_TABLE,PETSC_TRUE);
2221: PetscOptionsReal("-mat_use_hash_table","Use hash table factor","MatMPIBAIJSetHashTableFactor",fact,&fact,NULL);
2222: if (fact <= 1.0) fact = 1.39;
2223: MatMPIBAIJSetHashTableFactor(B,fact);
2224: PetscInfo1(B,"Hash table Factor used %5.2f\n",fact);
2225: }
2226: PetscOptionsEnd();
2227: return(0);
2228: }
2230: /*MC
2231: MATSBAIJ - MATSBAIJ = "sbaij" - A matrix type to be used for symmetric block sparse matrices.
2233: This matrix type is identical to MATSEQSBAIJ when constructed with a single process communicator,
2234: and MATMPISBAIJ otherwise.
2236: Options Database Keys:
2237: . -mat_type sbaij - sets the matrix type to "sbaij" during a call to MatSetFromOptions()
2239: Level: beginner
2241: .seealso: MatCreateSBAIJ, MATSEQSBAIJ, MATMPISBAIJ
2242: M*/
2244: /*@C
2245: MatMPISBAIJSetPreallocation - For good matrix assembly performance
2246: the user should preallocate the matrix storage by setting the parameters
2247: d_nz (or d_nnz) and o_nz (or o_nnz). By setting these parameters accurately,
2248: performance can be increased by more than a factor of 50.
2250: Collective on Mat
2252: Input Parameters:
2253: + B - the matrix
2254: . bs - size of block, the blocks are ALWAYS square. One can use MatSetBlockSizes() to set a different row and column blocksize but the row
2255: blocksize always defines the size of the blocks. The column blocksize sets the blocksize of the vectors obtained with MatCreateVecs()
2256: . d_nz - number of block nonzeros per block row in diagonal portion of local
2257: submatrix (same for all local rows)
2258: . d_nnz - array containing the number of block nonzeros in the various block rows
2259: in the upper triangular and diagonal part of the in diagonal portion of the local
2260: (possibly different for each block row) or NULL. If you plan to factor the matrix you must leave room
2261: for the diagonal entry and set a value even if it is zero.
2262: . o_nz - number of block nonzeros per block row in the off-diagonal portion of local
2263: submatrix (same for all local rows).
2264: - o_nnz - array containing the number of nonzeros in the various block rows of the
2265: off-diagonal portion of the local submatrix that is right of the diagonal
2266: (possibly different for each block row) or NULL.
2268: Options Database Keys:
2269: + -mat_no_unroll - uses code that does not unroll the loops in the
2270: block calculations (much slower)
2271: - -mat_block_size - size of the blocks to use
2273: Notes:
2275: If PETSC_DECIDE or PETSC_DETERMINE is used for a particular argument on one processor
2276: than it must be used on all processors that share the object for that argument.
2278: If the *_nnz parameter is given then the *_nz parameter is ignored
2280: Storage Information:
2281: For a square global matrix we define each processor's diagonal portion
2282: to be its local rows and the corresponding columns (a square submatrix);
2283: each processor's off-diagonal portion encompasses the remainder of the
2284: local matrix (a rectangular submatrix).
2286: The user can specify preallocated storage for the diagonal part of
2287: the local submatrix with either d_nz or d_nnz (not both). Set
2288: d_nz=PETSC_DEFAULT and d_nnz=NULL for PETSc to control dynamic
2289: memory allocation. Likewise, specify preallocated storage for the
2290: off-diagonal part of the local submatrix with o_nz or o_nnz (not both).
2292: You can call MatGetInfo() to get information on how effective the preallocation was;
2293: for example the fields mallocs,nz_allocated,nz_used,nz_unneeded;
2294: You can also run with the option -info and look for messages with the string
2295: malloc in them to see if additional memory allocation was needed.
2297: Consider a processor that owns rows 3, 4 and 5 of a parallel matrix. In
2298: the figure below we depict these three local rows and all columns (0-11).
2300: .vb
2301: 0 1 2 3 4 5 6 7 8 9 10 11
2302: --------------------------
2303: row 3 |. . . d d d o o o o o o
2304: row 4 |. . . d d d o o o o o o
2305: row 5 |. . . d d d o o o o o o
2306: --------------------------
2307: .ve
2309: Thus, any entries in the d locations are stored in the d (diagonal)
2310: submatrix, and any entries in the o locations are stored in the
2311: o (off-diagonal) submatrix. Note that the d matrix is stored in
2312: MatSeqSBAIJ format and the o submatrix in MATSEQBAIJ format.
2314: Now d_nz should indicate the number of block nonzeros per row in the upper triangular
2315: plus the diagonal part of the d matrix,
2316: and o_nz should indicate the number of block nonzeros per row in the o matrix
2318: In general, for PDE problems in which most nonzeros are near the diagonal,
2319: one expects d_nz >> o_nz. For large problems you MUST preallocate memory
2320: or you will get TERRIBLE performance; see the users' manual chapter on
2321: matrices.
2323: Level: intermediate
2325: .seealso: MatCreate(), MatCreateSeqSBAIJ(), MatSetValues(), MatCreateBAIJ(), PetscSplitOwnership()
2326: @*/
2327: PetscErrorCode MatMPISBAIJSetPreallocation(Mat B,PetscInt bs,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
2328: {
2335: PetscTryMethod(B,"MatMPISBAIJSetPreallocation_C",(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,bs,d_nz,d_nnz,o_nz,o_nnz));
2336: return(0);
2337: }
2339: /*@C
2340: MatCreateSBAIJ - Creates a sparse parallel matrix in symmetric block AIJ format
2341: (block compressed row). For good matrix assembly performance
2342: the user should preallocate the matrix storage by setting the parameters
2343: d_nz (or d_nnz) and o_nz (or o_nnz). By setting these parameters accurately,
2344: performance can be increased by more than a factor of 50.
2346: Collective
2348: Input Parameters:
2349: + comm - MPI communicator
2350: . bs - size of block, the blocks are ALWAYS square. One can use MatSetBlockSizes() to set a different row and column blocksize but the row
2351: blocksize always defines the size of the blocks. The column blocksize sets the blocksize of the vectors obtained with MatCreateVecs()
2352: . m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
2353: This value should be the same as the local size used in creating the
2354: y vector for the matrix-vector product y = Ax.
2355: . n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
2356: This value should be the same as the local size used in creating the
2357: x vector for the matrix-vector product y = Ax.
2358: . M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
2359: . N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
2360: . d_nz - number of block nonzeros per block row in diagonal portion of local
2361: submatrix (same for all local rows)
2362: . d_nnz - array containing the number of block nonzeros in the various block rows
2363: in the upper triangular portion of the in diagonal portion of the local
2364: (possibly different for each block block row) or NULL.
2365: If you plan to factor the matrix you must leave room for the diagonal entry and
2366: set its value even if it is zero.
2367: . o_nz - number of block nonzeros per block row in the off-diagonal portion of local
2368: submatrix (same for all local rows).
2369: - o_nnz - array containing the number of nonzeros in the various block rows of the
2370: off-diagonal portion of the local submatrix (possibly different for
2371: each block row) or NULL.
2373: Output Parameter:
2374: . A - the matrix
2376: Options Database Keys:
2377: + -mat_no_unroll - uses code that does not unroll the loops in the
2378: block calculations (much slower)
2379: . -mat_block_size - size of the blocks to use
2380: - -mat_mpi - use the parallel matrix data structures even on one processor
2381: (defaults to using SeqBAIJ format on one processor)
2383: It is recommended that one use the MatCreate(), MatSetType() and/or MatSetFromOptions(),
2384: MatXXXXSetPreallocation() paradigm instead of this routine directly.
2385: [MatXXXXSetPreallocation() is, for example, MatSeqAIJSetPreallocation]
2387: Notes:
2388: The number of rows and columns must be divisible by blocksize.
2389: This matrix type does not support complex Hermitian operation.
2391: The user MUST specify either the local or global matrix dimensions
2392: (possibly both).
2394: If PETSC_DECIDE or PETSC_DETERMINE is used for a particular argument on one processor
2395: than it must be used on all processors that share the object for that argument.
2397: If the *_nnz parameter is given then the *_nz parameter is ignored
2399: Storage Information:
2400: For a square global matrix we define each processor's diagonal portion
2401: to be its local rows and the corresponding columns (a square submatrix);
2402: each processor's off-diagonal portion encompasses the remainder of the
2403: local matrix (a rectangular submatrix).
2405: The user can specify preallocated storage for the diagonal part of
2406: the local submatrix with either d_nz or d_nnz (not both). Set
2407: d_nz=PETSC_DEFAULT and d_nnz=NULL for PETSc to control dynamic
2408: memory allocation. Likewise, specify preallocated storage for the
2409: off-diagonal part of the local submatrix with o_nz or o_nnz (not both).
2411: Consider a processor that owns rows 3, 4 and 5 of a parallel matrix. In
2412: the figure below we depict these three local rows and all columns (0-11).
2414: .vb
2415: 0 1 2 3 4 5 6 7 8 9 10 11
2416: --------------------------
2417: row 3 |. . . d d d o o o o o o
2418: row 4 |. . . d d d o o o o o o
2419: row 5 |. . . d d d o o o o o o
2420: --------------------------
2421: .ve
2423: Thus, any entries in the d locations are stored in the d (diagonal)
2424: submatrix, and any entries in the o locations are stored in the
2425: o (off-diagonal) submatrix. Note that the d matrix is stored in
2426: MatSeqSBAIJ format and the o submatrix in MATSEQBAIJ format.
2428: Now d_nz should indicate the number of block nonzeros per row in the upper triangular
2429: plus the diagonal part of the d matrix,
2430: and o_nz should indicate the number of block nonzeros per row in the o matrix.
2431: In general, for PDE problems in which most nonzeros are near the diagonal,
2432: one expects d_nz >> o_nz. For large problems you MUST preallocate memory
2433: or you will get TERRIBLE performance; see the users' manual chapter on
2434: matrices.
2436: Level: intermediate
2438: .seealso: MatCreate(), MatCreateSeqSBAIJ(), MatSetValues(), MatCreateBAIJ()
2439: @*/
2441: 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)
2442: {
2444: PetscMPIInt size;
2447: MatCreate(comm,A);
2448: MatSetSizes(*A,m,n,M,N);
2449: MPI_Comm_size(comm,&size);
2450: if (size > 1) {
2451: MatSetType(*A,MATMPISBAIJ);
2452: MatMPISBAIJSetPreallocation(*A,bs,d_nz,d_nnz,o_nz,o_nnz);
2453: } else {
2454: MatSetType(*A,MATSEQSBAIJ);
2455: MatSeqSBAIJSetPreallocation(*A,bs,d_nz,d_nnz);
2456: }
2457: return(0);
2458: }
2460: static PetscErrorCode MatDuplicate_MPISBAIJ(Mat matin,MatDuplicateOption cpvalues,Mat *newmat)
2461: {
2462: Mat mat;
2463: Mat_MPISBAIJ *a,*oldmat = (Mat_MPISBAIJ*)matin->data;
2465: PetscInt len=0,nt,bs=matin->rmap->bs,mbs=oldmat->mbs;
2466: PetscScalar *array;
2469: *newmat = NULL;
2471: MatCreate(PetscObjectComm((PetscObject)matin),&mat);
2472: MatSetSizes(mat,matin->rmap->n,matin->cmap->n,matin->rmap->N,matin->cmap->N);
2473: MatSetType(mat,((PetscObject)matin)->type_name);
2474: PetscLayoutReference(matin->rmap,&mat->rmap);
2475: PetscLayoutReference(matin->cmap,&mat->cmap);
2477: mat->factortype = matin->factortype;
2478: mat->preallocated = PETSC_TRUE;
2479: mat->assembled = PETSC_TRUE;
2480: mat->insertmode = NOT_SET_VALUES;
2482: a = (Mat_MPISBAIJ*)mat->data;
2483: a->bs2 = oldmat->bs2;
2484: a->mbs = oldmat->mbs;
2485: a->nbs = oldmat->nbs;
2486: a->Mbs = oldmat->Mbs;
2487: a->Nbs = oldmat->Nbs;
2489: a->size = oldmat->size;
2490: a->rank = oldmat->rank;
2491: a->donotstash = oldmat->donotstash;
2492: a->roworiented = oldmat->roworiented;
2493: a->rowindices = NULL;
2494: a->rowvalues = NULL;
2495: a->getrowactive = PETSC_FALSE;
2496: a->barray = NULL;
2497: a->rstartbs = oldmat->rstartbs;
2498: a->rendbs = oldmat->rendbs;
2499: a->cstartbs = oldmat->cstartbs;
2500: a->cendbs = oldmat->cendbs;
2502: /* hash table stuff */
2503: a->ht = NULL;
2504: a->hd = NULL;
2505: a->ht_size = 0;
2506: a->ht_flag = oldmat->ht_flag;
2507: a->ht_fact = oldmat->ht_fact;
2508: a->ht_total_ct = 0;
2509: a->ht_insert_ct = 0;
2511: PetscArraycpy(a->rangebs,oldmat->rangebs,a->size+2);
2512: if (oldmat->colmap) {
2513: #if defined(PETSC_USE_CTABLE)
2514: PetscTableCreateCopy(oldmat->colmap,&a->colmap);
2515: #else
2516: PetscMalloc1(a->Nbs,&a->colmap);
2517: PetscLogObjectMemory((PetscObject)mat,(a->Nbs)*sizeof(PetscInt));
2518: PetscArraycpy(a->colmap,oldmat->colmap,a->Nbs);
2519: #endif
2520: } else a->colmap = NULL;
2522: if (oldmat->garray && (len = ((Mat_SeqBAIJ*)(oldmat->B->data))->nbs)) {
2523: PetscMalloc1(len,&a->garray);
2524: PetscLogObjectMemory((PetscObject)mat,len*sizeof(PetscInt));
2525: PetscArraycpy(a->garray,oldmat->garray,len);
2526: } else a->garray = NULL;
2528: MatStashCreate_Private(PetscObjectComm((PetscObject)matin),matin->rmap->bs,&mat->bstash);
2529: VecDuplicate(oldmat->lvec,&a->lvec);
2530: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->lvec);
2531: VecScatterCopy(oldmat->Mvctx,&a->Mvctx);
2532: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->Mvctx);
2534: VecDuplicate(oldmat->slvec0,&a->slvec0);
2535: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec0);
2536: VecDuplicate(oldmat->slvec1,&a->slvec1);
2537: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec1);
2539: VecGetLocalSize(a->slvec1,&nt);
2540: VecGetArray(a->slvec1,&array);
2541: VecCreateSeqWithArray(PETSC_COMM_SELF,1,bs*mbs,array,&a->slvec1a);
2542: VecCreateSeqWithArray(PETSC_COMM_SELF,1,nt-bs*mbs,array+bs*mbs,&a->slvec1b);
2543: VecRestoreArray(a->slvec1,&array);
2544: VecGetArray(a->slvec0,&array);
2545: VecCreateSeqWithArray(PETSC_COMM_SELF,1,nt-bs*mbs,array+bs*mbs,&a->slvec0b);
2546: VecRestoreArray(a->slvec0,&array);
2547: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec0);
2548: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec1);
2549: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec0b);
2550: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec1a);
2551: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec1b);
2553: /* VecScatterCopy(oldmat->sMvctx,&a->sMvctx); - not written yet, replaced by the lazy trick: */
2554: PetscObjectReference((PetscObject)oldmat->sMvctx);
2555: a->sMvctx = oldmat->sMvctx;
2556: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->sMvctx);
2558: MatDuplicate(oldmat->A,cpvalues,&a->A);
2559: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->A);
2560: MatDuplicate(oldmat->B,cpvalues,&a->B);
2561: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->B);
2562: PetscFunctionListDuplicate(((PetscObject)matin)->qlist,&((PetscObject)mat)->qlist);
2563: *newmat = mat;
2564: return(0);
2565: }
2567: /* Used for both MPIBAIJ and MPISBAIJ matrices */
2568: #define MatLoad_MPISBAIJ_Binary MatLoad_MPIBAIJ_Binary
2570: PetscErrorCode MatLoad_MPISBAIJ(Mat mat,PetscViewer viewer)
2571: {
2573: PetscBool isbinary;
2576: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
2577: if (!isbinary) SETERRQ2(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer type %s not yet supported for reading %s matrices",((PetscObject)viewer)->type_name,((PetscObject)mat)->type_name);
2578: MatLoad_MPISBAIJ_Binary(mat,viewer);
2579: return(0);
2580: }
2582: /*XXXXX@
2583: MatMPISBAIJSetHashTableFactor - Sets the factor required to compute the size of the HashTable.
2585: Input Parameters:
2586: . mat - the matrix
2587: . fact - factor
2589: Not Collective on Mat, each process can have a different hash factor
2591: Level: advanced
2593: Notes:
2594: This can also be set by the command line option: -mat_use_hash_table fact
2596: .seealso: MatSetOption()
2597: @XXXXX*/
2599: PetscErrorCode MatGetRowMaxAbs_MPISBAIJ(Mat A,Vec v,PetscInt idx[])
2600: {
2601: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
2602: Mat_SeqBAIJ *b = (Mat_SeqBAIJ*)(a->B)->data;
2603: PetscReal atmp;
2604: PetscReal *work,*svalues,*rvalues;
2606: PetscInt i,bs,mbs,*bi,*bj,brow,j,ncols,krow,kcol,col,row,Mbs,bcol;
2607: PetscMPIInt rank,size;
2608: PetscInt *rowners_bs,dest,count,source;
2609: PetscScalar *va;
2610: MatScalar *ba;
2611: MPI_Status stat;
2614: if (idx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Send email to petsc-maint@mcs.anl.gov");
2615: MatGetRowMaxAbs(a->A,v,NULL);
2616: VecGetArray(v,&va);
2618: MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);
2619: MPI_Comm_rank(PetscObjectComm((PetscObject)A),&rank);
2621: bs = A->rmap->bs;
2622: mbs = a->mbs;
2623: Mbs = a->Mbs;
2624: ba = b->a;
2625: bi = b->i;
2626: bj = b->j;
2628: /* find ownerships */
2629: rowners_bs = A->rmap->range;
2631: /* each proc creates an array to be distributed */
2632: PetscCalloc1(bs*Mbs,&work);
2634: /* row_max for B */
2635: if (rank != size-1) {
2636: for (i=0; i<mbs; i++) {
2637: ncols = bi[1] - bi[0]; bi++;
2638: brow = bs*i;
2639: for (j=0; j<ncols; j++) {
2640: bcol = bs*(*bj);
2641: for (kcol=0; kcol<bs; kcol++) {
2642: col = bcol + kcol; /* local col index */
2643: col += rowners_bs[rank+1]; /* global col index */
2644: for (krow=0; krow<bs; krow++) {
2645: atmp = PetscAbsScalar(*ba); ba++;
2646: row = brow + krow; /* local row index */
2647: if (PetscRealPart(va[row]) < atmp) va[row] = atmp;
2648: if (work[col] < atmp) work[col] = atmp;
2649: }
2650: }
2651: bj++;
2652: }
2653: }
2655: /* send values to its owners */
2656: for (dest=rank+1; dest<size; dest++) {
2657: svalues = work + rowners_bs[dest];
2658: count = rowners_bs[dest+1]-rowners_bs[dest];
2659: MPI_Send(svalues,count,MPIU_REAL,dest,rank,PetscObjectComm((PetscObject)A));
2660: }
2661: }
2663: /* receive values */
2664: if (rank) {
2665: rvalues = work;
2666: count = rowners_bs[rank+1]-rowners_bs[rank];
2667: for (source=0; source<rank; source++) {
2668: MPI_Recv(rvalues,count,MPIU_REAL,MPI_ANY_SOURCE,MPI_ANY_TAG,PetscObjectComm((PetscObject)A),&stat);
2669: /* process values */
2670: for (i=0; i<count; i++) {
2671: if (PetscRealPart(va[i]) < rvalues[i]) va[i] = rvalues[i];
2672: }
2673: }
2674: }
2676: VecRestoreArray(v,&va);
2677: PetscFree(work);
2678: return(0);
2679: }
2681: PetscErrorCode MatSOR_MPISBAIJ(Mat matin,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
2682: {
2683: Mat_MPISBAIJ *mat = (Mat_MPISBAIJ*)matin->data;
2684: PetscErrorCode ierr;
2685: PetscInt mbs=mat->mbs,bs=matin->rmap->bs;
2686: PetscScalar *x,*ptr,*from;
2687: Vec bb1;
2688: const PetscScalar *b;
2691: 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);
2692: if (bs > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented");
2694: if (flag == SOR_APPLY_UPPER) {
2695: (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx);
2696: return(0);
2697: }
2699: if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP) {
2700: if (flag & SOR_ZERO_INITIAL_GUESS) {
2701: (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,lits,xx);
2702: its--;
2703: }
2705: VecDuplicate(bb,&bb1);
2706: while (its--) {
2708: /* lower triangular part: slvec0b = - B^T*xx */
2709: (*mat->B->ops->multtranspose)(mat->B,xx,mat->slvec0b);
2711: /* copy xx into slvec0a */
2712: VecGetArray(mat->slvec0,&ptr);
2713: VecGetArray(xx,&x);
2714: PetscArraycpy(ptr,x,bs*mbs);
2715: VecRestoreArray(mat->slvec0,&ptr);
2717: VecScale(mat->slvec0,-1.0);
2719: /* copy bb into slvec1a */
2720: VecGetArray(mat->slvec1,&ptr);
2721: VecGetArrayRead(bb,&b);
2722: PetscArraycpy(ptr,b,bs*mbs);
2723: VecRestoreArray(mat->slvec1,&ptr);
2725: /* set slvec1b = 0 */
2726: VecSet(mat->slvec1b,0.0);
2728: VecScatterBegin(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD);
2729: VecRestoreArray(xx,&x);
2730: VecRestoreArrayRead(bb,&b);
2731: VecScatterEnd(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD);
2733: /* upper triangular part: bb1 = bb1 - B*x */
2734: (*mat->B->ops->multadd)(mat->B,mat->slvec1b,mat->slvec1a,bb1);
2736: /* local diagonal sweep */
2737: (*mat->A->ops->sor)(mat->A,bb1,omega,SOR_SYMMETRIC_SWEEP,fshift,lits,lits,xx);
2738: }
2739: VecDestroy(&bb1);
2740: } else if ((flag & SOR_LOCAL_FORWARD_SWEEP) && (its == 1) && (flag & SOR_ZERO_INITIAL_GUESS)) {
2741: (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx);
2742: } else if ((flag & SOR_LOCAL_BACKWARD_SWEEP) && (its == 1) && (flag & SOR_ZERO_INITIAL_GUESS)) {
2743: (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx);
2744: } else if (flag & SOR_EISENSTAT) {
2745: Vec xx1;
2746: PetscBool hasop;
2747: const PetscScalar *diag;
2748: PetscScalar *sl,scale = (omega - 2.0)/omega;
2749: PetscInt i,n;
2751: if (!mat->xx1) {
2752: VecDuplicate(bb,&mat->xx1);
2753: VecDuplicate(bb,&mat->bb1);
2754: }
2755: xx1 = mat->xx1;
2756: bb1 = mat->bb1;
2758: (*mat->A->ops->sor)(mat->A,bb,omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_BACKWARD_SWEEP),fshift,lits,1,xx);
2760: if (!mat->diag) {
2761: /* this is wrong for same matrix with new nonzero values */
2762: MatCreateVecs(matin,&mat->diag,NULL);
2763: MatGetDiagonal(matin,mat->diag);
2764: }
2765: MatHasOperation(matin,MATOP_MULT_DIAGONAL_BLOCK,&hasop);
2767: if (hasop) {
2768: MatMultDiagonalBlock(matin,xx,bb1);
2769: VecAYPX(mat->slvec1a,scale,bb);
2770: } else {
2771: /*
2772: These two lines are replaced by code that may be a bit faster for a good compiler
2773: VecPointwiseMult(mat->slvec1a,mat->diag,xx);
2774: VecAYPX(mat->slvec1a,scale,bb);
2775: */
2776: VecGetArray(mat->slvec1a,&sl);
2777: VecGetArrayRead(mat->diag,&diag);
2778: VecGetArrayRead(bb,&b);
2779: VecGetArray(xx,&x);
2780: VecGetLocalSize(xx,&n);
2781: if (omega == 1.0) {
2782: for (i=0; i<n; i++) sl[i] = b[i] - diag[i]*x[i];
2783: PetscLogFlops(2.0*n);
2784: } else {
2785: for (i=0; i<n; i++) sl[i] = b[i] + scale*diag[i]*x[i];
2786: PetscLogFlops(3.0*n);
2787: }
2788: VecRestoreArray(mat->slvec1a,&sl);
2789: VecRestoreArrayRead(mat->diag,&diag);
2790: VecRestoreArrayRead(bb,&b);
2791: VecRestoreArray(xx,&x);
2792: }
2794: /* multiply off-diagonal portion of matrix */
2795: VecSet(mat->slvec1b,0.0);
2796: (*mat->B->ops->multtranspose)(mat->B,xx,mat->slvec0b);
2797: VecGetArray(mat->slvec0,&from);
2798: VecGetArray(xx,&x);
2799: PetscArraycpy(from,x,bs*mbs);
2800: VecRestoreArray(mat->slvec0,&from);
2801: VecRestoreArray(xx,&x);
2802: VecScatterBegin(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD);
2803: VecScatterEnd(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD);
2804: (*mat->B->ops->multadd)(mat->B,mat->slvec1b,mat->slvec1a,mat->slvec1a);
2806: /* local sweep */
2807: (*mat->A->ops->sor)(mat->A,mat->slvec1a,omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_FORWARD_SWEEP),fshift,lits,1,xx1);
2808: VecAXPY(xx,1.0,xx1);
2809: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatSORType is not supported for SBAIJ matrix format");
2810: return(0);
2811: }
2813: /*@
2814: MatCreateMPISBAIJWithArrays - creates a MPI SBAIJ matrix using arrays that contain in standard
2815: CSR format the local rows.
2817: Collective
2819: Input Parameters:
2820: + comm - MPI communicator
2821: . bs - the block size, only a block size of 1 is supported
2822: . m - number of local rows (Cannot be PETSC_DECIDE)
2823: . n - This value should be the same as the local size used in creating the
2824: x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have
2825: calculated if N is given) For square matrices n is almost always m.
2826: . M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
2827: . N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
2828: . i - row indices; that is i[0] = 0, i[row] = i[row-1] + number of block elements in that row block row of the matrix
2829: . j - column indices
2830: - a - matrix values
2832: Output Parameter:
2833: . mat - the matrix
2835: Level: intermediate
2837: Notes:
2838: The i, j, and a arrays ARE copied by this routine into the internal format used by PETSc;
2839: thus you CANNOT change the matrix entries by changing the values of a[] after you have
2840: called this routine. Use MatCreateMPIAIJWithSplitArrays() to avoid needing to copy the arrays.
2842: The i and j indices are 0 based, and i indices are indices corresponding to the local j array.
2844: .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatMPIAIJSetPreallocation(), MatMPIAIJSetPreallocationCSR(),
2845: MPIAIJ, MatCreateAIJ(), MatCreateMPIAIJWithSplitArrays()
2846: @*/
2847: 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)
2848: {
2852: if (i[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"i (row indices) must start with 0");
2853: if (m < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"local number of rows (m) cannot be PETSC_DECIDE, or negative");
2854: MatCreate(comm,mat);
2855: MatSetSizes(*mat,m,n,M,N);
2856: MatSetType(*mat,MATMPISBAIJ);
2857: MatMPISBAIJSetPreallocationCSR(*mat,bs,i,j,a);
2858: return(0);
2859: }
2861: /*@C
2862: MatMPISBAIJSetPreallocationCSR - Creates a sparse parallel matrix in SBAIJ format using the given nonzero structure and (optional) numerical values
2864: Collective
2866: Input Parameters:
2867: + B - the matrix
2868: . bs - the block size
2869: . i - the indices into j for the start of each local row (starts with zero)
2870: . j - the column indices for each local row (starts with zero) these must be sorted for each row
2871: - v - optional values in the matrix
2873: Level: advanced
2875: Notes:
2876: Though this routine has Preallocation() in the name it also sets the exact nonzero locations of the matrix entries
2877: and usually the numerical values as well
2879: Any entries below the diagonal are ignored
2881: .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatMPIBAIJSetPreallocation(), MatCreateAIJ(), MPIAIJ
2882: @*/
2883: PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat B,PetscInt bs,const PetscInt i[],const PetscInt j[], const PetscScalar v[])
2884: {
2888: PetscTryMethod(B,"MatMPISBAIJSetPreallocationCSR_C",(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]),(B,bs,i,j,v));
2889: return(0);
2890: }
2892: PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPISBAIJ(MPI_Comm comm,Mat inmat,PetscInt n,MatReuse scall,Mat *outmat)
2893: {
2895: PetscInt m,N,i,rstart,nnz,Ii,bs,cbs;
2896: PetscInt *indx;
2897: PetscScalar *values;
2900: MatGetSize(inmat,&m,&N);
2901: if (scall == MAT_INITIAL_MATRIX) { /* symbolic phase */
2902: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inmat->data;
2903: PetscInt *dnz,*onz,mbs,Nbs,nbs;
2904: PetscInt *bindx,rmax=a->rmax,j;
2905: PetscMPIInt rank,size;
2907: MatGetBlockSizes(inmat,&bs,&cbs);
2908: mbs = m/bs; Nbs = N/cbs;
2909: if (n == PETSC_DECIDE) {
2910: PetscSplitOwnershipBlock(comm,cbs,&n,&N);
2911: }
2912: nbs = n/cbs;
2914: PetscMalloc1(rmax,&bindx);
2915: MatPreallocateInitialize(comm,mbs,nbs,dnz,onz); /* inline function, output __end and __rstart are used below */
2917: MPI_Comm_rank(comm,&rank);
2918: MPI_Comm_rank(comm,&size);
2919: if (rank == size-1) {
2920: /* Check sum(nbs) = Nbs */
2921: if (__end != Nbs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Sum of local block columns %D != global block columns %D",__end,Nbs);
2922: }
2924: rstart = __rstart; /* block rstart of *outmat; see inline function MatPreallocateInitialize */
2925: MatSetOption(inmat,MAT_GETROW_UPPERTRIANGULAR,PETSC_TRUE);
2926: for (i=0; i<mbs; i++) {
2927: MatGetRow_SeqSBAIJ(inmat,i*bs,&nnz,&indx,NULL); /* non-blocked nnz and indx */
2928: nnz = nnz/bs;
2929: for (j=0; j<nnz; j++) bindx[j] = indx[j*bs]/bs;
2930: MatPreallocateSet(i+rstart,nnz,bindx,dnz,onz);
2931: MatRestoreRow_SeqSBAIJ(inmat,i*bs,&nnz,&indx,NULL);
2932: }
2933: MatSetOption(inmat,MAT_GETROW_UPPERTRIANGULAR,PETSC_FALSE);
2934: PetscFree(bindx);
2936: MatCreate(comm,outmat);
2937: MatSetSizes(*outmat,m,n,PETSC_DETERMINE,PETSC_DETERMINE);
2938: MatSetBlockSizes(*outmat,bs,cbs);
2939: MatSetType(*outmat,MATSBAIJ);
2940: MatSeqSBAIJSetPreallocation(*outmat,bs,0,dnz);
2941: MatMPISBAIJSetPreallocation(*outmat,bs,0,dnz,0,onz);
2942: MatPreallocateFinalize(dnz,onz);
2943: }
2945: /* numeric phase */
2946: MatGetBlockSizes(inmat,&bs,&cbs);
2947: MatGetOwnershipRange(*outmat,&rstart,NULL);
2949: MatSetOption(inmat,MAT_GETROW_UPPERTRIANGULAR,PETSC_TRUE);
2950: for (i=0; i<m; i++) {
2951: MatGetRow_SeqSBAIJ(inmat,i,&nnz,&indx,&values);
2952: Ii = i + rstart;
2953: MatSetValues(*outmat,1,&Ii,nnz,indx,values,INSERT_VALUES);
2954: MatRestoreRow_SeqSBAIJ(inmat,i,&nnz,&indx,&values);
2955: }
2956: MatSetOption(inmat,MAT_GETROW_UPPERTRIANGULAR,PETSC_FALSE);
2957: MatAssemblyBegin(*outmat,MAT_FINAL_ASSEMBLY);
2958: MatAssemblyEnd(*outmat,MAT_FINAL_ASSEMBLY);
2959: return(0);
2960: }