xref: /petsc/src/ksp/pc/impls/gamg/agg.c (revision 2e956fe4fc852fabc23b437482e1fb7b77fddb0d)
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;
12c8b0795cSMark F. Adams   PetscBool sym_graph;
139ab59c8bSMark Adams   PetscInt  square_graph;
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 /*@
52cab9ed1eSBarry Smith    PCGAMGSetSymGraph - 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:
61cab9ed1eSBarry Smith .  -pc_gamg_sym_graph <true,default=false> - symmetrize the graph before computing the aggregation
62c8b0795cSMark F. Adams 
63c8b0795cSMark F. Adams    Level: intermediate
64c8b0795cSMark F. Adams 
65db781477SPatrick Sanan .seealso: `PCGAMGSetSquareGraph()`
66c8b0795cSMark F. Adams @*/
67c8b0795cSMark F. Adams PetscErrorCode PCGAMGSetSymGraph(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);
72cac4c232SBarry Smith   PetscTryMethod(pc,"PCGAMGSetSymGraph_C",(PC,PetscBool),(pc,n));
73c8b0795cSMark F. Adams   PetscFunctionReturn(0);
74c8b0795cSMark F. Adams }
75c8b0795cSMark F. Adams 
76e0877f53SBarry Smith static PetscErrorCode PCGAMGSetSymGraph_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;
83c8b0795cSMark F. Adams   pc_gamg_agg->sym_graph = n;
842e68589bSMark F. Adams   PetscFunctionReturn(0);
852e68589bSMark F. Adams }
862e68589bSMark F. Adams 
87ef4ad70eSMark F. Adams /*@
88cab9ed1eSBarry Smith    PCGAMGSetSquareGraph -  Square the graph, ie. compute A'*A before aggregating it
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:
97cab9ed1eSBarry Smith .  -pc_gamg_square_graph <n,default = 1> - number of levels to square the graph on before aggregating it
98ef4ad70eSMark F. Adams 
99af3c827dSMark Adams    Notes:
100af3c827dSMark Adams    Squaring the graph increases the rate of coarsening (aggressive coarsening) and thereby reduces the complexity of the coarse grids, and generally results in slower solver converge rates. Reducing coarse grid complexity reduced the complexity of Galerkin coarse grid construction considerably.
101af3c827dSMark Adams 
102ef4ad70eSMark F. Adams    Level: intermediate
103ef4ad70eSMark F. Adams 
104db781477SPatrick Sanan .seealso: `PCGAMGSetSymGraph()`, `PCGAMGSetThreshold()`
105ef4ad70eSMark F. Adams @*/
1069ab59c8bSMark Adams PetscErrorCode PCGAMGSetSquareGraph(PC pc, PetscInt n)
107ef4ad70eSMark F. Adams {
108ef4ad70eSMark F. Adams   PetscFunctionBegin;
109ef4ad70eSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
110d5d25401SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
111cac4c232SBarry Smith   PetscTryMethod(pc,"PCGAMGSetSquareGraph_C",(PC,PetscInt),(pc,n));
112ef4ad70eSMark F. Adams   PetscFunctionReturn(0);
113ef4ad70eSMark F. Adams }
114ef4ad70eSMark F. Adams 
115e0877f53SBarry Smith static PetscErrorCode PCGAMGSetSquareGraph_AGG(PC pc, PetscInt n)
116ef4ad70eSMark F. Adams {
117ef4ad70eSMark F. Adams   PC_MG       *mg          = (PC_MG*)pc->data;
118ef4ad70eSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG*)mg->innerctx;
119ef4ad70eSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
120ef4ad70eSMark F. Adams 
121ef4ad70eSMark F. Adams   PetscFunctionBegin;
122ef4ad70eSMark F. Adams   pc_gamg_agg->square_graph = n;
123ef4ad70eSMark F. Adams   PetscFunctionReturn(0);
124ef4ad70eSMark F. Adams }
125ef4ad70eSMark F. Adams 
126e0877f53SBarry Smith static PetscErrorCode PCSetFromOptions_GAMG_AGG(PetscOptionItems *PetscOptionsObject,PC pc)
1272e68589bSMark F. Adams {
1282e68589bSMark F. Adams   PC_MG          *mg          = (PC_MG*)pc->data;
1292e68589bSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
130c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
1312e68589bSMark F. Adams 
1322e68589bSMark F. Adams   PetscFunctionBegin;
133d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject,"GAMG-AGG options");
1342e68589bSMark F. Adams   {
1359566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-pc_gamg_agg_nsmooths","smoothing steps for smoothed aggregation, usually 1","PCGAMGSetNSmooths",pc_gamg_agg->nsmooths,&pc_gamg_agg->nsmooths,NULL));
1369566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-pc_gamg_sym_graph","Set for asymmetric matrices","PCGAMGSetSymGraph",pc_gamg_agg->sym_graph,&pc_gamg_agg->sym_graph,NULL));
1379566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-pc_gamg_square_graph","Number of levels to square graph for faster coarsening and lower coarse grid complexity","PCGAMGSetSquareGraph",pc_gamg_agg->square_graph,&pc_gamg_agg->square_graph,NULL));
1382e68589bSMark F. Adams   }
139d0609cedSBarry Smith   PetscOptionsHeadEnd();
1402e68589bSMark F. Adams   PetscFunctionReturn(0);
1412e68589bSMark F. Adams }
1422e68589bSMark F. Adams 
1432e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
144e0877f53SBarry Smith static PetscErrorCode PCDestroy_GAMG_AGG(PC pc)
1452e68589bSMark F. Adams {
1462e68589bSMark F. Adams   PC_MG          *mg          = (PC_MG*)pc->data;
1472e68589bSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
1482e68589bSMark F. Adams 
1492e68589bSMark F. Adams   PetscFunctionBegin;
1509566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_gamg->subctx));
1519566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",NULL));
152*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNSmooths_C",NULL));
153*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetSymGraph_C",NULL));
154*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetSquareGraph_C",NULL));
1552e68589bSMark F. Adams   PetscFunctionReturn(0);
1562e68589bSMark F. Adams }
1572e68589bSMark F. Adams 
1582e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
1592e68589bSMark F. Adams /*
1602e68589bSMark F. Adams    PCSetCoordinates_AGG
161302f38e8SMark F. Adams      - collective
1622e68589bSMark F. Adams 
1632e68589bSMark F. Adams    Input Parameter:
1642e68589bSMark F. Adams    . pc - the preconditioner context
165a2f3521dSMark F. Adams    . ndm - dimesion of data (used for dof/vertex for Stokes)
166302f38e8SMark F. Adams    . a_nloc - number of vertices local
167302f38e8SMark F. Adams    . coords - [a_nloc][ndm] - interleaved coordinate data: {x_0, y_0, z_0, x_1, y_1, ...}
1682e68589bSMark F. Adams */
1691e6b0712SBarry Smith 
1701e6b0712SBarry Smith static PetscErrorCode PCSetCoordinates_AGG(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords)
1712e68589bSMark F. Adams {
1722e68589bSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
1732e68589bSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
17469344418SMark F. Adams   PetscInt       arrsz,kk,ii,jj,nloc,ndatarows,ndf;
175a2f3521dSMark F. Adams   Mat            mat = pc->pmat;
1762e68589bSMark F. Adams 
1772e68589bSMark F. Adams   PetscFunctionBegin;
178a2f3521dSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
179a2f3521dSMark F. Adams   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
180302f38e8SMark F. Adams   nloc = a_nloc;
1812e68589bSMark F. Adams 
1822e68589bSMark F. Adams   /* SA: null space vectors */
1839566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(mat, &ndf)); /* this does not work for Stokes */
18469344418SMark F. Adams   if (coords && ndf==1) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */
185a2f3521dSMark F. Adams   else if (coords) {
18663a3b9bcSJacob Faibussowitsch     PetscCheck(ndm <= ndf,PETSC_COMM_SELF,PETSC_ERR_PLIB,"degrees of motion %" PetscInt_FMT " > block size %" PetscInt_FMT,ndm,ndf);
18769344418SMark F. Adams     pc_gamg->data_cell_cols = (ndm==2 ? 3 : 6); /* displacement elasticity */
18869344418SMark F. Adams     if (ndm != ndf) {
18963a3b9bcSJacob 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);
190a2f3521dSMark F. Adams     }
191806fa848SBarry Smith   } else pc_gamg->data_cell_cols = ndf; /* no data, force SA with constant null space vectors */
19271959b99SBarry Smith   pc_gamg->data_cell_rows = ndatarows = ndf;
19363a3b9bcSJacob 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);
194c8b0795cSMark F. Adams   arrsz = nloc*pc_gamg->data_cell_rows*pc_gamg->data_cell_cols;
1952e68589bSMark F. Adams 
1967f66b68fSMark Adams   if (!pc_gamg->data || (pc_gamg->data_sz != arrsz)) {
1979566063dSJacob Faibussowitsch     PetscCall(PetscFree(pc_gamg->data));
1989566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(arrsz+1, &pc_gamg->data));
1992e68589bSMark F. Adams   }
2002e68589bSMark F. Adams   /* copy data in - column oriented */
2012e68589bSMark F. Adams   for (kk=0; kk<nloc; kk++) {
20269344418SMark F. Adams     const PetscInt M     = nloc*pc_gamg->data_cell_rows; /* stride into data */
20369344418SMark F. Adams     PetscReal      *data = &pc_gamg->data[kk*ndatarows]; /* start of cell */
204c8b0795cSMark F. Adams     if (pc_gamg->data_cell_cols==1) *data = 1.0;
2052e68589bSMark F. Adams     else {
20669344418SMark F. Adams       /* translational modes */
2072fa5cd67SKarl Rupp       for (ii=0;ii<ndatarows;ii++) {
2082fa5cd67SKarl Rupp         for (jj=0;jj<ndatarows;jj++) {
20969344418SMark F. Adams           if (ii==jj)data[ii*M + jj] = 1.0;
2102e68589bSMark F. Adams           else data[ii*M + jj] = 0.0;
2112fa5cd67SKarl Rupp         }
2122fa5cd67SKarl Rupp       }
21369344418SMark F. Adams 
21469344418SMark F. Adams       /* rotational modes */
2152e68589bSMark F. Adams       if (coords) {
21669344418SMark F. Adams         if (ndm == 2) {
2172e68589bSMark F. Adams           data   += 2*M;
2182e68589bSMark F. Adams           data[0] = -coords[2*kk+1];
2192e68589bSMark F. Adams           data[1] =  coords[2*kk];
220806fa848SBarry Smith         } else {
2212e68589bSMark F. Adams           data   += 3*M;
2222e68589bSMark F. Adams           data[0] = 0.0;             data[M+0] =  coords[3*kk+2]; data[2*M+0] = -coords[3*kk+1];
2232e68589bSMark F. Adams           data[1] = -coords[3*kk+2]; data[M+1] = 0.0;             data[2*M+1] =  coords[3*kk];
2242e68589bSMark F. Adams           data[2] =  coords[3*kk+1]; data[M+2] = -coords[3*kk];   data[2*M+2] = 0.0;
2252e68589bSMark F. Adams         }
2262e68589bSMark F. Adams       }
2272e68589bSMark F. Adams     }
2282e68589bSMark F. Adams   }
2292e68589bSMark F. Adams   pc_gamg->data_sz = arrsz;
2302e68589bSMark F. Adams   PetscFunctionReturn(0);
2312e68589bSMark F. Adams }
2322e68589bSMark F. Adams 
233b43b03e9SMark F. Adams typedef PetscInt NState;
234b43b03e9SMark F. Adams static const NState NOT_DONE=-2;
235b43b03e9SMark F. Adams static const NState DELETED =-1;
236b43b03e9SMark F. Adams static const NState REMOVED =-3;
237b43b03e9SMark F. Adams #define IS_SELECTED(s) (s!=DELETED && s!=NOT_DONE && s!=REMOVED)
238b43b03e9SMark F. Adams 
239c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
240c8b0795cSMark F. Adams /*
241b43b03e9SMark F. Adams    smoothAggs - greedy grab of with G1 (unsquared graph) -- AIJ specific
242b43b03e9SMark F. Adams      - AGG-MG specific: clears singletons out of 'selected_2'
243c8b0795cSMark F. Adams 
244c8b0795cSMark F. Adams    Input Parameter:
245fdb86571SRichard Tran Mills    . Gmat_2 - global matrix of graph (data not defined)   base (squared) graph
2462adfe9d3SBarry Smith    . Gmat_1 - base graph to grab with                 base graph
247c8b0795cSMark F. Adams    Input/Output Parameter:
2480cbbd2e1SMark F. Adams    . aggs_2 - linked list of aggs with gids)
249c8b0795cSMark F. Adams */
250e3df7ea3SBarry Smith static PetscErrorCode smoothAggs(PC pc,Mat Gmat_2, Mat Gmat_1,PetscCoarsenData *aggs_2)
251c8b0795cSMark F. Adams {
252c8b0795cSMark F. Adams   PetscBool      isMPI;
2530a545947SLisandro Dalcin   Mat_SeqAIJ     *matA_1, *matB_1=NULL;
2543b4367a7SBarry Smith   MPI_Comm       comm;
2550cbbd2e1SMark F. Adams   PetscInt       lid,*ii,*idx,ix,Iend,my0,kk,n,j;
2560a545947SLisandro Dalcin   Mat_MPIAIJ     *mpimat_2 = NULL, *mpimat_1=NULL;
257c8b0795cSMark F. Adams   const PetscInt nloc      = Gmat_2->rmap->n;
2580cbbd2e1SMark F. Adams   PetscScalar    *cpcol_1_state,*cpcol_2_state,*cpcol_2_par_orig,*lid_parent_gid;
2590cbbd2e1SMark F. Adams   PetscInt       *lid_cprowID_1;
260c8b0795cSMark F. Adams   NState         *lid_state;
2610cbbd2e1SMark F. Adams   Vec            ghost_par_orig2;
262c8b0795cSMark F. Adams 
263c8b0795cSMark F. Adams   PetscFunctionBegin;
2649566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Gmat_2,&comm));
2659566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Gmat_1,&my0,&Iend));
266c8b0795cSMark F. Adams 
267c8b0795cSMark F. Adams   /* get submatrices */
2689566063dSJacob Faibussowitsch   PetscCall(PetscStrbeginswith(((PetscObject)Gmat_1)->type_name,MATMPIAIJ,&isMPI));
269c8b0795cSMark F. Adams   if (isMPI) {
270c8b0795cSMark F. Adams     /* grab matrix objects */
271c8b0795cSMark F. Adams     mpimat_2 = (Mat_MPIAIJ*)Gmat_2->data;
272c8b0795cSMark F. Adams     mpimat_1 = (Mat_MPIAIJ*)Gmat_1->data;
273c8b0795cSMark F. Adams     matA_1   = (Mat_SeqAIJ*)mpimat_1->A->data;
274c8b0795cSMark F. Adams     matB_1   = (Mat_SeqAIJ*)mpimat_1->B->data;
275c8b0795cSMark F. Adams 
276c8b0795cSMark F. Adams     /* force compressed row storage for B matrix in AuxMat */
2779566063dSJacob Faibussowitsch     PetscCall(MatCheckCompressedRow(mpimat_1->B,matB_1->nonzerorowcnt,&matB_1->compressedrow,matB_1->i,Gmat_1->rmap->n,-1.0));
278c8b0795cSMark F. Adams 
2799566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nloc, &lid_cprowID_1));
2800cbbd2e1SMark F. Adams     for (lid = 0; lid < nloc; lid++) lid_cprowID_1[lid] = -1;
281c8b0795cSMark F. Adams     for (ix=0; ix<matB_1->compressedrow.nrows; ix++) {
282c8b0795cSMark F. Adams       PetscInt lid = matB_1->compressedrow.rindex[ix];
283c8b0795cSMark F. Adams       lid_cprowID_1[lid] = ix;
284c8b0795cSMark F. Adams     }
285806fa848SBarry Smith   } else {
28615687449SMark Adams     PetscBool isAIJ;
2879566063dSJacob Faibussowitsch     PetscCall(PetscStrbeginswith(((PetscObject)Gmat_1)->type_name,MATSEQAIJ,&isAIJ));
28828b400f6SJacob Faibussowitsch     PetscCheck(isAIJ,PETSC_COMM_SELF,PETSC_ERR_USER,"Require AIJ matrix.");
289c8b0795cSMark F. Adams     matA_1        = (Mat_SeqAIJ*)Gmat_1->data;
2900298fd71SBarry Smith     lid_cprowID_1 = NULL;
291c8b0795cSMark F. Adams   }
29278a438d6SMark Adams   if (nloc>0) {
29308401ef6SPierre Jolivet     PetscCheck(!matB_1 || matB_1->compressedrow.use,PETSC_COMM_SELF,PETSC_ERR_PLIB,"matB_1 && !matB_1->compressedrow.use: PETSc bug???");
29478a438d6SMark Adams   }
295c8b0795cSMark F. Adams   /* get state of locals and selected gid for deleted */
2969566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(nloc, &lid_state,nloc, &lid_parent_gid));
297c8b0795cSMark F. Adams   for (lid = 0; lid < nloc; lid++) {
2980cbbd2e1SMark F. Adams     lid_parent_gid[lid] = -1.0;
299c8b0795cSMark F. Adams     lid_state[lid]      = DELETED;
300c8b0795cSMark F. Adams   }
3010cbbd2e1SMark F. Adams 
3020cbbd2e1SMark F. Adams   /* set lid_state */
3030cbbd2e1SMark F. Adams   for (lid = 0; lid < nloc; lid++) {
304539c167fSBarry Smith     PetscCDIntNd *pos;
3059566063dSJacob Faibussowitsch     PetscCall(PetscCDGetHeadPos(aggs_2,lid,&pos));
306e78576d6SMark F. Adams     if (pos) {
307e78576d6SMark F. Adams       PetscInt gid1;
3082fa5cd67SKarl Rupp 
3099566063dSJacob Faibussowitsch       PetscCall(PetscCDIntNdGetID(pos, &gid1));
31063a3b9bcSJacob Faibussowitsch       PetscCheck(gid1 == lid+my0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"gid1 %" PetscInt_FMT " != lid %" PetscInt_FMT " + my0 %" PetscInt_FMT,gid1,lid,my0);
3110cbbd2e1SMark F. Adams       lid_state[lid] = gid1;
312b43b03e9SMark F. Adams     }
313b43b03e9SMark F. Adams   }
3140cbbd2e1SMark F. Adams 
3150cbbd2e1SMark F. Adams   /* map local to selected local, DELETED means a ghost owns it */
316c8b0795cSMark F. Adams   for (lid=kk=0; lid<nloc; lid++) {
317c8b0795cSMark F. Adams     NState state = lid_state[lid];
318c8b0795cSMark F. Adams     if (IS_SELECTED(state)) {
319539c167fSBarry Smith       PetscCDIntNd *pos;
3209566063dSJacob Faibussowitsch       PetscCall(PetscCDGetHeadPos(aggs_2,lid,&pos));
321e78576d6SMark F. Adams       while (pos) {
322e78576d6SMark F. Adams         PetscInt gid1;
3239566063dSJacob Faibussowitsch         PetscCall(PetscCDIntNdGetID(pos, &gid1));
3249566063dSJacob Faibussowitsch         PetscCall(PetscCDGetNextPos(aggs_2,lid,&pos));
3252fa5cd67SKarl Rupp         if (gid1 >= my0 && gid1 < Iend) lid_parent_gid[gid1-my0] = (PetscScalar)(lid + my0);
326c8b0795cSMark F. Adams       }
3270cbbd2e1SMark F. Adams     }
3280cbbd2e1SMark F. Adams   }
3290cbbd2e1SMark F. Adams   /* get 'cpcol_1/2_state' & cpcol_2_par_orig - uses mpimat_1/2->lvec for temp space */
330c8b0795cSMark F. Adams   if (isMPI) {
331c8b0795cSMark F. Adams     Vec tempVec;
332c8b0795cSMark F. Adams     /* get 'cpcol_1_state' */
3339566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(Gmat_1, &tempVec, NULL));
334c8b0795cSMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
335c8b0795cSMark F. Adams       PetscScalar v = (PetscScalar)lid_state[kk];
3369566063dSJacob Faibussowitsch       PetscCall(VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES));
337c8b0795cSMark F. Adams     }
3389566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(tempVec));
3399566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(tempVec));
3409566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD));
3419566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD));
3429566063dSJacob Faibussowitsch     PetscCall(VecGetArray(mpimat_1->lvec, &cpcol_1_state));
343c8b0795cSMark F. Adams     /* get 'cpcol_2_state' */
3449566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD));
3459566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD));
3469566063dSJacob Faibussowitsch     PetscCall(VecGetArray(mpimat_2->lvec, &cpcol_2_state));
3470cbbd2e1SMark F. Adams     /* get 'cpcol_2_par_orig' */
3480cbbd2e1SMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
3490cbbd2e1SMark F. Adams       PetscScalar v = (PetscScalar)lid_parent_gid[kk];
3509566063dSJacob Faibussowitsch       PetscCall(VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES));
3510cbbd2e1SMark F. Adams     }
3529566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(tempVec));
3539566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(tempVec));
3549566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(mpimat_2->lvec, &ghost_par_orig2));
3559566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD));
3569566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD));
3579566063dSJacob Faibussowitsch     PetscCall(VecGetArray(ghost_par_orig2, &cpcol_2_par_orig));
3580cbbd2e1SMark F. Adams 
3599566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&tempVec));
360c8b0795cSMark F. Adams   } /* ismpi */
361c8b0795cSMark F. Adams   for (lid=0; lid<nloc; lid++) {
362c8b0795cSMark F. Adams     NState state = lid_state[lid];
3630cbbd2e1SMark F. Adams     if (IS_SELECTED(state)) {
3640cbbd2e1SMark F. Adams       /* steal locals */
365c8b0795cSMark F. Adams       ii  = matA_1->i; n = ii[lid+1] - ii[lid];
366c8b0795cSMark F. Adams       idx = matA_1->j + ii[lid];
367c8b0795cSMark F. Adams       for (j=0; j<n; j++) {
3680cbbd2e1SMark F. Adams         PetscInt lidj   = idx[j], sgid;
369c8b0795cSMark F. Adams         NState   statej = lid_state[lidj];
3700cbbd2e1SMark F. Adams         if (statej==DELETED && (sgid=(PetscInt)PetscRealPart(lid_parent_gid[lidj])) != lid+my0) { /* steal local */
3710cbbd2e1SMark F. Adams           lid_parent_gid[lidj] = (PetscScalar)(lid+my0); /* send this if sgid is not local */
3720cbbd2e1SMark F. Adams           if (sgid >= my0 && sgid < Iend) {       /* I'm stealing this local from a local sgid */
3730cbbd2e1SMark F. Adams             PetscInt     hav=0,slid=sgid-my0,gidj=lidj+my0;
374539c167fSBarry Smith             PetscCDIntNd *pos,*last=NULL;
375c8b0795cSMark F. Adams             /* looking for local from local so id_llist_2 works */
3769566063dSJacob Faibussowitsch             PetscCall(PetscCDGetHeadPos(aggs_2,slid,&pos));
377e78576d6SMark F. Adams             while (pos) {
378e78576d6SMark F. Adams               PetscInt gid;
3799566063dSJacob Faibussowitsch               PetscCall(PetscCDIntNdGetID(pos, &gid));
3800cbbd2e1SMark F. Adams               if (gid == gidj) {
38128b400f6SJacob Faibussowitsch                 PetscCheck(last,PETSC_COMM_SELF,PETSC_ERR_PLIB,"last cannot be null");
3829566063dSJacob Faibussowitsch                 PetscCall(PetscCDRemoveNextNode(aggs_2, slid, last));
3839566063dSJacob Faibussowitsch                 PetscCall(PetscCDAppendNode(aggs_2, lid, pos));
3840cbbd2e1SMark F. Adams                 hav  = 1;
385c8b0795cSMark F. Adams                 break;
386806fa848SBarry Smith               } else last = pos;
3879566063dSJacob Faibussowitsch               PetscCall(PetscCDGetNextPos(aggs_2,slid,&pos));
388c8b0795cSMark F. Adams             }
389c8b0795cSMark F. Adams             if (hav != 1) {
39028b400f6SJacob Faibussowitsch               PetscCheck(hav,PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix");
39163a3b9bcSJacob Faibussowitsch               SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %" PetscInt_FMT " times???",hav);
392c8b0795cSMark F. Adams             }
393806fa848SBarry Smith           } else {            /* I'm stealing this local, owned by a ghost */
39408401ef6SPierre Jolivet             PetscCheck(sgid == -1,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Mat has an un-symmetric graph. Use '-%spc_gamg_sym_graph true' to symmetrize the graph or '-%spc_gamg_threshold -1' if the matrix is structurally symmetric.",((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "",((PetscObject)pc)->prefix ? ((PetscObject)pc)->prefix : "");
3959566063dSJacob Faibussowitsch             PetscCall(PetscCDAppendID(aggs_2, lid, lidj+my0));
396c8b0795cSMark F. Adams           }
397c8b0795cSMark F. Adams         }
3980cbbd2e1SMark F. Adams       } /* local neighbors */
399806fa848SBarry Smith     } else if (state == DELETED && lid_cprowID_1) {
4000cbbd2e1SMark F. Adams       PetscInt sgidold = (PetscInt)PetscRealPart(lid_parent_gid[lid]);
401c8b0795cSMark F. Adams       /* see if I have a selected ghost neighbor that will steal me */
402c8b0795cSMark F. Adams       if ((ix=lid_cprowID_1[lid]) != -1) {
403c8b0795cSMark F. Adams         ii  = matB_1->compressedrow.i; n = ii[ix+1] - ii[ix];
404c8b0795cSMark F. Adams         idx = matB_1->j + ii[ix];
405c8b0795cSMark F. Adams         for (j=0; j<n; j++) {
406c8b0795cSMark F. Adams           PetscInt cpid   = idx[j];
407c8b0795cSMark F. Adams           NState   statej = (NState)PetscRealPart(cpcol_1_state[cpid]);
408c8b0795cSMark F. Adams           if (IS_SELECTED(statej) && sgidold != (PetscInt)statej) { /* ghost will steal this, remove from my list */
4090cbbd2e1SMark F. Adams             lid_parent_gid[lid] = (PetscScalar)statej; /* send who selected */
4100cbbd2e1SMark F. Adams             if (sgidold>=my0 && sgidold<Iend) { /* this was mine */
4110cbbd2e1SMark F. Adams               PetscInt     hav=0,oldslidj=sgidold-my0;
412539c167fSBarry Smith               PetscCDIntNd *pos,*last=NULL;
4130cbbd2e1SMark F. Adams               /* remove from 'oldslidj' list */
4149566063dSJacob Faibussowitsch               PetscCall(PetscCDGetHeadPos(aggs_2,oldslidj,&pos));
415e78576d6SMark F. Adams               while (pos) {
416e78576d6SMark F. Adams                 PetscInt gid;
4179566063dSJacob Faibussowitsch                 PetscCall(PetscCDIntNdGetID(pos, &gid));
4180cbbd2e1SMark F. Adams                 if (lid+my0 == gid) {
4190cbbd2e1SMark F. Adams                   /* id_llist_2[lastid] = id_llist_2[flid];   /\* remove lid from oldslidj list *\/ */
42028b400f6SJacob Faibussowitsch                   PetscCheck(last,PETSC_COMM_SELF,PETSC_ERR_PLIB,"last cannot be null");
4219566063dSJacob Faibussowitsch                   PetscCall(PetscCDRemoveNextNode(aggs_2, oldslidj, last));
4220cbbd2e1SMark F. Adams                   /* ghost (PetscScalar)statej will add this later */
4230cbbd2e1SMark F. Adams                   hav = 1;
424c8b0795cSMark F. Adams                   break;
425806fa848SBarry Smith                 } else last = pos;
4269566063dSJacob Faibussowitsch                 PetscCall(PetscCDGetNextPos(aggs_2,oldslidj,&pos));
427c8b0795cSMark F. Adams               }
428c8b0795cSMark F. Adams               if (hav != 1) {
42928b400f6SJacob Faibussowitsch                 PetscCheck(hav,PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix");
43063a3b9bcSJacob Faibussowitsch                 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %" PetscInt_FMT " times???",hav);
431c8b0795cSMark F. Adams               }
432806fa848SBarry Smith             } else {
433d5d25401SBarry Smith               /* TODO: ghosts remove this later */
4340cbbd2e1SMark F. Adams             }
435c8b0795cSMark F. Adams           }
436c8b0795cSMark F. Adams         }
437c8b0795cSMark F. Adams       }
438c8b0795cSMark F. Adams     } /* selected/deleted */
439c8b0795cSMark F. Adams   } /* node loop */
440c8b0795cSMark F. Adams 
441c8b0795cSMark F. Adams   if (isMPI) {
4420cbbd2e1SMark F. Adams     PetscScalar     *cpcol_2_parent,*cpcol_2_gid;
4430cbbd2e1SMark F. Adams     Vec             tempVec,ghostgids2,ghostparents2;
4440cbbd2e1SMark F. Adams     PetscInt        cpid,nghost_2;
4451943db53SBarry Smith     PCGAMGHashTable gid_cpid;
446c8b0795cSMark F. Adams 
4479566063dSJacob Faibussowitsch     PetscCall(VecGetSize(mpimat_2->lvec, &nghost_2));
4489566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(Gmat_2, &tempVec, NULL));
4490cbbd2e1SMark F. Adams 
4500cbbd2e1SMark F. Adams     /* get 'cpcol_2_parent' */
451c8b0795cSMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
4529566063dSJacob Faibussowitsch       PetscCall(VecSetValues(tempVec, 1, &j, &lid_parent_gid[kk], INSERT_VALUES));
453c8b0795cSMark F. Adams     }
4549566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(tempVec));
4559566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(tempVec));
4569566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(mpimat_2->lvec, &ghostparents2));
4579566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD));
4589566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD));
4599566063dSJacob Faibussowitsch     PetscCall(VecGetArray(ghostparents2, &cpcol_2_parent));
4600cbbd2e1SMark F. Adams 
4610cbbd2e1SMark F. Adams     /* get 'cpcol_2_gid' */
4620cbbd2e1SMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
4630cbbd2e1SMark F. Adams       PetscScalar v = (PetscScalar)j;
4649566063dSJacob Faibussowitsch       PetscCall(VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES));
4650cbbd2e1SMark F. Adams     }
4669566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(tempVec));
4679566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(tempVec));
4689566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(mpimat_2->lvec, &ghostgids2));
4699566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD));
4709566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD));
4719566063dSJacob Faibussowitsch     PetscCall(VecGetArray(ghostgids2, &cpcol_2_gid));
4729566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&tempVec));
473c8b0795cSMark F. Adams 
4740cbbd2e1SMark F. Adams     /* look for deleted ghosts and add to table */
4759566063dSJacob Faibussowitsch     PetscCall(PCGAMGHashTableCreate(2*nghost_2+1, &gid_cpid));
4760cbbd2e1SMark F. Adams     for (cpid = 0; cpid < nghost_2; cpid++) {
4770cbbd2e1SMark F. Adams       NState state = (NState)PetscRealPart(cpcol_2_state[cpid]);
4780cbbd2e1SMark F. Adams       if (state==DELETED) {
4790cbbd2e1SMark F. Adams         PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]);
4800cbbd2e1SMark F. Adams         PetscInt sgid_old = (PetscInt)PetscRealPart(cpcol_2_par_orig[cpid]);
4810cbbd2e1SMark F. Adams         if (sgid_old == -1 && sgid_new != -1) {
4820cbbd2e1SMark F. Adams           PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]);
4839566063dSJacob Faibussowitsch           PetscCall(PCGAMGHashTableAdd(&gid_cpid, gid, cpid));
4840cbbd2e1SMark F. Adams         }
4850cbbd2e1SMark F. Adams       }
4860cbbd2e1SMark F. Adams     }
487c8b0795cSMark F. Adams 
4880cbbd2e1SMark F. Adams     /* look for deleted ghosts and see if they moved - remove it */
489c8b0795cSMark F. Adams     for (lid=0; lid<nloc; lid++) {
490c8b0795cSMark F. Adams       NState state = lid_state[lid];
491c8b0795cSMark F. Adams       if (IS_SELECTED(state)) {
492539c167fSBarry Smith         PetscCDIntNd *pos,*last=NULL;
493c8b0795cSMark F. Adams         /* look for deleted ghosts and see if they moved */
4949566063dSJacob Faibussowitsch         PetscCall(PetscCDGetHeadPos(aggs_2,lid,&pos));
495e78576d6SMark F. Adams         while (pos) {
496e78576d6SMark F. Adams           PetscInt gid;
4979566063dSJacob Faibussowitsch           PetscCall(PetscCDIntNdGetID(pos, &gid));
498e78576d6SMark F. Adams 
4990cbbd2e1SMark F. Adams           if (gid < my0 || gid >= Iend) {
5009566063dSJacob Faibussowitsch             PetscCall(PCGAMGHashTableFind(&gid_cpid, gid, &cpid));
5010cbbd2e1SMark F. Adams             if (cpid != -1) {
5020cbbd2e1SMark F. Adams               /* a moved ghost - */
5030cbbd2e1SMark F. Adams               /* id_llist_2[lastid] = id_llist_2[flid];    /\* remove 'flid' from list *\/ */
5049566063dSJacob Faibussowitsch               PetscCall(PetscCDRemoveNextNode(aggs_2, lid, last));
505806fa848SBarry Smith             } else last = pos;
506806fa848SBarry Smith           } else last = pos;
507e78576d6SMark F. Adams 
5089566063dSJacob Faibussowitsch           PetscCall(PetscCDGetNextPos(aggs_2,lid,&pos));
509c8b0795cSMark F. Adams         } /* loop over list of deleted */
510c8b0795cSMark F. Adams       } /* selected */
511c8b0795cSMark F. Adams     }
5129566063dSJacob Faibussowitsch     PetscCall(PCGAMGHashTableDestroy(&gid_cpid));
513c8b0795cSMark F. Adams 
5140cbbd2e1SMark F. Adams     /* look at ghosts, see if they changed - and it */
5150cbbd2e1SMark F. Adams     for (cpid = 0; cpid < nghost_2; cpid++) {
5160cbbd2e1SMark F. Adams       PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]);
5170cbbd2e1SMark F. Adams       if (sgid_new >= my0 && sgid_new < Iend) { /* this is mine */
5180cbbd2e1SMark F. Adams         PetscInt     gid     = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]);
5190cbbd2e1SMark F. Adams         PetscInt     slid_new=sgid_new-my0,hav=0;
520539c167fSBarry Smith         PetscCDIntNd *pos;
521539c167fSBarry Smith 
5220cbbd2e1SMark F. Adams         /* search for this gid to see if I have it */
5239566063dSJacob Faibussowitsch         PetscCall(PetscCDGetHeadPos(aggs_2,slid_new,&pos));
524e78576d6SMark F. Adams         while (pos) {
525e78576d6SMark F. Adams           PetscInt gidj;
5269566063dSJacob Faibussowitsch           PetscCall(PetscCDIntNdGetID(pos, &gidj));
5279566063dSJacob Faibussowitsch           PetscCall(PetscCDGetNextPos(aggs_2,slid_new,&pos));
528e78576d6SMark F. Adams 
5290cbbd2e1SMark F. Adams           if (gidj == gid) { hav = 1; break; }
530c8b0795cSMark F. Adams         }
531c8b0795cSMark F. Adams         if (hav != 1) {
532ffc955d6SMark F. Adams           /* insert 'flidj' into head of llist */
5339566063dSJacob Faibussowitsch           PetscCall(PetscCDAppendID(aggs_2, slid_new, gid));
534c8b0795cSMark F. Adams         }
535c8b0795cSMark F. Adams       }
536c8b0795cSMark F. Adams     }
537c8b0795cSMark F. Adams 
5389566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(mpimat_1->lvec, &cpcol_1_state));
5399566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(mpimat_2->lvec, &cpcol_2_state));
5409566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(ghostparents2, &cpcol_2_parent));
5419566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(ghostgids2, &cpcol_2_gid));
5429566063dSJacob Faibussowitsch     PetscCall(PetscFree(lid_cprowID_1));
5439566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&ghostgids2));
5449566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&ghostparents2));
5459566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&ghost_par_orig2));
546c8b0795cSMark F. Adams   }
547c8b0795cSMark F. Adams 
5489566063dSJacob Faibussowitsch   PetscCall(PetscFree2(lid_state,lid_parent_gid));
549c8b0795cSMark F. Adams   PetscFunctionReturn(0);
550c8b0795cSMark F. Adams }
5512e68589bSMark F. Adams 
5522e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
5532e68589bSMark F. Adams /*
554a2f3521dSMark F. Adams    PCSetData_AGG - called if data is not set with PCSetCoordinates.
555a2f3521dSMark F. Adams       Looks in Mat for near null space.
556a2f3521dSMark F. Adams       Does not work for Stokes
5572e68589bSMark F. Adams 
5582e68589bSMark F. Adams   Input Parameter:
5592e68589bSMark F. Adams    . pc -
560a2f3521dSMark F. Adams    . a_A - matrix to get (near) null space out of.
5612e68589bSMark F. Adams */
562e0877f53SBarry Smith static PetscErrorCode PCSetData_AGG(PC pc, Mat a_A)
5632e68589bSMark F. Adams {
564b8cd405aSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
565b8cd405aSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
566b8cd405aSMark F. Adams   MatNullSpace   mnull;
56766f2ef4bSMark Adams 
568b817416eSBarry Smith   PetscFunctionBegin;
5699566063dSJacob Faibussowitsch   PetscCall(MatGetNearNullSpace(a_A, &mnull));
570b8cd405aSMark F. Adams   if (!mnull) {
57166f2ef4bSMark Adams     DM dm;
5729566063dSJacob Faibussowitsch     PetscCall(PCGetDM(pc, &dm));
57366f2ef4bSMark Adams     if (!dm) {
5749566063dSJacob Faibussowitsch       PetscCall(MatGetDM(a_A, &dm));
57566f2ef4bSMark Adams     }
57666f2ef4bSMark Adams     if (dm) {
57766f2ef4bSMark Adams       PetscObject deformation;
578b0d30dd6SMatthew G. Knepley       PetscInt    Nf;
579b0d30dd6SMatthew G. Knepley 
5809566063dSJacob Faibussowitsch       PetscCall(DMGetNumFields(dm, &Nf));
581b0d30dd6SMatthew G. Knepley       if (Nf) {
5829566063dSJacob Faibussowitsch         PetscCall(DMGetField(dm, 0, NULL, &deformation));
5839566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery((PetscObject)deformation,"nearnullspace",(PetscObject*)&mnull));
58466f2ef4bSMark Adams         if (!mnull) {
5859566063dSJacob Faibussowitsch           PetscCall(PetscObjectQuery((PetscObject)deformation,"nullspace",(PetscObject*)&mnull));
58666f2ef4bSMark Adams         }
58766f2ef4bSMark Adams       }
58866f2ef4bSMark Adams     }
589b0d30dd6SMatthew G. Knepley   }
59066f2ef4bSMark Adams 
59166f2ef4bSMark Adams   if (!mnull) {
592a2f3521dSMark F. Adams     PetscInt bs,NN,MM;
5939566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSize(a_A, &bs));
5949566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(a_A, &MM, &NN));
59563a3b9bcSJacob Faibussowitsch     PetscCheck(MM % bs == 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MM %" PetscInt_FMT " must be divisible by bs %" PetscInt_FMT,MM,bs);
5969566063dSJacob Faibussowitsch     PetscCall(PCSetCoordinates_AGG(pc, bs, MM/bs, NULL));
597806fa848SBarry Smith   } else {
598b8cd405aSMark F. Adams     PetscReal         *nullvec;
599b8cd405aSMark F. Adams     PetscBool         has_const;
600b8cd405aSMark F. Adams     PetscInt          i,j,mlocal,nvec,bs;
601d5d25401SBarry Smith     const Vec         *vecs;
602d5d25401SBarry Smith     const PetscScalar *v;
603b817416eSBarry Smith 
6049566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(a_A,&mlocal,NULL));
6059566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs));
606a0dea87dSMark Adams     pc_gamg->data_sz = (nvec+!!has_const)*mlocal;
6079566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1((nvec+!!has_const)*mlocal,&nullvec));
608b8cd405aSMark F. Adams     if (has_const) for (i=0; i<mlocal; i++) nullvec[i] = 1.0;
609b8cd405aSMark F. Adams     for (i=0; i<nvec; i++) {
6109566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(vecs[i],&v));
611b8cd405aSMark F. Adams       for (j=0; j<mlocal; j++) nullvec[(i+!!has_const)*mlocal + j] = PetscRealPart(v[j]);
6129566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(vecs[i],&v));
613b8cd405aSMark F. Adams     }
614b8cd405aSMark F. Adams     pc_gamg->data           = nullvec;
615b8cd405aSMark F. Adams     pc_gamg->data_cell_cols = (nvec+!!has_const);
6169566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSize(a_A, &bs));
617b8cd405aSMark F. Adams     pc_gamg->data_cell_rows = bs;
618b8cd405aSMark F. Adams   }
6192e68589bSMark F. Adams   PetscFunctionReturn(0);
6202e68589bSMark F. Adams }
6212e68589bSMark F. Adams 
6222e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
6232e68589bSMark F. Adams /*
6242e68589bSMark F. Adams  formProl0
6252e68589bSMark F. Adams 
6262e68589bSMark F. Adams    Input Parameter:
6272adfe9d3SBarry Smith    . agg_llists - list of arrays with aggregates -- list from selected vertices of aggregate unselected vertices
6282adfe9d3SBarry Smith    . bs - row block size
6290cbbd2e1SMark F. Adams    . nSAvec - column bs of new P
6300cbbd2e1SMark F. Adams    . my0crs - global index of start of locals
6312adfe9d3SBarry Smith    . data_stride - bs*(nloc nodes + ghost nodes) [data_stride][nSAvec]
6322e68589bSMark F. Adams    . data_in[data_stride*nSAvec] - local data on fine grid
6332e68589bSMark F. Adams    . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist'
6342e68589bSMark F. Adams   Output Parameter:
6352e68589bSMark F. Adams    . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data
6362e68589bSMark F. Adams    . a_Prol - prolongation operator
6372e68589bSMark F. Adams */
6382adfe9d3SBarry 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)
6392e68589bSMark F. Adams {
640bd026e97SJed Brown   PetscInt        Istart,my0,Iend,nloc,clid,flid = 0,aggID,kk,jj,ii,mm,nSelected,minsz,nghosts,out_data_stride;
6413b4367a7SBarry Smith   MPI_Comm        comm;
6422e68589bSMark F. Adams   PetscReal       *out_data;
643539c167fSBarry Smith   PetscCDIntNd    *pos;
6441943db53SBarry Smith   PCGAMGHashTable fgid_flid;
6450cbbd2e1SMark F. Adams 
6462e68589bSMark F. Adams   PetscFunctionBegin;
6479566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)a_Prol,&comm));
6489566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(a_Prol, &Istart, &Iend));
64971959b99SBarry Smith   nloc = (Iend-Istart)/bs; my0 = Istart/bs;
65063a3b9bcSJacob 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);
6510cbbd2e1SMark F. Adams   Iend   /= bs;
6520cbbd2e1SMark F. Adams   nghosts = data_stride/bs - nloc;
6532e68589bSMark F. Adams 
6549566063dSJacob Faibussowitsch   PetscCall(PCGAMGHashTableCreate(2*nghosts+1, &fgid_flid));
6550cbbd2e1SMark F. Adams   for (kk=0; kk<nghosts; kk++) {
6569566063dSJacob Faibussowitsch     PetscCall(PCGAMGHashTableAdd(&fgid_flid, flid_fgid[nloc+kk], nloc+kk));
6572e68589bSMark F. Adams   }
6582e68589bSMark F. Adams 
6590cbbd2e1SMark F. Adams   /* count selected -- same as number of cols of P */
6600cbbd2e1SMark F. Adams   for (nSelected=mm=0; mm<nloc; mm++) {
661e78576d6SMark F. Adams     PetscBool ise;
6629566063dSJacob Faibussowitsch     PetscCall(PetscCDEmptyAt(agg_llists, mm, &ise));
663e78576d6SMark F. Adams     if (!ise) nSelected++;
6640cbbd2e1SMark F. Adams   }
6659566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRangeColumn(a_Prol, &ii, &jj));
66663a3b9bcSJacob Faibussowitsch   PetscCheck((ii/nSAvec) == my0crs,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %" PetscInt_FMT " /nSAvec %" PetscInt_FMT "  != my0crs %" PetscInt_FMT,ii,nSAvec,my0crs);
66763a3b9bcSJacob 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);
6680cbbd2e1SMark F. Adams 
6692e68589bSMark F. Adams   /* aloc space for coarse point data (output) */
6700cbbd2e1SMark F. Adams   out_data_stride = nSelected*nSAvec;
6712fa5cd67SKarl Rupp 
6729566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(out_data_stride*nSAvec, &out_data));
67333812677SSatish Balay   for (ii=0;ii<out_data_stride*nSAvec;ii++) out_data[ii]=PETSC_MAX_REAL;
6740cbbd2e1SMark F. Adams   *a_data_out = out_data; /* output - stride nSelected*nSAvec */
6752e68589bSMark F. Adams 
6762e68589bSMark F. Adams   /* find points and set prolongation */
677c8b0795cSMark F. Adams   minsz = 100;
6780cbbd2e1SMark F. Adams   for (mm = clid = 0; mm < nloc; mm++) {
6799566063dSJacob Faibussowitsch     PetscCall(PetscCDSizeAt(agg_llists, mm, &jj));
680e78576d6SMark F. Adams     if (jj > 0) {
6810cbbd2e1SMark F. Adams       const PetscInt lid = mm, cgid = my0crs + clid;
6820cbbd2e1SMark F. Adams       PetscInt       cids[100]; /* max bs */
6830cbbd2e1SMark F. Adams       PetscBLASInt   asz  =jj,M=asz*bs,N=nSAvec,INFO;
6842e68589bSMark F. Adams       PetscBLASInt   Mdata=M+((N-M>0) ? N-M : 0),LDA=Mdata,LWORK=N*bs;
6852e68589bSMark F. Adams       PetscScalar    *qqc,*qqr,*TAU,*WORK;
6862e68589bSMark F. Adams       PetscInt       *fids;
68765d7b583SSatish Balay       PetscReal      *data;
688b817416eSBarry Smith 
6890cbbd2e1SMark F. Adams       /* count agg */
6900cbbd2e1SMark F. Adams       if (asz<minsz) minsz = asz;
6910cbbd2e1SMark F. Adams 
6920cbbd2e1SMark F. Adams       /* get block */
6939566063dSJacob Faibussowitsch       PetscCall(PetscMalloc5(Mdata*N, &qqc,M*N, &qqr,N, &TAU,LWORK, &WORK,M, &fids));
6942e68589bSMark F. Adams 
6952e68589bSMark F. Adams       aggID = 0;
6969566063dSJacob Faibussowitsch       PetscCall(PetscCDGetHeadPos(agg_llists,lid,&pos));
697e78576d6SMark F. Adams       while (pos) {
698e78576d6SMark F. Adams         PetscInt gid1;
6999566063dSJacob Faibussowitsch         PetscCall(PetscCDIntNdGetID(pos, &gid1));
7009566063dSJacob Faibussowitsch         PetscCall(PetscCDGetNextPos(agg_llists,lid,&pos));
701e78576d6SMark F. Adams 
7020cbbd2e1SMark F. Adams         if (gid1 >= my0 && gid1 < Iend) flid = gid1 - my0;
7030cbbd2e1SMark F. Adams         else {
7049566063dSJacob Faibussowitsch           PetscCall(PCGAMGHashTableFind(&fgid_flid, gid1, &flid));
70508401ef6SPierre Jolivet           PetscCheck(flid >= 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot find gid1 in table");
7060cbbd2e1SMark F. Adams         }
7072e68589bSMark F. Adams         /* copy in B_i matrix - column oriented */
70865d7b583SSatish Balay         data = &data_in[flid*bs];
7093d3eaba7SBarry Smith         for (ii = 0; ii < bs; ii++) {
7102e68589bSMark F. Adams           for (jj = 0; jj < N; jj++) {
7110cbbd2e1SMark F. Adams             PetscReal d = data[jj*data_stride + ii];
7120cbbd2e1SMark F. Adams             qqc[jj*Mdata + aggID*bs + ii] = d;
7132e68589bSMark F. Adams           }
7142e68589bSMark F. Adams         }
7152e68589bSMark F. Adams         /* set fine IDs */
7162e68589bSMark F. Adams         for (kk=0; kk<bs; kk++) fids[aggID*bs + kk] = flid_fgid[flid]*bs + kk;
7172e68589bSMark F. Adams         aggID++;
7180cbbd2e1SMark F. Adams       }
7192e68589bSMark F. Adams 
7202e68589bSMark F. Adams       /* pad with zeros */
7212e68589bSMark F. Adams       for (ii = asz*bs; ii < Mdata; ii++) {
7222e68589bSMark F. Adams         for (jj = 0; jj < N; jj++, kk++) {
7232e68589bSMark F. Adams           qqc[jj*Mdata + ii] = .0;
7242e68589bSMark F. Adams         }
7252e68589bSMark F. Adams       }
7262e68589bSMark F. Adams 
7272e68589bSMark F. Adams       /* QR */
7289566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
7298b83055fSJed Brown       PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
7309566063dSJacob Faibussowitsch       PetscCall(PetscFPTrapPop());
73108401ef6SPierre Jolivet       PetscCheck(INFO == 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"xGEQRF error");
7322e68589bSMark F. Adams       /* get R - column oriented - output B_{i+1} */
7332e68589bSMark F. Adams       {
7342e68589bSMark F. Adams         PetscReal *data = &out_data[clid*nSAvec];
7352e68589bSMark F. Adams         for (jj = 0; jj < nSAvec; jj++) {
7362e68589bSMark F. Adams           for (ii = 0; ii < nSAvec; ii++) {
73708401ef6SPierre 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);
7380cbbd2e1SMark F. Adams            if (ii <= jj) data[jj*out_data_stride + ii] = PetscRealPart(qqc[jj*Mdata + ii]);
7390cbbd2e1SMark F. Adams            else data[jj*out_data_stride + ii] = 0.;
7402e68589bSMark F. Adams           }
7412e68589bSMark F. Adams         }
7422e68589bSMark F. Adams       }
7432e68589bSMark F. Adams 
7442e68589bSMark F. Adams       /* get Q - row oriented */
745c964aadfSJose E. Roman       PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
74663a3b9bcSJacob Faibussowitsch       PetscCheck(INFO == 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"xORGQR error arg %" PetscBLASInt_FMT,-INFO);
7472e68589bSMark F. Adams 
7482e68589bSMark F. Adams       for (ii = 0; ii < M; ii++) {
7492e68589bSMark F. Adams         for (jj = 0; jj < N; jj++) {
7502e68589bSMark F. Adams           qqr[N*ii + jj] = qqc[jj*Mdata + ii];
7512e68589bSMark F. Adams         }
7522e68589bSMark F. Adams       }
7532e68589bSMark F. Adams 
7542e68589bSMark F. Adams       /* add diagonal block of P0 */
755c8b0795cSMark F. Adams       for (kk=0; kk<N; kk++) {
756c8b0795cSMark F. Adams         cids[kk] = N*cgid + kk; /* global col IDs in P0 */
757c8b0795cSMark F. Adams       }
7589566063dSJacob Faibussowitsch       PetscCall(MatSetValues(a_Prol,M,fids,N,cids,qqr,INSERT_VALUES));
7599566063dSJacob Faibussowitsch       PetscCall(PetscFree5(qqc,qqr,TAU,WORK,fids));
760b43b03e9SMark F. Adams       clid++;
7610cbbd2e1SMark F. Adams     } /* coarse agg */
7620cbbd2e1SMark F. Adams   } /* for all fine nodes */
7639566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY));
7649566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY));
7659566063dSJacob Faibussowitsch   PetscCall(PCGAMGHashTableDestroy(&fgid_flid));
7662e68589bSMark F. Adams   PetscFunctionReturn(0);
7672e68589bSMark F. Adams }
7682e68589bSMark F. Adams 
7695adeb434SBarry Smith static PetscErrorCode PCView_GAMG_AGG(PC pc,PetscViewer viewer)
7705adeb434SBarry Smith {
7715adeb434SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
7725adeb434SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
7735adeb434SBarry Smith   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
7745adeb434SBarry Smith 
7755adeb434SBarry Smith   PetscFunctionBegin;
7769566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"      AGG specific options\n"));
7779566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"        Symmetric graph %s\n",pc_gamg_agg->sym_graph ? "true" : "false"));
77863a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"        Number of levels to square graph %" PetscInt_FMT "\n",pc_gamg_agg->square_graph));
77963a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"        Number smoothing steps %" PetscInt_FMT "\n",pc_gamg_agg->nsmooths));
7805adeb434SBarry Smith   PetscFunctionReturn(0);
7815adeb434SBarry Smith }
7825adeb434SBarry Smith 
7832e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
7842e68589bSMark F. Adams /*
785fd1112cbSBarry Smith    PCGAMGGraph_AGG
7862e68589bSMark F. Adams 
7872e68589bSMark F. Adams   Input Parameter:
7882e68589bSMark F. Adams    . pc - this
7892e68589bSMark F. Adams    . Amat - matrix on this fine level
7902e68589bSMark F. Adams   Output Parameter:
791c8b0795cSMark F. Adams    . a_Gmat -
7922e68589bSMark F. Adams */
793e0877f53SBarry Smith static PetscErrorCode PCGAMGGraph_AGG(PC pc,Mat Amat,Mat *a_Gmat)
794c8b0795cSMark F. Adams {
795c8b0795cSMark F. Adams   PC_MG                     *mg          = (PC_MG*)pc->data;
796c8b0795cSMark F. Adams   PC_GAMG                   *pc_gamg     = (PC_GAMG*)mg->innerctx;
797c1eae691SMark Adams   const PetscReal           vfilter      = pc_gamg->threshold[pc_gamg->current_level];
798c8b0795cSMark F. Adams   PC_GAMG_AGG               *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
799e0940f08SMark F. Adams   Mat                       Gmat;
8003b4367a7SBarry Smith   MPI_Comm                  comm;
80167744fedSMark F. Adams   PetscBool /* set,flg , */ symm;
802c8b0795cSMark F. Adams 
803c8b0795cSMark F. Adams   PetscFunctionBegin;
8049566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Amat,&comm));
805c8b0795cSMark F. Adams 
8069566063dSJacob Faibussowitsch   /* PetscCall(MatIsSymmetricKnown(Amat, &set, &flg)); || !(set && flg) -- this causes lot of symm calls */
807c666adbfSMark F. Adams   symm = (PetscBool)(pc_gamg_agg->sym_graph); /* && !pc_gamg_agg->square_graph; */
8080cbbd2e1SMark F. Adams 
809849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_GRAPH],0,0,0,0));
8105d5a09d9SMark Adams   PetscCall(PCGAMGCreateGraph(Amat, &Gmat, symm));
811849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_GRAPH],0,0,0,0));
812849bee69SMark Adams 
813849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_FILTER],0,0,0,0));
8145d5a09d9SMark Adams   PetscCall(PCGAMGFilterGraph(&Gmat, vfilter));
815849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_FILTER],0,0,0,0));
816849bee69SMark Adams 
817e0940f08SMark F. Adams   *a_Gmat = Gmat;
818849bee69SMark Adams 
819c8b0795cSMark F. Adams   PetscFunctionReturn(0);
820c8b0795cSMark F. Adams }
821c8b0795cSMark F. Adams 
822c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
823c8b0795cSMark F. Adams /*
824b43b03e9SMark F. Adams    PCGAMGCoarsen_AGG
825c8b0795cSMark F. Adams 
826c8b0795cSMark F. Adams   Input Parameter:
827e0940f08SMark F. Adams    . a_pc - this
828e0940f08SMark F. Adams   Input/Output Parameter:
8290cbbd2e1SMark F. Adams    . a_Gmat1 - graph on this fine level - coarsening can change this (squares it)
830c8b0795cSMark F. Adams   Output Parameter:
8310cbbd2e1SMark F. Adams    . agg_lists - list of aggregates
832c8b0795cSMark F. Adams */
833e0877f53SBarry Smith static PetscErrorCode PCGAMGCoarsen_AGG(PC a_pc,Mat *a_Gmat1,PetscCoarsenData **agg_lists)
834c8b0795cSMark F. Adams {
835e0940f08SMark F. Adams   PC_MG          *mg          = (PC_MG*)a_pc->data;
836c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
837c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
8380cbbd2e1SMark F. Adams   Mat            mat,Gmat2, Gmat1 = *a_Gmat1;  /* squared graph */
8390cbbd2e1SMark F. Adams   IS             perm;
8405cfd4bd9SMark Adams   PetscInt       Istart,Iend,Ii,nloc,bs,n,m;
841c8b0795cSMark F. Adams   PetscInt       *permute;
842c8b0795cSMark F. Adams   PetscBool      *bIndexSet;
843b43b03e9SMark F. Adams   MatCoarsen     crs;
8443b4367a7SBarry Smith   MPI_Comm       comm;
845aea53286SMark Adams   PetscReal      hashfact;
846e2c00dcbSBarry Smith   PetscInt       iSwapIndex;
8473b50add6SMark Adams   PetscRandom    random;
848c8b0795cSMark F. Adams 
849c8b0795cSMark F. Adams   PetscFunctionBegin;
850849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_COARSEN],0,0,0,0));
8519566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Gmat1,&comm));
8529566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(Gmat1, &n, &m));
8539566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(Gmat1, &bs));
85463a3b9bcSJacob Faibussowitsch   PetscCheck(bs == 1,PETSC_COMM_SELF,PETSC_ERR_PLIB,"bs %" PetscInt_FMT " must be 1",bs);
855c8b0795cSMark F. Adams   nloc = n/bs;
856c8b0795cSMark F. Adams 
8579ab59c8bSMark Adams   if (pc_gamg->current_level < pc_gamg_agg->square_graph) {
858849bee69SMark Adams     PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_SQUARE],0,0,0,0));
8599566063dSJacob Faibussowitsch     PetscCall(PCGAMGSquareGraph_GAMG(a_pc,Gmat1,&Gmat2));
860849bee69SMark Adams     PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_SQUARE],0,0,0,0));
861806fa848SBarry Smith   } else Gmat2 = Gmat1;
862c8b0795cSMark F. Adams 
8635cfd4bd9SMark Adams   /* get MIS aggs - randomize */
8649566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nloc, &permute));
8659566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(nloc, &bIndexSet));
8666e09b0e3SStefano Zampini   for (Ii = 0; Ii < nloc; Ii++) permute[Ii] = Ii;
8679566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_SELF,&random));
8689566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Gmat1, &Istart, &Iend));
869c8b0795cSMark F. Adams   for (Ii = 0; Ii < nloc; Ii++) {
8709566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValueReal(random,&hashfact));
871aea53286SMark Adams     iSwapIndex = (PetscInt) (hashfact*nloc)%nloc;
872c8b0795cSMark F. Adams     if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) {
873c8b0795cSMark F. Adams       PetscInt iTemp = permute[iSwapIndex];
874c8b0795cSMark F. Adams       permute[iSwapIndex]   = permute[Ii];
875c8b0795cSMark F. Adams       permute[Ii]           = iTemp;
876c8b0795cSMark F. Adams       bIndexSet[iSwapIndex] = PETSC_TRUE;
877c8b0795cSMark F. Adams     }
878c8b0795cSMark F. Adams   }
8799566063dSJacob Faibussowitsch   PetscCall(PetscFree(bIndexSet));
8809566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&random));
8819566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm));
882849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_MIS],0,0,0,0));
8839566063dSJacob Faibussowitsch   PetscCall(MatCoarsenCreate(comm, &crs));
8849566063dSJacob Faibussowitsch   PetscCall(MatCoarsenSetFromOptions(crs));
8859566063dSJacob Faibussowitsch   PetscCall(MatCoarsenSetGreedyOrdering(crs, perm));
8869566063dSJacob Faibussowitsch   PetscCall(MatCoarsenSetAdjacency(crs, Gmat2));
8879566063dSJacob Faibussowitsch   PetscCall(MatCoarsenSetStrictAggs(crs, PETSC_TRUE));
8889566063dSJacob Faibussowitsch   PetscCall(MatCoarsenApply(crs));
8899566063dSJacob Faibussowitsch   PetscCall(MatCoarsenGetData(crs, agg_lists)); /* output */
8909566063dSJacob Faibussowitsch   PetscCall(MatCoarsenDestroy(&crs));
891b43b03e9SMark F. Adams 
8929566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
8939566063dSJacob Faibussowitsch   PetscCall(PetscFree(permute));
894849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_MIS],0,0,0,0));
8959f3f12c8SMark F. Adams 
896c8b0795cSMark F. Adams   /* smooth aggs */
897e0940f08SMark F. Adams   if (Gmat2 != Gmat1) {
8980cbbd2e1SMark F. Adams     const PetscCoarsenData *llist = *agg_lists;
8999566063dSJacob Faibussowitsch     PetscCall(smoothAggs(a_pc,Gmat2, Gmat1, *agg_lists));
9009566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Gmat1));
901e0940f08SMark F. Adams     *a_Gmat1 = Gmat2; /* output */
9029566063dSJacob Faibussowitsch     PetscCall(PetscCDGetMat(llist, &mat));
90328b400f6SJacob Faibussowitsch     PetscCheck(!mat,comm,PETSC_ERR_ARG_WRONG, "Auxilary matrix with squared graph????");
904806fa848SBarry Smith   } else {
9050cbbd2e1SMark F. Adams     const PetscCoarsenData *llist = *agg_lists;
9069ab59c8bSMark Adams     /* see if we have a matrix that takes precedence (returned from MatCoarsenApply) */
9079566063dSJacob Faibussowitsch     PetscCall(PetscCDGetMat(llist, &mat));
9080cbbd2e1SMark F. Adams     if (mat) {
9099566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&Gmat1));
9100cbbd2e1SMark F. Adams       *a_Gmat1 = mat; /* output */
9110cbbd2e1SMark F. Adams     }
9120cbbd2e1SMark F. Adams   }
913849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_COARSEN],0,0,0,0));
914c8b0795cSMark F. Adams   PetscFunctionReturn(0);
915c8b0795cSMark F. Adams }
916c8b0795cSMark F. Adams 
917c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
918c8b0795cSMark F. Adams /*
9190cbbd2e1SMark F. Adams  PCGAMGProlongator_AGG
920c8b0795cSMark F. Adams 
921c8b0795cSMark F. Adams  Input Parameter:
922c8b0795cSMark F. Adams  . pc - this
923c8b0795cSMark F. Adams  . Amat - matrix on this fine level
924c8b0795cSMark F. Adams  . Graph - used to get ghost data for nodes in
9250cbbd2e1SMark F. Adams  . agg_lists - list of aggregates
926c8b0795cSMark F. Adams  Output Parameter:
927c8b0795cSMark F. Adams  . a_P_out - prolongation operator to the next level
928c8b0795cSMark F. Adams  */
929e0877f53SBarry Smith static PetscErrorCode PCGAMGProlongator_AGG(PC pc,Mat Amat,Mat Gmat,PetscCoarsenData *agg_lists,Mat *a_P_out)
9302e68589bSMark F. Adams {
9312e68589bSMark F. Adams   PC_MG          *mg       = (PC_MG*)pc->data;
9322e68589bSMark F. Adams   PC_GAMG        *pc_gamg  = (PC_GAMG*)mg->innerctx;
9334a2f8832SBarry Smith   const PetscInt col_bs = pc_gamg->data_cell_cols;
934c8b0795cSMark F. Adams   PetscInt       Istart,Iend,nloc,ii,jj,kk,my0,nLocalSelected,bs;
935c8b0795cSMark F. Adams   Mat            Prol;
936d5d25401SBarry Smith   PetscMPIInt    size;
9373b4367a7SBarry Smith   MPI_Comm       comm;
938c8b0795cSMark F. Adams   PetscReal      *data_w_ghost;
939c8b0795cSMark F. Adams   PetscInt       myCrs0, nbnodes=0, *flid_fgid;
940d9558ea9SBarry Smith   MatType        mtype;
9412e68589bSMark F. Adams 
9422e68589bSMark F. Adams   PetscFunctionBegin;
9439566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Amat,&comm));
94408401ef6SPierre Jolivet   PetscCheck(col_bs >= 1,comm,PETSC_ERR_PLIB,"Column bs cannot be less than 1");
945849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROL],0,0,0,0));
9469566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
9479566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Amat, &Istart, &Iend));
9489566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(Amat, &bs));
94971959b99SBarry Smith   nloc = (Iend-Istart)/bs; my0 = Istart/bs;
95063a3b9bcSJacob 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);
9512e68589bSMark F. Adams 
9522e68589bSMark F. Adams   /* get 'nLocalSelected' */
9530cbbd2e1SMark F. Adams   for (ii=0, nLocalSelected = 0; ii < nloc; ii++) {
954e78576d6SMark F. Adams     PetscBool ise;
9550cbbd2e1SMark F. Adams     /* filter out singletons 0 or 1? */
9569566063dSJacob Faibussowitsch     PetscCall(PetscCDEmptyAt(agg_lists, ii, &ise));
957e78576d6SMark F. Adams     if (!ise) nLocalSelected++;
9582e68589bSMark F. Adams   }
9592e68589bSMark F. Adams 
9602e68589bSMark F. Adams   /* create prolongator, create P matrix */
9619566063dSJacob Faibussowitsch   PetscCall(MatGetType(Amat,&mtype));
9629566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &Prol));
9639566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Prol,nloc*bs,nLocalSelected*col_bs,PETSC_DETERMINE,PETSC_DETERMINE));
9649566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSizes(Prol, bs, col_bs));
9659566063dSJacob Faibussowitsch   PetscCall(MatSetType(Prol, mtype));
9669566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(Prol,col_bs, NULL));
9679566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(Prol,col_bs, NULL, col_bs, NULL));
9682e68589bSMark F. Adams 
9692e68589bSMark F. Adams   /* can get all points "removed" */
9709566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Prol, &kk, &ii));
9717f66b68fSMark Adams   if (!ii) {
97263a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(pc,"%s: No selected points on coarse grid\n",((PetscObject)pc)->prefix));
9739566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Prol));
9740298fd71SBarry Smith     *a_P_out = NULL;  /* out */
975849bee69SMark Adams     PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROL],0,0,0,0));
9762e68589bSMark F. Adams     PetscFunctionReturn(0);
9772e68589bSMark F. Adams   }
97863a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(pc,"%s: New grid %" PetscInt_FMT " nodes\n",((PetscObject)pc)->prefix,ii/col_bs));
9799566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRangeColumn(Prol, &myCrs0, &kk));
9800cbbd2e1SMark F. Adams 
98163a3b9bcSJacob 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);
982c8b0795cSMark F. Adams   myCrs0 = myCrs0/col_bs;
98363a3b9bcSJacob 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);
9842e68589bSMark F. Adams 
9852e68589bSMark F. Adams   /* create global vector of data in 'data_w_ghost' */
986849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROLA],0,0,0,0));
987c5df96a5SBarry Smith   if (size > 1) { /*  */
9882e68589bSMark F. Adams     PetscReal *tmp_gdata,*tmp_ldata,*tp2;
9899566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nloc, &tmp_ldata));
9904a2f8832SBarry Smith     for (jj = 0; jj < col_bs; jj++) {
991c8b0795cSMark F. Adams       for (kk = 0; kk < bs; kk++) {
992a2f3521dSMark F. Adams         PetscInt        ii,stride;
993c8b0795cSMark F. Adams         const PetscReal *tp = pc_gamg->data + jj*bs*nloc + kk;
9942fa5cd67SKarl Rupp         for (ii = 0; ii < nloc; ii++, tp += bs) tmp_ldata[ii] = *tp;
9952fa5cd67SKarl Rupp 
9969566063dSJacob Faibussowitsch         PetscCall(PCGAMGGetDataWithGhosts(Gmat, 1, tmp_ldata, &stride, &tmp_gdata));
997a2f3521dSMark F. Adams 
9987f66b68fSMark Adams         if (!jj && !kk) { /* now I know how many todal nodes - allocate */
9999566063dSJacob Faibussowitsch           PetscCall(PetscMalloc1(stride*bs*col_bs, &data_w_ghost));
1000a2f3521dSMark F. Adams           nbnodes = bs*stride;
10012e68589bSMark F. Adams         }
1002a2f3521dSMark F. Adams         tp2 = data_w_ghost + jj*bs*stride + kk;
10032fa5cd67SKarl Rupp         for (ii = 0; ii < stride; ii++, tp2 += bs) *tp2 = tmp_gdata[ii];
10049566063dSJacob Faibussowitsch         PetscCall(PetscFree(tmp_gdata));
10052e68589bSMark F. Adams       }
10062e68589bSMark F. Adams     }
10079566063dSJacob Faibussowitsch     PetscCall(PetscFree(tmp_ldata));
1008806fa848SBarry Smith   } else {
1009c8b0795cSMark F. Adams     nbnodes      = bs*nloc;
1010c8b0795cSMark F. Adams     data_w_ghost = (PetscReal*)pc_gamg->data;
10112e68589bSMark F. Adams   }
10122e68589bSMark F. Adams 
10132e68589bSMark F. Adams   /* get P0 */
1014c5df96a5SBarry Smith   if (size > 1) {
10152e68589bSMark F. Adams     PetscReal *fid_glid_loc,*fiddata;
1016a2f3521dSMark F. Adams     PetscInt  stride;
10172e68589bSMark F. Adams 
10189566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nloc, &fid_glid_loc));
10192e68589bSMark F. Adams     for (kk=0; kk<nloc; kk++) fid_glid_loc[kk] = (PetscReal)(my0+kk);
10209566063dSJacob Faibussowitsch     PetscCall(PCGAMGGetDataWithGhosts(Gmat, 1, fid_glid_loc, &stride, &fiddata));
10219566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(stride, &flid_fgid));
1022a2f3521dSMark F. Adams     for (kk=0; kk<stride; kk++) flid_fgid[kk] = (PetscInt)fiddata[kk];
10239566063dSJacob Faibussowitsch     PetscCall(PetscFree(fiddata));
1024a2f3521dSMark F. Adams 
102563a3b9bcSJacob Faibussowitsch     PetscCheck(stride == nbnodes/bs,PETSC_COMM_SELF,PETSC_ERR_PLIB,"stride %" PetscInt_FMT " != nbnodes %" PetscInt_FMT "/bs %" PetscInt_FMT,stride,nbnodes,bs);
10269566063dSJacob Faibussowitsch     PetscCall(PetscFree(fid_glid_loc));
1027806fa848SBarry Smith   } else {
10289566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nloc, &flid_fgid));
10292e68589bSMark F. Adams     for (kk=0; kk<nloc; kk++) flid_fgid[kk] = my0 + kk;
10302e68589bSMark F. Adams   }
1031849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROLA],0,0,0,0));
1032849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_PROLB],0,0,0,0));
1033c8b0795cSMark F. Adams   {
10340298fd71SBarry Smith     PetscReal *data_out = NULL;
10359566063dSJacob Faibussowitsch     PetscCall(formProl0(agg_lists, bs, col_bs, myCrs0, nbnodes,data_w_ghost, flid_fgid, &data_out, Prol));
10369566063dSJacob Faibussowitsch     PetscCall(PetscFree(pc_gamg->data));
10372fa5cd67SKarl Rupp 
1038c8b0795cSMark F. Adams     pc_gamg->data           = data_out;
10394a2f8832SBarry Smith     pc_gamg->data_cell_rows = col_bs;
10404a2f8832SBarry Smith     pc_gamg->data_sz        = col_bs*col_bs*nLocalSelected;
1041c8b0795cSMark F. Adams   }
1042849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROLB],0,0,0,0));
1043849bee69SMark Adams   if (size > 1) {PetscCall(PetscFree(data_w_ghost));}
10449566063dSJacob Faibussowitsch   PetscCall(PetscFree(flid_fgid));
10452e68589bSMark F. Adams 
1046c8b0795cSMark F. Adams   *a_P_out = Prol;  /* out */
1047ed4e82eaSPeter Brune 
1048849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_PROL],0,0,0,0));
1049c8b0795cSMark F. Adams   PetscFunctionReturn(0);
1050c8b0795cSMark F. Adams }
1051c8b0795cSMark F. Adams 
1052c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
1053c8b0795cSMark F. Adams /*
1054fd1112cbSBarry Smith    PCGAMGOptProlongator_AGG
1055c8b0795cSMark F. Adams 
1056c8b0795cSMark F. Adams   Input Parameter:
1057c8b0795cSMark F. Adams    . pc - this
1058c8b0795cSMark F. Adams    . Amat - matrix on this fine level
1059c8b0795cSMark F. Adams  In/Output Parameter:
10602a808120SBarry Smith    . a_P - prolongation operator to the next level
1061c8b0795cSMark F. Adams */
1062e0877f53SBarry Smith static PetscErrorCode PCGAMGOptProlongator_AGG(PC pc,Mat Amat,Mat *a_P)
1063c8b0795cSMark F. Adams {
1064c8b0795cSMark F. Adams   PC_MG          *mg          = (PC_MG*)pc->data;
1065c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
1066c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
1067c8b0795cSMark F. Adams   PetscInt       jj;
1068c8b0795cSMark F. Adams   Mat            Prol  = *a_P;
10693b4367a7SBarry Smith   MPI_Comm       comm;
10702a808120SBarry Smith   KSP            eksp;
10712a808120SBarry Smith   Vec            bb, xx;
10722a808120SBarry Smith   PC             epc;
10732a808120SBarry Smith   PetscReal      alpha, emax, emin;
1074c8b0795cSMark F. Adams 
1075c8b0795cSMark F. Adams   PetscFunctionBegin;
10769566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Amat,&comm));
1077849bee69SMark Adams   PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_OPT],0,0,0,0));
1078c8b0795cSMark F. Adams 
1079d5d25401SBarry Smith   /* compute maximum singular value of operator to be used in smoother */
10802a808120SBarry Smith   if (0 < pc_gamg_agg->nsmooths) {
108118c3aa7eSMark     /* get eigen estimates */
108218c3aa7eSMark     if (pc_gamg->emax > 0) {
108318c3aa7eSMark       emin = pc_gamg->emin;
108418c3aa7eSMark       emax = pc_gamg->emax;
108518c3aa7eSMark     } else {
10860ed2132dSStefano Zampini       const char *prefix;
10870ed2132dSStefano Zampini 
10889566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(Amat, &bb, NULL));
10899566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(Amat, &xx, NULL));
10909566063dSJacob Faibussowitsch       PetscCall(KSPSetNoisy_Private(bb));
1091e696ed0bSMark F. Adams 
10929566063dSJacob Faibussowitsch       PetscCall(KSPCreate(comm,&eksp));
10939566063dSJacob Faibussowitsch       PetscCall(PCGetOptionsPrefix(pc,&prefix));
10949566063dSJacob Faibussowitsch       PetscCall(KSPSetOptionsPrefix(eksp,prefix));
10959566063dSJacob Faibussowitsch       PetscCall(KSPAppendOptionsPrefix(eksp,"pc_gamg_esteig_"));
109673f7197eSJed Brown       {
1097acbb5b45SJed Brown         PetscBool sflg;
10989566063dSJacob Faibussowitsch         PetscCall(MatGetOption(Amat, MAT_SPD, &sflg));
109990db8557SMark Adams         if (sflg) {
11009566063dSJacob Faibussowitsch           PetscCall(KSPSetType(eksp, KSPCG));
1101d24ecf33SMark         }
1102d24ecf33SMark       }
11039566063dSJacob Faibussowitsch       PetscCall(KSPSetErrorIfNotConverged(eksp,pc->erroriffailure));
11049566063dSJacob Faibussowitsch       PetscCall(KSPSetNormType(eksp, KSP_NORM_NONE));
1105c2436cacSMark F. Adams 
11069566063dSJacob Faibussowitsch       PetscCall(KSPSetInitialGuessNonzero(eksp, PETSC_FALSE));
11079566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(eksp, Amat, Amat));
11082e68589bSMark F. Adams 
11099566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(eksp, &epc));
11109566063dSJacob Faibussowitsch       PetscCall(PCSetType(epc, PCJACOBI));  /* smoother in smoothed agg. */
1111c2436cacSMark F. Adams 
11129566063dSJacob 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
1113074aec55SMark Adams 
11149566063dSJacob Faibussowitsch       PetscCall(KSPSetFromOptions(eksp));
11159566063dSJacob Faibussowitsch       PetscCall(KSPSetComputeSingularValues(eksp,PETSC_TRUE));
11169566063dSJacob Faibussowitsch       PetscCall(KSPSolve(eksp, bb, xx));
11179566063dSJacob Faibussowitsch       PetscCall(KSPCheckSolve(eksp,pc,xx));
11182e68589bSMark F. Adams 
11199566063dSJacob Faibussowitsch       PetscCall(KSPComputeExtremeSingularValues(eksp, &emax, &emin));
11209566063dSJacob Faibussowitsch       PetscCall(PetscInfo(pc,"%s: Smooth P0: max eigen=%e min=%e PC=%s\n",((PetscObject)pc)->prefix,(double)emax,(double)emin,PCJACOBI));
11219566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&xx));
11229566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&bb));
11239566063dSJacob Faibussowitsch       PetscCall(KSPDestroy(&eksp));
11242e68589bSMark F. Adams     }
112518c3aa7eSMark     if (pc_gamg->use_sa_esteig) {
112618c3aa7eSMark       mg->min_eigen_DinvA[pc_gamg->current_level] = emin;
112718c3aa7eSMark       mg->max_eigen_DinvA[pc_gamg->current_level] = emax;
112863a3b9bcSJacob 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));
112918c3aa7eSMark     } else {
113018c3aa7eSMark       mg->min_eigen_DinvA[pc_gamg->current_level] = 0;
113118c3aa7eSMark       mg->max_eigen_DinvA[pc_gamg->current_level] = 0;
113218c3aa7eSMark     }
113318c3aa7eSMark   } else {
113418c3aa7eSMark     mg->min_eigen_DinvA[pc_gamg->current_level] = 0;
113518c3aa7eSMark     mg->max_eigen_DinvA[pc_gamg->current_level] = 0;
113618c3aa7eSMark   }
11372e68589bSMark F. Adams 
11382a808120SBarry Smith   /* smooth P0 */
11392a808120SBarry Smith   for (jj = 0; jj < pc_gamg_agg->nsmooths; jj++) {
11402a808120SBarry Smith     Mat tMat;
11412a808120SBarry Smith     Vec diag;
11422a808120SBarry Smith 
1143849bee69SMark Adams     PetscCall(PetscLogEventBegin(petsc_gamg_setup_events[GAMG_OPTSM],0,0,0,0));
11442a808120SBarry Smith 
11452e68589bSMark F. Adams     /* smooth P1 := (I - omega/lam D^{-1}A)P0 */
11469566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2],0,0,0,0));
11479566063dSJacob Faibussowitsch     PetscCall(MatMatMult(Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat));
11489566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2],0,0,0,0));
11499566063dSJacob Faibussowitsch     PetscCall(MatProductClear(tMat));
11509566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(Amat, &diag, NULL));
11519566063dSJacob Faibussowitsch     PetscCall(MatGetDiagonal(Amat, diag)); /* effectively PCJACOBI */
11529566063dSJacob Faibussowitsch     PetscCall(VecReciprocal(diag));
11539566063dSJacob Faibussowitsch     PetscCall(MatDiagonalScale(tMat, diag, NULL));
11549566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&diag));
1155b4da3a1bSStefano Zampini 
1156d5d25401SBarry Smith     /* TODO: Set a PCFailedReason and exit the building of the AMG preconditioner */
115708401ef6SPierre Jolivet     PetscCheck(emax != 0.0,PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"Computed maximum singular value as zero");
1158d5d25401SBarry Smith     /* TODO: Document the 1.4 and don't hardwire it in this routine */
1159b4ec6429SMark F. Adams     alpha = -1.4/emax;
1160b4da3a1bSStefano Zampini 
11619566063dSJacob Faibussowitsch     PetscCall(MatAYPX(tMat, alpha, Prol, SUBSET_NONZERO_PATTERN));
11629566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Prol));
11632e68589bSMark F. Adams     Prol = tMat;
1164849bee69SMark Adams     PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_OPTSM],0,0,0,0));
11652e68589bSMark F. Adams   }
1166849bee69SMark Adams   PetscCall(PetscLogEventEnd(petsc_gamg_setup_events[GAMG_OPT],0,0,0,0));
1167c8b0795cSMark F. Adams   *a_P = Prol;
1168c8b0795cSMark F. Adams   PetscFunctionReturn(0);
11692e68589bSMark F. Adams }
11702e68589bSMark F. Adams 
1171c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
1172c8b0795cSMark F. Adams /*
1173c8b0795cSMark F. Adams    PCCreateGAMG_AGG
11742e68589bSMark F. Adams 
1175c8b0795cSMark F. Adams   Input Parameter:
1176c8b0795cSMark F. Adams    . pc -
1177c8b0795cSMark F. Adams */
1178c8b0795cSMark F. Adams PetscErrorCode  PCCreateGAMG_AGG(PC pc)
1179c8b0795cSMark F. Adams {
1180c8b0795cSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
1181c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1182c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg;
11832e68589bSMark F. Adams 
1184c8b0795cSMark F. Adams   PetscFunctionBegin;
1185c8b0795cSMark F. Adams   /* create sub context for SA */
11869566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(pc,&pc_gamg_agg));
1187c8b0795cSMark F. Adams   pc_gamg->subctx = pc_gamg_agg;
1188c8b0795cSMark F. Adams 
11891ab5ffc9SJed Brown   pc_gamg->ops->setfromoptions = PCSetFromOptions_GAMG_AGG;
11909b8ffb57SJed Brown   pc_gamg->ops->destroy        = PCDestroy_GAMG_AGG;
1191c8b0795cSMark F. Adams   /* reset does not do anything; setup not virtual */
1192c8b0795cSMark F. Adams 
1193c8b0795cSMark F. Adams   /* set internal function pointers */
1194fd1112cbSBarry Smith   pc_gamg->ops->graph             = PCGAMGGraph_AGG;
11951ab5ffc9SJed Brown   pc_gamg->ops->coarsen           = PCGAMGCoarsen_AGG;
11961ab5ffc9SJed Brown   pc_gamg->ops->prolongator       = PCGAMGProlongator_AGG;
1197fd1112cbSBarry Smith   pc_gamg->ops->optprolongator    = PCGAMGOptProlongator_AGG;
11981ab5ffc9SJed Brown   pc_gamg->ops->createdefaultdata = PCSetData_AGG;
11995adeb434SBarry Smith   pc_gamg->ops->view              = PCView_GAMG_AGG;
1200c8b0795cSMark F. Adams 
1201cab9ed1eSBarry Smith   pc_gamg_agg->square_graph = 1;
1202cab9ed1eSBarry Smith   pc_gamg_agg->sym_graph    = PETSC_FALSE;
1203cab9ed1eSBarry Smith   pc_gamg_agg->nsmooths     = 1;
1204cab9ed1eSBarry Smith 
12059566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNSmooths_C",PCGAMGSetNSmooths_AGG));
12069566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetSymGraph_C",PCGAMGSetSymGraph_AGG));
12079566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetSquareGraph_C",PCGAMGSetSquareGraph_AGG));
12089566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_AGG));
12092e68589bSMark F. Adams   PetscFunctionReturn(0);
12102e68589bSMark F. Adams }
1211