xref: /petsc/src/ksp/pc/impls/gamg/classical.c (revision 6d8694c4fbab79f9439f1ad13c0386ba7ee1ca4b)
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 
13cc4c1da9SBarry Smith /*@
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 
27562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`, `PCGAMGClassicalType`, `PCGAMGClassicalGetType()`
288eab0cc1SPeter Brune @*/
PCGAMGClassicalSetType(PC pc,PCGAMGClassicalType type)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 
37cc4c1da9SBarry Smith /*@
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 
50562efe2eSBarry Smith .seealso: [](ch_ksp), `PCGAMG`, `PCGAMGClassicalType`, `PCGAMGClassicalSetType()`
51c60c7ad4SBarry Smith @*/
PCGAMGClassicalGetType(PC pc,PCGAMGClassicalType * type)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 
PCGAMGClassicalSetType_GAMG(PC pc,PCGAMGClassicalType type)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 
PCGAMGClassicalGetType_GAMG(PC pc,PCGAMGClassicalType * type)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 
PCGAMGCreateGraph_Classical(PC pc,Mat A,Mat * G)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 
PCGAMGCoarsen_Classical(PC pc,Mat * G,PetscCoarsenData ** agg_lists)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;
167e0b7e82fSBarry Smith   const char *prefix;
1688e6d0c30SPeter Brune 
1698e6d0c30SPeter Brune   PetscFunctionBegin;
17028b400f6SJacob Faibussowitsch   PetscCheck(G, fcomm, PETSC_ERR_ARG_WRONGSTATE, "Must set Graph in PC in PCGAMG before coarsening");
1718e6d0c30SPeter Brune 
1729566063dSJacob Faibussowitsch   PetscCall(MatCoarsenCreate(fcomm, &crs));
173e0b7e82fSBarry Smith   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc, &prefix));
174e0b7e82fSBarry Smith   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)crs, prefix));
175e0b7e82fSBarry Smith   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)crs, "pc_gamg_"));
1769566063dSJacob Faibussowitsch   PetscCall(MatCoarsenSetFromOptions(crs));
1779566063dSJacob Faibussowitsch   PetscCall(MatCoarsenSetAdjacency(crs, *G));
1789566063dSJacob Faibussowitsch   PetscCall(MatCoarsenSetStrictAggs(crs, PETSC_TRUE));
1799566063dSJacob Faibussowitsch   PetscCall(MatCoarsenApply(crs));
1809566063dSJacob Faibussowitsch   PetscCall(MatCoarsenGetData(crs, agg_lists));
1819566063dSJacob Faibussowitsch   PetscCall(MatCoarsenDestroy(&crs));
1823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1838e6d0c30SPeter Brune }
1848e6d0c30SPeter Brune 
PCGAMGProlongator_Classical_Direct(PC pc,Mat A,PetscCoarsenData * agg_lists,Mat * P)1858926f930SMark Adams static PetscErrorCode PCGAMGProlongator_Classical_Direct(PC pc, Mat A, PetscCoarsenData *agg_lists, Mat *P)
186d71ae5a4SJacob Faibussowitsch {
1871ce39c63SPeter Brune   PC_MG             *mg   = (PC_MG *)pc->data;
1881ce39c63SPeter Brune   PC_GAMG           *gamg = (PC_GAMG *)mg->innerctx;
18963b0c588SPeter Brune   PetscBool          iscoarse, isMPIAIJ, isSEQAIJ;
19063b0c588SPeter Brune   PetscInt           fn, cn, fs, fe, cs, ce, i, j, ncols, col, row_f, row_c, cmax = 0, idx, noff;
19163b0c588SPeter Brune   PetscInt          *lcid, *gcid, *lsparse, *gsparse, *colmap, *pcols;
19263b0c588SPeter Brune   const PetscInt    *rcol;
19363b0c588SPeter Brune   PetscReal         *Amax_pos, *Amax_neg;
19463b0c588SPeter Brune   PetscScalar        g_pos, g_neg, a_pos, a_neg, diag, invdiag, alpha, beta, pij;
19563b0c588SPeter Brune   PetscScalar       *pvals;
19663b0c588SPeter Brune   const PetscScalar *rval;
19763b0c588SPeter Brune   Mat                lA, gA = NULL;
19863b0c588SPeter Brune   MatType            mtype;
19963b0c588SPeter Brune   Vec                C, lvec;
20087b9b13cSPeter Brune   PetscLayout        clayout;
20187b9b13cSPeter Brune   PetscSF            sf;
20287b9b13cSPeter Brune   Mat_MPIAIJ        *mpiaij;
2038e6d0c30SPeter Brune 
2048e6d0c30SPeter Brune   PetscFunctionBegin;
2059566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &fs, &fe));
20687b9b13cSPeter Brune   fn = fe - fs;
2079566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPIAIJ, &isMPIAIJ));
2089566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQAIJ, &isSEQAIJ));
2097827d75bSBarry Smith   PetscCheck(isMPIAIJ || isSEQAIJ, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Classical AMG requires MPIAIJ matrix");
21087b9b13cSPeter Brune   if (isMPIAIJ) {
21187b9b13cSPeter Brune     mpiaij = (Mat_MPIAIJ *)A->data;
21287b9b13cSPeter Brune     lA     = mpiaij->A;
21387b9b13cSPeter Brune     gA     = mpiaij->B;
21487b9b13cSPeter Brune     lvec   = mpiaij->lvec;
2159566063dSJacob Faibussowitsch     PetscCall(VecGetSize(lvec, &noff));
21687b9b13cSPeter Brune     colmap = mpiaij->garray;
2179566063dSJacob Faibussowitsch     PetscCall(MatGetLayouts(A, NULL, &clayout));
2189566063dSJacob Faibussowitsch     PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)A), &sf));
2199566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraphLayout(sf, clayout, noff, NULL, PETSC_COPY_VALUES, colmap));
2209566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(noff, &gcid));
22187b9b13cSPeter Brune   } else {
22287b9b13cSPeter Brune     lA = A;
22387b9b13cSPeter Brune   }
2249566063dSJacob Faibussowitsch   PetscCall(PetscMalloc5(fn, &lsparse, fn, &gsparse, fn, &lcid, fn, &Amax_pos, fn, &Amax_neg));
2258e6d0c30SPeter Brune 
2268e6d0c30SPeter Brune   /* count the number of coarse unknowns */
2278e6d0c30SPeter Brune   cn = 0;
2288e6d0c30SPeter Brune   for (i = 0; i < fn; i++) {
2298e6d0c30SPeter Brune     /* filter out singletons */
2305e62d202SMark Adams     PetscCall(PetscCDIsEmptyAt(agg_lists, i, &iscoarse));
2318e6d0c30SPeter Brune     lcid[i] = -1;
232ad540459SPierre Jolivet     if (!iscoarse) cn++;
2338e6d0c30SPeter Brune   }
2348e6d0c30SPeter Brune 
2358e6d0c30SPeter Brune   /* create the coarse vector */
2369566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)A), cn, PETSC_DECIDE, &C));
2379566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(C, &cs, &ce));
2388e6d0c30SPeter Brune 
2398e6d0c30SPeter Brune   cn = 0;
2408e6d0c30SPeter Brune   for (i = 0; i < fn; i++) {
2415e62d202SMark Adams     PetscCall(PetscCDIsEmptyAt(agg_lists, i, &iscoarse));
2428e6d0c30SPeter Brune     if (!iscoarse) {
2438e6d0c30SPeter Brune       lcid[i] = cs + cn;
2448e6d0c30SPeter Brune       cn++;
2458e6d0c30SPeter Brune     } else {
2468e6d0c30SPeter Brune       lcid[i] = -1;
2478e6d0c30SPeter Brune     }
2488e6d0c30SPeter Brune   }
2498e6d0c30SPeter Brune 
25087b9b13cSPeter Brune   if (gA) {
2519566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(sf, MPIU_INT, lcid, gcid, MPI_REPLACE));
2529566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(sf, MPIU_INT, lcid, gcid, MPI_REPLACE));
25387b9b13cSPeter Brune   }
2548e6d0c30SPeter Brune 
255b817416eSBarry Smith   /* determine the largest off-diagonal entries in each row */
2561ce39c63SPeter Brune   for (i = fs; i < fe; i++) {
2571ce39c63SPeter Brune     Amax_pos[i - fs] = 0.;
2581ce39c63SPeter Brune     Amax_neg[i - fs] = 0.;
2599566063dSJacob Faibussowitsch     PetscCall(MatGetRow(A, i, &ncols, &rcol, &rval));
2601ce39c63SPeter Brune     for (j = 0; j < ncols; j++) {
2611ce39c63SPeter Brune       if ((PetscRealPart(-rval[j]) > Amax_neg[i - fs]) && i != rcol[j]) Amax_neg[i - fs] = PetscAbsScalar(rval[j]);
2621ce39c63SPeter Brune       if ((PetscRealPart(rval[j]) > Amax_pos[i - fs]) && i != rcol[j]) Amax_pos[i - fs] = PetscAbsScalar(rval[j]);
2631ce39c63SPeter Brune     }
2641ce39c63SPeter Brune     if (ncols > cmax) cmax = ncols;
2659566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(A, i, &ncols, &rcol, &rval));
2661ce39c63SPeter Brune   }
2679566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(cmax, &pcols, cmax, &pvals));
2689566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&C));
2698e6d0c30SPeter Brune 
2708e6d0c30SPeter Brune   /* count the on and off processor sparsity patterns for the prolongator */
2718e6d0c30SPeter Brune   for (i = 0; i < fn; i++) {
2728e6d0c30SPeter Brune     /* on */
2738e6d0c30SPeter Brune     lsparse[i] = 0;
274e5a0faa4SPeter Brune     gsparse[i] = 0;
2758e6d0c30SPeter Brune     if (lcid[i] >= 0) {
2768e6d0c30SPeter Brune       lsparse[i] = 1;
2778e6d0c30SPeter Brune       gsparse[i] = 0;
2788e6d0c30SPeter Brune     } else {
2799566063dSJacob Faibussowitsch       PetscCall(MatGetRow(lA, i, &ncols, &rcol, &rval));
2808e6d0c30SPeter Brune       for (j = 0; j < ncols; j++) {
2811ce39c63SPeter Brune         col = rcol[j];
282ad540459SPierre 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;
2838e6d0c30SPeter Brune       }
2849566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(lA, i, &ncols, &rcol, &rval));
2858e6d0c30SPeter Brune       /* off */
2861ce39c63SPeter Brune       if (gA) {
2879566063dSJacob Faibussowitsch         PetscCall(MatGetRow(gA, i, &ncols, &rcol, &rval));
2888e6d0c30SPeter Brune         for (j = 0; j < ncols; j++) {
2891ce39c63SPeter Brune           col = rcol[j];
290ad540459SPierre 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;
2918e6d0c30SPeter Brune         }
2929566063dSJacob Faibussowitsch         PetscCall(MatRestoreRow(gA, i, &ncols, &rcol, &rval));
2938e6d0c30SPeter Brune       }
2948e6d0c30SPeter Brune     }
2951ce39c63SPeter Brune   }
2968e6d0c30SPeter Brune 
2978e6d0c30SPeter Brune   /* preallocate and create the prolongator */
2989566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), P));
2998926f930SMark Adams   PetscCall(MatGetType(A, &mtype));
3009566063dSJacob Faibussowitsch   PetscCall(MatSetType(*P, mtype));
3019566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*P, fn, cn, PETSC_DETERMINE, PETSC_DETERMINE));
3029566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(*P, 0, lsparse, 0, gsparse));
3039566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(*P, 0, lsparse));
3048e6d0c30SPeter Brune 
3058e6d0c30SPeter Brune   /* loop over local fine nodes -- get the diagonal, the sum of positive and negative strong and weak weights, and set up the row */
3068e6d0c30SPeter Brune   for (i = 0; i < fn; i++) {
3078e6d0c30SPeter Brune     /* determine on or off */
3088e6d0c30SPeter Brune     row_f = i + fs;
3098e6d0c30SPeter Brune     row_c = lcid[i];
3108e6d0c30SPeter Brune     if (row_c >= 0) {
3118e6d0c30SPeter Brune       pij = 1.;
3129566063dSJacob Faibussowitsch       PetscCall(MatSetValues(*P, 1, &row_f, 1, &row_c, &pij, INSERT_VALUES));
3138e6d0c30SPeter Brune     } else {
3148e6d0c30SPeter Brune       g_pos = 0.;
3158e6d0c30SPeter Brune       g_neg = 0.;
3168e6d0c30SPeter Brune       a_pos = 0.;
3178e6d0c30SPeter Brune       a_neg = 0.;
3188e6d0c30SPeter Brune       diag  = 0.;
3198e6d0c30SPeter Brune 
3201ce39c63SPeter Brune       /* local connections */
3219566063dSJacob Faibussowitsch       PetscCall(MatGetRow(lA, i, &ncols, &rcol, &rval));
3221ce39c63SPeter Brune       for (j = 0; j < ncols; j++) {
3231ce39c63SPeter Brune         col = rcol[j];
324c1eae691SMark Adams         if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0] * Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0] * Amax_neg[i])) {
3251ce39c63SPeter Brune           if (PetscRealPart(rval[j]) > 0.) {
3261ce39c63SPeter Brune             g_pos += rval[j];
3278e6d0c30SPeter Brune           } else {
3281ce39c63SPeter Brune             g_neg += rval[j];
3298e6d0c30SPeter Brune           }
3301ce39c63SPeter Brune         }
3311ce39c63SPeter Brune         if (col != i) {
3321ce39c63SPeter Brune           if (PetscRealPart(rval[j]) > 0.) {
3331ce39c63SPeter Brune             a_pos += rval[j];
3341ce39c63SPeter Brune           } else {
3351ce39c63SPeter Brune             a_neg += rval[j];
3361ce39c63SPeter Brune           }
3371ce39c63SPeter Brune         } else {
3381ce39c63SPeter Brune           diag = rval[j];
3391ce39c63SPeter Brune         }
3408e6d0c30SPeter Brune       }
3419566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(lA, i, &ncols, &rcol, &rval));
3428e6d0c30SPeter Brune 
3431ce39c63SPeter Brune       /* ghosted connections */
3448e6d0c30SPeter Brune       if (gA) {
3459566063dSJacob Faibussowitsch         PetscCall(MatGetRow(gA, i, &ncols, &rcol, &rval));
3461ce39c63SPeter Brune         for (j = 0; j < ncols; j++) {
3471ce39c63SPeter Brune           col = rcol[j];
348c1eae691SMark Adams           if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0] * Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0] * Amax_neg[i])) {
3491ce39c63SPeter Brune             if (PetscRealPart(rval[j]) > 0.) {
3501ce39c63SPeter Brune               g_pos += rval[j];
3518e6d0c30SPeter Brune             } else {
3521ce39c63SPeter Brune               g_neg += rval[j];
3538e6d0c30SPeter Brune             }
3541ce39c63SPeter Brune           }
3551ce39c63SPeter Brune           if (PetscRealPart(rval[j]) > 0.) {
3561ce39c63SPeter Brune             a_pos += rval[j];
3571ce39c63SPeter Brune           } else {
3581ce39c63SPeter Brune             a_neg += rval[j];
3591ce39c63SPeter Brune           }
3608e6d0c30SPeter Brune         }
3619566063dSJacob Faibussowitsch         PetscCall(MatRestoreRow(gA, i, &ncols, &rcol, &rval));
3628e6d0c30SPeter Brune       }
3638e6d0c30SPeter Brune 
3648e6d0c30SPeter Brune       if (g_neg == 0.) {
3658e6d0c30SPeter Brune         alpha = 0.;
3668e6d0c30SPeter Brune       } else {
3678e6d0c30SPeter Brune         alpha = -a_neg / g_neg;
3688e6d0c30SPeter Brune       }
3698e6d0c30SPeter Brune 
3708e6d0c30SPeter Brune       if (g_pos == 0.) {
3718e6d0c30SPeter Brune         diag += a_pos;
3728e6d0c30SPeter Brune         beta = 0.;
3738e6d0c30SPeter Brune       } else {
3748e6d0c30SPeter Brune         beta = -a_pos / g_pos;
3758e6d0c30SPeter Brune       }
376e5a0faa4SPeter Brune       if (diag == 0.) {
377e5a0faa4SPeter Brune         invdiag = 0.;
378e5a0faa4SPeter Brune       } else invdiag = 1. / diag;
3798e6d0c30SPeter Brune       /* on */
3809566063dSJacob Faibussowitsch       PetscCall(MatGetRow(lA, i, &ncols, &rcol, &rval));
3818e6d0c30SPeter Brune       idx = 0;
3828e6d0c30SPeter Brune       for (j = 0; j < ncols; j++) {
3831ce39c63SPeter Brune         col = rcol[j];
384c1eae691SMark Adams         if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0] * Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0] * Amax_neg[i])) {
3858e6d0c30SPeter Brune           row_f = i + fs;
3868e6d0c30SPeter Brune           row_c = lcid[col];
3878e6d0c30SPeter Brune           /* set the values for on-processor ones */
3881ce39c63SPeter Brune           if (PetscRealPart(rval[j]) < 0.) {
3891ce39c63SPeter Brune             pij = rval[j] * alpha * invdiag;
3908e6d0c30SPeter Brune           } else {
3911ce39c63SPeter Brune             pij = rval[j] * beta * invdiag;
3928e6d0c30SPeter Brune           }
3938e6d0c30SPeter Brune           if (PetscAbsScalar(pij) != 0.) {
3948e6d0c30SPeter Brune             pvals[idx] = pij;
3958e6d0c30SPeter Brune             pcols[idx] = row_c;
3968e6d0c30SPeter Brune             idx++;
3978e6d0c30SPeter Brune           }
3988e6d0c30SPeter Brune         }
3998e6d0c30SPeter Brune       }
4009566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(lA, i, &ncols, &rcol, &rval));
4018e6d0c30SPeter Brune       /* off */
4021ce39c63SPeter Brune       if (gA) {
4039566063dSJacob Faibussowitsch         PetscCall(MatGetRow(gA, i, &ncols, &rcol, &rval));
4048e6d0c30SPeter Brune         for (j = 0; j < ncols; j++) {
4051ce39c63SPeter Brune           col = rcol[j];
406c1eae691SMark Adams           if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0] * Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0] * Amax_neg[i])) {
4078e6d0c30SPeter Brune             row_f = i + fs;
4088e6d0c30SPeter Brune             row_c = gcid[col];
4098e6d0c30SPeter Brune             /* set the values for on-processor ones */
4101ce39c63SPeter Brune             if (PetscRealPart(rval[j]) < 0.) {
4111ce39c63SPeter Brune               pij = rval[j] * alpha * invdiag;
4128e6d0c30SPeter Brune             } else {
4131ce39c63SPeter Brune               pij = rval[j] * beta * invdiag;
4148e6d0c30SPeter Brune             }
4158e6d0c30SPeter Brune             if (PetscAbsScalar(pij) != 0.) {
4168e6d0c30SPeter Brune               pvals[idx] = pij;
4178e6d0c30SPeter Brune               pcols[idx] = row_c;
4188e6d0c30SPeter Brune               idx++;
4198e6d0c30SPeter Brune             }
4208e6d0c30SPeter Brune           }
4218e6d0c30SPeter Brune         }
4229566063dSJacob Faibussowitsch         PetscCall(MatRestoreRow(gA, i, &ncols, &rcol, &rval));
4233c9ab2c3SPeter Brune       }
4249566063dSJacob Faibussowitsch       PetscCall(MatSetValues(*P, 1, &row_f, idx, pcols, pvals, INSERT_VALUES));
4258e6d0c30SPeter Brune     }
4268e6d0c30SPeter Brune   }
4273c9ab2c3SPeter Brune 
4289566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY));
4299566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY));
4308e6d0c30SPeter Brune 
4319566063dSJacob Faibussowitsch   PetscCall(PetscFree5(lsparse, gsparse, lcid, Amax_pos, Amax_neg));
432e632b94dSBarry Smith 
4339566063dSJacob Faibussowitsch   PetscCall(PetscFree2(pcols, pvals));
43487b9b13cSPeter Brune   if (gA) {
4359566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&sf));
4369566063dSJacob Faibussowitsch     PetscCall(PetscFree(gcid));
43787b9b13cSPeter Brune   }
4383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4398e6d0c30SPeter Brune }
4408e6d0c30SPeter Brune 
PCGAMGTruncateProlongator_Private(PC pc,Mat * P)441ff6a9541SJacob Faibussowitsch static PetscErrorCode PCGAMGTruncateProlongator_Private(PC pc, Mat *P)
442d71ae5a4SJacob Faibussowitsch {
443586a8384SPeter Brune   PetscInt           j, i, ps, pf, pn, pcs, pcf, pcn, idx, cmax;
444586a8384SPeter Brune   const PetscScalar *pval;
445586a8384SPeter Brune   const PetscInt    *pcol;
446586a8384SPeter Brune   PetscScalar       *pnval;
447586a8384SPeter Brune   PetscInt          *pncol;
448586a8384SPeter Brune   PetscInt           ncols;
449586a8384SPeter Brune   Mat                Pnew;
450586a8384SPeter Brune   PetscInt          *lsparse, *gsparse;
451586a8384SPeter Brune   PetscReal          pmax_pos, pmax_neg, ptot_pos, ptot_neg, pthresh_pos, pthresh_neg;
452586a8384SPeter Brune   PC_MG             *mg      = (PC_MG *)pc->data;
453586a8384SPeter Brune   PC_GAMG           *pc_gamg = (PC_GAMG *)mg->innerctx;
454586a8384SPeter Brune   PC_GAMG_Classical *cls     = (PC_GAMG_Classical *)pc_gamg->subctx;
455d9558ea9SBarry Smith   MatType            mtype;
456586a8384SPeter Brune 
457586a8384SPeter Brune   PetscFunctionBegin;
458586a8384SPeter Brune   /* trim and rescale with reallocation */
4599566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(*P, &ps, &pf));
4609566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRangeColumn(*P, &pcs, &pcf));
461586a8384SPeter Brune   pn  = pf - ps;
462586a8384SPeter Brune   pcn = pcf - pcs;
4639566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(pn, &lsparse, pn, &gsparse));
464586a8384SPeter Brune   /* allocate */
465586a8384SPeter Brune   cmax = 0;
466586a8384SPeter Brune   for (i = ps; i < pf; i++) {
467b4dc3ebdSPeter Brune     lsparse[i - ps] = 0;
468b4dc3ebdSPeter Brune     gsparse[i - ps] = 0;
4699566063dSJacob Faibussowitsch     PetscCall(MatGetRow(*P, i, &ncols, &pcol, &pval));
470ad540459SPierre Jolivet     if (ncols > cmax) cmax = ncols;
471586a8384SPeter Brune     pmax_pos = 0.;
472586a8384SPeter Brune     pmax_neg = 0.;
473586a8384SPeter Brune     for (j = 0; j < ncols; j++) {
474586a8384SPeter Brune       if (PetscRealPart(pval[j]) > pmax_pos) {
475586a8384SPeter Brune         pmax_pos = PetscRealPart(pval[j]);
476586a8384SPeter Brune       } else if (PetscRealPart(pval[j]) < pmax_neg) {
477586a8384SPeter Brune         pmax_neg = PetscRealPart(pval[j]);
478586a8384SPeter Brune       }
479586a8384SPeter Brune     }
480586a8384SPeter Brune     for (j = 0; j < ncols; j++) {
4818eab0cc1SPeter Brune       if (PetscRealPart(pval[j]) >= pmax_pos * cls->interp_threshold || PetscRealPart(pval[j]) <= pmax_neg * cls->interp_threshold) {
482b4dc3ebdSPeter Brune         if (pcol[j] >= pcs && pcol[j] < pcf) {
483b4dc3ebdSPeter Brune           lsparse[i - ps]++;
484586a8384SPeter Brune         } else {
485b4dc3ebdSPeter Brune           gsparse[i - ps]++;
486586a8384SPeter Brune         }
487586a8384SPeter Brune       }
488586a8384SPeter Brune     }
4899566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(*P, i, &ncols, &pcol, &pval));
490586a8384SPeter Brune   }
491586a8384SPeter Brune 
4929566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(cmax, &pnval, cmax, &pncol));
493586a8384SPeter Brune 
4949566063dSJacob Faibussowitsch   PetscCall(MatGetType(*P, &mtype));
4959566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)*P), &Pnew));
4969566063dSJacob Faibussowitsch   PetscCall(MatSetType(Pnew, mtype));
4979566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Pnew, pn, pcn, PETSC_DETERMINE, PETSC_DETERMINE));
4989566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(Pnew, 0, lsparse));
4999566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(Pnew, 0, lsparse, 0, gsparse));
500586a8384SPeter Brune 
501586a8384SPeter Brune   for (i = ps; i < pf; i++) {
5029566063dSJacob Faibussowitsch     PetscCall(MatGetRow(*P, i, &ncols, &pcol, &pval));
503586a8384SPeter Brune     pmax_pos = 0.;
504586a8384SPeter Brune     pmax_neg = 0.;
505586a8384SPeter Brune     for (j = 0; j < ncols; j++) {
506586a8384SPeter Brune       if (PetscRealPart(pval[j]) > pmax_pos) {
507586a8384SPeter Brune         pmax_pos = PetscRealPart(pval[j]);
508586a8384SPeter Brune       } else if (PetscRealPart(pval[j]) < pmax_neg) {
509586a8384SPeter Brune         pmax_neg = PetscRealPart(pval[j]);
510586a8384SPeter Brune       }
511586a8384SPeter Brune     }
512586a8384SPeter Brune     pthresh_pos = 0.;
513586a8384SPeter Brune     pthresh_neg = 0.;
514586a8384SPeter Brune     ptot_pos    = 0.;
515586a8384SPeter Brune     ptot_neg    = 0.;
516586a8384SPeter Brune     for (j = 0; j < ncols; j++) {
5178eab0cc1SPeter Brune       if (PetscRealPart(pval[j]) >= cls->interp_threshold * pmax_pos) {
518586a8384SPeter Brune         pthresh_pos += PetscRealPart(pval[j]);
5198eab0cc1SPeter Brune       } else if (PetscRealPart(pval[j]) <= cls->interp_threshold * pmax_neg) {
520586a8384SPeter Brune         pthresh_neg += PetscRealPart(pval[j]);
521586a8384SPeter Brune       }
522586a8384SPeter Brune       if (PetscRealPart(pval[j]) > 0.) {
523586a8384SPeter Brune         ptot_pos += PetscRealPart(pval[j]);
524586a8384SPeter Brune       } else {
525586a8384SPeter Brune         ptot_neg += PetscRealPart(pval[j]);
526586a8384SPeter Brune       }
527586a8384SPeter Brune     }
5286bd8ea90SPeter Brune     if (PetscAbsReal(pthresh_pos) > 0.) ptot_pos /= pthresh_pos;
5296bd8ea90SPeter Brune     if (PetscAbsReal(pthresh_neg) > 0.) ptot_neg /= pthresh_neg;
530586a8384SPeter Brune     idx = 0;
531586a8384SPeter Brune     for (j = 0; j < ncols; j++) {
5328eab0cc1SPeter Brune       if (PetscRealPart(pval[j]) >= pmax_pos * cls->interp_threshold) {
533586a8384SPeter Brune         pnval[idx] = ptot_pos * pval[j];
534586a8384SPeter Brune         pncol[idx] = pcol[j];
535586a8384SPeter Brune         idx++;
5368eab0cc1SPeter Brune       } else if (PetscRealPart(pval[j]) <= pmax_neg * cls->interp_threshold) {
537586a8384SPeter Brune         pnval[idx] = ptot_neg * pval[j];
538586a8384SPeter Brune         pncol[idx] = pcol[j];
539586a8384SPeter Brune         idx++;
540586a8384SPeter Brune       }
541586a8384SPeter Brune     }
5429566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(*P, i, &ncols, &pcol, &pval));
5439566063dSJacob Faibussowitsch     PetscCall(MatSetValues(Pnew, 1, &i, idx, pncol, pnval, INSERT_VALUES));
544586a8384SPeter Brune   }
545586a8384SPeter Brune 
5469566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Pnew, MAT_FINAL_ASSEMBLY));
5479566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Pnew, MAT_FINAL_ASSEMBLY));
5489566063dSJacob Faibussowitsch   PetscCall(MatDestroy(P));
549586a8384SPeter Brune 
550586a8384SPeter Brune   *P = Pnew;
5519566063dSJacob Faibussowitsch   PetscCall(PetscFree2(lsparse, gsparse));
5529566063dSJacob Faibussowitsch   PetscCall(PetscFree2(pnval, pncol));
5533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
554586a8384SPeter Brune }
555586a8384SPeter Brune 
PCGAMGProlongator_Classical_Standard(PC pc,Mat A,PetscCoarsenData * agg_lists,Mat * P)5568926f930SMark Adams static PetscErrorCode PCGAMGProlongator_Classical_Standard(PC pc, Mat A, PetscCoarsenData *agg_lists, Mat *P)
557d71ae5a4SJacob Faibussowitsch {
558c634539dSPeter Brune   Mat                lA, *lAs;
559f9a65ec8SPeter Brune   MatType            mtype;
56063b0c588SPeter Brune   Vec                cv;
56163b0c588SPeter Brune   PetscInt          *gcid, *lcid, *lsparse, *gsparse, *picol;
562420cd23fSPeter Brune   PetscInt           fs, fe, cs, ce, nl, i, j, k, li, lni, ci, ncols, maxcols, fn, cn, cid;
563420cd23fSPeter Brune   PetscMPIInt        size;
56463b0c588SPeter Brune   const PetscInt    *lidx, *icol, *gidx;
56563b0c588SPeter Brune   PetscBool          iscoarse;
56663b0c588SPeter Brune   PetscScalar        vi, pentry, pjentry;
56763b0c588SPeter Brune   PetscScalar       *pcontrib, *pvcol;
56863b0c588SPeter Brune   const PetscScalar *vcol;
569ed4e82eaSPeter Brune   PetscReal          diag, jdiag, jwttotal;
570f9a65ec8SPeter Brune   PetscInt           pncols;
571c634539dSPeter Brune   PetscSF            sf;
572c634539dSPeter Brune   PetscLayout        clayout;
57363b0c588SPeter Brune   IS                 lis;
574f9a65ec8SPeter Brune 
575f9a65ec8SPeter Brune   PetscFunctionBegin;
5769566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
5779566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &fs, &fe));
578f9a65ec8SPeter Brune   fn = fe - fs;
5799566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF, fe - fs, fs, 1, &lis));
580c634539dSPeter Brune   if (size > 1) {
5819566063dSJacob Faibussowitsch     PetscCall(MatGetLayouts(A, NULL, &clayout));
582f9a65ec8SPeter Brune     /* increase the overlap by two to get neighbors of neighbors */
5839566063dSJacob Faibussowitsch     PetscCall(MatIncreaseOverlap(A, 1, &lis, 2));
5849566063dSJacob Faibussowitsch     PetscCall(ISSort(lis));
585f9a65ec8SPeter Brune     /* get the local part of A */
5869566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatrices(A, 1, &lis, &lis, MAT_INITIAL_MATRIX, &lAs));
587c634539dSPeter Brune     lA = lAs[0];
588c634539dSPeter Brune     /* build an SF out of it */
5899566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(lis, &nl));
5909566063dSJacob Faibussowitsch     PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)A), &sf));
5919566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(lis, &lidx));
5929566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraphLayout(sf, clayout, nl, NULL, PETSC_COPY_VALUES, lidx));
5939566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(lis, &lidx));
594c634539dSPeter Brune   } else {
595c634539dSPeter Brune     lA = A;
596c634539dSPeter Brune     nl = fn;
597c634539dSPeter Brune   }
598c634539dSPeter Brune   /* create a communication structure for the overlapped portion and transmit coarse indices */
5999566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(fn, &lsparse, fn, &gsparse, nl, &pcontrib));
600f9a65ec8SPeter Brune   /* create coarse vector */
601f9a65ec8SPeter Brune   cn = 0;
602f9a65ec8SPeter Brune   for (i = 0; i < fn; i++) {
6035e62d202SMark Adams     PetscCall(PetscCDIsEmptyAt(agg_lists, i, &iscoarse));
604ad540459SPierre Jolivet     if (!iscoarse) cn++;
605f9a65ec8SPeter Brune   }
6069566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(fn, &gcid));
6079566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)A), cn, PETSC_DECIDE, &cv));
6089566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(cv, &cs, &ce));
609f9a65ec8SPeter Brune   cn = 0;
610f9a65ec8SPeter Brune   for (i = 0; i < fn; i++) {
6115e62d202SMark Adams     PetscCall(PetscCDIsEmptyAt(agg_lists, i, &iscoarse));
612f9a65ec8SPeter Brune     if (!iscoarse) {
613c634539dSPeter Brune       gcid[i] = cs + cn;
614f9a65ec8SPeter Brune       cn++;
615f9a65ec8SPeter Brune     } else {
616c634539dSPeter Brune       gcid[i] = -1;
617f9a65ec8SPeter Brune     }
618f9a65ec8SPeter Brune   }
619c634539dSPeter Brune   if (size > 1) {
6209566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nl, &lcid));
6219566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(sf, MPIU_INT, gcid, lcid, MPI_REPLACE));
6229566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(sf, MPIU_INT, gcid, lcid, MPI_REPLACE));
623c634539dSPeter Brune   } else {
624c634539dSPeter Brune     lcid = gcid;
625c634539dSPeter Brune   }
626f9a65ec8SPeter Brune   /* count to preallocate the prolongator */
6279566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(lis, &gidx));
628f9a65ec8SPeter Brune   maxcols = 0;
629f9a65ec8SPeter Brune   /* count the number of unique contributing coarse cells for each fine */
630f9a65ec8SPeter Brune   for (i = 0; i < nl; i++) {
631ed4e82eaSPeter Brune     pcontrib[i] = 0.;
6329566063dSJacob Faibussowitsch     PetscCall(MatGetRow(lA, i, &ncols, &icol, NULL));
633f9a65ec8SPeter Brune     if (gidx[i] >= fs && gidx[i] < fe) {
634f9a65ec8SPeter Brune       li          = gidx[i] - fs;
635f9a65ec8SPeter Brune       lsparse[li] = 0;
636f9a65ec8SPeter Brune       gsparse[li] = 0;
637c634539dSPeter Brune       cid         = lcid[i];
638f9a65ec8SPeter Brune       if (cid >= 0) {
639f9a65ec8SPeter Brune         lsparse[li] = 1;
640f9a65ec8SPeter Brune       } else {
641f9a65ec8SPeter Brune         for (j = 0; j < ncols; j++) {
642c634539dSPeter Brune           if (lcid[icol[j]] >= 0) {
643f9a65ec8SPeter Brune             pcontrib[icol[j]] = 1.;
644f9a65ec8SPeter Brune           } else {
645f9a65ec8SPeter Brune             ci = icol[j];
6469566063dSJacob Faibussowitsch             PetscCall(MatRestoreRow(lA, i, &ncols, &icol, NULL));
6479566063dSJacob Faibussowitsch             PetscCall(MatGetRow(lA, ci, &ncols, &icol, NULL));
648f9a65ec8SPeter Brune             for (k = 0; k < ncols; k++) {
649ad540459SPierre Jolivet               if (lcid[icol[k]] >= 0) pcontrib[icol[k]] = 1.;
650f9a65ec8SPeter Brune             }
6519566063dSJacob Faibussowitsch             PetscCall(MatRestoreRow(lA, ci, &ncols, &icol, NULL));
6529566063dSJacob Faibussowitsch             PetscCall(MatGetRow(lA, i, &ncols, &icol, NULL));
653f9a65ec8SPeter Brune           }
654f9a65ec8SPeter Brune         }
655f9a65ec8SPeter Brune         for (j = 0; j < ncols; j++) {
656c634539dSPeter Brune           if (lcid[icol[j]] >= 0 && pcontrib[icol[j]] != 0.) {
657c634539dSPeter Brune             lni = lcid[icol[j]];
658f9a65ec8SPeter Brune             if (lni >= cs && lni < ce) {
659f9a65ec8SPeter Brune               lsparse[li]++;
660f9a65ec8SPeter Brune             } else {
661f9a65ec8SPeter Brune               gsparse[li]++;
662f9a65ec8SPeter Brune             }
663f9a65ec8SPeter Brune             pcontrib[icol[j]] = 0.;
664f9a65ec8SPeter Brune           } else {
665f9a65ec8SPeter Brune             ci = icol[j];
6669566063dSJacob Faibussowitsch             PetscCall(MatRestoreRow(lA, i, &ncols, &icol, NULL));
6679566063dSJacob Faibussowitsch             PetscCall(MatGetRow(lA, ci, &ncols, &icol, NULL));
668f9a65ec8SPeter Brune             for (k = 0; k < ncols; k++) {
669c634539dSPeter Brune               if (lcid[icol[k]] >= 0 && pcontrib[icol[k]] != 0.) {
670c634539dSPeter Brune                 lni = lcid[icol[k]];
671f9a65ec8SPeter Brune                 if (lni >= cs && lni < ce) {
672f9a65ec8SPeter Brune                   lsparse[li]++;
673f9a65ec8SPeter Brune                 } else {
674f9a65ec8SPeter Brune                   gsparse[li]++;
675f9a65ec8SPeter Brune                 }
676f9a65ec8SPeter Brune                 pcontrib[icol[k]] = 0.;
677f9a65ec8SPeter Brune               }
678f9a65ec8SPeter Brune             }
6799566063dSJacob Faibussowitsch             PetscCall(MatRestoreRow(lA, ci, &ncols, &icol, NULL));
6809566063dSJacob Faibussowitsch             PetscCall(MatGetRow(lA, i, &ncols, &icol, NULL));
681f9a65ec8SPeter Brune           }
682f9a65ec8SPeter Brune         }
683ed4e82eaSPeter Brune       }
684f9a65ec8SPeter Brune       if (lsparse[li] + gsparse[li] > maxcols) maxcols = lsparse[li] + gsparse[li];
685ed4e82eaSPeter Brune     }
6869566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(lA, i, &ncols, &icol, &vcol));
687f9a65ec8SPeter Brune   }
6889566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(maxcols, &picol, maxcols, &pvcol));
6899566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), P));
6909566063dSJacob Faibussowitsch   PetscCall(MatGetType(A, &mtype));
6919566063dSJacob Faibussowitsch   PetscCall(MatSetType(*P, mtype));
6929566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*P, fn, cn, PETSC_DETERMINE, PETSC_DETERMINE));
6939566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(*P, 0, lsparse, 0, gsparse));
6949566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(*P, 0, lsparse));
695f9a65ec8SPeter Brune   for (i = 0; i < nl; i++) {
696ed4e82eaSPeter Brune     diag = 0.;
697f9a65ec8SPeter Brune     if (gidx[i] >= fs && gidx[i] < fe) {
698f9a65ec8SPeter Brune       pncols = 0;
699c634539dSPeter Brune       cid    = lcid[i];
700f9a65ec8SPeter Brune       if (cid >= 0) {
701f9a65ec8SPeter Brune         pncols   = 1;
702f9a65ec8SPeter Brune         picol[0] = cid;
703f9a65ec8SPeter Brune         pvcol[0] = 1.;
704f9a65ec8SPeter Brune       } else {
7059566063dSJacob Faibussowitsch         PetscCall(MatGetRow(lA, i, &ncols, &icol, &vcol));
706f9a65ec8SPeter Brune         for (j = 0; j < ncols; j++) {
707ed4e82eaSPeter Brune           pentry = vcol[j];
708c634539dSPeter Brune           if (lcid[icol[j]] >= 0) {
709f9a65ec8SPeter Brune             /* coarse neighbor */
710ed4e82eaSPeter Brune             pcontrib[icol[j]] += pentry;
711ed4e82eaSPeter Brune           } else if (icol[j] != i) {
712f9a65ec8SPeter Brune             /* the neighbor is a strongly connected fine node */
713f9a65ec8SPeter Brune             ci = icol[j];
714f9a65ec8SPeter Brune             vi = vcol[j];
7159566063dSJacob Faibussowitsch             PetscCall(MatRestoreRow(lA, i, &ncols, &icol, &vcol));
7169566063dSJacob Faibussowitsch             PetscCall(MatGetRow(lA, ci, &ncols, &icol, &vcol));
717ed4e82eaSPeter Brune             jwttotal = 0.;
718f9a65ec8SPeter Brune             jdiag    = 0.;
719f9a65ec8SPeter Brune             for (k = 0; k < ncols; k++) {
720ad540459SPierre Jolivet               if (ci == icol[k]) jdiag = PetscRealPart(vcol[k]);
721f9a65ec8SPeter Brune             }
722f9a65ec8SPeter Brune             for (k = 0; k < ncols; k++) {
723c634539dSPeter Brune               if (lcid[icol[k]] >= 0 && jdiag * PetscRealPart(vcol[k]) < 0.) {
724ed4e82eaSPeter Brune                 pjentry = vcol[k];
7257779008dSPeter Brune                 jwttotal += PetscRealPart(pjentry);
726f9a65ec8SPeter Brune               }
727f9a65ec8SPeter Brune             }
728ed4e82eaSPeter Brune             if (jwttotal != 0.) {
7297779008dSPeter Brune               jwttotal = PetscRealPart(vi) / jwttotal;
730ed4e82eaSPeter Brune               for (k = 0; k < ncols; k++) {
731c634539dSPeter Brune                 if (lcid[icol[k]] >= 0 && jdiag * PetscRealPart(vcol[k]) < 0.) {
732586a8384SPeter Brune                   pjentry = vcol[k] * jwttotal;
733ed4e82eaSPeter Brune                   pcontrib[icol[k]] += pjentry;
734ed4e82eaSPeter Brune                 }
735ed4e82eaSPeter Brune               }
736ed4e82eaSPeter Brune             } else {
737ed4e82eaSPeter Brune               diag += PetscRealPart(vi);
738ed4e82eaSPeter Brune             }
7399566063dSJacob Faibussowitsch             PetscCall(MatRestoreRow(lA, ci, &ncols, &icol, &vcol));
7409566063dSJacob Faibussowitsch             PetscCall(MatGetRow(lA, i, &ncols, &icol, &vcol));
741ed4e82eaSPeter Brune           } else {
742ed4e82eaSPeter Brune             diag += PetscRealPart(vcol[j]);
743f9a65ec8SPeter Brune           }
744f9a65ec8SPeter Brune         }
745586a8384SPeter Brune         if (diag != 0.) {
746586a8384SPeter Brune           diag = 1. / diag;
747f9a65ec8SPeter Brune           for (j = 0; j < ncols; j++) {
748c634539dSPeter Brune             if (lcid[icol[j]] >= 0 && pcontrib[icol[j]] != 0.) {
749f9a65ec8SPeter Brune               /* the neighbor is a coarse node */
750ed4e82eaSPeter Brune               if (PetscAbsScalar(pcontrib[icol[j]]) > 0.0) {
751c634539dSPeter Brune                 lni           = lcid[icol[j]];
752586a8384SPeter Brune                 pvcol[pncols] = -pcontrib[icol[j]] * diag;
753f9a65ec8SPeter Brune                 picol[pncols] = lni;
754f9a65ec8SPeter Brune                 pncols++;
755ed4e82eaSPeter Brune               }
756ed4e82eaSPeter Brune               pcontrib[icol[j]] = 0.;
757f9a65ec8SPeter Brune             } else {
758f9a65ec8SPeter Brune               /* the neighbor is a strongly connected fine node */
759f9a65ec8SPeter Brune               ci = icol[j];
7609566063dSJacob Faibussowitsch               PetscCall(MatRestoreRow(lA, i, &ncols, &icol, &vcol));
7619566063dSJacob Faibussowitsch               PetscCall(MatGetRow(lA, ci, &ncols, &icol, &vcol));
762f9a65ec8SPeter Brune               for (k = 0; k < ncols; k++) {
763c634539dSPeter Brune                 if (lcid[icol[k]] >= 0 && pcontrib[icol[k]] != 0.) {
764ed4e82eaSPeter Brune                   if (PetscAbsScalar(pcontrib[icol[k]]) > 0.0) {
765c634539dSPeter Brune                     lni           = lcid[icol[k]];
766586a8384SPeter Brune                     pvcol[pncols] = -pcontrib[icol[k]] * diag;
767f9a65ec8SPeter Brune                     picol[pncols] = lni;
768f9a65ec8SPeter Brune                     pncols++;
769f9a65ec8SPeter Brune                   }
770ed4e82eaSPeter Brune                   pcontrib[icol[k]] = 0.;
771ed4e82eaSPeter Brune                 }
772f9a65ec8SPeter Brune               }
7739566063dSJacob Faibussowitsch               PetscCall(MatRestoreRow(lA, ci, &ncols, &icol, &vcol));
7749566063dSJacob Faibussowitsch               PetscCall(MatGetRow(lA, i, &ncols, &icol, &vcol));
775f9a65ec8SPeter Brune             }
776ed4e82eaSPeter Brune             pcontrib[icol[j]] = 0.;
777f9a65ec8SPeter Brune           }
7789566063dSJacob Faibussowitsch           PetscCall(MatRestoreRow(lA, i, &ncols, &icol, &vcol));
779f9a65ec8SPeter Brune         }
780586a8384SPeter Brune       }
781f9a65ec8SPeter Brune       ci = gidx[i];
78248a46eb9SPierre Jolivet       if (pncols > 0) PetscCall(MatSetValues(*P, 1, &ci, pncols, picol, pvcol, INSERT_VALUES));
783f9a65ec8SPeter Brune     }
784f9a65ec8SPeter Brune   }
7859566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(lis, &gidx));
7869566063dSJacob Faibussowitsch   PetscCall(PetscFree2(picol, pvcol));
7879566063dSJacob Faibussowitsch   PetscCall(PetscFree3(lsparse, gsparse, pcontrib));
7889566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&lis));
7899566063dSJacob Faibussowitsch   PetscCall(PetscFree(gcid));
790c634539dSPeter Brune   if (size > 1) {
7919566063dSJacob Faibussowitsch     PetscCall(PetscFree(lcid));
7929566063dSJacob Faibussowitsch     PetscCall(MatDestroyMatrices(1, &lAs));
7939566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&sf));
794c634539dSPeter Brune   }
7959566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&cv));
7969566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY));
7979566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY));
7983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7998eab0cc1SPeter Brune }
800f9a65ec8SPeter Brune 
PCGAMGOptProlongator_Classical_Jacobi(PC pc,Mat A,Mat * P)801ff6a9541SJacob Faibussowitsch static PetscErrorCode PCGAMGOptProlongator_Classical_Jacobi(PC pc, Mat A, Mat *P)
802d71ae5a4SJacob Faibussowitsch {
803b4dc3ebdSPeter Brune   PetscInt           f, s, n, cf, cs, i, idx;
804b4dc3ebdSPeter Brune   PetscInt          *coarserows;
805b4dc3ebdSPeter Brune   PetscInt           ncols;
806b4dc3ebdSPeter Brune   const PetscInt    *pcols;
807b4dc3ebdSPeter Brune   const PetscScalar *pvals;
808b4dc3ebdSPeter Brune   Mat                Pnew;
809b4dc3ebdSPeter Brune   Vec                diag;
810b4dc3ebdSPeter Brune   PC_MG             *mg      = (PC_MG *)pc->data;
811b4dc3ebdSPeter Brune   PC_GAMG           *pc_gamg = (PC_GAMG *)mg->innerctx;
812b4dc3ebdSPeter Brune   PC_GAMG_Classical *cls     = (PC_GAMG_Classical *)pc_gamg->subctx;
813b4dc3ebdSPeter Brune 
814b4dc3ebdSPeter Brune   PetscFunctionBegin;
815b4dc3ebdSPeter Brune   if (cls->nsmooths == 0) {
8169566063dSJacob Faibussowitsch     PetscCall(PCGAMGTruncateProlongator_Private(pc, P));
8173ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
818b4dc3ebdSPeter Brune   }
8199566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(*P, &s, &f));
820b4dc3ebdSPeter Brune   n = f - s;
8219566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRangeColumn(*P, &cs, &cf));
8229566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &coarserows));
823b4dc3ebdSPeter Brune   /* identify the rows corresponding to coarse unknowns */
824b4dc3ebdSPeter Brune   idx = 0;
825b4dc3ebdSPeter Brune   for (i = s; i < f; i++) {
8269566063dSJacob Faibussowitsch     PetscCall(MatGetRow(*P, i, &ncols, &pcols, &pvals));
827b4dc3ebdSPeter Brune     /* assume, for now, that it's a coarse unknown if it has a single unit entry */
828b4dc3ebdSPeter Brune     if (ncols == 1) {
829b4dc3ebdSPeter Brune       if (pvals[0] == 1.) {
830b4dc3ebdSPeter Brune         coarserows[idx] = i;
831b4dc3ebdSPeter Brune         idx++;
832b4dc3ebdSPeter Brune       }
833b4dc3ebdSPeter Brune     }
8349566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(*P, i, &ncols, &pcols, &pvals));
835b4dc3ebdSPeter Brune   }
8369566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &diag, NULL));
8379566063dSJacob Faibussowitsch   PetscCall(MatGetDiagonal(A, diag));
8389566063dSJacob Faibussowitsch   PetscCall(VecReciprocal(diag));
839b4dc3ebdSPeter Brune   for (i = 0; i < cls->nsmooths; i++) {
840fb842aefSJose E. Roman     PetscCall(MatMatMult(A, *P, MAT_INITIAL_MATRIX, PETSC_CURRENT, &Pnew));
8419566063dSJacob Faibussowitsch     PetscCall(MatZeroRows(Pnew, idx, coarserows, 0., NULL, NULL));
8429566063dSJacob Faibussowitsch     PetscCall(MatDiagonalScale(Pnew, diag, NULL));
8439566063dSJacob Faibussowitsch     PetscCall(MatAYPX(Pnew, -1.0, *P, DIFFERENT_NONZERO_PATTERN));
8449566063dSJacob Faibussowitsch     PetscCall(MatDestroy(P));
845b4dc3ebdSPeter Brune     *P   = Pnew;
846b4dc3ebdSPeter Brune     Pnew = NULL;
847b4dc3ebdSPeter Brune   }
8489566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&diag));
8499566063dSJacob Faibussowitsch   PetscCall(PetscFree(coarserows));
8509566063dSJacob Faibussowitsch   PetscCall(PCGAMGTruncateProlongator_Private(pc, P));
8513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
852b4dc3ebdSPeter Brune }
853b4dc3ebdSPeter Brune 
PCGAMGProlongator_Classical(PC pc,Mat A,PetscCoarsenData * agg_lists,Mat * P)8548926f930SMark Adams static PetscErrorCode PCGAMGProlongator_Classical(PC pc, Mat A, PetscCoarsenData *agg_lists, Mat *P)
855d71ae5a4SJacob Faibussowitsch {
8568926f930SMark Adams   PetscErrorCode (*f)(PC, Mat, PetscCoarsenData *, Mat *);
8578eab0cc1SPeter Brune   PC_MG             *mg      = (PC_MG *)pc->data;
8588eab0cc1SPeter Brune   PC_GAMG           *pc_gamg = (PC_GAMG *)mg->innerctx;
8598eab0cc1SPeter Brune   PC_GAMG_Classical *cls     = (PC_GAMG_Classical *)pc_gamg->subctx;
8608eab0cc1SPeter Brune 
8618eab0cc1SPeter Brune   PetscFunctionBegin;
8629566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(PCGAMGClassicalProlongatorList, cls->prolongtype, &f));
86328b400f6SJacob Faibussowitsch   PetscCheck(f, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Cannot find PCGAMG Classical prolongator type");
8648926f930SMark Adams   PetscCall((*f)(pc, A, agg_lists, P));
8653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
866f9a65ec8SPeter Brune }
867f9a65ec8SPeter Brune 
PCGAMGDestroy_Classical(PC pc)868d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGDestroy_Classical(PC pc)
869d71ae5a4SJacob Faibussowitsch {
8708e6d0c30SPeter Brune   PC_MG   *mg      = (PC_MG *)pc->data;
8718e6d0c30SPeter Brune   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
8728e6d0c30SPeter Brune 
8738e6d0c30SPeter Brune   PetscFunctionBegin;
8749566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg->subctx));
8759566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGClassicalSetType_C", NULL));
8769566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGClassicalGetType_C", NULL));
8773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8788e6d0c30SPeter Brune }
8798e6d0c30SPeter Brune 
PCGAMGSetFromOptions_Classical(PC pc,PetscOptionItems PetscOptionsObject)880*ce78bad3SBarry Smith static PetscErrorCode PCGAMGSetFromOptions_Classical(PC pc, PetscOptionItems PetscOptionsObject)
881d71ae5a4SJacob Faibussowitsch {
882586a8384SPeter Brune   PC_MG             *mg      = (PC_MG *)pc->data;
883586a8384SPeter Brune   PC_GAMG           *pc_gamg = (PC_GAMG *)mg->innerctx;
884586a8384SPeter Brune   PC_GAMG_Classical *cls     = (PC_GAMG_Classical *)pc_gamg->subctx;
8858eab0cc1SPeter Brune   char               tname[256];
8868eab0cc1SPeter Brune   PetscBool          flg;
887586a8384SPeter Brune 
8888e6d0c30SPeter Brune   PetscFunctionBegin;
889d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "GAMG-Classical options");
8909566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-pc_gamg_classical_type", "Type of Classical AMG prolongation", "PCGAMGClassicalSetType", PCGAMGClassicalProlongatorList, cls->prolongtype, tname, sizeof(tname), &flg));
8911baa6e33SBarry Smith   if (flg) PetscCall(PCGAMGClassicalSetType(pc, tname));
8929566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_gamg_classical_interp_threshold", "Threshold for classical interpolator entries", "", cls->interp_threshold, &cls->interp_threshold, NULL));
8939566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_gamg_classical_nsmooths", "Threshold for classical interpolator entries", "", cls->nsmooths, &cls->nsmooths, NULL));
894d0609cedSBarry Smith   PetscOptionsHeadEnd();
8953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8968e6d0c30SPeter Brune }
8978e6d0c30SPeter Brune 
PCGAMGSetData_Classical(PC pc,Mat A)898d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetData_Classical(PC pc, Mat A)
899d71ae5a4SJacob Faibussowitsch {
9008e6d0c30SPeter Brune   PC_MG   *mg      = (PC_MG *)pc->data;
9018e6d0c30SPeter Brune   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
9028e6d0c30SPeter Brune 
9038e6d0c30SPeter Brune   PetscFunctionBegin;
9048e6d0c30SPeter Brune   /* no data for classical AMG */
9058e6d0c30SPeter Brune   pc_gamg->data           = NULL;
906d2050638SMark Adams   pc_gamg->data_cell_cols = 0;
907d2050638SMark Adams   pc_gamg->data_cell_rows = 0;
9088e6d0c30SPeter Brune   pc_gamg->data_sz        = 0;
9093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9108e6d0c30SPeter Brune }
9118e6d0c30SPeter Brune 
PCGAMGClassicalFinalizePackage(void)912ff6a9541SJacob Faibussowitsch static PetscErrorCode PCGAMGClassicalFinalizePackage(void)
913d71ae5a4SJacob Faibussowitsch {
9148eab0cc1SPeter Brune   PetscFunctionBegin;
9158eab0cc1SPeter Brune   PCGAMGClassicalPackageInitialized = PETSC_FALSE;
9169566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDestroy(&PCGAMGClassicalProlongatorList));
9173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9188eab0cc1SPeter Brune }
9198eab0cc1SPeter Brune 
PCGAMGClassicalInitializePackage(void)920ff6a9541SJacob Faibussowitsch static PetscErrorCode PCGAMGClassicalInitializePackage(void)
921d71ae5a4SJacob Faibussowitsch {
9228eab0cc1SPeter Brune   PetscFunctionBegin;
9233ba16761SJacob Faibussowitsch   if (PCGAMGClassicalPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS);
9249566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&PCGAMGClassicalProlongatorList, PCGAMGCLASSICALDIRECT, PCGAMGProlongator_Classical_Direct));
9259566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&PCGAMGClassicalProlongatorList, PCGAMGCLASSICALSTANDARD, PCGAMGProlongator_Classical_Standard));
9269566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalize(PCGAMGClassicalFinalizePackage));
9273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9288eab0cc1SPeter Brune }
9298eab0cc1SPeter Brune 
PCCreateGAMG_Classical(PC pc)930d71ae5a4SJacob Faibussowitsch PetscErrorCode PCCreateGAMG_Classical(PC pc)
931d71ae5a4SJacob Faibussowitsch {
9328e6d0c30SPeter Brune   PC_MG             *mg      = (PC_MG *)pc->data;
9338e6d0c30SPeter Brune   PC_GAMG           *pc_gamg = (PC_GAMG *)mg->innerctx;
9348e6d0c30SPeter Brune   PC_GAMG_Classical *pc_gamg_classical;
9358e6d0c30SPeter Brune 
9368e6d0c30SPeter Brune   PetscFunctionBegin;
9379566063dSJacob Faibussowitsch   PetscCall(PCGAMGClassicalInitializePackage());
9388e6d0c30SPeter Brune   if (pc_gamg->subctx) {
9398e6d0c30SPeter Brune     /* call base class */
9409566063dSJacob Faibussowitsch     PetscCall(PCDestroy_GAMG(pc));
9418e6d0c30SPeter Brune   }
9428e6d0c30SPeter Brune 
9438e6d0c30SPeter Brune   /* create sub context for SA */
9444dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&pc_gamg_classical));
9458e6d0c30SPeter Brune   pc_gamg->subctx         = pc_gamg_classical;
9468e6d0c30SPeter Brune   pc->ops->setfromoptions = PCGAMGSetFromOptions_Classical;
9478e6d0c30SPeter Brune   /* reset does not do anything; setup not virtual */
9488e6d0c30SPeter Brune 
9498e6d0c30SPeter Brune   /* set internal function pointers */
9508e6d0c30SPeter Brune   pc_gamg->ops->destroy        = PCGAMGDestroy_Classical;
9512d776b49SBarry Smith   pc_gamg->ops->creategraph    = PCGAMGCreateGraph_Classical;
9528e6d0c30SPeter Brune   pc_gamg->ops->coarsen        = PCGAMGCoarsen_Classical;
9538eab0cc1SPeter Brune   pc_gamg->ops->prolongator    = PCGAMGProlongator_Classical;
954fd1112cbSBarry Smith   pc_gamg->ops->optprolongator = PCGAMGOptProlongator_Classical_Jacobi;
955586a8384SPeter Brune   pc_gamg->ops->setfromoptions = PCGAMGSetFromOptions_Classical;
9568e6d0c30SPeter Brune 
9578e6d0c30SPeter Brune   pc_gamg->ops->createdefaultdata     = PCGAMGSetData_Classical;
958586a8384SPeter Brune   pc_gamg_classical->interp_threshold = 0.2;
959b4dc3ebdSPeter Brune   pc_gamg_classical->nsmooths         = 0;
9609566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGClassicalSetType_C", PCGAMGClassicalSetType_GAMG));
9619566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGClassicalGetType_C", PCGAMGClassicalGetType_GAMG));
9629566063dSJacob Faibussowitsch   PetscCall(PCGAMGClassicalSetType(pc, PCGAMGCLASSICALSTANDARD));
9633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9648e6d0c30SPeter Brune }
965