1 /* 2 GAMG geometric-algebric multigrid PC - Mark Adams 2011 3 */ 4 5 #include <../src/ksp/pc/impls/gamg/gamg.h> /*I "petscpc.h" I*/ 6 #include <petscblaslapack.h> 7 #include <petscdm.h> 8 #include <petsc/private/kspimpl.h> 9 10 typedef struct { 11 PetscInt nsmooths; 12 PetscInt aggressive_coarsening_levels; // number of aggressive coarsening levels (square or MISk) 13 MatCoarsen crs; 14 } PC_GAMG_AGG; 15 16 /*@ 17 PCGAMGSetNSmooths - Set number of smoothing steps (1 is typical) used for multigrid on all the levels 18 19 Logically Collective on pc 20 21 Input Parameters: 22 . pc - the preconditioner context 23 24 Options Database Key: 25 . -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation 26 27 Level: intermediate 28 29 .seealso: `PCMG`, `PCGAMG` 30 @*/ 31 PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n) 32 { 33 PetscFunctionBegin; 34 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 35 PetscValidLogicalCollectiveInt(pc, n, 2); 36 PetscTryMethod(pc, "PCGAMGSetNSmooths_C", (PC, PetscInt), (pc, n)); 37 PetscFunctionReturn(0); 38 } 39 40 static PetscErrorCode PCGAMGSetNSmooths_AGG(PC pc, PetscInt n) 41 { 42 PC_MG *mg = (PC_MG *)pc->data; 43 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 44 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 45 46 PetscFunctionBegin; 47 pc_gamg_agg->nsmooths = n; 48 PetscFunctionReturn(0); 49 } 50 51 /*@ 52 PCGAMGSetAggressiveLevels - Use aggressive coarsening on first n levels 53 54 Logically Collective on pc 55 56 Input Parameters: 57 + pc - the preconditioner context 58 - n - 0, 1 or more 59 60 Options Database Key: 61 . -pc_gamg_aggressive_coarsening <n,default = 1> - Number of levels to square the graph on before aggregating it 62 63 Level: intermediate 64 65 .seealso: `PCGAMG`, `PCGAMGSetThreshold()` 66 @*/ 67 PetscErrorCode PCGAMGSetAggressiveLevels(PC pc, PetscInt n) 68 { 69 PetscFunctionBegin; 70 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 71 PetscValidLogicalCollectiveInt(pc, n, 2); 72 PetscTryMethod(pc, "PCGAMGSetAggressiveLevels_C", (PC, PetscInt), (pc, n)); 73 PetscFunctionReturn(0); 74 } 75 76 static PetscErrorCode PCGAMGSetAggressiveLevels_AGG(PC pc, PetscInt n) 77 { 78 PC_MG *mg = (PC_MG *)pc->data; 79 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 80 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 81 82 PetscFunctionBegin; 83 pc_gamg_agg->aggressive_coarsening_levels = n; 84 PetscFunctionReturn(0); 85 } 86 87 static PetscErrorCode PCSetFromOptions_GAMG_AGG(PC pc, PetscOptionItems *PetscOptionsObject) 88 { 89 PC_MG *mg = (PC_MG *)pc->data; 90 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 91 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 92 93 PetscFunctionBegin; 94 PetscOptionsHeadBegin(PetscOptionsObject, "GAMG-AGG options"); 95 { 96 PetscBool flg; 97 PetscCall(PetscOptionsInt("-pc_gamg_agg_nsmooths", "smoothing steps for smoothed aggregation, usually 1", "PCGAMGSetNSmooths", pc_gamg_agg->nsmooths, &pc_gamg_agg->nsmooths, NULL)); 98 pc_gamg_agg->aggressive_coarsening_levels = 1; 99 PetscCall( 100 PetscOptionsInt("-pc_gamg_square_graph", "Number of aggressive coarsening (MIS-2) levels from finest (alias for -pc_gamg_aggressive_coarsening, deprecated)", "PCGAMGSetAggressiveLevels", pc_gamg_agg->aggressive_coarsening_levels, &pc_gamg_agg->aggressive_coarsening_levels, &flg)); 101 if (!flg) { 102 PetscCall(PetscOptionsInt("-pc_gamg_aggressive_coarsening", "Number of aggressive coarsening (MIS-2) levels from finest", "PCGAMGSetAggressiveLevels", pc_gamg_agg->aggressive_coarsening_levels, &pc_gamg_agg->aggressive_coarsening_levels, NULL)); 103 } else { 104 PetscCall(PetscOptionsInt("-pc_gamg_aggressive_coarsening", "Number of aggressive coarsening (MIS-2) levels from finest", "PCGAMGSetAggressiveLevels", pc_gamg_agg->aggressive_coarsening_levels, &pc_gamg_agg->aggressive_coarsening_levels, &flg)); 105 if (flg) PetscCall(PetscInfo(pc, "Warning: both -pc_gamg_square_graph and -pc_gamg_aggressive_coarsening are used. -pc_gamg_square_graph is deprecated, Number of aggressive levels is %d\n", (int)pc_gamg_agg->aggressive_coarsening_levels)); 106 } 107 } 108 PetscOptionsHeadEnd(); 109 PetscFunctionReturn(0); 110 } 111 112 static PetscErrorCode PCDestroy_GAMG_AGG(PC pc) 113 { 114 PC_MG *mg = (PC_MG *)pc->data; 115 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 116 117 PetscFunctionBegin; 118 PetscCall(PetscFree(pc_gamg->subctx)); 119 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNSmooths_C", NULL)); 120 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveLevels_C", NULL)); 121 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", NULL)); 122 PetscFunctionReturn(0); 123 } 124 125 /* 126 PCSetCoordinates_AGG 127 - collective 128 129 Input Parameter: 130 . pc - the preconditioner context 131 . ndm - dimesion of data (used for dof/vertex for Stokes) 132 . a_nloc - number of vertices local 133 . coords - [a_nloc][ndm] - interleaved coordinate data: {x_0, y_0, z_0, x_1, y_1, ...} 134 */ 135 136 static PetscErrorCode PCSetCoordinates_AGG(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords) 137 { 138 PC_MG *mg = (PC_MG *)pc->data; 139 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 140 PetscInt arrsz, kk, ii, jj, nloc, ndatarows, ndf; 141 Mat mat = pc->pmat; 142 143 PetscFunctionBegin; 144 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 145 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 146 nloc = a_nloc; 147 148 /* SA: null space vectors */ 149 PetscCall(MatGetBlockSize(mat, &ndf)); /* this does not work for Stokes */ 150 if (coords && ndf == 1) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */ 151 else if (coords) { 152 PetscCheck(ndm <= ndf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "degrees of motion %" PetscInt_FMT " > block size %" PetscInt_FMT, ndm, ndf); 153 pc_gamg->data_cell_cols = (ndm == 2 ? 3 : 6); /* displacement elasticity */ 154 if (ndm != ndf) PetscCheck(pc_gamg->data_cell_cols == ndf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Don't know how to create null space for ndm=%" PetscInt_FMT ", ndf=%" PetscInt_FMT ". Use MatSetNearNullSpace().", ndm, ndf); 155 } else pc_gamg->data_cell_cols = ndf; /* no data, force SA with constant null space vectors */ 156 pc_gamg->data_cell_rows = ndatarows = ndf; 157 PetscCheck(pc_gamg->data_cell_cols > 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "pc_gamg->data_cell_cols %" PetscInt_FMT " <= 0", pc_gamg->data_cell_cols); 158 arrsz = nloc * pc_gamg->data_cell_rows * pc_gamg->data_cell_cols; 159 160 if (!pc_gamg->data || (pc_gamg->data_sz != arrsz)) { 161 PetscCall(PetscFree(pc_gamg->data)); 162 PetscCall(PetscMalloc1(arrsz + 1, &pc_gamg->data)); 163 } 164 /* copy data in - column oriented */ 165 for (kk = 0; kk < nloc; kk++) { 166 const PetscInt M = nloc * pc_gamg->data_cell_rows; /* stride into data */ 167 PetscReal *data = &pc_gamg->data[kk * ndatarows]; /* start of cell */ 168 if (pc_gamg->data_cell_cols == 1) *data = 1.0; 169 else { 170 /* translational modes */ 171 for (ii = 0; ii < ndatarows; ii++) { 172 for (jj = 0; jj < ndatarows; jj++) { 173 if (ii == jj) data[ii * M + jj] = 1.0; 174 else data[ii * M + jj] = 0.0; 175 } 176 } 177 178 /* rotational modes */ 179 if (coords) { 180 if (ndm == 2) { 181 data += 2 * M; 182 data[0] = -coords[2 * kk + 1]; 183 data[1] = coords[2 * kk]; 184 } else { 185 data += 3 * M; 186 data[0] = 0.0; 187 data[M + 0] = coords[3 * kk + 2]; 188 data[2 * M + 0] = -coords[3 * kk + 1]; 189 data[1] = -coords[3 * kk + 2]; 190 data[M + 1] = 0.0; 191 data[2 * M + 1] = coords[3 * kk]; 192 data[2] = coords[3 * kk + 1]; 193 data[M + 2] = -coords[3 * kk]; 194 data[2 * M + 2] = 0.0; 195 } 196 } 197 } 198 } 199 pc_gamg->data_sz = arrsz; 200 PetscFunctionReturn(0); 201 } 202 203 /* 204 PCSetData_AGG - called if data is not set with PCSetCoordinates. 205 Looks in Mat for near null space. 206 Does not work for Stokes 207 208 Input Parameter: 209 . pc - 210 . a_A - matrix to get (near) null space out of. 211 */ 212 static PetscErrorCode PCSetData_AGG(PC pc, Mat a_A) 213 { 214 PC_MG *mg = (PC_MG *)pc->data; 215 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 216 MatNullSpace mnull; 217 218 PetscFunctionBegin; 219 PetscCall(MatGetNearNullSpace(a_A, &mnull)); 220 if (!mnull) { 221 DM dm; 222 PetscCall(PCGetDM(pc, &dm)); 223 if (!dm) PetscCall(MatGetDM(a_A, &dm)); 224 if (dm) { 225 PetscObject deformation; 226 PetscInt Nf; 227 228 PetscCall(DMGetNumFields(dm, &Nf)); 229 if (Nf) { 230 PetscCall(DMGetField(dm, 0, NULL, &deformation)); 231 PetscCall(PetscObjectQuery((PetscObject)deformation, "nearnullspace", (PetscObject *)&mnull)); 232 if (!mnull) PetscCall(PetscObjectQuery((PetscObject)deformation, "nullspace", (PetscObject *)&mnull)); 233 } 234 } 235 } 236 237 if (!mnull) { 238 PetscInt bs, NN, MM; 239 PetscCall(MatGetBlockSize(a_A, &bs)); 240 PetscCall(MatGetLocalSize(a_A, &MM, &NN)); 241 PetscCheck(MM % bs == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MM %" PetscInt_FMT " must be divisible by bs %" PetscInt_FMT, MM, bs); 242 PetscCall(PCSetCoordinates_AGG(pc, bs, MM / bs, NULL)); 243 } else { 244 PetscReal *nullvec; 245 PetscBool has_const; 246 PetscInt i, j, mlocal, nvec, bs; 247 const Vec *vecs; 248 const PetscScalar *v; 249 250 PetscCall(MatGetLocalSize(a_A, &mlocal, NULL)); 251 PetscCall(MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs)); 252 pc_gamg->data_sz = (nvec + !!has_const) * mlocal; 253 PetscCall(PetscMalloc1((nvec + !!has_const) * mlocal, &nullvec)); 254 if (has_const) 255 for (i = 0; i < mlocal; i++) nullvec[i] = 1.0; 256 for (i = 0; i < nvec; i++) { 257 PetscCall(VecGetArrayRead(vecs[i], &v)); 258 for (j = 0; j < mlocal; j++) nullvec[(i + !!has_const) * mlocal + j] = PetscRealPart(v[j]); 259 PetscCall(VecRestoreArrayRead(vecs[i], &v)); 260 } 261 pc_gamg->data = nullvec; 262 pc_gamg->data_cell_cols = (nvec + !!has_const); 263 PetscCall(MatGetBlockSize(a_A, &bs)); 264 pc_gamg->data_cell_rows = bs; 265 } 266 PetscFunctionReturn(0); 267 } 268 269 /* 270 formProl0 - collect null space data for each aggregate, do QR, put R in coarse grid data and Q in P_0 271 272 Input Parameter: 273 . agg_llists - list of arrays with aggregates -- list from selected vertices of aggregate unselected vertices 274 . bs - row block size 275 . nSAvec - column bs of new P 276 . my0crs - global index of start of locals 277 . data_stride - bs*(nloc nodes + ghost nodes) [data_stride][nSAvec] 278 . data_in[data_stride*nSAvec] - local data on fine grid 279 . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist' 280 281 Output Parameter: 282 . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data 283 . a_Prol - prolongation operator 284 */ 285 static PetscErrorCode formProl0(PetscCoarsenData *agg_llists, PetscInt bs, PetscInt nSAvec, PetscInt my0crs, PetscInt data_stride, PetscReal data_in[], const PetscInt flid_fgid[], PetscReal **a_data_out, Mat a_Prol) 286 { 287 PetscInt Istart, my0, Iend, nloc, clid, flid = 0, aggID, kk, jj, ii, mm, nSelected, minsz, nghosts, out_data_stride; 288 MPI_Comm comm; 289 PetscReal *out_data; 290 PetscCDIntNd *pos; 291 PCGAMGHashTable fgid_flid; 292 293 PetscFunctionBegin; 294 PetscCall(PetscObjectGetComm((PetscObject)a_Prol, &comm)); 295 PetscCall(MatGetOwnershipRange(a_Prol, &Istart, &Iend)); 296 nloc = (Iend - Istart) / bs; 297 my0 = Istart / bs; 298 PetscCheck((Iend - Istart) % bs == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Iend %" PetscInt_FMT " - Istart %" PetscInt_FMT " must be divisible by bs %" PetscInt_FMT, Iend, Istart, bs); 299 Iend /= bs; 300 nghosts = data_stride / bs - nloc; 301 302 PetscCall(PCGAMGHashTableCreate(2 * nghosts + 1, &fgid_flid)); 303 for (kk = 0; kk < nghosts; kk++) PetscCall(PCGAMGHashTableAdd(&fgid_flid, flid_fgid[nloc + kk], nloc + kk)); 304 305 /* count selected -- same as number of cols of P */ 306 for (nSelected = mm = 0; mm < nloc; mm++) { 307 PetscBool ise; 308 PetscCall(PetscCDEmptyAt(agg_llists, mm, &ise)); 309 if (!ise) nSelected++; 310 } 311 PetscCall(MatGetOwnershipRangeColumn(a_Prol, &ii, &jj)); 312 PetscCheck((ii / nSAvec) == my0crs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ii %" PetscInt_FMT " /nSAvec %" PetscInt_FMT " != my0crs %" PetscInt_FMT, ii, nSAvec, my0crs); 313 PetscCheck(nSelected == (jj - ii) / nSAvec, PETSC_COMM_SELF, PETSC_ERR_PLIB, "nSelected %" PetscInt_FMT " != (jj %" PetscInt_FMT " - ii %" PetscInt_FMT ")/nSAvec %" PetscInt_FMT, nSelected, jj, ii, nSAvec); 314 315 /* aloc space for coarse point data (output) */ 316 out_data_stride = nSelected * nSAvec; 317 318 PetscCall(PetscMalloc1(out_data_stride * nSAvec, &out_data)); 319 for (ii = 0; ii < out_data_stride * nSAvec; ii++) out_data[ii] = PETSC_MAX_REAL; 320 *a_data_out = out_data; /* output - stride nSelected*nSAvec */ 321 322 /* find points and set prolongation */ 323 minsz = 100; 324 for (mm = clid = 0; mm < nloc; mm++) { 325 PetscCall(PetscCDSizeAt(agg_llists, mm, &jj)); 326 if (jj > 0) { 327 const PetscInt lid = mm, cgid = my0crs + clid; 328 PetscInt cids[100]; /* max bs */ 329 PetscBLASInt asz = jj, M = asz * bs, N = nSAvec, INFO; 330 PetscBLASInt Mdata = M + ((N - M > 0) ? N - M : 0), LDA = Mdata, LWORK = N * bs; 331 PetscScalar *qqc, *qqr, *TAU, *WORK; 332 PetscInt *fids; 333 PetscReal *data; 334 335 /* count agg */ 336 if (asz < minsz) minsz = asz; 337 338 /* get block */ 339 PetscCall(PetscMalloc5(Mdata * N, &qqc, M * N, &qqr, N, &TAU, LWORK, &WORK, M, &fids)); 340 341 aggID = 0; 342 PetscCall(PetscCDGetHeadPos(agg_llists, lid, &pos)); 343 while (pos) { 344 PetscInt gid1; 345 PetscCall(PetscCDIntNdGetID(pos, &gid1)); 346 PetscCall(PetscCDGetNextPos(agg_llists, lid, &pos)); 347 348 if (gid1 >= my0 && gid1 < Iend) flid = gid1 - my0; 349 else { 350 PetscCall(PCGAMGHashTableFind(&fgid_flid, gid1, &flid)); 351 PetscCheck(flid >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cannot find gid1 in table"); 352 } 353 /* copy in B_i matrix - column oriented */ 354 data = &data_in[flid * bs]; 355 for (ii = 0; ii < bs; ii++) { 356 for (jj = 0; jj < N; jj++) { 357 PetscReal d = data[jj * data_stride + ii]; 358 qqc[jj * Mdata + aggID * bs + ii] = d; 359 } 360 } 361 /* set fine IDs */ 362 for (kk = 0; kk < bs; kk++) fids[aggID * bs + kk] = flid_fgid[flid] * bs + kk; 363 aggID++; 364 } 365 366 /* pad with zeros */ 367 for (ii = asz * bs; ii < Mdata; ii++) { 368 for (jj = 0; jj < N; jj++, kk++) qqc[jj * Mdata + ii] = .0; 369 } 370 371 /* QR */ 372 PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 373 PetscCallBLAS("LAPACKgeqrf", LAPACKgeqrf_(&Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO)); 374 PetscCall(PetscFPTrapPop()); 375 PetscCheck(INFO == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xGEQRF error"); 376 /* get R - column oriented - output B_{i+1} */ 377 { 378 PetscReal *data = &out_data[clid * nSAvec]; 379 for (jj = 0; jj < nSAvec; jj++) { 380 for (ii = 0; ii < nSAvec; ii++) { 381 PetscCheck(data[jj * out_data_stride + ii] == PETSC_MAX_REAL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "data[jj*out_data_stride + ii] != %e", (double)PETSC_MAX_REAL); 382 if (ii <= jj) data[jj * out_data_stride + ii] = PetscRealPart(qqc[jj * Mdata + ii]); 383 else data[jj * out_data_stride + ii] = 0.; 384 } 385 } 386 } 387 388 /* get Q - row oriented */ 389 PetscCallBLAS("LAPACKorgqr", LAPACKorgqr_(&Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO)); 390 PetscCheck(INFO == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xORGQR error arg %" PetscBLASInt_FMT, -INFO); 391 392 for (ii = 0; ii < M; ii++) { 393 for (jj = 0; jj < N; jj++) qqr[N * ii + jj] = qqc[jj * Mdata + ii]; 394 } 395 396 /* add diagonal block of P0 */ 397 for (kk = 0; kk < N; kk++) { cids[kk] = N * cgid + kk; /* global col IDs in P0 */ } 398 PetscCall(MatSetValues(a_Prol, M, fids, N, cids, qqr, INSERT_VALUES)); 399 PetscCall(PetscFree5(qqc, qqr, TAU, WORK, fids)); 400 clid++; 401 } /* coarse agg */ 402 } /* for all fine nodes */ 403 PetscCall(MatAssemblyBegin(a_Prol, MAT_FINAL_ASSEMBLY)); 404 PetscCall(MatAssemblyEnd(a_Prol, MAT_FINAL_ASSEMBLY)); 405 PetscCall(PCGAMGHashTableDestroy(&fgid_flid)); 406 PetscFunctionReturn(0); 407 } 408 409 static PetscErrorCode PCView_GAMG_AGG(PC pc, PetscViewer viewer) 410 { 411 PC_MG *mg = (PC_MG *)pc->data; 412 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 413 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 414 415 PetscFunctionBegin; 416 PetscCall(PetscViewerASCIIPrintf(viewer, " AGG specific options\n")); 417 PetscCall(PetscViewerASCIIPrintf(viewer, " Number of levels to square graph %d\n", (int)pc_gamg_agg->aggressive_coarsening_levels)); 418 PetscCall(PetscViewerASCIIPrintf(viewer, " Number smoothing steps %d\n", (int)pc_gamg_agg->nsmooths)); 419 PetscFunctionReturn(0); 420 } 421 422 static PetscErrorCode PCGAMGCreateGraph_AGG(PC pc, Mat Amat, Mat *a_Gmat) 423 { 424 PC_MG *mg = (PC_MG *)pc->data; 425 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 426 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 427 const PetscReal vfilter = pc_gamg->threshold[pc_gamg->current_level]; 428 PetscBool ishem; 429 const char *prefix; 430 431 PetscFunctionBegin; 432 PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_GRAPH], 0, 0, 0, 0)); 433 /* Note: depending on the algorithm that will be used for computing the coarse grid points this should pass PETSC_TRUE or PETSC_FALSE as the first argument */ 434 435 /* MATCOARSENHEM requires numerical weights for edges so ensure they are computed */ 436 PetscCall(MatCoarsenCreate(PetscObjectComm((PetscObject)pc), &pc_gamg_agg->crs)); 437 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc, &prefix)); 438 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc_gamg_agg->crs, prefix)); 439 PetscCall(MatCoarsenSetFromOptions(pc_gamg_agg->crs)); 440 PetscCall(PetscObjectTypeCompare((PetscObject)pc_gamg_agg->crs, MATCOARSENHEM, &ishem)); 441 442 PetscCall(MatCreateGraph(Amat, PETSC_TRUE, (vfilter >= 0 || ishem) ? PETSC_TRUE : PETSC_FALSE, vfilter, a_Gmat)); 443 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_GRAPH], 0, 0, 0, 0)); 444 PetscFunctionReturn(0); 445 } 446 447 /* 448 PCGAMGCoarsen_AGG - supports squaring the graph (deprecated) and new graph for 449 communication of QR data used with HEM and MISk coarsening 450 451 Input Parameter: 452 . a_pc - this 453 454 Input/Output Parameter: 455 . a_Gmat1 - graph to coarsen (in), graph off processor edges for QR gather scatter (out) 456 457 Output Parameter: 458 . agg_lists - list of aggregates 459 460 */ 461 static PetscErrorCode PCGAMGCoarsen_AGG(PC a_pc, Mat *a_Gmat1, PetscCoarsenData **agg_lists) 462 { 463 PC_MG *mg = (PC_MG *)a_pc->data; 464 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 465 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 466 Mat mat, Gmat1 = *a_Gmat1; /* aggressive graph */ 467 IS perm; 468 PetscInt Istart, Iend, Ii, nloc, bs, nn; 469 PetscInt *permute, *degree; 470 PetscBool *bIndexSet; 471 PetscReal hashfact; 472 PetscInt iSwapIndex; 473 PetscRandom random; 474 475 PetscFunctionBegin; 476 PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0)); 477 PetscCall(MatGetLocalSize(Gmat1, &nn, NULL)); 478 PetscCall(MatGetBlockSize(Gmat1, &bs)); 479 PetscCheck(bs == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "bs %" PetscInt_FMT " must be 1", bs); 480 nloc = nn / bs; 481 482 /* get MIS aggs - randomize */ 483 PetscCall(PetscMalloc2(nloc, &permute, nloc, °ree)); 484 PetscCall(PetscCalloc1(nloc, &bIndexSet)); 485 for (Ii = 0; Ii < nloc; Ii++) permute[Ii] = Ii; 486 PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &random)); 487 PetscCall(MatGetOwnershipRange(Gmat1, &Istart, &Iend)); 488 for (Ii = 0; Ii < nloc; Ii++) { 489 PetscInt nc; 490 PetscCall(MatGetRow(Gmat1, Istart + Ii, &nc, NULL, NULL)); 491 degree[Ii] = nc; 492 PetscCall(MatRestoreRow(Gmat1, Istart + Ii, &nc, NULL, NULL)); 493 } 494 for (Ii = 0; Ii < nloc; Ii++) { 495 PetscCall(PetscRandomGetValueReal(random, &hashfact)); 496 iSwapIndex = (PetscInt)(hashfact * nloc) % nloc; 497 if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) { 498 PetscInt iTemp = permute[iSwapIndex]; 499 permute[iSwapIndex] = permute[Ii]; 500 permute[Ii] = iTemp; 501 iTemp = degree[iSwapIndex]; 502 degree[iSwapIndex] = degree[Ii]; 503 degree[Ii] = iTemp; 504 bIndexSet[iSwapIndex] = PETSC_TRUE; 505 } 506 } 507 // create minimum degree ordering 508 PetscCall(PetscSortIntWithArray(nloc, degree, permute)); 509 510 PetscCall(PetscFree(bIndexSet)); 511 PetscCall(PetscRandomDestroy(&random)); 512 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm)); 513 PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_MIS], 0, 0, 0, 0)); 514 PetscCall(MatCoarsenSetGreedyOrdering(pc_gamg_agg->crs, perm)); 515 PetscCall(MatCoarsenSetAdjacency(pc_gamg_agg->crs, Gmat1)); 516 PetscCall(MatCoarsenSetStrictAggs(pc_gamg_agg->crs, PETSC_TRUE)); 517 if (pc_gamg->current_level < pc_gamg_agg->aggressive_coarsening_levels) PetscCall(MatCoarsenMISKSetDistance(pc_gamg_agg->crs, 2)); // hardwire to MIS-2 518 else PetscCall(MatCoarsenMISKSetDistance(pc_gamg_agg->crs, 1)); // MIS 519 PetscCall(MatCoarsenApply(pc_gamg_agg->crs)); 520 PetscCall(MatCoarsenViewFromOptions(pc_gamg_agg->crs, NULL, "-mat_coarsen_view")); 521 PetscCall(MatCoarsenGetData(pc_gamg_agg->crs, agg_lists)); /* output */ 522 PetscCall(MatCoarsenDestroy(&pc_gamg_agg->crs)); 523 524 PetscCall(ISDestroy(&perm)); 525 PetscCall(PetscFree2(permute, degree)); 526 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_MIS], 0, 0, 0, 0)); 527 528 { 529 PetscCoarsenData *llist = *agg_lists; 530 /* see if we have a matrix that takes precedence (returned from MatCoarsenApply) */ 531 PetscCall(PetscCDGetMat(llist, &mat)); 532 if (mat) { 533 PetscCall(MatDestroy(&Gmat1)); 534 *a_Gmat1 = mat; /* output */ 535 } 536 } 537 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0)); 538 PetscFunctionReturn(0); 539 } 540 541 /* 542 PCGAMGProlongator_AGG 543 544 Input Parameter: 545 . pc - this 546 . Amat - matrix on this fine level 547 . Graph - used to get ghost data for nodes in 548 . agg_lists - list of aggregates 549 Output Parameter: 550 . a_P_out - prolongation operator to the next level 551 */ 552 static PetscErrorCode PCGAMGProlongator_AGG(PC pc, Mat Amat, Mat Gmat, PetscCoarsenData *agg_lists, Mat *a_P_out) 553 { 554 PC_MG *mg = (PC_MG *)pc->data; 555 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 556 const PetscInt col_bs = pc_gamg->data_cell_cols; 557 PetscInt Istart, Iend, nloc, ii, jj, kk, my0, nLocalSelected, bs; 558 Mat Prol; 559 PetscMPIInt size; 560 MPI_Comm comm; 561 PetscReal *data_w_ghost; 562 PetscInt myCrs0, nbnodes = 0, *flid_fgid; 563 MatType mtype; 564 565 PetscFunctionBegin; 566 PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm)); 567 PetscCheck(col_bs >= 1, comm, PETSC_ERR_PLIB, "Column bs cannot be less than 1"); 568 PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0)); 569 PetscCallMPI(MPI_Comm_size(comm, &size)); 570 PetscCall(MatGetOwnershipRange(Amat, &Istart, &Iend)); 571 PetscCall(MatGetBlockSize(Amat, &bs)); 572 nloc = (Iend - Istart) / bs; 573 my0 = Istart / bs; 574 PetscCheck((Iend - Istart) % bs == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "(Iend %" PetscInt_FMT " - Istart %" PetscInt_FMT ") not divisible by bs %" PetscInt_FMT, Iend, Istart, bs); 575 576 /* get 'nLocalSelected' */ 577 for (ii = 0, nLocalSelected = 0; ii < nloc; ii++) { 578 PetscBool ise; 579 /* filter out singletons 0 or 1? */ 580 PetscCall(PetscCDEmptyAt(agg_lists, ii, &ise)); 581 if (!ise) nLocalSelected++; 582 } 583 584 /* create prolongator, create P matrix */ 585 PetscCall(MatGetType(Amat, &mtype)); 586 PetscCall(MatCreate(comm, &Prol)); 587 PetscCall(MatSetSizes(Prol, nloc * bs, nLocalSelected * col_bs, PETSC_DETERMINE, PETSC_DETERMINE)); 588 PetscCall(MatSetBlockSizes(Prol, bs, col_bs)); 589 PetscCall(MatSetType(Prol, mtype)); 590 PetscCall(MatSeqAIJSetPreallocation(Prol, col_bs, NULL)); 591 PetscCall(MatMPIAIJSetPreallocation(Prol, col_bs, NULL, col_bs, NULL)); 592 593 /* can get all points "removed" */ 594 PetscCall(MatGetSize(Prol, &kk, &ii)); 595 if (!ii) { 596 PetscCall(PetscInfo(pc, "%s: No selected points on coarse grid\n", ((PetscObject)pc)->prefix)); 597 PetscCall(MatDestroy(&Prol)); 598 *a_P_out = NULL; /* out */ 599 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0)); 600 PetscFunctionReturn(0); 601 } 602 PetscCall(PetscInfo(pc, "%s: New grid %" PetscInt_FMT " nodes\n", ((PetscObject)pc)->prefix, ii / col_bs)); 603 PetscCall(MatGetOwnershipRangeColumn(Prol, &myCrs0, &kk)); 604 605 PetscCheck((kk - myCrs0) % col_bs == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "(kk %" PetscInt_FMT " -myCrs0 %" PetscInt_FMT ") not divisible by col_bs %" PetscInt_FMT, kk, myCrs0, col_bs); 606 myCrs0 = myCrs0 / col_bs; 607 PetscCheck((kk / col_bs - myCrs0) == nLocalSelected, PETSC_COMM_SELF, PETSC_ERR_PLIB, "(kk %" PetscInt_FMT "/col_bs %" PetscInt_FMT " - myCrs0 %" PetscInt_FMT ") != nLocalSelected %" PetscInt_FMT ")", kk, col_bs, myCrs0, nLocalSelected); 608 609 /* create global vector of data in 'data_w_ghost' */ 610 PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROLA], 0, 0, 0, 0)); 611 if (size > 1) { /* get ghost null space data */ 612 PetscReal *tmp_gdata, *tmp_ldata, *tp2; 613 PetscCall(PetscMalloc1(nloc, &tmp_ldata)); 614 for (jj = 0; jj < col_bs; jj++) { 615 for (kk = 0; kk < bs; kk++) { 616 PetscInt ii, stride; 617 const PetscReal *tp = pc_gamg->data + jj * bs * nloc + kk; 618 for (ii = 0; ii < nloc; ii++, tp += bs) tmp_ldata[ii] = *tp; 619 620 PetscCall(PCGAMGGetDataWithGhosts(Gmat, 1, tmp_ldata, &stride, &tmp_gdata)); 621 622 if (!jj && !kk) { /* now I know how many total nodes - allocate TODO: move below and do in one 'col_bs' call */ 623 PetscCall(PetscMalloc1(stride * bs * col_bs, &data_w_ghost)); 624 nbnodes = bs * stride; 625 } 626 tp2 = data_w_ghost + jj * bs * stride + kk; 627 for (ii = 0; ii < stride; ii++, tp2 += bs) *tp2 = tmp_gdata[ii]; 628 PetscCall(PetscFree(tmp_gdata)); 629 } 630 } 631 PetscCall(PetscFree(tmp_ldata)); 632 } else { 633 nbnodes = bs * nloc; 634 data_w_ghost = (PetscReal *)pc_gamg->data; 635 } 636 637 /* get 'flid_fgid' TODO - move up to get 'stride' and do get null space data above in one step (jj loop) */ 638 if (size > 1) { 639 PetscReal *fid_glid_loc, *fiddata; 640 PetscInt stride; 641 642 PetscCall(PetscMalloc1(nloc, &fid_glid_loc)); 643 for (kk = 0; kk < nloc; kk++) fid_glid_loc[kk] = (PetscReal)(my0 + kk); 644 PetscCall(PCGAMGGetDataWithGhosts(Gmat, 1, fid_glid_loc, &stride, &fiddata)); 645 PetscCall(PetscMalloc1(stride, &flid_fgid)); /* copy real data to in */ 646 for (kk = 0; kk < stride; kk++) flid_fgid[kk] = (PetscInt)fiddata[kk]; 647 PetscCall(PetscFree(fiddata)); 648 649 PetscCheck(stride == nbnodes / bs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "stride %" PetscInt_FMT " != nbnodes %" PetscInt_FMT "/bs %" PetscInt_FMT, stride, nbnodes, bs); 650 PetscCall(PetscFree(fid_glid_loc)); 651 } else { 652 PetscCall(PetscMalloc1(nloc, &flid_fgid)); 653 for (kk = 0; kk < nloc; kk++) flid_fgid[kk] = my0 + kk; 654 } 655 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROLA], 0, 0, 0, 0)); 656 /* get P0 */ 657 PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROLB], 0, 0, 0, 0)); 658 { 659 PetscReal *data_out = NULL; 660 PetscCall(formProl0(agg_lists, bs, col_bs, myCrs0, nbnodes, data_w_ghost, flid_fgid, &data_out, Prol)); 661 PetscCall(PetscFree(pc_gamg->data)); 662 663 pc_gamg->data = data_out; 664 pc_gamg->data_cell_rows = col_bs; 665 pc_gamg->data_sz = col_bs * col_bs * nLocalSelected; 666 } 667 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROLB], 0, 0, 0, 0)); 668 if (size > 1) PetscCall(PetscFree(data_w_ghost)); 669 PetscCall(PetscFree(flid_fgid)); 670 671 *a_P_out = Prol; /* out */ 672 673 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0)); 674 PetscFunctionReturn(0); 675 } 676 677 /* 678 PCGAMGOptProlongator_AGG 679 680 Input Parameter: 681 . pc - this 682 . Amat - matrix on this fine level 683 In/Output Parameter: 684 . a_P - prolongation operator to the next level 685 */ 686 static PetscErrorCode PCGAMGOptProlongator_AGG(PC pc, Mat Amat, Mat *a_P) 687 { 688 PC_MG *mg = (PC_MG *)pc->data; 689 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 690 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 691 PetscInt jj; 692 Mat Prol = *a_P; 693 MPI_Comm comm; 694 KSP eksp; 695 Vec bb, xx; 696 PC epc; 697 PetscReal alpha, emax, emin; 698 699 PetscFunctionBegin; 700 PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm)); 701 PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_OPT], 0, 0, 0, 0)); 702 703 /* compute maximum singular value of operator to be used in smoother */ 704 if (0 < pc_gamg_agg->nsmooths) { 705 /* get eigen estimates */ 706 if (pc_gamg->emax > 0) { 707 emin = pc_gamg->emin; 708 emax = pc_gamg->emax; 709 } else { 710 const char *prefix; 711 712 PetscCall(MatCreateVecs(Amat, &bb, NULL)); 713 PetscCall(MatCreateVecs(Amat, &xx, NULL)); 714 PetscCall(KSPSetNoisy_Private(bb)); 715 716 PetscCall(KSPCreate(comm, &eksp)); 717 PetscCall(PCGetOptionsPrefix(pc, &prefix)); 718 PetscCall(KSPSetOptionsPrefix(eksp, prefix)); 719 PetscCall(KSPAppendOptionsPrefix(eksp, "pc_gamg_esteig_")); 720 { 721 PetscBool isset, sflg; 722 PetscCall(MatIsSPDKnown(Amat, &isset, &sflg)); 723 if (isset && sflg) PetscCall(KSPSetType(eksp, KSPCG)); 724 } 725 PetscCall(KSPSetErrorIfNotConverged(eksp, pc->erroriffailure)); 726 PetscCall(KSPSetNormType(eksp, KSP_NORM_NONE)); 727 728 PetscCall(KSPSetInitialGuessNonzero(eksp, PETSC_FALSE)); 729 PetscCall(KSPSetOperators(eksp, Amat, Amat)); 730 731 PetscCall(KSPGetPC(eksp, &epc)); 732 PetscCall(PCSetType(epc, PCJACOBI)); /* smoother in smoothed agg. */ 733 734 PetscCall(KSPSetTolerances(eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10)); // 10 is safer, but 5 is often fine, can override with -pc_gamg_esteig_ksp_max_it -mg_levels_ksp_chebyshev_esteig 0,0.25,0,1.2 735 736 PetscCall(KSPSetFromOptions(eksp)); 737 PetscCall(KSPSetComputeSingularValues(eksp, PETSC_TRUE)); 738 PetscCall(KSPSolve(eksp, bb, xx)); 739 PetscCall(KSPCheckSolve(eksp, pc, xx)); 740 741 PetscCall(KSPComputeExtremeSingularValues(eksp, &emax, &emin)); 742 PetscCall(PetscInfo(pc, "%s: Smooth P0: max eigen=%e min=%e PC=%s\n", ((PetscObject)pc)->prefix, (double)emax, (double)emin, PCJACOBI)); 743 PetscCall(VecDestroy(&xx)); 744 PetscCall(VecDestroy(&bb)); 745 PetscCall(KSPDestroy(&eksp)); 746 } 747 if (pc_gamg->use_sa_esteig) { 748 mg->min_eigen_DinvA[pc_gamg->current_level] = emin; 749 mg->max_eigen_DinvA[pc_gamg->current_level] = emax; 750 PetscCall(PetscInfo(pc, "%s: Smooth P0: level %" PetscInt_FMT ", cache spectra %g %g\n", ((PetscObject)pc)->prefix, pc_gamg->current_level, (double)emin, (double)emax)); 751 } else { 752 mg->min_eigen_DinvA[pc_gamg->current_level] = 0; 753 mg->max_eigen_DinvA[pc_gamg->current_level] = 0; 754 } 755 } else { 756 mg->min_eigen_DinvA[pc_gamg->current_level] = 0; 757 mg->max_eigen_DinvA[pc_gamg->current_level] = 0; 758 } 759 760 /* smooth P0 */ 761 for (jj = 0; jj < pc_gamg_agg->nsmooths; jj++) { 762 Mat tMat; 763 Vec diag; 764 765 PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_OPTSM], 0, 0, 0, 0)); 766 767 /* smooth P1 := (I - omega/lam D^{-1}A)P0 */ 768 PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2], 0, 0, 0, 0)); 769 PetscCall(MatMatMult(Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat)); 770 PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2], 0, 0, 0, 0)); 771 PetscCall(MatProductClear(tMat)); 772 PetscCall(MatCreateVecs(Amat, &diag, NULL)); 773 PetscCall(MatGetDiagonal(Amat, diag)); /* effectively PCJACOBI */ 774 PetscCall(VecReciprocal(diag)); 775 PetscCall(MatDiagonalScale(tMat, diag, NULL)); 776 PetscCall(VecDestroy(&diag)); 777 778 /* TODO: Set a PCFailedReason and exit the building of the AMG preconditioner */ 779 PetscCheck(emax != 0.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_PLIB, "Computed maximum singular value as zero"); 780 /* TODO: Document the 1.4 and don't hardwire it in this routine */ 781 alpha = -1.4 / emax; 782 783 PetscCall(MatAYPX(tMat, alpha, Prol, SUBSET_NONZERO_PATTERN)); 784 PetscCall(MatDestroy(&Prol)); 785 Prol = tMat; 786 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_OPTSM], 0, 0, 0, 0)); 787 } 788 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_OPT], 0, 0, 0, 0)); 789 *a_P = Prol; 790 PetscFunctionReturn(0); 791 } 792 793 /* 794 PCCreateGAMG_AGG 795 796 Input Parameter: 797 . pc - 798 */ 799 PetscErrorCode PCCreateGAMG_AGG(PC pc) 800 { 801 PC_MG *mg = (PC_MG *)pc->data; 802 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 803 PC_GAMG_AGG *pc_gamg_agg; 804 805 PetscFunctionBegin; 806 /* create sub context for SA */ 807 PetscCall(PetscNew(&pc_gamg_agg)); 808 pc_gamg->subctx = pc_gamg_agg; 809 810 pc_gamg->ops->setfromoptions = PCSetFromOptions_GAMG_AGG; 811 pc_gamg->ops->destroy = PCDestroy_GAMG_AGG; 812 /* reset does not do anything; setup not virtual */ 813 814 /* set internal function pointers */ 815 pc_gamg->ops->creategraph = PCGAMGCreateGraph_AGG; 816 pc_gamg->ops->coarsen = PCGAMGCoarsen_AGG; 817 pc_gamg->ops->prolongator = PCGAMGProlongator_AGG; 818 pc_gamg->ops->optprolongator = PCGAMGOptProlongator_AGG; 819 pc_gamg->ops->createdefaultdata = PCSetData_AGG; 820 pc_gamg->ops->view = PCView_GAMG_AGG; 821 822 pc_gamg_agg->aggressive_coarsening_levels = 0; 823 pc_gamg_agg->nsmooths = 1; 824 825 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNSmooths_C", PCGAMGSetNSmooths_AGG)); 826 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveLevels_C", PCGAMGSetAggressiveLevels_AGG)); 827 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", PCSetCoordinates_AGG)); 828 PetscFunctionReturn(0); 829 } 830