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