Actual source code: mmsbaij.c

petsc-3.4.2 2013-07-02
  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->rmap->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->rangebs,*ec_owner,k;
 23:   const PetscInt *sowners;
 24:   PetscScalar    *ptr;

 27:   VecScatterDestroy(&sbaij->sMvctx);

 29:   /* For the first stab we make an array as long as the number of columns */
 30:   /* mark those columns that are in sbaij->B */
 31:   PetscMalloc(Nbs*sizeof(PetscInt),&indices);
 32:   PetscMemzero(indices,Nbs*sizeof(PetscInt));
 33:   for (i=0; i<mbs; i++) {
 34:     for (j=0; j<B->ilen[i]; j++) {
 35:       if (!indices[aj[B->i[i] + j]]) ec++;
 36:       indices[aj[B->i[i] + j]] = 1;
 37:     }
 38:   }

 40:   /* form arrays of columns we need */
 41:   PetscMalloc(ec*sizeof(PetscInt),&garray);
 42:   PetscMalloc2(2*ec,PetscInt,&sgarray,ec,PetscInt,&ec_owner);

 44:   ec = 0;
 45:   for (j=0; j<size; j++) {
 46:     for (i=owners[j]; i<owners[j+1]; i++) {
 47:       if (indices[i]) {
 48:         garray[ec]   = i;
 49:         ec_owner[ec] = j;
 50:         ec++;
 51:       }
 52:     }
 53:   }

 55:   /* make indices now point into garray */
 56:   for (i=0; i<ec; i++) indices[garray[i]] = i;

 58:   /* compact out the extra columns in B */
 59:   for (i=0; i<mbs; i++) {
 60:     for (j=0; j<B->ilen[i]; j++) aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
 61:   }
 62:   B->nbs = ec;

 64:   sbaij->B->cmap->n = sbaij->B->cmap->N = ec*mat->rmap->bs;

 66:   PetscLayoutSetUp((sbaij->B->cmap));
 67:   PetscFree(indices);

 69:   /* create local vector that is used to scatter into */
 70:   VecCreateSeq(PETSC_COMM_SELF,ec*bs,&sbaij->lvec);

 72:   /* create two temporary index sets for building scatter-gather */
 73:   PetscMalloc(2*ec*sizeof(PetscInt),&stmp);
 74:   ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,PETSC_COPY_VALUES,&from);
 75:   for (i=0; i<ec; i++) stmp[i] = i;
 76:   ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,PETSC_COPY_VALUES,&to);

 78:   /* generate the scatter context
 79:      -- Mvctx and lvec are not used by MatMult_MPISBAIJ(), but usefule for some applications */
 80:   VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&gvec);
 81:   VecScatterCreate(gvec,from,sbaij->lvec,to,&sbaij->Mvctx);
 82:   VecDestroy(&gvec);

 84:   sbaij->garray = garray;

 86:   PetscLogObjectParent(mat,sbaij->Mvctx);
 87:   PetscLogObjectParent(mat,sbaij->lvec);
 88:   PetscLogObjectParent(mat,from);
 89:   PetscLogObjectParent(mat,to);

 91:   ISDestroy(&from);
 92:   ISDestroy(&to);

 94:   /* create parallel vector that is used by SBAIJ matrix to scatter from/into */
 95:   lsize = (mbs + ec)*bs;
 96:   VecCreateMPI(PetscObjectComm((PetscObject)mat),lsize,PETSC_DETERMINE,&sbaij->slvec0);
 97:   VecDuplicate(sbaij->slvec0,&sbaij->slvec1);
 98:   VecGetSize(sbaij->slvec0,&vec_size);

100:   VecGetOwnershipRanges(sbaij->slvec0,&sowners);

102:   /* x index in the IS sfrom */
103:   for (i=0; i<ec; i++) {
104:     j          = ec_owner[i];
105:     sgarray[i] = garray[i] + (sowners[j]/bs - owners[j]);
106:   }
107:   /* b index in the IS sfrom */
108:   k = sowners[rank]/bs + mbs;
109:   for (i=ec,j=0; i< 2*ec; i++,j++) sgarray[i] = k + j;
110:   ISCreateBlock(PETSC_COMM_SELF,bs,2*ec,sgarray,PETSC_COPY_VALUES,&from);

112:   /* x index in the IS sto */
113:   k = sowners[rank]/bs + mbs;
114:   for (i=0; i<ec; i++) stmp[i] = (k + i);
115:   /* b index in the IS sto */
116:   for (i=ec; i<2*ec; i++) stmp[i] = sgarray[i-ec];

118:   ISCreateBlock(PETSC_COMM_SELF,bs,2*ec,stmp,PETSC_COPY_VALUES,&to);

120:   VecScatterCreate(sbaij->slvec0,from,sbaij->slvec1,to,&sbaij->sMvctx);

122:   VecGetLocalSize(sbaij->slvec1,&nt);
123:   VecGetArray(sbaij->slvec1,&ptr);
124:   VecCreateSeqWithArray(PETSC_COMM_SELF,1,bs*mbs,ptr,&sbaij->slvec1a);
125:   VecCreateSeqWithArray(PETSC_COMM_SELF,1,nt-bs*mbs,ptr+bs*mbs,&sbaij->slvec1b);
126:   VecRestoreArray(sbaij->slvec1,&ptr);

128:   VecGetArray(sbaij->slvec0,&ptr);
129:   VecCreateSeqWithArray(PETSC_COMM_SELF,1,nt-bs*mbs,ptr+bs*mbs,&sbaij->slvec0b);
130:   VecRestoreArray(sbaij->slvec0,&ptr);

132:   PetscFree(stmp);
133:   MPI_Barrier(PetscObjectComm((PetscObject)mat));

135:   PetscLogObjectParent(mat,sbaij->sMvctx);
136:   PetscLogObjectParent(mat,sbaij->slvec0);
137:   PetscLogObjectParent(mat,sbaij->slvec1);
138:   PetscLogObjectParent(mat,sbaij->slvec0b);
139:   PetscLogObjectParent(mat,sbaij->slvec1a);
140:   PetscLogObjectParent(mat,sbaij->slvec1b);
141:   PetscLogObjectParent(mat,from);
142:   PetscLogObjectParent(mat,to);

144:   PetscLogObjectMemory(mat,(ec+1)*sizeof(PetscInt));
145:   ISDestroy(&from);
146:   ISDestroy(&to);
147:   PetscFree2(sgarray,ec_owner);
148:   return(0);
149: }

153: PetscErrorCode MatSetUpMultiply_MPISBAIJ_2comm(Mat mat)
154: {
155:   Mat_MPISBAIJ   *baij = (Mat_MPISBAIJ*)mat->data;
156:   Mat_SeqBAIJ    *B    = (Mat_SeqBAIJ*)(baij->B->data);
158:   PetscInt       i,j,*aj = B->j,ec = 0,*garray;
159:   PetscInt       bs = mat->rmap->bs,*stmp;
160:   IS             from,to;
161:   Vec            gvec;
162: #if defined(PETSC_USE_CTABLE)
163:   PetscTable         gid1_lid1;
164:   PetscTablePosition tpos;
165:   PetscInt           gid,lid;
166: #else
167:   PetscInt Nbs = baij->Nbs,*indices;
168: #endif

171: #if defined(PETSC_USE_CTABLE)
172:   /* use a table - Mark Adams */
173:   PetscTableCreate(B->mbs,baij->Nbs+1,&gid1_lid1);
174:   for (i=0; i<B->mbs; i++) {
175:     for (j=0; j<B->ilen[i]; j++) {
176:       PetscInt data,gid1 = aj[B->i[i]+j] + 1;
177:       PetscTableFind(gid1_lid1,gid1,&data);
178:       if (!data) {
179:         /* one based table */
180:         PetscTableAdd(gid1_lid1,gid1,++ec,INSERT_VALUES);
181:       }
182:     }
183:   }
184:   /* form array of columns we need */
185:   PetscMalloc(ec*sizeof(PetscInt),&garray);
186:   PetscTableGetHeadPosition(gid1_lid1,&tpos);
187:   while (tpos) {
188:     PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);
189:     gid--; lid--;
190:     garray[lid] = gid;
191:   }
192:   PetscSortInt(ec,garray);
193:   PetscTableRemoveAll(gid1_lid1);
194:   for (i=0; i<ec; i++) {
195:     PetscTableAdd(gid1_lid1,garray[i]+1,i+1,INSERT_VALUES);
196:   }
197:   /* compact out the extra columns in B */
198:   for (i=0; i<B->mbs; i++) {
199:     for (j=0; j<B->ilen[i]; j++) {
200:       PetscInt gid1 = aj[B->i[i] + j] + 1;
201:       PetscTableFind(gid1_lid1,gid1,&lid);
202:       lid--;
203:       aj[B->i[i]+j] = lid;
204:     }
205:   }
206:   B->nbs = ec;

208:   baij->B->cmap->n = baij->B->cmap->N = ec*mat->rmap->bs;

210:   PetscLayoutSetUp((baij->B->cmap));
211:   PetscTableDestroy(&gid1_lid1);
212: #else
213:   /* For the first stab we make an array as long as the number of columns */
214:   /* mark those columns that are in baij->B */
215:   PetscMalloc(Nbs*sizeof(PetscInt),&indices);
216:   PetscMemzero(indices,Nbs*sizeof(PetscInt));
217:   for (i=0; i<B->mbs; i++) {
218:     for (j=0; j<B->ilen[i]; j++) {
219:       if (!indices[aj[B->i[i] + j]]) ec++;
220:       indices[aj[B->i[i] + j]] = 1;
221:     }
222:   }

224:   /* form array of columns we need */
225:   PetscMalloc(ec*sizeof(PetscInt),&garray);
226:   ec = 0;
227:   for (i=0; i<Nbs; i++) {
228:     if (indices[i]) {
229:       garray[ec++] = i;
230:     }
231:   }

233:   /* make indices now point into garray */
234:   for (i=0; i<ec; i++) indices[garray[i]] = i;

236:   /* compact out the extra columns in B */
237:   for (i=0; i<B->mbs; i++) {
238:     for (j=0; j<B->ilen[i]; j++) {
239:       aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
240:     }
241:   }
242:   B->nbs = ec;

244:   baij->B->cmap->n = ec*mat->rmap->bs;

246:   PetscFree(indices);
247: #endif

249:   /* create local vector that is used to scatter into */
250:   VecCreateSeq(PETSC_COMM_SELF,ec*bs,&baij->lvec);

252:   /* create two temporary index sets for building scatter-gather */
253:   ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,PETSC_COPY_VALUES,&from);

255:   PetscMalloc(ec*sizeof(PetscInt),&stmp);
256:   for (i=0; i<ec; i++) stmp[i] = i;
257:   ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,PETSC_OWN_POINTER,&to);

259:   /* create temporary global vector to generate scatter context */
260:   VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&gvec);
261:   VecScatterCreate(gvec,from,baij->lvec,to,&baij->Mvctx);
262:   VecDestroy(&gvec);

264:   PetscLogObjectParent(mat,baij->Mvctx);
265:   PetscLogObjectParent(mat,baij->lvec);
266:   PetscLogObjectParent(mat,from);
267:   PetscLogObjectParent(mat,to);

269:   baij->garray = garray;

271:   PetscLogObjectMemory(mat,(ec+1)*sizeof(PetscInt));
272:   ISDestroy(&from);
273:   ISDestroy(&to);
274:   return(0);
275: }

277: /*
278:      Takes the local part of an already assembled MPISBAIJ matrix
279:    and disassembles it. This is to allow new nonzeros into the matrix
280:    that require more communication in the matrix vector multiply.
281:    Thus certain data-structures must be rebuilt.

283:    Kind of slow! But that's what application programmers get when
284:    they are sloppy.
285: */
288: PetscErrorCode MatDisAssemble_MPISBAIJ(Mat A)
289: {
290:   Mat_MPISBAIJ   *baij  = (Mat_MPISBAIJ*)A->data;
291:   Mat            B      = baij->B,Bnew;
292:   Mat_SeqBAIJ    *Bbaij = (Mat_SeqBAIJ*)B->data;
294:   PetscInt       i,j,mbs=Bbaij->mbs,n = A->cmap->N,col,*garray=baij->garray;
295:   PetscInt       k,bs=A->rmap->bs,bs2=baij->bs2,*rvals,*nz,ec,m=A->rmap->n;
296:   MatScalar      *a = Bbaij->a;
297:   PetscScalar    *atmp;
298: #if defined(PETSC_USE_REAL_MAT_SINGLE)
299:   PetscInt l;
300: #endif

303: #if defined(PETSC_USE_REAL_MAT_SINGLE)
304:   PetscMalloc(A->rmap->bs*sizeof(PetscScalar),&atmp);
305: #endif
306:   /* free stuff related to matrix-vec multiply */
307:   VecGetSize(baij->lvec,&ec); /* needed for PetscLogObjectMemory below */
308:   VecDestroy(&baij->lvec);
309:   VecScatterDestroy(&baij->Mvctx);

311:   VecDestroy(&baij->slvec0);
312:   VecDestroy(&baij->slvec0b);
313:   VecDestroy(&baij->slvec1);
314:   VecDestroy(&baij->slvec1a);
315:   VecDestroy(&baij->slvec1b);

317:   if (baij->colmap) {
318: #if defined(PETSC_USE_CTABLE)
319:     PetscTableDestroy(&baij->colmap);
320: #else
321:     PetscFree(baij->colmap);
322:     PetscLogObjectMemory(A,-Bbaij->nbs*sizeof(PetscInt));
323: #endif
324:   }

326:   /* make sure that B is assembled so we can access its values */
327:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
328:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

330:   /* invent new B and copy stuff over */
331:   PetscMalloc(mbs*sizeof(PetscInt),&nz);
332:   for (i=0; i<mbs; i++) {
333:     nz[i] = Bbaij->i[i+1]-Bbaij->i[i];
334:   }
335:   MatCreate(PETSC_COMM_SELF,&Bnew);
336:   MatSetSizes(Bnew,m,n,m,n);
337:   MatSetType(Bnew,((PetscObject)B)->type_name);
338:   MatSeqBAIJSetPreallocation(Bnew,B->rmap->bs,0,nz);

340:   ((Mat_SeqSBAIJ*)Bnew->data)->nonew = Bbaij->nonew; /* Inherit insertion error options. */

342:   PetscFree(nz);

344:   PetscMalloc(bs*sizeof(PetscInt),&rvals);
345:   for (i=0; i<mbs; i++) {
346:     rvals[0] = bs*i;
347:     for (j=1; j<bs; j++) rvals[j] = rvals[j-1] + 1;
348:     for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) {
349:       col = garray[Bbaij->j[j]]*bs;
350:       for (k=0; k<bs; k++) {
351: #if defined(PETSC_USE_REAL_MAT_SINGLE)
352:         for (l=0; l<bs; l++) atmp[l] = a[j*bs2+l];
353: #else
354:         atmp = a+j*bs2 + k*bs;
355: #endif
356:         MatSetValues_SeqSBAIJ(Bnew,bs,rvals,1,&col,atmp,B->insertmode);
357:         col++;
358:       }
359:     }
360:   }
361: #if defined(PETSC_USE_REAL_MAT_SINGLE)
362:   PetscFree(atmp);
363: #endif
364:   PetscFree(baij->garray);

366:   baij->garray = 0;

368:   PetscFree(rvals);
369:   PetscLogObjectMemory(A,-ec*sizeof(PetscInt));
370:   MatDestroy(&B);
371:   PetscLogObjectParent(A,Bnew);

373:   baij->B = Bnew;

375:   A->was_assembled = PETSC_FALSE;
376:   return(0);
377: }