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
8: extern int MatSetValues_SeqSBAIJ(Mat,int,const int [],int,const int [],const PetscScalar [],InsertMode);
13: int MatSetUpMultiply_MPISBAIJ(Mat mat)
14: {
15: Mat_MPISBAIJ *sbaij = (Mat_MPISBAIJ*)mat->data;
16: Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)(sbaij->B->data);
17: int Nbs = sbaij->Nbs,i,j,*indices,*aj = B->j,ierr,ec = 0,*garray,*sgarray;
18: int bs = sbaij->bs,*stmp,mbs=sbaij->mbs, vec_size,nt;
19: IS from,to;
20: Vec gvec;
21: int rank=sbaij->rank,lsize,size=sbaij->size;
22: int *owners=sbaij->rowners,*sowners,*ec_owner,k;
23: PetscMap vecmap;
24: PetscScalar *ptr;
27: if (sbaij->lvec) {
28: VecDestroy(sbaij->lvec);
29: sbaij->lvec = 0;
30: }
31: if (sbaij->Mvctx) {
32: VecScatterDestroy(sbaij->Mvctx);
33: sbaij->Mvctx = 0;
34: }
36: /* For the first stab we make an array as long as the number of columns */
37: /* mark those columns that are in sbaij->B */
38: PetscMalloc((Nbs+1)*sizeof(int),&indices);
39: PetscMemzero(indices,Nbs*sizeof(int));
40: for (i=0; i<mbs; i++) {
41: for (j=0; j<B->ilen[i]; j++) {
42: if (!indices[aj[B->i[i] + j]]) ec++;
43: indices[aj[B->i[i] + j] ] = 1;
44: }
45: }
47: /* form arrays of columns we need */
48: PetscMalloc((ec+1)*sizeof(int),&garray);
49: PetscMalloc((3*ec+1)*sizeof(int),&sgarray);
50: ec_owner = sgarray + 2*ec;
51:
52: ec = 0;
53: for (j=0; j<size; j++){
54: for (i=owners[j]; i<owners[j+1]; i++){
55: if (indices[i]) {
56: garray[ec] = i;
57: ec_owner[ec] = j;
58: ec++;
59: }
60: }
61: }
63: /* make indices now point into garray */
64: for (i=0; i<ec; i++) indices[garray[i]] = i;
66: /* compact out the extra columns in B */
67: for (i=0; i<mbs; i++) {
68: for (j=0; j<B->ilen[i]; j++) aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
69: }
70: B->nbs = ec;
71: sbaij->B->n = ec*B->bs;
72: PetscFree(indices);
74: /* create local vector that is used to scatter into */
75: VecCreateSeq(PETSC_COMM_SELF,ec*bs,&sbaij->lvec);
77: /* create two temporary index sets for building scatter-gather */
78: PetscMalloc((2*ec+1)*sizeof(int),&stmp);
79: for (i=0; i<ec; i++) stmp[i] = bs*garray[i];
80: ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,&from);
81:
82: for (i=0; i<ec; i++) { stmp[i] = bs*i; }
83: ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,&to);
85: /* generate the scatter context
86: -- Mvctx and lvec are not used by MatMult_MPISBAIJ(), but usefule for some applications */
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 i,j,*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: #else
185: int Nbs = baij->Nbs,*indices;
186: #endif
190: #if defined (PETSC_USE_CTABLE)
191: /* use a table - Mark Adams */
192: PetscTableCreate(B->mbs,&gid1_lid1);
193: for (i=0; i<B->mbs; i++) {
194: for (j=0; j<B->ilen[i]; j++) {
195: int data,gid1 = aj[B->i[i]+j] + 1;
196: PetscTableFind(gid1_lid1,gid1,&data) ;
197: if (!data) {
198: /* one based table */
199: PetscTableAdd(gid1_lid1,gid1,++ec);
200: }
201: }
202: }
203: /* form array of columns we need */
204: PetscMalloc((ec+1)*sizeof(int),&garray);
205: PetscTableGetHeadPosition(gid1_lid1,&tpos);
206: while (tpos) {
207: PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);
208: gid--; lid--;
209: garray[lid] = gid;
210: }
211: PetscSortInt(ec,garray);
212: PetscTableRemoveAll(gid1_lid1);
213: for (i=0; i<ec; i++) {
214: PetscTableAdd(gid1_lid1,garray[i]+1,i+1);
215: }
216: /* compact out the extra columns in B */
217: for (i=0; i<B->mbs; i++) {
218: for (j=0; j<B->ilen[i]; j++) {
219: int gid1 = aj[B->i[i] + j] + 1;
220: PetscTableFind(gid1_lid1,gid1,&lid);
221: lid --;
222: aj[B->i[i]+j] = lid;
223: }
224: }
225: B->nbs = ec;
226: baij->B->n = ec*B->bs;
227: PetscTableDelete(gid1_lid1);
228: /* Mark Adams */
229: #else
230: /* For the first stab we make an array as long as the number of columns */
231: /* mark those columns that are in baij->B */
232: PetscMalloc((Nbs+1)*sizeof(int),&indices);
233: PetscMemzero(indices,Nbs*sizeof(int));
234: for (i=0; i<B->mbs; i++) {
235: for (j=0; j<B->ilen[i]; j++) {
236: if (!indices[aj[B->i[i] + j]]) ec++;
237: indices[aj[B->i[i] + j] ] = 1;
238: }
239: }
241: /* form array of columns we need */
242: PetscMalloc((ec+1)*sizeof(int),&garray);
243: ec = 0;
244: for (i=0; i<Nbs; i++) {
245: if (indices[i]) {
246: garray[ec++] = i;
247: }
248: }
250: /* make indices now point into garray */
251: for (i=0; i<ec; i++) {
252: indices[garray[i]] = i;
253: }
255: /* compact out the extra columns in B */
256: for (i=0; i<B->mbs; i++) {
257: for (j=0; j<B->ilen[i]; j++) {
258: aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
259: }
260: }
261: B->nbs = ec;
262: baij->B->n = ec*B->bs;
263: PetscFree(indices);
264: #endif
265:
266: /* create local vector that is used to scatter into */
267: VecCreateSeq(PETSC_COMM_SELF,ec*bs,&baij->lvec);
269: /* create two temporary index sets for building scatter-gather */
270: for (i=0; i<ec; i++) {
271: garray[i] = bs*garray[i];
272: }
273: ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,&from);
274: for (i=0; i<ec; i++) {
275: garray[i] = garray[i]/bs;
276: }
278: PetscMalloc((ec+1)*sizeof(int),&stmp);
279: for (i=0; i<ec; i++) { stmp[i] = bs*i; }
280: ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,&to);
281: PetscFree(stmp);
283: /* create temporary global vector to generate scatter context */
284: /* this is inefficient, but otherwise we must do either
285: 1) save garray until the first actual scatter when the vector is known or
286: 2) have another way of generating a scatter context without a vector.*/
287: VecCreateMPI(mat->comm,mat->n,mat->N,&gvec);
289: /* gnerate the scatter context */
290: VecScatterCreate(gvec,from,baij->lvec,to,&baij->Mvctx);
292: /*
293: Post the receives for the first matrix vector product. We sync-chronize after
294: this on the chance that the user immediately calls MatMult() after assemblying
295: the matrix.
296: */
297: VecScatterPostRecvs(gvec,baij->lvec,INSERT_VALUES,SCATTER_FORWARD,baij->Mvctx);
298: MPI_Barrier(mat->comm);
300: PetscLogObjectParent(mat,baij->Mvctx);
301: PetscLogObjectParent(mat,baij->lvec);
302: PetscLogObjectParent(mat,from);
303: PetscLogObjectParent(mat,to);
304: baij->garray = garray;
305: PetscLogObjectMemory(mat,(ec+1)*sizeof(int));
306: ISDestroy(from);
307: ISDestroy(to);
308: VecDestroy(gvec);
309: return(0);
310: }
313: /*
314: Takes the local part of an already assembled MPIBAIJ matrix
315: and disassembles it. This is to allow new nonzeros into the matrix
316: that require more communication in the matrix vector multiply.
317: Thus certain data-structures must be rebuilt.
319: Kind of slow! But that's what application programmers get when
320: they are sloppy.
321: */
324: int DisAssemble_MPISBAIJ(Mat A)
325: {
326: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)A->data;
327: Mat B = baij->B,Bnew;
328: Mat_SeqBAIJ *Bbaij = (Mat_SeqBAIJ*)B->data;
329: int ierr,i,j,mbs=Bbaij->mbs,n = A->N,col,*garray=baij->garray;
330: int k,bs=baij->bs,bs2=baij->bs2,*rvals,*nz,ec,m=A->m;
331: MatScalar *a = Bbaij->a;
332: PetscScalar *atmp;
333: #if defined(PETSC_USE_MAT_SINGLE)
334: int l;
335: #endif
338: #if defined(PETSC_USE_MAT_SINGLE)
339: PetscMalloc(baij->bs*sizeof(PetscScalar),&atmp);
340: #endif
341: /* free stuff related to matrix-vec multiply */
342: VecGetSize(baij->lvec,&ec); /* needed for PetscLogObjectMemory below */
343: VecDestroy(baij->lvec); baij->lvec = 0;
344: VecScatterDestroy(baij->Mvctx); baij->Mvctx = 0;
345: if (baij->colmap) {
346: #if defined (PETSC_USE_CTABLE)
347: PetscTableDelete(baij->colmap); baij->colmap = 0;
348: #else
349: PetscFree(baij->colmap);
350: baij->colmap = 0;
351: PetscLogObjectMemory(A,-Bbaij->nbs*sizeof(int));
352: #endif
353: }
355: /* make sure that B is assembled so we can access its values */
356: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
357: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
359: /* invent new B and copy stuff over */
360: PetscMalloc(mbs*sizeof(int),&nz);
361: for (i=0; i<mbs; i++) {
362: nz[i] = Bbaij->i[i+1]-Bbaij->i[i];
363: }
364: MatCreate(PETSC_COMM_SELF,m,n,m,n,&Bnew);
365: MatSetType(Bnew,B->type_name);
366: MatSeqBAIJSetPreallocation(Bnew,baij->bs,0,nz);
367: PetscFree(nz);
368:
369: PetscMalloc(bs*sizeof(int),&rvals);
370: for (i=0; i<mbs; i++) {
371: rvals[0] = bs*i;
372: for (j=1; j<bs; j++) { rvals[j] = rvals[j-1] + 1; }
373: for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) {
374: col = garray[Bbaij->j[j]]*bs;
375: for (k=0; k<bs; k++) {
376: #if defined(PETSC_USE_MAT_SINGLE)
377: for (l=0; l<bs; l++) atmp[l] = a[j*bs2+l];
378: #else
379: atmp = a+j*bs2 + k*bs;
380: #endif
381: MatSetValues_SeqSBAIJ(Bnew,bs,rvals,1,&col,atmp,B->insertmode);
382: col++;
383: }
384: }
385: }
386: #if defined(PETSC_USE_MAT_SINGLE)
387: PetscFree(atmp);
388: #endif
389: PetscFree(baij->garray);
390: baij->garray = 0;
391: PetscFree(rvals);
392: PetscLogObjectMemory(A,-ec*sizeof(int));
393: MatDestroy(B);
394: PetscLogObjectParent(A,Bnew);
395: baij->B = Bnew;
396: A->was_assembled = PETSC_FALSE;
397: return(0);
398: }