Actual source code: mpisbaij.c
2: #include <../src/mat/impls/baij/mpi/mpibaij.h> /*I "petscmat.h" I*/
3: #include <../src/mat/impls/sbaij/mpi/mpisbaij.h>
4: #include <../src/mat/impls/sbaij/seq/sbaij.h>
5: #include <petscblaslapack.h>
7: extern PetscErrorCode MatSetUpMultiply_MPISBAIJ(Mat);
8: extern PetscErrorCode MatSetUpMultiply_MPISBAIJ_2comm(Mat);
9: extern PetscErrorCode DisAssemble_MPISBAIJ(Mat);
10: extern PetscErrorCode MatIncreaseOverlap_MPISBAIJ(Mat,PetscInt,IS[],PetscInt);
11: extern PetscErrorCode MatGetValues_SeqSBAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar []);
12: extern PetscErrorCode MatGetValues_SeqBAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar []);
13: extern PetscErrorCode MatSetValues_SeqSBAIJ(Mat,PetscInt,const PetscInt [],PetscInt,const PetscInt [],const PetscScalar [],InsertMode);
14: extern PetscErrorCode MatSetValuesBlocked_SeqSBAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
15: extern PetscErrorCode MatSetValuesBlocked_SeqBAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
16: extern PetscErrorCode MatGetRow_SeqSBAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
17: extern PetscErrorCode MatRestoreRow_SeqSBAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
18: extern PetscErrorCode MatZeroRows_SeqSBAIJ(Mat,IS,PetscScalar*,Vec,Vec);
19: extern PetscErrorCode MatZeroRows_SeqBAIJ(Mat,IS,PetscScalar *,Vec,Vec);
20: extern PetscErrorCode MatGetRowMaxAbs_MPISBAIJ(Mat,Vec,PetscInt[]);
21: extern PetscErrorCode MatSOR_MPISBAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
23: EXTERN_C_BEGIN
26: PetscErrorCode MatStoreValues_MPISBAIJ(Mat mat)
27: {
28: Mat_MPISBAIJ *aij = (Mat_MPISBAIJ *)mat->data;
32: MatStoreValues(aij->A);
33: MatStoreValues(aij->B);
34: return(0);
35: }
36: EXTERN_C_END
38: EXTERN_C_BEGIN
41: PetscErrorCode MatRetrieveValues_MPISBAIJ(Mat mat)
42: {
43: Mat_MPISBAIJ *aij = (Mat_MPISBAIJ *)mat->data;
47: MatRetrieveValues(aij->A);
48: MatRetrieveValues(aij->B);
49: return(0);
50: }
51: EXTERN_C_END
54: #define MatSetValues_SeqSBAIJ_A_Private(row,col,value,addv) \
55: { \
56: \
57: brow = row/bs; \
58: rp = aj + ai[brow]; ap = aa + bs2*ai[brow]; \
59: rmax = aimax[brow]; nrow = ailen[brow]; \
60: bcol = col/bs; \
61: ridx = row % bs; cidx = col % bs; \
62: low = 0; high = nrow; \
63: while (high-low > 3) { \
64: t = (low+high)/2; \
65: if (rp[t] > bcol) high = t; \
66: else low = t; \
67: } \
68: for (_i=low; _i<high; _i++) { \
69: if (rp[_i] > bcol) break; \
70: if (rp[_i] == bcol) { \
71: bap = ap + bs2*_i + bs*cidx + ridx; \
72: if (addv == ADD_VALUES) *bap += value; \
73: else *bap = value; \
74: goto a_noinsert; \
75: } \
76: } \
77: if (a->nonew == 1) goto a_noinsert; \
78: if (a->nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) into matrix", row, col); \
79: MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,aimax,a->nonew,MatScalar); \
80: N = nrow++ - 1; \
81: /* shift up all the later entries in this row */ \
82: for (ii=N; ii>=_i; ii--) { \
83: rp[ii+1] = rp[ii]; \
84: PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar)); \
85: } \
86: if (N>=_i) { PetscMemzero(ap+bs2*_i,bs2*sizeof(MatScalar)); } \
87: rp[_i] = bcol; \
88: ap[bs2*_i + bs*cidx + ridx] = value; \
89: a_noinsert:; \
90: ailen[brow] = nrow; \
91: }
93: #define MatSetValues_SeqSBAIJ_B_Private(row,col,value,addv) \
94: { \
95: brow = row/bs; \
96: rp = bj + bi[brow]; ap = ba + bs2*bi[brow]; \
97: rmax = bimax[brow]; nrow = bilen[brow]; \
98: bcol = col/bs; \
99: ridx = row % bs; cidx = col % bs; \
100: low = 0; high = nrow; \
101: while (high-low > 3) { \
102: t = (low+high)/2; \
103: if (rp[t] > bcol) high = t; \
104: else low = t; \
105: } \
106: for (_i=low; _i<high; _i++) { \
107: if (rp[_i] > bcol) break; \
108: if (rp[_i] == bcol) { \
109: bap = ap + bs2*_i + bs*cidx + ridx; \
110: if (addv == ADD_VALUES) *bap += value; \
111: else *bap = value; \
112: goto b_noinsert; \
113: } \
114: } \
115: if (b->nonew == 1) goto b_noinsert; \
116: if (b->nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) into matrix", row, col); \
117: MatSeqXAIJReallocateAIJ(B,b->mbs,bs2,nrow,brow,bcol,rmax,ba,bi,bj,rp,ap,bimax,b->nonew,MatScalar); \
118: N = nrow++ - 1; \
119: /* shift up all the later entries in this row */ \
120: for (ii=N; ii>=_i; ii--) { \
121: rp[ii+1] = rp[ii]; \
122: PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar)); \
123: } \
124: if (N>=_i) { PetscMemzero(ap+bs2*_i,bs2*sizeof(MatScalar));} \
125: rp[_i] = bcol; \
126: ap[bs2*_i + bs*cidx + ridx] = value; \
127: b_noinsert:; \
128: bilen[brow] = nrow; \
129: }
131: /* Only add/insert a(i,j) with i<=j (blocks).
132: Any a(i,j) with i>j input by user is ingored.
133: */
136: PetscErrorCode MatSetValues_MPISBAIJ(Mat mat,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode addv)
137: {
138: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
139: MatScalar value;
140: PetscBool roworiented = baij->roworiented;
142: PetscInt i,j,row,col;
143: PetscInt rstart_orig=mat->rmap->rstart;
144: PetscInt rend_orig=mat->rmap->rend,cstart_orig=mat->cmap->rstart;
145: PetscInt cend_orig=mat->cmap->rend,bs=mat->rmap->bs;
147: /* Some Variables required in the macro */
148: Mat A = baij->A;
149: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)(A)->data;
150: PetscInt *aimax=a->imax,*ai=a->i,*ailen=a->ilen,*aj=a->j;
151: MatScalar *aa=a->a;
153: Mat B = baij->B;
154: Mat_SeqBAIJ *b = (Mat_SeqBAIJ*)(B)->data;
155: PetscInt *bimax=b->imax,*bi=b->i,*bilen=b->ilen,*bj=b->j;
156: MatScalar *ba=b->a;
158: PetscInt *rp,ii,nrow,_i,rmax,N,brow,bcol;
159: PetscInt low,high,t,ridx,cidx,bs2=a->bs2;
160: MatScalar *ap,*bap;
162: /* for stash */
163: PetscInt n_loc, *in_loc = PETSC_NULL;
164: MatScalar *v_loc = PETSC_NULL;
168: if (!baij->donotstash){
169: if (n > baij->n_loc) {
170: PetscFree(baij->in_loc);
171: PetscFree(baij->v_loc);
172: PetscMalloc(n*sizeof(PetscInt),&baij->in_loc);
173: PetscMalloc(n*sizeof(MatScalar),&baij->v_loc);
174: baij->n_loc = n;
175: }
176: in_loc = baij->in_loc;
177: v_loc = baij->v_loc;
178: }
180: for (i=0; i<m; i++) {
181: if (im[i] < 0) continue;
182: #if defined(PETSC_USE_DEBUG)
183: if (im[i] >= mat->rmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",im[i],mat->rmap->N-1);
184: #endif
185: if (im[i] >= rstart_orig && im[i] < rend_orig) { /* this processor entry */
186: row = im[i] - rstart_orig; /* local row index */
187: for (j=0; j<n; j++) {
188: if (im[i]/bs > in[j]/bs){
189: if (a->ignore_ltriangular){
190: continue; /* ignore lower triangular blocks */
191: } else {
192: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Lower triangular value cannot be set for sbaij format. Ignoring these values, run with -mat_ignore_lower_triangular or call MatSetOption(mat,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE)");
193: }
194: }
195: if (in[j] >= cstart_orig && in[j] < cend_orig){ /* diag entry (A) */
196: col = in[j] - cstart_orig; /* local col index */
197: brow = row/bs; bcol = col/bs;
198: if (brow > bcol) continue; /* ignore lower triangular blocks of A */
199: if (roworiented) value = v[i*n+j]; else value = v[i+j*m];
200: MatSetValues_SeqSBAIJ_A_Private(row,col,value,addv);
201: /* MatSetValues_SeqBAIJ(baij->A,1,&row,1,&col,&value,addv); */
202: } else if (in[j] < 0) continue;
203: #if defined(PETSC_USE_DEBUG)
204: else if (in[j] >= mat->cmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[j],mat->cmap->N-1);
205: #endif
206: else { /* off-diag entry (B) */
207: if (mat->was_assembled) {
208: if (!baij->colmap) {
209: CreateColmap_MPIBAIJ_Private(mat);
210: }
211: #if defined (PETSC_USE_CTABLE)
212: PetscTableFind(baij->colmap,in[j]/bs + 1,&col);
213: col = col - 1;
214: #else
215: col = baij->colmap[in[j]/bs] - 1;
216: #endif
217: if (col < 0 && !((Mat_SeqSBAIJ*)(baij->A->data))->nonew) {
218: DisAssemble_MPISBAIJ(mat);
219: col = in[j];
220: /* Reinitialize the variables required by MatSetValues_SeqBAIJ_B_Private() */
221: B = baij->B;
222: b = (Mat_SeqBAIJ*)(B)->data;
223: bimax=b->imax;bi=b->i;bilen=b->ilen;bj=b->j;
224: ba=b->a;
225: } else col += in[j]%bs;
226: } else col = in[j];
227: if (roworiented) value = v[i*n+j]; else value = v[i+j*m];
228: MatSetValues_SeqSBAIJ_B_Private(row,col,value,addv);
229: /* MatSetValues_SeqBAIJ(baij->B,1,&row,1,&col,&value,addv); */
230: }
231: }
232: } else { /* off processor entry */
233: if (mat->nooffprocentries) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Setting off process row %D even though MatSetOption(,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE) was set",im[i]);
234: if (!baij->donotstash) {
235: mat->assembled = PETSC_FALSE;
236: n_loc = 0;
237: for (j=0; j<n; j++){
238: if (im[i]/bs > in[j]/bs) continue; /* ignore lower triangular blocks */
239: in_loc[n_loc] = in[j];
240: if (roworiented) {
241: v_loc[n_loc] = v[i*n+j];
242: } else {
243: v_loc[n_loc] = v[j*m+i];
244: }
245: n_loc++;
246: }
247: MatStashValuesRow_Private(&mat->stash,im[i],n_loc,in_loc,v_loc,PETSC_FALSE);
248: }
249: }
250: }
251: return(0);
252: }
256: PetscErrorCode MatSetValuesBlocked_MPISBAIJ(Mat mat,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const MatScalar v[],InsertMode addv)
257: {
258: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
259: const MatScalar *value;
260: MatScalar *barray=baij->barray;
261: PetscBool roworiented = baij->roworiented,ignore_ltriangular = ((Mat_SeqSBAIJ*)baij->A->data)->ignore_ltriangular;
262: PetscErrorCode ierr;
263: PetscInt i,j,ii,jj,row,col,rstart=baij->rstartbs;
264: PetscInt rend=baij->rendbs,cstart=baij->rstartbs,stepval;
265: PetscInt cend=baij->rendbs,bs=mat->rmap->bs,bs2=baij->bs2;
268: if(!barray) {
269: PetscMalloc(bs2*sizeof(MatScalar),&barray);
270: baij->barray = barray;
271: }
273: if (roworiented) {
274: stepval = (n-1)*bs;
275: } else {
276: stepval = (m-1)*bs;
277: }
278: for (i=0; i<m; i++) {
279: if (im[i] < 0) continue;
280: #if defined(PETSC_USE_DEBUG)
281: if (im[i] >= baij->Mbs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large, row %D max %D",im[i],baij->Mbs-1);
282: #endif
283: if (im[i] >= rstart && im[i] < rend) {
284: row = im[i] - rstart;
285: for (j=0; j<n; j++) {
286: if (im[i] > in[j]) {
287: if (ignore_ltriangular) continue; /* ignore lower triangular blocks */
288: else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Lower triangular value cannot be set for sbaij format. Ignoring these values, run with -mat_ignore_lower_triangular or call MatSetOption(mat,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE)");
289: }
290: /* If NumCol = 1 then a copy is not required */
291: if ((roworiented) && (n == 1)) {
292: barray = (MatScalar*) v + i*bs2;
293: } else if((!roworiented) && (m == 1)) {
294: barray = (MatScalar*) v + j*bs2;
295: } else { /* Here a copy is required */
296: if (roworiented) {
297: value = v + i*(stepval+bs)*bs + j*bs;
298: } else {
299: value = v + j*(stepval+bs)*bs + i*bs;
300: }
301: for (ii=0; ii<bs; ii++,value+=stepval) {
302: for (jj=0; jj<bs; jj++) {
303: *barray++ = *value++;
304: }
305: }
306: barray -=bs2;
307: }
308:
309: if (in[j] >= cstart && in[j] < cend){
310: col = in[j] - cstart;
311: MatSetValuesBlocked_SeqSBAIJ(baij->A,1,&row,1,&col,barray,addv);
312: }
313: else if (in[j] < 0) continue;
314: #if defined(PETSC_USE_DEBUG)
315: else if (in[j] >= baij->Nbs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large, col %D max %D",in[j],baij->Nbs-1);
316: #endif
317: else {
318: if (mat->was_assembled) {
319: if (!baij->colmap) {
320: CreateColmap_MPIBAIJ_Private(mat);
321: }
323: #if defined(PETSC_USE_DEBUG)
324: #if defined (PETSC_USE_CTABLE)
325: { PetscInt data;
326: PetscTableFind(baij->colmap,in[j]+1,&data);
327: if ((data - 1) % bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect colmap");
328: }
329: #else
330: if ((baij->colmap[in[j]] - 1) % bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect colmap");
331: #endif
332: #endif
333: #if defined (PETSC_USE_CTABLE)
334: PetscTableFind(baij->colmap,in[j]+1,&col);
335: col = (col - 1)/bs;
336: #else
337: col = (baij->colmap[in[j]] - 1)/bs;
338: #endif
339: if (col < 0 && !((Mat_SeqBAIJ*)(baij->A->data))->nonew) {
340: DisAssemble_MPISBAIJ(mat);
341: col = in[j];
342: }
343: }
344: else col = in[j];
345: MatSetValuesBlocked_SeqBAIJ(baij->B,1,&row,1,&col,barray,addv);
346: }
347: }
348: } else {
349: if (mat->nooffprocentries) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Setting off process row %D even though MatSetOption(,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE) was set",im[i]);
350: if (!baij->donotstash) {
351: if (roworiented) {
352: MatStashValuesRowBlocked_Private(&mat->bstash,im[i],n,in,v,m,n,i);
353: } else {
354: MatStashValuesColBlocked_Private(&mat->bstash,im[i],n,in,v,m,n,i);
355: }
356: }
357: }
358: }
359: return(0);
360: }
364: PetscErrorCode MatGetValues_MPISBAIJ(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[])
365: {
366: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
368: PetscInt bs=mat->rmap->bs,i,j,bsrstart = mat->rmap->rstart,bsrend = mat->rmap->rend;
369: PetscInt bscstart = mat->cmap->rstart,bscend = mat->cmap->rend,row,col,data;
372: for (i=0; i<m; i++) {
373: if (idxm[i] < 0) continue; /* SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",idxm[i]); */
374: if (idxm[i] >= mat->rmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm[i],mat->rmap->N-1);
375: if (idxm[i] >= bsrstart && idxm[i] < bsrend) {
376: row = idxm[i] - bsrstart;
377: for (j=0; j<n; j++) {
378: if (idxn[j] < 0) continue; /* SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column %D",idxn[j]); */
379: if (idxn[j] >= mat->cmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",idxn[j],mat->cmap->N-1);
380: if (idxn[j] >= bscstart && idxn[j] < bscend){
381: col = idxn[j] - bscstart;
382: MatGetValues_SeqSBAIJ(baij->A,1,&row,1,&col,v+i*n+j);
383: } else {
384: if (!baij->colmap) {
385: CreateColmap_MPIBAIJ_Private(mat);
386: }
387: #if defined (PETSC_USE_CTABLE)
388: PetscTableFind(baij->colmap,idxn[j]/bs+1,&data);
389: data --;
390: #else
391: data = baij->colmap[idxn[j]/bs]-1;
392: #endif
393: if((data < 0) || (baij->garray[data/bs] != idxn[j]/bs)) *(v+i*n+j) = 0.0;
394: else {
395: col = data + idxn[j]%bs;
396: MatGetValues_SeqBAIJ(baij->B,1,&row,1,&col,v+i*n+j);
397: }
398: }
399: }
400: } else {
401: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only local values currently supported");
402: }
403: }
404: return(0);
405: }
409: PetscErrorCode MatNorm_MPISBAIJ(Mat mat,NormType type,PetscReal *norm)
410: {
411: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
413: PetscReal sum[2],*lnorm2;
416: if (baij->size == 1) {
417: MatNorm(baij->A,type,norm);
418: } else {
419: if (type == NORM_FROBENIUS) {
420: PetscMalloc(2*sizeof(PetscReal),&lnorm2);
421: MatNorm(baij->A,type,lnorm2);
422: *lnorm2 = (*lnorm2)*(*lnorm2); lnorm2++; /* squar power of norm(A) */
423: MatNorm(baij->B,type,lnorm2);
424: *lnorm2 = (*lnorm2)*(*lnorm2); lnorm2--; /* squar power of norm(B) */
425: MPI_Allreduce(lnorm2,&sum,2,MPIU_REAL,MPIU_SUM,((PetscObject)mat)->comm);
426: *norm = PetscSqrtReal(sum[0] + 2*sum[1]);
427: PetscFree(lnorm2);
428: } else if (type == NORM_INFINITY || type == NORM_1) { /* max row/column sum */
429: Mat_SeqSBAIJ *amat=(Mat_SeqSBAIJ*)baij->A->data;
430: Mat_SeqBAIJ *bmat=(Mat_SeqBAIJ*)baij->B->data;
431: PetscReal *rsum,*rsum2,vabs;
432: PetscInt *jj,*garray=baij->garray,rstart=baij->rstartbs,nz;
433: PetscInt brow,bcol,col,bs=baij->A->rmap->bs,row,grow,gcol,mbs=amat->mbs;
434: MatScalar *v;
436: PetscMalloc2(mat->cmap->N,PetscReal,&rsum,mat->cmap->N,PetscReal,&rsum2);
437: PetscMemzero(rsum,mat->cmap->N*sizeof(PetscReal));
438: /* Amat */
439: v = amat->a; jj = amat->j;
440: for (brow=0; brow<mbs; brow++) {
441: grow = bs*(rstart + brow);
442: nz = amat->i[brow+1] - amat->i[brow];
443: for (bcol=0; bcol<nz; bcol++){
444: gcol = bs*(rstart + *jj); jj++;
445: for (col=0; col<bs; col++){
446: for (row=0; row<bs; row++){
447: vabs = PetscAbsScalar(*v); v++;
448: rsum[gcol+col] += vabs;
449: /* non-diagonal block */
450: if (bcol > 0 && vabs > 0.0) rsum[grow+row] += vabs;
451: }
452: }
453: }
454: }
455: /* Bmat */
456: v = bmat->a; jj = bmat->j;
457: for (brow=0; brow<mbs; brow++) {
458: grow = bs*(rstart + brow);
459: nz = bmat->i[brow+1] - bmat->i[brow];
460: for (bcol=0; bcol<nz; bcol++){
461: gcol = bs*garray[*jj]; jj++;
462: for (col=0; col<bs; col++){
463: for (row=0; row<bs; row++){
464: vabs = PetscAbsScalar(*v); v++;
465: rsum[gcol+col] += vabs;
466: rsum[grow+row] += vabs;
467: }
468: }
469: }
470: }
471: MPI_Allreduce(rsum,rsum2,mat->cmap->N,MPIU_REAL,MPIU_SUM,((PetscObject)mat)->comm);
472: *norm = 0.0;
473: for (col=0; col<mat->cmap->N; col++) {
474: if (rsum2[col] > *norm) *norm = rsum2[col];
475: }
476: PetscFree2(rsum,rsum2);
477: } else {
478: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for this norm yet");
479: }
480: }
481: return(0);
482: }
486: PetscErrorCode MatAssemblyBegin_MPISBAIJ(Mat mat,MatAssemblyType mode)
487: {
488: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
490: PetscInt nstash,reallocs;
491: InsertMode addv;
494: if (baij->donotstash || mat->nooffprocentries) {
495: return(0);
496: }
498: /* make sure all processors are either in INSERTMODE or ADDMODE */
499: MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,((PetscObject)mat)->comm);
500: if (addv == (ADD_VALUES|INSERT_VALUES)) SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_ARG_WRONGSTATE,"Some processors inserted others added");
501: mat->insertmode = addv; /* in case this processor had no cache */
503: MatStashScatterBegin_Private(mat,&mat->stash,mat->rmap->range);
504: MatStashScatterBegin_Private(mat,&mat->bstash,baij->rangebs);
505: MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);
506: PetscInfo2(mat,"Stash has %D entries,uses %D mallocs.\n",nstash,reallocs);
507: MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);
508: PetscInfo2(mat,"Block-Stash has %D entries, uses %D mallocs.\n",nstash,reallocs);
509: return(0);
510: }
514: PetscErrorCode MatAssemblyEnd_MPISBAIJ(Mat mat,MatAssemblyType mode)
515: {
516: Mat_MPISBAIJ *baij=(Mat_MPISBAIJ*)mat->data;
517: Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)baij->A->data;
519: PetscInt i,j,rstart,ncols,flg,bs2=baij->bs2;
520: PetscInt *row,*col;
521: PetscBool other_disassembled;
522: PetscMPIInt n;
523: PetscBool r1,r2,r3;
524: MatScalar *val;
525: InsertMode addv = mat->insertmode;
527: /* do not use 'b=(Mat_SeqBAIJ*)baij->B->data' as B can be reset in disassembly */
530: if (!baij->donotstash && !mat->nooffprocentries) {
531: while (1) {
532: MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);
533: if (!flg) break;
535: for (i=0; i<n;) {
536: /* Now identify the consecutive vals belonging to the same row */
537: for (j=i,rstart=row[j]; j<n; j++) { if (row[j] != rstart) break; }
538: if (j < n) ncols = j-i;
539: else ncols = n-i;
540: /* Now assemble all these values with a single function call */
541: MatSetValues_MPISBAIJ(mat,1,row+i,ncols,col+i,val+i,addv);
542: i = j;
543: }
544: }
545: MatStashScatterEnd_Private(&mat->stash);
546: /* Now process the block-stash. Since the values are stashed column-oriented,
547: set the roworiented flag to column oriented, and after MatSetValues()
548: restore the original flags */
549: r1 = baij->roworiented;
550: r2 = a->roworiented;
551: r3 = ((Mat_SeqBAIJ*)baij->B->data)->roworiented;
552: baij->roworiented = PETSC_FALSE;
553: a->roworiented = PETSC_FALSE;
554: ((Mat_SeqBAIJ*)baij->B->data)->roworiented = PETSC_FALSE; /* b->roworinted */
555: while (1) {
556: MatStashScatterGetMesg_Private(&mat->bstash,&n,&row,&col,&val,&flg);
557: if (!flg) break;
558:
559: for (i=0; i<n;) {
560: /* Now identify the consecutive vals belonging to the same row */
561: for (j=i,rstart=row[j]; j<n; j++) { if (row[j] != rstart) break; }
562: if (j < n) ncols = j-i;
563: else ncols = n-i;
564: MatSetValuesBlocked_MPISBAIJ(mat,1,row+i,ncols,col+i,val+i*bs2,addv);
565: i = j;
566: }
567: }
568: MatStashScatterEnd_Private(&mat->bstash);
569: baij->roworiented = r1;
570: a->roworiented = r2;
571: ((Mat_SeqBAIJ*)baij->B->data)->roworiented = r3; /* b->roworinted */
572: }
574: MatAssemblyBegin(baij->A,mode);
575: MatAssemblyEnd(baij->A,mode);
577: /* determine if any processor has disassembled, if so we must
578: also disassemble ourselfs, in order that we may reassemble. */
579: /*
580: if nonzero structure of submatrix B cannot change then we know that
581: no processor disassembled thus we can skip this stuff
582: */
583: if (!((Mat_SeqBAIJ*)baij->B->data)->nonew) {
584: MPI_Allreduce(&mat->was_assembled,&other_disassembled,1,MPI_INT,MPI_PROD,((PetscObject)mat)->comm);
585: if (mat->was_assembled && !other_disassembled) {
586: DisAssemble_MPISBAIJ(mat);
587: }
588: }
590: if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
591: MatSetUpMultiply_MPISBAIJ(mat); /* setup Mvctx and sMvctx */
592: }
593: MatSetOption(baij->B,MAT_CHECK_COMPRESSED_ROW,PETSC_TRUE);
594: MatAssemblyBegin(baij->B,mode);
595: MatAssemblyEnd(baij->B,mode);
596:
597: PetscFree2(baij->rowvalues,baij->rowindices);
598: baij->rowvalues = 0;
600: return(0);
601: }
603: extern PetscErrorCode MatSetValues_MPIBAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
606: static PetscErrorCode MatView_MPISBAIJ_ASCIIorDraworSocket(Mat mat,PetscViewer viewer)
607: {
608: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
609: PetscErrorCode ierr;
610: PetscInt bs = mat->rmap->bs;
611: PetscMPIInt size = baij->size,rank = baij->rank;
612: PetscBool iascii,isdraw;
613: PetscViewer sviewer;
614: PetscViewerFormat format;
617: PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
618: PetscTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
619: if (iascii) {
620: PetscViewerGetFormat(viewer,&format);
621: if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
622: MatInfo info;
623: MPI_Comm_rank(((PetscObject)mat)->comm,&rank);
624: MatGetInfo(mat,MAT_LOCAL,&info);
625: PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);
626: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Local rows %D nz %D nz alloced %D bs %D mem %D\n",rank,mat->rmap->n,(PetscInt)info.nz_used,(PetscInt)info.nz_allocated,mat->rmap->bs,(PetscInt)info.memory);
627: MatGetInfo(baij->A,MAT_LOCAL,&info);
628: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] on-diagonal part: nz %D \n",rank,(PetscInt)info.nz_used);
629: MatGetInfo(baij->B,MAT_LOCAL,&info);
630: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] off-diagonal part: nz %D \n",rank,(PetscInt)info.nz_used);
631: PetscViewerFlush(viewer);
632: PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);
633: PetscViewerASCIIPrintf(viewer,"Information on VecScatter used in matrix-vector product: \n");
634: VecScatterView(baij->Mvctx,viewer);
635: return(0);
636: } else if (format == PETSC_VIEWER_ASCII_INFO) {
637: PetscViewerASCIIPrintf(viewer," block size is %D\n",bs);
638: return(0);
639: } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
640: return(0);
641: }
642: }
644: if (isdraw) {
645: PetscDraw draw;
646: PetscBool isnull;
647: PetscViewerDrawGetDraw(viewer,0,&draw);
648: PetscDrawIsNull(draw,&isnull); if (isnull) return(0);
649: }
651: if (size == 1) {
652: PetscObjectSetName((PetscObject)baij->A,((PetscObject)mat)->name);
653: MatView(baij->A,viewer);
654: } else {
655: /* assemble the entire matrix onto first processor. */
656: Mat A;
657: Mat_SeqSBAIJ *Aloc;
658: Mat_SeqBAIJ *Bloc;
659: PetscInt M = mat->rmap->N,N = mat->cmap->N,*ai,*aj,col,i,j,k,*rvals,mbs = baij->mbs;
660: MatScalar *a;
662: /* Should this be the same type as mat? */
663: MatCreate(((PetscObject)mat)->comm,&A);
664: if (!rank) {
665: MatSetSizes(A,M,N,M,N);
666: } else {
667: MatSetSizes(A,0,0,M,N);
668: }
669: MatSetType(A,MATMPISBAIJ);
670: MatMPISBAIJSetPreallocation(A,mat->rmap->bs,0,PETSC_NULL,0,PETSC_NULL);
671: MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);
672: PetscLogObjectParent(mat,A);
674: /* copy over the A part */
675: Aloc = (Mat_SeqSBAIJ*)baij->A->data;
676: ai = Aloc->i; aj = Aloc->j; a = Aloc->a;
677: PetscMalloc(bs*sizeof(PetscInt),&rvals);
679: for (i=0; i<mbs; i++) {
680: rvals[0] = bs*(baij->rstartbs + i);
681: for (j=1; j<bs; j++) { rvals[j] = rvals[j-1] + 1; }
682: for (j=ai[i]; j<ai[i+1]; j++) {
683: col = (baij->cstartbs+aj[j])*bs;
684: for (k=0; k<bs; k++) {
685: MatSetValues_MPISBAIJ(A,bs,rvals,1,&col,a,INSERT_VALUES);
686: col++; a += bs;
687: }
688: }
689: }
690: /* copy over the B part */
691: Bloc = (Mat_SeqBAIJ*)baij->B->data;
692: ai = Bloc->i; aj = Bloc->j; a = Bloc->a;
693: for (i=0; i<mbs; i++) {
694:
695: rvals[0] = bs*(baij->rstartbs + i);
696: for (j=1; j<bs; j++) { rvals[j] = rvals[j-1] + 1; }
697: for (j=ai[i]; j<ai[i+1]; j++) {
698: col = baij->garray[aj[j]]*bs;
699: for (k=0; k<bs; k++) {
700: MatSetValues_MPIBAIJ(A,bs,rvals,1,&col,a,INSERT_VALUES);
701: col++; a += bs;
702: }
703: }
704: }
705: PetscFree(rvals);
706: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
707: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
708: /*
709: Everyone has to call to draw the matrix since the graphics waits are
710: synchronized across all processors that share the PetscDraw object
711: */
712: PetscViewerGetSingleton(viewer,&sviewer);
713: if (!rank) {
714: PetscObjectSetName((PetscObject)((Mat_MPISBAIJ*)(A->data))->A,((PetscObject)mat)->name);
715: /* Set the type name to MATMPISBAIJ so that the correct type can be printed out by PetscObjectPrintClassNamePrefixType() in MatView_SeqSBAIJ_ASCII()*/
716: PetscStrcpy(((PetscObject)((Mat_MPISBAIJ*)(A->data))->A)->type_name,MATMPISBAIJ);
717: MatView(((Mat_MPISBAIJ*)(A->data))->A,sviewer);
718: }
719: PetscViewerRestoreSingleton(viewer,&sviewer);
720: MatDestroy(&A);
721: }
722: return(0);
723: }
727: PetscErrorCode MatView_MPISBAIJ(Mat mat,PetscViewer viewer)
728: {
730: PetscBool iascii,isdraw,issocket,isbinary;
733: PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
734: PetscTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
735: PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSOCKET,&issocket);
736: PetscTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
737: if (iascii || isdraw || issocket || isbinary) {
738: MatView_MPISBAIJ_ASCIIorDraworSocket(mat,viewer);
739: } else {
740: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported by MPISBAIJ matrices",((PetscObject)viewer)->type_name);
741: }
742: return(0);
743: }
747: PetscErrorCode MatDestroy_MPISBAIJ(Mat mat)
748: {
749: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
753: #if defined(PETSC_USE_LOG)
754: PetscLogObjectState((PetscObject)mat,"Rows=%D,Cols=%D",mat->rmap->N,mat->cmap->N);
755: #endif
756: MatStashDestroy_Private(&mat->stash);
757: MatStashDestroy_Private(&mat->bstash);
758: MatDestroy(&baij->A);
759: MatDestroy(&baij->B);
760: #if defined (PETSC_USE_CTABLE)
761: PetscTableDestroy(&baij->colmap);
762: #else
763: PetscFree(baij->colmap);
764: #endif
765: PetscFree(baij->garray);
766: VecDestroy(&baij->lvec);
767: VecScatterDestroy(&baij->Mvctx);
768: VecDestroy(&baij->slvec0);
769: VecDestroy(&baij->slvec0b);
770: VecDestroy(&baij->slvec1);
771: VecDestroy(&baij->slvec1a);
772: VecDestroy(&baij->slvec1b);
773: VecScatterDestroy(&baij->sMvctx);
774: PetscFree2(baij->rowvalues,baij->rowindices);
775: PetscFree(baij->barray);
776: PetscFree(baij->hd);
777: VecDestroy(&baij->diag);
778: VecDestroy(&baij->bb1);
779: VecDestroy(&baij->xx1);
780: #if defined(PETSC_USE_REAL_MAT_SINGLE)
781: PetscFree(baij->setvaluescopy);
782: #endif
783: PetscFree(baij->in_loc);
784: PetscFree(baij->v_loc);
785: PetscFree(baij->rangebs);
786: PetscFree(mat->data);
788: PetscObjectChangeTypeName((PetscObject)mat,0);
789: PetscObjectComposeFunction((PetscObject)mat,"MatStoreValues_C","",PETSC_NULL);
790: PetscObjectComposeFunction((PetscObject)mat,"MatRetrieveValues_C","",PETSC_NULL);
791: PetscObjectComposeFunction((PetscObject)mat,"MatGetDiagonalBlock_C","",PETSC_NULL);
792: PetscObjectComposeFunction((PetscObject)mat,"MatMPISBAIJSetPreallocation_C","",PETSC_NULL);
793: PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpisbaij_mpisbstrm_C","",PETSC_NULL);
794: return(0);
795: }
799: PetscErrorCode MatMult_MPISBAIJ_Hermitian(Mat A,Vec xx,Vec yy)
800: {
801: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
803: PetscInt nt,mbs=a->mbs,bs=A->rmap->bs;
804: PetscScalar *x,*from;
805:
807: VecGetLocalSize(xx,&nt);
808: if (nt != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible partition of A and xx");
810: /* diagonal part */
811: (*a->A->ops->mult)(a->A,xx,a->slvec1a);
812: VecSet(a->slvec1b,0.0);
814: /* subdiagonal part */
815: (*a->B->ops->multhermitiantranspose)(a->B,xx,a->slvec0b);
817: /* copy x into the vec slvec0 */
818: VecGetArray(a->slvec0,&from);
819: VecGetArray(xx,&x);
821: PetscMemcpy(from,x,bs*mbs*sizeof(MatScalar));
822: VecRestoreArray(a->slvec0,&from);
823: VecRestoreArray(xx,&x);
824:
825: VecScatterBegin(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
826: VecScatterEnd(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
827: /* supperdiagonal part */
828: (*a->B->ops->multadd)(a->B,a->slvec1b,a->slvec1a,yy);
829: return(0);
830: }
834: PetscErrorCode MatMult_MPISBAIJ(Mat A,Vec xx,Vec yy)
835: {
836: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
838: PetscInt nt,mbs=a->mbs,bs=A->rmap->bs;
839: PetscScalar *x,*from;
840:
842: VecGetLocalSize(xx,&nt);
843: if (nt != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible partition of A and xx");
845: /* diagonal part */
846: (*a->A->ops->mult)(a->A,xx,a->slvec1a);
847: VecSet(a->slvec1b,0.0);
849: /* subdiagonal part */
850: (*a->B->ops->multtranspose)(a->B,xx,a->slvec0b);
852: /* copy x into the vec slvec0 */
853: VecGetArray(a->slvec0,&from);
854: VecGetArray(xx,&x);
856: PetscMemcpy(from,x,bs*mbs*sizeof(MatScalar));
857: VecRestoreArray(a->slvec0,&from);
858: VecRestoreArray(xx,&x);
859:
860: VecScatterBegin(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
861: VecScatterEnd(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
862: /* supperdiagonal part */
863: (*a->B->ops->multadd)(a->B,a->slvec1b,a->slvec1a,yy);
864: return(0);
865: }
869: PetscErrorCode MatMult_MPISBAIJ_2comm(Mat A,Vec xx,Vec yy)
870: {
871: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
873: PetscInt nt;
876: VecGetLocalSize(xx,&nt);
877: if (nt != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible partition of A and xx");
879: VecGetLocalSize(yy,&nt);
880: if (nt != A->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible parition of A and yy");
882: VecScatterBegin(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);
883: /* do diagonal part */
884: (*a->A->ops->mult)(a->A,xx,yy);
885: /* do supperdiagonal part */
886: VecScatterEnd(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);
887: (*a->B->ops->multadd)(a->B,a->lvec,yy,yy);
888: /* do subdiagonal part */
889: (*a->B->ops->multtranspose)(a->B,xx,a->lvec);
890: VecScatterBegin(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);
891: VecScatterEnd(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);
893: return(0);
894: }
898: PetscErrorCode MatMultAdd_MPISBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
899: {
900: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
902: PetscInt mbs=a->mbs,bs=A->rmap->bs;
903: PetscScalar *x,*from,zero=0.0;
904:
906: /*
907: PetscSynchronizedPrintf(((PetscObject)A)->comm," MatMultAdd is called ...\n");
908: PetscSynchronizedFlush(((PetscObject)A)->comm);
909: */
910: /* diagonal part */
911: (*a->A->ops->multadd)(a->A,xx,yy,a->slvec1a);
912: VecSet(a->slvec1b,zero);
914: /* subdiagonal part */
915: (*a->B->ops->multtranspose)(a->B,xx,a->slvec0b);
917: /* copy x into the vec slvec0 */
918: VecGetArray(a->slvec0,&from);
919: VecGetArray(xx,&x);
920: PetscMemcpy(from,x,bs*mbs*sizeof(MatScalar));
921: VecRestoreArray(a->slvec0,&from);
922:
923: VecScatterBegin(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
924: VecRestoreArray(xx,&x);
925: VecScatterEnd(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
926:
927: /* supperdiagonal part */
928: (*a->B->ops->multadd)(a->B,a->slvec1b,a->slvec1a,zz);
929:
930: return(0);
931: }
935: PetscErrorCode MatMultAdd_MPISBAIJ_2comm(Mat A,Vec xx,Vec yy,Vec zz)
936: {
937: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
941: VecScatterBegin(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);
942: /* do diagonal part */
943: (*a->A->ops->multadd)(a->A,xx,yy,zz);
944: /* do supperdiagonal part */
945: VecScatterEnd(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);
946: (*a->B->ops->multadd)(a->B,a->lvec,zz,zz);
948: /* do subdiagonal part */
949: (*a->B->ops->multtranspose)(a->B,xx,a->lvec);
950: VecScatterBegin(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);
951: VecScatterEnd(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);
953: return(0);
954: }
956: /*
957: This only works correctly for square matrices where the subblock A->A is the
958: diagonal block
959: */
962: PetscErrorCode MatGetDiagonal_MPISBAIJ(Mat A,Vec v)
963: {
964: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
968: /* if (a->rmap->N != a->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Supports only square matrix where A->A is diag block"); */
969: MatGetDiagonal(a->A,v);
970: return(0);
971: }
975: PetscErrorCode MatScale_MPISBAIJ(Mat A,PetscScalar aa)
976: {
977: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
981: MatScale(a->A,aa);
982: MatScale(a->B,aa);
983: return(0);
984: }
988: PetscErrorCode MatGetRow_MPISBAIJ(Mat matin,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
989: {
990: Mat_MPISBAIJ *mat = (Mat_MPISBAIJ*)matin->data;
991: PetscScalar *vworkA,*vworkB,**pvA,**pvB,*v_p;
993: PetscInt bs = matin->rmap->bs,bs2 = mat->bs2,i,*cworkA,*cworkB,**pcA,**pcB;
994: PetscInt nztot,nzA,nzB,lrow,brstart = matin->rmap->rstart,brend = matin->rmap->rend;
995: PetscInt *cmap,*idx_p,cstart = mat->rstartbs;
998: if (mat->getrowactive) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Already active");
999: mat->getrowactive = PETSC_TRUE;
1001: if (!mat->rowvalues && (idx || v)) {
1002: /*
1003: allocate enough space to hold information from the longest row.
1004: */
1005: Mat_SeqSBAIJ *Aa = (Mat_SeqSBAIJ*)mat->A->data;
1006: Mat_SeqBAIJ *Ba = (Mat_SeqBAIJ*)mat->B->data;
1007: PetscInt max = 1,mbs = mat->mbs,tmp;
1008: for (i=0; i<mbs; i++) {
1009: tmp = Aa->i[i+1] - Aa->i[i] + Ba->i[i+1] - Ba->i[i]; /* row length */
1010: if (max < tmp) { max = tmp; }
1011: }
1012: PetscMalloc2(max*bs2,PetscScalar,&mat->rowvalues,max*bs2,PetscInt,&mat->rowindices);
1013: }
1014:
1015: if (row < brstart || row >= brend) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only local rows");
1016: lrow = row - brstart; /* local row index */
1018: pvA = &vworkA; pcA = &cworkA; pvB = &vworkB; pcB = &cworkB;
1019: if (!v) {pvA = 0; pvB = 0;}
1020: if (!idx) {pcA = 0; if (!v) pcB = 0;}
1021: (*mat->A->ops->getrow)(mat->A,lrow,&nzA,pcA,pvA);
1022: (*mat->B->ops->getrow)(mat->B,lrow,&nzB,pcB,pvB);
1023: nztot = nzA + nzB;
1025: cmap = mat->garray;
1026: if (v || idx) {
1027: if (nztot) {
1028: /* Sort by increasing column numbers, assuming A and B already sorted */
1029: PetscInt imark = -1;
1030: if (v) {
1031: *v = v_p = mat->rowvalues;
1032: for (i=0; i<nzB; i++) {
1033: if (cmap[cworkB[i]/bs] < cstart) v_p[i] = vworkB[i];
1034: else break;
1035: }
1036: imark = i;
1037: for (i=0; i<nzA; i++) v_p[imark+i] = vworkA[i];
1038: for (i=imark; i<nzB; i++) v_p[nzA+i] = vworkB[i];
1039: }
1040: if (idx) {
1041: *idx = idx_p = mat->rowindices;
1042: if (imark > -1) {
1043: for (i=0; i<imark; i++) {
1044: idx_p[i] = cmap[cworkB[i]/bs]*bs + cworkB[i]%bs;
1045: }
1046: } else {
1047: for (i=0; i<nzB; i++) {
1048: if (cmap[cworkB[i]/bs] < cstart)
1049: idx_p[i] = cmap[cworkB[i]/bs]*bs + cworkB[i]%bs ;
1050: else break;
1051: }
1052: imark = i;
1053: }
1054: for (i=0; i<nzA; i++) idx_p[imark+i] = cstart*bs + cworkA[i];
1055: for (i=imark; i<nzB; i++) idx_p[nzA+i] = cmap[cworkB[i]/bs]*bs + cworkB[i]%bs ;
1056: }
1057: } else {
1058: if (idx) *idx = 0;
1059: if (v) *v = 0;
1060: }
1061: }
1062: *nz = nztot;
1063: (*mat->A->ops->restorerow)(mat->A,lrow,&nzA,pcA,pvA);
1064: (*mat->B->ops->restorerow)(mat->B,lrow,&nzB,pcB,pvB);
1065: return(0);
1066: }
1070: PetscErrorCode MatRestoreRow_MPISBAIJ(Mat mat,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
1071: {
1072: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
1075: if (!baij->getrowactive) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"MatGetRow() must be called first");
1076: baij->getrowactive = PETSC_FALSE;
1077: return(0);
1078: }
1082: PetscErrorCode MatGetRowUpperTriangular_MPISBAIJ(Mat A)
1083: {
1084: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1085: Mat_SeqSBAIJ *aA = (Mat_SeqSBAIJ*)a->A->data;
1088: aA->getrow_utriangular = PETSC_TRUE;
1089: return(0);
1090: }
1093: PetscErrorCode MatRestoreRowUpperTriangular_MPISBAIJ(Mat A)
1094: {
1095: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1096: Mat_SeqSBAIJ *aA = (Mat_SeqSBAIJ*)a->A->data;
1099: aA->getrow_utriangular = PETSC_FALSE;
1100: return(0);
1101: }
1105: PetscErrorCode MatRealPart_MPISBAIJ(Mat A)
1106: {
1107: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1111: MatRealPart(a->A);
1112: MatRealPart(a->B);
1113: return(0);
1114: }
1118: PetscErrorCode MatImaginaryPart_MPISBAIJ(Mat A)
1119: {
1120: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1124: MatImaginaryPart(a->A);
1125: MatImaginaryPart(a->B);
1126: return(0);
1127: }
1131: PetscErrorCode MatZeroEntries_MPISBAIJ(Mat A)
1132: {
1133: Mat_MPISBAIJ *l = (Mat_MPISBAIJ*)A->data;
1137: MatZeroEntries(l->A);
1138: MatZeroEntries(l->B);
1139: return(0);
1140: }
1144: PetscErrorCode MatGetInfo_MPISBAIJ(Mat matin,MatInfoType flag,MatInfo *info)
1145: {
1146: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)matin->data;
1147: Mat A = a->A,B = a->B;
1149: PetscReal isend[5],irecv[5];
1152: info->block_size = (PetscReal)matin->rmap->bs;
1153: MatGetInfo(A,MAT_LOCAL,info);
1154: isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
1155: isend[3] = info->memory; isend[4] = info->mallocs;
1156: MatGetInfo(B,MAT_LOCAL,info);
1157: isend[0] += info->nz_used; isend[1] += info->nz_allocated; isend[2] += info->nz_unneeded;
1158: isend[3] += info->memory; isend[4] += info->mallocs;
1159: if (flag == MAT_LOCAL) {
1160: info->nz_used = isend[0];
1161: info->nz_allocated = isend[1];
1162: info->nz_unneeded = isend[2];
1163: info->memory = isend[3];
1164: info->mallocs = isend[4];
1165: } else if (flag == MAT_GLOBAL_MAX) {
1166: MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_MAX,((PetscObject)matin)->comm);
1167: info->nz_used = irecv[0];
1168: info->nz_allocated = irecv[1];
1169: info->nz_unneeded = irecv[2];
1170: info->memory = irecv[3];
1171: info->mallocs = irecv[4];
1172: } else if (flag == MAT_GLOBAL_SUM) {
1173: MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_SUM,((PetscObject)matin)->comm);
1174: info->nz_used = irecv[0];
1175: info->nz_allocated = irecv[1];
1176: info->nz_unneeded = irecv[2];
1177: info->memory = irecv[3];
1178: info->mallocs = irecv[4];
1179: } else {
1180: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown MatInfoType argument %d",(int)flag);
1181: }
1182: info->fill_ratio_given = 0; /* no parallel LU/ILU/Cholesky */
1183: info->fill_ratio_needed = 0;
1184: info->factor_mallocs = 0;
1185: return(0);
1186: }
1190: PetscErrorCode MatSetOption_MPISBAIJ(Mat A,MatOption op,PetscBool flg)
1191: {
1192: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1193: Mat_SeqSBAIJ *aA = (Mat_SeqSBAIJ*)a->A->data;
1197: switch (op) {
1198: case MAT_NEW_NONZERO_LOCATIONS:
1199: case MAT_NEW_NONZERO_ALLOCATION_ERR:
1200: case MAT_UNUSED_NONZERO_LOCATION_ERR:
1201: case MAT_KEEP_NONZERO_PATTERN:
1202: case MAT_NEW_NONZERO_LOCATION_ERR:
1203: MatSetOption(a->A,op,flg);
1204: MatSetOption(a->B,op,flg);
1205: break;
1206: case MAT_ROW_ORIENTED:
1207: a->roworiented = flg;
1208: MatSetOption(a->A,op,flg);
1209: MatSetOption(a->B,op,flg);
1210: break;
1211: case MAT_NEW_DIAGONALS:
1212: PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
1213: break;
1214: case MAT_IGNORE_OFF_PROC_ENTRIES:
1215: a->donotstash = flg;
1216: break;
1217: case MAT_USE_HASH_TABLE:
1218: a->ht_flag = flg;
1219: break;
1220: case MAT_HERMITIAN:
1221: if (!A->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatAssemblyEnd() first");
1222: MatSetOption(a->A,op,flg);
1223: A->ops->mult = MatMult_MPISBAIJ_Hermitian;
1224: break;
1225: case MAT_SPD:
1226: A->spd_set = PETSC_TRUE;
1227: A->spd = flg;
1228: if (flg) {
1229: A->symmetric = PETSC_TRUE;
1230: A->structurally_symmetric = PETSC_TRUE;
1231: A->symmetric_set = PETSC_TRUE;
1232: A->structurally_symmetric_set = PETSC_TRUE;
1233: }
1234: break;
1235: case MAT_SYMMETRIC:
1236: MatSetOption(a->A,op,flg);
1237: break;
1238: case MAT_STRUCTURALLY_SYMMETRIC:
1239: MatSetOption(a->A,op,flg);
1240: break;
1241: case MAT_SYMMETRY_ETERNAL:
1242: if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix must be symmetric");
1243: PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
1244: break;
1245: case MAT_IGNORE_LOWER_TRIANGULAR:
1246: aA->ignore_ltriangular = flg;
1247: break;
1248: case MAT_ERROR_LOWER_TRIANGULAR:
1249: aA->ignore_ltriangular = flg;
1250: break;
1251: case MAT_GETROW_UPPERTRIANGULAR:
1252: aA->getrow_utriangular = flg;
1253: break;
1254: default:
1255: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %d",op);
1256: }
1257: return(0);
1258: }
1262: PetscErrorCode MatTranspose_MPISBAIJ(Mat A,MatReuse reuse,Mat *B)
1263: {
1266: if (MAT_INITIAL_MATRIX || *B != A) {
1267: MatDuplicate(A,MAT_COPY_VALUES,B);
1268: }
1269: return(0);
1270: }
1274: PetscErrorCode MatDiagonalScale_MPISBAIJ(Mat mat,Vec ll,Vec rr)
1275: {
1276: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
1277: Mat a=baij->A, b=baij->B;
1279: PetscInt nv,m,n;
1280: PetscBool flg;
1283: if (ll != rr){
1284: VecEqual(ll,rr,&flg);
1285: if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"For symmetric format, left and right scaling vectors must be same\n");
1286: }
1287: if (!ll) return(0);
1289: MatGetLocalSize(mat,&m,&n);
1290: if (m != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"For symmetric format, local size %d %d must be same",m,n);
1291:
1292: VecGetLocalSize(rr,&nv);
1293: if (nv!=n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left and right vector non-conforming local size");
1295: VecScatterBegin(baij->Mvctx,rr,baij->lvec,INSERT_VALUES,SCATTER_FORWARD);
1296:
1297: /* left diagonalscale the off-diagonal part */
1298: (*b->ops->diagonalscale)(b,ll,PETSC_NULL);
1299:
1300: /* scale the diagonal part */
1301: (*a->ops->diagonalscale)(a,ll,rr);
1303: /* right diagonalscale the off-diagonal part */
1304: VecScatterEnd(baij->Mvctx,rr,baij->lvec,INSERT_VALUES,SCATTER_FORWARD);
1305: (*b->ops->diagonalscale)(b,PETSC_NULL,baij->lvec);
1306: return(0);
1307: }
1311: PetscErrorCode MatSetUnfactored_MPISBAIJ(Mat A)
1312: {
1313: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1317: MatSetUnfactored(a->A);
1318: return(0);
1319: }
1321: static PetscErrorCode MatDuplicate_MPISBAIJ(Mat,MatDuplicateOption,Mat *);
1325: PetscErrorCode MatEqual_MPISBAIJ(Mat A,Mat B,PetscBool *flag)
1326: {
1327: Mat_MPISBAIJ *matB = (Mat_MPISBAIJ*)B->data,*matA = (Mat_MPISBAIJ*)A->data;
1328: Mat a,b,c,d;
1329: PetscBool flg;
1333: a = matA->A; b = matA->B;
1334: c = matB->A; d = matB->B;
1336: MatEqual(a,c,&flg);
1337: if (flg) {
1338: MatEqual(b,d,&flg);
1339: }
1340: MPI_Allreduce(&flg,flag,1,MPI_INT,MPI_LAND,((PetscObject)A)->comm);
1341: return(0);
1342: }
1346: PetscErrorCode MatCopy_MPISBAIJ(Mat A,Mat B,MatStructure str)
1347: {
1349: Mat_MPISBAIJ *a = (Mat_MPISBAIJ *)A->data;
1350: Mat_MPISBAIJ *b = (Mat_MPISBAIJ *)B->data;
1353: /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
1354: if ((str != SAME_NONZERO_PATTERN) || (A->ops->copy != B->ops->copy)) {
1355: MatGetRowUpperTriangular(A);
1356: MatCopy_Basic(A,B,str);
1357: MatRestoreRowUpperTriangular(A);
1358: } else {
1359: MatCopy(a->A,b->A,str);
1360: MatCopy(a->B,b->B,str);
1361: }
1362: return(0);
1363: }
1367: PetscErrorCode MatSetUp_MPISBAIJ(Mat A)
1368: {
1372: MatMPISBAIJSetPreallocation(A,-PetscMax(A->rmap->bs,1),PETSC_DEFAULT,0,PETSC_DEFAULT,0);
1373: return(0);
1374: }
1378: PetscErrorCode MatAXPY_MPISBAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str)
1379: {
1381: Mat_MPISBAIJ *xx=(Mat_MPISBAIJ *)X->data,*yy=(Mat_MPISBAIJ *)Y->data;
1382: PetscBLASInt bnz,one=1;
1383: Mat_SeqSBAIJ *xa,*ya;
1384: Mat_SeqBAIJ *xb,*yb;
1387: if (str == SAME_NONZERO_PATTERN) {
1388: PetscScalar alpha = a;
1389: xa = (Mat_SeqSBAIJ *)xx->A->data;
1390: ya = (Mat_SeqSBAIJ *)yy->A->data;
1391: bnz = PetscBLASIntCast(xa->nz);
1392: BLASaxpy_(&bnz,&alpha,xa->a,&one,ya->a,&one);
1393: xb = (Mat_SeqBAIJ *)xx->B->data;
1394: yb = (Mat_SeqBAIJ *)yy->B->data;
1395: bnz = PetscBLASIntCast(xb->nz);
1396: BLASaxpy_(&bnz,&alpha,xb->a,&one,yb->a,&one);
1397: } else {
1398: MatGetRowUpperTriangular(X);
1399: MatAXPY_Basic(Y,a,X,str);
1400: MatRestoreRowUpperTriangular(X);
1401: }
1402: return(0);
1403: }
1407: PetscErrorCode MatSetBlockSize_MPISBAIJ(Mat A,PetscInt bs)
1408: {
1409: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1410: PetscInt rbs,cbs;
1411: PetscErrorCode ierr;
1414: MatSetBlockSize(a->A,bs);
1415: MatSetBlockSize(a->B,bs);
1416: PetscLayoutGetBlockSize(A->rmap,&rbs);
1417: PetscLayoutGetBlockSize(A->cmap,&cbs);
1418: if (rbs != bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Attempt to set block size %d with SBAIJ %d",bs,rbs);
1419: if (cbs != bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Attempt to set block size %d with SBAIJ %d",bs,cbs);
1420: return(0);
1421: }
1425: PetscErrorCode MatGetSubMatrices_MPISBAIJ(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1426: {
1428: PetscInt i;
1429: PetscBool flg;
1432: for (i=0; i<n; i++) {
1433: ISEqual(irow[i],icol[i],&flg);
1434: if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Can only get symmetric submatrix for MPISBAIJ matrices");
1435: }
1436: MatGetSubMatrices_MPIBAIJ(A,n,irow,icol,scall,B);
1437: return(0);
1438: }
1439:
1441: /* -------------------------------------------------------------------*/
1442: static struct _MatOps MatOps_Values = {
1443: MatSetValues_MPISBAIJ,
1444: MatGetRow_MPISBAIJ,
1445: MatRestoreRow_MPISBAIJ,
1446: MatMult_MPISBAIJ,
1447: /* 4*/ MatMultAdd_MPISBAIJ,
1448: MatMult_MPISBAIJ, /* transpose versions are same as non-transpose */
1449: MatMultAdd_MPISBAIJ,
1450: 0,
1451: 0,
1452: 0,
1453: /*10*/ 0,
1454: 0,
1455: 0,
1456: MatSOR_MPISBAIJ,
1457: MatTranspose_MPISBAIJ,
1458: /*15*/ MatGetInfo_MPISBAIJ,
1459: MatEqual_MPISBAIJ,
1460: MatGetDiagonal_MPISBAIJ,
1461: MatDiagonalScale_MPISBAIJ,
1462: MatNorm_MPISBAIJ,
1463: /*20*/ MatAssemblyBegin_MPISBAIJ,
1464: MatAssemblyEnd_MPISBAIJ,
1465: MatSetOption_MPISBAIJ,
1466: MatZeroEntries_MPISBAIJ,
1467: /*24*/ 0,
1468: 0,
1469: 0,
1470: 0,
1471: 0,
1472: /*29*/ MatSetUp_MPISBAIJ,
1473: 0,
1474: 0,
1475: 0,
1476: 0,
1477: /*34*/ MatDuplicate_MPISBAIJ,
1478: 0,
1479: 0,
1480: 0,
1481: 0,
1482: /*39*/ MatAXPY_MPISBAIJ,
1483: MatGetSubMatrices_MPISBAIJ,
1484: MatIncreaseOverlap_MPISBAIJ,
1485: MatGetValues_MPISBAIJ,
1486: MatCopy_MPISBAIJ,
1487: /*44*/ 0,
1488: MatScale_MPISBAIJ,
1489: 0,
1490: 0,
1491: 0,
1492: /*49*/ MatSetBlockSize_MPISBAIJ,
1493: 0,
1494: 0,
1495: 0,
1496: 0,
1497: /*54*/ 0,
1498: 0,
1499: MatSetUnfactored_MPISBAIJ,
1500: 0,
1501: MatSetValuesBlocked_MPISBAIJ,
1502: /*59*/ 0,
1503: 0,
1504: 0,
1505: 0,
1506: 0,
1507: /*64*/ 0,
1508: 0,
1509: 0,
1510: 0,
1511: 0,
1512: /*69*/ MatGetRowMaxAbs_MPISBAIJ,
1513: 0,
1514: 0,
1515: 0,
1516: 0,
1517: /*74*/ 0,
1518: 0,
1519: 0,
1520: 0,
1521: 0,
1522: /*79*/ 0,
1523: 0,
1524: 0,
1525: 0,
1526: MatLoad_MPISBAIJ,
1527: /*84*/ 0,
1528: 0,
1529: 0,
1530: 0,
1531: 0,
1532: /*89*/ 0,
1533: 0,
1534: 0,
1535: 0,
1536: 0,
1537: /*94*/ 0,
1538: 0,
1539: 0,
1540: 0,
1541: 0,
1542: /*99*/ 0,
1543: 0,
1544: 0,
1545: 0,
1546: 0,
1547: /*104*/0,
1548: MatRealPart_MPISBAIJ,
1549: MatImaginaryPart_MPISBAIJ,
1550: MatGetRowUpperTriangular_MPISBAIJ,
1551: MatRestoreRowUpperTriangular_MPISBAIJ,
1552: /*109*/0,
1553: 0,
1554: 0,
1555: 0,
1556: 0,
1557: /*114*/0,
1558: 0,
1559: 0,
1560: 0,
1561: 0,
1562: /*119*/0,
1563: 0,
1564: 0,
1565: 0
1566: };
1569: EXTERN_C_BEGIN
1572: PetscErrorCode MatGetDiagonalBlock_MPISBAIJ(Mat A,Mat *a)
1573: {
1575: *a = ((Mat_MPISBAIJ *)A->data)->A;
1576: return(0);
1577: }
1578: EXTERN_C_END
1580: EXTERN_C_BEGIN
1583: PetscErrorCode MatMPISBAIJSetPreallocation_MPISBAIJ(Mat B,PetscInt bs,PetscInt d_nz,PetscInt *d_nnz,PetscInt o_nz,PetscInt *o_nnz)
1584: {
1585: Mat_MPISBAIJ *b;
1587: PetscInt i,mbs,Mbs,newbs = PetscAbs(bs);
1590: if (bs < 0){
1591: PetscOptionsBegin(((PetscObject)B)->comm,((PetscObject)B)->prefix,"Options for MPISBAIJ matrix","Mat");
1592: PetscOptionsInt("-mat_block_size","Set the blocksize used to store the matrix","MatMPIBAIJSetPreallocation",newbs,&newbs,PETSC_NULL);
1593: PetscOptionsEnd();
1594: bs = PetscAbs(bs);
1595: }
1596: if ((d_nnz || o_nnz) && newbs != bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot change blocksize from command line if setting d_nnz or o_nnz");
1597: bs = newbs;
1599: if (d_nz == PETSC_DECIDE || d_nz == PETSC_DEFAULT) d_nz = 3;
1600: if (o_nz == PETSC_DECIDE || o_nz == PETSC_DEFAULT) o_nz = 1;
1601: if (d_nz < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"d_nz cannot be less than 0: value %D",d_nz);
1602: if (o_nz < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"o_nz cannot be less than 0: value %D",o_nz);
1604: B->rmap->bs = B->cmap->bs = bs;
1605: PetscLayoutSetUp(B->rmap);
1606: PetscLayoutSetUp(B->cmap);
1608: if (d_nnz) {
1609: for (i=0; i<B->rmap->n/bs; i++) {
1610: if (d_nnz[i] < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"d_nnz cannot be less than -1: local row %D value %D",i,d_nnz[i]);
1611: }
1612: }
1613: if (o_nnz) {
1614: for (i=0; i<B->rmap->n/bs; i++) {
1615: if (o_nnz[i] < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"o_nnz cannot be less than -1: local row %D value %D",i,o_nnz[i]);
1616: }
1617: }
1619: b = (Mat_MPISBAIJ*)B->data;
1620: mbs = B->rmap->n/bs;
1621: Mbs = B->rmap->N/bs;
1622: if (mbs*bs != B->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"No of local rows %D must be divisible by blocksize %D",B->rmap->N,bs);
1624: B->rmap->bs = bs;
1625: b->bs2 = bs*bs;
1626: b->mbs = mbs;
1627: b->nbs = mbs;
1628: b->Mbs = Mbs;
1629: b->Nbs = Mbs;
1631: for (i=0; i<=b->size; i++) {
1632: b->rangebs[i] = B->rmap->range[i]/bs;
1633: }
1634: b->rstartbs = B->rmap->rstart/bs;
1635: b->rendbs = B->rmap->rend/bs;
1636:
1637: b->cstartbs = B->cmap->rstart/bs;
1638: b->cendbs = B->cmap->rend/bs;
1640: if (!B->preallocated) {
1641: MatCreate(PETSC_COMM_SELF,&b->A);
1642: MatSetSizes(b->A,B->rmap->n,B->cmap->n,B->rmap->n,B->cmap->n);
1643: MatSetType(b->A,MATSEQSBAIJ);
1644: PetscLogObjectParent(B,b->A);
1645: MatCreate(PETSC_COMM_SELF,&b->B);
1646: MatSetSizes(b->B,B->rmap->n,B->cmap->N,B->rmap->n,B->cmap->N);
1647: MatSetType(b->B,MATSEQBAIJ);
1648: PetscLogObjectParent(B,b->B);
1649: MatStashCreate_Private(((PetscObject)B)->comm,bs,&B->bstash);
1650: }
1652: MatSeqSBAIJSetPreallocation(b->A,bs,d_nz,d_nnz);
1653: MatSeqBAIJSetPreallocation(b->B,bs,o_nz,o_nnz);
1654: B->preallocated = PETSC_TRUE;
1655: return(0);
1656: }
1657: EXTERN_C_END
1659: EXTERN_C_BEGIN
1662: PetscErrorCode MatMPISBAIJSetPreallocationCSR_MPISBAIJ(Mat B,PetscInt bs,const PetscInt ii[],const PetscInt jj[],const PetscScalar V[])
1663: {
1664: PetscInt m,rstart,cstart,cend;
1665: PetscInt i,j,d,nz,nz_max=0,*d_nnz=0,*o_nnz=0;
1666: const PetscInt *JJ=0;
1667: PetscScalar *values=0;
1672: if (bs < 1) SETERRQ1(((PetscObject)B)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid block size specified, must be positive but it is %D",bs);
1673: PetscLayoutSetBlockSize(B->rmap,bs);
1674: PetscLayoutSetBlockSize(B->cmap,bs);
1675: PetscLayoutSetUp(B->rmap);
1676: PetscLayoutSetUp(B->cmap);
1677: m = B->rmap->n/bs;
1678: rstart = B->rmap->rstart/bs;
1679: cstart = B->cmap->rstart/bs;
1680: cend = B->cmap->rend/bs;
1682: if (ii[0]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"ii[0] must be 0 but it is %D",ii[0]);
1683: PetscMalloc2(m,PetscInt,&d_nnz,m,PetscInt,&o_nnz);
1684: for (i=0; i<m; i++) {
1685: nz = ii[i+1] - ii[i];
1686: if (nz < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local row %D has a negative number of columns %D",i,nz);
1687: nz_max = PetscMax(nz_max,nz);
1688: JJ = jj + ii[i];
1689: for (j=0; j<nz; j++) {
1690: if (*JJ >= cstart) break;
1691: JJ++;
1692: }
1693: d = 0;
1694: for (; j<nz; j++) {
1695: if (*JJ++ >= cend) break;
1696: d++;
1697: }
1698: d_nnz[i] = d;
1699: o_nnz[i] = nz - d;
1700: }
1701: MatMPISBAIJSetPreallocation(B,bs,0,d_nnz,0,o_nnz);
1702: PetscFree2(d_nnz,o_nnz);
1704: values = (PetscScalar*)V;
1705: if (!values) {
1706: PetscMalloc(bs*bs*nz_max*sizeof(PetscScalar),&values);
1707: PetscMemzero(values,bs*bs*nz_max*sizeof(PetscScalar));
1708: }
1709: for (i=0; i<m; i++) {
1710: PetscInt row = i + rstart;
1711: PetscInt ncols = ii[i+1] - ii[i];
1712: const PetscInt *icols = jj + ii[i];
1713: const PetscScalar *svals = values + (V ? (bs*bs*ii[i]) : 0);
1714: MatSetValuesBlocked_MPISBAIJ(B,1,&row,ncols,icols,svals,INSERT_VALUES);
1715: }
1717: if (!V) { PetscFree(values); }
1718: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
1719: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
1720: MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
1721: return(0);
1722: }
1723: EXTERN_C_END
1725: EXTERN_C_BEGIN
1726: #if defined(PETSC_HAVE_MUMPS)
1727: extern PetscErrorCode MatGetFactor_sbaij_mumps(Mat,MatFactorType,Mat*);
1728: #endif
1729: #if defined(PETSC_HAVE_SPOOLES)
1730: extern PetscErrorCode MatGetFactor_mpisbaij_spooles(Mat,MatFactorType,Mat*);
1731: #endif
1732: #if defined(PETSC_HAVE_PASTIX)
1733: extern PetscErrorCode MatGetFactor_mpisbaij_pastix(Mat,MatFactorType,Mat*);
1734: #endif
1735: EXTERN_C_END
1737: /*MC
1738: MATMPISBAIJ - MATMPISBAIJ = "mpisbaij" - A matrix type to be used for distributed symmetric sparse block matrices,
1739: based on block compressed sparse row format. Only the upper triangular portion of the "diagonal" portion of
1740: the matrix is stored.
1742: For complex numbers by default this matrix is symmetric, NOT Hermitian symmetric. To make it Hermitian symmetric you
1743: can call MatSetOption(Mat, MAT_HERMITIAN);
1745: Options Database Keys:
1746: . -mat_type mpisbaij - sets the matrix type to "mpisbaij" during a call to MatSetFromOptions()
1748: Level: beginner
1750: .seealso: MatCreateMPISBAIJ
1751: M*/
1753: EXTERN_C_BEGIN
1754: extern PetscErrorCode MatConvert_MPISBAIJ_MPISBSTRM(Mat,const MatType,MatReuse,Mat*);
1755: EXTERN_C_END
1757: EXTERN_C_BEGIN
1760: PetscErrorCode MatCreate_MPISBAIJ(Mat B)
1761: {
1762: Mat_MPISBAIJ *b;
1764: PetscBool flg;
1768: PetscNewLog(B,Mat_MPISBAIJ,&b);
1769: B->data = (void*)b;
1770: PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
1772: B->ops->destroy = MatDestroy_MPISBAIJ;
1773: B->ops->view = MatView_MPISBAIJ;
1774: B->assembled = PETSC_FALSE;
1776: B->insertmode = NOT_SET_VALUES;
1777: MPI_Comm_rank(((PetscObject)B)->comm,&b->rank);
1778: MPI_Comm_size(((PetscObject)B)->comm,&b->size);
1780: /* build local table of row and column ownerships */
1781: PetscMalloc((b->size+2)*sizeof(PetscInt),&b->rangebs);
1783: /* build cache for off array entries formed */
1784: MatStashCreate_Private(((PetscObject)B)->comm,1,&B->stash);
1785: b->donotstash = PETSC_FALSE;
1786: b->colmap = PETSC_NULL;
1787: b->garray = PETSC_NULL;
1788: b->roworiented = PETSC_TRUE;
1790: /* stuff used in block assembly */
1791: b->barray = 0;
1793: /* stuff used for matrix vector multiply */
1794: b->lvec = 0;
1795: b->Mvctx = 0;
1796: b->slvec0 = 0;
1797: b->slvec0b = 0;
1798: b->slvec1 = 0;
1799: b->slvec1a = 0;
1800: b->slvec1b = 0;
1801: b->sMvctx = 0;
1803: /* stuff for MatGetRow() */
1804: b->rowindices = 0;
1805: b->rowvalues = 0;
1806: b->getrowactive = PETSC_FALSE;
1808: /* hash table stuff */
1809: b->ht = 0;
1810: b->hd = 0;
1811: b->ht_size = 0;
1812: b->ht_flag = PETSC_FALSE;
1813: b->ht_fact = 0;
1814: b->ht_total_ct = 0;
1815: b->ht_insert_ct = 0;
1817: /* stuff for MatGetSubMatrices_MPIBAIJ_local() */
1818: b->ijonly = PETSC_FALSE;
1820: b->in_loc = 0;
1821: b->v_loc = 0;
1822: b->n_loc = 0;
1823: PetscOptionsBegin(((PetscObject)B)->comm,PETSC_NULL,"Options for loading MPISBAIJ matrix 1","Mat");
1824: PetscOptionsBool("-mat_use_hash_table","Use hash table to save memory in constructing matrix","MatSetOption",PETSC_FALSE,&flg,PETSC_NULL);
1825: if (flg) {
1826: PetscReal fact = 1.39;
1827: MatSetOption(B,MAT_USE_HASH_TABLE,PETSC_TRUE);
1828: PetscOptionsReal("-mat_use_hash_table","Use hash table factor","MatMPIBAIJSetHashTableFactor",fact,&fact,PETSC_NULL);
1829: if (fact <= 1.0) fact = 1.39;
1830: MatMPIBAIJSetHashTableFactor(B,fact);
1831: PetscInfo1(B,"Hash table Factor used %5.2f\n",fact);
1832: }
1833: PetscOptionsEnd();
1835: #if defined(PETSC_HAVE_PASTIX)
1836: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_pastix_C",
1837: "MatGetFactor_mpisbaij_pastix",
1838: MatGetFactor_mpisbaij_pastix);
1839: #endif
1840: #if defined(PETSC_HAVE_MUMPS)
1841: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_mumps_C",
1842: "MatGetFactor_sbaij_mumps",
1843: MatGetFactor_sbaij_mumps);
1844: #endif
1845: #if defined(PETSC_HAVE_SPOOLES)
1846: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_spooles_C",
1847: "MatGetFactor_mpisbaij_spooles",
1848: MatGetFactor_mpisbaij_spooles);
1849: #endif
1850: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
1851: "MatStoreValues_MPISBAIJ",
1852: MatStoreValues_MPISBAIJ);
1853: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
1854: "MatRetrieveValues_MPISBAIJ",
1855: MatRetrieveValues_MPISBAIJ);
1856: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetDiagonalBlock_C",
1857: "MatGetDiagonalBlock_MPISBAIJ",
1858: MatGetDiagonalBlock_MPISBAIJ);
1859: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPISBAIJSetPreallocation_C",
1860: "MatMPISBAIJSetPreallocation_MPISBAIJ",
1861: MatMPISBAIJSetPreallocation_MPISBAIJ);
1862: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPISBAIJSetPreallocationCSR_C",
1863: "MatMPISBAIJSetPreallocationCSR_MPISBAIJ",
1864: MatMPISBAIJSetPreallocationCSR_MPISBAIJ);
1865: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpisbaij_mpisbstrm_C",
1866: "MatConvert_MPISBAIJ_MPISBSTRM",
1867: MatConvert_MPISBAIJ_MPISBSTRM);
1869: B->symmetric = PETSC_TRUE;
1870: B->structurally_symmetric = PETSC_TRUE;
1871: B->symmetric_set = PETSC_TRUE;
1872: B->structurally_symmetric_set = PETSC_TRUE;
1873: PetscObjectChangeTypeName((PetscObject)B,MATMPISBAIJ);
1874: return(0);
1875: }
1876: EXTERN_C_END
1878: /*MC
1879: MATSBAIJ - MATSBAIJ = "sbaij" - A matrix type to be used for symmetric block sparse matrices.
1881: This matrix type is identical to MATSEQSBAIJ when constructed with a single process communicator,
1882: and MATMPISBAIJ otherwise.
1884: Options Database Keys:
1885: . -mat_type sbaij - sets the matrix type to "sbaij" during a call to MatSetFromOptions()
1887: Level: beginner
1889: .seealso: MatCreateMPISBAIJ,MATSEQSBAIJ,MATMPISBAIJ
1890: M*/
1894: /*@C
1895: MatMPISBAIJSetPreallocation - For good matrix assembly performance
1896: the user should preallocate the matrix storage by setting the parameters
1897: d_nz (or d_nnz) and o_nz (or o_nnz). By setting these parameters accurately,
1898: performance can be increased by more than a factor of 50.
1900: Collective on Mat
1902: Input Parameters:
1903: + A - the matrix
1904: . bs - size of blockk
1905: . d_nz - number of block nonzeros per block row in diagonal portion of local
1906: submatrix (same for all local rows)
1907: . d_nnz - array containing the number of block nonzeros in the various block rows
1908: in the upper triangular and diagonal part of the in diagonal portion of the local
1909: (possibly different for each block row) or PETSC_NULL. If you plan to factor the matrix you must leave room
1910: for the diagonal entry and set a value even if it is zero.
1911: . o_nz - number of block nonzeros per block row in the off-diagonal portion of local
1912: submatrix (same for all local rows).
1913: - o_nnz - array containing the number of nonzeros in the various block rows of the
1914: off-diagonal portion of the local submatrix that is right of the diagonal
1915: (possibly different for each block row) or PETSC_NULL.
1918: Options Database Keys:
1919: . -mat_no_unroll - uses code that does not unroll the loops in the
1920: block calculations (much slower)
1921: . -mat_block_size - size of the blocks to use
1923: Notes:
1925: If PETSC_DECIDE or PETSC_DETERMINE is used for a particular argument on one processor
1926: than it must be used on all processors that share the object for that argument.
1928: If the *_nnz parameter is given then the *_nz parameter is ignored
1930: Storage Information:
1931: For a square global matrix we define each processor's diagonal portion
1932: to be its local rows and the corresponding columns (a square submatrix);
1933: each processor's off-diagonal portion encompasses the remainder of the
1934: local matrix (a rectangular submatrix).
1936: The user can specify preallocated storage for the diagonal part of
1937: the local submatrix with either d_nz or d_nnz (not both). Set
1938: d_nz=PETSC_DEFAULT and d_nnz=PETSC_NULL for PETSc to control dynamic
1939: memory allocation. Likewise, specify preallocated storage for the
1940: off-diagonal part of the local submatrix with o_nz or o_nnz (not both).
1942: You can call MatGetInfo() to get information on how effective the preallocation was;
1943: for example the fields mallocs,nz_allocated,nz_used,nz_unneeded;
1944: You can also run with the option -info and look for messages with the string
1945: malloc in them to see if additional memory allocation was needed.
1947: Consider a processor that owns rows 3, 4 and 5 of a parallel matrix. In
1948: the figure below we depict these three local rows and all columns (0-11).
1950: .vb
1951: 0 1 2 3 4 5 6 7 8 9 10 11
1952: -------------------
1953: row 3 | . . . d d d o o o o o o
1954: row 4 | . . . d d d o o o o o o
1955: row 5 | . . . d d d o o o o o o
1956: -------------------
1957: .ve
1958:
1959: Thus, any entries in the d locations are stored in the d (diagonal)
1960: submatrix, and any entries in the o locations are stored in the
1961: o (off-diagonal) submatrix. Note that the d matrix is stored in
1962: MatSeqSBAIJ format and the o submatrix in MATSEQBAIJ format.
1964: Now d_nz should indicate the number of block nonzeros per row in the upper triangular
1965: plus the diagonal part of the d matrix,
1966: and o_nz should indicate the number of block nonzeros per row in the o matrix
1968: In general, for PDE problems in which most nonzeros are near the diagonal,
1969: one expects d_nz >> o_nz. For large problems you MUST preallocate memory
1970: or you will get TERRIBLE performance; see the users' manual chapter on
1971: matrices.
1973: Level: intermediate
1975: .keywords: matrix, block, aij, compressed row, sparse, parallel
1977: .seealso: MatCreate(), MatCreateSeqSBAIJ(), MatSetValues(), MatCreateMPIBAIJ()
1978: @*/
1979: PetscErrorCode MatMPISBAIJSetPreallocation(Mat B,PetscInt bs,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
1980: {
1987: PetscTryMethod(B,"MatMPISBAIJSetPreallocation_C",(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,bs,d_nz,d_nnz,o_nz,o_nnz));
1988: return(0);
1989: }
1993: /*@C
1994: MatCreateMPISBAIJ - Creates a sparse parallel matrix in symmetric block AIJ format
1995: (block compressed row). For good matrix assembly performance
1996: the user should preallocate the matrix storage by setting the parameters
1997: d_nz (or d_nnz) and o_nz (or o_nnz). By setting these parameters accurately,
1998: performance can be increased by more than a factor of 50.
2000: Collective on MPI_Comm
2002: Input Parameters:
2003: + comm - MPI communicator
2004: . bs - size of blockk
2005: . m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
2006: This value should be the same as the local size used in creating the
2007: y vector for the matrix-vector product y = Ax.
2008: . n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
2009: This value should be the same as the local size used in creating the
2010: x vector for the matrix-vector product y = Ax.
2011: . M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
2012: . N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
2013: . d_nz - number of block nonzeros per block row in diagonal portion of local
2014: submatrix (same for all local rows)
2015: . d_nnz - array containing the number of block nonzeros in the various block rows
2016: in the upper triangular portion of the in diagonal portion of the local
2017: (possibly different for each block block row) or PETSC_NULL.
2018: If you plan to factor the matrix you must leave room for the diagonal entry and
2019: set its value even if it is zero.
2020: . o_nz - number of block nonzeros per block row in the off-diagonal portion of local
2021: submatrix (same for all local rows).
2022: - o_nnz - array containing the number of nonzeros in the various block rows of the
2023: off-diagonal portion of the local submatrix (possibly different for
2024: each block row) or PETSC_NULL.
2026: Output Parameter:
2027: . A - the matrix
2029: Options Database Keys:
2030: . -mat_no_unroll - uses code that does not unroll the loops in the
2031: block calculations (much slower)
2032: . -mat_block_size - size of the blocks to use
2033: . -mat_mpi - use the parallel matrix data structures even on one processor
2034: (defaults to using SeqBAIJ format on one processor)
2036: It is recommended that one use the MatCreate(), MatSetType() and/or MatSetFromOptions(),
2037: MatXXXXSetPreallocation() paradgm instead of this routine directly.
2038: [MatXXXXSetPreallocation() is, for example, MatSeqAIJSetPreallocation]
2040: Notes:
2041: The number of rows and columns must be divisible by blocksize.
2042: This matrix type does not support complex Hermitian operation.
2044: The user MUST specify either the local or global matrix dimensions
2045: (possibly both).
2047: If PETSC_DECIDE or PETSC_DETERMINE is used for a particular argument on one processor
2048: than it must be used on all processors that share the object for that argument.
2050: If the *_nnz parameter is given then the *_nz parameter is ignored
2052: Storage Information:
2053: For a square global matrix we define each processor's diagonal portion
2054: to be its local rows and the corresponding columns (a square submatrix);
2055: each processor's off-diagonal portion encompasses the remainder of the
2056: local matrix (a rectangular submatrix).
2058: The user can specify preallocated storage for the diagonal part of
2059: the local submatrix with either d_nz or d_nnz (not both). Set
2060: d_nz=PETSC_DEFAULT and d_nnz=PETSC_NULL for PETSc to control dynamic
2061: memory allocation. Likewise, specify preallocated storage for the
2062: off-diagonal part of the local submatrix with o_nz or o_nnz (not both).
2064: Consider a processor that owns rows 3, 4 and 5 of a parallel matrix. In
2065: the figure below we depict these three local rows and all columns (0-11).
2067: .vb
2068: 0 1 2 3 4 5 6 7 8 9 10 11
2069: -------------------
2070: row 3 | . . . d d d o o o o o o
2071: row 4 | . . . d d d o o o o o o
2072: row 5 | . . . d d d o o o o o o
2073: -------------------
2074: .ve
2075:
2076: Thus, any entries in the d locations are stored in the d (diagonal)
2077: submatrix, and any entries in the o locations are stored in the
2078: o (off-diagonal) submatrix. Note that the d matrix is stored in
2079: MatSeqSBAIJ format and the o submatrix in MATSEQBAIJ format.
2081: Now d_nz should indicate the number of block nonzeros per row in the upper triangular
2082: plus the diagonal part of the d matrix,
2083: and o_nz should indicate the number of block nonzeros per row in the o matrix.
2084: In general, for PDE problems in which most nonzeros are near the diagonal,
2085: one expects d_nz >> o_nz. For large problems you MUST preallocate memory
2086: or you will get TERRIBLE performance; see the users' manual chapter on
2087: matrices.
2089: Level: intermediate
2091: .keywords: matrix, block, aij, compressed row, sparse, parallel
2093: .seealso: MatCreate(), MatCreateSeqSBAIJ(), MatSetValues(), MatCreateMPIBAIJ()
2094: @*/
2096: PetscErrorCode MatCreateMPISBAIJ(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[],Mat *A)
2097: {
2099: PetscMPIInt size;
2102: MatCreate(comm,A);
2103: MatSetSizes(*A,m,n,M,N);
2104: MPI_Comm_size(comm,&size);
2105: if (size > 1) {
2106: MatSetType(*A,MATMPISBAIJ);
2107: MatMPISBAIJSetPreallocation(*A,bs,d_nz,d_nnz,o_nz,o_nnz);
2108: } else {
2109: MatSetType(*A,MATSEQSBAIJ);
2110: MatSeqSBAIJSetPreallocation(*A,bs,d_nz,d_nnz);
2111: }
2112: return(0);
2113: }
2118: static PetscErrorCode MatDuplicate_MPISBAIJ(Mat matin,MatDuplicateOption cpvalues,Mat *newmat)
2119: {
2120: Mat mat;
2121: Mat_MPISBAIJ *a,*oldmat = (Mat_MPISBAIJ*)matin->data;
2123: PetscInt len=0,nt,bs=matin->rmap->bs,mbs=oldmat->mbs;
2124: PetscScalar *array;
2127: *newmat = 0;
2128: MatCreate(((PetscObject)matin)->comm,&mat);
2129: MatSetSizes(mat,matin->rmap->n,matin->cmap->n,matin->rmap->N,matin->cmap->N);
2130: MatSetType(mat,((PetscObject)matin)->type_name);
2131: PetscMemcpy(mat->ops,matin->ops,sizeof(struct _MatOps));
2132: PetscLayoutReference(matin->rmap,&mat->rmap);
2133: PetscLayoutReference(matin->cmap,&mat->cmap);
2134:
2135: mat->factortype = matin->factortype;
2136: mat->preallocated = PETSC_TRUE;
2137: mat->assembled = PETSC_TRUE;
2138: mat->insertmode = NOT_SET_VALUES;
2140: a = (Mat_MPISBAIJ*)mat->data;
2141: a->bs2 = oldmat->bs2;
2142: a->mbs = oldmat->mbs;
2143: a->nbs = oldmat->nbs;
2144: a->Mbs = oldmat->Mbs;
2145: a->Nbs = oldmat->Nbs;
2148: a->size = oldmat->size;
2149: a->rank = oldmat->rank;
2150: a->donotstash = oldmat->donotstash;
2151: a->roworiented = oldmat->roworiented;
2152: a->rowindices = 0;
2153: a->rowvalues = 0;
2154: a->getrowactive = PETSC_FALSE;
2155: a->barray = 0;
2156: a->rstartbs = oldmat->rstartbs;
2157: a->rendbs = oldmat->rendbs;
2158: a->cstartbs = oldmat->cstartbs;
2159: a->cendbs = oldmat->cendbs;
2161: /* hash table stuff */
2162: a->ht = 0;
2163: a->hd = 0;
2164: a->ht_size = 0;
2165: a->ht_flag = oldmat->ht_flag;
2166: a->ht_fact = oldmat->ht_fact;
2167: a->ht_total_ct = 0;
2168: a->ht_insert_ct = 0;
2169:
2170: PetscMemcpy(a->rangebs,oldmat->rangebs,(a->size+2)*sizeof(PetscInt));
2171: if (oldmat->colmap) {
2172: #if defined (PETSC_USE_CTABLE)
2173: PetscTableCreateCopy(oldmat->colmap,&a->colmap);
2174: #else
2175: PetscMalloc((a->Nbs)*sizeof(PetscInt),&a->colmap);
2176: PetscLogObjectMemory(mat,(a->Nbs)*sizeof(PetscInt));
2177: PetscMemcpy(a->colmap,oldmat->colmap,(a->Nbs)*sizeof(PetscInt));
2178: #endif
2179: } else a->colmap = 0;
2181: if (oldmat->garray && (len = ((Mat_SeqBAIJ*)(oldmat->B->data))->nbs)) {
2182: PetscMalloc(len*sizeof(PetscInt),&a->garray);
2183: PetscLogObjectMemory(mat,len*sizeof(PetscInt));
2184: PetscMemcpy(a->garray,oldmat->garray,len*sizeof(PetscInt));
2185: } else a->garray = 0;
2186:
2187: MatStashCreate_Private(((PetscObject)matin)->comm,matin->rmap->bs,&mat->bstash);
2188: VecDuplicate(oldmat->lvec,&a->lvec);
2189: PetscLogObjectParent(mat,a->lvec);
2190: VecScatterCopy(oldmat->Mvctx,&a->Mvctx);
2191: PetscLogObjectParent(mat,a->Mvctx);
2193: VecDuplicate(oldmat->slvec0,&a->slvec0);
2194: PetscLogObjectParent(mat,a->slvec0);
2195: VecDuplicate(oldmat->slvec1,&a->slvec1);
2196: PetscLogObjectParent(mat,a->slvec1);
2198: VecGetLocalSize(a->slvec1,&nt);
2199: VecGetArray(a->slvec1,&array);
2200: VecCreateSeqWithArray(PETSC_COMM_SELF,bs*mbs,array,&a->slvec1a);
2201: VecCreateSeqWithArray(PETSC_COMM_SELF,nt-bs*mbs,array+bs*mbs,&a->slvec1b);
2202: VecRestoreArray(a->slvec1,&array);
2203: VecGetArray(a->slvec0,&array);
2204: VecCreateSeqWithArray(PETSC_COMM_SELF,nt-bs*mbs,array+bs*mbs,&a->slvec0b);
2205: VecRestoreArray(a->slvec0,&array);
2206: PetscLogObjectParent(mat,a->slvec0);
2207: PetscLogObjectParent(mat,a->slvec1);
2208: PetscLogObjectParent(mat,a->slvec0b);
2209: PetscLogObjectParent(mat,a->slvec1a);
2210: PetscLogObjectParent(mat,a->slvec1b);
2212: /* VecScatterCopy(oldmat->sMvctx,&a->sMvctx); - not written yet, replaced by the lazy trick: */
2213: PetscObjectReference((PetscObject)oldmat->sMvctx);
2214: a->sMvctx = oldmat->sMvctx;
2215: PetscLogObjectParent(mat,a->sMvctx);
2217: MatDuplicate(oldmat->A,cpvalues,&a->A);
2218: PetscLogObjectParent(mat,a->A);
2219: MatDuplicate(oldmat->B,cpvalues,&a->B);
2220: PetscLogObjectParent(mat,a->B);
2221: PetscFListDuplicate(((PetscObject)matin)->qlist,&((PetscObject)mat)->qlist);
2222: *newmat = mat;
2223: return(0);
2224: }
2228: PetscErrorCode MatLoad_MPISBAIJ(Mat newmat,PetscViewer viewer)
2229: {
2231: PetscInt i,nz,j,rstart,rend;
2232: PetscScalar *vals,*buf;
2233: MPI_Comm comm = ((PetscObject)viewer)->comm;
2234: MPI_Status status;
2235: PetscMPIInt rank,size,tag = ((PetscObject)viewer)->tag,*sndcounts = 0,*browners,maxnz,*rowners,*locrowlens,mmbs;
2236: PetscInt header[4],*rowlengths = 0,M,N,m,*cols;
2237: PetscInt *procsnz = 0,jj,*mycols,*ibuf;
2238: PetscInt bs=1,Mbs,mbs,extra_rows;
2239: PetscInt *dlens,*odlens,*mask,*masked1,*masked2,rowcount,odcount;
2240: PetscInt dcount,kmax,k,nzcount,tmp,sizesset=1,grows,gcols;
2241: int fd;
2242:
2244: PetscOptionsBegin(comm,PETSC_NULL,"Options for loading MPISBAIJ matrix 2","Mat");
2245: PetscOptionsInt("-matload_block_size","Set the blocksize used to store the matrix","MatLoad",bs,&bs,PETSC_NULL);
2246: PetscOptionsEnd();
2248: MPI_Comm_size(comm,&size);
2249: MPI_Comm_rank(comm,&rank);
2250: if (!rank) {
2251: PetscViewerBinaryGetDescriptor(viewer,&fd);
2252: PetscBinaryRead(fd,(char *)header,4,PETSC_INT);
2253: if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"not matrix object");
2254: if (header[3] < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as MPISBAIJ");
2255: }
2257: if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) sizesset = 0;
2259: MPI_Bcast(header+1,3,MPIU_INT,0,comm);
2260: M = header[1]; N = header[2];
2262: /* If global rows/cols are set to PETSC_DECIDE, set it to the sizes given in the file */
2263: if (sizesset && newmat->rmap->N < 0) newmat->rmap->N = M;
2264: if (sizesset && newmat->cmap->N < 0) newmat->cmap->N = N;
2265:
2266: /* If global sizes are set, check if they are consistent with that given in the file */
2267: if (sizesset) {
2268: MatGetSize(newmat,&grows,&gcols);
2269: }
2270: if (sizesset && newmat->rmap->N != grows) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Inconsistent # of rows:Matrix in file has (%d) and input matrix has (%d)",M,grows);
2271: if (sizesset && newmat->cmap->N != gcols) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Inconsistent # of cols:Matrix in file has (%d) and input matrix has (%d)",N,gcols);
2273: if (M != N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Can only do square matrices");
2275: /*
2276: This code adds extra rows to make sure the number of rows is
2277: divisible by the blocksize
2278: */
2279: Mbs = M/bs;
2280: extra_rows = bs - M + bs*(Mbs);
2281: if (extra_rows == bs) extra_rows = 0;
2282: else Mbs++;
2283: if (extra_rows &&!rank) {
2284: PetscInfo(viewer,"Padding loaded matrix to match blocksize\n");
2285: }
2287: /* determine ownership of all rows */
2288: if (newmat->rmap->n < 0) { /* PETSC_DECIDE */
2289: mbs = Mbs/size + ((Mbs % size) > rank);
2290: m = mbs*bs;
2291: } else { /* User Set */
2292: m = newmat->rmap->n;
2293: mbs = m/bs;
2294: }
2295: PetscMalloc2(size+1,PetscMPIInt,&rowners,size+1,PetscMPIInt,&browners);
2296: mmbs = PetscMPIIntCast(mbs);
2297: MPI_Allgather(&mmbs,1,MPI_INT,rowners+1,1,MPI_INT,comm);
2298: rowners[0] = 0;
2299: for (i=2; i<=size; i++) rowners[i] += rowners[i-1];
2300: for (i=0; i<=size; i++) browners[i] = rowners[i]*bs;
2301: rstart = rowners[rank];
2302: rend = rowners[rank+1];
2303:
2304: /* distribute row lengths to all processors */
2305: PetscMalloc((rend-rstart)*bs*sizeof(PetscMPIInt),&locrowlens);
2306: if (!rank) {
2307: PetscMalloc((M+extra_rows)*sizeof(PetscInt),&rowlengths);
2308: PetscBinaryRead(fd,rowlengths,M,PETSC_INT);
2309: for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;
2310: PetscMalloc(size*sizeof(PetscMPIInt),&sndcounts);
2311: for (i=0; i<size; i++) sndcounts[i] = browners[i+1] - browners[i];
2312: MPI_Scatterv(rowlengths,sndcounts,browners,MPIU_INT,locrowlens,(rend-rstart)*bs,MPIU_INT,0,comm);
2313: PetscFree(sndcounts);
2314: } else {
2315: MPI_Scatterv(0,0,0,MPIU_INT,locrowlens,(rend-rstart)*bs,MPIU_INT,0,comm);
2316: }
2317:
2318: if (!rank) { /* procs[0] */
2319: /* calculate the number of nonzeros on each processor */
2320: PetscMalloc(size*sizeof(PetscInt),&procsnz);
2321: PetscMemzero(procsnz,size*sizeof(PetscInt));
2322: for (i=0; i<size; i++) {
2323: for (j=rowners[i]*bs; j< rowners[i+1]*bs; j++) {
2324: procsnz[i] += rowlengths[j];
2325: }
2326: }
2327: PetscFree(rowlengths);
2328:
2329: /* determine max buffer needed and allocate it */
2330: maxnz = 0;
2331: for (i=0; i<size; i++) {
2332: maxnz = PetscMax(maxnz,procsnz[i]);
2333: }
2334: PetscMalloc(maxnz*sizeof(PetscInt),&cols);
2336: /* read in my part of the matrix column indices */
2337: nz = procsnz[0];
2338: PetscMalloc(nz*sizeof(PetscInt),&ibuf);
2339: mycols = ibuf;
2340: if (size == 1) nz -= extra_rows;
2341: PetscBinaryRead(fd,mycols,nz,PETSC_INT);
2342: if (size == 1) for (i=0; i< extra_rows; i++) { mycols[nz+i] = M+i; }
2344: /* read in every ones (except the last) and ship off */
2345: for (i=1; i<size-1; i++) {
2346: nz = procsnz[i];
2347: PetscBinaryRead(fd,cols,nz,PETSC_INT);
2348: MPI_Send(cols,nz,MPIU_INT,i,tag,comm);
2349: }
2350: /* read in the stuff for the last proc */
2351: if (size != 1) {
2352: nz = procsnz[size-1] - extra_rows; /* the extra rows are not on the disk */
2353: PetscBinaryRead(fd,cols,nz,PETSC_INT);
2354: for (i=0; i<extra_rows; i++) cols[nz+i] = M+i;
2355: MPI_Send(cols,nz+extra_rows,MPIU_INT,size-1,tag,comm);
2356: }
2357: PetscFree(cols);
2358: } else { /* procs[i], i>0 */
2359: /* determine buffer space needed for message */
2360: nz = 0;
2361: for (i=0; i<m; i++) {
2362: nz += locrowlens[i];
2363: }
2364: PetscMalloc(nz*sizeof(PetscInt),&ibuf);
2365: mycols = ibuf;
2366: /* receive message of column indices*/
2367: MPI_Recv(mycols,nz,MPIU_INT,0,tag,comm,&status);
2368: MPI_Get_count(&status,MPIU_INT,&maxnz);
2369: if (maxnz != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
2370: }
2372: /* loop over local rows, determining number of off diagonal entries */
2373: PetscMalloc2(rend-rstart,PetscInt,&dlens,rend-rstart,PetscInt,&odlens);
2374: PetscMalloc3(Mbs,PetscInt,&mask,Mbs,PetscInt,&masked1,Mbs,PetscInt,&masked2);
2375: PetscMemzero(mask,Mbs*sizeof(PetscInt));
2376: PetscMemzero(masked1,Mbs*sizeof(PetscInt));
2377: PetscMemzero(masked2,Mbs*sizeof(PetscInt));
2378: rowcount = 0;
2379: nzcount = 0;
2380: for (i=0; i<mbs; i++) {
2381: dcount = 0;
2382: odcount = 0;
2383: for (j=0; j<bs; j++) {
2384: kmax = locrowlens[rowcount];
2385: for (k=0; k<kmax; k++) {
2386: tmp = mycols[nzcount++]/bs; /* block col. index */
2387: if (!mask[tmp]) {
2388: mask[tmp] = 1;
2389: if (tmp < rstart || tmp >= rend) masked2[odcount++] = tmp; /* entry in off-diag portion */
2390: else masked1[dcount++] = tmp; /* entry in diag portion */
2391: }
2392: }
2393: rowcount++;
2394: }
2395:
2396: dlens[i] = dcount; /* d_nzz[i] */
2397: odlens[i] = odcount; /* o_nzz[i] */
2399: /* zero out the mask elements we set */
2400: for (j=0; j<dcount; j++) mask[masked1[j]] = 0;
2401: for (j=0; j<odcount; j++) mask[masked2[j]] = 0;
2402: }
2403: if (!sizesset) {
2404: MatSetSizes(newmat,m,m,M+extra_rows,N+extra_rows);
2405: }
2406: MatSetOption(newmat,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);
2407: MatMPISBAIJSetPreallocation(newmat,bs,0,dlens,0,odlens);
2408:
2409: if (!rank) {
2410: PetscMalloc(maxnz*sizeof(PetscScalar),&buf);
2411: /* read in my part of the matrix numerical values */
2412: nz = procsnz[0];
2413: vals = buf;
2414: mycols = ibuf;
2415: if (size == 1) nz -= extra_rows;
2416: PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
2417: if (size == 1) for (i=0; i< extra_rows; i++) { vals[nz+i] = 1.0; }
2419: /* insert into matrix */
2420: jj = rstart*bs;
2421: for (i=0; i<m; i++) {
2422: MatSetValues(newmat,1,&jj,locrowlens[i],mycols,vals,INSERT_VALUES);
2423: mycols += locrowlens[i];
2424: vals += locrowlens[i];
2425: jj++;
2426: }
2428: /* read in other processors (except the last one) and ship out */
2429: for (i=1; i<size-1; i++) {
2430: nz = procsnz[i];
2431: vals = buf;
2432: PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
2433: MPI_Send(vals,nz,MPIU_SCALAR,i,((PetscObject)newmat)->tag,comm);
2434: }
2435: /* the last proc */
2436: if (size != 1){
2437: nz = procsnz[i] - extra_rows;
2438: vals = buf;
2439: PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
2440: for (i=0; i<extra_rows; i++) vals[nz+i] = 1.0;
2441: MPI_Send(vals,nz+extra_rows,MPIU_SCALAR,size-1,((PetscObject)newmat)->tag,comm);
2442: }
2443: PetscFree(procsnz);
2445: } else {
2446: /* receive numeric values */
2447: PetscMalloc(nz*sizeof(PetscScalar),&buf);
2449: /* receive message of values*/
2450: vals = buf;
2451: mycols = ibuf;
2452: MPI_Recv(vals,nz,MPIU_SCALAR,0,((PetscObject)newmat)->tag,comm,&status);
2453: MPI_Get_count(&status,MPIU_SCALAR,&maxnz);
2454: if (maxnz != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
2456: /* insert into matrix */
2457: jj = rstart*bs;
2458: for (i=0; i<m; i++) {
2459: MatSetValues_MPISBAIJ(newmat,1,&jj,locrowlens[i],mycols,vals,INSERT_VALUES);
2460: mycols += locrowlens[i];
2461: vals += locrowlens[i];
2462: jj++;
2463: }
2464: }
2466: PetscFree(locrowlens);
2467: PetscFree(buf);
2468: PetscFree(ibuf);
2469: PetscFree2(rowners,browners);
2470: PetscFree2(dlens,odlens);
2471: PetscFree3(mask,masked1,masked2);
2472: MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);
2473: MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);
2474: return(0);
2475: }
2479: /*XXXXX@
2480: MatMPISBAIJSetHashTableFactor - Sets the factor required to compute the size of the HashTable.
2482: Input Parameters:
2483: . mat - the matrix
2484: . fact - factor
2486: Not Collective on Mat, each process can have a different hash factor
2488: Level: advanced
2490: Notes:
2491: This can also be set by the command line option: -mat_use_hash_table fact
2493: .keywords: matrix, hashtable, factor, HT
2495: .seealso: MatSetOption()
2496: @XXXXX*/
2501: PetscErrorCode MatGetRowMaxAbs_MPISBAIJ(Mat A,Vec v,PetscInt idx[])
2502: {
2503: Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
2504: Mat_SeqBAIJ *b = (Mat_SeqBAIJ*)(a->B)->data;
2505: PetscReal atmp;
2506: PetscReal *work,*svalues,*rvalues;
2508: PetscInt i,bs,mbs,*bi,*bj,brow,j,ncols,krow,kcol,col,row,Mbs,bcol;
2509: PetscMPIInt rank,size;
2510: PetscInt *rowners_bs,dest,count,source;
2511: PetscScalar *va;
2512: MatScalar *ba;
2513: MPI_Status stat;
2516: if (idx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Send email to petsc-maint@mcs.anl.gov");
2517: MatGetRowMaxAbs(a->A,v,PETSC_NULL);
2518: VecGetArray(v,&va);
2520: MPI_Comm_size(((PetscObject)A)->comm,&size);
2521: MPI_Comm_rank(((PetscObject)A)->comm,&rank);
2523: bs = A->rmap->bs;
2524: mbs = a->mbs;
2525: Mbs = a->Mbs;
2526: ba = b->a;
2527: bi = b->i;
2528: bj = b->j;
2530: /* find ownerships */
2531: rowners_bs = A->rmap->range;
2533: /* each proc creates an array to be distributed */
2534: PetscMalloc(bs*Mbs*sizeof(PetscReal),&work);
2535: PetscMemzero(work,bs*Mbs*sizeof(PetscReal));
2537: /* row_max for B */
2538: if (rank != size-1){
2539: for (i=0; i<mbs; i++) {
2540: ncols = bi[1] - bi[0]; bi++;
2541: brow = bs*i;
2542: for (j=0; j<ncols; j++){
2543: bcol = bs*(*bj);
2544: for (kcol=0; kcol<bs; kcol++){
2545: col = bcol + kcol; /* local col index */
2546: col += rowners_bs[rank+1]; /* global col index */
2547: for (krow=0; krow<bs; krow++){
2548: atmp = PetscAbsScalar(*ba); ba++;
2549: row = brow + krow; /* local row index */
2550: if (PetscRealPart(va[row]) < atmp) va[row] = atmp;
2551: if (work[col] < atmp) work[col] = atmp;
2552: }
2553: }
2554: bj++;
2555: }
2556: }
2558: /* send values to its owners */
2559: for (dest=rank+1; dest<size; dest++){
2560: svalues = work + rowners_bs[dest];
2561: count = rowners_bs[dest+1]-rowners_bs[dest];
2562: MPI_Send(svalues,count,MPIU_REAL,dest,rank,((PetscObject)A)->comm);
2563: }
2564: }
2565:
2566: /* receive values */
2567: if (rank){
2568: rvalues = work;
2569: count = rowners_bs[rank+1]-rowners_bs[rank];
2570: for (source=0; source<rank; source++){
2571: MPI_Recv(rvalues,count,MPIU_REAL,MPI_ANY_SOURCE,MPI_ANY_TAG,((PetscObject)A)->comm,&stat);
2572: /* process values */
2573: for (i=0; i<count; i++){
2574: if (PetscRealPart(va[i]) < rvalues[i]) va[i] = rvalues[i];
2575: }
2576: }
2577: }
2579: VecRestoreArray(v,&va);
2580: PetscFree(work);
2581: return(0);
2582: }
2586: PetscErrorCode MatSOR_MPISBAIJ(Mat matin,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
2587: {
2588: Mat_MPISBAIJ *mat = (Mat_MPISBAIJ*)matin->data;
2589: PetscErrorCode ierr;
2590: PetscInt mbs=mat->mbs,bs=matin->rmap->bs;
2591: PetscScalar *x,*ptr,*from;
2592: Vec bb1;
2593: const PetscScalar *b;
2594:
2596: if (its <= 0 || lits <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
2597: if (bs > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented");
2599: if (flag == SOR_APPLY_UPPER) {
2600: (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx);
2601: return(0);
2602: }
2604: if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP){
2605: if ( flag & SOR_ZERO_INITIAL_GUESS ) {
2606: (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,lits,xx);
2607: its--;
2608: }
2610: VecDuplicate(bb,&bb1);
2611: while (its--){
2612:
2613: /* lower triangular part: slvec0b = - B^T*xx */
2614: (*mat->B->ops->multtranspose)(mat->B,xx,mat->slvec0b);
2615:
2616: /* copy xx into slvec0a */
2617: VecGetArray(mat->slvec0,&ptr);
2618: VecGetArray(xx,&x);
2619: PetscMemcpy(ptr,x,bs*mbs*sizeof(MatScalar));
2620: VecRestoreArray(mat->slvec0,&ptr);
2622: VecScale(mat->slvec0,-1.0);
2624: /* copy bb into slvec1a */
2625: VecGetArray(mat->slvec1,&ptr);
2626: VecGetArrayRead(bb,&b);
2627: PetscMemcpy(ptr,b,bs*mbs*sizeof(MatScalar));
2628: VecRestoreArray(mat->slvec1,&ptr);
2630: /* set slvec1b = 0 */
2631: VecSet(mat->slvec1b,0.0);
2633: VecScatterBegin(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD);
2634: VecRestoreArray(xx,&x);
2635: VecRestoreArrayRead(bb,&b);
2636: VecScatterEnd(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD);
2638: /* upper triangular part: bb1 = bb1 - B*x */
2639: (*mat->B->ops->multadd)(mat->B,mat->slvec1b,mat->slvec1a,bb1);
2640:
2641: /* local diagonal sweep */
2642: (*mat->A->ops->sor)(mat->A,bb1,omega,SOR_SYMMETRIC_SWEEP,fshift,lits,lits,xx);
2643: }
2644: VecDestroy(&bb1);
2645: } else if ((flag & SOR_LOCAL_FORWARD_SWEEP) && (its == 1) && (flag & SOR_ZERO_INITIAL_GUESS)){
2646: (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx);
2647: } else if ((flag & SOR_LOCAL_BACKWARD_SWEEP) && (its == 1) && (flag & SOR_ZERO_INITIAL_GUESS)){
2648: (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx);
2649: } else if (flag & SOR_EISENSTAT) {
2650: Vec xx1;
2651: PetscBool hasop;
2652: const PetscScalar *diag;
2653: PetscScalar *sl,scale = (omega - 2.0)/omega;
2654: PetscInt i,n;
2656: if (!mat->xx1) {
2657: VecDuplicate(bb,&mat->xx1);
2658: VecDuplicate(bb,&mat->bb1);
2659: }
2660: xx1 = mat->xx1;
2661: bb1 = mat->bb1;
2663: (*mat->A->ops->sor)(mat->A,bb,omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_BACKWARD_SWEEP),fshift,lits,1,xx);
2665: if (!mat->diag) {
2666: /* this is wrong for same matrix with new nonzero values */
2667: MatGetVecs(matin,&mat->diag,PETSC_NULL);
2668: MatGetDiagonal(matin,mat->diag);
2669: }
2670: MatHasOperation(matin,MATOP_MULT_DIAGONAL_BLOCK,&hasop);
2672: if (hasop) {
2673: MatMultDiagonalBlock(matin,xx,bb1);
2674: VecAYPX(mat->slvec1a,scale,bb);
2675: } else {
2676: /*
2677: These two lines are replaced by code that may be a bit faster for a good compiler
2678: VecPointwiseMult(mat->slvec1a,mat->diag,xx);
2679: VecAYPX(mat->slvec1a,scale,bb);
2680: */
2681: VecGetArray(mat->slvec1a,&sl);
2682: VecGetArrayRead(mat->diag,&diag);
2683: VecGetArrayRead(bb,&b);
2684: VecGetArray(xx,&x);
2685: VecGetLocalSize(xx,&n);
2686: if (omega == 1.0) {
2687: for (i=0; i<n; i++) {
2688: sl[i] = b[i] - diag[i]*x[i];
2689: }
2690: PetscLogFlops(2.0*n);
2691: } else {
2692: for (i=0; i<n; i++) {
2693: sl[i] = b[i] + scale*diag[i]*x[i];
2694: }
2695: PetscLogFlops(3.0*n);
2696: }
2697: VecRestoreArray(mat->slvec1a,&sl);
2698: VecRestoreArrayRead(mat->diag,&diag);
2699: VecRestoreArrayRead(bb,&b);
2700: VecRestoreArray(xx,&x);
2701: }
2703: /* multiply off-diagonal portion of matrix */
2704: VecSet(mat->slvec1b,0.0);
2705: (*mat->B->ops->multtranspose)(mat->B,xx,mat->slvec0b);
2706: VecGetArray(mat->slvec0,&from);
2707: VecGetArray(xx,&x);
2708: PetscMemcpy(from,x,bs*mbs*sizeof(MatScalar));
2709: VecRestoreArray(mat->slvec0,&from);
2710: VecRestoreArray(xx,&x);
2711: VecScatterBegin(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD);
2712: VecScatterEnd(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD);
2713: (*mat->B->ops->multadd)(mat->B,mat->slvec1b,mat->slvec1a,mat->slvec1a);
2715: /* local sweep */
2716: (*mat->A->ops->sor)(mat->A,mat->slvec1a,omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_FORWARD_SWEEP),fshift,lits,1,xx1);
2717: VecAXPY(xx,1.0,xx1);
2718: } else {
2719: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatSORType is not supported for SBAIJ matrix format");
2720: }
2721: return(0);
2722: }
2726: PetscErrorCode MatSOR_MPISBAIJ_2comm(Mat matin,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
2727: {
2728: Mat_MPISBAIJ *mat = (Mat_MPISBAIJ*)matin->data;
2730: Vec lvec1,bb1;
2731:
2733: if (its <= 0 || lits <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
2734: if (matin->rmap->bs > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented");
2736: if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP){
2737: if ( flag & SOR_ZERO_INITIAL_GUESS ) {
2738: (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,lits,xx);
2739: its--;
2740: }
2742: VecDuplicate(mat->lvec,&lvec1);
2743: VecDuplicate(bb,&bb1);
2744: while (its--){
2745: VecScatterBegin(mat->Mvctx,xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD);
2746:
2747: /* lower diagonal part: bb1 = bb - B^T*xx */
2748: (*mat->B->ops->multtranspose)(mat->B,xx,lvec1);
2749: VecScale(lvec1,-1.0);
2751: VecScatterEnd(mat->Mvctx,xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD);
2752: VecCopy(bb,bb1);
2753: VecScatterBegin(mat->Mvctx,lvec1,bb1,ADD_VALUES,SCATTER_REVERSE);
2755: /* upper diagonal part: bb1 = bb1 - B*x */
2756: VecScale(mat->lvec,-1.0);
2757: (*mat->B->ops->multadd)(mat->B,mat->lvec,bb1,bb1);
2759: VecScatterEnd(mat->Mvctx,lvec1,bb1,ADD_VALUES,SCATTER_REVERSE);
2760:
2761: /* diagonal sweep */
2762: (*mat->A->ops->sor)(mat->A,bb1,omega,SOR_SYMMETRIC_SWEEP,fshift,lits,lits,xx);
2763: }
2764: VecDestroy(&lvec1);
2765: VecDestroy(&bb1);
2766: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatSORType is not supported for SBAIJ matrix format");
2767: return(0);
2768: }
2772: /*@
2773: MatCreateMPISBAIJWithArrays - creates a MPI SBAIJ matrix using arrays that contain in standard
2774: CSR format the local rows.
2776: Collective on MPI_Comm
2778: Input Parameters:
2779: + comm - MPI communicator
2780: . bs - the block size, only a block size of 1 is supported
2781: . m - number of local rows (Cannot be PETSC_DECIDE)
2782: . n - This value should be the same as the local size used in creating the
2783: x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have
2784: calculated if N is given) For square matrices n is almost always m.
2785: . M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
2786: . N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
2787: . i - row indices
2788: . j - column indices
2789: - a - matrix values
2791: Output Parameter:
2792: . mat - the matrix
2794: Level: intermediate
2796: Notes:
2797: The i, j, and a arrays ARE copied by this routine into the internal format used by PETSc;
2798: thus you CANNOT change the matrix entries by changing the values of a[] after you have
2799: called this routine. Use MatCreateMPIAIJWithSplitArrays() to avoid needing to copy the arrays.
2801: The i and j indices are 0 based, and i indices are indices corresponding to the local j array.
2803: .keywords: matrix, aij, compressed row, sparse, parallel
2805: .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatMPIAIJSetPreallocation(), MatMPIAIJSetPreallocationCSR(),
2806: MPIAIJ, MatCreateMPIAIJ(), MatCreateMPIAIJWithSplitArrays()
2807: @*/
2808: PetscErrorCode MatCreateMPISBAIJWithArrays(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt i[],const PetscInt j[],const PetscScalar a[],Mat *mat)
2809: {
2814: if (i[0]) {
2815: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"i (row indices) must start with 0");
2816: }
2817: if (m < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"local number of rows (m) cannot be PETSC_DECIDE, or negative");
2818: MatCreate(comm,mat);
2819: MatSetSizes(*mat,m,n,M,N);
2820: MatSetType(*mat,MATMPISBAIJ);
2821: MatMPISBAIJSetPreallocationCSR(*mat,bs,i,j,a);
2822: return(0);
2823: }
2828: /*@C
2829: MatMPISBAIJSetPreallocationCSR - Allocates memory for a sparse parallel matrix in BAIJ format
2830: (the default parallel PETSc format).
2832: Collective on MPI_Comm
2834: Input Parameters:
2835: + A - the matrix
2836: . bs - the block size
2837: . i - the indices into j for the start of each local row (starts with zero)
2838: . j - the column indices for each local row (starts with zero) these must be sorted for each row
2839: - v - optional values in the matrix
2841: Level: developer
2843: .keywords: matrix, aij, compressed row, sparse, parallel
2845: .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatMPIBAIJSetPreallocation(), MatCreateMPIAIJ(), MPIAIJ
2846: @*/
2847: PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat B,PetscInt bs,const PetscInt i[],const PetscInt j[], const PetscScalar v[])
2848: {
2852: PetscTryMethod(B,"MatMPISBAIJSetPreallocationCSR_C",(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]),(B,bs,i,j,v));
2853: return(0);
2854: }