Actual source code: mmaij.c


  2: /*
  3:    Support for the parallel AIJ matrix vector multiply
  4: */
  5: #include <../src/mat/impls/aij/mpi/mpiaij.h>
  6: #include <petsc/private/vecimpl.h>
  7: #include <petsc/private/isimpl.h>

  9: PetscErrorCode MatSetUpMultiply_MPIAIJ(Mat mat)
 10: {
 11:   Mat_MPIAIJ         *aij = (Mat_MPIAIJ*)mat->data;
 12:   Mat_SeqAIJ         *B   = (Mat_SeqAIJ*)(aij->B->data);
 13:   PetscInt           i,j,*aj = B->j,*garray;
 14:   PetscInt           ec = 0; /* Number of nonzero external columns */
 15:   IS                 from,to;
 16:   Vec                gvec;
 17: #if defined(PETSC_USE_CTABLE)
 18:   PetscTable         gid1_lid1;
 19:   PetscTablePosition tpos;
 20:   PetscInt           gid,lid;
 21: #else
 22:   PetscInt           N = mat->cmap->N,*indices;
 23: #endif

 25:   if (!aij->garray) {
 27: #if defined(PETSC_USE_CTABLE)
 28:     /* use a table */
 29:     PetscTableCreate(aij->B->rmap->n,mat->cmap->N+1,&gid1_lid1);
 30:     for (i=0; i<aij->B->rmap->n; i++) {
 31:       for (j=0; j<B->ilen[i]; j++) {
 32:         PetscInt data,gid1 = aj[B->i[i] + j] + 1;
 33:         PetscTableFind(gid1_lid1,gid1,&data);
 34:         if (!data) {
 35:           /* one based table */
 36:           PetscTableAdd(gid1_lid1,gid1,++ec,INSERT_VALUES);
 37:         }
 38:       }
 39:     }
 40:     /* form array of columns we need */
 41:     PetscMalloc1(ec,&garray);
 42:     PetscTableGetHeadPosition(gid1_lid1,&tpos);
 43:     while (tpos) {
 44:       PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);
 45:       gid--;
 46:       lid--;
 47:       garray[lid] = gid;
 48:     }
 49:     PetscSortInt(ec,garray); /* sort, and rebuild */
 50:     PetscTableRemoveAll(gid1_lid1);
 51:     for (i=0; i<ec; i++) {
 52:       PetscTableAdd(gid1_lid1,garray[i]+1,i+1,INSERT_VALUES);
 53:     }
 54:     /* compact out the extra columns in B */
 55:     for (i=0; i<aij->B->rmap->n; i++) {
 56:       for (j=0; j<B->ilen[i]; j++) {
 57:         PetscInt gid1 = aj[B->i[i] + j] + 1;
 58:         PetscTableFind(gid1_lid1,gid1,&lid);
 59:         lid--;
 60:         aj[B->i[i] + j] = lid;
 61:       }
 62:     }
 63:     PetscLayoutDestroy(&aij->B->cmap);
 64:     PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)aij->B),ec,ec,1,&aij->B->cmap);
 65:     PetscTableDestroy(&gid1_lid1);
 66: #else
 67:     /* Make an array as long as the number of columns */
 68:     /* mark those columns that are in aij->B */
 69:     PetscCalloc1(N,&indices);
 70:     for (i=0; i<aij->B->rmap->n; i++) {
 71:       for (j=0; j<B->ilen[i]; j++) {
 72:         if (!indices[aj[B->i[i] + j]]) ec++;
 73:         indices[aj[B->i[i] + j]] = 1;
 74:       }
 75:     }

 77:     /* form array of columns we need */
 78:     PetscMalloc1(ec,&garray);
 79:     ec   = 0;
 80:     for (i=0; i<N; i++) {
 81:       if (indices[i]) garray[ec++] = i;
 82:     }

 84:     /* make indices now point into garray */
 85:     for (i=0; i<ec; i++) {
 86:       indices[garray[i]] = i;
 87:     }

 89:     /* compact out the extra columns in B */
 90:     for (i=0; i<aij->B->rmap->n; i++) {
 91:       for (j=0; j<B->ilen[i]; j++) {
 92:         aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
 93:       }
 94:     }
 95:     PetscLayoutDestroy(&aij->B->cmap);
 96:     PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)aij->B),ec,ec,1,&aij->B->cmap);
 97:     PetscFree(indices);
 98: #endif
 99:   } else {
100:     garray = aij->garray;
101:   }

103:   if (!aij->lvec) {
105:     MatCreateVecs(aij->B,&aij->lvec,NULL);
106:   }
107:   VecGetSize(aij->lvec,&ec);

109:   /* create two temporary Index sets for build scatter gather */
110:   ISCreateGeneral(PETSC_COMM_SELF,ec,garray,PETSC_COPY_VALUES,&from);
111:   ISCreateStride(PETSC_COMM_SELF,ec,0,1,&to);

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

117:   /* generate the scatter context */
118:   VecScatterDestroy(&aij->Mvctx);
119:   VecScatterCreate(gvec,from,aij->lvec,to,&aij->Mvctx);
120:   VecScatterViewFromOptions(aij->Mvctx,(PetscObject)mat,"-matmult_vecscatter_view");
121:   PetscLogObjectParent((PetscObject)mat,(PetscObject)aij->Mvctx);
122:   PetscLogObjectParent((PetscObject)mat,(PetscObject)aij->lvec);
123:   PetscLogObjectMemory((PetscObject)mat,ec*sizeof(PetscInt));
124:   aij->garray = garray;

126:   PetscLogObjectParent((PetscObject)mat,(PetscObject)from);
127:   PetscLogObjectParent((PetscObject)mat,(PetscObject)to);

129:   ISDestroy(&from);
130:   ISDestroy(&to);
131:   VecDestroy(&gvec);
132:   return 0;
133: }

135: /* Disassemble the off-diag portion of the MPIAIJXxx matrix.
136:    In other words, change the B from reduced format using local col ids
137:    to expanded format using global col ids, which will make it easier to
138:    insert new nonzeros (with global col ids) into the matrix.
139:    The off-diag B determines communication in the matrix vector multiply.
140: */
141: PetscErrorCode MatDisAssemble_MPIAIJ(Mat A)
142: {
143:   Mat_MPIAIJ        *aij  = (Mat_MPIAIJ*)A->data;
144:   Mat               B     = aij->B,Bnew;
145:   Mat_SeqAIJ        *Baij = (Mat_SeqAIJ*)B->data;
146:   PetscInt          i,j,m = B->rmap->n,n = A->cmap->N,col,ct = 0,*garray = aij->garray,*nz,ec;
147:   PetscScalar       v;
148:   const PetscScalar *ba;

150:   /* free stuff related to matrix-vec multiply */
151:   VecGetSize(aij->lvec,&ec); /* needed for PetscLogObjectMemory below */
152:   VecDestroy(&aij->lvec);
153:   if (aij->colmap) {
154: #if defined(PETSC_USE_CTABLE)
155:     PetscTableDestroy(&aij->colmap);
156: #else
157:     PetscFree(aij->colmap);
158:     PetscLogObjectMemory((PetscObject)A,-aij->B->cmap->n*sizeof(PetscInt));
159: #endif
160:   }

162:   /* make sure that B is assembled so we can access its values */
163:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
164:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

166:   /* invent new B and copy stuff over */
167:   PetscMalloc1(m+1,&nz);
168:   for (i=0; i<m; i++) {
169:     nz[i] = Baij->i[i+1] - Baij->i[i];
170:   }
171:   MatCreate(PETSC_COMM_SELF,&Bnew);
172:   MatSetSizes(Bnew,m,n,m,n); /* Bnew now uses A->cmap->N as its col size */
173:   MatSetBlockSizesFromMats(Bnew,A,A);
174:   MatSetType(Bnew,((PetscObject)B)->type_name);
175:   MatSeqAIJSetPreallocation(Bnew,0,nz);

177:   if (Baij->nonew >= 0) { /* Inherit insertion error options (if positive). */
178:     ((Mat_SeqAIJ*)Bnew->data)->nonew = Baij->nonew;
179:   }

181:   /*
182:    Ensure that B's nonzerostate is monotonically increasing.
183:    Or should this follow the MatSetValues() loop to preserve B's nonzerstate across a MatDisAssemble() call?
184:    */
185:   Bnew->nonzerostate = B->nonzerostate;

187:   PetscFree(nz);
188:   MatSeqAIJGetArrayRead(B,&ba);
189:   for (i=0; i<m; i++) {
190:     for (j=Baij->i[i]; j<Baij->i[i+1]; j++) {
191:       col  = garray[Baij->j[ct]];
192:       v    = ba[ct++];
193:       MatSetValues(Bnew,1,&i,1,&col,&v,B->insertmode);
194:     }
195:   }
196:   MatSeqAIJRestoreArrayRead(B,&ba);

198:   PetscFree(aij->garray);
199:   PetscLogObjectMemory((PetscObject)A,-ec*sizeof(PetscInt));
200:   MatDestroy(&B);
201:   PetscLogObjectParent((PetscObject)A,(PetscObject)Bnew);

203:   aij->B           = Bnew;
204:   A->was_assembled = PETSC_FALSE;
205:   return 0;
206: }

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

210: static PetscInt *auglyrmapd = NULL,*auglyrmapo = NULL; /* mapping from the local ordering to the "diagonal" and "off-diagonal" parts of the local matrix */
211: static Vec auglydd          = NULL,auglyoo     = NULL; /* work vectors used to scale the two parts of the local matrix */

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

219:   MatGetOwnershipRange(inA,&cstart,&cend);
220:   MatGetSize(ina->A,NULL,&n);
221:   PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapd);
222:   nt   = 0;
223:   for (i=0; i<inA->rmap->mapping->n; i++) {
224:     if (inA->rmap->mapping->indices[i] >= cstart && inA->rmap->mapping->indices[i] < cend) {
225:       nt++;
226:       r_rmapd[i] = inA->rmap->mapping->indices[i] + 1;
227:     }
228:   }
230:   PetscMalloc1(n+1,&auglyrmapd);
231:   for (i=0; i<inA->rmap->mapping->n; i++) {
232:     if (r_rmapd[i]) {
233:       auglyrmapd[(r_rmapd[i]-1)-cstart] = i;
234:     }
235:   }
236:   PetscFree(r_rmapd);
237:   VecCreateSeq(PETSC_COMM_SELF,n,&auglydd);

239:   PetscCalloc1(inA->cmap->N+1,&lindices);
240:   for (i=0; i<ina->B->cmap->n; i++) {
241:     lindices[garray[i]] = i+1;
242:   }
243:   no   = inA->rmap->mapping->n - nt;
244:   PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapo);
245:   nt   = 0;
246:   for (i=0; i<inA->rmap->mapping->n; i++) {
247:     if (lindices[inA->rmap->mapping->indices[i]]) {
248:       nt++;
249:       r_rmapo[i] = lindices[inA->rmap->mapping->indices[i]];
250:     }
251:   }
253:   PetscFree(lindices);
254:   PetscMalloc1(nt+1,&auglyrmapo);
255:   for (i=0; i<inA->rmap->mapping->n; i++) {
256:     if (r_rmapo[i]) {
257:       auglyrmapo[(r_rmapo[i]-1)] = i;
258:     }
259:   }
260:   PetscFree(r_rmapo);
261:   VecCreateSeq(PETSC_COMM_SELF,nt,&auglyoo);
262:   return 0;
263: }

265: PetscErrorCode MatMPIAIJDiagonalScaleLocal(Mat A,Vec scale)
266: {
267:   /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */

269:   PetscTryMethod(A,"MatDiagonalScaleLocal_C",(Mat,Vec),(A,scale));
270:   return 0;
271: }

273: PetscErrorCode  MatDiagonalScaleLocal_MPIAIJ(Mat A,Vec scale)
274: {
275:   Mat_MPIAIJ        *a = (Mat_MPIAIJ*) A->data; /*access private part of matrix */
276:   PetscInt          n,i;
277:   PetscScalar       *d,*o;
278:   const PetscScalar *s;

280:   if (!auglyrmapd) {
281:     MatMPIAIJDiagonalScaleLocalSetUp(A,scale);
282:   }

284:   VecGetArrayRead(scale,&s);

286:   VecGetLocalSize(auglydd,&n);
287:   VecGetArray(auglydd,&d);
288:   for (i=0; i<n; i++) {
289:     d[i] = s[auglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */
290:   }
291:   VecRestoreArray(auglydd,&d);
292:   /* column scale "diagonal" portion of local matrix */
293:   MatDiagonalScale(a->A,NULL,auglydd);

295:   VecGetLocalSize(auglyoo,&n);
296:   VecGetArray(auglyoo,&o);
297:   for (i=0; i<n; i++) {
298:     o[i] = s[auglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */
299:   }
300:   VecRestoreArrayRead(scale,&s);
301:   VecRestoreArray(auglyoo,&o);
302:   /* column scale "off-diagonal" portion of local matrix */
303:   MatDiagonalScale(a->B,NULL,auglyoo);
304:   return 0;
305: }