Actual source code: mmaij.c

petsc-3.14.6 2021-03-30
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:     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+1,&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+1,&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) {
104:     /* create local vector that is used to scatter into */
105:     VecCreateSeq(PETSC_COMM_SELF,ec,&aij->lvec);
106:   } else {
107:     VecGetSize(aij->lvec,&ec);
108:   }

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

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

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

133:   PetscLogObjectParent((PetscObject)mat,(PetscObject)from);
134:   PetscLogObjectParent((PetscObject)mat,(PetscObject)to);

136:   ISDestroy(&from);
137:   ISDestroy(&to);
138:   VecDestroy(&gvec);
139:   return(0);
140: }

142: /*
143:      Takes the local part of an already assembled MPIAIJ matrix
144:    and disassembles it. This is to allow new nonzeros into the matrix
145:    that require more communication in the matrix vector multiply.
146:    Thus certain data-structures must be rebuilt.

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

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

173:   /* make sure that B is assembled so we can access its values */
174:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
175:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

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

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

192:   /*
193:    Ensure that B's nonzerostate is monotonically increasing.
194:    Or should this follow the MatSetValues() loop to preserve B's nonzerstate across a MatDisAssemble() call?
195:    */
196:   Bnew->nonzerostate = B->nonzerostate;

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

211:   aij->B           = Bnew;
212:   A->was_assembled = PETSC_FALSE;
213:   return(0);
214: }

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

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


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

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

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

276: PetscErrorCode MatMPIAIJDiagonalScaleLocal(Mat A,Vec scale)
277: {
278:   /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */

282:   PetscTryMethod(A,"MatDiagonalScaleLocal_C",(Mat,Vec),(A,scale));
283:   return(0);
284: }

286: PetscErrorCode  MatDiagonalScaleLocal_MPIAIJ(Mat A,Vec scale)
287: {
288:   Mat_MPIAIJ        *a = (Mat_MPIAIJ*) A->data; /*access private part of matrix */
289:   PetscErrorCode    ierr;
290:   PetscInt          n,i;
291:   PetscScalar       *d,*o;
292:   const PetscScalar *s;

295:   if (!auglyrmapd) {
296:     MatMPIAIJDiagonalScaleLocalSetUp(A,scale);
297:   }

299:   VecGetArrayRead(scale,&s);

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

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