xref: /petsc/src/ksp/pc/impls/gamg/agg.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
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 /*@
17f1580f4eSBarry Smith    PCGAMGSetNSmooths - Set number of smoothing steps (1 is typical) used for multigrid on all the levels
182e68589bSMark F. Adams 
19f1580f4eSBarry 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 
29f1580f4eSBarry Smith .seealso: `PCMG`, `PCGAMG`
302e68589bSMark F. Adams @*/
31*d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n)
32*d71ae5a4SJacob Faibussowitsch {
332e68589bSMark F. Adams   PetscFunctionBegin;
342e68589bSMark F. Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
35d5d25401SBarry Smith   PetscValidLogicalCollectiveInt(pc, n, 2);
36cac4c232SBarry Smith   PetscTryMethod(pc, "PCGAMGSetNSmooths_C", (PC, PetscInt), (pc, n));
372e68589bSMark F. Adams   PetscFunctionReturn(0);
382e68589bSMark F. Adams }
392e68589bSMark F. Adams 
40*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetNSmooths_AGG(PC pc, PetscInt n)
41*d71ae5a4SJacob Faibussowitsch {
422e68589bSMark F. Adams   PC_MG       *mg          = (PC_MG *)pc->data;
432e68589bSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
44c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
452e68589bSMark F. Adams 
462e68589bSMark F. Adams   PetscFunctionBegin;
47c8b0795cSMark F. Adams   pc_gamg_agg->nsmooths = n;
48c8b0795cSMark F. Adams   PetscFunctionReturn(0);
49c8b0795cSMark F. Adams }
50c8b0795cSMark F. Adams 
51c8b0795cSMark F. Adams /*@
52bae903cbSmarkadams4    PCGAMGSetSymmetrizeGraph - Symmetrize the graph before computing the aggregation. Some algorithms require the graph be symmetric
53c8b0795cSMark F. Adams 
54f1580f4eSBarry Smith    Logically Collective on pc
55c8b0795cSMark F. Adams 
56c8b0795cSMark F. Adams    Input Parameters:
57cab9ed1eSBarry Smith +  pc - the preconditioner context
58f1580f4eSBarry Smith -  n - `PETSC_TRUE` or `PETSC_FALSE`
59c8b0795cSMark F. Adams 
60c8b0795cSMark F. Adams    Options Database Key:
61bae903cbSmarkadams4 .  -pc_gamg_symmetrize_graph <true,default=false> - symmetrize the graph before computing the aggregation
62c8b0795cSMark F. Adams 
63c8b0795cSMark F. Adams    Level: intermediate
64c8b0795cSMark F. Adams 
65f1580f4eSBarry Smith    Developer Note:
66f1580f4eSBarry Smith    If the aggregation can hang with a nonsymmetric graph then the graph should always be symmetrized before calling the aggregation,
67f1580f4eSBarry Smith    it should not be up to the user.
68f1580f4eSBarry Smith 
69f1580f4eSBarry Smith .seealso: `PCGAMG`, `PCGAMGSetAggressiveLevels()`
70c8b0795cSMark F. Adams @*/
71*d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetSymmetrizeGraph(PC pc, PetscBool n)
72*d71ae5a4SJacob Faibussowitsch {
73c8b0795cSMark F. Adams   PetscFunctionBegin;
74c8b0795cSMark F. Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
75d5d25401SBarry Smith   PetscValidLogicalCollectiveBool(pc, n, 2);
76bae903cbSmarkadams4   PetscTryMethod(pc, "PCGAMGSetSymmetrizeGraph_C", (PC, PetscBool), (pc, n));
77c8b0795cSMark F. Adams   PetscFunctionReturn(0);
78c8b0795cSMark F. Adams }
79c8b0795cSMark F. Adams 
80*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetSymmetrizeGraph_AGG(PC pc, PetscBool n)
81*d71ae5a4SJacob Faibussowitsch {
82c8b0795cSMark F. Adams   PC_MG       *mg          = (PC_MG *)pc->data;
83c8b0795cSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
84c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
85c8b0795cSMark F. Adams 
86c8b0795cSMark F. Adams   PetscFunctionBegin;
87bae903cbSmarkadams4   pc_gamg_agg->symmetrize_graph = n;
882e68589bSMark F. Adams   PetscFunctionReturn(0);
892e68589bSMark F. Adams }
902e68589bSMark F. Adams 
91ef4ad70eSMark F. Adams /*@
92f1580f4eSBarry Smith    PCGAMGSetAggressiveLevels -  Use aggressive coarsening on first n levels
93ef4ad70eSMark F. Adams 
94f1580f4eSBarry Smith    Logically Collective on pc
95ef4ad70eSMark F. Adams 
96ef4ad70eSMark F. Adams    Input Parameters:
97cab9ed1eSBarry Smith +  pc - the preconditioner context
98d5d25401SBarry Smith -  n - 0, 1 or more
99ef4ad70eSMark F. Adams 
100ef4ad70eSMark F. Adams    Options Database Key:
101bae903cbSmarkadams4 .  -pc_gamg_aggressive_coarsening <n,default = 1> - Number of levels to square the graph on before aggregating it
102af3c827dSMark Adams 
103ef4ad70eSMark F. Adams    Level: intermediate
104ef4ad70eSMark F. Adams 
105f1580f4eSBarry Smith .seealso: `PCGAMG`, `PCGAMGSetSymmetrizeGraph()`, `PCGAMGSetThreshold()`
106ef4ad70eSMark F. Adams @*/
107*d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGAMGSetAggressiveLevels(PC pc, PetscInt n)
108*d71ae5a4SJacob Faibussowitsch {
109ef4ad70eSMark F. Adams   PetscFunctionBegin;
110ef4ad70eSMark F. Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
111d5d25401SBarry Smith   PetscValidLogicalCollectiveInt(pc, n, 2);
112bae903cbSmarkadams4   PetscTryMethod(pc, "PCGAMGSetAggressiveLevels_C", (PC, PetscInt), (pc, n));
113ef4ad70eSMark F. Adams   PetscFunctionReturn(0);
114ef4ad70eSMark F. Adams }
115ef4ad70eSMark F. Adams 
116*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGSetAggressiveLevels_AGG(PC pc, PetscInt n)
117*d71ae5a4SJacob Faibussowitsch {
118ef4ad70eSMark F. Adams   PC_MG       *mg          = (PC_MG *)pc->data;
119ef4ad70eSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
120ef4ad70eSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
121ef4ad70eSMark F. Adams 
122ef4ad70eSMark F. Adams   PetscFunctionBegin;
123bae903cbSmarkadams4   pc_gamg_agg->aggressive_coarsening_levels = n;
124ef4ad70eSMark F. Adams   PetscFunctionReturn(0);
125ef4ad70eSMark F. Adams }
126ef4ad70eSMark F. Adams 
127*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_GAMG_AGG(PC pc, PetscOptionItems *PetscOptionsObject)
128*d71ae5a4SJacob Faibussowitsch {
1292e68589bSMark F. Adams   PC_MG       *mg          = (PC_MG *)pc->data;
1302e68589bSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
131c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
1322e68589bSMark F. Adams 
1332e68589bSMark F. Adams   PetscFunctionBegin;
134d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "GAMG-AGG options");
1352e68589bSMark F. Adams   {
136bae903cbSmarkadams4     PetscBool flg;
1379566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-pc_gamg_agg_nsmooths", "smoothing steps for smoothed aggregation, usually 1", "PCGAMGSetNSmooths", pc_gamg_agg->nsmooths, &pc_gamg_agg->nsmooths, NULL));
138bae903cbSmarkadams4     PetscCall(PetscOptionsBool("-pc_gamg_symmetrize_graph", "Set for asymmetric matrices", "PCGAMGSetSymmetrizeGraph", pc_gamg_agg->symmetrize_graph, &pc_gamg_agg->symmetrize_graph, NULL));
139bae903cbSmarkadams4     pc_gamg_agg->aggressive_coarsening_levels = 1;
1409371c9d4SSatish Balay     PetscCall(
1419371c9d4SSatish 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));
142bae903cbSmarkadams4     if (!flg) {
143bae903cbSmarkadams4       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));
144bae903cbSmarkadams4     } else {
145bae903cbSmarkadams4       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));
146bae903cbSmarkadams4       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));
147bae903cbSmarkadams4     }
1482e68589bSMark F. Adams   }
149d0609cedSBarry Smith   PetscOptionsHeadEnd();
1502e68589bSMark F. Adams   PetscFunctionReturn(0);
1512e68589bSMark F. Adams }
1522e68589bSMark F. Adams 
153*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_GAMG_AGG(PC pc)
154*d71ae5a4SJacob Faibussowitsch {
1552e68589bSMark F. Adams   PC_MG   *mg      = (PC_MG *)pc->data;
1562e68589bSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
1572e68589bSMark F. Adams 
1582e68589bSMark F. Adams   PetscFunctionBegin;
1599566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg->subctx));
1602e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNSmooths_C", NULL));
161bae903cbSmarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetSymmetrizeGraph_C", NULL));
162bae903cbSmarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveLevels_C", NULL));
163bae903cbSmarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", NULL));
1642e68589bSMark F. Adams   PetscFunctionReturn(0);
1652e68589bSMark F. Adams }
1662e68589bSMark F. Adams 
1672e68589bSMark F. Adams /*
1682e68589bSMark F. Adams    PCSetCoordinates_AGG
169302f38e8SMark F. Adams      - collective
1702e68589bSMark F. Adams 
1712e68589bSMark F. Adams    Input Parameter:
1722e68589bSMark F. Adams    . pc - the preconditioner context
173a2f3521dSMark F. Adams    . ndm - dimesion of data (used for dof/vertex for Stokes)
174302f38e8SMark F. Adams    . a_nloc - number of vertices local
175302f38e8SMark F. Adams    . coords - [a_nloc][ndm] - interleaved coordinate data: {x_0, y_0, z_0, x_1, y_1, ...}
1762e68589bSMark F. Adams */
1771e6b0712SBarry Smith 
178*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetCoordinates_AGG(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords)
179*d71ae5a4SJacob Faibussowitsch {
1802e68589bSMark F. Adams   PC_MG   *mg      = (PC_MG *)pc->data;
1812e68589bSMark F. Adams   PC_GAMG *pc_gamg = (PC_GAMG *)mg->innerctx;
18269344418SMark F. Adams   PetscInt arrsz, kk, ii, jj, nloc, ndatarows, ndf;
183a2f3521dSMark F. Adams   Mat      mat = pc->pmat;
1842e68589bSMark F. Adams 
1852e68589bSMark F. Adams   PetscFunctionBegin;
186a2f3521dSMark F. Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
187a2f3521dSMark F. Adams   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
188302f38e8SMark F. Adams   nloc = a_nloc;
1892e68589bSMark F. Adams 
1902e68589bSMark F. Adams   /* SA: null space vectors */
1919566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(mat, &ndf));               /* this does not work for Stokes */
19269344418SMark F. Adams   if (coords && ndf == 1) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */
193a2f3521dSMark F. Adams   else if (coords) {
19463a3b9bcSJacob Faibussowitsch     PetscCheck(ndm <= ndf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "degrees of motion %" PetscInt_FMT " > block size %" PetscInt_FMT, ndm, ndf);
19569344418SMark F. Adams     pc_gamg->data_cell_cols = (ndm == 2 ? 3 : 6); /* displacement elasticity */
196ad540459SPierre 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);
197806fa848SBarry Smith   } else pc_gamg->data_cell_cols = ndf; /* no data, force SA with constant null space vectors */
19871959b99SBarry Smith   pc_gamg->data_cell_rows = ndatarows = ndf;
19963a3b9bcSJacob 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);
200c8b0795cSMark F. Adams   arrsz = nloc * pc_gamg->data_cell_rows * pc_gamg->data_cell_cols;
2012e68589bSMark F. Adams 
2027f66b68fSMark Adams   if (!pc_gamg->data || (pc_gamg->data_sz != arrsz)) {
2039566063dSJacob Faibussowitsch     PetscCall(PetscFree(pc_gamg->data));
2049566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(arrsz + 1, &pc_gamg->data));
2052e68589bSMark F. Adams   }
2062e68589bSMark F. Adams   /* copy data in - column oriented */
2072e68589bSMark F. Adams   for (kk = 0; kk < nloc; kk++) {
20869344418SMark F. Adams     const PetscInt M    = nloc * pc_gamg->data_cell_rows; /* stride into data */
20969344418SMark F. Adams     PetscReal     *data = &pc_gamg->data[kk * ndatarows]; /* start of cell */
210c8b0795cSMark F. Adams     if (pc_gamg->data_cell_cols == 1) *data = 1.0;
2112e68589bSMark F. Adams     else {
21269344418SMark F. Adams       /* translational modes */
2132fa5cd67SKarl Rupp       for (ii = 0; ii < ndatarows; ii++) {
2142fa5cd67SKarl Rupp         for (jj = 0; jj < ndatarows; jj++) {
21569344418SMark F. Adams           if (ii == jj) data[ii * M + jj] = 1.0;
2162e68589bSMark F. Adams           else data[ii * M + jj] = 0.0;
2172fa5cd67SKarl Rupp         }
2182fa5cd67SKarl Rupp       }
21969344418SMark F. Adams 
22069344418SMark F. Adams       /* rotational modes */
2212e68589bSMark F. Adams       if (coords) {
22269344418SMark F. Adams         if (ndm == 2) {
2232e68589bSMark F. Adams           data += 2 * M;
2242e68589bSMark F. Adams           data[0] = -coords[2 * kk + 1];
2252e68589bSMark F. Adams           data[1] = coords[2 * kk];
226806fa848SBarry Smith         } else {
2272e68589bSMark F. Adams           data += 3 * M;
2289371c9d4SSatish Balay           data[0]         = 0.0;
2299371c9d4SSatish Balay           data[M + 0]     = coords[3 * kk + 2];
2309371c9d4SSatish Balay           data[2 * M + 0] = -coords[3 * kk + 1];
2319371c9d4SSatish Balay           data[1]         = -coords[3 * kk + 2];
2329371c9d4SSatish Balay           data[M + 1]     = 0.0;
2339371c9d4SSatish Balay           data[2 * M + 1] = coords[3 * kk];
2349371c9d4SSatish Balay           data[2]         = coords[3 * kk + 1];
2359371c9d4SSatish Balay           data[M + 2]     = -coords[3 * kk];
2369371c9d4SSatish Balay           data[2 * M + 2] = 0.0;
2372e68589bSMark F. Adams         }
2382e68589bSMark F. Adams       }
2392e68589bSMark F. Adams     }
2402e68589bSMark F. Adams   }
2412e68589bSMark F. Adams   pc_gamg->data_sz = arrsz;
2422e68589bSMark F. Adams   PetscFunctionReturn(0);
2432e68589bSMark F. Adams }
2442e68589bSMark F. Adams 
2452e68589bSMark F. Adams /*
246a2f3521dSMark F. Adams    PCSetData_AGG - called if data is not set with PCSetCoordinates.
247a2f3521dSMark F. Adams       Looks in Mat for near null space.
248a2f3521dSMark F. Adams       Does not work for Stokes
2492e68589bSMark F. Adams 
2502e68589bSMark F. Adams   Input Parameter:
2512e68589bSMark F. Adams    . pc -
252a2f3521dSMark F. Adams    . a_A - matrix to get (near) null space out of.
2532e68589bSMark F. Adams */
254*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetData_AGG(PC pc, Mat a_A)
255*d71ae5a4SJacob Faibussowitsch {
256b8cd405aSMark F. Adams   PC_MG       *mg      = (PC_MG *)pc->data;
257b8cd405aSMark F. Adams   PC_GAMG     *pc_gamg = (PC_GAMG *)mg->innerctx;
258b8cd405aSMark F. Adams   MatNullSpace mnull;
25966f2ef4bSMark Adams 
260b817416eSBarry Smith   PetscFunctionBegin;
2619566063dSJacob Faibussowitsch   PetscCall(MatGetNearNullSpace(a_A, &mnull));
262b8cd405aSMark F. Adams   if (!mnull) {
26366f2ef4bSMark Adams     DM dm;
2649566063dSJacob Faibussowitsch     PetscCall(PCGetDM(pc, &dm));
26548a46eb9SPierre Jolivet     if (!dm) PetscCall(MatGetDM(a_A, &dm));
26666f2ef4bSMark Adams     if (dm) {
26766f2ef4bSMark Adams       PetscObject deformation;
268b0d30dd6SMatthew G. Knepley       PetscInt    Nf;
269b0d30dd6SMatthew G. Knepley 
2709566063dSJacob Faibussowitsch       PetscCall(DMGetNumFields(dm, &Nf));
271b0d30dd6SMatthew G. Knepley       if (Nf) {
2729566063dSJacob Faibussowitsch         PetscCall(DMGetField(dm, 0, NULL, &deformation));
2739566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery((PetscObject)deformation, "nearnullspace", (PetscObject *)&mnull));
27448a46eb9SPierre Jolivet         if (!mnull) PetscCall(PetscObjectQuery((PetscObject)deformation, "nullspace", (PetscObject *)&mnull));
27566f2ef4bSMark Adams       }
27666f2ef4bSMark Adams     }
277b0d30dd6SMatthew G. Knepley   }
27866f2ef4bSMark Adams 
27966f2ef4bSMark Adams   if (!mnull) {
280a2f3521dSMark F. Adams     PetscInt bs, NN, MM;
2819566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSize(a_A, &bs));
2829566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(a_A, &MM, &NN));
28363a3b9bcSJacob Faibussowitsch     PetscCheck(MM % bs == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MM %" PetscInt_FMT " must be divisible by bs %" PetscInt_FMT, MM, bs);
2849566063dSJacob Faibussowitsch     PetscCall(PCSetCoordinates_AGG(pc, bs, MM / bs, NULL));
285806fa848SBarry Smith   } else {
286b8cd405aSMark F. Adams     PetscReal         *nullvec;
287b8cd405aSMark F. Adams     PetscBool          has_const;
288b8cd405aSMark F. Adams     PetscInt           i, j, mlocal, nvec, bs;
289d5d25401SBarry Smith     const Vec         *vecs;
290d5d25401SBarry Smith     const PetscScalar *v;
291b817416eSBarry Smith 
2929566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(a_A, &mlocal, NULL));
2939566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs));
294a0dea87dSMark Adams     pc_gamg->data_sz = (nvec + !!has_const) * mlocal;
2959566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1((nvec + !!has_const) * mlocal, &nullvec));
2969371c9d4SSatish Balay     if (has_const)
2979371c9d4SSatish Balay       for (i = 0; i < mlocal; i++) nullvec[i] = 1.0;
298b8cd405aSMark F. Adams     for (i = 0; i < nvec; i++) {
2999566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(vecs[i], &v));
300b8cd405aSMark F. Adams       for (j = 0; j < mlocal; j++) nullvec[(i + !!has_const) * mlocal + j] = PetscRealPart(v[j]);
3019566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(vecs[i], &v));
302b8cd405aSMark F. Adams     }
303b8cd405aSMark F. Adams     pc_gamg->data           = nullvec;
304b8cd405aSMark F. Adams     pc_gamg->data_cell_cols = (nvec + !!has_const);
3059566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSize(a_A, &bs));
306b8cd405aSMark F. Adams     pc_gamg->data_cell_rows = bs;
307b8cd405aSMark F. Adams   }
3082e68589bSMark F. Adams   PetscFunctionReturn(0);
3092e68589bSMark F. Adams }
3102e68589bSMark F. Adams 
3112e68589bSMark F. Adams /*
312bae903cbSmarkadams4   formProl0 - collect null space data for each aggregate, do QR, put R in coarse grid data and Q in P_0
3132e68589bSMark F. Adams 
3142e68589bSMark F. Adams   Input Parameter:
3152adfe9d3SBarry Smith    . agg_llists - list of arrays with aggregates -- list from selected vertices of aggregate unselected vertices
3162adfe9d3SBarry Smith    . bs - row block size
3170cbbd2e1SMark F. Adams    . nSAvec - column bs of new P
3180cbbd2e1SMark F. Adams    . my0crs - global index of start of locals
3192adfe9d3SBarry Smith    . data_stride - bs*(nloc nodes + ghost nodes) [data_stride][nSAvec]
3202e68589bSMark F. Adams    . data_in[data_stride*nSAvec] - local data on fine grid
3212e68589bSMark F. Adams    . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist'
322bae903cbSmarkadams4 
3232e68589bSMark F. Adams   Output Parameter:
3242e68589bSMark F. Adams    . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data
3252e68589bSMark F. Adams    . a_Prol - prolongation operator
3262e68589bSMark F. Adams */
327*d71ae5a4SJacob Faibussowitsch static PetscErrorCode formProl0(PetscCoarsenData *agg_llists, PetscInt bs, PetscInt nSAvec, PetscInt my0crs, PetscInt data_stride, PetscReal data_in[], const PetscInt flid_fgid[], PetscReal **a_data_out, Mat a_Prol)
328*d71ae5a4SJacob Faibussowitsch {
329bd026e97SJed Brown   PetscInt        Istart, my0, Iend, nloc, clid, flid = 0, aggID, kk, jj, ii, mm, nSelected, minsz, nghosts, out_data_stride;
3303b4367a7SBarry Smith   MPI_Comm        comm;
3312e68589bSMark F. Adams   PetscReal      *out_data;
332539c167fSBarry Smith   PetscCDIntNd   *pos;
3331943db53SBarry Smith   PCGAMGHashTable fgid_flid;
3340cbbd2e1SMark F. Adams 
3352e68589bSMark F. Adams   PetscFunctionBegin;
3369566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)a_Prol, &comm));
3379566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(a_Prol, &Istart, &Iend));
3389371c9d4SSatish Balay   nloc = (Iend - Istart) / bs;
3399371c9d4SSatish Balay   my0  = Istart / bs;
34063a3b9bcSJacob 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);
3410cbbd2e1SMark F. Adams   Iend /= bs;
3420cbbd2e1SMark F. Adams   nghosts = data_stride / bs - nloc;
3432e68589bSMark F. Adams 
3449566063dSJacob Faibussowitsch   PetscCall(PCGAMGHashTableCreate(2 * nghosts + 1, &fgid_flid));
34548a46eb9SPierre Jolivet   for (kk = 0; kk < nghosts; kk++) PetscCall(PCGAMGHashTableAdd(&fgid_flid, flid_fgid[nloc + kk], nloc + kk));
3462e68589bSMark F. Adams 
3470cbbd2e1SMark F. Adams   /* count selected -- same as number of cols of P */
3480cbbd2e1SMark F. Adams   for (nSelected = mm = 0; mm < nloc; mm++) {
349e78576d6SMark F. Adams     PetscBool ise;
3509566063dSJacob Faibussowitsch     PetscCall(PetscCDEmptyAt(agg_llists, mm, &ise));
351e78576d6SMark F. Adams     if (!ise) nSelected++;
3520cbbd2e1SMark F. Adams   }
3539566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRangeColumn(a_Prol, &ii, &jj));
35463a3b9bcSJacob Faibussowitsch   PetscCheck((ii / nSAvec) == my0crs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ii %" PetscInt_FMT " /nSAvec %" PetscInt_FMT "  != my0crs %" PetscInt_FMT, ii, nSAvec, my0crs);
35563a3b9bcSJacob 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);
3560cbbd2e1SMark F. Adams 
3572e68589bSMark F. Adams   /* aloc space for coarse point data (output) */
3580cbbd2e1SMark F. Adams   out_data_stride = nSelected * nSAvec;
3592fa5cd67SKarl Rupp 
3609566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(out_data_stride * nSAvec, &out_data));
36133812677SSatish Balay   for (ii = 0; ii < out_data_stride * nSAvec; ii++) out_data[ii] = PETSC_MAX_REAL;
3620cbbd2e1SMark F. Adams   *a_data_out = out_data; /* output - stride nSelected*nSAvec */
3632e68589bSMark F. Adams 
3642e68589bSMark F. Adams   /* find points and set prolongation */
365c8b0795cSMark F. Adams   minsz = 100;
3660cbbd2e1SMark F. Adams   for (mm = clid = 0; mm < nloc; mm++) {
3679566063dSJacob Faibussowitsch     PetscCall(PetscCDSizeAt(agg_llists, mm, &jj));
368e78576d6SMark F. Adams     if (jj > 0) {
3690cbbd2e1SMark F. Adams       const PetscInt lid = mm, cgid = my0crs + clid;
3700cbbd2e1SMark F. Adams       PetscInt       cids[100]; /* max bs */
3710cbbd2e1SMark F. Adams       PetscBLASInt   asz = jj, M = asz * bs, N = nSAvec, INFO;
3722e68589bSMark F. Adams       PetscBLASInt   Mdata = M + ((N - M > 0) ? N - M : 0), LDA = Mdata, LWORK = N * bs;
3732e68589bSMark F. Adams       PetscScalar   *qqc, *qqr, *TAU, *WORK;
3742e68589bSMark F. Adams       PetscInt      *fids;
37565d7b583SSatish Balay       PetscReal     *data;
376b817416eSBarry Smith 
3770cbbd2e1SMark F. Adams       /* count agg */
3780cbbd2e1SMark F. Adams       if (asz < minsz) minsz = asz;
3790cbbd2e1SMark F. Adams 
3800cbbd2e1SMark F. Adams       /* get block */
3819566063dSJacob Faibussowitsch       PetscCall(PetscMalloc5(Mdata * N, &qqc, M * N, &qqr, N, &TAU, LWORK, &WORK, M, &fids));
3822e68589bSMark F. Adams 
3832e68589bSMark F. Adams       aggID = 0;
3849566063dSJacob Faibussowitsch       PetscCall(PetscCDGetHeadPos(agg_llists, lid, &pos));
385e78576d6SMark F. Adams       while (pos) {
386e78576d6SMark F. Adams         PetscInt gid1;
3879566063dSJacob Faibussowitsch         PetscCall(PetscCDIntNdGetID(pos, &gid1));
3889566063dSJacob Faibussowitsch         PetscCall(PetscCDGetNextPos(agg_llists, lid, &pos));
389e78576d6SMark F. Adams 
3900cbbd2e1SMark F. Adams         if (gid1 >= my0 && gid1 < Iend) flid = gid1 - my0;
3910cbbd2e1SMark F. Adams         else {
3929566063dSJacob Faibussowitsch           PetscCall(PCGAMGHashTableFind(&fgid_flid, gid1, &flid));
39308401ef6SPierre Jolivet           PetscCheck(flid >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cannot find gid1 in table");
3940cbbd2e1SMark F. Adams         }
3952e68589bSMark F. Adams         /* copy in B_i matrix - column oriented */
39665d7b583SSatish Balay         data = &data_in[flid * bs];
3973d3eaba7SBarry Smith         for (ii = 0; ii < bs; ii++) {
3982e68589bSMark F. Adams           for (jj = 0; jj < N; jj++) {
3990cbbd2e1SMark F. Adams             PetscReal d                       = data[jj * data_stride + ii];
4000cbbd2e1SMark F. Adams             qqc[jj * Mdata + aggID * bs + ii] = d;
4012e68589bSMark F. Adams           }
4022e68589bSMark F. Adams         }
4032e68589bSMark F. Adams         /* set fine IDs */
4042e68589bSMark F. Adams         for (kk = 0; kk < bs; kk++) fids[aggID * bs + kk] = flid_fgid[flid] * bs + kk;
4052e68589bSMark F. Adams         aggID++;
4060cbbd2e1SMark F. Adams       }
4072e68589bSMark F. Adams 
4082e68589bSMark F. Adams       /* pad with zeros */
4092e68589bSMark F. Adams       for (ii = asz * bs; ii < Mdata; ii++) {
410ad540459SPierre Jolivet         for (jj = 0; jj < N; jj++, kk++) qqc[jj * Mdata + ii] = .0;
4112e68589bSMark F. Adams       }
4122e68589bSMark F. Adams 
4132e68589bSMark F. Adams       /* QR */
4149566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
415792fecdfSBarry Smith       PetscCallBLAS("LAPACKgeqrf", LAPACKgeqrf_(&Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
4169566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPop());
41708401ef6SPierre Jolivet       PetscCheck(INFO == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xGEQRF error");
4182e68589bSMark F. Adams       /* get R - column oriented - output B_{i+1} */
4192e68589bSMark F. Adams       {
4202e68589bSMark F. Adams         PetscReal *data = &out_data[clid * nSAvec];
4212e68589bSMark F. Adams         for (jj = 0; jj < nSAvec; jj++) {
4222e68589bSMark F. Adams           for (ii = 0; ii < nSAvec; ii++) {
42308401ef6SPierre 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);
4240cbbd2e1SMark F. Adams             if (ii <= jj) data[jj * out_data_stride + ii] = PetscRealPart(qqc[jj * Mdata + ii]);
4250cbbd2e1SMark F. Adams             else data[jj * out_data_stride + ii] = 0.;
4262e68589bSMark F. Adams           }
4272e68589bSMark F. Adams         }
4282e68589bSMark F. Adams       }
4292e68589bSMark F. Adams 
4302e68589bSMark F. Adams       /* get Q - row oriented */
431792fecdfSBarry Smith       PetscCallBLAS("LAPACKorgqr", LAPACKorgqr_(&Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
43263a3b9bcSJacob Faibussowitsch       PetscCheck(INFO == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "xORGQR error arg %" PetscBLASInt_FMT, -INFO);
4332e68589bSMark F. Adams 
4342e68589bSMark F. Adams       for (ii = 0; ii < M; ii++) {
435ad540459SPierre Jolivet         for (jj = 0; jj < N; jj++) qqr[N * ii + jj] = qqc[jj * Mdata + ii];
4362e68589bSMark F. Adams       }
4372e68589bSMark F. Adams 
4382e68589bSMark F. Adams       /* add diagonal block of P0 */
4399371c9d4SSatish Balay       for (kk = 0; kk < N; kk++) { cids[kk] = N * cgid + kk; /* global col IDs in P0 */ }
4409566063dSJacob Faibussowitsch       PetscCall(MatSetValues(a_Prol, M, fids, N, cids, qqr, INSERT_VALUES));
4419566063dSJacob Faibussowitsch       PetscCall(PetscFree5(qqc, qqr, TAU, WORK, fids));
442b43b03e9SMark F. Adams       clid++;
4430cbbd2e1SMark F. Adams     } /* coarse agg */
4440cbbd2e1SMark F. Adams   }   /* for all fine nodes */
4459566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(a_Prol, MAT_FINAL_ASSEMBLY));
4469566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(a_Prol, MAT_FINAL_ASSEMBLY));
4479566063dSJacob Faibussowitsch   PetscCall(PCGAMGHashTableDestroy(&fgid_flid));
4482e68589bSMark F. Adams   PetscFunctionReturn(0);
4492e68589bSMark F. Adams }
4502e68589bSMark F. Adams 
451*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_GAMG_AGG(PC pc, PetscViewer viewer)
452*d71ae5a4SJacob Faibussowitsch {
4535adeb434SBarry Smith   PC_MG       *mg          = (PC_MG *)pc->data;
4545adeb434SBarry Smith   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
4555adeb434SBarry Smith   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
4565adeb434SBarry Smith 
4575adeb434SBarry Smith   PetscFunctionBegin;
4589566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "      AGG specific options\n"));
459bae903cbSmarkadams4   PetscCall(PetscViewerASCIIPrintf(viewer, "        Symmetric graph %s\n", pc_gamg_agg->symmetrize_graph ? "true" : "false"));
460bae903cbSmarkadams4   PetscCall(PetscViewerASCIIPrintf(viewer, "        Number of levels to square graph %d\n", (int)pc_gamg_agg->aggressive_coarsening_levels));
461bae903cbSmarkadams4   PetscCall(PetscViewerASCIIPrintf(viewer, "        Number smoothing steps %d\n", (int)pc_gamg_agg->nsmooths));
4625adeb434SBarry Smith   PetscFunctionReturn(0);
4635adeb434SBarry Smith }
4645adeb434SBarry Smith 
4652e68589bSMark F. Adams /*
466fd1112cbSBarry Smith    PCGAMGGraph_AGG
4672e68589bSMark F. Adams 
4682e68589bSMark F. Adams   Input Parameter:
4692e68589bSMark F. Adams    . pc - this
4702e68589bSMark F. Adams    . Amat - matrix on this fine level
471bae903cbSmarkadams4 
4722e68589bSMark F. Adams   Output Parameter:
473c8b0795cSMark F. Adams    . a_Gmat -
4742e68589bSMark F. Adams */
475*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGGraph_AGG(PC pc, Mat Amat, Mat *a_Gmat)
476*d71ae5a4SJacob Faibussowitsch {
477c8b0795cSMark F. Adams   PC_MG                    *mg          = (PC_MG *)pc->data;
478c8b0795cSMark F. Adams   PC_GAMG                  *pc_gamg     = (PC_GAMG *)mg->innerctx;
479c1eae691SMark Adams   const PetscReal           vfilter     = pc_gamg->threshold[pc_gamg->current_level];
480c8b0795cSMark F. Adams   PC_GAMG_AGG              *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
48172833a62Smarkadams4   Mat                       Gmat, F = NULL;
4823b4367a7SBarry Smith   MPI_Comm                  comm;
48367744fedSMark F. Adams   PetscBool /* set,flg , */ symm;
484c8b0795cSMark F. Adams 
485c8b0795cSMark F. Adams   PetscFunctionBegin;
4869566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm));
487c8b0795cSMark F. Adams 
4889566063dSJacob Faibussowitsch   /* PetscCall(MatIsSymmetricKnown(Amat, &set, &flg)); || !(set && flg) -- this causes lot of symm calls */
489bae903cbSmarkadams4   symm = (PetscBool)(pc_gamg_agg->symmetrize_graph); /* && !pc_gamg_agg->square_graph; */
4900cbbd2e1SMark F. Adams 
491849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_GRAPH], 0, 0, 0, 0));
49272833a62Smarkadams4   PetscCall(MatCreateGraph(Amat, symm, PETSC_TRUE, &Gmat));
493849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_GRAPH], 0, 0, 0, 0));
494849bee69SMark Adams 
495849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_FILTER], 0, 0, 0, 0));
49672833a62Smarkadams4   PetscCall(MatFilter(Gmat, vfilter, &F));
49772833a62Smarkadams4   if (F) {
49872833a62Smarkadams4     PetscCall(MatDestroy(&Gmat));
49972833a62Smarkadams4     Gmat = F;
50072833a62Smarkadams4   }
501849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_FILTER], 0, 0, 0, 0));
502849bee69SMark Adams 
503e0940f08SMark F. Adams   *a_Gmat = Gmat;
504849bee69SMark Adams 
505c8b0795cSMark F. Adams   PetscFunctionReturn(0);
506c8b0795cSMark F. Adams }
507c8b0795cSMark F. Adams 
508c8b0795cSMark F. Adams /*
509bae903cbSmarkadams4    PCGAMGCoarsen_AGG - supports squaring the graph (deprecated) and new graph for
510bae903cbSmarkadams4      communication of QR data used with HEM and MISk coarsening
511c8b0795cSMark F. Adams 
512c8b0795cSMark F. Adams   Input Parameter:
513e0940f08SMark F. Adams    . a_pc - this
514bae903cbSmarkadams4 
515e0940f08SMark F. Adams   Input/Output Parameter:
516bae903cbSmarkadams4    . a_Gmat1 - graph to coarsen (in), graph off processor edges for QR gather scatter (out)
517bae903cbSmarkadams4 
518c8b0795cSMark F. Adams   Output Parameter:
5190cbbd2e1SMark F. Adams    . agg_lists - list of aggregates
520bae903cbSmarkadams4 
521c8b0795cSMark F. Adams */
522*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGCoarsen_AGG(PC a_pc, Mat *a_Gmat1, PetscCoarsenData **agg_lists)
523*d71ae5a4SJacob Faibussowitsch {
524e0940f08SMark F. Adams   PC_MG       *mg          = (PC_MG *)a_pc->data;
525c8b0795cSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
526c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
527bae903cbSmarkadams4   Mat          mat, Gmat1 = *a_Gmat1; /* aggressive graph */
5280cbbd2e1SMark F. Adams   IS           perm;
529bae903cbSmarkadams4   PetscInt     Istart, Iend, Ii, nloc, bs, nn;
530bae903cbSmarkadams4   PetscInt    *permute, *degree;
531c8b0795cSMark F. Adams   PetscBool   *bIndexSet;
532b43b03e9SMark F. Adams   MatCoarsen   crs;
5333b4367a7SBarry Smith   MPI_Comm     comm;
534aea53286SMark Adams   PetscReal    hashfact;
535e2c00dcbSBarry Smith   PetscInt     iSwapIndex;
5363b50add6SMark Adams   PetscRandom  random;
537c8b0795cSMark F. Adams 
538c8b0795cSMark F. Adams   PetscFunctionBegin;
539849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0));
5409566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Gmat1, &comm));
541bae903cbSmarkadams4   PetscCall(MatGetLocalSize(Gmat1, &nn, NULL));
5429566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(Gmat1, &bs));
54363a3b9bcSJacob Faibussowitsch   PetscCheck(bs == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "bs %" PetscInt_FMT " must be 1", bs);
544bae903cbSmarkadams4   nloc = nn / bs;
545c8b0795cSMark F. Adams 
5465cfd4bd9SMark Adams   /* get MIS aggs - randomize */
547bae903cbSmarkadams4   PetscCall(PetscMalloc2(nloc, &permute, nloc, &degree));
5489566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(nloc, &bIndexSet));
5496e09b0e3SStefano Zampini   for (Ii = 0; Ii < nloc; Ii++) permute[Ii] = Ii;
5509566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &random));
5519566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Gmat1, &Istart, &Iend));
552c8b0795cSMark F. Adams   for (Ii = 0; Ii < nloc; Ii++) {
553bae903cbSmarkadams4     PetscInt nc;
554bae903cbSmarkadams4     PetscCall(MatGetRow(Gmat1, Istart + Ii, &nc, NULL, NULL));
555bae903cbSmarkadams4     degree[Ii] = nc;
556bae903cbSmarkadams4     PetscCall(MatRestoreRow(Gmat1, Istart + Ii, &nc, NULL, NULL));
557bae903cbSmarkadams4   }
558bae903cbSmarkadams4   for (Ii = 0; Ii < nloc; Ii++) {
5599566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValueReal(random, &hashfact));
560aea53286SMark Adams     iSwapIndex = (PetscInt)(hashfact * nloc) % nloc;
561c8b0795cSMark F. Adams     if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) {
562c8b0795cSMark F. Adams       PetscInt iTemp        = permute[iSwapIndex];
563c8b0795cSMark F. Adams       permute[iSwapIndex]   = permute[Ii];
564c8b0795cSMark F. Adams       permute[Ii]           = iTemp;
565bae903cbSmarkadams4       iTemp                 = degree[iSwapIndex];
566bae903cbSmarkadams4       degree[iSwapIndex]    = degree[Ii];
567bae903cbSmarkadams4       degree[Ii]            = iTemp;
568c8b0795cSMark F. Adams       bIndexSet[iSwapIndex] = PETSC_TRUE;
569c8b0795cSMark F. Adams     }
570c8b0795cSMark F. Adams   }
571bae903cbSmarkadams4   // create minimum degree ordering
572bae903cbSmarkadams4   PetscCall(PetscSortIntWithArray(nloc, degree, permute));
573bae903cbSmarkadams4 
5749566063dSJacob Faibussowitsch   PetscCall(PetscFree(bIndexSet));
5759566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&random));
5769566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm));
577849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_MIS], 0, 0, 0, 0));
5789566063dSJacob Faibussowitsch   PetscCall(MatCoarsenCreate(comm, &crs));
5799566063dSJacob Faibussowitsch   PetscCall(MatCoarsenSetFromOptions(crs));
5809566063dSJacob Faibussowitsch   PetscCall(MatCoarsenSetGreedyOrdering(crs, perm));
581bae903cbSmarkadams4   PetscCall(MatCoarsenSetAdjacency(crs, Gmat1));
5829566063dSJacob Faibussowitsch   PetscCall(MatCoarsenSetStrictAggs(crs, PETSC_TRUE));
583bae903cbSmarkadams4   if (pc_gamg->current_level < pc_gamg_agg->aggressive_coarsening_levels) PetscCall(MatCoarsenMISKSetDistance(crs, 2)); // hardwire to MIS-2
584bae903cbSmarkadams4   else PetscCall(MatCoarsenMISKSetDistance(crs, 1));                                                                    // MIS
5859566063dSJacob Faibussowitsch   PetscCall(MatCoarsenApply(crs));
586bae903cbSmarkadams4   PetscCall(MatCoarsenViewFromOptions(crs, NULL, "-mat_coarsen_view"));
5879566063dSJacob Faibussowitsch   PetscCall(MatCoarsenGetData(crs, agg_lists)); /* output */
5889566063dSJacob Faibussowitsch   PetscCall(MatCoarsenDestroy(&crs));
589b43b03e9SMark F. Adams 
5909566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
591bae903cbSmarkadams4   PetscCall(PetscFree2(permute, degree));
592849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_MIS], 0, 0, 0, 0));
5939f3f12c8SMark F. Adams 
594bae903cbSmarkadams4   {
595bae903cbSmarkadams4     PetscCoarsenData *llist = *agg_lists;
5969ab59c8bSMark Adams     /* see if we have a matrix that takes precedence (returned from MatCoarsenApply) */
5979566063dSJacob Faibussowitsch     PetscCall(PetscCDGetMat(llist, &mat));
5980cbbd2e1SMark F. Adams     if (mat) {
5999566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&Gmat1));
6000cbbd2e1SMark F. Adams       *a_Gmat1 = mat; /* output */
6010cbbd2e1SMark F. Adams     }
6020cbbd2e1SMark F. Adams   }
603849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_COARSEN], 0, 0, 0, 0));
604c8b0795cSMark F. Adams   PetscFunctionReturn(0);
605c8b0795cSMark F. Adams }
606c8b0795cSMark F. Adams 
607c8b0795cSMark F. Adams /*
6080cbbd2e1SMark F. Adams  PCGAMGProlongator_AGG
609c8b0795cSMark F. Adams 
610c8b0795cSMark F. Adams  Input Parameter:
611c8b0795cSMark F. Adams  . pc - this
612c8b0795cSMark F. Adams  . Amat - matrix on this fine level
613c8b0795cSMark F. Adams  . Graph - used to get ghost data for nodes in
6140cbbd2e1SMark F. Adams  . agg_lists - list of aggregates
615c8b0795cSMark F. Adams  Output Parameter:
616c8b0795cSMark F. Adams  . a_P_out - prolongation operator to the next level
617c8b0795cSMark F. Adams  */
618*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGProlongator_AGG(PC pc, Mat Amat, Mat Gmat, PetscCoarsenData *agg_lists, Mat *a_P_out)
619*d71ae5a4SJacob Faibussowitsch {
6202e68589bSMark F. Adams   PC_MG         *mg      = (PC_MG *)pc->data;
6212e68589bSMark F. Adams   PC_GAMG       *pc_gamg = (PC_GAMG *)mg->innerctx;
6224a2f8832SBarry Smith   const PetscInt col_bs  = pc_gamg->data_cell_cols;
623c8b0795cSMark F. Adams   PetscInt       Istart, Iend, nloc, ii, jj, kk, my0, nLocalSelected, bs;
624c8b0795cSMark F. Adams   Mat            Prol;
625d5d25401SBarry Smith   PetscMPIInt    size;
6263b4367a7SBarry Smith   MPI_Comm       comm;
627c8b0795cSMark F. Adams   PetscReal     *data_w_ghost;
628c8b0795cSMark F. Adams   PetscInt       myCrs0, nbnodes = 0, *flid_fgid;
629d9558ea9SBarry Smith   MatType        mtype;
6302e68589bSMark F. Adams 
6312e68589bSMark F. Adams   PetscFunctionBegin;
6329566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm));
63308401ef6SPierre Jolivet   PetscCheck(col_bs >= 1, comm, PETSC_ERR_PLIB, "Column bs cannot be less than 1");
634849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0));
6359566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
6369566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Amat, &Istart, &Iend));
6379566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(Amat, &bs));
6389371c9d4SSatish Balay   nloc = (Iend - Istart) / bs;
6399371c9d4SSatish Balay   my0  = Istart / bs;
64063a3b9bcSJacob 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);
6412e68589bSMark F. Adams 
6422e68589bSMark F. Adams   /* get 'nLocalSelected' */
6430cbbd2e1SMark F. Adams   for (ii = 0, nLocalSelected = 0; ii < nloc; ii++) {
644e78576d6SMark F. Adams     PetscBool ise;
6450cbbd2e1SMark F. Adams     /* filter out singletons 0 or 1? */
6469566063dSJacob Faibussowitsch     PetscCall(PetscCDEmptyAt(agg_lists, ii, &ise));
647e78576d6SMark F. Adams     if (!ise) nLocalSelected++;
6482e68589bSMark F. Adams   }
6492e68589bSMark F. Adams 
6502e68589bSMark F. Adams   /* create prolongator, create P matrix */
6519566063dSJacob Faibussowitsch   PetscCall(MatGetType(Amat, &mtype));
6529566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &Prol));
6539566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Prol, nloc * bs, nLocalSelected * col_bs, PETSC_DETERMINE, PETSC_DETERMINE));
6549566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSizes(Prol, bs, col_bs));
6559566063dSJacob Faibussowitsch   PetscCall(MatSetType(Prol, mtype));
6569566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(Prol, col_bs, NULL));
6579566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(Prol, col_bs, NULL, col_bs, NULL));
6582e68589bSMark F. Adams 
6592e68589bSMark F. Adams   /* can get all points "removed" */
6609566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Prol, &kk, &ii));
6617f66b68fSMark Adams   if (!ii) {
66263a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(pc, "%s: No selected points on coarse grid\n", ((PetscObject)pc)->prefix));
6639566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Prol));
6640298fd71SBarry Smith     *a_P_out = NULL; /* out */
665849bee69SMark Adams     PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0));
6662e68589bSMark F. Adams     PetscFunctionReturn(0);
6672e68589bSMark F. Adams   }
66863a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(pc, "%s: New grid %" PetscInt_FMT " nodes\n", ((PetscObject)pc)->prefix, ii / col_bs));
6699566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRangeColumn(Prol, &myCrs0, &kk));
6700cbbd2e1SMark F. Adams 
67163a3b9bcSJacob 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);
672c8b0795cSMark F. Adams   myCrs0 = myCrs0 / col_bs;
67363a3b9bcSJacob 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);
6742e68589bSMark F. Adams 
6752e68589bSMark F. Adams   /* create global vector of data in 'data_w_ghost' */
676849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROLA], 0, 0, 0, 0));
677bae903cbSmarkadams4   if (size > 1) { /* get ghost null space data */
6782e68589bSMark F. Adams     PetscReal *tmp_gdata, *tmp_ldata, *tp2;
6799566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nloc, &tmp_ldata));
6804a2f8832SBarry Smith     for (jj = 0; jj < col_bs; jj++) {
681c8b0795cSMark F. Adams       for (kk = 0; kk < bs; kk++) {
682a2f3521dSMark F. Adams         PetscInt         ii, stride;
683c8b0795cSMark F. Adams         const PetscReal *tp = pc_gamg->data + jj * bs * nloc + kk;
6842fa5cd67SKarl Rupp         for (ii = 0; ii < nloc; ii++, tp += bs) tmp_ldata[ii] = *tp;
6852fa5cd67SKarl Rupp 
6869566063dSJacob Faibussowitsch         PetscCall(PCGAMGGetDataWithGhosts(Gmat, 1, tmp_ldata, &stride, &tmp_gdata));
687a2f3521dSMark F. Adams 
688bae903cbSmarkadams4         if (!jj && !kk) { /* now I know how many total nodes - allocate TODO: move below and do in one 'col_bs' call */
6899566063dSJacob Faibussowitsch           PetscCall(PetscMalloc1(stride * bs * col_bs, &data_w_ghost));
690a2f3521dSMark F. Adams           nbnodes = bs * stride;
6912e68589bSMark F. Adams         }
692a2f3521dSMark F. Adams         tp2 = data_w_ghost + jj * bs * stride + kk;
6932fa5cd67SKarl Rupp         for (ii = 0; ii < stride; ii++, tp2 += bs) *tp2 = tmp_gdata[ii];
6949566063dSJacob Faibussowitsch         PetscCall(PetscFree(tmp_gdata));
6952e68589bSMark F. Adams       }
6962e68589bSMark F. Adams     }
6979566063dSJacob Faibussowitsch     PetscCall(PetscFree(tmp_ldata));
698806fa848SBarry Smith   } else {
699c8b0795cSMark F. Adams     nbnodes      = bs * nloc;
700c8b0795cSMark F. Adams     data_w_ghost = (PetscReal *)pc_gamg->data;
7012e68589bSMark F. Adams   }
7022e68589bSMark F. Adams 
703bae903cbSmarkadams4   /* get 'flid_fgid' TODO - move up to get 'stride' and do get null space data above in one step (jj loop) */
704c5df96a5SBarry Smith   if (size > 1) {
7052e68589bSMark F. Adams     PetscReal *fid_glid_loc, *fiddata;
706a2f3521dSMark F. Adams     PetscInt   stride;
7072e68589bSMark F. Adams 
7089566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nloc, &fid_glid_loc));
7092e68589bSMark F. Adams     for (kk = 0; kk < nloc; kk++) fid_glid_loc[kk] = (PetscReal)(my0 + kk);
7109566063dSJacob Faibussowitsch     PetscCall(PCGAMGGetDataWithGhosts(Gmat, 1, fid_glid_loc, &stride, &fiddata));
711bae903cbSmarkadams4     PetscCall(PetscMalloc1(stride, &flid_fgid)); /* copy real data to in */
712a2f3521dSMark F. Adams     for (kk = 0; kk < stride; kk++) flid_fgid[kk] = (PetscInt)fiddata[kk];
7139566063dSJacob Faibussowitsch     PetscCall(PetscFree(fiddata));
714a2f3521dSMark F. Adams 
71563a3b9bcSJacob Faibussowitsch     PetscCheck(stride == nbnodes / bs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "stride %" PetscInt_FMT " != nbnodes %" PetscInt_FMT "/bs %" PetscInt_FMT, stride, nbnodes, bs);
7169566063dSJacob Faibussowitsch     PetscCall(PetscFree(fid_glid_loc));
717806fa848SBarry Smith   } else {
7189566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nloc, &flid_fgid));
7192e68589bSMark F. Adams     for (kk = 0; kk < nloc; kk++) flid_fgid[kk] = my0 + kk;
7202e68589bSMark F. Adams   }
721849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROLA], 0, 0, 0, 0));
722bae903cbSmarkadams4   /* get P0 */
723849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROLB], 0, 0, 0, 0));
724c8b0795cSMark F. Adams   {
7250298fd71SBarry Smith     PetscReal *data_out = NULL;
7269566063dSJacob Faibussowitsch     PetscCall(formProl0(agg_lists, bs, col_bs, myCrs0, nbnodes, data_w_ghost, flid_fgid, &data_out, Prol));
7279566063dSJacob Faibussowitsch     PetscCall(PetscFree(pc_gamg->data));
7282fa5cd67SKarl Rupp 
729c8b0795cSMark F. Adams     pc_gamg->data           = data_out;
7304a2f8832SBarry Smith     pc_gamg->data_cell_rows = col_bs;
7314a2f8832SBarry Smith     pc_gamg->data_sz        = col_bs * col_bs * nLocalSelected;
732c8b0795cSMark F. Adams   }
733849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROLB], 0, 0, 0, 0));
73448a46eb9SPierre Jolivet   if (size > 1) PetscCall(PetscFree(data_w_ghost));
7359566063dSJacob Faibussowitsch   PetscCall(PetscFree(flid_fgid));
7362e68589bSMark F. Adams 
737c8b0795cSMark F. Adams   *a_P_out = Prol; /* out */
738ed4e82eaSPeter Brune 
739849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROL], 0, 0, 0, 0));
740c8b0795cSMark F. Adams   PetscFunctionReturn(0);
741c8b0795cSMark F. Adams }
742c8b0795cSMark F. Adams 
743c8b0795cSMark F. Adams /*
744fd1112cbSBarry Smith    PCGAMGOptProlongator_AGG
745c8b0795cSMark F. Adams 
746c8b0795cSMark F. Adams   Input Parameter:
747c8b0795cSMark F. Adams    . pc - this
748c8b0795cSMark F. Adams    . Amat - matrix on this fine level
749c8b0795cSMark F. Adams  In/Output Parameter:
7502a808120SBarry Smith    . a_P - prolongation operator to the next level
751c8b0795cSMark F. Adams */
752*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGAMGOptProlongator_AGG(PC pc, Mat Amat, Mat *a_P)
753*d71ae5a4SJacob Faibussowitsch {
754c8b0795cSMark F. Adams   PC_MG       *mg          = (PC_MG *)pc->data;
755c8b0795cSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG *)mg->innerctx;
756c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG *)pc_gamg->subctx;
757c8b0795cSMark F. Adams   PetscInt     jj;
758c8b0795cSMark F. Adams   Mat          Prol = *a_P;
7593b4367a7SBarry Smith   MPI_Comm     comm;
7602a808120SBarry Smith   KSP          eksp;
7612a808120SBarry Smith   Vec          bb, xx;
7622a808120SBarry Smith   PC           epc;
7632a808120SBarry Smith   PetscReal    alpha, emax, emin;
764c8b0795cSMark F. Adams 
765c8b0795cSMark F. Adams   PetscFunctionBegin;
7669566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm));
767849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_OPT], 0, 0, 0, 0));
768c8b0795cSMark F. Adams 
769d5d25401SBarry Smith   /* compute maximum singular value of operator to be used in smoother */
7702a808120SBarry Smith   if (0 < pc_gamg_agg->nsmooths) {
77118c3aa7eSMark     /* get eigen estimates */
77218c3aa7eSMark     if (pc_gamg->emax > 0) {
77318c3aa7eSMark       emin = pc_gamg->emin;
77418c3aa7eSMark       emax = pc_gamg->emax;
77518c3aa7eSMark     } else {
7760ed2132dSStefano Zampini       const char *prefix;
7770ed2132dSStefano Zampini 
7789566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(Amat, &bb, NULL));
7799566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(Amat, &xx, NULL));
7809566063dSJacob Faibussowitsch       PetscCall(KSPSetNoisy_Private(bb));
781e696ed0bSMark F. Adams 
7829566063dSJacob Faibussowitsch       PetscCall(KSPCreate(comm, &eksp));
7839566063dSJacob Faibussowitsch       PetscCall(PCGetOptionsPrefix(pc, &prefix));
7849566063dSJacob Faibussowitsch       PetscCall(KSPSetOptionsPrefix(eksp, prefix));
7859566063dSJacob Faibussowitsch       PetscCall(KSPAppendOptionsPrefix(eksp, "pc_gamg_esteig_"));
78673f7197eSJed Brown       {
787b94d7dedSBarry Smith         PetscBool isset, sflg;
788b94d7dedSBarry Smith         PetscCall(MatIsSPDKnown(Amat, &isset, &sflg));
789b94d7dedSBarry Smith         if (isset && sflg) PetscCall(KSPSetType(eksp, KSPCG));
790d24ecf33SMark       }
7919566063dSJacob Faibussowitsch       PetscCall(KSPSetErrorIfNotConverged(eksp, pc->erroriffailure));
7929566063dSJacob Faibussowitsch       PetscCall(KSPSetNormType(eksp, KSP_NORM_NONE));
793c2436cacSMark F. Adams 
7949566063dSJacob Faibussowitsch       PetscCall(KSPSetInitialGuessNonzero(eksp, PETSC_FALSE));
7959566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(eksp, Amat, Amat));
7962e68589bSMark F. Adams 
7979566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(eksp, &epc));
7989566063dSJacob Faibussowitsch       PetscCall(PCSetType(epc, PCJACOBI)); /* smoother in smoothed agg. */
799c2436cacSMark F. Adams 
8009566063dSJacob 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
801074aec55SMark Adams 
8029566063dSJacob Faibussowitsch       PetscCall(KSPSetFromOptions(eksp));
8039566063dSJacob Faibussowitsch       PetscCall(KSPSetComputeSingularValues(eksp, PETSC_TRUE));
8049566063dSJacob Faibussowitsch       PetscCall(KSPSolve(eksp, bb, xx));
8059566063dSJacob Faibussowitsch       PetscCall(KSPCheckSolve(eksp, pc, xx));
8062e68589bSMark F. Adams 
8079566063dSJacob Faibussowitsch       PetscCall(KSPComputeExtremeSingularValues(eksp, &emax, &emin));
8089566063dSJacob Faibussowitsch       PetscCall(PetscInfo(pc, "%s: Smooth P0: max eigen=%e min=%e PC=%s\n", ((PetscObject)pc)->prefix, (double)emax, (double)emin, PCJACOBI));
8099566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&xx));
8109566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&bb));
8119566063dSJacob Faibussowitsch       PetscCall(KSPDestroy(&eksp));
8122e68589bSMark F. Adams     }
81318c3aa7eSMark     if (pc_gamg->use_sa_esteig) {
81418c3aa7eSMark       mg->min_eigen_DinvA[pc_gamg->current_level] = emin;
81518c3aa7eSMark       mg->max_eigen_DinvA[pc_gamg->current_level] = emax;
81663a3b9bcSJacob 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));
81718c3aa7eSMark     } else {
81818c3aa7eSMark       mg->min_eigen_DinvA[pc_gamg->current_level] = 0;
81918c3aa7eSMark       mg->max_eigen_DinvA[pc_gamg->current_level] = 0;
82018c3aa7eSMark     }
82118c3aa7eSMark   } else {
82218c3aa7eSMark     mg->min_eigen_DinvA[pc_gamg->current_level] = 0;
82318c3aa7eSMark     mg->max_eigen_DinvA[pc_gamg->current_level] = 0;
82418c3aa7eSMark   }
8252e68589bSMark F. Adams 
8262a808120SBarry Smith   /* smooth P0 */
8272a808120SBarry Smith   for (jj = 0; jj < pc_gamg_agg->nsmooths; jj++) {
8282a808120SBarry Smith     Mat tMat;
8292a808120SBarry Smith     Vec diag;
8302a808120SBarry Smith 
831849bee69SMark Adams     PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_OPTSM], 0, 0, 0, 0));
8322a808120SBarry Smith 
8332e68589bSMark F. Adams     /* smooth P1 := (I - omega/lam D^{-1}A)P0 */
8349566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2], 0, 0, 0, 0));
8359566063dSJacob Faibussowitsch     PetscCall(MatMatMult(Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat));
8369566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2], 0, 0, 0, 0));
8379566063dSJacob Faibussowitsch     PetscCall(MatProductClear(tMat));
8389566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(Amat, &diag, NULL));
8399566063dSJacob Faibussowitsch     PetscCall(MatGetDiagonal(Amat, diag)); /* effectively PCJACOBI */
8409566063dSJacob Faibussowitsch     PetscCall(VecReciprocal(diag));
8419566063dSJacob Faibussowitsch     PetscCall(MatDiagonalScale(tMat, diag, NULL));
8429566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&diag));
843b4da3a1bSStefano Zampini 
844d5d25401SBarry Smith     /* TODO: Set a PCFailedReason and exit the building of the AMG preconditioner */
84508401ef6SPierre Jolivet     PetscCheck(emax != 0.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_PLIB, "Computed maximum singular value as zero");
846d5d25401SBarry Smith     /* TODO: Document the 1.4 and don't hardwire it in this routine */
847b4ec6429SMark F. Adams     alpha = -1.4 / emax;
848b4da3a1bSStefano Zampini 
8499566063dSJacob Faibussowitsch     PetscCall(MatAYPX(tMat, alpha, Prol, SUBSET_NONZERO_PATTERN));
8509566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Prol));
8512e68589bSMark F. Adams     Prol = tMat;
852849bee69SMark Adams     PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_OPTSM], 0, 0, 0, 0));
8532e68589bSMark F. Adams   }
854849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_OPT], 0, 0, 0, 0));
855c8b0795cSMark F. Adams   *a_P = Prol;
856c8b0795cSMark F. Adams   PetscFunctionReturn(0);
8572e68589bSMark F. Adams }
8582e68589bSMark F. Adams 
859c8b0795cSMark F. Adams /*
860c8b0795cSMark F. Adams    PCCreateGAMG_AGG
8612e68589bSMark F. Adams 
862c8b0795cSMark F. Adams   Input Parameter:
863c8b0795cSMark F. Adams    . pc -
864c8b0795cSMark F. Adams */
865*d71ae5a4SJacob Faibussowitsch PetscErrorCode PCCreateGAMG_AGG(PC pc)
866*d71ae5a4SJacob Faibussowitsch {
867c8b0795cSMark F. Adams   PC_MG       *mg      = (PC_MG *)pc->data;
868c8b0795cSMark F. Adams   PC_GAMG     *pc_gamg = (PC_GAMG *)mg->innerctx;
869c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg;
8702e68589bSMark F. Adams 
871c8b0795cSMark F. Adams   PetscFunctionBegin;
872c8b0795cSMark F. Adams   /* create sub context for SA */
8734dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&pc_gamg_agg));
874c8b0795cSMark F. Adams   pc_gamg->subctx = pc_gamg_agg;
875c8b0795cSMark F. Adams 
8761ab5ffc9SJed Brown   pc_gamg->ops->setfromoptions = PCSetFromOptions_GAMG_AGG;
8779b8ffb57SJed Brown   pc_gamg->ops->destroy        = PCDestroy_GAMG_AGG;
878c8b0795cSMark F. Adams   /* reset does not do anything; setup not virtual */
879c8b0795cSMark F. Adams 
880c8b0795cSMark F. Adams   /* set internal function pointers */
881fd1112cbSBarry Smith   pc_gamg->ops->graph             = PCGAMGGraph_AGG;
8821ab5ffc9SJed Brown   pc_gamg->ops->coarsen           = PCGAMGCoarsen_AGG;
8831ab5ffc9SJed Brown   pc_gamg->ops->prolongator       = PCGAMGProlongator_AGG;
884fd1112cbSBarry Smith   pc_gamg->ops->optprolongator    = PCGAMGOptProlongator_AGG;
8851ab5ffc9SJed Brown   pc_gamg->ops->createdefaultdata = PCSetData_AGG;
8865adeb434SBarry Smith   pc_gamg->ops->view              = PCView_GAMG_AGG;
887c8b0795cSMark F. Adams 
888bae903cbSmarkadams4   pc_gamg_agg->aggressive_coarsening_levels = 0;
889bae903cbSmarkadams4   pc_gamg_agg->symmetrize_graph             = PETSC_FALSE;
890cab9ed1eSBarry Smith   pc_gamg_agg->nsmooths                     = 1;
891cab9ed1eSBarry Smith 
8929566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetNSmooths_C", PCGAMGSetNSmooths_AGG));
893bae903cbSmarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetSymmetrizeGraph_C", PCGAMGSetSymmetrizeGraph_AGG));
894bae903cbSmarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGAMGSetAggressiveLevels_C", PCGAMGSetAggressiveLevels_AGG));
8959566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", PCSetCoordinates_AGG));
8962e68589bSMark F. Adams   PetscFunctionReturn(0);
8972e68589bSMark F. Adams }
898