1 #include <../src/ksp/pc/impls/gamg/gamg.h> /*I "petscpc.h" I*/ 2 #include <petscsf.h> 3 4 PetscFunctionList PCGAMGClassicalProlongatorList = NULL; 5 PetscBool PCGAMGClassicalPackageInitialized = PETSC_FALSE; 6 7 typedef struct { 8 PetscReal interp_threshold; /* interpolation threshold */ 9 char prolongtype[256]; 10 PetscInt nsmooths; /* number of jacobi smoothings on the prolongator */ 11 } PC_GAMG_Classical; 12 13 /*@C 14 PCGAMGClassicalSetType - Sets the type of classical interpolation to use 15 16 Collective on PC 17 18 Input Parameters: 19 . pc - the preconditioner context 20 21 Options Database Key: 22 . -pc_gamg_classical_type <direct,standard> - set type of Classical AMG prolongation 23 24 Level: intermediate 25 26 .seealso: `()` 27 @*/ 28 PetscErrorCode PCGAMGClassicalSetType(PC pc, PCGAMGClassicalType type) { 29 PetscFunctionBegin; 30 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 31 PetscTryMethod(pc, "PCGAMGClassicalSetType_C", (PC, PCGAMGClassicalType), (pc, type)); 32 PetscFunctionReturn(0); 33 } 34 35 /*@C 36 PCGAMGClassicalGetType - Gets the type of classical interpolation to use 37 38 Collective on PC 39 40 Input Parameter: 41 . pc - the preconditioner context 42 43 Output Parameter: 44 . type - the type used 45 46 Level: intermediate 47 48 .seealso: `()` 49 @*/ 50 PetscErrorCode PCGAMGClassicalGetType(PC pc, PCGAMGClassicalType *type) { 51 PetscFunctionBegin; 52 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 53 PetscUseMethod(pc, "PCGAMGClassicalGetType_C", (PC, PCGAMGClassicalType *), (pc, type)); 54 PetscFunctionReturn(0); 55 } 56 57 static PetscErrorCode PCGAMGClassicalSetType_GAMG(PC pc, PCGAMGClassicalType type) { 58 PC_MG *mg = (PC_MG *)pc->data; 59 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 60 PC_GAMG_Classical *cls = (PC_GAMG_Classical *)pc_gamg->subctx; 61 62 PetscFunctionBegin; 63 PetscCall(PetscStrcpy(cls->prolongtype, type)); 64 PetscFunctionReturn(0); 65 } 66 67 static PetscErrorCode PCGAMGClassicalGetType_GAMG(PC pc, PCGAMGClassicalType *type) { 68 PC_MG *mg = (PC_MG *)pc->data; 69 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 70 PC_GAMG_Classical *cls = (PC_GAMG_Classical *)pc_gamg->subctx; 71 72 PetscFunctionBegin; 73 *type = cls->prolongtype; 74 PetscFunctionReturn(0); 75 } 76 77 PetscErrorCode PCGAMGGraph_Classical(PC pc, Mat A, Mat *G) { 78 PetscInt s, f, n, idx, lidx, gidx; 79 PetscInt r, c, ncols; 80 const PetscInt *rcol; 81 const PetscScalar *rval; 82 PetscInt *gcol; 83 PetscScalar *gval; 84 PetscReal rmax; 85 PetscInt cmax = 0; 86 PC_MG *mg = (PC_MG *)pc->data; 87 PC_GAMG *gamg = (PC_GAMG *)mg->innerctx; 88 PetscInt *gsparse, *lsparse; 89 PetscScalar *Amax; 90 MatType mtype; 91 92 PetscFunctionBegin; 93 PetscCall(MatGetOwnershipRange(A, &s, &f)); 94 n = f - s; 95 PetscCall(PetscMalloc3(n, &lsparse, n, &gsparse, n, &Amax)); 96 97 for (r = 0; r < n; r++) { 98 lsparse[r] = 0; 99 gsparse[r] = 0; 100 } 101 102 for (r = s; r < f; r++) { 103 /* determine the maximum off-diagonal in each row */ 104 rmax = 0.; 105 PetscCall(MatGetRow(A, r, &ncols, &rcol, &rval)); 106 for (c = 0; c < ncols; c++) { 107 if (PetscRealPart(-rval[c]) > rmax && rcol[c] != r) rmax = PetscRealPart(-rval[c]); 108 } 109 Amax[r - s] = rmax; 110 if (ncols > cmax) cmax = ncols; 111 lidx = 0; 112 gidx = 0; 113 /* create the local and global sparsity patterns */ 114 for (c = 0; c < ncols; c++) { 115 if (PetscRealPart(-rval[c]) > gamg->threshold[0] * PetscRealPart(Amax[r - s]) || rcol[c] == r) { 116 if (rcol[c] < f && rcol[c] >= s) { 117 lidx++; 118 } else { 119 gidx++; 120 } 121 } 122 } 123 PetscCall(MatRestoreRow(A, r, &ncols, &rcol, &rval)); 124 lsparse[r - s] = lidx; 125 gsparse[r - s] = gidx; 126 } 127 PetscCall(PetscMalloc2(cmax, &gval, cmax, &gcol)); 128 129 PetscCall(MatCreate(PetscObjectComm((PetscObject)A), G)); 130 PetscCall(MatGetType(A, &mtype)); 131 PetscCall(MatSetType(*G, mtype)); 132 PetscCall(MatSetSizes(*G, n, n, PETSC_DETERMINE, PETSC_DETERMINE)); 133 PetscCall(MatMPIAIJSetPreallocation(*G, 0, lsparse, 0, gsparse)); 134 PetscCall(MatSeqAIJSetPreallocation(*G, 0, lsparse)); 135 for (r = s; r < f; r++) { 136 PetscCall(MatGetRow(A, r, &ncols, &rcol, &rval)); 137 idx = 0; 138 for (c = 0; c < ncols; c++) { 139 /* classical strength of connection */ 140 if (PetscRealPart(-rval[c]) > gamg->threshold[0] * PetscRealPart(Amax[r - s]) || rcol[c] == r) { 141 gcol[idx] = rcol[c]; 142 gval[idx] = rval[c]; 143 idx++; 144 } 145 } 146 PetscCall(MatSetValues(*G, 1, &r, idx, gcol, gval, INSERT_VALUES)); 147 PetscCall(MatRestoreRow(A, r, &ncols, &rcol, &rval)); 148 } 149 PetscCall(MatAssemblyBegin(*G, MAT_FINAL_ASSEMBLY)); 150 PetscCall(MatAssemblyEnd(*G, MAT_FINAL_ASSEMBLY)); 151 152 PetscCall(PetscFree2(gval, gcol)); 153 PetscCall(PetscFree3(lsparse, gsparse, Amax)); 154 PetscFunctionReturn(0); 155 } 156 157 PetscErrorCode PCGAMGCoarsen_Classical(PC pc, Mat *G, PetscCoarsenData **agg_lists) { 158 MatCoarsen crs; 159 MPI_Comm fcomm = ((PetscObject)pc)->comm; 160 161 PetscFunctionBegin; 162 PetscCheck(G, fcomm, PETSC_ERR_ARG_WRONGSTATE, "Must set Graph in PC in PCGAMG before coarsening"); 163 164 PetscCall(MatCoarsenCreate(fcomm, &crs)); 165 PetscCall(MatCoarsenSetFromOptions(crs)); 166 PetscCall(MatCoarsenSetAdjacency(crs, *G)); 167 PetscCall(MatCoarsenSetStrictAggs(crs, PETSC_TRUE)); 168 PetscCall(MatCoarsenApply(crs)); 169 PetscCall(MatCoarsenGetData(crs, agg_lists)); 170 PetscCall(MatCoarsenDestroy(&crs)); 171 PetscFunctionReturn(0); 172 } 173 174 PetscErrorCode PCGAMGProlongator_Classical_Direct(PC pc, Mat A, Mat G, PetscCoarsenData *agg_lists, Mat *P) { 175 PC_MG *mg = (PC_MG *)pc->data; 176 PC_GAMG *gamg = (PC_GAMG *)mg->innerctx; 177 PetscBool iscoarse, isMPIAIJ, isSEQAIJ; 178 PetscInt fn, cn, fs, fe, cs, ce, i, j, ncols, col, row_f, row_c, cmax = 0, idx, noff; 179 PetscInt *lcid, *gcid, *lsparse, *gsparse, *colmap, *pcols; 180 const PetscInt *rcol; 181 PetscReal *Amax_pos, *Amax_neg; 182 PetscScalar g_pos, g_neg, a_pos, a_neg, diag, invdiag, alpha, beta, pij; 183 PetscScalar *pvals; 184 const PetscScalar *rval; 185 Mat lA, gA = NULL; 186 MatType mtype; 187 Vec C, lvec; 188 PetscLayout clayout; 189 PetscSF sf; 190 Mat_MPIAIJ *mpiaij; 191 192 PetscFunctionBegin; 193 PetscCall(MatGetOwnershipRange(A, &fs, &fe)); 194 fn = fe - fs; 195 PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPIAIJ, &isMPIAIJ)); 196 PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQAIJ, &isSEQAIJ)); 197 PetscCheck(isMPIAIJ || isSEQAIJ, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Classical AMG requires MPIAIJ matrix"); 198 if (isMPIAIJ) { 199 mpiaij = (Mat_MPIAIJ *)A->data; 200 lA = mpiaij->A; 201 gA = mpiaij->B; 202 lvec = mpiaij->lvec; 203 PetscCall(VecGetSize(lvec, &noff)); 204 colmap = mpiaij->garray; 205 PetscCall(MatGetLayouts(A, NULL, &clayout)); 206 PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)A), &sf)); 207 PetscCall(PetscSFSetGraphLayout(sf, clayout, noff, NULL, PETSC_COPY_VALUES, colmap)); 208 PetscCall(PetscMalloc1(noff, &gcid)); 209 } else { 210 lA = A; 211 } 212 PetscCall(PetscMalloc5(fn, &lsparse, fn, &gsparse, fn, &lcid, fn, &Amax_pos, fn, &Amax_neg)); 213 214 /* count the number of coarse unknowns */ 215 cn = 0; 216 for (i = 0; i < fn; i++) { 217 /* filter out singletons */ 218 PetscCall(PetscCDEmptyAt(agg_lists, i, &iscoarse)); 219 lcid[i] = -1; 220 if (!iscoarse) cn++; 221 } 222 223 /* create the coarse vector */ 224 PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)A), cn, PETSC_DECIDE, &C)); 225 PetscCall(VecGetOwnershipRange(C, &cs, &ce)); 226 227 cn = 0; 228 for (i = 0; i < fn; i++) { 229 PetscCall(PetscCDEmptyAt(agg_lists, i, &iscoarse)); 230 if (!iscoarse) { 231 lcid[i] = cs + cn; 232 cn++; 233 } else { 234 lcid[i] = -1; 235 } 236 } 237 238 if (gA) { 239 PetscCall(PetscSFBcastBegin(sf, MPIU_INT, lcid, gcid, MPI_REPLACE)); 240 PetscCall(PetscSFBcastEnd(sf, MPIU_INT, lcid, gcid, MPI_REPLACE)); 241 } 242 243 /* determine the largest off-diagonal entries in each row */ 244 for (i = fs; i < fe; i++) { 245 Amax_pos[i - fs] = 0.; 246 Amax_neg[i - fs] = 0.; 247 PetscCall(MatGetRow(A, i, &ncols, &rcol, &rval)); 248 for (j = 0; j < ncols; j++) { 249 if ((PetscRealPart(-rval[j]) > Amax_neg[i - fs]) && i != rcol[j]) Amax_neg[i - fs] = PetscAbsScalar(rval[j]); 250 if ((PetscRealPart(rval[j]) > Amax_pos[i - fs]) && i != rcol[j]) Amax_pos[i - fs] = PetscAbsScalar(rval[j]); 251 } 252 if (ncols > cmax) cmax = ncols; 253 PetscCall(MatRestoreRow(A, i, &ncols, &rcol, &rval)); 254 } 255 PetscCall(PetscMalloc2(cmax, &pcols, cmax, &pvals)); 256 PetscCall(VecDestroy(&C)); 257 258 /* count the on and off processor sparsity patterns for the prolongator */ 259 for (i = 0; i < fn; i++) { 260 /* on */ 261 lsparse[i] = 0; 262 gsparse[i] = 0; 263 if (lcid[i] >= 0) { 264 lsparse[i] = 1; 265 gsparse[i] = 0; 266 } else { 267 PetscCall(MatGetRow(lA, i, &ncols, &rcol, &rval)); 268 for (j = 0; j < ncols; j++) { 269 col = rcol[j]; 270 if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0] * Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0] * Amax_neg[i])) lsparse[i] += 1; 271 } 272 PetscCall(MatRestoreRow(lA, i, &ncols, &rcol, &rval)); 273 /* off */ 274 if (gA) { 275 PetscCall(MatGetRow(gA, i, &ncols, &rcol, &rval)); 276 for (j = 0; j < ncols; j++) { 277 col = rcol[j]; 278 if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0] * Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0] * Amax_neg[i])) gsparse[i] += 1; 279 } 280 PetscCall(MatRestoreRow(gA, i, &ncols, &rcol, &rval)); 281 } 282 } 283 } 284 285 /* preallocate and create the prolongator */ 286 PetscCall(MatCreate(PetscObjectComm((PetscObject)A), P)); 287 PetscCall(MatGetType(G, &mtype)); 288 PetscCall(MatSetType(*P, mtype)); 289 PetscCall(MatSetSizes(*P, fn, cn, PETSC_DETERMINE, PETSC_DETERMINE)); 290 PetscCall(MatMPIAIJSetPreallocation(*P, 0, lsparse, 0, gsparse)); 291 PetscCall(MatSeqAIJSetPreallocation(*P, 0, lsparse)); 292 293 /* loop over local fine nodes -- get the diagonal, the sum of positive and negative strong and weak weights, and set up the row */ 294 for (i = 0; i < fn; i++) { 295 /* determine on or off */ 296 row_f = i + fs; 297 row_c = lcid[i]; 298 if (row_c >= 0) { 299 pij = 1.; 300 PetscCall(MatSetValues(*P, 1, &row_f, 1, &row_c, &pij, INSERT_VALUES)); 301 } else { 302 g_pos = 0.; 303 g_neg = 0.; 304 a_pos = 0.; 305 a_neg = 0.; 306 diag = 0.; 307 308 /* local connections */ 309 PetscCall(MatGetRow(lA, i, &ncols, &rcol, &rval)); 310 for (j = 0; j < ncols; j++) { 311 col = rcol[j]; 312 if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0] * Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0] * Amax_neg[i])) { 313 if (PetscRealPart(rval[j]) > 0.) { 314 g_pos += rval[j]; 315 } else { 316 g_neg += rval[j]; 317 } 318 } 319 if (col != i) { 320 if (PetscRealPart(rval[j]) > 0.) { 321 a_pos += rval[j]; 322 } else { 323 a_neg += rval[j]; 324 } 325 } else { 326 diag = rval[j]; 327 } 328 } 329 PetscCall(MatRestoreRow(lA, i, &ncols, &rcol, &rval)); 330 331 /* ghosted connections */ 332 if (gA) { 333 PetscCall(MatGetRow(gA, i, &ncols, &rcol, &rval)); 334 for (j = 0; j < ncols; j++) { 335 col = rcol[j]; 336 if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0] * Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0] * Amax_neg[i])) { 337 if (PetscRealPart(rval[j]) > 0.) { 338 g_pos += rval[j]; 339 } else { 340 g_neg += rval[j]; 341 } 342 } 343 if (PetscRealPart(rval[j]) > 0.) { 344 a_pos += rval[j]; 345 } else { 346 a_neg += rval[j]; 347 } 348 } 349 PetscCall(MatRestoreRow(gA, i, &ncols, &rcol, &rval)); 350 } 351 352 if (g_neg == 0.) { 353 alpha = 0.; 354 } else { 355 alpha = -a_neg / g_neg; 356 } 357 358 if (g_pos == 0.) { 359 diag += a_pos; 360 beta = 0.; 361 } else { 362 beta = -a_pos / g_pos; 363 } 364 if (diag == 0.) { 365 invdiag = 0.; 366 } else invdiag = 1. / diag; 367 /* on */ 368 PetscCall(MatGetRow(lA, i, &ncols, &rcol, &rval)); 369 idx = 0; 370 for (j = 0; j < ncols; j++) { 371 col = rcol[j]; 372 if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0] * Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0] * Amax_neg[i])) { 373 row_f = i + fs; 374 row_c = lcid[col]; 375 /* set the values for on-processor ones */ 376 if (PetscRealPart(rval[j]) < 0.) { 377 pij = rval[j] * alpha * invdiag; 378 } else { 379 pij = rval[j] * beta * invdiag; 380 } 381 if (PetscAbsScalar(pij) != 0.) { 382 pvals[idx] = pij; 383 pcols[idx] = row_c; 384 idx++; 385 } 386 } 387 } 388 PetscCall(MatRestoreRow(lA, i, &ncols, &rcol, &rval)); 389 /* off */ 390 if (gA) { 391 PetscCall(MatGetRow(gA, i, &ncols, &rcol, &rval)); 392 for (j = 0; j < ncols; j++) { 393 col = rcol[j]; 394 if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0] * Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0] * Amax_neg[i])) { 395 row_f = i + fs; 396 row_c = gcid[col]; 397 /* set the values for on-processor ones */ 398 if (PetscRealPart(rval[j]) < 0.) { 399 pij = rval[j] * alpha * invdiag; 400 } else { 401 pij = rval[j] * beta * invdiag; 402 } 403 if (PetscAbsScalar(pij) != 0.) { 404 pvals[idx] = pij; 405 pcols[idx] = row_c; 406 idx++; 407 } 408 } 409 } 410 PetscCall(MatRestoreRow(gA, i, &ncols, &rcol, &rval)); 411 } 412 PetscCall(MatSetValues(*P, 1, &row_f, idx, pcols, pvals, INSERT_VALUES)); 413 } 414 } 415 416 PetscCall(MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY)); 417 PetscCall(MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY)); 418 419 PetscCall(PetscFree5(lsparse, gsparse, lcid, Amax_pos, Amax_neg)); 420 421 PetscCall(PetscFree2(pcols, pvals)); 422 if (gA) { 423 PetscCall(PetscSFDestroy(&sf)); 424 PetscCall(PetscFree(gcid)); 425 } 426 PetscFunctionReturn(0); 427 } 428 429 PetscErrorCode PCGAMGTruncateProlongator_Private(PC pc, Mat *P) { 430 PetscInt j, i, ps, pf, pn, pcs, pcf, pcn, idx, cmax; 431 const PetscScalar *pval; 432 const PetscInt *pcol; 433 PetscScalar *pnval; 434 PetscInt *pncol; 435 PetscInt ncols; 436 Mat Pnew; 437 PetscInt *lsparse, *gsparse; 438 PetscReal pmax_pos, pmax_neg, ptot_pos, ptot_neg, pthresh_pos, pthresh_neg; 439 PC_MG *mg = (PC_MG *)pc->data; 440 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 441 PC_GAMG_Classical *cls = (PC_GAMG_Classical *)pc_gamg->subctx; 442 MatType mtype; 443 444 PetscFunctionBegin; 445 /* trim and rescale with reallocation */ 446 PetscCall(MatGetOwnershipRange(*P, &ps, &pf)); 447 PetscCall(MatGetOwnershipRangeColumn(*P, &pcs, &pcf)); 448 pn = pf - ps; 449 pcn = pcf - pcs; 450 PetscCall(PetscMalloc2(pn, &lsparse, pn, &gsparse)); 451 /* allocate */ 452 cmax = 0; 453 for (i = ps; i < pf; i++) { 454 lsparse[i - ps] = 0; 455 gsparse[i - ps] = 0; 456 PetscCall(MatGetRow(*P, i, &ncols, &pcol, &pval)); 457 if (ncols > cmax) cmax = ncols; 458 pmax_pos = 0.; 459 pmax_neg = 0.; 460 for (j = 0; j < ncols; j++) { 461 if (PetscRealPart(pval[j]) > pmax_pos) { 462 pmax_pos = PetscRealPart(pval[j]); 463 } else if (PetscRealPart(pval[j]) < pmax_neg) { 464 pmax_neg = PetscRealPart(pval[j]); 465 } 466 } 467 for (j = 0; j < ncols; j++) { 468 if (PetscRealPart(pval[j]) >= pmax_pos * cls->interp_threshold || PetscRealPart(pval[j]) <= pmax_neg * cls->interp_threshold) { 469 if (pcol[j] >= pcs && pcol[j] < pcf) { 470 lsparse[i - ps]++; 471 } else { 472 gsparse[i - ps]++; 473 } 474 } 475 } 476 PetscCall(MatRestoreRow(*P, i, &ncols, &pcol, &pval)); 477 } 478 479 PetscCall(PetscMalloc2(cmax, &pnval, cmax, &pncol)); 480 481 PetscCall(MatGetType(*P, &mtype)); 482 PetscCall(MatCreate(PetscObjectComm((PetscObject)*P), &Pnew)); 483 PetscCall(MatSetType(Pnew, mtype)); 484 PetscCall(MatSetSizes(Pnew, pn, pcn, PETSC_DETERMINE, PETSC_DETERMINE)); 485 PetscCall(MatSeqAIJSetPreallocation(Pnew, 0, lsparse)); 486 PetscCall(MatMPIAIJSetPreallocation(Pnew, 0, lsparse, 0, gsparse)); 487 488 for (i = ps; i < pf; i++) { 489 PetscCall(MatGetRow(*P, i, &ncols, &pcol, &pval)); 490 pmax_pos = 0.; 491 pmax_neg = 0.; 492 for (j = 0; j < ncols; j++) { 493 if (PetscRealPart(pval[j]) > pmax_pos) { 494 pmax_pos = PetscRealPart(pval[j]); 495 } else if (PetscRealPart(pval[j]) < pmax_neg) { 496 pmax_neg = PetscRealPart(pval[j]); 497 } 498 } 499 pthresh_pos = 0.; 500 pthresh_neg = 0.; 501 ptot_pos = 0.; 502 ptot_neg = 0.; 503 for (j = 0; j < ncols; j++) { 504 if (PetscRealPart(pval[j]) >= cls->interp_threshold * pmax_pos) { 505 pthresh_pos += PetscRealPart(pval[j]); 506 } else if (PetscRealPart(pval[j]) <= cls->interp_threshold * pmax_neg) { 507 pthresh_neg += PetscRealPart(pval[j]); 508 } 509 if (PetscRealPart(pval[j]) > 0.) { 510 ptot_pos += PetscRealPart(pval[j]); 511 } else { 512 ptot_neg += PetscRealPart(pval[j]); 513 } 514 } 515 if (PetscAbsReal(pthresh_pos) > 0.) ptot_pos /= pthresh_pos; 516 if (PetscAbsReal(pthresh_neg) > 0.) ptot_neg /= pthresh_neg; 517 idx = 0; 518 for (j = 0; j < ncols; j++) { 519 if (PetscRealPart(pval[j]) >= pmax_pos * cls->interp_threshold) { 520 pnval[idx] = ptot_pos * pval[j]; 521 pncol[idx] = pcol[j]; 522 idx++; 523 } else if (PetscRealPart(pval[j]) <= pmax_neg * cls->interp_threshold) { 524 pnval[idx] = ptot_neg * pval[j]; 525 pncol[idx] = pcol[j]; 526 idx++; 527 } 528 } 529 PetscCall(MatRestoreRow(*P, i, &ncols, &pcol, &pval)); 530 PetscCall(MatSetValues(Pnew, 1, &i, idx, pncol, pnval, INSERT_VALUES)); 531 } 532 533 PetscCall(MatAssemblyBegin(Pnew, MAT_FINAL_ASSEMBLY)); 534 PetscCall(MatAssemblyEnd(Pnew, MAT_FINAL_ASSEMBLY)); 535 PetscCall(MatDestroy(P)); 536 537 *P = Pnew; 538 PetscCall(PetscFree2(lsparse, gsparse)); 539 PetscCall(PetscFree2(pnval, pncol)); 540 PetscFunctionReturn(0); 541 } 542 543 PetscErrorCode PCGAMGProlongator_Classical_Standard(PC pc, Mat A, Mat G, PetscCoarsenData *agg_lists, Mat *P) { 544 Mat lA, *lAs; 545 MatType mtype; 546 Vec cv; 547 PetscInt *gcid, *lcid, *lsparse, *gsparse, *picol; 548 PetscInt fs, fe, cs, ce, nl, i, j, k, li, lni, ci, ncols, maxcols, fn, cn, cid; 549 PetscMPIInt size; 550 const PetscInt *lidx, *icol, *gidx; 551 PetscBool iscoarse; 552 PetscScalar vi, pentry, pjentry; 553 PetscScalar *pcontrib, *pvcol; 554 const PetscScalar *vcol; 555 PetscReal diag, jdiag, jwttotal; 556 PetscInt pncols; 557 PetscSF sf; 558 PetscLayout clayout; 559 IS lis; 560 561 PetscFunctionBegin; 562 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size)); 563 PetscCall(MatGetOwnershipRange(A, &fs, &fe)); 564 fn = fe - fs; 565 PetscCall(ISCreateStride(PETSC_COMM_SELF, fe - fs, fs, 1, &lis)); 566 if (size > 1) { 567 PetscCall(MatGetLayouts(A, NULL, &clayout)); 568 /* increase the overlap by two to get neighbors of neighbors */ 569 PetscCall(MatIncreaseOverlap(A, 1, &lis, 2)); 570 PetscCall(ISSort(lis)); 571 /* get the local part of A */ 572 PetscCall(MatCreateSubMatrices(A, 1, &lis, &lis, MAT_INITIAL_MATRIX, &lAs)); 573 lA = lAs[0]; 574 /* build an SF out of it */ 575 PetscCall(ISGetLocalSize(lis, &nl)); 576 PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)A), &sf)); 577 PetscCall(ISGetIndices(lis, &lidx)); 578 PetscCall(PetscSFSetGraphLayout(sf, clayout, nl, NULL, PETSC_COPY_VALUES, lidx)); 579 PetscCall(ISRestoreIndices(lis, &lidx)); 580 } else { 581 lA = A; 582 nl = fn; 583 } 584 /* create a communication structure for the overlapped portion and transmit coarse indices */ 585 PetscCall(PetscMalloc3(fn, &lsparse, fn, &gsparse, nl, &pcontrib)); 586 /* create coarse vector */ 587 cn = 0; 588 for (i = 0; i < fn; i++) { 589 PetscCall(PetscCDEmptyAt(agg_lists, i, &iscoarse)); 590 if (!iscoarse) cn++; 591 } 592 PetscCall(PetscMalloc1(fn, &gcid)); 593 PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)A), cn, PETSC_DECIDE, &cv)); 594 PetscCall(VecGetOwnershipRange(cv, &cs, &ce)); 595 cn = 0; 596 for (i = 0; i < fn; i++) { 597 PetscCall(PetscCDEmptyAt(agg_lists, i, &iscoarse)); 598 if (!iscoarse) { 599 gcid[i] = cs + cn; 600 cn++; 601 } else { 602 gcid[i] = -1; 603 } 604 } 605 if (size > 1) { 606 PetscCall(PetscMalloc1(nl, &lcid)); 607 PetscCall(PetscSFBcastBegin(sf, MPIU_INT, gcid, lcid, MPI_REPLACE)); 608 PetscCall(PetscSFBcastEnd(sf, MPIU_INT, gcid, lcid, MPI_REPLACE)); 609 } else { 610 lcid = gcid; 611 } 612 /* count to preallocate the prolongator */ 613 PetscCall(ISGetIndices(lis, &gidx)); 614 maxcols = 0; 615 /* count the number of unique contributing coarse cells for each fine */ 616 for (i = 0; i < nl; i++) { 617 pcontrib[i] = 0.; 618 PetscCall(MatGetRow(lA, i, &ncols, &icol, NULL)); 619 if (gidx[i] >= fs && gidx[i] < fe) { 620 li = gidx[i] - fs; 621 lsparse[li] = 0; 622 gsparse[li] = 0; 623 cid = lcid[i]; 624 if (cid >= 0) { 625 lsparse[li] = 1; 626 } else { 627 for (j = 0; j < ncols; j++) { 628 if (lcid[icol[j]] >= 0) { 629 pcontrib[icol[j]] = 1.; 630 } else { 631 ci = icol[j]; 632 PetscCall(MatRestoreRow(lA, i, &ncols, &icol, NULL)); 633 PetscCall(MatGetRow(lA, ci, &ncols, &icol, NULL)); 634 for (k = 0; k < ncols; k++) { 635 if (lcid[icol[k]] >= 0) pcontrib[icol[k]] = 1.; 636 } 637 PetscCall(MatRestoreRow(lA, ci, &ncols, &icol, NULL)); 638 PetscCall(MatGetRow(lA, i, &ncols, &icol, NULL)); 639 } 640 } 641 for (j = 0; j < ncols; j++) { 642 if (lcid[icol[j]] >= 0 && pcontrib[icol[j]] != 0.) { 643 lni = lcid[icol[j]]; 644 if (lni >= cs && lni < ce) { 645 lsparse[li]++; 646 } else { 647 gsparse[li]++; 648 } 649 pcontrib[icol[j]] = 0.; 650 } else { 651 ci = icol[j]; 652 PetscCall(MatRestoreRow(lA, i, &ncols, &icol, NULL)); 653 PetscCall(MatGetRow(lA, ci, &ncols, &icol, NULL)); 654 for (k = 0; k < ncols; k++) { 655 if (lcid[icol[k]] >= 0 && pcontrib[icol[k]] != 0.) { 656 lni = lcid[icol[k]]; 657 if (lni >= cs && lni < ce) { 658 lsparse[li]++; 659 } else { 660 gsparse[li]++; 661 } 662 pcontrib[icol[k]] = 0.; 663 } 664 } 665 PetscCall(MatRestoreRow(lA, ci, &ncols, &icol, NULL)); 666 PetscCall(MatGetRow(lA, i, &ncols, &icol, NULL)); 667 } 668 } 669 } 670 if (lsparse[li] + gsparse[li] > maxcols) maxcols = lsparse[li] + gsparse[li]; 671 } 672 PetscCall(MatRestoreRow(lA, i, &ncols, &icol, &vcol)); 673 } 674 PetscCall(PetscMalloc2(maxcols, &picol, maxcols, &pvcol)); 675 PetscCall(MatCreate(PetscObjectComm((PetscObject)A), P)); 676 PetscCall(MatGetType(A, &mtype)); 677 PetscCall(MatSetType(*P, mtype)); 678 PetscCall(MatSetSizes(*P, fn, cn, PETSC_DETERMINE, PETSC_DETERMINE)); 679 PetscCall(MatMPIAIJSetPreallocation(*P, 0, lsparse, 0, gsparse)); 680 PetscCall(MatSeqAIJSetPreallocation(*P, 0, lsparse)); 681 for (i = 0; i < nl; i++) { 682 diag = 0.; 683 if (gidx[i] >= fs && gidx[i] < fe) { 684 pncols = 0; 685 cid = lcid[i]; 686 if (cid >= 0) { 687 pncols = 1; 688 picol[0] = cid; 689 pvcol[0] = 1.; 690 } else { 691 PetscCall(MatGetRow(lA, i, &ncols, &icol, &vcol)); 692 for (j = 0; j < ncols; j++) { 693 pentry = vcol[j]; 694 if (lcid[icol[j]] >= 0) { 695 /* coarse neighbor */ 696 pcontrib[icol[j]] += pentry; 697 } else if (icol[j] != i) { 698 /* the neighbor is a strongly connected fine node */ 699 ci = icol[j]; 700 vi = vcol[j]; 701 PetscCall(MatRestoreRow(lA, i, &ncols, &icol, &vcol)); 702 PetscCall(MatGetRow(lA, ci, &ncols, &icol, &vcol)); 703 jwttotal = 0.; 704 jdiag = 0.; 705 for (k = 0; k < ncols; k++) { 706 if (ci == icol[k]) jdiag = PetscRealPart(vcol[k]); 707 } 708 for (k = 0; k < ncols; k++) { 709 if (lcid[icol[k]] >= 0 && jdiag * PetscRealPart(vcol[k]) < 0.) { 710 pjentry = vcol[k]; 711 jwttotal += PetscRealPart(pjentry); 712 } 713 } 714 if (jwttotal != 0.) { 715 jwttotal = PetscRealPart(vi) / jwttotal; 716 for (k = 0; k < ncols; k++) { 717 if (lcid[icol[k]] >= 0 && jdiag * PetscRealPart(vcol[k]) < 0.) { 718 pjentry = vcol[k] * jwttotal; 719 pcontrib[icol[k]] += pjentry; 720 } 721 } 722 } else { 723 diag += PetscRealPart(vi); 724 } 725 PetscCall(MatRestoreRow(lA, ci, &ncols, &icol, &vcol)); 726 PetscCall(MatGetRow(lA, i, &ncols, &icol, &vcol)); 727 } else { 728 diag += PetscRealPart(vcol[j]); 729 } 730 } 731 if (diag != 0.) { 732 diag = 1. / diag; 733 for (j = 0; j < ncols; j++) { 734 if (lcid[icol[j]] >= 0 && pcontrib[icol[j]] != 0.) { 735 /* the neighbor is a coarse node */ 736 if (PetscAbsScalar(pcontrib[icol[j]]) > 0.0) { 737 lni = lcid[icol[j]]; 738 pvcol[pncols] = -pcontrib[icol[j]] * diag; 739 picol[pncols] = lni; 740 pncols++; 741 } 742 pcontrib[icol[j]] = 0.; 743 } else { 744 /* the neighbor is a strongly connected fine node */ 745 ci = icol[j]; 746 PetscCall(MatRestoreRow(lA, i, &ncols, &icol, &vcol)); 747 PetscCall(MatGetRow(lA, ci, &ncols, &icol, &vcol)); 748 for (k = 0; k < ncols; k++) { 749 if (lcid[icol[k]] >= 0 && pcontrib[icol[k]] != 0.) { 750 if (PetscAbsScalar(pcontrib[icol[k]]) > 0.0) { 751 lni = lcid[icol[k]]; 752 pvcol[pncols] = -pcontrib[icol[k]] * diag; 753 picol[pncols] = lni; 754 pncols++; 755 } 756 pcontrib[icol[k]] = 0.; 757 } 758 } 759 PetscCall(MatRestoreRow(lA, ci, &ncols, &icol, &vcol)); 760 PetscCall(MatGetRow(lA, i, &ncols, &icol, &vcol)); 761 } 762 pcontrib[icol[j]] = 0.; 763 } 764 PetscCall(MatRestoreRow(lA, i, &ncols, &icol, &vcol)); 765 } 766 } 767 ci = gidx[i]; 768 if (pncols > 0) PetscCall(MatSetValues(*P, 1, &ci, pncols, picol, pvcol, INSERT_VALUES)); 769 } 770 } 771 PetscCall(ISRestoreIndices(lis, &gidx)); 772 PetscCall(PetscFree2(picol, pvcol)); 773 PetscCall(PetscFree3(lsparse, gsparse, pcontrib)); 774 PetscCall(ISDestroy(&lis)); 775 PetscCall(PetscFree(gcid)); 776 if (size > 1) { 777 PetscCall(PetscFree(lcid)); 778 PetscCall(MatDestroyMatrices(1, &lAs)); 779 PetscCall(PetscSFDestroy(&sf)); 780 } 781 PetscCall(VecDestroy(&cv)); 782 PetscCall(MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY)); 783 PetscCall(MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY)); 784 PetscFunctionReturn(0); 785 } 786 787 PetscErrorCode PCGAMGOptProlongator_Classical_Jacobi(PC pc, Mat A, Mat *P) { 788 PetscInt f, s, n, cf, cs, i, idx; 789 PetscInt *coarserows; 790 PetscInt ncols; 791 const PetscInt *pcols; 792 const PetscScalar *pvals; 793 Mat Pnew; 794 Vec diag; 795 PC_MG *mg = (PC_MG *)pc->data; 796 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 797 PC_GAMG_Classical *cls = (PC_GAMG_Classical *)pc_gamg->subctx; 798 799 PetscFunctionBegin; 800 if (cls->nsmooths == 0) { 801 PetscCall(PCGAMGTruncateProlongator_Private(pc, P)); 802 PetscFunctionReturn(0); 803 } 804 PetscCall(MatGetOwnershipRange(*P, &s, &f)); 805 n = f - s; 806 PetscCall(MatGetOwnershipRangeColumn(*P, &cs, &cf)); 807 PetscCall(PetscMalloc1(n, &coarserows)); 808 /* identify the rows corresponding to coarse unknowns */ 809 idx = 0; 810 for (i = s; i < f; i++) { 811 PetscCall(MatGetRow(*P, i, &ncols, &pcols, &pvals)); 812 /* assume, for now, that it's a coarse unknown if it has a single unit entry */ 813 if (ncols == 1) { 814 if (pvals[0] == 1.) { 815 coarserows[idx] = i; 816 idx++; 817 } 818 } 819 PetscCall(MatRestoreRow(*P, i, &ncols, &pcols, &pvals)); 820 } 821 PetscCall(MatCreateVecs(A, &diag, NULL)); 822 PetscCall(MatGetDiagonal(A, diag)); 823 PetscCall(VecReciprocal(diag)); 824 for (i = 0; i < cls->nsmooths; i++) { 825 PetscCall(MatMatMult(A, *P, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Pnew)); 826 PetscCall(MatZeroRows(Pnew, idx, coarserows, 0., NULL, NULL)); 827 PetscCall(MatDiagonalScale(Pnew, diag, NULL)); 828 PetscCall(MatAYPX(Pnew, -1.0, *P, DIFFERENT_NONZERO_PATTERN)); 829 PetscCall(MatDestroy(P)); 830 *P = Pnew; 831 Pnew = NULL; 832 } 833 PetscCall(VecDestroy(&diag)); 834 PetscCall(PetscFree(coarserows)); 835 PetscCall(PCGAMGTruncateProlongator_Private(pc, P)); 836 PetscFunctionReturn(0); 837 } 838 839 static PetscErrorCode PCGAMGProlongator_Classical(PC pc, Mat A, Mat G, PetscCoarsenData *agg_lists, Mat *P) { 840 PetscErrorCode (*f)(PC, Mat, Mat, PetscCoarsenData *, Mat *); 841 PC_MG *mg = (PC_MG *)pc->data; 842 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 843 PC_GAMG_Classical *cls = (PC_GAMG_Classical *)pc_gamg->subctx; 844 845 PetscFunctionBegin; 846 PetscCall(PetscFunctionListFind(PCGAMGClassicalProlongatorList, cls->prolongtype, &f)); 847 PetscCheck(f, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Cannot find PCGAMG Classical prolongator type"); 848 PetscCall((*f)(pc, A, G, agg_lists, P)); 849 PetscFunctionReturn(0); 850 } 851 852 static PetscErrorCode PCGAMGDestroy_Classical(PC pc) { 853 PC_MG *mg = (PC_MG *)pc->data; 854 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 855 856 PetscFunctionBegin; 857 PetscCall(PetscFree(pc_gamg->subctx)); 858 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGClassicalSetType_C", NULL)); 859 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGClassicalGetType_C", NULL)); 860 PetscFunctionReturn(0); 861 } 862 863 PetscErrorCode PCGAMGSetFromOptions_Classical(PC pc, PetscOptionItems *PetscOptionsObject) { 864 PC_MG *mg = (PC_MG *)pc->data; 865 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 866 PC_GAMG_Classical *cls = (PC_GAMG_Classical *)pc_gamg->subctx; 867 char tname[256]; 868 PetscBool flg; 869 870 PetscFunctionBegin; 871 PetscOptionsHeadBegin(PetscOptionsObject, "GAMG-Classical options"); 872 PetscCall(PetscOptionsFList("-pc_gamg_classical_type", "Type of Classical AMG prolongation", "PCGAMGClassicalSetType", PCGAMGClassicalProlongatorList, cls->prolongtype, tname, sizeof(tname), &flg)); 873 if (flg) PetscCall(PCGAMGClassicalSetType(pc, tname)); 874 PetscCall(PetscOptionsReal("-pc_gamg_classical_interp_threshold", "Threshold for classical interpolator entries", "", cls->interp_threshold, &cls->interp_threshold, NULL)); 875 PetscCall(PetscOptionsInt("-pc_gamg_classical_nsmooths", "Threshold for classical interpolator entries", "", cls->nsmooths, &cls->nsmooths, NULL)); 876 PetscOptionsHeadEnd(); 877 PetscFunctionReturn(0); 878 } 879 880 static PetscErrorCode PCGAMGSetData_Classical(PC pc, Mat A) { 881 PC_MG *mg = (PC_MG *)pc->data; 882 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 883 884 PetscFunctionBegin; 885 /* no data for classical AMG */ 886 pc_gamg->data = NULL; 887 pc_gamg->data_cell_cols = 0; 888 pc_gamg->data_cell_rows = 0; 889 pc_gamg->data_sz = 0; 890 PetscFunctionReturn(0); 891 } 892 893 PetscErrorCode PCGAMGClassicalFinalizePackage(void) { 894 PetscFunctionBegin; 895 PCGAMGClassicalPackageInitialized = PETSC_FALSE; 896 PetscCall(PetscFunctionListDestroy(&PCGAMGClassicalProlongatorList)); 897 PetscFunctionReturn(0); 898 } 899 900 PetscErrorCode PCGAMGClassicalInitializePackage(void) { 901 PetscFunctionBegin; 902 if (PCGAMGClassicalPackageInitialized) PetscFunctionReturn(0); 903 PetscCall(PetscFunctionListAdd(&PCGAMGClassicalProlongatorList, PCGAMGCLASSICALDIRECT, PCGAMGProlongator_Classical_Direct)); 904 PetscCall(PetscFunctionListAdd(&PCGAMGClassicalProlongatorList, PCGAMGCLASSICALSTANDARD, PCGAMGProlongator_Classical_Standard)); 905 PetscCall(PetscRegisterFinalize(PCGAMGClassicalFinalizePackage)); 906 PetscFunctionReturn(0); 907 } 908 909 /* -------------------------------------------------------------------------- */ 910 /* 911 PCCreateGAMG_Classical 912 913 */ 914 PetscErrorCode PCCreateGAMG_Classical(PC pc) { 915 PC_MG *mg = (PC_MG *)pc->data; 916 PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 917 PC_GAMG_Classical *pc_gamg_classical; 918 919 PetscFunctionBegin; 920 PetscCall(PCGAMGClassicalInitializePackage()); 921 if (pc_gamg->subctx) { 922 /* call base class */ 923 PetscCall(PCDestroy_GAMG(pc)); 924 } 925 926 /* create sub context for SA */ 927 PetscCall(PetscNewLog(pc, &pc_gamg_classical)); 928 pc_gamg->subctx = pc_gamg_classical; 929 pc->ops->setfromoptions = PCGAMGSetFromOptions_Classical; 930 /* reset does not do anything; setup not virtual */ 931 932 /* set internal function pointers */ 933 pc_gamg->ops->destroy = PCGAMGDestroy_Classical; 934 pc_gamg->ops->graph = PCGAMGGraph_Classical; 935 pc_gamg->ops->coarsen = PCGAMGCoarsen_Classical; 936 pc_gamg->ops->prolongator = PCGAMGProlongator_Classical; 937 pc_gamg->ops->optprolongator = PCGAMGOptProlongator_Classical_Jacobi; 938 pc_gamg->ops->setfromoptions = PCGAMGSetFromOptions_Classical; 939 940 pc_gamg->ops->createdefaultdata = PCGAMGSetData_Classical; 941 pc_gamg_classical->interp_threshold = 0.2; 942 pc_gamg_classical->nsmooths = 0; 943 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGClassicalSetType_C", PCGAMGClassicalSetType_GAMG)); 944 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGClassicalGetType_C", PCGAMGClassicalGetType_GAMG)); 945 PetscCall(PCGAMGClassicalSetType(pc, PCGAMGCLASSICALSTANDARD)); 946 PetscFunctionReturn(0); 947 } 948