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: PetscMapSetUp((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: VecCreateMPIWithArray(((PetscObject)mat)->comm,mat->cmap->n,mat->cmap->N,PETSC_NULL,&gvec);
86: VecScatterCreate(gvec,from,sbaij->lvec,to,&sbaij->Mvctx);
87: VecDestroy(gvec);
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(((PetscObject)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: VecScatterCreate(sbaij->slvec0,from,sbaij->slvec1,to,&sbaij->sMvctx);
127:
128: VecGetLocalSize(sbaij->slvec1,&nt);
129: VecGetArray(sbaij->slvec1,&ptr);
130: VecCreateSeqWithArray(PETSC_COMM_SELF,bs*mbs,ptr,&sbaij->slvec1a);
131: VecCreateSeqWithArray(PETSC_COMM_SELF,nt-bs*mbs,ptr+bs*mbs,&sbaij->slvec1b);
132: VecRestoreArray(sbaij->slvec1,&ptr);
134: VecGetArray(sbaij->slvec0,&ptr);
135: VecCreateSeqWithArray(PETSC_COMM_SELF,nt-bs*mbs,ptr+bs*mbs,&sbaij->slvec0b);
136: VecRestoreArray(sbaij->slvec0,&ptr);
138: PetscFree(stmp);
139: MPI_Barrier(((PetscObject)mat)->comm);
140:
141: PetscLogObjectParent(mat,sbaij->sMvctx);
142: PetscLogObjectParent(mat,sbaij->slvec0);
143: PetscLogObjectParent(mat,sbaij->slvec1);
144: PetscLogObjectParent(mat,sbaij->slvec0b);
145: PetscLogObjectParent(mat,sbaij->slvec1a);
146: PetscLogObjectParent(mat,sbaij->slvec1b);
147: PetscLogObjectParent(mat,from);
148: PetscLogObjectParent(mat,to);
149:
150: PetscLogObjectMemory(mat,(ec+1)*sizeof(PetscInt));
151: ISDestroy(from);
152: ISDestroy(to);
153: PetscFree(sgarray);
154: return(0);
155: }
159: PetscErrorCode MatSetUpMultiply_MPISBAIJ_2comm(Mat mat)
160: {
161: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
162: Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)(baij->B->data);
163: PetscErrorCode ierr;
164: PetscInt i,j,*aj = B->j,ec = 0,*garray;
165: PetscInt bs = mat->rmap->bs,*stmp;
166: IS from,to;
167: Vec gvec;
168: #if defined (PETSC_USE_CTABLE)
169: PetscTable gid1_lid1;
170: PetscTablePosition tpos;
171: PetscInt gid,lid;
172: #else
173: PetscInt Nbs = baij->Nbs,*indices;
174: #endif
177: #if defined (PETSC_USE_CTABLE)
178: /* use a table - Mark Adams */
179: PetscTableCreate(B->mbs,&gid1_lid1);
180: for (i=0; i<B->mbs; i++) {
181: for (j=0; j<B->ilen[i]; j++) {
182: PetscInt data,gid1 = aj[B->i[i]+j] + 1;
183: PetscTableFind(gid1_lid1,gid1,&data);
184: if (!data) {
185: /* one based table */
186: PetscTableAdd(gid1_lid1,gid1,++ec);
187: }
188: }
189: }
190: /* form array of columns we need */
191: PetscMalloc((ec+1)*sizeof(PetscInt),&garray);
192: PetscTableGetHeadPosition(gid1_lid1,&tpos);
193: while (tpos) {
194: PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);
195: gid--; lid--;
196: garray[lid] = gid;
197: }
198: PetscSortInt(ec,garray);
199: PetscTableRemoveAll(gid1_lid1);
200: for (i=0; i<ec; i++) {
201: PetscTableAdd(gid1_lid1,garray[i]+1,i+1);
202: }
203: /* compact out the extra columns in B */
204: for (i=0; i<B->mbs; i++) {
205: for (j=0; j<B->ilen[i]; j++) {
206: PetscInt gid1 = aj[B->i[i] + j] + 1;
207: PetscTableFind(gid1_lid1,gid1,&lid);
208: lid --;
209: aj[B->i[i]+j] = lid;
210: }
211: }
212: B->nbs = ec;
213: baij->B->cmap->n = baij->B->cmap->N = ec*mat->rmap->bs;
214: PetscMapSetUp((baij->B->cmap));
215: PetscTableDestroy(gid1_lid1);
216: /* Mark Adams */
217: #else
218: /* For the first stab we make an array as long as the number of columns */
219: /* mark those columns that are in baij->B */
220: PetscMalloc((Nbs+1)*sizeof(PetscInt),&indices);
221: PetscMemzero(indices,Nbs*sizeof(PetscInt));
222: for (i=0; i<B->mbs; i++) {
223: for (j=0; j<B->ilen[i]; j++) {
224: if (!indices[aj[B->i[i] + j]]) ec++;
225: indices[aj[B->i[i] + j] ] = 1;
226: }
227: }
229: /* form array of columns we need */
230: PetscMalloc((ec+1)*sizeof(PetscInt),&garray);
231: ec = 0;
232: for (i=0; i<Nbs; i++) {
233: if (indices[i]) {
234: garray[ec++] = i;
235: }
236: }
238: /* make indices now point into garray */
239: for (i=0; i<ec; i++) {
240: indices[garray[i]] = i;
241: }
243: /* compact out the extra columns in B */
244: for (i=0; i<B->mbs; i++) {
245: for (j=0; j<B->ilen[i]; j++) {
246: aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
247: }
248: }
249: B->nbs = ec;
250: baij->B->cmap->n = ec*mat->rmap->bs;
251: PetscFree(indices);
252: #endif
253:
254: /* create local vector that is used to scatter into */
255: VecCreateSeq(PETSC_COMM_SELF,ec*bs,&baij->lvec);
257: /* create two temporary index sets for building scatter-gather */
258: for (i=0; i<ec; i++) {
259: garray[i] = bs*garray[i];
260: }
261: ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,&from);
262: for (i=0; i<ec; i++) {
263: garray[i] = garray[i]/bs;
264: }
266: PetscMalloc((ec+1)*sizeof(PetscInt),&stmp);
267: for (i=0; i<ec; i++) { stmp[i] = bs*i; }
268: ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,&to);
269: PetscFree(stmp);
271: /* create temporary global vector to generate scatter context */
272: VecCreateMPIWithArray(((PetscObject)mat)->comm,mat->cmap->n,mat->cmap->N,PETSC_NULL,&gvec);
273: VecScatterCreate(gvec,from,baij->lvec,to,&baij->Mvctx);
274: VecDestroy(gvec);
276: PetscLogObjectParent(mat,baij->Mvctx);
277: PetscLogObjectParent(mat,baij->lvec);
278: PetscLogObjectParent(mat,from);
279: PetscLogObjectParent(mat,to);
280: baij->garray = garray;
281: PetscLogObjectMemory(mat,(ec+1)*sizeof(PetscInt));
282: ISDestroy(from);
283: ISDestroy(to);
284: return(0);
285: }
287: /*
288: Takes the local part of an already assembled MPISBAIJ matrix
289: and disassembles it. This is to allow new nonzeros into the matrix
290: that require more communication in the matrix vector multiply.
291: Thus certain data-structures must be rebuilt.
293: Kind of slow! But that's what application programmers get when
294: they are sloppy.
295: */
298: PetscErrorCode DisAssemble_MPISBAIJ(Mat A)
299: {
300: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)A->data;
301: Mat B = baij->B,Bnew;
302: Mat_SeqBAIJ *Bbaij = (Mat_SeqBAIJ*)B->data;
304: PetscInt i,j,mbs=Bbaij->mbs,n = A->cmap->N,col,*garray=baij->garray;
305: PetscInt k,bs=A->rmap->bs,bs2=baij->bs2,*rvals,*nz,ec,m=A->rmap->n;
306: MatScalar *a = Bbaij->a;
307: PetscScalar *atmp;
308: #if defined(PETSC_USE_MAT_SINGLE)
309: PetscInt l;
310: #endif
313: #if defined(PETSC_USE_MAT_SINGLE)
314: PetscMalloc(A->rmap->bs*sizeof(PetscScalar),&atmp);
315: #endif
316: /* free stuff related to matrix-vec multiply */
317: VecGetSize(baij->lvec,&ec); /* needed for PetscLogObjectMemory below */
318: VecDestroy(baij->lvec);
319: baij->lvec = 0;
320: VecScatterDestroy(baij->Mvctx);
321: baij->Mvctx = 0;
323: VecDestroy(baij->slvec0);
324: VecDestroy(baij->slvec0b);
325: baij->slvec0 = 0;
326: VecDestroy(baij->slvec1);
327: VecDestroy(baij->slvec1a);
328: VecDestroy(baij->slvec1b);
329: baij->slvec1 = 0;
331: if (baij->colmap) {
332: #if defined (PETSC_USE_CTABLE)
333: PetscTableDestroy(baij->colmap); baij->colmap = 0;
334: #else
335: PetscFree(baij->colmap);
336: baij->colmap = 0;
337: PetscLogObjectMemory(A,-Bbaij->nbs*sizeof(PetscInt));
338: #endif
339: }
341: /* make sure that B is assembled so we can access its values */
342: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
343: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
345: /* invent new B and copy stuff over */
346: PetscMalloc(mbs*sizeof(PetscInt),&nz);
347: for (i=0; i<mbs; i++) {
348: nz[i] = Bbaij->i[i+1]-Bbaij->i[i];
349: }
350: MatCreate(PETSC_COMM_SELF,&Bnew);
351: MatSetSizes(Bnew,m,n,m,n);
352: MatSetType(Bnew,((PetscObject)B)->type_name);
353: MatSeqBAIJSetPreallocation(Bnew,B->rmap->bs,0,nz);
354: PetscFree(nz);
355:
356: PetscMalloc(bs*sizeof(PetscInt),&rvals);
357: for (i=0; i<mbs; i++) {
358: rvals[0] = bs*i;
359: for (j=1; j<bs; j++) { rvals[j] = rvals[j-1] + 1; }
360: for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) {
361: col = garray[Bbaij->j[j]]*bs;
362: for (k=0; k<bs; k++) {
363: #if defined(PETSC_USE_MAT_SINGLE)
364: for (l=0; l<bs; l++) atmp[l] = a[j*bs2+l];
365: #else
366: atmp = a+j*bs2 + k*bs;
367: #endif
368: MatSetValues_SeqSBAIJ(Bnew,bs,rvals,1,&col,atmp,B->insertmode);
369: col++;
370: }
371: }
372: }
373: #if defined(PETSC_USE_MAT_SINGLE)
374: PetscFree(atmp);
375: #endif
376: PetscFree(baij->garray);
377: baij->garray = 0;
378: PetscFree(rvals);
379: PetscLogObjectMemory(A,-ec*sizeof(PetscInt));
380: MatDestroy(B);
381: PetscLogObjectParent(A,Bnew);
382: baij->B = Bnew;
383: A->was_assembled = PETSC_FALSE;
384: return(0);
385: }