Actual source code: mmbaij.c

  2: /*
  3:    Support for the parallel BAIJ matrix vector multiply
  4: */
  5: #include <../src/mat/impls/baij/mpi/mpibaij.h>


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


 30: #if defined (PETSC_USE_CTABLE)
 31:   /* use a table - Mark Adams */
 32:   PetscTableCreate(B->mbs,&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);
 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);
 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;
 67:   PetscLayoutSetUp((baij->B->cmap));
 68:   PetscTableDestroy(&gid1_lid1);
 69: #else
 70:   /* Make an array as long as the number of columns */
 71:   /* mark those columns that are in baij->B */
 72:   PetscMalloc((Nbs+1)*sizeof(PetscInt),&indices);
 73:   PetscMemzero(indices,Nbs*sizeof(PetscInt));
 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:   PetscMalloc((ec+1)*sizeof(PetscInt),&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;
103:   PetscLayoutSetUp((baij->B->cmap));
104:   PetscFree(indices);
105: #endif  

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

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

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

117:   /* create temporary global vector to generate scatter context */
118:   VecCreateMPIWithArray(((PetscObject)mat)->comm,mat->cmap->n,mat->cmap->N,PETSC_NULL,&gvec);

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

122:   PetscLogObjectParent(mat,baij->Mvctx);
123:   PetscLogObjectParent(mat,baij->lvec);
124:   PetscLogObjectParent(mat,from);
125:   PetscLogObjectParent(mat,to);
126:   baij->garray = garray;
127:   PetscLogObjectMemory(mat,(ec+1)*sizeof(PetscInt));
128:   ISDestroy(&from);
129:   ISDestroy(&to);
130:   VecDestroy(&gvec);
131:   return(0);
132: }

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

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


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

171:   /* make sure that B is assembled so we can access its values */
172:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
173:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

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

186:   for (i=0; i<mbs; i++) {
187:     for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) {
188:       col  = garray[Bbaij->j[j]];
189:       atmp = a + j*bs2;
190:       MatSetValuesBlocked_SeqBAIJ(Bnew,1,&i,1,&col,atmp,B->insertmode);
191:     }
192:   }
193:   MatSetOption(Bnew,MAT_ROW_ORIENTED,PETSC_TRUE);

195:   PetscFree(nz);
196:   PetscFree(baij->garray);
197:   PetscLogObjectMemory(A,-ec*sizeof(PetscInt));
198:   MatDestroy(&B);
199:   PetscLogObjectParent(A,Bnew);
200:   baij->B = Bnew;
201:   A->was_assembled = PETSC_FALSE;
202:   A->assembled     = PETSC_FALSE;
203:   return(0);
204: }

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

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


215: PetscErrorCode MatMPIBAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale)
216: {
217:   Mat_MPIBAIJ    *ina = (Mat_MPIBAIJ*) inA->data; /*access private part of matrix */
218:   Mat_SeqBAIJ    *B = (Mat_SeqBAIJ*)ina->B->data;
220:   PetscInt       bs = inA->rmap->bs,i,n,nt,j,cstart,cend,no,*garray = ina->garray,*lindices;
221:   PetscInt       *r_rmapd,*r_rmapo;
222: 
224:   MatGetOwnershipRange(inA,&cstart,&cend);
225:   MatGetSize(ina->A,PETSC_NULL,&n);
226:   PetscMalloc((inA->rmap->bmapping->n+1)*sizeof(PetscInt),&r_rmapd);
227:   PetscMemzero(r_rmapd,inA->rmap->bmapping->n*sizeof(PetscInt));
228:   nt   = 0;
229:   for (i=0; i<inA->rmap->bmapping->n; i++) {
230:     if (inA->rmap->bmapping->indices[i]*bs >= cstart && inA->rmap->bmapping->indices[i]*bs < cend) {
231:       nt++;
232:       r_rmapd[i] = inA->rmap->bmapping->indices[i] + 1;
233:     }
234:   }
235:   if (nt*bs != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt*bs %D n %D",nt*bs,n);
236:   PetscMalloc((n+1)*sizeof(PetscInt),&uglyrmapd);
237:   for (i=0; i<inA->rmap->bmapping->n; i++) {
238:     if (r_rmapd[i]){
239:       for (j=0; j<bs; j++) {
240:         uglyrmapd[(r_rmapd[i]-1)*bs+j-cstart] = i*bs + j;
241:       }
242:     }
243:   }
244:   PetscFree(r_rmapd);
245:   VecCreateSeq(PETSC_COMM_SELF,n,&uglydd);

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

275:   return(0);
276: }

280: PetscErrorCode  MatMPIBAIJDiagonalScaleLocal(Mat A,Vec scale)
281: {
282:   /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */

286:   PetscTryMethod(A,"MatDiagonalScaleLocal_C",(Mat,Vec),(A,scale));
287:   return(0);
288: }

293: PetscErrorCode  MatDiagonalScaleLocal_MPIBAIJ(Mat A,Vec scale)
294: {
295:   Mat_MPIBAIJ    *a = (Mat_MPIBAIJ*) A->data; /*access private part of matrix */
297:   PetscInt       n,i;
298:   PetscScalar    *d,*o,*s;
299: 
301:   if (!uglyrmapd) {
302:     MatMPIBAIJDiagonalScaleLocalSetUp(A,scale);
303:   }

305:   VecGetArray(scale,&s);
306: 
307:   VecGetLocalSize(uglydd,&n);
308:   VecGetArray(uglydd,&d);
309:   for (i=0; i<n; i++) {
310:     d[i] = s[uglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */
311:   }
312:   VecRestoreArray(uglydd,&d);
313:   /* column scale "diagonal" portion of local matrix */
314:   MatDiagonalScale(a->A,PETSC_NULL,uglydd);

316:   VecGetLocalSize(uglyoo,&n);
317:   VecGetArray(uglyoo,&o);
318:   for (i=0; i<n; i++) {
319:     o[i] = s[uglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */
320:   }
321:   VecRestoreArray(scale,&s);
322:   VecRestoreArray(uglyoo,&o);
323:   /* column scale "off-diagonal" portion of local matrix */
324:   MatDiagonalScale(a->B,PETSC_NULL,uglyoo);

326:   return(0);
327: }