Actual source code: mmbaij.c

petsc-3.6.1 2015-08-06
Report Typos and Errors
  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:   PetscMalloc1(ec+1,&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:   PetscCalloc1(Nbs+1,&indices);
 74:   for (i=0; i<B->mbs; i++) {
 75:     for (j=0; j<B->ilen[i]; j++) {
 76:       if (!indices[aj[B->i[i] + j]]) ec++;
 77:       indices[aj[B->i[i] + j]] = 1;
 78:     }
 79:   }

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

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

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

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

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

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

114:   PetscMalloc1(ec+1,&stmp);
115:   for (i=0; i<ec; i++) stmp[i] = i;
116:   ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,PETSC_OWN_POINTER,&to);

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

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

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

128:   baij->garray = garray;

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

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

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


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

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

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

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

190:   MatSetOption(Bnew,MAT_ROW_ORIENTED,PETSC_FALSE);
191:   /*
192:    Ensure that B's nonzerostate is monotonically increasing.
193:    Or should this follow the MatSetValuesBlocked() loop to preserve B's nonzerstate across a MatDisAssemble() call?
194:    */
195:   Bnew->nonzerostate = B->nonzerostate;

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

206:   PetscFree(nz);
207:   PetscFree(baij->garray);
208:   PetscLogObjectMemory((PetscObject)A,-ec*sizeof(PetscInt));
209:   MatDestroy(&B);
210:   PetscLogObjectParent((PetscObject)A,(PetscObject)Bnew);

212:   baij->B          = Bnew;
213:   A->was_assembled = PETSC_FALSE;
214:   A->assembled     = PETSC_FALSE;
215:   return(0);
216: }

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

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


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

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

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

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

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

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

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

311:   VecGetArray(scale,&s);

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

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