xref: /petsc/src/ksp/pc/impls/gamg/agg.c (revision f1580f4e3ce5d5b2393648fd039d0d41b440385d)
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   PetscBool symmetrize_graph;
13bae903cbSmarkadams4   PetscInt  aggressive_coarsening_levels; // number of aggressive coarsening levels (square or MISk)
142e68589bSMark F. Adams } PC_GAMG_AGG;
152e68589bSMark F. Adams 
162e68589bSMark F. Adams /*@
17*f1580f4eSBarry Smith    PCGAMGSetNSmooths - Set number of smoothing steps (1 is typical) used for multigrid on all the levels
182e68589bSMark F. Adams 
19*f1580f4eSBarry Smith    Logically Collective on pc
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 
29*f1580f4eSBarry Smith .seealso: `PCMG`, `PCGAMG`
302e68589bSMark F. Adams @*/
319371c9d4SSatish Balay PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n) {
322e68589bSMark F. Adams   PetscFunctionBegin;
332e68589bSMark F. Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
34d5d25401SBarry Smith   PetscValidLogicalCollectiveInt(pc, n, 2);
35cac4c232SBarry Smith   PetscTryMethod(pc, "PCGAMGSetNSmooths_C", (PC, PetscInt), (pc, n));
362e68589bSMark F. Adams   PetscFunctionReturn(0);
372e68589bSMark F. Adams }
382e68589bSMark F. Adams 
399371c9d4SSatish Balay static PetscErrorCode PCGAMGSetNSmooths_AGG(PC pc, PetscInt n) {
402e68589bSMark F. Adams   PC_MG       *mg          = (PC_MG *)pc->data;
412e68589bSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
42c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
432e68589bSMark F. Adams 
442e68589bSMark F. Adams   PetscFunctionBegin;
45c8b0795cSMark F. Adams   pc_gamg_agg->nsmooths = n;
46c8b0795cSMark F. Adams   PetscFunctionReturn(0);
47c8b0795cSMark F. Adams }
48c8b0795cSMark F. Adams 
49c8b0795cSMark F. Adams /*@
50bae903cbSmarkadams4    PCGAMGSetSymmetrizeGraph - Symmetrize the graph before computing the aggregation. Some algorithms require the graph be symmetric
51c8b0795cSMark F. Adams 
52*f1580f4eSBarry Smith    Logically Collective on pc
53c8b0795cSMark F. Adams 
54c8b0795cSMark F. Adams    Input Parameters:
55cab9ed1eSBarry Smith +  pc - the preconditioner context
56*f1580f4eSBarry Smith -  n - `PETSC_TRUE` or `PETSC_FALSE`
57c8b0795cSMark F. Adams 
58c8b0795cSMark F. Adams    Options Database Key:
59bae903cbSmarkadams4 .  -pc_gamg_symmetrize_graph <true,default=false> - symmetrize the graph before computing the aggregation
60c8b0795cSMark F. Adams 
61c8b0795cSMark F. Adams    Level: intermediate
62c8b0795cSMark F. Adams 
63*f1580f4eSBarry Smith    Developer Note:
64*f1580f4eSBarry Smith    If the aggregation can hang with a nonsymmetric graph then the graph should always be symmetrized before calling the aggregation,
65*f1580f4eSBarry Smith    it should not be up to the user.
66*f1580f4eSBarry Smith 
67*f1580f4eSBarry Smith .seealso: `PCGAMG`, `PCGAMGSetAggressiveLevels()`
68c8b0795cSMark F. Adams @*/
699371c9d4SSatish Balay PetscErrorCode PCGAMGSetSymmetrizeGraph(PC pc, PetscBool n) {
70c8b0795cSMark F. Adams   PetscFunctionBegin;
71c8b0795cSMark F. Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
72d5d25401SBarry Smith   PetscValidLogicalCollectiveBool(pc, n, 2);
73bae903cbSmarkadams4   PetscTryMethod(pc, "PCGAMGSetSymmetrizeGraph_C", (PC, PetscBool), (pc, n));
74c8b0795cSMark F. Adams   PetscFunctionReturn(0);
75c8b0795cSMark F. Adams }
76c8b0795cSMark F. Adams 
779371c9d4SSatish Balay static PetscErrorCode PCGAMGSetSymmetrizeGraph_AGG(PC pc, PetscBool n) {
78c8b0795cSMark F. Adams   PC_MG       *mg          = (PC_MG *)pc->data;
79c8b0795cSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
80c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
81c8b0795cSMark F. Adams 
82c8b0795cSMark F. Adams   PetscFunctionBegin;
83bae903cbSmarkadams4   pc_gamg_agg->symmetrize_graph = n;
842e68589bSMark F. Adams   PetscFunctionReturn(0);
852e68589bSMark F. Adams }
862e68589bSMark F. Adams 
87ef4ad70eSMark F. Adams /*@
88*f1580f4eSBarry Smith    PCGAMGSetAggressiveLevels -  Use aggressive coarsening on first n levels
89ef4ad70eSMark F. Adams 
90*f1580f4eSBarry Smith    Logically Collective on pc
91ef4ad70eSMark F. Adams 
92ef4ad70eSMark F. Adams    Input Parameters:
93cab9ed1eSBarry Smith +  pc - the preconditioner context
94d5d25401SBarry Smith -  n - 0, 1 or more
95ef4ad70eSMark F. Adams 
96ef4ad70eSMark F. Adams    Options Database Key:
97bae903cbSmarkadams4 .  -pc_gamg_aggressive_coarsening <n,default = 1> - Number of levels to square the graph on before aggregating it
98af3c827dSMark Adams 
99ef4ad70eSMark F. Adams    Level: intermediate
100ef4ad70eSMark F. Adams 
101*f1580f4eSBarry Smith .seealso: `PCGAMG`, `PCGAMGSetSymmetrizeGraph()`, `PCGAMGSetThreshold()`
102ef4ad70eSMark F. Adams @*/
1039371c9d4SSatish Balay PetscErrorCode PCGAMGSetAggressiveLevels(PC pc, PetscInt n) {
104ef4ad70eSMark F. Adams   PetscFunctionBegin;
105ef4ad70eSMark F. Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
106d5d25401SBarry Smith   PetscValidLogicalCollectiveInt(pc, n, 2);
107bae903cbSmarkadams4   PetscTryMethod(pc, "PCGAMGSetAggressiveLevels_C", (PC, PetscInt), (pc, n));
108ef4ad70eSMark F. Adams   PetscFunctionReturn(0);
109ef4ad70eSMark F. Adams }
110ef4ad70eSMark F. Adams 
1119371c9d4SSatish Balay static PetscErrorCode PCGAMGSetAggressiveLevels_AGG(PC pc, PetscInt n) {
112ef4ad70eSMark F. Adams   PC_MG       *mg          = (PC_MG *)pc->data;
113ef4ad70eSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
114ef4ad70eSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
115ef4ad70eSMark F. Adams 
116ef4ad70eSMark F. Adams   PetscFunctionBegin;
117bae903cbSmarkadams4   pc_gamg_agg->aggressive_coarsening_levels = n;
118ef4ad70eSMark F. Adams   PetscFunctionReturn(0);
119ef4ad70eSMark F. Adams }
120ef4ad70eSMark F. Adams 
1219371c9d4SSatish Balay static PetscErrorCode PCSetFromOptions_GAMG_AGG(PC pc, PetscOptionItems *PetscOptionsObject) {
1222e68589bSMark F. Adams   PC_MG       *mg          = (PC_MG *)pc->data;
1232e68589bSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
124c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
1252e68589bSMark F. Adams 
1262e68589bSMark F. Adams   PetscFunctionBegin;
127d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "GAMG-AGG options");
1282e68589bSMark F. Adams   {
129bae903cbSmarkadams4     PetscBool flg;
1309566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-pc_gamg_agg_nsmooths", "smoothing steps for smoothed aggregation, usually 1", "PCGAMGSetNSmooths", pc_gamg_agg->nsmooths, &pc_gamg_agg->nsmooths, NULL));
131bae903cbSmarkadams4     PetscCall(PetscOptionsBool("-pc_gamg_symmetrize_graph", "Set for asymmetric matrices", "PCGAMGSetSymmetrizeGraph", pc_gamg_agg->symmetrize_graph, &pc_gamg_agg->symmetrize_graph, NULL));
132bae903cbSmarkadams4     pc_gamg_agg->aggressive_coarsening_levels = 1;
1339371c9d4SSatish Balay     PetscCall(
1349371c9d4SSatish 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));
135bae903cbSmarkadams4     if (!flg) {
136bae903cbSmarkadams4       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));
137bae903cbSmarkadams4     } else {
138bae903cbSmarkadams4       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));
139bae903cbSmarkadams4       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));
140bae903cbSmarkadams4     }
1412e68589bSMark F. Adams   }
142d0609cedSBarry Smith   PetscOptionsHeadEnd();
1432e68589bSMark F. Adams   PetscFunctionReturn(0);
1442e68589bSMark F. Adams }
1452e68589bSMark F. Adams 
1469371c9d4SSatish Balay static PetscErrorCode PCDestroy_GAMG_AGG(PC pc) {
1472e68589bSMark F. Adams   PC_MG   *mg      = (PC_MG *)pc->data;
1482e68589bSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1492e68589bSMark F. Adams 
1502e68589bSMark F. Adams   PetscFunctionBegin;
1519566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg->subctx));
1522e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNSmooths_C", NULL));
153bae903cbSmarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetSymmetrizeGraph_C", NULL));
154bae903cbSmarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveLevels_C", NULL));
155bae903cbSmarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", NULL));
1562e68589bSMark F. Adams   PetscFunctionReturn(0);
1572e68589bSMark F. Adams }
1582e68589bSMark F. Adams 
1592e68589bSMark F. Adams /*
1602e68589bSMark F. Adams    PCSetCoordinates_AGG
161302f38e8SMark F. Adams      - collective
1622e68589bSMark F. Adams 
1632e68589bSMark F. Adams    Input Parameter:
1642e68589bSMark F. Adams    . pc - the preconditioner context
165a2f3521dSMark F. Adams    . ndm - dimesion of data (used for dof/vertex for Stokes)
166302f38e8SMark F. Adams    . a_nloc - number of vertices local
167302f38e8SMark F. Adams    . coords - [a_nloc][ndm] - interleaved coordinate data: {x_0, y_0, z_0, x_1, y_1, ...}
1682e68589bSMark F. Adams */
1691e6b0712SBarry Smith 
1709371c9d4SSatish Balay static PetscErrorCode PCSetCoordinates_AGG(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords) {
1712e68589bSMark F. Adams   PC_MG   *mg      = (PC_MG *)pc->data;
1722e68589bSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
17369344418SMark F. Adams   PetscInt arrsz, kk, ii, jj, nloc, ndatarows, ndf;
174a2f3521dSMark F. Adams   Mat      mat = pc->pmat;
1752e68589bSMark F. Adams 
1762e68589bSMark F. Adams   PetscFunctionBegin;
177a2f3521dSMark F. Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
178a2f3521dSMark F. Adams   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
179302f38e8SMark F. Adams   nloc = a_nloc;
1802e68589bSMark F. Adams 
1812e68589bSMark F. Adams   /* SA: null space vectors */
1829566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(mat, &ndf));               /* this does not work for Stokes */
18369344418SMark F. Adams   if (coords && ndf == 1) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */
184a2f3521dSMark F. Adams   else if (coords) {
18563a3b9bcSJacob Faibussowitsch     PetscCheck(ndm <= ndf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "degrees of motion %" PetscInt_FMT " > block size %" PetscInt_FMT, ndm, ndf);
18669344418SMark F. Adams     pc_gamg->data_cell_cols = (ndm == 2 ? 3 : 6); /* displacement elasticity */
187ad540459SPierre 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);
188806fa848SBarry Smith   } else pc_gamg->data_cell_cols = ndf; /* no data, force SA with constant null space vectors */
18971959b99SBarry Smith   pc_gamg->data_cell_rows = ndatarows = ndf;
19063a3b9bcSJacob 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);
191c8b0795cSMark F. Adams   arrsz = nloc * pc_gamg->data_cell_rows * pc_gamg->data_cell_cols;
1922e68589bSMark F. Adams 
1937f66b68fSMark Adams   if (!pc_gamg->data || (pc_gamg->data_sz != arrsz)) {
1949566063dSJacob Faibussowitsch     PetscCall(PetscFree(pc_gamg->data));
1959566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(arrsz + 1, &pc_gamg->data));
1962e68589bSMark F. Adams   }
1972e68589bSMark F. Adams   /* copy data in - column oriented */
1982e68589bSMark F. Adams   for (kk = 0; kk < nloc; kk++) {
19969344418SMark F. Adams     const PetscInt M    = nloc * pc_gamg->data_cell_rows; /* stride into data */
20069344418SMark F. Adams     PetscReal     *data = &pc_gamg->data[kk * ndatarows]; /* start of cell */
201c8b0795cSMark F. Adams     if (pc_gamg->data_cell_cols == 1) *data = 1.0;
2022e68589bSMark F. Adams     else {
20369344418SMark F. Adams       /* translational modes */
2042fa5cd67SKarl Rupp       for (ii = 0; ii < ndatarows; ii++) {
2052fa5cd67SKarl Rupp         for (jj = 0; jj < ndatarows; jj++) {
20669344418SMark F. Adams           if (ii == jj) data[ii * M + jj] = 1.0;
2072e68589bSMark F. Adams           else data[ii * M + jj] = 0.0;
2082fa5cd67SKarl Rupp         }
2092fa5cd67SKarl Rupp       }
21069344418SMark F. Adams 
21169344418SMark F. Adams       /* rotational modes */
2122e68589bSMark F. Adams       if (coords) {
21369344418SMark F. Adams         if (ndm == 2) {
2142e68589bSMark F. Adams           data += 2 * M;
2152e68589bSMark F. Adams           data[0] = -coords[2 * kk + 1];
2162e68589bSMark F. Adams           data[1] = coords[2 * kk];
217806fa848SBarry Smith         } else {
2182e68589bSMark F. Adams           data += 3 * M;
2199371c9d4SSatish Balay           data[0]         = 0.0;
2209371c9d4SSatish Balay           data[M + 0]     = coords[3 * kk + 2];
2219371c9d4SSatish Balay           data[2 * M + 0] = -coords[3 * kk + 1];
2229371c9d4SSatish Balay           data[1]         = -coords[3 * kk + 2];
2239371c9d4SSatish Balay           data[M + 1]     = 0.0;
2249371c9d4SSatish Balay           data[2 * M + 1] = coords[3 * kk];
2259371c9d4SSatish Balay           data[2]         = coords[3 * kk + 1];
2269371c9d4SSatish Balay           data[M + 2]     = -coords[3 * kk];
2279371c9d4SSatish Balay           data[2 * M + 2] = 0.0;
2282e68589bSMark F. Adams         }
2292e68589bSMark F. Adams       }
2302e68589bSMark F. Adams     }
2312e68589bSMark F. Adams   }
2322e68589bSMark F. Adams   pc_gamg->data_sz = arrsz;
2332e68589bSMark F. Adams   PetscFunctionReturn(0);
2342e68589bSMark F. Adams }
2352e68589bSMark F. Adams 
2362e68589bSMark F. Adams /*
237a2f3521dSMark F. Adams    PCSetData_AGG - called if data is not set with PCSetCoordinates.
238a2f3521dSMark F. Adams       Looks in Mat for near null space.
239a2f3521dSMark F. Adams       Does not work for Stokes
2402e68589bSMark F. Adams 
2412e68589bSMark F. Adams   Input Parameter:
2422e68589bSMark F. Adams    . pc -
243a2f3521dSMark F. Adams    . a_A - matrix to get (near) null space out of.
2442e68589bSMark F. Adams */
2459371c9d4SSatish Balay static PetscErrorCode PCSetData_AGG(PC pc, Mat a_A) {
246b8cd405aSMark F. Adams   PC_MG       *mg      = (PC_MG *)pc->data;
247b8cd405aSMark F. Adams   PC_GAMG     *pc_gamg = (PC_GAMG *)mg->innerctx;
248b8cd405aSMark F. Adams   MatNullSpace mnull;
24966f2ef4bSMark Adams 
250b817416eSBarry Smith   PetscFunctionBegin;
2519566063dSJacob Faibussowitsch   PetscCall(MatGetNearNullSpace(a_A, &mnull));
252b8cd405aSMark F. Adams   if (!mnull) {
25366f2ef4bSMark Adams     DM dm;
2549566063dSJacob Faibussowitsch     PetscCall(PCGetDM(pc, &dm));
25548a46eb9SPierre Jolivet     if (!dm) PetscCall(MatGetDM(a_A, &dm));
25666f2ef4bSMark Adams     if (dm) {
25766f2ef4bSMark Adams       PetscObject deformation;
258b0d30dd6SMatthew G. Knepley       PetscInt    Nf;
259b0d30dd6SMatthew G. Knepley 
2609566063dSJacob Faibussowitsch       PetscCall(DMGetNumFields(dm, &Nf));
261b0d30dd6SMatthew G. Knepley       if (Nf) {
2629566063dSJacob Faibussowitsch         PetscCall(DMGetField(dm, 0, NULL, &deformation));
2639566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery((PetscObject)deformation, "nearnullspace", (PetscObject *)&mnull));
26448a46eb9SPierre Jolivet         if (!mnull) PetscCall(PetscObjectQuery((PetscObject)deformation, "nullspace", (PetscObject *)&mnull));
26566f2ef4bSMark Adams       }
26666f2ef4bSMark Adams     }
267b0d30dd6SMatthew G. Knepley   }
26866f2ef4bSMark Adams 
26966f2ef4bSMark Adams   if (!mnull) {
270a2f3521dSMark F. Adams     PetscInt bs, NN, MM;
2719566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSize(a_A, &bs));
2729566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(a_A, &MM, &NN));
27363a3b9bcSJacob Faibussowitsch     PetscCheck(MM % bs == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MM %" PetscInt_FMT " must be divisible by bs %" PetscInt_FMT, MM, bs);
2749566063dSJacob Faibussowitsch     PetscCall(PCSetCoordinates_AGG(pc, bs, MM / bs, NULL));
275806fa848SBarry Smith   } else {
276b8cd405aSMark F. Adams     PetscReal         *nullvec;
277b8cd405aSMark F. Adams     PetscBool          has_const;
278b8cd405aSMark F. Adams     PetscInt           i, j, mlocal, nvec, bs;
279d5d25401SBarry Smith     const Vec         *vecs;
280d5d25401SBarry Smith     const PetscScalar *v;
281b817416eSBarry Smith 
2829566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(a_A, &mlocal, NULL));
2839566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs));
284a0dea87dSMark Adams     pc_gamg->data_sz = (nvec + !!has_const) * mlocal;
2859566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1((nvec + !!has_const) * mlocal, &nullvec));
2869371c9d4SSatish Balay     if (has_const)
2879371c9d4SSatish Balay       for (i = 0; i < mlocal; i++) nullvec[i] = 1.0;
288b8cd405aSMark F. Adams     for (i = 0; i < nvec; i++) {
2899566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(vecs[i], &v));
290b8cd405aSMark F. Adams       for (j = 0; j < mlocal; j++) nullvec[(i + !!has_const) * mlocal + j] = PetscRealPart(v[j]);
2919566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(vecs[i], &v));
292b8cd405aSMark F. Adams     }
293b8cd405aSMark F. Adams     pc_gamg->data           = nullvec;
294b8cd405aSMark F. Adams     pc_gamg->data_cell_cols = (nvec + !!has_const);
2959566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSize(a_A, &bs));
296b8cd405aSMark F. Adams     pc_gamg->data_cell_rows = bs;
297b8cd405aSMark F. Adams   }
2982e68589bSMark F. Adams   PetscFunctionReturn(0);
2992e68589bSMark F. Adams }
3002e68589bSMark F. Adams 
3012e68589bSMark F. Adams /*
302bae903cbSmarkadams4   formProl0 - collect null space data for each aggregate, do QR, put R in coarse grid data and Q in P_0
3032e68589bSMark F. Adams 
3042e68589bSMark F. Adams   Input Parameter:
3052adfe9d3SBarry Smith    . agg_llists - list of arrays with aggregates -- list from selected vertices of aggregate unselected vertices
3062adfe9d3SBarry Smith    . bs - row block size
3070cbbd2e1SMark F. Adams    . nSAvec - column bs of new P
3080cbbd2e1SMark F. Adams    . my0crs - global index of start of locals
3092adfe9d3SBarry Smith    . data_stride - bs*(nloc nodes + ghost nodes) [data_stride][nSAvec]
3102e68589bSMark F. Adams    . data_in[data_stride*nSAvec] - local data on fine grid
3112e68589bSMark F. Adams    . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist'
312bae903cbSmarkadams4 
3132e68589bSMark F. Adams   Output Parameter:
3142e68589bSMark F. Adams    . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data
3152e68589bSMark F. Adams    . a_Prol - prolongation operator
3162e68589bSMark F. Adams */
3179371c9d4SSatish Balay 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) {
318bd026e97SJed Brown   PetscInt        Istart, my0, Iend, nloc, clid, flid = 0, aggID, kk, jj, ii, mm, nSelected, minsz, nghosts, out_data_stride;
3193b4367a7SBarry Smith   MPI_Comm        comm;
3202e68589bSMark F. Adams   PetscReal      *out_data;
321539c167fSBarry Smith   PetscCDIntNd   *pos;
3221943db53SBarry Smith   PCGAMGHashTable fgid_flid;
3230cbbd2e1SMark F. Adams 
3242e68589bSMark F. Adams   PetscFunctionBegin;
3259566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)a_Prol, &comm));
3269566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(a_Prol, &Istart, &Iend));
3279371c9d4SSatish Balay   nloc = (Iend - Istart) / bs;
3289371c9d4SSatish Balay   my0  = Istart / bs;
32963a3b9bcSJacob 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);
3300cbbd2e1SMark F. Adams   Iend /= bs;
3310cbbd2e1SMark F. Adams   nghosts = data_stride / bs - nloc;
3322e68589bSMark F. Adams 
3339566063dSJacob Faibussowitsch   PetscCall(PCGAMGHashTableCreate(2 * nghosts + 1, &fgid_flid));
33448a46eb9SPierre Jolivet   for (kk = 0; kk < nghosts; kk++) PetscCall(PCGAMGHashTableAdd(&fgid_flid, flid_fgid[nloc + kk], nloc + kk));
3352e68589bSMark F. Adams 
3360cbbd2e1SMark F. Adams   /* count selected -- same as number of cols of P */
3370cbbd2e1SMark F. Adams   for (nSelected = mm = 0; mm < nloc; mm++) {
338e78576d6SMark F. Adams     PetscBool ise;
3399566063dSJacob Faibussowitsch     PetscCall(PetscCDEmptyAt(agg_llists, mm, &ise));
340e78576d6SMark F. Adams     if (!ise) nSelected++;
3410cbbd2e1SMark F. Adams   }
3429566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRangeColumn(a_Prol, &ii, &jj));
34363a3b9bcSJacob Faibussowitsch   PetscCheck((ii / nSAvec) == my0crs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ii %" PetscInt_FMT " /nSAvec %" PetscInt_FMT "  != my0crs %" PetscInt_FMT, ii, nSAvec, my0crs);
34463a3b9bcSJacob 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);
3450cbbd2e1SMark F. Adams 
3462e68589bSMark F. Adams   /* aloc space for coarse point data (output) */
3470cbbd2e1SMark F. Adams   out_data_stride = nSelected * nSAvec;
3482fa5cd67SKarl Rupp 
3499566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(out_data_stride * nSAvec, &out_data));
35033812677SSatish Balay   for (ii = 0; ii < out_data_stride * nSAvec; ii++) out_data[ii] = PETSC_MAX_REAL;
3510cbbd2e1SMark F. Adams   *a_data_out = out_data; /* output - stride nSelected*nSAvec */
3522e68589bSMark F. Adams 
3532e68589bSMark F. Adams   /* find points and set prolongation */
354c8b0795cSMark F. Adams   minsz = 100;
3550cbbd2e1SMark F. Adams   for (mm = clid = 0; mm < nloc; mm++) {
3569566063dSJacob Faibussowitsch     PetscCall(PetscCDSizeAt(agg_llists, mm, &jj));
357e78576d6SMark F. Adams     if (jj > 0) {
3580cbbd2e1SMark F. Adams       const PetscInt lid = mm, cgid = my0crs + clid;
3590cbbd2e1SMark F. Adams       PetscInt       cids[100]; /* max bs */
3600cbbd2e1SMark F. Adams       PetscBLASInt   asz = jj, M = asz * bs, N = nSAvec, INFO;
3612e68589bSMark F. Adams       PetscBLASInt   Mdata = M + ((N - M > 0) ? N - M : 0), LDA = Mdata, LWORK = N * bs;
3622e68589bSMark F. Adams       PetscScalar   *qqc, *qqr, *TAU, *WORK;
3632e68589bSMark F. Adams       PetscInt      *fids;
36465d7b583SSatish Balay       PetscReal     *data;
365b817416eSBarry Smith 
3660cbbd2e1SMark F. Adams       /* count agg */
3670cbbd2e1SMark F. Adams       if (asz < minsz) minsz = asz;
3680cbbd2e1SMark F. Adams 
3690cbbd2e1SMark F. Adams       /* get block */
3709566063dSJacob Faibussowitsch       PetscCall(PetscMalloc5(Mdata * N, &qqc, M * N, &qqr, N, &TAU, LWORK, &WORK, M, &fids));
3712e68589bSMark F. Adams 
3722e68589bSMark F. Adams       aggID = 0;
3739566063dSJacob Faibussowitsch       PetscCall(PetscCDGetHeadPos(agg_llists, lid, &pos));
374e78576d6SMark F. Adams       while (pos) {
375e78576d6SMark F. Adams         PetscInt gid1;
3769566063dSJacob Faibussowitsch         PetscCall(PetscCDIntNdGetID(pos, &gid1));
3779566063dSJacob Faibussowitsch         PetscCall(PetscCDGetNextPos(agg_llists, lid, &pos));
378e78576d6SMark F. Adams 
3790cbbd2e1SMark F. Adams         if (gid1 >= my0 && gid1 < Iend) flid = gid1 - my0;
3800cbbd2e1SMark F. Adams         else {
3819566063dSJacob Faibussowitsch           PetscCall(PCGAMGHashTableFind(&fgid_flid, gid1, &flid));
38208401ef6SPierre Jolivet           PetscCheck(flid >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cannot find gid1 in table");
3830cbbd2e1SMark F. Adams         }
3842e68589bSMark F. Adams         /* copy in B_i matrix - column oriented */
38565d7b583SSatish Balay         data = &data_in[flid * bs];
3863d3eaba7SBarry Smith         for (ii = 0; ii < bs; ii++) {
3872e68589bSMark F. Adams           for (jj = 0; jj < N; jj++) {
3880cbbd2e1SMark F. Adams             PetscReal d                       = data[jj * data_stride + ii];
3890cbbd2e1SMark F. Adams             qqc[jj * Mdata + aggID * bs + ii] = d;
3902e68589bSMark F. Adams           }
3912e68589bSMark F. Adams         }
3922e68589bSMark F. Adams         /* set fine IDs */
3932e68589bSMark F. Adams         for (kk = 0; kk < bs; kk++) fids[aggID * bs + kk] = flid_fgid[flid] * bs + kk;
3942e68589bSMark F. Adams         aggID++;
3950cbbd2e1SMark F. Adams       }
3962e68589bSMark F. Adams 
3972e68589bSMark F. Adams       /* pad with zeros */
3982e68589bSMark F. Adams       for (ii = asz * bs; ii < Mdata; ii++) {
399ad540459SPierre Jolivet         for (jj = 0; jj < N; jj++, kk++) qqc[jj * Mdata + ii] = .0;
4002e68589bSMark F. Adams       }
4012e68589bSMark F. Adams 
4022e68589bSMark F. Adams       /* QR */
4039566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
404792fecdfSBarry Smith       PetscCallBLAS("LAPACKgeqrf", LAPACKgeqrf_(&Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
4059566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPop());
40608401ef6SPierre Jolivet       PetscCheck(INFO == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xGEQRF error");
4072e68589bSMark F. Adams       /* get R - column oriented - output B_{i+1} */
4082e68589bSMark F. Adams       {
4092e68589bSMark F. Adams         PetscReal *data = &out_data[clid * nSAvec];
4102e68589bSMark F. Adams         for (jj = 0; jj < nSAvec; jj++) {
4112e68589bSMark F. Adams           for (ii = 0; ii < nSAvec; ii++) {
41208401ef6SPierre 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);
4130cbbd2e1SMark F. Adams             if (ii <= jj) data[jj * out_data_stride + ii] = PetscRealPart(qqc[jj * Mdata + ii]);
4140cbbd2e1SMark F. Adams             else data[jj * out_data_stride + ii] = 0.;
4152e68589bSMark F. Adams           }
4162e68589bSMark F. Adams         }
4172e68589bSMark F. Adams       }
4182e68589bSMark F. Adams 
4192e68589bSMark F. Adams       /* get Q - row oriented */
420792fecdfSBarry Smith       PetscCallBLAS("LAPACKorgqr", LAPACKorgqr_(&Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
42163a3b9bcSJacob Faibussowitsch       PetscCheck(INFO == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xORGQR error arg %" PetscBLASInt_FMT, -INFO);
4222e68589bSMark F. Adams 
4232e68589bSMark F. Adams       for (ii = 0; ii < M; ii++) {
424ad540459SPierre Jolivet         for (jj = 0; jj < N; jj++) qqr[N * ii + jj] = qqc[jj * Mdata + ii];
4252e68589bSMark F. Adams       }
4262e68589bSMark F. Adams 
4272e68589bSMark F. Adams       /* add diagonal block of P0 */
4289371c9d4SSatish Balay       for (kk = 0; kk < N; kk++) { cids[kk] = N * cgid + kk; /* global col IDs in P0 */ }
4299566063dSJacob Faibussowitsch       PetscCall(MatSetValues(a_Prol, M, fids, N, cids, qqr, INSERT_VALUES));
4309566063dSJacob Faibussowitsch       PetscCall(PetscFree5(qqc, qqr, TAU, WORK, fids));
431b43b03e9SMark F. Adams       clid++;
4320cbbd2e1SMark F. Adams     } /* coarse agg */
4330cbbd2e1SMark F. Adams   }   /* for all fine nodes */
4349566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(a_Prol, MAT_FINAL_ASSEMBLY));
4359566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(a_Prol, MAT_FINAL_ASSEMBLY));
4369566063dSJacob Faibussowitsch   PetscCall(PCGAMGHashTableDestroy(&fgid_flid));
4372e68589bSMark F. Adams   PetscFunctionReturn(0);
4382e68589bSMark F. Adams }
4392e68589bSMark F. Adams 
4409371c9d4SSatish Balay static PetscErrorCode PCView_GAMG_AGG(PC pc, PetscViewer viewer) {
4415adeb434SBarry Smith   PC_MG       *mg          = (PC_MG *)pc->data;
4425adeb434SBarry Smith   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
4435adeb434SBarry Smith   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
4445adeb434SBarry Smith 
4455adeb434SBarry Smith   PetscFunctionBegin;
4469566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "      AGG specific options\n"));
447bae903cbSmarkadams4   PetscCall(PetscViewerASCIIPrintf(viewer, "        Symmetric graph %s\n", pc_gamg_agg->symmetrize_graph ? "true" : "false"));
448bae903cbSmarkadams4   PetscCall(PetscViewerASCIIPrintf(viewer, "        Number of levels to square graph %d\n", (int)pc_gamg_agg->aggressive_coarsening_levels));
449bae903cbSmarkadams4   PetscCall(PetscViewerASCIIPrintf(viewer, "        Number smoothing steps %d\n", (int)pc_gamg_agg->nsmooths));
4505adeb434SBarry Smith   PetscFunctionReturn(0);
4515adeb434SBarry Smith }
4525adeb434SBarry Smith 
4532e68589bSMark F. Adams /*
454fd1112cbSBarry Smith    PCGAMGGraph_AGG
4552e68589bSMark F. Adams 
4562e68589bSMark F. Adams   Input Parameter:
4572e68589bSMark F. Adams    . pc - this
4582e68589bSMark F. Adams    . Amat - matrix on this fine level
459bae903cbSmarkadams4 
4602e68589bSMark F. Adams   Output Parameter:
461c8b0795cSMark F. Adams    . a_Gmat -
4622e68589bSMark F. Adams */
4639371c9d4SSatish Balay static PetscErrorCode PCGAMGGraph_AGG(PC pc, Mat Amat, Mat *a_Gmat) {
464c8b0795cSMark F. Adams   PC_MG                    *mg          = (PC_MG *)pc->data;
465c8b0795cSMark F. Adams   PC_GAMG                  *pc_gamg     = (PC_GAMG *)mg->innerctx;
466c1eae691SMark Adams   const PetscReal           vfilter     = pc_gamg->threshold[pc_gamg->current_level];
467c8b0795cSMark F. Adams   PC_GAMG_AGG              *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
46872833a62Smarkadams4   Mat                       Gmat, F = NULL;
4693b4367a7SBarry Smith   MPI_Comm                  comm;
47067744fedSMark F. Adams   PetscBool /* set,flg , */ symm;
471c8b0795cSMark F. Adams 
472c8b0795cSMark F. Adams   PetscFunctionBegin;
4739566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm));
474c8b0795cSMark F. Adams 
4759566063dSJacob Faibussowitsch   /* PetscCall(MatIsSymmetricKnown(Amat, &set, &flg)); || !(set && flg) -- this causes lot of symm calls */
476bae903cbSmarkadams4   symm = (PetscBool)(pc_gamg_agg->symmetrize_graph); /* && !pc_gamg_agg->square_graph; */
4770cbbd2e1SMark F. Adams 
478849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_GRAPH], 0, 0, 0, 0));
47972833a62Smarkadams4   PetscCall(MatCreateGraph(Amat, symm, PETSC_TRUE, &Gmat));
480849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_GRAPH], 0, 0, 0, 0));
481849bee69SMark Adams 
482849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_FILTER], 0, 0, 0, 0));
48372833a62Smarkadams4   PetscCall(MatFilter(Gmat, vfilter, &F));
48472833a62Smarkadams4   if (F) {
48572833a62Smarkadams4     PetscCall(MatDestroy(&Gmat));
48672833a62Smarkadams4     Gmat = F;
48772833a62Smarkadams4   }
488849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_FILTER], 0, 0, 0, 0));
489849bee69SMark Adams 
490e0940f08SMark F. Adams   *a_Gmat = Gmat;
491849bee69SMark Adams 
492c8b0795cSMark F. Adams   PetscFunctionReturn(0);
493c8b0795cSMark F. Adams }
494c8b0795cSMark F. Adams 
495c8b0795cSMark F. Adams /*
496bae903cbSmarkadams4    PCGAMGCoarsen_AGG - supports squaring the graph (deprecated) and new graph for
497bae903cbSmarkadams4      communication of QR data used with HEM and MISk coarsening
498c8b0795cSMark F. Adams 
499c8b0795cSMark F. Adams   Input Parameter:
500e0940f08SMark F. Adams    . a_pc - this
501bae903cbSmarkadams4 
502e0940f08SMark F. Adams   Input/Output Parameter:
503bae903cbSmarkadams4    . a_Gmat1 - graph to coarsen (in), graph off processor edges for QR gather scatter (out)
504bae903cbSmarkadams4 
505c8b0795cSMark F. Adams   Output Parameter:
5060cbbd2e1SMark F. Adams    . agg_lists - list of aggregates
507bae903cbSmarkadams4 
508c8b0795cSMark F. Adams */
5099371c9d4SSatish Balay static PetscErrorCode PCGAMGCoarsen_AGG(PC a_pc, Mat *a_Gmat1, PetscCoarsenData **agg_lists) {
510e0940f08SMark F. Adams   PC_MG       *mg          = (PC_MG *)a_pc->data;
511c8b0795cSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
512c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
513bae903cbSmarkadams4   Mat          mat, Gmat1 = *a_Gmat1; /* aggressive graph */
5140cbbd2e1SMark F. Adams   IS           perm;
515bae903cbSmarkadams4   PetscInt     Istart, Iend, Ii, nloc, bs, nn;
516bae903cbSmarkadams4   PetscInt    *permute, *degree;
517c8b0795cSMark F. Adams   PetscBool   *bIndexSet;
518b43b03e9SMark F. Adams   MatCoarsen   crs;
5193b4367a7SBarry Smith   MPI_Comm     comm;
520aea53286SMark Adams   PetscReal    hashfact;
521e2c00dcbSBarry Smith   PetscInt     iSwapIndex;
5223b50add6SMark Adams   PetscRandom  random;
523c8b0795cSMark F. Adams 
524c8b0795cSMark F. Adams   PetscFunctionBegin;
525849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0));
5269566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Gmat1, &comm));
527bae903cbSmarkadams4   PetscCall(MatGetLocalSize(Gmat1, &nn, NULL));
5289566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(Gmat1, &bs));
52963a3b9bcSJacob Faibussowitsch   PetscCheck(bs == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "bs %" PetscInt_FMT " must be 1", bs);
530bae903cbSmarkadams4   nloc = nn / bs;
531c8b0795cSMark F. Adams 
5325cfd4bd9SMark Adams   /* get MIS aggs - randomize */
533bae903cbSmarkadams4   PetscCall(PetscMalloc2(nloc, &permute, nloc, &degree));
5349566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(nloc, &bIndexSet));
5356e09b0e3SStefano Zampini   for (Ii = 0; Ii < nloc; Ii++) permute[Ii] = Ii;
5369566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &random));
5379566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Gmat1, &Istart, &Iend));
538c8b0795cSMark F. Adams   for (Ii = 0; Ii < nloc; Ii++) {
539bae903cbSmarkadams4     PetscInt nc;
540bae903cbSmarkadams4     PetscCall(MatGetRow(Gmat1, Istart + Ii, &nc, NULL, NULL));
541bae903cbSmarkadams4     degree[Ii] = nc;
542bae903cbSmarkadams4     PetscCall(MatRestoreRow(Gmat1, Istart + Ii, &nc, NULL, NULL));
543bae903cbSmarkadams4   }
544bae903cbSmarkadams4   for (Ii = 0; Ii < nloc; Ii++) {
5459566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValueReal(random, &hashfact));
546aea53286SMark Adams     iSwapIndex = (PetscInt)(hashfact * nloc) % nloc;
547c8b0795cSMark F. Adams     if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) {
548c8b0795cSMark F. Adams       PetscInt iTemp        = permute[iSwapIndex];
549c8b0795cSMark F. Adams       permute[iSwapIndex]   = permute[Ii];
550c8b0795cSMark F. Adams       permute[Ii]           = iTemp;
551bae903cbSmarkadams4       iTemp                 = degree[iSwapIndex];
552bae903cbSmarkadams4       degree[iSwapIndex]    = degree[Ii];
553bae903cbSmarkadams4       degree[Ii]            = iTemp;
554c8b0795cSMark F. Adams       bIndexSet[iSwapIndex] = PETSC_TRUE;
555c8b0795cSMark F. Adams     }
556c8b0795cSMark F. Adams   }
557bae903cbSmarkadams4   // create minimum degree ordering
558bae903cbSmarkadams4   PetscCall(PetscSortIntWithArray(nloc, degree, permute));
559bae903cbSmarkadams4 
5609566063dSJacob Faibussowitsch   PetscCall(PetscFree(bIndexSet));
5619566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&random));
5629566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm));
563849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_MIS], 0, 0, 0, 0));
5649566063dSJacob Faibussowitsch   PetscCall(MatCoarsenCreate(comm, &crs));
5659566063dSJacob Faibussowitsch   PetscCall(MatCoarsenSetFromOptions(crs));
5669566063dSJacob Faibussowitsch   PetscCall(MatCoarsenSetGreedyOrdering(crs, perm));
567bae903cbSmarkadams4   PetscCall(MatCoarsenSetAdjacency(crs, Gmat1));
5689566063dSJacob Faibussowitsch   PetscCall(MatCoarsenSetStrictAggs(crs, PETSC_TRUE));
569bae903cbSmarkadams4   if (pc_gamg->current_level < pc_gamg_agg->aggressive_coarsening_levels) PetscCall(MatCoarsenMISKSetDistance(crs, 2)); // hardwire to MIS-2
570bae903cbSmarkadams4   else PetscCall(MatCoarsenMISKSetDistance(crs, 1));                                                                    // MIS
5719566063dSJacob Faibussowitsch   PetscCall(MatCoarsenApply(crs));
572bae903cbSmarkadams4   PetscCall(MatCoarsenViewFromOptions(crs, NULL, "-mat_coarsen_view"));
5739566063dSJacob Faibussowitsch   PetscCall(MatCoarsenGetData(crs, agg_lists)); /* output */
5749566063dSJacob Faibussowitsch   PetscCall(MatCoarsenDestroy(&crs));
575b43b03e9SMark F. Adams 
5769566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
577bae903cbSmarkadams4   PetscCall(PetscFree2(permute, degree));
578849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_MIS], 0, 0, 0, 0));
5799f3f12c8SMark F. Adams 
580bae903cbSmarkadams4   {
581bae903cbSmarkadams4     PetscCoarsenData *llist = *agg_lists;
5829ab59c8bSMark Adams     /* see if we have a matrix that takes precedence (returned from MatCoarsenApply) */
5839566063dSJacob Faibussowitsch     PetscCall(PetscCDGetMat(llist, &mat));
5840cbbd2e1SMark F. Adams     if (mat) {
5859566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&Gmat1));
5860cbbd2e1SMark F. Adams       *a_Gmat1 = mat; /* output */
5870cbbd2e1SMark F. Adams     }
5880cbbd2e1SMark F. Adams   }
589849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0));
590c8b0795cSMark F. Adams   PetscFunctionReturn(0);
591c8b0795cSMark F. Adams }
592c8b0795cSMark F. Adams 
593c8b0795cSMark F. Adams /*
5940cbbd2e1SMark F. Adams  PCGAMGProlongator_AGG
595c8b0795cSMark F. Adams 
596c8b0795cSMark F. Adams  Input Parameter:
597c8b0795cSMark F. Adams  . pc - this
598c8b0795cSMark F. Adams  . Amat - matrix on this fine level
599c8b0795cSMark F. Adams  . Graph - used to get ghost data for nodes in
6000cbbd2e1SMark F. Adams  . agg_lists - list of aggregates
601c8b0795cSMark F. Adams  Output Parameter:
602c8b0795cSMark F. Adams  . a_P_out - prolongation operator to the next level
603c8b0795cSMark F. Adams  */
6049371c9d4SSatish Balay static PetscErrorCode PCGAMGProlongator_AGG(PC pc, Mat Amat, Mat Gmat, PetscCoarsenData *agg_lists, Mat *a_P_out) {
6052e68589bSMark F. Adams   PC_MG         *mg      = (PC_MG *)pc->data;
6062e68589bSMark F. Adams   PC_GAMG       *pc_gamg = (PC_GAMG *)mg->innerctx;
6074a2f8832SBarry Smith   const PetscInt col_bs  = pc_gamg->data_cell_cols;
608c8b0795cSMark F. Adams   PetscInt       Istart, Iend, nloc, ii, jj, kk, my0, nLocalSelected, bs;
609c8b0795cSMark F. Adams   Mat            Prol;
610d5d25401SBarry Smith   PetscMPIInt    size;
6113b4367a7SBarry Smith   MPI_Comm       comm;
612c8b0795cSMark F. Adams   PetscReal     *data_w_ghost;
613c8b0795cSMark F. Adams   PetscInt       myCrs0, nbnodes = 0, *flid_fgid;
614d9558ea9SBarry Smith   MatType        mtype;
6152e68589bSMark F. Adams 
6162e68589bSMark F. Adams   PetscFunctionBegin;
6179566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm));
61808401ef6SPierre Jolivet   PetscCheck(col_bs >= 1, comm, PETSC_ERR_PLIB, "Column bs cannot be less than 1");
619849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0));
6209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
6219566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Amat, &Istart, &Iend));
6229566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(Amat, &bs));
6239371c9d4SSatish Balay   nloc = (Iend - Istart) / bs;
6249371c9d4SSatish Balay   my0  = Istart / bs;
62563a3b9bcSJacob 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);
6262e68589bSMark F. Adams 
6272e68589bSMark F. Adams   /* get 'nLocalSelected' */
6280cbbd2e1SMark F. Adams   for (ii = 0, nLocalSelected = 0; ii < nloc; ii++) {
629e78576d6SMark F. Adams     PetscBool ise;
6300cbbd2e1SMark F. Adams     /* filter out singletons 0 or 1? */
6319566063dSJacob Faibussowitsch     PetscCall(PetscCDEmptyAt(agg_lists, ii, &ise));
632e78576d6SMark F. Adams     if (!ise) nLocalSelected++;
6332e68589bSMark F. Adams   }
6342e68589bSMark F. Adams 
6352e68589bSMark F. Adams   /* create prolongator, create P matrix */
6369566063dSJacob Faibussowitsch   PetscCall(MatGetType(Amat, &mtype));
6379566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &Prol));
6389566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Prol, nloc * bs, nLocalSelected * col_bs, PETSC_DETERMINE, PETSC_DETERMINE));
6399566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSizes(Prol, bs, col_bs));
6409566063dSJacob Faibussowitsch   PetscCall(MatSetType(Prol, mtype));
6419566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(Prol, col_bs, NULL));
6429566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(Prol, col_bs, NULL, col_bs, NULL));
6432e68589bSMark F. Adams 
6442e68589bSMark F. Adams   /* can get all points "removed" */
6459566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Prol, &kk, &ii));
6467f66b68fSMark Adams   if (!ii) {
64763a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(pc, "%s: No selected points on coarse grid\n", ((PetscObject)pc)->prefix));
6489566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Prol));
6490298fd71SBarry Smith     *a_P_out = NULL; /* out */
650849bee69SMark Adams     PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0));
6512e68589bSMark F. Adams     PetscFunctionReturn(0);
6522e68589bSMark F. Adams   }
65363a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(pc, "%s: New grid %" PetscInt_FMT " nodes\n", ((PetscObject)pc)->prefix, ii / col_bs));
6549566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRangeColumn(Prol, &myCrs0, &kk));
6550cbbd2e1SMark F. Adams 
65663a3b9bcSJacob 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);
657c8b0795cSMark F. Adams   myCrs0 = myCrs0 / col_bs;
65863a3b9bcSJacob 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);
6592e68589bSMark F. Adams 
6602e68589bSMark F. Adams   /* create global vector of data in 'data_w_ghost' */
661849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROLA], 0, 0, 0, 0));
662bae903cbSmarkadams4   if (size > 1) { /* get ghost null space data */
6632e68589bSMark F. Adams     PetscReal *tmp_gdata, *tmp_ldata, *tp2;
6649566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nloc, &tmp_ldata));
6654a2f8832SBarry Smith     for (jj = 0; jj < col_bs; jj++) {
666c8b0795cSMark F. Adams       for (kk = 0; kk < bs; kk++) {
667a2f3521dSMark F. Adams         PetscInt         ii, stride;
668c8b0795cSMark F. Adams         const PetscReal *tp = pc_gamg->data + jj * bs * nloc + kk;
6692fa5cd67SKarl Rupp         for (ii = 0; ii < nloc; ii++, tp += bs) tmp_ldata[ii] = *tp;
6702fa5cd67SKarl Rupp 
6719566063dSJacob Faibussowitsch         PetscCall(PCGAMGGetDataWithGhosts(Gmat, 1, tmp_ldata, &stride, &tmp_gdata));
672a2f3521dSMark F. Adams 
673bae903cbSmarkadams4         if (!jj && !kk) { /* now I know how many total nodes - allocate TODO: move below and do in one 'col_bs' call */
6749566063dSJacob Faibussowitsch           PetscCall(PetscMalloc1(stride * bs * col_bs, &data_w_ghost));
675a2f3521dSMark F. Adams           nbnodes = bs * stride;
6762e68589bSMark F. Adams         }
677a2f3521dSMark F. Adams         tp2 = data_w_ghost + jj * bs * stride + kk;
6782fa5cd67SKarl Rupp         for (ii = 0; ii < stride; ii++, tp2 += bs) *tp2 = tmp_gdata[ii];
6799566063dSJacob Faibussowitsch         PetscCall(PetscFree(tmp_gdata));
6802e68589bSMark F. Adams       }
6812e68589bSMark F. Adams     }
6829566063dSJacob Faibussowitsch     PetscCall(PetscFree(tmp_ldata));
683806fa848SBarry Smith   } else {
684c8b0795cSMark F. Adams     nbnodes      = bs * nloc;
685c8b0795cSMark F. Adams     data_w_ghost = (PetscReal *)pc_gamg->data;
6862e68589bSMark F. Adams   }
6872e68589bSMark F. Adams 
688bae903cbSmarkadams4   /* get 'flid_fgid' TODO - move up to get 'stride' and do get null space data above in one step (jj loop) */
689c5df96a5SBarry Smith   if (size > 1) {
6902e68589bSMark F. Adams     PetscReal *fid_glid_loc, *fiddata;
691a2f3521dSMark F. Adams     PetscInt   stride;
6922e68589bSMark F. Adams 
6939566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nloc, &fid_glid_loc));
6942e68589bSMark F. Adams     for (kk = 0; kk < nloc; kk++) fid_glid_loc[kk] = (PetscReal)(my0 + kk);
6959566063dSJacob Faibussowitsch     PetscCall(PCGAMGGetDataWithGhosts(Gmat, 1, fid_glid_loc, &stride, &fiddata));
696bae903cbSmarkadams4     PetscCall(PetscMalloc1(stride, &flid_fgid)); /* copy real data to in */
697a2f3521dSMark F. Adams     for (kk = 0; kk < stride; kk++) flid_fgid[kk] = (PetscInt)fiddata[kk];
6989566063dSJacob Faibussowitsch     PetscCall(PetscFree(fiddata));
699a2f3521dSMark F. Adams 
70063a3b9bcSJacob Faibussowitsch     PetscCheck(stride == nbnodes / bs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "stride %" PetscInt_FMT " != nbnodes %" PetscInt_FMT "/bs %" PetscInt_FMT, stride, nbnodes, bs);
7019566063dSJacob Faibussowitsch     PetscCall(PetscFree(fid_glid_loc));
702806fa848SBarry Smith   } else {
7039566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nloc, &flid_fgid));
7042e68589bSMark F. Adams     for (kk = 0; kk < nloc; kk++) flid_fgid[kk] = my0 + kk;
7052e68589bSMark F. Adams   }
706849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROLA], 0, 0, 0, 0));
707bae903cbSmarkadams4   /* get P0 */
708849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROLB], 0, 0, 0, 0));
709c8b0795cSMark F. Adams   {
7100298fd71SBarry Smith     PetscReal *data_out = NULL;
7119566063dSJacob Faibussowitsch     PetscCall(formProl0(agg_lists, bs, col_bs, myCrs0, nbnodes, data_w_ghost, flid_fgid, &data_out, Prol));
7129566063dSJacob Faibussowitsch     PetscCall(PetscFree(pc_gamg->data));
7132fa5cd67SKarl Rupp 
714c8b0795cSMark F. Adams     pc_gamg->data           = data_out;
7154a2f8832SBarry Smith     pc_gamg->data_cell_rows = col_bs;
7164a2f8832SBarry Smith     pc_gamg->data_sz        = col_bs * col_bs * nLocalSelected;
717c8b0795cSMark F. Adams   }
718849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROLB], 0, 0, 0, 0));
71948a46eb9SPierre Jolivet   if (size > 1) PetscCall(PetscFree(data_w_ghost));
7209566063dSJacob Faibussowitsch   PetscCall(PetscFree(flid_fgid));
7212e68589bSMark F. Adams 
722c8b0795cSMark F. Adams   *a_P_out = Prol; /* out */
723ed4e82eaSPeter Brune 
724849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0));
725c8b0795cSMark F. Adams   PetscFunctionReturn(0);
726c8b0795cSMark F. Adams }
727c8b0795cSMark F. Adams 
728c8b0795cSMark F. Adams /*
729fd1112cbSBarry Smith    PCGAMGOptProlongator_AGG
730c8b0795cSMark F. Adams 
731c8b0795cSMark F. Adams   Input Parameter:
732c8b0795cSMark F. Adams    . pc - this
733c8b0795cSMark F. Adams    . Amat - matrix on this fine level
734c8b0795cSMark F. Adams  In/Output Parameter:
7352a808120SBarry Smith    . a_P - prolongation operator to the next level
736c8b0795cSMark F. Adams */
7379371c9d4SSatish Balay static PetscErrorCode PCGAMGOptProlongator_AGG(PC pc, Mat Amat, Mat *a_P) {
738c8b0795cSMark F. Adams   PC_MG       *mg          = (PC_MG *)pc->data;
739c8b0795cSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
740c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
741c8b0795cSMark F. Adams   PetscInt     jj;
742c8b0795cSMark F. Adams   Mat          Prol = *a_P;
7433b4367a7SBarry Smith   MPI_Comm     comm;
7442a808120SBarry Smith   KSP          eksp;
7452a808120SBarry Smith   Vec          bb, xx;
7462a808120SBarry Smith   PC           epc;
7472a808120SBarry Smith   PetscReal    alpha, emax, emin;
748c8b0795cSMark F. Adams 
749c8b0795cSMark F. Adams   PetscFunctionBegin;
7509566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm));
751849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_OPT], 0, 0, 0, 0));
752c8b0795cSMark F. Adams 
753d5d25401SBarry Smith   /* compute maximum singular value of operator to be used in smoother */
7542a808120SBarry Smith   if (0 < pc_gamg_agg->nsmooths) {
75518c3aa7eSMark     /* get eigen estimates */
75618c3aa7eSMark     if (pc_gamg->emax > 0) {
75718c3aa7eSMark       emin = pc_gamg->emin;
75818c3aa7eSMark       emax = pc_gamg->emax;
75918c3aa7eSMark     } else {
7600ed2132dSStefano Zampini       const char *prefix;
7610ed2132dSStefano Zampini 
7629566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(Amat, &bb, NULL));
7639566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(Amat, &xx, NULL));
7649566063dSJacob Faibussowitsch       PetscCall(KSPSetNoisy_Private(bb));
765e696ed0bSMark F. Adams 
7669566063dSJacob Faibussowitsch       PetscCall(KSPCreate(comm, &eksp));
7679566063dSJacob Faibussowitsch       PetscCall(PCGetOptionsPrefix(pc, &prefix));
7689566063dSJacob Faibussowitsch       PetscCall(KSPSetOptionsPrefix(eksp, prefix));
7699566063dSJacob Faibussowitsch       PetscCall(KSPAppendOptionsPrefix(eksp, "pc_gamg_esteig_"));
77073f7197eSJed Brown       {
771b94d7dedSBarry Smith         PetscBool isset, sflg;
772b94d7dedSBarry Smith         PetscCall(MatIsSPDKnown(Amat, &isset, &sflg));
773b94d7dedSBarry Smith         if (isset && sflg) PetscCall(KSPSetType(eksp, KSPCG));
774d24ecf33SMark       }
7759566063dSJacob Faibussowitsch       PetscCall(KSPSetErrorIfNotConverged(eksp, pc->erroriffailure));
7769566063dSJacob Faibussowitsch       PetscCall(KSPSetNormType(eksp, KSP_NORM_NONE));
777c2436cacSMark F. Adams 
7789566063dSJacob Faibussowitsch       PetscCall(KSPSetInitialGuessNonzero(eksp, PETSC_FALSE));
7799566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(eksp, Amat, Amat));
7802e68589bSMark F. Adams 
7819566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(eksp, &epc));
7829566063dSJacob Faibussowitsch       PetscCall(PCSetType(epc, PCJACOBI)); /* smoother in smoothed agg. */
783c2436cacSMark F. Adams 
7849566063dSJacob 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
785074aec55SMark Adams 
7869566063dSJacob Faibussowitsch       PetscCall(KSPSetFromOptions(eksp));
7879566063dSJacob Faibussowitsch       PetscCall(KSPSetComputeSingularValues(eksp, PETSC_TRUE));
7889566063dSJacob Faibussowitsch       PetscCall(KSPSolve(eksp, bb, xx));
7899566063dSJacob Faibussowitsch       PetscCall(KSPCheckSolve(eksp, pc, xx));
7902e68589bSMark F. Adams 
7919566063dSJacob Faibussowitsch       PetscCall(KSPComputeExtremeSingularValues(eksp, &emax, &emin));
7929566063dSJacob Faibussowitsch       PetscCall(PetscInfo(pc, "%s: Smooth P0: max eigen=%e min=%e PC=%s\n", ((PetscObject)pc)->prefix, (double)emax, (double)emin, PCJACOBI));
7939566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&xx));
7949566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&bb));
7959566063dSJacob Faibussowitsch       PetscCall(KSPDestroy(&eksp));
7962e68589bSMark F. Adams     }
79718c3aa7eSMark     if (pc_gamg->use_sa_esteig) {
79818c3aa7eSMark       mg->min_eigen_DinvA[pc_gamg->current_level] = emin;
79918c3aa7eSMark       mg->max_eigen_DinvA[pc_gamg->current_level] = emax;
80063a3b9bcSJacob 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));
80118c3aa7eSMark     } else {
80218c3aa7eSMark       mg->min_eigen_DinvA[pc_gamg->current_level] = 0;
80318c3aa7eSMark       mg->max_eigen_DinvA[pc_gamg->current_level] = 0;
80418c3aa7eSMark     }
80518c3aa7eSMark   } else {
80618c3aa7eSMark     mg->min_eigen_DinvA[pc_gamg->current_level] = 0;
80718c3aa7eSMark     mg->max_eigen_DinvA[pc_gamg->current_level] = 0;
80818c3aa7eSMark   }
8092e68589bSMark F. Adams 
8102a808120SBarry Smith   /* smooth P0 */
8112a808120SBarry Smith   for (jj = 0; jj < pc_gamg_agg->nsmooths; jj++) {
8122a808120SBarry Smith     Mat tMat;
8132a808120SBarry Smith     Vec diag;
8142a808120SBarry Smith 
815849bee69SMark Adams     PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_OPTSM], 0, 0, 0, 0));
8162a808120SBarry Smith 
8172e68589bSMark F. Adams     /* smooth P1 := (I - omega/lam D^{-1}A)P0 */
8189566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2], 0, 0, 0, 0));
8199566063dSJacob Faibussowitsch     PetscCall(MatMatMult(Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat));
8209566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2], 0, 0, 0, 0));
8219566063dSJacob Faibussowitsch     PetscCall(MatProductClear(tMat));
8229566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(Amat, &diag, NULL));
8239566063dSJacob Faibussowitsch     PetscCall(MatGetDiagonal(Amat, diag)); /* effectively PCJACOBI */
8249566063dSJacob Faibussowitsch     PetscCall(VecReciprocal(diag));
8259566063dSJacob Faibussowitsch     PetscCall(MatDiagonalScale(tMat, diag, NULL));
8269566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&diag));
827b4da3a1bSStefano Zampini 
828d5d25401SBarry Smith     /* TODO: Set a PCFailedReason and exit the building of the AMG preconditioner */
82908401ef6SPierre Jolivet     PetscCheck(emax != 0.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_PLIB, "Computed maximum singular value as zero");
830d5d25401SBarry Smith     /* TODO: Document the 1.4 and don't hardwire it in this routine */
831b4ec6429SMark F. Adams     alpha = -1.4 / emax;
832b4da3a1bSStefano Zampini 
8339566063dSJacob Faibussowitsch     PetscCall(MatAYPX(tMat, alpha, Prol, SUBSET_NONZERO_PATTERN));
8349566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Prol));
8352e68589bSMark F. Adams     Prol = tMat;
836849bee69SMark Adams     PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_OPTSM], 0, 0, 0, 0));
8372e68589bSMark F. Adams   }
838849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_OPT], 0, 0, 0, 0));
839c8b0795cSMark F. Adams   *a_P = Prol;
840c8b0795cSMark F. Adams   PetscFunctionReturn(0);
8412e68589bSMark F. Adams }
8422e68589bSMark F. Adams 
843c8b0795cSMark F. Adams /*
844c8b0795cSMark F. Adams    PCCreateGAMG_AGG
8452e68589bSMark F. Adams 
846c8b0795cSMark F. Adams   Input Parameter:
847c8b0795cSMark F. Adams    . pc -
848c8b0795cSMark F. Adams */
8499371c9d4SSatish Balay PetscErrorCode PCCreateGAMG_AGG(PC pc) {
850c8b0795cSMark F. Adams   PC_MG       *mg      = (PC_MG *)pc->data;
851c8b0795cSMark F. Adams   PC_GAMG     *pc_gamg = (PC_GAMG *)mg->innerctx;
852c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg;
8532e68589bSMark F. Adams 
854c8b0795cSMark F. Adams   PetscFunctionBegin;
855c8b0795cSMark F. Adams   /* create sub context for SA */
8569566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(pc, &pc_gamg_agg));
857c8b0795cSMark F. Adams   pc_gamg->subctx = pc_gamg_agg;
858c8b0795cSMark F. Adams 
8591ab5ffc9SJed Brown   pc_gamg->ops->setfromoptions = PCSetFromOptions_GAMG_AGG;
8609b8ffb57SJed Brown   pc_gamg->ops->destroy        = PCDestroy_GAMG_AGG;
861c8b0795cSMark F. Adams   /* reset does not do anything; setup not virtual */
862c8b0795cSMark F. Adams 
863c8b0795cSMark F. Adams   /* set internal function pointers */
864fd1112cbSBarry Smith   pc_gamg->ops->graph             = PCGAMGGraph_AGG;
8651ab5ffc9SJed Brown   pc_gamg->ops->coarsen           = PCGAMGCoarsen_AGG;
8661ab5ffc9SJed Brown   pc_gamg->ops->prolongator       = PCGAMGProlongator_AGG;
867fd1112cbSBarry Smith   pc_gamg->ops->optprolongator    = PCGAMGOptProlongator_AGG;
8681ab5ffc9SJed Brown   pc_gamg->ops->createdefaultdata = PCSetData_AGG;
8695adeb434SBarry Smith   pc_gamg->ops->view              = PCView_GAMG_AGG;
870c8b0795cSMark F. Adams 
871bae903cbSmarkadams4   pc_gamg_agg->aggressive_coarsening_levels = 0;
872bae903cbSmarkadams4   pc_gamg_agg->symmetrize_graph             = PETSC_FALSE;
873cab9ed1eSBarry Smith   pc_gamg_agg->nsmooths                     = 1;
874cab9ed1eSBarry Smith 
8759566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNSmooths_C", PCGAMGSetNSmooths_AGG));
876bae903cbSmarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetSymmetrizeGraph_C", PCGAMGSetSymmetrizeGraph_AGG));
877bae903cbSmarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveLevels_C", PCGAMGSetAggressiveLevels_AGG));
8789566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", PCSetCoordinates_AGG));
8792e68589bSMark F. Adams   PetscFunctionReturn(0);
8802e68589bSMark F. Adams }
881