xref: /petsc/src/ksp/pc/impls/gamg/agg.c (revision c3339decea92175325d9368fa13196bcd0e0e58b)
12e68589bSMark F. Adams /*
2b817416eSBarry Smith  GAMG geometric-algebric multigrid PC - Mark Adams 2011
32e68589bSMark F. Adams  */
42e68589bSMark F. Adams 
52e68589bSMark F. Adams #include <../src/ksp/pc/impls/gamg/gamg.h> /*I "petscpc.h" I*/
62e68589bSMark F. Adams #include <petscblaslapack.h>
7539c167fSBarry Smith #include <petscdm.h>
873f7197eSJed Brown #include <petsc/private/kspimpl.h>
92e68589bSMark F. Adams 
102e68589bSMark F. Adams typedef struct {
11c8b0795cSMark F. Adams   PetscInt   nsmooths;
12bae903cbSmarkadams4   PetscInt   aggressive_coarsening_levels; // number of aggressive coarsening levels (square or MISk)
132d776b49SBarry Smith   MatCoarsen crs;
142e68589bSMark F. Adams } PC_GAMG_AGG;
152e68589bSMark F. Adams 
162e68589bSMark F. Adams /*@
17f1580f4eSBarry Smith    PCGAMGSetNSmooths - Set number of smoothing steps (1 is typical) used for multigrid on all the levels
182e68589bSMark F. Adams 
19*c3339decSBarry Smith    Logically Collective
202e68589bSMark F. Adams 
212e68589bSMark F. Adams    Input Parameters:
222e68589bSMark F. Adams .  pc - the preconditioner context
232e68589bSMark F. Adams 
242e68589bSMark F. Adams    Options Database Key:
25cab9ed1eSBarry Smith .  -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation
262e68589bSMark F. Adams 
272e68589bSMark F. Adams    Level: intermediate
282e68589bSMark F. Adams 
29f1580f4eSBarry Smith .seealso: `PCMG`, `PCGAMG`
302e68589bSMark F. Adams @*/
31d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n)
32d71ae5a4SJacob Faibussowitsch {
332e68589bSMark F. Adams   PetscFunctionBegin;
342e68589bSMark F. Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
35d5d25401SBarry Smith   PetscValidLogicalCollectiveInt(pc, n, 2);
36cac4c232SBarry Smith   PetscTryMethod(pc, "PCGAMGSetNSmooths_C", (PC, PetscInt), (pc, n));
372e68589bSMark F. Adams   PetscFunctionReturn(0);
382e68589bSMark F. Adams }
392e68589bSMark F. Adams 
40d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetNSmooths_AGG(PC pc, PetscInt n)
41d71ae5a4SJacob Faibussowitsch {
422e68589bSMark F. Adams   PC_MG       *mg          = (PC_MG *)pc->data;
432e68589bSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
44c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
452e68589bSMark F. Adams 
462e68589bSMark F. Adams   PetscFunctionBegin;
47c8b0795cSMark F. Adams   pc_gamg_agg->nsmooths = n;
48c8b0795cSMark F. Adams   PetscFunctionReturn(0);
49c8b0795cSMark F. Adams }
50c8b0795cSMark F. Adams 
51c8b0795cSMark F. Adams /*@
52f1580f4eSBarry Smith    PCGAMGSetAggressiveLevels -  Use aggressive coarsening on first n levels
53ef4ad70eSMark F. Adams 
54*c3339decSBarry Smith    Logically Collective
55ef4ad70eSMark F. Adams 
56ef4ad70eSMark F. Adams    Input Parameters:
57cab9ed1eSBarry Smith +  pc - the preconditioner context
58d5d25401SBarry Smith -  n - 0, 1 or more
59ef4ad70eSMark F. Adams 
60ef4ad70eSMark F. Adams    Options Database Key:
61bae903cbSmarkadams4 .  -pc_gamg_aggressive_coarsening <n,default = 1> - Number of levels to square the graph on before aggregating it
62af3c827dSMark Adams 
63ef4ad70eSMark F. Adams    Level: intermediate
64ef4ad70eSMark F. Adams 
652d776b49SBarry Smith .seealso: `PCGAMG`, `PCGAMGSetThreshold()`
66ef4ad70eSMark F. Adams @*/
67d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetAggressiveLevels(PC pc, PetscInt n)
68d71ae5a4SJacob Faibussowitsch {
69ef4ad70eSMark F. Adams   PetscFunctionBegin;
70ef4ad70eSMark F. Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
71d5d25401SBarry Smith   PetscValidLogicalCollectiveInt(pc, n, 2);
72bae903cbSmarkadams4   PetscTryMethod(pc, "PCGAMGSetAggressiveLevels_C", (PC, PetscInt), (pc, n));
73ef4ad70eSMark F. Adams   PetscFunctionReturn(0);
74ef4ad70eSMark F. Adams }
75ef4ad70eSMark F. Adams 
76d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetAggressiveLevels_AGG(PC pc, PetscInt n)
77d71ae5a4SJacob Faibussowitsch {
78ef4ad70eSMark F. Adams   PC_MG       *mg          = (PC_MG *)pc->data;
79ef4ad70eSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
80ef4ad70eSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
81ef4ad70eSMark F. Adams 
82ef4ad70eSMark F. Adams   PetscFunctionBegin;
83bae903cbSmarkadams4   pc_gamg_agg->aggressive_coarsening_levels = n;
84ef4ad70eSMark F. Adams   PetscFunctionReturn(0);
85ef4ad70eSMark F. Adams }
86ef4ad70eSMark F. Adams 
87d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_GAMG_AGG(PC pc, PetscOptionItems *PetscOptionsObject)
88d71ae5a4SJacob Faibussowitsch {
892e68589bSMark F. Adams   PC_MG       *mg          = (PC_MG *)pc->data;
902e68589bSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
91c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
922e68589bSMark F. Adams 
932e68589bSMark F. Adams   PetscFunctionBegin;
94d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "GAMG-AGG options");
952e68589bSMark F. Adams   {
96bae903cbSmarkadams4     PetscBool flg;
979566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-pc_gamg_agg_nsmooths", "smoothing steps for smoothed aggregation, usually 1", "PCGAMGSetNSmooths", pc_gamg_agg->nsmooths, &pc_gamg_agg->nsmooths, NULL));
98bae903cbSmarkadams4     pc_gamg_agg->aggressive_coarsening_levels = 1;
999371c9d4SSatish Balay     PetscCall(
1009371c9d4SSatish Balay       PetscOptionsInt("-pc_gamg_square_graph", "Number of aggressive coarsening (MIS-2) levels from finest (alias for -pc_gamg_aggressive_coarsening, deprecated)", "PCGAMGSetAggressiveLevels", pc_gamg_agg->aggressive_coarsening_levels, &pc_gamg_agg->aggressive_coarsening_levels, &flg));
101bae903cbSmarkadams4     if (!flg) {
102bae903cbSmarkadams4       PetscCall(PetscOptionsInt("-pc_gamg_aggressive_coarsening", "Number of aggressive coarsening (MIS-2) levels from finest", "PCGAMGSetAggressiveLevels", pc_gamg_agg->aggressive_coarsening_levels, &pc_gamg_agg->aggressive_coarsening_levels, NULL));
103bae903cbSmarkadams4     } else {
104bae903cbSmarkadams4       PetscCall(PetscOptionsInt("-pc_gamg_aggressive_coarsening", "Number of aggressive coarsening (MIS-2) levels from finest", "PCGAMGSetAggressiveLevels", pc_gamg_agg->aggressive_coarsening_levels, &pc_gamg_agg->aggressive_coarsening_levels, &flg));
105bae903cbSmarkadams4       if (flg) PetscCall(PetscInfo(pc, "Warning: both -pc_gamg_square_graph and -pc_gamg_aggressive_coarsening are used. -pc_gamg_square_graph is deprecated, Number of aggressive levels is %d\n", (int)pc_gamg_agg->aggressive_coarsening_levels));
106bae903cbSmarkadams4     }
1072e68589bSMark F. Adams   }
108d0609cedSBarry Smith   PetscOptionsHeadEnd();
1092e68589bSMark F. Adams   PetscFunctionReturn(0);
1102e68589bSMark F. Adams }
1112e68589bSMark F. Adams 
112d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_GAMG_AGG(PC pc)
113d71ae5a4SJacob Faibussowitsch {
1142e68589bSMark F. Adams   PC_MG   *mg      = (PC_MG *)pc->data;
1152e68589bSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1162e68589bSMark F. Adams 
1172e68589bSMark F. Adams   PetscFunctionBegin;
1189566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg->subctx));
1192e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNSmooths_C", NULL));
120bae903cbSmarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveLevels_C", NULL));
121bae903cbSmarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", NULL));
1222e68589bSMark F. Adams   PetscFunctionReturn(0);
1232e68589bSMark F. Adams }
1242e68589bSMark F. Adams 
1252e68589bSMark F. Adams /*
1262e68589bSMark F. Adams    PCSetCoordinates_AGG
127302f38e8SMark F. Adams      - collective
1282e68589bSMark F. Adams 
1292e68589bSMark F. Adams    Input Parameter:
1302e68589bSMark F. Adams    . pc - the preconditioner context
131a2f3521dSMark F. Adams    . ndm - dimesion of data (used for dof/vertex for Stokes)
132302f38e8SMark F. Adams    . a_nloc - number of vertices local
133302f38e8SMark F. Adams    . coords - [a_nloc][ndm] - interleaved coordinate data: {x_0, y_0, z_0, x_1, y_1, ...}
1342e68589bSMark F. Adams */
1351e6b0712SBarry Smith 
136d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetCoordinates_AGG(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords)
137d71ae5a4SJacob Faibussowitsch {
1382e68589bSMark F. Adams   PC_MG   *mg      = (PC_MG *)pc->data;
1392e68589bSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
14069344418SMark F. Adams   PetscInt arrsz, kk, ii, jj, nloc, ndatarows, ndf;
141a2f3521dSMark F. Adams   Mat      mat = pc->pmat;
1422e68589bSMark F. Adams 
1432e68589bSMark F. Adams   PetscFunctionBegin;
144a2f3521dSMark F. Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
145a2f3521dSMark F. Adams   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
146302f38e8SMark F. Adams   nloc = a_nloc;
1472e68589bSMark F. Adams 
1482e68589bSMark F. Adams   /* SA: null space vectors */
1499566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(mat, &ndf));               /* this does not work for Stokes */
15069344418SMark F. Adams   if (coords && ndf == 1) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */
151a2f3521dSMark F. Adams   else if (coords) {
15263a3b9bcSJacob Faibussowitsch     PetscCheck(ndm <= ndf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "degrees of motion %" PetscInt_FMT " > block size %" PetscInt_FMT, ndm, ndf);
15369344418SMark F. Adams     pc_gamg->data_cell_cols = (ndm == 2 ? 3 : 6); /* displacement elasticity */
154ad540459SPierre Jolivet     if (ndm != ndf) PetscCheck(pc_gamg->data_cell_cols == ndf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Don't know how to create null space for ndm=%" PetscInt_FMT ", ndf=%" PetscInt_FMT ".  Use MatSetNearNullSpace().", ndm, ndf);
155806fa848SBarry Smith   } else pc_gamg->data_cell_cols = ndf; /* no data, force SA with constant null space vectors */
15671959b99SBarry Smith   pc_gamg->data_cell_rows = ndatarows = ndf;
15763a3b9bcSJacob Faibussowitsch   PetscCheck(pc_gamg->data_cell_cols > 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "pc_gamg->data_cell_cols %" PetscInt_FMT " <= 0", pc_gamg->data_cell_cols);
158c8b0795cSMark F. Adams   arrsz = nloc * pc_gamg->data_cell_rows * pc_gamg->data_cell_cols;
1592e68589bSMark F. Adams 
1607f66b68fSMark Adams   if (!pc_gamg->data || (pc_gamg->data_sz != arrsz)) {
1619566063dSJacob Faibussowitsch     PetscCall(PetscFree(pc_gamg->data));
1629566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(arrsz + 1, &pc_gamg->data));
1632e68589bSMark F. Adams   }
1642e68589bSMark F. Adams   /* copy data in - column oriented */
1652e68589bSMark F. Adams   for (kk = 0; kk < nloc; kk++) {
16669344418SMark F. Adams     const PetscInt M    = nloc * pc_gamg->data_cell_rows; /* stride into data */
16769344418SMark F. Adams     PetscReal     *data = &pc_gamg->data[kk * ndatarows]; /* start of cell */
168c8b0795cSMark F. Adams     if (pc_gamg->data_cell_cols == 1) *data = 1.0;
1692e68589bSMark F. Adams     else {
17069344418SMark F. Adams       /* translational modes */
1712fa5cd67SKarl Rupp       for (ii = 0; ii < ndatarows; ii++) {
1722fa5cd67SKarl Rupp         for (jj = 0; jj < ndatarows; jj++) {
17369344418SMark F. Adams           if (ii == jj) data[ii * M + jj] = 1.0;
1742e68589bSMark F. Adams           else data[ii * M + jj] = 0.0;
1752fa5cd67SKarl Rupp         }
1762fa5cd67SKarl Rupp       }
17769344418SMark F. Adams 
17869344418SMark F. Adams       /* rotational modes */
1792e68589bSMark F. Adams       if (coords) {
18069344418SMark F. Adams         if (ndm == 2) {
1812e68589bSMark F. Adams           data += 2 * M;
1822e68589bSMark F. Adams           data[0] = -coords[2 * kk + 1];
1832e68589bSMark F. Adams           data[1] = coords[2 * kk];
184806fa848SBarry Smith         } else {
1852e68589bSMark F. Adams           data += 3 * M;
1869371c9d4SSatish Balay           data[0]         = 0.0;
1879371c9d4SSatish Balay           data[M + 0]     = coords[3 * kk + 2];
1889371c9d4SSatish Balay           data[2 * M + 0] = -coords[3 * kk + 1];
1899371c9d4SSatish Balay           data[1]         = -coords[3 * kk + 2];
1909371c9d4SSatish Balay           data[M + 1]     = 0.0;
1919371c9d4SSatish Balay           data[2 * M + 1] = coords[3 * kk];
1929371c9d4SSatish Balay           data[2]         = coords[3 * kk + 1];
1939371c9d4SSatish Balay           data[M + 2]     = -coords[3 * kk];
1949371c9d4SSatish Balay           data[2 * M + 2] = 0.0;
1952e68589bSMark F. Adams         }
1962e68589bSMark F. Adams       }
1972e68589bSMark F. Adams     }
1982e68589bSMark F. Adams   }
1992e68589bSMark F. Adams   pc_gamg->data_sz = arrsz;
2002e68589bSMark F. Adams   PetscFunctionReturn(0);
2012e68589bSMark F. Adams }
2022e68589bSMark F. Adams 
2032e68589bSMark F. Adams /*
204a2f3521dSMark F. Adams    PCSetData_AGG - called if data is not set with PCSetCoordinates.
205a2f3521dSMark F. Adams       Looks in Mat for near null space.
206a2f3521dSMark F. Adams       Does not work for Stokes
2072e68589bSMark F. Adams 
2082e68589bSMark F. Adams   Input Parameter:
2092e68589bSMark F. Adams    . pc -
210a2f3521dSMark F. Adams    . a_A - matrix to get (near) null space out of.
2112e68589bSMark F. Adams */
212d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetData_AGG(PC pc, Mat a_A)
213d71ae5a4SJacob Faibussowitsch {
214b8cd405aSMark F. Adams   PC_MG       *mg      = (PC_MG *)pc->data;
215b8cd405aSMark F. Adams   PC_GAMG     *pc_gamg = (PC_GAMG *)mg->innerctx;
216b8cd405aSMark F. Adams   MatNullSpace mnull;
21766f2ef4bSMark Adams 
218b817416eSBarry Smith   PetscFunctionBegin;
2199566063dSJacob Faibussowitsch   PetscCall(MatGetNearNullSpace(a_A, &mnull));
220b8cd405aSMark F. Adams   if (!mnull) {
22166f2ef4bSMark Adams     DM dm;
2229566063dSJacob Faibussowitsch     PetscCall(PCGetDM(pc, &dm));
22348a46eb9SPierre Jolivet     if (!dm) PetscCall(MatGetDM(a_A, &dm));
22466f2ef4bSMark Adams     if (dm) {
22566f2ef4bSMark Adams       PetscObject deformation;
226b0d30dd6SMatthew G. Knepley       PetscInt    Nf;
227b0d30dd6SMatthew G. Knepley 
2289566063dSJacob Faibussowitsch       PetscCall(DMGetNumFields(dm, &Nf));
229b0d30dd6SMatthew G. Knepley       if (Nf) {
2309566063dSJacob Faibussowitsch         PetscCall(DMGetField(dm, 0, NULL, &deformation));
2319566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery((PetscObject)deformation, "nearnullspace", (PetscObject *)&mnull));
23248a46eb9SPierre Jolivet         if (!mnull) PetscCall(PetscObjectQuery((PetscObject)deformation, "nullspace", (PetscObject *)&mnull));
23366f2ef4bSMark Adams       }
23466f2ef4bSMark Adams     }
235b0d30dd6SMatthew G. Knepley   }
23666f2ef4bSMark Adams 
23766f2ef4bSMark Adams   if (!mnull) {
238a2f3521dSMark F. Adams     PetscInt bs, NN, MM;
2399566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSize(a_A, &bs));
2409566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(a_A, &MM, &NN));
24163a3b9bcSJacob Faibussowitsch     PetscCheck(MM % bs == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MM %" PetscInt_FMT " must be divisible by bs %" PetscInt_FMT, MM, bs);
2429566063dSJacob Faibussowitsch     PetscCall(PCSetCoordinates_AGG(pc, bs, MM / bs, NULL));
243806fa848SBarry Smith   } else {
244b8cd405aSMark F. Adams     PetscReal         *nullvec;
245b8cd405aSMark F. Adams     PetscBool          has_const;
246b8cd405aSMark F. Adams     PetscInt           i, j, mlocal, nvec, bs;
247d5d25401SBarry Smith     const Vec         *vecs;
248d5d25401SBarry Smith     const PetscScalar *v;
249b817416eSBarry Smith 
2509566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(a_A, &mlocal, NULL));
2519566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs));
252a0dea87dSMark Adams     pc_gamg->data_sz = (nvec + !!has_const) * mlocal;
2539566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1((nvec + !!has_const) * mlocal, &nullvec));
2549371c9d4SSatish Balay     if (has_const)
2559371c9d4SSatish Balay       for (i = 0; i < mlocal; i++) nullvec[i] = 1.0;
256b8cd405aSMark F. Adams     for (i = 0; i < nvec; i++) {
2579566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(vecs[i], &v));
258b8cd405aSMark F. Adams       for (j = 0; j < mlocal; j++) nullvec[(i + !!has_const) * mlocal + j] = PetscRealPart(v[j]);
2599566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(vecs[i], &v));
260b8cd405aSMark F. Adams     }
261b8cd405aSMark F. Adams     pc_gamg->data           = nullvec;
262b8cd405aSMark F. Adams     pc_gamg->data_cell_cols = (nvec + !!has_const);
2639566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSize(a_A, &bs));
264b8cd405aSMark F. Adams     pc_gamg->data_cell_rows = bs;
265b8cd405aSMark F. Adams   }
2662e68589bSMark F. Adams   PetscFunctionReturn(0);
2672e68589bSMark F. Adams }
2682e68589bSMark F. Adams 
2692e68589bSMark F. Adams /*
270bae903cbSmarkadams4   formProl0 - collect null space data for each aggregate, do QR, put R in coarse grid data and Q in P_0
2712e68589bSMark F. Adams 
2722e68589bSMark F. Adams   Input Parameter:
2732adfe9d3SBarry Smith    . agg_llists - list of arrays with aggregates -- list from selected vertices of aggregate unselected vertices
2742adfe9d3SBarry Smith    . bs - row block size
2750cbbd2e1SMark F. Adams    . nSAvec - column bs of new P
2760cbbd2e1SMark F. Adams    . my0crs - global index of start of locals
2772adfe9d3SBarry Smith    . data_stride - bs*(nloc nodes + ghost nodes) [data_stride][nSAvec]
2782e68589bSMark F. Adams    . data_in[data_stride*nSAvec] - local data on fine grid
2792e68589bSMark F. Adams    . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist'
280bae903cbSmarkadams4 
2812e68589bSMark F. Adams   Output Parameter:
2822e68589bSMark F. Adams    . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data
2832e68589bSMark F. Adams    . a_Prol - prolongation operator
2842e68589bSMark F. Adams */
285d71ae5a4SJacob Faibussowitsch static PetscErrorCode formProl0(PetscCoarsenData *agg_llists, PetscInt bs, PetscInt nSAvec, PetscInt my0crs, PetscInt data_stride, PetscReal data_in[], const PetscInt flid_fgid[], PetscReal **a_data_out, Mat a_Prol)
286d71ae5a4SJacob Faibussowitsch {
287bd026e97SJed Brown   PetscInt        Istart, my0, Iend, nloc, clid, flid = 0, aggID, kk, jj, ii, mm, nSelected, minsz, nghosts, out_data_stride;
2883b4367a7SBarry Smith   MPI_Comm        comm;
2892e68589bSMark F. Adams   PetscReal      *out_data;
290539c167fSBarry Smith   PetscCDIntNd   *pos;
2911943db53SBarry Smith   PCGAMGHashTable fgid_flid;
2920cbbd2e1SMark F. Adams 
2932e68589bSMark F. Adams   PetscFunctionBegin;
2949566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)a_Prol, &comm));
2959566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(a_Prol, &Istart, &Iend));
2969371c9d4SSatish Balay   nloc = (Iend - Istart) / bs;
2979371c9d4SSatish Balay   my0  = Istart / bs;
29863a3b9bcSJacob Faibussowitsch   PetscCheck((Iend - Istart) % bs == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Iend %" PetscInt_FMT " - Istart %" PetscInt_FMT " must be divisible by bs %" PetscInt_FMT, Iend, Istart, bs);
2990cbbd2e1SMark F. Adams   Iend /= bs;
3000cbbd2e1SMark F. Adams   nghosts = data_stride / bs - nloc;
3012e68589bSMark F. Adams 
3029566063dSJacob Faibussowitsch   PetscCall(PCGAMGHashTableCreate(2 * nghosts + 1, &fgid_flid));
30348a46eb9SPierre Jolivet   for (kk = 0; kk < nghosts; kk++) PetscCall(PCGAMGHashTableAdd(&fgid_flid, flid_fgid[nloc + kk], nloc + kk));
3042e68589bSMark F. Adams 
3050cbbd2e1SMark F. Adams   /* count selected -- same as number of cols of P */
3060cbbd2e1SMark F. Adams   for (nSelected = mm = 0; mm < nloc; mm++) {
307e78576d6SMark F. Adams     PetscBool ise;
3089566063dSJacob Faibussowitsch     PetscCall(PetscCDEmptyAt(agg_llists, mm, &ise));
309e78576d6SMark F. Adams     if (!ise) nSelected++;
3100cbbd2e1SMark F. Adams   }
3119566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRangeColumn(a_Prol, &ii, &jj));
31263a3b9bcSJacob Faibussowitsch   PetscCheck((ii / nSAvec) == my0crs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ii %" PetscInt_FMT " /nSAvec %" PetscInt_FMT "  != my0crs %" PetscInt_FMT, ii, nSAvec, my0crs);
31363a3b9bcSJacob Faibussowitsch   PetscCheck(nSelected == (jj - ii) / nSAvec, PETSC_COMM_SELF, PETSC_ERR_PLIB, "nSelected %" PetscInt_FMT " != (jj %" PetscInt_FMT " - ii %" PetscInt_FMT ")/nSAvec %" PetscInt_FMT, nSelected, jj, ii, nSAvec);
3140cbbd2e1SMark F. Adams 
3152e68589bSMark F. Adams   /* aloc space for coarse point data (output) */
3160cbbd2e1SMark F. Adams   out_data_stride = nSelected * nSAvec;
3172fa5cd67SKarl Rupp 
3189566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(out_data_stride * nSAvec, &out_data));
31933812677SSatish Balay   for (ii = 0; ii < out_data_stride * nSAvec; ii++) out_data[ii] = PETSC_MAX_REAL;
3200cbbd2e1SMark F. Adams   *a_data_out = out_data; /* output - stride nSelected*nSAvec */
3212e68589bSMark F. Adams 
3222e68589bSMark F. Adams   /* find points and set prolongation */
323c8b0795cSMark F. Adams   minsz = 100;
3240cbbd2e1SMark F. Adams   for (mm = clid = 0; mm < nloc; mm++) {
3259566063dSJacob Faibussowitsch     PetscCall(PetscCDSizeAt(agg_llists, mm, &jj));
326e78576d6SMark F. Adams     if (jj > 0) {
3270cbbd2e1SMark F. Adams       const PetscInt lid = mm, cgid = my0crs + clid;
3280cbbd2e1SMark F. Adams       PetscInt       cids[100]; /* max bs */
3290cbbd2e1SMark F. Adams       PetscBLASInt   asz = jj, M = asz * bs, N = nSAvec, INFO;
3302e68589bSMark F. Adams       PetscBLASInt   Mdata = M + ((N - M > 0) ? N - M : 0), LDA = Mdata, LWORK = N * bs;
3312e68589bSMark F. Adams       PetscScalar   *qqc, *qqr, *TAU, *WORK;
3322e68589bSMark F. Adams       PetscInt      *fids;
33365d7b583SSatish Balay       PetscReal     *data;
334b817416eSBarry Smith 
3350cbbd2e1SMark F. Adams       /* count agg */
3360cbbd2e1SMark F. Adams       if (asz < minsz) minsz = asz;
3370cbbd2e1SMark F. Adams 
3380cbbd2e1SMark F. Adams       /* get block */
3399566063dSJacob Faibussowitsch       PetscCall(PetscMalloc5(Mdata * N, &qqc, M * N, &qqr, N, &TAU, LWORK, &WORK, M, &fids));
3402e68589bSMark F. Adams 
3412e68589bSMark F. Adams       aggID = 0;
3429566063dSJacob Faibussowitsch       PetscCall(PetscCDGetHeadPos(agg_llists, lid, &pos));
343e78576d6SMark F. Adams       while (pos) {
344e78576d6SMark F. Adams         PetscInt gid1;
3459566063dSJacob Faibussowitsch         PetscCall(PetscCDIntNdGetID(pos, &gid1));
3469566063dSJacob Faibussowitsch         PetscCall(PetscCDGetNextPos(agg_llists, lid, &pos));
347e78576d6SMark F. Adams 
3480cbbd2e1SMark F. Adams         if (gid1 >= my0 && gid1 < Iend) flid = gid1 - my0;
3490cbbd2e1SMark F. Adams         else {
3509566063dSJacob Faibussowitsch           PetscCall(PCGAMGHashTableFind(&fgid_flid, gid1, &flid));
35108401ef6SPierre Jolivet           PetscCheck(flid >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cannot find gid1 in table");
3520cbbd2e1SMark F. Adams         }
3532e68589bSMark F. Adams         /* copy in B_i matrix - column oriented */
35465d7b583SSatish Balay         data = &data_in[flid * bs];
3553d3eaba7SBarry Smith         for (ii = 0; ii < bs; ii++) {
3562e68589bSMark F. Adams           for (jj = 0; jj < N; jj++) {
3570cbbd2e1SMark F. Adams             PetscReal d                       = data[jj * data_stride + ii];
3580cbbd2e1SMark F. Adams             qqc[jj * Mdata + aggID * bs + ii] = d;
3592e68589bSMark F. Adams           }
3602e68589bSMark F. Adams         }
3612e68589bSMark F. Adams         /* set fine IDs */
3622e68589bSMark F. Adams         for (kk = 0; kk < bs; kk++) fids[aggID * bs + kk] = flid_fgid[flid] * bs + kk;
3632e68589bSMark F. Adams         aggID++;
3640cbbd2e1SMark F. Adams       }
3652e68589bSMark F. Adams 
3662e68589bSMark F. Adams       /* pad with zeros */
3672e68589bSMark F. Adams       for (ii = asz * bs; ii < Mdata; ii++) {
368ad540459SPierre Jolivet         for (jj = 0; jj < N; jj++, kk++) qqc[jj * Mdata + ii] = .0;
3692e68589bSMark F. Adams       }
3702e68589bSMark F. Adams 
3712e68589bSMark F. Adams       /* QR */
3729566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
373792fecdfSBarry Smith       PetscCallBLAS("LAPACKgeqrf", LAPACKgeqrf_(&Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
3749566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPop());
37508401ef6SPierre Jolivet       PetscCheck(INFO == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xGEQRF error");
3762e68589bSMark F. Adams       /* get R - column oriented - output B_{i+1} */
3772e68589bSMark F. Adams       {
3782e68589bSMark F. Adams         PetscReal *data = &out_data[clid * nSAvec];
3792e68589bSMark F. Adams         for (jj = 0; jj < nSAvec; jj++) {
3802e68589bSMark F. Adams           for (ii = 0; ii < nSAvec; ii++) {
38108401ef6SPierre Jolivet             PetscCheck(data[jj * out_data_stride + ii] == PETSC_MAX_REAL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "data[jj*out_data_stride + ii] != %e", (double)PETSC_MAX_REAL);
3820cbbd2e1SMark F. Adams             if (ii <= jj) data[jj * out_data_stride + ii] = PetscRealPart(qqc[jj * Mdata + ii]);
3830cbbd2e1SMark F. Adams             else data[jj * out_data_stride + ii] = 0.;
3842e68589bSMark F. Adams           }
3852e68589bSMark F. Adams         }
3862e68589bSMark F. Adams       }
3872e68589bSMark F. Adams 
3882e68589bSMark F. Adams       /* get Q - row oriented */
389792fecdfSBarry Smith       PetscCallBLAS("LAPACKorgqr", LAPACKorgqr_(&Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
39063a3b9bcSJacob Faibussowitsch       PetscCheck(INFO == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xORGQR error arg %" PetscBLASInt_FMT, -INFO);
3912e68589bSMark F. Adams 
3922e68589bSMark F. Adams       for (ii = 0; ii < M; ii++) {
393ad540459SPierre Jolivet         for (jj = 0; jj < N; jj++) qqr[N * ii + jj] = qqc[jj * Mdata + ii];
3942e68589bSMark F. Adams       }
3952e68589bSMark F. Adams 
3962e68589bSMark F. Adams       /* add diagonal block of P0 */
3979371c9d4SSatish Balay       for (kk = 0; kk < N; kk++) { cids[kk] = N * cgid + kk; /* global col IDs in P0 */ }
3989566063dSJacob Faibussowitsch       PetscCall(MatSetValues(a_Prol, M, fids, N, cids, qqr, INSERT_VALUES));
3999566063dSJacob Faibussowitsch       PetscCall(PetscFree5(qqc, qqr, TAU, WORK, fids));
400b43b03e9SMark F. Adams       clid++;
4010cbbd2e1SMark F. Adams     } /* coarse agg */
4020cbbd2e1SMark F. Adams   }   /* for all fine nodes */
4039566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(a_Prol, MAT_FINAL_ASSEMBLY));
4049566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(a_Prol, MAT_FINAL_ASSEMBLY));
4059566063dSJacob Faibussowitsch   PetscCall(PCGAMGHashTableDestroy(&fgid_flid));
4062e68589bSMark F. Adams   PetscFunctionReturn(0);
4072e68589bSMark F. Adams }
4082e68589bSMark F. Adams 
409d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_GAMG_AGG(PC pc, PetscViewer viewer)
410d71ae5a4SJacob Faibussowitsch {
4115adeb434SBarry Smith   PC_MG       *mg          = (PC_MG *)pc->data;
4125adeb434SBarry Smith   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
4135adeb434SBarry Smith   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
4145adeb434SBarry Smith 
4155adeb434SBarry Smith   PetscFunctionBegin;
4169566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "      AGG specific options\n"));
417bae903cbSmarkadams4   PetscCall(PetscViewerASCIIPrintf(viewer, "        Number of levels to square graph %d\n", (int)pc_gamg_agg->aggressive_coarsening_levels));
418bae903cbSmarkadams4   PetscCall(PetscViewerASCIIPrintf(viewer, "        Number smoothing steps %d\n", (int)pc_gamg_agg->nsmooths));
4195adeb434SBarry Smith   PetscFunctionReturn(0);
4205adeb434SBarry Smith }
4215adeb434SBarry Smith 
4222d776b49SBarry Smith static PetscErrorCode PCGAMGCreateGraph_AGG(PC pc, Mat Amat, Mat *a_Gmat)
423d71ae5a4SJacob Faibussowitsch {
424c8b0795cSMark F. Adams   PC_MG          *mg          = (PC_MG *)pc->data;
425c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG *)mg->innerctx;
426c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
4272d776b49SBarry Smith   const PetscReal vfilter     = pc_gamg->threshold[pc_gamg->current_level];
4282d776b49SBarry Smith   PetscBool       ishem;
4292d776b49SBarry Smith   const char     *prefix;
430c8b0795cSMark F. Adams 
431c8b0795cSMark F. Adams   PetscFunctionBegin;
432849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_GRAPH], 0, 0, 0, 0));
4332d776b49SBarry Smith   /* Note: depending on the algorithm that will be used for computing the coarse grid points this should pass PETSC_TRUE or PETSC_FALSE as the first argument */
4342d776b49SBarry Smith 
4352d776b49SBarry Smith   /* MATCOARSENHEM requires numerical weights for edges so ensure they are computed */
4362d776b49SBarry Smith   PetscCall(MatCoarsenCreate(PetscObjectComm((PetscObject)pc), &pc_gamg_agg->crs));
4372d776b49SBarry Smith   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc, &prefix));
4382d776b49SBarry Smith   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc_gamg_agg->crs, prefix));
4392d776b49SBarry Smith   PetscCall(MatCoarsenSetFromOptions(pc_gamg_agg->crs));
4402d776b49SBarry Smith   PetscCall(PetscObjectTypeCompare((PetscObject)pc_gamg_agg->crs, MATCOARSENHEM, &ishem));
4412d776b49SBarry Smith 
4422d776b49SBarry Smith   PetscCall(MatCreateGraph(Amat, PETSC_TRUE, (vfilter >= 0 || ishem) ? PETSC_TRUE : PETSC_FALSE, vfilter, a_Gmat));
443849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_GRAPH], 0, 0, 0, 0));
444c8b0795cSMark F. Adams   PetscFunctionReturn(0);
445c8b0795cSMark F. Adams }
446c8b0795cSMark F. Adams 
447c8b0795cSMark F. Adams /*
448bae903cbSmarkadams4    PCGAMGCoarsen_AGG - supports squaring the graph (deprecated) and new graph for
449bae903cbSmarkadams4      communication of QR data used with HEM and MISk coarsening
450c8b0795cSMark F. Adams 
451c8b0795cSMark F. Adams   Input Parameter:
452e0940f08SMark F. Adams    . a_pc - this
453bae903cbSmarkadams4 
454e0940f08SMark F. Adams   Input/Output Parameter:
455bae903cbSmarkadams4    . a_Gmat1 - graph to coarsen (in), graph off processor edges for QR gather scatter (out)
456bae903cbSmarkadams4 
457c8b0795cSMark F. Adams   Output Parameter:
4580cbbd2e1SMark F. Adams    . agg_lists - list of aggregates
459bae903cbSmarkadams4 
460c8b0795cSMark F. Adams */
461d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGCoarsen_AGG(PC a_pc, Mat *a_Gmat1, PetscCoarsenData **agg_lists)
462d71ae5a4SJacob Faibussowitsch {
463e0940f08SMark F. Adams   PC_MG       *mg          = (PC_MG *)a_pc->data;
464c8b0795cSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
465c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
466bae903cbSmarkadams4   Mat          mat, Gmat1 = *a_Gmat1; /* aggressive graph */
4670cbbd2e1SMark F. Adams   IS           perm;
468bae903cbSmarkadams4   PetscInt     Istart, Iend, Ii, nloc, bs, nn;
469bae903cbSmarkadams4   PetscInt    *permute, *degree;
470c8b0795cSMark F. Adams   PetscBool   *bIndexSet;
471aea53286SMark Adams   PetscReal    hashfact;
472e2c00dcbSBarry Smith   PetscInt     iSwapIndex;
4733b50add6SMark Adams   PetscRandom  random;
474c8b0795cSMark F. Adams 
475c8b0795cSMark F. Adams   PetscFunctionBegin;
476849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0));
477bae903cbSmarkadams4   PetscCall(MatGetLocalSize(Gmat1, &nn, NULL));
4789566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(Gmat1, &bs));
47963a3b9bcSJacob Faibussowitsch   PetscCheck(bs == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "bs %" PetscInt_FMT " must be 1", bs);
480bae903cbSmarkadams4   nloc = nn / bs;
481c8b0795cSMark F. Adams 
4825cfd4bd9SMark Adams   /* get MIS aggs - randomize */
483bae903cbSmarkadams4   PetscCall(PetscMalloc2(nloc, &permute, nloc, &degree));
4849566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(nloc, &bIndexSet));
4856e09b0e3SStefano Zampini   for (Ii = 0; Ii < nloc; Ii++) permute[Ii] = Ii;
4869566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &random));
4879566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Gmat1, &Istart, &Iend));
488c8b0795cSMark F. Adams   for (Ii = 0; Ii < nloc; Ii++) {
489bae903cbSmarkadams4     PetscInt nc;
490bae903cbSmarkadams4     PetscCall(MatGetRow(Gmat1, Istart + Ii, &nc, NULL, NULL));
491bae903cbSmarkadams4     degree[Ii] = nc;
492bae903cbSmarkadams4     PetscCall(MatRestoreRow(Gmat1, Istart + Ii, &nc, NULL, NULL));
493bae903cbSmarkadams4   }
494bae903cbSmarkadams4   for (Ii = 0; Ii < nloc; Ii++) {
4959566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValueReal(random, &hashfact));
496aea53286SMark Adams     iSwapIndex = (PetscInt)(hashfact * nloc) % nloc;
497c8b0795cSMark F. Adams     if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) {
498c8b0795cSMark F. Adams       PetscInt iTemp        = permute[iSwapIndex];
499c8b0795cSMark F. Adams       permute[iSwapIndex]   = permute[Ii];
500c8b0795cSMark F. Adams       permute[Ii]           = iTemp;
501bae903cbSmarkadams4       iTemp                 = degree[iSwapIndex];
502bae903cbSmarkadams4       degree[iSwapIndex]    = degree[Ii];
503bae903cbSmarkadams4       degree[Ii]            = iTemp;
504c8b0795cSMark F. Adams       bIndexSet[iSwapIndex] = PETSC_TRUE;
505c8b0795cSMark F. Adams     }
506c8b0795cSMark F. Adams   }
507bae903cbSmarkadams4   // create minimum degree ordering
508bae903cbSmarkadams4   PetscCall(PetscSortIntWithArray(nloc, degree, permute));
509bae903cbSmarkadams4 
5109566063dSJacob Faibussowitsch   PetscCall(PetscFree(bIndexSet));
5119566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&random));
5129566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm));
513849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_MIS], 0, 0, 0, 0));
5142d776b49SBarry Smith   PetscCall(MatCoarsenSetGreedyOrdering(pc_gamg_agg->crs, perm));
5152d776b49SBarry Smith   PetscCall(MatCoarsenSetAdjacency(pc_gamg_agg->crs, Gmat1));
5162d776b49SBarry Smith   PetscCall(MatCoarsenSetStrictAggs(pc_gamg_agg->crs, PETSC_TRUE));
5172d776b49SBarry Smith   if (pc_gamg->current_level < pc_gamg_agg->aggressive_coarsening_levels) PetscCall(MatCoarsenMISKSetDistance(pc_gamg_agg->crs, 2)); // hardwire to MIS-2
5182d776b49SBarry Smith   else PetscCall(MatCoarsenMISKSetDistance(pc_gamg_agg->crs, 1));                                                                    // MIS
5192d776b49SBarry Smith   PetscCall(MatCoarsenApply(pc_gamg_agg->crs));
5202d776b49SBarry Smith   PetscCall(MatCoarsenViewFromOptions(pc_gamg_agg->crs, NULL, "-mat_coarsen_view"));
5212d776b49SBarry Smith   PetscCall(MatCoarsenGetData(pc_gamg_agg->crs, agg_lists)); /* output */
5222d776b49SBarry Smith   PetscCall(MatCoarsenDestroy(&pc_gamg_agg->crs));
523b43b03e9SMark F. Adams 
5249566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
525bae903cbSmarkadams4   PetscCall(PetscFree2(permute, degree));
526849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_MIS], 0, 0, 0, 0));
5279f3f12c8SMark F. Adams 
528bae903cbSmarkadams4   {
529bae903cbSmarkadams4     PetscCoarsenData *llist = *agg_lists;
5309ab59c8bSMark Adams     /* see if we have a matrix that takes precedence (returned from MatCoarsenApply) */
5319566063dSJacob Faibussowitsch     PetscCall(PetscCDGetMat(llist, &mat));
5320cbbd2e1SMark F. Adams     if (mat) {
5339566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&Gmat1));
5340cbbd2e1SMark F. Adams       *a_Gmat1 = mat; /* output */
5350cbbd2e1SMark F. Adams     }
5360cbbd2e1SMark F. Adams   }
537849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0));
538c8b0795cSMark F. Adams   PetscFunctionReturn(0);
539c8b0795cSMark F. Adams }
540c8b0795cSMark F. Adams 
541c8b0795cSMark F. Adams /*
5420cbbd2e1SMark F. Adams  PCGAMGProlongator_AGG
543c8b0795cSMark F. Adams 
544c8b0795cSMark F. Adams  Input Parameter:
545c8b0795cSMark F. Adams  . pc - this
546c8b0795cSMark F. Adams  . Amat - matrix on this fine level
547c8b0795cSMark F. Adams  . Graph - used to get ghost data for nodes in
5480cbbd2e1SMark F. Adams  . agg_lists - list of aggregates
549c8b0795cSMark F. Adams  Output Parameter:
550c8b0795cSMark F. Adams  . a_P_out - prolongation operator to the next level
551c8b0795cSMark F. Adams  */
552d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGProlongator_AGG(PC pc, Mat Amat, Mat Gmat, PetscCoarsenData *agg_lists, Mat *a_P_out)
553d71ae5a4SJacob Faibussowitsch {
5542e68589bSMark F. Adams   PC_MG         *mg      = (PC_MG *)pc->data;
5552e68589bSMark F. Adams   PC_GAMG       *pc_gamg = (PC_GAMG *)mg->innerctx;
5564a2f8832SBarry Smith   const PetscInt col_bs  = pc_gamg->data_cell_cols;
557c8b0795cSMark F. Adams   PetscInt       Istart, Iend, nloc, ii, jj, kk, my0, nLocalSelected, bs;
558c8b0795cSMark F. Adams   Mat            Prol;
559d5d25401SBarry Smith   PetscMPIInt    size;
5603b4367a7SBarry Smith   MPI_Comm       comm;
561c8b0795cSMark F. Adams   PetscReal     *data_w_ghost;
562c8b0795cSMark F. Adams   PetscInt       myCrs0, nbnodes = 0, *flid_fgid;
563d9558ea9SBarry Smith   MatType        mtype;
5642e68589bSMark F. Adams 
5652e68589bSMark F. Adams   PetscFunctionBegin;
5669566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm));
56708401ef6SPierre Jolivet   PetscCheck(col_bs >= 1, comm, PETSC_ERR_PLIB, "Column bs cannot be less than 1");
568849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0));
5699566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
5709566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Amat, &Istart, &Iend));
5719566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(Amat, &bs));
5729371c9d4SSatish Balay   nloc = (Iend - Istart) / bs;
5739371c9d4SSatish Balay   my0  = Istart / bs;
57463a3b9bcSJacob Faibussowitsch   PetscCheck((Iend - Istart) % bs == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "(Iend %" PetscInt_FMT " - Istart %" PetscInt_FMT ") not divisible by bs %" PetscInt_FMT, Iend, Istart, bs);
5752e68589bSMark F. Adams 
5762e68589bSMark F. Adams   /* get 'nLocalSelected' */
5770cbbd2e1SMark F. Adams   for (ii = 0, nLocalSelected = 0; ii < nloc; ii++) {
578e78576d6SMark F. Adams     PetscBool ise;
5790cbbd2e1SMark F. Adams     /* filter out singletons 0 or 1? */
5809566063dSJacob Faibussowitsch     PetscCall(PetscCDEmptyAt(agg_lists, ii, &ise));
581e78576d6SMark F. Adams     if (!ise) nLocalSelected++;
5822e68589bSMark F. Adams   }
5832e68589bSMark F. Adams 
5842e68589bSMark F. Adams   /* create prolongator, create P matrix */
5859566063dSJacob Faibussowitsch   PetscCall(MatGetType(Amat, &mtype));
5869566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &Prol));
5879566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Prol, nloc * bs, nLocalSelected * col_bs, PETSC_DETERMINE, PETSC_DETERMINE));
5889566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSizes(Prol, bs, col_bs));
5899566063dSJacob Faibussowitsch   PetscCall(MatSetType(Prol, mtype));
5909566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(Prol, col_bs, NULL));
5919566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(Prol, col_bs, NULL, col_bs, NULL));
5922e68589bSMark F. Adams 
5932e68589bSMark F. Adams   /* can get all points "removed" */
5949566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Prol, &kk, &ii));
5957f66b68fSMark Adams   if (!ii) {
59663a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(pc, "%s: No selected points on coarse grid\n", ((PetscObject)pc)->prefix));
5979566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Prol));
5980298fd71SBarry Smith     *a_P_out = NULL; /* out */
599849bee69SMark Adams     PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0));
6002e68589bSMark F. Adams     PetscFunctionReturn(0);
6012e68589bSMark F. Adams   }
60263a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(pc, "%s: New grid %" PetscInt_FMT " nodes\n", ((PetscObject)pc)->prefix, ii / col_bs));
6039566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRangeColumn(Prol, &myCrs0, &kk));
6040cbbd2e1SMark F. Adams 
60563a3b9bcSJacob Faibussowitsch   PetscCheck((kk - myCrs0) % col_bs == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "(kk %" PetscInt_FMT " -myCrs0 %" PetscInt_FMT ") not divisible by col_bs %" PetscInt_FMT, kk, myCrs0, col_bs);
606c8b0795cSMark F. Adams   myCrs0 = myCrs0 / col_bs;
60763a3b9bcSJacob Faibussowitsch   PetscCheck((kk / col_bs - myCrs0) == nLocalSelected, PETSC_COMM_SELF, PETSC_ERR_PLIB, "(kk %" PetscInt_FMT "/col_bs %" PetscInt_FMT " - myCrs0 %" PetscInt_FMT ") != nLocalSelected %" PetscInt_FMT ")", kk, col_bs, myCrs0, nLocalSelected);
6082e68589bSMark F. Adams 
6092e68589bSMark F. Adams   /* create global vector of data in 'data_w_ghost' */
610849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROLA], 0, 0, 0, 0));
611bae903cbSmarkadams4   if (size > 1) { /* get ghost null space data */
6122e68589bSMark F. Adams     PetscReal *tmp_gdata, *tmp_ldata, *tp2;
6139566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nloc, &tmp_ldata));
6144a2f8832SBarry Smith     for (jj = 0; jj < col_bs; jj++) {
615c8b0795cSMark F. Adams       for (kk = 0; kk < bs; kk++) {
616a2f3521dSMark F. Adams         PetscInt         ii, stride;
617c8b0795cSMark F. Adams         const PetscReal *tp = pc_gamg->data + jj * bs * nloc + kk;
6182fa5cd67SKarl Rupp         for (ii = 0; ii < nloc; ii++, tp += bs) tmp_ldata[ii] = *tp;
6192fa5cd67SKarl Rupp 
6209566063dSJacob Faibussowitsch         PetscCall(PCGAMGGetDataWithGhosts(Gmat, 1, tmp_ldata, &stride, &tmp_gdata));
621a2f3521dSMark F. Adams 
622bae903cbSmarkadams4         if (!jj && !kk) { /* now I know how many total nodes - allocate TODO: move below and do in one 'col_bs' call */
6239566063dSJacob Faibussowitsch           PetscCall(PetscMalloc1(stride * bs * col_bs, &data_w_ghost));
624a2f3521dSMark F. Adams           nbnodes = bs * stride;
6252e68589bSMark F. Adams         }
626a2f3521dSMark F. Adams         tp2 = data_w_ghost + jj * bs * stride + kk;
6272fa5cd67SKarl Rupp         for (ii = 0; ii < stride; ii++, tp2 += bs) *tp2 = tmp_gdata[ii];
6289566063dSJacob Faibussowitsch         PetscCall(PetscFree(tmp_gdata));
6292e68589bSMark F. Adams       }
6302e68589bSMark F. Adams     }
6319566063dSJacob Faibussowitsch     PetscCall(PetscFree(tmp_ldata));
632806fa848SBarry Smith   } else {
633c8b0795cSMark F. Adams     nbnodes      = bs * nloc;
634c8b0795cSMark F. Adams     data_w_ghost = (PetscReal *)pc_gamg->data;
6352e68589bSMark F. Adams   }
6362e68589bSMark F. Adams 
637bae903cbSmarkadams4   /* get 'flid_fgid' TODO - move up to get 'stride' and do get null space data above in one step (jj loop) */
638c5df96a5SBarry Smith   if (size > 1) {
6392e68589bSMark F. Adams     PetscReal *fid_glid_loc, *fiddata;
640a2f3521dSMark F. Adams     PetscInt   stride;
6412e68589bSMark F. Adams 
6429566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nloc, &fid_glid_loc));
6432e68589bSMark F. Adams     for (kk = 0; kk < nloc; kk++) fid_glid_loc[kk] = (PetscReal)(my0 + kk);
6449566063dSJacob Faibussowitsch     PetscCall(PCGAMGGetDataWithGhosts(Gmat, 1, fid_glid_loc, &stride, &fiddata));
645bae903cbSmarkadams4     PetscCall(PetscMalloc1(stride, &flid_fgid)); /* copy real data to in */
646a2f3521dSMark F. Adams     for (kk = 0; kk < stride; kk++) flid_fgid[kk] = (PetscInt)fiddata[kk];
6479566063dSJacob Faibussowitsch     PetscCall(PetscFree(fiddata));
648a2f3521dSMark F. Adams 
64963a3b9bcSJacob Faibussowitsch     PetscCheck(stride == nbnodes / bs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "stride %" PetscInt_FMT " != nbnodes %" PetscInt_FMT "/bs %" PetscInt_FMT, stride, nbnodes, bs);
6509566063dSJacob Faibussowitsch     PetscCall(PetscFree(fid_glid_loc));
651806fa848SBarry Smith   } else {
6529566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nloc, &flid_fgid));
6532e68589bSMark F. Adams     for (kk = 0; kk < nloc; kk++) flid_fgid[kk] = my0 + kk;
6542e68589bSMark F. Adams   }
655849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROLA], 0, 0, 0, 0));
656bae903cbSmarkadams4   /* get P0 */
657849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROLB], 0, 0, 0, 0));
658c8b0795cSMark F. Adams   {
6590298fd71SBarry Smith     PetscReal *data_out = NULL;
6609566063dSJacob Faibussowitsch     PetscCall(formProl0(agg_lists, bs, col_bs, myCrs0, nbnodes, data_w_ghost, flid_fgid, &data_out, Prol));
6619566063dSJacob Faibussowitsch     PetscCall(PetscFree(pc_gamg->data));
6622fa5cd67SKarl Rupp 
663c8b0795cSMark F. Adams     pc_gamg->data           = data_out;
6644a2f8832SBarry Smith     pc_gamg->data_cell_rows = col_bs;
6654a2f8832SBarry Smith     pc_gamg->data_sz        = col_bs * col_bs * nLocalSelected;
666c8b0795cSMark F. Adams   }
667849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROLB], 0, 0, 0, 0));
66848a46eb9SPierre Jolivet   if (size > 1) PetscCall(PetscFree(data_w_ghost));
6699566063dSJacob Faibussowitsch   PetscCall(PetscFree(flid_fgid));
6702e68589bSMark F. Adams 
671c8b0795cSMark F. Adams   *a_P_out = Prol; /* out */
672ed4e82eaSPeter Brune 
673849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0));
674c8b0795cSMark F. Adams   PetscFunctionReturn(0);
675c8b0795cSMark F. Adams }
676c8b0795cSMark F. Adams 
677c8b0795cSMark F. Adams /*
678fd1112cbSBarry Smith    PCGAMGOptProlongator_AGG
679c8b0795cSMark F. Adams 
680c8b0795cSMark F. Adams   Input Parameter:
681c8b0795cSMark F. Adams    . pc - this
682c8b0795cSMark F. Adams    . Amat - matrix on this fine level
683c8b0795cSMark F. Adams  In/Output Parameter:
6842a808120SBarry Smith    . a_P - prolongation operator to the next level
685c8b0795cSMark F. Adams */
686d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGOptProlongator_AGG(PC pc, Mat Amat, Mat *a_P)
687d71ae5a4SJacob Faibussowitsch {
688c8b0795cSMark F. Adams   PC_MG       *mg          = (PC_MG *)pc->data;
689c8b0795cSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
690c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
691c8b0795cSMark F. Adams   PetscInt     jj;
692c8b0795cSMark F. Adams   Mat          Prol = *a_P;
6933b4367a7SBarry Smith   MPI_Comm     comm;
6942a808120SBarry Smith   KSP          eksp;
6952a808120SBarry Smith   Vec          bb, xx;
6962a808120SBarry Smith   PC           epc;
6972a808120SBarry Smith   PetscReal    alpha, emax, emin;
698c8b0795cSMark F. Adams 
699c8b0795cSMark F. Adams   PetscFunctionBegin;
7009566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm));
701849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_OPT], 0, 0, 0, 0));
702c8b0795cSMark F. Adams 
703d5d25401SBarry Smith   /* compute maximum singular value of operator to be used in smoother */
7042a808120SBarry Smith   if (0 < pc_gamg_agg->nsmooths) {
70518c3aa7eSMark     /* get eigen estimates */
70618c3aa7eSMark     if (pc_gamg->emax > 0) {
70718c3aa7eSMark       emin = pc_gamg->emin;
70818c3aa7eSMark       emax = pc_gamg->emax;
70918c3aa7eSMark     } else {
7100ed2132dSStefano Zampini       const char *prefix;
7110ed2132dSStefano Zampini 
7129566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(Amat, &bb, NULL));
7139566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(Amat, &xx, NULL));
7149566063dSJacob Faibussowitsch       PetscCall(KSPSetNoisy_Private(bb));
715e696ed0bSMark F. Adams 
7169566063dSJacob Faibussowitsch       PetscCall(KSPCreate(comm, &eksp));
7179566063dSJacob Faibussowitsch       PetscCall(PCGetOptionsPrefix(pc, &prefix));
7189566063dSJacob Faibussowitsch       PetscCall(KSPSetOptionsPrefix(eksp, prefix));
7199566063dSJacob Faibussowitsch       PetscCall(KSPAppendOptionsPrefix(eksp, "pc_gamg_esteig_"));
72073f7197eSJed Brown       {
721b94d7dedSBarry Smith         PetscBool isset, sflg;
722b94d7dedSBarry Smith         PetscCall(MatIsSPDKnown(Amat, &isset, &sflg));
723b94d7dedSBarry Smith         if (isset && sflg) PetscCall(KSPSetType(eksp, KSPCG));
724d24ecf33SMark       }
7259566063dSJacob Faibussowitsch       PetscCall(KSPSetErrorIfNotConverged(eksp, pc->erroriffailure));
7269566063dSJacob Faibussowitsch       PetscCall(KSPSetNormType(eksp, KSP_NORM_NONE));
727c2436cacSMark F. Adams 
7289566063dSJacob Faibussowitsch       PetscCall(KSPSetInitialGuessNonzero(eksp, PETSC_FALSE));
7299566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(eksp, Amat, Amat));
7302e68589bSMark F. Adams 
7319566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(eksp, &epc));
7329566063dSJacob Faibussowitsch       PetscCall(PCSetType(epc, PCJACOBI)); /* smoother in smoothed agg. */
733c2436cacSMark F. Adams 
7349566063dSJacob Faibussowitsch       PetscCall(KSPSetTolerances(eksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 10)); // 10 is safer, but 5 is often fine, can override with -pc_gamg_esteig_ksp_max_it -mg_levels_ksp_chebyshev_esteig 0,0.25,0,1.2
735074aec55SMark Adams 
7369566063dSJacob Faibussowitsch       PetscCall(KSPSetFromOptions(eksp));
7379566063dSJacob Faibussowitsch       PetscCall(KSPSetComputeSingularValues(eksp, PETSC_TRUE));
7389566063dSJacob Faibussowitsch       PetscCall(KSPSolve(eksp, bb, xx));
7399566063dSJacob Faibussowitsch       PetscCall(KSPCheckSolve(eksp, pc, xx));
7402e68589bSMark F. Adams 
7419566063dSJacob Faibussowitsch       PetscCall(KSPComputeExtremeSingularValues(eksp, &emax, &emin));
7429566063dSJacob Faibussowitsch       PetscCall(PetscInfo(pc, "%s: Smooth P0: max eigen=%e min=%e PC=%s\n", ((PetscObject)pc)->prefix, (double)emax, (double)emin, PCJACOBI));
7439566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&xx));
7449566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&bb));
7459566063dSJacob Faibussowitsch       PetscCall(KSPDestroy(&eksp));
7462e68589bSMark F. Adams     }
74718c3aa7eSMark     if (pc_gamg->use_sa_esteig) {
74818c3aa7eSMark       mg->min_eigen_DinvA[pc_gamg->current_level] = emin;
74918c3aa7eSMark       mg->max_eigen_DinvA[pc_gamg->current_level] = emax;
75063a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(pc, "%s: Smooth P0: level %" PetscInt_FMT ", cache spectra %g %g\n", ((PetscObject)pc)->prefix, pc_gamg->current_level, (double)emin, (double)emax));
75118c3aa7eSMark     } else {
75218c3aa7eSMark       mg->min_eigen_DinvA[pc_gamg->current_level] = 0;
75318c3aa7eSMark       mg->max_eigen_DinvA[pc_gamg->current_level] = 0;
75418c3aa7eSMark     }
75518c3aa7eSMark   } else {
75618c3aa7eSMark     mg->min_eigen_DinvA[pc_gamg->current_level] = 0;
75718c3aa7eSMark     mg->max_eigen_DinvA[pc_gamg->current_level] = 0;
75818c3aa7eSMark   }
7592e68589bSMark F. Adams 
7602a808120SBarry Smith   /* smooth P0 */
7612a808120SBarry Smith   for (jj = 0; jj < pc_gamg_agg->nsmooths; jj++) {
7622a808120SBarry Smith     Mat tMat;
7632a808120SBarry Smith     Vec diag;
7642a808120SBarry Smith 
765849bee69SMark Adams     PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_OPTSM], 0, 0, 0, 0));
7662a808120SBarry Smith 
7672e68589bSMark F. Adams     /* smooth P1 := (I - omega/lam D^{-1}A)P0 */
7689566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2], 0, 0, 0, 0));
7699566063dSJacob Faibussowitsch     PetscCall(MatMatMult(Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat));
7709566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2], 0, 0, 0, 0));
7719566063dSJacob Faibussowitsch     PetscCall(MatProductClear(tMat));
7729566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(Amat, &diag, NULL));
7739566063dSJacob Faibussowitsch     PetscCall(MatGetDiagonal(Amat, diag)); /* effectively PCJACOBI */
7749566063dSJacob Faibussowitsch     PetscCall(VecReciprocal(diag));
7759566063dSJacob Faibussowitsch     PetscCall(MatDiagonalScale(tMat, diag, NULL));
7769566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&diag));
777b4da3a1bSStefano Zampini 
778d5d25401SBarry Smith     /* TODO: Set a PCFailedReason and exit the building of the AMG preconditioner */
77908401ef6SPierre Jolivet     PetscCheck(emax != 0.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_PLIB, "Computed maximum singular value as zero");
780d5d25401SBarry Smith     /* TODO: Document the 1.4 and don't hardwire it in this routine */
781b4ec6429SMark F. Adams     alpha = -1.4 / emax;
782b4da3a1bSStefano Zampini 
7839566063dSJacob Faibussowitsch     PetscCall(MatAYPX(tMat, alpha, Prol, SUBSET_NONZERO_PATTERN));
7849566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Prol));
7852e68589bSMark F. Adams     Prol = tMat;
786849bee69SMark Adams     PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_OPTSM], 0, 0, 0, 0));
7872e68589bSMark F. Adams   }
788849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_OPT], 0, 0, 0, 0));
789c8b0795cSMark F. Adams   *a_P = Prol;
790c8b0795cSMark F. Adams   PetscFunctionReturn(0);
7912e68589bSMark F. Adams }
7922e68589bSMark F. Adams 
793c8b0795cSMark F. Adams /*
794c8b0795cSMark F. Adams    PCCreateGAMG_AGG
7952e68589bSMark F. Adams 
796c8b0795cSMark F. Adams   Input Parameter:
797c8b0795cSMark F. Adams    . pc -
798c8b0795cSMark F. Adams */
799d71ae5a4SJacob Faibussowitsch PetscErrorCode PCCreateGAMG_AGG(PC pc)
800d71ae5a4SJacob Faibussowitsch {
801c8b0795cSMark F. Adams   PC_MG       *mg      = (PC_MG *)pc->data;
802c8b0795cSMark F. Adams   PC_GAMG     *pc_gamg = (PC_GAMG *)mg->innerctx;
803c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg;
8042e68589bSMark F. Adams 
805c8b0795cSMark F. Adams   PetscFunctionBegin;
806c8b0795cSMark F. Adams   /* create sub context for SA */
8074dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&pc_gamg_agg));
808c8b0795cSMark F. Adams   pc_gamg->subctx = pc_gamg_agg;
809c8b0795cSMark F. Adams 
8101ab5ffc9SJed Brown   pc_gamg->ops->setfromoptions = PCSetFromOptions_GAMG_AGG;
8119b8ffb57SJed Brown   pc_gamg->ops->destroy        = PCDestroy_GAMG_AGG;
812c8b0795cSMark F. Adams   /* reset does not do anything; setup not virtual */
813c8b0795cSMark F. Adams 
814c8b0795cSMark F. Adams   /* set internal function pointers */
8152d776b49SBarry Smith   pc_gamg->ops->creategraph       = PCGAMGCreateGraph_AGG;
8161ab5ffc9SJed Brown   pc_gamg->ops->coarsen           = PCGAMGCoarsen_AGG;
8171ab5ffc9SJed Brown   pc_gamg->ops->prolongator       = PCGAMGProlongator_AGG;
818fd1112cbSBarry Smith   pc_gamg->ops->optprolongator    = PCGAMGOptProlongator_AGG;
8191ab5ffc9SJed Brown   pc_gamg->ops->createdefaultdata = PCSetData_AGG;
8205adeb434SBarry Smith   pc_gamg->ops->view              = PCView_GAMG_AGG;
821c8b0795cSMark F. Adams 
822bae903cbSmarkadams4   pc_gamg_agg->aggressive_coarsening_levels = 0;
823cab9ed1eSBarry Smith   pc_gamg_agg->nsmooths                     = 1;
824cab9ed1eSBarry Smith 
8259566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNSmooths_C", PCGAMGSetNSmooths_AGG));
826bae903cbSmarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveLevels_C", PCGAMGSetAggressiveLevels_AGG));
8279566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", PCSetCoordinates_AGG));
8282e68589bSMark F. Adams   PetscFunctionReturn(0);
8292e68589bSMark F. Adams }
830