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; // number of smoothing steps to construct prolongation 12 PetscInt aggressive_coarsening_levels; // number of aggressive coarsening levels (square or MISk) 13 PetscInt aggressive_mis_k; // the k in MIS-k 14 PetscBool use_aggressive_square_graph; 15 PetscBool use_minimum_degree_ordering; 16 PetscBool use_low_mem_filter; 17 MatCoarsen crs; 18 } PC_GAMG_AGG; 19 20 /*@ 21 PCGAMGSetNSmooths - Set number of smoothing steps (1 is typical) used to construct the prolongation operator 22 23 Logically Collective 24 25 Input Parameters: 26 + pc - the preconditioner context 27 - n - the number of smooths 28 29 Options Database Key: 30 . -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use 31 32 Level: intermediate 33 34 Note: 35 This is a different concept from the number smoothing steps used during the linear solution process which 36 can be set with `-mg_levels_ksp_max_it` 37 38 Developer Note: 39 This should be named `PCGAMGAGGSetNSmooths()`. 40 41 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCMG`, `PCGAMG` 42 @*/ 43 PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n) 44 { 45 PetscFunctionBegin; 46 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 47 PetscValidLogicalCollectiveInt(pc, n, 2); 48 PetscTryMethod(pc, "PCGAMGSetNSmooths_C", (PC, PetscInt), (pc, n)); 49 PetscFunctionReturn(PETSC_SUCCESS); 50 } 51 52 static PetscErrorCode PCGAMGSetNSmooths_AGG(PC pc, PetscInt n) 53 { 54 PC_MG *mg = (PC_MG *)pc->data; 55 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 56 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 57 58 PetscFunctionBegin; 59 pc_gamg_agg->nsmooths = n; 60 PetscFunctionReturn(PETSC_SUCCESS); 61 } 62 63 /*@ 64 PCGAMGSetAggressiveLevels - Use aggressive coarsening on first n levels 65 66 Logically Collective 67 68 Input Parameters: 69 + pc - the preconditioner context 70 - n - 0, 1 or more 71 72 Options Database Key: 73 . -pc_gamg_aggressive_coarsening <n,default = 1> - Number of levels on which to square the graph on before aggregating it 74 75 Level: intermediate 76 77 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `PCGAMGSetThreshold()`, `PCGAMGMISkSetAggressive()`, `PCGAMGSetAggressiveSquareGraph()`, `PCGAMGMISkSetMinDegreeOrdering()`, `PCGAMGSetLowMemoryFilter()` 78 @*/ 79 PetscErrorCode PCGAMGSetAggressiveLevels(PC pc, PetscInt n) 80 { 81 PetscFunctionBegin; 82 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 83 PetscValidLogicalCollectiveInt(pc, n, 2); 84 PetscTryMethod(pc, "PCGAMGSetAggressiveLevels_C", (PC, PetscInt), (pc, n)); 85 PetscFunctionReturn(PETSC_SUCCESS); 86 } 87 88 /*@ 89 PCGAMGMISkSetAggressive - Number (k) distance in MIS coarsening (>2 is 'aggressive') 90 91 Logically Collective 92 93 Input Parameters: 94 + pc - the preconditioner context 95 - n - 1 or more (default = 2) 96 97 Options Database Key: 98 . -pc_gamg_aggressive_mis_k <n,default=2> - Number (k) distance in MIS coarsening (>2 is 'aggressive') 99 100 Level: intermediate 101 102 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `PCGAMGSetThreshold()`, `PCGAMGSetAggressiveLevels()`, `PCGAMGSetAggressiveSquareGraph()`, `PCGAMGMISkSetMinDegreeOrdering()`, `PCGAMGSetLowMemoryFilter()` 103 @*/ 104 PetscErrorCode PCGAMGMISkSetAggressive(PC pc, PetscInt n) 105 { 106 PetscFunctionBegin; 107 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 108 PetscValidLogicalCollectiveInt(pc, n, 2); 109 PetscTryMethod(pc, "PCGAMGMISkSetAggressive_C", (PC, PetscInt), (pc, n)); 110 PetscFunctionReturn(PETSC_SUCCESS); 111 } 112 113 /*@ 114 PCGAMGSetAggressiveSquareGraph - Use graph square A'A for aggressive coarsening, old method 115 116 Logically Collective 117 118 Input Parameters: 119 + pc - the preconditioner context 120 - b - default false - MIS-k is faster 121 122 Options Database Key: 123 . -pc_gamg_aggressive_square_graph <bool,default=false> - Use square graph (A'A) or MIS-k (k=2) for aggressive coarsening 124 125 Level: intermediate 126 127 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `PCGAMGSetThreshold()`, `PCGAMGSetAggressiveLevels()`, `PCGAMGMISkSetAggressive()`, `PCGAMGMISkSetMinDegreeOrdering()`, `PCGAMGSetLowMemoryFilter()` 128 @*/ 129 PetscErrorCode PCGAMGSetAggressiveSquareGraph(PC pc, PetscBool b) 130 { 131 PetscFunctionBegin; 132 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 133 PetscValidLogicalCollectiveBool(pc, b, 2); 134 PetscTryMethod(pc, "PCGAMGSetAggressiveSquareGraph_C", (PC, PetscBool), (pc, b)); 135 PetscFunctionReturn(PETSC_SUCCESS); 136 } 137 138 /*@ 139 PCGAMGMISkSetMinDegreeOrdering - Use minimum degree ordering in greedy MIS algorithm 140 141 Logically Collective 142 143 Input Parameters: 144 + pc - the preconditioner context 145 - b - default true 146 147 Options Database Key: 148 . -pc_gamg_mis_k_minimum_degree_ordering <bool,default=true> - Use minimum degree ordering in greedy MIS algorithm 149 150 Level: intermediate 151 152 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCGAMG`, `PCGAMGSetThreshold()`, `PCGAMGSetAggressiveLevels()`, `PCGAMGMISkSetAggressive()`, `PCGAMGSetAggressiveSquareGraph()`, `PCGAMGSetLowMemoryFilter()` 153 @*/ 154 PetscErrorCode PCGAMGMISkSetMinDegreeOrdering(PC pc, PetscBool b) 155 { 156 PetscFunctionBegin; 157 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 158 PetscValidLogicalCollectiveBool(pc, b, 2); 159 PetscTryMethod(pc, "PCGAMGMISkSetMinDegreeOrdering_C", (PC, PetscBool), (pc, b)); 160 PetscFunctionReturn(PETSC_SUCCESS); 161 } 162 163 /*@ 164 PCGAMGSetLowMemoryFilter - Use low memory graph/matrix filter 165 166 Logically Collective 167 168 Input Parameters: 169 + pc - the preconditioner context 170 - b - default false 171 172 Options Database Key: 173 . -pc_gamg_low_memory_threshold_filter <bool,default=false> - Use low memory graph/matrix filter 174 175 Level: intermediate 176 177 .seealso: [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), `PCGAMG`, `PCGAMGSetThreshold()`, `PCGAMGSetAggressiveLevels()`, 178 `PCGAMGMISkSetAggressive()`, `PCGAMGSetAggressiveSquareGraph()`, `PCGAMGMISkSetMinDegreeOrdering()` 179 @*/ 180 PetscErrorCode PCGAMGSetLowMemoryFilter(PC pc, PetscBool b) 181 { 182 PetscFunctionBegin; 183 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 184 PetscValidLogicalCollectiveBool(pc, b, 2); 185 PetscTryMethod(pc, "PCGAMGSetLowMemoryFilter_C", (PC, PetscBool), (pc, b)); 186 PetscFunctionReturn(PETSC_SUCCESS); 187 } 188 189 static PetscErrorCode PCGAMGSetAggressiveLevels_AGG(PC pc, PetscInt n) 190 { 191 PC_MG *mg = (PC_MG *)pc->data; 192 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 193 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 194 195 PetscFunctionBegin; 196 pc_gamg_agg->aggressive_coarsening_levels = n; 197 PetscFunctionReturn(PETSC_SUCCESS); 198 } 199 200 static PetscErrorCode PCGAMGMISkSetAggressive_AGG(PC pc, PetscInt n) 201 { 202 PC_MG *mg = (PC_MG *)pc->data; 203 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 204 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 205 206 PetscFunctionBegin; 207 pc_gamg_agg->aggressive_mis_k = n; 208 PetscFunctionReturn(PETSC_SUCCESS); 209 } 210 211 static PetscErrorCode PCGAMGSetAggressiveSquareGraph_AGG(PC pc, PetscBool b) 212 { 213 PC_MG *mg = (PC_MG *)pc->data; 214 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 215 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 216 217 PetscFunctionBegin; 218 pc_gamg_agg->use_aggressive_square_graph = b; 219 PetscFunctionReturn(PETSC_SUCCESS); 220 } 221 222 static PetscErrorCode PCGAMGSetLowMemoryFilter_AGG(PC pc, PetscBool b) 223 { 224 PC_MG *mg = (PC_MG *)pc->data; 225 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 226 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 227 228 PetscFunctionBegin; 229 pc_gamg_agg->use_low_mem_filter = b; 230 PetscFunctionReturn(PETSC_SUCCESS); 231 } 232 233 static PetscErrorCode PCGAMGMISkSetMinDegreeOrdering_AGG(PC pc, PetscBool b) 234 { 235 PC_MG *mg = (PC_MG *)pc->data; 236 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 237 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 238 239 PetscFunctionBegin; 240 pc_gamg_agg->use_minimum_degree_ordering = b; 241 PetscFunctionReturn(PETSC_SUCCESS); 242 } 243 244 static PetscErrorCode PCSetFromOptions_GAMG_AGG(PC pc, PetscOptionItems *PetscOptionsObject) 245 { 246 PC_MG *mg = (PC_MG *)pc->data; 247 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 248 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 249 PetscBool n_aggressive_flg, old_sq_provided = PETSC_FALSE, new_sq_provided = PETSC_FALSE, new_sqr_graph = pc_gamg_agg->use_aggressive_square_graph; 250 PetscInt nsq_graph_old = 0; 251 252 PetscFunctionBegin; 253 PetscOptionsHeadBegin(PetscOptionsObject, "GAMG-AGG options"); 254 PetscCall(PetscOptionsInt("-pc_gamg_agg_nsmooths", "number of smoothing steps to construct prolongation, usually 1", "PCGAMGSetNSmooths", pc_gamg_agg->nsmooths, &pc_gamg_agg->nsmooths, NULL)); 255 // aggressive coarsening logic with deprecated -pc_gamg_square_graph 256 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, &n_aggressive_flg)); 257 if (!n_aggressive_flg) 258 PetscCall(PetscOptionsInt("-pc_gamg_square_graph", "Number of aggressive coarsening (MIS-2) levels from finest (deprecated alias for -pc_gamg_aggressive_coarsening)", "PCGAMGSetAggressiveLevels", nsq_graph_old, &nsq_graph_old, &old_sq_provided)); 259 PetscCall(PetscOptionsBool("-pc_gamg_aggressive_square_graph", "Use square graph (A'A) or MIS-k (k=2) for aggressive coarsening", "PCGAMGSetAggressiveSquareGraph", new_sqr_graph, &pc_gamg_agg->use_aggressive_square_graph, &new_sq_provided)); 260 if (!new_sq_provided && old_sq_provided) { 261 pc_gamg_agg->aggressive_coarsening_levels = nsq_graph_old; // could be zero 262 pc_gamg_agg->use_aggressive_square_graph = PETSC_TRUE; 263 } 264 if (new_sq_provided && old_sq_provided) 265 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)); 266 PetscCall(PetscOptionsBool("-pc_gamg_mis_k_minimum_degree_ordering", "Use minimum degree ordering for greedy MIS", "PCGAMGMISkSetMinDegreeOrdering", pc_gamg_agg->use_minimum_degree_ordering, &pc_gamg_agg->use_minimum_degree_ordering, NULL)); 267 PetscCall(PetscOptionsBool("-pc_gamg_low_memory_threshold_filter", "Use the (built-in) low memory graph/matrix filter", "PCGAMGSetLowMemoryFilter", pc_gamg_agg->use_low_mem_filter, &pc_gamg_agg->use_low_mem_filter, NULL)); 268 PetscCall(PetscOptionsInt("-pc_gamg_aggressive_mis_k", "Number of levels of multigrid to use.", "PCGAMGMISkSetAggressive", pc_gamg_agg->aggressive_mis_k, &pc_gamg_agg->aggressive_mis_k, NULL)); 269 PetscOptionsHeadEnd(); 270 PetscFunctionReturn(PETSC_SUCCESS); 271 } 272 273 static PetscErrorCode PCDestroy_GAMG_AGG(PC pc) 274 { 275 PC_MG *mg = (PC_MG *)pc->data; 276 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 277 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 278 279 PetscFunctionBegin; 280 PetscCall(MatCoarsenDestroy(&pc_gamg_agg->crs)); 281 PetscCall(PetscFree(pc_gamg->subctx)); 282 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNSmooths_C", NULL)); 283 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveLevels_C", NULL)); 284 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetAggressive_C", NULL)); 285 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetMinDegreeOrdering_C", NULL)); 286 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetLowMemoryFilter_C", NULL)); 287 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveSquareGraph_C", NULL)); 288 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", NULL)); 289 PetscFunctionReturn(PETSC_SUCCESS); 290 } 291 292 /* 293 PCSetCoordinates_AGG 294 295 Collective 296 297 Input Parameter: 298 . pc - the preconditioner context 299 . ndm - dimension of data (used for dof/vertex for Stokes) 300 . a_nloc - number of vertices local 301 . coords - [a_nloc][ndm] - interleaved coordinate data: {x_0, y_0, z_0, x_1, y_1, ...} 302 */ 303 304 static PetscErrorCode PCSetCoordinates_AGG(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords) 305 { 306 PC_MG *mg = (PC_MG *)pc->data; 307 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 308 PetscInt arrsz, kk, ii, jj, nloc, ndatarows, ndf; 309 Mat mat = pc->pmat; 310 311 PetscFunctionBegin; 312 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 313 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 314 nloc = a_nloc; 315 316 /* SA: null space vectors */ 317 PetscCall(MatGetBlockSize(mat, &ndf)); /* this does not work for Stokes */ 318 if (coords && ndf == 1) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */ 319 else if (coords) { 320 PetscCheck(ndm <= ndf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "degrees of motion %" PetscInt_FMT " > block size %" PetscInt_FMT, ndm, ndf); 321 pc_gamg->data_cell_cols = (ndm == 2 ? 3 : 6); /* displacement elasticity */ 322 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); 323 } else pc_gamg->data_cell_cols = ndf; /* no data, force SA with constant null space vectors */ 324 pc_gamg->data_cell_rows = ndatarows = ndf; 325 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); 326 arrsz = nloc * pc_gamg->data_cell_rows * pc_gamg->data_cell_cols; 327 328 if (!pc_gamg->data || (pc_gamg->data_sz != arrsz)) { 329 PetscCall(PetscFree(pc_gamg->data)); 330 PetscCall(PetscMalloc1(arrsz + 1, &pc_gamg->data)); 331 } 332 /* copy data in - column-oriented */ 333 for (kk = 0; kk < nloc; kk++) { 334 const PetscInt M = nloc * pc_gamg->data_cell_rows; /* stride into data */ 335 PetscReal *data = &pc_gamg->data[kk * ndatarows]; /* start of cell */ 336 337 if (pc_gamg->data_cell_cols == 1) *data = 1.0; 338 else { 339 /* translational modes */ 340 for (ii = 0; ii < ndatarows; ii++) { 341 for (jj = 0; jj < ndatarows; jj++) { 342 if (ii == jj) data[ii * M + jj] = 1.0; 343 else data[ii * M + jj] = 0.0; 344 } 345 } 346 347 /* rotational modes */ 348 if (coords) { 349 if (ndm == 2) { 350 data += 2 * M; 351 data[0] = -coords[2 * kk + 1]; 352 data[1] = coords[2 * kk]; 353 } else { 354 data += 3 * M; 355 data[0] = 0.0; 356 data[M + 0] = coords[3 * kk + 2]; 357 data[2 * M + 0] = -coords[3 * kk + 1]; 358 data[1] = -coords[3 * kk + 2]; 359 data[M + 1] = 0.0; 360 data[2 * M + 1] = coords[3 * kk]; 361 data[2] = coords[3 * kk + 1]; 362 data[M + 2] = -coords[3 * kk]; 363 data[2 * M + 2] = 0.0; 364 } 365 } 366 } 367 } 368 pc_gamg->data_sz = arrsz; 369 PetscFunctionReturn(PETSC_SUCCESS); 370 } 371 372 /* 373 PCSetData_AGG - called if data is not set with PCSetCoordinates. 374 Looks in Mat for near null space. 375 Does not work for Stokes 376 377 Input Parameter: 378 . pc - 379 . a_A - matrix to get (near) null space out of. 380 */ 381 static PetscErrorCode PCSetData_AGG(PC pc, Mat a_A) 382 { 383 PC_MG *mg = (PC_MG *)pc->data; 384 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 385 MatNullSpace mnull; 386 387 PetscFunctionBegin; 388 PetscCall(MatGetNearNullSpace(a_A, &mnull)); 389 if (!mnull) { 390 DM dm; 391 392 PetscCall(PCGetDM(pc, &dm)); 393 if (!dm) PetscCall(MatGetDM(a_A, &dm)); 394 if (dm) { 395 PetscObject deformation; 396 PetscInt Nf; 397 398 PetscCall(DMGetNumFields(dm, &Nf)); 399 if (Nf) { 400 PetscCall(DMGetField(dm, 0, NULL, &deformation)); 401 PetscCall(PetscObjectQuery((PetscObject)deformation, "nearnullspace", (PetscObject *)&mnull)); 402 if (!mnull) PetscCall(PetscObjectQuery((PetscObject)deformation, "nullspace", (PetscObject *)&mnull)); 403 } 404 } 405 } 406 407 if (!mnull) { 408 PetscInt bs, NN, MM; 409 410 PetscCall(MatGetBlockSize(a_A, &bs)); 411 PetscCall(MatGetLocalSize(a_A, &MM, &NN)); 412 PetscCheck(MM % bs == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MM %" PetscInt_FMT " must be divisible by bs %" PetscInt_FMT, MM, bs); 413 PetscCall(PCSetCoordinates_AGG(pc, bs, MM / bs, NULL)); 414 } else { 415 PetscReal *nullvec; 416 PetscBool has_const; 417 PetscInt i, j, mlocal, nvec, bs; 418 const Vec *vecs; 419 const PetscScalar *v; 420 421 PetscCall(MatGetLocalSize(a_A, &mlocal, NULL)); 422 PetscCall(MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs)); 423 for (i = 0; i < nvec; i++) { 424 PetscCall(VecGetLocalSize(vecs[i], &j)); 425 PetscCheck(j == mlocal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Attached null space vector size %" PetscInt_FMT " != matrix size %" PetscInt_FMT, j, mlocal); 426 } 427 pc_gamg->data_sz = (nvec + !!has_const) * mlocal; 428 PetscCall(PetscMalloc1((nvec + !!has_const) * mlocal, &nullvec)); 429 if (has_const) 430 for (i = 0; i < mlocal; i++) nullvec[i] = 1.0; 431 for (i = 0; i < nvec; i++) { 432 PetscCall(VecGetArrayRead(vecs[i], &v)); 433 for (j = 0; j < mlocal; j++) nullvec[(i + !!has_const) * mlocal + j] = PetscRealPart(v[j]); 434 PetscCall(VecRestoreArrayRead(vecs[i], &v)); 435 } 436 pc_gamg->data = nullvec; 437 pc_gamg->data_cell_cols = (nvec + !!has_const); 438 PetscCall(MatGetBlockSize(a_A, &bs)); 439 pc_gamg->data_cell_rows = bs; 440 } 441 PetscFunctionReturn(PETSC_SUCCESS); 442 } 443 444 /* 445 formProl0 - collect null space data for each aggregate, do QR, put R in coarse grid data and Q in P_0 446 447 Input Parameter: 448 . agg_llists - list of arrays with aggregates -- list from selected vertices of aggregate unselected vertices 449 . bs - row block size 450 . nSAvec - column bs of new P 451 . my0crs - global index of start of locals 452 . data_stride - bs*(nloc nodes + ghost nodes) [data_stride][nSAvec] 453 . data_in[data_stride*nSAvec] - local data on fine grid 454 . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist' 455 456 Output Parameter: 457 . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data 458 . a_Prol - prolongation operator 459 */ 460 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) 461 { 462 PetscInt Istart, my0, Iend, nloc, clid, flid = 0, aggID, kk, jj, ii, mm, nSelected, minsz, nghosts, out_data_stride; 463 MPI_Comm comm; 464 PetscReal *out_data; 465 PetscCDIntNd *pos; 466 PCGAMGHashTable fgid_flid; 467 468 PetscFunctionBegin; 469 PetscCall(PetscObjectGetComm((PetscObject)a_Prol, &comm)); 470 PetscCall(MatGetOwnershipRange(a_Prol, &Istart, &Iend)); 471 nloc = (Iend - Istart) / bs; 472 my0 = Istart / bs; 473 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); 474 Iend /= bs; 475 nghosts = data_stride / bs - nloc; 476 477 PetscCall(PCGAMGHashTableCreate(2 * nghosts + 1, &fgid_flid)); 478 for (kk = 0; kk < nghosts; kk++) PetscCall(PCGAMGHashTableAdd(&fgid_flid, flid_fgid[nloc + kk], nloc + kk)); 479 480 /* count selected -- same as number of cols of P */ 481 for (nSelected = mm = 0; mm < nloc; mm++) { 482 PetscBool ise; 483 484 PetscCall(PetscCDIsEmptyAt(agg_llists, mm, &ise)); 485 if (!ise) nSelected++; 486 } 487 PetscCall(MatGetOwnershipRangeColumn(a_Prol, &ii, &jj)); 488 PetscCheck((ii / nSAvec) == my0crs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ii %" PetscInt_FMT " /nSAvec %" PetscInt_FMT " != my0crs %" PetscInt_FMT, ii, nSAvec, my0crs); 489 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); 490 491 /* aloc space for coarse point data (output) */ 492 out_data_stride = nSelected * nSAvec; 493 494 PetscCall(PetscMalloc1(out_data_stride * nSAvec, &out_data)); 495 for (ii = 0; ii < out_data_stride * nSAvec; ii++) out_data[ii] = PETSC_MAX_REAL; 496 *a_data_out = out_data; /* output - stride nSelected*nSAvec */ 497 498 /* find points and set prolongation */ 499 minsz = 100; 500 for (mm = clid = 0; mm < nloc; mm++) { 501 PetscCall(PetscCDCountAt(agg_llists, mm, &jj)); 502 if (jj > 0) { 503 const PetscInt lid = mm, cgid = my0crs + clid; 504 PetscInt cids[100]; /* max bs */ 505 PetscBLASInt asz = jj, M = asz * bs, N = nSAvec, INFO; 506 PetscBLASInt Mdata = M + ((N - M > 0) ? N - M : 0), LDA = Mdata, LWORK = N * bs; 507 PetscScalar *qqc, *qqr, *TAU, *WORK; 508 PetscInt *fids; 509 PetscReal *data; 510 511 /* count agg */ 512 if (asz < minsz) minsz = asz; 513 514 /* get block */ 515 PetscCall(PetscMalloc5(Mdata * N, &qqc, M * N, &qqr, N, &TAU, LWORK, &WORK, M, &fids)); 516 517 aggID = 0; 518 PetscCall(PetscCDGetHeadPos(agg_llists, lid, &pos)); 519 while (pos) { 520 PetscInt gid1; 521 522 PetscCall(PetscCDIntNdGetID(pos, &gid1)); 523 PetscCall(PetscCDGetNextPos(agg_llists, lid, &pos)); 524 525 if (gid1 >= my0 && gid1 < Iend) flid = gid1 - my0; 526 else { 527 PetscCall(PCGAMGHashTableFind(&fgid_flid, gid1, &flid)); 528 PetscCheck(flid >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cannot find gid1 in table"); 529 } 530 /* copy in B_i matrix - column-oriented */ 531 data = &data_in[flid * bs]; 532 for (ii = 0; ii < bs; ii++) { 533 for (jj = 0; jj < N; jj++) { 534 PetscReal d = data[jj * data_stride + ii]; 535 536 qqc[jj * Mdata + aggID * bs + ii] = d; 537 } 538 } 539 /* set fine IDs */ 540 for (kk = 0; kk < bs; kk++) fids[aggID * bs + kk] = flid_fgid[flid] * bs + kk; 541 aggID++; 542 } 543 544 /* pad with zeros */ 545 for (ii = asz * bs; ii < Mdata; ii++) { 546 for (jj = 0; jj < N; jj++, kk++) qqc[jj * Mdata + ii] = .0; 547 } 548 549 /* QR */ 550 PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 551 PetscCallBLAS("LAPACKgeqrf", LAPACKgeqrf_(&Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO)); 552 PetscCall(PetscFPTrapPop()); 553 PetscCheck(INFO == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xGEQRF error"); 554 /* get R - column-oriented - output B_{i+1} */ 555 { 556 PetscReal *data = &out_data[clid * nSAvec]; 557 558 for (jj = 0; jj < nSAvec; jj++) { 559 for (ii = 0; ii < nSAvec; ii++) { 560 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); 561 if (ii <= jj) data[jj * out_data_stride + ii] = PetscRealPart(qqc[jj * Mdata + ii]); 562 else data[jj * out_data_stride + ii] = 0.; 563 } 564 } 565 } 566 567 /* get Q - row-oriented */ 568 PetscCallBLAS("LAPACKorgqr", LAPACKorgqr_(&Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO)); 569 PetscCheck(INFO == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xORGQR error arg %" PetscBLASInt_FMT, -INFO); 570 571 for (ii = 0; ii < M; ii++) { 572 for (jj = 0; jj < N; jj++) qqr[N * ii + jj] = qqc[jj * Mdata + ii]; 573 } 574 575 /* add diagonal block of P0 */ 576 for (kk = 0; kk < N; kk++) { cids[kk] = N * cgid + kk; /* global col IDs in P0 */ } 577 PetscCall(MatSetValues(a_Prol, M, fids, N, cids, qqr, INSERT_VALUES)); 578 PetscCall(PetscFree5(qqc, qqr, TAU, WORK, fids)); 579 clid++; 580 } /* coarse agg */ 581 } /* for all fine nodes */ 582 PetscCall(MatAssemblyBegin(a_Prol, MAT_FINAL_ASSEMBLY)); 583 PetscCall(MatAssemblyEnd(a_Prol, MAT_FINAL_ASSEMBLY)); 584 PetscCall(PCGAMGHashTableDestroy(&fgid_flid)); 585 PetscFunctionReturn(PETSC_SUCCESS); 586 } 587 588 static PetscErrorCode PCView_GAMG_AGG(PC pc, PetscViewer viewer) 589 { 590 PC_MG *mg = (PC_MG *)pc->data; 591 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 592 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 593 594 PetscFunctionBegin; 595 PetscCall(PetscViewerASCIIPrintf(viewer, " AGG specific options\n")); 596 PetscCall(PetscViewerASCIIPrintf(viewer, " Number of levels of aggressive coarsening %d\n", (int)pc_gamg_agg->aggressive_coarsening_levels)); 597 if (pc_gamg_agg->aggressive_coarsening_levels > 0) { 598 PetscCall(PetscViewerASCIIPrintf(viewer, " %s aggressive coarsening\n", !pc_gamg_agg->use_aggressive_square_graph ? "MIS-k" : "Square graph")); 599 if (!pc_gamg_agg->use_aggressive_square_graph) PetscCall(PetscViewerASCIIPrintf(viewer, " MIS-%d coarsening on aggressive levels\n", (int)pc_gamg_agg->aggressive_mis_k)); 600 } 601 PetscCall(PetscViewerASCIIPushTab(viewer)); 602 PetscCall(PetscViewerASCIIPushTab(viewer)); 603 PetscCall(PetscViewerASCIIPushTab(viewer)); 604 PetscCall(PetscViewerASCIIPushTab(viewer)); 605 if (pc_gamg_agg->crs) PetscCall(MatCoarsenView(pc_gamg_agg->crs, viewer)); 606 else PetscCall(PetscViewerASCIIPrintf(viewer, "Coarsening algorithm not yet selected\n")); 607 PetscCall(PetscViewerASCIIPopTab(viewer)); 608 PetscCall(PetscViewerASCIIPopTab(viewer)); 609 PetscCall(PetscViewerASCIIPopTab(viewer)); 610 PetscCall(PetscViewerASCIIPopTab(viewer)); 611 PetscCall(PetscViewerASCIIPrintf(viewer, " Number smoothing steps to construct prolongation %d\n", (int)pc_gamg_agg->nsmooths)); 612 PetscFunctionReturn(PETSC_SUCCESS); 613 } 614 615 static PetscErrorCode PCGAMGCreateGraph_AGG(PC pc, Mat Amat, Mat *a_Gmat) 616 { 617 PC_MG *mg = (PC_MG *)pc->data; 618 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 619 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 620 const PetscReal vfilter = pc_gamg->threshold[pc_gamg->current_level]; 621 PetscBool ishem, ismis; 622 const char *prefix; 623 MatInfo info0, info1; 624 PetscInt bs; 625 626 PetscFunctionBegin; 627 PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0)); 628 /* 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 */ 629 /* MATCOARSENHEM requires numerical weights for edges so ensure they are computed */ 630 PetscCall(MatCoarsenDestroy(&pc_gamg_agg->crs)); 631 PetscCall(MatCoarsenCreate(PetscObjectComm((PetscObject)pc), &pc_gamg_agg->crs)); 632 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc, &prefix)); 633 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc_gamg_agg->crs, prefix)); 634 PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)pc_gamg_agg->crs, "pc_gamg_")); 635 PetscCall(MatCoarsenSetFromOptions(pc_gamg_agg->crs)); 636 PetscCall(MatGetBlockSize(Amat, &bs)); 637 // check for valid indices wrt bs 638 for (int ii = 0; ii < pc_gamg_agg->crs->strength_index_size; ii++) { 639 PetscCheck(pc_gamg_agg->crs->strength_index[ii] >= 0 && pc_gamg_agg->crs->strength_index[ii] < bs, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Indices (%d) must be non-negative and < block size (%d), NB, can not use -mat_coarsen_strength_index with -mat_coarsen_strength_index", 640 (int)pc_gamg_agg->crs->strength_index[ii], (int)bs); 641 } 642 PetscCall(PetscObjectTypeCompare((PetscObject)pc_gamg_agg->crs, MATCOARSENHEM, &ishem)); 643 if (ishem) { 644 if (pc_gamg_agg->aggressive_coarsening_levels) PetscCall(PetscInfo(pc, "HEM and aggressive coarsening ignored: HEM using %d iterations\n", (int)pc_gamg_agg->crs->max_it)); 645 pc_gamg_agg->aggressive_coarsening_levels = 0; // aggressive and HEM does not make sense 646 PetscCall(MatCoarsenSetMaximumIterations(pc_gamg_agg->crs, pc_gamg_agg->crs->max_it)); // for code coverage 647 PetscCall(MatCoarsenSetThreshold(pc_gamg_agg->crs, vfilter)); // for code coverage 648 } else { 649 PetscCall(PetscObjectTypeCompare((PetscObject)pc_gamg_agg->crs, MATCOARSENMIS, &ismis)); 650 if (ismis && pc_gamg_agg->aggressive_coarsening_levels && !pc_gamg_agg->use_aggressive_square_graph) { 651 PetscCall(PetscInfo(pc, "MIS and aggressive coarsening and no square graph: force square graph\n")); 652 pc_gamg_agg->use_aggressive_square_graph = PETSC_TRUE; 653 } 654 } 655 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0)); 656 PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_GRAPH], 0, 0, 0, 0)); 657 PetscCall(MatGetInfo(Amat, MAT_LOCAL, &info0)); /* global reduction */ 658 659 if (ishem || pc_gamg_agg->use_low_mem_filter) { 660 PetscCall(MatCreateGraph(Amat, PETSC_TRUE, (vfilter >= 0 || ishem) ? PETSC_TRUE : PETSC_FALSE, vfilter, pc_gamg_agg->crs->strength_index_size, pc_gamg_agg->crs->strength_index, a_Gmat)); 661 } else { 662 // make scalar graph, symetrize if not know to be symmetric, scale, but do not filter (expensive) 663 PetscCall(MatCreateGraph(Amat, PETSC_TRUE, PETSC_TRUE, -1, pc_gamg_agg->crs->strength_index_size, pc_gamg_agg->crs->strength_index, a_Gmat)); 664 if (vfilter >= 0) { 665 PetscInt Istart, Iend, ncols, nnz0, nnz1, NN, MM, nloc; 666 Mat tGmat, Gmat = *a_Gmat; 667 MPI_Comm comm; 668 const PetscScalar *vals; 669 const PetscInt *idx; 670 PetscInt *d_nnz, *o_nnz, kk, *garray = NULL, *AJ, maxcols = 0; 671 MatScalar *AA; // this is checked in graph 672 PetscBool isseqaij; 673 Mat a, b, c; 674 MatType jtype; 675 676 PetscCall(PetscObjectGetComm((PetscObject)Gmat, &comm)); 677 PetscCall(PetscObjectBaseTypeCompare((PetscObject)Gmat, MATSEQAIJ, &isseqaij)); 678 PetscCall(MatGetType(Gmat, &jtype)); 679 PetscCall(MatCreate(comm, &tGmat)); 680 PetscCall(MatSetType(tGmat, jtype)); 681 682 /* TODO GPU: this can be called when filter = 0 -> Probably provide MatAIJThresholdCompress that compresses the entries below a threshold? 683 Also, if the matrix is symmetric, can we skip this 684 operation? It can be very expensive on large matrices. */ 685 686 // global sizes 687 PetscCall(MatGetSize(Gmat, &MM, &NN)); 688 PetscCall(MatGetOwnershipRange(Gmat, &Istart, &Iend)); 689 nloc = Iend - Istart; 690 PetscCall(PetscMalloc2(nloc, &d_nnz, nloc, &o_nnz)); 691 if (isseqaij) { 692 a = Gmat; 693 b = NULL; 694 } else { 695 Mat_MPIAIJ *d = (Mat_MPIAIJ *)Gmat->data; 696 697 a = d->A; 698 b = d->B; 699 garray = d->garray; 700 } 701 /* Determine upper bound on non-zeros needed in new filtered matrix */ 702 for (PetscInt row = 0; row < nloc; row++) { 703 PetscCall(MatGetRow(a, row, &ncols, NULL, NULL)); 704 d_nnz[row] = ncols; 705 if (ncols > maxcols) maxcols = ncols; 706 PetscCall(MatRestoreRow(a, row, &ncols, NULL, NULL)); 707 } 708 if (b) { 709 for (PetscInt row = 0; row < nloc; row++) { 710 PetscCall(MatGetRow(b, row, &ncols, NULL, NULL)); 711 o_nnz[row] = ncols; 712 if (ncols > maxcols) maxcols = ncols; 713 PetscCall(MatRestoreRow(b, row, &ncols, NULL, NULL)); 714 } 715 } 716 PetscCall(MatSetSizes(tGmat, nloc, nloc, MM, MM)); 717 PetscCall(MatSetBlockSizes(tGmat, 1, 1)); 718 PetscCall(MatSeqAIJSetPreallocation(tGmat, 0, d_nnz)); 719 PetscCall(MatMPIAIJSetPreallocation(tGmat, 0, d_nnz, 0, o_nnz)); 720 PetscCall(MatSetOption(tGmat, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE)); 721 PetscCall(PetscFree2(d_nnz, o_nnz)); 722 PetscCall(PetscMalloc2(maxcols, &AA, maxcols, &AJ)); 723 nnz0 = nnz1 = 0; 724 for (c = a, kk = 0; c && kk < 2; c = b, kk++) { 725 for (PetscInt row = 0, grow = Istart, ncol_row, jj; row < nloc; row++, grow++) { 726 PetscCall(MatGetRow(c, row, &ncols, &idx, &vals)); 727 for (ncol_row = jj = 0; jj < ncols; jj++, nnz0++) { 728 PetscScalar sv = PetscAbs(PetscRealPart(vals[jj])); 729 if (PetscRealPart(sv) > vfilter) { 730 PetscInt cid = idx[jj] + Istart; //diag 731 732 nnz1++; 733 if (c != a) cid = garray[idx[jj]]; 734 AA[ncol_row] = vals[jj]; 735 AJ[ncol_row] = cid; 736 ncol_row++; 737 } 738 } 739 PetscCall(MatRestoreRow(c, row, &ncols, &idx, &vals)); 740 PetscCall(MatSetValues(tGmat, 1, &grow, ncol_row, AJ, AA, INSERT_VALUES)); 741 } 742 } 743 PetscCall(PetscFree2(AA, AJ)); 744 PetscCall(MatAssemblyBegin(tGmat, MAT_FINAL_ASSEMBLY)); 745 PetscCall(MatAssemblyEnd(tGmat, MAT_FINAL_ASSEMBLY)); 746 PetscCall(MatPropagateSymmetryOptions(Gmat, tGmat)); /* Normal Mat options are not relevant ? */ 747 PetscCall(PetscInfo(pc, "\t %g%% nnz after filtering, with threshold %g, %g nnz ave. (N=%" PetscInt_FMT ", max row size %" PetscInt_FMT "\n", (!nnz0) ? 1. : 100. * (double)nnz1 / (double)nnz0, (double)vfilter, (!nloc) ? 1. : (double)nnz0 / (double)nloc, MM, maxcols)); 748 PetscCall(MatViewFromOptions(tGmat, NULL, "-mat_filter_graph_view")); 749 PetscCall(MatDestroy(&Gmat)); 750 *a_Gmat = tGmat; 751 } 752 } 753 754 PetscCall(MatGetInfo(*a_Gmat, MAT_LOCAL, &info1)); /* global reduction */ 755 if (info0.nz_used > 0) PetscCall(PetscInfo(pc, "Filtering left %g %% edges in graph (%e %e)\n", 100.0 * info1.nz_used * (double)(bs * bs) / info0.nz_used, info0.nz_used, info1.nz_used)); 756 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_GRAPH], 0, 0, 0, 0)); 757 PetscFunctionReturn(PETSC_SUCCESS); 758 } 759 760 typedef PetscInt NState; 761 static const NState NOT_DONE = -2; 762 static const NState DELETED = -1; 763 static const NState REMOVED = -3; 764 #define IS_SELECTED(s) (s != DELETED && s != NOT_DONE && s != REMOVED) 765 766 /* 767 fixAggregatesWithSquare - greedy grab of with G1 (unsquared graph) -- AIJ specific -- change to fixAggregatesWithSquare -- TODD 768 - AGG-MG specific: clears singletons out of 'selected_2' 769 770 Input Parameter: 771 . Gmat_2 - global matrix of squared graph (data not defined) 772 . Gmat_1 - base graph to grab with base graph 773 Input/Output Parameter: 774 . aggs_2 - linked list of aggs with gids) 775 */ 776 static PetscErrorCode fixAggregatesWithSquare(PC pc, Mat Gmat_2, Mat Gmat_1, PetscCoarsenData *aggs_2) 777 { 778 PetscBool isMPI; 779 Mat_SeqAIJ *matA_1, *matB_1 = NULL; 780 MPI_Comm comm; 781 PetscInt lid, *ii, *idx, ix, Iend, my0, kk, n, j; 782 Mat_MPIAIJ *mpimat_2 = NULL, *mpimat_1 = NULL; 783 const PetscInt nloc = Gmat_2->rmap->n; 784 PetscScalar *cpcol_1_state, *cpcol_2_state, *cpcol_2_par_orig, *lid_parent_gid; 785 PetscInt *lid_cprowID_1 = NULL; 786 NState *lid_state; 787 Vec ghost_par_orig2; 788 PetscMPIInt rank; 789 790 PetscFunctionBegin; 791 PetscCall(PetscObjectGetComm((PetscObject)Gmat_2, &comm)); 792 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 793 PetscCall(MatGetOwnershipRange(Gmat_1, &my0, &Iend)); 794 795 /* get submatrices */ 796 PetscCall(PetscStrbeginswith(((PetscObject)Gmat_1)->type_name, MATMPIAIJ, &isMPI)); 797 PetscCall(PetscInfo(pc, "isMPI = %s\n", isMPI ? "yes" : "no")); 798 PetscCall(PetscMalloc3(nloc, &lid_state, nloc, &lid_parent_gid, nloc, &lid_cprowID_1)); 799 for (lid = 0; lid < nloc; lid++) lid_cprowID_1[lid] = -1; 800 if (isMPI) { 801 /* grab matrix objects */ 802 mpimat_2 = (Mat_MPIAIJ *)Gmat_2->data; 803 mpimat_1 = (Mat_MPIAIJ *)Gmat_1->data; 804 matA_1 = (Mat_SeqAIJ *)mpimat_1->A->data; 805 matB_1 = (Mat_SeqAIJ *)mpimat_1->B->data; 806 807 /* force compressed row storage for B matrix in AuxMat */ 808 PetscCall(MatCheckCompressedRow(mpimat_1->B, matB_1->nonzerorowcnt, &matB_1->compressedrow, matB_1->i, Gmat_1->rmap->n, -1.0)); 809 for (ix = 0; ix < matB_1->compressedrow.nrows; ix++) { 810 PetscInt lid = matB_1->compressedrow.rindex[ix]; 811 812 PetscCheck(lid <= nloc && lid >= -1, PETSC_COMM_SELF, PETSC_ERR_USER, "lid %d out of range. nloc = %d", (int)lid, (int)nloc); 813 if (lid != -1) lid_cprowID_1[lid] = ix; 814 } 815 } else { 816 PetscBool isAIJ; 817 818 PetscCall(PetscStrbeginswith(((PetscObject)Gmat_1)->type_name, MATSEQAIJ, &isAIJ)); 819 PetscCheck(isAIJ, PETSC_COMM_SELF, PETSC_ERR_USER, "Require AIJ matrix."); 820 matA_1 = (Mat_SeqAIJ *)Gmat_1->data; 821 } 822 if (nloc > 0) { PetscCheck(!matB_1 || matB_1->compressedrow.use, PETSC_COMM_SELF, PETSC_ERR_PLIB, "matB_1 && !matB_1->compressedrow.use: PETSc bug???"); } 823 /* get state of locals and selected gid for deleted */ 824 for (lid = 0; lid < nloc; lid++) { 825 lid_parent_gid[lid] = -1.0; 826 lid_state[lid] = DELETED; 827 } 828 829 /* set lid_state */ 830 for (lid = 0; lid < nloc; lid++) { 831 PetscCDIntNd *pos; 832 833 PetscCall(PetscCDGetHeadPos(aggs_2, lid, &pos)); 834 if (pos) { 835 PetscInt gid1; 836 837 PetscCall(PetscCDIntNdGetID(pos, &gid1)); 838 PetscCheck(gid1 == lid + my0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "gid1 %d != lid %d + my0 %d", (int)gid1, (int)lid, (int)my0); 839 lid_state[lid] = gid1; 840 } 841 } 842 843 /* map local to selected local, DELETED means a ghost owns it */ 844 for (lid = kk = 0; lid < nloc; lid++) { 845 NState state = lid_state[lid]; 846 if (IS_SELECTED(state)) { 847 PetscCDIntNd *pos; 848 849 PetscCall(PetscCDGetHeadPos(aggs_2, lid, &pos)); 850 while (pos) { 851 PetscInt gid1; 852 853 PetscCall(PetscCDIntNdGetID(pos, &gid1)); 854 PetscCall(PetscCDGetNextPos(aggs_2, lid, &pos)); 855 if (gid1 >= my0 && gid1 < Iend) lid_parent_gid[gid1 - my0] = (PetscScalar)(lid + my0); 856 } 857 } 858 } 859 /* get 'cpcol_1/2_state' & cpcol_2_par_orig - uses mpimat_1/2->lvec for temp space */ 860 if (isMPI) { 861 Vec tempVec; 862 /* get 'cpcol_1_state' */ 863 PetscCall(MatCreateVecs(Gmat_1, &tempVec, NULL)); 864 for (kk = 0, j = my0; kk < nloc; kk++, j++) { 865 PetscScalar v = (PetscScalar)lid_state[kk]; 866 867 PetscCall(VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES)); 868 } 869 PetscCall(VecAssemblyBegin(tempVec)); 870 PetscCall(VecAssemblyEnd(tempVec)); 871 PetscCall(VecScatterBegin(mpimat_1->Mvctx, tempVec, mpimat_1->lvec, INSERT_VALUES, SCATTER_FORWARD)); 872 PetscCall(VecScatterEnd(mpimat_1->Mvctx, tempVec, mpimat_1->lvec, INSERT_VALUES, SCATTER_FORWARD)); 873 PetscCall(VecGetArray(mpimat_1->lvec, &cpcol_1_state)); 874 /* get 'cpcol_2_state' */ 875 PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, mpimat_2->lvec, INSERT_VALUES, SCATTER_FORWARD)); 876 PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, mpimat_2->lvec, INSERT_VALUES, SCATTER_FORWARD)); 877 PetscCall(VecGetArray(mpimat_2->lvec, &cpcol_2_state)); 878 /* get 'cpcol_2_par_orig' */ 879 for (kk = 0, j = my0; kk < nloc; kk++, j++) { 880 PetscScalar v = (PetscScalar)lid_parent_gid[kk]; 881 882 PetscCall(VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES)); 883 } 884 PetscCall(VecAssemblyBegin(tempVec)); 885 PetscCall(VecAssemblyEnd(tempVec)); 886 PetscCall(VecDuplicate(mpimat_2->lvec, &ghost_par_orig2)); 887 PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, ghost_par_orig2, INSERT_VALUES, SCATTER_FORWARD)); 888 PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, ghost_par_orig2, INSERT_VALUES, SCATTER_FORWARD)); 889 PetscCall(VecGetArray(ghost_par_orig2, &cpcol_2_par_orig)); 890 891 PetscCall(VecDestroy(&tempVec)); 892 } /* ismpi */ 893 for (lid = 0; lid < nloc; lid++) { 894 NState state = lid_state[lid]; 895 896 if (IS_SELECTED(state)) { 897 /* steal locals */ 898 ii = matA_1->i; 899 n = ii[lid + 1] - ii[lid]; 900 idx = matA_1->j + ii[lid]; 901 for (j = 0; j < n; j++) { 902 PetscInt lidj = idx[j], sgid; 903 NState statej = lid_state[lidj]; 904 905 if (statej == DELETED && (sgid = (PetscInt)PetscRealPart(lid_parent_gid[lidj])) != lid + my0) { /* steal local */ 906 lid_parent_gid[lidj] = (PetscScalar)(lid + my0); /* send this if sgid is not local */ 907 if (sgid >= my0 && sgid < Iend) { /* I'm stealing this local from a local sgid */ 908 PetscInt hav = 0, slid = sgid - my0, gidj = lidj + my0; 909 PetscCDIntNd *pos, *last = NULL; 910 911 /* looking for local from local so id_llist_2 works */ 912 PetscCall(PetscCDGetHeadPos(aggs_2, slid, &pos)); 913 while (pos) { 914 PetscInt gid; 915 PetscCall(PetscCDIntNdGetID(pos, &gid)); 916 if (gid == gidj) { 917 PetscCheck(last, PETSC_COMM_SELF, PETSC_ERR_PLIB, "last cannot be null"); 918 PetscCall(PetscCDRemoveNextNode(aggs_2, slid, last)); 919 PetscCall(PetscCDAppendNode(aggs_2, lid, pos)); 920 hav = 1; 921 break; 922 } else last = pos; 923 PetscCall(PetscCDGetNextPos(aggs_2, slid, &pos)); 924 } 925 if (hav != 1) { 926 PetscCheck(hav, PETSC_COMM_SELF, PETSC_ERR_PLIB, "failed to find adj in 'selected' lists - structurally unsymmetric matrix"); 927 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "found node %d times???", (int)hav); 928 } 929 } else { /* I'm stealing this local, owned by a ghost */ 930 PetscCheck(sgid == -1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Mat has an un-symmetric graph. Use '-%spc_gamg_sym_graph true' to symmetrize the graph or '-%spc_gamg_threshold -1' if the matrix is structurally symmetric.", 931 ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "", ((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : ""); 932 PetscCall(PetscCDAppendID(aggs_2, lid, lidj + my0)); 933 } 934 } 935 } /* local neighbors */ 936 } else if (state == DELETED /* && lid_cprowID_1 */) { 937 PetscInt sgidold = (PetscInt)PetscRealPart(lid_parent_gid[lid]); 938 939 /* see if I have a selected ghost neighbor that will steal me */ 940 if ((ix = lid_cprowID_1[lid]) != -1) { 941 ii = matB_1->compressedrow.i; 942 n = ii[ix + 1] - ii[ix]; 943 idx = matB_1->j + ii[ix]; 944 for (j = 0; j < n; j++) { 945 PetscInt cpid = idx[j]; 946 NState statej = (NState)PetscRealPart(cpcol_1_state[cpid]); 947 948 if (IS_SELECTED(statej) && sgidold != (PetscInt)statej) { /* ghost will steal this, remove from my list */ 949 lid_parent_gid[lid] = (PetscScalar)statej; /* send who selected */ 950 if (sgidold >= my0 && sgidold < Iend) { /* this was mine */ 951 PetscInt hav = 0, oldslidj = sgidold - my0; 952 PetscCDIntNd *pos, *last = NULL; 953 954 /* remove from 'oldslidj' list */ 955 PetscCall(PetscCDGetHeadPos(aggs_2, oldslidj, &pos)); 956 while (pos) { 957 PetscInt gid; 958 959 PetscCall(PetscCDIntNdGetID(pos, &gid)); 960 if (lid + my0 == gid) { 961 /* id_llist_2[lastid] = id_llist_2[flid]; /\* remove lid from oldslidj list *\/ */ 962 PetscCheck(last, PETSC_COMM_SELF, PETSC_ERR_PLIB, "last cannot be null"); 963 PetscCall(PetscCDRemoveNextNode(aggs_2, oldslidj, last)); 964 /* ghost (PetscScalar)statej will add this later */ 965 hav = 1; 966 break; 967 } else last = pos; 968 PetscCall(PetscCDGetNextPos(aggs_2, oldslidj, &pos)); 969 } 970 if (hav != 1) { 971 PetscCheck(hav, PETSC_COMM_SELF, PETSC_ERR_PLIB, "failed to find (hav=%d) adj in 'selected' lists - structurally unsymmetric matrix", (int)hav); 972 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "found node %d times???", (int)hav); 973 } 974 } else { 975 /* TODO: ghosts remove this later */ 976 } 977 } 978 } 979 } 980 } /* selected/deleted */ 981 } /* node loop */ 982 983 if (isMPI) { 984 PetscScalar *cpcol_2_parent, *cpcol_2_gid; 985 Vec tempVec, ghostgids2, ghostparents2; 986 PetscInt cpid, nghost_2; 987 PCGAMGHashTable gid_cpid; 988 989 PetscCall(VecGetSize(mpimat_2->lvec, &nghost_2)); 990 PetscCall(MatCreateVecs(Gmat_2, &tempVec, NULL)); 991 992 /* get 'cpcol_2_parent' */ 993 for (kk = 0, j = my0; kk < nloc; kk++, j++) { PetscCall(VecSetValues(tempVec, 1, &j, &lid_parent_gid[kk], INSERT_VALUES)); } 994 PetscCall(VecAssemblyBegin(tempVec)); 995 PetscCall(VecAssemblyEnd(tempVec)); 996 PetscCall(VecDuplicate(mpimat_2->lvec, &ghostparents2)); 997 PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, ghostparents2, INSERT_VALUES, SCATTER_FORWARD)); 998 PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, ghostparents2, INSERT_VALUES, SCATTER_FORWARD)); 999 PetscCall(VecGetArray(ghostparents2, &cpcol_2_parent)); 1000 1001 /* get 'cpcol_2_gid' */ 1002 for (kk = 0, j = my0; kk < nloc; kk++, j++) { 1003 PetscScalar v = (PetscScalar)j; 1004 1005 PetscCall(VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES)); 1006 } 1007 PetscCall(VecAssemblyBegin(tempVec)); 1008 PetscCall(VecAssemblyEnd(tempVec)); 1009 PetscCall(VecDuplicate(mpimat_2->lvec, &ghostgids2)); 1010 PetscCall(VecScatterBegin(mpimat_2->Mvctx, tempVec, ghostgids2, INSERT_VALUES, SCATTER_FORWARD)); 1011 PetscCall(VecScatterEnd(mpimat_2->Mvctx, tempVec, ghostgids2, INSERT_VALUES, SCATTER_FORWARD)); 1012 PetscCall(VecGetArray(ghostgids2, &cpcol_2_gid)); 1013 PetscCall(VecDestroy(&tempVec)); 1014 1015 /* look for deleted ghosts and add to table */ 1016 PetscCall(PCGAMGHashTableCreate(2 * nghost_2 + 1, &gid_cpid)); 1017 for (cpid = 0; cpid < nghost_2; cpid++) { 1018 NState state = (NState)PetscRealPart(cpcol_2_state[cpid]); 1019 1020 if (state == DELETED) { 1021 PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]); 1022 PetscInt sgid_old = (PetscInt)PetscRealPart(cpcol_2_par_orig[cpid]); 1023 1024 if (sgid_old == -1 && sgid_new != -1) { 1025 PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]); 1026 1027 PetscCall(PCGAMGHashTableAdd(&gid_cpid, gid, cpid)); 1028 } 1029 } 1030 } 1031 1032 /* look for deleted ghosts and see if they moved - remove it */ 1033 for (lid = 0; lid < nloc; lid++) { 1034 NState state = lid_state[lid]; 1035 if (IS_SELECTED(state)) { 1036 PetscCDIntNd *pos, *last = NULL; 1037 /* look for deleted ghosts and see if they moved */ 1038 PetscCall(PetscCDGetHeadPos(aggs_2, lid, &pos)); 1039 while (pos) { 1040 PetscInt gid; 1041 PetscCall(PetscCDIntNdGetID(pos, &gid)); 1042 1043 if (gid < my0 || gid >= Iend) { 1044 PetscCall(PCGAMGHashTableFind(&gid_cpid, gid, &cpid)); 1045 if (cpid != -1) { 1046 /* a moved ghost - */ 1047 /* id_llist_2[lastid] = id_llist_2[flid]; /\* remove 'flid' from list *\/ */ 1048 PetscCall(PetscCDRemoveNextNode(aggs_2, lid, last)); 1049 } else last = pos; 1050 } else last = pos; 1051 1052 PetscCall(PetscCDGetNextPos(aggs_2, lid, &pos)); 1053 } /* loop over list of deleted */ 1054 } /* selected */ 1055 } 1056 PetscCall(PCGAMGHashTableDestroy(&gid_cpid)); 1057 1058 /* look at ghosts, see if they changed - and it */ 1059 for (cpid = 0; cpid < nghost_2; cpid++) { 1060 PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]); 1061 1062 if (sgid_new >= my0 && sgid_new < Iend) { /* this is mine */ 1063 PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]); 1064 PetscInt slid_new = sgid_new - my0, hav = 0; 1065 PetscCDIntNd *pos; 1066 1067 /* search for this gid to see if I have it */ 1068 PetscCall(PetscCDGetHeadPos(aggs_2, slid_new, &pos)); 1069 while (pos) { 1070 PetscInt gidj; 1071 1072 PetscCall(PetscCDIntNdGetID(pos, &gidj)); 1073 PetscCall(PetscCDGetNextPos(aggs_2, slid_new, &pos)); 1074 1075 if (gidj == gid) { 1076 hav = 1; 1077 break; 1078 } 1079 } 1080 if (hav != 1) { 1081 /* insert 'flidj' into head of llist */ 1082 PetscCall(PetscCDAppendID(aggs_2, slid_new, gid)); 1083 } 1084 } 1085 } 1086 PetscCall(VecRestoreArray(mpimat_1->lvec, &cpcol_1_state)); 1087 PetscCall(VecRestoreArray(mpimat_2->lvec, &cpcol_2_state)); 1088 PetscCall(VecRestoreArray(ghostparents2, &cpcol_2_parent)); 1089 PetscCall(VecRestoreArray(ghostgids2, &cpcol_2_gid)); 1090 PetscCall(VecDestroy(&ghostgids2)); 1091 PetscCall(VecDestroy(&ghostparents2)); 1092 PetscCall(VecDestroy(&ghost_par_orig2)); 1093 } 1094 PetscCall(PetscFree3(lid_state, lid_parent_gid, lid_cprowID_1)); 1095 PetscFunctionReturn(PETSC_SUCCESS); 1096 } 1097 1098 /* 1099 PCGAMGCoarsen_AGG - supports squaring the graph (deprecated) and new graph for 1100 communication of QR data used with HEM and MISk coarsening 1101 1102 Input Parameter: 1103 . a_pc - this 1104 1105 Input/Output Parameter: 1106 . a_Gmat1 - graph to coarsen (in), graph off processor edges for QR gather scatter (out) 1107 1108 Output Parameter: 1109 . agg_lists - list of aggregates 1110 1111 */ 1112 static PetscErrorCode PCGAMGCoarsen_AGG(PC a_pc, Mat *a_Gmat1, PetscCoarsenData **agg_lists) 1113 { 1114 PC_MG *mg = (PC_MG *)a_pc->data; 1115 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1116 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 1117 Mat Gmat2, Gmat1 = *a_Gmat1; /* aggressive graph */ 1118 IS perm; 1119 PetscInt Istart, Iend, Ii, nloc, bs, nn; 1120 PetscInt *permute, *degree; 1121 PetscBool *bIndexSet; 1122 PetscReal hashfact; 1123 PetscInt iSwapIndex; 1124 PetscRandom random; 1125 MPI_Comm comm; 1126 1127 PetscFunctionBegin; 1128 PetscCall(PetscObjectGetComm((PetscObject)Gmat1, &comm)); 1129 PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0)); 1130 PetscCall(MatGetLocalSize(Gmat1, &nn, NULL)); 1131 PetscCall(MatGetBlockSize(Gmat1, &bs)); 1132 PetscCheck(bs == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "bs %" PetscInt_FMT " must be 1", bs); 1133 nloc = nn / bs; 1134 /* get MIS aggs - randomize */ 1135 PetscCall(PetscMalloc2(nloc, &permute, nloc, °ree)); 1136 PetscCall(PetscCalloc1(nloc, &bIndexSet)); 1137 for (Ii = 0; Ii < nloc; Ii++) permute[Ii] = Ii; 1138 PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &random)); 1139 PetscCall(MatGetOwnershipRange(Gmat1, &Istart, &Iend)); 1140 for (Ii = 0; Ii < nloc; Ii++) { 1141 PetscInt nc; 1142 1143 PetscCall(MatGetRow(Gmat1, Istart + Ii, &nc, NULL, NULL)); 1144 degree[Ii] = nc; 1145 PetscCall(MatRestoreRow(Gmat1, Istart + Ii, &nc, NULL, NULL)); 1146 } 1147 for (Ii = 0; Ii < nloc; Ii++) { 1148 PetscCall(PetscRandomGetValueReal(random, &hashfact)); 1149 iSwapIndex = (PetscInt)(hashfact * nloc) % nloc; 1150 if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) { 1151 PetscInt iTemp = permute[iSwapIndex]; 1152 1153 permute[iSwapIndex] = permute[Ii]; 1154 permute[Ii] = iTemp; 1155 iTemp = degree[iSwapIndex]; 1156 degree[iSwapIndex] = degree[Ii]; 1157 degree[Ii] = iTemp; 1158 bIndexSet[iSwapIndex] = PETSC_TRUE; 1159 } 1160 } 1161 // apply minimum degree ordering -- NEW 1162 if (pc_gamg_agg->use_minimum_degree_ordering) { PetscCall(PetscSortIntWithArray(nloc, degree, permute)); } 1163 PetscCall(PetscFree(bIndexSet)); 1164 PetscCall(PetscRandomDestroy(&random)); 1165 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm)); 1166 PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_MIS], 0, 0, 0, 0)); 1167 // square graph 1168 if (pc_gamg->current_level < pc_gamg_agg->aggressive_coarsening_levels && pc_gamg_agg->use_aggressive_square_graph) { 1169 PetscCall(PCGAMGSquareGraph_GAMG(a_pc, Gmat1, &Gmat2)); 1170 } else Gmat2 = Gmat1; 1171 // switch to old MIS-1 for square graph 1172 if (pc_gamg->current_level < pc_gamg_agg->aggressive_coarsening_levels) { 1173 if (!pc_gamg_agg->use_aggressive_square_graph) PetscCall(MatCoarsenMISKSetDistance(pc_gamg_agg->crs, pc_gamg_agg->aggressive_mis_k)); // hardwire to MIS-2 1174 else PetscCall(MatCoarsenSetType(pc_gamg_agg->crs, MATCOARSENMIS)); // old MIS -- side effect 1175 } else if (pc_gamg_agg->use_aggressive_square_graph && pc_gamg_agg->aggressive_coarsening_levels > 0) { // we reset the MIS 1176 const char *prefix; 1177 1178 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)a_pc, &prefix)); 1179 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc_gamg_agg->crs, prefix)); 1180 PetscCall(MatCoarsenSetFromOptions(pc_gamg_agg->crs)); // get the default back on non-aggressive levels when square graph switched to old MIS 1181 } 1182 PetscCall(MatCoarsenSetAdjacency(pc_gamg_agg->crs, Gmat2)); 1183 PetscCall(MatCoarsenSetStrictAggs(pc_gamg_agg->crs, PETSC_TRUE)); 1184 PetscCall(MatCoarsenSetGreedyOrdering(pc_gamg_agg->crs, perm)); 1185 PetscCall(MatCoarsenApply(pc_gamg_agg->crs)); 1186 PetscCall(MatCoarsenGetData(pc_gamg_agg->crs, agg_lists)); /* output */ 1187 1188 PetscCall(ISDestroy(&perm)); 1189 PetscCall(PetscFree2(permute, degree)); 1190 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_MIS], 0, 0, 0, 0)); 1191 1192 if (Gmat2 != Gmat1) { // square graph, we need ghosts for selected 1193 PetscCoarsenData *llist = *agg_lists; 1194 1195 PetscCall(fixAggregatesWithSquare(a_pc, Gmat2, Gmat1, *agg_lists)); 1196 PetscCall(MatDestroy(&Gmat1)); 1197 *a_Gmat1 = Gmat2; /* output */ 1198 PetscCall(PetscCDSetMat(llist, *a_Gmat1)); /* Need a graph with ghosts here */ 1199 } 1200 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0)); 1201 PetscFunctionReturn(PETSC_SUCCESS); 1202 } 1203 1204 /* 1205 PCGAMGConstructProlongator_AGG 1206 1207 Input Parameter: 1208 . pc - this 1209 . Amat - matrix on this fine level 1210 . Graph - used to get ghost data for nodes in 1211 . agg_lists - list of aggregates 1212 Output Parameter: 1213 . a_P_out - prolongation operator to the next level 1214 */ 1215 static PetscErrorCode PCGAMGConstructProlongator_AGG(PC pc, Mat Amat, PetscCoarsenData *agg_lists, Mat *a_P_out) 1216 { 1217 PC_MG *mg = (PC_MG *)pc->data; 1218 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1219 const PetscInt col_bs = pc_gamg->data_cell_cols; 1220 PetscInt Istart, Iend, nloc, ii, jj, kk, my0, nLocalSelected, bs; 1221 Mat Gmat, Prol; 1222 PetscMPIInt size; 1223 MPI_Comm comm; 1224 PetscReal *data_w_ghost; 1225 PetscInt myCrs0, nbnodes = 0, *flid_fgid; 1226 MatType mtype; 1227 1228 PetscFunctionBegin; 1229 PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm)); 1230 PetscCheck(col_bs >= 1, comm, PETSC_ERR_PLIB, "Column bs cannot be less than 1"); 1231 PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0)); 1232 PetscCallMPI(MPI_Comm_size(comm, &size)); 1233 PetscCall(MatGetOwnershipRange(Amat, &Istart, &Iend)); 1234 PetscCall(MatGetBlockSize(Amat, &bs)); 1235 nloc = (Iend - Istart) / bs; 1236 my0 = Istart / bs; 1237 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); 1238 PetscCall(PetscCDGetMat(agg_lists, &Gmat)); // get auxiliary matrix for ghost edges for size > 1 1239 1240 /* get 'nLocalSelected' */ 1241 for (ii = 0, nLocalSelected = 0; ii < nloc; ii++) { 1242 PetscBool ise; 1243 1244 /* filter out singletons 0 or 1? */ 1245 PetscCall(PetscCDIsEmptyAt(agg_lists, ii, &ise)); 1246 if (!ise) nLocalSelected++; 1247 } 1248 1249 /* create prolongator, create P matrix */ 1250 PetscCall(MatGetType(Amat, &mtype)); 1251 PetscCall(MatCreate(comm, &Prol)); 1252 PetscCall(MatSetSizes(Prol, nloc * bs, nLocalSelected * col_bs, PETSC_DETERMINE, PETSC_DETERMINE)); 1253 PetscCall(MatSetBlockSizes(Prol, bs, col_bs)); // should this be before MatSetSizes? 1254 PetscCall(MatSetType(Prol, mtype)); 1255 #if PetscDefined(HAVE_DEVICE) 1256 PetscBool flg; 1257 PetscCall(MatBoundToCPU(Amat, &flg)); 1258 PetscCall(MatBindToCPU(Prol, flg)); 1259 if (flg) PetscCall(MatSetBindingPropagates(Prol, PETSC_TRUE)); 1260 #endif 1261 PetscCall(MatSeqAIJSetPreallocation(Prol, col_bs, NULL)); 1262 PetscCall(MatMPIAIJSetPreallocation(Prol, col_bs, NULL, col_bs, NULL)); 1263 1264 /* can get all points "removed" */ 1265 PetscCall(MatGetSize(Prol, &kk, &ii)); 1266 if (!ii) { 1267 PetscCall(PetscInfo(pc, "%s: No selected points on coarse grid\n", ((PetscObject)pc)->prefix)); 1268 PetscCall(MatDestroy(&Prol)); 1269 *a_P_out = NULL; /* out */ 1270 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0)); 1271 PetscFunctionReturn(PETSC_SUCCESS); 1272 } 1273 PetscCall(PetscInfo(pc, "%s: New grid %" PetscInt_FMT " nodes\n", ((PetscObject)pc)->prefix, ii / col_bs)); 1274 PetscCall(MatGetOwnershipRangeColumn(Prol, &myCrs0, &kk)); 1275 1276 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); 1277 myCrs0 = myCrs0 / col_bs; 1278 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); 1279 1280 /* create global vector of data in 'data_w_ghost' */ 1281 PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROLA], 0, 0, 0, 0)); 1282 if (size > 1) { /* get ghost null space data */ 1283 PetscReal *tmp_gdata, *tmp_ldata, *tp2; 1284 1285 PetscCall(PetscMalloc1(nloc, &tmp_ldata)); 1286 for (jj = 0; jj < col_bs; jj++) { 1287 for (kk = 0; kk < bs; kk++) { 1288 PetscInt ii, stride; 1289 const PetscReal *tp = PetscSafePointerPlusOffset(pc_gamg->data, jj * bs * nloc + kk); 1290 1291 for (ii = 0; ii < nloc; ii++, tp += bs) tmp_ldata[ii] = *tp; 1292 1293 PetscCall(PCGAMGGetDataWithGhosts(Gmat, 1, tmp_ldata, &stride, &tmp_gdata)); 1294 1295 if (!jj && !kk) { /* now I know how many total nodes - allocate TODO: move below and do in one 'col_bs' call */ 1296 PetscCall(PetscMalloc1(stride * bs * col_bs, &data_w_ghost)); 1297 nbnodes = bs * stride; 1298 } 1299 tp2 = PetscSafePointerPlusOffset(data_w_ghost, jj * bs * stride + kk); 1300 for (ii = 0; ii < stride; ii++, tp2 += bs) *tp2 = tmp_gdata[ii]; 1301 PetscCall(PetscFree(tmp_gdata)); 1302 } 1303 } 1304 PetscCall(PetscFree(tmp_ldata)); 1305 } else { 1306 nbnodes = bs * nloc; 1307 data_w_ghost = (PetscReal *)pc_gamg->data; 1308 } 1309 1310 /* get 'flid_fgid' TODO - move up to get 'stride' and do get null space data above in one step (jj loop) */ 1311 if (size > 1) { 1312 PetscReal *fid_glid_loc, *fiddata; 1313 PetscInt stride; 1314 1315 PetscCall(PetscMalloc1(nloc, &fid_glid_loc)); 1316 for (kk = 0; kk < nloc; kk++) fid_glid_loc[kk] = (PetscReal)(my0 + kk); 1317 PetscCall(PCGAMGGetDataWithGhosts(Gmat, 1, fid_glid_loc, &stride, &fiddata)); 1318 PetscCall(PetscMalloc1(stride, &flid_fgid)); /* copy real data to in */ 1319 for (kk = 0; kk < stride; kk++) flid_fgid[kk] = (PetscInt)fiddata[kk]; 1320 PetscCall(PetscFree(fiddata)); 1321 1322 PetscCheck(stride == nbnodes / bs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "stride %" PetscInt_FMT " != nbnodes %" PetscInt_FMT "/bs %" PetscInt_FMT, stride, nbnodes, bs); 1323 PetscCall(PetscFree(fid_glid_loc)); 1324 } else { 1325 PetscCall(PetscMalloc1(nloc, &flid_fgid)); 1326 for (kk = 0; kk < nloc; kk++) flid_fgid[kk] = my0 + kk; 1327 } 1328 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROLA], 0, 0, 0, 0)); 1329 /* get P0 */ 1330 PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROLB], 0, 0, 0, 0)); 1331 { 1332 PetscReal *data_out = NULL; 1333 1334 PetscCall(formProl0(agg_lists, bs, col_bs, myCrs0, nbnodes, data_w_ghost, flid_fgid, &data_out, Prol)); 1335 PetscCall(PetscFree(pc_gamg->data)); 1336 1337 pc_gamg->data = data_out; 1338 pc_gamg->data_cell_rows = col_bs; 1339 pc_gamg->data_sz = col_bs * col_bs * nLocalSelected; 1340 } 1341 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROLB], 0, 0, 0, 0)); 1342 if (size > 1) PetscCall(PetscFree(data_w_ghost)); 1343 PetscCall(PetscFree(flid_fgid)); 1344 1345 *a_P_out = Prol; /* out */ 1346 PetscCall(MatViewFromOptions(Prol, NULL, "-pc_gamg_agg_view_initial_prolongation")); 1347 1348 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0)); 1349 PetscFunctionReturn(PETSC_SUCCESS); 1350 } 1351 1352 /* 1353 PCGAMGOptimizeProlongator_AGG - given the initial prolongator optimizes it by smoothed aggregation pc_gamg_agg->nsmooths times 1354 1355 Input Parameter: 1356 . pc - this 1357 . Amat - matrix on this fine level 1358 In/Output Parameter: 1359 . a_P - prolongation operator to the next level 1360 */ 1361 static PetscErrorCode PCGAMGOptimizeProlongator_AGG(PC pc, Mat Amat, Mat *a_P) 1362 { 1363 PC_MG *mg = (PC_MG *)pc->data; 1364 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1365 PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx; 1366 PetscInt jj; 1367 Mat Prol = *a_P; 1368 MPI_Comm comm; 1369 KSP eksp; 1370 Vec bb, xx; 1371 PC epc; 1372 PetscReal alpha, emax, emin; 1373 1374 PetscFunctionBegin; 1375 PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm)); 1376 PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_OPT], 0, 0, 0, 0)); 1377 1378 /* compute maximum singular value of operator to be used in smoother */ 1379 if (0 < pc_gamg_agg->nsmooths) { 1380 /* get eigen estimates */ 1381 if (pc_gamg->emax > 0) { 1382 emin = pc_gamg->emin; 1383 emax = pc_gamg->emax; 1384 } else { 1385 const char *prefix; 1386 1387 PetscCall(MatCreateVecs(Amat, &bb, NULL)); 1388 PetscCall(MatCreateVecs(Amat, &xx, NULL)); 1389 PetscCall(KSPSetNoisy_Private(bb)); 1390 1391 PetscCall(KSPCreate(comm, &eksp)); 1392 PetscCall(KSPSetNestLevel(eksp, pc->kspnestlevel)); 1393 PetscCall(PCGetOptionsPrefix(pc, &prefix)); 1394 PetscCall(KSPSetOptionsPrefix(eksp, prefix)); 1395 PetscCall(KSPAppendOptionsPrefix(eksp, "pc_gamg_esteig_")); 1396 { 1397 PetscBool isset, sflg; 1398 1399 PetscCall(MatIsSPDKnown(Amat, &isset, &sflg)); 1400 if (isset && sflg) PetscCall(KSPSetType(eksp, KSPCG)); 1401 } 1402 PetscCall(KSPSetErrorIfNotConverged(eksp, pc->erroriffailure)); 1403 PetscCall(KSPSetNormType(eksp, KSP_NORM_NONE)); 1404 1405 PetscCall(KSPSetInitialGuessNonzero(eksp, PETSC_FALSE)); 1406 PetscCall(KSPSetOperators(eksp, Amat, Amat)); 1407 1408 PetscCall(KSPGetPC(eksp, &epc)); 1409 PetscCall(PCSetType(epc, PCJACOBI)); /* smoother in smoothed agg. */ 1410 1411 PetscCall(KSPSetTolerances(eksp, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT, 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 1412 1413 PetscCall(KSPSetFromOptions(eksp)); 1414 PetscCall(KSPSetComputeSingularValues(eksp, PETSC_TRUE)); 1415 PetscCall(KSPSolve(eksp, bb, xx)); 1416 PetscCall(KSPCheckSolve(eksp, pc, xx)); 1417 1418 PetscCall(KSPComputeExtremeSingularValues(eksp, &emax, &emin)); 1419 PetscCall(PetscInfo(pc, "%s: Smooth P0: max eigen=%e min=%e PC=%s\n", ((PetscObject)pc)->prefix, (double)emax, (double)emin, PCJACOBI)); 1420 PetscCall(VecDestroy(&xx)); 1421 PetscCall(VecDestroy(&bb)); 1422 PetscCall(KSPDestroy(&eksp)); 1423 } 1424 if (pc_gamg->use_sa_esteig) { 1425 mg->min_eigen_DinvA[pc_gamg->current_level] = emin; 1426 mg->max_eigen_DinvA[pc_gamg->current_level] = emax; 1427 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)); 1428 } else { 1429 mg->min_eigen_DinvA[pc_gamg->current_level] = 0; 1430 mg->max_eigen_DinvA[pc_gamg->current_level] = 0; 1431 } 1432 } else { 1433 mg->min_eigen_DinvA[pc_gamg->current_level] = 0; 1434 mg->max_eigen_DinvA[pc_gamg->current_level] = 0; 1435 } 1436 1437 /* smooth P0 */ 1438 if (pc_gamg_agg->nsmooths > 0) { 1439 Vec diag; 1440 1441 /* TODO: Set a PCFailedReason and exit the building of the AMG preconditioner */ 1442 PetscCheck(emax != 0.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_PLIB, "Computed maximum singular value as zero"); 1443 1444 PetscCall(MatCreateVecs(Amat, &diag, NULL)); 1445 PetscCall(MatGetDiagonal(Amat, diag)); /* effectively PCJACOBI */ 1446 PetscCall(VecReciprocal(diag)); 1447 1448 for (jj = 0; jj < pc_gamg_agg->nsmooths; jj++) { 1449 Mat tMat; 1450 1451 PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_OPTSM], 0, 0, 0, 0)); 1452 /* 1453 Smooth aggregation on the prolongator 1454 1455 P_{i} := (I - 1.4/emax D^{-1}A) P_i\{i-1} 1456 */ 1457 PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2], 0, 0, 0, 0)); 1458 PetscCall(MatMatMult(Amat, Prol, MAT_INITIAL_MATRIX, PETSC_CURRENT, &tMat)); 1459 PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2], 0, 0, 0, 0)); 1460 PetscCall(MatProductClear(tMat)); 1461 PetscCall(MatDiagonalScale(tMat, diag, NULL)); 1462 1463 /* TODO: Document the 1.4 and don't hardwire it in this routine */ 1464 alpha = -1.4 / emax; 1465 PetscCall(MatAYPX(tMat, alpha, Prol, SUBSET_NONZERO_PATTERN)); 1466 PetscCall(MatDestroy(&Prol)); 1467 Prol = tMat; 1468 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_OPTSM], 0, 0, 0, 0)); 1469 } 1470 PetscCall(VecDestroy(&diag)); 1471 } 1472 PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_OPT], 0, 0, 0, 0)); 1473 PetscCall(MatViewFromOptions(Prol, NULL, "-pc_gamg_agg_view_prolongation")); 1474 *a_P = Prol; 1475 PetscFunctionReturn(PETSC_SUCCESS); 1476 } 1477 1478 /*MC 1479 PCGAMGAGG - Smooth aggregation, {cite}`vanek1996algebraic`, {cite}`vanek2001convergence`, variant of PETSc's algebraic multigrid (`PCGAMG`) preconditioner 1480 1481 Options Database Keys: 1482 + -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation to construct prolongation 1483 . -pc_gamg_aggressive_coarsening <n,default=1> - number of aggressive coarsening (MIS-2) levels from finest. 1484 . -pc_gamg_aggressive_square_graph <bool,default=false> - Use square graph (A'A) or MIS-k (k=2) for aggressive coarsening 1485 . -pc_gamg_mis_k_minimum_degree_ordering <bool,default=true> - Use minimum degree ordering in greedy MIS algorithm 1486 . -pc_gamg_pc_gamg_asm_hem_aggs <n,default=0> - Number of HEM aggregation steps for ASM smoother 1487 - -pc_gamg_aggressive_mis_k <n,default=2> - Number (k) distance in MIS coarsening (>2 is 'aggressive') 1488 1489 Level: intermediate 1490 1491 Notes: 1492 To obtain good performance for `PCGAMG` for vector valued problems you must 1493 call `MatSetBlockSize()` to indicate the number of degrees of freedom per grid point. 1494 Call `MatSetNearNullSpace()` (or `PCSetCoordinates()` if solving the equations of elasticity) to indicate the near null space of the operator 1495 1496 The many options for `PCMG` and `PCGAMG` such as controlling the smoothers on each level etc. also work for `PCGAMGAGG` 1497 1498 .seealso: `PCGAMG`, [the Users Manual section on PCGAMG](sec_amg), [the Users Manual section on PCMG](sec_mg), [](ch_ksp), `PCCreate()`, `PCSetType()`, 1499 `MatSetBlockSize()`, `PCMGType`, `PCSetCoordinates()`, `MatSetNearNullSpace()`, `PCGAMGSetType()`, 1500 `PCGAMGAGG`, `PCGAMGGEO`, `PCGAMGCLASSICAL`, `PCGAMGSetProcEqLim()`, `PCGAMGSetCoarseEqLim()`, `PCGAMGSetRepartition()`, `PCGAMGRegister()`, 1501 `PCGAMGSetReuseInterpolation()`, `PCGAMGASMSetUseAggs()`, `PCGAMGSetParallelCoarseGridSolve()`, `PCGAMGSetNlevels()`, `PCGAMGSetThreshold()`, 1502 `PCGAMGGetType()`, `PCGAMGSetUseSAEstEig()` 1503 M*/ 1504 PetscErrorCode PCCreateGAMG_AGG(PC pc) 1505 { 1506 PC_MG *mg = (PC_MG *)pc->data; 1507 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 1508 PC_GAMG_AGG *pc_gamg_agg; 1509 1510 PetscFunctionBegin; 1511 /* create sub context for SA */ 1512 PetscCall(PetscNew(&pc_gamg_agg)); 1513 pc_gamg->subctx = pc_gamg_agg; 1514 1515 pc_gamg->ops->setfromoptions = PCSetFromOptions_GAMG_AGG; 1516 pc_gamg->ops->destroy = PCDestroy_GAMG_AGG; 1517 /* reset does not do anything; setup not virtual */ 1518 1519 /* set internal function pointers */ 1520 pc_gamg->ops->creategraph = PCGAMGCreateGraph_AGG; 1521 pc_gamg->ops->coarsen = PCGAMGCoarsen_AGG; 1522 pc_gamg->ops->prolongator = PCGAMGConstructProlongator_AGG; 1523 pc_gamg->ops->optprolongator = PCGAMGOptimizeProlongator_AGG; 1524 pc_gamg->ops->createdefaultdata = PCSetData_AGG; 1525 pc_gamg->ops->view = PCView_GAMG_AGG; 1526 1527 pc_gamg_agg->nsmooths = 1; 1528 pc_gamg_agg->aggressive_coarsening_levels = 1; 1529 pc_gamg_agg->use_aggressive_square_graph = PETSC_TRUE; 1530 pc_gamg_agg->use_minimum_degree_ordering = PETSC_FALSE; 1531 pc_gamg_agg->use_low_mem_filter = PETSC_FALSE; 1532 pc_gamg_agg->aggressive_mis_k = 2; 1533 1534 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNSmooths_C", PCGAMGSetNSmooths_AGG)); 1535 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveLevels_C", PCGAMGSetAggressiveLevels_AGG)); 1536 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveSquareGraph_C", PCGAMGSetAggressiveSquareGraph_AGG)); 1537 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetMinDegreeOrdering_C", PCGAMGMISkSetMinDegreeOrdering_AGG)); 1538 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetLowMemoryFilter_C", PCGAMGSetLowMemoryFilter_AGG)); 1539 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGMISkSetAggressive_C", PCGAMGMISkSetAggressive_AGG)); 1540 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", PCSetCoordinates_AGG)); 1541 PetscFunctionReturn(PETSC_SUCCESS); 1542 } 1543