Actual source code: mmsbaij.c

  2: /*
  3:    Support for the parallel SBAIJ matrix vector multiply
  4: */
 5:  #include src/mat/impls/sbaij/mpi/mpisbaij.h

  7: EXTERN PetscErrorCode MatSetValues_SeqSBAIJ(Mat,PetscInt,const PetscInt [],PetscInt,const PetscInt [],const PetscScalar [],InsertMode);


 12: PetscErrorCode MatSetUpMultiply_MPISBAIJ(Mat mat)
 13: {
 14:   Mat_MPISBAIJ   *sbaij = (Mat_MPISBAIJ*)mat->data;
 15:   Mat_SeqBAIJ    *B = (Mat_SeqBAIJ*)(sbaij->B->data);
 17:   PetscInt       Nbs = sbaij->Nbs,i,j,*indices,*aj = B->j,ec = 0,*garray,*sgarray;
 18:   PetscInt       bs = mat->bs,*stmp,mbs=sbaij->mbs, vec_size,nt;
 19:   IS             from,to;
 20:   Vec            gvec;
 21:   PetscMPIInt    rank=sbaij->rank,lsize,size=sbaij->size;
 22:   PetscInt       *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(PetscInt),&indices);
 39:   PetscMemzero(indices,Nbs*sizeof(PetscInt));
 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(PetscInt),&garray);
 49:   PetscMalloc((3*ec+1)*sizeof(PetscInt),&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*mat->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(PetscInt),&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(PetscInt));
162:   ISDestroy(from);
163:   ISDestroy(to);
164:   VecDestroy(gvec);
165:   PetscFree(sgarray);

167:   return(0);
168: }

172: PetscErrorCode MatSetUpMultiply_MPISBAIJ_2comm(Mat mat)
173: {
174:   Mat_MPISBAIJ       *baij = (Mat_MPISBAIJ*)mat->data;
175:   Mat_SeqBAIJ        *B = (Mat_SeqBAIJ*)(baij->B->data);
176:   PetscErrorCode     ierr;
177:   PetscInt           i,j,*aj = B->j,ec = 0,*garray;
178:   PetscInt           bs = mat->bs,*stmp;
179:   IS                 from,to;
180:   Vec                gvec;
181: #if defined (PETSC_USE_CTABLE)
182:   PetscTable         gid1_lid1;
183:   PetscTablePosition tpos;
184:   PetscInt           gid,lid;
185: #else
186:   PetscInt           Nbs = baij->Nbs,*indices;
187: #endif  


191: #if defined (PETSC_USE_CTABLE)
192:   /* use a table - Mark Adams */
193:   PetscTableCreate(B->mbs,&gid1_lid1);
194:   for (i=0; i<B->mbs; i++) {
195:     for (j=0; j<B->ilen[i]; j++) {
196:       PetscInt data,gid1 = aj[B->i[i]+j] + 1;
197:       PetscTableFind(gid1_lid1,gid1,&data);
198:       if (!data) {
199:         /* one based table */
200:         PetscTableAdd(gid1_lid1,gid1,++ec);
201:       }
202:     }
203:   }
204:   /* form array of columns we need */
205:   PetscMalloc((ec+1)*sizeof(PetscInt),&garray);
206:   PetscTableGetHeadPosition(gid1_lid1,&tpos);
207:   while (tpos) {
208:     PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);
209:     gid--; lid--;
210:     garray[lid] = gid;
211:   }
212:   PetscSortInt(ec,garray);
213:   PetscTableRemoveAll(gid1_lid1);
214:   for (i=0; i<ec; i++) {
215:     PetscTableAdd(gid1_lid1,garray[i]+1,i+1);
216:   }
217:   /* compact out the extra columns in B */
218:   for (i=0; i<B->mbs; i++) {
219:     for (j=0; j<B->ilen[i]; j++) {
220:       PetscInt gid1 = aj[B->i[i] + j] + 1;
221:       PetscTableFind(gid1_lid1,gid1,&lid);
222:       lid --;
223:       aj[B->i[i]+j] = lid;
224:     }
225:   }
226:   B->nbs     = ec;
227:   baij->B->n = ec*mat->bs;
228:   PetscTableDelete(gid1_lid1);
229:   /* Mark Adams */
230: #else
231:   /* For the first stab we make an array as long as the number of columns */
232:   /* mark those columns that are in baij->B */
233:   PetscMalloc((Nbs+1)*sizeof(PetscInt),&indices);
234:   PetscMemzero(indices,Nbs*sizeof(PetscInt));
235:   for (i=0; i<B->mbs; i++) {
236:     for (j=0; j<B->ilen[i]; j++) {
237:       if (!indices[aj[B->i[i] + j]]) ec++;
238:       indices[aj[B->i[i] + j] ] = 1;
239:     }
240:   }

242:   /* form array of columns we need */
243:   PetscMalloc((ec+1)*sizeof(PetscInt),&garray);
244:   ec = 0;
245:   for (i=0; i<Nbs; i++) {
246:     if (indices[i]) {
247:       garray[ec++] = i;
248:     }
249:   }

251:   /* make indices now point into garray */
252:   for (i=0; i<ec; i++) {
253:     indices[garray[i]] = i;
254:   }

256:   /* compact out the extra columns in B */
257:   for (i=0; i<B->mbs; i++) {
258:     for (j=0; j<B->ilen[i]; j++) {
259:       aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
260:     }
261:   }
262:   B->nbs       = ec;
263:   baij->B->n   = ec*mat->bs;
264:   PetscFree(indices);
265: #endif  
266: 
267:   /* create local vector that is used to scatter into */
268:   VecCreateSeq(PETSC_COMM_SELF,ec*bs,&baij->lvec);

270:   /* create two temporary index sets for building scatter-gather */
271:   for (i=0; i<ec; i++) {
272:     garray[i] = bs*garray[i];
273:   }
274:   ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,&from);
275:   for (i=0; i<ec; i++) {
276:     garray[i] = garray[i]/bs;
277:   }

279:   PetscMalloc((ec+1)*sizeof(PetscInt),&stmp);
280:   for (i=0; i<ec; i++) { stmp[i] = bs*i; }
281:   ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,&to);
282:   PetscFree(stmp);

284:   /* create temporary global vector to generate scatter context */
285:   /* this is inefficient, but otherwise we must do either 
286:      1) save garray until the first actual scatter when the vector is known or
287:      2) have another way of generating a scatter context without a vector.*/
288:   VecCreateMPI(mat->comm,mat->n,mat->N,&gvec);

290:   /* gnerate the scatter context */
291:   VecScatterCreate(gvec,from,baij->lvec,to,&baij->Mvctx);

293:   /*
294:       Post the receives for the first matrix vector product. We sync-chronize after
295:     this on the chance that the user immediately calls MatMult() after assemblying 
296:     the matrix.
297:   */
298:   VecScatterPostRecvs(gvec,baij->lvec,INSERT_VALUES,SCATTER_FORWARD,baij->Mvctx);
299:   MPI_Barrier(mat->comm);

301:   PetscLogObjectParent(mat,baij->Mvctx);
302:   PetscLogObjectParent(mat,baij->lvec);
303:   PetscLogObjectParent(mat,from);
304:   PetscLogObjectParent(mat,to);
305:   baij->garray = garray;
306:   PetscLogObjectMemory(mat,(ec+1)*sizeof(PetscInt));
307:   ISDestroy(from);
308:   ISDestroy(to);
309:   VecDestroy(gvec);
310:   return(0);
311: }


314: /*
315:      Takes the local part of an already assembled MPIBAIJ matrix
316:    and disassembles it. This is to allow new nonzeros into the matrix
317:    that require more communication in the matrix vector multiply. 
318:    Thus certain data-structures must be rebuilt.

320:    Kind of slow! But that's what application programmers get when 
321:    they are sloppy.
322: */
325: PetscErrorCode DisAssemble_MPISBAIJ(Mat A)
326: {
327:   Mat_MPISBAIJ   *baij = (Mat_MPISBAIJ*)A->data;
328:   Mat           B = baij->B,Bnew;
329:   Mat_SeqBAIJ   *Bbaij = (Mat_SeqBAIJ*)B->data;
331:   PetscInt i,j,mbs=Bbaij->mbs,n = A->N,col,*garray=baij->garray;
332:   PetscInt           k,bs=A->bs,bs2=baij->bs2,*rvals,*nz,ec,m=A->m;
333:   MatScalar     *a = Bbaij->a;
334:   PetscScalar   *atmp;
335: #if defined(PETSC_USE_MAT_SINGLE)
336:   PetscInt           l;
337: #endif

340: #if defined(PETSC_USE_MAT_SINGLE)
341:   PetscMalloc(baij->bs*sizeof(PetscScalar),&atmp);
342: #endif
343:   /* free stuff related to matrix-vec multiply */
344:   VecGetSize(baij->lvec,&ec); /* needed for PetscLogObjectMemory below */
345:   VecDestroy(baij->lvec); baij->lvec = 0;
346:   VecScatterDestroy(baij->Mvctx); baij->Mvctx = 0;
347:   if (baij->colmap) {
348: #if defined (PETSC_USE_CTABLE)
349:     PetscTableDelete(baij->colmap); baij->colmap = 0;
350: #else
351:     PetscFree(baij->colmap);
352:     baij->colmap = 0;
353:     PetscLogObjectMemory(A,-Bbaij->nbs*sizeof(PetscInt));
354: #endif
355:   }

357:   /* make sure that B is assembled so we can access its values */
358:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
359:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

361:   /* invent new B and copy stuff over */
362:   PetscMalloc(mbs*sizeof(PetscInt),&nz);
363:   for (i=0; i<mbs; i++) {
364:     nz[i] = Bbaij->i[i+1]-Bbaij->i[i];
365:   }
366:   MatCreate(PETSC_COMM_SELF,m,n,m,n,&Bnew);
367:   MatSetType(Bnew,B->type_name);
368:   MatSeqBAIJSetPreallocation(Bnew,B->bs,0,nz);
369:   PetscFree(nz);
370: 
371:   PetscMalloc(bs*sizeof(PetscInt),&rvals);
372:   for (i=0; i<mbs; i++) {
373:     rvals[0] = bs*i;
374:     for (j=1; j<bs; j++) { rvals[j] = rvals[j-1] + 1; }
375:     for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) {
376:       col = garray[Bbaij->j[j]]*bs;
377:       for (k=0; k<bs; k++) {
378: #if defined(PETSC_USE_MAT_SINGLE)
379:         for (l=0; l<bs; l++) atmp[l] = a[j*bs2+l];
380: #else
381:         atmp = a+j*bs2 + k*bs;
382: #endif
383:         MatSetValues_SeqSBAIJ(Bnew,bs,rvals,1,&col,atmp,B->insertmode);
384:         col++;
385:       }
386:     }
387:   }
388: #if defined(PETSC_USE_MAT_SINGLE)
389:   PetscFree(atmp);
390: #endif
391:   PetscFree(baij->garray);
392:   baij->garray = 0;
393:   PetscFree(rvals);
394:   PetscLogObjectMemory(A,-ec*sizeof(PetscInt));
395:   MatDestroy(B);
396:   PetscLogObjectParent(A,Bnew);
397:   baij->B = Bnew;
398:   A->was_assembled = PETSC_FALSE;
399:   return(0);
400: }