Actual source code: mpisbaij.c
1: #define PETSCMAT_DLL
3: #include src/mat/impls/baij/mpi/mpibaij.h
4: #include mpisbaij.h
5: #include src/mat/impls/sbaij/seq/sbaij.h
7: EXTERN PetscErrorCode MatSetUpMultiply_MPISBAIJ(Mat);
8: EXTERN PetscErrorCode MatSetUpMultiply_MPISBAIJ_2comm(Mat);
9: EXTERN PetscErrorCode DisAssemble_MPISBAIJ(Mat);
10: EXTERN PetscErrorCode MatIncreaseOverlap_MPISBAIJ(Mat,PetscInt,IS[],PetscInt);
11: EXTERN PetscErrorCode MatGetValues_SeqSBAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar []);
12: EXTERN PetscErrorCode MatGetValues_SeqBAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar []);
13: EXTERN PetscErrorCode MatSetValues_SeqSBAIJ(Mat,PetscInt,const PetscInt [],PetscInt,const PetscInt [],const PetscScalar [],InsertMode);
14: EXTERN PetscErrorCode MatSetValuesBlocked_SeqSBAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
15: EXTERN PetscErrorCode MatSetValuesBlocked_SeqBAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
16: EXTERN PetscErrorCode MatGetRow_SeqSBAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
17: EXTERN PetscErrorCode MatRestoreRow_SeqSBAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
18: EXTERN PetscErrorCode MatZeroRows_SeqSBAIJ(Mat,IS,PetscScalar*);
19: EXTERN PetscErrorCode MatZeroRows_SeqBAIJ(Mat,IS,PetscScalar *);
20: EXTERN PetscErrorCode MatGetRowMax_MPISBAIJ(Mat,Vec);
21: EXTERN PetscErrorCode MatRelax_MPISBAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
23: /* UGLY, ugly, ugly
24: When MatScalar == PetscScalar the function MatSetValuesBlocked_MPIBAIJ_MatScalar() does
25: not exist. Otherwise ..._MatScalar() takes matrix elements in single precision and
26: inserts them into the single precision data structure. The function MatSetValuesBlocked_MPIBAIJ()
27: converts the entries into single precision and then calls ..._MatScalar() to put them
28: into the single precision data structures.
29: */
30: #if defined(PETSC_USE_MAT_SINGLE)
31: EXTERN PetscErrorCode MatSetValuesBlocked_SeqSBAIJ_MatScalar(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const MatScalar[],InsertMode);
32: EXTERN PetscErrorCode MatSetValues_MPISBAIJ_MatScalar(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const MatScalar[],InsertMode);
33: EXTERN PetscErrorCode MatSetValuesBlocked_MPISBAIJ_MatScalar(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const MatScalar[],InsertMode);
34: EXTERN PetscErrorCode MatSetValues_MPISBAIJ_HT_MatScalar(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const MatScalar[],InsertMode);
35: EXTERN PetscErrorCode MatSetValuesBlocked_MPISBAIJ_HT_MatScalar(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const MatScalar[],InsertMode);
36: #else
37: #define MatSetValuesBlocked_SeqSBAIJ_MatScalar MatSetValuesBlocked_SeqSBAIJ
38: #define MatSetValues_MPISBAIJ_MatScalar MatSetValues_MPISBAIJ
39: #define MatSetValuesBlocked_MPISBAIJ_MatScalar MatSetValuesBlocked_MPISBAIJ
40: #define MatSetValues_MPISBAIJ_HT_MatScalar MatSetValues_MPISBAIJ_HT
41: #define MatSetValuesBlocked_MPISBAIJ_HT_MatScalar MatSetValuesBlocked_MPISBAIJ_HT
42: #endif
47: PetscErrorCode MatStoreValues_MPISBAIJ(Mat mat)
48: {
49: Mat_MPISBAIJ *aij = (Mat_MPISBAIJ *)mat->data;
53: MatStoreValues(aij->A);
54: MatStoreValues(aij->B);
55: return(0);
56: }
62: PetscErrorCode MatRetrieveValues_MPISBAIJ(Mat mat)
63: {
64: Mat_MPISBAIJ *aij = (Mat_MPISBAIJ *)mat->data;
68: MatRetrieveValues(aij->A);
69: MatRetrieveValues(aij->B);
70: return(0);
71: }
75: #define CHUNKSIZE 10
77: #define MatSetValues_SeqSBAIJ_A_Private(row,col,value,addv) \
78: { \
79: \
80: brow = row/bs; \
81: rp = aj + ai[brow]; ap = aa + bs2*ai[brow]; \
82: rmax = aimax[brow]; nrow = ailen[brow]; \
83: bcol = col/bs; \
84: ridx = row % bs; cidx = col % bs; \
85: low = 0; high = nrow; \
86: while (high-low > 3) { \
87: t = (low+high)/2; \
88: if (rp[t] > bcol) high = t; \
89: else low = t; \
90: } \
91: for (_i=low; _i<high; _i++) { \
92: if (rp[_i] > bcol) break; \
93: if (rp[_i] == bcol) { \
94: bap = ap + bs2*_i + bs*cidx + ridx; \
95: if (addv == ADD_VALUES) *bap += value; \
96: else *bap = value; \
97: goto a_noinsert; \
98: } \
99: } \
100: if (a->nonew == 1) goto a_noinsert; \
101: if (a->nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) into matrix", row, col); \
102: MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,aimax,a->nonew,MatScalar); \
103: N = nrow++ - 1; \
104: /* shift up all the later entries in this row */ \
105: for (ii=N; ii>=_i; ii--) { \
106: rp[ii+1] = rp[ii]; \
107: PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar)); \
108: } \
109: if (N>=_i) { PetscMemzero(ap+bs2*_i,bs2*sizeof(MatScalar)); } \
110: rp[_i] = bcol; \
111: ap[bs2*_i + bs*cidx + ridx] = value; \
112: a_noinsert:; \
113: ailen[brow] = nrow; \
114: }
115: #ifndef MatSetValues_SeqBAIJ_B_Private
116: #define MatSetValues_SeqSBAIJ_B_Private(row,col,value,addv) \
117: { \
118: brow = row/bs; \
119: rp = bj + bi[brow]; ap = ba + bs2*bi[brow]; \
120: rmax = bimax[brow]; nrow = bilen[brow]; \
121: bcol = col/bs; \
122: ridx = row % bs; cidx = col % bs; \
123: low = 0; high = nrow; \
124: while (high-low > 3) { \
125: t = (low+high)/2; \
126: if (rp[t] > bcol) high = t; \
127: else low = t; \
128: } \
129: for (_i=low; _i<high; _i++) { \
130: if (rp[_i] > bcol) break; \
131: if (rp[_i] == bcol) { \
132: bap = ap + bs2*_i + bs*cidx + ridx; \
133: if (addv == ADD_VALUES) *bap += value; \
134: else *bap = value; \
135: goto b_noinsert; \
136: } \
137: } \
138: if (b->nonew == 1) goto b_noinsert; \
139: if (b->nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) into matrix", row, col); \
140: MatSeqXAIJReallocateAIJ(B,b->mbs,bs2,nrow,brow,bcol,rmax,ba,bi,bj,rp,ap,bimax,b->nonew,MatScalar); \
141: N = nrow++ - 1; \
142: /* shift up all the later entries in this row */ \
143: for (ii=N; ii>=_i; ii--) { \
144: rp[ii+1] = rp[ii]; \
145: PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar)); \
146: } \
147: if (N>=_i) { PetscMemzero(ap+bs2*_i,bs2*sizeof(MatScalar));} \
148: rp[_i] = bcol; \
149: ap[bs2*_i + bs*cidx + ridx] = value; \
150: b_noinsert:; \
151: bilen[brow] = nrow; \
152: }
153: #endif
155: #if defined(PETSC_USE_MAT_SINGLE)
158: PetscErrorCode MatSetValues_MPISBAIJ(Mat mat,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode addv)
159: {
160: Mat_MPISBAIJ *b = (Mat_MPISBAIJ*)mat->data;
162: PetscInt i,N = m*n;
163: MatScalar *vsingle;
166: if (N > b->setvalueslen) {
167: PetscFree(b->setvaluescopy);
168: PetscMalloc(N*sizeof(MatScalar),&b->setvaluescopy);
169: b->setvalueslen = N;
170: }
171: vsingle = b->setvaluescopy;
173: for (i=0; i<N; i++) {
174: vsingle[i] = v[i];
175: }
176: MatSetValues_MPISBAIJ_MatScalar(mat,m,im,n,in,vsingle,addv);
177: return(0);
178: }
182: PetscErrorCode MatSetValuesBlocked_MPISBAIJ(Mat mat,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode addv)
183: {
184: Mat_MPIBAIJ *b = (Mat_MPIBAIJ*)mat->data;
186: PetscInt i,N = m*n*b->bs2;
187: MatScalar *vsingle;
190: if (N > b->setvalueslen) {
191: PetscFree(b->setvaluescopy);
192: PetscMalloc(N*sizeof(MatScalar),&b->setvaluescopy);
193: b->setvalueslen = N;
194: }
195: vsingle = b->setvaluescopy;
196: for (i=0; i<N; i++) {
197: vsingle[i] = v[i];
198: }
199: MatSetValuesBlocked_MPISBAIJ_MatScalar(mat,m,im,n,in,vsingle,addv);
200: return(0);
201: }
202: #endif
204: /* Only add/insert a(i,j) with i<=j (blocks).
205: Any a(i,j) with i>j input by user is ingored.
206: */
209: PetscErrorCode MatSetValues_MPISBAIJ_MatScalar(Mat mat,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const MatScalar v[],InsertMode addv)
210: {
211: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
212: MatScalar value;
213: PetscTruth roworiented = baij->roworiented;
215: PetscInt i,j,row,col;
216: PetscInt rstart_orig=mat->rmap.rstart;
217: PetscInt rend_orig=mat->rmap.rend,cstart_orig=mat->cmap.rstart;
218: PetscInt cend_orig=mat->cmap.rend,bs=mat->rmap.bs;
220: /* Some Variables required in the macro */
221: Mat A = baij->A;
222: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)(A)->data;
223: PetscInt *aimax=a->imax,*ai=a->i,*ailen=a->ilen,*aj=a->j;
224: MatScalar *aa=a->a;
226: Mat B = baij->B;
227: Mat_SeqBAIJ *b = (Mat_SeqBAIJ*)(B)->data;
228: PetscInt *bimax=b->imax,*bi=b->i,*bilen=b->ilen,*bj=b->j;
229: MatScalar *ba=b->a;
231: PetscInt *rp,ii,nrow,_i,rmax,N,brow,bcol;
232: PetscInt low,high,t,ridx,cidx,bs2=a->bs2;
233: MatScalar *ap,*bap;
235: /* for stash */
236: PetscInt n_loc, *in_loc = PETSC_NULL;
237: MatScalar *v_loc = PETSC_NULL;
241: if (!baij->donotstash){
242: if (n > baij->n_loc) {
243: PetscFree(baij->in_loc);
244: PetscFree(baij->v_loc);
245: PetscMalloc(n*sizeof(PetscInt),&baij->in_loc);
246: PetscMalloc(n*sizeof(MatScalar),&baij->v_loc);
247: baij->n_loc = n;
248: }
249: in_loc = baij->in_loc;
250: v_loc = baij->v_loc;
251: }
253: for (i=0; i<m; i++) {
254: if (im[i] < 0) continue;
255: #if defined(PETSC_USE_DEBUG)
256: if (im[i] >= mat->rmap.N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",im[i],mat->rmap.N-1);
257: #endif
258: if (im[i] >= rstart_orig && im[i] < rend_orig) { /* this processor entry */
259: row = im[i] - rstart_orig; /* local row index */
260: for (j=0; j<n; j++) {
261: if (im[i]/bs > in[j]/bs){
262: if (a->ignore_ltriangular){
263: continue; /* ignore lower triangular blocks */
264: } else {
265: SETERRQ(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)");
266: }
267: }
268: if (in[j] >= cstart_orig && in[j] < cend_orig){ /* diag entry (A) */
269: col = in[j] - cstart_orig; /* local col index */
270: brow = row/bs; bcol = col/bs;
271: if (brow > bcol) continue; /* ignore lower triangular blocks of A */
272: if (roworiented) value = v[i*n+j]; else value = v[i+j*m];
273: MatSetValues_SeqSBAIJ_A_Private(row,col,value,addv);
274: /* MatSetValues_SeqBAIJ(baij->A,1,&row,1,&col,&value,addv); */
275: } else if (in[j] < 0) continue;
276: #if defined(PETSC_USE_DEBUG)
277: else if (in[j] >= mat->cmap.N) {SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[j],mat->cmap.N-1);}
278: #endif
279: else { /* off-diag entry (B) */
280: if (mat->was_assembled) {
281: if (!baij->colmap) {
282: CreateColmap_MPIBAIJ_Private(mat);
283: }
284: #if defined (PETSC_USE_CTABLE)
285: PetscTableFind(baij->colmap,in[j]/bs + 1,&col);
286: col = col - 1;
287: #else
288: col = baij->colmap[in[j]/bs] - 1;
289: #endif
290: if (col < 0 && !((Mat_SeqSBAIJ*)(baij->A->data))->nonew) {
291: DisAssemble_MPISBAIJ(mat);
292: col = in[j];
293: /* Reinitialize the variables required by MatSetValues_SeqBAIJ_B_Private() */
294: B = baij->B;
295: b = (Mat_SeqBAIJ*)(B)->data;
296: bimax=b->imax;bi=b->i;bilen=b->ilen;bj=b->j;
297: ba=b->a;
298: } else col += in[j]%bs;
299: } else col = in[j];
300: if (roworiented) value = v[i*n+j]; else value = v[i+j*m];
301: MatSetValues_SeqSBAIJ_B_Private(row,col,value,addv);
302: /* MatSetValues_SeqBAIJ(baij->B,1,&row,1,&col,&value,addv); */
303: }
304: }
305: } else { /* off processor entry */
306: if (!baij->donotstash) {
307: n_loc = 0;
308: for (j=0; j<n; j++){
309: if (im[i]/bs > in[j]/bs) continue; /* ignore lower triangular blocks */
310: in_loc[n_loc] = in[j];
311: if (roworiented) {
312: v_loc[n_loc] = v[i*n+j];
313: } else {
314: v_loc[n_loc] = v[j*m+i];
315: }
316: n_loc++;
317: }
318: MatStashValuesRow_Private(&mat->stash,im[i],n_loc,in_loc,v_loc);
319: }
320: }
321: }
322: return(0);
323: }
327: PetscErrorCode MatSetValuesBlocked_MPISBAIJ_MatScalar(Mat mat,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const MatScalar v[],InsertMode addv)
328: {
329: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
330: const MatScalar *value;
331: MatScalar *barray=baij->barray;
332: PetscTruth roworiented = baij->roworiented;
333: PetscErrorCode ierr;
334: PetscInt i,j,ii,jj,row,col,rstart=baij->rstartbs;
335: PetscInt rend=baij->rendbs,cstart=baij->rstartbs,stepval;
336: PetscInt cend=baij->rendbs,bs=mat->rmap.bs,bs2=baij->bs2;
339: if(!barray) {
340: PetscMalloc(bs2*sizeof(MatScalar),&barray);
341: baij->barray = barray;
342: }
344: if (roworiented) {
345: stepval = (n-1)*bs;
346: } else {
347: stepval = (m-1)*bs;
348: }
349: for (i=0; i<m; i++) {
350: if (im[i] < 0) continue;
351: #if defined(PETSC_USE_DEBUG)
352: if (im[i] >= baij->Mbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large, row %D max %D",im[i],baij->Mbs-1);
353: #endif
354: if (im[i] >= rstart && im[i] < rend) {
355: row = im[i] - rstart;
356: for (j=0; j<n; j++) {
357: /* If NumCol = 1 then a copy is not required */
358: if ((roworiented) && (n == 1)) {
359: barray = (MatScalar*) v + i*bs2;
360: } else if((!roworiented) && (m == 1)) {
361: barray = (MatScalar*) v + j*bs2;
362: } else { /* Here a copy is required */
363: if (roworiented) {
364: value = v + i*(stepval+bs)*bs + j*bs;
365: } else {
366: value = v + j*(stepval+bs)*bs + i*bs;
367: }
368: for (ii=0; ii<bs; ii++,value+=stepval) {
369: for (jj=0; jj<bs; jj++) {
370: *barray++ = *value++;
371: }
372: }
373: barray -=bs2;
374: }
375:
376: if (in[j] >= cstart && in[j] < cend){
377: col = in[j] - cstart;
378: MatSetValuesBlocked_SeqSBAIJ(baij->A,1,&row,1,&col,barray,addv);
379: }
380: else if (in[j] < 0) continue;
381: #if defined(PETSC_USE_DEBUG)
382: else if (in[j] >= baij->Nbs) {SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large, col %D max %D",in[j],baij->Nbs-1);}
383: #endif
384: else {
385: if (mat->was_assembled) {
386: if (!baij->colmap) {
387: CreateColmap_MPIBAIJ_Private(mat);
388: }
390: #if defined(PETSC_USE_DEBUG)
391: #if defined (PETSC_USE_CTABLE)
392: { PetscInt data;
393: PetscTableFind(baij->colmap,in[j]+1,&data);
394: if ((data - 1) % bs) SETERRQ(PETSC_ERR_PLIB,"Incorrect colmap");
395: }
396: #else
397: if ((baij->colmap[in[j]] - 1) % bs) SETERRQ(PETSC_ERR_PLIB,"Incorrect colmap");
398: #endif
399: #endif
400: #if defined (PETSC_USE_CTABLE)
401: PetscTableFind(baij->colmap,in[j]+1,&col);
402: col = (col - 1)/bs;
403: #else
404: col = (baij->colmap[in[j]] - 1)/bs;
405: #endif
406: if (col < 0 && !((Mat_SeqBAIJ*)(baij->A->data))->nonew) {
407: DisAssemble_MPISBAIJ(mat);
408: col = in[j];
409: }
410: }
411: else col = in[j];
412: MatSetValuesBlocked_SeqBAIJ(baij->B,1,&row,1,&col,barray,addv);
413: }
414: }
415: } else {
416: if (!baij->donotstash) {
417: if (roworiented) {
418: MatStashValuesRowBlocked_Private(&mat->bstash,im[i],n,in,v,m,n,i);
419: } else {
420: MatStashValuesColBlocked_Private(&mat->bstash,im[i],n,in,v,m,n,i);
421: }
422: }
423: }
424: }
425: return(0);
426: }
430: PetscErrorCode MatGetValues_MPISBAIJ(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[])
431: {
432: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
434: PetscInt bs=mat->rmap.bs,i,j,bsrstart = mat->rmap.rstart,bsrend = mat->rmap.rend;
435: PetscInt bscstart = mat->cmap.rstart,bscend = mat->cmap.rend,row,col,data;
438: for (i=0; i<m; i++) {
439: if (idxm[i] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",idxm[i]);
440: if (idxm[i] >= mat->rmap.N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm[i],mat->rmap.N-1);
441: if (idxm[i] >= bsrstart && idxm[i] < bsrend) {
442: row = idxm[i] - bsrstart;
443: for (j=0; j<n; j++) {
444: if (idxn[j] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative column %D",idxn[j]);
445: if (idxn[j] >= mat->cmap.N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",idxn[j],mat->cmap.N-1);
446: if (idxn[j] >= bscstart && idxn[j] < bscend){
447: col = idxn[j] - bscstart;
448: MatGetValues_SeqSBAIJ(baij->A,1,&row,1,&col,v+i*n+j);
449: } else {
450: if (!baij->colmap) {
451: CreateColmap_MPIBAIJ_Private(mat);
452: }
453: #if defined (PETSC_USE_CTABLE)
454: PetscTableFind(baij->colmap,idxn[j]/bs+1,&data);
455: data --;
456: #else
457: data = baij->colmap[idxn[j]/bs]-1;
458: #endif
459: if((data < 0) || (baij->garray[data/bs] != idxn[j]/bs)) *(v+i*n+j) = 0.0;
460: else {
461: col = data + idxn[j]%bs;
462: MatGetValues_SeqBAIJ(baij->B,1,&row,1,&col,v+i*n+j);
463: }
464: }
465: }
466: } else {
467: SETERRQ(PETSC_ERR_SUP,"Only local values currently supported");
468: }
469: }
470: return(0);
471: }
475: PetscErrorCode MatNorm_MPISBAIJ(Mat mat,NormType type,PetscReal *norm)
476: {
477: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
479: PetscReal sum[2],*lnorm2;
482: if (baij->size == 1) {
483: MatNorm(baij->A,type,norm);
484: } else {
485: if (type == NORM_FROBENIUS) {
486: PetscMalloc(2*sizeof(PetscReal),&lnorm2);
487: MatNorm(baij->A,type,lnorm2);
488: *lnorm2 = (*lnorm2)*(*lnorm2); lnorm2++; /* squar power of norm(A) */
489: MatNorm(baij->B,type,lnorm2);
490: *lnorm2 = (*lnorm2)*(*lnorm2); lnorm2--; /* squar power of norm(B) */
491: MPI_Allreduce(lnorm2,&sum,2,MPIU_REAL,MPI_SUM,mat->comm);
492: *norm = sqrt(sum[0] + 2*sum[1]);
493: PetscFree(lnorm2);
494: } else if (type == NORM_INFINITY || type == NORM_1) { /* max row/column sum */
495: Mat_SeqSBAIJ *amat=(Mat_SeqSBAIJ*)baij->A->data;
496: Mat_SeqBAIJ *bmat=(Mat_SeqBAIJ*)baij->B->data;
497: PetscReal *rsum,*rsum2,vabs;
498: PetscInt *jj,*garray=baij->garray,rstart=baij->rstartbs,nz;
499: PetscInt brow,bcol,col,bs=baij->A->rmap.bs,row,grow,gcol,mbs=amat->mbs;
500: MatScalar *v;
502: PetscMalloc((2*mat->cmap.N+1)*sizeof(PetscReal),&rsum);
503: rsum2 = rsum + mat->cmap.N;
504: PetscMemzero(rsum,mat->cmap.N*sizeof(PetscReal));
505: /* Amat */
506: v = amat->a; jj = amat->j;
507: for (brow=0; brow<mbs; brow++) {
508: grow = bs*(rstart + brow);
509: nz = amat->i[brow+1] - amat->i[brow];
510: for (bcol=0; bcol<nz; bcol++){
511: gcol = bs*(rstart + *jj); jj++;
512: for (col=0; col<bs; col++){
513: for (row=0; row<bs; row++){
514: vabs = PetscAbsScalar(*v); v++;
515: rsum[gcol+col] += vabs;
516: /* non-diagonal block */
517: if (bcol > 0 && vabs > 0.0) rsum[grow+row] += vabs;
518: }
519: }
520: }
521: }
522: /* Bmat */
523: v = bmat->a; jj = bmat->j;
524: for (brow=0; brow<mbs; brow++) {
525: grow = bs*(rstart + brow);
526: nz = bmat->i[brow+1] - bmat->i[brow];
527: for (bcol=0; bcol<nz; bcol++){
528: gcol = bs*garray[*jj]; jj++;
529: for (col=0; col<bs; col++){
530: for (row=0; row<bs; row++){
531: vabs = PetscAbsScalar(*v); v++;
532: rsum[gcol+col] += vabs;
533: rsum[grow+row] += vabs;
534: }
535: }
536: }
537: }
538: MPI_Allreduce(rsum,rsum2,mat->cmap.N,MPIU_REAL,MPI_SUM,mat->comm);
539: *norm = 0.0;
540: for (col=0; col<mat->cmap.N; col++) {
541: if (rsum2[col] > *norm) *norm = rsum2[col];
542: }
543: PetscFree(rsum);
544: } else {
545: SETERRQ(PETSC_ERR_SUP,"No support for this norm yet");
546: }
547: }
548: return(0);
549: }
553: PetscErrorCode MatAssemblyBegin_MPISBAIJ(Mat mat,MatAssemblyType mode)
554: {
555: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
557: PetscInt nstash,reallocs;
558: InsertMode addv;
561: if (baij->donotstash) {
562: return(0);
563: }
565: /* make sure all processors are either in INSERTMODE or ADDMODE */
566: MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,mat->comm);
567: if (addv == (ADD_VALUES|INSERT_VALUES)) {
568: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Some processors inserted others added");
569: }
570: mat->insertmode = addv; /* in case this processor had no cache */
572: MatStashScatterBegin_Private(&mat->stash,mat->rmap.range);
573: MatStashScatterBegin_Private(&mat->bstash,baij->rangebs);
574: MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);
575: PetscInfo2(0,"Stash has %D entries,uses %D mallocs.\n",nstash,reallocs);
576: MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);
577: PetscInfo2(0,"Block-Stash has %D entries, uses %D mallocs.\n",nstash,reallocs);
578: return(0);
579: }
583: PetscErrorCode MatAssemblyEnd_MPISBAIJ(Mat mat,MatAssemblyType mode)
584: {
585: Mat_MPISBAIJ *baij=(Mat_MPISBAIJ*)mat->data;
586: Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)baij->A->data;
588: PetscInt i,j,rstart,ncols,flg,bs2=baij->bs2;
589: PetscInt *row,*col,other_disassembled;
590: PetscMPIInt n;
591: PetscTruth r1,r2,r3;
592: MatScalar *val;
593: InsertMode addv = mat->insertmode;
595: /* do not use 'b=(Mat_SeqBAIJ*)baij->B->data' as B can be reset in disassembly */
598: if (!baij->donotstash) {
599: while (1) {
600: MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);
601: if (!flg) break;
603: for (i=0; i<n;) {
604: /* Now identify the consecutive vals belonging to the same row */
605: for (j=i,rstart=row[j]; j<n; j++) { if (row[j] != rstart) break; }
606: if (j < n) ncols = j-i;
607: else ncols = n-i;
608: /* Now assemble all these values with a single function call */
609: MatSetValues_MPISBAIJ_MatScalar(mat,1,row+i,ncols,col+i,val+i,addv);
610: i = j;
611: }
612: }
613: MatStashScatterEnd_Private(&mat->stash);
614: /* Now process the block-stash. Since the values are stashed column-oriented,
615: set the roworiented flag to column oriented, and after MatSetValues()
616: restore the original flags */
617: r1 = baij->roworiented;
618: r2 = a->roworiented;
619: r3 = ((Mat_SeqBAIJ*)baij->B->data)->roworiented;
620: baij->roworiented = PETSC_FALSE;
621: a->roworiented = PETSC_FALSE;
622: ((Mat_SeqBAIJ*)baij->B->data)->roworiented = PETSC_FALSE; /* b->roworinted */
623: while (1) {
624: MatStashScatterGetMesg_Private(&mat->bstash,&n,&row,&col,&val,&flg);
625: if (!flg) break;
626:
627: for (i=0; i<n;) {
628: /* Now identify the consecutive vals belonging to the same row */
629: for (j=i,rstart=row[j]; j<n; j++) { if (row[j] != rstart) break; }
630: if (j < n) ncols = j-i;
631: else ncols = n-i;
632: MatSetValuesBlocked_MPISBAIJ_MatScalar(mat,1,row+i,ncols,col+i,val+i*bs2,addv);
633: i = j;
634: }
635: }
636: MatStashScatterEnd_Private(&mat->bstash);
637: baij->roworiented = r1;
638: a->roworiented = r2;
639: ((Mat_SeqBAIJ*)baij->B->data)->roworiented = r3; /* b->roworinted */
640: }
642: MatAssemblyBegin(baij->A,mode);
643: MatAssemblyEnd(baij->A,mode);
645: /* determine if any processor has disassembled, if so we must
646: also disassemble ourselfs, in order that we may reassemble. */
647: /*
648: if nonzero structure of submatrix B cannot change then we know that
649: no processor disassembled thus we can skip this stuff
650: */
651: if (!((Mat_SeqBAIJ*)baij->B->data)->nonew) {
652: MPI_Allreduce(&mat->was_assembled,&other_disassembled,1,MPI_INT,MPI_PROD,mat->comm);
653: if (mat->was_assembled && !other_disassembled) {
654: DisAssemble_MPISBAIJ(mat);
655: }
656: }
658: if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
659: MatSetUpMultiply_MPISBAIJ(mat); /* setup Mvctx and sMvctx */
660: }
661: ((Mat_SeqBAIJ*)baij->B->data)->compressedrow.use = PETSC_TRUE; /* b->compressedrow.use */
662: MatAssemblyBegin(baij->B,mode);
663: MatAssemblyEnd(baij->B,mode);
664:
665: PetscFree(baij->rowvalues);
666: baij->rowvalues = 0;
668: return(0);
669: }
674: static PetscErrorCode MatView_MPISBAIJ_ASCIIorDraworSocket(Mat mat,PetscViewer viewer)
675: {
676: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
677: PetscErrorCode ierr;
678: PetscInt bs = mat->rmap.bs;
679: PetscMPIInt size = baij->size,rank = baij->rank;
680: PetscTruth iascii,isdraw;
681: PetscViewer sviewer;
682: PetscViewerFormat format;
685: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
686: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
687: if (iascii) {
688: PetscViewerGetFormat(viewer,&format);
689: if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
690: MatInfo info;
691: MPI_Comm_rank(mat->comm,&rank);
692: MatGetInfo(mat,MAT_LOCAL,&info);
693: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Local rows %D nz %D nz alloced %D bs %D mem %D\n",
694: rank,mat->rmap.N,(PetscInt)info.nz_used*bs,(PetscInt)info.nz_allocated*bs,
695: mat->rmap.bs,(PetscInt)info.memory);
696: MatGetInfo(baij->A,MAT_LOCAL,&info);
697: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] on-diagonal part: nz %D \n",rank,(PetscInt)info.nz_used*bs);
698: MatGetInfo(baij->B,MAT_LOCAL,&info);
699: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] off-diagonal part: nz %D \n",rank,(PetscInt)info.nz_used*bs);
700: PetscViewerFlush(viewer);
701: VecScatterView(baij->Mvctx,viewer);
702: return(0);
703: } else if (format == PETSC_VIEWER_ASCII_INFO) {
704: PetscViewerASCIIPrintf(viewer," block size is %D\n",bs);
705: return(0);
706: } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
707: return(0);
708: }
709: }
711: if (isdraw) {
712: PetscDraw draw;
713: PetscTruth isnull;
714: PetscViewerDrawGetDraw(viewer,0,&draw);
715: PetscDrawIsNull(draw,&isnull); if (isnull) return(0);
716: }
718: if (size == 1) {
719: PetscObjectSetName((PetscObject)baij->A,mat->name);
720: MatView(baij->A,viewer);
721: } else {
722: /* assemble the entire matrix onto first processor. */
723: Mat A;
724: Mat_SeqSBAIJ *Aloc;
725: Mat_SeqBAIJ *Bloc;
726: PetscInt M = mat->rmap.N,N = mat->cmap.N,*ai,*aj,col,i,j,k,*rvals,mbs = baij->mbs;
727: MatScalar *a;
729: /* Should this be the same type as mat? */
730: MatCreate(mat->comm,&A);
731: if (!rank) {
732: MatSetSizes(A,M,N,M,N);
733: } else {
734: MatSetSizes(A,0,0,M,N);
735: }
736: MatSetType(A,MATMPISBAIJ);
737: MatMPISBAIJSetPreallocation(A,mat->rmap.bs,0,PETSC_NULL,0,PETSC_NULL);
738: PetscLogObjectParent(mat,A);
740: /* copy over the A part */
741: Aloc = (Mat_SeqSBAIJ*)baij->A->data;
742: ai = Aloc->i; aj = Aloc->j; a = Aloc->a;
743: PetscMalloc(bs*sizeof(PetscInt),&rvals);
745: for (i=0; i<mbs; i++) {
746: rvals[0] = bs*(baij->rstartbs + i);
747: for (j=1; j<bs; j++) { rvals[j] = rvals[j-1] + 1; }
748: for (j=ai[i]; j<ai[i+1]; j++) {
749: col = (baij->cstartbs+aj[j])*bs;
750: for (k=0; k<bs; k++) {
751: MatSetValues_MPISBAIJ_MatScalar(A,bs,rvals,1,&col,a,INSERT_VALUES);
752: col++; a += bs;
753: }
754: }
755: }
756: /* copy over the B part */
757: Bloc = (Mat_SeqBAIJ*)baij->B->data;
758: ai = Bloc->i; aj = Bloc->j; a = Bloc->a;
759: for (i=0; i<mbs; i++) {
760:
761: rvals[0] = bs*(baij->rstartbs + i);
762: for (j=1; j<bs; j++) { rvals[j] = rvals[j-1] + 1; }
763: for (j=ai[i]; j<ai[i+1]; j++) {
764: col = baij->garray[aj[j]]*bs;
765: for (k=0; k<bs; k++) {
766: MatSetValues_MPIBAIJ(A,bs,rvals,1,&col,a,INSERT_VALUES);
767: col++; a += bs;
768: }
769: }
770: }
771: PetscFree(rvals);
772: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
773: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
774: /*
775: Everyone has to call to draw the matrix since the graphics waits are
776: synchronized across all processors that share the PetscDraw object
777: */
778: PetscViewerGetSingleton(viewer,&sviewer);
779: if (!rank) {
780: PetscObjectSetName((PetscObject)((Mat_MPISBAIJ*)(A->data))->A,mat->name);
781: MatView(((Mat_MPISBAIJ*)(A->data))->A,sviewer);
782: }
783: PetscViewerRestoreSingleton(viewer,&sviewer);
784: MatDestroy(A);
785: }
786: return(0);
787: }
791: PetscErrorCode MatView_MPISBAIJ(Mat mat,PetscViewer viewer)
792: {
794: PetscTruth iascii,isdraw,issocket,isbinary;
797: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
798: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
799: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);
800: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
801: if (iascii || isdraw || issocket || isbinary) {
802: MatView_MPISBAIJ_ASCIIorDraworSocket(mat,viewer);
803: } else {
804: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by MPISBAIJ matrices",((PetscObject)viewer)->type_name);
805: }
806: return(0);
807: }
811: PetscErrorCode MatDestroy_MPISBAIJ(Mat mat)
812: {
813: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
817: #if defined(PETSC_USE_LOG)
818: PetscLogObjectState((PetscObject)mat,"Rows=%D,Cols=%D",mat->rmap.N,mat->cmap.N);
819: #endif
820: MatStashDestroy_Private(&mat->stash);
821: MatStashDestroy_Private(&mat->bstash);
822: MatDestroy(baij->A);
823: MatDestroy(baij->B);
824: #if defined (PETSC_USE_CTABLE)
825: if (baij->colmap) {PetscTableDelete(baij->colmap);}
826: #else
827: PetscFree(baij->colmap);
828: #endif
829: PetscFree(baij->garray);
830: if (baij->lvec) {VecDestroy(baij->lvec);}
831: if (baij->Mvctx) {VecScatterDestroy(baij->Mvctx);}
832: if (baij->slvec0) {
833: VecDestroy(baij->slvec0);
834: VecDestroy(baij->slvec0b);
835: }
836: if (baij->slvec1) {
837: VecDestroy(baij->slvec1);
838: VecDestroy(baij->slvec1a);
839: VecDestroy(baij->slvec1b);
840: }
841: if (baij->sMvctx) {VecScatterDestroy(baij->sMvctx);}
842: PetscFree(baij->rowvalues);
843: PetscFree(baij->barray);
844: PetscFree(baij->hd);
845: #if defined(PETSC_USE_MAT_SINGLE)
846: PetscFree(baij->setvaluescopy);
847: #endif
848: PetscFree(baij->in_loc);
849: PetscFree(baij->v_loc);
850: PetscFree(baij->rangebs);
851: PetscFree(baij);
853: PetscObjectChangeTypeName((PetscObject)mat,0);
854: PetscObjectComposeFunction((PetscObject)mat,"MatStoreValues_C","",PETSC_NULL);
855: PetscObjectComposeFunction((PetscObject)mat,"MatRetrieveValues_C","",PETSC_NULL);
856: PetscObjectComposeFunction((PetscObject)mat,"MatGetDiagonalBlock_C","",PETSC_NULL);
857: PetscObjectComposeFunction((PetscObject)mat,"MatMPISBAIJSetPreallocation_C","",PETSC_NULL);
858: return(0);
859: }
863: PetscErrorCode MatMult_MPISBAIJ(Mat A,Vec xx,Vec yy)
864: {
865: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
867: PetscInt nt,mbs=a->mbs,bs=A->rmap.bs;
868: PetscScalar *x,*from,zero=0.0;
869:
871: VecGetLocalSize(xx,&nt);
872: if (nt != A->cmap.n) {
873: SETERRQ(PETSC_ERR_ARG_SIZ,"Incompatible partition of A and xx");
874: }
876: /* diagonal part */
877: (*a->A->ops->mult)(a->A,xx,a->slvec1a);
878: VecSet(a->slvec1b,zero);
880: /* subdiagonal part */
881: (*a->B->ops->multtranspose)(a->B,xx,a->slvec0b);
882: CHKMEMQ;
883: /* copy x into the vec slvec0 */
884: VecGetArray(a->slvec0,&from);
885: VecGetArray(xx,&x);
886: CHKMEMQ;
887: PetscMemcpy(from,x,bs*mbs*sizeof(MatScalar));
888: CHKMEMQ;
889: VecRestoreArray(a->slvec0,&from);
890:
891: CHKMEMQ;
892: VecScatterBegin(a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD,a->sMvctx);
893: CHKMEMQ;
894: VecRestoreArray(xx,&x);
895: CHKMEMQ;
896: VecScatterEnd(a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD,a->sMvctx);
897: CHKMEMQ;
898: /* supperdiagonal part */
899: (*a->B->ops->multadd)(a->B,a->slvec1b,a->slvec1a,yy);
900: CHKMEMQ;
901: return(0);
902: }
906: PetscErrorCode MatMult_MPISBAIJ_2comm(Mat A,Vec xx,Vec yy)
907: {
908: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
910: PetscInt nt;
913: VecGetLocalSize(xx,&nt);
914: if (nt != A->cmap.n) {
915: SETERRQ(PETSC_ERR_ARG_SIZ,"Incompatible partition of A and xx");
916: }
917: VecGetLocalSize(yy,&nt);
918: if (nt != A->rmap.N) {
919: SETERRQ(PETSC_ERR_ARG_SIZ,"Incompatible parition of A and yy");
920: }
922: VecScatterBegin(xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);
923: /* do diagonal part */
924: (*a->A->ops->mult)(a->A,xx,yy);
925: /* do supperdiagonal part */
926: VecScatterEnd(xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);
927: (*a->B->ops->multadd)(a->B,a->lvec,yy,yy);
928: /* do subdiagonal part */
929: (*a->B->ops->multtranspose)(a->B,xx,a->lvec);
930: VecScatterBegin(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);
931: VecScatterEnd(a->lvec,yy,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);
933: return(0);
934: }
938: PetscErrorCode MatMultAdd_MPISBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
939: {
940: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
942: PetscInt mbs=a->mbs,bs=A->rmap.bs;
943: PetscScalar *x,*from,zero=0.0;
944:
946: /*
947: PetscSynchronizedPrintf(A->comm," MatMultAdd is called ...\n");
948: PetscSynchronizedFlush(A->comm);
949: */
950: /* diagonal part */
951: (*a->A->ops->multadd)(a->A,xx,yy,a->slvec1a);
952: VecSet(a->slvec1b,zero);
954: /* subdiagonal part */
955: (*a->B->ops->multtranspose)(a->B,xx,a->slvec0b);
957: /* copy x into the vec slvec0 */
958: VecGetArray(a->slvec0,&from);
959: VecGetArray(xx,&x);
960: PetscMemcpy(from,x,bs*mbs*sizeof(MatScalar));
961: VecRestoreArray(a->slvec0,&from);
962:
963: VecScatterBegin(a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD,a->sMvctx);
964: VecRestoreArray(xx,&x);
965: VecScatterEnd(a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD,a->sMvctx);
966:
967: /* supperdiagonal part */
968: (*a->B->ops->multadd)(a->B,a->slvec1b,a->slvec1a,zz);
969:
970: return(0);
971: }
975: PetscErrorCode MatMultAdd_MPISBAIJ_2comm(Mat A,Vec xx,Vec yy,Vec zz)
976: {
977: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
981: VecScatterBegin(xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);
982: /* do diagonal part */
983: (*a->A->ops->multadd)(a->A,xx,yy,zz);
984: /* do supperdiagonal part */
985: VecScatterEnd(xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);
986: (*a->B->ops->multadd)(a->B,a->lvec,zz,zz);
988: /* do subdiagonal part */
989: (*a->B->ops->multtranspose)(a->B,xx,a->lvec);
990: VecScatterBegin(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);
991: VecScatterEnd(a->lvec,zz,ADD_VALUES,SCATTER_REVERSE,a->Mvctx);
993: return(0);
994: }
996: /*
997: This only works correctly for square matrices where the subblock A->A is the
998: diagonal block
999: */
1002: PetscErrorCode MatGetDiagonal_MPISBAIJ(Mat A,Vec v)
1003: {
1004: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1008: /* if (a->rmap.N != a->cmap.N) SETERRQ(PETSC_ERR_SUP,"Supports only square matrix where A->A is diag block"); */
1009: MatGetDiagonal(a->A,v);
1010: return(0);
1011: }
1015: PetscErrorCode MatScale_MPISBAIJ(Mat A,PetscScalar aa)
1016: {
1017: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1021: MatScale(a->A,aa);
1022: MatScale(a->B,aa);
1023: return(0);
1024: }
1028: PetscErrorCode MatGetRow_MPISBAIJ(Mat matin,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
1029: {
1030: Mat_MPISBAIJ *mat = (Mat_MPISBAIJ*)matin->data;
1031: PetscScalar *vworkA,*vworkB,**pvA,**pvB,*v_p;
1033: PetscInt bs = matin->rmap.bs,bs2 = mat->bs2,i,*cworkA,*cworkB,**pcA,**pcB;
1034: PetscInt nztot,nzA,nzB,lrow,brstart = matin->rmap.rstart,brend = matin->rmap.rend;
1035: PetscInt *cmap,*idx_p,cstart = mat->rstartbs;
1038: if (mat->getrowactive == PETSC_TRUE) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Already active");
1039: mat->getrowactive = PETSC_TRUE;
1041: if (!mat->rowvalues && (idx || v)) {
1042: /*
1043: allocate enough space to hold information from the longest row.
1044: */
1045: Mat_SeqSBAIJ *Aa = (Mat_SeqSBAIJ*)mat->A->data;
1046: Mat_SeqBAIJ *Ba = (Mat_SeqBAIJ*)mat->B->data;
1047: PetscInt max = 1,mbs = mat->mbs,tmp;
1048: for (i=0; i<mbs; i++) {
1049: tmp = Aa->i[i+1] - Aa->i[i] + Ba->i[i+1] - Ba->i[i]; /* row length */
1050: if (max < tmp) { max = tmp; }
1051: }
1052: PetscMalloc(max*bs2*(sizeof(PetscInt)+sizeof(PetscScalar)),&mat->rowvalues);
1053: mat->rowindices = (PetscInt*)(mat->rowvalues + max*bs2);
1054: }
1055:
1056: if (row < brstart || row >= brend) SETERRQ(PETSC_ERR_SUP,"Only local rows")
1057: lrow = row - brstart; /* local row index */
1059: pvA = &vworkA; pcA = &cworkA; pvB = &vworkB; pcB = &cworkB;
1060: if (!v) {pvA = 0; pvB = 0;}
1061: if (!idx) {pcA = 0; if (!v) pcB = 0;}
1062: (*mat->A->ops->getrow)(mat->A,lrow,&nzA,pcA,pvA);
1063: (*mat->B->ops->getrow)(mat->B,lrow,&nzB,pcB,pvB);
1064: nztot = nzA + nzB;
1066: cmap = mat->garray;
1067: if (v || idx) {
1068: if (nztot) {
1069: /* Sort by increasing column numbers, assuming A and B already sorted */
1070: PetscInt imark = -1;
1071: if (v) {
1072: *v = v_p = mat->rowvalues;
1073: for (i=0; i<nzB; i++) {
1074: if (cmap[cworkB[i]/bs] < cstart) v_p[i] = vworkB[i];
1075: else break;
1076: }
1077: imark = i;
1078: for (i=0; i<nzA; i++) v_p[imark+i] = vworkA[i];
1079: for (i=imark; i<nzB; i++) v_p[nzA+i] = vworkB[i];
1080: }
1081: if (idx) {
1082: *idx = idx_p = mat->rowindices;
1083: if (imark > -1) {
1084: for (i=0; i<imark; i++) {
1085: idx_p[i] = cmap[cworkB[i]/bs]*bs + cworkB[i]%bs;
1086: }
1087: } else {
1088: for (i=0; i<nzB; i++) {
1089: if (cmap[cworkB[i]/bs] < cstart)
1090: idx_p[i] = cmap[cworkB[i]/bs]*bs + cworkB[i]%bs ;
1091: else break;
1092: }
1093: imark = i;
1094: }
1095: for (i=0; i<nzA; i++) idx_p[imark+i] = cstart*bs + cworkA[i];
1096: for (i=imark; i<nzB; i++) idx_p[nzA+i] = cmap[cworkB[i]/bs]*bs + cworkB[i]%bs ;
1097: }
1098: } else {
1099: if (idx) *idx = 0;
1100: if (v) *v = 0;
1101: }
1102: }
1103: *nz = nztot;
1104: (*mat->A->ops->restorerow)(mat->A,lrow,&nzA,pcA,pvA);
1105: (*mat->B->ops->restorerow)(mat->B,lrow,&nzB,pcB,pvB);
1106: return(0);
1107: }
1111: PetscErrorCode MatRestoreRow_MPISBAIJ(Mat mat,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
1112: {
1113: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
1116: if (!baij->getrowactive) {
1117: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"MatGetRow() must be called first");
1118: }
1119: baij->getrowactive = PETSC_FALSE;
1120: return(0);
1121: }
1125: PetscErrorCode MatGetRowUpperTriangular_MPISBAIJ(Mat A)
1126: {
1127: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1128: Mat_SeqSBAIJ *aA = (Mat_SeqSBAIJ*)a->A->data;
1131: aA->getrow_utriangular = PETSC_TRUE;
1132: return(0);
1133: }
1136: PetscErrorCode MatRestoreRowUpperTriangular_MPISBAIJ(Mat A)
1137: {
1138: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1139: Mat_SeqSBAIJ *aA = (Mat_SeqSBAIJ*)a->A->data;
1142: aA->getrow_utriangular = PETSC_FALSE;
1143: return(0);
1144: }
1148: PetscErrorCode MatRealPart_MPISBAIJ(Mat A)
1149: {
1150: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1154: MatRealPart(a->A);
1155: MatRealPart(a->B);
1156: return(0);
1157: }
1161: PetscErrorCode MatImaginaryPart_MPISBAIJ(Mat A)
1162: {
1163: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1167: MatImaginaryPart(a->A);
1168: MatImaginaryPart(a->B);
1169: return(0);
1170: }
1174: PetscErrorCode MatZeroEntries_MPISBAIJ(Mat A)
1175: {
1176: Mat_MPISBAIJ *l = (Mat_MPISBAIJ*)A->data;
1180: MatZeroEntries(l->A);
1181: MatZeroEntries(l->B);
1182: return(0);
1183: }
1187: PetscErrorCode MatGetInfo_MPISBAIJ(Mat matin,MatInfoType flag,MatInfo *info)
1188: {
1189: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)matin->data;
1190: Mat A = a->A,B = a->B;
1192: PetscReal isend[5],irecv[5];
1195: info->block_size = (PetscReal)matin->rmap.bs;
1196: MatGetInfo(A,MAT_LOCAL,info);
1197: isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
1198: isend[3] = info->memory; isend[4] = info->mallocs;
1199: MatGetInfo(B,MAT_LOCAL,info);
1200: isend[0] += info->nz_used; isend[1] += info->nz_allocated; isend[2] += info->nz_unneeded;
1201: isend[3] += info->memory; isend[4] += info->mallocs;
1202: if (flag == MAT_LOCAL) {
1203: info->nz_used = isend[0];
1204: info->nz_allocated = isend[1];
1205: info->nz_unneeded = isend[2];
1206: info->memory = isend[3];
1207: info->mallocs = isend[4];
1208: } else if (flag == MAT_GLOBAL_MAX) {
1209: MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_MAX,matin->comm);
1210: info->nz_used = irecv[0];
1211: info->nz_allocated = irecv[1];
1212: info->nz_unneeded = irecv[2];
1213: info->memory = irecv[3];
1214: info->mallocs = irecv[4];
1215: } else if (flag == MAT_GLOBAL_SUM) {
1216: MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_SUM,matin->comm);
1217: info->nz_used = irecv[0];
1218: info->nz_allocated = irecv[1];
1219: info->nz_unneeded = irecv[2];
1220: info->memory = irecv[3];
1221: info->mallocs = irecv[4];
1222: } else {
1223: SETERRQ1(PETSC_ERR_ARG_WRONG,"Unknown MatInfoType argument %d",(int)flag);
1224: }
1225: info->rows_global = (PetscReal)A->rmap.N;
1226: info->columns_global = (PetscReal)A->cmap.N;
1227: info->rows_local = (PetscReal)A->rmap.N;
1228: info->columns_local = (PetscReal)A->cmap.N;
1229: info->fill_ratio_given = 0; /* no parallel LU/ILU/Cholesky */
1230: info->fill_ratio_needed = 0;
1231: info->factor_mallocs = 0;
1232: return(0);
1233: }
1237: PetscErrorCode MatSetOption_MPISBAIJ(Mat A,MatOption op)
1238: {
1239: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1240: Mat_SeqSBAIJ *aA = (Mat_SeqSBAIJ*)a->A->data;
1244: switch (op) {
1245: case MAT_NO_NEW_NONZERO_LOCATIONS:
1246: case MAT_YES_NEW_NONZERO_LOCATIONS:
1247: case MAT_COLUMNS_UNSORTED:
1248: case MAT_COLUMNS_SORTED:
1249: case MAT_NEW_NONZERO_ALLOCATION_ERR:
1250: case MAT_KEEP_ZEROED_ROWS:
1251: case MAT_NEW_NONZERO_LOCATION_ERR:
1252: MatSetOption(a->A,op);
1253: MatSetOption(a->B,op);
1254: break;
1255: case MAT_ROW_ORIENTED:
1256: a->roworiented = PETSC_TRUE;
1257: MatSetOption(a->A,op);
1258: MatSetOption(a->B,op);
1259: break;
1260: case MAT_ROWS_SORTED:
1261: case MAT_ROWS_UNSORTED:
1262: case MAT_YES_NEW_DIAGONALS:
1263: PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
1264: break;
1265: case MAT_COLUMN_ORIENTED:
1266: a->roworiented = PETSC_FALSE;
1267: MatSetOption(a->A,op);
1268: MatSetOption(a->B,op);
1269: break;
1270: case MAT_IGNORE_OFF_PROC_ENTRIES:
1271: a->donotstash = PETSC_TRUE;
1272: break;
1273: case MAT_NO_NEW_DIAGONALS:
1274: SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
1275: case MAT_USE_HASH_TABLE:
1276: a->ht_flag = PETSC_TRUE;
1277: break;
1278: case MAT_NOT_SYMMETRIC:
1279: case MAT_NOT_STRUCTURALLY_SYMMETRIC:
1280: case MAT_HERMITIAN:
1281: SETERRQ(PETSC_ERR_SUP,"Matrix must be symmetric");
1282: case MAT_SYMMETRIC:
1283: case MAT_STRUCTURALLY_SYMMETRIC:
1284: case MAT_NOT_HERMITIAN:
1285: case MAT_SYMMETRY_ETERNAL:
1286: case MAT_NOT_SYMMETRY_ETERNAL:
1287: PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
1288: break;
1289: case MAT_IGNORE_LOWER_TRIANGULAR:
1290: aA->ignore_ltriangular = PETSC_TRUE;
1291: break;
1292: case MAT_ERROR_LOWER_TRIANGULAR:
1293: aA->ignore_ltriangular = PETSC_FALSE;
1294: break;
1295: case MAT_GETROW_UPPERTRIANGULAR:
1296: aA->getrow_utriangular = PETSC_TRUE;
1297: break;
1298: default:
1299: SETERRQ1(PETSC_ERR_SUP,"unknown option %d",op);
1300: }
1301: return(0);
1302: }
1306: PetscErrorCode MatTranspose_MPISBAIJ(Mat A,Mat *B)
1307: {
1310: MatDuplicate(A,MAT_COPY_VALUES,B);
1311: return(0);
1312: }
1316: PetscErrorCode MatDiagonalScale_MPISBAIJ(Mat mat,Vec ll,Vec rr)
1317: {
1318: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
1319: Mat a=baij->A, b=baij->B;
1321: PetscInt nv,m,n;
1322: PetscTruth flg;
1325: if (ll != rr){
1326: VecEqual(ll,rr,&flg);
1327: if (!flg)
1328: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"For symmetric format, left and right scaling vectors must be same\n");
1329: }
1330: if (!ll) return(0);
1332: MatGetLocalSize(mat,&m,&n);
1333: if (m != n) SETERRQ2(PETSC_ERR_ARG_SIZ,"For symmetric format, local size %d %d must be same",m,n);
1334:
1335: VecGetLocalSize(rr,&nv);
1336: if (nv!=n) SETERRQ(PETSC_ERR_ARG_SIZ,"Left and right vector non-conforming local size");
1338: VecScatterBegin(rr,baij->lvec,INSERT_VALUES,SCATTER_FORWARD,baij->Mvctx);
1339:
1340: /* left diagonalscale the off-diagonal part */
1341: (*b->ops->diagonalscale)(b,ll,PETSC_NULL);
1342:
1343: /* scale the diagonal part */
1344: (*a->ops->diagonalscale)(a,ll,rr);
1346: /* right diagonalscale the off-diagonal part */
1347: VecScatterEnd(rr,baij->lvec,INSERT_VALUES,SCATTER_FORWARD,baij->Mvctx);
1348: (*b->ops->diagonalscale)(b,PETSC_NULL,baij->lvec);
1349: return(0);
1350: }
1354: PetscErrorCode MatSetUnfactored_MPISBAIJ(Mat A)
1355: {
1356: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1360: MatSetUnfactored(a->A);
1361: return(0);
1362: }
1364: static PetscErrorCode MatDuplicate_MPISBAIJ(Mat,MatDuplicateOption,Mat *);
1368: PetscErrorCode MatEqual_MPISBAIJ(Mat A,Mat B,PetscTruth *flag)
1369: {
1370: Mat_MPISBAIJ *matB = (Mat_MPISBAIJ*)B->data,*matA = (Mat_MPISBAIJ*)A->data;
1371: Mat a,b,c,d;
1372: PetscTruth flg;
1376: a = matA->A; b = matA->B;
1377: c = matB->A; d = matB->B;
1379: MatEqual(a,c,&flg);
1380: if (flg) {
1381: MatEqual(b,d,&flg);
1382: }
1383: MPI_Allreduce(&flg,flag,1,MPI_INT,MPI_LAND,A->comm);
1384: return(0);
1385: }
1389: PetscErrorCode MatCopy_MPISBAIJ(Mat A,Mat B,MatStructure str)
1390: {
1392: Mat_MPISBAIJ *a = (Mat_MPISBAIJ *)A->data;
1393: Mat_MPISBAIJ *b = (Mat_MPISBAIJ *)B->data;
1396: /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
1397: if ((str != SAME_NONZERO_PATTERN) || (A->ops->copy != B->ops->copy)) {
1398: MatGetRowUpperTriangular(A);
1399: MatCopy_Basic(A,B,str);
1400: MatRestoreRowUpperTriangular(A);
1401: } else {
1402: MatCopy(a->A,b->A,str);
1403: MatCopy(a->B,b->B,str);
1404: }
1405: return(0);
1406: }
1410: PetscErrorCode MatSetUpPreallocation_MPISBAIJ(Mat A)
1411: {
1415: MatMPISBAIJSetPreallocation(A,PetscMax(A->rmap.bs,1),PETSC_DEFAULT,0,PETSC_DEFAULT,0);
1416: return(0);
1417: }
1419: #include petscblaslapack.h
1422: PetscErrorCode MatAXPY_MPISBAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str)
1423: {
1425: Mat_MPISBAIJ *xx=(Mat_MPISBAIJ *)X->data,*yy=(Mat_MPISBAIJ *)Y->data;
1426: PetscBLASInt bnz,one=1;
1427: Mat_SeqSBAIJ *xa,*ya;
1428: Mat_SeqBAIJ *xb,*yb;
1431: if (str == SAME_NONZERO_PATTERN) {
1432: PetscScalar alpha = a;
1433: xa = (Mat_SeqSBAIJ *)xx->A->data;
1434: ya = (Mat_SeqSBAIJ *)yy->A->data;
1435: bnz = (PetscBLASInt)xa->nz;
1436: BLASaxpy_(&bnz,&alpha,xa->a,&one,ya->a,&one);
1437: xb = (Mat_SeqBAIJ *)xx->B->data;
1438: yb = (Mat_SeqBAIJ *)yy->B->data;
1439: bnz = (PetscBLASInt)xb->nz;
1440: BLASaxpy_(&bnz,&alpha,xb->a,&one,yb->a,&one);
1441: } else {
1442: MatGetRowUpperTriangular(X);
1443: MatAXPY_Basic(Y,a,X,str);
1444: MatRestoreRowUpperTriangular(X);
1445: }
1446: return(0);
1447: }
1451: PetscErrorCode MatGetSubMatrices_MPISBAIJ(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1452: {
1454: PetscInt i;
1455: PetscTruth flg;
1458: for (i=0; i<n; i++) {
1459: ISEqual(irow[i],icol[i],&flg);
1460: if (!flg) {
1461: SETERRQ(PETSC_ERR_SUP,"Can only get symmetric submatrix for MPISBAIJ matrices");
1462: }
1463: }
1464: MatGetSubMatrices_MPIBAIJ(A,n,irow,icol,scall,B);
1465: return(0);
1466: }
1467:
1469: /* -------------------------------------------------------------------*/
1470: static struct _MatOps MatOps_Values = {
1471: MatSetValues_MPISBAIJ,
1472: MatGetRow_MPISBAIJ,
1473: MatRestoreRow_MPISBAIJ,
1474: MatMult_MPISBAIJ,
1475: /* 4*/ MatMultAdd_MPISBAIJ,
1476: MatMult_MPISBAIJ, /* transpose versions are same as non-transpose */
1477: MatMultAdd_MPISBAIJ,
1478: 0,
1479: 0,
1480: 0,
1481: /*10*/ 0,
1482: 0,
1483: 0,
1484: MatRelax_MPISBAIJ,
1485: MatTranspose_MPISBAIJ,
1486: /*15*/ MatGetInfo_MPISBAIJ,
1487: MatEqual_MPISBAIJ,
1488: MatGetDiagonal_MPISBAIJ,
1489: MatDiagonalScale_MPISBAIJ,
1490: MatNorm_MPISBAIJ,
1491: /*20*/ MatAssemblyBegin_MPISBAIJ,
1492: MatAssemblyEnd_MPISBAIJ,
1493: 0,
1494: MatSetOption_MPISBAIJ,
1495: MatZeroEntries_MPISBAIJ,
1496: /*25*/ 0,
1497: 0,
1498: 0,
1499: 0,
1500: 0,
1501: /*30*/ MatSetUpPreallocation_MPISBAIJ,
1502: 0,
1503: 0,
1504: 0,
1505: 0,
1506: /*35*/ MatDuplicate_MPISBAIJ,
1507: 0,
1508: 0,
1509: 0,
1510: 0,
1511: /*40*/ MatAXPY_MPISBAIJ,
1512: MatGetSubMatrices_MPISBAIJ,
1513: MatIncreaseOverlap_MPISBAIJ,
1514: MatGetValues_MPISBAIJ,
1515: MatCopy_MPISBAIJ,
1516: /*45*/ 0,
1517: MatScale_MPISBAIJ,
1518: 0,
1519: 0,
1520: 0,
1521: /*50*/ 0,
1522: 0,
1523: 0,
1524: 0,
1525: 0,
1526: /*55*/ 0,
1527: 0,
1528: MatSetUnfactored_MPISBAIJ,
1529: 0,
1530: MatSetValuesBlocked_MPISBAIJ,
1531: /*60*/ 0,
1532: 0,
1533: 0,
1534: 0,
1535: 0,
1536: /*65*/ 0,
1537: 0,
1538: 0,
1539: 0,
1540: 0,
1541: /*70*/ MatGetRowMax_MPISBAIJ,
1542: 0,
1543: 0,
1544: 0,
1545: 0,
1546: /*75*/ 0,
1547: 0,
1548: 0,
1549: 0,
1550: 0,
1551: /*80*/ 0,
1552: 0,
1553: 0,
1554: 0,
1555: MatLoad_MPISBAIJ,
1556: /*85*/ 0,
1557: 0,
1558: 0,
1559: 0,
1560: 0,
1561: /*90*/ 0,
1562: 0,
1563: 0,
1564: 0,
1565: 0,
1566: /*95*/ 0,
1567: 0,
1568: 0,
1569: 0,
1570: 0,
1571: /*100*/0,
1572: 0,
1573: 0,
1574: 0,
1575: 0,
1576: /*105*/0,
1577: MatRealPart_MPISBAIJ,
1578: MatImaginaryPart_MPISBAIJ,
1579: MatGetRowUpperTriangular_MPISBAIJ,
1580: MatRestoreRowUpperTriangular_MPISBAIJ
1581: };
1587: PetscErrorCode MatGetDiagonalBlock_MPISBAIJ(Mat A,PetscTruth *iscopy,MatReuse reuse,Mat *a)
1588: {
1590: *a = ((Mat_MPISBAIJ *)A->data)->A;
1591: *iscopy = PETSC_FALSE;
1592: return(0);
1593: }
1599: PetscErrorCode MatMPISBAIJSetPreallocation_MPISBAIJ(Mat B,PetscInt bs,PetscInt d_nz,PetscInt *d_nnz,PetscInt o_nz,PetscInt *o_nnz)
1600: {
1601: Mat_MPISBAIJ *b;
1603: PetscInt i,mbs,Mbs;
1606: PetscOptionsBegin(B->comm,B->prefix,"Options for MPISBAIJ matrix","Mat");
1607: PetscOptionsInt("-mat_block_size","Set the blocksize used to store the matrix","MatMPIBAIJSetPreallocation",bs,&bs,PETSC_NULL);
1608: PetscOptionsEnd();
1610: if (bs < 1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid block size specified, must be positive");
1611: if (d_nz == PETSC_DECIDE || d_nz == PETSC_DEFAULT) d_nz = 3;
1612: if (o_nz == PETSC_DECIDE || o_nz == PETSC_DEFAULT) o_nz = 1;
1613: if (d_nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"d_nz cannot be less than 0: value %D",d_nz);
1614: if (o_nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"o_nz cannot be less than 0: value %D",o_nz);
1616: B->rmap.bs = B->cmap.bs = bs;
1617: PetscMapInitialize(B->comm,&B->rmap);
1618: PetscMapInitialize(B->comm,&B->cmap);
1620: if (d_nnz) {
1621: for (i=0; i<B->rmap.n/bs; i++) {
1622: if (d_nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"d_nnz cannot be less than -1: local row %D value %D",i,d_nnz[i]);
1623: }
1624: }
1625: if (o_nnz) {
1626: for (i=0; i<B->rmap.n/bs; i++) {
1627: if (o_nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"o_nnz cannot be less than -1: local row %D value %D",i,o_nnz[i]);
1628: }
1629: }
1630: B->preallocated = PETSC_TRUE;
1632: b = (Mat_MPISBAIJ*)B->data;
1633: mbs = B->rmap.n/bs;
1634: Mbs = B->rmap.N/bs;
1635: if (mbs*bs != B->rmap.n) {
1636: SETERRQ2(PETSC_ERR_ARG_SIZ,"No of local rows %D must be divisible by blocksize %D",B->rmap.N,bs);
1637: }
1639: B->rmap.bs = bs;
1640: b->bs2 = bs*bs;
1641: b->mbs = mbs;
1642: b->nbs = mbs;
1643: b->Mbs = Mbs;
1644: b->Nbs = Mbs;
1646: for (i=0; i<=b->size; i++) {
1647: b->rangebs[i] = B->rmap.range[i]/bs;
1648: }
1649: b->rstartbs = B->rmap.rstart/bs;
1650: b->rendbs = B->rmap.rend/bs;
1651:
1652: b->cstartbs = B->cmap.rstart/bs;
1653: b->cendbs = B->cmap.rend/bs;
1654:
1655: MatCreate(PETSC_COMM_SELF,&b->A);
1656: MatSetSizes(b->A,B->rmap.n,B->cmap.n,B->rmap.n,B->cmap.n);
1657: MatSetType(b->A,MATSEQSBAIJ);
1658: MatSeqSBAIJSetPreallocation(b->A,bs,d_nz,d_nnz);
1659: PetscLogObjectParent(B,b->A);
1661: MatCreate(PETSC_COMM_SELF,&b->B);
1662: MatSetSizes(b->B,B->rmap.n,B->cmap.N,B->rmap.n,B->cmap.N);
1663: MatSetType(b->B,MATSEQBAIJ);
1664: MatSeqBAIJSetPreallocation(b->B,bs,o_nz,o_nnz);
1665: PetscLogObjectParent(B,b->B);
1667: /* build cache for off array entries formed */
1668: MatStashCreate_Private(B->comm,bs,&B->bstash);
1670: return(0);
1671: }
1674: /*MC
1675: MATMPISBAIJ - MATMPISBAIJ = "mpisbaij" - A matrix type to be used for distributed symmetric sparse block matrices,
1676: based on block compressed sparse row format. Only the upper triangular portion of the matrix is stored.
1678: Options Database Keys:
1679: . -mat_type mpisbaij - sets the matrix type to "mpisbaij" during a call to MatSetFromOptions()
1681: Level: beginner
1683: .seealso: MatCreateMPISBAIJ
1684: M*/
1689: PetscErrorCode MatCreate_MPISBAIJ(Mat B)
1690: {
1691: Mat_MPISBAIJ *b;
1693: PetscTruth flg;
1697: PetscNew(Mat_MPISBAIJ,&b);
1698: B->data = (void*)b;
1699: PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
1701: B->ops->destroy = MatDestroy_MPISBAIJ;
1702: B->ops->view = MatView_MPISBAIJ;
1703: B->mapping = 0;
1704: B->factor = 0;
1705: B->assembled = PETSC_FALSE;
1707: B->insertmode = NOT_SET_VALUES;
1708: MPI_Comm_rank(B->comm,&b->rank);
1709: MPI_Comm_size(B->comm,&b->size);
1711: /* build local table of row and column ownerships */
1712: PetscMalloc((b->size+2)*sizeof(PetscInt),&b->rangebs);
1714: /* build cache for off array entries formed */
1715: MatStashCreate_Private(B->comm,1,&B->stash);
1716: b->donotstash = PETSC_FALSE;
1717: b->colmap = PETSC_NULL;
1718: b->garray = PETSC_NULL;
1719: b->roworiented = PETSC_TRUE;
1721: #if defined(PETSC_USE_MAT_SINGLE)
1722: /* stuff for MatSetValues_XXX in single precision */
1723: b->setvalueslen = 0;
1724: b->setvaluescopy = PETSC_NULL;
1725: #endif
1727: /* stuff used in block assembly */
1728: b->barray = 0;
1730: /* stuff used for matrix vector multiply */
1731: b->lvec = 0;
1732: b->Mvctx = 0;
1733: b->slvec0 = 0;
1734: b->slvec0b = 0;
1735: b->slvec1 = 0;
1736: b->slvec1a = 0;
1737: b->slvec1b = 0;
1738: b->sMvctx = 0;
1740: /* stuff for MatGetRow() */
1741: b->rowindices = 0;
1742: b->rowvalues = 0;
1743: b->getrowactive = PETSC_FALSE;
1745: /* hash table stuff */
1746: b->ht = 0;
1747: b->hd = 0;
1748: b->ht_size = 0;
1749: b->ht_flag = PETSC_FALSE;
1750: b->ht_fact = 0;
1751: b->ht_total_ct = 0;
1752: b->ht_insert_ct = 0;
1754: b->in_loc = 0;
1755: b->v_loc = 0;
1756: b->n_loc = 0;
1757: PetscOptionsBegin(B->comm,PETSC_NULL,"Options for loading MPISBAIJ matrix 1","Mat");
1758: PetscOptionsTruth("-mat_use_hash_table","Use hash table to save memory in constructing matrix","MatSetOption",PETSC_FALSE,&flg,PETSC_NULL);
1759: if (flg) {
1760: PetscReal fact = 1.39;
1761: MatSetOption(B,MAT_USE_HASH_TABLE);
1762: PetscOptionsReal("-mat_use_hash_table","Use hash table factor","MatMPIBAIJSetHashTableFactor",fact,&fact,PETSC_NULL);
1763: if (fact <= 1.0) fact = 1.39;
1764: MatMPIBAIJSetHashTableFactor(B,fact);
1765: PetscInfo1(0,"Hash table Factor used %5.2f\n",fact);
1766: }
1767: PetscOptionsEnd();
1769: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
1770: "MatStoreValues_MPISBAIJ",
1771: MatStoreValues_MPISBAIJ);
1772: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
1773: "MatRetrieveValues_MPISBAIJ",
1774: MatRetrieveValues_MPISBAIJ);
1775: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetDiagonalBlock_C",
1776: "MatGetDiagonalBlock_MPISBAIJ",
1777: MatGetDiagonalBlock_MPISBAIJ);
1778: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPISBAIJSetPreallocation_C",
1779: "MatMPISBAIJSetPreallocation_MPISBAIJ",
1780: MatMPISBAIJSetPreallocation_MPISBAIJ);
1781: B->symmetric = PETSC_TRUE;
1782: B->structurally_symmetric = PETSC_TRUE;
1783: B->symmetric_set = PETSC_TRUE;
1784: B->structurally_symmetric_set = PETSC_TRUE;
1785: PetscObjectChangeTypeName((PetscObject)B,MATMPISBAIJ);
1786: return(0);
1787: }
1790: /*MC
1791: MATSBAIJ - MATSBAIJ = "sbaij" - A matrix type to be used for symmetric block sparse matrices.
1793: This matrix type is identical to MATSEQSBAIJ when constructed with a single process communicator,
1794: and MATMPISBAIJ otherwise.
1796: Options Database Keys:
1797: . -mat_type sbaij - sets the matrix type to "sbaij" during a call to MatSetFromOptions()
1799: Level: beginner
1801: .seealso: MatCreateMPISBAIJ,MATSEQSBAIJ,MATMPISBAIJ
1802: M*/
1807: PetscErrorCode MatCreate_SBAIJ(Mat A)
1808: {
1810: PetscMPIInt size;
1813: MPI_Comm_size(A->comm,&size);
1814: if (size == 1) {
1815: MatSetType(A,MATSEQSBAIJ);
1816: } else {
1817: MatSetType(A,MATMPISBAIJ);
1818: }
1819: return(0);
1820: }
1825: /*@C
1826: MatMPISBAIJSetPreallocation - For good matrix assembly performance
1827: the user should preallocate the matrix storage by setting the parameters
1828: d_nz (or d_nnz) and o_nz (or o_nnz). By setting these parameters accurately,
1829: performance can be increased by more than a factor of 50.
1831: Collective on Mat
1833: Input Parameters:
1834: + A - the matrix
1835: . bs - size of blockk
1836: . d_nz - number of block nonzeros per block row in diagonal portion of local
1837: submatrix (same for all local rows)
1838: . d_nnz - array containing the number of block nonzeros in the various block rows
1839: in the upper triangular and diagonal part of the in diagonal portion of the local
1840: (possibly different for each block row) or PETSC_NULL. You must leave room
1841: for the diagonal entry even if it is zero.
1842: . o_nz - number of block nonzeros per block row in the off-diagonal portion of local
1843: submatrix (same for all local rows).
1844: - o_nnz - array containing the number of nonzeros in the various block rows of the
1845: off-diagonal portion of the local submatrix (possibly different for
1846: each block row) or PETSC_NULL.
1849: Options Database Keys:
1850: . -mat_no_unroll - uses code that does not unroll the loops in the
1851: block calculations (much slower)
1852: . -mat_block_size - size of the blocks to use
1854: Notes:
1856: If PETSC_DECIDE or PETSC_DETERMINE is used for a particular argument on one processor
1857: than it must be used on all processors that share the object for that argument.
1859: If the *_nnz parameter is given then the *_nz parameter is ignored
1861: Storage Information:
1862: For a square global matrix we define each processor's diagonal portion
1863: to be its local rows and the corresponding columns (a square submatrix);
1864: each processor's off-diagonal portion encompasses the remainder of the
1865: local matrix (a rectangular submatrix).
1867: The user can specify preallocated storage for the diagonal part of
1868: the local submatrix with either d_nz or d_nnz (not both). Set
1869: d_nz=PETSC_DEFAULT and d_nnz=PETSC_NULL for PETSc to control dynamic
1870: memory allocation. Likewise, specify preallocated storage for the
1871: off-diagonal part of the local submatrix with o_nz or o_nnz (not both).
1873: Consider a processor that owns rows 3, 4 and 5 of a parallel matrix. In
1874: the figure below we depict these three local rows and all columns (0-11).
1876: .vb
1877: 0 1 2 3 4 5 6 7 8 9 10 11
1878: -------------------
1879: row 3 | o o o d d d o o o o o o
1880: row 4 | o o o d d d o o o o o o
1881: row 5 | o o o d d d o o o o o o
1882: -------------------
1883: .ve
1884:
1885: Thus, any entries in the d locations are stored in the d (diagonal)
1886: submatrix, and any entries in the o locations are stored in the
1887: o (off-diagonal) submatrix. Note that the d matrix is stored in
1888: MatSeqSBAIJ format and the o submatrix in MATSEQBAIJ format.
1890: Now d_nz should indicate the number of block nonzeros per row in the upper triangular
1891: plus the diagonal part of the d matrix,
1892: and o_nz should indicate the number of block nonzeros per row in the o matrix.
1893: In general, for PDE problems in which most nonzeros are near the diagonal,
1894: one expects d_nz >> o_nz. For large problems you MUST preallocate memory
1895: or you will get TERRIBLE performance; see the users' manual chapter on
1896: matrices.
1898: Level: intermediate
1900: .keywords: matrix, block, aij, compressed row, sparse, parallel
1902: .seealso: MatCreate(), MatCreateSeqSBAIJ(), MatSetValues(), MatCreateMPIBAIJ()
1903: @*/
1904: PetscErrorCode MatMPISBAIJSetPreallocation(Mat B,PetscInt bs,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
1905: {
1906: PetscErrorCode ierr,(*f)(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1909: PetscObjectQueryFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",(void (**)(void))&f);
1910: if (f) {
1911: (*f)(B,bs,d_nz,d_nnz,o_nz,o_nnz);
1912: }
1913: return(0);
1914: }
1918: /*@C
1919: MatCreateMPISBAIJ - Creates a sparse parallel matrix in symmetric block AIJ format
1920: (block compressed row). For good matrix assembly performance
1921: the user should preallocate the matrix storage by setting the parameters
1922: d_nz (or d_nnz) and o_nz (or o_nnz). By setting these parameters accurately,
1923: performance can be increased by more than a factor of 50.
1925: Collective on MPI_Comm
1927: Input Parameters:
1928: + comm - MPI communicator
1929: . bs - size of blockk
1930: . m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
1931: This value should be the same as the local size used in creating the
1932: y vector for the matrix-vector product y = Ax.
1933: . n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
1934: This value should be the same as the local size used in creating the
1935: x vector for the matrix-vector product y = Ax.
1936: . M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
1937: . N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
1938: . d_nz - number of block nonzeros per block row in diagonal portion of local
1939: submatrix (same for all local rows)
1940: . d_nnz - array containing the number of block nonzeros in the various block rows
1941: in the upper triangular portion of the in diagonal portion of the local
1942: (possibly different for each block block row) or PETSC_NULL.
1943: You must leave room for the diagonal entry even if it is zero.
1944: . o_nz - number of block nonzeros per block row in the off-diagonal portion of local
1945: submatrix (same for all local rows).
1946: - o_nnz - array containing the number of nonzeros in the various block rows of the
1947: off-diagonal portion of the local submatrix (possibly different for
1948: each block row) or PETSC_NULL.
1950: Output Parameter:
1951: . A - the matrix
1953: Options Database Keys:
1954: . -mat_no_unroll - uses code that does not unroll the loops in the
1955: block calculations (much slower)
1956: . -mat_block_size - size of the blocks to use
1957: . -mat_mpi - use the parallel matrix data structures even on one processor
1958: (defaults to using SeqBAIJ format on one processor)
1960: Notes:
1961: The number of rows and columns must be divisible by blocksize.
1963: The user MUST specify either the local or global matrix dimensions
1964: (possibly both).
1966: If PETSC_DECIDE or PETSC_DETERMINE is used for a particular argument on one processor
1967: than it must be used on all processors that share the object for that argument.
1969: If the *_nnz parameter is given then the *_nz parameter is ignored
1971: Storage Information:
1972: For a square global matrix we define each processor's diagonal portion
1973: to be its local rows and the corresponding columns (a square submatrix);
1974: each processor's off-diagonal portion encompasses the remainder of the
1975: local matrix (a rectangular submatrix).
1977: The user can specify preallocated storage for the diagonal part of
1978: the local submatrix with either d_nz or d_nnz (not both). Set
1979: d_nz=PETSC_DEFAULT and d_nnz=PETSC_NULL for PETSc to control dynamic
1980: memory allocation. Likewise, specify preallocated storage for the
1981: off-diagonal part of the local submatrix with o_nz or o_nnz (not both).
1983: Consider a processor that owns rows 3, 4 and 5 of a parallel matrix. In
1984: the figure below we depict these three local rows and all columns (0-11).
1986: .vb
1987: 0 1 2 3 4 5 6 7 8 9 10 11
1988: -------------------
1989: row 3 | o o o d d d o o o o o o
1990: row 4 | o o o d d d o o o o o o
1991: row 5 | o o o d d d o o o o o o
1992: -------------------
1993: .ve
1994:
1995: Thus, any entries in the d locations are stored in the d (diagonal)
1996: submatrix, and any entries in the o locations are stored in the
1997: o (off-diagonal) submatrix. Note that the d matrix is stored in
1998: MatSeqSBAIJ format and the o submatrix in MATSEQBAIJ format.
2000: Now d_nz should indicate the number of block nonzeros per row in the upper triangular
2001: plus the diagonal part of the d matrix,
2002: and o_nz should indicate the number of block nonzeros per row in the o matrix.
2003: In general, for PDE problems in which most nonzeros are near the diagonal,
2004: one expects d_nz >> o_nz. For large problems you MUST preallocate memory
2005: or you will get TERRIBLE performance; see the users' manual chapter on
2006: matrices.
2008: Level: intermediate
2010: .keywords: matrix, block, aij, compressed row, sparse, parallel
2012: .seealso: MatCreate(), MatCreateSeqSBAIJ(), MatSetValues(), MatCreateMPIBAIJ()
2013: @*/
2015: PetscErrorCode MatCreateMPISBAIJ(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)
2016: {
2018: PetscMPIInt size;
2021: MatCreate(comm,A);
2022: MatSetSizes(*A,m,n,M,N);
2023: MPI_Comm_size(comm,&size);
2024: if (size > 1) {
2025: MatSetType(*A,MATMPISBAIJ);
2026: MatMPISBAIJSetPreallocation(*A,bs,d_nz,d_nnz,o_nz,o_nnz);
2027: } else {
2028: MatSetType(*A,MATSEQSBAIJ);
2029: MatSeqSBAIJSetPreallocation(*A,bs,d_nz,d_nnz);
2030: }
2031: return(0);
2032: }
2037: static PetscErrorCode MatDuplicate_MPISBAIJ(Mat matin,MatDuplicateOption cpvalues,Mat *newmat)
2038: {
2039: Mat mat;
2040: Mat_MPISBAIJ *a,*oldmat = (Mat_MPISBAIJ*)matin->data;
2042: PetscInt len=0,nt,bs=matin->rmap.bs,mbs=oldmat->mbs;
2043: PetscScalar *array;
2046: *newmat = 0;
2047: MatCreate(matin->comm,&mat);
2048: MatSetSizes(mat,matin->rmap.n,matin->cmap.n,matin->rmap.N,matin->cmap.N);
2049: MatSetType(mat,matin->type_name);
2050: PetscMemcpy(mat->ops,matin->ops,sizeof(struct _MatOps));
2051: PetscMapCopy(matin->comm,&matin->rmap,&mat->rmap);
2052: PetscMapCopy(matin->comm,&matin->cmap,&mat->cmap);
2053:
2054: mat->factor = matin->factor;
2055: mat->preallocated = PETSC_TRUE;
2056: mat->assembled = PETSC_TRUE;
2057: mat->insertmode = NOT_SET_VALUES;
2059: a = (Mat_MPISBAIJ*)mat->data;
2060: a->bs2 = oldmat->bs2;
2061: a->mbs = oldmat->mbs;
2062: a->nbs = oldmat->nbs;
2063: a->Mbs = oldmat->Mbs;
2064: a->Nbs = oldmat->Nbs;
2067: a->size = oldmat->size;
2068: a->rank = oldmat->rank;
2069: a->donotstash = oldmat->donotstash;
2070: a->roworiented = oldmat->roworiented;
2071: a->rowindices = 0;
2072: a->rowvalues = 0;
2073: a->getrowactive = PETSC_FALSE;
2074: a->barray = 0;
2075: a->rstartbs = oldmat->rstartbs;
2076: a->rendbs = oldmat->rendbs;
2077: a->cstartbs = oldmat->cstartbs;
2078: a->cendbs = oldmat->cendbs;
2080: /* hash table stuff */
2081: a->ht = 0;
2082: a->hd = 0;
2083: a->ht_size = 0;
2084: a->ht_flag = oldmat->ht_flag;
2085: a->ht_fact = oldmat->ht_fact;
2086: a->ht_total_ct = 0;
2087: a->ht_insert_ct = 0;
2088:
2089: PetscMemcpy(a->rangebs,oldmat->rangebs,(a->size+2)*sizeof(PetscInt));
2090: MatStashCreate_Private(matin->comm,1,&mat->stash);
2091: MatStashCreate_Private(matin->comm,matin->rmap.bs,&mat->bstash);
2092: if (oldmat->colmap) {
2093: #if defined (PETSC_USE_CTABLE)
2094: PetscTableCreateCopy(oldmat->colmap,&a->colmap);
2095: #else
2096: PetscMalloc((a->Nbs)*sizeof(PetscInt),&a->colmap);
2097: PetscLogObjectMemory(mat,(a->Nbs)*sizeof(PetscInt));
2098: PetscMemcpy(a->colmap,oldmat->colmap,(a->Nbs)*sizeof(PetscInt));
2099: #endif
2100: } else a->colmap = 0;
2102: if (oldmat->garray && (len = ((Mat_SeqBAIJ*)(oldmat->B->data))->nbs)) {
2103: PetscMalloc(len*sizeof(PetscInt),&a->garray);
2104: PetscLogObjectMemory(mat,len*sizeof(PetscInt));
2105: PetscMemcpy(a->garray,oldmat->garray,len*sizeof(PetscInt));
2106: } else a->garray = 0;
2107:
2108: VecDuplicate(oldmat->lvec,&a->lvec);
2109: PetscLogObjectParent(mat,a->lvec);
2110: VecScatterCopy(oldmat->Mvctx,&a->Mvctx);
2111: PetscLogObjectParent(mat,a->Mvctx);
2113: VecDuplicate(oldmat->slvec0,&a->slvec0);
2114: PetscLogObjectParent(mat,a->slvec0);
2115: VecDuplicate(oldmat->slvec1,&a->slvec1);
2116: PetscLogObjectParent(mat,a->slvec1);
2118: VecGetLocalSize(a->slvec1,&nt);
2119: VecGetArray(a->slvec1,&array);
2120: VecCreateSeqWithArray(PETSC_COMM_SELF,bs*mbs,array,&a->slvec1a);
2121: VecCreateSeqWithArray(PETSC_COMM_SELF,nt-bs*mbs,array+bs*mbs,&a->slvec1b);
2122: VecRestoreArray(a->slvec1,&array);
2123: VecGetArray(a->slvec0,&array);
2124: VecCreateSeqWithArray(PETSC_COMM_SELF,nt-bs*mbs,array+bs*mbs,&a->slvec0b);
2125: VecRestoreArray(a->slvec0,&array);
2126: PetscLogObjectParent(mat,a->slvec0);
2127: PetscLogObjectParent(mat,a->slvec1);
2128: PetscLogObjectParent(mat,a->slvec0b);
2129: PetscLogObjectParent(mat,a->slvec1a);
2130: PetscLogObjectParent(mat,a->slvec1b);
2132: /* VecScatterCopy(oldmat->sMvctx,&a->sMvctx); - not written yet, replaced by the lazy trick: */
2133: PetscObjectReference((PetscObject)oldmat->sMvctx);
2134: a->sMvctx = oldmat->sMvctx;
2135: PetscLogObjectParent(mat,a->sMvctx);
2137: MatDuplicate(oldmat->A,cpvalues,&a->A);
2138: PetscLogObjectParent(mat,a->A);
2139: MatDuplicate(oldmat->B,cpvalues,&a->B);
2140: PetscLogObjectParent(mat,a->B);
2141: PetscFListDuplicate(mat->qlist,&matin->qlist);
2142: *newmat = mat;
2143: return(0);
2144: }
2146: #include petscsys.h
2150: PetscErrorCode MatLoad_MPISBAIJ(PetscViewer viewer, MatType type,Mat *newmat)
2151: {
2152: Mat A;
2154: PetscInt i,nz,j,rstart,rend;
2155: PetscScalar *vals,*buf;
2156: MPI_Comm comm = ((PetscObject)viewer)->comm;
2157: MPI_Status status;
2158: PetscMPIInt rank,size,tag = ((PetscObject)viewer)->tag,*sndcounts = 0,*browners,maxnz,*rowners,*locrowlens;
2159: PetscInt header[4],*rowlengths = 0,M,N,m,*cols;
2160: PetscInt *procsnz = 0,jj,*mycols,*ibuf;
2161: PetscInt bs=1,Mbs,mbs,extra_rows;
2162: PetscInt *dlens,*odlens,*mask,*masked1,*masked2,rowcount,odcount;
2163: PetscInt dcount,kmax,k,nzcount,tmp;
2164: int fd;
2165:
2167: PetscOptionsBegin(comm,PETSC_NULL,"Options for loading MPISBAIJ matrix 2","Mat");
2168: PetscOptionsInt("-matload_block_size","Set the blocksize used to store the matrix","MatLoad",bs,&bs,PETSC_NULL);
2169: PetscOptionsEnd();
2171: MPI_Comm_size(comm,&size);
2172: MPI_Comm_rank(comm,&rank);
2173: if (!rank) {
2174: PetscViewerBinaryGetDescriptor(viewer,&fd);
2175: PetscBinaryRead(fd,(char *)header,4,PETSC_INT);
2176: if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object");
2177: if (header[3] < 0) {
2178: SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as MPISBAIJ");
2179: }
2180: }
2182: MPI_Bcast(header+1,3,MPIU_INT,0,comm);
2183: M = header[1]; N = header[2];
2185: if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices");
2187: /*
2188: This code adds extra rows to make sure the number of rows is
2189: divisible by the blocksize
2190: */
2191: Mbs = M/bs;
2192: extra_rows = bs - M + bs*(Mbs);
2193: if (extra_rows == bs) extra_rows = 0;
2194: else Mbs++;
2195: if (extra_rows &&!rank) {
2196: PetscInfo(0,"Padding loaded matrix to match blocksize\n");
2197: }
2199: /* determine ownership of all rows */
2200: mbs = Mbs/size + ((Mbs % size) > rank);
2201: m = mbs*bs;
2202: PetscMalloc(2*(size+2)*sizeof(PetscMPIInt),&rowners);
2203: browners = rowners + size + 1;
2204: MPI_Allgather(&mbs,1,MPI_INT,rowners+1,1,MPI_INT,comm);
2205: rowners[0] = 0;
2206: for (i=2; i<=size; i++) rowners[i] += rowners[i-1];
2207: for (i=0; i<=size; i++) browners[i] = rowners[i]*bs;
2208: rstart = rowners[rank];
2209: rend = rowners[rank+1];
2210:
2211: /* distribute row lengths to all processors */
2212: PetscMalloc((rend-rstart)*bs*sizeof(PetscMPIInt),&locrowlens);
2213: if (!rank) {
2214: PetscMalloc((M+extra_rows)*sizeof(PetscInt),&rowlengths);
2215: PetscBinaryRead(fd,rowlengths,M,PETSC_INT);
2216: for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;
2217: PetscMalloc(size*sizeof(PetscMPIInt),&sndcounts);
2218: for (i=0; i<size; i++) sndcounts[i] = browners[i+1] - browners[i];
2219: MPI_Scatterv(rowlengths,sndcounts,browners,MPIU_INT,locrowlens,(rend-rstart)*bs,MPIU_INT,0,comm);
2220: PetscFree(sndcounts);
2221: } else {
2222: MPI_Scatterv(0,0,0,MPIU_INT,locrowlens,(rend-rstart)*bs,MPIU_INT,0,comm);
2223: }
2224:
2225: if (!rank) { /* procs[0] */
2226: /* calculate the number of nonzeros on each processor */
2227: PetscMalloc(size*sizeof(PetscInt),&procsnz);
2228: PetscMemzero(procsnz,size*sizeof(PetscInt));
2229: for (i=0; i<size; i++) {
2230: for (j=rowners[i]*bs; j< rowners[i+1]*bs; j++) {
2231: procsnz[i] += rowlengths[j];
2232: }
2233: }
2234: PetscFree(rowlengths);
2235:
2236: /* determine max buffer needed and allocate it */
2237: maxnz = 0;
2238: for (i=0; i<size; i++) {
2239: maxnz = PetscMax(maxnz,procsnz[i]);
2240: }
2241: PetscMalloc(maxnz*sizeof(PetscInt),&cols);
2243: /* read in my part of the matrix column indices */
2244: nz = procsnz[0];
2245: PetscMalloc(nz*sizeof(PetscInt),&ibuf);
2246: mycols = ibuf;
2247: if (size == 1) nz -= extra_rows;
2248: PetscBinaryRead(fd,mycols,nz,PETSC_INT);
2249: if (size == 1) for (i=0; i< extra_rows; i++) { mycols[nz+i] = M+i; }
2251: /* read in every ones (except the last) and ship off */
2252: for (i=1; i<size-1; i++) {
2253: nz = procsnz[i];
2254: PetscBinaryRead(fd,cols,nz,PETSC_INT);
2255: MPI_Send(cols,nz,MPIU_INT,i,tag,comm);
2256: }
2257: /* read in the stuff for the last proc */
2258: if (size != 1) {
2259: nz = procsnz[size-1] - extra_rows; /* the extra rows are not on the disk */
2260: PetscBinaryRead(fd,cols,nz,PETSC_INT);
2261: for (i=0; i<extra_rows; i++) cols[nz+i] = M+i;
2262: MPI_Send(cols,nz+extra_rows,MPIU_INT,size-1,tag,comm);
2263: }
2264: PetscFree(cols);
2265: } else { /* procs[i], i>0 */
2266: /* determine buffer space needed for message */
2267: nz = 0;
2268: for (i=0; i<m; i++) {
2269: nz += locrowlens[i];
2270: }
2271: PetscMalloc(nz*sizeof(PetscInt),&ibuf);
2272: mycols = ibuf;
2273: /* receive message of column indices*/
2274: MPI_Recv(mycols,nz,MPIU_INT,0,tag,comm,&status);
2275: MPI_Get_count(&status,MPIU_INT,&maxnz);
2276: if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
2277: }
2279: /* loop over local rows, determining number of off diagonal entries */
2280: PetscMalloc(2*(rend-rstart+1)*sizeof(PetscInt),&dlens);
2281: odlens = dlens + (rend-rstart);
2282: PetscMalloc(3*Mbs*sizeof(PetscInt),&mask);
2283: PetscMemzero(mask,3*Mbs*sizeof(PetscInt));
2284: masked1 = mask + Mbs;
2285: masked2 = masked1 + Mbs;
2286: rowcount = 0; nzcount = 0;
2287: for (i=0; i<mbs; i++) {
2288: dcount = 0;
2289: odcount = 0;
2290: for (j=0; j<bs; j++) {
2291: kmax = locrowlens[rowcount];
2292: for (k=0; k<kmax; k++) {
2293: tmp = mycols[nzcount++]/bs; /* block col. index */
2294: if (!mask[tmp]) {
2295: mask[tmp] = 1;
2296: if (tmp < rstart || tmp >= rend) masked2[odcount++] = tmp; /* entry in off-diag portion */
2297: else masked1[dcount++] = tmp; /* entry in diag portion */
2298: }
2299: }
2300: rowcount++;
2301: }
2302:
2303: dlens[i] = dcount; /* d_nzz[i] */
2304: odlens[i] = odcount; /* o_nzz[i] */
2306: /* zero out the mask elements we set */
2307: for (j=0; j<dcount; j++) mask[masked1[j]] = 0;
2308: for (j=0; j<odcount; j++) mask[masked2[j]] = 0;
2309: }
2310:
2311: /* create our matrix */
2312: MatCreate(comm,&A);
2313: MatSetSizes(A,m,m,PETSC_DETERMINE,PETSC_DETERMINE);
2314: MatSetType(A,type);
2315: MatMPISBAIJSetPreallocation(A,bs,0,dlens,0,odlens);
2316: MatSetOption(A,MAT_COLUMNS_SORTED);
2317:
2318: if (!rank) {
2319: PetscMalloc(maxnz*sizeof(PetscScalar),&buf);
2320: /* read in my part of the matrix numerical values */
2321: nz = procsnz[0];
2322: vals = buf;
2323: mycols = ibuf;
2324: if (size == 1) nz -= extra_rows;
2325: PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
2326: if (size == 1) for (i=0; i< extra_rows; i++) { vals[nz+i] = 1.0; }
2328: /* insert into matrix */
2329: jj = rstart*bs;
2330: for (i=0; i<m; i++) {
2331: MatSetValues(A,1,&jj,locrowlens[i],mycols,vals,INSERT_VALUES);
2332: mycols += locrowlens[i];
2333: vals += locrowlens[i];
2334: jj++;
2335: }
2337: /* read in other processors (except the last one) and ship out */
2338: for (i=1; i<size-1; i++) {
2339: nz = procsnz[i];
2340: vals = buf;
2341: PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
2342: MPI_Send(vals,nz,MPIU_SCALAR,i,A->tag,comm);
2343: }
2344: /* the last proc */
2345: if (size != 1){
2346: nz = procsnz[i] - extra_rows;
2347: vals = buf;
2348: PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
2349: for (i=0; i<extra_rows; i++) vals[nz+i] = 1.0;
2350: MPI_Send(vals,nz+extra_rows,MPIU_SCALAR,size-1,A->tag,comm);
2351: }
2352: PetscFree(procsnz);
2354: } else {
2355: /* receive numeric values */
2356: PetscMalloc(nz*sizeof(PetscScalar),&buf);
2358: /* receive message of values*/
2359: vals = buf;
2360: mycols = ibuf;
2361: MPI_Recv(vals,nz,MPIU_SCALAR,0,A->tag,comm,&status);
2362: MPI_Get_count(&status,MPIU_SCALAR,&maxnz);
2363: if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
2365: /* insert into matrix */
2366: jj = rstart*bs;
2367: for (i=0; i<m; i++) {
2368: MatSetValues_MPISBAIJ(A,1,&jj,locrowlens[i],mycols,vals,INSERT_VALUES);
2369: mycols += locrowlens[i];
2370: vals += locrowlens[i];
2371: jj++;
2372: }
2373: }
2375: PetscFree(locrowlens);
2376: PetscFree(buf);
2377: PetscFree(ibuf);
2378: PetscFree(rowners);
2379: PetscFree(dlens);
2380: PetscFree(mask);
2381: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2382: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2383: *newmat = A;
2384: return(0);
2385: }
2389: /*XXXXX@
2390: MatMPISBAIJSetHashTableFactor - Sets the factor required to compute the size of the HashTable.
2392: Input Parameters:
2393: . mat - the matrix
2394: . fact - factor
2396: Collective on Mat
2398: Level: advanced
2400: Notes:
2401: This can also be set by the command line option: -mat_use_hash_table fact
2403: .keywords: matrix, hashtable, factor, HT
2405: .seealso: MatSetOption()
2406: @XXXXX*/
2411: PetscErrorCode MatGetRowMax_MPISBAIJ(Mat A,Vec v)
2412: {
2413: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
2414: Mat_SeqBAIJ *b = (Mat_SeqBAIJ*)(a->B)->data;
2415: PetscReal atmp;
2416: PetscReal *work,*svalues,*rvalues;
2418: PetscInt i,bs,mbs,*bi,*bj,brow,j,ncols,krow,kcol,col,row,Mbs,bcol;
2419: PetscMPIInt rank,size;
2420: PetscInt *rowners_bs,dest,count,source;
2421: PetscScalar *va;
2422: MatScalar *ba;
2423: MPI_Status stat;
2426: MatGetRowMax(a->A,v);
2427: VecGetArray(v,&va);
2429: MPI_Comm_size(A->comm,&size);
2430: MPI_Comm_rank(A->comm,&rank);
2432: bs = A->rmap.bs;
2433: mbs = a->mbs;
2434: Mbs = a->Mbs;
2435: ba = b->a;
2436: bi = b->i;
2437: bj = b->j;
2439: /* find ownerships */
2440: rowners_bs = A->rmap.range;
2442: /* each proc creates an array to be distributed */
2443: PetscMalloc(bs*Mbs*sizeof(PetscReal),&work);
2444: PetscMemzero(work,bs*Mbs*sizeof(PetscReal));
2446: /* row_max for B */
2447: if (rank != size-1){
2448: for (i=0; i<mbs; i++) {
2449: ncols = bi[1] - bi[0]; bi++;
2450: brow = bs*i;
2451: for (j=0; j<ncols; j++){
2452: bcol = bs*(*bj);
2453: for (kcol=0; kcol<bs; kcol++){
2454: col = bcol + kcol; /* local col index */
2455: col += rowners_bs[rank+1]; /* global col index */
2456: for (krow=0; krow<bs; krow++){
2457: atmp = PetscAbsScalar(*ba); ba++;
2458: row = brow + krow; /* local row index */
2459: if (PetscRealPart(va[row]) < atmp) va[row] = atmp;
2460: if (work[col] < atmp) work[col] = atmp;
2461: }
2462: }
2463: bj++;
2464: }
2465: }
2467: /* send values to its owners */
2468: for (dest=rank+1; dest<size; dest++){
2469: svalues = work + rowners_bs[dest];
2470: count = rowners_bs[dest+1]-rowners_bs[dest];
2471: MPI_Send(svalues,count,MPIU_REAL,dest,rank,A->comm);
2472: }
2473: }
2474:
2475: /* receive values */
2476: if (rank){
2477: rvalues = work;
2478: count = rowners_bs[rank+1]-rowners_bs[rank];
2479: for (source=0; source<rank; source++){
2480: MPI_Recv(rvalues,count,MPIU_REAL,MPI_ANY_SOURCE,MPI_ANY_TAG,A->comm,&stat);
2481: /* process values */
2482: for (i=0; i<count; i++){
2483: if (PetscRealPart(va[i]) < rvalues[i]) va[i] = rvalues[i];
2484: }
2485: }
2486: }
2488: VecRestoreArray(v,&va);
2489: PetscFree(work);
2490: return(0);
2491: }
2495: PetscErrorCode MatRelax_MPISBAIJ(Mat matin,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
2496: {
2497: Mat_MPISBAIJ *mat = (Mat_MPISBAIJ*)matin->data;
2499: PetscInt mbs=mat->mbs,bs=matin->rmap.bs;
2500: PetscScalar *x,*b,*ptr,zero=0.0;
2501: Vec bb1;
2502:
2504: if (its <= 0 || lits <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
2505: if (bs > 1)
2506: SETERRQ(PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented");
2508: if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP){
2509: if ( flag & SOR_ZERO_INITIAL_GUESS ) {
2510: (*mat->A->ops->relax)(mat->A,bb,omega,flag,fshift,lits,lits,xx);
2511: its--;
2512: }
2514: VecDuplicate(bb,&bb1);
2515: while (its--){
2516:
2517: /* lower triangular part: slvec0b = - B^T*xx */
2518: (*mat->B->ops->multtranspose)(mat->B,xx,mat->slvec0b);
2519:
2520: /* copy xx into slvec0a */
2521: VecGetArray(mat->slvec0,&ptr);
2522: VecGetArray(xx,&x);
2523: PetscMemcpy(ptr,x,bs*mbs*sizeof(MatScalar));
2524: VecRestoreArray(mat->slvec0,&ptr);
2526: VecScale(mat->slvec0,-1.0);
2528: /* copy bb into slvec1a */
2529: VecGetArray(mat->slvec1,&ptr);
2530: VecGetArray(bb,&b);
2531: PetscMemcpy(ptr,b,bs*mbs*sizeof(MatScalar));
2532: VecRestoreArray(mat->slvec1,&ptr);
2534: /* set slvec1b = 0 */
2535: VecSet(mat->slvec1b,zero);
2537: VecScatterBegin(mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD,mat->sMvctx);
2538: VecRestoreArray(xx,&x);
2539: VecRestoreArray(bb,&b);
2540: VecScatterEnd(mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD,mat->sMvctx);
2542: /* upper triangular part: bb1 = bb1 - B*x */
2543: (*mat->B->ops->multadd)(mat->B,mat->slvec1b,mat->slvec1a,bb1);
2544:
2545: /* local diagonal sweep */
2546: (*mat->A->ops->relax)(mat->A,bb1,omega,SOR_SYMMETRIC_SWEEP,fshift,lits,lits,xx);
2547: }
2548: VecDestroy(bb1);
2549: } else {
2550: SETERRQ(PETSC_ERR_SUP,"MatSORType is not supported for SBAIJ matrix format");
2551: }
2552: return(0);
2553: }
2557: PetscErrorCode MatRelax_MPISBAIJ_2comm(Mat matin,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
2558: {
2559: Mat_MPISBAIJ *mat = (Mat_MPISBAIJ*)matin->data;
2561: Vec lvec1,bb1;
2562:
2564: if (its <= 0 || lits <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
2565: if (matin->rmap.bs > 1)
2566: SETERRQ(PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented");
2568: if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP){
2569: if ( flag & SOR_ZERO_INITIAL_GUESS ) {
2570: (*mat->A->ops->relax)(mat->A,bb,omega,flag,fshift,lits,lits,xx);
2571: its--;
2572: }
2574: VecDuplicate(mat->lvec,&lvec1);
2575: VecDuplicate(bb,&bb1);
2576: while (its--){
2577: VecScatterBegin(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,mat->Mvctx);
2578:
2579: /* lower diagonal part: bb1 = bb - B^T*xx */
2580: (*mat->B->ops->multtranspose)(mat->B,xx,lvec1);
2581: VecScale(lvec1,-1.0);
2583: VecScatterEnd(xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD,mat->Mvctx);
2584: VecCopy(bb,bb1);
2585: VecScatterBegin(lvec1,bb1,ADD_VALUES,SCATTER_REVERSE,mat->Mvctx);
2587: /* upper diagonal part: bb1 = bb1 - B*x */
2588: VecScale(mat->lvec,-1.0);
2589: (*mat->B->ops->multadd)(mat->B,mat->lvec,bb1,bb1);
2591: VecScatterEnd(lvec1,bb1,ADD_VALUES,SCATTER_REVERSE,mat->Mvctx);
2592:
2593: /* diagonal sweep */
2594: (*mat->A->ops->relax)(mat->A,bb1,omega,SOR_SYMMETRIC_SWEEP,fshift,lits,lits,xx);
2595: }
2596: VecDestroy(lvec1);
2597: VecDestroy(bb1);
2598: } else {
2599: SETERRQ(PETSC_ERR_SUP,"MatSORType is not supported for SBAIJ matrix format");
2600: }
2601: return(0);
2602: }