Actual source code: mpibaij.c
1: #include <../src/mat/impls/baij/mpi/mpibaij.h>
3: #include <petsc/private/hashseti.h>
4: #include <petscblaslapack.h>
5: #include <petscsf.h>
7: #if defined(PETSC_HAVE_HYPRE)
8: PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat,MatType,MatReuse,Mat*);
9: #endif
11: PetscErrorCode MatGetRowMaxAbs_MPIBAIJ(Mat A,Vec v,PetscInt idx[])
12: {
13: Mat_MPIBAIJ *a = (Mat_MPIBAIJ*)A->data;
14: PetscInt i,*idxb = NULL,m = A->rmap->n,bs = A->cmap->bs;
15: PetscScalar *va,*vv;
16: Vec vB,vA;
17: const PetscScalar *vb;
19: VecCreateSeq(PETSC_COMM_SELF,m,&vA);
20: MatGetRowMaxAbs(a->A,vA,idx);
22: VecGetArrayWrite(vA,&va);
23: if (idx) {
24: for (i=0; i<m; i++) {
25: if (PetscAbsScalar(va[i])) idx[i] += A->cmap->rstart;
26: }
27: }
29: VecCreateSeq(PETSC_COMM_SELF,m,&vB);
30: PetscMalloc1(m,&idxb);
31: MatGetRowMaxAbs(a->B,vB,idxb);
33: VecGetArrayWrite(v,&vv);
34: VecGetArrayRead(vB,&vb);
35: for (i=0; i<m; i++) {
36: if (PetscAbsScalar(va[i]) < PetscAbsScalar(vb[i])) {
37: vv[i] = vb[i];
38: if (idx) idx[i] = bs*a->garray[idxb[i]/bs] + (idxb[i] % bs);
39: } else {
40: vv[i] = va[i];
41: if (idx && PetscAbsScalar(va[i]) == PetscAbsScalar(vb[i]) && idxb[i] != -1 && idx[i] > bs*a->garray[idxb[i]/bs] + (idxb[i] % bs))
42: idx[i] = bs*a->garray[idxb[i]/bs] + (idxb[i] % bs);
43: }
44: }
45: VecRestoreArrayWrite(vA,&vv);
46: VecRestoreArrayWrite(vA,&va);
47: VecRestoreArrayRead(vB,&vb);
48: PetscFree(idxb);
49: VecDestroy(&vA);
50: VecDestroy(&vB);
51: return 0;
52: }
54: PetscErrorCode MatStoreValues_MPIBAIJ(Mat mat)
55: {
56: Mat_MPIBAIJ *aij = (Mat_MPIBAIJ*)mat->data;
58: MatStoreValues(aij->A);
59: MatStoreValues(aij->B);
60: return 0;
61: }
63: PetscErrorCode MatRetrieveValues_MPIBAIJ(Mat mat)
64: {
65: Mat_MPIBAIJ *aij = (Mat_MPIBAIJ*)mat->data;
67: MatRetrieveValues(aij->A);
68: MatRetrieveValues(aij->B);
69: return 0;
70: }
72: /*
73: Local utility routine that creates a mapping from the global column
74: number to the local number in the off-diagonal part of the local
75: storage of the matrix. This is done in a non scalable way since the
76: length of colmap equals the global matrix length.
77: */
78: PetscErrorCode MatCreateColmap_MPIBAIJ_Private(Mat mat)
79: {
80: Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)mat->data;
81: Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)baij->B->data;
82: PetscInt nbs = B->nbs,i,bs=mat->rmap->bs;
84: #if defined(PETSC_USE_CTABLE)
85: PetscTableCreate(baij->nbs,baij->Nbs+1,&baij->colmap);
86: for (i=0; i<nbs; i++) {
87: PetscTableAdd(baij->colmap,baij->garray[i]+1,i*bs+1,INSERT_VALUES);
88: }
89: #else
90: PetscCalloc1(baij->Nbs+1,&baij->colmap);
91: PetscLogObjectMemory((PetscObject)mat,baij->Nbs*sizeof(PetscInt));
92: for (i=0; i<nbs; i++) baij->colmap[baij->garray[i]] = i*bs+1;
93: #endif
94: return 0;
95: }
97: #define MatSetValues_SeqBAIJ_A_Private(row,col,value,addv,orow,ocol) \
98: { \
99: brow = row/bs; \
100: rp = aj + ai[brow]; ap = aa + bs2*ai[brow]; \
101: rmax = aimax[brow]; nrow = ailen[brow]; \
102: bcol = col/bs; \
103: ridx = row % bs; cidx = col % bs; \
104: low = 0; high = nrow; \
105: while (high-low > 3) { \
106: t = (low+high)/2; \
107: if (rp[t] > bcol) high = t; \
108: else low = t; \
109: } \
110: for (_i=low; _i<high; _i++) { \
111: if (rp[_i] > bcol) break; \
112: if (rp[_i] == bcol) { \
113: bap = ap + bs2*_i + bs*cidx + ridx; \
114: if (addv == ADD_VALUES) *bap += value; \
115: else *bap = value; \
116: goto a_noinsert; \
117: } \
118: } \
119: if (a->nonew == 1) goto a_noinsert; \
121: MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,aimax,a->nonew,MatScalar); \
122: N = nrow++ - 1; \
123: /* shift up all the later entries in this row */ \
124: PetscArraymove(rp+_i+1,rp+_i,N-_i+1);\
125: PetscArraymove(ap+bs2*(_i+1),ap+bs2*_i,bs2*(N-_i+1)); \
126: PetscArrayzero(ap+bs2*_i,bs2); \
127: rp[_i] = bcol; \
128: ap[bs2*_i + bs*cidx + ridx] = value; \
129: a_noinsert:; \
130: ailen[brow] = nrow; \
131: }
133: #define MatSetValues_SeqBAIJ_B_Private(row,col,value,addv,orow,ocol) \
134: { \
135: brow = row/bs; \
136: rp = bj + bi[brow]; ap = ba + bs2*bi[brow]; \
137: rmax = bimax[brow]; nrow = bilen[brow]; \
138: bcol = col/bs; \
139: ridx = row % bs; cidx = col % bs; \
140: low = 0; high = nrow; \
141: while (high-low > 3) { \
142: t = (low+high)/2; \
143: if (rp[t] > bcol) high = t; \
144: else low = t; \
145: } \
146: for (_i=low; _i<high; _i++) { \
147: if (rp[_i] > bcol) break; \
148: if (rp[_i] == bcol) { \
149: bap = ap + bs2*_i + bs*cidx + ridx; \
150: if (addv == ADD_VALUES) *bap += value; \
151: else *bap = value; \
152: goto b_noinsert; \
153: } \
154: } \
155: if (b->nonew == 1) goto b_noinsert; \
157: MatSeqXAIJReallocateAIJ(B,b->mbs,bs2,nrow,brow,bcol,rmax,ba,bi,bj,rp,ap,bimax,b->nonew,MatScalar); \
158: N = nrow++ - 1; \
159: /* shift up all the later entries in this row */ \
160: PetscArraymove(rp+_i+1,rp+_i,N-_i+1);\
161: PetscArraymove(ap+bs2*(_i+1),ap+bs2*_i,bs2*(N-_i+1));\
162: PetscArrayzero(ap+bs2*_i,bs2); \
163: rp[_i] = bcol; \
164: ap[bs2*_i + bs*cidx + ridx] = value; \
165: b_noinsert:; \
166: bilen[brow] = nrow; \
167: }
169: PetscErrorCode MatSetValues_MPIBAIJ(Mat mat,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode addv)
170: {
171: Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)mat->data;
172: MatScalar value;
173: PetscBool roworiented = baij->roworiented;
174: PetscInt i,j,row,col;
175: PetscInt rstart_orig=mat->rmap->rstart;
176: PetscInt rend_orig =mat->rmap->rend,cstart_orig=mat->cmap->rstart;
177: PetscInt cend_orig =mat->cmap->rend,bs=mat->rmap->bs;
179: /* Some Variables required in the macro */
180: Mat A = baij->A;
181: Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)(A)->data;
182: PetscInt *aimax=a->imax,*ai=a->i,*ailen=a->ilen,*aj=a->j;
183: MatScalar *aa =a->a;
185: Mat B = baij->B;
186: Mat_SeqBAIJ *b = (Mat_SeqBAIJ*)(B)->data;
187: PetscInt *bimax=b->imax,*bi=b->i,*bilen=b->ilen,*bj=b->j;
188: MatScalar *ba =b->a;
190: PetscInt *rp,ii,nrow,_i,rmax,N,brow,bcol;
191: PetscInt low,high,t,ridx,cidx,bs2=a->bs2;
192: MatScalar *ap,*bap;
194: for (i=0; i<m; i++) {
195: if (im[i] < 0) continue;
197: if (im[i] >= rstart_orig && im[i] < rend_orig) {
198: row = im[i] - rstart_orig;
199: for (j=0; j<n; j++) {
200: if (in[j] >= cstart_orig && in[j] < cend_orig) {
201: col = in[j] - cstart_orig;
202: if (roworiented) value = v[i*n+j];
203: else value = v[i+j*m];
204: MatSetValues_SeqBAIJ_A_Private(row,col,value,addv,im[i],in[j]);
205: } else if (in[j] < 0) continue;
207: else {
208: if (mat->was_assembled) {
209: if (!baij->colmap) {
210: MatCreateColmap_MPIBAIJ_Private(mat);
211: }
212: #if defined(PETSC_USE_CTABLE)
213: PetscTableFind(baij->colmap,in[j]/bs + 1,&col);
214: col = col - 1;
215: #else
216: col = baij->colmap[in[j]/bs] - 1;
217: #endif
218: if (col < 0 && !((Mat_SeqBAIJ*)(baij->B->data))->nonew) {
219: MatDisAssemble_MPIBAIJ(mat);
220: col = in[j];
221: /* Reinitialize the variables required by MatSetValues_SeqBAIJ_B_Private() */
222: B = baij->B;
223: b = (Mat_SeqBAIJ*)(B)->data;
224: bimax=b->imax;bi=b->i;bilen=b->ilen;bj=b->j;
225: ba =b->a;
227: else col += in[j]%bs;
228: } else col = in[j];
229: if (roworiented) value = v[i*n+j];
230: else value = v[i+j*m];
231: MatSetValues_SeqBAIJ_B_Private(row,col,value,addv,im[i],in[j]);
232: /* MatSetValues_SeqBAIJ(baij->B,1,&row,1,&col,&value,addv); */
233: }
234: }
235: } else {
237: if (!baij->donotstash) {
238: mat->assembled = PETSC_FALSE;
239: if (roworiented) {
240: MatStashValuesRow_Private(&mat->stash,im[i],n,in,v+i*n,PETSC_FALSE);
241: } else {
242: MatStashValuesCol_Private(&mat->stash,im[i],n,in,v+i,m,PETSC_FALSE);
243: }
244: }
245: }
246: }
247: return 0;
248: }
250: static inline PetscErrorCode MatSetValuesBlocked_SeqBAIJ_Inlined(Mat A,PetscInt row,PetscInt col,const PetscScalar v[],InsertMode is,PetscInt orow,PetscInt ocol)
251: {
252: Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
253: PetscInt *rp,low,high,t,ii,jj,nrow,i,rmax,N;
254: PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen;
255: PetscInt *aj =a->j,nonew=a->nonew,bs2=a->bs2,bs=A->rmap->bs;
256: PetscBool roworiented=a->roworiented;
257: const PetscScalar *value = v;
258: MatScalar *ap,*aa = a->a,*bap;
260: rp = aj + ai[row];
261: ap = aa + bs2*ai[row];
262: rmax = imax[row];
263: nrow = ailen[row];
264: value = v;
265: low = 0;
266: high = nrow;
267: while (high-low > 7) {
268: t = (low+high)/2;
269: if (rp[t] > col) high = t;
270: else low = t;
271: }
272: for (i=low; i<high; i++) {
273: if (rp[i] > col) break;
274: if (rp[i] == col) {
275: bap = ap + bs2*i;
276: if (roworiented) {
277: if (is == ADD_VALUES) {
278: for (ii=0; ii<bs; ii++) {
279: for (jj=ii; jj<bs2; jj+=bs) {
280: bap[jj] += *value++;
281: }
282: }
283: } else {
284: for (ii=0; ii<bs; ii++) {
285: for (jj=ii; jj<bs2; jj+=bs) {
286: bap[jj] = *value++;
287: }
288: }
289: }
290: } else {
291: if (is == ADD_VALUES) {
292: for (ii=0; ii<bs; ii++,value+=bs) {
293: for (jj=0; jj<bs; jj++) {
294: bap[jj] += value[jj];
295: }
296: bap += bs;
297: }
298: } else {
299: for (ii=0; ii<bs; ii++,value+=bs) {
300: for (jj=0; jj<bs; jj++) {
301: bap[jj] = value[jj];
302: }
303: bap += bs;
304: }
305: }
306: }
307: goto noinsert2;
308: }
309: }
310: if (nonew == 1) goto noinsert2;
312: MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar);
313: N = nrow++ - 1; high++;
314: /* shift up all the later entries in this row */
315: PetscArraymove(rp+i+1,rp+i,N-i+1);
316: PetscArraymove(ap+bs2*(i+1),ap+bs2*i,bs2*(N-i+1));
317: rp[i] = col;
318: bap = ap + bs2*i;
319: if (roworiented) {
320: for (ii=0; ii<bs; ii++) {
321: for (jj=ii; jj<bs2; jj+=bs) {
322: bap[jj] = *value++;
323: }
324: }
325: } else {
326: for (ii=0; ii<bs; ii++) {
327: for (jj=0; jj<bs; jj++) {
328: *bap++ = *value++;
329: }
330: }
331: }
332: noinsert2:;
333: ailen[row] = nrow;
334: return 0;
335: }
337: /*
338: This routine should be optimized so that the block copy at ** Here a copy is required ** below is not needed
339: by passing additional stride information into the MatSetValuesBlocked_SeqBAIJ_Inlined() routine
340: */
341: PetscErrorCode MatSetValuesBlocked_MPIBAIJ(Mat mat,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode addv)
342: {
343: Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)mat->data;
344: const PetscScalar *value;
345: MatScalar *barray = baij->barray;
346: PetscBool roworiented = baij->roworiented;
347: PetscInt i,j,ii,jj,row,col,rstart=baij->rstartbs;
348: PetscInt rend=baij->rendbs,cstart=baij->cstartbs,stepval;
349: PetscInt cend=baij->cendbs,bs=mat->rmap->bs,bs2=baij->bs2;
351: if (!barray) {
352: PetscMalloc1(bs2,&barray);
353: baij->barray = barray;
354: }
356: if (roworiented) stepval = (n-1)*bs;
357: else stepval = (m-1)*bs;
359: for (i=0; i<m; i++) {
360: if (im[i] < 0) continue;
362: if (im[i] >= rstart && im[i] < rend) {
363: row = im[i] - rstart;
364: for (j=0; j<n; j++) {
365: /* If NumCol = 1 then a copy is not required */
366: if ((roworiented) && (n == 1)) {
367: barray = (MatScalar*)v + i*bs2;
368: } else if ((!roworiented) && (m == 1)) {
369: barray = (MatScalar*)v + j*bs2;
370: } else { /* Here a copy is required */
371: if (roworiented) {
372: value = v + (i*(stepval+bs) + j)*bs;
373: } else {
374: value = v + (j*(stepval+bs) + i)*bs;
375: }
376: for (ii=0; ii<bs; ii++,value+=bs+stepval) {
377: for (jj=0; jj<bs; jj++) barray[jj] = value[jj];
378: barray += bs;
379: }
380: barray -= bs2;
381: }
383: if (in[j] >= cstart && in[j] < cend) {
384: col = in[j] - cstart;
385: MatSetValuesBlocked_SeqBAIJ_Inlined(baij->A,row,col,barray,addv,im[i],in[j]);
386: } else if (in[j] < 0) continue;
388: else {
389: if (mat->was_assembled) {
390: if (!baij->colmap) {
391: MatCreateColmap_MPIBAIJ_Private(mat);
392: }
394: #if defined(PETSC_USE_DEBUG)
395: #if defined(PETSC_USE_CTABLE)
396: { PetscInt data;
397: PetscTableFind(baij->colmap,in[j]+1,&data);
399: }
400: #else
402: #endif
403: #endif
404: #if defined(PETSC_USE_CTABLE)
405: PetscTableFind(baij->colmap,in[j]+1,&col);
406: col = (col - 1)/bs;
407: #else
408: col = (baij->colmap[in[j]] - 1)/bs;
409: #endif
410: if (col < 0 && !((Mat_SeqBAIJ*)(baij->B->data))->nonew) {
411: MatDisAssemble_MPIBAIJ(mat);
412: col = in[j];
414: } else col = in[j];
415: MatSetValuesBlocked_SeqBAIJ_Inlined(baij->B,row,col,barray,addv,im[i],in[j]);
416: }
417: }
418: } else {
420: if (!baij->donotstash) {
421: if (roworiented) {
422: MatStashValuesRowBlocked_Private(&mat->bstash,im[i],n,in,v,m,n,i);
423: } else {
424: MatStashValuesColBlocked_Private(&mat->bstash,im[i],n,in,v,m,n,i);
425: }
426: }
427: }
428: }
429: return 0;
430: }
432: #define HASH_KEY 0.6180339887
433: #define HASH(size,key,tmp) (tmp = (key)*HASH_KEY,(PetscInt)((size)*(tmp-(PetscInt)tmp)))
434: /* #define HASH(size,key) ((PetscInt)((size)*fmod(((key)*HASH_KEY),1))) */
435: /* #define HASH(size,key,tmp) ((PetscInt)((size)*fmod(((key)*HASH_KEY),1))) */
436: PetscErrorCode MatSetValues_MPIBAIJ_HT(Mat mat,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode addv)
437: {
438: Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)mat->data;
439: PetscBool roworiented = baij->roworiented;
440: PetscInt i,j,row,col;
441: PetscInt rstart_orig=mat->rmap->rstart;
442: PetscInt rend_orig =mat->rmap->rend,Nbs=baij->Nbs;
443: PetscInt h1,key,size=baij->ht_size,bs=mat->rmap->bs,*HT=baij->ht,idx;
444: PetscReal tmp;
445: MatScalar **HD = baij->hd,value;
446: PetscInt total_ct=baij->ht_total_ct,insert_ct=baij->ht_insert_ct;
448: for (i=0; i<m; i++) {
449: if (PetscDefined(USE_DEBUG)) {
452: }
453: row = im[i];
454: if (row >= rstart_orig && row < rend_orig) {
455: for (j=0; j<n; j++) {
456: col = in[j];
457: if (roworiented) value = v[i*n+j];
458: else value = v[i+j*m];
459: /* Look up PetscInto the Hash Table */
460: key = (row/bs)*Nbs+(col/bs)+1;
461: h1 = HASH(size,key,tmp);
463: idx = h1;
464: if (PetscDefined(USE_DEBUG)) {
465: insert_ct++;
466: total_ct++;
467: if (HT[idx] != key) {
468: for (idx=h1; (idx<size) && (HT[idx]!=key); idx++,total_ct++) ;
469: if (idx == size) {
470: for (idx=0; (idx<h1) && (HT[idx]!=key); idx++,total_ct++) ;
472: }
473: }
474: } else if (HT[idx] != key) {
475: for (idx=h1; (idx<size) && (HT[idx]!=key); idx++) ;
476: if (idx == size) {
477: for (idx=0; (idx<h1) && (HT[idx]!=key); idx++) ;
479: }
480: }
481: /* A HASH table entry is found, so insert the values at the correct address */
482: if (addv == ADD_VALUES) *(HD[idx]+ (col % bs)*bs + (row % bs)) += value;
483: else *(HD[idx]+ (col % bs)*bs + (row % bs)) = value;
484: }
485: } else if (!baij->donotstash) {
486: if (roworiented) {
487: MatStashValuesRow_Private(&mat->stash,im[i],n,in,v+i*n,PETSC_FALSE);
488: } else {
489: MatStashValuesCol_Private(&mat->stash,im[i],n,in,v+i,m,PETSC_FALSE);
490: }
491: }
492: }
493: if (PetscDefined(USE_DEBUG)) {
494: baij->ht_total_ct += total_ct;
495: baij->ht_insert_ct += insert_ct;
496: }
497: return 0;
498: }
500: PetscErrorCode MatSetValuesBlocked_MPIBAIJ_HT(Mat mat,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode addv)
501: {
502: Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)mat->data;
503: PetscBool roworiented = baij->roworiented;
504: PetscInt i,j,ii,jj,row,col;
505: PetscInt rstart=baij->rstartbs;
506: PetscInt rend =mat->rmap->rend,stepval,bs=mat->rmap->bs,bs2=baij->bs2,nbs2=n*bs2;
507: PetscInt h1,key,size=baij->ht_size,idx,*HT=baij->ht,Nbs=baij->Nbs;
508: PetscReal tmp;
509: MatScalar **HD = baij->hd,*baij_a;
510: const PetscScalar *v_t,*value;
511: PetscInt total_ct=baij->ht_total_ct,insert_ct=baij->ht_insert_ct;
513: if (roworiented) stepval = (n-1)*bs;
514: else stepval = (m-1)*bs;
516: for (i=0; i<m; i++) {
517: if (PetscDefined(USE_DEBUG)) {
520: }
521: row = im[i];
522: v_t = v + i*nbs2;
523: if (row >= rstart && row < rend) {
524: for (j=0; j<n; j++) {
525: col = in[j];
527: /* Look up into the Hash Table */
528: key = row*Nbs+col+1;
529: h1 = HASH(size,key,tmp);
531: idx = h1;
532: if (PetscDefined(USE_DEBUG)) {
533: total_ct++;
534: insert_ct++;
535: if (HT[idx] != key) {
536: for (idx=h1; (idx<size) && (HT[idx]!=key); idx++,total_ct++) ;
537: if (idx == size) {
538: for (idx=0; (idx<h1) && (HT[idx]!=key); idx++,total_ct++) ;
540: }
541: }
542: } else if (HT[idx] != key) {
543: for (idx=h1; (idx<size) && (HT[idx]!=key); idx++) ;
544: if (idx == size) {
545: for (idx=0; (idx<h1) && (HT[idx]!=key); idx++) ;
547: }
548: }
549: baij_a = HD[idx];
550: if (roworiented) {
551: /*value = v + i*(stepval+bs)*bs + j*bs;*/
552: /* value = v + (i*(stepval+bs)+j)*bs; */
553: value = v_t;
554: v_t += bs;
555: if (addv == ADD_VALUES) {
556: for (ii=0; ii<bs; ii++,value+=stepval) {
557: for (jj=ii; jj<bs2; jj+=bs) {
558: baij_a[jj] += *value++;
559: }
560: }
561: } else {
562: for (ii=0; ii<bs; ii++,value+=stepval) {
563: for (jj=ii; jj<bs2; jj+=bs) {
564: baij_a[jj] = *value++;
565: }
566: }
567: }
568: } else {
569: value = v + j*(stepval+bs)*bs + i*bs;
570: if (addv == ADD_VALUES) {
571: for (ii=0; ii<bs; ii++,value+=stepval,baij_a+=bs) {
572: for (jj=0; jj<bs; jj++) {
573: baij_a[jj] += *value++;
574: }
575: }
576: } else {
577: for (ii=0; ii<bs; ii++,value+=stepval,baij_a+=bs) {
578: for (jj=0; jj<bs; jj++) {
579: baij_a[jj] = *value++;
580: }
581: }
582: }
583: }
584: }
585: } else {
586: if (!baij->donotstash) {
587: if (roworiented) {
588: MatStashValuesRowBlocked_Private(&mat->bstash,im[i],n,in,v,m,n,i);
589: } else {
590: MatStashValuesColBlocked_Private(&mat->bstash,im[i],n,in,v,m,n,i);
591: }
592: }
593: }
594: }
595: if (PetscDefined(USE_DEBUG)) {
596: baij->ht_total_ct += total_ct;
597: baij->ht_insert_ct += insert_ct;
598: }
599: return 0;
600: }
602: PetscErrorCode MatGetValues_MPIBAIJ(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[])
603: {
604: Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)mat->data;
605: PetscInt bs = mat->rmap->bs,i,j,bsrstart = mat->rmap->rstart,bsrend = mat->rmap->rend;
606: PetscInt bscstart = mat->cmap->rstart,bscend = mat->cmap->rend,row,col,data;
608: for (i=0; i<m; i++) {
609: if (idxm[i] < 0) continue; /* negative row */
611: if (idxm[i] >= bsrstart && idxm[i] < bsrend) {
612: row = idxm[i] - bsrstart;
613: for (j=0; j<n; j++) {
614: if (idxn[j] < 0) continue; /* negative column */
616: if (idxn[j] >= bscstart && idxn[j] < bscend) {
617: col = idxn[j] - bscstart;
618: MatGetValues_SeqBAIJ(baij->A,1,&row,1,&col,v+i*n+j);
619: } else {
620: if (!baij->colmap) {
621: MatCreateColmap_MPIBAIJ_Private(mat);
622: }
623: #if defined(PETSC_USE_CTABLE)
624: PetscTableFind(baij->colmap,idxn[j]/bs+1,&data);
625: data--;
626: #else
627: data = baij->colmap[idxn[j]/bs]-1;
628: #endif
629: if ((data < 0) || (baij->garray[data/bs] != idxn[j]/bs)) *(v+i*n+j) = 0.0;
630: else {
631: col = data + idxn[j]%bs;
632: MatGetValues_SeqBAIJ(baij->B,1,&row,1,&col,v+i*n+j);
633: }
634: }
635: }
636: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only local values currently supported");
637: }
638: return 0;
639: }
641: PetscErrorCode MatNorm_MPIBAIJ(Mat mat,NormType type,PetscReal *nrm)
642: {
643: Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)mat->data;
644: Mat_SeqBAIJ *amat = (Mat_SeqBAIJ*)baij->A->data,*bmat = (Mat_SeqBAIJ*)baij->B->data;
645: PetscInt i,j,bs2=baij->bs2,bs=baij->A->rmap->bs,nz,row,col;
646: PetscReal sum = 0.0;
647: MatScalar *v;
649: if (baij->size == 1) {
650: MatNorm(baij->A,type,nrm);
651: } else {
652: if (type == NORM_FROBENIUS) {
653: v = amat->a;
654: nz = amat->nz*bs2;
655: for (i=0; i<nz; i++) {
656: sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
657: }
658: v = bmat->a;
659: nz = bmat->nz*bs2;
660: for (i=0; i<nz; i++) {
661: sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
662: }
663: MPIU_Allreduce(&sum,nrm,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)mat));
664: *nrm = PetscSqrtReal(*nrm);
665: } else if (type == NORM_1) { /* max column sum */
666: PetscReal *tmp,*tmp2;
667: PetscInt *jj,*garray=baij->garray,cstart=baij->rstartbs;
668: PetscCalloc1(mat->cmap->N,&tmp);
669: PetscMalloc1(mat->cmap->N,&tmp2);
670: v = amat->a; jj = amat->j;
671: for (i=0; i<amat->nz; i++) {
672: for (j=0; j<bs; j++) {
673: col = bs*(cstart + *jj) + j; /* column index */
674: for (row=0; row<bs; row++) {
675: tmp[col] += PetscAbsScalar(*v); v++;
676: }
677: }
678: jj++;
679: }
680: v = bmat->a; jj = bmat->j;
681: for (i=0; i<bmat->nz; i++) {
682: for (j=0; j<bs; j++) {
683: col = bs*garray[*jj] + j;
684: for (row=0; row<bs; row++) {
685: tmp[col] += PetscAbsScalar(*v); v++;
686: }
687: }
688: jj++;
689: }
690: MPIU_Allreduce(tmp,tmp2,mat->cmap->N,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)mat));
691: *nrm = 0.0;
692: for (j=0; j<mat->cmap->N; j++) {
693: if (tmp2[j] > *nrm) *nrm = tmp2[j];
694: }
695: PetscFree(tmp);
696: PetscFree(tmp2);
697: } else if (type == NORM_INFINITY) { /* max row sum */
698: PetscReal *sums;
699: PetscMalloc1(bs,&sums);
700: sum = 0.0;
701: for (j=0; j<amat->mbs; j++) {
702: for (row=0; row<bs; row++) sums[row] = 0.0;
703: v = amat->a + bs2*amat->i[j];
704: nz = amat->i[j+1]-amat->i[j];
705: for (i=0; i<nz; i++) {
706: for (col=0; col<bs; col++) {
707: for (row=0; row<bs; row++) {
708: sums[row] += PetscAbsScalar(*v); v++;
709: }
710: }
711: }
712: v = bmat->a + bs2*bmat->i[j];
713: nz = bmat->i[j+1]-bmat->i[j];
714: for (i=0; i<nz; i++) {
715: for (col=0; col<bs; col++) {
716: for (row=0; row<bs; row++) {
717: sums[row] += PetscAbsScalar(*v); v++;
718: }
719: }
720: }
721: for (row=0; row<bs; row++) {
722: if (sums[row] > sum) sum = sums[row];
723: }
724: }
725: MPIU_Allreduce(&sum,nrm,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)mat));
726: PetscFree(sums);
727: } else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"No support for this norm yet");
728: }
729: return 0;
730: }
732: /*
733: Creates the hash table, and sets the table
734: This table is created only once.
735: If new entried need to be added to the matrix
736: then the hash table has to be destroyed and
737: recreated.
738: */
739: PetscErrorCode MatCreateHashTable_MPIBAIJ_Private(Mat mat,PetscReal factor)
740: {
741: Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)mat->data;
742: Mat A = baij->A,B=baij->B;
743: Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)B->data;
744: PetscInt i,j,k,nz=a->nz+b->nz,h1,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
745: PetscInt ht_size,bs2=baij->bs2,rstart=baij->rstartbs;
746: PetscInt cstart=baij->cstartbs,*garray=baij->garray,row,col,Nbs=baij->Nbs;
747: PetscInt *HT,key;
748: MatScalar **HD;
749: PetscReal tmp;
750: #if defined(PETSC_USE_INFO)
751: PetscInt ct=0,max=0;
752: #endif
754: if (baij->ht) return 0;
756: baij->ht_size = (PetscInt)(factor*nz);
757: ht_size = baij->ht_size;
759: /* Allocate Memory for Hash Table */
760: PetscCalloc2(ht_size,&baij->hd,ht_size,&baij->ht);
761: HD = baij->hd;
762: HT = baij->ht;
764: /* Loop Over A */
765: for (i=0; i<a->mbs; i++) {
766: for (j=ai[i]; j<ai[i+1]; j++) {
767: row = i+rstart;
768: col = aj[j]+cstart;
770: key = row*Nbs + col + 1;
771: h1 = HASH(ht_size,key,tmp);
772: for (k=0; k<ht_size; k++) {
773: if (!HT[(h1+k)%ht_size]) {
774: HT[(h1+k)%ht_size] = key;
775: HD[(h1+k)%ht_size] = a->a + j*bs2;
776: break;
777: #if defined(PETSC_USE_INFO)
778: } else {
779: ct++;
780: #endif
781: }
782: }
783: #if defined(PETSC_USE_INFO)
784: if (k> max) max = k;
785: #endif
786: }
787: }
788: /* Loop Over B */
789: for (i=0; i<b->mbs; i++) {
790: for (j=bi[i]; j<bi[i+1]; j++) {
791: row = i+rstart;
792: col = garray[bj[j]];
793: key = row*Nbs + col + 1;
794: h1 = HASH(ht_size,key,tmp);
795: for (k=0; k<ht_size; k++) {
796: if (!HT[(h1+k)%ht_size]) {
797: HT[(h1+k)%ht_size] = key;
798: HD[(h1+k)%ht_size] = b->a + j*bs2;
799: break;
800: #if defined(PETSC_USE_INFO)
801: } else {
802: ct++;
803: #endif
804: }
805: }
806: #if defined(PETSC_USE_INFO)
807: if (k> max) max = k;
808: #endif
809: }
810: }
812: /* Print Summary */
813: #if defined(PETSC_USE_INFO)
814: for (i=0,j=0; i<ht_size; i++) {
815: if (HT[i]) j++;
816: }
817: PetscInfo(mat,"Average Search = %5.2g,max search = %" PetscInt_FMT "\n",(!j) ? (double)0.0:(double)(((PetscReal)(ct+j))/(double)j),max);
818: #endif
819: return 0;
820: }
822: PetscErrorCode MatAssemblyBegin_MPIBAIJ(Mat mat,MatAssemblyType mode)
823: {
824: Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)mat->data;
825: PetscInt nstash,reallocs;
827: if (baij->donotstash || mat->nooffprocentries) return 0;
829: MatStashScatterBegin_Private(mat,&mat->stash,mat->rmap->range);
830: MatStashScatterBegin_Private(mat,&mat->bstash,baij->rangebs);
831: MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);
832: PetscInfo(mat,"Stash has %" PetscInt_FMT " entries,uses %" PetscInt_FMT " mallocs.\n",nstash,reallocs);
833: MatStashGetInfo_Private(&mat->bstash,&nstash,&reallocs);
834: PetscInfo(mat,"Block-Stash has %" PetscInt_FMT " entries, uses %" PetscInt_FMT " mallocs.\n",nstash,reallocs);
835: return 0;
836: }
838: PetscErrorCode MatAssemblyEnd_MPIBAIJ(Mat mat,MatAssemblyType mode)
839: {
840: Mat_MPIBAIJ *baij=(Mat_MPIBAIJ*)mat->data;
841: Mat_SeqBAIJ *a =(Mat_SeqBAIJ*)baij->A->data;
842: PetscInt i,j,rstart,ncols,flg,bs2=baij->bs2;
843: PetscInt *row,*col;
844: PetscBool r1,r2,r3,other_disassembled;
845: MatScalar *val;
846: PetscMPIInt n;
848: /* do not use 'b=(Mat_SeqBAIJ*)baij->B->data' as B can be reset in disassembly */
849: if (!baij->donotstash && !mat->nooffprocentries) {
850: while (1) {
851: MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);
852: if (!flg) break;
854: for (i=0; i<n;) {
855: /* Now identify the consecutive vals belonging to the same row */
856: for (j=i,rstart=row[j]; j<n; j++) {
857: if (row[j] != rstart) break;
858: }
859: if (j < n) ncols = j-i;
860: else ncols = n-i;
861: /* Now assemble all these values with a single function call */
862: MatSetValues_MPIBAIJ(mat,1,row+i,ncols,col+i,val+i,mat->insertmode);
863: i = j;
864: }
865: }
866: MatStashScatterEnd_Private(&mat->stash);
867: /* Now process the block-stash. Since the values are stashed column-oriented,
868: set the roworiented flag to column oriented, and after MatSetValues()
869: restore the original flags */
870: r1 = baij->roworiented;
871: r2 = a->roworiented;
872: r3 = ((Mat_SeqBAIJ*)baij->B->data)->roworiented;
874: baij->roworiented = PETSC_FALSE;
875: a->roworiented = PETSC_FALSE;
877: (((Mat_SeqBAIJ*)baij->B->data))->roworiented = PETSC_FALSE; /* b->roworiented */
878: while (1) {
879: MatStashScatterGetMesg_Private(&mat->bstash,&n,&row,&col,&val,&flg);
880: if (!flg) break;
882: for (i=0; i<n;) {
883: /* Now identify the consecutive vals belonging to the same row */
884: for (j=i,rstart=row[j]; j<n; j++) {
885: if (row[j] != rstart) break;
886: }
887: if (j < n) ncols = j-i;
888: else ncols = n-i;
889: MatSetValuesBlocked_MPIBAIJ(mat,1,row+i,ncols,col+i,val+i*bs2,mat->insertmode);
890: i = j;
891: }
892: }
893: MatStashScatterEnd_Private(&mat->bstash);
895: baij->roworiented = r1;
896: a->roworiented = r2;
898: ((Mat_SeqBAIJ*)baij->B->data)->roworiented = r3; /* b->roworiented */
899: }
901: MatAssemblyBegin(baij->A,mode);
902: MatAssemblyEnd(baij->A,mode);
904: /* determine if any processor has disassembled, if so we must
905: also disassemble ourselves, in order that we may reassemble. */
906: /*
907: if nonzero structure of submatrix B cannot change then we know that
908: no processor disassembled thus we can skip this stuff
909: */
910: if (!((Mat_SeqBAIJ*)baij->B->data)->nonew) {
911: MPIU_Allreduce(&mat->was_assembled,&other_disassembled,1,MPIU_BOOL,MPI_PROD,PetscObjectComm((PetscObject)mat));
912: if (mat->was_assembled && !other_disassembled) {
913: MatDisAssemble_MPIBAIJ(mat);
914: }
915: }
917: if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
918: MatSetUpMultiply_MPIBAIJ(mat);
919: }
920: MatAssemblyBegin(baij->B,mode);
921: MatAssemblyEnd(baij->B,mode);
923: #if defined(PETSC_USE_INFO)
924: if (baij->ht && mode== MAT_FINAL_ASSEMBLY) {
925: PetscInfo(mat,"Average Hash Table Search in MatSetValues = %5.2f\n",(double)((PetscReal)baij->ht_total_ct)/baij->ht_insert_ct);
927: baij->ht_total_ct = 0;
928: baij->ht_insert_ct = 0;
929: }
930: #endif
931: if (baij->ht_flag && !baij->ht && mode == MAT_FINAL_ASSEMBLY) {
932: MatCreateHashTable_MPIBAIJ_Private(mat,baij->ht_fact);
934: mat->ops->setvalues = MatSetValues_MPIBAIJ_HT;
935: mat->ops->setvaluesblocked = MatSetValuesBlocked_MPIBAIJ_HT;
936: }
938: PetscFree2(baij->rowvalues,baij->rowindices);
940: baij->rowvalues = NULL;
942: /* if no new nonzero locations are allowed in matrix then only set the matrix state the first time through */
943: if ((!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) || !((Mat_SeqBAIJ*)(baij->A->data))->nonew) {
944: PetscObjectState state = baij->A->nonzerostate + baij->B->nonzerostate;
945: MPIU_Allreduce(&state,&mat->nonzerostate,1,MPIU_INT64,MPI_SUM,PetscObjectComm((PetscObject)mat));
946: }
947: return 0;
948: }
950: extern PetscErrorCode MatView_SeqBAIJ(Mat,PetscViewer);
951: #include <petscdraw.h>
952: static PetscErrorCode MatView_MPIBAIJ_ASCIIorDraworSocket(Mat mat,PetscViewer viewer)
953: {
954: Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)mat->data;
955: PetscMPIInt rank = baij->rank;
956: PetscInt bs = mat->rmap->bs;
957: PetscBool iascii,isdraw;
958: PetscViewer sviewer;
959: PetscViewerFormat format;
961: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
962: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
963: if (iascii) {
964: PetscViewerGetFormat(viewer,&format);
965: if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
966: MatInfo info;
967: MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&rank);
968: MatGetInfo(mat,MAT_LOCAL,&info);
969: PetscViewerASCIIPushSynchronized(viewer);
970: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Local rows %" PetscInt_FMT " nz %" PetscInt_FMT " nz alloced %" PetscInt_FMT " bs %" PetscInt_FMT " mem %g\n",
971: rank,mat->rmap->n,(PetscInt)info.nz_used,(PetscInt)info.nz_allocated,mat->rmap->bs,(double)info.memory));
972: MatGetInfo(baij->A,MAT_LOCAL,&info);
973: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] on-diagonal part: nz %" PetscInt_FMT " \n",rank,(PetscInt)info.nz_used);
974: MatGetInfo(baij->B,MAT_LOCAL,&info);
975: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] off-diagonal part: nz %" PetscInt_FMT " \n",rank,(PetscInt)info.nz_used);
976: PetscViewerFlush(viewer);
977: PetscViewerASCIIPopSynchronized(viewer);
978: PetscViewerASCIIPrintf(viewer,"Information on VecScatter used in matrix-vector product: \n");
979: VecScatterView(baij->Mvctx,viewer);
980: return 0;
981: } else if (format == PETSC_VIEWER_ASCII_INFO) {
982: PetscViewerASCIIPrintf(viewer," block size is %" PetscInt_FMT "\n",bs);
983: return 0;
984: } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
985: return 0;
986: }
987: }
989: if (isdraw) {
990: PetscDraw draw;
991: PetscBool isnull;
992: PetscViewerDrawGetDraw(viewer,0,&draw);
993: PetscDrawIsNull(draw,&isnull);
994: if (isnull) return 0;
995: }
997: {
998: /* assemble the entire matrix onto first processor. */
999: Mat A;
1000: Mat_SeqBAIJ *Aloc;
1001: PetscInt M = mat->rmap->N,N = mat->cmap->N,*ai,*aj,col,i,j,k,*rvals,mbs = baij->mbs;
1002: MatScalar *a;
1003: const char *matname;
1005: /* Here we are creating a temporary matrix, so will assume MPIBAIJ is acceptable */
1006: /* Perhaps this should be the type of mat? */
1007: MatCreate(PetscObjectComm((PetscObject)mat),&A);
1008: if (rank == 0) {
1009: MatSetSizes(A,M,N,M,N);
1010: } else {
1011: MatSetSizes(A,0,0,M,N);
1012: }
1013: MatSetType(A,MATMPIBAIJ);
1014: MatMPIBAIJSetPreallocation(A,mat->rmap->bs,0,NULL,0,NULL);
1015: MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);
1016: PetscLogObjectParent((PetscObject)mat,(PetscObject)A);
1018: /* copy over the A part */
1019: Aloc = (Mat_SeqBAIJ*)baij->A->data;
1020: ai = Aloc->i; aj = Aloc->j; a = Aloc->a;
1021: PetscMalloc1(bs,&rvals);
1023: for (i=0; i<mbs; i++) {
1024: rvals[0] = bs*(baij->rstartbs + i);
1025: for (j=1; j<bs; j++) rvals[j] = rvals[j-1] + 1;
1026: for (j=ai[i]; j<ai[i+1]; j++) {
1027: col = (baij->cstartbs+aj[j])*bs;
1028: for (k=0; k<bs; k++) {
1029: MatSetValues_MPIBAIJ(A,bs,rvals,1,&col,a,INSERT_VALUES);
1030: col++; a += bs;
1031: }
1032: }
1033: }
1034: /* copy over the B part */
1035: Aloc = (Mat_SeqBAIJ*)baij->B->data;
1036: ai = Aloc->i; aj = Aloc->j; a = Aloc->a;
1037: for (i=0; i<mbs; i++) {
1038: rvals[0] = bs*(baij->rstartbs + i);
1039: for (j=1; j<bs; j++) rvals[j] = rvals[j-1] + 1;
1040: for (j=ai[i]; j<ai[i+1]; j++) {
1041: col = baij->garray[aj[j]]*bs;
1042: for (k=0; k<bs; k++) {
1043: MatSetValues_MPIBAIJ(A,bs,rvals,1,&col,a,INSERT_VALUES);
1044: col++; a += bs;
1045: }
1046: }
1047: }
1048: PetscFree(rvals);
1049: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
1050: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
1051: /*
1052: Everyone has to call to draw the matrix since the graphics waits are
1053: synchronized across all processors that share the PetscDraw object
1054: */
1055: PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);
1056: PetscObjectGetName((PetscObject)mat,&matname);
1057: if (rank == 0) {
1058: PetscObjectSetName((PetscObject)((Mat_MPIBAIJ*)(A->data))->A,matname);
1059: MatView_SeqBAIJ(((Mat_MPIBAIJ*)(A->data))->A,sviewer);
1060: }
1061: PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);
1062: PetscViewerFlush(viewer);
1063: MatDestroy(&A);
1064: }
1065: return 0;
1066: }
1068: /* Used for both MPIBAIJ and MPISBAIJ matrices */
1069: PetscErrorCode MatView_MPIBAIJ_Binary(Mat mat,PetscViewer viewer)
1070: {
1071: Mat_MPIBAIJ *aij = (Mat_MPIBAIJ*)mat->data;
1072: Mat_SeqBAIJ *A = (Mat_SeqBAIJ*)aij->A->data;
1073: Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)aij->B->data;
1074: const PetscInt *garray = aij->garray;
1075: PetscInt header[4],M,N,m,rs,cs,bs,nz,cnt,i,j,ja,jb,k,l;
1076: PetscInt *rowlens,*colidxs;
1077: PetscScalar *matvals;
1079: PetscViewerSetUp(viewer);
1081: M = mat->rmap->N;
1082: N = mat->cmap->N;
1083: m = mat->rmap->n;
1084: rs = mat->rmap->rstart;
1085: cs = mat->cmap->rstart;
1086: bs = mat->rmap->bs;
1087: nz = bs*bs*(A->nz + B->nz);
1089: /* write matrix header */
1090: header[0] = MAT_FILE_CLASSID;
1091: header[1] = M; header[2] = N; header[3] = nz;
1092: MPI_Reduce(&nz,&header[3],1,MPIU_INT,MPI_SUM,0,PetscObjectComm((PetscObject)mat));
1093: PetscViewerBinaryWrite(viewer,header,4,PETSC_INT);
1095: /* fill in and store row lengths */
1096: PetscMalloc1(m,&rowlens);
1097: for (cnt=0, i=0; i<A->mbs; i++)
1098: for (j=0; j<bs; j++)
1099: rowlens[cnt++] = bs*(A->i[i+1] - A->i[i] + B->i[i+1] - B->i[i]);
1100: PetscViewerBinaryWriteAll(viewer,rowlens,m,rs,M,PETSC_INT);
1101: PetscFree(rowlens);
1103: /* fill in and store column indices */
1104: PetscMalloc1(nz,&colidxs);
1105: for (cnt=0, i=0; i<A->mbs; i++) {
1106: for (k=0; k<bs; k++) {
1107: for (jb=B->i[i]; jb<B->i[i+1]; jb++) {
1108: if (garray[B->j[jb]] > cs/bs) break;
1109: for (l=0; l<bs; l++)
1110: colidxs[cnt++] = bs*garray[B->j[jb]] + l;
1111: }
1112: for (ja=A->i[i]; ja<A->i[i+1]; ja++)
1113: for (l=0; l<bs; l++)
1114: colidxs[cnt++] = bs*A->j[ja] + l + cs;
1115: for (; jb<B->i[i+1]; jb++)
1116: for (l=0; l<bs; l++)
1117: colidxs[cnt++] = bs*garray[B->j[jb]] + l;
1118: }
1119: }
1121: PetscViewerBinaryWriteAll(viewer,colidxs,nz,PETSC_DECIDE,PETSC_DECIDE,PETSC_INT);
1122: PetscFree(colidxs);
1124: /* fill in and store nonzero values */
1125: PetscMalloc1(nz,&matvals);
1126: for (cnt=0, i=0; i<A->mbs; i++) {
1127: for (k=0; k<bs; k++) {
1128: for (jb=B->i[i]; jb<B->i[i+1]; jb++) {
1129: if (garray[B->j[jb]] > cs/bs) break;
1130: for (l=0; l<bs; l++)
1131: matvals[cnt++] = B->a[bs*(bs*jb + l) + k];
1132: }
1133: for (ja=A->i[i]; ja<A->i[i+1]; ja++)
1134: for (l=0; l<bs; l++)
1135: matvals[cnt++] = A->a[bs*(bs*ja + l) + k];
1136: for (; jb<B->i[i+1]; jb++)
1137: for (l=0; l<bs; l++)
1138: matvals[cnt++] = B->a[bs*(bs*jb + l) + k];
1139: }
1140: }
1141: PetscViewerBinaryWriteAll(viewer,matvals,nz,PETSC_DECIDE,PETSC_DECIDE,PETSC_SCALAR);
1142: PetscFree(matvals);
1144: /* write block size option to the viewer's .info file */
1145: MatView_Binary_BlockSizes(mat,viewer);
1146: return 0;
1147: }
1149: PetscErrorCode MatView_MPIBAIJ(Mat mat,PetscViewer viewer)
1150: {
1151: PetscBool iascii,isdraw,issocket,isbinary;
1153: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1154: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
1155: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSOCKET,&issocket);
1156: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
1157: if (iascii || isdraw || issocket) {
1158: MatView_MPIBAIJ_ASCIIorDraworSocket(mat,viewer);
1159: } else if (isbinary) {
1160: MatView_MPIBAIJ_Binary(mat,viewer);
1161: }
1162: return 0;
1163: }
1165: PetscErrorCode MatDestroy_MPIBAIJ(Mat mat)
1166: {
1167: Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)mat->data;
1169: #if defined(PETSC_USE_LOG)
1170: PetscLogObjectState((PetscObject)mat,"Rows=%" PetscInt_FMT ",Cols=%" PetscInt_FMT,mat->rmap->N,mat->cmap->N);
1171: #endif
1172: MatStashDestroy_Private(&mat->stash);
1173: MatStashDestroy_Private(&mat->bstash);
1174: MatDestroy(&baij->A);
1175: MatDestroy(&baij->B);
1176: #if defined(PETSC_USE_CTABLE)
1177: PetscTableDestroy(&baij->colmap);
1178: #else
1179: PetscFree(baij->colmap);
1180: #endif
1181: PetscFree(baij->garray);
1182: VecDestroy(&baij->lvec);
1183: VecScatterDestroy(&baij->Mvctx);
1184: PetscFree2(baij->rowvalues,baij->rowindices);
1185: PetscFree(baij->barray);
1186: PetscFree2(baij->hd,baij->ht);
1187: PetscFree(baij->rangebs);
1188: PetscFree(mat->data);
1190: PetscObjectChangeTypeName((PetscObject)mat,NULL);
1191: PetscObjectComposeFunction((PetscObject)mat,"MatStoreValues_C",NULL);
1192: PetscObjectComposeFunction((PetscObject)mat,"MatRetrieveValues_C",NULL);
1193: PetscObjectComposeFunction((PetscObject)mat,"MatMPIBAIJSetPreallocation_C",NULL);
1194: PetscObjectComposeFunction((PetscObject)mat,"MatMPIBAIJSetPreallocationCSR_C",NULL);
1195: PetscObjectComposeFunction((PetscObject)mat,"MatDiagonalScaleLocal_C",NULL);
1196: PetscObjectComposeFunction((PetscObject)mat,"MatSetHashTableFactor_C",NULL);
1197: PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpibaij_mpisbaij_C",NULL);
1198: PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpibaij_mpibstrm_C",NULL);
1199: #if defined(PETSC_HAVE_HYPRE)
1200: PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpibaij_hypre_C",NULL);
1201: #endif
1202: PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpibaij_is_C",NULL);
1203: return 0;
1204: }
1206: PetscErrorCode MatMult_MPIBAIJ(Mat A,Vec xx,Vec yy)
1207: {
1208: Mat_MPIBAIJ *a = (Mat_MPIBAIJ*)A->data;
1209: PetscInt nt;
1211: VecGetLocalSize(xx,&nt);
1213: VecGetLocalSize(yy,&nt);
1215: VecScatterBegin(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);
1216: (*a->A->ops->mult)(a->A,xx,yy);
1217: VecScatterEnd(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);
1218: (*a->B->ops->multadd)(a->B,a->lvec,yy,yy);
1219: return 0;
1220: }
1222: PetscErrorCode MatMultAdd_MPIBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
1223: {
1224: Mat_MPIBAIJ *a = (Mat_MPIBAIJ*)A->data;
1226: VecScatterBegin(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);
1227: (*a->A->ops->multadd)(a->A,xx,yy,zz);
1228: VecScatterEnd(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);
1229: (*a->B->ops->multadd)(a->B,a->lvec,zz,zz);
1230: return 0;
1231: }
1233: PetscErrorCode MatMultTranspose_MPIBAIJ(Mat A,Vec xx,Vec yy)
1234: {
1235: Mat_MPIBAIJ *a = (Mat_MPIBAIJ*)A->data;
1237: /* do nondiagonal part */
1238: (*a->B->ops->multtranspose)(a->B,xx,a->lvec);
1239: /* do local part */
1240: (*a->A->ops->multtranspose)(a->A,xx,yy);
1241: /* add partial results together */
1242: VecScatterBegin(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);
1243: VecScatterEnd(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);
1244: return 0;
1245: }
1247: PetscErrorCode MatMultTransposeAdd_MPIBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
1248: {
1249: Mat_MPIBAIJ *a = (Mat_MPIBAIJ*)A->data;
1251: /* do nondiagonal part */
1252: (*a->B->ops->multtranspose)(a->B,xx,a->lvec);
1253: /* do local part */
1254: (*a->A->ops->multtransposeadd)(a->A,xx,yy,zz);
1255: /* add partial results together */
1256: VecScatterBegin(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);
1257: VecScatterEnd(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);
1258: return 0;
1259: }
1261: /*
1262: This only works correctly for square matrices where the subblock A->A is the
1263: diagonal block
1264: */
1265: PetscErrorCode MatGetDiagonal_MPIBAIJ(Mat A,Vec v)
1266: {
1268: MatGetDiagonal(((Mat_MPIBAIJ*)A->data)->A,v);
1269: return 0;
1270: }
1272: PetscErrorCode MatScale_MPIBAIJ(Mat A,PetscScalar aa)
1273: {
1274: Mat_MPIBAIJ *a = (Mat_MPIBAIJ*)A->data;
1276: MatScale(a->A,aa);
1277: MatScale(a->B,aa);
1278: return 0;
1279: }
1281: PetscErrorCode MatGetRow_MPIBAIJ(Mat matin,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
1282: {
1283: Mat_MPIBAIJ *mat = (Mat_MPIBAIJ*)matin->data;
1284: PetscScalar *vworkA,*vworkB,**pvA,**pvB,*v_p;
1285: PetscInt bs = matin->rmap->bs,bs2 = mat->bs2,i,*cworkA,*cworkB,**pcA,**pcB;
1286: PetscInt nztot,nzA,nzB,lrow,brstart = matin->rmap->rstart,brend = matin->rmap->rend;
1287: PetscInt *cmap,*idx_p,cstart = mat->cstartbs;
1291: mat->getrowactive = PETSC_TRUE;
1293: if (!mat->rowvalues && (idx || v)) {
1294: /*
1295: allocate enough space to hold information from the longest row.
1296: */
1297: Mat_SeqBAIJ *Aa = (Mat_SeqBAIJ*)mat->A->data,*Ba = (Mat_SeqBAIJ*)mat->B->data;
1298: PetscInt max = 1,mbs = mat->mbs,tmp;
1299: for (i=0; i<mbs; i++) {
1300: tmp = Aa->i[i+1] - Aa->i[i] + Ba->i[i+1] - Ba->i[i];
1301: if (max < tmp) max = tmp;
1302: }
1303: PetscMalloc2(max*bs2,&mat->rowvalues,max*bs2,&mat->rowindices);
1304: }
1305: lrow = row - brstart;
1307: pvA = &vworkA; pcA = &cworkA; pvB = &vworkB; pcB = &cworkB;
1308: if (!v) {pvA = NULL; pvB = NULL;}
1309: if (!idx) {pcA = NULL; if (!v) pcB = NULL;}
1310: (*mat->A->ops->getrow)(mat->A,lrow,&nzA,pcA,pvA);
1311: (*mat->B->ops->getrow)(mat->B,lrow,&nzB,pcB,pvB);
1312: nztot = nzA + nzB;
1314: cmap = mat->garray;
1315: if (v || idx) {
1316: if (nztot) {
1317: /* Sort by increasing column numbers, assuming A and B already sorted */
1318: PetscInt imark = -1;
1319: if (v) {
1320: *v = v_p = mat->rowvalues;
1321: for (i=0; i<nzB; i++) {
1322: if (cmap[cworkB[i]/bs] < cstart) v_p[i] = vworkB[i];
1323: else break;
1324: }
1325: imark = i;
1326: for (i=0; i<nzA; i++) v_p[imark+i] = vworkA[i];
1327: for (i=imark; i<nzB; i++) v_p[nzA+i] = vworkB[i];
1328: }
1329: if (idx) {
1330: *idx = idx_p = mat->rowindices;
1331: if (imark > -1) {
1332: for (i=0; i<imark; i++) {
1333: idx_p[i] = cmap[cworkB[i]/bs]*bs + cworkB[i]%bs;
1334: }
1335: } else {
1336: for (i=0; i<nzB; i++) {
1337: if (cmap[cworkB[i]/bs] < cstart) idx_p[i] = cmap[cworkB[i]/bs]*bs + cworkB[i]%bs;
1338: else break;
1339: }
1340: imark = i;
1341: }
1342: for (i=0; i<nzA; i++) idx_p[imark+i] = cstart*bs + cworkA[i];
1343: for (i=imark; i<nzB; i++) idx_p[nzA+i] = cmap[cworkB[i]/bs]*bs + cworkB[i]%bs ;
1344: }
1345: } else {
1346: if (idx) *idx = NULL;
1347: if (v) *v = NULL;
1348: }
1349: }
1350: *nz = nztot;
1351: (*mat->A->ops->restorerow)(mat->A,lrow,&nzA,pcA,pvA);
1352: (*mat->B->ops->restorerow)(mat->B,lrow,&nzB,pcB,pvB);
1353: return 0;
1354: }
1356: PetscErrorCode MatRestoreRow_MPIBAIJ(Mat mat,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
1357: {
1358: Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)mat->data;
1361: baij->getrowactive = PETSC_FALSE;
1362: return 0;
1363: }
1365: PetscErrorCode MatZeroEntries_MPIBAIJ(Mat A)
1366: {
1367: Mat_MPIBAIJ *l = (Mat_MPIBAIJ*)A->data;
1369: MatZeroEntries(l->A);
1370: MatZeroEntries(l->B);
1371: return 0;
1372: }
1374: PetscErrorCode MatGetInfo_MPIBAIJ(Mat matin,MatInfoType flag,MatInfo *info)
1375: {
1376: Mat_MPIBAIJ *a = (Mat_MPIBAIJ*)matin->data;
1377: Mat A = a->A,B = a->B;
1378: PetscLogDouble isend[5],irecv[5];
1380: info->block_size = (PetscReal)matin->rmap->bs;
1382: MatGetInfo(A,MAT_LOCAL,info);
1384: isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
1385: isend[3] = info->memory; isend[4] = info->mallocs;
1387: MatGetInfo(B,MAT_LOCAL,info);
1389: isend[0] += info->nz_used; isend[1] += info->nz_allocated; isend[2] += info->nz_unneeded;
1390: isend[3] += info->memory; isend[4] += info->mallocs;
1392: if (flag == MAT_LOCAL) {
1393: info->nz_used = isend[0];
1394: info->nz_allocated = isend[1];
1395: info->nz_unneeded = isend[2];
1396: info->memory = isend[3];
1397: info->mallocs = isend[4];
1398: } else if (flag == MAT_GLOBAL_MAX) {
1399: MPIU_Allreduce(isend,irecv,5,MPIU_PETSCLOGDOUBLE,MPI_MAX,PetscObjectComm((PetscObject)matin));
1401: info->nz_used = irecv[0];
1402: info->nz_allocated = irecv[1];
1403: info->nz_unneeded = irecv[2];
1404: info->memory = irecv[3];
1405: info->mallocs = irecv[4];
1406: } else if (flag == MAT_GLOBAL_SUM) {
1407: MPIU_Allreduce(isend,irecv,5,MPIU_PETSCLOGDOUBLE,MPI_SUM,PetscObjectComm((PetscObject)matin));
1409: info->nz_used = irecv[0];
1410: info->nz_allocated = irecv[1];
1411: info->nz_unneeded = irecv[2];
1412: info->memory = irecv[3];
1413: info->mallocs = irecv[4];
1414: } else SETERRQ(PetscObjectComm((PetscObject)matin),PETSC_ERR_ARG_WRONG,"Unknown MatInfoType argument %d",(int)flag);
1415: info->fill_ratio_given = 0; /* no parallel LU/ILU/Cholesky */
1416: info->fill_ratio_needed = 0;
1417: info->factor_mallocs = 0;
1418: return 0;
1419: }
1421: PetscErrorCode MatSetOption_MPIBAIJ(Mat A,MatOption op,PetscBool flg)
1422: {
1423: Mat_MPIBAIJ *a = (Mat_MPIBAIJ*)A->data;
1425: switch (op) {
1426: case MAT_NEW_NONZERO_LOCATIONS:
1427: case MAT_NEW_NONZERO_ALLOCATION_ERR:
1428: case MAT_UNUSED_NONZERO_LOCATION_ERR:
1429: case MAT_KEEP_NONZERO_PATTERN:
1430: case MAT_NEW_NONZERO_LOCATION_ERR:
1431: MatCheckPreallocated(A,1);
1432: MatSetOption(a->A,op,flg);
1433: MatSetOption(a->B,op,flg);
1434: break;
1435: case MAT_ROW_ORIENTED:
1436: MatCheckPreallocated(A,1);
1437: a->roworiented = flg;
1439: MatSetOption(a->A,op,flg);
1440: MatSetOption(a->B,op,flg);
1441: break;
1442: case MAT_FORCE_DIAGONAL_ENTRIES:
1443: case MAT_SORTED_FULL:
1444: PetscInfo(A,"Option %s ignored\n",MatOptions[op]);
1445: break;
1446: case MAT_IGNORE_OFF_PROC_ENTRIES:
1447: a->donotstash = flg;
1448: break;
1449: case MAT_USE_HASH_TABLE:
1450: a->ht_flag = flg;
1451: a->ht_fact = 1.39;
1452: break;
1453: case MAT_SYMMETRIC:
1454: case MAT_STRUCTURALLY_SYMMETRIC:
1455: case MAT_HERMITIAN:
1456: case MAT_SUBMAT_SINGLEIS:
1457: case MAT_SYMMETRY_ETERNAL:
1458: MatCheckPreallocated(A,1);
1459: MatSetOption(a->A,op,flg);
1460: break;
1461: default:
1462: SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"unknown option %d",op);
1463: }
1464: return 0;
1465: }
1467: PetscErrorCode MatTranspose_MPIBAIJ(Mat A,MatReuse reuse,Mat *matout)
1468: {
1469: Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)A->data;
1470: Mat_SeqBAIJ *Aloc;
1471: Mat B;
1472: PetscInt M =A->rmap->N,N=A->cmap->N,*ai,*aj,i,*rvals,j,k,col;
1473: PetscInt bs=A->rmap->bs,mbs=baij->mbs;
1474: MatScalar *a;
1476: if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_INPLACE_MATRIX) {
1477: MatCreate(PetscObjectComm((PetscObject)A),&B);
1478: MatSetSizes(B,A->cmap->n,A->rmap->n,N,M);
1479: MatSetType(B,((PetscObject)A)->type_name);
1480: /* Do not know preallocation information, but must set block size */
1481: MatMPIBAIJSetPreallocation(B,A->rmap->bs,PETSC_DECIDE,NULL,PETSC_DECIDE,NULL);
1482: } else {
1483: B = *matout;
1484: }
1486: /* copy over the A part */
1487: Aloc = (Mat_SeqBAIJ*)baij->A->data;
1488: ai = Aloc->i; aj = Aloc->j; a = Aloc->a;
1489: PetscMalloc1(bs,&rvals);
1491: for (i=0; i<mbs; i++) {
1492: rvals[0] = bs*(baij->rstartbs + i);
1493: for (j=1; j<bs; j++) rvals[j] = rvals[j-1] + 1;
1494: for (j=ai[i]; j<ai[i+1]; j++) {
1495: col = (baij->cstartbs+aj[j])*bs;
1496: for (k=0; k<bs; k++) {
1497: MatSetValues_MPIBAIJ(B,1,&col,bs,rvals,a,INSERT_VALUES);
1499: col++; a += bs;
1500: }
1501: }
1502: }
1503: /* copy over the B part */
1504: Aloc = (Mat_SeqBAIJ*)baij->B->data;
1505: ai = Aloc->i; aj = Aloc->j; a = Aloc->a;
1506: for (i=0; i<mbs; i++) {
1507: rvals[0] = bs*(baij->rstartbs + i);
1508: for (j=1; j<bs; j++) rvals[j] = rvals[j-1] + 1;
1509: for (j=ai[i]; j<ai[i+1]; j++) {
1510: col = baij->garray[aj[j]]*bs;
1511: for (k=0; k<bs; k++) {
1512: MatSetValues_MPIBAIJ(B,1,&col,bs,rvals,a,INSERT_VALUES);
1513: col++;
1514: a += bs;
1515: }
1516: }
1517: }
1518: PetscFree(rvals);
1519: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
1520: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
1522: if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_REUSE_MATRIX) *matout = B;
1523: else {
1524: MatHeaderMerge(A,&B);
1525: }
1526: return 0;
1527: }
1529: PetscErrorCode MatDiagonalScale_MPIBAIJ(Mat mat,Vec ll,Vec rr)
1530: {
1531: Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)mat->data;
1532: Mat a = baij->A,b = baij->B;
1533: PetscInt s1,s2,s3;
1535: MatGetLocalSize(mat,&s2,&s3);
1536: if (rr) {
1537: VecGetLocalSize(rr,&s1);
1539: /* Overlap communication with computation. */
1540: VecScatterBegin(baij->Mvctx,rr,baij->lvec,INSERT_VALUES,SCATTER_FORWARD);
1541: }
1542: if (ll) {
1543: VecGetLocalSize(ll,&s1);
1545: (*b->ops->diagonalscale)(b,ll,NULL);
1546: }
1547: /* scale the diagonal block */
1548: (*a->ops->diagonalscale)(a,ll,rr);
1550: if (rr) {
1551: /* Do a scatter end and then right scale the off-diagonal block */
1552: VecScatterEnd(baij->Mvctx,rr,baij->lvec,INSERT_VALUES,SCATTER_FORWARD);
1553: (*b->ops->diagonalscale)(b,NULL,baij->lvec);
1554: }
1555: return 0;
1556: }
1558: PetscErrorCode MatZeroRows_MPIBAIJ(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
1559: {
1560: Mat_MPIBAIJ *l = (Mat_MPIBAIJ *) A->data;
1561: PetscInt *lrows;
1562: PetscInt r, len;
1563: PetscBool cong;
1565: /* get locally owned rows */
1566: MatZeroRowsMapLocal_Private(A,N,rows,&len,&lrows);
1567: /* fix right hand side if needed */
1568: if (x && b) {
1569: const PetscScalar *xx;
1570: PetscScalar *bb;
1572: VecGetArrayRead(x,&xx);
1573: VecGetArray(b,&bb);
1574: for (r = 0; r < len; ++r) bb[lrows[r]] = diag*xx[lrows[r]];
1575: VecRestoreArrayRead(x,&xx);
1576: VecRestoreArray(b,&bb);
1577: }
1579: /* actually zap the local rows */
1580: /*
1581: Zero the required rows. If the "diagonal block" of the matrix
1582: is square and the user wishes to set the diagonal we use separate
1583: code so that MatSetValues() is not called for each diagonal allocating
1584: new memory, thus calling lots of mallocs and slowing things down.
1586: */
1587: /* must zero l->B before l->A because the (diag) case below may put values into l->B*/
1588: MatZeroRows_SeqBAIJ(l->B,len,lrows,0.0,NULL,NULL);
1589: MatHasCongruentLayouts(A,&cong);
1590: if ((diag != 0.0) && cong) {
1591: MatZeroRows_SeqBAIJ(l->A,len,lrows,diag,NULL,NULL);
1592: } else if (diag != 0.0) {
1593: MatZeroRows_SeqBAIJ(l->A,len,lrows,0.0,NULL,NULL);
1595: MAT_NEW_NONZERO_LOCATIONS,MAT_NEW_NONZERO_LOCATION_ERR,MAT_NEW_NONZERO_ALLOCATION_ERR");
1596: for (r = 0; r < len; ++r) {
1597: const PetscInt row = lrows[r] + A->rmap->rstart;
1598: MatSetValues(A,1,&row,1,&row,&diag,INSERT_VALUES);
1599: }
1600: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
1601: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
1602: } else {
1603: MatZeroRows_SeqBAIJ(l->A,len,lrows,0.0,NULL,NULL);
1604: }
1605: PetscFree(lrows);
1607: /* only change matrix nonzero state if pattern was allowed to be changed */
1608: if (!((Mat_SeqBAIJ*)(l->A->data))->keepnonzeropattern) {
1609: PetscObjectState state = l->A->nonzerostate + l->B->nonzerostate;
1610: MPIU_Allreduce(&state,&A->nonzerostate,1,MPIU_INT64,MPI_SUM,PetscObjectComm((PetscObject)A));
1611: }
1612: return 0;
1613: }
1615: PetscErrorCode MatZeroRowsColumns_MPIBAIJ(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
1616: {
1617: Mat_MPIBAIJ *l = (Mat_MPIBAIJ*)A->data;
1618: PetscMPIInt n = A->rmap->n,p = 0;
1619: PetscInt i,j,k,r,len = 0,row,col,count;
1620: PetscInt *lrows,*owners = A->rmap->range;
1621: PetscSFNode *rrows;
1622: PetscSF sf;
1623: const PetscScalar *xx;
1624: PetscScalar *bb,*mask;
1625: Vec xmask,lmask;
1626: Mat_SeqBAIJ *baij = (Mat_SeqBAIJ*)l->B->data;
1627: PetscInt bs = A->rmap->bs, bs2 = baij->bs2;
1628: PetscScalar *aa;
1630: /* Create SF where leaves are input rows and roots are owned rows */
1631: PetscMalloc1(n, &lrows);
1632: for (r = 0; r < n; ++r) lrows[r] = -1;
1633: PetscMalloc1(N, &rrows);
1634: for (r = 0; r < N; ++r) {
1635: const PetscInt idx = rows[r];
1637: if (idx < owners[p] || owners[p+1] <= idx) { /* short-circuit the search if the last p owns this row too */
1638: PetscLayoutFindOwner(A->rmap,idx,&p);
1639: }
1640: rrows[r].rank = p;
1641: rrows[r].index = rows[r] - owners[p];
1642: }
1643: PetscSFCreate(PetscObjectComm((PetscObject) A), &sf);
1644: PetscSFSetGraph(sf, n, N, NULL, PETSC_OWN_POINTER, rrows, PETSC_OWN_POINTER);
1645: /* Collect flags for rows to be zeroed */
1646: PetscSFReduceBegin(sf, MPIU_INT, (PetscInt *) rows, lrows, MPI_LOR);
1647: PetscSFReduceEnd(sf, MPIU_INT, (PetscInt *) rows, lrows, MPI_LOR);
1648: PetscSFDestroy(&sf);
1649: /* Compress and put in row numbers */
1650: for (r = 0; r < n; ++r) if (lrows[r] >= 0) lrows[len++] = r;
1651: /* zero diagonal part of matrix */
1652: MatZeroRowsColumns(l->A,len,lrows,diag,x,b);
1653: /* handle off diagonal part of matrix */
1654: MatCreateVecs(A,&xmask,NULL);
1655: VecDuplicate(l->lvec,&lmask);
1656: VecGetArray(xmask,&bb);
1657: for (i=0; i<len; i++) bb[lrows[i]] = 1;
1658: VecRestoreArray(xmask,&bb);
1659: VecScatterBegin(l->Mvctx,xmask,lmask,ADD_VALUES,SCATTER_FORWARD);
1660: VecScatterEnd(l->Mvctx,xmask,lmask,ADD_VALUES,SCATTER_FORWARD);
1661: VecDestroy(&xmask);
1662: if (x) {
1663: VecScatterBegin(l->Mvctx,x,l->lvec,INSERT_VALUES,SCATTER_FORWARD);
1664: VecScatterEnd(l->Mvctx,x,l->lvec,INSERT_VALUES,SCATTER_FORWARD);
1665: VecGetArrayRead(l->lvec,&xx);
1666: VecGetArray(b,&bb);
1667: }
1668: VecGetArray(lmask,&mask);
1669: /* remove zeroed rows of off diagonal matrix */
1670: for (i = 0; i < len; ++i) {
1671: row = lrows[i];
1672: count = (baij->i[row/bs +1] - baij->i[row/bs])*bs;
1673: aa = ((MatScalar*)(baij->a)) + baij->i[row/bs]*bs2 + (row%bs);
1674: for (k = 0; k < count; ++k) {
1675: aa[0] = 0.0;
1676: aa += bs;
1677: }
1678: }
1679: /* loop over all elements of off process part of matrix zeroing removed columns*/
1680: for (i = 0; i < l->B->rmap->N; ++i) {
1681: row = i/bs;
1682: for (j = baij->i[row]; j < baij->i[row+1]; ++j) {
1683: for (k = 0; k < bs; ++k) {
1684: col = bs*baij->j[j] + k;
1685: if (PetscAbsScalar(mask[col])) {
1686: aa = ((MatScalar*)(baij->a)) + j*bs2 + (i%bs) + bs*k;
1687: if (x) bb[i] -= aa[0]*xx[col];
1688: aa[0] = 0.0;
1689: }
1690: }
1691: }
1692: }
1693: if (x) {
1694: VecRestoreArray(b,&bb);
1695: VecRestoreArrayRead(l->lvec,&xx);
1696: }
1697: VecRestoreArray(lmask,&mask);
1698: VecDestroy(&lmask);
1699: PetscFree(lrows);
1701: /* only change matrix nonzero state if pattern was allowed to be changed */
1702: if (!((Mat_SeqBAIJ*)(l->A->data))->keepnonzeropattern) {
1703: PetscObjectState state = l->A->nonzerostate + l->B->nonzerostate;
1704: MPIU_Allreduce(&state,&A->nonzerostate,1,MPIU_INT64,MPI_SUM,PetscObjectComm((PetscObject)A));
1705: }
1706: return 0;
1707: }
1709: PetscErrorCode MatSetUnfactored_MPIBAIJ(Mat A)
1710: {
1711: Mat_MPIBAIJ *a = (Mat_MPIBAIJ*)A->data;
1713: MatSetUnfactored(a->A);
1714: return 0;
1715: }
1717: static PetscErrorCode MatDuplicate_MPIBAIJ(Mat,MatDuplicateOption,Mat*);
1719: PetscErrorCode MatEqual_MPIBAIJ(Mat A,Mat B,PetscBool *flag)
1720: {
1721: Mat_MPIBAIJ *matB = (Mat_MPIBAIJ*)B->data,*matA = (Mat_MPIBAIJ*)A->data;
1722: Mat a,b,c,d;
1723: PetscBool flg;
1725: a = matA->A; b = matA->B;
1726: c = matB->A; d = matB->B;
1728: MatEqual(a,c,&flg);
1729: if (flg) {
1730: MatEqual(b,d,&flg);
1731: }
1732: MPIU_Allreduce(&flg,flag,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));
1733: return 0;
1734: }
1736: PetscErrorCode MatCopy_MPIBAIJ(Mat A,Mat B,MatStructure str)
1737: {
1738: Mat_MPIBAIJ *a = (Mat_MPIBAIJ*)A->data;
1739: Mat_MPIBAIJ *b = (Mat_MPIBAIJ*)B->data;
1741: /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
1742: if ((str != SAME_NONZERO_PATTERN) || (A->ops->copy != B->ops->copy)) {
1743: MatCopy_Basic(A,B,str);
1744: } else {
1745: MatCopy(a->A,b->A,str);
1746: MatCopy(a->B,b->B,str);
1747: }
1748: PetscObjectStateIncrease((PetscObject)B);
1749: return 0;
1750: }
1752: PetscErrorCode MatSetUp_MPIBAIJ(Mat A)
1753: {
1754: MatMPIBAIJSetPreallocation(A,A->rmap->bs,PETSC_DEFAULT,NULL,PETSC_DEFAULT,NULL);
1755: return 0;
1756: }
1758: PetscErrorCode MatAXPYGetPreallocation_MPIBAIJ(Mat Y,const PetscInt *yltog,Mat X,const PetscInt *xltog,PetscInt *nnz)
1759: {
1760: PetscInt bs = Y->rmap->bs,m = Y->rmap->N/bs;
1761: Mat_SeqBAIJ *x = (Mat_SeqBAIJ*)X->data;
1762: Mat_SeqBAIJ *y = (Mat_SeqBAIJ*)Y->data;
1764: MatAXPYGetPreallocation_MPIX_private(m,x->i,x->j,xltog,y->i,y->j,yltog,nnz);
1765: return 0;
1766: }
1768: PetscErrorCode MatAXPY_MPIBAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str)
1769: {
1770: Mat_MPIBAIJ *xx=(Mat_MPIBAIJ*)X->data,*yy=(Mat_MPIBAIJ*)Y->data;
1771: PetscBLASInt bnz,one=1;
1772: Mat_SeqBAIJ *x,*y;
1773: PetscInt bs2 = Y->rmap->bs*Y->rmap->bs;
1775: if (str == SAME_NONZERO_PATTERN) {
1776: PetscScalar alpha = a;
1777: x = (Mat_SeqBAIJ*)xx->A->data;
1778: y = (Mat_SeqBAIJ*)yy->A->data;
1779: PetscBLASIntCast(x->nz*bs2,&bnz);
1780: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&bnz,&alpha,x->a,&one,y->a,&one));
1781: x = (Mat_SeqBAIJ*)xx->B->data;
1782: y = (Mat_SeqBAIJ*)yy->B->data;
1783: PetscBLASIntCast(x->nz*bs2,&bnz);
1784: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&bnz,&alpha,x->a,&one,y->a,&one));
1785: PetscObjectStateIncrease((PetscObject)Y);
1786: } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
1787: MatAXPY_Basic(Y,a,X,str);
1788: } else {
1789: Mat B;
1790: PetscInt *nnz_d,*nnz_o,bs=Y->rmap->bs;
1791: PetscMalloc1(yy->A->rmap->N,&nnz_d);
1792: PetscMalloc1(yy->B->rmap->N,&nnz_o);
1793: MatCreate(PetscObjectComm((PetscObject)Y),&B);
1794: PetscObjectSetName((PetscObject)B,((PetscObject)Y)->name);
1795: MatSetSizes(B,Y->rmap->n,Y->cmap->n,Y->rmap->N,Y->cmap->N);
1796: MatSetBlockSizesFromMats(B,Y,Y);
1797: MatSetType(B,MATMPIBAIJ);
1798: MatAXPYGetPreallocation_SeqBAIJ(yy->A,xx->A,nnz_d);
1799: MatAXPYGetPreallocation_MPIBAIJ(yy->B,yy->garray,xx->B,xx->garray,nnz_o);
1800: MatMPIBAIJSetPreallocation(B,bs,0,nnz_d,0,nnz_o);
1801: /* MatAXPY_BasicWithPreallocation() for BAIJ matrix is much slower than AIJ, even for bs=1 ! */
1802: MatAXPY_BasicWithPreallocation(B,Y,a,X,str);
1803: MatHeaderMerge(Y,&B);
1804: PetscFree(nnz_d);
1805: PetscFree(nnz_o);
1806: }
1807: return 0;
1808: }
1810: PetscErrorCode MatConjugate_MPIBAIJ(Mat mat)
1811: {
1812: if (PetscDefined(USE_COMPLEX)) {
1813: Mat_MPIBAIJ *a = (Mat_MPIBAIJ*)mat->data;
1815: MatConjugate_SeqBAIJ(a->A);
1816: MatConjugate_SeqBAIJ(a->B);
1817: }
1818: return 0;
1819: }
1821: PetscErrorCode MatRealPart_MPIBAIJ(Mat A)
1822: {
1823: Mat_MPIBAIJ *a = (Mat_MPIBAIJ*)A->data;
1825: MatRealPart(a->A);
1826: MatRealPart(a->B);
1827: return 0;
1828: }
1830: PetscErrorCode MatImaginaryPart_MPIBAIJ(Mat A)
1831: {
1832: Mat_MPIBAIJ *a = (Mat_MPIBAIJ*)A->data;
1834: MatImaginaryPart(a->A);
1835: MatImaginaryPart(a->B);
1836: return 0;
1837: }
1839: PetscErrorCode MatCreateSubMatrix_MPIBAIJ(Mat mat,IS isrow,IS iscol,MatReuse call,Mat *newmat)
1840: {
1841: IS iscol_local;
1842: PetscInt csize;
1844: ISGetLocalSize(iscol,&csize);
1845: if (call == MAT_REUSE_MATRIX) {
1846: PetscObjectQuery((PetscObject)*newmat,"ISAllGather",(PetscObject*)&iscol_local);
1848: } else {
1849: ISAllGather(iscol,&iscol_local);
1850: }
1851: MatCreateSubMatrix_MPIBAIJ_Private(mat,isrow,iscol_local,csize,call,newmat);
1852: if (call == MAT_INITIAL_MATRIX) {
1853: PetscObjectCompose((PetscObject)*newmat,"ISAllGather",(PetscObject)iscol_local);
1854: ISDestroy(&iscol_local);
1855: }
1856: return 0;
1857: }
1859: /*
1860: Not great since it makes two copies of the submatrix, first an SeqBAIJ
1861: in local and then by concatenating the local matrices the end result.
1862: Writing it directly would be much like MatCreateSubMatrices_MPIBAIJ().
1863: This routine is used for BAIJ and SBAIJ matrices (unfortunate dependency).
1864: */
1865: PetscErrorCode MatCreateSubMatrix_MPIBAIJ_Private(Mat mat,IS isrow,IS iscol,PetscInt csize,MatReuse call,Mat *newmat)
1866: {
1867: PetscMPIInt rank,size;
1868: PetscInt i,m,n,rstart,row,rend,nz,*cwork,j,bs;
1869: PetscInt *ii,*jj,nlocal,*dlens,*olens,dlen,olen,jend,mglobal;
1870: Mat M,Mreuse;
1871: MatScalar *vwork,*aa;
1872: MPI_Comm comm;
1873: IS isrow_new, iscol_new;
1874: Mat_SeqBAIJ *aij;
1876: PetscObjectGetComm((PetscObject)mat,&comm);
1877: MPI_Comm_rank(comm,&rank);
1878: MPI_Comm_size(comm,&size);
1879: /* The compression and expansion should be avoided. Doesn't point
1880: out errors, might change the indices, hence buggey */
1881: ISCompressIndicesGeneral(mat->rmap->N,mat->rmap->n,mat->rmap->bs,1,&isrow,&isrow_new);
1882: ISCompressIndicesGeneral(mat->cmap->N,mat->cmap->n,mat->cmap->bs,1,&iscol,&iscol_new);
1884: if (call == MAT_REUSE_MATRIX) {
1885: PetscObjectQuery((PetscObject)*newmat,"SubMatrix",(PetscObject*)&Mreuse);
1887: MatCreateSubMatrices_MPIBAIJ_local(mat,1,&isrow_new,&iscol_new,MAT_REUSE_MATRIX,&Mreuse);
1888: } else {
1889: MatCreateSubMatrices_MPIBAIJ_local(mat,1,&isrow_new,&iscol_new,MAT_INITIAL_MATRIX,&Mreuse);
1890: }
1891: ISDestroy(&isrow_new);
1892: ISDestroy(&iscol_new);
1893: /*
1894: m - number of local rows
1895: n - number of columns (same on all processors)
1896: rstart - first row in new global matrix generated
1897: */
1898: MatGetBlockSize(mat,&bs);
1899: MatGetSize(Mreuse,&m,&n);
1900: m = m/bs;
1901: n = n/bs;
1903: if (call == MAT_INITIAL_MATRIX) {
1904: aij = (Mat_SeqBAIJ*)(Mreuse)->data;
1905: ii = aij->i;
1906: jj = aij->j;
1908: /*
1909: Determine the number of non-zeros in the diagonal and off-diagonal
1910: portions of the matrix in order to do correct preallocation
1911: */
1913: /* first get start and end of "diagonal" columns */
1914: if (csize == PETSC_DECIDE) {
1915: ISGetSize(isrow,&mglobal);
1916: if (mglobal == n*bs) { /* square matrix */
1917: nlocal = m;
1918: } else {
1919: nlocal = n/size + ((n % size) > rank);
1920: }
1921: } else {
1922: nlocal = csize/bs;
1923: }
1924: MPI_Scan(&nlocal,&rend,1,MPIU_INT,MPI_SUM,comm);
1925: rstart = rend - nlocal;
1928: /* next, compute all the lengths */
1929: PetscMalloc2(m+1,&dlens,m+1,&olens);
1930: for (i=0; i<m; i++) {
1931: jend = ii[i+1] - ii[i];
1932: olen = 0;
1933: dlen = 0;
1934: for (j=0; j<jend; j++) {
1935: if (*jj < rstart || *jj >= rend) olen++;
1936: else dlen++;
1937: jj++;
1938: }
1939: olens[i] = olen;
1940: dlens[i] = dlen;
1941: }
1942: MatCreate(comm,&M);
1943: MatSetSizes(M,bs*m,bs*nlocal,PETSC_DECIDE,bs*n);
1944: MatSetType(M,((PetscObject)mat)->type_name);
1945: MatMPIBAIJSetPreallocation(M,bs,0,dlens,0,olens);
1946: MatMPISBAIJSetPreallocation(M,bs,0,dlens,0,olens);
1947: PetscFree2(dlens,olens);
1948: } else {
1949: PetscInt ml,nl;
1951: M = *newmat;
1952: MatGetLocalSize(M,&ml,&nl);
1954: MatZeroEntries(M);
1955: /*
1956: The next two lines are needed so we may call MatSetValues_MPIAIJ() below directly,
1957: rather than the slower MatSetValues().
1958: */
1959: M->was_assembled = PETSC_TRUE;
1960: M->assembled = PETSC_FALSE;
1961: }
1962: MatSetOption(M,MAT_ROW_ORIENTED,PETSC_FALSE);
1963: MatGetOwnershipRange(M,&rstart,&rend);
1964: aij = (Mat_SeqBAIJ*)(Mreuse)->data;
1965: ii = aij->i;
1966: jj = aij->j;
1967: aa = aij->a;
1968: for (i=0; i<m; i++) {
1969: row = rstart/bs + i;
1970: nz = ii[i+1] - ii[i];
1971: cwork = jj; jj += nz;
1972: vwork = aa; aa += nz*bs*bs;
1973: MatSetValuesBlocked_MPIBAIJ(M,1,&row,nz,cwork,vwork,INSERT_VALUES);
1974: }
1976: MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
1977: MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
1978: *newmat = M;
1980: /* save submatrix used in processor for next request */
1981: if (call == MAT_INITIAL_MATRIX) {
1982: PetscObjectCompose((PetscObject)M,"SubMatrix",(PetscObject)Mreuse);
1983: PetscObjectDereference((PetscObject)Mreuse);
1984: }
1985: return 0;
1986: }
1988: PetscErrorCode MatPermute_MPIBAIJ(Mat A,IS rowp,IS colp,Mat *B)
1989: {
1990: MPI_Comm comm,pcomm;
1991: PetscInt clocal_size,nrows;
1992: const PetscInt *rows;
1993: PetscMPIInt size;
1994: IS crowp,lcolp;
1996: PetscObjectGetComm((PetscObject)A,&comm);
1997: /* make a collective version of 'rowp' */
1998: PetscObjectGetComm((PetscObject)rowp,&pcomm);
1999: if (pcomm==comm) {
2000: crowp = rowp;
2001: } else {
2002: ISGetSize(rowp,&nrows);
2003: ISGetIndices(rowp,&rows);
2004: ISCreateGeneral(comm,nrows,rows,PETSC_COPY_VALUES,&crowp);
2005: ISRestoreIndices(rowp,&rows);
2006: }
2007: ISSetPermutation(crowp);
2008: /* make a local version of 'colp' */
2009: PetscObjectGetComm((PetscObject)colp,&pcomm);
2010: MPI_Comm_size(pcomm,&size);
2011: if (size==1) {
2012: lcolp = colp;
2013: } else {
2014: ISAllGather(colp,&lcolp);
2015: }
2016: ISSetPermutation(lcolp);
2017: /* now we just get the submatrix */
2018: MatGetLocalSize(A,NULL,&clocal_size);
2019: MatCreateSubMatrix_MPIBAIJ_Private(A,crowp,lcolp,clocal_size,MAT_INITIAL_MATRIX,B);
2020: /* clean up */
2021: if (pcomm!=comm) {
2022: ISDestroy(&crowp);
2023: }
2024: if (size>1) {
2025: ISDestroy(&lcolp);
2026: }
2027: return 0;
2028: }
2030: PetscErrorCode MatGetGhosts_MPIBAIJ(Mat mat,PetscInt *nghosts,const PetscInt *ghosts[])
2031: {
2032: Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*) mat->data;
2033: Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)baij->B->data;
2035: if (nghosts) *nghosts = B->nbs;
2036: if (ghosts) *ghosts = baij->garray;
2037: return 0;
2038: }
2040: PetscErrorCode MatGetSeqNonzeroStructure_MPIBAIJ(Mat A,Mat *newmat)
2041: {
2042: Mat B;
2043: Mat_MPIBAIJ *a = (Mat_MPIBAIJ*)A->data;
2044: Mat_SeqBAIJ *ad = (Mat_SeqBAIJ*)a->A->data,*bd = (Mat_SeqBAIJ*)a->B->data;
2045: Mat_SeqAIJ *b;
2046: PetscMPIInt size,rank,*recvcounts = NULL,*displs = NULL;
2047: PetscInt sendcount,i,*rstarts = A->rmap->range,n,cnt,j,bs = A->rmap->bs;
2048: PetscInt m,*garray = a->garray,*lens,*jsendbuf,*a_jsendbuf,*b_jsendbuf;
2050: MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);
2051: MPI_Comm_rank(PetscObjectComm((PetscObject)A),&rank);
2053: /* ----------------------------------------------------------------
2054: Tell every processor the number of nonzeros per row
2055: */
2056: PetscMalloc1(A->rmap->N/bs,&lens);
2057: for (i=A->rmap->rstart/bs; i<A->rmap->rend/bs; i++) {
2058: lens[i] = ad->i[i-A->rmap->rstart/bs+1] - ad->i[i-A->rmap->rstart/bs] + bd->i[i-A->rmap->rstart/bs+1] - bd->i[i-A->rmap->rstart/bs];
2059: }
2060: PetscMalloc1(2*size,&recvcounts);
2061: displs = recvcounts + size;
2062: for (i=0; i<size; i++) {
2063: recvcounts[i] = A->rmap->range[i+1]/bs - A->rmap->range[i]/bs;
2064: displs[i] = A->rmap->range[i]/bs;
2065: }
2066: MPI_Allgatherv(MPI_IN_PLACE,0,MPI_DATATYPE_NULL,lens,recvcounts,displs,MPIU_INT,PetscObjectComm((PetscObject)A));
2067: /* ---------------------------------------------------------------
2068: Create the sequential matrix of the same type as the local block diagonal
2069: */
2070: MatCreate(PETSC_COMM_SELF,&B);
2071: MatSetSizes(B,A->rmap->N/bs,A->cmap->N/bs,PETSC_DETERMINE,PETSC_DETERMINE);
2072: MatSetType(B,MATSEQAIJ);
2073: MatSeqAIJSetPreallocation(B,0,lens);
2074: b = (Mat_SeqAIJ*)B->data;
2076: /*--------------------------------------------------------------------
2077: Copy my part of matrix column indices over
2078: */
2079: sendcount = ad->nz + bd->nz;
2080: jsendbuf = b->j + b->i[rstarts[rank]/bs];
2081: a_jsendbuf = ad->j;
2082: b_jsendbuf = bd->j;
2083: n = A->rmap->rend/bs - A->rmap->rstart/bs;
2084: cnt = 0;
2085: for (i=0; i<n; i++) {
2087: /* put in lower diagonal portion */
2088: m = bd->i[i+1] - bd->i[i];
2089: while (m > 0) {
2090: /* is it above diagonal (in bd (compressed) numbering) */
2091: if (garray[*b_jsendbuf] > A->rmap->rstart/bs + i) break;
2092: jsendbuf[cnt++] = garray[*b_jsendbuf++];
2093: m--;
2094: }
2096: /* put in diagonal portion */
2097: for (j=ad->i[i]; j<ad->i[i+1]; j++) {
2098: jsendbuf[cnt++] = A->rmap->rstart/bs + *a_jsendbuf++;
2099: }
2101: /* put in upper diagonal portion */
2102: while (m-- > 0) {
2103: jsendbuf[cnt++] = garray[*b_jsendbuf++];
2104: }
2105: }
2108: /*--------------------------------------------------------------------
2109: Gather all column indices to all processors
2110: */
2111: for (i=0; i<size; i++) {
2112: recvcounts[i] = 0;
2113: for (j=A->rmap->range[i]/bs; j<A->rmap->range[i+1]/bs; j++) {
2114: recvcounts[i] += lens[j];
2115: }
2116: }
2117: displs[0] = 0;
2118: for (i=1; i<size; i++) {
2119: displs[i] = displs[i-1] + recvcounts[i-1];
2120: }
2121: MPI_Allgatherv(MPI_IN_PLACE,0,MPI_DATATYPE_NULL,b->j,recvcounts,displs,MPIU_INT,PetscObjectComm((PetscObject)A));
2122: /*--------------------------------------------------------------------
2123: Assemble the matrix into useable form (note numerical values not yet set)
2124: */
2125: /* set the b->ilen (length of each row) values */
2126: PetscArraycpy(b->ilen,lens,A->rmap->N/bs);
2127: /* set the b->i indices */
2128: b->i[0] = 0;
2129: for (i=1; i<=A->rmap->N/bs; i++) {
2130: b->i[i] = b->i[i-1] + lens[i-1];
2131: }
2132: PetscFree(lens);
2133: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
2134: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
2135: PetscFree(recvcounts);
2137: if (A->symmetric) {
2138: MatSetOption(B,MAT_SYMMETRIC,PETSC_TRUE);
2139: } else if (A->hermitian) {
2140: MatSetOption(B,MAT_HERMITIAN,PETSC_TRUE);
2141: } else if (A->structurally_symmetric) {
2142: MatSetOption(B,MAT_STRUCTURALLY_SYMMETRIC,PETSC_TRUE);
2143: }
2144: *newmat = B;
2145: return 0;
2146: }
2148: PetscErrorCode MatSOR_MPIBAIJ(Mat matin,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
2149: {
2150: Mat_MPIBAIJ *mat = (Mat_MPIBAIJ*)matin->data;
2151: Vec bb1 = NULL;
2153: if (flag == SOR_APPLY_UPPER) {
2154: (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx);
2155: return 0;
2156: }
2158: if (its > 1 || ~flag & SOR_ZERO_INITIAL_GUESS) {
2159: VecDuplicate(bb,&bb1);
2160: }
2162: if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP) {
2163: if (flag & SOR_ZERO_INITIAL_GUESS) {
2164: (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx);
2165: its--;
2166: }
2168: while (its--) {
2169: VecScatterBegin(mat->Mvctx,xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD);
2170: VecScatterEnd(mat->Mvctx,xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD);
2172: /* update rhs: bb1 = bb - B*x */
2173: VecScale(mat->lvec,-1.0);
2174: (*mat->B->ops->multadd)(mat->B,mat->lvec,bb,bb1);
2176: /* local sweep */
2177: (*mat->A->ops->sor)(mat->A,bb1,omega,SOR_SYMMETRIC_SWEEP,fshift,lits,1,xx);
2178: }
2179: } else if (flag & SOR_LOCAL_FORWARD_SWEEP) {
2180: if (flag & SOR_ZERO_INITIAL_GUESS) {
2181: (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx);
2182: its--;
2183: }
2184: while (its--) {
2185: VecScatterBegin(mat->Mvctx,xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD);
2186: VecScatterEnd(mat->Mvctx,xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD);
2188: /* update rhs: bb1 = bb - B*x */
2189: VecScale(mat->lvec,-1.0);
2190: (*mat->B->ops->multadd)(mat->B,mat->lvec,bb,bb1);
2192: /* local sweep */
2193: (*mat->A->ops->sor)(mat->A,bb1,omega,SOR_FORWARD_SWEEP,fshift,lits,1,xx);
2194: }
2195: } else if (flag & SOR_LOCAL_BACKWARD_SWEEP) {
2196: if (flag & SOR_ZERO_INITIAL_GUESS) {
2197: (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx);
2198: its--;
2199: }
2200: while (its--) {
2201: VecScatterBegin(mat->Mvctx,xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD);
2202: VecScatterEnd(mat->Mvctx,xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD);
2204: /* update rhs: bb1 = bb - B*x */
2205: VecScale(mat->lvec,-1.0);
2206: (*mat->B->ops->multadd)(mat->B,mat->lvec,bb,bb1);
2208: /* local sweep */
2209: (*mat->A->ops->sor)(mat->A,bb1,omega,SOR_BACKWARD_SWEEP,fshift,lits,1,xx);
2210: }
2211: } else SETERRQ(PetscObjectComm((PetscObject)matin),PETSC_ERR_SUP,"Parallel version of SOR requested not supported");
2213: VecDestroy(&bb1);
2214: return 0;
2215: }
2217: PetscErrorCode MatGetColumnReductions_MPIBAIJ(Mat A,PetscInt type,PetscReal *reductions)
2218: {
2219: Mat_MPIBAIJ *aij = (Mat_MPIBAIJ*)A->data;
2220: PetscInt m,N,i,*garray = aij->garray;
2221: PetscInt ib,jb,bs = A->rmap->bs;
2222: Mat_SeqBAIJ *a_aij = (Mat_SeqBAIJ*) aij->A->data;
2223: MatScalar *a_val = a_aij->a;
2224: Mat_SeqBAIJ *b_aij = (Mat_SeqBAIJ*) aij->B->data;
2225: MatScalar *b_val = b_aij->a;
2226: PetscReal *work;
2228: MatGetSize(A,&m,&N);
2229: PetscCalloc1(N,&work);
2230: if (type == NORM_2) {
2231: for (i=a_aij->i[0]; i<a_aij->i[aij->A->rmap->n/bs]; i++) {
2232: for (jb=0; jb<bs; jb++) {
2233: for (ib=0; ib<bs; ib++) {
2234: work[A->cmap->rstart + a_aij->j[i] * bs + jb] += PetscAbsScalar(*a_val * *a_val);
2235: a_val++;
2236: }
2237: }
2238: }
2239: for (i=b_aij->i[0]; i<b_aij->i[aij->B->rmap->n/bs]; i++) {
2240: for (jb=0; jb<bs; jb++) {
2241: for (ib=0; ib<bs; ib++) {
2242: work[garray[b_aij->j[i]] * bs + jb] += PetscAbsScalar(*b_val * *b_val);
2243: b_val++;
2244: }
2245: }
2246: }
2247: } else if (type == NORM_1) {
2248: for (i=a_aij->i[0]; i<a_aij->i[aij->A->rmap->n/bs]; i++) {
2249: for (jb=0; jb<bs; jb++) {
2250: for (ib=0; ib<bs; ib++) {
2251: work[A->cmap->rstart + a_aij->j[i] * bs + jb] += PetscAbsScalar(*a_val);
2252: a_val++;
2253: }
2254: }
2255: }
2256: for (i=b_aij->i[0]; i<b_aij->i[aij->B->rmap->n/bs]; i++) {
2257: for (jb=0; jb<bs; jb++) {
2258: for (ib=0; ib<bs; ib++) {
2259: work[garray[b_aij->j[i]] * bs + jb] += PetscAbsScalar(*b_val);
2260: b_val++;
2261: }
2262: }
2263: }
2264: } else if (type == NORM_INFINITY) {
2265: for (i=a_aij->i[0]; i<a_aij->i[aij->A->rmap->n/bs]; i++) {
2266: for (jb=0; jb<bs; jb++) {
2267: for (ib=0; ib<bs; ib++) {
2268: int col = A->cmap->rstart + a_aij->j[i] * bs + jb;
2269: work[col] = PetscMax(PetscAbsScalar(*a_val), work[col]);
2270: a_val++;
2271: }
2272: }
2273: }
2274: for (i=b_aij->i[0]; i<b_aij->i[aij->B->rmap->n/bs]; i++) {
2275: for (jb=0; jb<bs; jb++) {
2276: for (ib=0; ib<bs; ib++) {
2277: int col = garray[b_aij->j[i]] * bs + jb;
2278: work[col] = PetscMax(PetscAbsScalar(*b_val), work[col]);
2279: b_val++;
2280: }
2281: }
2282: }
2283: } else if (type == REDUCTION_SUM_REALPART || type == REDUCTION_MEAN_REALPART) {
2284: for (i=a_aij->i[0]; i<a_aij->i[aij->A->rmap->n/bs]; i++) {
2285: for (jb=0; jb<bs; jb++) {
2286: for (ib=0; ib<bs; ib++) {
2287: work[A->cmap->rstart + a_aij->j[i] * bs + jb] += PetscRealPart(*a_val);
2288: a_val++;
2289: }
2290: }
2291: }
2292: for (i=b_aij->i[0]; i<b_aij->i[aij->B->rmap->n/bs]; i++) {
2293: for (jb=0; jb<bs; jb++) {
2294: for (ib=0; ib<bs; ib++) {
2295: work[garray[b_aij->j[i]] * bs + jb] += PetscRealPart(*b_val);
2296: b_val++;
2297: }
2298: }
2299: }
2300: } else if (type == REDUCTION_SUM_IMAGINARYPART || type == REDUCTION_MEAN_IMAGINARYPART) {
2301: for (i=a_aij->i[0]; i<a_aij->i[aij->A->rmap->n/bs]; i++) {
2302: for (jb=0; jb<bs; jb++) {
2303: for (ib=0; ib<bs; ib++) {
2304: work[A->cmap->rstart + a_aij->j[i] * bs + jb] += PetscImaginaryPart(*a_val);
2305: a_val++;
2306: }
2307: }
2308: }
2309: for (i=b_aij->i[0]; i<b_aij->i[aij->B->rmap->n/bs]; i++) {
2310: for (jb=0; jb<bs; jb++) {
2311: for (ib=0; ib<bs; ib++) {
2312: work[garray[b_aij->j[i]] * bs + jb] += PetscImaginaryPart(*b_val);
2313: b_val++;
2314: }
2315: }
2316: }
2317: } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Unknown reduction type");
2318: if (type == NORM_INFINITY) {
2319: MPIU_Allreduce(work,reductions,N,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)A));
2320: } else {
2321: MPIU_Allreduce(work,reductions,N,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)A));
2322: }
2323: PetscFree(work);
2324: if (type == NORM_2) {
2325: for (i=0; i<N; i++) reductions[i] = PetscSqrtReal(reductions[i]);
2326: } else if (type == REDUCTION_MEAN_REALPART || type == REDUCTION_MEAN_IMAGINARYPART) {
2327: for (i=0; i<N; i++) reductions[i] /= m;
2328: }
2329: return 0;
2330: }
2332: PetscErrorCode MatInvertBlockDiagonal_MPIBAIJ(Mat A,const PetscScalar **values)
2333: {
2334: Mat_MPIBAIJ *a = (Mat_MPIBAIJ*) A->data;
2336: MatInvertBlockDiagonal(a->A,values);
2337: A->factorerrortype = a->A->factorerrortype;
2338: A->factorerror_zeropivot_value = a->A->factorerror_zeropivot_value;
2339: A->factorerror_zeropivot_row = a->A->factorerror_zeropivot_row;
2340: return 0;
2341: }
2343: PetscErrorCode MatShift_MPIBAIJ(Mat Y,PetscScalar a)
2344: {
2345: Mat_MPIBAIJ *maij = (Mat_MPIBAIJ*)Y->data;
2346: Mat_SeqBAIJ *aij = (Mat_SeqBAIJ*)maij->A->data;
2348: if (!Y->preallocated) {
2349: MatMPIBAIJSetPreallocation(Y,Y->rmap->bs,1,NULL,0,NULL);
2350: } else if (!aij->nz) {
2351: PetscInt nonew = aij->nonew;
2352: MatSeqBAIJSetPreallocation(maij->A,Y->rmap->bs,1,NULL);
2353: aij->nonew = nonew;
2354: }
2355: MatShift_Basic(Y,a);
2356: return 0;
2357: }
2359: PetscErrorCode MatMissingDiagonal_MPIBAIJ(Mat A,PetscBool *missing,PetscInt *d)
2360: {
2361: Mat_MPIBAIJ *a = (Mat_MPIBAIJ*)A->data;
2364: MatMissingDiagonal(a->A,missing,d);
2365: if (d) {
2366: PetscInt rstart;
2367: MatGetOwnershipRange(A,&rstart,NULL);
2368: *d += rstart/A->rmap->bs;
2370: }
2371: return 0;
2372: }
2374: PetscErrorCode MatGetDiagonalBlock_MPIBAIJ(Mat A,Mat *a)
2375: {
2376: *a = ((Mat_MPIBAIJ*)A->data)->A;
2377: return 0;
2378: }
2380: /* -------------------------------------------------------------------*/
2381: static struct _MatOps MatOps_Values = {MatSetValues_MPIBAIJ,
2382: MatGetRow_MPIBAIJ,
2383: MatRestoreRow_MPIBAIJ,
2384: MatMult_MPIBAIJ,
2385: /* 4*/ MatMultAdd_MPIBAIJ,
2386: MatMultTranspose_MPIBAIJ,
2387: MatMultTransposeAdd_MPIBAIJ,
2388: NULL,
2389: NULL,
2390: NULL,
2391: /*10*/ NULL,
2392: NULL,
2393: NULL,
2394: MatSOR_MPIBAIJ,
2395: MatTranspose_MPIBAIJ,
2396: /*15*/ MatGetInfo_MPIBAIJ,
2397: MatEqual_MPIBAIJ,
2398: MatGetDiagonal_MPIBAIJ,
2399: MatDiagonalScale_MPIBAIJ,
2400: MatNorm_MPIBAIJ,
2401: /*20*/ MatAssemblyBegin_MPIBAIJ,
2402: MatAssemblyEnd_MPIBAIJ,
2403: MatSetOption_MPIBAIJ,
2404: MatZeroEntries_MPIBAIJ,
2405: /*24*/ MatZeroRows_MPIBAIJ,
2406: NULL,
2407: NULL,
2408: NULL,
2409: NULL,
2410: /*29*/ MatSetUp_MPIBAIJ,
2411: NULL,
2412: NULL,
2413: MatGetDiagonalBlock_MPIBAIJ,
2414: NULL,
2415: /*34*/ MatDuplicate_MPIBAIJ,
2416: NULL,
2417: NULL,
2418: NULL,
2419: NULL,
2420: /*39*/ MatAXPY_MPIBAIJ,
2421: MatCreateSubMatrices_MPIBAIJ,
2422: MatIncreaseOverlap_MPIBAIJ,
2423: MatGetValues_MPIBAIJ,
2424: MatCopy_MPIBAIJ,
2425: /*44*/ NULL,
2426: MatScale_MPIBAIJ,
2427: MatShift_MPIBAIJ,
2428: NULL,
2429: MatZeroRowsColumns_MPIBAIJ,
2430: /*49*/ NULL,
2431: NULL,
2432: NULL,
2433: NULL,
2434: NULL,
2435: /*54*/ MatFDColoringCreate_MPIXAIJ,
2436: NULL,
2437: MatSetUnfactored_MPIBAIJ,
2438: MatPermute_MPIBAIJ,
2439: MatSetValuesBlocked_MPIBAIJ,
2440: /*59*/ MatCreateSubMatrix_MPIBAIJ,
2441: MatDestroy_MPIBAIJ,
2442: MatView_MPIBAIJ,
2443: NULL,
2444: NULL,
2445: /*64*/ NULL,
2446: NULL,
2447: NULL,
2448: NULL,
2449: NULL,
2450: /*69*/ MatGetRowMaxAbs_MPIBAIJ,
2451: NULL,
2452: NULL,
2453: NULL,
2454: NULL,
2455: /*74*/ NULL,
2456: MatFDColoringApply_BAIJ,
2457: NULL,
2458: NULL,
2459: NULL,
2460: /*79*/ NULL,
2461: NULL,
2462: NULL,
2463: NULL,
2464: MatLoad_MPIBAIJ,
2465: /*84*/ NULL,
2466: NULL,
2467: NULL,
2468: NULL,
2469: NULL,
2470: /*89*/ NULL,
2471: NULL,
2472: NULL,
2473: NULL,
2474: NULL,
2475: /*94*/ NULL,
2476: NULL,
2477: NULL,
2478: NULL,
2479: NULL,
2480: /*99*/ NULL,
2481: NULL,
2482: NULL,
2483: MatConjugate_MPIBAIJ,
2484: NULL,
2485: /*104*/NULL,
2486: MatRealPart_MPIBAIJ,
2487: MatImaginaryPart_MPIBAIJ,
2488: NULL,
2489: NULL,
2490: /*109*/NULL,
2491: NULL,
2492: NULL,
2493: NULL,
2494: MatMissingDiagonal_MPIBAIJ,
2495: /*114*/MatGetSeqNonzeroStructure_MPIBAIJ,
2496: NULL,
2497: MatGetGhosts_MPIBAIJ,
2498: NULL,
2499: NULL,
2500: /*119*/NULL,
2501: NULL,
2502: NULL,
2503: NULL,
2504: MatGetMultiProcBlock_MPIBAIJ,
2505: /*124*/NULL,
2506: MatGetColumnReductions_MPIBAIJ,
2507: MatInvertBlockDiagonal_MPIBAIJ,
2508: NULL,
2509: NULL,
2510: /*129*/ NULL,
2511: NULL,
2512: NULL,
2513: NULL,
2514: NULL,
2515: /*134*/ NULL,
2516: NULL,
2517: NULL,
2518: NULL,
2519: NULL,
2520: /*139*/ MatSetBlockSizes_Default,
2521: NULL,
2522: NULL,
2523: MatFDColoringSetUp_MPIXAIJ,
2524: NULL,
2525: /*144*/MatCreateMPIMatConcatenateSeqMat_MPIBAIJ,
2526: NULL,
2527: NULL,
2528: NULL
2529: };
2531: PETSC_INTERN PetscErrorCode MatConvert_MPIBAIJ_MPISBAIJ(Mat,MatType,MatReuse,Mat*);
2532: PETSC_INTERN PetscErrorCode MatConvert_XAIJ_IS(Mat,MatType,MatReuse,Mat*);
2534: PetscErrorCode MatMPIBAIJSetPreallocationCSR_MPIBAIJ(Mat B,PetscInt bs,const PetscInt ii[],const PetscInt jj[],const PetscScalar V[])
2535: {
2536: PetscInt m,rstart,cstart,cend;
2537: PetscInt i,j,dlen,olen,nz,nz_max=0,*d_nnz=NULL,*o_nnz=NULL;
2538: const PetscInt *JJ =NULL;
2539: PetscScalar *values=NULL;
2540: PetscBool roworiented = ((Mat_MPIBAIJ*)B->data)->roworiented;
2541: PetscBool nooffprocentries;
2543: PetscLayoutSetBlockSize(B->rmap,bs);
2544: PetscLayoutSetBlockSize(B->cmap,bs);
2545: PetscLayoutSetUp(B->rmap);
2546: PetscLayoutSetUp(B->cmap);
2547: PetscLayoutGetBlockSize(B->rmap,&bs);
2548: m = B->rmap->n/bs;
2549: rstart = B->rmap->rstart/bs;
2550: cstart = B->cmap->rstart/bs;
2551: cend = B->cmap->rend/bs;
2554: PetscMalloc2(m,&d_nnz,m,&o_nnz);
2555: for (i=0; i<m; i++) {
2556: nz = ii[i+1] - ii[i];
2558: nz_max = PetscMax(nz_max,nz);
2559: dlen = 0;
2560: olen = 0;
2561: JJ = jj + ii[i];
2562: for (j=0; j<nz; j++) {
2563: if (*JJ < cstart || *JJ >= cend) olen++;
2564: else dlen++;
2565: JJ++;
2566: }
2567: d_nnz[i] = dlen;
2568: o_nnz[i] = olen;
2569: }
2570: MatMPIBAIJSetPreallocation(B,bs,0,d_nnz,0,o_nnz);
2571: PetscFree2(d_nnz,o_nnz);
2573: values = (PetscScalar*)V;
2574: if (!values) {
2575: PetscCalloc1(bs*bs*nz_max,&values);
2576: }
2577: for (i=0; i<m; i++) {
2578: PetscInt row = i + rstart;
2579: PetscInt ncols = ii[i+1] - ii[i];
2580: const PetscInt *icols = jj + ii[i];
2581: if (bs == 1 || !roworiented) { /* block ordering matches the non-nested layout of MatSetValues so we can insert entire rows */
2582: const PetscScalar *svals = values + (V ? (bs*bs*ii[i]) : 0);
2583: MatSetValuesBlocked_MPIBAIJ(B,1,&row,ncols,icols,svals,INSERT_VALUES);
2584: } else { /* block ordering does not match so we can only insert one block at a time. */
2585: PetscInt j;
2586: for (j=0; j<ncols; j++) {
2587: const PetscScalar *svals = values + (V ? (bs*bs*(ii[i]+j)) : 0);
2588: MatSetValuesBlocked_MPIBAIJ(B,1,&row,1,&icols[j],svals,INSERT_VALUES);
2589: }
2590: }
2591: }
2593: if (!V) PetscFree(values);
2594: nooffprocentries = B->nooffprocentries;
2595: B->nooffprocentries = PETSC_TRUE;
2596: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
2597: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
2598: B->nooffprocentries = nooffprocentries;
2600: MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
2601: return 0;
2602: }
2604: /*@C
2605: MatMPIBAIJSetPreallocationCSR - Creates a sparse parallel matrix in BAIJ format using the given nonzero structure and (optional) numerical values
2607: Collective
2609: Input Parameters:
2610: + B - the matrix
2611: . bs - the block size
2612: . i - the indices into j for the start of each local row (starts with zero)
2613: . j - the column indices for each local row (starts with zero) these must be sorted for each row
2614: - v - optional values in the matrix
2616: Level: advanced
2618: Notes:
2619: The order of the entries in values is specified by the MatOption MAT_ROW_ORIENTED. For example, C programs
2620: may want to use the default MAT_ROW_ORIENTED=PETSC_TRUE and use an array v[nnz][bs][bs] where the second index is
2621: over rows within a block and the last index is over columns within a block row. Fortran programs will likely set
2622: MAT_ROW_ORIENTED=PETSC_FALSE and use a Fortran array v(bs,bs,nnz) in which the first index is over rows within a
2623: block column and the second index is over columns within a block.
2625: Though this routine has Preallocation() in the name it also sets the exact nonzero locations of the matrix entries and usually the numerical values as well
2627: .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatMPIBAIJSetPreallocation(), MatCreateAIJ(), MPIAIJ, MatCreateMPIBAIJWithArrays(), MPIBAIJ
2628: @*/
2629: PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat B,PetscInt bs,const PetscInt i[],const PetscInt j[], const PetscScalar v[])
2630: {
2634: PetscTryMethod(B,"MatMPIBAIJSetPreallocationCSR_C",(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]),(B,bs,i,j,v));
2635: return 0;
2636: }
2638: PetscErrorCode MatMPIBAIJSetPreallocation_MPIBAIJ(Mat B,PetscInt bs,PetscInt d_nz,const PetscInt *d_nnz,PetscInt o_nz,const PetscInt *o_nnz)
2639: {
2640: Mat_MPIBAIJ *b;
2641: PetscInt i;
2642: PetscMPIInt size;
2644: MatSetBlockSize(B,PetscAbs(bs));
2645: PetscLayoutSetUp(B->rmap);
2646: PetscLayoutSetUp(B->cmap);
2647: PetscLayoutGetBlockSize(B->rmap,&bs);
2649: if (d_nnz) {
2650: for (i=0; i<B->rmap->n/bs; i++) {
2652: }
2653: }
2654: if (o_nnz) {
2655: for (i=0; i<B->rmap->n/bs; i++) {
2657: }
2658: }
2660: b = (Mat_MPIBAIJ*)B->data;
2661: b->bs2 = bs*bs;
2662: b->mbs = B->rmap->n/bs;
2663: b->nbs = B->cmap->n/bs;
2664: b->Mbs = B->rmap->N/bs;
2665: b->Nbs = B->cmap->N/bs;
2667: for (i=0; i<=b->size; i++) {
2668: b->rangebs[i] = B->rmap->range[i]/bs;
2669: }
2670: b->rstartbs = B->rmap->rstart/bs;
2671: b->rendbs = B->rmap->rend/bs;
2672: b->cstartbs = B->cmap->rstart/bs;
2673: b->cendbs = B->cmap->rend/bs;
2675: #if defined(PETSC_USE_CTABLE)
2676: PetscTableDestroy(&b->colmap);
2677: #else
2678: PetscFree(b->colmap);
2679: #endif
2680: PetscFree(b->garray);
2681: VecDestroy(&b->lvec);
2682: VecScatterDestroy(&b->Mvctx);
2684: /* Because the B will have been resized we simply destroy it and create a new one each time */
2685: MPI_Comm_size(PetscObjectComm((PetscObject)B),&size);
2686: MatDestroy(&b->B);
2687: MatCreate(PETSC_COMM_SELF,&b->B);
2688: MatSetSizes(b->B,B->rmap->n,size > 1 ? B->cmap->N : 0,B->rmap->n,size > 1 ? B->cmap->N : 0);
2689: MatSetType(b->B,MATSEQBAIJ);
2690: PetscLogObjectParent((PetscObject)B,(PetscObject)b->B);
2692: if (!B->preallocated) {
2693: MatCreate(PETSC_COMM_SELF,&b->A);
2694: MatSetSizes(b->A,B->rmap->n,B->cmap->n,B->rmap->n,B->cmap->n);
2695: MatSetType(b->A,MATSEQBAIJ);
2696: PetscLogObjectParent((PetscObject)B,(PetscObject)b->A);
2697: MatStashCreate_Private(PetscObjectComm((PetscObject)B),bs,&B->bstash);
2698: }
2700: MatSeqBAIJSetPreallocation(b->A,bs,d_nz,d_nnz);
2701: MatSeqBAIJSetPreallocation(b->B,bs,o_nz,o_nnz);
2702: B->preallocated = PETSC_TRUE;
2703: B->was_assembled = PETSC_FALSE;
2704: B->assembled = PETSC_FALSE;
2705: return 0;
2706: }
2708: extern PetscErrorCode MatDiagonalScaleLocal_MPIBAIJ(Mat,Vec);
2709: extern PetscErrorCode MatSetHashTableFactor_MPIBAIJ(Mat,PetscReal);
2711: PETSC_INTERN PetscErrorCode MatConvert_MPIBAIJ_MPIAdj(Mat B, MatType newtype,MatReuse reuse,Mat *adj)
2712: {
2713: Mat_MPIBAIJ *b = (Mat_MPIBAIJ*)B->data;
2714: Mat_SeqBAIJ *d = (Mat_SeqBAIJ*) b->A->data,*o = (Mat_SeqBAIJ*) b->B->data;
2715: PetscInt M = B->rmap->n/B->rmap->bs,i,*ii,*jj,cnt,j,k,rstart = B->rmap->rstart/B->rmap->bs;
2716: const PetscInt *id = d->i, *jd = d->j, *io = o->i, *jo = o->j, *garray = b->garray;
2718: PetscMalloc1(M+1,&ii);
2719: ii[0] = 0;
2720: for (i=0; i<M; i++) {
2723: ii[i+1] = ii[i] + id[i+1] - id[i] + io[i+1] - io[i];
2724: /* remove one from count of matrix has diagonal */
2725: for (j=id[i]; j<id[i+1]; j++) {
2726: if (jd[j] == i) {ii[i+1]--;break;}
2727: }
2728: }
2729: PetscMalloc1(ii[M],&jj);
2730: cnt = 0;
2731: for (i=0; i<M; i++) {
2732: for (j=io[i]; j<io[i+1]; j++) {
2733: if (garray[jo[j]] > rstart) break;
2734: jj[cnt++] = garray[jo[j]];
2735: }
2736: for (k=id[i]; k<id[i+1]; k++) {
2737: if (jd[k] != i) {
2738: jj[cnt++] = rstart + jd[k];
2739: }
2740: }
2741: for (; j<io[i+1]; j++) {
2742: jj[cnt++] = garray[jo[j]];
2743: }
2744: }
2745: MatCreateMPIAdj(PetscObjectComm((PetscObject)B),M,B->cmap->N/B->rmap->bs,ii,jj,NULL,adj);
2746: return 0;
2747: }
2749: #include <../src/mat/impls/aij/mpi/mpiaij.h>
2751: PETSC_INTERN PetscErrorCode MatConvert_SeqBAIJ_SeqAIJ(Mat,MatType,MatReuse,Mat*);
2753: PETSC_INTERN PetscErrorCode MatConvert_MPIBAIJ_MPIAIJ(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
2754: {
2755: Mat_MPIBAIJ *a = (Mat_MPIBAIJ*)A->data;
2756: Mat_MPIAIJ *b;
2757: Mat B;
2761: if (reuse == MAT_REUSE_MATRIX) {
2762: B = *newmat;
2763: } else {
2764: MatCreate(PetscObjectComm((PetscObject)A),&B);
2765: MatSetType(B,MATMPIAIJ);
2766: MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
2767: MatSetBlockSizes(B,A->rmap->bs,A->cmap->bs);
2768: MatSeqAIJSetPreallocation(B,0,NULL);
2769: MatMPIAIJSetPreallocation(B,0,NULL,0,NULL);
2770: }
2771: b = (Mat_MPIAIJ*) B->data;
2773: if (reuse == MAT_REUSE_MATRIX) {
2774: MatConvert_SeqBAIJ_SeqAIJ(a->A, MATSEQAIJ, MAT_REUSE_MATRIX, &b->A);
2775: MatConvert_SeqBAIJ_SeqAIJ(a->B, MATSEQAIJ, MAT_REUSE_MATRIX, &b->B);
2776: } else {
2777: MatDestroy(&b->A);
2778: MatDestroy(&b->B);
2779: MatDisAssemble_MPIBAIJ(A);
2780: MatConvert_SeqBAIJ_SeqAIJ(a->A, MATSEQAIJ, MAT_INITIAL_MATRIX, &b->A);
2781: MatConvert_SeqBAIJ_SeqAIJ(a->B, MATSEQAIJ, MAT_INITIAL_MATRIX, &b->B);
2782: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2783: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2784: }
2785: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
2786: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
2788: if (reuse == MAT_INPLACE_MATRIX) {
2789: MatHeaderReplace(A,&B);
2790: } else {
2791: *newmat = B;
2792: }
2793: return 0;
2794: }
2796: /*MC
2797: MATMPIBAIJ - MATMPIBAIJ = "mpibaij" - A matrix type to be used for distributed block sparse matrices.
2799: Options Database Keys:
2800: + -mat_type mpibaij - sets the matrix type to "mpibaij" during a call to MatSetFromOptions()
2801: . -mat_block_size <bs> - set the blocksize used to store the matrix
2802: . -mat_baij_mult_version version - indicate the version of the matrix-vector product to use (0 often indicates using BLAS)
2803: - -mat_use_hash_table <fact> - set hash table factor
2805: Level: beginner
2807: Notes:
2808: MatSetOptions(,MAT_STRUCTURE_ONLY,PETSC_TRUE) may be called for this matrix type. In this no
2809: space is allocated for the nonzero entries and any entries passed with MatSetValues() are ignored
2811: .seealso: MatCreateBAIJ
2812: M*/
2814: PETSC_INTERN PetscErrorCode MatConvert_MPIBAIJ_MPIBSTRM(Mat,MatType,MatReuse,Mat*);
2816: PETSC_EXTERN PetscErrorCode MatCreate_MPIBAIJ(Mat B)
2817: {
2818: Mat_MPIBAIJ *b;
2820: PetscBool flg = PETSC_FALSE;
2822: PetscNewLog(B,&b);
2823: B->data = (void*)b;
2825: PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
2826: B->assembled = PETSC_FALSE;
2828: B->insertmode = NOT_SET_VALUES;
2829: MPI_Comm_rank(PetscObjectComm((PetscObject)B),&b->rank);
2830: MPI_Comm_size(PetscObjectComm((PetscObject)B),&b->size);
2832: /* build local table of row and column ownerships */
2833: PetscMalloc1(b->size+1,&b->rangebs);
2835: /* build cache for off array entries formed */
2836: MatStashCreate_Private(PetscObjectComm((PetscObject)B),1,&B->stash);
2838: b->donotstash = PETSC_FALSE;
2839: b->colmap = NULL;
2840: b->garray = NULL;
2841: b->roworiented = PETSC_TRUE;
2843: /* stuff used in block assembly */
2844: b->barray = NULL;
2846: /* stuff used for matrix vector multiply */
2847: b->lvec = NULL;
2848: b->Mvctx = NULL;
2850: /* stuff for MatGetRow() */
2851: b->rowindices = NULL;
2852: b->rowvalues = NULL;
2853: b->getrowactive = PETSC_FALSE;
2855: /* hash table stuff */
2856: b->ht = NULL;
2857: b->hd = NULL;
2858: b->ht_size = 0;
2859: b->ht_flag = PETSC_FALSE;
2860: b->ht_fact = 0;
2861: b->ht_total_ct = 0;
2862: b->ht_insert_ct = 0;
2864: /* stuff for MatCreateSubMatrices_MPIBAIJ_local() */
2865: b->ijonly = PETSC_FALSE;
2867: PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpibaij_mpiadj_C",MatConvert_MPIBAIJ_MPIAdj);
2868: PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpibaij_mpiaij_C",MatConvert_MPIBAIJ_MPIAIJ);
2869: PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpibaij_mpisbaij_C",MatConvert_MPIBAIJ_MPISBAIJ);
2870: #if defined(PETSC_HAVE_HYPRE)
2871: PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpibaij_hypre_C",MatConvert_AIJ_HYPRE);
2872: #endif
2873: PetscObjectComposeFunction((PetscObject)B,"MatStoreValues_C",MatStoreValues_MPIBAIJ);
2874: PetscObjectComposeFunction((PetscObject)B,"MatRetrieveValues_C",MatRetrieveValues_MPIBAIJ);
2875: PetscObjectComposeFunction((PetscObject)B,"MatMPIBAIJSetPreallocation_C",MatMPIBAIJSetPreallocation_MPIBAIJ);
2876: PetscObjectComposeFunction((PetscObject)B,"MatMPIBAIJSetPreallocationCSR_C",MatMPIBAIJSetPreallocationCSR_MPIBAIJ);
2877: PetscObjectComposeFunction((PetscObject)B,"MatDiagonalScaleLocal_C",MatDiagonalScaleLocal_MPIBAIJ);
2878: PetscObjectComposeFunction((PetscObject)B,"MatSetHashTableFactor_C",MatSetHashTableFactor_MPIBAIJ);
2879: PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpibaij_is_C",MatConvert_XAIJ_IS);
2880: PetscObjectChangeTypeName((PetscObject)B,MATMPIBAIJ);
2882: PetscOptionsBegin(PetscObjectComm((PetscObject)B),NULL,"Options for loading MPIBAIJ matrix 1","Mat");
2883: PetscOptionsName("-mat_use_hash_table","Use hash table to save time in constructing matrix","MatSetOption",&flg);
2884: if (flg) {
2885: PetscReal fact = 1.39;
2886: MatSetOption(B,MAT_USE_HASH_TABLE,PETSC_TRUE);
2887: PetscOptionsReal("-mat_use_hash_table","Use hash table factor","MatMPIBAIJSetHashTableFactor",fact,&fact,NULL);
2888: if (fact <= 1.0) fact = 1.39;
2889: MatMPIBAIJSetHashTableFactor(B,fact);
2890: PetscInfo(B,"Hash table Factor used %5.2g\n",(double)fact);
2891: }
2892: PetscOptionsEnd();
2893: return 0;
2894: }
2896: /*MC
2897: MATBAIJ - MATBAIJ = "baij" - A matrix type to be used for block sparse matrices.
2899: This matrix type is identical to MATSEQBAIJ when constructed with a single process communicator,
2900: and MATMPIBAIJ otherwise.
2902: Options Database Keys:
2903: . -mat_type baij - sets the matrix type to "baij" during a call to MatSetFromOptions()
2905: Level: beginner
2907: .seealso: MatCreateBAIJ(),MATSEQBAIJ,MATMPIBAIJ, MatMPIBAIJSetPreallocation(), MatMPIBAIJSetPreallocationCSR()
2908: M*/
2910: /*@C
2911: MatMPIBAIJSetPreallocation - Allocates memory for a sparse parallel matrix in block AIJ format
2912: (block compressed row). For good matrix assembly performance
2913: the user should preallocate the matrix storage by setting the parameters
2914: d_nz (or d_nnz) and o_nz (or o_nnz). By setting these parameters accurately,
2915: performance can be increased by more than a factor of 50.
2917: Collective on Mat
2919: Input Parameters:
2920: + B - the matrix
2921: . bs - size of block, the blocks are ALWAYS square. One can use MatSetBlockSizes() to set a different row and column blocksize but the row
2922: blocksize always defines the size of the blocks. The column blocksize sets the blocksize of the vectors obtained with MatCreateVecs()
2923: . d_nz - number of block nonzeros per block row in diagonal portion of local
2924: submatrix (same for all local rows)
2925: . d_nnz - array containing the number of block nonzeros in the various block rows
2926: of the in diagonal portion of the local (possibly different for each block
2927: row) or NULL. If you plan to factor the matrix you must leave room for the diagonal entry and
2928: set it even if it is zero.
2929: . o_nz - number of block nonzeros per block row in the off-diagonal portion of local
2930: submatrix (same for all local rows).
2931: - o_nnz - array containing the number of nonzeros in the various block rows of the
2932: off-diagonal portion of the local submatrix (possibly different for
2933: each block row) or NULL.
2935: If the *_nnz parameter is given then the *_nz parameter is ignored
2937: Options Database Keys:
2938: + -mat_block_size - size of the blocks to use
2939: - -mat_use_hash_table <fact> - set hash table factor
2941: Notes:
2942: If PETSC_DECIDE or PETSC_DETERMINE is used for a particular argument on one processor
2943: than it must be used on all processors that share the object for that argument.
2945: Storage Information:
2946: For a square global matrix we define each processor's diagonal portion
2947: to be its local rows and the corresponding columns (a square submatrix);
2948: each processor's off-diagonal portion encompasses the remainder of the
2949: local matrix (a rectangular submatrix).
2951: The user can specify preallocated storage for the diagonal part of
2952: the local submatrix with either d_nz or d_nnz (not both). Set
2953: d_nz=PETSC_DEFAULT and d_nnz=NULL for PETSc to control dynamic
2954: memory allocation. Likewise, specify preallocated storage for the
2955: off-diagonal part of the local submatrix with o_nz or o_nnz (not both).
2957: Consider a processor that owns rows 3, 4 and 5 of a parallel matrix. In
2958: the figure below we depict these three local rows and all columns (0-11).
2960: .vb
2961: 0 1 2 3 4 5 6 7 8 9 10 11
2962: --------------------------
2963: row 3 |o o o d d d o o o o o o
2964: row 4 |o o o d d d o o o o o o
2965: row 5 |o o o d d d o o o o o o
2966: --------------------------
2967: .ve
2969: Thus, any entries in the d locations are stored in the d (diagonal)
2970: submatrix, and any entries in the o locations are stored in the
2971: o (off-diagonal) submatrix. Note that the d and the o submatrices are
2972: stored simply in the MATSEQBAIJ format for compressed row storage.
2974: Now d_nz should indicate the number of block nonzeros per row in the d matrix,
2975: and o_nz should indicate the number of block nonzeros per row in the o matrix.
2976: In general, for PDE problems in which most nonzeros are near the diagonal,
2977: one expects d_nz >> o_nz. For large problems you MUST preallocate memory
2978: or you will get TERRIBLE performance; see the users' manual chapter on
2979: matrices.
2981: You can call MatGetInfo() to get information on how effective the preallocation was;
2982: for example the fields mallocs,nz_allocated,nz_used,nz_unneeded;
2983: You can also run with the option -info and look for messages with the string
2984: malloc in them to see if additional memory allocation was needed.
2986: Level: intermediate
2988: .seealso: MatCreate(), MatCreateSeqBAIJ(), MatSetValues(), MatCreateBAIJ(), MatMPIBAIJSetPreallocationCSR(), PetscSplitOwnership()
2989: @*/
2990: PetscErrorCode MatMPIBAIJSetPreallocation(Mat B,PetscInt bs,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
2991: {
2995: PetscTryMethod(B,"MatMPIBAIJSetPreallocation_C",(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,bs,d_nz,d_nnz,o_nz,o_nnz));
2996: return 0;
2997: }
2999: /*@C
3000: MatCreateBAIJ - Creates a sparse parallel matrix in block AIJ format
3001: (block compressed row). For good matrix assembly performance
3002: the user should preallocate the matrix storage by setting the parameters
3003: d_nz (or d_nnz) and o_nz (or o_nnz). By setting these parameters accurately,
3004: performance can be increased by more than a factor of 50.
3006: Collective
3008: Input Parameters:
3009: + comm - MPI communicator
3010: . bs - size of block, the blocks are ALWAYS square. One can use MatSetBlockSizes() to set a different row and column blocksize but the row
3011: blocksize always defines the size of the blocks. The column blocksize sets the blocksize of the vectors obtained with MatCreateVecs()
3012: . m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
3013: This value should be the same as the local size used in creating the
3014: y vector for the matrix-vector product y = Ax.
3015: . n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
3016: This value should be the same as the local size used in creating the
3017: x vector for the matrix-vector product y = Ax.
3018: . M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
3019: . N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
3020: . d_nz - number of nonzero blocks per block row in diagonal portion of local
3021: submatrix (same for all local rows)
3022: . d_nnz - array containing the number of nonzero blocks in the various block rows
3023: of the in diagonal portion of the local (possibly different for each block
3024: row) or NULL. If you plan to factor the matrix you must leave room for the diagonal entry
3025: and set it even if it is zero.
3026: . o_nz - number of nonzero blocks per block row in the off-diagonal portion of local
3027: submatrix (same for all local rows).
3028: - o_nnz - array containing the number of nonzero blocks in the various block rows of the
3029: off-diagonal portion of the local submatrix (possibly different for
3030: each block row) or NULL.
3032: Output Parameter:
3033: . A - the matrix
3035: Options Database Keys:
3036: + -mat_block_size - size of the blocks to use
3037: - -mat_use_hash_table <fact> - set hash table factor
3039: It is recommended that one use the MatCreate(), MatSetType() and/or MatSetFromOptions(),
3040: MatXXXXSetPreallocation() paradigm instead of this routine directly.
3041: [MatXXXXSetPreallocation() is, for example, MatSeqAIJSetPreallocation]
3043: Notes:
3044: If the *_nnz parameter is given then the *_nz parameter is ignored
3046: A nonzero block is any block that as 1 or more nonzeros in it
3048: The user MUST specify either the local or global matrix dimensions
3049: (possibly both).
3051: If PETSC_DECIDE or PETSC_DETERMINE is used for a particular argument on one processor
3052: than it must be used on all processors that share the object for that argument.
3054: Storage Information:
3055: For a square global matrix we define each processor's diagonal portion
3056: to be its local rows and the corresponding columns (a square submatrix);
3057: each processor's off-diagonal portion encompasses the remainder of the
3058: local matrix (a rectangular submatrix).
3060: The user can specify preallocated storage for the diagonal part of
3061: the local submatrix with either d_nz or d_nnz (not both). Set
3062: d_nz=PETSC_DEFAULT and d_nnz=NULL for PETSc to control dynamic
3063: memory allocation. Likewise, specify preallocated storage for the
3064: off-diagonal part of the local submatrix with o_nz or o_nnz (not both).
3066: Consider a processor that owns rows 3, 4 and 5 of a parallel matrix. In
3067: the figure below we depict these three local rows and all columns (0-11).
3069: .vb
3070: 0 1 2 3 4 5 6 7 8 9 10 11
3071: --------------------------
3072: row 3 |o o o d d d o o o o o o
3073: row 4 |o o o d d d o o o o o o
3074: row 5 |o o o d d d o o o o o o
3075: --------------------------
3076: .ve
3078: Thus, any entries in the d locations are stored in the d (diagonal)
3079: submatrix, and any entries in the o locations are stored in the
3080: o (off-diagonal) submatrix. Note that the d and the o submatrices are
3081: stored simply in the MATSEQBAIJ format for compressed row storage.
3083: Now d_nz should indicate the number of block nonzeros per row in the d matrix,
3084: and o_nz should indicate the number of block nonzeros per row in the o matrix.
3085: In general, for PDE problems in which most nonzeros are near the diagonal,
3086: one expects d_nz >> o_nz. For large problems you MUST preallocate memory
3087: or you will get TERRIBLE performance; see the users' manual chapter on
3088: matrices.
3090: Level: intermediate
3092: .seealso: MatCreate(), MatCreateSeqBAIJ(), MatSetValues(), MatCreateBAIJ(), MatMPIBAIJSetPreallocation(), MatMPIBAIJSetPreallocationCSR()
3093: @*/
3094: PetscErrorCode MatCreateBAIJ(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)
3095: {
3096: PetscMPIInt size;
3098: MatCreate(comm,A);
3099: MatSetSizes(*A,m,n,M,N);
3100: MPI_Comm_size(comm,&size);
3101: if (size > 1) {
3102: MatSetType(*A,MATMPIBAIJ);
3103: MatMPIBAIJSetPreallocation(*A,bs,d_nz,d_nnz,o_nz,o_nnz);
3104: } else {
3105: MatSetType(*A,MATSEQBAIJ);
3106: MatSeqBAIJSetPreallocation(*A,bs,d_nz,d_nnz);
3107: }
3108: return 0;
3109: }
3111: static PetscErrorCode MatDuplicate_MPIBAIJ(Mat matin,MatDuplicateOption cpvalues,Mat *newmat)
3112: {
3113: Mat mat;
3114: Mat_MPIBAIJ *a,*oldmat = (Mat_MPIBAIJ*)matin->data;
3115: PetscInt len=0;
3117: *newmat = NULL;
3118: MatCreate(PetscObjectComm((PetscObject)matin),&mat);
3119: MatSetSizes(mat,matin->rmap->n,matin->cmap->n,matin->rmap->N,matin->cmap->N);
3120: MatSetType(mat,((PetscObject)matin)->type_name);
3122: mat->factortype = matin->factortype;
3123: mat->preallocated = PETSC_TRUE;
3124: mat->assembled = PETSC_TRUE;
3125: mat->insertmode = NOT_SET_VALUES;
3127: a = (Mat_MPIBAIJ*)mat->data;
3128: mat->rmap->bs = matin->rmap->bs;
3129: a->bs2 = oldmat->bs2;
3130: a->mbs = oldmat->mbs;
3131: a->nbs = oldmat->nbs;
3132: a->Mbs = oldmat->Mbs;
3133: a->Nbs = oldmat->Nbs;
3135: PetscLayoutReference(matin->rmap,&mat->rmap);
3136: PetscLayoutReference(matin->cmap,&mat->cmap);
3138: a->size = oldmat->size;
3139: a->rank = oldmat->rank;
3140: a->donotstash = oldmat->donotstash;
3141: a->roworiented = oldmat->roworiented;
3142: a->rowindices = NULL;
3143: a->rowvalues = NULL;
3144: a->getrowactive = PETSC_FALSE;
3145: a->barray = NULL;
3146: a->rstartbs = oldmat->rstartbs;
3147: a->rendbs = oldmat->rendbs;
3148: a->cstartbs = oldmat->cstartbs;
3149: a->cendbs = oldmat->cendbs;
3151: /* hash table stuff */
3152: a->ht = NULL;
3153: a->hd = NULL;
3154: a->ht_size = 0;
3155: a->ht_flag = oldmat->ht_flag;
3156: a->ht_fact = oldmat->ht_fact;
3157: a->ht_total_ct = 0;
3158: a->ht_insert_ct = 0;
3160: PetscArraycpy(a->rangebs,oldmat->rangebs,a->size+1);
3161: if (oldmat->colmap) {
3162: #if defined(PETSC_USE_CTABLE)
3163: PetscTableCreateCopy(oldmat->colmap,&a->colmap);
3164: #else
3165: PetscMalloc1(a->Nbs,&a->colmap);
3166: PetscLogObjectMemory((PetscObject)mat,(a->Nbs)*sizeof(PetscInt));
3167: PetscArraycpy(a->colmap,oldmat->colmap,a->Nbs);
3168: #endif
3169: } else a->colmap = NULL;
3171: if (oldmat->garray && (len = ((Mat_SeqBAIJ*)(oldmat->B->data))->nbs)) {
3172: PetscMalloc1(len,&a->garray);
3173: PetscLogObjectMemory((PetscObject)mat,len*sizeof(PetscInt));
3174: PetscArraycpy(a->garray,oldmat->garray,len);
3175: } else a->garray = NULL;
3177: MatStashCreate_Private(PetscObjectComm((PetscObject)matin),matin->rmap->bs,&mat->bstash);
3178: VecDuplicate(oldmat->lvec,&a->lvec);
3179: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->lvec);
3180: VecScatterCopy(oldmat->Mvctx,&a->Mvctx);
3181: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->Mvctx);
3183: MatDuplicate(oldmat->A,cpvalues,&a->A);
3184: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->A);
3185: MatDuplicate(oldmat->B,cpvalues,&a->B);
3186: PetscLogObjectParent((PetscObject)mat,(PetscObject)a->B);
3187: PetscFunctionListDuplicate(((PetscObject)matin)->qlist,&((PetscObject)mat)->qlist);
3188: *newmat = mat;
3189: return 0;
3190: }
3192: /* Used for both MPIBAIJ and MPISBAIJ matrices */
3193: PetscErrorCode MatLoad_MPIBAIJ_Binary(Mat mat,PetscViewer viewer)
3194: {
3195: PetscInt header[4],M,N,nz,bs,m,n,mbs,nbs,rows,cols,sum,i,j,k;
3196: PetscInt *rowidxs,*colidxs,rs,cs,ce;
3197: PetscScalar *matvals;
3199: PetscViewerSetUp(viewer);
3201: /* read in matrix header */
3202: PetscViewerBinaryRead(viewer,header,4,NULL,PETSC_INT);
3204: M = header[1]; N = header[2]; nz = header[3];
3209: /* set block sizes from the viewer's .info file */
3210: MatLoad_Binary_BlockSizes(mat,viewer);
3211: /* set local sizes if not set already */
3212: if (mat->rmap->n < 0 && M == N) mat->rmap->n = mat->cmap->n;
3213: if (mat->cmap->n < 0 && M == N) mat->cmap->n = mat->rmap->n;
3214: /* set global sizes if not set already */
3215: if (mat->rmap->N < 0) mat->rmap->N = M;
3216: if (mat->cmap->N < 0) mat->cmap->N = N;
3217: PetscLayoutSetUp(mat->rmap);
3218: PetscLayoutSetUp(mat->cmap);
3220: /* check if the matrix sizes are correct */
3221: MatGetSize(mat,&rows,&cols);
3223: MatGetBlockSize(mat,&bs);
3224: MatGetLocalSize(mat,&m,&n);
3225: PetscLayoutGetRange(mat->rmap,&rs,NULL);
3226: PetscLayoutGetRange(mat->cmap,&cs,&ce);
3227: mbs = m/bs; nbs = n/bs;
3229: /* read in row lengths and build row indices */
3230: PetscMalloc1(m+1,&rowidxs);
3231: PetscViewerBinaryReadAll(viewer,rowidxs+1,m,PETSC_DECIDE,M,PETSC_INT);
3232: rowidxs[0] = 0; for (i=0; i<m; i++) rowidxs[i+1] += rowidxs[i];
3233: MPIU_Allreduce(&rowidxs[m],&sum,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)viewer));
3236: /* read in column indices and matrix values */
3237: PetscMalloc2(rowidxs[m],&colidxs,rowidxs[m],&matvals);
3238: PetscViewerBinaryReadAll(viewer,colidxs,rowidxs[m],PETSC_DETERMINE,PETSC_DETERMINE,PETSC_INT);
3239: PetscViewerBinaryReadAll(viewer,matvals,rowidxs[m],PETSC_DETERMINE,PETSC_DETERMINE,PETSC_SCALAR);
3241: { /* preallocate matrix storage */
3242: PetscBT bt; /* helper bit set to count diagonal nonzeros */
3243: PetscHSetI ht; /* helper hash set to count off-diagonal nonzeros */
3244: PetscBool sbaij,done;
3245: PetscInt *d_nnz,*o_nnz;
3247: PetscBTCreate(nbs,&bt);
3248: PetscHSetICreate(&ht);
3249: PetscCalloc2(mbs,&d_nnz,mbs,&o_nnz);
3250: PetscObjectTypeCompare((PetscObject)mat,MATMPISBAIJ,&sbaij);
3251: for (i=0; i<mbs; i++) {
3252: PetscBTMemzero(nbs,bt);
3253: PetscHSetIClear(ht);
3254: for (k=0; k<bs; k++) {
3255: PetscInt row = bs*i + k;
3256: for (j=rowidxs[row]; j<rowidxs[row+1]; j++) {
3257: PetscInt col = colidxs[j];
3258: if (!sbaij || col >= row) {
3259: if (col >= cs && col < ce) {
3260: if (!PetscBTLookupSet(bt,(col-cs)/bs)) d_nnz[i]++;
3261: } else {
3262: PetscHSetIQueryAdd(ht,col/bs,&done);
3263: if (done) o_nnz[i]++;
3264: }
3265: }
3266: }
3267: }
3268: }
3269: PetscBTDestroy(&bt);
3270: PetscHSetIDestroy(&ht);
3271: MatMPIBAIJSetPreallocation(mat,bs,0,d_nnz,0,o_nnz);
3272: MatMPISBAIJSetPreallocation(mat,bs,0,d_nnz,0,o_nnz);
3273: PetscFree2(d_nnz,o_nnz);
3274: }
3276: /* store matrix values */
3277: for (i=0; i<m; i++) {
3278: PetscInt row = rs + i, s = rowidxs[i], e = rowidxs[i+1];
3279: (*mat->ops->setvalues)(mat,1,&row,e-s,colidxs+s,matvals+s,INSERT_VALUES);
3280: }
3282: PetscFree(rowidxs);
3283: PetscFree2(colidxs,matvals);
3284: MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
3285: MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);
3286: return 0;
3287: }
3289: PetscErrorCode MatLoad_MPIBAIJ(Mat mat,PetscViewer viewer)
3290: {
3291: PetscBool isbinary;
3293: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
3295: MatLoad_MPIBAIJ_Binary(mat,viewer);
3296: return 0;
3297: }
3299: /*@
3300: MatMPIBAIJSetHashTableFactor - Sets the factor required to compute the size of the HashTable.
3302: Input Parameters:
3303: + mat - the matrix
3304: - fact - factor
3306: Not Collective, each process can use a different factor
3308: Level: advanced
3310: Notes:
3311: This can also be set by the command line option: -mat_use_hash_table <fact>
3313: .seealso: MatSetOption()
3314: @*/
3315: PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat mat,PetscReal fact)
3316: {
3317: PetscTryMethod(mat,"MatSetHashTableFactor_C",(Mat,PetscReal),(mat,fact));
3318: return 0;
3319: }
3321: PetscErrorCode MatSetHashTableFactor_MPIBAIJ(Mat mat,PetscReal fact)
3322: {
3323: Mat_MPIBAIJ *baij;
3325: baij = (Mat_MPIBAIJ*)mat->data;
3326: baij->ht_fact = fact;
3327: return 0;
3328: }
3330: PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat A,Mat *Ad,Mat *Ao,const PetscInt *colmap[])
3331: {
3332: Mat_MPIBAIJ *a = (Mat_MPIBAIJ*)A->data;
3333: PetscBool flg;
3335: PetscObjectTypeCompare((PetscObject)A,MATMPIBAIJ,&flg);
3337: if (Ad) *Ad = a->A;
3338: if (Ao) *Ao = a->B;
3339: if (colmap) *colmap = a->garray;
3340: return 0;
3341: }
3343: /*
3344: Special version for direct calls from Fortran (to eliminate two function call overheads
3345: */
3346: #if defined(PETSC_HAVE_FORTRAN_CAPS)
3347: #define matmpibaijsetvaluesblocked_ MATMPIBAIJSETVALUESBLOCKED
3348: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
3349: #define matmpibaijsetvaluesblocked_ matmpibaijsetvaluesblocked
3350: #endif
3352: /*@C
3353: MatMPIBAIJSetValuesBlocked - Direct Fortran call to replace call to MatSetValuesBlocked()
3355: Collective on Mat
3357: Input Parameters:
3358: + mat - the matrix
3359: . min - number of input rows
3360: . im - input rows
3361: . nin - number of input columns
3362: . in - input columns
3363: . v - numerical values input
3364: - addvin - INSERT_VALUES or ADD_VALUES
3366: Notes:
3367: This has a complete copy of MatSetValuesBlocked_MPIBAIJ() which is terrible code un-reuse.
3369: Level: advanced
3371: .seealso: MatSetValuesBlocked()
3372: @*/
3373: PetscErrorCode matmpibaijsetvaluesblocked_(Mat *matin,PetscInt *min,const PetscInt im[],PetscInt *nin,const PetscInt in[],const MatScalar v[],InsertMode *addvin)
3374: {
3375: /* convert input arguments to C version */
3376: Mat mat = *matin;
3377: PetscInt m = *min, n = *nin;
3378: InsertMode addv = *addvin;
3380: Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)mat->data;
3381: const MatScalar *value;
3382: MatScalar *barray = baij->barray;
3383: PetscBool roworiented = baij->roworiented;
3384: PetscInt i,j,ii,jj,row,col,rstart=baij->rstartbs;
3385: PetscInt rend=baij->rendbs,cstart=baij->cstartbs,stepval;
3386: PetscInt cend=baij->cendbs,bs=mat->rmap->bs,bs2=baij->bs2;
3388: /* tasks normally handled by MatSetValuesBlocked() */
3389: if (mat->insertmode == NOT_SET_VALUES) mat->insertmode = addv;
3392: if (mat->assembled) {
3393: mat->was_assembled = PETSC_TRUE;
3394: mat->assembled = PETSC_FALSE;
3395: }
3396: PetscLogEventBegin(MAT_SetValues,mat,0,0,0);
3398: if (!barray) {
3399: PetscMalloc1(bs2,&barray);
3400: baij->barray = barray;
3401: }
3403: if (roworiented) stepval = (n-1)*bs;
3404: else stepval = (m-1)*bs;
3406: for (i=0; i<m; i++) {
3407: if (im[i] < 0) continue;
3409: if (im[i] >= rstart && im[i] < rend) {
3410: row = im[i] - rstart;
3411: for (j=0; j<n; j++) {
3412: /* If NumCol = 1 then a copy is not required */
3413: if ((roworiented) && (n == 1)) {
3414: barray = (MatScalar*)v + i*bs2;
3415: } else if ((!roworiented) && (m == 1)) {
3416: barray = (MatScalar*)v + j*bs2;
3417: } else { /* Here a copy is required */
3418: if (roworiented) {
3419: value = v + i*(stepval+bs)*bs + j*bs;
3420: } else {
3421: value = v + j*(stepval+bs)*bs + i*bs;
3422: }
3423: for (ii=0; ii<bs; ii++,value+=stepval) {
3424: for (jj=0; jj<bs; jj++) {
3425: *barray++ = *value++;
3426: }
3427: }
3428: barray -=bs2;
3429: }
3431: if (in[j] >= cstart && in[j] < cend) {
3432: col = in[j] - cstart;
3433: MatSetValuesBlocked_SeqBAIJ_Inlined(baij->A,row,col,barray,addv,im[i],in[j]);
3434: } else if (in[j] < 0) continue;
3436: else {
3437: if (mat->was_assembled) {
3438: if (!baij->colmap) {
3439: MatCreateColmap_MPIBAIJ_Private(mat);
3440: }
3442: #if defined(PETSC_USE_DEBUG)
3443: #if defined(PETSC_USE_CTABLE)
3444: { PetscInt data;
3445: PetscTableFind(baij->colmap,in[j]+1,&data);
3447: }
3448: #else
3450: #endif
3451: #endif
3452: #if defined(PETSC_USE_CTABLE)
3453: PetscTableFind(baij->colmap,in[j]+1,&col);
3454: col = (col - 1)/bs;
3455: #else
3456: col = (baij->colmap[in[j]] - 1)/bs;
3457: #endif
3458: if (col < 0 && !((Mat_SeqBAIJ*)(baij->A->data))->nonew) {
3459: MatDisAssemble_MPIBAIJ(mat);
3460: col = in[j];
3461: }
3462: } else col = in[j];
3463: MatSetValuesBlocked_SeqBAIJ_Inlined(baij->B,row,col,barray,addv,im[i],in[j]);
3464: }
3465: }
3466: } else {
3467: if (!baij->donotstash) {
3468: if (roworiented) {
3469: MatStashValuesRowBlocked_Private(&mat->bstash,im[i],n,in,v,m,n,i);
3470: } else {
3471: MatStashValuesColBlocked_Private(&mat->bstash,im[i],n,in,v,m,n,i);
3472: }
3473: }
3474: }
3475: }
3477: /* task normally handled by MatSetValuesBlocked() */
3478: PetscLogEventEnd(MAT_SetValues,mat,0,0,0);
3479: return 0;
3480: }
3482: /*@
3483: MatCreateMPIBAIJWithArrays - creates a MPI BAIJ matrix using arrays that contain in standard block
3484: CSR format the local rows.
3486: Collective
3488: Input Parameters:
3489: + comm - MPI communicator
3490: . bs - the block size, only a block size of 1 is supported
3491: . m - number of local rows (Cannot be PETSC_DECIDE)
3492: . n - This value should be the same as the local size used in creating the
3493: x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have
3494: calculated if N is given) For square matrices n is almost always m.
3495: . M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
3496: . N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
3497: . i - row indices; that is i[0] = 0, i[row] = i[row-1] + number of block elements in that rowth block row of the matrix
3498: . j - column indices
3499: - a - matrix values
3501: Output Parameter:
3502: . mat - the matrix
3504: Level: intermediate
3506: Notes:
3507: The i, j, and a arrays ARE copied by this routine into the internal format used by PETSc;
3508: thus you CANNOT change the matrix entries by changing the values of a[] after you have
3509: called this routine. Use MatCreateMPIAIJWithSplitArrays() to avoid needing to copy the arrays.
3511: The order of the entries in values is the same as the block compressed sparse row storage format; that is, it is
3512: the same as a three dimensional array in Fortran values(bs,bs,nnz) that contains the first column of the first
3513: block, followed by the second column of the first block etc etc. That is, the blocks are contiguous in memory
3514: with column-major ordering within blocks.
3516: The i and j indices are 0 based, and i indices are indices corresponding to the local j array.
3518: .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatMPIAIJSetPreallocation(), MatMPIAIJSetPreallocationCSR(),
3519: MPIAIJ, MatCreateAIJ(), MatCreateMPIAIJWithSplitArrays()
3520: @*/
3521: PetscErrorCode MatCreateMPIBAIJWithArrays(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt i[],const PetscInt j[],const PetscScalar a[],Mat *mat)
3522: {
3525: MatCreate(comm,mat);
3526: MatSetSizes(*mat,m,n,M,N);
3527: MatSetType(*mat,MATMPIBAIJ);
3528: MatSetBlockSize(*mat,bs);
3529: MatSetUp(*mat);
3530: MatSetOption(*mat,MAT_ROW_ORIENTED,PETSC_FALSE);
3531: MatMPIBAIJSetPreallocationCSR(*mat,bs,i,j,a);
3532: MatSetOption(*mat,MAT_ROW_ORIENTED,PETSC_TRUE);
3533: return 0;
3534: }
3536: PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPIBAIJ(MPI_Comm comm,Mat inmat,PetscInt n,MatReuse scall,Mat *outmat)
3537: {
3539: PetscInt m,N,i,rstart,nnz,Ii,bs,cbs;
3540: PetscInt *indx;
3541: PetscScalar *values;
3543: MatGetSize(inmat,&m,&N);
3544: if (scall == MAT_INITIAL_MATRIX) { /* symbolic phase */
3545: Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)inmat->data;
3546: PetscInt *dnz,*onz,mbs,Nbs,nbs;
3547: PetscInt *bindx,rmax=a->rmax,j;
3548: PetscMPIInt rank,size;
3550: MatGetBlockSizes(inmat,&bs,&cbs);
3551: mbs = m/bs; Nbs = N/cbs;
3552: if (n == PETSC_DECIDE) {
3553: PetscSplitOwnershipBlock(comm,cbs,&n,&N);
3554: }
3555: nbs = n/cbs;
3557: PetscMalloc1(rmax,&bindx);
3558: MatPreallocateInitialize(comm,mbs,nbs,dnz,onz); /* inline function, output __end and __rstart are used below */
3560: MPI_Comm_rank(comm,&rank);
3561: MPI_Comm_rank(comm,&size);
3562: if (rank == size-1) {
3563: /* Check sum(nbs) = Nbs */
3565: }
3567: rstart = __rstart; /* block rstart of *outmat; see inline function MatPreallocateInitialize */
3568: for (i=0; i<mbs; i++) {
3569: MatGetRow_SeqBAIJ(inmat,i*bs,&nnz,&indx,NULL); /* non-blocked nnz and indx */
3570: nnz = nnz/bs;
3571: for (j=0; j<nnz; j++) bindx[j] = indx[j*bs]/bs;
3572: MatPreallocateSet(i+rstart,nnz,bindx,dnz,onz);
3573: MatRestoreRow_SeqBAIJ(inmat,i*bs,&nnz,&indx,NULL);
3574: }
3575: PetscFree(bindx);
3577: MatCreate(comm,outmat);
3578: MatSetSizes(*outmat,m,n,PETSC_DETERMINE,PETSC_DETERMINE);
3579: MatSetBlockSizes(*outmat,bs,cbs);
3580: MatSetType(*outmat,MATBAIJ);
3581: MatSeqBAIJSetPreallocation(*outmat,bs,0,dnz);
3582: MatMPIBAIJSetPreallocation(*outmat,bs,0,dnz,0,onz);
3583: MatPreallocateFinalize(dnz,onz);
3584: MatSetOption(*outmat,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE);
3585: }
3587: /* numeric phase */
3588: MatGetBlockSizes(inmat,&bs,&cbs);
3589: MatGetOwnershipRange(*outmat,&rstart,NULL);
3591: for (i=0; i<m; i++) {
3592: MatGetRow_SeqBAIJ(inmat,i,&nnz,&indx,&values);
3593: Ii = i + rstart;
3594: MatSetValues(*outmat,1,&Ii,nnz,indx,values,INSERT_VALUES);
3595: MatRestoreRow_SeqBAIJ(inmat,i,&nnz,&indx,&values);
3596: }
3597: MatAssemblyBegin(*outmat,MAT_FINAL_ASSEMBLY);
3598: MatAssemblyEnd(*outmat,MAT_FINAL_ASSEMBLY);
3599: return 0;
3600: }