xref: /petsc/src/ksp/pc/impls/gamg/agg.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 /*@
172e68589bSMark F. Adams    PCGAMGSetNSmooths - Set number of smoothing steps (1 is typical)
182e68589bSMark F. Adams 
19d5d25401SBarry 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 
29db781477SPatrick Sanan .seealso: `()`
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 
52d5d25401SBarry Smith    Logically Collective on PC
53c8b0795cSMark F. Adams 
54c8b0795cSMark F. Adams    Input Parameters:
55cab9ed1eSBarry Smith +  pc - the preconditioner context
56a2b725a8SWilliam Gropp -  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 
63bae903cbSmarkadams4 .seealso: `PCGAMGSetAggressiveLevels()`
64c8b0795cSMark F. Adams @*/
659371c9d4SSatish Balay PetscErrorCode PCGAMGSetSymmetrizeGraph(PC pc, PetscBool n) {
66c8b0795cSMark F. Adams   PetscFunctionBegin;
67c8b0795cSMark F. Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
68d5d25401SBarry Smith   PetscValidLogicalCollectiveBool(pc, n, 2);
69bae903cbSmarkadams4   PetscTryMethod(pc, "PCGAMGSetSymmetrizeGraph_C", (PC, PetscBool), (pc, n));
70c8b0795cSMark F. Adams   PetscFunctionReturn(0);
71c8b0795cSMark F. Adams }
72c8b0795cSMark F. Adams 
739371c9d4SSatish Balay static PetscErrorCode PCGAMGSetSymmetrizeGraph_AGG(PC pc, PetscBool n) {
74c8b0795cSMark F. Adams   PC_MG       *mg          = (PC_MG *)pc->data;
75c8b0795cSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
76c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
77c8b0795cSMark F. Adams 
78c8b0795cSMark F. Adams   PetscFunctionBegin;
79bae903cbSmarkadams4   pc_gamg_agg->symmetrize_graph = n;
802e68589bSMark F. Adams   PetscFunctionReturn(0);
812e68589bSMark F. Adams }
822e68589bSMark F. Adams 
83ef4ad70eSMark F. Adams /*@
84bae903cbSmarkadams4    PCGAMGSetAggressiveLevels -  Aggressive coarsening on first n levels
85ef4ad70eSMark F. Adams 
86d5d25401SBarry Smith    Logically Collective on PC
87ef4ad70eSMark F. Adams 
88ef4ad70eSMark F. Adams    Input Parameters:
89cab9ed1eSBarry Smith +  pc - the preconditioner context
90d5d25401SBarry Smith -  n - 0, 1 or more
91ef4ad70eSMark F. Adams 
92ef4ad70eSMark F. Adams    Options Database Key:
93bae903cbSmarkadams4 .  -pc_gamg_aggressive_coarsening <n,default = 1> - Number of levels to square the graph on before aggregating it
94af3c827dSMark Adams 
95ef4ad70eSMark F. Adams    Level: intermediate
96ef4ad70eSMark F. Adams 
97bae903cbSmarkadams4 .seealso: `PCGAMGSetSymmetrizeGraph()`, `PCGAMGSetThreshold()`
98ef4ad70eSMark F. Adams @*/
999371c9d4SSatish Balay PetscErrorCode PCGAMGSetAggressiveLevels(PC pc, PetscInt n) {
100ef4ad70eSMark F. Adams   PetscFunctionBegin;
101ef4ad70eSMark F. Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
102d5d25401SBarry Smith   PetscValidLogicalCollectiveInt(pc, n, 2);
103bae903cbSmarkadams4   PetscTryMethod(pc, "PCGAMGSetAggressiveLevels_C", (PC, PetscInt), (pc, n));
104ef4ad70eSMark F. Adams   PetscFunctionReturn(0);
105ef4ad70eSMark F. Adams }
106ef4ad70eSMark F. Adams 
1079371c9d4SSatish Balay static PetscErrorCode PCGAMGSetAggressiveLevels_AGG(PC pc, PetscInt n) {
108ef4ad70eSMark F. Adams   PC_MG       *mg          = (PC_MG *)pc->data;
109ef4ad70eSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
110ef4ad70eSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
111ef4ad70eSMark F. Adams 
112ef4ad70eSMark F. Adams   PetscFunctionBegin;
113bae903cbSmarkadams4   pc_gamg_agg->aggressive_coarsening_levels = n;
114ef4ad70eSMark F. Adams   PetscFunctionReturn(0);
115ef4ad70eSMark F. Adams }
116ef4ad70eSMark F. Adams 
1179371c9d4SSatish Balay static PetscErrorCode PCSetFromOptions_GAMG_AGG(PC pc, PetscOptionItems *PetscOptionsObject) {
1182e68589bSMark F. Adams   PC_MG       *mg          = (PC_MG *)pc->data;
1192e68589bSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
120c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
1212e68589bSMark F. Adams 
1222e68589bSMark F. Adams   PetscFunctionBegin;
123d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "GAMG-AGG options");
1242e68589bSMark F. Adams   {
125bae903cbSmarkadams4     PetscBool flg;
1269566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-pc_gamg_agg_nsmooths", "smoothing steps for smoothed aggregation, usually 1", "PCGAMGSetNSmooths", pc_gamg_agg->nsmooths, &pc_gamg_agg->nsmooths, NULL));
127bae903cbSmarkadams4     PetscCall(PetscOptionsBool("-pc_gamg_symmetrize_graph", "Set for asymmetric matrices", "PCGAMGSetSymmetrizeGraph", pc_gamg_agg->symmetrize_graph, &pc_gamg_agg->symmetrize_graph, NULL));
128bae903cbSmarkadams4     pc_gamg_agg->aggressive_coarsening_levels = 1;
1299371c9d4SSatish Balay     PetscCall(
1309371c9d4SSatish 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));
131bae903cbSmarkadams4     if (!flg) {
132bae903cbSmarkadams4       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));
133bae903cbSmarkadams4     } else {
134bae903cbSmarkadams4       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));
135bae903cbSmarkadams4       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));
136bae903cbSmarkadams4     }
1372e68589bSMark F. Adams   }
138d0609cedSBarry Smith   PetscOptionsHeadEnd();
1392e68589bSMark F. Adams   PetscFunctionReturn(0);
1402e68589bSMark F. Adams }
1412e68589bSMark F. Adams 
1422e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
1439371c9d4SSatish Balay static PetscErrorCode PCDestroy_GAMG_AGG(PC pc) {
1442e68589bSMark F. Adams   PC_MG   *mg      = (PC_MG *)pc->data;
1452e68589bSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1462e68589bSMark F. Adams 
1472e68589bSMark F. Adams   PetscFunctionBegin;
1489566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg->subctx));
1492e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNSmooths_C", NULL));
150bae903cbSmarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetSymmetrizeGraph_C", NULL));
151bae903cbSmarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveLevels_C", NULL));
152bae903cbSmarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", NULL));
1532e68589bSMark F. Adams   PetscFunctionReturn(0);
1542e68589bSMark F. Adams }
1552e68589bSMark F. Adams 
1562e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
1572e68589bSMark F. Adams /*
1582e68589bSMark F. Adams    PCSetCoordinates_AGG
159302f38e8SMark F. Adams      - collective
1602e68589bSMark F. Adams 
1612e68589bSMark F. Adams    Input Parameter:
1622e68589bSMark F. Adams    . pc - the preconditioner context
163a2f3521dSMark F. Adams    . ndm - dimesion of data (used for dof/vertex for Stokes)
164302f38e8SMark F. Adams    . a_nloc - number of vertices local
165302f38e8SMark F. Adams    . coords - [a_nloc][ndm] - interleaved coordinate data: {x_0, y_0, z_0, x_1, y_1, ...}
1662e68589bSMark F. Adams */
1671e6b0712SBarry Smith 
1689371c9d4SSatish Balay static PetscErrorCode PCSetCoordinates_AGG(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords) {
1692e68589bSMark F. Adams   PC_MG   *mg      = (PC_MG *)pc->data;
1702e68589bSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
17169344418SMark F. Adams   PetscInt arrsz, kk, ii, jj, nloc, ndatarows, ndf;
172a2f3521dSMark F. Adams   Mat      mat = pc->pmat;
1732e68589bSMark F. Adams 
1742e68589bSMark F. Adams   PetscFunctionBegin;
175a2f3521dSMark F. Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
176a2f3521dSMark F. Adams   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
177302f38e8SMark F. Adams   nloc = a_nloc;
1782e68589bSMark F. Adams 
1792e68589bSMark F. Adams   /* SA: null space vectors */
1809566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(mat, &ndf));               /* this does not work for Stokes */
18169344418SMark F. Adams   if (coords && ndf == 1) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */
182a2f3521dSMark F. Adams   else if (coords) {
18363a3b9bcSJacob Faibussowitsch     PetscCheck(ndm <= ndf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "degrees of motion %" PetscInt_FMT " > block size %" PetscInt_FMT, ndm, ndf);
18469344418SMark F. Adams     pc_gamg->data_cell_cols = (ndm == 2 ? 3 : 6); /* displacement elasticity */
1859371c9d4SSatish Balay     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); }
186806fa848SBarry Smith   } else pc_gamg->data_cell_cols = ndf; /* no data, force SA with constant null space vectors */
18771959b99SBarry Smith   pc_gamg->data_cell_rows = ndatarows = ndf;
18863a3b9bcSJacob 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);
189c8b0795cSMark F. Adams   arrsz = nloc * pc_gamg->data_cell_rows * pc_gamg->data_cell_cols;
1902e68589bSMark F. Adams 
1917f66b68fSMark Adams   if (!pc_gamg->data || (pc_gamg->data_sz != arrsz)) {
1929566063dSJacob Faibussowitsch     PetscCall(PetscFree(pc_gamg->data));
1939566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(arrsz + 1, &pc_gamg->data));
1942e68589bSMark F. Adams   }
1952e68589bSMark F. Adams   /* copy data in - column oriented */
1962e68589bSMark F. Adams   for (kk = 0; kk < nloc; kk++) {
19769344418SMark F. Adams     const PetscInt M    = nloc * pc_gamg->data_cell_rows; /* stride into data */
19869344418SMark F. Adams     PetscReal     *data = &pc_gamg->data[kk * ndatarows]; /* start of cell */
199c8b0795cSMark F. Adams     if (pc_gamg->data_cell_cols == 1) *data = 1.0;
2002e68589bSMark F. Adams     else {
20169344418SMark F. Adams       /* translational modes */
2022fa5cd67SKarl Rupp       for (ii = 0; ii < ndatarows; ii++) {
2032fa5cd67SKarl Rupp         for (jj = 0; jj < ndatarows; jj++) {
20469344418SMark F. Adams           if (ii == jj) data[ii * M + jj] = 1.0;
2052e68589bSMark F. Adams           else data[ii * M + jj] = 0.0;
2062fa5cd67SKarl Rupp         }
2072fa5cd67SKarl Rupp       }
20869344418SMark F. Adams 
20969344418SMark F. Adams       /* rotational modes */
2102e68589bSMark F. Adams       if (coords) {
21169344418SMark F. Adams         if (ndm == 2) {
2122e68589bSMark F. Adams           data += 2 * M;
2132e68589bSMark F. Adams           data[0] = -coords[2 * kk + 1];
2142e68589bSMark F. Adams           data[1] = coords[2 * kk];
215806fa848SBarry Smith         } else {
2162e68589bSMark F. Adams           data += 3 * M;
2179371c9d4SSatish Balay           data[0]         = 0.0;
2189371c9d4SSatish Balay           data[M + 0]     = coords[3 * kk + 2];
2199371c9d4SSatish Balay           data[2 * M + 0] = -coords[3 * kk + 1];
2209371c9d4SSatish Balay           data[1]         = -coords[3 * kk + 2];
2219371c9d4SSatish Balay           data[M + 1]     = 0.0;
2229371c9d4SSatish Balay           data[2 * M + 1] = coords[3 * kk];
2239371c9d4SSatish Balay           data[2]         = coords[3 * kk + 1];
2249371c9d4SSatish Balay           data[M + 2]     = -coords[3 * kk];
2259371c9d4SSatish Balay           data[2 * M + 2] = 0.0;
2262e68589bSMark F. Adams         }
2272e68589bSMark F. Adams       }
2282e68589bSMark F. Adams     }
2292e68589bSMark F. Adams   }
2302e68589bSMark F. Adams   pc_gamg->data_sz = arrsz;
2312e68589bSMark F. Adams   PetscFunctionReturn(0);
2322e68589bSMark F. Adams }
2332e68589bSMark F. Adams 
2342e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
2352e68589bSMark F. Adams /*
236a2f3521dSMark F. Adams    PCSetData_AGG - called if data is not set with PCSetCoordinates.
237a2f3521dSMark F. Adams       Looks in Mat for near null space.
238a2f3521dSMark F. Adams       Does not work for Stokes
2392e68589bSMark F. Adams 
2402e68589bSMark F. Adams   Input Parameter:
2412e68589bSMark F. Adams    . pc -
242a2f3521dSMark F. Adams    . a_A - matrix to get (near) null space out of.
2432e68589bSMark F. Adams */
2449371c9d4SSatish Balay static PetscErrorCode PCSetData_AGG(PC pc, Mat a_A) {
245b8cd405aSMark F. Adams   PC_MG       *mg      = (PC_MG *)pc->data;
246b8cd405aSMark F. Adams   PC_GAMG     *pc_gamg = (PC_GAMG *)mg->innerctx;
247b8cd405aSMark F. Adams   MatNullSpace mnull;
24866f2ef4bSMark Adams 
249b817416eSBarry Smith   PetscFunctionBegin;
2509566063dSJacob Faibussowitsch   PetscCall(MatGetNearNullSpace(a_A, &mnull));
251b8cd405aSMark F. Adams   if (!mnull) {
25266f2ef4bSMark Adams     DM dm;
2539566063dSJacob Faibussowitsch     PetscCall(PCGetDM(pc, &dm));
254*48a46eb9SPierre Jolivet     if (!dm) PetscCall(MatGetDM(a_A, &dm));
25566f2ef4bSMark Adams     if (dm) {
25666f2ef4bSMark Adams       PetscObject deformation;
257b0d30dd6SMatthew G. Knepley       PetscInt    Nf;
258b0d30dd6SMatthew G. Knepley 
2599566063dSJacob Faibussowitsch       PetscCall(DMGetNumFields(dm, &Nf));
260b0d30dd6SMatthew G. Knepley       if (Nf) {
2619566063dSJacob Faibussowitsch         PetscCall(DMGetField(dm, 0, NULL, &deformation));
2629566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery((PetscObject)deformation, "nearnullspace", (PetscObject *)&mnull));
263*48a46eb9SPierre Jolivet         if (!mnull) PetscCall(PetscObjectQuery((PetscObject)deformation, "nullspace", (PetscObject *)&mnull));
26466f2ef4bSMark Adams       }
26566f2ef4bSMark Adams     }
266b0d30dd6SMatthew G. Knepley   }
26766f2ef4bSMark Adams 
26866f2ef4bSMark Adams   if (!mnull) {
269a2f3521dSMark F. Adams     PetscInt bs, NN, MM;
2709566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSize(a_A, &bs));
2719566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(a_A, &MM, &NN));
27263a3b9bcSJacob Faibussowitsch     PetscCheck(MM % bs == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MM %" PetscInt_FMT " must be divisible by bs %" PetscInt_FMT, MM, bs);
2739566063dSJacob Faibussowitsch     PetscCall(PCSetCoordinates_AGG(pc, bs, MM / bs, NULL));
274806fa848SBarry Smith   } else {
275b8cd405aSMark F. Adams     PetscReal         *nullvec;
276b8cd405aSMark F. Adams     PetscBool          has_const;
277b8cd405aSMark F. Adams     PetscInt           i, j, mlocal, nvec, bs;
278d5d25401SBarry Smith     const Vec         *vecs;
279d5d25401SBarry Smith     const PetscScalar *v;
280b817416eSBarry Smith 
2819566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(a_A, &mlocal, NULL));
2829566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs));
283a0dea87dSMark Adams     pc_gamg->data_sz = (nvec + !!has_const) * mlocal;
2849566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1((nvec + !!has_const) * mlocal, &nullvec));
2859371c9d4SSatish Balay     if (has_const)
2869371c9d4SSatish Balay       for (i = 0; i < mlocal; i++) nullvec[i] = 1.0;
287b8cd405aSMark F. Adams     for (i = 0; i < nvec; i++) {
2889566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(vecs[i], &v));
289b8cd405aSMark F. Adams       for (j = 0; j < mlocal; j++) nullvec[(i + !!has_const) * mlocal + j] = PetscRealPart(v[j]);
2909566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(vecs[i], &v));
291b8cd405aSMark F. Adams     }
292b8cd405aSMark F. Adams     pc_gamg->data           = nullvec;
293b8cd405aSMark F. Adams     pc_gamg->data_cell_cols = (nvec + !!has_const);
2949566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSize(a_A, &bs));
295b8cd405aSMark F. Adams     pc_gamg->data_cell_rows = bs;
296b8cd405aSMark F. Adams   }
2972e68589bSMark F. Adams   PetscFunctionReturn(0);
2982e68589bSMark F. Adams }
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
316bae903cbSmarkadams4 
3172e68589bSMark F. Adams */
3189371c9d4SSatish 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) {
319bd026e97SJed Brown   PetscInt        Istart, my0, Iend, nloc, clid, flid = 0, aggID, kk, jj, ii, mm, nSelected, minsz, nghosts, out_data_stride;
3203b4367a7SBarry Smith   MPI_Comm        comm;
3212e68589bSMark F. Adams   PetscReal      *out_data;
322539c167fSBarry Smith   PetscCDIntNd   *pos;
3231943db53SBarry Smith   PCGAMGHashTable fgid_flid;
3240cbbd2e1SMark F. Adams 
3252e68589bSMark F. Adams   PetscFunctionBegin;
3269566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)a_Prol, &comm));
3279566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(a_Prol, &Istart, &Iend));
3289371c9d4SSatish Balay   nloc = (Iend - Istart) / bs;
3299371c9d4SSatish Balay   my0  = Istart / bs;
33063a3b9bcSJacob 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);
3310cbbd2e1SMark F. Adams   Iend /= bs;
3320cbbd2e1SMark F. Adams   nghosts = data_stride / bs - nloc;
3332e68589bSMark F. Adams 
3349566063dSJacob Faibussowitsch   PetscCall(PCGAMGHashTableCreate(2 * nghosts + 1, &fgid_flid));
335*48a46eb9SPierre Jolivet   for (kk = 0; kk < nghosts; kk++) PetscCall(PCGAMGHashTableAdd(&fgid_flid, flid_fgid[nloc + kk], nloc + kk));
3362e68589bSMark F. Adams 
3370cbbd2e1SMark F. Adams   /* count selected -- same as number of cols of P */
3380cbbd2e1SMark F. Adams   for (nSelected = mm = 0; mm < nloc; mm++) {
339e78576d6SMark F. Adams     PetscBool ise;
3409566063dSJacob Faibussowitsch     PetscCall(PetscCDEmptyAt(agg_llists, mm, &ise));
341e78576d6SMark F. Adams     if (!ise) nSelected++;
3420cbbd2e1SMark F. Adams   }
3439566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRangeColumn(a_Prol, &ii, &jj));
34463a3b9bcSJacob Faibussowitsch   PetscCheck((ii / nSAvec) == my0crs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ii %" PetscInt_FMT " /nSAvec %" PetscInt_FMT "  != my0crs %" PetscInt_FMT, ii, nSAvec, my0crs);
34563a3b9bcSJacob 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);
3460cbbd2e1SMark F. Adams 
3472e68589bSMark F. Adams   /* aloc space for coarse point data (output) */
3480cbbd2e1SMark F. Adams   out_data_stride = nSelected * nSAvec;
3492fa5cd67SKarl Rupp 
3509566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(out_data_stride * nSAvec, &out_data));
35133812677SSatish Balay   for (ii = 0; ii < out_data_stride * nSAvec; ii++) out_data[ii] = PETSC_MAX_REAL;
3520cbbd2e1SMark F. Adams   *a_data_out = out_data; /* output - stride nSelected*nSAvec */
3532e68589bSMark F. Adams 
3542e68589bSMark F. Adams   /* find points and set prolongation */
355c8b0795cSMark F. Adams   minsz = 100;
3560cbbd2e1SMark F. Adams   for (mm = clid = 0; mm < nloc; mm++) {
3579566063dSJacob Faibussowitsch     PetscCall(PetscCDSizeAt(agg_llists, mm, &jj));
358e78576d6SMark F. Adams     if (jj > 0) {
3590cbbd2e1SMark F. Adams       const PetscInt lid = mm, cgid = my0crs + clid;
3600cbbd2e1SMark F. Adams       PetscInt       cids[100]; /* max bs */
3610cbbd2e1SMark F. Adams       PetscBLASInt   asz = jj, M = asz * bs, N = nSAvec, INFO;
3622e68589bSMark F. Adams       PetscBLASInt   Mdata = M + ((N - M > 0) ? N - M : 0), LDA = Mdata, LWORK = N * bs;
3632e68589bSMark F. Adams       PetscScalar   *qqc, *qqr, *TAU, *WORK;
3642e68589bSMark F. Adams       PetscInt      *fids;
36565d7b583SSatish Balay       PetscReal     *data;
366b817416eSBarry Smith 
3670cbbd2e1SMark F. Adams       /* count agg */
3680cbbd2e1SMark F. Adams       if (asz < minsz) minsz = asz;
3690cbbd2e1SMark F. Adams 
3700cbbd2e1SMark F. Adams       /* get block */
3719566063dSJacob Faibussowitsch       PetscCall(PetscMalloc5(Mdata * N, &qqc, M * N, &qqr, N, &TAU, LWORK, &WORK, M, &fids));
3722e68589bSMark F. Adams 
3732e68589bSMark F. Adams       aggID = 0;
3749566063dSJacob Faibussowitsch       PetscCall(PetscCDGetHeadPos(agg_llists, lid, &pos));
375e78576d6SMark F. Adams       while (pos) {
376e78576d6SMark F. Adams         PetscInt gid1;
3779566063dSJacob Faibussowitsch         PetscCall(PetscCDIntNdGetID(pos, &gid1));
3789566063dSJacob Faibussowitsch         PetscCall(PetscCDGetNextPos(agg_llists, lid, &pos));
379e78576d6SMark F. Adams 
3800cbbd2e1SMark F. Adams         if (gid1 >= my0 && gid1 < Iend) flid = gid1 - my0;
3810cbbd2e1SMark F. Adams         else {
3829566063dSJacob Faibussowitsch           PetscCall(PCGAMGHashTableFind(&fgid_flid, gid1, &flid));
38308401ef6SPierre Jolivet           PetscCheck(flid >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cannot find gid1 in table");
3840cbbd2e1SMark F. Adams         }
3852e68589bSMark F. Adams         /* copy in B_i matrix - column oriented */
38665d7b583SSatish Balay         data = &data_in[flid * bs];
3873d3eaba7SBarry Smith         for (ii = 0; ii < bs; ii++) {
3882e68589bSMark F. Adams           for (jj = 0; jj < N; jj++) {
3890cbbd2e1SMark F. Adams             PetscReal d                       = data[jj * data_stride + ii];
3900cbbd2e1SMark F. Adams             qqc[jj * Mdata + aggID * bs + ii] = d;
3912e68589bSMark F. Adams           }
3922e68589bSMark F. Adams         }
3932e68589bSMark F. Adams         /* set fine IDs */
3942e68589bSMark F. Adams         for (kk = 0; kk < bs; kk++) fids[aggID * bs + kk] = flid_fgid[flid] * bs + kk;
3952e68589bSMark F. Adams         aggID++;
3960cbbd2e1SMark F. Adams       }
3972e68589bSMark F. Adams 
3982e68589bSMark F. Adams       /* pad with zeros */
3992e68589bSMark F. Adams       for (ii = asz * bs; ii < Mdata; ii++) {
4009371c9d4SSatish Balay         for (jj = 0; jj < N; jj++, kk++) { qqc[jj * Mdata + ii] = .0; }
4012e68589bSMark F. Adams       }
4022e68589bSMark F. Adams 
4032e68589bSMark F. Adams       /* QR */
4049566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
405792fecdfSBarry Smith       PetscCallBLAS("LAPACKgeqrf", LAPACKgeqrf_(&Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
4069566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPop());
40708401ef6SPierre Jolivet       PetscCheck(INFO == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xGEQRF error");
4082e68589bSMark F. Adams       /* get R - column oriented - output B_{i+1} */
4092e68589bSMark F. Adams       {
4102e68589bSMark F. Adams         PetscReal *data = &out_data[clid * nSAvec];
4112e68589bSMark F. Adams         for (jj = 0; jj < nSAvec; jj++) {
4122e68589bSMark F. Adams           for (ii = 0; ii < nSAvec; ii++) {
41308401ef6SPierre 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);
4140cbbd2e1SMark F. Adams             if (ii <= jj) data[jj * out_data_stride + ii] = PetscRealPart(qqc[jj * Mdata + ii]);
4150cbbd2e1SMark F. Adams             else data[jj * out_data_stride + ii] = 0.;
4162e68589bSMark F. Adams           }
4172e68589bSMark F. Adams         }
4182e68589bSMark F. Adams       }
4192e68589bSMark F. Adams 
4202e68589bSMark F. Adams       /* get Q - row oriented */
421792fecdfSBarry Smith       PetscCallBLAS("LAPACKorgqr", LAPACKorgqr_(&Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
42263a3b9bcSJacob Faibussowitsch       PetscCheck(INFO == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xORGQR error arg %" PetscBLASInt_FMT, -INFO);
4232e68589bSMark F. Adams 
4242e68589bSMark F. Adams       for (ii = 0; ii < M; ii++) {
4259371c9d4SSatish Balay         for (jj = 0; jj < N; jj++) { qqr[N * ii + jj] = qqc[jj * Mdata + ii]; }
4262e68589bSMark F. Adams       }
4272e68589bSMark F. Adams 
4282e68589bSMark F. Adams       /* add diagonal block of P0 */
4299371c9d4SSatish Balay       for (kk = 0; kk < N; kk++) { cids[kk] = N * cgid + kk; /* global col IDs in P0 */ }
4309566063dSJacob Faibussowitsch       PetscCall(MatSetValues(a_Prol, M, fids, N, cids, qqr, INSERT_VALUES));
4319566063dSJacob Faibussowitsch       PetscCall(PetscFree5(qqc, qqr, TAU, WORK, fids));
432b43b03e9SMark F. Adams       clid++;
4330cbbd2e1SMark F. Adams     } /* coarse agg */
4340cbbd2e1SMark F. Adams   }   /* for all fine nodes */
4359566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(a_Prol, MAT_FINAL_ASSEMBLY));
4369566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(a_Prol, MAT_FINAL_ASSEMBLY));
4379566063dSJacob Faibussowitsch   PetscCall(PCGAMGHashTableDestroy(&fgid_flid));
4382e68589bSMark F. Adams   PetscFunctionReturn(0);
4392e68589bSMark F. Adams }
4402e68589bSMark F. Adams 
4419371c9d4SSatish Balay static PetscErrorCode PCView_GAMG_AGG(PC pc, PetscViewer viewer) {
4425adeb434SBarry Smith   PC_MG       *mg          = (PC_MG *)pc->data;
4435adeb434SBarry Smith   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
4445adeb434SBarry Smith   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
4455adeb434SBarry Smith 
4465adeb434SBarry Smith   PetscFunctionBegin;
4479566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "      AGG specific options\n"));
448bae903cbSmarkadams4   PetscCall(PetscViewerASCIIPrintf(viewer, "        Symmetric graph %s\n", pc_gamg_agg->symmetrize_graph ? "true" : "false"));
449bae903cbSmarkadams4   PetscCall(PetscViewerASCIIPrintf(viewer, "        Number of levels to square graph %d\n", (int)pc_gamg_agg->aggressive_coarsening_levels));
450bae903cbSmarkadams4   PetscCall(PetscViewerASCIIPrintf(viewer, "        Number smoothing steps %d\n", (int)pc_gamg_agg->nsmooths));
4515adeb434SBarry Smith   PetscFunctionReturn(0);
4525adeb434SBarry Smith }
4535adeb434SBarry Smith 
4542e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
4552e68589bSMark F. Adams /*
456fd1112cbSBarry Smith    PCGAMGGraph_AGG
4572e68589bSMark F. Adams 
4582e68589bSMark F. Adams   Input Parameter:
4592e68589bSMark F. Adams    . pc - this
4602e68589bSMark F. Adams    . Amat - matrix on this fine level
461bae903cbSmarkadams4 
4622e68589bSMark F. Adams   Output Parameter:
463c8b0795cSMark F. Adams    . a_Gmat -
4642e68589bSMark F. Adams */
4659371c9d4SSatish Balay static PetscErrorCode PCGAMGGraph_AGG(PC pc, Mat Amat, Mat *a_Gmat) {
466c8b0795cSMark F. Adams   PC_MG                    *mg          = (PC_MG *)pc->data;
467c8b0795cSMark F. Adams   PC_GAMG                  *pc_gamg     = (PC_GAMG *)mg->innerctx;
468c1eae691SMark Adams   const PetscReal           vfilter     = pc_gamg->threshold[pc_gamg->current_level];
469c8b0795cSMark F. Adams   PC_GAMG_AGG              *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
47072833a62Smarkadams4   Mat                       Gmat, F = NULL;
4713b4367a7SBarry Smith   MPI_Comm                  comm;
47267744fedSMark F. Adams   PetscBool /* set,flg , */ symm;
473c8b0795cSMark F. Adams 
474c8b0795cSMark F. Adams   PetscFunctionBegin;
4759566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm));
476c8b0795cSMark F. Adams 
4779566063dSJacob Faibussowitsch   /* PetscCall(MatIsSymmetricKnown(Amat, &set, &flg)); || !(set && flg) -- this causes lot of symm calls */
478bae903cbSmarkadams4   symm = (PetscBool)(pc_gamg_agg->symmetrize_graph); /* && !pc_gamg_agg->square_graph; */
4790cbbd2e1SMark F. Adams 
480849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_GRAPH], 0, 0, 0, 0));
48172833a62Smarkadams4   PetscCall(MatCreateGraph(Amat, symm, PETSC_TRUE, &Gmat));
482849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_GRAPH], 0, 0, 0, 0));
483849bee69SMark Adams 
484849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_FILTER], 0, 0, 0, 0));
48572833a62Smarkadams4   PetscCall(MatFilter(Gmat, vfilter, &F));
48672833a62Smarkadams4   if (F) {
48772833a62Smarkadams4     PetscCall(MatDestroy(&Gmat));
48872833a62Smarkadams4     Gmat = F;
48972833a62Smarkadams4   }
490849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_FILTER], 0, 0, 0, 0));
491849bee69SMark Adams 
492e0940f08SMark F. Adams   *a_Gmat = Gmat;
493849bee69SMark Adams 
494c8b0795cSMark F. Adams   PetscFunctionReturn(0);
495c8b0795cSMark F. Adams }
496c8b0795cSMark F. Adams 
497c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
498c8b0795cSMark F. Adams /*
499bae903cbSmarkadams4    PCGAMGCoarsen_AGG - supports squaring the graph (deprecated) and new graph for
500bae903cbSmarkadams4      communication of QR data used with HEM and MISk coarsening
501c8b0795cSMark F. Adams 
502c8b0795cSMark F. Adams   Input Parameter:
503e0940f08SMark F. Adams    . a_pc - this
504bae903cbSmarkadams4 
505e0940f08SMark F. Adams   Input/Output Parameter:
506bae903cbSmarkadams4    . a_Gmat1 - graph to coarsen (in), graph off processor edges for QR gather scatter (out)
507bae903cbSmarkadams4 
508c8b0795cSMark F. Adams   Output Parameter:
5090cbbd2e1SMark F. Adams    . agg_lists - list of aggregates
510bae903cbSmarkadams4 
511c8b0795cSMark F. Adams */
5129371c9d4SSatish Balay static PetscErrorCode PCGAMGCoarsen_AGG(PC a_pc, Mat *a_Gmat1, PetscCoarsenData **agg_lists) {
513e0940f08SMark F. Adams   PC_MG       *mg          = (PC_MG *)a_pc->data;
514c8b0795cSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
515c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
516bae903cbSmarkadams4   Mat          mat, Gmat1 = *a_Gmat1; /* aggressive graph */
5170cbbd2e1SMark F. Adams   IS           perm;
518bae903cbSmarkadams4   PetscInt     Istart, Iend, Ii, nloc, bs, nn;
519bae903cbSmarkadams4   PetscInt    *permute, *degree;
520c8b0795cSMark F. Adams   PetscBool   *bIndexSet;
521b43b03e9SMark F. Adams   MatCoarsen   crs;
5223b4367a7SBarry Smith   MPI_Comm     comm;
523aea53286SMark Adams   PetscReal    hashfact;
524e2c00dcbSBarry Smith   PetscInt     iSwapIndex;
5253b50add6SMark Adams   PetscRandom  random;
526c8b0795cSMark F. Adams 
527c8b0795cSMark F. Adams   PetscFunctionBegin;
528849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0));
5299566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Gmat1, &comm));
530bae903cbSmarkadams4   PetscCall(MatGetLocalSize(Gmat1, &nn, NULL));
5319566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(Gmat1, &bs));
53263a3b9bcSJacob Faibussowitsch   PetscCheck(bs == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "bs %" PetscInt_FMT " must be 1", bs);
533bae903cbSmarkadams4   nloc = nn / bs;
534c8b0795cSMark F. Adams 
5355cfd4bd9SMark Adams   /* get MIS aggs - randomize */
536bae903cbSmarkadams4   PetscCall(PetscMalloc2(nloc, &permute, nloc, &degree));
5379566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(nloc, &bIndexSet));
5386e09b0e3SStefano Zampini   for (Ii = 0; Ii < nloc; Ii++) permute[Ii] = Ii;
5399566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &random));
5409566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Gmat1, &Istart, &Iend));
541c8b0795cSMark F. Adams   for (Ii = 0; Ii < nloc; Ii++) {
542bae903cbSmarkadams4     PetscInt nc;
543bae903cbSmarkadams4     PetscCall(MatGetRow(Gmat1, Istart + Ii, &nc, NULL, NULL));
544bae903cbSmarkadams4     degree[Ii] = nc;
545bae903cbSmarkadams4     PetscCall(MatRestoreRow(Gmat1, Istart + Ii, &nc, NULL, NULL));
546bae903cbSmarkadams4   }
547bae903cbSmarkadams4   for (Ii = 0; Ii < nloc; Ii++) {
5489566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValueReal(random, &hashfact));
549aea53286SMark Adams     iSwapIndex = (PetscInt)(hashfact * nloc) % nloc;
550c8b0795cSMark F. Adams     if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) {
551c8b0795cSMark F. Adams       PetscInt iTemp        = permute[iSwapIndex];
552c8b0795cSMark F. Adams       permute[iSwapIndex]   = permute[Ii];
553c8b0795cSMark F. Adams       permute[Ii]           = iTemp;
554bae903cbSmarkadams4       iTemp                 = degree[iSwapIndex];
555bae903cbSmarkadams4       degree[iSwapIndex]    = degree[Ii];
556bae903cbSmarkadams4       degree[Ii]            = iTemp;
557c8b0795cSMark F. Adams       bIndexSet[iSwapIndex] = PETSC_TRUE;
558c8b0795cSMark F. Adams     }
559c8b0795cSMark F. Adams   }
560bae903cbSmarkadams4   // create minimum degree ordering
561bae903cbSmarkadams4   PetscCall(PetscSortIntWithArray(nloc, degree, permute));
562bae903cbSmarkadams4 
5639566063dSJacob Faibussowitsch   PetscCall(PetscFree(bIndexSet));
5649566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&random));
5659566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm));
566849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_MIS], 0, 0, 0, 0));
5679566063dSJacob Faibussowitsch   PetscCall(MatCoarsenCreate(comm, &crs));
5689566063dSJacob Faibussowitsch   PetscCall(MatCoarsenSetFromOptions(crs));
5699566063dSJacob Faibussowitsch   PetscCall(MatCoarsenSetGreedyOrdering(crs, perm));
570bae903cbSmarkadams4   PetscCall(MatCoarsenSetAdjacency(crs, Gmat1));
5719566063dSJacob Faibussowitsch   PetscCall(MatCoarsenSetStrictAggs(crs, PETSC_TRUE));
572bae903cbSmarkadams4   if (pc_gamg->current_level < pc_gamg_agg->aggressive_coarsening_levels) PetscCall(MatCoarsenMISKSetDistance(crs, 2)); // hardwire to MIS-2
573bae903cbSmarkadams4   else PetscCall(MatCoarsenMISKSetDistance(crs, 1));                                                                    // MIS
5749566063dSJacob Faibussowitsch   PetscCall(MatCoarsenApply(crs));
575bae903cbSmarkadams4   PetscCall(MatCoarsenViewFromOptions(crs, NULL, "-mat_coarsen_view"));
5769566063dSJacob Faibussowitsch   PetscCall(MatCoarsenGetData(crs, agg_lists)); /* output */
5779566063dSJacob Faibussowitsch   PetscCall(MatCoarsenDestroy(&crs));
578b43b03e9SMark F. Adams 
5799566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
580bae903cbSmarkadams4   PetscCall(PetscFree2(permute, degree));
581849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_MIS], 0, 0, 0, 0));
5829f3f12c8SMark F. Adams 
583bae903cbSmarkadams4   {
584bae903cbSmarkadams4     PetscCoarsenData *llist = *agg_lists;
5859ab59c8bSMark Adams     /* see if we have a matrix that takes precedence (returned from MatCoarsenApply) */
5869566063dSJacob Faibussowitsch     PetscCall(PetscCDGetMat(llist, &mat));
5870cbbd2e1SMark F. Adams     if (mat) {
5889566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&Gmat1));
5890cbbd2e1SMark F. Adams       *a_Gmat1 = mat; /* output */
5900cbbd2e1SMark F. Adams     }
5910cbbd2e1SMark F. Adams   }
592849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0));
593c8b0795cSMark F. Adams   PetscFunctionReturn(0);
594c8b0795cSMark F. Adams }
595c8b0795cSMark F. Adams 
596c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
597c8b0795cSMark F. Adams /*
5980cbbd2e1SMark F. Adams  PCGAMGProlongator_AGG
599c8b0795cSMark F. Adams 
600c8b0795cSMark F. Adams  Input Parameter:
601c8b0795cSMark F. Adams  . pc - this
602c8b0795cSMark F. Adams  . Amat - matrix on this fine level
603c8b0795cSMark F. Adams  . Graph - used to get ghost data for nodes in
6040cbbd2e1SMark F. Adams  . agg_lists - list of aggregates
605c8b0795cSMark F. Adams  Output Parameter:
606c8b0795cSMark F. Adams  . a_P_out - prolongation operator to the next level
607c8b0795cSMark F. Adams  */
6089371c9d4SSatish Balay static PetscErrorCode PCGAMGProlongator_AGG(PC pc, Mat Amat, Mat Gmat, PetscCoarsenData *agg_lists, Mat *a_P_out) {
6092e68589bSMark F. Adams   PC_MG         *mg      = (PC_MG *)pc->data;
6102e68589bSMark F. Adams   PC_GAMG       *pc_gamg = (PC_GAMG *)mg->innerctx;
6114a2f8832SBarry Smith   const PetscInt col_bs  = pc_gamg->data_cell_cols;
612c8b0795cSMark F. Adams   PetscInt       Istart, Iend, nloc, ii, jj, kk, my0, nLocalSelected, bs;
613c8b0795cSMark F. Adams   Mat            Prol;
614d5d25401SBarry Smith   PetscMPIInt    size;
6153b4367a7SBarry Smith   MPI_Comm       comm;
616c8b0795cSMark F. Adams   PetscReal     *data_w_ghost;
617c8b0795cSMark F. Adams   PetscInt       myCrs0, nbnodes = 0, *flid_fgid;
618d9558ea9SBarry Smith   MatType        mtype;
6192e68589bSMark F. Adams 
6202e68589bSMark F. Adams   PetscFunctionBegin;
6219566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm));
62208401ef6SPierre Jolivet   PetscCheck(col_bs >= 1, comm, PETSC_ERR_PLIB, "Column bs cannot be less than 1");
623849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0));
6249566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
6259566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Amat, &Istart, &Iend));
6269566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(Amat, &bs));
6279371c9d4SSatish Balay   nloc = (Iend - Istart) / bs;
6289371c9d4SSatish Balay   my0  = Istart / bs;
62963a3b9bcSJacob 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);
6302e68589bSMark F. Adams 
6312e68589bSMark F. Adams   /* get 'nLocalSelected' */
6320cbbd2e1SMark F. Adams   for (ii = 0, nLocalSelected = 0; ii < nloc; ii++) {
633e78576d6SMark F. Adams     PetscBool ise;
6340cbbd2e1SMark F. Adams     /* filter out singletons 0 or 1? */
6359566063dSJacob Faibussowitsch     PetscCall(PetscCDEmptyAt(agg_lists, ii, &ise));
636e78576d6SMark F. Adams     if (!ise) nLocalSelected++;
6372e68589bSMark F. Adams   }
6382e68589bSMark F. Adams 
6392e68589bSMark F. Adams   /* create prolongator, create P matrix */
6409566063dSJacob Faibussowitsch   PetscCall(MatGetType(Amat, &mtype));
6419566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &Prol));
6429566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Prol, nloc * bs, nLocalSelected * col_bs, PETSC_DETERMINE, PETSC_DETERMINE));
6439566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSizes(Prol, bs, col_bs));
6449566063dSJacob Faibussowitsch   PetscCall(MatSetType(Prol, mtype));
6459566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(Prol, col_bs, NULL));
6469566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(Prol, col_bs, NULL, col_bs, NULL));
6472e68589bSMark F. Adams 
6482e68589bSMark F. Adams   /* can get all points "removed" */
6499566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Prol, &kk, &ii));
6507f66b68fSMark Adams   if (!ii) {
65163a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(pc, "%s: No selected points on coarse grid\n", ((PetscObject)pc)->prefix));
6529566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Prol));
6530298fd71SBarry Smith     *a_P_out = NULL; /* out */
654849bee69SMark Adams     PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0));
6552e68589bSMark F. Adams     PetscFunctionReturn(0);
6562e68589bSMark F. Adams   }
65763a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(pc, "%s: New grid %" PetscInt_FMT " nodes\n", ((PetscObject)pc)->prefix, ii / col_bs));
6589566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRangeColumn(Prol, &myCrs0, &kk));
6590cbbd2e1SMark F. Adams 
66063a3b9bcSJacob 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);
661c8b0795cSMark F. Adams   myCrs0 = myCrs0 / col_bs;
66263a3b9bcSJacob 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);
6632e68589bSMark F. Adams 
6642e68589bSMark F. Adams   /* create global vector of data in 'data_w_ghost' */
665849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROLA], 0, 0, 0, 0));
666bae903cbSmarkadams4   if (size > 1) { /* get ghost null space data */
6672e68589bSMark F. Adams     PetscReal *tmp_gdata, *tmp_ldata, *tp2;
6689566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nloc, &tmp_ldata));
6694a2f8832SBarry Smith     for (jj = 0; jj < col_bs; jj++) {
670c8b0795cSMark F. Adams       for (kk = 0; kk < bs; kk++) {
671a2f3521dSMark F. Adams         PetscInt         ii, stride;
672c8b0795cSMark F. Adams         const PetscReal *tp = pc_gamg->data + jj * bs * nloc + kk;
6732fa5cd67SKarl Rupp         for (ii = 0; ii < nloc; ii++, tp += bs) tmp_ldata[ii] = *tp;
6742fa5cd67SKarl Rupp 
6759566063dSJacob Faibussowitsch         PetscCall(PCGAMGGetDataWithGhosts(Gmat, 1, tmp_ldata, &stride, &tmp_gdata));
676a2f3521dSMark F. Adams 
677bae903cbSmarkadams4         if (!jj && !kk) { /* now I know how many total nodes - allocate TODO: move below and do in one 'col_bs' call */
6789566063dSJacob Faibussowitsch           PetscCall(PetscMalloc1(stride * bs * col_bs, &data_w_ghost));
679a2f3521dSMark F. Adams           nbnodes = bs * stride;
6802e68589bSMark F. Adams         }
681a2f3521dSMark F. Adams         tp2 = data_w_ghost + jj * bs * stride + kk;
6822fa5cd67SKarl Rupp         for (ii = 0; ii < stride; ii++, tp2 += bs) *tp2 = tmp_gdata[ii];
6839566063dSJacob Faibussowitsch         PetscCall(PetscFree(tmp_gdata));
6842e68589bSMark F. Adams       }
6852e68589bSMark F. Adams     }
6869566063dSJacob Faibussowitsch     PetscCall(PetscFree(tmp_ldata));
687806fa848SBarry Smith   } else {
688c8b0795cSMark F. Adams     nbnodes      = bs * nloc;
689c8b0795cSMark F. Adams     data_w_ghost = (PetscReal *)pc_gamg->data;
6902e68589bSMark F. Adams   }
6912e68589bSMark F. Adams 
692bae903cbSmarkadams4   /* get 'flid_fgid' TODO - move up to get 'stride' and do get null space data above in one step (jj loop) */
693c5df96a5SBarry Smith   if (size > 1) {
6942e68589bSMark F. Adams     PetscReal *fid_glid_loc, *fiddata;
695a2f3521dSMark F. Adams     PetscInt   stride;
6962e68589bSMark F. Adams 
6979566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nloc, &fid_glid_loc));
6982e68589bSMark F. Adams     for (kk = 0; kk < nloc; kk++) fid_glid_loc[kk] = (PetscReal)(my0 + kk);
6999566063dSJacob Faibussowitsch     PetscCall(PCGAMGGetDataWithGhosts(Gmat, 1, fid_glid_loc, &stride, &fiddata));
700bae903cbSmarkadams4     PetscCall(PetscMalloc1(stride, &flid_fgid)); /* copy real data to in */
701a2f3521dSMark F. Adams     for (kk = 0; kk < stride; kk++) flid_fgid[kk] = (PetscInt)fiddata[kk];
7029566063dSJacob Faibussowitsch     PetscCall(PetscFree(fiddata));
703a2f3521dSMark F. Adams 
70463a3b9bcSJacob Faibussowitsch     PetscCheck(stride == nbnodes / bs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "stride %" PetscInt_FMT " != nbnodes %" PetscInt_FMT "/bs %" PetscInt_FMT, stride, nbnodes, bs);
7059566063dSJacob Faibussowitsch     PetscCall(PetscFree(fid_glid_loc));
706806fa848SBarry Smith   } else {
7079566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nloc, &flid_fgid));
7082e68589bSMark F. Adams     for (kk = 0; kk < nloc; kk++) flid_fgid[kk] = my0 + kk;
7092e68589bSMark F. Adams   }
710849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROLA], 0, 0, 0, 0));
711bae903cbSmarkadams4   /* get P0 */
712849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROLB], 0, 0, 0, 0));
713c8b0795cSMark F. Adams   {
7140298fd71SBarry Smith     PetscReal *data_out = NULL;
7159566063dSJacob Faibussowitsch     PetscCall(formProl0(agg_lists, bs, col_bs, myCrs0, nbnodes, data_w_ghost, flid_fgid, &data_out, Prol));
7169566063dSJacob Faibussowitsch     PetscCall(PetscFree(pc_gamg->data));
7172fa5cd67SKarl Rupp 
718c8b0795cSMark F. Adams     pc_gamg->data           = data_out;
7194a2f8832SBarry Smith     pc_gamg->data_cell_rows = col_bs;
7204a2f8832SBarry Smith     pc_gamg->data_sz        = col_bs * col_bs * nLocalSelected;
721c8b0795cSMark F. Adams   }
722849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROLB], 0, 0, 0, 0));
723*48a46eb9SPierre Jolivet   if (size > 1) PetscCall(PetscFree(data_w_ghost));
7249566063dSJacob Faibussowitsch   PetscCall(PetscFree(flid_fgid));
7252e68589bSMark F. Adams 
726c8b0795cSMark F. Adams   *a_P_out = Prol; /* out */
727ed4e82eaSPeter Brune 
728849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0));
729c8b0795cSMark F. Adams   PetscFunctionReturn(0);
730c8b0795cSMark F. Adams }
731c8b0795cSMark F. Adams 
732c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
733c8b0795cSMark F. Adams /*
734fd1112cbSBarry Smith    PCGAMGOptProlongator_AGG
735c8b0795cSMark F. Adams 
736c8b0795cSMark F. Adams   Input Parameter:
737c8b0795cSMark F. Adams    . pc - this
738c8b0795cSMark F. Adams    . Amat - matrix on this fine level
739c8b0795cSMark F. Adams  In/Output Parameter:
7402a808120SBarry Smith    . a_P - prolongation operator to the next level
741c8b0795cSMark F. Adams */
7429371c9d4SSatish Balay static PetscErrorCode PCGAMGOptProlongator_AGG(PC pc, Mat Amat, Mat *a_P) {
743c8b0795cSMark F. Adams   PC_MG       *mg          = (PC_MG *)pc->data;
744c8b0795cSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
745c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
746c8b0795cSMark F. Adams   PetscInt     jj;
747c8b0795cSMark F. Adams   Mat          Prol = *a_P;
7483b4367a7SBarry Smith   MPI_Comm     comm;
7492a808120SBarry Smith   KSP          eksp;
7502a808120SBarry Smith   Vec          bb, xx;
7512a808120SBarry Smith   PC           epc;
7522a808120SBarry Smith   PetscReal    alpha, emax, emin;
753c8b0795cSMark F. Adams 
754c8b0795cSMark F. Adams   PetscFunctionBegin;
7559566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm));
756849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_OPT], 0, 0, 0, 0));
757c8b0795cSMark F. Adams 
758d5d25401SBarry Smith   /* compute maximum singular value of operator to be used in smoother */
7592a808120SBarry Smith   if (0 < pc_gamg_agg->nsmooths) {
76018c3aa7eSMark     /* get eigen estimates */
76118c3aa7eSMark     if (pc_gamg->emax > 0) {
76218c3aa7eSMark       emin = pc_gamg->emin;
76318c3aa7eSMark       emax = pc_gamg->emax;
76418c3aa7eSMark     } else {
7650ed2132dSStefano Zampini       const char *prefix;
7660ed2132dSStefano Zampini 
7679566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(Amat, &bb, NULL));
7689566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(Amat, &xx, NULL));
7699566063dSJacob Faibussowitsch       PetscCall(KSPSetNoisy_Private(bb));
770e696ed0bSMark F. Adams 
7719566063dSJacob Faibussowitsch       PetscCall(KSPCreate(comm, &eksp));
7729566063dSJacob Faibussowitsch       PetscCall(PCGetOptionsPrefix(pc, &prefix));
7739566063dSJacob Faibussowitsch       PetscCall(KSPSetOptionsPrefix(eksp, prefix));
7749566063dSJacob Faibussowitsch       PetscCall(KSPAppendOptionsPrefix(eksp, "pc_gamg_esteig_"));
77573f7197eSJed Brown       {
776b94d7dedSBarry Smith         PetscBool isset, sflg;
777b94d7dedSBarry Smith         PetscCall(MatIsSPDKnown(Amat, &isset, &sflg));
778b94d7dedSBarry Smith         if (isset && sflg) PetscCall(KSPSetType(eksp, KSPCG));
779d24ecf33SMark       }
7809566063dSJacob Faibussowitsch       PetscCall(KSPSetErrorIfNotConverged(eksp, pc->erroriffailure));
7819566063dSJacob Faibussowitsch       PetscCall(KSPSetNormType(eksp, KSP_NORM_NONE));
782c2436cacSMark F. Adams 
7839566063dSJacob Faibussowitsch       PetscCall(KSPSetInitialGuessNonzero(eksp, PETSC_FALSE));
7849566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(eksp, Amat, Amat));
7852e68589bSMark F. Adams 
7869566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(eksp, &epc));
7879566063dSJacob Faibussowitsch       PetscCall(PCSetType(epc, PCJACOBI)); /* smoother in smoothed agg. */
788c2436cacSMark F. Adams 
7899566063dSJacob 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
790074aec55SMark Adams 
7919566063dSJacob Faibussowitsch       PetscCall(KSPSetFromOptions(eksp));
7929566063dSJacob Faibussowitsch       PetscCall(KSPSetComputeSingularValues(eksp, PETSC_TRUE));
7939566063dSJacob Faibussowitsch       PetscCall(KSPSolve(eksp, bb, xx));
7949566063dSJacob Faibussowitsch       PetscCall(KSPCheckSolve(eksp, pc, xx));
7952e68589bSMark F. Adams 
7969566063dSJacob Faibussowitsch       PetscCall(KSPComputeExtremeSingularValues(eksp, &emax, &emin));
7979566063dSJacob Faibussowitsch       PetscCall(PetscInfo(pc, "%s: Smooth P0: max eigen=%e min=%e PC=%s\n", ((PetscObject)pc)->prefix, (double)emax, (double)emin, PCJACOBI));
7989566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&xx));
7999566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&bb));
8009566063dSJacob Faibussowitsch       PetscCall(KSPDestroy(&eksp));
8012e68589bSMark F. Adams     }
80218c3aa7eSMark     if (pc_gamg->use_sa_esteig) {
80318c3aa7eSMark       mg->min_eigen_DinvA[pc_gamg->current_level] = emin;
80418c3aa7eSMark       mg->max_eigen_DinvA[pc_gamg->current_level] = emax;
80563a3b9bcSJacob 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));
80618c3aa7eSMark     } else {
80718c3aa7eSMark       mg->min_eigen_DinvA[pc_gamg->current_level] = 0;
80818c3aa7eSMark       mg->max_eigen_DinvA[pc_gamg->current_level] = 0;
80918c3aa7eSMark     }
81018c3aa7eSMark   } else {
81118c3aa7eSMark     mg->min_eigen_DinvA[pc_gamg->current_level] = 0;
81218c3aa7eSMark     mg->max_eigen_DinvA[pc_gamg->current_level] = 0;
81318c3aa7eSMark   }
8142e68589bSMark F. Adams 
8152a808120SBarry Smith   /* smooth P0 */
8162a808120SBarry Smith   for (jj = 0; jj < pc_gamg_agg->nsmooths; jj++) {
8172a808120SBarry Smith     Mat tMat;
8182a808120SBarry Smith     Vec diag;
8192a808120SBarry Smith 
820849bee69SMark Adams     PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_OPTSM], 0, 0, 0, 0));
8212a808120SBarry Smith 
8222e68589bSMark F. Adams     /* smooth P1 := (I - omega/lam D^{-1}A)P0 */
8239566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2], 0, 0, 0, 0));
8249566063dSJacob Faibussowitsch     PetscCall(MatMatMult(Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat));
8259566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2], 0, 0, 0, 0));
8269566063dSJacob Faibussowitsch     PetscCall(MatProductClear(tMat));
8279566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(Amat, &diag, NULL));
8289566063dSJacob Faibussowitsch     PetscCall(MatGetDiagonal(Amat, diag)); /* effectively PCJACOBI */
8299566063dSJacob Faibussowitsch     PetscCall(VecReciprocal(diag));
8309566063dSJacob Faibussowitsch     PetscCall(MatDiagonalScale(tMat, diag, NULL));
8319566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&diag));
832b4da3a1bSStefano Zampini 
833d5d25401SBarry Smith     /* TODO: Set a PCFailedReason and exit the building of the AMG preconditioner */
83408401ef6SPierre Jolivet     PetscCheck(emax != 0.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_PLIB, "Computed maximum singular value as zero");
835d5d25401SBarry Smith     /* TODO: Document the 1.4 and don't hardwire it in this routine */
836b4ec6429SMark F. Adams     alpha = -1.4 / emax;
837b4da3a1bSStefano Zampini 
8389566063dSJacob Faibussowitsch     PetscCall(MatAYPX(tMat, alpha, Prol, SUBSET_NONZERO_PATTERN));
8399566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Prol));
8402e68589bSMark F. Adams     Prol = tMat;
841849bee69SMark Adams     PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_OPTSM], 0, 0, 0, 0));
8422e68589bSMark F. Adams   }
843849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_OPT], 0, 0, 0, 0));
844c8b0795cSMark F. Adams   *a_P = Prol;
845c8b0795cSMark F. Adams   PetscFunctionReturn(0);
8462e68589bSMark F. Adams }
8472e68589bSMark F. Adams 
848c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
849c8b0795cSMark F. Adams /*
850c8b0795cSMark F. Adams    PCCreateGAMG_AGG
8512e68589bSMark F. Adams 
852c8b0795cSMark F. Adams   Input Parameter:
853c8b0795cSMark F. Adams    . pc -
854c8b0795cSMark F. Adams */
8559371c9d4SSatish Balay PetscErrorCode PCCreateGAMG_AGG(PC pc) {
856c8b0795cSMark F. Adams   PC_MG       *mg      = (PC_MG *)pc->data;
857c8b0795cSMark F. Adams   PC_GAMG     *pc_gamg = (PC_GAMG *)mg->innerctx;
858c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg;
8592e68589bSMark F. Adams 
860c8b0795cSMark F. Adams   PetscFunctionBegin;
861c8b0795cSMark F. Adams   /* create sub context for SA */
8629566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(pc, &pc_gamg_agg));
863c8b0795cSMark F. Adams   pc_gamg->subctx = pc_gamg_agg;
864c8b0795cSMark F. Adams 
8651ab5ffc9SJed Brown   pc_gamg->ops->setfromoptions = PCSetFromOptions_GAMG_AGG;
8669b8ffb57SJed Brown   pc_gamg->ops->destroy        = PCDestroy_GAMG_AGG;
867c8b0795cSMark F. Adams   /* reset does not do anything; setup not virtual */
868c8b0795cSMark F. Adams 
869c8b0795cSMark F. Adams   /* set internal function pointers */
870fd1112cbSBarry Smith   pc_gamg->ops->graph             = PCGAMGGraph_AGG;
8711ab5ffc9SJed Brown   pc_gamg->ops->coarsen           = PCGAMGCoarsen_AGG;
8721ab5ffc9SJed Brown   pc_gamg->ops->prolongator       = PCGAMGProlongator_AGG;
873fd1112cbSBarry Smith   pc_gamg->ops->optprolongator    = PCGAMGOptProlongator_AGG;
8741ab5ffc9SJed Brown   pc_gamg->ops->createdefaultdata = PCSetData_AGG;
8755adeb434SBarry Smith   pc_gamg->ops->view              = PCView_GAMG_AGG;
876c8b0795cSMark F. Adams 
877bae903cbSmarkadams4   pc_gamg_agg->aggressive_coarsening_levels = 0;
878bae903cbSmarkadams4   pc_gamg_agg->symmetrize_graph             = PETSC_FALSE;
879cab9ed1eSBarry Smith   pc_gamg_agg->nsmooths                     = 1;
880cab9ed1eSBarry Smith 
8819566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNSmooths_C", PCGAMGSetNSmooths_AGG));
882bae903cbSmarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetSymmetrizeGraph_C", PCGAMGSetSymmetrizeGraph_AGG));
883bae903cbSmarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveLevels_C", PCGAMGSetAggressiveLevels_AGG));
8849566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", PCSetCoordinates_AGG));
8852e68589bSMark F. Adams   PetscFunctionReturn(0);
8862e68589bSMark F. Adams }
887