Actual source code: mmsbaij.c

petsc-3.14.6 2021-03-30
Report Typos and Errors

  2: /*
  3:    Support for the parallel SBAIJ matrix vector multiply
  4: */
  5: #include <../src/mat/impls/sbaij/mpi/mpisbaij.h>

  7: PetscErrorCode MatSetUpMultiply_MPISBAIJ(Mat mat)
  8: {
  9:   Mat_MPISBAIJ   *sbaij = (Mat_MPISBAIJ*)mat->data;
 10:   Mat_SeqBAIJ    *B     = (Mat_SeqBAIJ*)(sbaij->B->data);
 12:   PetscInt       Nbs = sbaij->Nbs,i,j,*aj = B->j,ec = 0,*garray,*sgarray;
 13:   PetscInt       bs  = mat->rmap->bs,*stmp,mbs=sbaij->mbs, vec_size,nt;
 14:   IS             from,to;
 15:   Vec            gvec;
 16:   PetscMPIInt    rank = sbaij->rank,lsize;
 17:   PetscInt       *owners = sbaij->rangebs,*ec_owner,k;
 18:   const PetscInt *sowners;
 19:   PetscScalar    *ptr;
 20: #if defined(PETSC_USE_CTABLE)
 21:   PetscTable         gid1_lid1; /* one-based gid to lid table */
 22:   PetscTablePosition tpos;
 23:   PetscInt           gid,lid;
 24: #else
 25:   PetscInt           *indices;
 26: #endif

 29: #if defined(PETSC_USE_CTABLE)
 30:   PetscTableCreate(mbs,Nbs+1,&gid1_lid1);
 31:   for (i=0; i<B->mbs; i++) {
 32:     for (j=0; j<B->ilen[i]; j++) {
 33:       PetscInt data,gid1 = aj[B->i[i]+j] + 1;
 34:       PetscTableFind(gid1_lid1,gid1,&data);
 35:       if (!data) {PetscTableAdd(gid1_lid1,gid1,++ec,INSERT_VALUES);}
 36:     }
 37:   }
 38:   /* form array of columns we need */
 39:   PetscMalloc1(ec,&garray);
 40:   PetscTableGetHeadPosition(gid1_lid1,&tpos);
 41:   while (tpos) {
 42:     PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);
 43:     gid--; lid--;
 44:     garray[lid] = gid;
 45:   }
 46:   PetscSortInt(ec,garray);
 47:   PetscTableRemoveAll(gid1_lid1);
 48:   for (i=0; i<ec; i++) {PetscTableAdd(gid1_lid1,garray[i]+1,i+1,INSERT_VALUES);}
 49:   /* compact out the extra columns in B */
 50:   for (i=0; i<B->mbs; i++) {
 51:     for (j=0; j<B->ilen[i]; j++) {
 52:       PetscInt gid1 = aj[B->i[i] + j] + 1;
 53:       PetscTableFind(gid1_lid1,gid1,&lid);
 54:       lid--;
 55:       aj[B->i[i]+j] = lid;
 56:     }
 57:   }
 58:   PetscTableDestroy(&gid1_lid1);
 59:   PetscMalloc2(2*ec,&sgarray,ec,&ec_owner);
 60:   for (i=j=0; i<ec; i++) {
 61:     while (garray[i]>=owners[j+1]) j++;
 62:     ec_owner[i] = j;
 63:   }
 64: #else
 65:   /* For the first stab we make an array as long as the number of columns */
 66:   /* mark those columns that are in sbaij->B */
 67:   PetscCalloc1(Nbs,&indices);
 68:   for (i=0; i<mbs; i++) {
 69:     for (j=0; j<B->ilen[i]; j++) {
 70:       if (!indices[aj[B->i[i] + j]]) ec++;
 71:       indices[aj[B->i[i] + j]] = 1;
 72:     }
 73:   }

 75:   /* form arrays of columns we need */
 76:   PetscMalloc1(ec,&garray);
 77:   PetscMalloc2(2*ec,&sgarray,ec,&ec_owner);

 79:   ec = 0;
 80:   for (j=0; j<sbaij->size; j++) {
 81:     for (i=owners[j]; i<owners[j+1]; i++) {
 82:       if (indices[i]) {
 83:         garray[ec]   = i;
 84:         ec_owner[ec] = j;
 85:         ec++;
 86:       }
 87:     }
 88:   }

 90:   /* make indices now point into garray */
 91:   for (i=0; i<ec; i++) indices[garray[i]] = i;

 93:   /* compact out the extra columns in B */
 94:   for (i=0; i<mbs; i++) {
 95:     for (j=0; j<B->ilen[i]; j++) aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
 96:   }
 97:   PetscFree(indices);
 98: #endif
 99:   B->nbs = ec;
100:   PetscLayoutDestroy(&sbaij->B->cmap);
101:   PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)sbaij->B),ec*mat->rmap->bs,ec*mat->rmap->bs,mat->rmap->bs,&sbaij->B->cmap);

103:   VecScatterDestroy(&sbaij->sMvctx);
104:   /* create local vector that is used to scatter into */
105:   VecCreateSeq(PETSC_COMM_SELF,ec*bs,&sbaij->lvec);

107:   /* create two temporary index sets for building scatter-gather */
108:   PetscMalloc1(2*ec,&stmp);
109:   ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,PETSC_COPY_VALUES,&from);
110:   for (i=0; i<ec; i++) stmp[i] = i;
111:   ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,PETSC_COPY_VALUES,&to);

113:   /* generate the scatter context
114:      -- Mvctx and lvec are not used by MatMult_MPISBAIJ(), but have other uses, such as in MatDiagonalScale_MPISBAIJ */
115:   VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&gvec);
116:   VecScatterCreate(gvec,from,sbaij->lvec,to,&sbaij->Mvctx);
117:   VecDestroy(&gvec);

119:   sbaij->garray = garray;

121:   PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->Mvctx);
122:   PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->lvec);
123:   PetscLogObjectParent((PetscObject)mat,(PetscObject)from);
124:   PetscLogObjectParent((PetscObject)mat,(PetscObject)to);

126:   ISDestroy(&from);
127:   ISDestroy(&to);

129:   /* create parallel vector that is used by SBAIJ matrix to scatter from/into */
130:   lsize = (mbs + ec)*bs;
131:   VecCreateMPI(PetscObjectComm((PetscObject)mat),lsize,PETSC_DETERMINE,&sbaij->slvec0);
132:   VecDuplicate(sbaij->slvec0,&sbaij->slvec1);
133:   VecGetSize(sbaij->slvec0,&vec_size);

135:   VecGetOwnershipRanges(sbaij->slvec0,&sowners);

137:   /* x index in the IS sfrom */
138:   for (i=0; i<ec; i++) {
139:     j          = ec_owner[i];
140:     sgarray[i] = garray[i] + (sowners[j]/bs - owners[j]);
141:   }
142:   /* b index in the IS sfrom */
143:   k = sowners[rank]/bs + mbs;
144:   for (i=ec,j=0; i< 2*ec; i++,j++) sgarray[i] = k + j;
145:   ISCreateBlock(PETSC_COMM_SELF,bs,2*ec,sgarray,PETSC_COPY_VALUES,&from);

147:   /* x index in the IS sto */
148:   k = sowners[rank]/bs + mbs;
149:   for (i=0; i<ec; i++) stmp[i] = (k + i);
150:   /* b index in the IS sto */
151:   for (i=ec; i<2*ec; i++) stmp[i] = sgarray[i-ec];

153:   ISCreateBlock(PETSC_COMM_SELF,bs,2*ec,stmp,PETSC_COPY_VALUES,&to);

155:   VecScatterCreate(sbaij->slvec0,from,sbaij->slvec1,to,&sbaij->sMvctx);

157:   VecGetLocalSize(sbaij->slvec1,&nt);
158:   VecGetArray(sbaij->slvec1,&ptr);
159:   VecCreateSeqWithArray(PETSC_COMM_SELF,1,bs*mbs,ptr,&sbaij->slvec1a);
160:   VecCreateSeqWithArray(PETSC_COMM_SELF,1,nt-bs*mbs,ptr+bs*mbs,&sbaij->slvec1b);
161:   VecRestoreArray(sbaij->slvec1,&ptr);

163:   VecGetArray(sbaij->slvec0,&ptr);
164:   VecCreateSeqWithArray(PETSC_COMM_SELF,1,nt-bs*mbs,ptr+bs*mbs,&sbaij->slvec0b);
165:   VecRestoreArray(sbaij->slvec0,&ptr);

167:   PetscFree(stmp);
168:   MPI_Barrier(PetscObjectComm((PetscObject)mat));

170:   PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->sMvctx);
171:   PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec0);
172:   PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec1);
173:   PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec0b);
174:   PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec1a);
175:   PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec1b);
176:   PetscLogObjectParent((PetscObject)mat,(PetscObject)from);
177:   PetscLogObjectParent((PetscObject)mat,(PetscObject)to);

179:   PetscLogObjectMemory((PetscObject)mat,(ec+1)*sizeof(PetscInt));
180:   ISDestroy(&from);
181:   ISDestroy(&to);
182:   PetscFree2(sgarray,ec_owner);
183:   return(0);
184: }

186: /*
187:      Takes the local part of an already assembled MPISBAIJ matrix
188:    and disassembles it. This is to allow new nonzeros into the matrix
189:    that require more communication in the matrix vector multiply.
190:    Thus certain data-structures must be rebuilt.

192:    Kind of slow! But that's what application programmers get when
193:    they are sloppy.
194: */
195: PetscErrorCode MatDisAssemble_MPISBAIJ(Mat A)
196: {
197:   Mat_MPISBAIJ   *baij  = (Mat_MPISBAIJ*)A->data;
198:   Mat            B      = baij->B,Bnew;
199:   Mat_SeqBAIJ    *Bbaij = (Mat_SeqBAIJ*)B->data;
201:   PetscInt       i,j,mbs=Bbaij->mbs,n = A->cmap->N,col,*garray=baij->garray;
202:   PetscInt       k,bs=A->rmap->bs,bs2=baij->bs2,*rvals,*nz,ec,m=A->rmap->n;
203:   MatScalar      *a = Bbaij->a;
204:   PetscScalar    *atmp;
205: #if defined(PETSC_USE_REAL_MAT_SINGLE)
206:   PetscInt l;
207: #endif

210: #if defined(PETSC_USE_REAL_MAT_SINGLE)
211:   PetscMalloc1(A->rmap->bs,&atmp);
212: #endif
213:   /* free stuff related to matrix-vec multiply */
214:   VecGetSize(baij->lvec,&ec); /* needed for PetscLogObjectMemory below */
215:   VecDestroy(&baij->lvec);
216:   VecScatterDestroy(&baij->Mvctx);

218:   VecDestroy(&baij->slvec0);
219:   VecDestroy(&baij->slvec0b);
220:   VecDestroy(&baij->slvec1);
221:   VecDestroy(&baij->slvec1a);
222:   VecDestroy(&baij->slvec1b);

224:   if (baij->colmap) {
225: #if defined(PETSC_USE_CTABLE)
226:     PetscTableDestroy(&baij->colmap);
227: #else
228:     PetscFree(baij->colmap);
229:     PetscLogObjectMemory((PetscObject)A,-Bbaij->nbs*sizeof(PetscInt));
230: #endif
231:   }

233:   /* make sure that B is assembled so we can access its values */
234:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
235:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

237:   /* invent new B and copy stuff over */
238:   PetscMalloc1(mbs,&nz);
239:   for (i=0; i<mbs; i++) {
240:     nz[i] = Bbaij->i[i+1]-Bbaij->i[i];
241:   }
242:   MatCreate(PETSC_COMM_SELF,&Bnew);
243:   MatSetSizes(Bnew,m,n,m,n);
244:   MatSetType(Bnew,((PetscObject)B)->type_name);
245:   MatSeqBAIJSetPreallocation(Bnew,B->rmap->bs,0,nz);
246:   PetscFree(nz);

248:   if (Bbaij->nonew >= 0) { /* Inherit insertion error options (if positive). */
249:     ((Mat_SeqSBAIJ*)Bnew->data)->nonew = Bbaij->nonew;
250:   }

252:   /*
253:    Ensure that B's nonzerostate is monotonically increasing.
254:    Or should this follow the MatSetValues() loop to preserve B's nonzerstate across a MatDisAssemble() call?
255:    */
256:   Bnew->nonzerostate = B->nonzerostate;
257:   PetscMalloc1(bs,&rvals);
258:   for (i=0; i<mbs; i++) {
259:     rvals[0] = bs*i;
260:     for (j=1; j<bs; j++) rvals[j] = rvals[j-1] + 1;
261:     for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) {
262:       col = garray[Bbaij->j[j]]*bs;
263:       for (k=0; k<bs; k++) {
264: #if defined(PETSC_USE_REAL_MAT_SINGLE)
265:         for (l=0; l<bs; l++) atmp[l] = a[j*bs2+l];
266: #else
267:         atmp = a+j*bs2 + k*bs;
268: #endif
269:         MatSetValues_SeqSBAIJ(Bnew,bs,rvals,1,&col,atmp,B->insertmode);
270:         col++;
271:       }
272:     }
273:   }
274: #if defined(PETSC_USE_REAL_MAT_SINGLE)
275:   PetscFree(atmp);
276: #endif
277:   PetscFree(baij->garray);

279:   baij->garray = NULL;

281:   PetscFree(rvals);
282:   PetscLogObjectMemory((PetscObject)A,-ec*sizeof(PetscInt));
283:   MatDestroy(&B);
284:   PetscLogObjectParent((PetscObject)A,(PetscObject)Bnew);

286:   baij->B = Bnew;

288:   A->was_assembled = PETSC_FALSE;
289:   return(0);
290: }