18e6d0c30SPeter Brune #include <../src/ksp/pc/impls/gamg/gamg.h> /*I "petscpc.h" I*/ 287b9b13cSPeter Brune #include <petscsf.h> 38e6d0c30SPeter Brune 48eab0cc1SPeter Brune PetscFunctionList PCGAMGClassicalProlongatorList = NULL; 58eab0cc1SPeter Brune 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 14*f1580f4eSBarry Smith PCGAMGClassicalSetType - Sets the type of classical interpolation to use with `PCGAMG` 158eab0cc1SPeter Brune 16*f1580f4eSBarry Smith Collective on pc 178eab0cc1SPeter Brune 188eab0cc1SPeter Brune Input Parameters: 198eab0cc1SPeter Brune . pc - the preconditioner context 208eab0cc1SPeter Brune 218eab0cc1SPeter Brune Options Database Key: 22*f1580f4eSBarry Smith . -pc_gamg_classical_type <direct,standard> - set type of classical AMG prolongation 238eab0cc1SPeter Brune 248eab0cc1SPeter Brune Level: intermediate 258eab0cc1SPeter Brune 26*f1580f4eSBarry Smith .seealso: `PCGAMG` 278eab0cc1SPeter Brune @*/ 289371c9d4SSatish Balay PetscErrorCode PCGAMGClassicalSetType(PC pc, PCGAMGClassicalType type) { 298eab0cc1SPeter Brune PetscFunctionBegin; 308eab0cc1SPeter Brune PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 31cac4c232SBarry Smith PetscTryMethod(pc, "PCGAMGClassicalSetType_C", (PC, PCGAMGClassicalType), (pc, type)); 328eab0cc1SPeter Brune PetscFunctionReturn(0); 338eab0cc1SPeter Brune } 348eab0cc1SPeter Brune 35c60c7ad4SBarry Smith /*@C 36*f1580f4eSBarry Smith PCGAMGClassicalGetType - Gets the type of classical interpolation to use with `PCGAMG` 37c60c7ad4SBarry Smith 38*f1580f4eSBarry Smith Collective on pc 39c60c7ad4SBarry Smith 40c60c7ad4SBarry Smith Input Parameter: 41c60c7ad4SBarry Smith . pc - the preconditioner context 42c60c7ad4SBarry Smith 43c60c7ad4SBarry Smith Output Parameter: 44c60c7ad4SBarry Smith . type - the type used 45c60c7ad4SBarry Smith 46c60c7ad4SBarry Smith Level: intermediate 47c60c7ad4SBarry Smith 48*f1580f4eSBarry Smith .seealso: `PCGAMG` 49c60c7ad4SBarry Smith @*/ 509371c9d4SSatish Balay PetscErrorCode PCGAMGClassicalGetType(PC pc, PCGAMGClassicalType *type) { 51c60c7ad4SBarry Smith PetscFunctionBegin; 52c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 53cac4c232SBarry Smith PetscUseMethod(pc, "PCGAMGClassicalGetType_C", (PC, PCGAMGClassicalType *), (pc, type)); 54c60c7ad4SBarry Smith PetscFunctionReturn(0); 55c60c7ad4SBarry Smith } 56c60c7ad4SBarry Smith 579371c9d4SSatish Balay static PetscErrorCode PCGAMGClassicalSetType_GAMG(PC pc, PCGAMGClassicalType type) { 588eab0cc1SPeter Brune PC_MG *mg = (PC_MG *)pc->data; 598eab0cc1SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 608eab0cc1SPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical *)pc_gamg->subctx; 618eab0cc1SPeter Brune 628eab0cc1SPeter Brune PetscFunctionBegin; 639566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(cls->prolongtype, type)); 648eab0cc1SPeter Brune PetscFunctionReturn(0); 658eab0cc1SPeter Brune } 668e6d0c30SPeter Brune 679371c9d4SSatish Balay static PetscErrorCode PCGAMGClassicalGetType_GAMG(PC pc, PCGAMGClassicalType *type) { 68c60c7ad4SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 69c60c7ad4SBarry Smith PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 70c60c7ad4SBarry Smith PC_GAMG_Classical *cls = (PC_GAMG_Classical *)pc_gamg->subctx; 71c60c7ad4SBarry Smith 72c60c7ad4SBarry Smith PetscFunctionBegin; 73c60c7ad4SBarry Smith *type = cls->prolongtype; 74c60c7ad4SBarry Smith PetscFunctionReturn(0); 75c60c7ad4SBarry Smith } 76c60c7ad4SBarry Smith 779371c9d4SSatish Balay PetscErrorCode PCGAMGGraph_Classical(PC pc, Mat A, Mat *G) { 78550383edSPeter Brune PetscInt s, f, n, idx, lidx, gidx; 79e5a0faa4SPeter Brune PetscInt r, c, ncols; 808e6d0c30SPeter Brune const PetscInt *rcol; 818e6d0c30SPeter Brune const PetscScalar *rval; 82e5a0faa4SPeter Brune PetscInt *gcol; 838e6d0c30SPeter Brune PetscScalar *gval; 84e5a0faa4SPeter Brune PetscReal rmax; 85550383edSPeter Brune PetscInt cmax = 0; 86b817416eSBarry Smith PC_MG *mg = (PC_MG *)pc->data; 87b817416eSBarry Smith PC_GAMG *gamg = (PC_GAMG *)mg->innerctx; 888e6d0c30SPeter Brune PetscInt *gsparse, *lsparse; 89e5a0faa4SPeter Brune PetscScalar *Amax; 908e6d0c30SPeter Brune MatType mtype; 918e6d0c30SPeter Brune 928e6d0c30SPeter Brune PetscFunctionBegin; 939566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &s, &f)); 94550383edSPeter Brune n = f - s; 959566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(n, &lsparse, n, &gsparse, n, &Amax)); 968e6d0c30SPeter Brune 97550383edSPeter Brune for (r = 0; r < n; r++) { 988e6d0c30SPeter Brune lsparse[r] = 0; 99550383edSPeter Brune gsparse[r] = 0; 1008e6d0c30SPeter Brune } 1018e6d0c30SPeter Brune 102550383edSPeter Brune for (r = s; r < f; r++) { 103e5a0faa4SPeter Brune /* determine the maximum off-diagonal in each row */ 104e5a0faa4SPeter Brune rmax = 0.; 1059566063dSJacob Faibussowitsch PetscCall(MatGetRow(A, r, &ncols, &rcol, &rval)); 106e5a0faa4SPeter Brune for (c = 0; c < ncols; c++) { 107ad540459SPierre Jolivet if (PetscRealPart(-rval[c]) > rmax && rcol[c] != r) rmax = PetscRealPart(-rval[c]); 108e5a0faa4SPeter Brune } 109550383edSPeter Brune Amax[r - s] = rmax; 110550383edSPeter Brune if (ncols > cmax) cmax = ncols; 111550383edSPeter Brune lidx = 0; 112550383edSPeter Brune gidx = 0; 113e5a0faa4SPeter Brune /* create the local and global sparsity patterns */ 1148e6d0c30SPeter Brune for (c = 0; c < ncols; c++) { 115c1eae691SMark Adams if (PetscRealPart(-rval[c]) > gamg->threshold[0] * PetscRealPart(Amax[r - s]) || rcol[c] == r) { 116550383edSPeter Brune if (rcol[c] < f && rcol[c] >= s) { 117550383edSPeter Brune lidx++; 118550383edSPeter Brune } else { 119550383edSPeter Brune gidx++; 1208e6d0c30SPeter Brune } 1218e6d0c30SPeter Brune } 1228e6d0c30SPeter Brune } 1239566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A, r, &ncols, &rcol, &rval)); 124550383edSPeter Brune lsparse[r - s] = lidx; 125550383edSPeter Brune gsparse[r - s] = gidx; 1268e6d0c30SPeter Brune } 1279566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(cmax, &gval, cmax, &gcol)); 128e5a0faa4SPeter Brune 1299566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), G)); 1309566063dSJacob Faibussowitsch PetscCall(MatGetType(A, &mtype)); 1319566063dSJacob Faibussowitsch PetscCall(MatSetType(*G, mtype)); 1329566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*G, n, n, PETSC_DETERMINE, PETSC_DETERMINE)); 1339566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(*G, 0, lsparse, 0, gsparse)); 1349566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(*G, 0, lsparse)); 1358e6d0c30SPeter Brune for (r = s; r < f; r++) { 1369566063dSJacob Faibussowitsch PetscCall(MatGetRow(A, r, &ncols, &rcol, &rval)); 1378e6d0c30SPeter Brune idx = 0; 1388e6d0c30SPeter Brune for (c = 0; c < ncols; c++) { 1398e6d0c30SPeter Brune /* classical strength of connection */ 140c1eae691SMark Adams if (PetscRealPart(-rval[c]) > gamg->threshold[0] * PetscRealPart(Amax[r - s]) || rcol[c] == r) { 1418e6d0c30SPeter Brune gcol[idx] = rcol[c]; 1428e6d0c30SPeter Brune gval[idx] = rval[c]; 1438e6d0c30SPeter Brune idx++; 1448e6d0c30SPeter Brune } 1458e6d0c30SPeter Brune } 1469566063dSJacob Faibussowitsch PetscCall(MatSetValues(*G, 1, &r, idx, gcol, gval, INSERT_VALUES)); 1479566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A, r, &ncols, &rcol, &rval)); 1488e6d0c30SPeter Brune } 1499566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*G, MAT_FINAL_ASSEMBLY)); 1509566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*G, MAT_FINAL_ASSEMBLY)); 1518e6d0c30SPeter Brune 1529566063dSJacob Faibussowitsch PetscCall(PetscFree2(gval, gcol)); 1539566063dSJacob Faibussowitsch PetscCall(PetscFree3(lsparse, gsparse, Amax)); 1548e6d0c30SPeter Brune PetscFunctionReturn(0); 1558e6d0c30SPeter Brune } 1568e6d0c30SPeter Brune 1579371c9d4SSatish Balay PetscErrorCode PCGAMGCoarsen_Classical(PC pc, Mat *G, PetscCoarsenData **agg_lists) { 1588e6d0c30SPeter Brune MatCoarsen crs; 1598e6d0c30SPeter Brune MPI_Comm fcomm = ((PetscObject)pc)->comm; 1608e6d0c30SPeter Brune 1618e6d0c30SPeter Brune PetscFunctionBegin; 16228b400f6SJacob Faibussowitsch PetscCheck(G, fcomm, PETSC_ERR_ARG_WRONGSTATE, "Must set Graph in PC in PCGAMG before coarsening"); 1638e6d0c30SPeter Brune 1649566063dSJacob Faibussowitsch PetscCall(MatCoarsenCreate(fcomm, &crs)); 1659566063dSJacob Faibussowitsch PetscCall(MatCoarsenSetFromOptions(crs)); 1669566063dSJacob Faibussowitsch PetscCall(MatCoarsenSetAdjacency(crs, *G)); 1679566063dSJacob Faibussowitsch PetscCall(MatCoarsenSetStrictAggs(crs, PETSC_TRUE)); 1689566063dSJacob Faibussowitsch PetscCall(MatCoarsenApply(crs)); 1699566063dSJacob Faibussowitsch PetscCall(MatCoarsenGetData(crs, agg_lists)); 1709566063dSJacob Faibussowitsch PetscCall(MatCoarsenDestroy(&crs)); 1718e6d0c30SPeter Brune PetscFunctionReturn(0); 1728e6d0c30SPeter Brune } 1738e6d0c30SPeter Brune 1749371c9d4SSatish Balay PetscErrorCode PCGAMGProlongator_Classical_Direct(PC pc, Mat A, Mat G, PetscCoarsenData *agg_lists, Mat *P) { 1751ce39c63SPeter Brune PC_MG *mg = (PC_MG *)pc->data; 1761ce39c63SPeter Brune PC_GAMG *gamg = (PC_GAMG *)mg->innerctx; 17763b0c588SPeter Brune PetscBool iscoarse, isMPIAIJ, isSEQAIJ; 17863b0c588SPeter Brune PetscInt fn, cn, fs, fe, cs, ce, i, j, ncols, col, row_f, row_c, cmax = 0, idx, noff; 17963b0c588SPeter Brune PetscInt *lcid, *gcid, *lsparse, *gsparse, *colmap, *pcols; 18063b0c588SPeter Brune const PetscInt *rcol; 18163b0c588SPeter Brune PetscReal *Amax_pos, *Amax_neg; 18263b0c588SPeter Brune PetscScalar g_pos, g_neg, a_pos, a_neg, diag, invdiag, alpha, beta, pij; 18363b0c588SPeter Brune PetscScalar *pvals; 18463b0c588SPeter Brune const PetscScalar *rval; 18563b0c588SPeter Brune Mat lA, gA = NULL; 18663b0c588SPeter Brune MatType mtype; 18763b0c588SPeter Brune Vec C, lvec; 18887b9b13cSPeter Brune PetscLayout clayout; 18987b9b13cSPeter Brune PetscSF sf; 19087b9b13cSPeter Brune Mat_MPIAIJ *mpiaij; 1918e6d0c30SPeter Brune 1928e6d0c30SPeter Brune PetscFunctionBegin; 1939566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &fs, &fe)); 19487b9b13cSPeter Brune fn = fe - fs; 1959566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPIAIJ, &isMPIAIJ)); 1969566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQAIJ, &isSEQAIJ)); 1977827d75bSBarry Smith PetscCheck(isMPIAIJ || isSEQAIJ, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Classical AMG requires MPIAIJ matrix"); 19887b9b13cSPeter Brune if (isMPIAIJ) { 19987b9b13cSPeter Brune mpiaij = (Mat_MPIAIJ *)A->data; 20087b9b13cSPeter Brune lA = mpiaij->A; 20187b9b13cSPeter Brune gA = mpiaij->B; 20287b9b13cSPeter Brune lvec = mpiaij->lvec; 2039566063dSJacob Faibussowitsch PetscCall(VecGetSize(lvec, &noff)); 20487b9b13cSPeter Brune colmap = mpiaij->garray; 2059566063dSJacob Faibussowitsch PetscCall(MatGetLayouts(A, NULL, &clayout)); 2069566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)A), &sf)); 2079566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraphLayout(sf, clayout, noff, NULL, PETSC_COPY_VALUES, colmap)); 2089566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(noff, &gcid)); 20987b9b13cSPeter Brune } else { 21087b9b13cSPeter Brune lA = A; 21187b9b13cSPeter Brune } 2129566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(fn, &lsparse, fn, &gsparse, fn, &lcid, fn, &Amax_pos, fn, &Amax_neg)); 2138e6d0c30SPeter Brune 2148e6d0c30SPeter Brune /* count the number of coarse unknowns */ 2158e6d0c30SPeter Brune cn = 0; 2168e6d0c30SPeter Brune for (i = 0; i < fn; i++) { 2178e6d0c30SPeter Brune /* filter out singletons */ 2189566063dSJacob Faibussowitsch PetscCall(PetscCDEmptyAt(agg_lists, i, &iscoarse)); 2198e6d0c30SPeter Brune lcid[i] = -1; 220ad540459SPierre Jolivet if (!iscoarse) cn++; 2218e6d0c30SPeter Brune } 2228e6d0c30SPeter Brune 2238e6d0c30SPeter Brune /* create the coarse vector */ 2249566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)A), cn, PETSC_DECIDE, &C)); 2259566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(C, &cs, &ce)); 2268e6d0c30SPeter Brune 2278e6d0c30SPeter Brune cn = 0; 2288e6d0c30SPeter Brune for (i = 0; i < fn; i++) { 2299566063dSJacob Faibussowitsch PetscCall(PetscCDEmptyAt(agg_lists, i, &iscoarse)); 2308e6d0c30SPeter Brune if (!iscoarse) { 2318e6d0c30SPeter Brune lcid[i] = cs + cn; 2328e6d0c30SPeter Brune cn++; 2338e6d0c30SPeter Brune } else { 2348e6d0c30SPeter Brune lcid[i] = -1; 2358e6d0c30SPeter Brune } 2368e6d0c30SPeter Brune } 2378e6d0c30SPeter Brune 23887b9b13cSPeter Brune if (gA) { 2399566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sf, MPIU_INT, lcid, gcid, MPI_REPLACE)); 2409566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sf, MPIU_INT, lcid, gcid, MPI_REPLACE)); 24187b9b13cSPeter Brune } 2428e6d0c30SPeter Brune 243b817416eSBarry Smith /* determine the largest off-diagonal entries in each row */ 2441ce39c63SPeter Brune for (i = fs; i < fe; i++) { 2451ce39c63SPeter Brune Amax_pos[i - fs] = 0.; 2461ce39c63SPeter Brune Amax_neg[i - fs] = 0.; 2479566063dSJacob Faibussowitsch PetscCall(MatGetRow(A, i, &ncols, &rcol, &rval)); 2481ce39c63SPeter Brune for (j = 0; j < ncols; j++) { 2491ce39c63SPeter Brune if ((PetscRealPart(-rval[j]) > Amax_neg[i - fs]) && i != rcol[j]) Amax_neg[i - fs] = PetscAbsScalar(rval[j]); 2501ce39c63SPeter Brune if ((PetscRealPart(rval[j]) > Amax_pos[i - fs]) && i != rcol[j]) Amax_pos[i - fs] = PetscAbsScalar(rval[j]); 2511ce39c63SPeter Brune } 2521ce39c63SPeter Brune if (ncols > cmax) cmax = ncols; 2539566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A, i, &ncols, &rcol, &rval)); 2541ce39c63SPeter Brune } 2559566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(cmax, &pcols, cmax, &pvals)); 2569566063dSJacob Faibussowitsch PetscCall(VecDestroy(&C)); 2578e6d0c30SPeter Brune 2588e6d0c30SPeter Brune /* count the on and off processor sparsity patterns for the prolongator */ 2598e6d0c30SPeter Brune for (i = 0; i < fn; i++) { 2608e6d0c30SPeter Brune /* on */ 2618e6d0c30SPeter Brune lsparse[i] = 0; 262e5a0faa4SPeter Brune gsparse[i] = 0; 2638e6d0c30SPeter Brune if (lcid[i] >= 0) { 2648e6d0c30SPeter Brune lsparse[i] = 1; 2658e6d0c30SPeter Brune gsparse[i] = 0; 2668e6d0c30SPeter Brune } else { 2679566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA, i, &ncols, &rcol, &rval)); 2688e6d0c30SPeter Brune for (j = 0; j < ncols; j++) { 2691ce39c63SPeter Brune col = rcol[j]; 270ad540459SPierre 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; 2718e6d0c30SPeter Brune } 2729566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA, i, &ncols, &rcol, &rval)); 2738e6d0c30SPeter Brune /* off */ 2741ce39c63SPeter Brune if (gA) { 2759566063dSJacob Faibussowitsch PetscCall(MatGetRow(gA, i, &ncols, &rcol, &rval)); 2768e6d0c30SPeter Brune for (j = 0; j < ncols; j++) { 2771ce39c63SPeter Brune col = rcol[j]; 278ad540459SPierre 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; 2798e6d0c30SPeter Brune } 2809566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(gA, i, &ncols, &rcol, &rval)); 2818e6d0c30SPeter Brune } 2828e6d0c30SPeter Brune } 2831ce39c63SPeter Brune } 2848e6d0c30SPeter Brune 2858e6d0c30SPeter Brune /* preallocate and create the prolongator */ 2869566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), P)); 2879566063dSJacob Faibussowitsch PetscCall(MatGetType(G, &mtype)); 2889566063dSJacob Faibussowitsch PetscCall(MatSetType(*P, mtype)); 2899566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*P, fn, cn, PETSC_DETERMINE, PETSC_DETERMINE)); 2909566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(*P, 0, lsparse, 0, gsparse)); 2919566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(*P, 0, lsparse)); 2928e6d0c30SPeter Brune 2938e6d0c30SPeter Brune /* loop over local fine nodes -- get the diagonal, the sum of positive and negative strong and weak weights, and set up the row */ 2948e6d0c30SPeter Brune for (i = 0; i < fn; i++) { 2958e6d0c30SPeter Brune /* determine on or off */ 2968e6d0c30SPeter Brune row_f = i + fs; 2978e6d0c30SPeter Brune row_c = lcid[i]; 2988e6d0c30SPeter Brune if (row_c >= 0) { 2998e6d0c30SPeter Brune pij = 1.; 3009566063dSJacob Faibussowitsch PetscCall(MatSetValues(*P, 1, &row_f, 1, &row_c, &pij, INSERT_VALUES)); 3018e6d0c30SPeter Brune } else { 3028e6d0c30SPeter Brune g_pos = 0.; 3038e6d0c30SPeter Brune g_neg = 0.; 3048e6d0c30SPeter Brune a_pos = 0.; 3058e6d0c30SPeter Brune a_neg = 0.; 3068e6d0c30SPeter Brune diag = 0.; 3078e6d0c30SPeter Brune 3081ce39c63SPeter Brune /* local connections */ 3099566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA, i, &ncols, &rcol, &rval)); 3101ce39c63SPeter Brune for (j = 0; j < ncols; j++) { 3111ce39c63SPeter Brune col = rcol[j]; 312c1eae691SMark Adams if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0] * Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0] * Amax_neg[i])) { 3131ce39c63SPeter Brune if (PetscRealPart(rval[j]) > 0.) { 3141ce39c63SPeter Brune g_pos += rval[j]; 3158e6d0c30SPeter Brune } else { 3161ce39c63SPeter Brune g_neg += rval[j]; 3178e6d0c30SPeter Brune } 3181ce39c63SPeter Brune } 3191ce39c63SPeter Brune if (col != i) { 3201ce39c63SPeter Brune if (PetscRealPart(rval[j]) > 0.) { 3211ce39c63SPeter Brune a_pos += rval[j]; 3221ce39c63SPeter Brune } else { 3231ce39c63SPeter Brune a_neg += rval[j]; 3241ce39c63SPeter Brune } 3251ce39c63SPeter Brune } else { 3261ce39c63SPeter Brune diag = rval[j]; 3271ce39c63SPeter Brune } 3288e6d0c30SPeter Brune } 3299566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA, i, &ncols, &rcol, &rval)); 3308e6d0c30SPeter Brune 3311ce39c63SPeter Brune /* ghosted connections */ 3328e6d0c30SPeter Brune if (gA) { 3339566063dSJacob Faibussowitsch PetscCall(MatGetRow(gA, i, &ncols, &rcol, &rval)); 3341ce39c63SPeter Brune for (j = 0; j < ncols; j++) { 3351ce39c63SPeter Brune col = rcol[j]; 336c1eae691SMark Adams if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0] * Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0] * Amax_neg[i])) { 3371ce39c63SPeter Brune if (PetscRealPart(rval[j]) > 0.) { 3381ce39c63SPeter Brune g_pos += rval[j]; 3398e6d0c30SPeter Brune } else { 3401ce39c63SPeter Brune g_neg += rval[j]; 3418e6d0c30SPeter Brune } 3421ce39c63SPeter Brune } 3431ce39c63SPeter Brune if (PetscRealPart(rval[j]) > 0.) { 3441ce39c63SPeter Brune a_pos += rval[j]; 3451ce39c63SPeter Brune } else { 3461ce39c63SPeter Brune a_neg += rval[j]; 3471ce39c63SPeter Brune } 3488e6d0c30SPeter Brune } 3499566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(gA, i, &ncols, &rcol, &rval)); 3508e6d0c30SPeter Brune } 3518e6d0c30SPeter Brune 3528e6d0c30SPeter Brune if (g_neg == 0.) { 3538e6d0c30SPeter Brune alpha = 0.; 3548e6d0c30SPeter Brune } else { 3558e6d0c30SPeter Brune alpha = -a_neg / g_neg; 3568e6d0c30SPeter Brune } 3578e6d0c30SPeter Brune 3588e6d0c30SPeter Brune if (g_pos == 0.) { 3598e6d0c30SPeter Brune diag += a_pos; 3608e6d0c30SPeter Brune beta = 0.; 3618e6d0c30SPeter Brune } else { 3628e6d0c30SPeter Brune beta = -a_pos / g_pos; 3638e6d0c30SPeter Brune } 364e5a0faa4SPeter Brune if (diag == 0.) { 365e5a0faa4SPeter Brune invdiag = 0.; 366e5a0faa4SPeter Brune } else invdiag = 1. / diag; 3678e6d0c30SPeter Brune /* on */ 3689566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA, i, &ncols, &rcol, &rval)); 3698e6d0c30SPeter Brune idx = 0; 3708e6d0c30SPeter Brune for (j = 0; j < ncols; j++) { 3711ce39c63SPeter Brune col = rcol[j]; 372c1eae691SMark Adams if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0] * Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0] * Amax_neg[i])) { 3738e6d0c30SPeter Brune row_f = i + fs; 3748e6d0c30SPeter Brune row_c = lcid[col]; 3758e6d0c30SPeter Brune /* set the values for on-processor ones */ 3761ce39c63SPeter Brune if (PetscRealPart(rval[j]) < 0.) { 3771ce39c63SPeter Brune pij = rval[j] * alpha * invdiag; 3788e6d0c30SPeter Brune } else { 3791ce39c63SPeter Brune pij = rval[j] * beta * invdiag; 3808e6d0c30SPeter Brune } 3818e6d0c30SPeter Brune if (PetscAbsScalar(pij) != 0.) { 3828e6d0c30SPeter Brune pvals[idx] = pij; 3838e6d0c30SPeter Brune pcols[idx] = row_c; 3848e6d0c30SPeter Brune idx++; 3858e6d0c30SPeter Brune } 3868e6d0c30SPeter Brune } 3878e6d0c30SPeter Brune } 3889566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA, i, &ncols, &rcol, &rval)); 3898e6d0c30SPeter Brune /* off */ 3901ce39c63SPeter Brune if (gA) { 3919566063dSJacob Faibussowitsch PetscCall(MatGetRow(gA, i, &ncols, &rcol, &rval)); 3928e6d0c30SPeter Brune for (j = 0; j < ncols; j++) { 3931ce39c63SPeter Brune col = rcol[j]; 394c1eae691SMark Adams if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0] * Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0] * Amax_neg[i])) { 3958e6d0c30SPeter Brune row_f = i + fs; 3968e6d0c30SPeter Brune row_c = gcid[col]; 3978e6d0c30SPeter Brune /* set the values for on-processor ones */ 3981ce39c63SPeter Brune if (PetscRealPart(rval[j]) < 0.) { 3991ce39c63SPeter Brune pij = rval[j] * alpha * invdiag; 4008e6d0c30SPeter Brune } else { 4011ce39c63SPeter Brune pij = rval[j] * beta * invdiag; 4028e6d0c30SPeter Brune } 4038e6d0c30SPeter Brune if (PetscAbsScalar(pij) != 0.) { 4048e6d0c30SPeter Brune pvals[idx] = pij; 4058e6d0c30SPeter Brune pcols[idx] = row_c; 4068e6d0c30SPeter Brune idx++; 4078e6d0c30SPeter Brune } 4088e6d0c30SPeter Brune } 4098e6d0c30SPeter Brune } 4109566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(gA, i, &ncols, &rcol, &rval)); 4113c9ab2c3SPeter Brune } 4129566063dSJacob Faibussowitsch PetscCall(MatSetValues(*P, 1, &row_f, idx, pcols, pvals, INSERT_VALUES)); 4138e6d0c30SPeter Brune } 4148e6d0c30SPeter Brune } 4153c9ab2c3SPeter Brune 4169566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY)); 4179566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY)); 4188e6d0c30SPeter Brune 4199566063dSJacob Faibussowitsch PetscCall(PetscFree5(lsparse, gsparse, lcid, Amax_pos, Amax_neg)); 420e632b94dSBarry Smith 4219566063dSJacob Faibussowitsch PetscCall(PetscFree2(pcols, pvals)); 42287b9b13cSPeter Brune if (gA) { 4239566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sf)); 4249566063dSJacob Faibussowitsch PetscCall(PetscFree(gcid)); 42587b9b13cSPeter Brune } 4268e6d0c30SPeter Brune PetscFunctionReturn(0); 4278e6d0c30SPeter Brune } 4288e6d0c30SPeter Brune 4299371c9d4SSatish Balay PetscErrorCode PCGAMGTruncateProlongator_Private(PC pc, Mat *P) { 430586a8384SPeter Brune PetscInt j, i, ps, pf, pn, pcs, pcf, pcn, idx, cmax; 431586a8384SPeter Brune const PetscScalar *pval; 432586a8384SPeter Brune const PetscInt *pcol; 433586a8384SPeter Brune PetscScalar *pnval; 434586a8384SPeter Brune PetscInt *pncol; 435586a8384SPeter Brune PetscInt ncols; 436586a8384SPeter Brune Mat Pnew; 437586a8384SPeter Brune PetscInt *lsparse, *gsparse; 438586a8384SPeter Brune PetscReal pmax_pos, pmax_neg, ptot_pos, ptot_neg, pthresh_pos, pthresh_neg; 439586a8384SPeter Brune PC_MG *mg = (PC_MG *)pc->data; 440586a8384SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 441586a8384SPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical *)pc_gamg->subctx; 442d9558ea9SBarry Smith MatType mtype; 443586a8384SPeter Brune 444586a8384SPeter Brune PetscFunctionBegin; 445586a8384SPeter Brune /* trim and rescale with reallocation */ 4469566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(*P, &ps, &pf)); 4479566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(*P, &pcs, &pcf)); 448586a8384SPeter Brune pn = pf - ps; 449586a8384SPeter Brune pcn = pcf - pcs; 4509566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(pn, &lsparse, pn, &gsparse)); 451586a8384SPeter Brune /* allocate */ 452586a8384SPeter Brune cmax = 0; 453586a8384SPeter Brune for (i = ps; i < pf; i++) { 454b4dc3ebdSPeter Brune lsparse[i - ps] = 0; 455b4dc3ebdSPeter Brune gsparse[i - ps] = 0; 4569566063dSJacob Faibussowitsch PetscCall(MatGetRow(*P, i, &ncols, &pcol, &pval)); 457ad540459SPierre Jolivet if (ncols > cmax) cmax = ncols; 458586a8384SPeter Brune pmax_pos = 0.; 459586a8384SPeter Brune pmax_neg = 0.; 460586a8384SPeter Brune for (j = 0; j < ncols; j++) { 461586a8384SPeter Brune if (PetscRealPart(pval[j]) > pmax_pos) { 462586a8384SPeter Brune pmax_pos = PetscRealPart(pval[j]); 463586a8384SPeter Brune } else if (PetscRealPart(pval[j]) < pmax_neg) { 464586a8384SPeter Brune pmax_neg = PetscRealPart(pval[j]); 465586a8384SPeter Brune } 466586a8384SPeter Brune } 467586a8384SPeter Brune for (j = 0; j < ncols; j++) { 4688eab0cc1SPeter Brune if (PetscRealPart(pval[j]) >= pmax_pos * cls->interp_threshold || PetscRealPart(pval[j]) <= pmax_neg * cls->interp_threshold) { 469b4dc3ebdSPeter Brune if (pcol[j] >= pcs && pcol[j] < pcf) { 470b4dc3ebdSPeter Brune lsparse[i - ps]++; 471586a8384SPeter Brune } else { 472b4dc3ebdSPeter Brune gsparse[i - ps]++; 473586a8384SPeter Brune } 474586a8384SPeter Brune } 475586a8384SPeter Brune } 4769566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(*P, i, &ncols, &pcol, &pval)); 477586a8384SPeter Brune } 478586a8384SPeter Brune 4799566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(cmax, &pnval, cmax, &pncol)); 480586a8384SPeter Brune 4819566063dSJacob Faibussowitsch PetscCall(MatGetType(*P, &mtype)); 4829566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)*P), &Pnew)); 4839566063dSJacob Faibussowitsch PetscCall(MatSetType(Pnew, mtype)); 4849566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Pnew, pn, pcn, PETSC_DETERMINE, PETSC_DETERMINE)); 4859566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(Pnew, 0, lsparse)); 4869566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(Pnew, 0, lsparse, 0, gsparse)); 487586a8384SPeter Brune 488586a8384SPeter Brune for (i = ps; i < pf; i++) { 4899566063dSJacob Faibussowitsch PetscCall(MatGetRow(*P, i, &ncols, &pcol, &pval)); 490586a8384SPeter Brune pmax_pos = 0.; 491586a8384SPeter Brune pmax_neg = 0.; 492586a8384SPeter Brune for (j = 0; j < ncols; j++) { 493586a8384SPeter Brune if (PetscRealPart(pval[j]) > pmax_pos) { 494586a8384SPeter Brune pmax_pos = PetscRealPart(pval[j]); 495586a8384SPeter Brune } else if (PetscRealPart(pval[j]) < pmax_neg) { 496586a8384SPeter Brune pmax_neg = PetscRealPart(pval[j]); 497586a8384SPeter Brune } 498586a8384SPeter Brune } 499586a8384SPeter Brune pthresh_pos = 0.; 500586a8384SPeter Brune pthresh_neg = 0.; 501586a8384SPeter Brune ptot_pos = 0.; 502586a8384SPeter Brune ptot_neg = 0.; 503586a8384SPeter Brune for (j = 0; j < ncols; j++) { 5048eab0cc1SPeter Brune if (PetscRealPart(pval[j]) >= cls->interp_threshold * pmax_pos) { 505586a8384SPeter Brune pthresh_pos += PetscRealPart(pval[j]); 5068eab0cc1SPeter Brune } else if (PetscRealPart(pval[j]) <= cls->interp_threshold * pmax_neg) { 507586a8384SPeter Brune pthresh_neg += PetscRealPart(pval[j]); 508586a8384SPeter Brune } 509586a8384SPeter Brune if (PetscRealPart(pval[j]) > 0.) { 510586a8384SPeter Brune ptot_pos += PetscRealPart(pval[j]); 511586a8384SPeter Brune } else { 512586a8384SPeter Brune ptot_neg += PetscRealPart(pval[j]); 513586a8384SPeter Brune } 514586a8384SPeter Brune } 5156bd8ea90SPeter Brune if (PetscAbsReal(pthresh_pos) > 0.) ptot_pos /= pthresh_pos; 5166bd8ea90SPeter Brune if (PetscAbsReal(pthresh_neg) > 0.) ptot_neg /= pthresh_neg; 517586a8384SPeter Brune idx = 0; 518586a8384SPeter Brune for (j = 0; j < ncols; j++) { 5198eab0cc1SPeter Brune if (PetscRealPart(pval[j]) >= pmax_pos * cls->interp_threshold) { 520586a8384SPeter Brune pnval[idx] = ptot_pos * pval[j]; 521586a8384SPeter Brune pncol[idx] = pcol[j]; 522586a8384SPeter Brune idx++; 5238eab0cc1SPeter Brune } else if (PetscRealPart(pval[j]) <= pmax_neg * cls->interp_threshold) { 524586a8384SPeter Brune pnval[idx] = ptot_neg * pval[j]; 525586a8384SPeter Brune pncol[idx] = pcol[j]; 526586a8384SPeter Brune idx++; 527586a8384SPeter Brune } 528586a8384SPeter Brune } 5299566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(*P, i, &ncols, &pcol, &pval)); 5309566063dSJacob Faibussowitsch PetscCall(MatSetValues(Pnew, 1, &i, idx, pncol, pnval, INSERT_VALUES)); 531586a8384SPeter Brune } 532586a8384SPeter Brune 5339566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Pnew, MAT_FINAL_ASSEMBLY)); 5349566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Pnew, MAT_FINAL_ASSEMBLY)); 5359566063dSJacob Faibussowitsch PetscCall(MatDestroy(P)); 536586a8384SPeter Brune 537586a8384SPeter Brune *P = Pnew; 5389566063dSJacob Faibussowitsch PetscCall(PetscFree2(lsparse, gsparse)); 5399566063dSJacob Faibussowitsch PetscCall(PetscFree2(pnval, pncol)); 540586a8384SPeter Brune PetscFunctionReturn(0); 541586a8384SPeter Brune } 542586a8384SPeter Brune 5439371c9d4SSatish Balay PetscErrorCode PCGAMGProlongator_Classical_Standard(PC pc, Mat A, Mat G, PetscCoarsenData *agg_lists, Mat *P) { 544c634539dSPeter Brune Mat lA, *lAs; 545f9a65ec8SPeter Brune MatType mtype; 54663b0c588SPeter Brune Vec cv; 54763b0c588SPeter Brune PetscInt *gcid, *lcid, *lsparse, *gsparse, *picol; 548420cd23fSPeter Brune PetscInt fs, fe, cs, ce, nl, i, j, k, li, lni, ci, ncols, maxcols, fn, cn, cid; 549420cd23fSPeter Brune PetscMPIInt size; 55063b0c588SPeter Brune const PetscInt *lidx, *icol, *gidx; 55163b0c588SPeter Brune PetscBool iscoarse; 55263b0c588SPeter Brune PetscScalar vi, pentry, pjentry; 55363b0c588SPeter Brune PetscScalar *pcontrib, *pvcol; 55463b0c588SPeter Brune const PetscScalar *vcol; 555ed4e82eaSPeter Brune PetscReal diag, jdiag, jwttotal; 556f9a65ec8SPeter Brune PetscInt pncols; 557c634539dSPeter Brune PetscSF sf; 558c634539dSPeter Brune PetscLayout clayout; 55963b0c588SPeter Brune IS lis; 560f9a65ec8SPeter Brune 561f9a65ec8SPeter Brune PetscFunctionBegin; 5629566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size)); 5639566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &fs, &fe)); 564f9a65ec8SPeter Brune fn = fe - fs; 5659566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, fe - fs, fs, 1, &lis)); 566c634539dSPeter Brune if (size > 1) { 5679566063dSJacob Faibussowitsch PetscCall(MatGetLayouts(A, NULL, &clayout)); 568f9a65ec8SPeter Brune /* increase the overlap by two to get neighbors of neighbors */ 5699566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(A, 1, &lis, 2)); 5709566063dSJacob Faibussowitsch PetscCall(ISSort(lis)); 571f9a65ec8SPeter Brune /* get the local part of A */ 5729566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(A, 1, &lis, &lis, MAT_INITIAL_MATRIX, &lAs)); 573c634539dSPeter Brune lA = lAs[0]; 574c634539dSPeter Brune /* build an SF out of it */ 5759566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(lis, &nl)); 5769566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)A), &sf)); 5779566063dSJacob Faibussowitsch PetscCall(ISGetIndices(lis, &lidx)); 5789566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraphLayout(sf, clayout, nl, NULL, PETSC_COPY_VALUES, lidx)); 5799566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(lis, &lidx)); 580c634539dSPeter Brune } else { 581c634539dSPeter Brune lA = A; 582c634539dSPeter Brune nl = fn; 583c634539dSPeter Brune } 584c634539dSPeter Brune /* create a communication structure for the overlapped portion and transmit coarse indices */ 5859566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(fn, &lsparse, fn, &gsparse, nl, &pcontrib)); 586f9a65ec8SPeter Brune /* create coarse vector */ 587f9a65ec8SPeter Brune cn = 0; 588f9a65ec8SPeter Brune for (i = 0; i < fn; i++) { 5899566063dSJacob Faibussowitsch PetscCall(PetscCDEmptyAt(agg_lists, i, &iscoarse)); 590ad540459SPierre Jolivet if (!iscoarse) cn++; 591f9a65ec8SPeter Brune } 5929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(fn, &gcid)); 5939566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)A), cn, PETSC_DECIDE, &cv)); 5949566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(cv, &cs, &ce)); 595f9a65ec8SPeter Brune cn = 0; 596f9a65ec8SPeter Brune for (i = 0; i < fn; i++) { 5979566063dSJacob Faibussowitsch PetscCall(PetscCDEmptyAt(agg_lists, i, &iscoarse)); 598f9a65ec8SPeter Brune if (!iscoarse) { 599c634539dSPeter Brune gcid[i] = cs + cn; 600f9a65ec8SPeter Brune cn++; 601f9a65ec8SPeter Brune } else { 602c634539dSPeter Brune gcid[i] = -1; 603f9a65ec8SPeter Brune } 604f9a65ec8SPeter Brune } 605c634539dSPeter Brune if (size > 1) { 6069566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nl, &lcid)); 6079566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sf, MPIU_INT, gcid, lcid, MPI_REPLACE)); 6089566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sf, MPIU_INT, gcid, lcid, MPI_REPLACE)); 609c634539dSPeter Brune } else { 610c634539dSPeter Brune lcid = gcid; 611c634539dSPeter Brune } 612f9a65ec8SPeter Brune /* count to preallocate the prolongator */ 6139566063dSJacob Faibussowitsch PetscCall(ISGetIndices(lis, &gidx)); 614f9a65ec8SPeter Brune maxcols = 0; 615f9a65ec8SPeter Brune /* count the number of unique contributing coarse cells for each fine */ 616f9a65ec8SPeter Brune for (i = 0; i < nl; i++) { 617ed4e82eaSPeter Brune pcontrib[i] = 0.; 6189566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA, i, &ncols, &icol, NULL)); 619f9a65ec8SPeter Brune if (gidx[i] >= fs && gidx[i] < fe) { 620f9a65ec8SPeter Brune li = gidx[i] - fs; 621f9a65ec8SPeter Brune lsparse[li] = 0; 622f9a65ec8SPeter Brune gsparse[li] = 0; 623c634539dSPeter Brune cid = lcid[i]; 624f9a65ec8SPeter Brune if (cid >= 0) { 625f9a65ec8SPeter Brune lsparse[li] = 1; 626f9a65ec8SPeter Brune } else { 627f9a65ec8SPeter Brune for (j = 0; j < ncols; j++) { 628c634539dSPeter Brune if (lcid[icol[j]] >= 0) { 629f9a65ec8SPeter Brune pcontrib[icol[j]] = 1.; 630f9a65ec8SPeter Brune } else { 631f9a65ec8SPeter Brune ci = icol[j]; 6329566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA, i, &ncols, &icol, NULL)); 6339566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA, ci, &ncols, &icol, NULL)); 634f9a65ec8SPeter Brune for (k = 0; k < ncols; k++) { 635ad540459SPierre Jolivet if (lcid[icol[k]] >= 0) pcontrib[icol[k]] = 1.; 636f9a65ec8SPeter Brune } 6379566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA, ci, &ncols, &icol, NULL)); 6389566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA, i, &ncols, &icol, NULL)); 639f9a65ec8SPeter Brune } 640f9a65ec8SPeter Brune } 641f9a65ec8SPeter Brune for (j = 0; j < ncols; j++) { 642c634539dSPeter Brune if (lcid[icol[j]] >= 0 && pcontrib[icol[j]] != 0.) { 643c634539dSPeter Brune lni = lcid[icol[j]]; 644f9a65ec8SPeter Brune if (lni >= cs && lni < ce) { 645f9a65ec8SPeter Brune lsparse[li]++; 646f9a65ec8SPeter Brune } else { 647f9a65ec8SPeter Brune gsparse[li]++; 648f9a65ec8SPeter Brune } 649f9a65ec8SPeter Brune pcontrib[icol[j]] = 0.; 650f9a65ec8SPeter Brune } else { 651f9a65ec8SPeter Brune ci = icol[j]; 6529566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA, i, &ncols, &icol, NULL)); 6539566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA, ci, &ncols, &icol, NULL)); 654f9a65ec8SPeter Brune for (k = 0; k < ncols; k++) { 655c634539dSPeter Brune if (lcid[icol[k]] >= 0 && pcontrib[icol[k]] != 0.) { 656c634539dSPeter Brune lni = lcid[icol[k]]; 657f9a65ec8SPeter Brune if (lni >= cs && lni < ce) { 658f9a65ec8SPeter Brune lsparse[li]++; 659f9a65ec8SPeter Brune } else { 660f9a65ec8SPeter Brune gsparse[li]++; 661f9a65ec8SPeter Brune } 662f9a65ec8SPeter Brune pcontrib[icol[k]] = 0.; 663f9a65ec8SPeter Brune } 664f9a65ec8SPeter Brune } 6659566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA, ci, &ncols, &icol, NULL)); 6669566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA, i, &ncols, &icol, NULL)); 667f9a65ec8SPeter Brune } 668f9a65ec8SPeter Brune } 669ed4e82eaSPeter Brune } 670f9a65ec8SPeter Brune if (lsparse[li] + gsparse[li] > maxcols) maxcols = lsparse[li] + gsparse[li]; 671ed4e82eaSPeter Brune } 6729566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA, i, &ncols, &icol, &vcol)); 673f9a65ec8SPeter Brune } 6749566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(maxcols, &picol, maxcols, &pvcol)); 6759566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), P)); 6769566063dSJacob Faibussowitsch PetscCall(MatGetType(A, &mtype)); 6779566063dSJacob Faibussowitsch PetscCall(MatSetType(*P, mtype)); 6789566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*P, fn, cn, PETSC_DETERMINE, PETSC_DETERMINE)); 6799566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(*P, 0, lsparse, 0, gsparse)); 6809566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(*P, 0, lsparse)); 681f9a65ec8SPeter Brune for (i = 0; i < nl; i++) { 682ed4e82eaSPeter Brune diag = 0.; 683f9a65ec8SPeter Brune if (gidx[i] >= fs && gidx[i] < fe) { 684f9a65ec8SPeter Brune pncols = 0; 685c634539dSPeter Brune cid = lcid[i]; 686f9a65ec8SPeter Brune if (cid >= 0) { 687f9a65ec8SPeter Brune pncols = 1; 688f9a65ec8SPeter Brune picol[0] = cid; 689f9a65ec8SPeter Brune pvcol[0] = 1.; 690f9a65ec8SPeter Brune } else { 6919566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA, i, &ncols, &icol, &vcol)); 692f9a65ec8SPeter Brune for (j = 0; j < ncols; j++) { 693ed4e82eaSPeter Brune pentry = vcol[j]; 694c634539dSPeter Brune if (lcid[icol[j]] >= 0) { 695f9a65ec8SPeter Brune /* coarse neighbor */ 696ed4e82eaSPeter Brune pcontrib[icol[j]] += pentry; 697ed4e82eaSPeter Brune } else if (icol[j] != i) { 698f9a65ec8SPeter Brune /* the neighbor is a strongly connected fine node */ 699f9a65ec8SPeter Brune ci = icol[j]; 700f9a65ec8SPeter Brune vi = vcol[j]; 7019566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA, i, &ncols, &icol, &vcol)); 7029566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA, ci, &ncols, &icol, &vcol)); 703ed4e82eaSPeter Brune jwttotal = 0.; 704f9a65ec8SPeter Brune jdiag = 0.; 705f9a65ec8SPeter Brune for (k = 0; k < ncols; k++) { 706ad540459SPierre Jolivet if (ci == icol[k]) jdiag = PetscRealPart(vcol[k]); 707f9a65ec8SPeter Brune } 708f9a65ec8SPeter Brune for (k = 0; k < ncols; k++) { 709c634539dSPeter Brune if (lcid[icol[k]] >= 0 && jdiag * PetscRealPart(vcol[k]) < 0.) { 710ed4e82eaSPeter Brune pjentry = vcol[k]; 7117779008dSPeter Brune jwttotal += PetscRealPart(pjentry); 712f9a65ec8SPeter Brune } 713f9a65ec8SPeter Brune } 714ed4e82eaSPeter Brune if (jwttotal != 0.) { 7157779008dSPeter Brune jwttotal = PetscRealPart(vi) / jwttotal; 716ed4e82eaSPeter Brune for (k = 0; k < ncols; k++) { 717c634539dSPeter Brune if (lcid[icol[k]] >= 0 && jdiag * PetscRealPart(vcol[k]) < 0.) { 718586a8384SPeter Brune pjentry = vcol[k] * jwttotal; 719ed4e82eaSPeter Brune pcontrib[icol[k]] += pjentry; 720ed4e82eaSPeter Brune } 721ed4e82eaSPeter Brune } 722ed4e82eaSPeter Brune } else { 723ed4e82eaSPeter Brune diag += PetscRealPart(vi); 724ed4e82eaSPeter Brune } 7259566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA, ci, &ncols, &icol, &vcol)); 7269566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA, i, &ncols, &icol, &vcol)); 727ed4e82eaSPeter Brune } else { 728ed4e82eaSPeter Brune diag += PetscRealPart(vcol[j]); 729f9a65ec8SPeter Brune } 730f9a65ec8SPeter Brune } 731586a8384SPeter Brune if (diag != 0.) { 732586a8384SPeter Brune diag = 1. / diag; 733f9a65ec8SPeter Brune for (j = 0; j < ncols; j++) { 734c634539dSPeter Brune if (lcid[icol[j]] >= 0 && pcontrib[icol[j]] != 0.) { 735f9a65ec8SPeter Brune /* the neighbor is a coarse node */ 736ed4e82eaSPeter Brune if (PetscAbsScalar(pcontrib[icol[j]]) > 0.0) { 737c634539dSPeter Brune lni = lcid[icol[j]]; 738586a8384SPeter Brune pvcol[pncols] = -pcontrib[icol[j]] * diag; 739f9a65ec8SPeter Brune picol[pncols] = lni; 740f9a65ec8SPeter Brune pncols++; 741ed4e82eaSPeter Brune } 742ed4e82eaSPeter Brune pcontrib[icol[j]] = 0.; 743f9a65ec8SPeter Brune } else { 744f9a65ec8SPeter Brune /* the neighbor is a strongly connected fine node */ 745f9a65ec8SPeter Brune ci = icol[j]; 7469566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA, i, &ncols, &icol, &vcol)); 7479566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA, ci, &ncols, &icol, &vcol)); 748f9a65ec8SPeter Brune for (k = 0; k < ncols; k++) { 749c634539dSPeter Brune if (lcid[icol[k]] >= 0 && pcontrib[icol[k]] != 0.) { 750ed4e82eaSPeter Brune if (PetscAbsScalar(pcontrib[icol[k]]) > 0.0) { 751c634539dSPeter Brune lni = lcid[icol[k]]; 752586a8384SPeter Brune pvcol[pncols] = -pcontrib[icol[k]] * diag; 753f9a65ec8SPeter Brune picol[pncols] = lni; 754f9a65ec8SPeter Brune pncols++; 755f9a65ec8SPeter Brune } 756ed4e82eaSPeter Brune pcontrib[icol[k]] = 0.; 757ed4e82eaSPeter Brune } 758f9a65ec8SPeter Brune } 7599566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA, ci, &ncols, &icol, &vcol)); 7609566063dSJacob Faibussowitsch PetscCall(MatGetRow(lA, i, &ncols, &icol, &vcol)); 761f9a65ec8SPeter Brune } 762ed4e82eaSPeter Brune pcontrib[icol[j]] = 0.; 763f9a65ec8SPeter Brune } 7649566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(lA, i, &ncols, &icol, &vcol)); 765f9a65ec8SPeter Brune } 766586a8384SPeter Brune } 767f9a65ec8SPeter Brune ci = gidx[i]; 76848a46eb9SPierre Jolivet if (pncols > 0) PetscCall(MatSetValues(*P, 1, &ci, pncols, picol, pvcol, INSERT_VALUES)); 769f9a65ec8SPeter Brune } 770f9a65ec8SPeter Brune } 7719566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(lis, &gidx)); 7729566063dSJacob Faibussowitsch PetscCall(PetscFree2(picol, pvcol)); 7739566063dSJacob Faibussowitsch PetscCall(PetscFree3(lsparse, gsparse, pcontrib)); 7749566063dSJacob Faibussowitsch PetscCall(ISDestroy(&lis)); 7759566063dSJacob Faibussowitsch PetscCall(PetscFree(gcid)); 776c634539dSPeter Brune if (size > 1) { 7779566063dSJacob Faibussowitsch PetscCall(PetscFree(lcid)); 7789566063dSJacob Faibussowitsch PetscCall(MatDestroyMatrices(1, &lAs)); 7799566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sf)); 780c634539dSPeter Brune } 7819566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cv)); 7829566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY)); 7839566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY)); 7848eab0cc1SPeter Brune PetscFunctionReturn(0); 7858eab0cc1SPeter Brune } 786f9a65ec8SPeter Brune 7879371c9d4SSatish Balay PetscErrorCode PCGAMGOptProlongator_Classical_Jacobi(PC pc, Mat A, Mat *P) { 788b4dc3ebdSPeter Brune PetscInt f, s, n, cf, cs, i, idx; 789b4dc3ebdSPeter Brune PetscInt *coarserows; 790b4dc3ebdSPeter Brune PetscInt ncols; 791b4dc3ebdSPeter Brune const PetscInt *pcols; 792b4dc3ebdSPeter Brune const PetscScalar *pvals; 793b4dc3ebdSPeter Brune Mat Pnew; 794b4dc3ebdSPeter Brune Vec diag; 795b4dc3ebdSPeter Brune PC_MG *mg = (PC_MG *)pc->data; 796b4dc3ebdSPeter Brune PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 797b4dc3ebdSPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical *)pc_gamg->subctx; 798b4dc3ebdSPeter Brune 799b4dc3ebdSPeter Brune PetscFunctionBegin; 800b4dc3ebdSPeter Brune if (cls->nsmooths == 0) { 8019566063dSJacob Faibussowitsch PetscCall(PCGAMGTruncateProlongator_Private(pc, P)); 802b4dc3ebdSPeter Brune PetscFunctionReturn(0); 803b4dc3ebdSPeter Brune } 8049566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(*P, &s, &f)); 805b4dc3ebdSPeter Brune n = f - s; 8069566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(*P, &cs, &cf)); 8079566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &coarserows)); 808b4dc3ebdSPeter Brune /* identify the rows corresponding to coarse unknowns */ 809b4dc3ebdSPeter Brune idx = 0; 810b4dc3ebdSPeter Brune for (i = s; i < f; i++) { 8119566063dSJacob Faibussowitsch PetscCall(MatGetRow(*P, i, &ncols, &pcols, &pvals)); 812b4dc3ebdSPeter Brune /* assume, for now, that it's a coarse unknown if it has a single unit entry */ 813b4dc3ebdSPeter Brune if (ncols == 1) { 814b4dc3ebdSPeter Brune if (pvals[0] == 1.) { 815b4dc3ebdSPeter Brune coarserows[idx] = i; 816b4dc3ebdSPeter Brune idx++; 817b4dc3ebdSPeter Brune } 818b4dc3ebdSPeter Brune } 8199566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(*P, i, &ncols, &pcols, &pvals)); 820b4dc3ebdSPeter Brune } 8219566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &diag, NULL)); 8229566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(A, diag)); 8239566063dSJacob Faibussowitsch PetscCall(VecReciprocal(diag)); 824b4dc3ebdSPeter Brune for (i = 0; i < cls->nsmooths; i++) { 8259566063dSJacob Faibussowitsch PetscCall(MatMatMult(A, *P, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Pnew)); 8269566063dSJacob Faibussowitsch PetscCall(MatZeroRows(Pnew, idx, coarserows, 0., NULL, NULL)); 8279566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(Pnew, diag, NULL)); 8289566063dSJacob Faibussowitsch PetscCall(MatAYPX(Pnew, -1.0, *P, DIFFERENT_NONZERO_PATTERN)); 8299566063dSJacob Faibussowitsch PetscCall(MatDestroy(P)); 830b4dc3ebdSPeter Brune *P = Pnew; 831b4dc3ebdSPeter Brune Pnew = NULL; 832b4dc3ebdSPeter Brune } 8339566063dSJacob Faibussowitsch PetscCall(VecDestroy(&diag)); 8349566063dSJacob Faibussowitsch PetscCall(PetscFree(coarserows)); 8359566063dSJacob Faibussowitsch PetscCall(PCGAMGTruncateProlongator_Private(pc, P)); 836b4dc3ebdSPeter Brune PetscFunctionReturn(0); 837b4dc3ebdSPeter Brune } 838b4dc3ebdSPeter Brune 8399371c9d4SSatish Balay static PetscErrorCode PCGAMGProlongator_Classical(PC pc, Mat A, Mat G, PetscCoarsenData *agg_lists, Mat *P) { 8408eab0cc1SPeter Brune PetscErrorCode (*f)(PC, Mat, Mat, PetscCoarsenData *, Mat *); 8418eab0cc1SPeter Brune PC_MG *mg = (PC_MG *)pc->data; 8428eab0cc1SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 8438eab0cc1SPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical *)pc_gamg->subctx; 8448eab0cc1SPeter Brune 8458eab0cc1SPeter Brune PetscFunctionBegin; 8469566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(PCGAMGClassicalProlongatorList, cls->prolongtype, &f)); 84728b400f6SJacob Faibussowitsch PetscCheck(f, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Cannot find PCGAMG Classical prolongator type"); 8489566063dSJacob Faibussowitsch PetscCall((*f)(pc, A, G, agg_lists, P)); 849f9a65ec8SPeter Brune PetscFunctionReturn(0); 850f9a65ec8SPeter Brune } 851f9a65ec8SPeter Brune 8529371c9d4SSatish Balay static PetscErrorCode PCGAMGDestroy_Classical(PC pc) { 8538e6d0c30SPeter Brune PC_MG *mg = (PC_MG *)pc->data; 8548e6d0c30SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 8558e6d0c30SPeter Brune 8568e6d0c30SPeter Brune PetscFunctionBegin; 8579566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_gamg->subctx)); 8589566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGClassicalSetType_C", NULL)); 8599566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGClassicalGetType_C", NULL)); 8608e6d0c30SPeter Brune PetscFunctionReturn(0); 8618e6d0c30SPeter Brune } 8628e6d0c30SPeter Brune 8639371c9d4SSatish Balay PetscErrorCode PCGAMGSetFromOptions_Classical(PC pc, PetscOptionItems *PetscOptionsObject) { 864586a8384SPeter Brune PC_MG *mg = (PC_MG *)pc->data; 865586a8384SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 866586a8384SPeter Brune PC_GAMG_Classical *cls = (PC_GAMG_Classical *)pc_gamg->subctx; 8678eab0cc1SPeter Brune char tname[256]; 8688eab0cc1SPeter Brune PetscBool flg; 869586a8384SPeter Brune 8708e6d0c30SPeter Brune PetscFunctionBegin; 871d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "GAMG-Classical options"); 8729566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-pc_gamg_classical_type", "Type of Classical AMG prolongation", "PCGAMGClassicalSetType", PCGAMGClassicalProlongatorList, cls->prolongtype, tname, sizeof(tname), &flg)); 8731baa6e33SBarry Smith if (flg) PetscCall(PCGAMGClassicalSetType(pc, tname)); 8749566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_gamg_classical_interp_threshold", "Threshold for classical interpolator entries", "", cls->interp_threshold, &cls->interp_threshold, NULL)); 8759566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_gamg_classical_nsmooths", "Threshold for classical interpolator entries", "", cls->nsmooths, &cls->nsmooths, NULL)); 876d0609cedSBarry Smith PetscOptionsHeadEnd(); 8778e6d0c30SPeter Brune PetscFunctionReturn(0); 8788e6d0c30SPeter Brune } 8798e6d0c30SPeter Brune 8809371c9d4SSatish Balay static PetscErrorCode PCGAMGSetData_Classical(PC pc, Mat A) { 8818e6d0c30SPeter Brune PC_MG *mg = (PC_MG *)pc->data; 8828e6d0c30SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 8838e6d0c30SPeter Brune 8848e6d0c30SPeter Brune PetscFunctionBegin; 8858e6d0c30SPeter Brune /* no data for classical AMG */ 8868e6d0c30SPeter Brune pc_gamg->data = NULL; 887d2050638SMark Adams pc_gamg->data_cell_cols = 0; 888d2050638SMark Adams pc_gamg->data_cell_rows = 0; 8898e6d0c30SPeter Brune pc_gamg->data_sz = 0; 8908e6d0c30SPeter Brune PetscFunctionReturn(0); 8918e6d0c30SPeter Brune } 8928e6d0c30SPeter Brune 8939371c9d4SSatish Balay PetscErrorCode PCGAMGClassicalFinalizePackage(void) { 8948eab0cc1SPeter Brune PetscFunctionBegin; 8958eab0cc1SPeter Brune PCGAMGClassicalPackageInitialized = PETSC_FALSE; 8969566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&PCGAMGClassicalProlongatorList)); 8978eab0cc1SPeter Brune PetscFunctionReturn(0); 8988eab0cc1SPeter Brune } 8998eab0cc1SPeter Brune 9009371c9d4SSatish Balay PetscErrorCode PCGAMGClassicalInitializePackage(void) { 9018eab0cc1SPeter Brune PetscFunctionBegin; 9028eab0cc1SPeter Brune if (PCGAMGClassicalPackageInitialized) PetscFunctionReturn(0); 9039566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&PCGAMGClassicalProlongatorList, PCGAMGCLASSICALDIRECT, PCGAMGProlongator_Classical_Direct)); 9049566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&PCGAMGClassicalProlongatorList, PCGAMGCLASSICALSTANDARD, PCGAMGProlongator_Classical_Standard)); 9059566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(PCGAMGClassicalFinalizePackage)); 9068eab0cc1SPeter Brune PetscFunctionReturn(0); 9078eab0cc1SPeter Brune } 9088eab0cc1SPeter Brune 9098e6d0c30SPeter Brune /* 9108e6d0c30SPeter Brune PCCreateGAMG_Classical 9118e6d0c30SPeter Brune 9128e6d0c30SPeter Brune */ 9139371c9d4SSatish Balay PetscErrorCode PCCreateGAMG_Classical(PC pc) { 9148e6d0c30SPeter Brune PC_MG *mg = (PC_MG *)pc->data; 9158e6d0c30SPeter Brune PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx; 9168e6d0c30SPeter Brune PC_GAMG_Classical *pc_gamg_classical; 9178e6d0c30SPeter Brune 9188e6d0c30SPeter Brune PetscFunctionBegin; 9199566063dSJacob Faibussowitsch PetscCall(PCGAMGClassicalInitializePackage()); 9208e6d0c30SPeter Brune if (pc_gamg->subctx) { 9218e6d0c30SPeter Brune /* call base class */ 9229566063dSJacob Faibussowitsch PetscCall(PCDestroy_GAMG(pc)); 9238e6d0c30SPeter Brune } 9248e6d0c30SPeter Brune 9258e6d0c30SPeter Brune /* create sub context for SA */ 9269566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc, &pc_gamg_classical)); 9278e6d0c30SPeter Brune pc_gamg->subctx = pc_gamg_classical; 9288e6d0c30SPeter Brune pc->ops->setfromoptions = PCGAMGSetFromOptions_Classical; 9298e6d0c30SPeter Brune /* reset does not do anything; setup not virtual */ 9308e6d0c30SPeter Brune 9318e6d0c30SPeter Brune /* set internal function pointers */ 9328e6d0c30SPeter Brune pc_gamg->ops->destroy = PCGAMGDestroy_Classical; 9338e6d0c30SPeter Brune pc_gamg->ops->graph = PCGAMGGraph_Classical; 9348e6d0c30SPeter Brune pc_gamg->ops->coarsen = PCGAMGCoarsen_Classical; 9358eab0cc1SPeter Brune pc_gamg->ops->prolongator = PCGAMGProlongator_Classical; 936fd1112cbSBarry Smith pc_gamg->ops->optprolongator = PCGAMGOptProlongator_Classical_Jacobi; 937586a8384SPeter Brune pc_gamg->ops->setfromoptions = PCGAMGSetFromOptions_Classical; 9388e6d0c30SPeter Brune 9398e6d0c30SPeter Brune pc_gamg->ops->createdefaultdata = PCGAMGSetData_Classical; 940586a8384SPeter Brune pc_gamg_classical->interp_threshold = 0.2; 941b4dc3ebdSPeter Brune pc_gamg_classical->nsmooths = 0; 9429566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGClassicalSetType_C", PCGAMGClassicalSetType_GAMG)); 9439566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGClassicalGetType_C", PCGAMGClassicalGetType_GAMG)); 9449566063dSJacob Faibussowitsch PetscCall(PCGAMGClassicalSetType(pc, PCGAMGCLASSICALSTANDARD)); 9458e6d0c30SPeter Brune PetscFunctionReturn(0); 9468e6d0c30SPeter Brune } 947