Actual source code: util.c

petsc-3.14.6 2021-03-30
Report Typos and Errors
  1: /*
  2:  GAMG geometric-algebric multigrid PC - Mark Adams 2011
  3:  */
  4: #include <petsc/private/matimpl.h>
  5: #include <../src/ksp/pc/impls/gamg/gamg.h>

  7: /*
  8:    Produces a set of block column indices of the matrix row, one for each block represented in the original row

 10:    n - the number of block indices in cc[]
 11:    cc - the block indices (must be large enough to contain the indices)
 12: */
 13: PETSC_STATIC_INLINE PetscErrorCode MatCollapseRow(Mat Amat,PetscInt row,PetscInt bs,PetscInt *n,PetscInt *cc)
 14: {
 15:   PetscInt       cnt = -1,nidx,j;
 16:   const PetscInt *idx;

 20:   MatGetRow(Amat,row,&nidx,&idx,NULL);
 21:   if (nidx) {
 22:     cnt = 0;
 23:     cc[cnt] = idx[0]/bs;
 24:     for (j=1; j<nidx; j++) {
 25:       if (cc[cnt] < idx[j]/bs) cc[++cnt] = idx[j]/bs;
 26:     }
 27:   }
 28:   MatRestoreRow(Amat,row,&nidx,&idx,NULL);
 29:   *n = cnt+1;
 30:   return(0);
 31: }

 33: /*
 34:     Produces a set of block column indices of the matrix block row, one for each block represented in the original set of rows

 36:     ncollapsed - the number of block indices
 37:     collapsed - the block indices (must be large enough to contain the indices)
 38: */
 39: PETSC_STATIC_INLINE PetscErrorCode MatCollapseRows(Mat Amat,PetscInt start,PetscInt bs,PetscInt *w0,PetscInt *w1,PetscInt *w2,PetscInt *ncollapsed,PetscInt **collapsed)
 40: {
 41:   PetscInt       i,nprev,*cprev = w0,ncur = 0,*ccur = w1,*merged = w2,*cprevtmp;

 45:   MatCollapseRow(Amat,start,bs,&nprev,cprev);
 46:   for (i=start+1; i<start+bs; i++) {
 47:     MatCollapseRow(Amat,i,bs,&ncur,ccur);
 48:     PetscMergeIntArray(nprev,cprev,ncur,ccur,&nprev,&merged);
 49:     cprevtmp = cprev; cprev = merged; merged = cprevtmp;
 50:   }
 51:   *ncollapsed = nprev;
 52:   if (collapsed) *collapsed  = cprev;
 53:   return(0);
 54: }


 57: /* -------------------------------------------------------------------------- */
 58: /*
 59:    PCGAMGCreateGraph - create simple scaled scalar graph from matrix

 61:  Input Parameter:
 62:  . Amat - matrix
 63:  Output Parameter:
 64:  . a_Gmaat - eoutput scalar graph (symmetric?)
 65:  */
 66: PetscErrorCode PCGAMGCreateGraph(Mat Amat, Mat *a_Gmat)
 67: {
 69:   PetscInt       Istart,Iend,Ii,jj,kk,ncols,nloc,NN,MM,bs;
 70:   MPI_Comm       comm;
 71:   Mat            Gmat;

 74:   PetscObjectGetComm((PetscObject)Amat,&comm);
 75:   MatGetOwnershipRange(Amat, &Istart, &Iend);
 76:   MatGetSize(Amat, &MM, &NN);
 77:   MatGetBlockSize(Amat, &bs);
 78:   nloc = (Iend-Istart)/bs;

 80: #if defined PETSC_GAMG_USE_LOG
 81:   PetscLogEventBegin(petsc_gamg_setup_events[GRAPH],0,0,0,0);
 82: #endif

 84:   if (bs > 1) {
 85:     const PetscScalar *vals;
 86:     const PetscInt    *idx;
 87:     PetscInt          *d_nnz, *o_nnz,*w0,*w1,*w2;
 88:     PetscBool         ismpiaij,isseqaij;

 90:     /*
 91:        Determine the preallocation needed for the scalar matrix derived from the vector matrix.
 92:     */

 94:     PetscObjectBaseTypeCompare((PetscObject)Amat,MATSEQAIJ,&isseqaij);
 95:     PetscObjectBaseTypeCompare((PetscObject)Amat,MATMPIAIJ,&ismpiaij);
 96:     PetscMalloc2(nloc, &d_nnz,isseqaij ? 0 : nloc, &o_nnz);

 98:     if (isseqaij) {
 99:       PetscInt max_d_nnz;

101:       /*
102:           Determine exact preallocation count for (sequential) scalar matrix
103:       */
104:       MatSeqAIJGetMaxRowNonzeros(Amat,&max_d_nnz);
105:       max_d_nnz = PetscMin(nloc,bs*max_d_nnz);
106:       PetscMalloc3(max_d_nnz, &w0,max_d_nnz, &w1,max_d_nnz, &w2);
107:       for (Ii = 0, jj = 0; Ii < Iend; Ii += bs, jj++) {
108:         MatCollapseRows(Amat,Ii,bs,w0,w1,w2,&d_nnz[jj],NULL);
109:       }
110:       PetscFree3(w0,w1,w2);

112:     } else if (ismpiaij) {
113:       Mat            Daij,Oaij;
114:       const PetscInt *garray;
115:       PetscInt       max_d_nnz;

117:       MatMPIAIJGetSeqAIJ(Amat,&Daij,&Oaij,&garray);

119:       /*
120:           Determine exact preallocation count for diagonal block portion of scalar matrix
121:       */
122:       MatSeqAIJGetMaxRowNonzeros(Daij,&max_d_nnz);
123:       max_d_nnz = PetscMin(nloc,bs*max_d_nnz);
124:       PetscMalloc3(max_d_nnz, &w0,max_d_nnz, &w1,max_d_nnz, &w2);
125:       for (Ii = 0, jj = 0; Ii < Iend - Istart; Ii += bs, jj++) {
126:         MatCollapseRows(Daij,Ii,bs,w0,w1,w2,&d_nnz[jj],NULL);
127:       }
128:       PetscFree3(w0,w1,w2);

130:       /*
131:          Over estimate (usually grossly over), preallocation count for off-diagonal portion of scalar matrix
132:       */
133:       for (Ii = 0, jj = 0; Ii < Iend - Istart; Ii += bs, jj++) {
134:         o_nnz[jj] = 0;
135:         for (kk=0; kk<bs; kk++) { /* rows that get collapsed to a single row */
136:           MatGetRow(Oaij,Ii+kk,&ncols,NULL,NULL);
137:           o_nnz[jj] += ncols;
138:           MatRestoreRow(Oaij,Ii+kk,&ncols,NULL,NULL);
139:         }
140:         if (o_nnz[jj] > (NN/bs-nloc)) o_nnz[jj] = NN/bs-nloc;
141:       }

143:     } else SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Require AIJ matrix type");

145:     /* get scalar copy (norms) of matrix */
146:     MatCreate(comm, &Gmat);
147:     MatSetSizes(Gmat,nloc,nloc,PETSC_DETERMINE,PETSC_DETERMINE);
148:     MatSetBlockSizes(Gmat, 1, 1);
149:     MatSetType(Gmat, MATAIJ);
150:     MatSeqAIJSetPreallocation(Gmat,0,d_nnz);
151:     MatMPIAIJSetPreallocation(Gmat,0,d_nnz,0,o_nnz);
152:     PetscFree2(d_nnz,o_nnz);

154:     for (Ii = Istart; Ii < Iend; Ii++) {
155:       PetscInt dest_row = Ii/bs;
156:       MatGetRow(Amat,Ii,&ncols,&idx,&vals);
157:       for (jj=0; jj<ncols; jj++) {
158:         PetscInt    dest_col = idx[jj]/bs;
159:         PetscScalar sv       = PetscAbs(PetscRealPart(vals[jj]));
160:         MatSetValues(Gmat,1,&dest_row,1,&dest_col,&sv,ADD_VALUES);
161:       }
162:       MatRestoreRow(Amat,Ii,&ncols,&idx,&vals);
163:     }
164:     MatAssemblyBegin(Gmat,MAT_FINAL_ASSEMBLY);
165:     MatAssemblyEnd(Gmat,MAT_FINAL_ASSEMBLY);
166:   } else {
167:     /* just copy scalar matrix - abs() not taken here but scaled later */
168:     MatDuplicate(Amat, MAT_COPY_VALUES, &Gmat);
169:   }
170:   MatPropagateSymmetryOptions(Amat, Gmat);

172: #if defined PETSC_GAMG_USE_LOG
173:   PetscLogEventEnd(petsc_gamg_setup_events[GRAPH],0,0,0,0);
174: #endif

176:   *a_Gmat = Gmat;
177:   return(0);
178: }

180: /* -------------------------------------------------------------------------- */
181: /*@C
182:    PCGAMGFilterGraph - filter (remove zero and possibly small values from the) graph and make it symmetric if requested

184:    Collective on Mat

186:    Input Parameter:
187: +   a_Gmat - the graph
188: .   vfilter - threshold parameter [0,1)
189: -   symm - make the result symmetric

191:    Level: developer

193:    Notes:
194:     This is called before graph coarsers are called.

196: .seealso: PCGAMGSetThreshold()
197: @*/
198: PetscErrorCode PCGAMGFilterGraph(Mat *a_Gmat,PetscReal vfilter,PetscBool symm)
199: {
200:   PetscErrorCode    ierr;
201:   PetscInt          Istart,Iend,Ii,jj,ncols,nnz0,nnz1, NN, MM, nloc;
202:   PetscMPIInt       rank;
203:   Mat               Gmat  = *a_Gmat, tGmat, matTrans;
204:   MPI_Comm          comm;
205:   const PetscScalar *vals;
206:   const PetscInt    *idx;
207:   PetscInt          *d_nnz, *o_nnz;
208:   Vec               diag;

211: #if defined PETSC_GAMG_USE_LOG
212:   PetscLogEventBegin(petsc_gamg_setup_events[GRAPH],0,0,0,0);
213: #endif
214:   /* scale Gmat for all values between -1 and 1 */
215:   MatCreateVecs(Gmat, &diag, NULL);
216:   MatGetDiagonal(Gmat, diag);
217:   VecReciprocal(diag);
218:   VecSqrtAbs(diag);
219:   MatDiagonalScale(Gmat, diag, diag);
220:   VecDestroy(&diag);

222:   if (vfilter < 0.0 && !symm) {
223:     /* Just use the provided matrix as the graph but make all values positive */
224:     MatInfo     info;
225:     PetscScalar *avals;
226:     PetscBool isaij,ismpiaij;
227:     PetscObjectBaseTypeCompare((PetscObject)Gmat,MATSEQAIJ,&isaij);
228:     PetscObjectBaseTypeCompare((PetscObject)Gmat,MATMPIAIJ,&ismpiaij);
229:     if (!isaij && !ismpiaij) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Require (MPI)AIJ matrix type");
230:     if (isaij) {
231:       MatGetInfo(Gmat,MAT_LOCAL,&info);
232:       MatSeqAIJGetArray(Gmat,&avals);
233:       for (jj = 0; jj<info.nz_used; jj++) avals[jj] = PetscAbsScalar(avals[jj]);
234:       MatSeqAIJRestoreArray(Gmat,&avals);
235:     } else {
236:       Mat_MPIAIJ  *aij = (Mat_MPIAIJ*)Gmat->data;
237:       MatGetInfo(aij->A,MAT_LOCAL,&info);
238:       MatSeqAIJGetArray(aij->A,&avals);
239:       for (jj = 0; jj<info.nz_used; jj++) avals[jj] = PetscAbsScalar(avals[jj]);
240:       MatSeqAIJRestoreArray(aij->A,&avals);
241:       MatGetInfo(aij->B,MAT_LOCAL,&info);
242:       MatSeqAIJGetArray(aij->B,&avals);
243:       for (jj = 0; jj<info.nz_used; jj++) avals[jj] = PetscAbsScalar(avals[jj]);
244:       MatSeqAIJRestoreArray(aij->B,&avals);
245:     }
246: #if defined PETSC_GAMG_USE_LOG
247:     PetscLogEventEnd(petsc_gamg_setup_events[GRAPH],0,0,0,0);
248: #endif
249:     return(0);
250:   }

252:   PetscObjectGetComm((PetscObject)Gmat,&comm);
253:   MPI_Comm_rank(comm,&rank);
254:   MatGetOwnershipRange(Gmat, &Istart, &Iend);
255:   nloc = Iend - Istart;
256:   MatGetSize(Gmat, &MM, &NN);

258:   if (symm) {
259:     MatTranspose(Gmat, MAT_INITIAL_MATRIX, &matTrans);
260:   }

262:   /* Determine upper bound on nonzeros needed in new filtered matrix */
263:   PetscMalloc2(nloc, &d_nnz,nloc, &o_nnz);
264:   for (Ii = Istart, jj = 0; Ii < Iend; Ii++, jj++) {
265:     MatGetRow(Gmat,Ii,&ncols,NULL,NULL);
266:     d_nnz[jj] = ncols;
267:     o_nnz[jj] = ncols;
268:     MatRestoreRow(Gmat,Ii,&ncols,NULL,NULL);
269:     if (symm) {
270:       MatGetRow(matTrans,Ii,&ncols,NULL,NULL);
271:       d_nnz[jj] += ncols;
272:       o_nnz[jj] += ncols;
273:       MatRestoreRow(matTrans,Ii,&ncols,NULL,NULL);
274:     }
275:     if (d_nnz[jj] > nloc) d_nnz[jj] = nloc;
276:     if (o_nnz[jj] > (MM-nloc)) o_nnz[jj] = MM - nloc;
277:   }
278:   MatCreate(comm, &tGmat);
279:   MatSetSizes(tGmat,nloc,nloc,MM,MM);
280:   MatSetBlockSizes(tGmat, 1, 1);
281:   MatSetType(tGmat, MATAIJ);
282:   MatSeqAIJSetPreallocation(tGmat,0,d_nnz);
283:   MatMPIAIJSetPreallocation(tGmat,0,d_nnz,0,o_nnz);
284:   PetscFree2(d_nnz,o_nnz);
285:   if (symm) {
286:     MatDestroy(&matTrans);
287:   } else {
288:     /* all entries are generated locally so MatAssembly will be slightly faster for large process counts */
289:     MatSetOption(tGmat,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE);
290:   }

292:   for (Ii = Istart, nnz0 = nnz1 = 0; Ii < Iend; Ii++) {
293:     MatGetRow(Gmat,Ii,&ncols,&idx,&vals);
294:     for (jj=0; jj<ncols; jj++,nnz0++) {
295:       PetscScalar sv = PetscAbs(PetscRealPart(vals[jj]));
296:       if (PetscRealPart(sv) > vfilter) {
297:         nnz1++;
298:         if (symm) {
299:           sv  *= 0.5;
300:           MatSetValues(tGmat,1,&Ii,1,&idx[jj],&sv,ADD_VALUES);
301:           MatSetValues(tGmat,1,&idx[jj],1,&Ii,&sv,ADD_VALUES);
302:         } else {
303:           MatSetValues(tGmat,1,&Ii,1,&idx[jj],&sv,ADD_VALUES);
304:         }
305:       }
306:     }
307:     MatRestoreRow(Gmat,Ii,&ncols,&idx,&vals);
308:   }
309:   MatAssemblyBegin(tGmat,MAT_FINAL_ASSEMBLY);
310:   MatAssemblyEnd(tGmat,MAT_FINAL_ASSEMBLY);
311:   if (symm) {
312:     MatSetOption(tGmat,MAT_SYMMETRIC,PETSC_TRUE);
313:   } else {
314:     MatPropagateSymmetryOptions(Gmat,tGmat);
315:   }
316: #if defined PETSC_GAMG_USE_LOG
317:   PetscLogEventEnd(petsc_gamg_setup_events[GRAPH],0,0,0,0);
318: #endif

320: #if defined(PETSC_USE_INFO)
321:   {
322:     double t1 = (!nnz0) ? 1. : 100.*(double)nnz1/(double)nnz0, t2 = (!nloc) ? 1. : (double)nnz0/(double)nloc;
323:     PetscInfo4(*a_Gmat,"\t %g%% nnz after filtering, with threshold %g, %g nnz ave. (N=%D)\n",t1,vfilter,t2,MM);
324:   }
325: #endif
326:   MatDestroy(&Gmat);
327:   *a_Gmat = tGmat;
328:   return(0);
329: }

331: /* -------------------------------------------------------------------------- */
332: /*
333:    PCGAMGGetDataWithGhosts - hacks into Mat MPIAIJ so this must have size > 1

335:    Input Parameter:
336:    . Gmat - MPIAIJ matrix for scattters
337:    . data_sz - number of data terms per node (# cols in output)
338:    . data_in[nloc*data_sz] - column oriented data
339:    Output Parameter:
340:    . a_stride - numbrt of rows of output
341:    . a_data_out[stride*data_sz] - output data with ghosts
342: */
343: PetscErrorCode PCGAMGGetDataWithGhosts(Mat Gmat,PetscInt data_sz,PetscReal data_in[],PetscInt *a_stride,PetscReal **a_data_out)
344: {
346:   Vec            tmp_crds;
347:   Mat_MPIAIJ     *mpimat = (Mat_MPIAIJ*)Gmat->data;
348:   PetscInt       nnodes,num_ghosts,dir,kk,jj,my0,Iend,nloc;
349:   PetscScalar    *data_arr;
350:   PetscReal      *datas;
351:   PetscBool      isMPIAIJ;

354:   PetscObjectBaseTypeCompare((PetscObject)Gmat, MATMPIAIJ, &isMPIAIJ);
355:   MatGetOwnershipRange(Gmat, &my0, &Iend);
356:   nloc      = Iend - my0;
357:   VecGetLocalSize(mpimat->lvec, &num_ghosts);
358:   nnodes    = num_ghosts + nloc;
359:   *a_stride = nnodes;
360:   MatCreateVecs(Gmat, &tmp_crds, NULL);

362:   PetscMalloc1(data_sz*nnodes, &datas);
363:   for (dir=0; dir<data_sz; dir++) {
364:     /* set local, and global */
365:     for (kk=0; kk<nloc; kk++) {
366:       PetscInt    gid = my0 + kk;
367:       PetscScalar crd = (PetscScalar)data_in[dir*nloc + kk]; /* col oriented */
368:       datas[dir*nnodes + kk] = PetscRealPart(crd);

370:       VecSetValues(tmp_crds, 1, &gid, &crd, INSERT_VALUES);
371:     }
372:     VecAssemblyBegin(tmp_crds);
373:     VecAssemblyEnd(tmp_crds);
374:     /* get ghost datas */
375:     VecScatterBegin(mpimat->Mvctx,tmp_crds,mpimat->lvec,INSERT_VALUES,SCATTER_FORWARD);
376:     VecScatterEnd(mpimat->Mvctx,tmp_crds,mpimat->lvec,INSERT_VALUES,SCATTER_FORWARD);
377:     VecGetArray(mpimat->lvec, &data_arr);
378:     for (kk=nloc,jj=0;jj<num_ghosts;kk++,jj++) datas[dir*nnodes + kk] = PetscRealPart(data_arr[jj]);
379:     VecRestoreArray(mpimat->lvec, &data_arr);
380:   }
381:   VecDestroy(&tmp_crds);
382:   *a_data_out = datas;
383:   return(0);
384: }

386: PetscErrorCode PCGAMGHashTableCreate(PetscInt a_size, PCGAMGHashTable *a_tab)
387: {
389:   PetscInt       kk;

392:   a_tab->size = a_size;
393:   PetscMalloc2(a_size, &a_tab->table,a_size, &a_tab->data);
394:   for (kk=0; kk<a_size; kk++) a_tab->table[kk] = -1;
395:   return(0);
396: }

398: PetscErrorCode PCGAMGHashTableDestroy(PCGAMGHashTable *a_tab)
399: {

403:   PetscFree2(a_tab->table,a_tab->data);
404:   return(0);
405: }

407: PetscErrorCode PCGAMGHashTableAdd(PCGAMGHashTable *a_tab, PetscInt a_key, PetscInt a_data)
408: {
409:   PetscInt kk,idx;

412:   if (a_key<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Negative key %D.",a_key);
413:   for (kk = 0, idx = GAMG_HASH(a_key); kk < a_tab->size; kk++, idx = (idx==(a_tab->size-1)) ? 0 : idx + 1) {
414:     if (a_tab->table[idx] == a_key) {
415:       /* exists */
416:       a_tab->data[idx] = a_data;
417:       break;
418:     } else if (a_tab->table[idx] == -1) {
419:       /* add */
420:       a_tab->table[idx] = a_key;
421:       a_tab->data[idx]  = a_data;
422:       break;
423:     }
424:   }
425:   if (kk==a_tab->size) {
426:     /* this is not to efficient, waiting until completely full */
427:     PetscInt       oldsize = a_tab->size, new_size = 2*a_tab->size + 5, *oldtable = a_tab->table, *olddata = a_tab->data;

430:     a_tab->size = new_size;
431:     PetscMalloc2(a_tab->size, &a_tab->table,a_tab->size, &a_tab->data);
432:     for (kk=0;kk<a_tab->size;kk++) a_tab->table[kk] = -1;
433:     for (kk=0;kk<oldsize;kk++) {
434:       if (oldtable[kk] != -1) {
435:         PCGAMGHashTableAdd(a_tab, oldtable[kk], olddata[kk]);
436:        }
437:     }
438:     PetscFree2(oldtable,olddata);
439:     PCGAMGHashTableAdd(a_tab, a_key, a_data);
440:   }
441:   return(0);
442: }