Actual source code: mmaij.c

petsc-3.7.3 2016-08-01
Report Typos and Errors
  2: /*
  3:    Support for the parallel AIJ matrix vector multiply
  4: */
  5: #include <../src/mat/impls/aij/mpi/mpiaij.h>
  6: #include <petsc/private/isimpl.h>    /* needed because accesses data structure of ISLocalToGlobalMapping directly */

 10: PetscErrorCode MatSetUpMultiply_MPIAIJ(Mat mat)
 11: {
 12:   Mat_MPIAIJ     *aij = (Mat_MPIAIJ*)mat->data;
 13:   Mat_SeqAIJ     *B   = (Mat_SeqAIJ*)(aij->B->data);
 15:   PetscInt       i,j,*aj = B->j,ec = 0,*garray;
 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 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+1,&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:   aij->B->cmap->n = aij->B->cmap->N = ec;
 64:   aij->B->cmap->bs = 1;

 66:   PetscLayoutSetUp((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+1,&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+1,&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:   aij->B->cmap->n = aij->B->cmap->N = ec;
 98:   aij->B->cmap->bs = 1;

100:   PetscLayoutSetUp((aij->B->cmap));
101:   PetscFree(indices);
102: #endif
103:   /* create local vector that is used to scatter into */
104:   VecCreateSeq(PETSC_COMM_SELF,ec,&aij->lvec);

106:   /* create two temporary Index sets for build scatter gather */
107:   ISCreateGeneral(((PetscObject)mat)->comm,ec,garray,PETSC_COPY_VALUES,&from);

109:   ISCreateStride(PETSC_COMM_SELF,ec,0,1,&to);

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

115:   /* generate the scatter context */
116:   VecScatterCreate(gvec,from,aij->lvec,to,&aij->Mvctx);
117:   PetscLogObjectParent((PetscObject)mat,(PetscObject)aij->Mvctx);
118:   PetscLogObjectParent((PetscObject)mat,(PetscObject)aij->lvec);
119:   PetscLogObjectParent((PetscObject)mat,(PetscObject)from);
120:   PetscLogObjectParent((PetscObject)mat,(PetscObject)to);

122:   aij->garray = garray;

124:   PetscLogObjectMemory((PetscObject)mat,(ec+1)*sizeof(PetscInt));
125:   ISDestroy(&from);
126:   ISDestroy(&to);
127:   VecDestroy(&gvec);
128:   return(0);
129: }


134: /*
135:      Takes the local part of an already assembled MPIAIJ 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: */
143: PetscErrorCode MatDisAssemble_MPIAIJ(Mat A)
144: {
145:   Mat_MPIAIJ     *aij  = (Mat_MPIAIJ*)A->data;
146:   Mat            B     = aij->B,Bnew;
147:   Mat_SeqAIJ     *Baij = (Mat_SeqAIJ*)B->data;
149:   PetscInt       i,j,m = B->rmap->n,n = A->cmap->N,col,ct = 0,*garray = aij->garray,*nz,ec;
150:   PetscScalar    v;

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

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

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

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

188:   PetscFree(nz);
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    = Baij->a[ct++];
193:       MatSetValues(Bnew,1,&i,1,&col,&v,B->insertmode);
194:     }
195:   }
196:   PetscFree(aij->garray);
197:   PetscLogObjectMemory((PetscObject)A,-ec*sizeof(PetscInt));
198:   MatDestroy(&B);
199:   PetscLogObjectParent((PetscObject)A,(PetscObject)Bnew);

201:   aij->B           = Bnew;
202:   A->was_assembled = PETSC_FALSE;
203:   return(0);
204: }

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

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


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

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

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

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

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

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

290:   if (!auglyrmapd) {
291:     MatMPIAIJDiagonalScaleLocalSetUp(A,scale);
292:   }

294:   VecGetArray(scale,&s);

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

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