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