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);
 14:   PetscInt       i,j,*aj = B->j,*garray;
 15:   PetscInt       ec = 0; /* Number of nonzero external columns */
 16:   IS             from,to;
 17:   Vec            gvec;
 18: #if defined(PETSC_USE_CTABLE)
 19:   PetscTable         gid1_lid1;
 20:   PetscTablePosition tpos;
 21:   PetscInt           gid,lid;
 22: #else
 23:   PetscInt N = mat->cmap->N,*indices;
 24: #endif

 27:   if (!aij->garray) {
 28:     if (!aij->B) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing B mat");
 29: #if defined(PETSC_USE_CTABLE)
 30:     /* use a table */
 31:     PetscTableCreate(aij->B->rmap->n,mat->cmap->N+1,&gid1_lid1);
 32:     for (i=0; i<aij->B->rmap->n; i++) {
 33:       for (j=0; j<B->ilen[i]; j++) {
 34:         PetscInt data,gid1 = aj[B->i[i] + j] + 1;
 35:         PetscTableFind(gid1_lid1,gid1,&data);
 36:         if (!data) {
 37:           /* one based table */
 38:           PetscTableAdd(gid1_lid1,gid1,++ec,INSERT_VALUES);
 39:         }
 40:       }
 41:     }
 42:     /* form array of columns we need */
 43:     PetscMalloc1(ec,&garray);
 44:     PetscTableGetHeadPosition(gid1_lid1,&tpos);
 45:     while (tpos) {
 46:       PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);
 47:       gid--;
 48:       lid--;
 49:       garray[lid] = gid;
 50:     }
 51:     PetscSortInt(ec,garray); /* sort, and rebuild */
 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<aij->B->rmap->n; 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:     PetscLayoutDestroy(&aij->B->cmap);
 66:     PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)aij->B),ec,ec,1,&aij->B->cmap);
 67:     PetscTableDestroy(&gid1_lid1);
 68: #else
 69:     /* Make an array as long as the number of columns */
 70:     /* mark those columns that are in aij->B */
 71:     PetscCalloc1(N,&indices);
 72:     for (i=0; i<aij->B->rmap->n; i++) {
 73:       for (j=0; j<B->ilen[i]; j++) {
 74:         if (!indices[aj[B->i[i] + j]]) ec++;
 75:         indices[aj[B->i[i] + j]] = 1;
 76:       }
 77:     }

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

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

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

105:   if (!aij->lvec) {
106:     if (!aij->B) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing B mat");
107:     MatCreateVecs(aij->B,&aij->lvec,NULL);
108:   }
109:   VecGetSize(aij->lvec,&ec);

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

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

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

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

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

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

142:    Kind of slow! But that's what application programmers get when
143:    they are sloppy.
144: */
145: PetscErrorCode MatDisAssemble_MPIAIJ(Mat A)
146: {
147:   Mat_MPIAIJ     *aij  = (Mat_MPIAIJ*)A->data;
148:   Mat            B     = aij->B,Bnew;
149:   Mat_SeqAIJ     *Baij = (Mat_SeqAIJ*)B->data;
151:   PetscInt       i,j,m = B->rmap->n,n = A->cmap->N,col,ct = 0,*garray = aij->garray,*nz,ec;
152:   PetscScalar    v;

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

167:   /* make sure that B is assembled so we can access its values */
168:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
169:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

171:   /* invent new B and copy stuff over */
172:   PetscMalloc1(m+1,&nz);
173:   for (i=0; i<m; i++) {
174:     nz[i] = Baij->i[i+1] - Baij->i[i];
175:   }
176:   MatCreate(PETSC_COMM_SELF,&Bnew);
177:   MatSetSizes(Bnew,m,n,m,n);
178:   MatSetBlockSizesFromMats(Bnew,A,A);
179:   MatSetType(Bnew,((PetscObject)B)->type_name);
180:   MatSeqAIJSetPreallocation(Bnew,0,nz);

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

186:   /*
187:    Ensure that B's nonzerostate is monotonically increasing.
188:    Or should this follow the MatSetValues() loop to preserve B's nonzerstate across a MatDisAssemble() call?
189:    */
190:   Bnew->nonzerostate = B->nonzerostate;

192:   PetscFree(nz);
193:   for (i=0; i<m; i++) {
194:     for (j=Baij->i[i]; j<Baij->i[i+1]; j++) {
195:       col  = garray[Baij->j[ct]];
196:       v    = Baij->a[ct++];
197:       MatSetValues(Bnew,1,&i,1,&col,&v,B->insertmode);
198:     }
199:   }
200:   PetscFree(aij->garray);
201:   PetscLogObjectMemory((PetscObject)A,-ec*sizeof(PetscInt));
202:   MatDestroy(&B);
203:   PetscLogObjectParent((PetscObject)A,(PetscObject)Bnew);

205:   aij->B           = Bnew;
206:   A->was_assembled = PETSC_FALSE;
207:   return(0);
208: }

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

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

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

223:   MatGetOwnershipRange(inA,&cstart,&cend);
224:   MatGetSize(ina->A,NULL,&n);
225:   PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapd);
226:   nt   = 0;
227:   for (i=0; i<inA->rmap->mapping->n; i++) {
228:     if (inA->rmap->mapping->indices[i] >= cstart && inA->rmap->mapping->indices[i] < cend) {
229:       nt++;
230:       r_rmapd[i] = inA->rmap->mapping->indices[i] + 1;
231:     }
232:   }
233:   if (nt != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D n %D",nt,n);
234:   PetscMalloc1(n+1,&auglyrmapd);
235:   for (i=0; i<inA->rmap->mapping->n; i++) {
236:     if (r_rmapd[i]) {
237:       auglyrmapd[(r_rmapd[i]-1)-cstart] = i;
238:     }
239:   }
240:   PetscFree(r_rmapd);
241:   VecCreateSeq(PETSC_COMM_SELF,n,&auglydd);

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

269: PetscErrorCode MatMPIAIJDiagonalScaleLocal(Mat A,Vec scale)
270: {
271:   /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */

275:   PetscTryMethod(A,"MatDiagonalScaleLocal_C",(Mat,Vec),(A,scale));
276:   return(0);
277: }

279: PetscErrorCode  MatDiagonalScaleLocal_MPIAIJ(Mat A,Vec scale)
280: {
281:   Mat_MPIAIJ        *a = (Mat_MPIAIJ*) A->data; /*access private part of matrix */
282:   PetscErrorCode    ierr;
283:   PetscInt          n,i;
284:   PetscScalar       *d,*o;
285:   const PetscScalar *s;

288:   if (!auglyrmapd) {
289:     MatMPIAIJDiagonalScaleLocalSetUp(A,scale);
290:   }

292:   VecGetArrayRead(scale,&s);

294:   VecGetLocalSize(auglydd,&n);
295:   VecGetArray(auglydd,&d);
296:   for (i=0; i<n; i++) {
297:     d[i] = s[auglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */
298:   }
299:   VecRestoreArray(auglydd,&d);
300:   /* column scale "diagonal" portion of local matrix */
301:   MatDiagonalScale(a->A,NULL,auglydd);

303:   VecGetLocalSize(auglyoo,&n);
304:   VecGetArray(auglyoo,&o);
305:   for (i=0; i<n; i++) {
306:     o[i] = s[auglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */
307:   }
308:   VecRestoreArrayRead(scale,&s);
309:   VecRestoreArray(auglyoo,&o);
310:   /* column scale "off-diagonal" portion of local matrix */
311:   MatDiagonalScale(a->B,NULL,auglyoo);
312:   return(0);
313: }