18e6d0c30SPeter Brune #include <../src/ksp/pc/impls/gamg/gamg.h> /*I "petscpc.h" I*/ 287b9b13cSPeter Brune #include <petscsf.h> 38e6d0c30SPeter Brune 4ff6a9541SJacob Faibussowitsch static PetscFunctionList PCGAMGClassicalProlongatorList = NULL; 5ff6a9541SJacob Faibussowitsch static PetscBool PCGAMGClassicalPackageInitialized = PETSC_FALSE; 68eab0cc1SPeter Brune 78e6d0c30SPeter Brune typedef struct { 8586a8384SPeter Brune PetscReal interp_threshold; /* interpolation threshold */ 98eab0cc1SPeter Brune char prolongtype[256]; 10b4dc3ebdSPeter Brune PetscInt nsmooths; /* number of jacobi smoothings on the prolongator */ 118e6d0c30SPeter Brune } PC_GAMG_Classical; 128e6d0c30SPeter Brune 137779008dSPeter Brune /*@C 14f1580f4eSBarry Smith PCGAMGClassicalSetType - Sets the type of classical interpolation to use with `PCGAMG` 158eab0cc1SPeter Brune 16c3339decSBarry Smith Collective 178eab0cc1SPeter Brune 1804c3f3b8SBarry Smith Input Parameters: 1904c3f3b8SBarry Smith + pc - the preconditioner context 2004c3f3b8SBarry Smith - type - the interpolation to use, see `PCGAMGClassicalType()` 218eab0cc1SPeter Brune 228eab0cc1SPeter Brune Options Database Key: 23f1580f4eSBarry Smith . -pc_gamg_classical_type <direct,standard> - set type of classical AMG prolongation 248eab0cc1SPeter Brune 258eab0cc1SPeter Brune Level: intermediate 268eab0cc1SPeter Brune 27*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`, `PCGAMGClassicalType`, `PCGAMGClassicalGetType()` 288eab0cc1SPeter Brune @*/ 29d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGClassicalSetType(PC pc, PCGAMGClassicalType type) 30d71ae5a4SJacob Faibussowitsch { 318eab0cc1SPeter Brune PetscFunctionBegin; 328eab0cc1SPeter Brune PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 33cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGClassicalSetType_C", (PC, PCGAMGClassicalType), (pc, type)); 343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 358eab0cc1SPeter Brune } 368eab0cc1SPeter Brune 37c60c7ad4SBarry Smith /*@C 38f1580f4eSBarry Smith PCGAMGClassicalGetType - Gets the type of classical interpolation to use with `PCGAMG` 39c60c7ad4SBarry Smith 40c3339decSBarry Smith Collective 41c60c7ad4SBarry Smith 42c60c7ad4SBarry Smith Input Parameter: 43c60c7ad4SBarry Smith . pc - the preconditioner context 44c60c7ad4SBarry Smith 45c60c7ad4SBarry Smith Output Parameter: 4604c3f3b8SBarry Smith . type - the type used, see `PCGAMGClassicalType()` 47c60c7ad4SBarry Smith 48c60c7ad4SBarry Smith Level: intermediate 49c60c7ad4SBarry Smith 50*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`, `PCGAMGClassicalType`, `PCGAMGClassicalSetType()` 51c60c7ad4SBarry Smith @*/ 52d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGClassicalGetType(PC pc, PCGAMGClassicalType *type) 53d71ae5a4SJacob Faibussowitsch { 54c60c7ad4SBarry Smith PetscFunctionBegin; 55c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 56cac4c232SBarry Smith PetscUseMethod(pc, "PCGAMGClassicalGetType_C", (PC, PCGAMGClassicalType *), (pc, type)); 573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 58c60c7ad4SBarry Smith } 59c60c7ad4SBarry Smith 60d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGClassicalSetType_GAMG(PC pc, PCGAMGClassicalType type) 61d71ae5a4SJacob Faibussowitsch { 628eab0cc1SPeter Brune PC_MG *mg = (PC_MG *)pc->data; 638eab0cc1SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 648eab0cc1SPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical *)pc_gamg->subctx; 658eab0cc1SPeter Brune 668eab0cc1SPeter Brune PetscFunctionBegin; 67c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(cls->prolongtype, type, sizeof(cls->prolongtype))); 683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 698eab0cc1SPeter Brune } 708e6d0c30SPeter Brune 71d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGClassicalGetType_GAMG(PC pc, PCGAMGClassicalType *type) 72d71ae5a4SJacob Faibussowitsch { 73c60c7ad4SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 74c60c7ad4SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 75c60c7ad4SBarry Smith PC_GAMG_Classical *cls = (PC_GAMG_Classical *)pc_gamg->subctx; 76c60c7ad4SBarry Smith 77c60c7ad4SBarry Smith PetscFunctionBegin; 78c60c7ad4SBarry Smith *type = cls->prolongtype; 793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 80c60c7ad4SBarry Smith } 81c60c7ad4SBarry Smith 82ff6a9541SJacob Faibussowitsch static PetscErrorCode PCGAMGCreateGraph_Classical(PC pc, Mat A, Mat *G) 83d71ae5a4SJacob Faibussowitsch { 84550383edSPeter Brune PetscInt s, f, n, idx, lidx, gidx; 85e5a0faa4SPeter Brune PetscInt r, c, ncols; 868e6d0c30SPeter Brune const PetscInt *rcol; 878e6d0c30SPeter Brune const PetscScalar *rval; 88e5a0faa4SPeter Brune PetscInt *gcol; 898e6d0c30SPeter Brune PetscScalar *gval; 90e5a0faa4SPeter Brune PetscReal rmax; 91550383edSPeter Brune PetscInt cmax = 0; 92b817416eSBarry Smith PC_MG *mg = (PC_MG *)pc->data; 93b817416eSBarry Smith PC_GAMG *gamg = (PC_GAMG *)mg->innerctx; 948e6d0c30SPeter Brune PetscInt *gsparse, *lsparse; 95e5a0faa4SPeter Brune PetscScalar *Amax; 968e6d0c30SPeter Brune MatType mtype; 978e6d0c30SPeter Brune 988e6d0c30SPeter Brune PetscFunctionBegin; 999566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &s, &f)); 100550383edSPeter Brune n = f - s; 1019566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(n, &lsparse, n, &gsparse, n, &Amax)); 1028e6d0c30SPeter Brune 103550383edSPeter Brune for (r = 0; r < n; r++) { 1048e6d0c30SPeter Brune lsparse[r] = 0; 105550383edSPeter Brune gsparse[r] = 0; 1068e6d0c30SPeter Brune } 1078e6d0c30SPeter Brune 108550383edSPeter Brune for (r = s; r < f; r++) { 109e5a0faa4SPeter Brune /* determine the maximum off-diagonal in each row */ 110e5a0faa4SPeter Brune rmax = 0.; 1119566063dSJacob Faibussowitsch PetscCall(MatGetRow(A, r, &ncols, &rcol, &rval)); 112e5a0faa4SPeter Brune for (c = 0; c < ncols; c++) { 113ad540459SPierre Jolivet if (PetscRealPart(-rval[c]) > rmax && rcol[c] != r) rmax = PetscRealPart(-rval[c]); 114e5a0faa4SPeter Brune } 115550383edSPeter Brune Amax[r - s] = rmax; 116550383edSPeter Brune if (ncols > cmax) cmax = ncols; 117550383edSPeter Brune lidx = 0; 118550383edSPeter Brune gidx = 0; 119e5a0faa4SPeter Brune /* create the local and global sparsity patterns */ 1208e6d0c30SPeter Brune for (c = 0; c < ncols; c++) { 121c1eae691SMark Adams if (PetscRealPart(-rval[c]) > gamg->threshold[0] * PetscRealPart(Amax[r - s]) || rcol[c] == r) { 122550383edSPeter Brune if (rcol[c] < f && rcol[c] >= s) { 123550383edSPeter Brune lidx++; 124550383edSPeter Brune } else { 125550383edSPeter Brune gidx++; 1268e6d0c30SPeter Brune } 1278e6d0c30SPeter Brune } 1288e6d0c30SPeter Brune } 1299566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A, r, &ncols, &rcol, &rval)); 130550383edSPeter Brune lsparse[r - s] = lidx; 131550383edSPeter Brune gsparse[r - s] = gidx; 1328e6d0c30SPeter Brune } 1339566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(cmax, &gval, cmax, &gcol)); 134e5a0faa4SPeter Brune 1359566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), G)); 1369566063dSJacob Faibussowitsch PetscCall(MatGetType(A, &mtype)); 1379566063dSJacob Faibussowitsch PetscCall(MatSetType(*G, mtype)); 1389566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*G, n, n, PETSC_DETERMINE, PETSC_DETERMINE)); 1399566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(*G, 0, lsparse, 0, gsparse)); 1409566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(*G, 0, lsparse)); 1418e6d0c30SPeter Brune for (r = s; r < f; r++) { 1429566063dSJacob Faibussowitsch PetscCall(MatGetRow(A, r, &ncols, &rcol, &rval)); 1438e6d0c30SPeter Brune idx = 0; 1448e6d0c30SPeter Brune for (c = 0; c < ncols; c++) { 1458e6d0c30SPeter Brune /* classical strength of connection */ 146c1eae691SMark Adams if (PetscRealPart(-rval[c]) > gamg->threshold[0] * PetscRealPart(Amax[r - s]) || rcol[c] == r) { 1478e6d0c30SPeter Brune gcol[idx] = rcol[c]; 1488e6d0c30SPeter Brune gval[idx] = rval[c]; 1498e6d0c30SPeter Brune idx++; 1508e6d0c30SPeter Brune } 1518e6d0c30SPeter Brune } 1529566063dSJacob Faibussowitsch PetscCall(MatSetValues(*G, 1, &r, idx, gcol, gval, INSERT_VALUES)); 1539566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A, r, &ncols, &rcol, &rval)); 1548e6d0c30SPeter Brune } 1559566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*G, MAT_FINAL_ASSEMBLY)); 1569566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*G, MAT_FINAL_ASSEMBLY)); 1578e6d0c30SPeter Brune 1589566063dSJacob Faibussowitsch PetscCall(PetscFree2(gval, gcol)); 1599566063dSJacob Faibussowitsch PetscCall(PetscFree3(lsparse, gsparse, Amax)); 1603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1618e6d0c30SPeter Brune } 1628e6d0c30SPeter Brune 163ff6a9541SJacob Faibussowitsch static PetscErrorCode PCGAMGCoarsen_Classical(PC pc, Mat *G, PetscCoarsenData **agg_lists) 164d71ae5a4SJacob Faibussowitsch { 1658e6d0c30SPeter Brune MatCoarsen crs; 1668e6d0c30SPeter Brune MPI_Comm fcomm = ((PetscObject)pc)->comm; 1678e6d0c30SPeter Brune 1688e6d0c30SPeter Brune PetscFunctionBegin; 16928b400f6SJacob Faibussowitsch PetscCheck(G, fcomm, PETSC_ERR_ARG_WRONGSTATE, "Must set Graph in PC in PCGAMG before coarsening"); 1708e6d0c30SPeter Brune 1719566063dSJacob Faibussowitsch PetscCall(MatCoarsenCreate(fcomm, &crs)); 1729566063dSJacob Faibussowitsch PetscCall(MatCoarsenSetFromOptions(crs)); 1739566063dSJacob Faibussowitsch PetscCall(MatCoarsenSetAdjacency(crs, *G)); 1749566063dSJacob Faibussowitsch PetscCall(MatCoarsenSetStrictAggs(crs, PETSC_TRUE)); 1759566063dSJacob Faibussowitsch PetscCall(MatCoarsenApply(crs)); 1769566063dSJacob Faibussowitsch PetscCall(MatCoarsenGetData(crs, agg_lists)); 1779566063dSJacob Faibussowitsch PetscCall(MatCoarsenDestroy(&crs)); 1783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1798e6d0c30SPeter Brune } 1808e6d0c30SPeter Brune 181ff6a9541SJacob Faibussowitsch static PetscErrorCode PCGAMGProlongator_Classical_Direct(PC pc, Mat A, Mat G, PetscCoarsenData *agg_lists, Mat *P) 182d71ae5a4SJacob Faibussowitsch { 1831ce39c63SPeter Brune PC_MG *mg = (PC_MG *)pc->data; 1841ce39c63SPeter Brune PC_GAMG *gamg = (PC_GAMG *)mg->innerctx; 18563b0c588SPeter Brune PetscBool iscoarse, isMPIAIJ, isSEQAIJ; 18663b0c588SPeter Brune PetscInt fn, cn, fs, fe, cs, ce, i, j, ncols, col, row_f, row_c, cmax = 0, idx, noff; 18763b0c588SPeter Brune PetscInt *lcid, *gcid, *lsparse, *gsparse, *colmap, *pcols; 18863b0c588SPeter Brune const PetscInt *rcol; 18963b0c588SPeter Brune PetscReal *Amax_pos, *Amax_neg; 19063b0c588SPeter Brune PetscScalar g_pos, g_neg, a_pos, a_neg, diag, invdiag, alpha, beta, pij; 19163b0c588SPeter Brune PetscScalar *pvals; 19263b0c588SPeter Brune const PetscScalar *rval; 19363b0c588SPeter Brune Mat lA, gA = NULL; 19463b0c588SPeter Brune MatType mtype; 19563b0c588SPeter Brune Vec C, lvec; 19687b9b13cSPeter Brune PetscLayout clayout; 19787b9b13cSPeter Brune PetscSF sf; 19887b9b13cSPeter Brune Mat_MPIAIJ *mpiaij; 1998e6d0c30SPeter Brune 2008e6d0c30SPeter Brune PetscFunctionBegin; 2019566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &fs, &fe)); 20287b9b13cSPeter Brune fn = fe - fs; 2039566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPIAIJ, &isMPIAIJ)); 2049566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQAIJ, &isSEQAIJ)); 2057827d75bSBarry Smith PetscCheck(isMPIAIJ || isSEQAIJ, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Classical AMG requires MPIAIJ matrix"); 20687b9b13cSPeter Brune if (isMPIAIJ) { 20787b9b13cSPeter Brune mpiaij = (Mat_MPIAIJ *)A->data; 20887b9b13cSPeter Brune lA = mpiaij->A; 20987b9b13cSPeter Brune gA = mpiaij->B; 21087b9b13cSPeter Brune lvec = mpiaij->lvec; 2119566063dSJacob Faibussowitsch PetscCall(VecGetSize(lvec, &noff)); 21287b9b13cSPeter Brune colmap = mpiaij->garray; 2139566063dSJacob Faibussowitsch PetscCall(MatGetLayouts(A, NULL, &clayout)); 2149566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)A), &sf)); 2159566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraphLayout(sf, clayout, noff, NULL, PETSC_COPY_VALUES, colmap)); 2169566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(noff, &gcid)); 21787b9b13cSPeter Brune } else { 21887b9b13cSPeter Brune lA = A; 21987b9b13cSPeter Brune } 2209566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(fn, &lsparse, fn, &gsparse, fn, &lcid, fn, &Amax_pos, fn, &Amax_neg)); 2218e6d0c30SPeter Brune 2228e6d0c30SPeter Brune /* count the number of coarse unknowns */ 2238e6d0c30SPeter Brune cn = 0; 2248e6d0c30SPeter Brune for (i = 0; i < fn; i++) { 2258e6d0c30SPeter Brune /* filter out singletons */ 2269566063dSJacob Faibussowitsch PetscCall(PetscCDEmptyAt(agg_lists, i, &iscoarse)); 2278e6d0c30SPeter Brune lcid[i] = -1; 228ad540459SPierre Jolivet if (!iscoarse) cn++; 2298e6d0c30SPeter Brune } 2308e6d0c30SPeter Brune 2318e6d0c30SPeter Brune /* create the coarse vector */ 2329566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)A), cn, PETSC_DECIDE, &C)); 2339566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(C, &cs, &ce)); 2348e6d0c30SPeter Brune 2358e6d0c30SPeter Brune cn = 0; 2368e6d0c30SPeter Brune for (i = 0; i < fn; i++) { 2379566063dSJacob Faibussowitsch PetscCall(PetscCDEmptyAt(agg_lists, i, &iscoarse)); 2388e6d0c30SPeter Brune if (!iscoarse) { 2398e6d0c30SPeter Brune lcid[i] = cs + cn; 2408e6d0c30SPeter Brune cn++; 2418e6d0c30SPeter Brune } else { 2428e6d0c30SPeter Brune lcid[i] = -1; 2438e6d0c30SPeter Brune } 2448e6d0c30SPeter Brune } 2458e6d0c30SPeter Brune 24687b9b13cSPeter Brune if (gA) { 2479566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sf, MPIU_INT, lcid, gcid, MPI_REPLACE)); 2489566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sf, MPIU_INT, lcid, gcid, MPI_REPLACE)); 24987b9b13cSPeter Brune } 2508e6d0c30SPeter Brune 251b817416eSBarry Smith /* determine the largest off-diagonal entries in each row */ 2521ce39c63SPeter Brune for (i = fs; i < fe; i++) { 2531ce39c63SPeter Brune Amax_pos[i - fs] = 0.; 2541ce39c63SPeter Brune Amax_neg[i - fs] = 0.; 2559566063dSJacob Faibussowitsch PetscCall(MatGetRow(A, i, &ncols, &rcol, &rval)); 2561ce39c63SPeter Brune for (j = 0; j < ncols; j++) { 2571ce39c63SPeter Brune if ((PetscRealPart(-rval[j]) > Amax_neg[i - fs]) && i != rcol[j]) Amax_neg[i - fs] = PetscAbsScalar(rval[j]); 2581ce39c63SPeter Brune if ((PetscRealPart(rval[j]) > Amax_pos[i - fs]) && i != rcol[j]) Amax_pos[i - fs] = PetscAbsScalar(rval[j]); 2591ce39c63SPeter Brune } 2601ce39c63SPeter Brune if (ncols > cmax) cmax = ncols; 2619566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A, i, &ncols, &rcol, &rval)); 2621ce39c63SPeter Brune } 2639566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(cmax, &pcols, cmax, &pvals)); 2649566063dSJacob Faibussowitsch PetscCall(VecDestroy(&C)); 2658e6d0c30SPeter Brune 2668e6d0c30SPeter Brune /* count the on and off processor sparsity patterns for the prolongator */ 2678e6d0c30SPeter Brune for (i = 0; i < fn; i++) { 2688e6d0c30SPeter Brune /* on */ 2698e6d0c30SPeter Brune lsparse[i] = 0; 270e5a0faa4SPeter Brune gsparse[i] = 0; 2718e6d0c30SPeter Brune if (lcid[i] >= 0) { 2728e6d0c30SPeter Brune lsparse[i] = 1; 2738e6d0c30SPeter Brune gsparse[i] = 0; 2748e6d0c30SPeter Brune } else { 2759566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA, i, &ncols, &rcol, &rval)); 2768e6d0c30SPeter Brune for (j = 0; j < ncols; j++) { 2771ce39c63SPeter Brune col = rcol[j]; 278ad540459SPierre Jolivet 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; 2798e6d0c30SPeter Brune } 2809566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA, i, &ncols, &rcol, &rval)); 2818e6d0c30SPeter Brune /* off */ 2821ce39c63SPeter Brune if (gA) { 2839566063dSJacob Faibussowitsch PetscCall(MatGetRow(gA, i, &ncols, &rcol, &rval)); 2848e6d0c30SPeter Brune for (j = 0; j < ncols; j++) { 2851ce39c63SPeter Brune col = rcol[j]; 286ad540459SPierre Jolivet 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; 2878e6d0c30SPeter Brune } 2889566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(gA, i, &ncols, &rcol, &rval)); 2898e6d0c30SPeter Brune } 2908e6d0c30SPeter Brune } 2911ce39c63SPeter Brune } 2928e6d0c30SPeter Brune 2938e6d0c30SPeter Brune /* preallocate and create the prolongator */ 2949566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), P)); 2959566063dSJacob Faibussowitsch PetscCall(MatGetType(G, &mtype)); 2969566063dSJacob Faibussowitsch PetscCall(MatSetType(*P, mtype)); 2979566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*P, fn, cn, PETSC_DETERMINE, PETSC_DETERMINE)); 2989566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(*P, 0, lsparse, 0, gsparse)); 2999566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(*P, 0, lsparse)); 3008e6d0c30SPeter Brune 3018e6d0c30SPeter Brune /* loop over local fine nodes -- get the diagonal, the sum of positive and negative strong and weak weights, and set up the row */ 3028e6d0c30SPeter Brune for (i = 0; i < fn; i++) { 3038e6d0c30SPeter Brune /* determine on or off */ 3048e6d0c30SPeter Brune row_f = i + fs; 3058e6d0c30SPeter Brune row_c = lcid[i]; 3068e6d0c30SPeter Brune if (row_c >= 0) { 3078e6d0c30SPeter Brune pij = 1.; 3089566063dSJacob Faibussowitsch PetscCall(MatSetValues(*P, 1, &row_f, 1, &row_c, &pij, INSERT_VALUES)); 3098e6d0c30SPeter Brune } else { 3108e6d0c30SPeter Brune g_pos = 0.; 3118e6d0c30SPeter Brune g_neg = 0.; 3128e6d0c30SPeter Brune a_pos = 0.; 3138e6d0c30SPeter Brune a_neg = 0.; 3148e6d0c30SPeter Brune diag = 0.; 3158e6d0c30SPeter Brune 3161ce39c63SPeter Brune /* local connections */ 3179566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA, i, &ncols, &rcol, &rval)); 3181ce39c63SPeter Brune for (j = 0; j < ncols; j++) { 3191ce39c63SPeter Brune col = rcol[j]; 320c1eae691SMark Adams if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0] * Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0] * Amax_neg[i])) { 3211ce39c63SPeter Brune if (PetscRealPart(rval[j]) > 0.) { 3221ce39c63SPeter Brune g_pos += rval[j]; 3238e6d0c30SPeter Brune } else { 3241ce39c63SPeter Brune g_neg += rval[j]; 3258e6d0c30SPeter Brune } 3261ce39c63SPeter Brune } 3271ce39c63SPeter Brune if (col != i) { 3281ce39c63SPeter Brune if (PetscRealPart(rval[j]) > 0.) { 3291ce39c63SPeter Brune a_pos += rval[j]; 3301ce39c63SPeter Brune } else { 3311ce39c63SPeter Brune a_neg += rval[j]; 3321ce39c63SPeter Brune } 3331ce39c63SPeter Brune } else { 3341ce39c63SPeter Brune diag = rval[j]; 3351ce39c63SPeter Brune } 3368e6d0c30SPeter Brune } 3379566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA, i, &ncols, &rcol, &rval)); 3388e6d0c30SPeter Brune 3391ce39c63SPeter Brune /* ghosted connections */ 3408e6d0c30SPeter Brune if (gA) { 3419566063dSJacob Faibussowitsch PetscCall(MatGetRow(gA, i, &ncols, &rcol, &rval)); 3421ce39c63SPeter Brune for (j = 0; j < ncols; j++) { 3431ce39c63SPeter Brune col = rcol[j]; 344c1eae691SMark Adams if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0] * Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0] * Amax_neg[i])) { 3451ce39c63SPeter Brune if (PetscRealPart(rval[j]) > 0.) { 3461ce39c63SPeter Brune g_pos += rval[j]; 3478e6d0c30SPeter Brune } else { 3481ce39c63SPeter Brune g_neg += rval[j]; 3498e6d0c30SPeter Brune } 3501ce39c63SPeter Brune } 3511ce39c63SPeter Brune if (PetscRealPart(rval[j]) > 0.) { 3521ce39c63SPeter Brune a_pos += rval[j]; 3531ce39c63SPeter Brune } else { 3541ce39c63SPeter Brune a_neg += rval[j]; 3551ce39c63SPeter Brune } 3568e6d0c30SPeter Brune } 3579566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(gA, i, &ncols, &rcol, &rval)); 3588e6d0c30SPeter Brune } 3598e6d0c30SPeter Brune 3608e6d0c30SPeter Brune if (g_neg == 0.) { 3618e6d0c30SPeter Brune alpha = 0.; 3628e6d0c30SPeter Brune } else { 3638e6d0c30SPeter Brune alpha = -a_neg / g_neg; 3648e6d0c30SPeter Brune } 3658e6d0c30SPeter Brune 3668e6d0c30SPeter Brune if (g_pos == 0.) { 3678e6d0c30SPeter Brune diag += a_pos; 3688e6d0c30SPeter Brune beta = 0.; 3698e6d0c30SPeter Brune } else { 3708e6d0c30SPeter Brune beta = -a_pos / g_pos; 3718e6d0c30SPeter Brune } 372e5a0faa4SPeter Brune if (diag == 0.) { 373e5a0faa4SPeter Brune invdiag = 0.; 374e5a0faa4SPeter Brune } else invdiag = 1. / diag; 3758e6d0c30SPeter Brune /* on */ 3769566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA, i, &ncols, &rcol, &rval)); 3778e6d0c30SPeter Brune idx = 0; 3788e6d0c30SPeter Brune for (j = 0; j < ncols; j++) { 3791ce39c63SPeter Brune col = rcol[j]; 380c1eae691SMark Adams if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0] * Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0] * Amax_neg[i])) { 3818e6d0c30SPeter Brune row_f = i + fs; 3828e6d0c30SPeter Brune row_c = lcid[col]; 3838e6d0c30SPeter Brune /* set the values for on-processor ones */ 3841ce39c63SPeter Brune if (PetscRealPart(rval[j]) < 0.) { 3851ce39c63SPeter Brune pij = rval[j] * alpha * invdiag; 3868e6d0c30SPeter Brune } else { 3871ce39c63SPeter Brune pij = rval[j] * beta * invdiag; 3888e6d0c30SPeter Brune } 3898e6d0c30SPeter Brune if (PetscAbsScalar(pij) != 0.) { 3908e6d0c30SPeter Brune pvals[idx] = pij; 3918e6d0c30SPeter Brune pcols[idx] = row_c; 3928e6d0c30SPeter Brune idx++; 3938e6d0c30SPeter Brune } 3948e6d0c30SPeter Brune } 3958e6d0c30SPeter Brune } 3969566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA, i, &ncols, &rcol, &rval)); 3978e6d0c30SPeter Brune /* off */ 3981ce39c63SPeter Brune if (gA) { 3999566063dSJacob Faibussowitsch PetscCall(MatGetRow(gA, i, &ncols, &rcol, &rval)); 4008e6d0c30SPeter Brune for (j = 0; j < ncols; j++) { 4011ce39c63SPeter Brune col = rcol[j]; 402c1eae691SMark Adams if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0] * Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0] * Amax_neg[i])) { 4038e6d0c30SPeter Brune row_f = i + fs; 4048e6d0c30SPeter Brune row_c = gcid[col]; 4058e6d0c30SPeter Brune /* set the values for on-processor ones */ 4061ce39c63SPeter Brune if (PetscRealPart(rval[j]) < 0.) { 4071ce39c63SPeter Brune pij = rval[j] * alpha * invdiag; 4088e6d0c30SPeter Brune } else { 4091ce39c63SPeter Brune pij = rval[j] * beta * invdiag; 4108e6d0c30SPeter Brune } 4118e6d0c30SPeter Brune if (PetscAbsScalar(pij) != 0.) { 4128e6d0c30SPeter Brune pvals[idx] = pij; 4138e6d0c30SPeter Brune pcols[idx] = row_c; 4148e6d0c30SPeter Brune idx++; 4158e6d0c30SPeter Brune } 4168e6d0c30SPeter Brune } 4178e6d0c30SPeter Brune } 4189566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(gA, i, &ncols, &rcol, &rval)); 4193c9ab2c3SPeter Brune } 4209566063dSJacob Faibussowitsch PetscCall(MatSetValues(*P, 1, &row_f, idx, pcols, pvals, INSERT_VALUES)); 4218e6d0c30SPeter Brune } 4228e6d0c30SPeter Brune } 4233c9ab2c3SPeter Brune 4249566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY)); 4259566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY)); 4268e6d0c30SPeter Brune 4279566063dSJacob Faibussowitsch PetscCall(PetscFree5(lsparse, gsparse, lcid, Amax_pos, Amax_neg)); 428e632b94dSBarry Smith 4299566063dSJacob Faibussowitsch PetscCall(PetscFree2(pcols, pvals)); 43087b9b13cSPeter Brune if (gA) { 4319566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sf)); 4329566063dSJacob Faibussowitsch PetscCall(PetscFree(gcid)); 43387b9b13cSPeter Brune } 4343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4358e6d0c30SPeter Brune } 4368e6d0c30SPeter Brune 437ff6a9541SJacob Faibussowitsch static PetscErrorCode PCGAMGTruncateProlongator_Private(PC pc, Mat *P) 438d71ae5a4SJacob Faibussowitsch { 439586a8384SPeter Brune PetscInt j, i, ps, pf, pn, pcs, pcf, pcn, idx, cmax; 440586a8384SPeter Brune const PetscScalar *pval; 441586a8384SPeter Brune const PetscInt *pcol; 442586a8384SPeter Brune PetscScalar *pnval; 443586a8384SPeter Brune PetscInt *pncol; 444586a8384SPeter Brune PetscInt ncols; 445586a8384SPeter Brune Mat Pnew; 446586a8384SPeter Brune PetscInt *lsparse, *gsparse; 447586a8384SPeter Brune PetscReal pmax_pos, pmax_neg, ptot_pos, ptot_neg, pthresh_pos, pthresh_neg; 448586a8384SPeter Brune PC_MG *mg = (PC_MG *)pc->data; 449586a8384SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 450586a8384SPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical *)pc_gamg->subctx; 451d9558ea9SBarry Smith MatType mtype; 452586a8384SPeter Brune 453586a8384SPeter Brune PetscFunctionBegin; 454586a8384SPeter Brune /* trim and rescale with reallocation */ 4559566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(*P, &ps, &pf)); 4569566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(*P, &pcs, &pcf)); 457586a8384SPeter Brune pn = pf - ps; 458586a8384SPeter Brune pcn = pcf - pcs; 4599566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(pn, &lsparse, pn, &gsparse)); 460586a8384SPeter Brune /* allocate */ 461586a8384SPeter Brune cmax = 0; 462586a8384SPeter Brune for (i = ps; i < pf; i++) { 463b4dc3ebdSPeter Brune lsparse[i - ps] = 0; 464b4dc3ebdSPeter Brune gsparse[i - ps] = 0; 4659566063dSJacob Faibussowitsch PetscCall(MatGetRow(*P, i, &ncols, &pcol, &pval)); 466ad540459SPierre Jolivet if (ncols > cmax) cmax = ncols; 467586a8384SPeter Brune pmax_pos = 0.; 468586a8384SPeter Brune pmax_neg = 0.; 469586a8384SPeter Brune for (j = 0; j < ncols; j++) { 470586a8384SPeter Brune if (PetscRealPart(pval[j]) > pmax_pos) { 471586a8384SPeter Brune pmax_pos = PetscRealPart(pval[j]); 472586a8384SPeter Brune } else if (PetscRealPart(pval[j]) < pmax_neg) { 473586a8384SPeter Brune pmax_neg = PetscRealPart(pval[j]); 474586a8384SPeter Brune } 475586a8384SPeter Brune } 476586a8384SPeter Brune for (j = 0; j < ncols; j++) { 4778eab0cc1SPeter Brune if (PetscRealPart(pval[j]) >= pmax_pos * cls->interp_threshold || PetscRealPart(pval[j]) <= pmax_neg * cls->interp_threshold) { 478b4dc3ebdSPeter Brune if (pcol[j] >= pcs && pcol[j] < pcf) { 479b4dc3ebdSPeter Brune lsparse[i - ps]++; 480586a8384SPeter Brune } else { 481b4dc3ebdSPeter Brune gsparse[i - ps]++; 482586a8384SPeter Brune } 483586a8384SPeter Brune } 484586a8384SPeter Brune } 4859566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(*P, i, &ncols, &pcol, &pval)); 486586a8384SPeter Brune } 487586a8384SPeter Brune 4889566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(cmax, &pnval, cmax, &pncol)); 489586a8384SPeter Brune 4909566063dSJacob Faibussowitsch PetscCall(MatGetType(*P, &mtype)); 4919566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)*P), &Pnew)); 4929566063dSJacob Faibussowitsch PetscCall(MatSetType(Pnew, mtype)); 4939566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Pnew, pn, pcn, PETSC_DETERMINE, PETSC_DETERMINE)); 4949566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(Pnew, 0, lsparse)); 4959566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(Pnew, 0, lsparse, 0, gsparse)); 496586a8384SPeter Brune 497586a8384SPeter Brune for (i = ps; i < pf; i++) { 4989566063dSJacob Faibussowitsch PetscCall(MatGetRow(*P, i, &ncols, &pcol, &pval)); 499586a8384SPeter Brune pmax_pos = 0.; 500586a8384SPeter Brune pmax_neg = 0.; 501586a8384SPeter Brune for (j = 0; j < ncols; j++) { 502586a8384SPeter Brune if (PetscRealPart(pval[j]) > pmax_pos) { 503586a8384SPeter Brune pmax_pos = PetscRealPart(pval[j]); 504586a8384SPeter Brune } else if (PetscRealPart(pval[j]) < pmax_neg) { 505586a8384SPeter Brune pmax_neg = PetscRealPart(pval[j]); 506586a8384SPeter Brune } 507586a8384SPeter Brune } 508586a8384SPeter Brune pthresh_pos = 0.; 509586a8384SPeter Brune pthresh_neg = 0.; 510586a8384SPeter Brune ptot_pos = 0.; 511586a8384SPeter Brune ptot_neg = 0.; 512586a8384SPeter Brune for (j = 0; j < ncols; j++) { 5138eab0cc1SPeter Brune if (PetscRealPart(pval[j]) >= cls->interp_threshold * pmax_pos) { 514586a8384SPeter Brune pthresh_pos += PetscRealPart(pval[j]); 5158eab0cc1SPeter Brune } else if (PetscRealPart(pval[j]) <= cls->interp_threshold * pmax_neg) { 516586a8384SPeter Brune pthresh_neg += PetscRealPart(pval[j]); 517586a8384SPeter Brune } 518586a8384SPeter Brune if (PetscRealPart(pval[j]) > 0.) { 519586a8384SPeter Brune ptot_pos += PetscRealPart(pval[j]); 520586a8384SPeter Brune } else { 521586a8384SPeter Brune ptot_neg += PetscRealPart(pval[j]); 522586a8384SPeter Brune } 523586a8384SPeter Brune } 5246bd8ea90SPeter Brune if (PetscAbsReal(pthresh_pos) > 0.) ptot_pos /= pthresh_pos; 5256bd8ea90SPeter Brune if (PetscAbsReal(pthresh_neg) > 0.) ptot_neg /= pthresh_neg; 526586a8384SPeter Brune idx = 0; 527586a8384SPeter Brune for (j = 0; j < ncols; j++) { 5288eab0cc1SPeter Brune if (PetscRealPart(pval[j]) >= pmax_pos * cls->interp_threshold) { 529586a8384SPeter Brune pnval[idx] = ptot_pos * pval[j]; 530586a8384SPeter Brune pncol[idx] = pcol[j]; 531586a8384SPeter Brune idx++; 5328eab0cc1SPeter Brune } else if (PetscRealPart(pval[j]) <= pmax_neg * cls->interp_threshold) { 533586a8384SPeter Brune pnval[idx] = ptot_neg * pval[j]; 534586a8384SPeter Brune pncol[idx] = pcol[j]; 535586a8384SPeter Brune idx++; 536586a8384SPeter Brune } 537586a8384SPeter Brune } 5389566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(*P, i, &ncols, &pcol, &pval)); 5399566063dSJacob Faibussowitsch PetscCall(MatSetValues(Pnew, 1, &i, idx, pncol, pnval, INSERT_VALUES)); 540586a8384SPeter Brune } 541586a8384SPeter Brune 5429566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Pnew, MAT_FINAL_ASSEMBLY)); 5439566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Pnew, MAT_FINAL_ASSEMBLY)); 5449566063dSJacob Faibussowitsch PetscCall(MatDestroy(P)); 545586a8384SPeter Brune 546586a8384SPeter Brune *P = Pnew; 5479566063dSJacob Faibussowitsch PetscCall(PetscFree2(lsparse, gsparse)); 5489566063dSJacob Faibussowitsch PetscCall(PetscFree2(pnval, pncol)); 5493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 550586a8384SPeter Brune } 551586a8384SPeter Brune 552ff6a9541SJacob Faibussowitsch static PetscErrorCode PCGAMGProlongator_Classical_Standard(PC pc, Mat A, Mat G, PetscCoarsenData *agg_lists, Mat *P) 553d71ae5a4SJacob Faibussowitsch { 554c634539dSPeter Brune Mat lA, *lAs; 555f9a65ec8SPeter Brune MatType mtype; 55663b0c588SPeter Brune Vec cv; 55763b0c588SPeter Brune PetscInt *gcid, *lcid, *lsparse, *gsparse, *picol; 558420cd23fSPeter Brune PetscInt fs, fe, cs, ce, nl, i, j, k, li, lni, ci, ncols, maxcols, fn, cn, cid; 559420cd23fSPeter Brune PetscMPIInt size; 56063b0c588SPeter Brune const PetscInt *lidx, *icol, *gidx; 56163b0c588SPeter Brune PetscBool iscoarse; 56263b0c588SPeter Brune PetscScalar vi, pentry, pjentry; 56363b0c588SPeter Brune PetscScalar *pcontrib, *pvcol; 56463b0c588SPeter Brune const PetscScalar *vcol; 565ed4e82eaSPeter Brune PetscReal diag, jdiag, jwttotal; 566f9a65ec8SPeter Brune PetscInt pncols; 567c634539dSPeter Brune PetscSF sf; 568c634539dSPeter Brune PetscLayout clayout; 56963b0c588SPeter Brune IS lis; 570f9a65ec8SPeter Brune 571f9a65ec8SPeter Brune PetscFunctionBegin; 5729566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size)); 5739566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &fs, &fe)); 574f9a65ec8SPeter Brune fn = fe - fs; 5759566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, fe - fs, fs, 1, &lis)); 576c634539dSPeter Brune if (size > 1) { 5779566063dSJacob Faibussowitsch PetscCall(MatGetLayouts(A, NULL, &clayout)); 578f9a65ec8SPeter Brune /* increase the overlap by two to get neighbors of neighbors */ 5799566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(A, 1, &lis, 2)); 5809566063dSJacob Faibussowitsch PetscCall(ISSort(lis)); 581f9a65ec8SPeter Brune /* get the local part of A */ 5829566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(A, 1, &lis, &lis, MAT_INITIAL_MATRIX, &lAs)); 583c634539dSPeter Brune lA = lAs[0]; 584c634539dSPeter Brune /* build an SF out of it */ 5859566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(lis, &nl)); 5869566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)A), &sf)); 5879566063dSJacob Faibussowitsch PetscCall(ISGetIndices(lis, &lidx)); 5889566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraphLayout(sf, clayout, nl, NULL, PETSC_COPY_VALUES, lidx)); 5899566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(lis, &lidx)); 590c634539dSPeter Brune } else { 591c634539dSPeter Brune lA = A; 592c634539dSPeter Brune nl = fn; 593c634539dSPeter Brune } 594c634539dSPeter Brune /* create a communication structure for the overlapped portion and transmit coarse indices */ 5959566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(fn, &lsparse, fn, &gsparse, nl, &pcontrib)); 596f9a65ec8SPeter Brune /* create coarse vector */ 597f9a65ec8SPeter Brune cn = 0; 598f9a65ec8SPeter Brune for (i = 0; i < fn; i++) { 5999566063dSJacob Faibussowitsch PetscCall(PetscCDEmptyAt(agg_lists, i, &iscoarse)); 600ad540459SPierre Jolivet if (!iscoarse) cn++; 601f9a65ec8SPeter Brune } 6029566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(fn, &gcid)); 6039566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)A), cn, PETSC_DECIDE, &cv)); 6049566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(cv, &cs, &ce)); 605f9a65ec8SPeter Brune cn = 0; 606f9a65ec8SPeter Brune for (i = 0; i < fn; i++) { 6079566063dSJacob Faibussowitsch PetscCall(PetscCDEmptyAt(agg_lists, i, &iscoarse)); 608f9a65ec8SPeter Brune if (!iscoarse) { 609c634539dSPeter Brune gcid[i] = cs + cn; 610f9a65ec8SPeter Brune cn++; 611f9a65ec8SPeter Brune } else { 612c634539dSPeter Brune gcid[i] = -1; 613f9a65ec8SPeter Brune } 614f9a65ec8SPeter Brune } 615c634539dSPeter Brune if (size > 1) { 6169566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nl, &lcid)); 6179566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sf, MPIU_INT, gcid, lcid, MPI_REPLACE)); 6189566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sf, MPIU_INT, gcid, lcid, MPI_REPLACE)); 619c634539dSPeter Brune } else { 620c634539dSPeter Brune lcid = gcid; 621c634539dSPeter Brune } 622f9a65ec8SPeter Brune /* count to preallocate the prolongator */ 6239566063dSJacob Faibussowitsch PetscCall(ISGetIndices(lis, &gidx)); 624f9a65ec8SPeter Brune maxcols = 0; 625f9a65ec8SPeter Brune /* count the number of unique contributing coarse cells for each fine */ 626f9a65ec8SPeter Brune for (i = 0; i < nl; i++) { 627ed4e82eaSPeter Brune pcontrib[i] = 0.; 6289566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA, i, &ncols, &icol, NULL)); 629f9a65ec8SPeter Brune if (gidx[i] >= fs && gidx[i] < fe) { 630f9a65ec8SPeter Brune li = gidx[i] - fs; 631f9a65ec8SPeter Brune lsparse[li] = 0; 632f9a65ec8SPeter Brune gsparse[li] = 0; 633c634539dSPeter Brune cid = lcid[i]; 634f9a65ec8SPeter Brune if (cid >= 0) { 635f9a65ec8SPeter Brune lsparse[li] = 1; 636f9a65ec8SPeter Brune } else { 637f9a65ec8SPeter Brune for (j = 0; j < ncols; j++) { 638c634539dSPeter Brune if (lcid[icol[j]] >= 0) { 639f9a65ec8SPeter Brune pcontrib[icol[j]] = 1.; 640f9a65ec8SPeter Brune } else { 641f9a65ec8SPeter Brune ci = icol[j]; 6429566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA, i, &ncols, &icol, NULL)); 6439566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA, ci, &ncols, &icol, NULL)); 644f9a65ec8SPeter Brune for (k = 0; k < ncols; k++) { 645ad540459SPierre Jolivet if (lcid[icol[k]] >= 0) pcontrib[icol[k]] = 1.; 646f9a65ec8SPeter Brune } 6479566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA, ci, &ncols, &icol, NULL)); 6489566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA, i, &ncols, &icol, NULL)); 649f9a65ec8SPeter Brune } 650f9a65ec8SPeter Brune } 651f9a65ec8SPeter Brune for (j = 0; j < ncols; j++) { 652c634539dSPeter Brune if (lcid[icol[j]] >= 0 && pcontrib[icol[j]] != 0.) { 653c634539dSPeter Brune lni = lcid[icol[j]]; 654f9a65ec8SPeter Brune if (lni >= cs && lni < ce) { 655f9a65ec8SPeter Brune lsparse[li]++; 656f9a65ec8SPeter Brune } else { 657f9a65ec8SPeter Brune gsparse[li]++; 658f9a65ec8SPeter Brune } 659f9a65ec8SPeter Brune pcontrib[icol[j]] = 0.; 660f9a65ec8SPeter Brune } else { 661f9a65ec8SPeter Brune ci = icol[j]; 6629566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA, i, &ncols, &icol, NULL)); 6639566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA, ci, &ncols, &icol, NULL)); 664f9a65ec8SPeter Brune for (k = 0; k < ncols; k++) { 665c634539dSPeter Brune if (lcid[icol[k]] >= 0 && pcontrib[icol[k]] != 0.) { 666c634539dSPeter Brune lni = lcid[icol[k]]; 667f9a65ec8SPeter Brune if (lni >= cs && lni < ce) { 668f9a65ec8SPeter Brune lsparse[li]++; 669f9a65ec8SPeter Brune } else { 670f9a65ec8SPeter Brune gsparse[li]++; 671f9a65ec8SPeter Brune } 672f9a65ec8SPeter Brune pcontrib[icol[k]] = 0.; 673f9a65ec8SPeter Brune } 674f9a65ec8SPeter Brune } 6759566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA, ci, &ncols, &icol, NULL)); 6769566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA, i, &ncols, &icol, NULL)); 677f9a65ec8SPeter Brune } 678f9a65ec8SPeter Brune } 679ed4e82eaSPeter Brune } 680f9a65ec8SPeter Brune if (lsparse[li] + gsparse[li] > maxcols) maxcols = lsparse[li] + gsparse[li]; 681ed4e82eaSPeter Brune } 6829566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA, i, &ncols, &icol, &vcol)); 683f9a65ec8SPeter Brune } 6849566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(maxcols, &picol, maxcols, &pvcol)); 6859566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), P)); 6869566063dSJacob Faibussowitsch PetscCall(MatGetType(A, &mtype)); 6879566063dSJacob Faibussowitsch PetscCall(MatSetType(*P, mtype)); 6889566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*P, fn, cn, PETSC_DETERMINE, PETSC_DETERMINE)); 6899566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(*P, 0, lsparse, 0, gsparse)); 6909566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(*P, 0, lsparse)); 691f9a65ec8SPeter Brune for (i = 0; i < nl; i++) { 692ed4e82eaSPeter Brune diag = 0.; 693f9a65ec8SPeter Brune if (gidx[i] >= fs && gidx[i] < fe) { 694f9a65ec8SPeter Brune pncols = 0; 695c634539dSPeter Brune cid = lcid[i]; 696f9a65ec8SPeter Brune if (cid >= 0) { 697f9a65ec8SPeter Brune pncols = 1; 698f9a65ec8SPeter Brune picol[0] = cid; 699f9a65ec8SPeter Brune pvcol[0] = 1.; 700f9a65ec8SPeter Brune } else { 7019566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA, i, &ncols, &icol, &vcol)); 702f9a65ec8SPeter Brune for (j = 0; j < ncols; j++) { 703ed4e82eaSPeter Brune pentry = vcol[j]; 704c634539dSPeter Brune if (lcid[icol[j]] >= 0) { 705f9a65ec8SPeter Brune /* coarse neighbor */ 706ed4e82eaSPeter Brune pcontrib[icol[j]] += pentry; 707ed4e82eaSPeter Brune } else if (icol[j] != i) { 708f9a65ec8SPeter Brune /* the neighbor is a strongly connected fine node */ 709f9a65ec8SPeter Brune ci = icol[j]; 710f9a65ec8SPeter Brune vi = vcol[j]; 7119566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA, i, &ncols, &icol, &vcol)); 7129566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA, ci, &ncols, &icol, &vcol)); 713ed4e82eaSPeter Brune jwttotal = 0.; 714f9a65ec8SPeter Brune jdiag = 0.; 715f9a65ec8SPeter Brune for (k = 0; k < ncols; k++) { 716ad540459SPierre Jolivet if (ci == icol[k]) jdiag = PetscRealPart(vcol[k]); 717f9a65ec8SPeter Brune } 718f9a65ec8SPeter Brune for (k = 0; k < ncols; k++) { 719c634539dSPeter Brune if (lcid[icol[k]] >= 0 && jdiag * PetscRealPart(vcol[k]) < 0.) { 720ed4e82eaSPeter Brune pjentry = vcol[k]; 7217779008dSPeter Brune jwttotal += PetscRealPart(pjentry); 722f9a65ec8SPeter Brune } 723f9a65ec8SPeter Brune } 724ed4e82eaSPeter Brune if (jwttotal != 0.) { 7257779008dSPeter Brune jwttotal = PetscRealPart(vi) / jwttotal; 726ed4e82eaSPeter Brune for (k = 0; k < ncols; k++) { 727c634539dSPeter Brune if (lcid[icol[k]] >= 0 && jdiag * PetscRealPart(vcol[k]) < 0.) { 728586a8384SPeter Brune pjentry = vcol[k] * jwttotal; 729ed4e82eaSPeter Brune pcontrib[icol[k]] += pjentry; 730ed4e82eaSPeter Brune } 731ed4e82eaSPeter Brune } 732ed4e82eaSPeter Brune } else { 733ed4e82eaSPeter Brune diag += PetscRealPart(vi); 734ed4e82eaSPeter Brune } 7359566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA, ci, &ncols, &icol, &vcol)); 7369566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA, i, &ncols, &icol, &vcol)); 737ed4e82eaSPeter Brune } else { 738ed4e82eaSPeter Brune diag += PetscRealPart(vcol[j]); 739f9a65ec8SPeter Brune } 740f9a65ec8SPeter Brune } 741586a8384SPeter Brune if (diag != 0.) { 742586a8384SPeter Brune diag = 1. / diag; 743f9a65ec8SPeter Brune for (j = 0; j < ncols; j++) { 744c634539dSPeter Brune if (lcid[icol[j]] >= 0 && pcontrib[icol[j]] != 0.) { 745f9a65ec8SPeter Brune /* the neighbor is a coarse node */ 746ed4e82eaSPeter Brune if (PetscAbsScalar(pcontrib[icol[j]]) > 0.0) { 747c634539dSPeter Brune lni = lcid[icol[j]]; 748586a8384SPeter Brune pvcol[pncols] = -pcontrib[icol[j]] * diag; 749f9a65ec8SPeter Brune picol[pncols] = lni; 750f9a65ec8SPeter Brune pncols++; 751ed4e82eaSPeter Brune } 752ed4e82eaSPeter Brune pcontrib[icol[j]] = 0.; 753f9a65ec8SPeter Brune } else { 754f9a65ec8SPeter Brune /* the neighbor is a strongly connected fine node */ 755f9a65ec8SPeter Brune ci = icol[j]; 7569566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA, i, &ncols, &icol, &vcol)); 7579566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA, ci, &ncols, &icol, &vcol)); 758f9a65ec8SPeter Brune for (k = 0; k < ncols; k++) { 759c634539dSPeter Brune if (lcid[icol[k]] >= 0 && pcontrib[icol[k]] != 0.) { 760ed4e82eaSPeter Brune if (PetscAbsScalar(pcontrib[icol[k]]) > 0.0) { 761c634539dSPeter Brune lni = lcid[icol[k]]; 762586a8384SPeter Brune pvcol[pncols] = -pcontrib[icol[k]] * diag; 763f9a65ec8SPeter Brune picol[pncols] = lni; 764f9a65ec8SPeter Brune pncols++; 765f9a65ec8SPeter Brune } 766ed4e82eaSPeter Brune pcontrib[icol[k]] = 0.; 767ed4e82eaSPeter Brune } 768f9a65ec8SPeter Brune } 7699566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA, ci, &ncols, &icol, &vcol)); 7709566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA, i, &ncols, &icol, &vcol)); 771f9a65ec8SPeter Brune } 772ed4e82eaSPeter Brune pcontrib[icol[j]] = 0.; 773f9a65ec8SPeter Brune } 7749566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA, i, &ncols, &icol, &vcol)); 775f9a65ec8SPeter Brune } 776586a8384SPeter Brune } 777f9a65ec8SPeter Brune ci = gidx[i]; 77848a46eb9SPierre Jolivet if (pncols > 0) PetscCall(MatSetValues(*P, 1, &ci, pncols, picol, pvcol, INSERT_VALUES)); 779f9a65ec8SPeter Brune } 780f9a65ec8SPeter Brune } 7819566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(lis, &gidx)); 7829566063dSJacob Faibussowitsch PetscCall(PetscFree2(picol, pvcol)); 7839566063dSJacob Faibussowitsch PetscCall(PetscFree3(lsparse, gsparse, pcontrib)); 7849566063dSJacob Faibussowitsch PetscCall(ISDestroy(&lis)); 7859566063dSJacob Faibussowitsch PetscCall(PetscFree(gcid)); 786c634539dSPeter Brune if (size > 1) { 7879566063dSJacob Faibussowitsch PetscCall(PetscFree(lcid)); 7889566063dSJacob Faibussowitsch PetscCall(MatDestroyMatrices(1, &lAs)); 7899566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sf)); 790c634539dSPeter Brune } 7919566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cv)); 7929566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY)); 7939566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY)); 7943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7958eab0cc1SPeter Brune } 796f9a65ec8SPeter Brune 797ff6a9541SJacob Faibussowitsch static PetscErrorCode PCGAMGOptProlongator_Classical_Jacobi(PC pc, Mat A, Mat *P) 798d71ae5a4SJacob Faibussowitsch { 799b4dc3ebdSPeter Brune PetscInt f, s, n, cf, cs, i, idx; 800b4dc3ebdSPeter Brune PetscInt *coarserows; 801b4dc3ebdSPeter Brune PetscInt ncols; 802b4dc3ebdSPeter Brune const PetscInt *pcols; 803b4dc3ebdSPeter Brune const PetscScalar *pvals; 804b4dc3ebdSPeter Brune Mat Pnew; 805b4dc3ebdSPeter Brune Vec diag; 806b4dc3ebdSPeter Brune PC_MG *mg = (PC_MG *)pc->data; 807b4dc3ebdSPeter Brune PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 808b4dc3ebdSPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical *)pc_gamg->subctx; 809b4dc3ebdSPeter Brune 810b4dc3ebdSPeter Brune PetscFunctionBegin; 811b4dc3ebdSPeter Brune if (cls->nsmooths == 0) { 8129566063dSJacob Faibussowitsch PetscCall(PCGAMGTruncateProlongator_Private(pc, P)); 8133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 814b4dc3ebdSPeter Brune } 8159566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(*P, &s, &f)); 816b4dc3ebdSPeter Brune n = f - s; 8179566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(*P, &cs, &cf)); 8189566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &coarserows)); 819b4dc3ebdSPeter Brune /* identify the rows corresponding to coarse unknowns */ 820b4dc3ebdSPeter Brune idx = 0; 821b4dc3ebdSPeter Brune for (i = s; i < f; i++) { 8229566063dSJacob Faibussowitsch PetscCall(MatGetRow(*P, i, &ncols, &pcols, &pvals)); 823b4dc3ebdSPeter Brune /* assume, for now, that it's a coarse unknown if it has a single unit entry */ 824b4dc3ebdSPeter Brune if (ncols == 1) { 825b4dc3ebdSPeter Brune if (pvals[0] == 1.) { 826b4dc3ebdSPeter Brune coarserows[idx] = i; 827b4dc3ebdSPeter Brune idx++; 828b4dc3ebdSPeter Brune } 829b4dc3ebdSPeter Brune } 8309566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(*P, i, &ncols, &pcols, &pvals)); 831b4dc3ebdSPeter Brune } 8329566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &diag, NULL)); 8339566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(A, diag)); 8349566063dSJacob Faibussowitsch PetscCall(VecReciprocal(diag)); 835b4dc3ebdSPeter Brune for (i = 0; i < cls->nsmooths; i++) { 8369566063dSJacob Faibussowitsch PetscCall(MatMatMult(A, *P, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Pnew)); 8379566063dSJacob Faibussowitsch PetscCall(MatZeroRows(Pnew, idx, coarserows, 0., NULL, NULL)); 8389566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(Pnew, diag, NULL)); 8399566063dSJacob Faibussowitsch PetscCall(MatAYPX(Pnew, -1.0, *P, DIFFERENT_NONZERO_PATTERN)); 8409566063dSJacob Faibussowitsch PetscCall(MatDestroy(P)); 841b4dc3ebdSPeter Brune *P = Pnew; 842b4dc3ebdSPeter Brune Pnew = NULL; 843b4dc3ebdSPeter Brune } 8449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&diag)); 8459566063dSJacob Faibussowitsch PetscCall(PetscFree(coarserows)); 8469566063dSJacob Faibussowitsch PetscCall(PCGAMGTruncateProlongator_Private(pc, P)); 8473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 848b4dc3ebdSPeter Brune } 849b4dc3ebdSPeter Brune 850d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGProlongator_Classical(PC pc, Mat A, Mat G, PetscCoarsenData *agg_lists, Mat *P) 851d71ae5a4SJacob Faibussowitsch { 8528eab0cc1SPeter Brune PetscErrorCode (*f)(PC, Mat, Mat, PetscCoarsenData *, Mat *); 8538eab0cc1SPeter Brune PC_MG *mg = (PC_MG *)pc->data; 8548eab0cc1SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 8558eab0cc1SPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical *)pc_gamg->subctx; 8568eab0cc1SPeter Brune 8578eab0cc1SPeter Brune PetscFunctionBegin; 8589566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(PCGAMGClassicalProlongatorList, cls->prolongtype, &f)); 85928b400f6SJacob Faibussowitsch PetscCheck(f, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Cannot find PCGAMG Classical prolongator type"); 8609566063dSJacob Faibussowitsch PetscCall((*f)(pc, A, G, agg_lists, P)); 8613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 862f9a65ec8SPeter Brune } 863f9a65ec8SPeter Brune 864d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGDestroy_Classical(PC pc) 865d71ae5a4SJacob Faibussowitsch { 8668e6d0c30SPeter Brune PC_MG *mg = (PC_MG *)pc->data; 8678e6d0c30SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 8688e6d0c30SPeter Brune 8698e6d0c30SPeter Brune PetscFunctionBegin; 8709566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->subctx)); 8719566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGClassicalSetType_C", NULL)); 8729566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGClassicalGetType_C", NULL)); 8733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8748e6d0c30SPeter Brune } 8758e6d0c30SPeter Brune 876ff6a9541SJacob Faibussowitsch static PetscErrorCode PCGAMGSetFromOptions_Classical(PC pc, PetscOptionItems *PetscOptionsObject) 877d71ae5a4SJacob Faibussowitsch { 878586a8384SPeter Brune PC_MG *mg = (PC_MG *)pc->data; 879586a8384SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 880586a8384SPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical *)pc_gamg->subctx; 8818eab0cc1SPeter Brune char tname[256]; 8828eab0cc1SPeter Brune PetscBool flg; 883586a8384SPeter Brune 8848e6d0c30SPeter Brune PetscFunctionBegin; 885d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "GAMG-Classical options"); 8869566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-pc_gamg_classical_type", "Type of Classical AMG prolongation", "PCGAMGClassicalSetType", PCGAMGClassicalProlongatorList, cls->prolongtype, tname, sizeof(tname), &flg)); 8871baa6e33SBarry Smith if (flg) PetscCall(PCGAMGClassicalSetType(pc, tname)); 8889566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_gamg_classical_interp_threshold", "Threshold for classical interpolator entries", "", cls->interp_threshold, &cls->interp_threshold, NULL)); 8899566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_gamg_classical_nsmooths", "Threshold for classical interpolator entries", "", cls->nsmooths, &cls->nsmooths, NULL)); 890d0609cedSBarry Smith PetscOptionsHeadEnd(); 8913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8928e6d0c30SPeter Brune } 8938e6d0c30SPeter Brune 894d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetData_Classical(PC pc, Mat A) 895d71ae5a4SJacob Faibussowitsch { 8968e6d0c30SPeter Brune PC_MG *mg = (PC_MG *)pc->data; 8978e6d0c30SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 8988e6d0c30SPeter Brune 8998e6d0c30SPeter Brune PetscFunctionBegin; 9008e6d0c30SPeter Brune /* no data for classical AMG */ 9018e6d0c30SPeter Brune pc_gamg->data = NULL; 902d2050638SMark Adams pc_gamg->data_cell_cols = 0; 903d2050638SMark Adams pc_gamg->data_cell_rows = 0; 9048e6d0c30SPeter Brune pc_gamg->data_sz = 0; 9053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9068e6d0c30SPeter Brune } 9078e6d0c30SPeter Brune 908ff6a9541SJacob Faibussowitsch static PetscErrorCode PCGAMGClassicalFinalizePackage(void) 909d71ae5a4SJacob Faibussowitsch { 9108eab0cc1SPeter Brune PetscFunctionBegin; 9118eab0cc1SPeter Brune PCGAMGClassicalPackageInitialized = PETSC_FALSE; 9129566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&PCGAMGClassicalProlongatorList)); 9133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9148eab0cc1SPeter Brune } 9158eab0cc1SPeter Brune 916ff6a9541SJacob Faibussowitsch static PetscErrorCode PCGAMGClassicalInitializePackage(void) 917d71ae5a4SJacob Faibussowitsch { 9188eab0cc1SPeter Brune PetscFunctionBegin; 9193ba16761SJacob Faibussowitsch if (PCGAMGClassicalPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS); 9209566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&PCGAMGClassicalProlongatorList, PCGAMGCLASSICALDIRECT, PCGAMGProlongator_Classical_Direct)); 9219566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&PCGAMGClassicalProlongatorList, PCGAMGCLASSICALSTANDARD, PCGAMGProlongator_Classical_Standard)); 9229566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(PCGAMGClassicalFinalizePackage)); 9233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9248eab0cc1SPeter Brune } 9258eab0cc1SPeter Brune 926d71ae5a4SJacob Faibussowitsch PetscErrorCode PCCreateGAMG_Classical(PC pc) 927d71ae5a4SJacob Faibussowitsch { 9288e6d0c30SPeter Brune PC_MG *mg = (PC_MG *)pc->data; 9298e6d0c30SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 9308e6d0c30SPeter Brune PC_GAMG_Classical *pc_gamg_classical; 9318e6d0c30SPeter Brune 9328e6d0c30SPeter Brune PetscFunctionBegin; 9339566063dSJacob Faibussowitsch PetscCall(PCGAMGClassicalInitializePackage()); 9348e6d0c30SPeter Brune if (pc_gamg->subctx) { 9358e6d0c30SPeter Brune /* call base class */ 9369566063dSJacob Faibussowitsch PetscCall(PCDestroy_GAMG(pc)); 9378e6d0c30SPeter Brune } 9388e6d0c30SPeter Brune 9398e6d0c30SPeter Brune /* create sub context for SA */ 9404dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&pc_gamg_classical)); 9418e6d0c30SPeter Brune pc_gamg->subctx = pc_gamg_classical; 9428e6d0c30SPeter Brune pc->ops->setfromoptions = PCGAMGSetFromOptions_Classical; 9438e6d0c30SPeter Brune /* reset does not do anything; setup not virtual */ 9448e6d0c30SPeter Brune 9458e6d0c30SPeter Brune /* set internal function pointers */ 9468e6d0c30SPeter Brune pc_gamg->ops->destroy = PCGAMGDestroy_Classical; 9472d776b49SBarry Smith pc_gamg->ops->creategraph = PCGAMGCreateGraph_Classical; 9488e6d0c30SPeter Brune pc_gamg->ops->coarsen = PCGAMGCoarsen_Classical; 9498eab0cc1SPeter Brune pc_gamg->ops->prolongator = PCGAMGProlongator_Classical; 950fd1112cbSBarry Smith pc_gamg->ops->optprolongator = PCGAMGOptProlongator_Classical_Jacobi; 951586a8384SPeter Brune pc_gamg->ops->setfromoptions = PCGAMGSetFromOptions_Classical; 9528e6d0c30SPeter Brune 9538e6d0c30SPeter Brune pc_gamg->ops->createdefaultdata = PCGAMGSetData_Classical; 954586a8384SPeter Brune pc_gamg_classical->interp_threshold = 0.2; 955b4dc3ebdSPeter Brune pc_gamg_classical->nsmooths = 0; 9569566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGClassicalSetType_C", PCGAMGClassicalSetType_GAMG)); 9579566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGClassicalGetType_C", PCGAMGClassicalGetType_GAMG)); 9589566063dSJacob Faibussowitsch PetscCall(PCGAMGClassicalSetType(pc, PCGAMGCLASSICALSTANDARD)); 9593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9608e6d0c30SPeter Brune } 961