Actual source code: mmaij.c

petsc-3.11.4 2019-09-28
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/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,ec = 0,*garray;
 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

 26:   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+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:   } else {
104:     garray = aij->garray;
105:   }

107:   if (!aij->lvec) {
108:     /* create local vector that is used to scatter into */
109:     VecCreateSeq(PETSC_COMM_SELF,ec,&aij->lvec);
110:   } else {
111:     VecGetSize(aij->lvec,&ec);
112:   }

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

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

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

123:   /* generate the scatter context */
124:   if (aij->Mvctx_mpi1_flg) {
125:     VecScatterDestroy(&aij->Mvctx_mpi1);
126:     VecScatterCreate(gvec,from,aij->lvec,to,&aij->Mvctx_mpi1);
127:     VecScatterSetType(aij->Mvctx_mpi1,VECSCATTERMPI1);
128:     PetscLogObjectParent((PetscObject)mat,(PetscObject)aij->Mvctx_mpi1);
129:   } else {
130:     VecScatterDestroy(&aij->Mvctx);
131:     VecScatterCreate(gvec,from,aij->lvec,to,&aij->Mvctx);
132:     PetscLogObjectParent((PetscObject)mat,(PetscObject)aij->Mvctx);
133:     PetscLogObjectParent((PetscObject)mat,(PetscObject)aij->lvec);
134:     PetscLogObjectMemory((PetscObject)mat,(ec+1)*sizeof(PetscInt));
135:   }
136:   aij->garray = garray;

138:   PetscLogObjectParent((PetscObject)mat,(PetscObject)from);
139:   PetscLogObjectParent((PetscObject)mat,(PetscObject)to);

141:   ISDestroy(&from);
142:   ISDestroy(&to);
143:   VecDestroy(&gvec);
144:   return(0);
145: }

147: /*
148:      Takes the local part of an already assembled MPIAIJ matrix
149:    and disassembles it. This is to allow new nonzeros into the matrix
150:    that require more communication in the matrix vector multiply.
151:    Thus certain data-structures must be rebuilt.

153:    Kind of slow! But that's what application programmers get when
154:    they are sloppy.
155: */
156: PetscErrorCode MatDisAssemble_MPIAIJ(Mat A)
157: {
158:   Mat_MPIAIJ     *aij  = (Mat_MPIAIJ*)A->data;
159:   Mat            B     = aij->B,Bnew;
160:   Mat_SeqAIJ     *Baij = (Mat_SeqAIJ*)B->data;
162:   PetscInt       i,j,m = B->rmap->n,n = A->cmap->N,col,ct = 0,*garray = aij->garray,*nz,ec;
163:   PetscScalar    v;

166:   /* free stuff related to matrix-vec multiply */
167:   VecGetSize(aij->lvec,&ec); /* needed for PetscLogObjectMemory below */
168:   VecDestroy(&aij->lvec);
169:   if (aij->colmap) {
170: #if defined(PETSC_USE_CTABLE)
171:     PetscTableDestroy(&aij->colmap);
172: #else
173:     PetscFree(aij->colmap);
174:     PetscLogObjectMemory((PetscObject)A,-aij->B->cmap->n*sizeof(PetscInt));
175: #endif
176:   }

178:   /* make sure that B is assembled so we can access its values */
179:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
180:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

182:   /* invent new B and copy stuff over */
183:   PetscMalloc1(m+1,&nz);
184:   for (i=0; i<m; i++) {
185:     nz[i] = Baij->i[i+1] - Baij->i[i];
186:   }
187:   MatCreate(PETSC_COMM_SELF,&Bnew);
188:   MatSetSizes(Bnew,m,n,m,n);
189:   MatSetBlockSizesFromMats(Bnew,A,A);
190:   MatSetType(Bnew,((PetscObject)B)->type_name);
191:   MatSeqAIJSetPreallocation(Bnew,0,nz);

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

197:   /*
198:    Ensure that B's nonzerostate is monotonically increasing.
199:    Or should this follow the MatSetValues() loop to preserve B's nonzerstate across a MatDisAssemble() call?
200:    */
201:   Bnew->nonzerostate = B->nonzerostate;

203:   PetscFree(nz);
204:   for (i=0; i<m; i++) {
205:     for (j=Baij->i[i]; j<Baij->i[i+1]; j++) {
206:       col  = garray[Baij->j[ct]];
207:       v    = Baij->a[ct++];
208:       MatSetValues(Bnew,1,&i,1,&col,&v,B->insertmode);
209:     }
210:   }
211:   PetscFree(aij->garray);
212:   PetscLogObjectMemory((PetscObject)A,-ec*sizeof(PetscInt));
213:   MatDestroy(&B);
214:   PetscLogObjectParent((PetscObject)A,(PetscObject)Bnew);

216:   aij->B           = Bnew;
217:   A->was_assembled = PETSC_FALSE;
218:   return(0);
219: }

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

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


227: PetscErrorCode MatMPIAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale)
228: {
229:   Mat_MPIAIJ     *ina = (Mat_MPIAIJ*) inA->data; /*access private part of matrix */
231:   PetscInt       i,n,nt,cstart,cend,no,*garray = ina->garray,*lindices;
232:   PetscInt       *r_rmapd,*r_rmapo;

235:   MatGetOwnershipRange(inA,&cstart,&cend);
236:   MatGetSize(ina->A,NULL,&n);
237:   PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapd);
238:   nt   = 0;
239:   for (i=0; i<inA->rmap->mapping->n; i++) {
240:     if (inA->rmap->mapping->indices[i] >= cstart && inA->rmap->mapping->indices[i] < cend) {
241:       nt++;
242:       r_rmapd[i] = inA->rmap->mapping->indices[i] + 1;
243:     }
244:   }
245:   if (nt != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D n %D",nt,n);
246:   PetscMalloc1(n+1,&auglyrmapd);
247:   for (i=0; i<inA->rmap->mapping->n; i++) {
248:     if (r_rmapd[i]) {
249:       auglyrmapd[(r_rmapd[i]-1)-cstart] = i;
250:     }
251:   }
252:   PetscFree(r_rmapd);
253:   VecCreateSeq(PETSC_COMM_SELF,n,&auglydd);

255:   PetscCalloc1(inA->cmap->N+1,&lindices);
256:   for (i=0; i<ina->B->cmap->n; i++) {
257:     lindices[garray[i]] = i+1;
258:   }
259:   no   = inA->rmap->mapping->n - nt;
260:   PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapo);
261:   nt   = 0;
262:   for (i=0; i<inA->rmap->mapping->n; i++) {
263:     if (lindices[inA->rmap->mapping->indices[i]]) {
264:       nt++;
265:       r_rmapo[i] = lindices[inA->rmap->mapping->indices[i]];
266:     }
267:   }
268:   if (nt > no) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D no %D",nt,n);
269:   PetscFree(lindices);
270:   PetscMalloc1(nt+1,&auglyrmapo);
271:   for (i=0; i<inA->rmap->mapping->n; i++) {
272:     if (r_rmapo[i]) {
273:       auglyrmapo[(r_rmapo[i]-1)] = i;
274:     }
275:   }
276:   PetscFree(r_rmapo);
277:   VecCreateSeq(PETSC_COMM_SELF,nt,&auglyoo);
278:   return(0);
279: }

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

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

291: PetscErrorCode  MatDiagonalScaleLocal_MPIAIJ(Mat A,Vec scale)
292: {
293:   Mat_MPIAIJ        *a = (Mat_MPIAIJ*) A->data; /*access private part of matrix */
294:   PetscErrorCode    ierr;
295:   PetscInt          n,i;
296:   PetscScalar       *d,*o;
297:   const PetscScalar *s;

300:   if (!auglyrmapd) {
301:     MatMPIAIJDiagonalScaleLocalSetUp(A,scale);
302:   }

304:   VecGetArrayRead(scale,&s);

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

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