Actual source code: mmsbaij.c
1: /*$Id: mmsbaij.c,v 1.10 2001/08/07 03:03:05 balay Exp $*/
3: /*
4: Support for the parallel SBAIJ matrix vector multiply
5: */
6: #include src/mat/impls/sbaij/mpi/mpisbaij.h
7: #include src/vec/vecimpl.h
9: extern int MatSetValues_SeqSBAIJ(Mat,int,const int [],int,const int [],const PetscScalar [],InsertMode);
14: int MatSetUpMultiply_MPISBAIJ(Mat mat)
15: {
16: Mat_MPISBAIJ *sbaij = (Mat_MPISBAIJ*)mat->data;
17: Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)(sbaij->B->data);
18: int Nbs = sbaij->Nbs,i,j,*indices,*aj = B->j,ierr,ec = 0,*garray,*sgarray;
19: int bs = sbaij->bs,*stmp,mbs=sbaij->mbs, vec_size,nt;
20: IS from,to;
21: Vec gvec;
22: int rank=sbaij->rank,lsize,size=sbaij->size;
23: int *owners=sbaij->rowners,*sowners,*ec_owner,k;
24: PetscMap vecmap;
25: PetscScalar *ptr;
28: if (sbaij->lvec) {
29: VecDestroy(sbaij->lvec);
30: sbaij->lvec = 0;
31: }
32: if (sbaij->Mvctx) {
33: VecScatterDestroy(sbaij->Mvctx);
34: sbaij->Mvctx = 0;
35: }
37: /* For the first stab we make an array as long as the number of columns */
38: /* mark those columns that are in sbaij->B */
39: PetscMalloc((Nbs+1)*sizeof(int),&indices);
40: PetscMemzero(indices,Nbs*sizeof(int));
41: for (i=0; i<mbs; i++) {
42: for (j=0; j<B->ilen[i]; j++) {
43: if (!indices[aj[B->i[i] + j]]) ec++;
44: indices[aj[B->i[i] + j] ] = 1;
45: }
46: }
48: /* form arrays of columns we need */
49: PetscMalloc((ec+1)*sizeof(int),&garray);
50: PetscMalloc((3*ec+1)*sizeof(int),&sgarray);
51: ec_owner = sgarray + 2*ec;
52:
53: ec = 0;
54: for (j=0; j<size; j++){
55: for (i=owners[j]; i<owners[j+1]; i++){
56: if (indices[i]) {
57: garray[ec] = i;
58: ec_owner[ec] = j;
59: ec++;
60: }
61: }
62: }
64: /* make indices now point into garray */
65: for (i=0; i<ec; i++) indices[garray[i]] = i;
67: /* compact out the extra columns in B */
68: for (i=0; i<mbs; i++) {
69: for (j=0; j<B->ilen[i]; j++) aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
70: }
71: B->nbs = ec;
72: sbaij->B->n = ec*B->bs;
73: PetscFree(indices);
75: /* create local vector that is used to scatter into */
76: VecCreateSeq(PETSC_COMM_SELF,ec*bs,&sbaij->lvec);
78: /* create two temporary index sets for building scatter-gather */
79: PetscMalloc((2*ec+1)*sizeof(int),&stmp);
80: for (i=0; i<ec; i++) stmp[i] = bs*garray[i];
81: ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,&from);
82:
83: for (i=0; i<ec; i++) { stmp[i] = bs*i; }
84: ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,&to);
86: /* generate the scatter context */
87: VecCreateMPI(mat->comm,mat->n,mat->N,&gvec);
88: VecScatterCreate(gvec,from,sbaij->lvec,to,&sbaij->Mvctx);
89: VecScatterPostRecvs(gvec,sbaij->lvec,INSERT_VALUES,SCATTER_FORWARD,sbaij->Mvctx);
91: sbaij->garray = garray;
92: PetscLogObjectParent(mat,sbaij->Mvctx);
93: PetscLogObjectParent(mat,sbaij->lvec);
94: PetscLogObjectParent(mat,from);
95: PetscLogObjectParent(mat,to);
97: ISDestroy(from);
98: ISDestroy(to);
100: /* create parallel vector that is used by SBAIJ matrix to scatter from/into */
101: lsize = (mbs + ec)*bs;
102: VecCreateMPI(mat->comm,lsize,PETSC_DETERMINE,&sbaij->slvec0);
103: VecDuplicate(sbaij->slvec0,&sbaij->slvec1);
104: VecGetSize(sbaij->slvec0,&vec_size);
106: VecGetPetscMap(sbaij->slvec0,&vecmap);
107: PetscMapGetGlobalRange(vecmap,&sowners);
108:
109: /* x index in the IS sfrom */
110: for (i=0; i<ec; i++) {
111: j = ec_owner[i];
112: sgarray[i] = garray[i] + (sowners[j]/bs - owners[j]);
113: }
114: /* b index in the IS sfrom */
115: k = sowners[rank]/bs + mbs;
116: for (i=ec,j=0; i< 2*ec; i++,j++) sgarray[i] = k + j;
117:
118: for (i=0; i<2*ec; i++) stmp[i] = bs*sgarray[i];
119: ISCreateBlock(PETSC_COMM_SELF,bs,2*ec,stmp,&from);
120:
121: /* x index in the IS sto */
122: k = sowners[rank]/bs + mbs;
123: for (i=0; i<ec; i++) stmp[i] = bs*(k + i);
124: /* b index in the IS sto */
125: for (i=ec; i<2*ec; i++) stmp[i] = bs*sgarray[i-ec];
127: ISCreateBlock(PETSC_COMM_SELF,bs,2*ec,stmp,&to);
129: /* gnerate the SBAIJ scatter context */
130: VecScatterCreate(sbaij->slvec0,from,sbaij->slvec1,to,&sbaij->sMvctx);
131:
132: /*
133: Post the receives for the first matrix vector product. We sync-chronize after
134: this on the chance that the user immediately calls MatMult() after assemblying
135: the matrix.
136: */
137: VecScatterPostRecvs(sbaij->slvec0,sbaij->slvec1,INSERT_VALUES,SCATTER_FORWARD,sbaij->sMvctx);
139: VecGetLocalSize(sbaij->slvec1,&nt);
140: VecGetArray(sbaij->slvec1,&ptr);
141: VecCreateSeqWithArray(PETSC_COMM_SELF,bs*mbs,ptr,&sbaij->slvec1a);
142: VecCreateSeqWithArray(PETSC_COMM_SELF,nt-bs*mbs,ptr+bs*mbs,&sbaij->slvec1b);
143: VecRestoreArray(sbaij->slvec1,&ptr);
145: VecGetArray(sbaij->slvec0,&ptr);
146: VecCreateSeqWithArray(PETSC_COMM_SELF,nt-bs*mbs,ptr+bs*mbs,&sbaij->slvec0b);
147: VecRestoreArray(sbaij->slvec0,&ptr);
149: PetscFree(stmp);
150: MPI_Barrier(mat->comm);
151:
152: PetscLogObjectParent(mat,sbaij->sMvctx);
153: PetscLogObjectParent(mat,sbaij->slvec0);
154: PetscLogObjectParent(mat,sbaij->slvec1);
155: PetscLogObjectParent(mat,sbaij->slvec0b);
156: PetscLogObjectParent(mat,sbaij->slvec1a);
157: PetscLogObjectParent(mat,sbaij->slvec1b);
158: PetscLogObjectParent(mat,from);
159: PetscLogObjectParent(mat,to);
160:
161: PetscLogObjectMemory(mat,(ec+1)*sizeof(int));
162: ISDestroy(from);
163: ISDestroy(to);
164: VecDestroy(gvec);
165: PetscFree(sgarray);
167: return(0);
168: }
172: int MatSetUpMultiply_MPISBAIJ_2comm(Mat mat)
173: {
174: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
175: Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)(baij->B->data);
176: int Nbs = baij->Nbs,i,j,*indices,*aj = B->j,ierr,ec = 0,*garray;
177: int bs = baij->bs,*stmp;
178: IS from,to;
179: Vec gvec;
180: #if defined (PETSC_USE_CTABLE)
181: PetscTable gid1_lid1;
182: PetscTablePosition tpos;
183: int gid,lid;
184: #endif
188: #if defined (PETSC_USE_CTABLE)
189: /* use a table - Mark Adams */
190: PetscTableCreate(B->mbs,&gid1_lid1);
191: for (i=0; i<B->mbs; i++) {
192: for (j=0; j<B->ilen[i]; j++) {
193: int data,gid1 = aj[B->i[i]+j] + 1;
194: PetscTableFind(gid1_lid1,gid1,&data) ;
195: if (!data) {
196: /* one based table */
197: PetscTableAdd(gid1_lid1,gid1,++ec);
198: }
199: }
200: }
201: /* form array of columns we need */
202: PetscMalloc((ec+1)*sizeof(int),&garray);
203: PetscTableGetHeadPosition(gid1_lid1,&tpos);
204: while (tpos) {
205: PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);
206: gid--; lid--;
207: garray[lid] = gid;
208: }
209: PetscSortInt(ec,garray);
210: /* qsort(garray, ec, sizeof(int), intcomparcarc); */
211: PetscTableRemoveAll(gid1_lid1);
212: for (i=0; i<ec; i++) {
213: PetscTableAdd(gid1_lid1,garray[i]+1,i+1);
214: }
215: /* compact out the extra columns in B */
216: for (i=0; i<B->mbs; i++) {
217: for (j=0; j<B->ilen[i]; j++) {
218: int gid1 = aj[B->i[i] + j] + 1;
219: PetscTableFind(gid1_lid1,gid1,&lid);
220: lid --;
221: aj[B->i[i]+j] = lid;
222: }
223: }
224: B->nbs = ec;
225: baij->B->n = ec*B->bs;
226: PetscTableDelete(gid1_lid1);
227: /* Mark Adams */
228: #else
229: /* For the first stab we make an array as long as the number of columns */
230: /* mark those columns that are in baij->B */
231: PetscMalloc((Nbs+1)*sizeof(int),&indices);
232: PetscMemzero(indices,Nbs*sizeof(int));
233: for (i=0; i<B->mbs; i++) {
234: for (j=0; j<B->ilen[i]; j++) {
235: if (!indices[aj[B->i[i] + j]]) ec++;
236: indices[aj[B->i[i] + j] ] = 1;
237: }
238: }
240: /* form array of columns we need */
241: PetscMalloc((ec+1)*sizeof(int),&garray);
242: ec = 0;
243: for (i=0; i<Nbs; i++) {
244: if (indices[i]) {
245: garray[ec++] = i;
246: }
247: }
249: /* make indices now point into garray */
250: for (i=0; i<ec; i++) {
251: indices[garray[i]] = i;
252: }
254: /* compact out the extra columns in B */
255: for (i=0; i<B->mbs; i++) {
256: for (j=0; j<B->ilen[i]; j++) {
257: aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
258: }
259: }
260: B->nbs = ec;
261: baij->B->n = ec*B->bs;
262: PetscFree(indices);
263: #endif
264:
265: /* create local vector that is used to scatter into */
266: VecCreateSeq(PETSC_COMM_SELF,ec*bs,&baij->lvec);
268: /* create two temporary index sets for building scatter-gather */
269: for (i=0; i<ec; i++) {
270: garray[i] = bs*garray[i];
271: }
272: ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,&from);
273: for (i=0; i<ec; i++) {
274: garray[i] = garray[i]/bs;
275: }
277: PetscMalloc((ec+1)*sizeof(int),&stmp);
278: for (i=0; i<ec; i++) { stmp[i] = bs*i; }
279: ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,&to);
280: PetscFree(stmp);
282: /* create temporary global vector to generate scatter context */
283: /* this is inefficient, but otherwise we must do either
284: 1) save garray until the first actual scatter when the vector is known or
285: 2) have another way of generating a scatter context without a vector.*/
286: VecCreateMPI(mat->comm,mat->n,mat->N,&gvec);
288: /* gnerate the scatter context */
289: VecScatterCreate(gvec,from,baij->lvec,to,&baij->Mvctx);
291: /*
292: Post the receives for the first matrix vector product. We sync-chronize after
293: this on the chance that the user immediately calls MatMult() after assemblying
294: the matrix.
295: */
296: VecScatterPostRecvs(gvec,baij->lvec,INSERT_VALUES,SCATTER_FORWARD,baij->Mvctx);
297: MPI_Barrier(mat->comm);
299: PetscLogObjectParent(mat,baij->Mvctx);
300: PetscLogObjectParent(mat,baij->lvec);
301: PetscLogObjectParent(mat,from);
302: PetscLogObjectParent(mat,to);
303: baij->garray = garray;
304: PetscLogObjectMemory(mat,(ec+1)*sizeof(int));
305: ISDestroy(from);
306: ISDestroy(to);
307: VecDestroy(gvec);
308: return(0);
309: }
312: /*
313: Takes the local part of an already assembled MPIBAIJ matrix
314: and disassembles it. This is to allow new nonzeros into the matrix
315: that require more communication in the matrix vector multiply.
316: Thus certain data-structures must be rebuilt.
318: Kind of slow! But that's what application programmers get when
319: they are sloppy.
320: */
323: int DisAssemble_MPISBAIJ(Mat A)
324: {
325: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)A->data;
326: Mat B = baij->B,Bnew;
327: Mat_SeqBAIJ *Bbaij = (Mat_SeqBAIJ*)B->data;
328: int ierr,i,j,mbs=Bbaij->mbs,n = A->N,col,*garray=baij->garray;
329: int k,bs=baij->bs,bs2=baij->bs2,*rvals,*nz,ec,m=A->m;
330: MatScalar *a = Bbaij->a;
331: PetscScalar *atmp;
332: #if defined(PETSC_USE_MAT_SINGLE)
333: int l;
334: #endif
337: #if defined(PETSC_USE_MAT_SINGLE)
338: PetscMalloc(baij->bs*sizeof(PetscScalar),&atmp);
339: #endif
340: /* free stuff related to matrix-vec multiply */
341: VecGetSize(baij->lvec,&ec); /* needed for PetscLogObjectMemory below */
342: VecDestroy(baij->lvec); baij->lvec = 0;
343: VecScatterDestroy(baij->Mvctx); baij->Mvctx = 0;
344: if (baij->colmap) {
345: #if defined (PETSC_USE_CTABLE)
346: PetscTableDelete(baij->colmap); baij->colmap = 0;
347: #else
348: PetscFree(baij->colmap);
349: baij->colmap = 0;
350: PetscLogObjectMemory(A,-Bbaij->nbs*sizeof(int));
351: #endif
352: }
354: /* make sure that B is assembled so we can access its values */
355: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
356: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
358: /* invent new B and copy stuff over */
359: PetscMalloc(mbs*sizeof(int),&nz);
360: for (i=0; i<mbs; i++) {
361: nz[i] = Bbaij->i[i+1]-Bbaij->i[i];
362: }
363: MatCreateSeqBAIJ(PETSC_COMM_SELF,baij->bs,m,n,0,nz,&Bnew);
364: PetscFree(nz);
365:
366: PetscMalloc(bs*sizeof(int),&rvals);
367: for (i=0; i<mbs; i++) {
368: rvals[0] = bs*i;
369: for (j=1; j<bs; j++) { rvals[j] = rvals[j-1] + 1; }
370: for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) {
371: col = garray[Bbaij->j[j]]*bs;
372: for (k=0; k<bs; k++) {
373: #if defined(PETSC_USE_MAT_SINGLE)
374: for (l=0; l<bs; l++) atmp[l] = a[j*bs2+l];
375: #else
376: atmp = a+j*bs2 + k*bs;
377: #endif
378: MatSetValues_SeqSBAIJ(Bnew,bs,rvals,1,&col,atmp,B->insertmode);
379: col++;
380: }
381: }
382: }
383: #if defined(PETSC_USE_MAT_SINGLE)
384: PetscFree(atmp);
385: #endif
386: PetscFree(baij->garray);
387: baij->garray = 0;
388: PetscFree(rvals);
389: PetscLogObjectMemory(A,-ec*sizeof(int));
390: MatDestroy(B);
391: PetscLogObjectParent(A,Bnew);
392: baij->B = Bnew;
393: A->was_assembled = PETSC_FALSE;
394: return(0);
395: }