Actual source code: mmsbaij.c
1: #define PETSCMAT_DLL
3: /*
4: Support for the parallel SBAIJ matrix vector multiply
5: */
6: #include src/mat/impls/sbaij/mpi/mpisbaij.h
8: EXTERN PetscErrorCode MatSetValues_SeqSBAIJ(Mat,PetscInt,const PetscInt [],PetscInt,const PetscInt [],const PetscScalar [],InsertMode);
13: PetscErrorCode MatSetUpMultiply_MPISBAIJ(Mat mat)
14: {
15: Mat_MPISBAIJ *sbaij = (Mat_MPISBAIJ*)mat->data;
16: Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)(sbaij->B->data);
18: PetscInt Nbs = sbaij->Nbs,i,j,*indices,*aj = B->j,ec = 0,*garray,*sgarray;
19: PetscInt bs = mat->rmap.bs,*stmp,mbs=sbaij->mbs, vec_size,nt;
20: IS from,to;
21: Vec gvec;
22: PetscMPIInt rank=sbaij->rank,lsize,size=sbaij->size;
23: PetscInt *owners=sbaij->rangebs,*sowners,*ec_owner,k;
24: PetscScalar *ptr;
27: if (sbaij->sMvctx) {
28: /* This two lines should be in DisAssemble_MPISBAIJ(). Don't know why it causes crash there? */
29: VecScatterDestroy(sbaij->sMvctx);
30: sbaij->sMvctx = 0;
31: }
32:
33: /* For the first stab we make an array as long as the number of columns */
34: /* mark those columns that are in sbaij->B */
35: PetscMalloc((Nbs+1)*sizeof(PetscInt),&indices);
36: PetscMemzero(indices,Nbs*sizeof(PetscInt));
37: for (i=0; i<mbs; i++) {
38: for (j=0; j<B->ilen[i]; j++) {
39: if (!indices[aj[B->i[i] + j]]) ec++;
40: indices[aj[B->i[i] + j] ] = 1;
41: }
42: }
44: /* form arrays of columns we need */
45: PetscMalloc((ec+1)*sizeof(PetscInt),&garray);
46: PetscMalloc((3*ec+1)*sizeof(PetscInt),&sgarray);
47: ec_owner = sgarray + 2*ec;
48:
49: ec = 0;
50: for (j=0; j<size; j++){
51: for (i=owners[j]; i<owners[j+1]; i++){
52: if (indices[i]) {
53: garray[ec] = i;
54: ec_owner[ec] = j;
55: ec++;
56: }
57: }
58: }
60: /* make indices now point into garray */
61: for (i=0; i<ec; i++) indices[garray[i]] = i;
63: /* compact out the extra columns in B */
64: for (i=0; i<mbs; i++) {
65: for (j=0; j<B->ilen[i]; j++) aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
66: }
67: B->nbs = ec;
68: sbaij->B->cmap.n = sbaij->B->cmap.N = ec*mat->rmap.bs;
69: PetscMapInitialize(sbaij->B->comm,&(sbaij->B->cmap));
70: PetscFree(indices);
72: /* create local vector that is used to scatter into */
73: VecCreateSeq(PETSC_COMM_SELF,ec*bs,&sbaij->lvec);
75: /* create two temporary index sets for building scatter-gather */
76: PetscMalloc((2*ec+1)*sizeof(PetscInt),&stmp);
77: for (i=0; i<ec; i++) stmp[i] = bs*garray[i];
78: ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,&from);
79:
80: for (i=0; i<ec; i++) { stmp[i] = bs*i; }
81: ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,&to);
83: /* generate the scatter context
84: -- Mvctx and lvec are not used by MatMult_MPISBAIJ(), but usefule for some applications */
85: VecCreateMPI(mat->comm,mat->cmap.n,mat->cmap.N,&gvec);
86: VecScatterCreate(gvec,from,sbaij->lvec,to,&sbaij->Mvctx);
87: VecScatterPostRecvs(gvec,sbaij->lvec,INSERT_VALUES,SCATTER_FORWARD,sbaij->Mvctx);
89: sbaij->garray = garray;
90: PetscLogObjectParent(mat,sbaij->Mvctx);
91: PetscLogObjectParent(mat,sbaij->lvec);
92: PetscLogObjectParent(mat,from);
93: PetscLogObjectParent(mat,to);
95: ISDestroy(from);
96: ISDestroy(to);
98: /* create parallel vector that is used by SBAIJ matrix to scatter from/into */
99: lsize = (mbs + ec)*bs;
100: VecCreateMPI(mat->comm,lsize,PETSC_DETERMINE,&sbaij->slvec0);
101: VecDuplicate(sbaij->slvec0,&sbaij->slvec1);
102: VecGetSize(sbaij->slvec0,&vec_size);
104: sowners = sbaij->slvec0->map.range;
105:
106: /* x index in the IS sfrom */
107: for (i=0; i<ec; i++) {
108: j = ec_owner[i];
109: sgarray[i] = garray[i] + (sowners[j]/bs - owners[j]);
110: }
111: /* b index in the IS sfrom */
112: k = sowners[rank]/bs + mbs;
113: for (i=ec,j=0; i< 2*ec; i++,j++) sgarray[i] = k + j;
114:
115: for (i=0; i<2*ec; i++) stmp[i] = bs*sgarray[i];
116: ISCreateBlock(PETSC_COMM_SELF,bs,2*ec,stmp,&from);
117:
118: /* x index in the IS sto */
119: k = sowners[rank]/bs + mbs;
120: for (i=0; i<ec; i++) stmp[i] = bs*(k + i);
121: /* b index in the IS sto */
122: for (i=ec; i<2*ec; i++) stmp[i] = bs*sgarray[i-ec];
124: ISCreateBlock(PETSC_COMM_SELF,bs,2*ec,stmp,&to);
126: /* gnerate the SBAIJ scatter context */
127: VecScatterCreate(sbaij->slvec0,from,sbaij->slvec1,to,&sbaij->sMvctx);
128:
129: /*
130: Post the receives for the first matrix vector product. We sync-chronize after
131: this on the chance that the user immediately calls MatMult() after assemblying
132: the matrix.
133: */
134: VecScatterPostRecvs(sbaij->slvec0,sbaij->slvec1,INSERT_VALUES,SCATTER_FORWARD,sbaij->sMvctx);
136: VecGetLocalSize(sbaij->slvec1,&nt);
137: VecGetArray(sbaij->slvec1,&ptr);
138: VecCreateSeqWithArray(PETSC_COMM_SELF,bs*mbs,ptr,&sbaij->slvec1a);
139: VecCreateSeqWithArray(PETSC_COMM_SELF,nt-bs*mbs,ptr+bs*mbs,&sbaij->slvec1b);
140: VecRestoreArray(sbaij->slvec1,&ptr);
142: VecGetArray(sbaij->slvec0,&ptr);
143: VecCreateSeqWithArray(PETSC_COMM_SELF,nt-bs*mbs,ptr+bs*mbs,&sbaij->slvec0b);
144: VecRestoreArray(sbaij->slvec0,&ptr);
146: PetscFree(stmp);
147: MPI_Barrier(mat->comm);
148:
149: PetscLogObjectParent(mat,sbaij->sMvctx);
150: PetscLogObjectParent(mat,sbaij->slvec0);
151: PetscLogObjectParent(mat,sbaij->slvec1);
152: PetscLogObjectParent(mat,sbaij->slvec0b);
153: PetscLogObjectParent(mat,sbaij->slvec1a);
154: PetscLogObjectParent(mat,sbaij->slvec1b);
155: PetscLogObjectParent(mat,from);
156: PetscLogObjectParent(mat,to);
157:
158: PetscLogObjectMemory(mat,(ec+1)*sizeof(PetscInt));
159: ISDestroy(from);
160: ISDestroy(to);
161: VecDestroy(gvec);
162: PetscFree(sgarray);
163: return(0);
164: }
168: PetscErrorCode MatSetUpMultiply_MPISBAIJ_2comm(Mat mat)
169: {
170: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
171: Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)(baij->B->data);
172: PetscErrorCode ierr;
173: PetscInt i,j,*aj = B->j,ec = 0,*garray;
174: PetscInt bs = mat->rmap.bs,*stmp;
175: IS from,to;
176: Vec gvec;
177: #if defined (PETSC_USE_CTABLE)
178: PetscTable gid1_lid1;
179: PetscTablePosition tpos;
180: PetscInt gid,lid;
181: #else
182: PetscInt Nbs = baij->Nbs,*indices;
183: #endif
186: #if defined (PETSC_USE_CTABLE)
187: /* use a table - Mark Adams */
188: PetscTableCreate(B->mbs,&gid1_lid1);
189: for (i=0; i<B->mbs; i++) {
190: for (j=0; j<B->ilen[i]; j++) {
191: PetscInt data,gid1 = aj[B->i[i]+j] + 1;
192: PetscTableFind(gid1_lid1,gid1,&data);
193: if (!data) {
194: /* one based table */
195: PetscTableAdd(gid1_lid1,gid1,++ec);
196: }
197: }
198: }
199: /* form array of columns we need */
200: PetscMalloc((ec+1)*sizeof(PetscInt),&garray);
201: PetscTableGetHeadPosition(gid1_lid1,&tpos);
202: while (tpos) {
203: PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);
204: gid--; lid--;
205: garray[lid] = gid;
206: }
207: PetscSortInt(ec,garray);
208: PetscTableRemoveAll(gid1_lid1);
209: for (i=0; i<ec; i++) {
210: PetscTableAdd(gid1_lid1,garray[i]+1,i+1);
211: }
212: /* compact out the extra columns in B */
213: for (i=0; i<B->mbs; i++) {
214: for (j=0; j<B->ilen[i]; j++) {
215: PetscInt gid1 = aj[B->i[i] + j] + 1;
216: PetscTableFind(gid1_lid1,gid1,&lid);
217: lid --;
218: aj[B->i[i]+j] = lid;
219: }
220: }
221: B->nbs = ec;
222: baij->B->cmap.n = baij->B->cmap.N = ec*mat->rmap.bs;
223: PetscMapInitialize(baij->B->comm,&(baij->B->cmap));
224: PetscTableDelete(gid1_lid1);
225: /* Mark Adams */
226: #else
227: /* For the first stab we make an array as long as the number of columns */
228: /* mark those columns that are in baij->B */
229: PetscMalloc((Nbs+1)*sizeof(PetscInt),&indices);
230: PetscMemzero(indices,Nbs*sizeof(PetscInt));
231: for (i=0; i<B->mbs; i++) {
232: for (j=0; j<B->ilen[i]; j++) {
233: if (!indices[aj[B->i[i] + j]]) ec++;
234: indices[aj[B->i[i] + j] ] = 1;
235: }
236: }
238: /* form array of columns we need */
239: PetscMalloc((ec+1)*sizeof(PetscInt),&garray);
240: ec = 0;
241: for (i=0; i<Nbs; i++) {
242: if (indices[i]) {
243: garray[ec++] = i;
244: }
245: }
247: /* make indices now point into garray */
248: for (i=0; i<ec; i++) {
249: indices[garray[i]] = i;
250: }
252: /* compact out the extra columns in B */
253: for (i=0; i<B->mbs; i++) {
254: for (j=0; j<B->ilen[i]; j++) {
255: aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
256: }
257: }
258: B->nbs = ec;
259: baij->B->cmap.n = ec*mat->rmap.bs;
260: PetscFree(indices);
261: #endif
262:
263: /* create local vector that is used to scatter into */
264: VecCreateSeq(PETSC_COMM_SELF,ec*bs,&baij->lvec);
266: /* create two temporary index sets for building scatter-gather */
267: for (i=0; i<ec; i++) {
268: garray[i] = bs*garray[i];
269: }
270: ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,&from);
271: for (i=0; i<ec; i++) {
272: garray[i] = garray[i]/bs;
273: }
275: PetscMalloc((ec+1)*sizeof(PetscInt),&stmp);
276: for (i=0; i<ec; i++) { stmp[i] = bs*i; }
277: ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,&to);
278: PetscFree(stmp);
280: /* create temporary global vector to generate scatter context */
281: /* this is inefficient, but otherwise we must do either
282: 1) save garray until the first actual scatter when the vector is known or
283: 2) have another way of generating a scatter context without a vector.*/
284: VecCreateMPI(mat->comm,mat->cmap.n,mat->cmap.N,&gvec);
286: /* gnerate the scatter context */
287: VecScatterCreate(gvec,from,baij->lvec,to,&baij->Mvctx);
289: /*
290: Post the receives for the first matrix vector product. We sync-chronize after
291: this on the chance that the user immediately calls MatMult() after assemblying
292: the matrix.
293: */
294: VecScatterPostRecvs(gvec,baij->lvec,INSERT_VALUES,SCATTER_FORWARD,baij->Mvctx);
295: MPI_Barrier(mat->comm);
297: PetscLogObjectParent(mat,baij->Mvctx);
298: PetscLogObjectParent(mat,baij->lvec);
299: PetscLogObjectParent(mat,from);
300: PetscLogObjectParent(mat,to);
301: baij->garray = garray;
302: PetscLogObjectMemory(mat,(ec+1)*sizeof(PetscInt));
303: ISDestroy(from);
304: ISDestroy(to);
305: VecDestroy(gvec);
306: return(0);
307: }
309: /*
310: Takes the local part of an already assembled MPISBAIJ matrix
311: and disassembles it. This is to allow new nonzeros into the matrix
312: that require more communication in the matrix vector multiply.
313: Thus certain data-structures must be rebuilt.
315: Kind of slow! But that's what application programmers get when
316: they are sloppy.
317: */
320: PetscErrorCode DisAssemble_MPISBAIJ(Mat A)
321: {
322: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)A->data;
323: Mat B = baij->B,Bnew;
324: Mat_SeqBAIJ *Bbaij = (Mat_SeqBAIJ*)B->data;
326: PetscInt i,j,mbs=Bbaij->mbs,n = A->cmap.n,col,*garray=baij->garray;
327: PetscInt k,bs=A->rmap.bs,bs2=baij->bs2,*rvals,*nz,ec,m=A->rmap.n;
328: MatScalar *a = Bbaij->a;
329: PetscScalar *atmp;
330: #if defined(PETSC_USE_MAT_SINGLE)
331: PetscInt l;
332: #endif
335: #if defined(PETSC_USE_MAT_SINGLE)
336: PetscMalloc(A->rmap.bs*sizeof(PetscScalar),&atmp);
337: #endif
338: /* free stuff related to matrix-vec multiply */
339: VecGetSize(baij->lvec,&ec); /* needed for PetscLogObjectMemory below */
340: VecDestroy(baij->lvec);
341: baij->lvec = 0;
342: VecScatterDestroy(baij->Mvctx);
343: baij->Mvctx = 0;
345: VecDestroy(baij->slvec0);
346: VecDestroy(baij->slvec0b);
347: baij->slvec0 = 0;
348: VecDestroy(baij->slvec1);
349: VecDestroy(baij->slvec1a);
350: VecDestroy(baij->slvec1b);
351: baij->slvec1 = 0;
353: if (baij->colmap) {
354: #if defined (PETSC_USE_CTABLE)
355: PetscTableDelete(baij->colmap); baij->colmap = 0;
356: #else
357: PetscFree(baij->colmap);
358: baij->colmap = 0;
359: PetscLogObjectMemory(A,-Bbaij->nbs*sizeof(PetscInt));
360: #endif
361: }
363: /* make sure that B is assembled so we can access its values */
364: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
365: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
367: /* invent new B and copy stuff over */
368: PetscMalloc(mbs*sizeof(PetscInt),&nz);
369: for (i=0; i<mbs; i++) {
370: nz[i] = Bbaij->i[i+1]-Bbaij->i[i];
371: }
372: MatCreate(PETSC_COMM_SELF,&Bnew);
373: MatSetSizes(Bnew,m,n,m,n);
374: MatSetType(Bnew,B->type_name);
375: MatSeqBAIJSetPreallocation(Bnew,B->rmap.bs,0,nz);
376: PetscFree(nz);
377:
378: PetscMalloc(bs*sizeof(PetscInt),&rvals);
379: for (i=0; i<mbs; i++) {
380: rvals[0] = bs*i;
381: for (j=1; j<bs; j++) { rvals[j] = rvals[j-1] + 1; }
382: for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) {
383: col = garray[Bbaij->j[j]]*bs;
384: for (k=0; k<bs; k++) {
385: #if defined(PETSC_USE_MAT_SINGLE)
386: for (l=0; l<bs; l++) atmp[l] = a[j*bs2+l];
387: #else
388: atmp = a+j*bs2 + k*bs;
389: #endif
390: MatSetValues_SeqSBAIJ(Bnew,bs,rvals,1,&col,atmp,B->insertmode);
391: col++;
392: }
393: }
394: }
395: #if defined(PETSC_USE_MAT_SINGLE)
396: PetscFree(atmp);
397: #endif
398: PetscFree(baij->garray);
399: baij->garray = 0;
400: PetscFree(rvals);
401: PetscLogObjectMemory(A,-ec*sizeof(PetscInt));
402: MatDestroy(B);
403: PetscLogObjectParent(A,Bnew);
404: baij->B = Bnew;
405: A->was_assembled = PETSC_FALSE;
406: return(0);
407: }