xref: /petsc/src/ksp/pc/impls/gamg/agg.c (revision b94d7ded0a05f1bbd5e48daa6f92b28259c75b44)
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 @*/
312e68589bSMark F. Adams PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n)
322e68589bSMark F. Adams {
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 
40e0877f53SBarry Smith static PetscErrorCode PCGAMGSetNSmooths_AGG(PC pc, PetscInt n)
412e68589bSMark F. Adams {
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 
54d5d25401SBarry Smith    Logically Collective on PC
55c8b0795cSMark F. Adams 
56c8b0795cSMark F. Adams    Input Parameters:
57cab9ed1eSBarry Smith +  pc - the preconditioner context
58a2b725a8SWilliam Gropp -  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 
65bae903cbSmarkadams4 .seealso: `PCGAMGSetAggressiveLevels()`
66c8b0795cSMark F. Adams @*/
67bae903cbSmarkadams4 PetscErrorCode PCGAMGSetSymmetrizeGraph(PC pc, PetscBool n)
68c8b0795cSMark F. Adams {
69c8b0795cSMark F. Adams   PetscFunctionBegin;
70c8b0795cSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
71d5d25401SBarry Smith   PetscValidLogicalCollectiveBool(pc,n,2);
72bae903cbSmarkadams4   PetscTryMethod(pc,"PCGAMGSetSymmetrizeGraph_C",(PC,PetscBool),(pc,n));
73c8b0795cSMark F. Adams   PetscFunctionReturn(0);
74c8b0795cSMark F. Adams }
75c8b0795cSMark F. Adams 
76bae903cbSmarkadams4 static PetscErrorCode PCGAMGSetSymmetrizeGraph_AGG(PC pc, PetscBool n)
77c8b0795cSMark F. Adams {
78c8b0795cSMark F. Adams   PC_MG       *mg          = (PC_MG*)pc->data;
79c8b0795cSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG*)mg->innerctx;
80c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
81c8b0795cSMark F. Adams 
82c8b0795cSMark F. Adams   PetscFunctionBegin;
83bae903cbSmarkadams4   pc_gamg_agg->symmetrize_graph = n;
842e68589bSMark F. Adams   PetscFunctionReturn(0);
852e68589bSMark F. Adams }
862e68589bSMark F. Adams 
87ef4ad70eSMark F. Adams /*@
88bae903cbSmarkadams4    PCGAMGSetAggressiveLevels -  Aggressive coarsening on first n levels
89ef4ad70eSMark F. Adams 
90d5d25401SBarry Smith    Logically Collective on PC
91ef4ad70eSMark F. Adams 
92ef4ad70eSMark F. Adams    Input Parameters:
93cab9ed1eSBarry Smith +  pc - the preconditioner context
94d5d25401SBarry Smith -  n - 0, 1 or more
95ef4ad70eSMark F. Adams 
96ef4ad70eSMark F. Adams    Options Database Key:
97bae903cbSmarkadams4 .  -pc_gamg_aggressive_coarsening <n,default = 1> - Number of levels to square the graph on before aggregating it
98af3c827dSMark Adams 
99ef4ad70eSMark F. Adams    Level: intermediate
100ef4ad70eSMark F. Adams 
101bae903cbSmarkadams4 .seealso: `PCGAMGSetSymmetrizeGraph()`, `PCGAMGSetThreshold()`
102ef4ad70eSMark F. Adams @*/
103bae903cbSmarkadams4 PetscErrorCode PCGAMGSetAggressiveLevels(PC pc, PetscInt n)
104ef4ad70eSMark F. Adams {
105ef4ad70eSMark F. Adams   PetscFunctionBegin;
106ef4ad70eSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
107d5d25401SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
108bae903cbSmarkadams4   PetscTryMethod(pc,"PCGAMGSetAggressiveLevels_C",(PC,PetscInt),(pc,n));
109ef4ad70eSMark F. Adams   PetscFunctionReturn(0);
110ef4ad70eSMark F. Adams }
111ef4ad70eSMark F. Adams 
112bae903cbSmarkadams4 static PetscErrorCode PCGAMGSetAggressiveLevels_AGG(PC pc, PetscInt n)
113ef4ad70eSMark F. Adams {
114ef4ad70eSMark F. Adams   PC_MG       *mg          = (PC_MG*)pc->data;
115ef4ad70eSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG*)mg->innerctx;
116ef4ad70eSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
117ef4ad70eSMark F. Adams 
118ef4ad70eSMark F. Adams   PetscFunctionBegin;
119bae903cbSmarkadams4   pc_gamg_agg->aggressive_coarsening_levels = n;
120ef4ad70eSMark F. Adams   PetscFunctionReturn(0);
121ef4ad70eSMark F. Adams }
122ef4ad70eSMark F. Adams 
123e0877f53SBarry Smith static PetscErrorCode PCSetFromOptions_GAMG_AGG(PetscOptionItems *PetscOptionsObject,PC pc)
1242e68589bSMark F. Adams {
1252e68589bSMark F. Adams   PC_MG          *mg          = (PC_MG*)pc->data;
1262e68589bSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
127c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
1282e68589bSMark F. Adams 
1292e68589bSMark F. Adams   PetscFunctionBegin;
130d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject,"GAMG-AGG options");
1312e68589bSMark F. Adams   {
132bae903cbSmarkadams4     PetscBool flg;
1339566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-pc_gamg_agg_nsmooths","smoothing steps for smoothed aggregation, usually 1","PCGAMGSetNSmooths",pc_gamg_agg->nsmooths,&pc_gamg_agg->nsmooths,NULL));
134bae903cbSmarkadams4     PetscCall(PetscOptionsBool("-pc_gamg_symmetrize_graph","Set for asymmetric matrices","PCGAMGSetSymmetrizeGraph",pc_gamg_agg->symmetrize_graph,&pc_gamg_agg->symmetrize_graph,NULL));
135bae903cbSmarkadams4     pc_gamg_agg->aggressive_coarsening_levels = 1;
136bae903cbSmarkadams4     PetscCall(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));
137bae903cbSmarkadams4     if (!flg) {
138bae903cbSmarkadams4       PetscCall(PetscOptionsInt("-pc_gamg_aggressive_coarsening","Number of aggressive coarsening (MIS-2) levels from finest","PCGAMGSetAggressiveLevels",pc_gamg_agg->aggressive_coarsening_levels,&pc_gamg_agg->aggressive_coarsening_levels,NULL));
139bae903cbSmarkadams4     } else {
140bae903cbSmarkadams4       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));
141bae903cbSmarkadams4       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));
142bae903cbSmarkadams4     }
1432e68589bSMark F. Adams   }
144d0609cedSBarry Smith   PetscOptionsHeadEnd();
1452e68589bSMark F. Adams   PetscFunctionReturn(0);
1462e68589bSMark F. Adams }
1472e68589bSMark F. Adams 
1482e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
149e0877f53SBarry Smith static PetscErrorCode PCDestroy_GAMG_AGG(PC pc)
1502e68589bSMark F. Adams {
1512e68589bSMark F. Adams   PC_MG          *mg          = (PC_MG*)pc->data;
1522e68589bSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
1532e68589bSMark F. Adams 
1542e68589bSMark F. Adams   PetscFunctionBegin;
1559566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg->subctx));
1562e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNSmooths_C",NULL));
157bae903cbSmarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetSymmetrizeGraph_C",NULL));
158bae903cbSmarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetAggressiveLevels_C",NULL));
159bae903cbSmarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",NULL));
1602e68589bSMark F. Adams   PetscFunctionReturn(0);
1612e68589bSMark F. Adams }
1622e68589bSMark F. Adams 
1632e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
1642e68589bSMark F. Adams /*
1652e68589bSMark F. Adams    PCSetCoordinates_AGG
166302f38e8SMark F. Adams      - collective
1672e68589bSMark F. Adams 
1682e68589bSMark F. Adams    Input Parameter:
1692e68589bSMark F. Adams    . pc - the preconditioner context
170a2f3521dSMark F. Adams    . ndm - dimesion of data (used for dof/vertex for Stokes)
171302f38e8SMark F. Adams    . a_nloc - number of vertices local
172302f38e8SMark F. Adams    . coords - [a_nloc][ndm] - interleaved coordinate data: {x_0, y_0, z_0, x_1, y_1, ...}
1732e68589bSMark F. Adams */
1741e6b0712SBarry Smith 
1751e6b0712SBarry Smith static PetscErrorCode PCSetCoordinates_AGG(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords)
1762e68589bSMark F. Adams {
1772e68589bSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
1782e68589bSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
17969344418SMark F. Adams   PetscInt       arrsz,kk,ii,jj,nloc,ndatarows,ndf;
180a2f3521dSMark F. Adams   Mat            mat = pc->pmat;
1812e68589bSMark F. Adams 
1822e68589bSMark F. Adams   PetscFunctionBegin;
183a2f3521dSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
184a2f3521dSMark F. Adams   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
185302f38e8SMark F. Adams   nloc = a_nloc;
1862e68589bSMark F. Adams 
1872e68589bSMark F. Adams   /* SA: null space vectors */
1889566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(mat, &ndf)); /* this does not work for Stokes */
18969344418SMark F. Adams   if (coords && ndf==1) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */
190a2f3521dSMark F. Adams   else if (coords) {
19163a3b9bcSJacob Faibussowitsch     PetscCheck(ndm <= ndf,PETSC_COMM_SELF,PETSC_ERR_PLIB,"degrees of motion %" PetscInt_FMT " > block size %" PetscInt_FMT,ndm,ndf);
19269344418SMark F. Adams     pc_gamg->data_cell_cols = (ndm==2 ? 3 : 6); /* displacement elasticity */
19369344418SMark F. Adams     if (ndm != ndf) {
19463a3b9bcSJacob Faibussowitsch       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);
195a2f3521dSMark F. Adams     }
196806fa848SBarry Smith   } else pc_gamg->data_cell_cols = ndf; /* no data, force SA with constant null space vectors */
19771959b99SBarry Smith   pc_gamg->data_cell_rows = ndatarows = ndf;
19863a3b9bcSJacob 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);
199c8b0795cSMark F. Adams   arrsz = nloc*pc_gamg->data_cell_rows*pc_gamg->data_cell_cols;
2002e68589bSMark F. Adams 
2017f66b68fSMark Adams   if (!pc_gamg->data || (pc_gamg->data_sz != arrsz)) {
2029566063dSJacob Faibussowitsch     PetscCall(PetscFree(pc_gamg->data));
2039566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(arrsz+1, &pc_gamg->data));
2042e68589bSMark F. Adams   }
2052e68589bSMark F. Adams   /* copy data in - column oriented */
2062e68589bSMark F. Adams   for (kk=0; kk<nloc; kk++) {
20769344418SMark F. Adams     const PetscInt M     = nloc*pc_gamg->data_cell_rows; /* stride into data */
20869344418SMark F. Adams     PetscReal      *data = &pc_gamg->data[kk*ndatarows]; /* start of cell */
209c8b0795cSMark F. Adams     if (pc_gamg->data_cell_cols==1) *data = 1.0;
2102e68589bSMark F. Adams     else {
21169344418SMark F. Adams       /* translational modes */
2122fa5cd67SKarl Rupp       for (ii=0;ii<ndatarows;ii++) {
2132fa5cd67SKarl Rupp         for (jj=0;jj<ndatarows;jj++) {
21469344418SMark F. Adams           if (ii==jj)data[ii*M + jj] = 1.0;
2152e68589bSMark F. Adams           else data[ii*M + jj] = 0.0;
2162fa5cd67SKarl Rupp         }
2172fa5cd67SKarl Rupp       }
21869344418SMark F. Adams 
21969344418SMark F. Adams       /* rotational modes */
2202e68589bSMark F. Adams       if (coords) {
22169344418SMark F. Adams         if (ndm == 2) {
2222e68589bSMark F. Adams           data   += 2*M;
2232e68589bSMark F. Adams           data[0] = -coords[2*kk+1];
2242e68589bSMark F. Adams           data[1] =  coords[2*kk];
225806fa848SBarry Smith         } else {
2262e68589bSMark F. Adams           data   += 3*M;
2272e68589bSMark F. Adams           data[0] = 0.0;             data[M+0] =  coords[3*kk+2]; data[2*M+0] = -coords[3*kk+1];
2282e68589bSMark F. Adams           data[1] = -coords[3*kk+2]; data[M+1] = 0.0;             data[2*M+1] =  coords[3*kk];
2292e68589bSMark F. Adams           data[2] =  coords[3*kk+1]; data[M+2] = -coords[3*kk];   data[2*M+2] = 0.0;
2302e68589bSMark F. Adams         }
2312e68589bSMark F. Adams       }
2322e68589bSMark F. Adams     }
2332e68589bSMark F. Adams   }
2342e68589bSMark F. Adams   pc_gamg->data_sz = arrsz;
2352e68589bSMark F. Adams   PetscFunctionReturn(0);
2362e68589bSMark F. Adams }
2372e68589bSMark F. Adams 
2382e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
2392e68589bSMark F. Adams /*
240a2f3521dSMark F. Adams    PCSetData_AGG - called if data is not set with PCSetCoordinates.
241a2f3521dSMark F. Adams       Looks in Mat for near null space.
242a2f3521dSMark F. Adams       Does not work for Stokes
2432e68589bSMark F. Adams 
2442e68589bSMark F. Adams   Input Parameter:
2452e68589bSMark F. Adams    . pc -
246a2f3521dSMark F. Adams    . a_A - matrix to get (near) null space out of.
2472e68589bSMark F. Adams */
248e0877f53SBarry Smith static PetscErrorCode PCSetData_AGG(PC pc, Mat a_A)
2492e68589bSMark F. Adams {
250b8cd405aSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
251b8cd405aSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
252b8cd405aSMark F. Adams   MatNullSpace   mnull;
25366f2ef4bSMark Adams 
254b817416eSBarry Smith   PetscFunctionBegin;
2559566063dSJacob Faibussowitsch   PetscCall(MatGetNearNullSpace(a_A, &mnull));
256b8cd405aSMark F. Adams   if (!mnull) {
25766f2ef4bSMark Adams     DM dm;
2589566063dSJacob Faibussowitsch     PetscCall(PCGetDM(pc, &dm));
25966f2ef4bSMark Adams     if (!dm) {
2609566063dSJacob Faibussowitsch       PetscCall(MatGetDM(a_A, &dm));
26166f2ef4bSMark Adams     }
26266f2ef4bSMark Adams     if (dm) {
26366f2ef4bSMark Adams       PetscObject deformation;
264b0d30dd6SMatthew G. Knepley       PetscInt    Nf;
265b0d30dd6SMatthew G. Knepley 
2669566063dSJacob Faibussowitsch       PetscCall(DMGetNumFields(dm, &Nf));
267b0d30dd6SMatthew G. Knepley       if (Nf) {
2689566063dSJacob Faibussowitsch         PetscCall(DMGetField(dm, 0, NULL, &deformation));
2699566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery((PetscObject)deformation,"nearnullspace",(PetscObject*)&mnull));
27066f2ef4bSMark Adams         if (!mnull) {
2719566063dSJacob Faibussowitsch           PetscCall(PetscObjectQuery((PetscObject)deformation,"nullspace",(PetscObject*)&mnull));
27266f2ef4bSMark Adams         }
27366f2ef4bSMark Adams       }
27466f2ef4bSMark Adams     }
275b0d30dd6SMatthew G. Knepley   }
27666f2ef4bSMark Adams 
27766f2ef4bSMark Adams   if (!mnull) {
278a2f3521dSMark F. Adams     PetscInt bs,NN,MM;
2799566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSize(a_A, &bs));
2809566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(a_A, &MM, &NN));
28163a3b9bcSJacob Faibussowitsch     PetscCheck(MM % bs == 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MM %" PetscInt_FMT " must be divisible by bs %" PetscInt_FMT,MM,bs);
2829566063dSJacob Faibussowitsch     PetscCall(PCSetCoordinates_AGG(pc, bs, MM/bs, NULL));
283806fa848SBarry Smith   } else {
284b8cd405aSMark F. Adams     PetscReal         *nullvec;
285b8cd405aSMark F. Adams     PetscBool         has_const;
286b8cd405aSMark F. Adams     PetscInt          i,j,mlocal,nvec,bs;
287d5d25401SBarry Smith     const Vec         *vecs;
288d5d25401SBarry Smith     const PetscScalar *v;
289b817416eSBarry Smith 
2909566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(a_A,&mlocal,NULL));
2919566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs));
292a0dea87dSMark Adams     pc_gamg->data_sz = (nvec+!!has_const)*mlocal;
2939566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1((nvec+!!has_const)*mlocal,&nullvec));
294b8cd405aSMark F. Adams     if (has_const) for (i=0; i<mlocal; i++) nullvec[i] = 1.0;
295b8cd405aSMark F. Adams     for (i=0; i<nvec; i++) {
2969566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(vecs[i],&v));
297b8cd405aSMark F. Adams       for (j=0; j<mlocal; j++) nullvec[(i+!!has_const)*mlocal + j] = PetscRealPart(v[j]);
2989566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(vecs[i],&v));
299b8cd405aSMark F. Adams     }
300b8cd405aSMark F. Adams     pc_gamg->data           = nullvec;
301b8cd405aSMark F. Adams     pc_gamg->data_cell_cols = (nvec+!!has_const);
3029566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSize(a_A, &bs));
303b8cd405aSMark F. Adams     pc_gamg->data_cell_rows = bs;
304b8cd405aSMark F. Adams   }
3052e68589bSMark F. Adams   PetscFunctionReturn(0);
3062e68589bSMark F. Adams }
3072e68589bSMark F. Adams 
3082e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
3092e68589bSMark F. Adams /*
310bae903cbSmarkadams4   formProl0 - collect null space data for each aggregate, do QR, put R in coarse grid data and Q in P_0
3112e68589bSMark F. Adams 
3122e68589bSMark F. Adams   Input Parameter:
3132adfe9d3SBarry Smith    . agg_llists - list of arrays with aggregates -- list from selected vertices of aggregate unselected vertices
3142adfe9d3SBarry Smith    . bs - row block size
3150cbbd2e1SMark F. Adams    . nSAvec - column bs of new P
3160cbbd2e1SMark F. Adams    . my0crs - global index of start of locals
3172adfe9d3SBarry Smith    . data_stride - bs*(nloc nodes + ghost nodes) [data_stride][nSAvec]
3182e68589bSMark F. Adams    . data_in[data_stride*nSAvec] - local data on fine grid
3192e68589bSMark F. Adams    . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist'
320bae903cbSmarkadams4 
3212e68589bSMark F. Adams   Output Parameter:
3222e68589bSMark F. Adams    . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data
3232e68589bSMark F. Adams    . a_Prol - prolongation operator
324bae903cbSmarkadams4 
3252e68589bSMark F. Adams */
3262adfe9d3SBarry Smith 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)
3272e68589bSMark F. Adams {
328bd026e97SJed Brown   PetscInt        Istart,my0,Iend,nloc,clid,flid = 0,aggID,kk,jj,ii,mm,nSelected,minsz,nghosts,out_data_stride;
3293b4367a7SBarry Smith   MPI_Comm        comm;
3302e68589bSMark F. Adams   PetscReal       *out_data;
331539c167fSBarry Smith   PetscCDIntNd    *pos;
3321943db53SBarry Smith   PCGAMGHashTable fgid_flid;
3330cbbd2e1SMark F. Adams 
3342e68589bSMark F. Adams   PetscFunctionBegin;
3359566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)a_Prol,&comm));
3369566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(a_Prol, &Istart, &Iend));
33771959b99SBarry Smith   nloc = (Iend-Istart)/bs; my0 = Istart/bs;
33863a3b9bcSJacob 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);
3390cbbd2e1SMark F. Adams   Iend   /= bs;
3400cbbd2e1SMark F. Adams   nghosts = data_stride/bs - nloc;
3412e68589bSMark F. Adams 
3429566063dSJacob Faibussowitsch   PetscCall(PCGAMGHashTableCreate(2*nghosts+1, &fgid_flid));
3430cbbd2e1SMark F. Adams   for (kk=0; kk<nghosts; kk++) {
3449566063dSJacob Faibussowitsch     PetscCall(PCGAMGHashTableAdd(&fgid_flid, flid_fgid[nloc+kk], nloc+kk));
3452e68589bSMark F. Adams   }
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++) {
4102e68589bSMark F. Adams         for (jj = 0; jj < N; jj++, kk++) {
4112e68589bSMark F. Adams           qqc[jj*Mdata + ii] = .0;
4122e68589bSMark F. Adams         }
4132e68589bSMark F. Adams       }
4142e68589bSMark F. Adams 
4152e68589bSMark F. Adams       /* QR */
4169566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
4178b83055fSJed Brown       PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
4189566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPop());
41908401ef6SPierre Jolivet       PetscCheck(INFO == 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"xGEQRF error");
4202e68589bSMark F. Adams       /* get R - column oriented - output B_{i+1} */
4212e68589bSMark F. Adams       {
4222e68589bSMark F. Adams         PetscReal *data = &out_data[clid*nSAvec];
4232e68589bSMark F. Adams         for (jj = 0; jj < nSAvec; jj++) {
4242e68589bSMark F. Adams           for (ii = 0; ii < nSAvec; ii++) {
42508401ef6SPierre 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);
4260cbbd2e1SMark F. Adams            if (ii <= jj) data[jj*out_data_stride + ii] = PetscRealPart(qqc[jj*Mdata + ii]);
4270cbbd2e1SMark F. Adams            else data[jj*out_data_stride + ii] = 0.;
4282e68589bSMark F. Adams           }
4292e68589bSMark F. Adams         }
4302e68589bSMark F. Adams       }
4312e68589bSMark F. Adams 
4322e68589bSMark F. Adams       /* get Q - row oriented */
433c964aadfSJose E. Roman       PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
43463a3b9bcSJacob Faibussowitsch       PetscCheck(INFO == 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"xORGQR error arg %" PetscBLASInt_FMT,-INFO);
4352e68589bSMark F. Adams 
4362e68589bSMark F. Adams       for (ii = 0; ii < M; ii++) {
4372e68589bSMark F. Adams         for (jj = 0; jj < N; jj++) {
4382e68589bSMark F. Adams           qqr[N*ii + jj] = qqc[jj*Mdata + ii];
4392e68589bSMark F. Adams         }
4402e68589bSMark F. Adams       }
4412e68589bSMark F. Adams 
4422e68589bSMark F. Adams       /* add diagonal block of P0 */
443c8b0795cSMark F. Adams       for (kk=0; kk<N; kk++) {
444c8b0795cSMark F. Adams         cids[kk] = N*cgid + kk; /* global col IDs in P0 */
445c8b0795cSMark F. Adams       }
4469566063dSJacob Faibussowitsch       PetscCall(MatSetValues(a_Prol,M,fids,N,cids,qqr,INSERT_VALUES));
4479566063dSJacob Faibussowitsch       PetscCall(PetscFree5(qqc,qqr,TAU,WORK,fids));
448b43b03e9SMark F. Adams       clid++;
4490cbbd2e1SMark F. Adams     } /* coarse agg */
4500cbbd2e1SMark F. Adams   } /* for all fine nodes */
4519566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY));
4529566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY));
4539566063dSJacob Faibussowitsch   PetscCall(PCGAMGHashTableDestroy(&fgid_flid));
4542e68589bSMark F. Adams   PetscFunctionReturn(0);
4552e68589bSMark F. Adams }
4562e68589bSMark F. Adams 
4575adeb434SBarry Smith static PetscErrorCode PCView_GAMG_AGG(PC pc,PetscViewer viewer)
4585adeb434SBarry Smith {
4595adeb434SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
4605adeb434SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
4615adeb434SBarry Smith   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
4625adeb434SBarry Smith 
4635adeb434SBarry Smith   PetscFunctionBegin;
4649566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"      AGG specific options\n"));
465bae903cbSmarkadams4   PetscCall(PetscViewerASCIIPrintf(viewer,"        Symmetric graph %s\n",pc_gamg_agg->symmetrize_graph ? "true" : "false"));
466bae903cbSmarkadams4   PetscCall(PetscViewerASCIIPrintf(viewer,"        Number of levels to square graph %d\n",(int)pc_gamg_agg->aggressive_coarsening_levels));
467bae903cbSmarkadams4   PetscCall(PetscViewerASCIIPrintf(viewer,"        Number smoothing steps %d\n",(int)pc_gamg_agg->nsmooths));
4685adeb434SBarry Smith   PetscFunctionReturn(0);
4695adeb434SBarry Smith }
4705adeb434SBarry Smith 
4712e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
4722e68589bSMark F. Adams /*
473fd1112cbSBarry Smith    PCGAMGGraph_AGG
4742e68589bSMark F. Adams 
4752e68589bSMark F. Adams   Input Parameter:
4762e68589bSMark F. Adams    . pc - this
4772e68589bSMark F. Adams    . Amat - matrix on this fine level
478bae903cbSmarkadams4 
4792e68589bSMark F. Adams   Output Parameter:
480c8b0795cSMark F. Adams    . a_Gmat -
4812e68589bSMark F. Adams */
482e0877f53SBarry Smith static PetscErrorCode PCGAMGGraph_AGG(PC pc,Mat Amat,Mat *a_Gmat)
483c8b0795cSMark F. Adams {
484c8b0795cSMark F. Adams   PC_MG                     *mg          = (PC_MG*)pc->data;
485c8b0795cSMark F. Adams   PC_GAMG                   *pc_gamg     = (PC_GAMG*)mg->innerctx;
486c1eae691SMark Adams   const PetscReal           vfilter      = pc_gamg->threshold[pc_gamg->current_level];
487c8b0795cSMark F. Adams   PC_GAMG_AGG               *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
48872833a62Smarkadams4   Mat                       Gmat,F=NULL;
4893b4367a7SBarry Smith   MPI_Comm                  comm;
49067744fedSMark F. Adams   PetscBool /* set,flg , */ symm;
491c8b0795cSMark F. Adams 
492c8b0795cSMark F. Adams   PetscFunctionBegin;
4939566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Amat,&comm));
494c8b0795cSMark F. Adams 
4959566063dSJacob Faibussowitsch   /* PetscCall(MatIsSymmetricKnown(Amat, &set, &flg)); || !(set && flg) -- this causes lot of symm calls */
496bae903cbSmarkadams4   symm = (PetscBool)(pc_gamg_agg->symmetrize_graph); /* && !pc_gamg_agg->square_graph; */
4970cbbd2e1SMark F. Adams 
498849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_GRAPH],0,0,0,0));
49972833a62Smarkadams4   PetscCall(MatCreateGraph(Amat, symm, PETSC_TRUE, &Gmat));
500849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_GRAPH],0,0,0,0));
501849bee69SMark Adams 
502849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_FILTER],0,0,0,0));
50372833a62Smarkadams4   PetscCall(MatFilter(Gmat, vfilter,&F));
50472833a62Smarkadams4   if (F) {
50572833a62Smarkadams4     PetscCall(MatDestroy(&Gmat));
50672833a62Smarkadams4     Gmat = F;
50772833a62Smarkadams4   }
508849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_FILTER],0,0,0,0));
509849bee69SMark Adams 
510e0940f08SMark F. Adams   *a_Gmat = Gmat;
511849bee69SMark Adams 
512c8b0795cSMark F. Adams   PetscFunctionReturn(0);
513c8b0795cSMark F. Adams }
514c8b0795cSMark F. Adams 
515c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
516c8b0795cSMark F. Adams /*
517bae903cbSmarkadams4    PCGAMGCoarsen_AGG - supports squaring the graph (deprecated) and new graph for
518bae903cbSmarkadams4      communication of QR data used with HEM and MISk coarsening
519c8b0795cSMark F. Adams 
520c8b0795cSMark F. Adams   Input Parameter:
521e0940f08SMark F. Adams    . a_pc - this
522bae903cbSmarkadams4 
523e0940f08SMark F. Adams   Input/Output Parameter:
524bae903cbSmarkadams4    . a_Gmat1 - graph to coarsen (in), graph off processor edges for QR gather scatter (out)
525bae903cbSmarkadams4 
526c8b0795cSMark F. Adams   Output Parameter:
5270cbbd2e1SMark F. Adams    . agg_lists - list of aggregates
528bae903cbSmarkadams4 
529c8b0795cSMark F. Adams */
530e0877f53SBarry Smith static PetscErrorCode PCGAMGCoarsen_AGG(PC a_pc,Mat *a_Gmat1, PetscCoarsenData **agg_lists)
531c8b0795cSMark F. Adams {
532e0940f08SMark F. Adams   PC_MG          *mg          = (PC_MG*)a_pc->data;
533c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
534c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
535bae903cbSmarkadams4   Mat            mat, Gmat1 = *a_Gmat1;  /* aggressive graph */
5360cbbd2e1SMark F. Adams   IS             perm;
537bae903cbSmarkadams4   PetscInt       Istart,Iend,Ii,nloc,bs,nn;
538bae903cbSmarkadams4   PetscInt       *permute,*degree;
539c8b0795cSMark F. Adams   PetscBool      *bIndexSet;
540b43b03e9SMark F. Adams   MatCoarsen     crs;
5413b4367a7SBarry Smith   MPI_Comm       comm;
542aea53286SMark Adams   PetscReal      hashfact;
543e2c00dcbSBarry Smith   PetscInt       iSwapIndex;
5443b50add6SMark Adams   PetscRandom    random;
545c8b0795cSMark F. Adams 
546c8b0795cSMark F. Adams   PetscFunctionBegin;
547849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_COARSEN],0,0,0,0));
5489566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Gmat1,&comm));
549bae903cbSmarkadams4   PetscCall(MatGetLocalSize(Gmat1, &nn, NULL));
5509566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(Gmat1, &bs));
55163a3b9bcSJacob Faibussowitsch   PetscCheck(bs == 1,PETSC_COMM_SELF,PETSC_ERR_PLIB,"bs %" PetscInt_FMT " must be 1",bs);
552bae903cbSmarkadams4   nloc = nn/bs;
553c8b0795cSMark F. Adams 
5545cfd4bd9SMark Adams   /* get MIS aggs - randomize */
555bae903cbSmarkadams4   PetscCall(PetscMalloc2(nloc, &permute,nloc, &degree));
5569566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(nloc, &bIndexSet));
5576e09b0e3SStefano Zampini   for (Ii = 0; Ii < nloc; Ii++) permute[Ii] = Ii;
5589566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_SELF,&random));
5599566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Gmat1, &Istart, &Iend));
560c8b0795cSMark F. Adams   for (Ii = 0; Ii < nloc; Ii++) {
561bae903cbSmarkadams4     PetscInt nc;
562bae903cbSmarkadams4     PetscCall(MatGetRow(Gmat1,Istart+Ii,&nc,NULL,NULL));
563bae903cbSmarkadams4     degree[Ii] = nc;
564bae903cbSmarkadams4     PetscCall(MatRestoreRow(Gmat1,Istart+Ii,&nc,NULL,NULL));
565bae903cbSmarkadams4   }
566bae903cbSmarkadams4   for (Ii = 0; Ii < nloc; Ii++) {
5679566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValueReal(random,&hashfact));
568aea53286SMark Adams     iSwapIndex = (PetscInt) (hashfact*nloc)%nloc;
569c8b0795cSMark F. Adams     if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) {
570c8b0795cSMark F. Adams       PetscInt iTemp = permute[iSwapIndex];
571c8b0795cSMark F. Adams       permute[iSwapIndex]   = permute[Ii];
572c8b0795cSMark F. Adams       permute[Ii]           = iTemp;
573bae903cbSmarkadams4       iTemp = degree[iSwapIndex];
574bae903cbSmarkadams4       degree[iSwapIndex]   = degree[Ii];
575bae903cbSmarkadams4       degree[Ii]           = iTemp;
576c8b0795cSMark F. Adams       bIndexSet[iSwapIndex] = PETSC_TRUE;
577c8b0795cSMark F. Adams     }
578c8b0795cSMark F. Adams   }
579bae903cbSmarkadams4   // create minimum degree ordering
580bae903cbSmarkadams4   PetscCall(PetscSortIntWithArray(nloc,degree,permute));
581bae903cbSmarkadams4 
5829566063dSJacob Faibussowitsch   PetscCall(PetscFree(bIndexSet));
5839566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&random));
5849566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm));
585849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_MIS],0,0,0,0));
5869566063dSJacob Faibussowitsch   PetscCall(MatCoarsenCreate(comm, &crs));
5879566063dSJacob Faibussowitsch   PetscCall(MatCoarsenSetFromOptions(crs));
5889566063dSJacob Faibussowitsch   PetscCall(MatCoarsenSetGreedyOrdering(crs, perm));
589bae903cbSmarkadams4   PetscCall(MatCoarsenSetAdjacency(crs, Gmat1));
5909566063dSJacob Faibussowitsch   PetscCall(MatCoarsenSetStrictAggs(crs, PETSC_TRUE));
591bae903cbSmarkadams4   if (pc_gamg->current_level < pc_gamg_agg->aggressive_coarsening_levels) PetscCall(MatCoarsenMISKSetDistance(crs,2)); // hardwire to MIS-2
592bae903cbSmarkadams4   else PetscCall(MatCoarsenMISKSetDistance(crs,1)); // MIS
5939566063dSJacob Faibussowitsch   PetscCall(MatCoarsenApply(crs));
594bae903cbSmarkadams4   PetscCall(MatCoarsenViewFromOptions(crs,NULL,"-mat_coarsen_view"));
5959566063dSJacob Faibussowitsch   PetscCall(MatCoarsenGetData(crs, agg_lists)); /* output */
5969566063dSJacob Faibussowitsch   PetscCall(MatCoarsenDestroy(&crs));
597b43b03e9SMark F. Adams 
5989566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
599bae903cbSmarkadams4   PetscCall(PetscFree2(permute,degree));
600849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_MIS],0,0,0,0));
6019f3f12c8SMark F. Adams 
602bae903cbSmarkadams4   {
603bae903cbSmarkadams4     PetscCoarsenData *llist = *agg_lists;
6049ab59c8bSMark Adams     /* see if we have a matrix that takes precedence (returned from MatCoarsenApply) */
6059566063dSJacob Faibussowitsch     PetscCall(PetscCDGetMat(llist, &mat));
6060cbbd2e1SMark F. Adams     if (mat) {
6079566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&Gmat1));
6080cbbd2e1SMark F. Adams       *a_Gmat1 = mat; /* output */
6090cbbd2e1SMark F. Adams     }
6100cbbd2e1SMark F. Adams   }
611849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_COARSEN],0,0,0,0));
612c8b0795cSMark F. Adams   PetscFunctionReturn(0);
613c8b0795cSMark F. Adams }
614c8b0795cSMark F. Adams 
615c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
616c8b0795cSMark F. Adams /*
6170cbbd2e1SMark F. Adams  PCGAMGProlongator_AGG
618c8b0795cSMark F. Adams 
619c8b0795cSMark F. Adams  Input Parameter:
620c8b0795cSMark F. Adams  . pc - this
621c8b0795cSMark F. Adams  . Amat - matrix on this fine level
622c8b0795cSMark F. Adams  . Graph - used to get ghost data for nodes in
6230cbbd2e1SMark F. Adams  . agg_lists - list of aggregates
624c8b0795cSMark F. Adams  Output Parameter:
625c8b0795cSMark F. Adams  . a_P_out - prolongation operator to the next level
626c8b0795cSMark F. Adams  */
627e0877f53SBarry Smith static PetscErrorCode PCGAMGProlongator_AGG(PC pc,Mat Amat,Mat Gmat,PetscCoarsenData *agg_lists,Mat *a_P_out)
6282e68589bSMark F. Adams {
6292e68589bSMark F. Adams   PC_MG          *mg       = (PC_MG*)pc->data;
6302e68589bSMark F. Adams   PC_GAMG        *pc_gamg  = (PC_GAMG*)mg->innerctx;
6314a2f8832SBarry Smith   const PetscInt col_bs = pc_gamg->data_cell_cols;
632c8b0795cSMark F. Adams   PetscInt       Istart,Iend,nloc,ii,jj,kk,my0,nLocalSelected,bs;
633c8b0795cSMark F. Adams   Mat            Prol;
634d5d25401SBarry Smith   PetscMPIInt    size;
6353b4367a7SBarry Smith   MPI_Comm       comm;
636c8b0795cSMark F. Adams   PetscReal      *data_w_ghost;
637c8b0795cSMark F. Adams   PetscInt       myCrs0, nbnodes=0, *flid_fgid;
638d9558ea9SBarry Smith   MatType        mtype;
6392e68589bSMark F. Adams 
6402e68589bSMark F. Adams   PetscFunctionBegin;
6419566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Amat,&comm));
64208401ef6SPierre Jolivet   PetscCheck(col_bs >= 1,comm,PETSC_ERR_PLIB,"Column bs cannot be less than 1");
643849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROL],0,0,0,0));
6449566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
6459566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Amat, &Istart, &Iend));
6469566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(Amat, &bs));
64771959b99SBarry Smith   nloc = (Iend-Istart)/bs; my0 = Istart/bs;
64863a3b9bcSJacob 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);
6492e68589bSMark F. Adams 
6502e68589bSMark F. Adams   /* get 'nLocalSelected' */
6510cbbd2e1SMark F. Adams   for (ii=0, nLocalSelected = 0; ii < nloc; ii++) {
652e78576d6SMark F. Adams     PetscBool ise;
6530cbbd2e1SMark F. Adams     /* filter out singletons 0 or 1? */
6549566063dSJacob Faibussowitsch     PetscCall(PetscCDEmptyAt(agg_lists, ii, &ise));
655e78576d6SMark F. Adams     if (!ise) nLocalSelected++;
6562e68589bSMark F. Adams   }
6572e68589bSMark F. Adams 
6582e68589bSMark F. Adams   /* create prolongator, create P matrix */
6599566063dSJacob Faibussowitsch   PetscCall(MatGetType(Amat,&mtype));
6609566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &Prol));
6619566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Prol,nloc*bs,nLocalSelected*col_bs,PETSC_DETERMINE,PETSC_DETERMINE));
6629566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSizes(Prol, bs, col_bs));
6639566063dSJacob Faibussowitsch   PetscCall(MatSetType(Prol, mtype));
6649566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(Prol,col_bs, NULL));
6659566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(Prol,col_bs, NULL, col_bs, NULL));
6662e68589bSMark F. Adams 
6672e68589bSMark F. Adams   /* can get all points "removed" */
6689566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Prol, &kk, &ii));
6697f66b68fSMark Adams   if (!ii) {
67063a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(pc,"%s: No selected points on coarse grid\n",((PetscObject)pc)->prefix));
6719566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Prol));
6720298fd71SBarry Smith     *a_P_out = NULL;  /* out */
673849bee69SMark Adams     PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROL],0,0,0,0));
6742e68589bSMark F. Adams     PetscFunctionReturn(0);
6752e68589bSMark F. Adams   }
67663a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(pc,"%s: New grid %" PetscInt_FMT " nodes\n",((PetscObject)pc)->prefix,ii/col_bs));
6779566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRangeColumn(Prol, &myCrs0, &kk));
6780cbbd2e1SMark F. Adams 
67963a3b9bcSJacob 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);
680c8b0795cSMark F. Adams   myCrs0 = myCrs0/col_bs;
68163a3b9bcSJacob 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);
6822e68589bSMark F. Adams 
6832e68589bSMark F. Adams   /* create global vector of data in 'data_w_ghost' */
684849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROLA],0,0,0,0));
685bae903cbSmarkadams4   if (size > 1) { /* get ghost null space data */
6862e68589bSMark F. Adams     PetscReal *tmp_gdata,*tmp_ldata,*tp2;
6879566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nloc, &tmp_ldata));
6884a2f8832SBarry Smith     for (jj = 0; jj < col_bs; jj++) {
689c8b0795cSMark F. Adams       for (kk = 0; kk < bs; kk++) {
690a2f3521dSMark F. Adams         PetscInt        ii,stride;
691c8b0795cSMark F. Adams         const PetscReal *tp = pc_gamg->data + jj*bs*nloc + kk;
6922fa5cd67SKarl Rupp         for (ii = 0; ii < nloc; ii++, tp += bs) tmp_ldata[ii] = *tp;
6932fa5cd67SKarl Rupp 
6949566063dSJacob Faibussowitsch         PetscCall(PCGAMGGetDataWithGhosts(Gmat, 1, tmp_ldata, &stride, &tmp_gdata));
695a2f3521dSMark F. Adams 
696bae903cbSmarkadams4         if (!jj && !kk) { /* now I know how many total nodes - allocate TODO: move below and do in one 'col_bs' call */
6979566063dSJacob Faibussowitsch           PetscCall(PetscMalloc1(stride*bs*col_bs, &data_w_ghost));
698a2f3521dSMark F. Adams           nbnodes = bs*stride;
6992e68589bSMark F. Adams         }
700a2f3521dSMark F. Adams         tp2 = data_w_ghost + jj*bs*stride + kk;
7012fa5cd67SKarl Rupp         for (ii = 0; ii < stride; ii++, tp2 += bs) *tp2 = tmp_gdata[ii];
7029566063dSJacob Faibussowitsch         PetscCall(PetscFree(tmp_gdata));
7032e68589bSMark F. Adams       }
7042e68589bSMark F. Adams     }
7059566063dSJacob Faibussowitsch     PetscCall(PetscFree(tmp_ldata));
706806fa848SBarry Smith   } else {
707c8b0795cSMark F. Adams     nbnodes      = bs*nloc;
708c8b0795cSMark F. Adams     data_w_ghost = (PetscReal*)pc_gamg->data;
7092e68589bSMark F. Adams   }
7102e68589bSMark F. Adams 
711bae903cbSmarkadams4   /* get 'flid_fgid' TODO - move up to get 'stride' and do get null space data above in one step (jj loop) */
712c5df96a5SBarry Smith   if (size > 1) {
7132e68589bSMark F. Adams     PetscReal *fid_glid_loc,*fiddata;
714a2f3521dSMark F. Adams     PetscInt  stride;
7152e68589bSMark F. Adams 
7169566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nloc, &fid_glid_loc));
7172e68589bSMark F. Adams     for (kk=0; kk<nloc; kk++) fid_glid_loc[kk] = (PetscReal)(my0+kk);
7189566063dSJacob Faibussowitsch     PetscCall(PCGAMGGetDataWithGhosts(Gmat, 1, fid_glid_loc, &stride, &fiddata));
719bae903cbSmarkadams4     PetscCall(PetscMalloc1(stride, &flid_fgid)); /* copy real data to in */
720a2f3521dSMark F. Adams     for (kk=0; kk<stride; kk++) flid_fgid[kk] = (PetscInt)fiddata[kk];
7219566063dSJacob Faibussowitsch     PetscCall(PetscFree(fiddata));
722a2f3521dSMark F. Adams 
72363a3b9bcSJacob Faibussowitsch     PetscCheck(stride == nbnodes/bs,PETSC_COMM_SELF,PETSC_ERR_PLIB,"stride %" PetscInt_FMT " != nbnodes %" PetscInt_FMT "/bs %" PetscInt_FMT,stride,nbnodes,bs);
7249566063dSJacob Faibussowitsch     PetscCall(PetscFree(fid_glid_loc));
725806fa848SBarry Smith   } else {
7269566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nloc, &flid_fgid));
7272e68589bSMark F. Adams     for (kk=0; kk<nloc; kk++) flid_fgid[kk] = my0 + kk;
7282e68589bSMark F. Adams   }
729849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROLA],0,0,0,0));
730bae903cbSmarkadams4   /* get P0 */
731849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROLB],0,0,0,0));
732c8b0795cSMark F. Adams   {
7330298fd71SBarry Smith     PetscReal *data_out = NULL;
7349566063dSJacob Faibussowitsch     PetscCall(formProl0(agg_lists, bs, col_bs, myCrs0, nbnodes,data_w_ghost, flid_fgid, &data_out, Prol));
7359566063dSJacob Faibussowitsch     PetscCall(PetscFree(pc_gamg->data));
7362fa5cd67SKarl Rupp 
737c8b0795cSMark F. Adams     pc_gamg->data           = data_out;
7384a2f8832SBarry Smith     pc_gamg->data_cell_rows = col_bs;
7394a2f8832SBarry Smith     pc_gamg->data_sz        = col_bs*col_bs*nLocalSelected;
740c8b0795cSMark F. Adams   }
741849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROLB],0,0,0,0));
742849bee69SMark Adams   if (size > 1) {PetscCall(PetscFree(data_w_ghost));}
7439566063dSJacob Faibussowitsch   PetscCall(PetscFree(flid_fgid));
7442e68589bSMark F. Adams 
745c8b0795cSMark F. Adams   *a_P_out = Prol;  /* out */
746ed4e82eaSPeter Brune 
747849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROL],0,0,0,0));
748c8b0795cSMark F. Adams   PetscFunctionReturn(0);
749c8b0795cSMark F. Adams }
750c8b0795cSMark F. Adams 
751c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
752c8b0795cSMark F. Adams /*
753fd1112cbSBarry Smith    PCGAMGOptProlongator_AGG
754c8b0795cSMark F. Adams 
755c8b0795cSMark F. Adams   Input Parameter:
756c8b0795cSMark F. Adams    . pc - this
757c8b0795cSMark F. Adams    . Amat - matrix on this fine level
758c8b0795cSMark F. Adams  In/Output Parameter:
7592a808120SBarry Smith    . a_P - prolongation operator to the next level
760c8b0795cSMark F. Adams */
761e0877f53SBarry Smith static PetscErrorCode PCGAMGOptProlongator_AGG(PC pc,Mat Amat,Mat *a_P)
762c8b0795cSMark F. Adams {
763c8b0795cSMark F. Adams   PC_MG          *mg          = (PC_MG*)pc->data;
764c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
765c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
766c8b0795cSMark F. Adams   PetscInt       jj;
767c8b0795cSMark F. Adams   Mat            Prol  = *a_P;
7683b4367a7SBarry Smith   MPI_Comm       comm;
7692a808120SBarry Smith   KSP            eksp;
7702a808120SBarry Smith   Vec            bb, xx;
7712a808120SBarry Smith   PC             epc;
7722a808120SBarry Smith   PetscReal      alpha, emax, emin;
773c8b0795cSMark F. Adams 
774c8b0795cSMark F. Adams   PetscFunctionBegin;
7759566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Amat,&comm));
776849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_OPT],0,0,0,0));
777c8b0795cSMark F. Adams 
778d5d25401SBarry Smith   /* compute maximum singular value of operator to be used in smoother */
7792a808120SBarry Smith   if (0 < pc_gamg_agg->nsmooths) {
78018c3aa7eSMark     /* get eigen estimates */
78118c3aa7eSMark     if (pc_gamg->emax > 0) {
78218c3aa7eSMark       emin = pc_gamg->emin;
78318c3aa7eSMark       emax = pc_gamg->emax;
78418c3aa7eSMark     } else {
7850ed2132dSStefano Zampini       const char *prefix;
7860ed2132dSStefano Zampini 
7879566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(Amat, &bb, NULL));
7889566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(Amat, &xx, NULL));
7899566063dSJacob Faibussowitsch       PetscCall(KSPSetNoisy_Private(bb));
790e696ed0bSMark F. Adams 
7919566063dSJacob Faibussowitsch       PetscCall(KSPCreate(comm,&eksp));
7929566063dSJacob Faibussowitsch       PetscCall(PCGetOptionsPrefix(pc,&prefix));
7939566063dSJacob Faibussowitsch       PetscCall(KSPSetOptionsPrefix(eksp,prefix));
7949566063dSJacob Faibussowitsch       PetscCall(KSPAppendOptionsPrefix(eksp,"pc_gamg_esteig_"));
79573f7197eSJed Brown       {
796*b94d7dedSBarry Smith         PetscBool isset,sflg;
797*b94d7dedSBarry Smith         PetscCall(MatIsSPDKnown(Amat, &isset, &sflg));
798*b94d7dedSBarry Smith         if (isset && sflg) PetscCall(KSPSetType(eksp, KSPCG));
799d24ecf33SMark       }
8009566063dSJacob Faibussowitsch       PetscCall(KSPSetErrorIfNotConverged(eksp,pc->erroriffailure));
8019566063dSJacob Faibussowitsch       PetscCall(KSPSetNormType(eksp, KSP_NORM_NONE));
802c2436cacSMark F. Adams 
8039566063dSJacob Faibussowitsch       PetscCall(KSPSetInitialGuessNonzero(eksp, PETSC_FALSE));
8049566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(eksp, Amat, Amat));
8052e68589bSMark F. Adams 
8069566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(eksp, &epc));
8079566063dSJacob Faibussowitsch       PetscCall(PCSetType(epc, PCJACOBI));  /* smoother in smoothed agg. */
808c2436cacSMark F. Adams 
8099566063dSJacob 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
810074aec55SMark Adams 
8119566063dSJacob Faibussowitsch       PetscCall(KSPSetFromOptions(eksp));
8129566063dSJacob Faibussowitsch       PetscCall(KSPSetComputeSingularValues(eksp,PETSC_TRUE));
8139566063dSJacob Faibussowitsch       PetscCall(KSPSolve(eksp, bb, xx));
8149566063dSJacob Faibussowitsch       PetscCall(KSPCheckSolve(eksp,pc,xx));
8152e68589bSMark F. Adams 
8169566063dSJacob Faibussowitsch       PetscCall(KSPComputeExtremeSingularValues(eksp, &emax, &emin));
8179566063dSJacob Faibussowitsch       PetscCall(PetscInfo(pc,"%s: Smooth P0: max eigen=%e min=%e PC=%s\n",((PetscObject)pc)->prefix,(double)emax,(double)emin,PCJACOBI));
8189566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&xx));
8199566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&bb));
8209566063dSJacob Faibussowitsch       PetscCall(KSPDestroy(&eksp));
8212e68589bSMark F. Adams     }
82218c3aa7eSMark     if (pc_gamg->use_sa_esteig) {
82318c3aa7eSMark       mg->min_eigen_DinvA[pc_gamg->current_level] = emin;
82418c3aa7eSMark       mg->max_eigen_DinvA[pc_gamg->current_level] = emax;
82563a3b9bcSJacob 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));
82618c3aa7eSMark     } else {
82718c3aa7eSMark       mg->min_eigen_DinvA[pc_gamg->current_level] = 0;
82818c3aa7eSMark       mg->max_eigen_DinvA[pc_gamg->current_level] = 0;
82918c3aa7eSMark     }
83018c3aa7eSMark   } else {
83118c3aa7eSMark     mg->min_eigen_DinvA[pc_gamg->current_level] = 0;
83218c3aa7eSMark     mg->max_eigen_DinvA[pc_gamg->current_level] = 0;
83318c3aa7eSMark   }
8342e68589bSMark F. Adams 
8352a808120SBarry Smith   /* smooth P0 */
8362a808120SBarry Smith   for (jj = 0; jj < pc_gamg_agg->nsmooths; jj++) {
8372a808120SBarry Smith     Mat tMat;
8382a808120SBarry Smith     Vec diag;
8392a808120SBarry Smith 
840849bee69SMark Adams     PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_OPTSM],0,0,0,0));
8412a808120SBarry Smith 
8422e68589bSMark F. Adams     /* smooth P1 := (I - omega/lam D^{-1}A)P0 */
8439566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2],0,0,0,0));
8449566063dSJacob Faibussowitsch     PetscCall(MatMatMult(Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat));
8459566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2],0,0,0,0));
8469566063dSJacob Faibussowitsch     PetscCall(MatProductClear(tMat));
8479566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(Amat, &diag, NULL));
8489566063dSJacob Faibussowitsch     PetscCall(MatGetDiagonal(Amat, diag)); /* effectively PCJACOBI */
8499566063dSJacob Faibussowitsch     PetscCall(VecReciprocal(diag));
8509566063dSJacob Faibussowitsch     PetscCall(MatDiagonalScale(tMat, diag, NULL));
8519566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&diag));
852b4da3a1bSStefano Zampini 
853d5d25401SBarry Smith     /* TODO: Set a PCFailedReason and exit the building of the AMG preconditioner */
85408401ef6SPierre Jolivet     PetscCheck(emax != 0.0,PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"Computed maximum singular value as zero");
855d5d25401SBarry Smith     /* TODO: Document the 1.4 and don't hardwire it in this routine */
856b4ec6429SMark F. Adams     alpha = -1.4/emax;
857b4da3a1bSStefano Zampini 
8589566063dSJacob Faibussowitsch     PetscCall(MatAYPX(tMat, alpha, Prol, SUBSET_NONZERO_PATTERN));
8599566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Prol));
8602e68589bSMark F. Adams     Prol = tMat;
861849bee69SMark Adams     PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_OPTSM],0,0,0,0));
8622e68589bSMark F. Adams   }
863849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_OPT],0,0,0,0));
864c8b0795cSMark F. Adams   *a_P = Prol;
865c8b0795cSMark F. Adams   PetscFunctionReturn(0);
8662e68589bSMark F. Adams }
8672e68589bSMark F. Adams 
868c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
869c8b0795cSMark F. Adams /*
870c8b0795cSMark F. Adams    PCCreateGAMG_AGG
8712e68589bSMark F. Adams 
872c8b0795cSMark F. Adams   Input Parameter:
873c8b0795cSMark F. Adams    . pc -
874c8b0795cSMark F. Adams */
875c8b0795cSMark F. Adams PetscErrorCode  PCCreateGAMG_AGG(PC pc)
876c8b0795cSMark F. Adams {
877c8b0795cSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
878c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
879c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg;
8802e68589bSMark F. Adams 
881c8b0795cSMark F. Adams   PetscFunctionBegin;
882c8b0795cSMark F. Adams   /* create sub context for SA */
8839566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(pc,&pc_gamg_agg));
884c8b0795cSMark F. Adams   pc_gamg->subctx = pc_gamg_agg;
885c8b0795cSMark F. Adams 
8861ab5ffc9SJed Brown   pc_gamg->ops->setfromoptions = PCSetFromOptions_GAMG_AGG;
8879b8ffb57SJed Brown   pc_gamg->ops->destroy        = PCDestroy_GAMG_AGG;
888c8b0795cSMark F. Adams   /* reset does not do anything; setup not virtual */
889c8b0795cSMark F. Adams 
890c8b0795cSMark F. Adams   /* set internal function pointers */
891fd1112cbSBarry Smith   pc_gamg->ops->graph             = PCGAMGGraph_AGG;
8921ab5ffc9SJed Brown   pc_gamg->ops->coarsen           = PCGAMGCoarsen_AGG;
8931ab5ffc9SJed Brown   pc_gamg->ops->prolongator       = PCGAMGProlongator_AGG;
894fd1112cbSBarry Smith   pc_gamg->ops->optprolongator    = PCGAMGOptProlongator_AGG;
8951ab5ffc9SJed Brown   pc_gamg->ops->createdefaultdata = PCSetData_AGG;
8965adeb434SBarry Smith   pc_gamg->ops->view              = PCView_GAMG_AGG;
897c8b0795cSMark F. Adams 
898bae903cbSmarkadams4   pc_gamg_agg->aggressive_coarsening_levels = 0;
899bae903cbSmarkadams4   pc_gamg_agg->symmetrize_graph    = PETSC_FALSE;
900cab9ed1eSBarry Smith   pc_gamg_agg->nsmooths     = 1;
901cab9ed1eSBarry Smith 
9029566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNSmooths_C",PCGAMGSetNSmooths_AGG));
903bae903cbSmarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetSymmetrizeGraph_C",PCGAMGSetSymmetrizeGraph_AGG));
904bae903cbSmarkadams4   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetAggressiveLevels_C",PCGAMGSetAggressiveLevels_AGG));
9059566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_AGG));
9062e68589bSMark F. Adams   PetscFunctionReturn(0);
9072e68589bSMark F. Adams }
908