Actual source code: mmsell.c

  1: /*
  2:    Support for the parallel SELL matrix vector multiply
  3: */
  4: #include <../src/mat/impls/sell/mpi/mpisell.h>
  5: #include <petsc/private/isimpl.h>

  7: /*
  8:    Takes the local part of an already assembled MPISELL matrix
  9:    and disassembles it. This is to allow new nonzeros into the matrix
 10:    that require more communication in the matrix vector multiply.
 11:    Thus certain data-structures must be rebuilt.

 13:    Kind of slow! But that's what application programmers get when
 14:    they are sloppy.
 15: */
 16: PetscErrorCode MatDisAssemble_MPISELL(Mat A)
 17: {
 18:   Mat_MPISELL    *sell=(Mat_MPISELL*)A->data;
 19:   Mat            B=sell->B,Bnew;
 20:   Mat_SeqSELL    *Bsell=(Mat_SeqSELL*)B->data;
 21:   PetscInt       i,j,totalslices,N=A->cmap->N,ec,row;
 22:   PetscBool      isnonzero;

 24:   /* free stuff related to matrix-vec multiply */
 25:   VecGetSize(sell->lvec,&ec); /* needed for PetscLogObjectMemory below */
 26:   VecDestroy(&sell->lvec);
 27:   VecScatterDestroy(&sell->Mvctx);
 28:   if (sell->colmap) {
 29: #if defined(PETSC_USE_CTABLE)
 30:     PetscTableDestroy(&sell->colmap);
 31: #else
 32:     PetscFree(sell->colmap);
 33:     PetscLogObjectMemory((PetscObject)A,-sell->B->cmap->n*sizeof(PetscInt));
 34: #endif
 35:   }

 37:   /* make sure that B is assembled so we can access its values */
 38:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 39:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

 41:   /* invent new B and copy stuff over */
 42:   MatCreate(PETSC_COMM_SELF,&Bnew);
 43:   MatSetSizes(Bnew,B->rmap->n,N,B->rmap->n,N);
 44:   MatSetBlockSizesFromMats(Bnew,A,A);
 45:   MatSetType(Bnew,((PetscObject)B)->type_name);
 46:   MatSeqSELLSetPreallocation(Bnew,0,Bsell->rlen);
 47:   if (Bsell->nonew >= 0) { /* Inherit insertion error options (if positive). */
 48:     ((Mat_SeqSELL*)Bnew->data)->nonew = Bsell->nonew;
 49:   }

 51:   /*
 52:    Ensure that B's nonzerostate is monotonically increasing.
 53:    Or should this follow the MatSetValues() loop to preserve B's nonzerstate across a MatDisAssemble() call?
 54:    */
 55:   Bnew->nonzerostate = B->nonzerostate;

 57:   totalslices = B->rmap->n/8+((B->rmap->n & 0x07)?1:0); /* floor(n/8) */
 58:   for (i=0; i<totalslices; i++) { /* loop over slices */
 59:     for (j=Bsell->sliidx[i],row=0; j<Bsell->sliidx[i+1]; j++,row=((row+1)&0x07)) {
 60:       isnonzero = (PetscBool)((j-Bsell->sliidx[i])/8 < Bsell->rlen[8*i+row]);
 61:       if (isnonzero) {
 62:         MatSetValue(Bnew,8*i+row,sell->garray[Bsell->colidx[j]],Bsell->val[j],B->insertmode);
 63:       }
 64:     }
 65:   }

 67:   PetscFree(sell->garray);
 68:   PetscLogObjectMemory((PetscObject)A,-ec*sizeof(PetscInt));
 69:   MatDestroy(&B);
 70:   PetscLogObjectParent((PetscObject)A,(PetscObject)Bnew);

 72:   sell->B          = Bnew;
 73:   A->was_assembled = PETSC_FALSE;
 74:   A->assembled     = PETSC_FALSE;
 75:   return 0;
 76: }

 78: PetscErrorCode MatSetUpMultiply_MPISELL(Mat mat)
 79: {
 80:   Mat_MPISELL    *sell=(Mat_MPISELL*)mat->data;
 81:   Mat_SeqSELL    *B=(Mat_SeqSELL*)(sell->B->data);
 82:   PetscInt       i,j,*bcolidx=B->colidx,ec=0,*garray,totalslices;
 83:   IS             from,to;
 84:   Vec            gvec;
 85:   PetscBool      isnonzero;
 86: #if defined(PETSC_USE_CTABLE)
 87:   PetscTable         gid1_lid1;
 88:   PetscTablePosition tpos;
 89:   PetscInt           gid,lid;
 90: #else
 91:   PetscInt N = mat->cmap->N,*indices;
 92: #endif

 94:   totalslices = sell->B->rmap->n/8+((sell->B->rmap->n & 0x07)?1:0); /* floor(n/8) */

 96:   /* ec counts the number of columns that contain nonzeros */
 97: #if defined(PETSC_USE_CTABLE)
 98:   /* use a table */
 99:   PetscTableCreate(sell->B->rmap->n,mat->cmap->N+1,&gid1_lid1);
100:   for (i=0; i<totalslices; i++) { /* loop over slices */
101:     for (j=B->sliidx[i]; j<B->sliidx[i+1]; j++) {
102:       isnonzero = (PetscBool)((j-B->sliidx[i])/8 < B->rlen[(i<<3)+(j&0x07)]);
103:       if (isnonzero) { /* check the mask bit */
104:         PetscInt data,gid1 = bcolidx[j] + 1;
105:         PetscTableFind(gid1_lid1,gid1,&data);
106:         if (!data) {
107:           /* one based table */
108:           PetscTableAdd(gid1_lid1,gid1,++ec,INSERT_VALUES);
109:         }
110:       }
111:     }
112:   }

114:   /* form array of columns we need */
115:   PetscMalloc1(ec,&garray);
116:   PetscTableGetHeadPosition(gid1_lid1,&tpos);
117:   while (tpos) {
118:     PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);
119:     gid--;
120:     lid--;
121:     garray[lid] = gid;
122:   }
123:   PetscSortInt(ec,garray); /* sort, and rebuild */
124:   PetscTableRemoveAll(gid1_lid1);
125:   for (i=0; i<ec; i++) {
126:     PetscTableAdd(gid1_lid1,garray[i]+1,i+1,INSERT_VALUES);
127:   }

129:   /* compact out the extra columns in B */
130:   for (i=0; i<totalslices; i++) { /* loop over slices */
131:     for (j=B->sliidx[i]; j<B->sliidx[i+1]; j++) {
132:       isnonzero = (PetscBool)((j-B->sliidx[i])/8 < B->rlen[(i<<3)+(j&0x07)]);
133:       if (isnonzero) {
134:         PetscInt gid1 = bcolidx[j] + 1;
135:         PetscTableFind(gid1_lid1,gid1,&lid);
136:         lid--;
137:         bcolidx[j] = lid;
138:       }
139:     }
140:   }
141:   PetscLayoutDestroy(&sell->B->cmap);
142:   PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)sell->B),ec,ec,1,&sell->B->cmap);
143:   PetscTableDestroy(&gid1_lid1);
144: #else
145:   /* Make an array as long as the number of columns */
146:   PetscCalloc1(N,&indices);
147:   /* mark those columns that are in sell->B */
148:   for (i=0; i<totalslices; i++) { /* loop over slices */
149:     for (j=B->sliidx[i]; j<B->sliidx[i+1]; j++) {
150:       isnonzero = (PetscBool)((j-B->sliidx[i])/8 < B->rlen[(i<<3)+(j&0x07)]);
151:       if (isnonzero) {
152:         if (!indices[bcolidx[j]]) ec++;
153:         indices[bcolidx[j]] = 1;
154:       }
155:     }
156:   }

158:   /* form array of columns we need */
159:   PetscMalloc1(ec,&garray);
160:   ec   = 0;
161:   for (i=0; i<N; i++) {
162:     if (indices[i]) garray[ec++] = i;
163:   }

165:   /* make indices now point into garray */
166:   for (i=0; i<ec; i++) {
167:     indices[garray[i]] = i;
168:   }

170:   /* compact out the extra columns in B */
171:   for (i=0; i<totalslices; i++) { /* loop over slices */
172:     for (j=B->sliidx[i]; j<B->sliidx[i+1]; j++) {
173:       isnonzero = (PetscBool)((j-B->sliidx[i])/8 < B->rlen[(i<<3)+(j&0x07)]);
174:       if (isnonzero) bcolidx[j] = indices[bcolidx[j]];
175:     }
176:   }
177:   PetscLayoutDestroy(&sell->B->cmap);
178:   PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)sell->B),ec,ec,1,&sell->B->cmap);
179:   PetscFree(indices);
180: #endif
181:   /* create local vector that is used to scatter into */
182:   VecCreateSeq(PETSC_COMM_SELF,ec,&sell->lvec);
183:   /* create two temporary Index sets for build scatter gather */
184:   ISCreateGeneral(PETSC_COMM_SELF,ec,garray,PETSC_COPY_VALUES,&from);
185:   ISCreateStride(PETSC_COMM_SELF,ec,0,1,&to);

187:   /* create temporary global vector to generate scatter context */
188:   /* This does not allocate the array's memory so is efficient */
189:   VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&gvec);

191:   /* generate the scatter context */
192:   VecScatterCreate(gvec,from,sell->lvec,to,&sell->Mvctx);
193:   VecScatterViewFromOptions(sell->Mvctx,(PetscObject)mat,"-matmult_vecscatter_view");
194:   PetscLogObjectParent((PetscObject)mat,(PetscObject)sell->Mvctx);
195:   PetscLogObjectParent((PetscObject)mat,(PetscObject)sell->lvec);
196:   PetscLogObjectParent((PetscObject)mat,(PetscObject)from);
197:   PetscLogObjectParent((PetscObject)mat,(PetscObject)to);

199:   sell->garray = garray;

201:   PetscLogObjectMemory((PetscObject)mat,ec*sizeof(PetscInt));
202:   ISDestroy(&from);
203:   ISDestroy(&to);
204:   VecDestroy(&gvec);
205:   return 0;
206: }

208: /*      ugly stuff added for Glenn someday we should fix this up */
209: static PetscInt *auglyrmapd = NULL,*auglyrmapo = NULL; /* mapping from the local ordering to the "diagonal" and "off-diagonal" parts of the local matrix */
210: static Vec auglydd          = NULL,auglyoo     = NULL; /* work vectors used to scale the two parts of the local matrix */

212: PetscErrorCode MatMPISELLDiagonalScaleLocalSetUp(Mat inA,Vec scale)
213: {
214:   Mat_MPISELL    *ina=(Mat_MPISELL*)inA->data; /*access private part of matrix */
215:   PetscInt       i,n,nt,cstart,cend,no,*garray=ina->garray,*lindices;
216:   PetscInt       *r_rmapd,*r_rmapo;

218:   MatGetOwnershipRange(inA,&cstart,&cend);
219:   MatGetSize(ina->A,NULL,&n);
220:   PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapd);
221:   nt   = 0;
222:   for (i=0; i<inA->rmap->mapping->n; i++) {
223:     if (inA->rmap->mapping->indices[i] >= cstart && inA->rmap->mapping->indices[i] < cend) {
224:       nt++;
225:       r_rmapd[i] = inA->rmap->mapping->indices[i] + 1;
226:     }
227:   }
229:   PetscMalloc1(n+1,&auglyrmapd);
230:   for (i=0; i<inA->rmap->mapping->n; i++) {
231:     if (r_rmapd[i]) {
232:       auglyrmapd[(r_rmapd[i]-1)-cstart] = i;
233:     }
234:   }
235:   PetscFree(r_rmapd);
236:   VecCreateSeq(PETSC_COMM_SELF,n,&auglydd);
237:   PetscCalloc1(inA->cmap->N+1,&lindices);
238:   for (i=0; i<ina->B->cmap->n; i++) {
239:     lindices[garray[i]] = i+1;
240:   }
241:   no   = inA->rmap->mapping->n - nt;
242:   PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapo);
243:   nt   = 0;
244:   for (i=0; i<inA->rmap->mapping->n; i++) {
245:     if (lindices[inA->rmap->mapping->indices[i]]) {
246:       nt++;
247:       r_rmapo[i] = lindices[inA->rmap->mapping->indices[i]];
248:     }
249:   }
251:   PetscFree(lindices);
252:   PetscMalloc1(nt+1,&auglyrmapo);
253:   for (i=0; i<inA->rmap->mapping->n; i++) {
254:     if (r_rmapo[i]) {
255:       auglyrmapo[(r_rmapo[i]-1)] = i;
256:     }
257:   }
258:   PetscFree(r_rmapo);
259:   VecCreateSeq(PETSC_COMM_SELF,nt,&auglyoo);
260:   return 0;
261: }

263: PetscErrorCode MatDiagonalScaleLocal_MPISELL(Mat A,Vec scale)
264: {
265:   Mat_MPISELL       *a=(Mat_MPISELL*)A->data; /*access private part of matrix */
266:   PetscInt          n,i;
267:   PetscScalar       *d,*o;
268:   const PetscScalar *s;

270:   if (!auglyrmapd) {
271:     MatMPISELLDiagonalScaleLocalSetUp(A,scale);
272:   }
273:   VecGetArrayRead(scale,&s);
274:   VecGetLocalSize(auglydd,&n);
275:   VecGetArray(auglydd,&d);
276:   for (i=0; i<n; i++) {
277:     d[i] = s[auglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */
278:   }
279:   VecRestoreArray(auglydd,&d);
280:   /* column scale "diagonal" portion of local matrix */
281:   MatDiagonalScale(a->A,NULL,auglydd);
282:   VecGetLocalSize(auglyoo,&n);
283:   VecGetArray(auglyoo,&o);
284:   for (i=0; i<n; i++) {
285:     o[i] = s[auglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */
286:   }
287:   VecRestoreArrayRead(scale,&s);
288:   VecRestoreArray(auglyoo,&o);
289:   /* column scale "off-diagonal" portion of local matrix */
290:   MatDiagonalScale(a->B,NULL,auglyoo);
291:   return 0;
292: }