Actual source code: mmbaij.c

petsc-3.4.2 2013-07-02
  2: /*
  3:    Support for the parallel BAIJ matrix vector multiply
  4: */
  5: #include <../src/mat/impls/baij/mpi/mpibaij.h>
  6: #include <petsc-private/isimpl.h>    /* needed because accesses data structure of ISLocalToGlobalMapping directly */

  8: extern PetscErrorCode MatSetValuesBlocked_SeqBAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);

 12: PetscErrorCode MatSetUpMultiply_MPIBAIJ(Mat mat)
 13: {
 14:   Mat_MPIBAIJ    *baij = (Mat_MPIBAIJ*)mat->data;
 15:   Mat_SeqBAIJ    *B    = (Mat_SeqBAIJ*)(baij->B->data);
 17:   PetscInt       i,j,*aj = B->j,ec = 0,*garray;
 18:   PetscInt       bs = mat->rmap->bs,*stmp;
 19:   IS             from,to;
 20:   Vec            gvec;
 21: #if defined(PETSC_USE_CTABLE)
 22:   PetscTable         gid1_lid1;
 23:   PetscTablePosition tpos;
 24:   PetscInt           gid,lid;
 25: #else
 26:   PetscInt Nbs = baij->Nbs,*indices;
 27: #endif

 30: #if defined(PETSC_USE_CTABLE)
 31:   /* use a table - Mark Adams */
 32:   PetscTableCreate(B->mbs,baij->Nbs+1,&gid1_lid1);
 33:   for (i=0; i<B->mbs; i++) {
 34:     for (j=0; j<B->ilen[i]; j++) {
 35:       PetscInt data,gid1 = aj[B->i[i]+j] + 1;
 36:       PetscTableFind(gid1_lid1,gid1,&data);
 37:       if (!data) {
 38:         /* one based table */
 39:         PetscTableAdd(gid1_lid1,gid1,++ec,INSERT_VALUES);
 40:       }
 41:     }
 42:   }
 43:   /* form array of columns we need */
 44:   PetscMalloc((ec+1)*sizeof(PetscInt),&garray);
 45:   PetscTableGetHeadPosition(gid1_lid1,&tpos);
 46:   while (tpos) {
 47:     PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);
 48:     gid--; lid--;
 49:     garray[lid] = gid;
 50:   }
 51:   PetscSortInt(ec,garray);
 52:   PetscTableRemoveAll(gid1_lid1);
 53:   for (i=0; i<ec; i++) {
 54:     PetscTableAdd(gid1_lid1,garray[i]+1,i+1,INSERT_VALUES);
 55:   }
 56:   /* compact out the extra columns in B */
 57:   for (i=0; i<B->mbs; i++) {
 58:     for (j=0; j<B->ilen[i]; j++) {
 59:       PetscInt gid1 = aj[B->i[i] + j] + 1;
 60:       PetscTableFind(gid1_lid1,gid1,&lid);
 61:       lid--;
 62:       aj[B->i[i]+j] = lid;
 63:     }
 64:   }
 65:   B->nbs           = ec;
 66:   baij->B->cmap->n = baij->B->cmap->N = ec*mat->rmap->bs;

 68:   PetscLayoutSetUp((baij->B->cmap));
 69:   PetscTableDestroy(&gid1_lid1);
 70: #else
 71:   /* Make an array as long as the number of columns */
 72:   /* mark those columns that are in baij->B */
 73:   PetscMalloc((Nbs+1)*sizeof(PetscInt),&indices);
 74:   PetscMemzero(indices,Nbs*sizeof(PetscInt));
 75:   for (i=0; i<B->mbs; i++) {
 76:     for (j=0; j<B->ilen[i]; j++) {
 77:       if (!indices[aj[B->i[i] + j]]) ec++;
 78:       indices[aj[B->i[i] + j]] = 1;
 79:     }
 80:   }

 82:   /* form array of columns we need */
 83:   PetscMalloc((ec+1)*sizeof(PetscInt),&garray);
 84:   ec   = 0;
 85:   for (i=0; i<Nbs; i++) {
 86:     if (indices[i]) {
 87:       garray[ec++] = i;
 88:     }
 89:   }

 91:   /* make indices now point into garray */
 92:   for (i=0; i<ec; i++) {
 93:     indices[garray[i]] = i;
 94:   }

 96:   /* compact out the extra columns in B */
 97:   for (i=0; i<B->mbs; i++) {
 98:     for (j=0; j<B->ilen[i]; j++) {
 99:       aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
100:     }
101:   }
102:   B->nbs           = ec;
103:   baij->B->cmap->n = baij->B->cmap->N  = ec*mat->rmap->bs;

105:   PetscLayoutSetUp((baij->B->cmap));
106:   PetscFree(indices);
107: #endif

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

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

115:   PetscMalloc((ec+1)*sizeof(PetscInt),&stmp);
116:   for (i=0; i<ec; i++) stmp[i] = i;
117:   ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,PETSC_OWN_POINTER,&to);

119:   /* create temporary global vector to generate scatter context */
120:   VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&gvec);

122:   VecScatterCreate(gvec,from,baij->lvec,to,&baij->Mvctx);

124:   PetscLogObjectParent(mat,baij->Mvctx);
125:   PetscLogObjectParent(mat,baij->lvec);
126:   PetscLogObjectParent(mat,from);
127:   PetscLogObjectParent(mat,to);

129:   baij->garray = garray;

131:   PetscLogObjectMemory(mat,(ec+1)*sizeof(PetscInt));
132:   ISDestroy(&from);
133:   ISDestroy(&to);
134:   VecDestroy(&gvec);
135:   return(0);
136: }

138: /*
139:      Takes the local part of an already assembled MPIBAIJ matrix
140:    and disassembles it. This is to allow new nonzeros into the matrix
141:    that require more communication in the matrix vector multiply.
142:    Thus certain data-structures must be rebuilt.

144:    Kind of slow! But that's what application programmers get when
145:    they are sloppy.
146: */
149: PetscErrorCode MatDisAssemble_MPIBAIJ(Mat A)
150: {
151:   Mat_MPIBAIJ    *baij  = (Mat_MPIBAIJ*)A->data;
152:   Mat            B      = baij->B,Bnew;
153:   Mat_SeqBAIJ    *Bbaij = (Mat_SeqBAIJ*)B->data;
155:   PetscInt       i,j,mbs=Bbaij->mbs,n = A->cmap->N,col,*garray=baij->garray;
156:   PetscInt       bs2 = baij->bs2,*nz,ec,m = A->rmap->n;
157:   MatScalar      *a  = Bbaij->a;
158:   MatScalar      *atmp;


162:   /* free stuff related to matrix-vec multiply */
163:   VecGetSize(baij->lvec,&ec); /* needed for PetscLogObjectMemory below */
164:   VecDestroy(&baij->lvec); baij->lvec = 0;
165:   VecScatterDestroy(&baij->Mvctx); baij->Mvctx = 0;
166:   if (baij->colmap) {
167: #if defined(PETSC_USE_CTABLE)
168:     PetscTableDestroy(&baij->colmap);
169: #else
170:     PetscFree(baij->colmap);
171:     PetscLogObjectMemory(A,-Bbaij->nbs*sizeof(PetscInt));
172: #endif
173:   }

175:   /* make sure that B is assembled so we can access its values */
176:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
177:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

179:   /* invent new B and copy stuff over */
180:   PetscMalloc(mbs*sizeof(PetscInt),&nz);
181:   for (i=0; i<mbs; i++) {
182:     nz[i] = Bbaij->i[i+1]-Bbaij->i[i];
183:   }
184:   MatCreate(PetscObjectComm((PetscObject)B),&Bnew);
185:   MatSetSizes(Bnew,m,n,m,n);
186:   MatSetType(Bnew,((PetscObject)B)->type_name);
187:   MatSeqBAIJSetPreallocation(Bnew,B->rmap->bs,0,nz);

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

191:   MatSetOption(Bnew,MAT_ROW_ORIENTED,PETSC_FALSE);

193:   for (i=0; i<mbs; i++) {
194:     for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) {
195:       col  = garray[Bbaij->j[j]];
196:       atmp = a + j*bs2;
197:       MatSetValuesBlocked_SeqBAIJ(Bnew,1,&i,1,&col,atmp,B->insertmode);
198:     }
199:   }
200:   MatSetOption(Bnew,MAT_ROW_ORIENTED,PETSC_TRUE);

202:   PetscFree(nz);
203:   PetscFree(baij->garray);
204:   PetscLogObjectMemory(A,-ec*sizeof(PetscInt));
205:   MatDestroy(&B);
206:   PetscLogObjectParent(A,Bnew);

208:   baij->B          = Bnew;
209:   A->was_assembled = PETSC_FALSE;
210:   A->assembled     = PETSC_FALSE;
211:   return(0);
212: }

214: /*      ugly stuff added for Glenn someday we should fix this up */

216: static PetscInt *uglyrmapd = 0,*uglyrmapo = 0;  /* mapping from the local ordering to the "diagonal" and "off-diagonal" parts of the local matrix */
217: static Vec      uglydd     = 0,uglyoo     = 0;  /* work vectors used to scale the two parts of the local matrix */


222: PetscErrorCode MatMPIBAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale)
223: {
224:   Mat_MPIBAIJ    *ina = (Mat_MPIBAIJ*) inA->data; /*access private part of matrix */
225:   Mat_SeqBAIJ    *B   = (Mat_SeqBAIJ*)ina->B->data;
227:   PetscInt       bs = inA->rmap->bs,i,n,nt,j,cstart,cend,no,*garray = ina->garray,*lindices;
228:   PetscInt       *r_rmapd,*r_rmapo;

231:   MatGetOwnershipRange(inA,&cstart,&cend);
232:   MatGetSize(ina->A,NULL,&n);
233:   PetscMalloc((inA->rmap->bmapping->n+1)*sizeof(PetscInt),&r_rmapd);
234:   PetscMemzero(r_rmapd,inA->rmap->bmapping->n*sizeof(PetscInt));
235:   nt   = 0;
236:   for (i=0; i<inA->rmap->bmapping->n; i++) {
237:     if (inA->rmap->bmapping->indices[i]*bs >= cstart && inA->rmap->bmapping->indices[i]*bs < cend) {
238:       nt++;
239:       r_rmapd[i] = inA->rmap->bmapping->indices[i] + 1;
240:     }
241:   }
242:   if (nt*bs != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt*bs %D n %D",nt*bs,n);
243:   PetscMalloc((n+1)*sizeof(PetscInt),&uglyrmapd);
244:   for (i=0; i<inA->rmap->bmapping->n; i++) {
245:     if (r_rmapd[i]) {
246:       for (j=0; j<bs; j++) {
247:         uglyrmapd[(r_rmapd[i]-1)*bs+j-cstart] = i*bs + j;
248:       }
249:     }
250:   }
251:   PetscFree(r_rmapd);
252:   VecCreateSeq(PETSC_COMM_SELF,n,&uglydd);

254:   PetscMalloc((ina->Nbs+1)*sizeof(PetscInt),&lindices);
255:   PetscMemzero(lindices,ina->Nbs*sizeof(PetscInt));
256:   for (i=0; i<B->nbs; i++) {
257:     lindices[garray[i]] = i+1;
258:   }
259:   no   = inA->rmap->bmapping->n - nt;
260:   PetscMalloc((inA->rmap->bmapping->n+1)*sizeof(PetscInt),&r_rmapo);
261:   PetscMemzero(r_rmapo,inA->rmap->bmapping->n*sizeof(PetscInt));
262:   nt   = 0;
263:   for (i=0; i<inA->rmap->bmapping->n; i++) {
264:     if (lindices[inA->rmap->bmapping->indices[i]]) {
265:       nt++;
266:       r_rmapo[i] = lindices[inA->rmap->bmapping->indices[i]];
267:     }
268:   }
269:   if (nt > no) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D no %D",nt,n);
270:   PetscFree(lindices);
271:   PetscMalloc((nt*bs+1)*sizeof(PetscInt),&uglyrmapo);
272:   for (i=0; i<inA->rmap->bmapping->n; i++) {
273:     if (r_rmapo[i]) {
274:       for (j=0; j<bs; j++) {
275:         uglyrmapo[(r_rmapo[i]-1)*bs+j] = i*bs + j;
276:       }
277:     }
278:   }
279:   PetscFree(r_rmapo);
280:   VecCreateSeq(PETSC_COMM_SELF,nt*bs,&uglyoo);
281:   return(0);
282: }

286: PetscErrorCode  MatMPIBAIJDiagonalScaleLocal(Mat A,Vec scale)
287: {
288:   /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */

292:   PetscTryMethod(A,"MatDiagonalScaleLocal_C",(Mat,Vec),(A,scale));
293:   return(0);
294: }

298: PetscErrorCode  MatDiagonalScaleLocal_MPIBAIJ(Mat A,Vec scale)
299: {
300:   Mat_MPIBAIJ    *a = (Mat_MPIBAIJ*) A->data; /*access private part of matrix */
302:   PetscInt       n,i;
303:   PetscScalar    *d,*o,*s;

306:   if (!uglyrmapd) {
307:     MatMPIBAIJDiagonalScaleLocalSetUp(A,scale);
308:   }

310:   VecGetArray(scale,&s);

312:   VecGetLocalSize(uglydd,&n);
313:   VecGetArray(uglydd,&d);
314:   for (i=0; i<n; i++) {
315:     d[i] = s[uglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */
316:   }
317:   VecRestoreArray(uglydd,&d);
318:   /* column scale "diagonal" portion of local matrix */
319:   MatDiagonalScale(a->A,NULL,uglydd);

321:   VecGetLocalSize(uglyoo,&n);
322:   VecGetArray(uglyoo,&o);
323:   for (i=0; i<n; i++) {
324:     o[i] = s[uglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */
325:   }
326:   VecRestoreArray(scale,&s);
327:   VecRestoreArray(uglyoo,&o);
328:   /* column scale "off-diagonal" portion of local matrix */
329:   MatDiagonalScale(a->B,NULL,uglyoo);
330:   return(0);
331: }