xref: /petsc/src/ksp/pc/impls/gamg/classical.c (revision ff6a95418ff72e09afb68819ffbfc86bad111fa0)
18e6d0c30SPeter Brune #include <../src/ksp/pc/impls/gamg/gamg.h> /*I "petscpc.h" I*/
287b9b13cSPeter Brune #include <petscsf.h>
38e6d0c30SPeter Brune 
4*ff6a9541SJacob Faibussowitsch static PetscFunctionList PCGAMGClassicalProlongatorList    = NULL;
5*ff6a9541SJacob Faibussowitsch static PetscBool         PCGAMGClassicalPackageInitialized = PETSC_FALSE;
68eab0cc1SPeter Brune 
78e6d0c30SPeter Brune typedef struct {
8586a8384SPeter Brune   PetscReal interp_threshold; /* interpolation threshold */
98eab0cc1SPeter Brune   char      prolongtype[256];
10b4dc3ebdSPeter Brune   PetscInt  nsmooths; /* number of jacobi smoothings on the prolongator */
118e6d0c30SPeter Brune } PC_GAMG_Classical;
128e6d0c30SPeter Brune 
137779008dSPeter Brune /*@C
14f1580f4eSBarry Smith    PCGAMGClassicalSetType - Sets the type of classical interpolation to use with `PCGAMG`
158eab0cc1SPeter Brune 
16c3339decSBarry Smith    Collective
178eab0cc1SPeter Brune 
188eab0cc1SPeter Brune    Input Parameters:
198eab0cc1SPeter Brune .  pc - the preconditioner context
208eab0cc1SPeter Brune 
218eab0cc1SPeter Brune    Options Database Key:
22f1580f4eSBarry Smith .  -pc_gamg_classical_type <direct,standard> - set type of classical AMG prolongation
238eab0cc1SPeter Brune 
248eab0cc1SPeter Brune    Level: intermediate
258eab0cc1SPeter Brune 
26f1580f4eSBarry Smith .seealso: `PCGAMG`
278eab0cc1SPeter Brune @*/
28d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGClassicalSetType(PC pc, PCGAMGClassicalType type)
29d71ae5a4SJacob Faibussowitsch {
308eab0cc1SPeter Brune   PetscFunctionBegin;
318eab0cc1SPeter Brune   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
32cac4c232SBarry Smith   PetscTryMethod(pc, "PCGAMGClassicalSetType_C", (PC, PCGAMGClassicalType), (pc, type));
333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
348eab0cc1SPeter Brune }
358eab0cc1SPeter Brune 
36c60c7ad4SBarry Smith /*@C
37f1580f4eSBarry Smith    PCGAMGClassicalGetType - Gets the type of classical interpolation to use with `PCGAMG`
38c60c7ad4SBarry Smith 
39c3339decSBarry Smith    Collective
40c60c7ad4SBarry Smith 
41c60c7ad4SBarry Smith    Input Parameter:
42c60c7ad4SBarry Smith .  pc - the preconditioner context
43c60c7ad4SBarry Smith 
44c60c7ad4SBarry Smith    Output Parameter:
45c60c7ad4SBarry Smith .  type - the type used
46c60c7ad4SBarry Smith 
47c60c7ad4SBarry Smith    Level: intermediate
48c60c7ad4SBarry Smith 
49f1580f4eSBarry Smith .seealso: `PCGAMG`
50c60c7ad4SBarry Smith @*/
51d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGClassicalGetType(PC pc, PCGAMGClassicalType *type)
52d71ae5a4SJacob Faibussowitsch {
53c60c7ad4SBarry Smith   PetscFunctionBegin;
54c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
55cac4c232SBarry Smith   PetscUseMethod(pc, "PCGAMGClassicalGetType_C", (PC, PCGAMGClassicalType *), (pc, type));
563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
57c60c7ad4SBarry Smith }
58c60c7ad4SBarry Smith 
59d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGClassicalSetType_GAMG(PC pc, PCGAMGClassicalType type)
60d71ae5a4SJacob Faibussowitsch {
618eab0cc1SPeter Brune   PC_MG             *mg      = (PC_MG *)pc->data;
628eab0cc1SPeter Brune   PC_GAMG           *pc_gamg = (PC_GAMG *)mg->innerctx;
638eab0cc1SPeter Brune   PC_GAMG_Classical *cls     = (PC_GAMG_Classical *)pc_gamg->subctx;
648eab0cc1SPeter Brune 
658eab0cc1SPeter Brune   PetscFunctionBegin;
669566063dSJacob Faibussowitsch   PetscCall(PetscStrcpy(cls->prolongtype, type));
673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
688eab0cc1SPeter Brune }
698e6d0c30SPeter Brune 
70d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGClassicalGetType_GAMG(PC pc, PCGAMGClassicalType *type)
71d71ae5a4SJacob Faibussowitsch {
72c60c7ad4SBarry Smith   PC_MG             *mg      = (PC_MG *)pc->data;
73c60c7ad4SBarry Smith   PC_GAMG           *pc_gamg = (PC_GAMG *)mg->innerctx;
74c60c7ad4SBarry Smith   PC_GAMG_Classical *cls     = (PC_GAMG_Classical *)pc_gamg->subctx;
75c60c7ad4SBarry Smith 
76c60c7ad4SBarry Smith   PetscFunctionBegin;
77c60c7ad4SBarry Smith   *type = cls->prolongtype;
783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
79c60c7ad4SBarry Smith }
80c60c7ad4SBarry Smith 
81*ff6a9541SJacob Faibussowitsch static PetscErrorCode PCGAMGCreateGraph_Classical(PC pc, Mat A, Mat *G)
82d71ae5a4SJacob Faibussowitsch {
83550383edSPeter Brune   PetscInt           s, f, n, idx, lidx, gidx;
84e5a0faa4SPeter Brune   PetscInt           r, c, ncols;
858e6d0c30SPeter Brune   const PetscInt    *rcol;
868e6d0c30SPeter Brune   const PetscScalar *rval;
87e5a0faa4SPeter Brune   PetscInt          *gcol;
888e6d0c30SPeter Brune   PetscScalar       *gval;
89e5a0faa4SPeter Brune   PetscReal          rmax;
90550383edSPeter Brune   PetscInt           cmax = 0;
91b817416eSBarry Smith   PC_MG             *mg   = (PC_MG *)pc->data;
92b817416eSBarry Smith   PC_GAMG           *gamg = (PC_GAMG *)mg->innerctx;
938e6d0c30SPeter Brune   PetscInt          *gsparse, *lsparse;
94e5a0faa4SPeter Brune   PetscScalar       *Amax;
958e6d0c30SPeter Brune   MatType            mtype;
968e6d0c30SPeter Brune 
978e6d0c30SPeter Brune   PetscFunctionBegin;
989566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &s, &f));
99550383edSPeter Brune   n = f - s;
1009566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(n, &lsparse, n, &gsparse, n, &Amax));
1018e6d0c30SPeter Brune 
102550383edSPeter Brune   for (r = 0; r < n; r++) {
1038e6d0c30SPeter Brune     lsparse[r] = 0;
104550383edSPeter Brune     gsparse[r] = 0;
1058e6d0c30SPeter Brune   }
1068e6d0c30SPeter Brune 
107550383edSPeter Brune   for (r = s; r < f; r++) {
108e5a0faa4SPeter Brune     /* determine the maximum off-diagonal in each row */
109e5a0faa4SPeter Brune     rmax = 0.;
1109566063dSJacob Faibussowitsch     PetscCall(MatGetRow(A, r, &ncols, &rcol, &rval));
111e5a0faa4SPeter Brune     for (c = 0; c < ncols; c++) {
112ad540459SPierre Jolivet       if (PetscRealPart(-rval[c]) > rmax && rcol[c] != r) rmax = PetscRealPart(-rval[c]);
113e5a0faa4SPeter Brune     }
114550383edSPeter Brune     Amax[r - s] = rmax;
115550383edSPeter Brune     if (ncols > cmax) cmax = ncols;
116550383edSPeter Brune     lidx = 0;
117550383edSPeter Brune     gidx = 0;
118e5a0faa4SPeter Brune     /* create the local and global sparsity patterns */
1198e6d0c30SPeter Brune     for (c = 0; c < ncols; c++) {
120c1eae691SMark Adams       if (PetscRealPart(-rval[c]) > gamg->threshold[0] * PetscRealPart(Amax[r - s]) || rcol[c] == r) {
121550383edSPeter Brune         if (rcol[c] < f && rcol[c] >= s) {
122550383edSPeter Brune           lidx++;
123550383edSPeter Brune         } else {
124550383edSPeter Brune           gidx++;
1258e6d0c30SPeter Brune         }
1268e6d0c30SPeter Brune       }
1278e6d0c30SPeter Brune     }
1289566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(A, r, &ncols, &rcol, &rval));
129550383edSPeter Brune     lsparse[r - s] = lidx;
130550383edSPeter Brune     gsparse[r - s] = gidx;
1318e6d0c30SPeter Brune   }
1329566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(cmax, &gval, cmax, &gcol));
133e5a0faa4SPeter Brune 
1349566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), G));
1359566063dSJacob Faibussowitsch   PetscCall(MatGetType(A, &mtype));
1369566063dSJacob Faibussowitsch   PetscCall(MatSetType(*G, mtype));
1379566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*G, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
1389566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(*G, 0, lsparse, 0, gsparse));
1399566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(*G, 0, lsparse));
1408e6d0c30SPeter Brune   for (r = s; r < f; r++) {
1419566063dSJacob Faibussowitsch     PetscCall(MatGetRow(A, r, &ncols, &rcol, &rval));
1428e6d0c30SPeter Brune     idx = 0;
1438e6d0c30SPeter Brune     for (c = 0; c < ncols; c++) {
1448e6d0c30SPeter Brune       /* classical strength of connection */
145c1eae691SMark Adams       if (PetscRealPart(-rval[c]) > gamg->threshold[0] * PetscRealPart(Amax[r - s]) || rcol[c] == r) {
1468e6d0c30SPeter Brune         gcol[idx] = rcol[c];
1478e6d0c30SPeter Brune         gval[idx] = rval[c];
1488e6d0c30SPeter Brune         idx++;
1498e6d0c30SPeter Brune       }
1508e6d0c30SPeter Brune     }
1519566063dSJacob Faibussowitsch     PetscCall(MatSetValues(*G, 1, &r, idx, gcol, gval, INSERT_VALUES));
1529566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(A, r, &ncols, &rcol, &rval));
1538e6d0c30SPeter Brune   }
1549566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*G, MAT_FINAL_ASSEMBLY));
1559566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*G, MAT_FINAL_ASSEMBLY));
1568e6d0c30SPeter Brune 
1579566063dSJacob Faibussowitsch   PetscCall(PetscFree2(gval, gcol));
1589566063dSJacob Faibussowitsch   PetscCall(PetscFree3(lsparse, gsparse, Amax));
1593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1608e6d0c30SPeter Brune }
1618e6d0c30SPeter Brune 
162*ff6a9541SJacob Faibussowitsch static PetscErrorCode PCGAMGCoarsen_Classical(PC pc, Mat *G, PetscCoarsenData **agg_lists)
163d71ae5a4SJacob Faibussowitsch {
1648e6d0c30SPeter Brune   MatCoarsen crs;
1658e6d0c30SPeter Brune   MPI_Comm   fcomm = ((PetscObject)pc)->comm;
1668e6d0c30SPeter Brune 
1678e6d0c30SPeter Brune   PetscFunctionBegin;
16828b400f6SJacob Faibussowitsch   PetscCheck(G, fcomm, PETSC_ERR_ARG_WRONGSTATE, "Must set Graph in PC in PCGAMG before coarsening");
1698e6d0c30SPeter Brune 
1709566063dSJacob Faibussowitsch   PetscCall(MatCoarsenCreate(fcomm, &crs));
1719566063dSJacob Faibussowitsch   PetscCall(MatCoarsenSetFromOptions(crs));
1729566063dSJacob Faibussowitsch   PetscCall(MatCoarsenSetAdjacency(crs, *G));
1739566063dSJacob Faibussowitsch   PetscCall(MatCoarsenSetStrictAggs(crs, PETSC_TRUE));
1749566063dSJacob Faibussowitsch   PetscCall(MatCoarsenApply(crs));
1759566063dSJacob Faibussowitsch   PetscCall(MatCoarsenGetData(crs, agg_lists));
1769566063dSJacob Faibussowitsch   PetscCall(MatCoarsenDestroy(&crs));
1773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1788e6d0c30SPeter Brune }
1798e6d0c30SPeter Brune 
180*ff6a9541SJacob Faibussowitsch static PetscErrorCode PCGAMGProlongator_Classical_Direct(PC pc, Mat A, Mat G, PetscCoarsenData *agg_lists, Mat *P)
181d71ae5a4SJacob Faibussowitsch {
1821ce39c63SPeter Brune   PC_MG             *mg   = (PC_MG *)pc->data;
1831ce39c63SPeter Brune   PC_GAMG           *gamg = (PC_GAMG *)mg->innerctx;
18463b0c588SPeter Brune   PetscBool          iscoarse, isMPIAIJ, isSEQAIJ;
18563b0c588SPeter Brune   PetscInt           fn, cn, fs, fe, cs, ce, i, j, ncols, col, row_f, row_c, cmax = 0, idx, noff;
18663b0c588SPeter Brune   PetscInt          *lcid, *gcid, *lsparse, *gsparse, *colmap, *pcols;
18763b0c588SPeter Brune   const PetscInt    *rcol;
18863b0c588SPeter Brune   PetscReal         *Amax_pos, *Amax_neg;
18963b0c588SPeter Brune   PetscScalar        g_pos, g_neg, a_pos, a_neg, diag, invdiag, alpha, beta, pij;
19063b0c588SPeter Brune   PetscScalar       *pvals;
19163b0c588SPeter Brune   const PetscScalar *rval;
19263b0c588SPeter Brune   Mat                lA, gA = NULL;
19363b0c588SPeter Brune   MatType            mtype;
19463b0c588SPeter Brune   Vec                C, lvec;
19587b9b13cSPeter Brune   PetscLayout        clayout;
19687b9b13cSPeter Brune   PetscSF            sf;
19787b9b13cSPeter Brune   Mat_MPIAIJ        *mpiaij;
1988e6d0c30SPeter Brune 
1998e6d0c30SPeter Brune   PetscFunctionBegin;
2009566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &fs, &fe));
20187b9b13cSPeter Brune   fn = fe - fs;
2029566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPIAIJ, &isMPIAIJ));
2039566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQAIJ, &isSEQAIJ));
2047827d75bSBarry Smith   PetscCheck(isMPIAIJ || isSEQAIJ, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Classical AMG requires MPIAIJ matrix");
20587b9b13cSPeter Brune   if (isMPIAIJ) {
20687b9b13cSPeter Brune     mpiaij = (Mat_MPIAIJ *)A->data;
20787b9b13cSPeter Brune     lA     = mpiaij->A;
20887b9b13cSPeter Brune     gA     = mpiaij->B;
20987b9b13cSPeter Brune     lvec   = mpiaij->lvec;
2109566063dSJacob Faibussowitsch     PetscCall(VecGetSize(lvec, &noff));
21187b9b13cSPeter Brune     colmap = mpiaij->garray;
2129566063dSJacob Faibussowitsch     PetscCall(MatGetLayouts(A, NULL, &clayout));
2139566063dSJacob Faibussowitsch     PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)A), &sf));
2149566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraphLayout(sf, clayout, noff, NULL, PETSC_COPY_VALUES, colmap));
2159566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(noff, &gcid));
21687b9b13cSPeter Brune   } else {
21787b9b13cSPeter Brune     lA = A;
21887b9b13cSPeter Brune   }
2199566063dSJacob Faibussowitsch   PetscCall(PetscMalloc5(fn, &lsparse, fn, &gsparse, fn, &lcid, fn, &Amax_pos, fn, &Amax_neg));
2208e6d0c30SPeter Brune 
2218e6d0c30SPeter Brune   /* count the number of coarse unknowns */
2228e6d0c30SPeter Brune   cn = 0;
2238e6d0c30SPeter Brune   for (i = 0; i < fn; i++) {
2248e6d0c30SPeter Brune     /* filter out singletons */
2259566063dSJacob Faibussowitsch     PetscCall(PetscCDEmptyAt(agg_lists, i, &iscoarse));
2268e6d0c30SPeter Brune     lcid[i] = -1;
227ad540459SPierre Jolivet     if (!iscoarse) cn++;
2288e6d0c30SPeter Brune   }
2298e6d0c30SPeter Brune 
2308e6d0c30SPeter Brune   /* create the coarse vector */
2319566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)A), cn, PETSC_DECIDE, &C));
2329566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(C, &cs, &ce));
2338e6d0c30SPeter Brune 
2348e6d0c30SPeter Brune   cn = 0;
2358e6d0c30SPeter Brune   for (i = 0; i < fn; i++) {
2369566063dSJacob Faibussowitsch     PetscCall(PetscCDEmptyAt(agg_lists, i, &iscoarse));
2378e6d0c30SPeter Brune     if (!iscoarse) {
2388e6d0c30SPeter Brune       lcid[i] = cs + cn;
2398e6d0c30SPeter Brune       cn++;
2408e6d0c30SPeter Brune     } else {
2418e6d0c30SPeter Brune       lcid[i] = -1;
2428e6d0c30SPeter Brune     }
2438e6d0c30SPeter Brune   }
2448e6d0c30SPeter Brune 
24587b9b13cSPeter Brune   if (gA) {
2469566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(sf, MPIU_INT, lcid, gcid, MPI_REPLACE));
2479566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(sf, MPIU_INT, lcid, gcid, MPI_REPLACE));
24887b9b13cSPeter Brune   }
2498e6d0c30SPeter Brune 
250b817416eSBarry Smith   /* determine the largest off-diagonal entries in each row */
2511ce39c63SPeter Brune   for (i = fs; i < fe; i++) {
2521ce39c63SPeter Brune     Amax_pos[i - fs] = 0.;
2531ce39c63SPeter Brune     Amax_neg[i - fs] = 0.;
2549566063dSJacob Faibussowitsch     PetscCall(MatGetRow(A, i, &ncols, &rcol, &rval));
2551ce39c63SPeter Brune     for (j = 0; j < ncols; j++) {
2561ce39c63SPeter Brune       if ((PetscRealPart(-rval[j]) > Amax_neg[i - fs]) && i != rcol[j]) Amax_neg[i - fs] = PetscAbsScalar(rval[j]);
2571ce39c63SPeter Brune       if ((PetscRealPart(rval[j]) > Amax_pos[i - fs]) && i != rcol[j]) Amax_pos[i - fs] = PetscAbsScalar(rval[j]);
2581ce39c63SPeter Brune     }
2591ce39c63SPeter Brune     if (ncols > cmax) cmax = ncols;
2609566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(A, i, &ncols, &rcol, &rval));
2611ce39c63SPeter Brune   }
2629566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(cmax, &pcols, cmax, &pvals));
2639566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&C));
2648e6d0c30SPeter Brune 
2658e6d0c30SPeter Brune   /* count the on and off processor sparsity patterns for the prolongator */
2668e6d0c30SPeter Brune   for (i = 0; i < fn; i++) {
2678e6d0c30SPeter Brune     /* on */
2688e6d0c30SPeter Brune     lsparse[i] = 0;
269e5a0faa4SPeter Brune     gsparse[i] = 0;
2708e6d0c30SPeter Brune     if (lcid[i] >= 0) {
2718e6d0c30SPeter Brune       lsparse[i] = 1;
2728e6d0c30SPeter Brune       gsparse[i] = 0;
2738e6d0c30SPeter Brune     } else {
2749566063dSJacob Faibussowitsch       PetscCall(MatGetRow(lA, i, &ncols, &rcol, &rval));
2758e6d0c30SPeter Brune       for (j = 0; j < ncols; j++) {
2761ce39c63SPeter Brune         col = rcol[j];
277ad540459SPierre 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;
2788e6d0c30SPeter Brune       }
2799566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(lA, i, &ncols, &rcol, &rval));
2808e6d0c30SPeter Brune       /* off */
2811ce39c63SPeter Brune       if (gA) {
2829566063dSJacob Faibussowitsch         PetscCall(MatGetRow(gA, i, &ncols, &rcol, &rval));
2838e6d0c30SPeter Brune         for (j = 0; j < ncols; j++) {
2841ce39c63SPeter Brune           col = rcol[j];
285ad540459SPierre 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;
2868e6d0c30SPeter Brune         }
2879566063dSJacob Faibussowitsch         PetscCall(MatRestoreRow(gA, i, &ncols, &rcol, &rval));
2888e6d0c30SPeter Brune       }
2898e6d0c30SPeter Brune     }
2901ce39c63SPeter Brune   }
2918e6d0c30SPeter Brune 
2928e6d0c30SPeter Brune   /* preallocate and create the prolongator */
2939566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), P));
2949566063dSJacob Faibussowitsch   PetscCall(MatGetType(G, &mtype));
2959566063dSJacob Faibussowitsch   PetscCall(MatSetType(*P, mtype));
2969566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*P, fn, cn, PETSC_DETERMINE, PETSC_DETERMINE));
2979566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(*P, 0, lsparse, 0, gsparse));
2989566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(*P, 0, lsparse));
2998e6d0c30SPeter Brune 
3008e6d0c30SPeter Brune   /* loop over local fine nodes -- get the diagonal, the sum of positive and negative strong and weak weights, and set up the row */
3018e6d0c30SPeter Brune   for (i = 0; i < fn; i++) {
3028e6d0c30SPeter Brune     /* determine on or off */
3038e6d0c30SPeter Brune     row_f = i + fs;
3048e6d0c30SPeter Brune     row_c = lcid[i];
3058e6d0c30SPeter Brune     if (row_c >= 0) {
3068e6d0c30SPeter Brune       pij = 1.;
3079566063dSJacob Faibussowitsch       PetscCall(MatSetValues(*P, 1, &row_f, 1, &row_c, &pij, INSERT_VALUES));
3088e6d0c30SPeter Brune     } else {
3098e6d0c30SPeter Brune       g_pos = 0.;
3108e6d0c30SPeter Brune       g_neg = 0.;
3118e6d0c30SPeter Brune       a_pos = 0.;
3128e6d0c30SPeter Brune       a_neg = 0.;
3138e6d0c30SPeter Brune       diag  = 0.;
3148e6d0c30SPeter Brune 
3151ce39c63SPeter Brune       /* local connections */
3169566063dSJacob Faibussowitsch       PetscCall(MatGetRow(lA, i, &ncols, &rcol, &rval));
3171ce39c63SPeter Brune       for (j = 0; j < ncols; j++) {
3181ce39c63SPeter Brune         col = rcol[j];
319c1eae691SMark Adams         if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0] * Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0] * Amax_neg[i])) {
3201ce39c63SPeter Brune           if (PetscRealPart(rval[j]) > 0.) {
3211ce39c63SPeter Brune             g_pos += rval[j];
3228e6d0c30SPeter Brune           } else {
3231ce39c63SPeter Brune             g_neg += rval[j];
3248e6d0c30SPeter Brune           }
3251ce39c63SPeter Brune         }
3261ce39c63SPeter Brune         if (col != i) {
3271ce39c63SPeter Brune           if (PetscRealPart(rval[j]) > 0.) {
3281ce39c63SPeter Brune             a_pos += rval[j];
3291ce39c63SPeter Brune           } else {
3301ce39c63SPeter Brune             a_neg += rval[j];
3311ce39c63SPeter Brune           }
3321ce39c63SPeter Brune         } else {
3331ce39c63SPeter Brune           diag = rval[j];
3341ce39c63SPeter Brune         }
3358e6d0c30SPeter Brune       }
3369566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(lA, i, &ncols, &rcol, &rval));
3378e6d0c30SPeter Brune 
3381ce39c63SPeter Brune       /* ghosted connections */
3398e6d0c30SPeter Brune       if (gA) {
3409566063dSJacob Faibussowitsch         PetscCall(MatGetRow(gA, i, &ncols, &rcol, &rval));
3411ce39c63SPeter Brune         for (j = 0; j < ncols; j++) {
3421ce39c63SPeter Brune           col = rcol[j];
343c1eae691SMark Adams           if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0] * Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0] * Amax_neg[i])) {
3441ce39c63SPeter Brune             if (PetscRealPart(rval[j]) > 0.) {
3451ce39c63SPeter Brune               g_pos += rval[j];
3468e6d0c30SPeter Brune             } else {
3471ce39c63SPeter Brune               g_neg += rval[j];
3488e6d0c30SPeter Brune             }
3491ce39c63SPeter Brune           }
3501ce39c63SPeter Brune           if (PetscRealPart(rval[j]) > 0.) {
3511ce39c63SPeter Brune             a_pos += rval[j];
3521ce39c63SPeter Brune           } else {
3531ce39c63SPeter Brune             a_neg += rval[j];
3541ce39c63SPeter Brune           }
3558e6d0c30SPeter Brune         }
3569566063dSJacob Faibussowitsch         PetscCall(MatRestoreRow(gA, i, &ncols, &rcol, &rval));
3578e6d0c30SPeter Brune       }
3588e6d0c30SPeter Brune 
3598e6d0c30SPeter Brune       if (g_neg == 0.) {
3608e6d0c30SPeter Brune         alpha = 0.;
3618e6d0c30SPeter Brune       } else {
3628e6d0c30SPeter Brune         alpha = -a_neg / g_neg;
3638e6d0c30SPeter Brune       }
3648e6d0c30SPeter Brune 
3658e6d0c30SPeter Brune       if (g_pos == 0.) {
3668e6d0c30SPeter Brune         diag += a_pos;
3678e6d0c30SPeter Brune         beta = 0.;
3688e6d0c30SPeter Brune       } else {
3698e6d0c30SPeter Brune         beta = -a_pos / g_pos;
3708e6d0c30SPeter Brune       }
371e5a0faa4SPeter Brune       if (diag == 0.) {
372e5a0faa4SPeter Brune         invdiag = 0.;
373e5a0faa4SPeter Brune       } else invdiag = 1. / diag;
3748e6d0c30SPeter Brune       /* on */
3759566063dSJacob Faibussowitsch       PetscCall(MatGetRow(lA, i, &ncols, &rcol, &rval));
3768e6d0c30SPeter Brune       idx = 0;
3778e6d0c30SPeter Brune       for (j = 0; j < ncols; j++) {
3781ce39c63SPeter Brune         col = rcol[j];
379c1eae691SMark Adams         if (lcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0] * Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0] * Amax_neg[i])) {
3808e6d0c30SPeter Brune           row_f = i + fs;
3818e6d0c30SPeter Brune           row_c = lcid[col];
3828e6d0c30SPeter Brune           /* set the values for on-processor ones */
3831ce39c63SPeter Brune           if (PetscRealPart(rval[j]) < 0.) {
3841ce39c63SPeter Brune             pij = rval[j] * alpha * invdiag;
3858e6d0c30SPeter Brune           } else {
3861ce39c63SPeter Brune             pij = rval[j] * beta * invdiag;
3878e6d0c30SPeter Brune           }
3888e6d0c30SPeter Brune           if (PetscAbsScalar(pij) != 0.) {
3898e6d0c30SPeter Brune             pvals[idx] = pij;
3908e6d0c30SPeter Brune             pcols[idx] = row_c;
3918e6d0c30SPeter Brune             idx++;
3928e6d0c30SPeter Brune           }
3938e6d0c30SPeter Brune         }
3948e6d0c30SPeter Brune       }
3959566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(lA, i, &ncols, &rcol, &rval));
3968e6d0c30SPeter Brune       /* off */
3971ce39c63SPeter Brune       if (gA) {
3989566063dSJacob Faibussowitsch         PetscCall(MatGetRow(gA, i, &ncols, &rcol, &rval));
3998e6d0c30SPeter Brune         for (j = 0; j < ncols; j++) {
4001ce39c63SPeter Brune           col = rcol[j];
401c1eae691SMark Adams           if (gcid[col] >= 0 && (PetscRealPart(rval[j]) > gamg->threshold[0] * Amax_pos[i] || PetscRealPart(-rval[j]) > gamg->threshold[0] * Amax_neg[i])) {
4028e6d0c30SPeter Brune             row_f = i + fs;
4038e6d0c30SPeter Brune             row_c = gcid[col];
4048e6d0c30SPeter Brune             /* set the values for on-processor ones */
4051ce39c63SPeter Brune             if (PetscRealPart(rval[j]) < 0.) {
4061ce39c63SPeter Brune               pij = rval[j] * alpha * invdiag;
4078e6d0c30SPeter Brune             } else {
4081ce39c63SPeter Brune               pij = rval[j] * beta * invdiag;
4098e6d0c30SPeter Brune             }
4108e6d0c30SPeter Brune             if (PetscAbsScalar(pij) != 0.) {
4118e6d0c30SPeter Brune               pvals[idx] = pij;
4128e6d0c30SPeter Brune               pcols[idx] = row_c;
4138e6d0c30SPeter Brune               idx++;
4148e6d0c30SPeter Brune             }
4158e6d0c30SPeter Brune           }
4168e6d0c30SPeter Brune         }
4179566063dSJacob Faibussowitsch         PetscCall(MatRestoreRow(gA, i, &ncols, &rcol, &rval));
4183c9ab2c3SPeter Brune       }
4199566063dSJacob Faibussowitsch       PetscCall(MatSetValues(*P, 1, &row_f, idx, pcols, pvals, INSERT_VALUES));
4208e6d0c30SPeter Brune     }
4218e6d0c30SPeter Brune   }
4223c9ab2c3SPeter Brune 
4239566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY));
4249566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY));
4258e6d0c30SPeter Brune 
4269566063dSJacob Faibussowitsch   PetscCall(PetscFree5(lsparse, gsparse, lcid, Amax_pos, Amax_neg));
427e632b94dSBarry Smith 
4289566063dSJacob Faibussowitsch   PetscCall(PetscFree2(pcols, pvals));
42987b9b13cSPeter Brune   if (gA) {
4309566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&sf));
4319566063dSJacob Faibussowitsch     PetscCall(PetscFree(gcid));
43287b9b13cSPeter Brune   }
4333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4348e6d0c30SPeter Brune }
4358e6d0c30SPeter Brune 
436*ff6a9541SJacob Faibussowitsch static PetscErrorCode PCGAMGTruncateProlongator_Private(PC pc, Mat *P)
437d71ae5a4SJacob Faibussowitsch {
438586a8384SPeter Brune   PetscInt           j, i, ps, pf, pn, pcs, pcf, pcn, idx, cmax;
439586a8384SPeter Brune   const PetscScalar *pval;
440586a8384SPeter Brune   const PetscInt    *pcol;
441586a8384SPeter Brune   PetscScalar       *pnval;
442586a8384SPeter Brune   PetscInt          *pncol;
443586a8384SPeter Brune   PetscInt           ncols;
444586a8384SPeter Brune   Mat                Pnew;
445586a8384SPeter Brune   PetscInt          *lsparse, *gsparse;
446586a8384SPeter Brune   PetscReal          pmax_pos, pmax_neg, ptot_pos, ptot_neg, pthresh_pos, pthresh_neg;
447586a8384SPeter Brune   PC_MG             *mg      = (PC_MG *)pc->data;
448586a8384SPeter Brune   PC_GAMG           *pc_gamg = (PC_GAMG *)mg->innerctx;
449586a8384SPeter Brune   PC_GAMG_Classical *cls     = (PC_GAMG_Classical *)pc_gamg->subctx;
450d9558ea9SBarry Smith   MatType            mtype;
451586a8384SPeter Brune 
452586a8384SPeter Brune   PetscFunctionBegin;
453586a8384SPeter Brune   /* trim and rescale with reallocation */
4549566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(*P, &ps, &pf));
4559566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRangeColumn(*P, &pcs, &pcf));
456586a8384SPeter Brune   pn  = pf - ps;
457586a8384SPeter Brune   pcn = pcf - pcs;
4589566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(pn, &lsparse, pn, &gsparse));
459586a8384SPeter Brune   /* allocate */
460586a8384SPeter Brune   cmax = 0;
461586a8384SPeter Brune   for (i = ps; i < pf; i++) {
462b4dc3ebdSPeter Brune     lsparse[i - ps] = 0;
463b4dc3ebdSPeter Brune     gsparse[i - ps] = 0;
4649566063dSJacob Faibussowitsch     PetscCall(MatGetRow(*P, i, &ncols, &pcol, &pval));
465ad540459SPierre Jolivet     if (ncols > cmax) cmax = ncols;
466586a8384SPeter Brune     pmax_pos = 0.;
467586a8384SPeter Brune     pmax_neg = 0.;
468586a8384SPeter Brune     for (j = 0; j < ncols; j++) {
469586a8384SPeter Brune       if (PetscRealPart(pval[j]) > pmax_pos) {
470586a8384SPeter Brune         pmax_pos = PetscRealPart(pval[j]);
471586a8384SPeter Brune       } else if (PetscRealPart(pval[j]) < pmax_neg) {
472586a8384SPeter Brune         pmax_neg = PetscRealPart(pval[j]);
473586a8384SPeter Brune       }
474586a8384SPeter Brune     }
475586a8384SPeter Brune     for (j = 0; j < ncols; j++) {
4768eab0cc1SPeter Brune       if (PetscRealPart(pval[j]) >= pmax_pos * cls->interp_threshold || PetscRealPart(pval[j]) <= pmax_neg * cls->interp_threshold) {
477b4dc3ebdSPeter Brune         if (pcol[j] >= pcs && pcol[j] < pcf) {
478b4dc3ebdSPeter Brune           lsparse[i - ps]++;
479586a8384SPeter Brune         } else {
480b4dc3ebdSPeter Brune           gsparse[i - ps]++;
481586a8384SPeter Brune         }
482586a8384SPeter Brune       }
483586a8384SPeter Brune     }
4849566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(*P, i, &ncols, &pcol, &pval));
485586a8384SPeter Brune   }
486586a8384SPeter Brune 
4879566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(cmax, &pnval, cmax, &pncol));
488586a8384SPeter Brune 
4899566063dSJacob Faibussowitsch   PetscCall(MatGetType(*P, &mtype));
4909566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)*P), &Pnew));
4919566063dSJacob Faibussowitsch   PetscCall(MatSetType(Pnew, mtype));
4929566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Pnew, pn, pcn, PETSC_DETERMINE, PETSC_DETERMINE));
4939566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(Pnew, 0, lsparse));
4949566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(Pnew, 0, lsparse, 0, gsparse));
495586a8384SPeter Brune 
496586a8384SPeter Brune   for (i = ps; i < pf; i++) {
4979566063dSJacob Faibussowitsch     PetscCall(MatGetRow(*P, i, &ncols, &pcol, &pval));
498586a8384SPeter Brune     pmax_pos = 0.;
499586a8384SPeter Brune     pmax_neg = 0.;
500586a8384SPeter Brune     for (j = 0; j < ncols; j++) {
501586a8384SPeter Brune       if (PetscRealPart(pval[j]) > pmax_pos) {
502586a8384SPeter Brune         pmax_pos = PetscRealPart(pval[j]);
503586a8384SPeter Brune       } else if (PetscRealPart(pval[j]) < pmax_neg) {
504586a8384SPeter Brune         pmax_neg = PetscRealPart(pval[j]);
505586a8384SPeter Brune       }
506586a8384SPeter Brune     }
507586a8384SPeter Brune     pthresh_pos = 0.;
508586a8384SPeter Brune     pthresh_neg = 0.;
509586a8384SPeter Brune     ptot_pos    = 0.;
510586a8384SPeter Brune     ptot_neg    = 0.;
511586a8384SPeter Brune     for (j = 0; j < ncols; j++) {
5128eab0cc1SPeter Brune       if (PetscRealPart(pval[j]) >= cls->interp_threshold * pmax_pos) {
513586a8384SPeter Brune         pthresh_pos += PetscRealPart(pval[j]);
5148eab0cc1SPeter Brune       } else if (PetscRealPart(pval[j]) <= cls->interp_threshold * pmax_neg) {
515586a8384SPeter Brune         pthresh_neg += PetscRealPart(pval[j]);
516586a8384SPeter Brune       }
517586a8384SPeter Brune       if (PetscRealPart(pval[j]) > 0.) {
518586a8384SPeter Brune         ptot_pos += PetscRealPart(pval[j]);
519586a8384SPeter Brune       } else {
520586a8384SPeter Brune         ptot_neg += PetscRealPart(pval[j]);
521586a8384SPeter Brune       }
522586a8384SPeter Brune     }
5236bd8ea90SPeter Brune     if (PetscAbsReal(pthresh_pos) > 0.) ptot_pos /= pthresh_pos;
5246bd8ea90SPeter Brune     if (PetscAbsReal(pthresh_neg) > 0.) ptot_neg /= pthresh_neg;
525586a8384SPeter Brune     idx = 0;
526586a8384SPeter Brune     for (j = 0; j < ncols; j++) {
5278eab0cc1SPeter Brune       if (PetscRealPart(pval[j]) >= pmax_pos * cls->interp_threshold) {
528586a8384SPeter Brune         pnval[idx] = ptot_pos * pval[j];
529586a8384SPeter Brune         pncol[idx] = pcol[j];
530586a8384SPeter Brune         idx++;
5318eab0cc1SPeter Brune       } else if (PetscRealPart(pval[j]) <= pmax_neg * cls->interp_threshold) {
532586a8384SPeter Brune         pnval[idx] = ptot_neg * pval[j];
533586a8384SPeter Brune         pncol[idx] = pcol[j];
534586a8384SPeter Brune         idx++;
535586a8384SPeter Brune       }
536586a8384SPeter Brune     }
5379566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(*P, i, &ncols, &pcol, &pval));
5389566063dSJacob Faibussowitsch     PetscCall(MatSetValues(Pnew, 1, &i, idx, pncol, pnval, INSERT_VALUES));
539586a8384SPeter Brune   }
540586a8384SPeter Brune 
5419566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Pnew, MAT_FINAL_ASSEMBLY));
5429566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Pnew, MAT_FINAL_ASSEMBLY));
5439566063dSJacob Faibussowitsch   PetscCall(MatDestroy(P));
544586a8384SPeter Brune 
545586a8384SPeter Brune   *P = Pnew;
5469566063dSJacob Faibussowitsch   PetscCall(PetscFree2(lsparse, gsparse));
5479566063dSJacob Faibussowitsch   PetscCall(PetscFree2(pnval, pncol));
5483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
549586a8384SPeter Brune }
550586a8384SPeter Brune 
551*ff6a9541SJacob Faibussowitsch static PetscErrorCode PCGAMGProlongator_Classical_Standard(PC pc, Mat A, Mat G, PetscCoarsenData *agg_lists, Mat *P)
552d71ae5a4SJacob Faibussowitsch {
553c634539dSPeter Brune   Mat                lA, *lAs;
554f9a65ec8SPeter Brune   MatType            mtype;
55563b0c588SPeter Brune   Vec                cv;
55663b0c588SPeter Brune   PetscInt          *gcid, *lcid, *lsparse, *gsparse, *picol;
557420cd23fSPeter Brune   PetscInt           fs, fe, cs, ce, nl, i, j, k, li, lni, ci, ncols, maxcols, fn, cn, cid;
558420cd23fSPeter Brune   PetscMPIInt        size;
55963b0c588SPeter Brune   const PetscInt    *lidx, *icol, *gidx;
56063b0c588SPeter Brune   PetscBool          iscoarse;
56163b0c588SPeter Brune   PetscScalar        vi, pentry, pjentry;
56263b0c588SPeter Brune   PetscScalar       *pcontrib, *pvcol;
56363b0c588SPeter Brune   const PetscScalar *vcol;
564ed4e82eaSPeter Brune   PetscReal          diag, jdiag, jwttotal;
565f9a65ec8SPeter Brune   PetscInt           pncols;
566c634539dSPeter Brune   PetscSF            sf;
567c634539dSPeter Brune   PetscLayout        clayout;
56863b0c588SPeter Brune   IS                 lis;
569f9a65ec8SPeter Brune 
570f9a65ec8SPeter Brune   PetscFunctionBegin;
5719566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
5729566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &fs, &fe));
573f9a65ec8SPeter Brune   fn = fe - fs;
5749566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF, fe - fs, fs, 1, &lis));
575c634539dSPeter Brune   if (size > 1) {
5769566063dSJacob Faibussowitsch     PetscCall(MatGetLayouts(A, NULL, &clayout));
577f9a65ec8SPeter Brune     /* increase the overlap by two to get neighbors of neighbors */
5789566063dSJacob Faibussowitsch     PetscCall(MatIncreaseOverlap(A, 1, &lis, 2));
5799566063dSJacob Faibussowitsch     PetscCall(ISSort(lis));
580f9a65ec8SPeter Brune     /* get the local part of A */
5819566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatrices(A, 1, &lis, &lis, MAT_INITIAL_MATRIX, &lAs));
582c634539dSPeter Brune     lA = lAs[0];
583c634539dSPeter Brune     /* build an SF out of it */
5849566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(lis, &nl));
5859566063dSJacob Faibussowitsch     PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)A), &sf));
5869566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(lis, &lidx));
5879566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraphLayout(sf, clayout, nl, NULL, PETSC_COPY_VALUES, lidx));
5889566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(lis, &lidx));
589c634539dSPeter Brune   } else {
590c634539dSPeter Brune     lA = A;
591c634539dSPeter Brune     nl = fn;
592c634539dSPeter Brune   }
593c634539dSPeter Brune   /* create a communication structure for the overlapped portion and transmit coarse indices */
5949566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(fn, &lsparse, fn, &gsparse, nl, &pcontrib));
595f9a65ec8SPeter Brune   /* create coarse vector */
596f9a65ec8SPeter Brune   cn = 0;
597f9a65ec8SPeter Brune   for (i = 0; i < fn; i++) {
5989566063dSJacob Faibussowitsch     PetscCall(PetscCDEmptyAt(agg_lists, i, &iscoarse));
599ad540459SPierre Jolivet     if (!iscoarse) cn++;
600f9a65ec8SPeter Brune   }
6019566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(fn, &gcid));
6029566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)A), cn, PETSC_DECIDE, &cv));
6039566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(cv, &cs, &ce));
604f9a65ec8SPeter Brune   cn = 0;
605f9a65ec8SPeter Brune   for (i = 0; i < fn; i++) {
6069566063dSJacob Faibussowitsch     PetscCall(PetscCDEmptyAt(agg_lists, i, &iscoarse));
607f9a65ec8SPeter Brune     if (!iscoarse) {
608c634539dSPeter Brune       gcid[i] = cs + cn;
609f9a65ec8SPeter Brune       cn++;
610f9a65ec8SPeter Brune     } else {
611c634539dSPeter Brune       gcid[i] = -1;
612f9a65ec8SPeter Brune     }
613f9a65ec8SPeter Brune   }
614c634539dSPeter Brune   if (size > 1) {
6159566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nl, &lcid));
6169566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(sf, MPIU_INT, gcid, lcid, MPI_REPLACE));
6179566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(sf, MPIU_INT, gcid, lcid, MPI_REPLACE));
618c634539dSPeter Brune   } else {
619c634539dSPeter Brune     lcid = gcid;
620c634539dSPeter Brune   }
621f9a65ec8SPeter Brune   /* count to preallocate the prolongator */
6229566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(lis, &gidx));
623f9a65ec8SPeter Brune   maxcols = 0;
624f9a65ec8SPeter Brune   /* count the number of unique contributing coarse cells for each fine */
625f9a65ec8SPeter Brune   for (i = 0; i < nl; i++) {
626ed4e82eaSPeter Brune     pcontrib[i] = 0.;
6279566063dSJacob Faibussowitsch     PetscCall(MatGetRow(lA, i, &ncols, &icol, NULL));
628f9a65ec8SPeter Brune     if (gidx[i] >= fs && gidx[i] < fe) {
629f9a65ec8SPeter Brune       li          = gidx[i] - fs;
630f9a65ec8SPeter Brune       lsparse[li] = 0;
631f9a65ec8SPeter Brune       gsparse[li] = 0;
632c634539dSPeter Brune       cid         = lcid[i];
633f9a65ec8SPeter Brune       if (cid >= 0) {
634f9a65ec8SPeter Brune         lsparse[li] = 1;
635f9a65ec8SPeter Brune       } else {
636f9a65ec8SPeter Brune         for (j = 0; j < ncols; j++) {
637c634539dSPeter Brune           if (lcid[icol[j]] >= 0) {
638f9a65ec8SPeter Brune             pcontrib[icol[j]] = 1.;
639f9a65ec8SPeter Brune           } else {
640f9a65ec8SPeter Brune             ci = icol[j];
6419566063dSJacob Faibussowitsch             PetscCall(MatRestoreRow(lA, i, &ncols, &icol, NULL));
6429566063dSJacob Faibussowitsch             PetscCall(MatGetRow(lA, ci, &ncols, &icol, NULL));
643f9a65ec8SPeter Brune             for (k = 0; k < ncols; k++) {
644ad540459SPierre Jolivet               if (lcid[icol[k]] >= 0) pcontrib[icol[k]] = 1.;
645f9a65ec8SPeter Brune             }
6469566063dSJacob Faibussowitsch             PetscCall(MatRestoreRow(lA, ci, &ncols, &icol, NULL));
6479566063dSJacob Faibussowitsch             PetscCall(MatGetRow(lA, i, &ncols, &icol, NULL));
648f9a65ec8SPeter Brune           }
649f9a65ec8SPeter Brune         }
650f9a65ec8SPeter Brune         for (j = 0; j < ncols; j++) {
651c634539dSPeter Brune           if (lcid[icol[j]] >= 0 && pcontrib[icol[j]] != 0.) {
652c634539dSPeter Brune             lni = lcid[icol[j]];
653f9a65ec8SPeter Brune             if (lni >= cs && lni < ce) {
654f9a65ec8SPeter Brune               lsparse[li]++;
655f9a65ec8SPeter Brune             } else {
656f9a65ec8SPeter Brune               gsparse[li]++;
657f9a65ec8SPeter Brune             }
658f9a65ec8SPeter Brune             pcontrib[icol[j]] = 0.;
659f9a65ec8SPeter Brune           } else {
660f9a65ec8SPeter Brune             ci = icol[j];
6619566063dSJacob Faibussowitsch             PetscCall(MatRestoreRow(lA, i, &ncols, &icol, NULL));
6629566063dSJacob Faibussowitsch             PetscCall(MatGetRow(lA, ci, &ncols, &icol, NULL));
663f9a65ec8SPeter Brune             for (k = 0; k < ncols; k++) {
664c634539dSPeter Brune               if (lcid[icol[k]] >= 0 && pcontrib[icol[k]] != 0.) {
665c634539dSPeter Brune                 lni = lcid[icol[k]];
666f9a65ec8SPeter Brune                 if (lni >= cs && lni < ce) {
667f9a65ec8SPeter Brune                   lsparse[li]++;
668f9a65ec8SPeter Brune                 } else {
669f9a65ec8SPeter Brune                   gsparse[li]++;
670f9a65ec8SPeter Brune                 }
671f9a65ec8SPeter Brune                 pcontrib[icol[k]] = 0.;
672f9a65ec8SPeter Brune               }
673f9a65ec8SPeter Brune             }
6749566063dSJacob Faibussowitsch             PetscCall(MatRestoreRow(lA, ci, &ncols, &icol, NULL));
6759566063dSJacob Faibussowitsch             PetscCall(MatGetRow(lA, i, &ncols, &icol, NULL));
676f9a65ec8SPeter Brune           }
677f9a65ec8SPeter Brune         }
678ed4e82eaSPeter Brune       }
679f9a65ec8SPeter Brune       if (lsparse[li] + gsparse[li] > maxcols) maxcols = lsparse[li] + gsparse[li];
680ed4e82eaSPeter Brune     }
6819566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(lA, i, &ncols, &icol, &vcol));
682f9a65ec8SPeter Brune   }
6839566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(maxcols, &picol, maxcols, &pvcol));
6849566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), P));
6859566063dSJacob Faibussowitsch   PetscCall(MatGetType(A, &mtype));
6869566063dSJacob Faibussowitsch   PetscCall(MatSetType(*P, mtype));
6879566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*P, fn, cn, PETSC_DETERMINE, PETSC_DETERMINE));
6889566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(*P, 0, lsparse, 0, gsparse));
6899566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(*P, 0, lsparse));
690f9a65ec8SPeter Brune   for (i = 0; i < nl; i++) {
691ed4e82eaSPeter Brune     diag = 0.;
692f9a65ec8SPeter Brune     if (gidx[i] >= fs && gidx[i] < fe) {
693f9a65ec8SPeter Brune       pncols = 0;
694c634539dSPeter Brune       cid    = lcid[i];
695f9a65ec8SPeter Brune       if (cid >= 0) {
696f9a65ec8SPeter Brune         pncols   = 1;
697f9a65ec8SPeter Brune         picol[0] = cid;
698f9a65ec8SPeter Brune         pvcol[0] = 1.;
699f9a65ec8SPeter Brune       } else {
7009566063dSJacob Faibussowitsch         PetscCall(MatGetRow(lA, i, &ncols, &icol, &vcol));
701f9a65ec8SPeter Brune         for (j = 0; j < ncols; j++) {
702ed4e82eaSPeter Brune           pentry = vcol[j];
703c634539dSPeter Brune           if (lcid[icol[j]] >= 0) {
704f9a65ec8SPeter Brune             /* coarse neighbor */
705ed4e82eaSPeter Brune             pcontrib[icol[j]] += pentry;
706ed4e82eaSPeter Brune           } else if (icol[j] != i) {
707f9a65ec8SPeter Brune             /* the neighbor is a strongly connected fine node */
708f9a65ec8SPeter Brune             ci = icol[j];
709f9a65ec8SPeter Brune             vi = vcol[j];
7109566063dSJacob Faibussowitsch             PetscCall(MatRestoreRow(lA, i, &ncols, &icol, &vcol));
7119566063dSJacob Faibussowitsch             PetscCall(MatGetRow(lA, ci, &ncols, &icol, &vcol));
712ed4e82eaSPeter Brune             jwttotal = 0.;
713f9a65ec8SPeter Brune             jdiag    = 0.;
714f9a65ec8SPeter Brune             for (k = 0; k < ncols; k++) {
715ad540459SPierre Jolivet               if (ci == icol[k]) jdiag = PetscRealPart(vcol[k]);
716f9a65ec8SPeter Brune             }
717f9a65ec8SPeter Brune             for (k = 0; k < ncols; k++) {
718c634539dSPeter Brune               if (lcid[icol[k]] >= 0 && jdiag * PetscRealPart(vcol[k]) < 0.) {
719ed4e82eaSPeter Brune                 pjentry = vcol[k];
7207779008dSPeter Brune                 jwttotal += PetscRealPart(pjentry);
721f9a65ec8SPeter Brune               }
722f9a65ec8SPeter Brune             }
723ed4e82eaSPeter Brune             if (jwttotal != 0.) {
7247779008dSPeter Brune               jwttotal = PetscRealPart(vi) / jwttotal;
725ed4e82eaSPeter Brune               for (k = 0; k < ncols; k++) {
726c634539dSPeter Brune                 if (lcid[icol[k]] >= 0 && jdiag * PetscRealPart(vcol[k]) < 0.) {
727586a8384SPeter Brune                   pjentry = vcol[k] * jwttotal;
728ed4e82eaSPeter Brune                   pcontrib[icol[k]] += pjentry;
729ed4e82eaSPeter Brune                 }
730ed4e82eaSPeter Brune               }
731ed4e82eaSPeter Brune             } else {
732ed4e82eaSPeter Brune               diag += PetscRealPart(vi);
733ed4e82eaSPeter Brune             }
7349566063dSJacob Faibussowitsch             PetscCall(MatRestoreRow(lA, ci, &ncols, &icol, &vcol));
7359566063dSJacob Faibussowitsch             PetscCall(MatGetRow(lA, i, &ncols, &icol, &vcol));
736ed4e82eaSPeter Brune           } else {
737ed4e82eaSPeter Brune             diag += PetscRealPart(vcol[j]);
738f9a65ec8SPeter Brune           }
739f9a65ec8SPeter Brune         }
740586a8384SPeter Brune         if (diag != 0.) {
741586a8384SPeter Brune           diag = 1. / diag;
742f9a65ec8SPeter Brune           for (j = 0; j < ncols; j++) {
743c634539dSPeter Brune             if (lcid[icol[j]] >= 0 && pcontrib[icol[j]] != 0.) {
744f9a65ec8SPeter Brune               /* the neighbor is a coarse node */
745ed4e82eaSPeter Brune               if (PetscAbsScalar(pcontrib[icol[j]]) > 0.0) {
746c634539dSPeter Brune                 lni           = lcid[icol[j]];
747586a8384SPeter Brune                 pvcol[pncols] = -pcontrib[icol[j]] * diag;
748f9a65ec8SPeter Brune                 picol[pncols] = lni;
749f9a65ec8SPeter Brune                 pncols++;
750ed4e82eaSPeter Brune               }
751ed4e82eaSPeter Brune               pcontrib[icol[j]] = 0.;
752f9a65ec8SPeter Brune             } else {
753f9a65ec8SPeter Brune               /* the neighbor is a strongly connected fine node */
754f9a65ec8SPeter Brune               ci = icol[j];
7559566063dSJacob Faibussowitsch               PetscCall(MatRestoreRow(lA, i, &ncols, &icol, &vcol));
7569566063dSJacob Faibussowitsch               PetscCall(MatGetRow(lA, ci, &ncols, &icol, &vcol));
757f9a65ec8SPeter Brune               for (k = 0; k < ncols; k++) {
758c634539dSPeter Brune                 if (lcid[icol[k]] >= 0 && pcontrib[icol[k]] != 0.) {
759ed4e82eaSPeter Brune                   if (PetscAbsScalar(pcontrib[icol[k]]) > 0.0) {
760c634539dSPeter Brune                     lni           = lcid[icol[k]];
761586a8384SPeter Brune                     pvcol[pncols] = -pcontrib[icol[k]] * diag;
762f9a65ec8SPeter Brune                     picol[pncols] = lni;
763f9a65ec8SPeter Brune                     pncols++;
764f9a65ec8SPeter Brune                   }
765ed4e82eaSPeter Brune                   pcontrib[icol[k]] = 0.;
766ed4e82eaSPeter Brune                 }
767f9a65ec8SPeter Brune               }
7689566063dSJacob Faibussowitsch               PetscCall(MatRestoreRow(lA, ci, &ncols, &icol, &vcol));
7699566063dSJacob Faibussowitsch               PetscCall(MatGetRow(lA, i, &ncols, &icol, &vcol));
770f9a65ec8SPeter Brune             }
771ed4e82eaSPeter Brune             pcontrib[icol[j]] = 0.;
772f9a65ec8SPeter Brune           }
7739566063dSJacob Faibussowitsch           PetscCall(MatRestoreRow(lA, i, &ncols, &icol, &vcol));
774f9a65ec8SPeter Brune         }
775586a8384SPeter Brune       }
776f9a65ec8SPeter Brune       ci = gidx[i];
77748a46eb9SPierre Jolivet       if (pncols > 0) PetscCall(MatSetValues(*P, 1, &ci, pncols, picol, pvcol, INSERT_VALUES));
778f9a65ec8SPeter Brune     }
779f9a65ec8SPeter Brune   }
7809566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(lis, &gidx));
7819566063dSJacob Faibussowitsch   PetscCall(PetscFree2(picol, pvcol));
7829566063dSJacob Faibussowitsch   PetscCall(PetscFree3(lsparse, gsparse, pcontrib));
7839566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&lis));
7849566063dSJacob Faibussowitsch   PetscCall(PetscFree(gcid));
785c634539dSPeter Brune   if (size > 1) {
7869566063dSJacob Faibussowitsch     PetscCall(PetscFree(lcid));
7879566063dSJacob Faibussowitsch     PetscCall(MatDestroyMatrices(1, &lAs));
7889566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&sf));
789c634539dSPeter Brune   }
7909566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&cv));
7919566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY));
7929566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY));
7933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7948eab0cc1SPeter Brune }
795f9a65ec8SPeter Brune 
796*ff6a9541SJacob Faibussowitsch static PetscErrorCode PCGAMGOptProlongator_Classical_Jacobi(PC pc, Mat A, Mat *P)
797d71ae5a4SJacob Faibussowitsch {
798b4dc3ebdSPeter Brune   PetscInt           f, s, n, cf, cs, i, idx;
799b4dc3ebdSPeter Brune   PetscInt          *coarserows;
800b4dc3ebdSPeter Brune   PetscInt           ncols;
801b4dc3ebdSPeter Brune   const PetscInt    *pcols;
802b4dc3ebdSPeter Brune   const PetscScalar *pvals;
803b4dc3ebdSPeter Brune   Mat                Pnew;
804b4dc3ebdSPeter Brune   Vec                diag;
805b4dc3ebdSPeter Brune   PC_MG             *mg      = (PC_MG *)pc->data;
806b4dc3ebdSPeter Brune   PC_GAMG           *pc_gamg = (PC_GAMG *)mg->innerctx;
807b4dc3ebdSPeter Brune   PC_GAMG_Classical *cls     = (PC_GAMG_Classical *)pc_gamg->subctx;
808b4dc3ebdSPeter Brune 
809b4dc3ebdSPeter Brune   PetscFunctionBegin;
810b4dc3ebdSPeter Brune   if (cls->nsmooths == 0) {
8119566063dSJacob Faibussowitsch     PetscCall(PCGAMGTruncateProlongator_Private(pc, P));
8123ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
813b4dc3ebdSPeter Brune   }
8149566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(*P, &s, &f));
815b4dc3ebdSPeter Brune   n = f - s;
8169566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRangeColumn(*P, &cs, &cf));
8179566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &coarserows));
818b4dc3ebdSPeter Brune   /* identify the rows corresponding to coarse unknowns */
819b4dc3ebdSPeter Brune   idx = 0;
820b4dc3ebdSPeter Brune   for (i = s; i < f; i++) {
8219566063dSJacob Faibussowitsch     PetscCall(MatGetRow(*P, i, &ncols, &pcols, &pvals));
822b4dc3ebdSPeter Brune     /* assume, for now, that it's a coarse unknown if it has a single unit entry */
823b4dc3ebdSPeter Brune     if (ncols == 1) {
824b4dc3ebdSPeter Brune       if (pvals[0] == 1.) {
825b4dc3ebdSPeter Brune         coarserows[idx] = i;
826b4dc3ebdSPeter Brune         idx++;
827b4dc3ebdSPeter Brune       }
828b4dc3ebdSPeter Brune     }
8299566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(*P, i, &ncols, &pcols, &pvals));
830b4dc3ebdSPeter Brune   }
8319566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &diag, NULL));
8329566063dSJacob Faibussowitsch   PetscCall(MatGetDiagonal(A, diag));
8339566063dSJacob Faibussowitsch   PetscCall(VecReciprocal(diag));
834b4dc3ebdSPeter Brune   for (i = 0; i < cls->nsmooths; i++) {
8359566063dSJacob Faibussowitsch     PetscCall(MatMatMult(A, *P, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Pnew));
8369566063dSJacob Faibussowitsch     PetscCall(MatZeroRows(Pnew, idx, coarserows, 0., NULL, NULL));
8379566063dSJacob Faibussowitsch     PetscCall(MatDiagonalScale(Pnew, diag, NULL));
8389566063dSJacob Faibussowitsch     PetscCall(MatAYPX(Pnew, -1.0, *P, DIFFERENT_NONZERO_PATTERN));
8399566063dSJacob Faibussowitsch     PetscCall(MatDestroy(P));
840b4dc3ebdSPeter Brune     *P   = Pnew;
841b4dc3ebdSPeter Brune     Pnew = NULL;
842b4dc3ebdSPeter Brune   }
8439566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&diag));
8449566063dSJacob Faibussowitsch   PetscCall(PetscFree(coarserows));
8459566063dSJacob Faibussowitsch   PetscCall(PCGAMGTruncateProlongator_Private(pc, P));
8463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
847b4dc3ebdSPeter Brune }
848b4dc3ebdSPeter Brune 
849d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGProlongator_Classical(PC pc, Mat A, Mat G, PetscCoarsenData *agg_lists, Mat *P)
850d71ae5a4SJacob Faibussowitsch {
8518eab0cc1SPeter Brune   PetscErrorCode (*f)(PC, Mat, Mat, PetscCoarsenData *, Mat *);
8528eab0cc1SPeter Brune   PC_MG             *mg      = (PC_MG *)pc->data;
8538eab0cc1SPeter Brune   PC_GAMG           *pc_gamg = (PC_GAMG *)mg->innerctx;
8548eab0cc1SPeter Brune   PC_GAMG_Classical *cls     = (PC_GAMG_Classical *)pc_gamg->subctx;
8558eab0cc1SPeter Brune 
8568eab0cc1SPeter Brune   PetscFunctionBegin;
8579566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(PCGAMGClassicalProlongatorList, cls->prolongtype, &f));
85828b400f6SJacob Faibussowitsch   PetscCheck(f, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Cannot find PCGAMG Classical prolongator type");
8599566063dSJacob Faibussowitsch   PetscCall((*f)(pc, A, G, agg_lists, P));
8603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
861f9a65ec8SPeter Brune }
862f9a65ec8SPeter Brune 
863d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGDestroy_Classical(PC pc)
864d71ae5a4SJacob Faibussowitsch {
8658e6d0c30SPeter Brune   PC_MG   *mg      = (PC_MG *)pc->data;
8668e6d0c30SPeter Brune   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
8678e6d0c30SPeter Brune 
8688e6d0c30SPeter Brune   PetscFunctionBegin;
8699566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg->subctx));
8709566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGClassicalSetType_C", NULL));
8719566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGClassicalGetType_C", NULL));
8723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8738e6d0c30SPeter Brune }
8748e6d0c30SPeter Brune 
875*ff6a9541SJacob Faibussowitsch static PetscErrorCode PCGAMGSetFromOptions_Classical(PC pc, PetscOptionItems *PetscOptionsObject)
876d71ae5a4SJacob Faibussowitsch {
877586a8384SPeter Brune   PC_MG             *mg      = (PC_MG *)pc->data;
878586a8384SPeter Brune   PC_GAMG           *pc_gamg = (PC_GAMG *)mg->innerctx;
879586a8384SPeter Brune   PC_GAMG_Classical *cls     = (PC_GAMG_Classical *)pc_gamg->subctx;
8808eab0cc1SPeter Brune   char               tname[256];
8818eab0cc1SPeter Brune   PetscBool          flg;
882586a8384SPeter Brune 
8838e6d0c30SPeter Brune   PetscFunctionBegin;
884d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "GAMG-Classical options");
8859566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-pc_gamg_classical_type", "Type of Classical AMG prolongation", "PCGAMGClassicalSetType", PCGAMGClassicalProlongatorList, cls->prolongtype, tname, sizeof(tname), &flg));
8861baa6e33SBarry Smith   if (flg) PetscCall(PCGAMGClassicalSetType(pc, tname));
8879566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_gamg_classical_interp_threshold", "Threshold for classical interpolator entries", "", cls->interp_threshold, &cls->interp_threshold, NULL));
8889566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_gamg_classical_nsmooths", "Threshold for classical interpolator entries", "", cls->nsmooths, &cls->nsmooths, NULL));
889d0609cedSBarry Smith   PetscOptionsHeadEnd();
8903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8918e6d0c30SPeter Brune }
8928e6d0c30SPeter Brune 
893d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetData_Classical(PC pc, Mat A)
894d71ae5a4SJacob Faibussowitsch {
8958e6d0c30SPeter Brune   PC_MG   *mg      = (PC_MG *)pc->data;
8968e6d0c30SPeter Brune   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
8978e6d0c30SPeter Brune 
8988e6d0c30SPeter Brune   PetscFunctionBegin;
8998e6d0c30SPeter Brune   /* no data for classical AMG */
9008e6d0c30SPeter Brune   pc_gamg->data           = NULL;
901d2050638SMark Adams   pc_gamg->data_cell_cols = 0;
902d2050638SMark Adams   pc_gamg->data_cell_rows = 0;
9038e6d0c30SPeter Brune   pc_gamg->data_sz        = 0;
9043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9058e6d0c30SPeter Brune }
9068e6d0c30SPeter Brune 
907*ff6a9541SJacob Faibussowitsch static PetscErrorCode PCGAMGClassicalFinalizePackage(void)
908d71ae5a4SJacob Faibussowitsch {
9098eab0cc1SPeter Brune   PetscFunctionBegin;
9108eab0cc1SPeter Brune   PCGAMGClassicalPackageInitialized = PETSC_FALSE;
9119566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDestroy(&PCGAMGClassicalProlongatorList));
9123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9138eab0cc1SPeter Brune }
9148eab0cc1SPeter Brune 
915*ff6a9541SJacob Faibussowitsch static PetscErrorCode PCGAMGClassicalInitializePackage(void)
916d71ae5a4SJacob Faibussowitsch {
9178eab0cc1SPeter Brune   PetscFunctionBegin;
9183ba16761SJacob Faibussowitsch   if (PCGAMGClassicalPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS);
9199566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&PCGAMGClassicalProlongatorList, PCGAMGCLASSICALDIRECT, PCGAMGProlongator_Classical_Direct));
9209566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&PCGAMGClassicalProlongatorList, PCGAMGCLASSICALSTANDARD, PCGAMGProlongator_Classical_Standard));
9219566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalize(PCGAMGClassicalFinalizePackage));
9223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9238eab0cc1SPeter Brune }
9248eab0cc1SPeter Brune 
9258e6d0c30SPeter Brune /*
9268e6d0c30SPeter Brune    PCCreateGAMG_Classical
9278e6d0c30SPeter Brune 
9288e6d0c30SPeter Brune */
929d71ae5a4SJacob Faibussowitsch PetscErrorCode PCCreateGAMG_Classical(PC pc)
930d71ae5a4SJacob Faibussowitsch {
9318e6d0c30SPeter Brune   PC_MG             *mg      = (PC_MG *)pc->data;
9328e6d0c30SPeter Brune   PC_GAMG           *pc_gamg = (PC_GAMG *)mg->innerctx;
9338e6d0c30SPeter Brune   PC_GAMG_Classical *pc_gamg_classical;
9348e6d0c30SPeter Brune 
9358e6d0c30SPeter Brune   PetscFunctionBegin;
9369566063dSJacob Faibussowitsch   PetscCall(PCGAMGClassicalInitializePackage());
9378e6d0c30SPeter Brune   if (pc_gamg->subctx) {
9388e6d0c30SPeter Brune     /* call base class */
9399566063dSJacob Faibussowitsch     PetscCall(PCDestroy_GAMG(pc));
9408e6d0c30SPeter Brune   }
9418e6d0c30SPeter Brune 
9428e6d0c30SPeter Brune   /* create sub context for SA */
9434dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&pc_gamg_classical));
9448e6d0c30SPeter Brune   pc_gamg->subctx         = pc_gamg_classical;
9458e6d0c30SPeter Brune   pc->ops->setfromoptions = PCGAMGSetFromOptions_Classical;
9468e6d0c30SPeter Brune   /* reset does not do anything; setup not virtual */
9478e6d0c30SPeter Brune 
9488e6d0c30SPeter Brune   /* set internal function pointers */
9498e6d0c30SPeter Brune   pc_gamg->ops->destroy        = PCGAMGDestroy_Classical;
9502d776b49SBarry Smith   pc_gamg->ops->creategraph    = PCGAMGCreateGraph_Classical;
9518e6d0c30SPeter Brune   pc_gamg->ops->coarsen        = PCGAMGCoarsen_Classical;
9528eab0cc1SPeter Brune   pc_gamg->ops->prolongator    = PCGAMGProlongator_Classical;
953fd1112cbSBarry Smith   pc_gamg->ops->optprolongator = PCGAMGOptProlongator_Classical_Jacobi;
954586a8384SPeter Brune   pc_gamg->ops->setfromoptions = PCGAMGSetFromOptions_Classical;
9558e6d0c30SPeter Brune 
9568e6d0c30SPeter Brune   pc_gamg->ops->createdefaultdata     = PCGAMGSetData_Classical;
957586a8384SPeter Brune   pc_gamg_classical->interp_threshold = 0.2;
958b4dc3ebdSPeter Brune   pc_gamg_classical->nsmooths         = 0;
9599566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGClassicalSetType_C", PCGAMGClassicalSetType_GAMG));
9609566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGClassicalGetType_C", PCGAMGClassicalGetType_GAMG));
9619566063dSJacob Faibussowitsch   PetscCall(PCGAMGClassicalSetType(pc, PCGAMGCLASSICALSTANDARD));
9623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9638e6d0c30SPeter Brune }
964