xref: /petsc/src/ksp/pc/impls/gamg/agg.c (revision 074aec55f51aed287b387eef488abfd3536af0fd)
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 
292e68589bSMark F. Adams .seealso: ()
302e68589bSMark F. Adams @*/
312e68589bSMark F. Adams PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n)
322e68589bSMark F. Adams {
332e68589bSMark F. Adams   PetscErrorCode ierr;
342e68589bSMark F. Adams 
352e68589bSMark F. Adams   PetscFunctionBegin;
362e68589bSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
37d5d25401SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
382e68589bSMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetNSmooths_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
392e68589bSMark F. Adams   PetscFunctionReturn(0);
402e68589bSMark F. Adams }
412e68589bSMark F. Adams 
42e0877f53SBarry Smith static PetscErrorCode PCGAMGSetNSmooths_AGG(PC pc, PetscInt n)
432e68589bSMark F. Adams {
442e68589bSMark F. Adams   PC_MG       *mg          = (PC_MG*)pc->data;
452e68589bSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG*)mg->innerctx;
46c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
472e68589bSMark F. Adams 
482e68589bSMark F. Adams   PetscFunctionBegin;
49c8b0795cSMark F. Adams   pc_gamg_agg->nsmooths = n;
50c8b0795cSMark F. Adams   PetscFunctionReturn(0);
51c8b0795cSMark F. Adams }
52c8b0795cSMark F. Adams 
53c8b0795cSMark F. Adams /*@
54cab9ed1eSBarry Smith    PCGAMGSetSymGraph - Symmetrize the graph before computing the aggregation. Some algorithms require the graph be symmetric
55c8b0795cSMark F. Adams 
56d5d25401SBarry Smith    Logically Collective on PC
57c8b0795cSMark F. Adams 
58c8b0795cSMark F. Adams    Input Parameters:
59cab9ed1eSBarry Smith +  pc - the preconditioner context
60a2b725a8SWilliam Gropp -  n - PETSC_TRUE or PETSC_FALSE
61c8b0795cSMark F. Adams 
62c8b0795cSMark F. Adams    Options Database Key:
63cab9ed1eSBarry Smith .  -pc_gamg_sym_graph <true,default=false> - symmetrize the graph before computing the aggregation
64c8b0795cSMark F. Adams 
65c8b0795cSMark F. Adams    Level: intermediate
66c8b0795cSMark F. Adams 
67cab9ed1eSBarry Smith .seealso: PCGAMGSetSquareGraph()
68c8b0795cSMark F. Adams @*/
69c8b0795cSMark F. Adams PetscErrorCode PCGAMGSetSymGraph(PC pc, PetscBool n)
70c8b0795cSMark F. Adams {
71c8b0795cSMark F. Adams   PetscErrorCode ierr;
72c8b0795cSMark F. Adams 
73c8b0795cSMark F. Adams   PetscFunctionBegin;
74c8b0795cSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
75d5d25401SBarry Smith   PetscValidLogicalCollectiveBool(pc,n,2);
76c8b0795cSMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetSymGraph_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
77c8b0795cSMark F. Adams   PetscFunctionReturn(0);
78c8b0795cSMark F. Adams }
79c8b0795cSMark F. Adams 
80e0877f53SBarry Smith static PetscErrorCode PCGAMGSetSymGraph_AGG(PC pc, PetscBool n)
81c8b0795cSMark F. Adams {
82c8b0795cSMark F. Adams   PC_MG       *mg          = (PC_MG*)pc->data;
83c8b0795cSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG*)mg->innerctx;
84c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
85c8b0795cSMark F. Adams 
86c8b0795cSMark F. Adams   PetscFunctionBegin;
87c8b0795cSMark F. Adams   pc_gamg_agg->sym_graph = n;
882e68589bSMark F. Adams   PetscFunctionReturn(0);
892e68589bSMark F. Adams }
902e68589bSMark F. Adams 
91ef4ad70eSMark F. Adams /*@
92cab9ed1eSBarry Smith    PCGAMGSetSquareGraph -  Square the graph, ie. compute A'*A before aggregating it
93ef4ad70eSMark F. Adams 
94d5d25401SBarry Smith    Logically Collective on PC
95ef4ad70eSMark F. Adams 
96ef4ad70eSMark F. Adams    Input Parameters:
97cab9ed1eSBarry Smith +  pc - the preconditioner context
98d5d25401SBarry Smith -  n - 0, 1 or more
99ef4ad70eSMark F. Adams 
100ef4ad70eSMark F. Adams    Options Database Key:
101cab9ed1eSBarry Smith .  -pc_gamg_square_graph <n,default = 1> - number of levels to square the graph on before aggregating it
102ef4ad70eSMark F. Adams 
103af3c827dSMark Adams    Notes:
104af3c827dSMark 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.
105af3c827dSMark Adams 
106ef4ad70eSMark F. Adams    Level: intermediate
107ef4ad70eSMark F. Adams 
108af3c827dSMark Adams .seealso: PCGAMGSetSymGraph(), PCGAMGSetThreshold()
109ef4ad70eSMark F. Adams @*/
1109ab59c8bSMark Adams PetscErrorCode PCGAMGSetSquareGraph(PC pc, PetscInt n)
111ef4ad70eSMark F. Adams {
112ef4ad70eSMark F. Adams   PetscErrorCode ierr;
113ef4ad70eSMark F. Adams 
114ef4ad70eSMark F. Adams   PetscFunctionBegin;
115ef4ad70eSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
116d5d25401SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
1179ab59c8bSMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetSquareGraph_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
118ef4ad70eSMark F. Adams   PetscFunctionReturn(0);
119ef4ad70eSMark F. Adams }
120ef4ad70eSMark F. Adams 
121e0877f53SBarry Smith static PetscErrorCode PCGAMGSetSquareGraph_AGG(PC pc, PetscInt n)
122ef4ad70eSMark F. Adams {
123ef4ad70eSMark F. Adams   PC_MG       *mg          = (PC_MG*)pc->data;
124ef4ad70eSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG*)mg->innerctx;
125ef4ad70eSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
126ef4ad70eSMark F. Adams 
127ef4ad70eSMark F. Adams   PetscFunctionBegin;
128ef4ad70eSMark F. Adams   pc_gamg_agg->square_graph = n;
129ef4ad70eSMark F. Adams   PetscFunctionReturn(0);
130ef4ad70eSMark F. Adams }
131ef4ad70eSMark F. Adams 
132e0877f53SBarry Smith static PetscErrorCode PCSetFromOptions_GAMG_AGG(PetscOptionItems *PetscOptionsObject,PC pc)
1332e68589bSMark F. Adams {
1342e68589bSMark F. Adams   PetscErrorCode ierr;
1352e68589bSMark F. Adams   PC_MG          *mg          = (PC_MG*)pc->data;
1362e68589bSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
137c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
1382e68589bSMark F. Adams 
1392e68589bSMark F. Adams   PetscFunctionBegin;
140e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"GAMG-AGG options");CHKERRQ(ierr);
1412e68589bSMark F. Adams   {
1428afaa268SBarry Smith     ierr = PetscOptionsInt("-pc_gamg_agg_nsmooths","smoothing steps for smoothed aggregation, usually 1","PCGAMGSetNSmooths",pc_gamg_agg->nsmooths,&pc_gamg_agg->nsmooths,NULL);CHKERRQ(ierr);
1438afaa268SBarry Smith     ierr = PetscOptionsBool("-pc_gamg_sym_graph","Set for asymmetric matrices","PCGAMGSetSymGraph",pc_gamg_agg->sym_graph,&pc_gamg_agg->sym_graph,NULL);CHKERRQ(ierr);
1449ab59c8bSMark Adams     ierr = 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);CHKERRQ(ierr);
1452e68589bSMark F. Adams   }
1462e68589bSMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
1472e68589bSMark F. Adams   PetscFunctionReturn(0);
1482e68589bSMark F. Adams }
1492e68589bSMark F. Adams 
1502e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
151e0877f53SBarry Smith static PetscErrorCode PCDestroy_GAMG_AGG(PC pc)
1522e68589bSMark F. Adams {
1532e68589bSMark F. Adams   PetscErrorCode ierr;
1542e68589bSMark F. Adams   PC_MG          *mg          = (PC_MG*)pc->data;
1552e68589bSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
1562e68589bSMark F. Adams 
1572e68589bSMark F. Adams   PetscFunctionBegin;
1589b8ffb57SJed Brown   ierr = PetscFree(pc_gamg->subctx);CHKERRQ(ierr);
1593ae0bb68SMark Adams   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",NULL);CHKERRQ(ierr);
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;
1792e68589bSMark F. Adams   PetscErrorCode ierr;
18069344418SMark F. Adams   PetscInt       arrsz,kk,ii,jj,nloc,ndatarows,ndf;
181a2f3521dSMark F. Adams   Mat            mat = pc->pmat;
1822e68589bSMark F. Adams 
1832e68589bSMark F. Adams   PetscFunctionBegin;
184a2f3521dSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
185a2f3521dSMark F. Adams   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
186302f38e8SMark F. Adams   nloc = a_nloc;
1872e68589bSMark F. Adams 
1882e68589bSMark F. Adams   /* SA: null space vectors */
18969344418SMark F. Adams   ierr = MatGetBlockSize(mat, &ndf);CHKERRQ(ierr); /* this does not work for Stokes */
19069344418SMark F. Adams   if (coords && ndf==1) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */
191a2f3521dSMark F. Adams   else if (coords) {
1922c71b3e2SJacob Faibussowitsch     PetscCheckFalse(ndm > ndf,PETSC_COMM_SELF,PETSC_ERR_PLIB,"degrees of motion %D > block size %D",ndm,ndf);
19369344418SMark F. Adams     pc_gamg->data_cell_cols = (ndm==2 ? 3 : 6); /* displacement elasticity */
19469344418SMark F. Adams     if (ndm != ndf) {
1952c71b3e2SJacob Faibussowitsch       PetscCheckFalse(pc_gamg->data_cell_cols != ndf,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Don't know how to create null space for ndm=%D, ndf=%D.  Use MatSetNearNullSpace().",ndm,ndf);
196a2f3521dSMark F. Adams     }
197806fa848SBarry Smith   } else pc_gamg->data_cell_cols = ndf; /* no data, force SA with constant null space vectors */
19871959b99SBarry Smith   pc_gamg->data_cell_rows = ndatarows = ndf;
1992c71b3e2SJacob Faibussowitsch   PetscCheckFalse(pc_gamg->data_cell_cols <= 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"pc_gamg->data_cell_cols %D <= 0",pc_gamg->data_cell_cols);
200c8b0795cSMark F. Adams   arrsz = nloc*pc_gamg->data_cell_rows*pc_gamg->data_cell_cols;
2012e68589bSMark F. Adams 
2027f66b68fSMark Adams   if (!pc_gamg->data || (pc_gamg->data_sz != arrsz)) {
2032e68589bSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
204854ce69bSBarry Smith     ierr = PetscMalloc1(arrsz+1, &pc_gamg->data);CHKERRQ(ierr);
2052e68589bSMark F. Adams   }
2062e68589bSMark F. Adams   /* copy data in - column oriented */
2072e68589bSMark F. Adams   for (kk=0; kk<nloc; kk++) {
20869344418SMark F. Adams     const PetscInt M     = nloc*pc_gamg->data_cell_rows; /* stride into data */
20969344418SMark F. Adams     PetscReal      *data = &pc_gamg->data[kk*ndatarows]; /* start of cell */
210c8b0795cSMark F. Adams     if (pc_gamg->data_cell_cols==1) *data = 1.0;
2112e68589bSMark F. Adams     else {
21269344418SMark F. Adams       /* translational modes */
2132fa5cd67SKarl Rupp       for (ii=0;ii<ndatarows;ii++) {
2142fa5cd67SKarl Rupp         for (jj=0;jj<ndatarows;jj++) {
21569344418SMark F. Adams           if (ii==jj)data[ii*M + jj] = 1.0;
2162e68589bSMark F. Adams           else data[ii*M + jj] = 0.0;
2172fa5cd67SKarl Rupp         }
2182fa5cd67SKarl Rupp       }
21969344418SMark F. Adams 
22069344418SMark F. Adams       /* rotational modes */
2212e68589bSMark F. Adams       if (coords) {
22269344418SMark F. Adams         if (ndm == 2) {
2232e68589bSMark F. Adams           data   += 2*M;
2242e68589bSMark F. Adams           data[0] = -coords[2*kk+1];
2252e68589bSMark F. Adams           data[1] =  coords[2*kk];
226806fa848SBarry Smith         } else {
2272e68589bSMark F. Adams           data   += 3*M;
2282e68589bSMark F. Adams           data[0] = 0.0;             data[M+0] =  coords[3*kk+2]; data[2*M+0] = -coords[3*kk+1];
2292e68589bSMark F. Adams           data[1] = -coords[3*kk+2]; data[M+1] = 0.0;             data[2*M+1] =  coords[3*kk];
2302e68589bSMark F. Adams           data[2] =  coords[3*kk+1]; data[M+2] = -coords[3*kk];   data[2*M+2] = 0.0;
2312e68589bSMark F. Adams         }
2322e68589bSMark F. Adams       }
2332e68589bSMark F. Adams     }
2342e68589bSMark F. Adams   }
2352e68589bSMark F. Adams   pc_gamg->data_sz = arrsz;
2362e68589bSMark F. Adams   PetscFunctionReturn(0);
2372e68589bSMark F. Adams }
2382e68589bSMark F. Adams 
239b43b03e9SMark F. Adams typedef PetscInt NState;
240b43b03e9SMark F. Adams static const NState NOT_DONE=-2;
241b43b03e9SMark F. Adams static const NState DELETED =-1;
242b43b03e9SMark F. Adams static const NState REMOVED =-3;
243b43b03e9SMark F. Adams #define IS_SELECTED(s) (s!=DELETED && s!=NOT_DONE && s!=REMOVED)
244b43b03e9SMark F. Adams 
245c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
246c8b0795cSMark F. Adams /*
247b43b03e9SMark F. Adams    smoothAggs - greedy grab of with G1 (unsquared graph) -- AIJ specific
248b43b03e9SMark F. Adams      - AGG-MG specific: clears singletons out of 'selected_2'
249c8b0795cSMark F. Adams 
250c8b0795cSMark F. Adams    Input Parameter:
251fdb86571SRichard Tran Mills    . Gmat_2 - global matrix of graph (data not defined)   base (squared) graph
2522adfe9d3SBarry Smith    . Gmat_1 - base graph to grab with                 base graph
253c8b0795cSMark F. Adams    Input/Output Parameter:
2540cbbd2e1SMark F. Adams    . aggs_2 - linked list of aggs with gids)
255c8b0795cSMark F. Adams */
256e3df7ea3SBarry Smith static PetscErrorCode smoothAggs(PC pc,Mat Gmat_2, Mat Gmat_1,PetscCoarsenData *aggs_2)
257c8b0795cSMark F. Adams {
258c8b0795cSMark F. Adams   PetscErrorCode ierr;
259c8b0795cSMark F. Adams   PetscBool      isMPI;
2600a545947SLisandro Dalcin   Mat_SeqAIJ     *matA_1, *matB_1=NULL;
2613b4367a7SBarry Smith   MPI_Comm       comm;
2620cbbd2e1SMark F. Adams   PetscInt       lid,*ii,*idx,ix,Iend,my0,kk,n,j;
2630a545947SLisandro Dalcin   Mat_MPIAIJ     *mpimat_2 = NULL, *mpimat_1=NULL;
264c8b0795cSMark F. Adams   const PetscInt nloc      = Gmat_2->rmap->n;
2650cbbd2e1SMark F. Adams   PetscScalar    *cpcol_1_state,*cpcol_2_state,*cpcol_2_par_orig,*lid_parent_gid;
2660cbbd2e1SMark F. Adams   PetscInt       *lid_cprowID_1;
267c8b0795cSMark F. Adams   NState         *lid_state;
2680cbbd2e1SMark F. Adams   Vec            ghost_par_orig2;
269c8b0795cSMark F. Adams 
270c8b0795cSMark F. Adams   PetscFunctionBegin;
2713b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Gmat_2,&comm);CHKERRQ(ierr);
272c8b0795cSMark F. Adams   ierr = MatGetOwnershipRange(Gmat_1,&my0,&Iend);CHKERRQ(ierr);
273c8b0795cSMark F. Adams 
274c8b0795cSMark F. Adams   /* get submatrices */
275b92f168fSBarry Smith   ierr = PetscStrbeginswith(((PetscObject)Gmat_1)->type_name,MATMPIAIJ,&isMPI);CHKERRQ(ierr);
276c8b0795cSMark F. Adams   if (isMPI) {
277c8b0795cSMark F. Adams     /* grab matrix objects */
278c8b0795cSMark F. Adams     mpimat_2 = (Mat_MPIAIJ*)Gmat_2->data;
279c8b0795cSMark F. Adams     mpimat_1 = (Mat_MPIAIJ*)Gmat_1->data;
280c8b0795cSMark F. Adams     matA_1   = (Mat_SeqAIJ*)mpimat_1->A->data;
281c8b0795cSMark F. Adams     matB_1   = (Mat_SeqAIJ*)mpimat_1->B->data;
282c8b0795cSMark F. Adams 
283c8b0795cSMark F. Adams     /* force compressed row storage for B matrix in AuxMat */
28411e456e1SBarry Smith     ierr = MatCheckCompressedRow(mpimat_1->B,matB_1->nonzerorowcnt,&matB_1->compressedrow,matB_1->i,Gmat_1->rmap->n,-1.0);CHKERRQ(ierr);
285c8b0795cSMark F. Adams 
286785e854fSJed Brown     ierr = PetscMalloc1(nloc, &lid_cprowID_1);CHKERRQ(ierr);
2870cbbd2e1SMark F. Adams     for (lid = 0; lid < nloc; lid++) lid_cprowID_1[lid] = -1;
288c8b0795cSMark F. Adams     for (ix=0; ix<matB_1->compressedrow.nrows; ix++) {
289c8b0795cSMark F. Adams       PetscInt lid = matB_1->compressedrow.rindex[ix];
290c8b0795cSMark F. Adams       lid_cprowID_1[lid] = ix;
291c8b0795cSMark F. Adams     }
292806fa848SBarry Smith   } else {
29315687449SMark Adams     PetscBool isAIJ;
294b92f168fSBarry Smith     ierr = PetscStrbeginswith(((PetscObject)Gmat_1)->type_name,MATSEQAIJ,&isAIJ);CHKERRQ(ierr);
2952c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!isAIJ,PETSC_COMM_SELF,PETSC_ERR_USER,"Require AIJ matrix.");
296c8b0795cSMark F. Adams     matA_1        = (Mat_SeqAIJ*)Gmat_1->data;
2970298fd71SBarry Smith     lid_cprowID_1 = NULL;
298c8b0795cSMark F. Adams   }
29978a438d6SMark Adams   if (nloc>0) {
3002c71b3e2SJacob Faibussowitsch     PetscCheckFalse(matB_1 && !matB_1->compressedrow.use,PETSC_COMM_SELF,PETSC_ERR_PLIB,"matB_1 && !matB_1->compressedrow.use: PETSc bug???");
30178a438d6SMark Adams   }
302c8b0795cSMark F. Adams   /* get state of locals and selected gid for deleted */
303e632b94dSBarry Smith   ierr = PetscMalloc2(nloc, &lid_state,nloc, &lid_parent_gid);CHKERRQ(ierr);
304c8b0795cSMark F. Adams   for (lid = 0; lid < nloc; lid++) {
3050cbbd2e1SMark F. Adams     lid_parent_gid[lid] = -1.0;
306c8b0795cSMark F. Adams     lid_state[lid]      = DELETED;
307c8b0795cSMark F. Adams   }
3080cbbd2e1SMark F. Adams 
3090cbbd2e1SMark F. Adams   /* set lid_state */
3100cbbd2e1SMark F. Adams   for (lid = 0; lid < nloc; lid++) {
311539c167fSBarry Smith     PetscCDIntNd *pos;
312e78576d6SMark F. Adams     ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr);
313e78576d6SMark F. Adams     if (pos) {
314e78576d6SMark F. Adams       PetscInt gid1;
3152fa5cd67SKarl Rupp 
316484f0a72SBarry Smith       ierr = PetscCDIntNdGetID(pos, &gid1);CHKERRQ(ierr);
3172c71b3e2SJacob Faibussowitsch       PetscCheckFalse(gid1 != lid+my0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"gid1 %D != lid %D + my0 %D",gid1,lid,my0);
3180cbbd2e1SMark F. Adams       lid_state[lid] = gid1;
319b43b03e9SMark F. Adams     }
320b43b03e9SMark F. Adams   }
3210cbbd2e1SMark F. Adams 
3220cbbd2e1SMark F. Adams   /* map local to selected local, DELETED means a ghost owns it */
323c8b0795cSMark F. Adams   for (lid=kk=0; lid<nloc; lid++) {
324c8b0795cSMark F. Adams     NState state = lid_state[lid];
325c8b0795cSMark F. Adams     if (IS_SELECTED(state)) {
326539c167fSBarry Smith       PetscCDIntNd *pos;
327e78576d6SMark F. Adams       ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr);
328e78576d6SMark F. Adams       while (pos) {
329e78576d6SMark F. Adams         PetscInt gid1;
330484f0a72SBarry Smith         ierr = PetscCDIntNdGetID(pos, &gid1);CHKERRQ(ierr);
331e78576d6SMark F. Adams         ierr = PetscCDGetNextPos(aggs_2,lid,&pos);CHKERRQ(ierr);
3322fa5cd67SKarl Rupp         if (gid1 >= my0 && gid1 < Iend) lid_parent_gid[gid1-my0] = (PetscScalar)(lid + my0);
333c8b0795cSMark F. Adams       }
3340cbbd2e1SMark F. Adams     }
3350cbbd2e1SMark F. Adams   }
3360cbbd2e1SMark F. Adams   /* get 'cpcol_1/2_state' & cpcol_2_par_orig - uses mpimat_1/2->lvec for temp space */
337c8b0795cSMark F. Adams   if (isMPI) {
338c8b0795cSMark F. Adams     Vec tempVec;
339c8b0795cSMark F. Adams     /* get 'cpcol_1_state' */
3400a545947SLisandro Dalcin     ierr = MatCreateVecs(Gmat_1, &tempVec, NULL);CHKERRQ(ierr);
341c8b0795cSMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
342c8b0795cSMark F. Adams       PetscScalar v = (PetscScalar)lid_state[kk];
343c8b0795cSMark F. Adams       ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr);
344c8b0795cSMark F. Adams     }
345c8b0795cSMark F. Adams     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
346c8b0795cSMark F. Adams     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
347806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
348806fa848SBarry Smith     ierr = VecScatterEnd(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
349c8b0795cSMark F. Adams     ierr = VecGetArray(mpimat_1->lvec, &cpcol_1_state);CHKERRQ(ierr);
350c8b0795cSMark F. Adams     /* get 'cpcol_2_state' */
351806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
352806fa848SBarry Smith     ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
353c8b0795cSMark F. Adams     ierr = VecGetArray(mpimat_2->lvec, &cpcol_2_state);CHKERRQ(ierr);
3540cbbd2e1SMark F. Adams     /* get 'cpcol_2_par_orig' */
3550cbbd2e1SMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
3560cbbd2e1SMark F. Adams       PetscScalar v = (PetscScalar)lid_parent_gid[kk];
3570cbbd2e1SMark F. Adams       ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr);
3580cbbd2e1SMark F. Adams     }
3590cbbd2e1SMark F. Adams     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
3600cbbd2e1SMark F. Adams     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
3610cbbd2e1SMark F. Adams     ierr = VecDuplicate(mpimat_2->lvec, &ghost_par_orig2);CHKERRQ(ierr);
362806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
363806fa848SBarry Smith     ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
3640cbbd2e1SMark F. Adams     ierr = VecGetArray(ghost_par_orig2, &cpcol_2_par_orig);CHKERRQ(ierr);
3650cbbd2e1SMark F. Adams 
366c8b0795cSMark F. Adams     ierr = VecDestroy(&tempVec);CHKERRQ(ierr);
367c8b0795cSMark F. Adams   } /* ismpi */
368c8b0795cSMark F. Adams   for (lid=0; lid<nloc; lid++) {
369c8b0795cSMark F. Adams     NState state = lid_state[lid];
3700cbbd2e1SMark F. Adams     if (IS_SELECTED(state)) {
3710cbbd2e1SMark F. Adams       /* steal locals */
372c8b0795cSMark F. Adams       ii  = matA_1->i; n = ii[lid+1] - ii[lid];
373c8b0795cSMark F. Adams       idx = matA_1->j + ii[lid];
374c8b0795cSMark F. Adams       for (j=0; j<n; j++) {
3750cbbd2e1SMark F. Adams         PetscInt lidj   = idx[j], sgid;
376c8b0795cSMark F. Adams         NState   statej = lid_state[lidj];
3770cbbd2e1SMark F. Adams         if (statej==DELETED && (sgid=(PetscInt)PetscRealPart(lid_parent_gid[lidj])) != lid+my0) { /* steal local */
3780cbbd2e1SMark F. Adams           lid_parent_gid[lidj] = (PetscScalar)(lid+my0); /* send this if sgid is not local */
3790cbbd2e1SMark F. Adams           if (sgid >= my0 && sgid < Iend) {       /* I'm stealing this local from a local sgid */
3800cbbd2e1SMark F. Adams             PetscInt     hav=0,slid=sgid-my0,gidj=lidj+my0;
381539c167fSBarry Smith             PetscCDIntNd *pos,*last=NULL;
382c8b0795cSMark F. Adams             /* looking for local from local so id_llist_2 works */
383e78576d6SMark F. Adams             ierr = PetscCDGetHeadPos(aggs_2,slid,&pos);CHKERRQ(ierr);
384e78576d6SMark F. Adams             while (pos) {
385e78576d6SMark F. Adams               PetscInt gid;
386484f0a72SBarry Smith               ierr = PetscCDIntNdGetID(pos, &gid);CHKERRQ(ierr);
3870cbbd2e1SMark F. Adams               if (gid == gidj) {
3882c71b3e2SJacob Faibussowitsch                 PetscCheckFalse(!last,PETSC_COMM_SELF,PETSC_ERR_PLIB,"last cannot be null");
38941b27cdeSMark F. Adams                 ierr = PetscCDRemoveNextNode(aggs_2, slid, last);CHKERRQ(ierr);
39041b27cdeSMark F. Adams                 ierr = PetscCDAppendNode(aggs_2, lid, pos);CHKERRQ(ierr);
3910cbbd2e1SMark F. Adams                 hav  = 1;
392c8b0795cSMark F. Adams                 break;
393806fa848SBarry Smith               } else last = pos;
394e78576d6SMark F. Adams               ierr = PetscCDGetNextPos(aggs_2,slid,&pos);CHKERRQ(ierr);
395c8b0795cSMark F. Adams             }
396c8b0795cSMark F. Adams             if (hav != 1) {
3972c71b3e2SJacob Faibussowitsch               PetscCheckFalse(!hav,PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix");
39898921bdaSJacob Faibussowitsch               SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %D times???",hav);
399c8b0795cSMark F. Adams             }
400806fa848SBarry Smith           } else {            /* I'm stealing this local, owned by a ghost */
4012c71b3e2SJacob Faibussowitsch             PetscCheckFalse(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 : "");
40241b27cdeSMark F. Adams             ierr = PetscCDAppendID(aggs_2, lid, lidj+my0);CHKERRQ(ierr);
403c8b0795cSMark F. Adams           }
404c8b0795cSMark F. Adams         }
4050cbbd2e1SMark F. Adams       } /* local neighbors */
406806fa848SBarry Smith     } else if (state == DELETED && lid_cprowID_1) {
4070cbbd2e1SMark F. Adams       PetscInt sgidold = (PetscInt)PetscRealPart(lid_parent_gid[lid]);
408c8b0795cSMark F. Adams       /* see if I have a selected ghost neighbor that will steal me */
409c8b0795cSMark F. Adams       if ((ix=lid_cprowID_1[lid]) != -1) {
410c8b0795cSMark F. Adams         ii  = matB_1->compressedrow.i; n = ii[ix+1] - ii[ix];
411c8b0795cSMark F. Adams         idx = matB_1->j + ii[ix];
412c8b0795cSMark F. Adams         for (j=0; j<n; j++) {
413c8b0795cSMark F. Adams           PetscInt cpid   = idx[j];
414c8b0795cSMark F. Adams           NState   statej = (NState)PetscRealPart(cpcol_1_state[cpid]);
415c8b0795cSMark F. Adams           if (IS_SELECTED(statej) && sgidold != (PetscInt)statej) { /* ghost will steal this, remove from my list */
4160cbbd2e1SMark F. Adams             lid_parent_gid[lid] = (PetscScalar)statej; /* send who selected */
4170cbbd2e1SMark F. Adams             if (sgidold>=my0 && sgidold<Iend) { /* this was mine */
4180cbbd2e1SMark F. Adams               PetscInt     hav=0,oldslidj=sgidold-my0;
419539c167fSBarry Smith               PetscCDIntNd *pos,*last=NULL;
4200cbbd2e1SMark F. Adams               /* remove from 'oldslidj' list */
421e78576d6SMark F. Adams               ierr = PetscCDGetHeadPos(aggs_2,oldslidj,&pos);CHKERRQ(ierr);
422e78576d6SMark F. Adams               while (pos) {
423e78576d6SMark F. Adams                 PetscInt gid;
424484f0a72SBarry Smith                 ierr = PetscCDIntNdGetID(pos, &gid);CHKERRQ(ierr);
4250cbbd2e1SMark F. Adams                 if (lid+my0 == gid) {
4260cbbd2e1SMark F. Adams                   /* id_llist_2[lastid] = id_llist_2[flid];   /\* remove lid from oldslidj list *\/ */
4272c71b3e2SJacob Faibussowitsch                   PetscCheckFalse(!last,PETSC_COMM_SELF,PETSC_ERR_PLIB,"last cannot be null");
42841b27cdeSMark F. Adams                   ierr = PetscCDRemoveNextNode(aggs_2, oldslidj, last);CHKERRQ(ierr);
4290cbbd2e1SMark F. Adams                   /* ghost (PetscScalar)statej will add this later */
4300cbbd2e1SMark F. Adams                   hav = 1;
431c8b0795cSMark F. Adams                   break;
432806fa848SBarry Smith                 } else last = pos;
433e78576d6SMark F. Adams                 ierr = PetscCDGetNextPos(aggs_2,oldslidj,&pos);CHKERRQ(ierr);
434c8b0795cSMark F. Adams               }
435c8b0795cSMark F. Adams               if (hav != 1) {
4362c71b3e2SJacob Faibussowitsch                 PetscCheckFalse(!hav,PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix");
43798921bdaSJacob Faibussowitsch                 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %D times???",hav);
438c8b0795cSMark F. Adams               }
439806fa848SBarry Smith             } else {
440d5d25401SBarry Smith               /* TODO: ghosts remove this later */
4410cbbd2e1SMark F. Adams             }
442c8b0795cSMark F. Adams           }
443c8b0795cSMark F. Adams         }
444c8b0795cSMark F. Adams       }
445c8b0795cSMark F. Adams     } /* selected/deleted */
446c8b0795cSMark F. Adams   } /* node loop */
447c8b0795cSMark F. Adams 
448c8b0795cSMark F. Adams   if (isMPI) {
4490cbbd2e1SMark F. Adams     PetscScalar     *cpcol_2_parent,*cpcol_2_gid;
4500cbbd2e1SMark F. Adams     Vec             tempVec,ghostgids2,ghostparents2;
4510cbbd2e1SMark F. Adams     PetscInt        cpid,nghost_2;
4521943db53SBarry Smith     PCGAMGHashTable gid_cpid;
453c8b0795cSMark F. Adams 
4540cbbd2e1SMark F. Adams     ierr = VecGetSize(mpimat_2->lvec, &nghost_2);CHKERRQ(ierr);
4550a545947SLisandro Dalcin     ierr = MatCreateVecs(Gmat_2, &tempVec, NULL);CHKERRQ(ierr);
4560cbbd2e1SMark F. Adams 
4570cbbd2e1SMark F. Adams     /* get 'cpcol_2_parent' */
458c8b0795cSMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
4590cbbd2e1SMark F. Adams       ierr = VecSetValues(tempVec, 1, &j, &lid_parent_gid[kk], INSERT_VALUES);CHKERRQ(ierr);
460c8b0795cSMark F. Adams     }
461c8b0795cSMark F. Adams     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
462c8b0795cSMark F. Adams     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
4630cbbd2e1SMark F. Adams     ierr = VecDuplicate(mpimat_2->lvec, &ghostparents2);CHKERRQ(ierr);
464806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
465806fa848SBarry Smith     ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4660cbbd2e1SMark F. Adams     ierr = VecGetArray(ghostparents2, &cpcol_2_parent);CHKERRQ(ierr);
4670cbbd2e1SMark F. Adams 
4680cbbd2e1SMark F. Adams     /* get 'cpcol_2_gid' */
4690cbbd2e1SMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
4700cbbd2e1SMark F. Adams       PetscScalar v = (PetscScalar)j;
4710cbbd2e1SMark F. Adams       ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr);
4720cbbd2e1SMark F. Adams     }
4730cbbd2e1SMark F. Adams     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
4740cbbd2e1SMark F. Adams     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
4750cbbd2e1SMark F. Adams     ierr = VecDuplicate(mpimat_2->lvec, &ghostgids2);CHKERRQ(ierr);
476806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
477806fa848SBarry Smith     ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4780cbbd2e1SMark F. Adams     ierr = VecGetArray(ghostgids2, &cpcol_2_gid);CHKERRQ(ierr);
479c8b0795cSMark F. Adams     ierr = VecDestroy(&tempVec);CHKERRQ(ierr);
480c8b0795cSMark F. Adams 
4810cbbd2e1SMark F. Adams     /* look for deleted ghosts and add to table */
4821943db53SBarry Smith     ierr = PCGAMGHashTableCreate(2*nghost_2+1, &gid_cpid);CHKERRQ(ierr);
4830cbbd2e1SMark F. Adams     for (cpid = 0; cpid < nghost_2; cpid++) {
4840cbbd2e1SMark F. Adams       NState state = (NState)PetscRealPart(cpcol_2_state[cpid]);
4850cbbd2e1SMark F. Adams       if (state==DELETED) {
4860cbbd2e1SMark F. Adams         PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]);
4870cbbd2e1SMark F. Adams         PetscInt sgid_old = (PetscInt)PetscRealPart(cpcol_2_par_orig[cpid]);
4880cbbd2e1SMark F. Adams         if (sgid_old == -1 && sgid_new != -1) {
4890cbbd2e1SMark F. Adams           PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]);
4901943db53SBarry Smith           ierr = PCGAMGHashTableAdd(&gid_cpid, gid, cpid);CHKERRQ(ierr);
4910cbbd2e1SMark F. Adams         }
4920cbbd2e1SMark F. Adams       }
4930cbbd2e1SMark F. Adams     }
494c8b0795cSMark F. Adams 
4950cbbd2e1SMark F. Adams     /* look for deleted ghosts and see if they moved - remove it */
496c8b0795cSMark F. Adams     for (lid=0; lid<nloc; lid++) {
497c8b0795cSMark F. Adams       NState state = lid_state[lid];
498c8b0795cSMark F. Adams       if (IS_SELECTED(state)) {
499539c167fSBarry Smith         PetscCDIntNd *pos,*last=NULL;
500c8b0795cSMark F. Adams         /* look for deleted ghosts and see if they moved */
501e78576d6SMark F. Adams         ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr);
502e78576d6SMark F. Adams         while (pos) {
503e78576d6SMark F. Adams           PetscInt gid;
504484f0a72SBarry Smith           ierr = PetscCDIntNdGetID(pos, &gid);CHKERRQ(ierr);
505e78576d6SMark F. Adams 
5060cbbd2e1SMark F. Adams           if (gid < my0 || gid >= Iend) {
5071943db53SBarry Smith             ierr = PCGAMGHashTableFind(&gid_cpid, gid, &cpid);CHKERRQ(ierr);
5080cbbd2e1SMark F. Adams             if (cpid != -1) {
5090cbbd2e1SMark F. Adams               /* a moved ghost - */
5100cbbd2e1SMark F. Adams               /* id_llist_2[lastid] = id_llist_2[flid];    /\* remove 'flid' from list *\/ */
51141b27cdeSMark F. Adams               ierr = PetscCDRemoveNextNode(aggs_2, lid, last);CHKERRQ(ierr);
512806fa848SBarry Smith             } else last = pos;
513806fa848SBarry Smith           } else last = pos;
514e78576d6SMark F. Adams 
515e78576d6SMark F. Adams           ierr = PetscCDGetNextPos(aggs_2,lid,&pos);CHKERRQ(ierr);
516c8b0795cSMark F. Adams         } /* loop over list of deleted */
517c8b0795cSMark F. Adams       } /* selected */
518c8b0795cSMark F. Adams     }
5191943db53SBarry Smith     ierr = PCGAMGHashTableDestroy(&gid_cpid);CHKERRQ(ierr);
520c8b0795cSMark F. Adams 
5210cbbd2e1SMark F. Adams     /* look at ghosts, see if they changed - and it */
5220cbbd2e1SMark F. Adams     for (cpid = 0; cpid < nghost_2; cpid++) {
5230cbbd2e1SMark F. Adams       PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]);
5240cbbd2e1SMark F. Adams       if (sgid_new >= my0 && sgid_new < Iend) { /* this is mine */
5250cbbd2e1SMark F. Adams         PetscInt     gid     = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]);
5260cbbd2e1SMark F. Adams         PetscInt     slid_new=sgid_new-my0,hav=0;
527539c167fSBarry Smith         PetscCDIntNd *pos;
528539c167fSBarry Smith 
5290cbbd2e1SMark F. Adams         /* search for this gid to see if I have it */
530e78576d6SMark F. Adams         ierr = PetscCDGetHeadPos(aggs_2,slid_new,&pos);CHKERRQ(ierr);
531e78576d6SMark F. Adams         while (pos) {
532e78576d6SMark F. Adams           PetscInt gidj;
533484f0a72SBarry Smith           ierr = PetscCDIntNdGetID(pos, &gidj);CHKERRQ(ierr);
534e78576d6SMark F. Adams           ierr = PetscCDGetNextPos(aggs_2,slid_new,&pos);CHKERRQ(ierr);
535e78576d6SMark F. Adams 
5360cbbd2e1SMark F. Adams           if (gidj == gid) { hav = 1; break; }
537c8b0795cSMark F. Adams         }
538c8b0795cSMark F. Adams         if (hav != 1) {
539ffc955d6SMark F. Adams           /* insert 'flidj' into head of llist */
54041b27cdeSMark F. Adams           ierr = PetscCDAppendID(aggs_2, slid_new, gid);CHKERRQ(ierr);
541c8b0795cSMark F. Adams         }
542c8b0795cSMark F. Adams       }
543c8b0795cSMark F. Adams     }
544c8b0795cSMark F. Adams 
5450cbbd2e1SMark F. Adams     ierr = VecRestoreArray(mpimat_1->lvec, &cpcol_1_state);CHKERRQ(ierr);
5460cbbd2e1SMark F. Adams     ierr = VecRestoreArray(mpimat_2->lvec, &cpcol_2_state);CHKERRQ(ierr);
5470cbbd2e1SMark F. Adams     ierr = VecRestoreArray(ghostparents2, &cpcol_2_parent);CHKERRQ(ierr);
5480cbbd2e1SMark F. Adams     ierr = VecRestoreArray(ghostgids2, &cpcol_2_gid);CHKERRQ(ierr);
549c8b0795cSMark F. Adams     ierr = PetscFree(lid_cprowID_1);CHKERRQ(ierr);
5500cbbd2e1SMark F. Adams     ierr = VecDestroy(&ghostgids2);CHKERRQ(ierr);
5510cbbd2e1SMark F. Adams     ierr = VecDestroy(&ghostparents2);CHKERRQ(ierr);
5520cbbd2e1SMark F. Adams     ierr = VecDestroy(&ghost_par_orig2);CHKERRQ(ierr);
553c8b0795cSMark F. Adams   }
554c8b0795cSMark F. Adams 
555e632b94dSBarry Smith   ierr = PetscFree2(lid_state,lid_parent_gid);CHKERRQ(ierr);
556c8b0795cSMark F. Adams   PetscFunctionReturn(0);
557c8b0795cSMark F. Adams }
5582e68589bSMark F. Adams 
5592e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
5602e68589bSMark F. Adams /*
561a2f3521dSMark F. Adams    PCSetData_AGG - called if data is not set with PCSetCoordinates.
562a2f3521dSMark F. Adams       Looks in Mat for near null space.
563a2f3521dSMark F. Adams       Does not work for Stokes
5642e68589bSMark F. Adams 
5652e68589bSMark F. Adams   Input Parameter:
5662e68589bSMark F. Adams    . pc -
567a2f3521dSMark F. Adams    . a_A - matrix to get (near) null space out of.
5682e68589bSMark F. Adams */
569e0877f53SBarry Smith static PetscErrorCode PCSetData_AGG(PC pc, Mat a_A)
5702e68589bSMark F. Adams {
5712e68589bSMark F. Adams   PetscErrorCode ierr;
572b8cd405aSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
573b8cd405aSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
574b8cd405aSMark F. Adams   MatNullSpace   mnull;
57566f2ef4bSMark Adams 
576b817416eSBarry Smith   PetscFunctionBegin;
577b8cd405aSMark F. Adams   ierr = MatGetNearNullSpace(a_A, &mnull);CHKERRQ(ierr);
578b8cd405aSMark F. Adams   if (!mnull) {
57966f2ef4bSMark Adams     DM dm;
58066f2ef4bSMark Adams     ierr = PCGetDM(pc, &dm);CHKERRQ(ierr);
58166f2ef4bSMark Adams     if (!dm) {
58266f2ef4bSMark Adams       ierr = MatGetDM(a_A, &dm);CHKERRQ(ierr);
58366f2ef4bSMark Adams     }
58466f2ef4bSMark Adams     if (dm) {
58566f2ef4bSMark Adams       PetscObject deformation;
586b0d30dd6SMatthew G. Knepley       PetscInt    Nf;
587b0d30dd6SMatthew G. Knepley 
588b0d30dd6SMatthew G. Knepley       ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr);
589b0d30dd6SMatthew G. Knepley       if (Nf) {
59044a7f3ddSMatthew G. Knepley         ierr = DMGetField(dm, 0, NULL, &deformation);CHKERRQ(ierr);
59166f2ef4bSMark Adams         ierr = PetscObjectQuery((PetscObject)deformation,"nearnullspace",(PetscObject*)&mnull);CHKERRQ(ierr);
59266f2ef4bSMark Adams         if (!mnull) {
59366f2ef4bSMark Adams           ierr = PetscObjectQuery((PetscObject)deformation,"nullspace",(PetscObject*)&mnull);CHKERRQ(ierr);
59466f2ef4bSMark Adams         }
59566f2ef4bSMark Adams       }
59666f2ef4bSMark Adams     }
597b0d30dd6SMatthew G. Knepley   }
59866f2ef4bSMark Adams 
59966f2ef4bSMark Adams   if (!mnull) {
600a2f3521dSMark F. Adams     PetscInt bs,NN,MM;
6019d1b15c3SMark F. Adams     ierr = MatGetBlockSize(a_A, &bs);CHKERRQ(ierr);
60271959b99SBarry Smith     ierr = MatGetLocalSize(a_A, &MM, &NN);CHKERRQ(ierr);
6032c71b3e2SJacob Faibussowitsch     PetscCheckFalse(MM % bs,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MM %D must be divisible by bs %D",MM,bs);
6040298fd71SBarry Smith     ierr = PCSetCoordinates_AGG(pc, bs, MM/bs, NULL);CHKERRQ(ierr);
605806fa848SBarry Smith   } else {
606b8cd405aSMark F. Adams     PetscReal         *nullvec;
607b8cd405aSMark F. Adams     PetscBool         has_const;
608b8cd405aSMark F. Adams     PetscInt          i,j,mlocal,nvec,bs;
609d5d25401SBarry Smith     const Vec         *vecs;
610d5d25401SBarry Smith     const PetscScalar *v;
611b817416eSBarry Smith 
6120298fd71SBarry Smith     ierr = MatGetLocalSize(a_A,&mlocal,NULL);CHKERRQ(ierr);
613b8cd405aSMark F. Adams     ierr = MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs);CHKERRQ(ierr);
614a0dea87dSMark Adams     pc_gamg->data_sz = (nvec+!!has_const)*mlocal;
615785e854fSJed Brown     ierr = PetscMalloc1((nvec+!!has_const)*mlocal,&nullvec);CHKERRQ(ierr);
616b8cd405aSMark F. Adams     if (has_const) for (i=0; i<mlocal; i++) nullvec[i] = 1.0;
617b8cd405aSMark F. Adams     for (i=0; i<nvec; i++) {
618b8cd405aSMark F. Adams       ierr = VecGetArrayRead(vecs[i],&v);CHKERRQ(ierr);
619b8cd405aSMark F. Adams       for (j=0; j<mlocal; j++) nullvec[(i+!!has_const)*mlocal + j] = PetscRealPart(v[j]);
620b8cd405aSMark F. Adams       ierr = VecRestoreArrayRead(vecs[i],&v);CHKERRQ(ierr);
621b8cd405aSMark F. Adams     }
622b8cd405aSMark F. Adams     pc_gamg->data           = nullvec;
623b8cd405aSMark F. Adams     pc_gamg->data_cell_cols = (nvec+!!has_const);
62406e133e7SMark Adams     ierr = MatGetBlockSize(a_A, &bs);CHKERRQ(ierr);
625b8cd405aSMark F. Adams     pc_gamg->data_cell_rows = bs;
626b8cd405aSMark F. Adams   }
6272e68589bSMark F. Adams   PetscFunctionReturn(0);
6282e68589bSMark F. Adams }
6292e68589bSMark F. Adams 
6302e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
6312e68589bSMark F. Adams /*
6322e68589bSMark F. Adams  formProl0
6332e68589bSMark F. Adams 
6342e68589bSMark F. Adams    Input Parameter:
6352adfe9d3SBarry Smith    . agg_llists - list of arrays with aggregates -- list from selected vertices of aggregate unselected vertices
6362adfe9d3SBarry Smith    . bs - row block size
6370cbbd2e1SMark F. Adams    . nSAvec - column bs of new P
6380cbbd2e1SMark F. Adams    . my0crs - global index of start of locals
6392adfe9d3SBarry Smith    . data_stride - bs*(nloc nodes + ghost nodes) [data_stride][nSAvec]
6402e68589bSMark F. Adams    . data_in[data_stride*nSAvec] - local data on fine grid
6412e68589bSMark F. Adams    . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist'
6422e68589bSMark F. Adams   Output Parameter:
6432e68589bSMark F. Adams    . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data
6442e68589bSMark F. Adams    . a_Prol - prolongation operator
6452e68589bSMark F. Adams */
6462adfe9d3SBarry 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)
6472e68589bSMark F. Adams {
6482e68589bSMark F. Adams   PetscErrorCode  ierr;
649bd026e97SJed Brown   PetscInt        Istart,my0,Iend,nloc,clid,flid = 0,aggID,kk,jj,ii,mm,nSelected,minsz,nghosts,out_data_stride;
6503b4367a7SBarry Smith   MPI_Comm        comm;
6512e68589bSMark F. Adams   PetscReal       *out_data;
652539c167fSBarry Smith   PetscCDIntNd    *pos;
6531943db53SBarry Smith   PCGAMGHashTable fgid_flid;
6540cbbd2e1SMark F. Adams 
6552e68589bSMark F. Adams   PetscFunctionBegin;
6563b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)a_Prol,&comm);CHKERRQ(ierr);
6572e68589bSMark F. Adams   ierr = MatGetOwnershipRange(a_Prol, &Istart, &Iend);CHKERRQ(ierr);
65871959b99SBarry Smith   nloc = (Iend-Istart)/bs; my0 = Istart/bs;
6592c71b3e2SJacob Faibussowitsch   PetscCheckFalse((Iend-Istart) % bs,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Iend %D - Istart %D must be divisible by bs %D",Iend,Istart,bs);
6600cbbd2e1SMark F. Adams   Iend   /= bs;
6610cbbd2e1SMark F. Adams   nghosts = data_stride/bs - nloc;
6622e68589bSMark F. Adams 
6631943db53SBarry Smith   ierr = PCGAMGHashTableCreate(2*nghosts+1, &fgid_flid);CHKERRQ(ierr);
6640cbbd2e1SMark F. Adams   for (kk=0; kk<nghosts; kk++) {
6651943db53SBarry Smith     ierr = PCGAMGHashTableAdd(&fgid_flid, flid_fgid[nloc+kk], nloc+kk);CHKERRQ(ierr);
6662e68589bSMark F. Adams   }
6672e68589bSMark F. Adams 
6680cbbd2e1SMark F. Adams   /* count selected -- same as number of cols of P */
6690cbbd2e1SMark F. Adams   for (nSelected=mm=0; mm<nloc; mm++) {
670e78576d6SMark F. Adams     PetscBool ise;
671e78576d6SMark F. Adams     ierr = PetscCDEmptyAt(agg_llists, mm, &ise);CHKERRQ(ierr);
672e78576d6SMark F. Adams     if (!ise) nSelected++;
6730cbbd2e1SMark F. Adams   }
6740cbbd2e1SMark F. Adams   ierr = MatGetOwnershipRangeColumn(a_Prol, &ii, &jj);CHKERRQ(ierr);
6752c71b3e2SJacob Faibussowitsch   PetscCheckFalse((ii/nSAvec) != my0crs,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D /nSAvec %D  != my0crs %D",ii,nSAvec,my0crs);
6762c71b3e2SJacob Faibussowitsch   PetscCheckFalse(nSelected != (jj-ii)/nSAvec,PETSC_COMM_SELF,PETSC_ERR_PLIB,"nSelected %D != (jj %D - ii %D)/nSAvec %D",nSelected,jj,ii,nSAvec);
6770cbbd2e1SMark F. Adams 
6782e68589bSMark F. Adams   /* aloc space for coarse point data (output) */
6790cbbd2e1SMark F. Adams   out_data_stride = nSelected*nSAvec;
6802fa5cd67SKarl Rupp 
681785e854fSJed Brown   ierr = PetscMalloc1(out_data_stride*nSAvec, &out_data);CHKERRQ(ierr);
68233812677SSatish Balay   for (ii=0;ii<out_data_stride*nSAvec;ii++) out_data[ii]=PETSC_MAX_REAL;
6830cbbd2e1SMark F. Adams   *a_data_out = out_data; /* output - stride nSelected*nSAvec */
6842e68589bSMark F. Adams 
6852e68589bSMark F. Adams   /* find points and set prolongation */
686c8b0795cSMark F. Adams   minsz = 100;
6870cbbd2e1SMark F. Adams   for (mm = clid = 0; mm < nloc; mm++) {
688e78576d6SMark F. Adams     ierr = PetscCDSizeAt(agg_llists, mm, &jj);CHKERRQ(ierr);
689e78576d6SMark F. Adams     if (jj > 0) {
6900cbbd2e1SMark F. Adams       const PetscInt lid = mm, cgid = my0crs + clid;
6910cbbd2e1SMark F. Adams       PetscInt       cids[100]; /* max bs */
6920cbbd2e1SMark F. Adams       PetscBLASInt   asz  =jj,M=asz*bs,N=nSAvec,INFO;
6932e68589bSMark F. Adams       PetscBLASInt   Mdata=M+((N-M>0) ? N-M : 0),LDA=Mdata,LWORK=N*bs;
6942e68589bSMark F. Adams       PetscScalar    *qqc,*qqr,*TAU,*WORK;
6952e68589bSMark F. Adams       PetscInt       *fids;
69665d7b583SSatish Balay       PetscReal      *data;
697b817416eSBarry Smith 
6980cbbd2e1SMark F. Adams       /* count agg */
6990cbbd2e1SMark F. Adams       if (asz<minsz) minsz = asz;
7000cbbd2e1SMark F. Adams 
7010cbbd2e1SMark F. Adams       /* get block */
702e632b94dSBarry Smith       ierr = PetscMalloc5(Mdata*N, &qqc,M*N, &qqr,N, &TAU,LWORK, &WORK,M, &fids);CHKERRQ(ierr);
7032e68589bSMark F. Adams 
7042e68589bSMark F. Adams       aggID = 0;
705e78576d6SMark F. Adams       ierr  = PetscCDGetHeadPos(agg_llists,lid,&pos);CHKERRQ(ierr);
706e78576d6SMark F. Adams       while (pos) {
707e78576d6SMark F. Adams         PetscInt gid1;
708484f0a72SBarry Smith         ierr = PetscCDIntNdGetID(pos, &gid1);CHKERRQ(ierr);
709e78576d6SMark F. Adams         ierr = PetscCDGetNextPos(agg_llists,lid,&pos);CHKERRQ(ierr);
710e78576d6SMark F. Adams 
7110cbbd2e1SMark F. Adams         if (gid1 >= my0 && gid1 < Iend) flid = gid1 - my0;
7120cbbd2e1SMark F. Adams         else {
7131943db53SBarry Smith           ierr = PCGAMGHashTableFind(&fgid_flid, gid1, &flid);CHKERRQ(ierr);
7142c71b3e2SJacob Faibussowitsch           PetscCheckFalse(flid < 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot find gid1 in table");
7150cbbd2e1SMark F. Adams         }
7162e68589bSMark F. Adams         /* copy in B_i matrix - column oriented */
71765d7b583SSatish Balay         data = &data_in[flid*bs];
7183d3eaba7SBarry Smith         for (ii = 0; ii < bs; ii++) {
7192e68589bSMark F. Adams           for (jj = 0; jj < N; jj++) {
7200cbbd2e1SMark F. Adams             PetscReal d = data[jj*data_stride + ii];
7210cbbd2e1SMark F. Adams             qqc[jj*Mdata + aggID*bs + ii] = d;
7222e68589bSMark F. Adams           }
7232e68589bSMark F. Adams         }
7242e68589bSMark F. Adams         /* set fine IDs */
7252e68589bSMark F. Adams         for (kk=0; kk<bs; kk++) fids[aggID*bs + kk] = flid_fgid[flid]*bs + kk;
7262e68589bSMark F. Adams         aggID++;
7270cbbd2e1SMark F. Adams       }
7282e68589bSMark F. Adams 
7292e68589bSMark F. Adams       /* pad with zeros */
7302e68589bSMark F. Adams       for (ii = asz*bs; ii < Mdata; ii++) {
7312e68589bSMark F. Adams         for (jj = 0; jj < N; jj++, kk++) {
7322e68589bSMark F. Adams           qqc[jj*Mdata + ii] = .0;
7332e68589bSMark F. Adams         }
7342e68589bSMark F. Adams       }
7352e68589bSMark F. Adams 
7362e68589bSMark F. Adams       /* QR */
73784df3f90SSatish Balay       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
7388b83055fSJed Brown       PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
73984df3f90SSatish Balay       ierr = PetscFPTrapPop();CHKERRQ(ierr);
7402c71b3e2SJacob Faibussowitsch       PetscCheckFalse(INFO != 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"xGEQRF error");
7412e68589bSMark F. Adams       /* get R - column oriented - output B_{i+1} */
7422e68589bSMark F. Adams       {
7432e68589bSMark F. Adams         PetscReal *data = &out_data[clid*nSAvec];
7442e68589bSMark F. Adams         for (jj = 0; jj < nSAvec; jj++) {
7452e68589bSMark F. Adams           for (ii = 0; ii < nSAvec; ii++) {
7462c71b3e2SJacob Faibussowitsch             PetscCheckFalse(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);
7470cbbd2e1SMark F. Adams            if (ii <= jj) data[jj*out_data_stride + ii] = PetscRealPart(qqc[jj*Mdata + ii]);
7480cbbd2e1SMark F. Adams            else data[jj*out_data_stride + ii] = 0.;
7492e68589bSMark F. Adams           }
7502e68589bSMark F. Adams         }
7512e68589bSMark F. Adams       }
7522e68589bSMark F. Adams 
7532e68589bSMark F. Adams       /* get Q - row oriented */
754c964aadfSJose E. Roman       PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
7552c71b3e2SJacob Faibussowitsch       PetscCheckFalse(INFO != 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"xORGQR error arg %d",-INFO);
7562e68589bSMark F. Adams 
7572e68589bSMark F. Adams       for (ii = 0; ii < M; ii++) {
7582e68589bSMark F. Adams         for (jj = 0; jj < N; jj++) {
7592e68589bSMark F. Adams           qqr[N*ii + jj] = qqc[jj*Mdata + ii];
7602e68589bSMark F. Adams         }
7612e68589bSMark F. Adams       }
7622e68589bSMark F. Adams 
7632e68589bSMark F. Adams       /* add diagonal block of P0 */
764c8b0795cSMark F. Adams       for (kk=0; kk<N; kk++) {
765c8b0795cSMark F. Adams         cids[kk] = N*cgid + kk; /* global col IDs in P0 */
766c8b0795cSMark F. Adams       }
7672e68589bSMark F. Adams       ierr = MatSetValues(a_Prol,M,fids,N,cids,qqr,INSERT_VALUES);CHKERRQ(ierr);
768e632b94dSBarry Smith       ierr = PetscFree5(qqc,qqr,TAU,WORK,fids);CHKERRQ(ierr);
769b43b03e9SMark F. Adams       clid++;
7700cbbd2e1SMark F. Adams     } /* coarse agg */
7710cbbd2e1SMark F. Adams   } /* for all fine nodes */
7720cbbd2e1SMark F. Adams   ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7730cbbd2e1SMark F. Adams   ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7741943db53SBarry Smith   ierr = PCGAMGHashTableDestroy(&fgid_flid);CHKERRQ(ierr);
7752e68589bSMark F. Adams   PetscFunctionReturn(0);
7762e68589bSMark F. Adams }
7772e68589bSMark F. Adams 
7785adeb434SBarry Smith static PetscErrorCode PCView_GAMG_AGG(PC pc,PetscViewer viewer)
7795adeb434SBarry Smith {
7805adeb434SBarry Smith   PetscErrorCode ierr;
7815adeb434SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
7825adeb434SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
7835adeb434SBarry Smith   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
7845adeb434SBarry Smith 
7855adeb434SBarry Smith   PetscFunctionBegin;
7865adeb434SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"      AGG specific options\n");CHKERRQ(ierr);
7875adeb434SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"        Symmetric graph %s\n",pc_gamg_agg->sym_graph ? "true" : "false");CHKERRQ(ierr);
788cab9ed1eSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"        Number of levels to square graph %D\n",pc_gamg_agg->square_graph);CHKERRQ(ierr);
789cab9ed1eSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"        Number smoothing steps %D\n",pc_gamg_agg->nsmooths);CHKERRQ(ierr);
7905adeb434SBarry Smith   PetscFunctionReturn(0);
7915adeb434SBarry Smith }
7925adeb434SBarry Smith 
7932e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
7942e68589bSMark F. Adams /*
795fd1112cbSBarry Smith    PCGAMGGraph_AGG
7962e68589bSMark F. Adams 
7972e68589bSMark F. Adams   Input Parameter:
7982e68589bSMark F. Adams    . pc - this
7992e68589bSMark F. Adams    . Amat - matrix on this fine level
8002e68589bSMark F. Adams   Output Parameter:
801c8b0795cSMark F. Adams    . a_Gmat -
8022e68589bSMark F. Adams */
803e0877f53SBarry Smith static PetscErrorCode PCGAMGGraph_AGG(PC pc,Mat Amat,Mat *a_Gmat)
804c8b0795cSMark F. Adams {
805c8b0795cSMark F. Adams   PetscErrorCode            ierr;
806c8b0795cSMark F. Adams   PC_MG                     *mg          = (PC_MG*)pc->data;
807c8b0795cSMark F. Adams   PC_GAMG                   *pc_gamg     = (PC_GAMG*)mg->innerctx;
808c1eae691SMark Adams   const PetscReal           vfilter      = pc_gamg->threshold[pc_gamg->current_level];
809c8b0795cSMark F. Adams   PC_GAMG_AGG               *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
810e0940f08SMark F. Adams   Mat                       Gmat;
8113b4367a7SBarry Smith   MPI_Comm                  comm;
81267744fedSMark F. Adams   PetscBool /* set,flg , */ symm;
813c8b0795cSMark F. Adams 
814c8b0795cSMark F. Adams   PetscFunctionBegin;
8153b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr);
816fd1112cbSBarry Smith   ierr = PetscLogEventBegin(PC_GAMGGraph_AGG,0,0,0,0);CHKERRQ(ierr);
817c8b0795cSMark F. Adams 
81867744fedSMark F. Adams   /* ierr = MatIsSymmetricKnown(Amat, &set, &flg);CHKERRQ(ierr); || !(set && flg) -- this causes lot of symm calls */
819c666adbfSMark F. Adams   symm = (PetscBool)(pc_gamg_agg->sym_graph); /* && !pc_gamg_agg->square_graph; */
8200cbbd2e1SMark F. Adams 
8212d7fac45SMark F. Adams   ierr = PCGAMGCreateGraph(Amat, &Gmat);CHKERRQ(ierr);
822bf4339c2SBarry Smith   ierr = PCGAMGFilterGraph(&Gmat, vfilter, symm);CHKERRQ(ierr);
823e0940f08SMark F. Adams   *a_Gmat = Gmat;
824fd1112cbSBarry Smith   ierr = PetscLogEventEnd(PC_GAMGGraph_AGG,0,0,0,0);CHKERRQ(ierr);
825c8b0795cSMark F. Adams   PetscFunctionReturn(0);
826c8b0795cSMark F. Adams }
827c8b0795cSMark F. Adams 
828c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
829c8b0795cSMark F. Adams /*
830b43b03e9SMark F. Adams    PCGAMGCoarsen_AGG
831c8b0795cSMark F. Adams 
832c8b0795cSMark F. Adams   Input Parameter:
833e0940f08SMark F. Adams    . a_pc - this
834e0940f08SMark F. Adams   Input/Output Parameter:
8350cbbd2e1SMark F. Adams    . a_Gmat1 - graph on this fine level - coarsening can change this (squares it)
836c8b0795cSMark F. Adams   Output Parameter:
8370cbbd2e1SMark F. Adams    . agg_lists - list of aggregates
838c8b0795cSMark F. Adams */
839e0877f53SBarry Smith static PetscErrorCode PCGAMGCoarsen_AGG(PC a_pc,Mat *a_Gmat1,PetscCoarsenData **agg_lists)
840c8b0795cSMark F. Adams {
841c8b0795cSMark F. Adams   PetscErrorCode ierr;
842e0940f08SMark F. Adams   PC_MG          *mg          = (PC_MG*)a_pc->data;
843c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
844c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
8450cbbd2e1SMark F. Adams   Mat            mat,Gmat2, Gmat1 = *a_Gmat1;  /* squared graph */
8460cbbd2e1SMark F. Adams   IS             perm;
8475cfd4bd9SMark Adams   PetscInt       Istart,Iend,Ii,nloc,bs,n,m;
848c8b0795cSMark F. Adams   PetscInt       *permute;
849c8b0795cSMark F. Adams   PetscBool      *bIndexSet;
850b43b03e9SMark F. Adams   MatCoarsen     crs;
8513b4367a7SBarry Smith   MPI_Comm       comm;
852aea53286SMark Adams   PetscReal      hashfact;
853e2c00dcbSBarry Smith   PetscInt       iSwapIndex;
8543b50add6SMark Adams   PetscRandom    random;
855c8b0795cSMark F. Adams 
856c8b0795cSMark F. Adams   PetscFunctionBegin;
8570cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr);
8583b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Gmat1,&comm);CHKERRQ(ierr);
859e0940f08SMark F. Adams   ierr = MatGetLocalSize(Gmat1, &n, &m);CHKERRQ(ierr);
86071959b99SBarry Smith   ierr = MatGetBlockSize(Gmat1, &bs);CHKERRQ(ierr);
8612c71b3e2SJacob Faibussowitsch   PetscCheckFalse(bs != 1,PETSC_COMM_SELF,PETSC_ERR_PLIB,"bs %D must be 1",bs);
862c8b0795cSMark F. Adams   nloc = n/bs;
863c8b0795cSMark F. Adams 
8649ab59c8bSMark Adams   if (pc_gamg->current_level < pc_gamg_agg->square_graph) {
8654b1575e2SStefano Zampini     ierr = PCGAMGSquareGraph_GAMG(a_pc,Gmat1,&Gmat2);CHKERRQ(ierr);
866806fa848SBarry Smith   } else Gmat2 = Gmat1;
867c8b0795cSMark F. Adams 
8685cfd4bd9SMark Adams   /* get MIS aggs - randomize */
869785e854fSJed Brown   ierr = PetscMalloc1(nloc, &permute);CHKERRQ(ierr);
870e2c00dcbSBarry Smith   ierr = PetscCalloc1(nloc, &bIndexSet);CHKERRQ(ierr);
8716e09b0e3SStefano Zampini   for (Ii = 0; Ii < nloc; Ii++) permute[Ii] = Ii;
8723b50add6SMark Adams   ierr = PetscRandomCreate(PETSC_COMM_SELF,&random);CHKERRQ(ierr);
8735cfd4bd9SMark Adams   ierr = MatGetOwnershipRange(Gmat1, &Istart, &Iend);CHKERRQ(ierr);
874c8b0795cSMark F. Adams   for (Ii = 0; Ii < nloc; Ii++) {
8753b50add6SMark Adams     ierr = PetscRandomGetValueReal(random,&hashfact);CHKERRQ(ierr);
876aea53286SMark Adams     iSwapIndex = (PetscInt) (hashfact*nloc)%nloc;
877c8b0795cSMark F. Adams     if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) {
878c8b0795cSMark F. Adams       PetscInt iTemp = permute[iSwapIndex];
879c8b0795cSMark F. Adams       permute[iSwapIndex]   = permute[Ii];
880c8b0795cSMark F. Adams       permute[Ii]           = iTemp;
881c8b0795cSMark F. Adams       bIndexSet[iSwapIndex] = PETSC_TRUE;
882c8b0795cSMark F. Adams     }
883c8b0795cSMark F. Adams   }
884c8b0795cSMark F. Adams   ierr = PetscFree(bIndexSet);CHKERRQ(ierr);
8853b50add6SMark Adams   ierr = PetscRandomDestroy(&random);CHKERRQ(ierr);
886806fa848SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm);CHKERRQ(ierr);
8870cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr);
8883b4367a7SBarry Smith   ierr = MatCoarsenCreate(comm, &crs);CHKERRQ(ierr);
8899057884aSMark F. Adams   ierr = MatCoarsenSetFromOptions(crs);CHKERRQ(ierr);
890b43b03e9SMark F. Adams   ierr = MatCoarsenSetGreedyOrdering(crs, perm);CHKERRQ(ierr);
891b43b03e9SMark F. Adams   ierr = MatCoarsenSetAdjacency(crs, Gmat2);CHKERRQ(ierr);
892b43b03e9SMark F. Adams   ierr = MatCoarsenSetStrictAggs(crs, PETSC_TRUE);CHKERRQ(ierr);
893b43b03e9SMark F. Adams   ierr = MatCoarsenApply(crs);CHKERRQ(ierr);
8940cbbd2e1SMark F. Adams   ierr = MatCoarsenGetData(crs, agg_lists);CHKERRQ(ierr); /* output */
895b43b03e9SMark F. Adams   ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr);
896b43b03e9SMark F. Adams 
897c8b0795cSMark F. Adams   ierr = ISDestroy(&perm);CHKERRQ(ierr);
898c8b0795cSMark F. Adams   ierr = PetscFree(permute);CHKERRQ(ierr);
8990cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr);
9009f3f12c8SMark F. Adams 
901c8b0795cSMark F. Adams   /* smooth aggs */
902e0940f08SMark F. Adams   if (Gmat2 != Gmat1) {
9030cbbd2e1SMark F. Adams     const PetscCoarsenData *llist = *agg_lists;
904e3df7ea3SBarry Smith     ierr     = smoothAggs(a_pc,Gmat2, Gmat1, *agg_lists);CHKERRQ(ierr);
905c8b0795cSMark F. Adams     ierr     = MatDestroy(&Gmat1);CHKERRQ(ierr);
906e0940f08SMark F. Adams     *a_Gmat1 = Gmat2; /* output */
90741b27cdeSMark F. Adams     ierr     = PetscCDGetMat(llist, &mat);CHKERRQ(ierr);
9082c71b3e2SJacob Faibussowitsch     PetscCheckFalse(mat,comm,PETSC_ERR_ARG_WRONG, "Auxilary matrix with squared graph????");
909806fa848SBarry Smith   } else {
9100cbbd2e1SMark F. Adams     const PetscCoarsenData *llist = *agg_lists;
9119ab59c8bSMark Adams     /* see if we have a matrix that takes precedence (returned from MatCoarsenApply) */
91241b27cdeSMark F. Adams     ierr = PetscCDGetMat(llist, &mat);CHKERRQ(ierr);
9130cbbd2e1SMark F. Adams     if (mat) {
9140cbbd2e1SMark F. Adams       ierr     = MatDestroy(&Gmat1);CHKERRQ(ierr);
9150cbbd2e1SMark F. Adams       *a_Gmat1 = mat; /* output */
9160cbbd2e1SMark F. Adams     }
9170cbbd2e1SMark F. Adams   }
9180cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr);
919c8b0795cSMark F. Adams   PetscFunctionReturn(0);
920c8b0795cSMark F. Adams }
921c8b0795cSMark F. Adams 
922c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
923c8b0795cSMark F. Adams /*
9240cbbd2e1SMark F. Adams  PCGAMGProlongator_AGG
925c8b0795cSMark F. Adams 
926c8b0795cSMark F. Adams  Input Parameter:
927c8b0795cSMark F. Adams  . pc - this
928c8b0795cSMark F. Adams  . Amat - matrix on this fine level
929c8b0795cSMark F. Adams  . Graph - used to get ghost data for nodes in
9300cbbd2e1SMark F. Adams  . agg_lists - list of aggregates
931c8b0795cSMark F. Adams  Output Parameter:
932c8b0795cSMark F. Adams  . a_P_out - prolongation operator to the next level
933c8b0795cSMark F. Adams  */
934e0877f53SBarry Smith static PetscErrorCode PCGAMGProlongator_AGG(PC pc,Mat Amat,Mat Gmat,PetscCoarsenData *agg_lists,Mat *a_P_out)
9352e68589bSMark F. Adams {
9362e68589bSMark F. Adams   PC_MG          *mg       = (PC_MG*)pc->data;
9372e68589bSMark F. Adams   PC_GAMG        *pc_gamg  = (PC_GAMG*)mg->innerctx;
9384a2f8832SBarry Smith   const PetscInt col_bs = pc_gamg->data_cell_cols;
9392e68589bSMark F. Adams   PetscErrorCode ierr;
940c8b0795cSMark F. Adams   PetscInt       Istart,Iend,nloc,ii,jj,kk,my0,nLocalSelected,bs;
941c8b0795cSMark F. Adams   Mat            Prol;
942d5d25401SBarry Smith   PetscMPIInt    size;
9433b4367a7SBarry Smith   MPI_Comm       comm;
944c8b0795cSMark F. Adams   PetscReal      *data_w_ghost;
945c8b0795cSMark F. Adams   PetscInt       myCrs0, nbnodes=0, *flid_fgid;
946d9558ea9SBarry Smith   MatType        mtype;
9472e68589bSMark F. Adams 
9482e68589bSMark F. Adams   PetscFunctionBegin;
9493b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr);
9502c71b3e2SJacob Faibussowitsch   PetscCheckFalse(col_bs < 1,comm,PETSC_ERR_PLIB,"Column bs cannot be less than 1");
9510cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
952ffc4695bSBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr);
9532e68589bSMark F. Adams   ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr);
954c8b0795cSMark F. Adams   ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr);
95571959b99SBarry Smith   nloc = (Iend-Istart)/bs; my0 = Istart/bs;
9562c71b3e2SJacob Faibussowitsch   PetscCheckFalse((Iend-Istart) % bs,PETSC_COMM_SELF,PETSC_ERR_PLIB,"(Iend %D - Istart %D) not divisible by bs %D",Iend,Istart,bs);
9572e68589bSMark F. Adams 
9582e68589bSMark F. Adams   /* get 'nLocalSelected' */
9590cbbd2e1SMark F. Adams   for (ii=0, nLocalSelected = 0; ii < nloc; ii++) {
960e78576d6SMark F. Adams     PetscBool ise;
9610cbbd2e1SMark F. Adams     /* filter out singletons 0 or 1? */
962e78576d6SMark F. Adams     ierr = PetscCDEmptyAt(agg_lists, ii, &ise);CHKERRQ(ierr);
963e78576d6SMark F. Adams     if (!ise) nLocalSelected++;
9642e68589bSMark F. Adams   }
9652e68589bSMark F. Adams 
9662e68589bSMark F. Adams   /* create prolongator, create P matrix */
967d9558ea9SBarry Smith   ierr = MatGetType(Amat,&mtype);CHKERRQ(ierr);
9683b4367a7SBarry Smith   ierr = MatCreate(comm, &Prol);CHKERRQ(ierr);
969806fa848SBarry Smith   ierr = MatSetSizes(Prol,nloc*bs,nLocalSelected*col_bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
970a2f3521dSMark F. Adams   ierr = MatSetBlockSizes(Prol, bs, col_bs);CHKERRQ(ierr);
971d9558ea9SBarry Smith   ierr = MatSetType(Prol, mtype);CHKERRQ(ierr);
9724a2f8832SBarry Smith   ierr = MatSeqAIJSetPreallocation(Prol,col_bs, NULL);CHKERRQ(ierr);
9734a2f8832SBarry Smith   ierr = MatMPIAIJSetPreallocation(Prol,col_bs, NULL, col_bs, NULL);CHKERRQ(ierr);
9742e68589bSMark F. Adams 
9752e68589bSMark F. Adams   /* can get all points "removed" */
976c8b0795cSMark F. Adams   ierr =  MatGetSize(Prol, &kk, &ii);CHKERRQ(ierr);
9777f66b68fSMark Adams   if (!ii) {
978569f4572SMark Adams     ierr = PetscInfo(pc,"No selected points on coarse grid\n");CHKERRQ(ierr);
9792e68589bSMark F. Adams     ierr = MatDestroy(&Prol);CHKERRQ(ierr);
9800298fd71SBarry Smith     *a_P_out = NULL;  /* out */
981e87675ddSBarry Smith     ierr = PetscLogEventEnd(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
9822e68589bSMark F. Adams     PetscFunctionReturn(0);
9832e68589bSMark F. Adams   }
9847d3de750SJacob Faibussowitsch   ierr = PetscInfo(pc,"New grid %D nodes\n",ii/col_bs);CHKERRQ(ierr);
985c8b0795cSMark F. Adams   ierr = MatGetOwnershipRangeColumn(Prol, &myCrs0, &kk);CHKERRQ(ierr);
9860cbbd2e1SMark F. Adams 
9872c71b3e2SJacob Faibussowitsch   PetscCheckFalse((kk-myCrs0) % col_bs,PETSC_COMM_SELF,PETSC_ERR_PLIB,"(kk %D -myCrs0 %D) not divisible by col_bs %D",kk,myCrs0,col_bs);
988c8b0795cSMark F. Adams   myCrs0 = myCrs0/col_bs;
9892c71b3e2SJacob Faibussowitsch   PetscCheckFalse((kk/col_bs-myCrs0) != nLocalSelected,PETSC_COMM_SELF,PETSC_ERR_PLIB,"(kk %D/col_bs %D - myCrs0 %D) != nLocalSelected %D)",kk,col_bs,myCrs0,nLocalSelected);
9902e68589bSMark F. Adams 
9912e68589bSMark F. Adams   /* create global vector of data in 'data_w_ghost' */
9920cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr);
993c5df96a5SBarry Smith   if (size > 1) { /*  */
9942e68589bSMark F. Adams     PetscReal *tmp_gdata,*tmp_ldata,*tp2;
995785e854fSJed Brown     ierr = PetscMalloc1(nloc, &tmp_ldata);CHKERRQ(ierr);
9964a2f8832SBarry Smith     for (jj = 0; jj < col_bs; jj++) {
997c8b0795cSMark F. Adams       for (kk = 0; kk < bs; kk++) {
998a2f3521dSMark F. Adams         PetscInt        ii,stride;
999c8b0795cSMark F. Adams         const PetscReal *tp = pc_gamg->data + jj*bs*nloc + kk;
10002fa5cd67SKarl Rupp         for (ii = 0; ii < nloc; ii++, tp += bs) tmp_ldata[ii] = *tp;
10012fa5cd67SKarl Rupp 
1002806fa848SBarry Smith         ierr = PCGAMGGetDataWithGhosts(Gmat, 1, tmp_ldata, &stride, &tmp_gdata);CHKERRQ(ierr);
1003a2f3521dSMark F. Adams 
10047f66b68fSMark Adams         if (!jj && !kk) { /* now I know how many todal nodes - allocate */
10054a2f8832SBarry Smith           ierr    = PetscMalloc1(stride*bs*col_bs, &data_w_ghost);CHKERRQ(ierr);
1006a2f3521dSMark F. Adams           nbnodes = bs*stride;
10072e68589bSMark F. Adams         }
1008a2f3521dSMark F. Adams         tp2 = data_w_ghost + jj*bs*stride + kk;
10092fa5cd67SKarl Rupp         for (ii = 0; ii < stride; ii++, tp2 += bs) *tp2 = tmp_gdata[ii];
10102e68589bSMark F. Adams         ierr = PetscFree(tmp_gdata);CHKERRQ(ierr);
10112e68589bSMark F. Adams       }
10122e68589bSMark F. Adams     }
10132e68589bSMark F. Adams     ierr = PetscFree(tmp_ldata);CHKERRQ(ierr);
1014806fa848SBarry Smith   } else {
1015c8b0795cSMark F. Adams     nbnodes      = bs*nloc;
1016c8b0795cSMark F. Adams     data_w_ghost = (PetscReal*)pc_gamg->data;
10172e68589bSMark F. Adams   }
10182e68589bSMark F. Adams 
10192e68589bSMark F. Adams   /* get P0 */
1020c5df96a5SBarry Smith   if (size > 1) {
10212e68589bSMark F. Adams     PetscReal *fid_glid_loc,*fiddata;
1022a2f3521dSMark F. Adams     PetscInt  stride;
10232e68589bSMark F. Adams 
1024785e854fSJed Brown     ierr = PetscMalloc1(nloc, &fid_glid_loc);CHKERRQ(ierr);
10252e68589bSMark F. Adams     for (kk=0; kk<nloc; kk++) fid_glid_loc[kk] = (PetscReal)(my0+kk);
1026806fa848SBarry Smith     ierr = PCGAMGGetDataWithGhosts(Gmat, 1, fid_glid_loc, &stride, &fiddata);CHKERRQ(ierr);
1027785e854fSJed Brown     ierr = PetscMalloc1(stride, &flid_fgid);CHKERRQ(ierr);
1028a2f3521dSMark F. Adams     for (kk=0; kk<stride; kk++) flid_fgid[kk] = (PetscInt)fiddata[kk];
10292e68589bSMark F. Adams     ierr = PetscFree(fiddata);CHKERRQ(ierr);
1030a2f3521dSMark F. Adams 
10312c71b3e2SJacob Faibussowitsch     PetscCheckFalse(stride != nbnodes/bs,PETSC_COMM_SELF,PETSC_ERR_PLIB,"stride %D != nbnodes %D/bs %D",stride,nbnodes,bs);
10322e68589bSMark F. Adams     ierr = PetscFree(fid_glid_loc);CHKERRQ(ierr);
1033806fa848SBarry Smith   } else {
1034785e854fSJed Brown     ierr = PetscMalloc1(nloc, &flid_fgid);CHKERRQ(ierr);
10352e68589bSMark F. Adams     for (kk=0; kk<nloc; kk++) flid_fgid[kk] = my0 + kk;
10362e68589bSMark F. Adams   }
10370cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr);
10380cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr);
1039c8b0795cSMark F. Adams   {
10400298fd71SBarry Smith     PetscReal *data_out = NULL;
10414a2f8832SBarry Smith     ierr = formProl0(agg_lists, bs, col_bs, myCrs0, nbnodes,data_w_ghost, flid_fgid, &data_out, Prol);CHKERRQ(ierr);
1042c8b0795cSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
10432fa5cd67SKarl Rupp 
1044c8b0795cSMark F. Adams     pc_gamg->data           = data_out;
10454a2f8832SBarry Smith     pc_gamg->data_cell_rows = col_bs;
10464a2f8832SBarry Smith     pc_gamg->data_sz        = col_bs*col_bs*nLocalSelected;
1047c8b0795cSMark F. Adams   }
10480cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr);
104961ba4676SBarry Smith   if (size > 1) {ierr = PetscFree(data_w_ghost);CHKERRQ(ierr);}
10502e68589bSMark F. Adams   ierr = PetscFree(flid_fgid);CHKERRQ(ierr);
10512e68589bSMark F. Adams 
1052c8b0795cSMark F. Adams   *a_P_out = Prol;  /* out */
1053ed4e82eaSPeter Brune 
10540cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
1055c8b0795cSMark F. Adams   PetscFunctionReturn(0);
1056c8b0795cSMark F. Adams }
1057c8b0795cSMark F. Adams 
1058c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
1059c8b0795cSMark F. Adams /*
1060fd1112cbSBarry Smith    PCGAMGOptProlongator_AGG
1061c8b0795cSMark F. Adams 
1062c8b0795cSMark F. Adams   Input Parameter:
1063c8b0795cSMark F. Adams    . pc - this
1064c8b0795cSMark F. Adams    . Amat - matrix on this fine level
1065c8b0795cSMark F. Adams  In/Output Parameter:
10662a808120SBarry Smith    . a_P - prolongation operator to the next level
1067c8b0795cSMark F. Adams */
1068e0877f53SBarry Smith static PetscErrorCode PCGAMGOptProlongator_AGG(PC pc,Mat Amat,Mat *a_P)
1069c8b0795cSMark F. Adams {
1070c8b0795cSMark F. Adams   PetscErrorCode ierr;
1071c8b0795cSMark F. Adams   PC_MG          *mg          = (PC_MG*)pc->data;
1072c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
1073c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
1074c8b0795cSMark F. Adams   PetscInt       jj;
1075c8b0795cSMark F. Adams   Mat            Prol  = *a_P;
10763b4367a7SBarry Smith   MPI_Comm       comm;
10772a808120SBarry Smith   KSP            eksp;
10782a808120SBarry Smith   Vec            bb, xx;
10792a808120SBarry Smith   PC             epc;
10802a808120SBarry Smith   PetscReal      alpha, emax, emin;
1081c8b0795cSMark F. Adams 
1082c8b0795cSMark F. Adams   PetscFunctionBegin;
10833b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr);
1084fd1112cbSBarry Smith   ierr = PetscLogEventBegin(PC_GAMGOptProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
1085c8b0795cSMark F. Adams 
1086d5d25401SBarry Smith   /* compute maximum singular value of operator to be used in smoother */
10872a808120SBarry Smith   if (0 < pc_gamg_agg->nsmooths) {
108818c3aa7eSMark     /* get eigen estimates */
108918c3aa7eSMark     if (pc_gamg->emax > 0) {
109018c3aa7eSMark       emin = pc_gamg->emin;
109118c3aa7eSMark       emax = pc_gamg->emax;
109218c3aa7eSMark     } else {
10930ed2132dSStefano Zampini       const char *prefix;
10940ed2132dSStefano Zampini 
10950a545947SLisandro Dalcin       ierr = MatCreateVecs(Amat, &bb, NULL);CHKERRQ(ierr);
10960a545947SLisandro Dalcin       ierr = MatCreateVecs(Amat, &xx, NULL);CHKERRQ(ierr);
109773f7197eSJed Brown       ierr = KSPSetNoisy_Private(bb);CHKERRQ(ierr);
1098e696ed0bSMark F. Adams 
10993b4367a7SBarry Smith       ierr = KSPCreate(comm,&eksp);CHKERRQ(ierr);
11000ed2132dSStefano Zampini       ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
11010ed2132dSStefano Zampini       ierr = KSPSetOptionsPrefix(eksp,prefix);CHKERRQ(ierr);
110273f7197eSJed Brown       ierr = KSPAppendOptionsPrefix(eksp,"pc_gamg_esteig_");CHKERRQ(ierr);
110373f7197eSJed Brown       {
1104acbb5b45SJed Brown         PetscBool sflg;
110590db8557SMark Adams         ierr = MatGetOption(Amat, MAT_SPD, &sflg);CHKERRQ(ierr);
110690db8557SMark Adams         if (sflg) {
1107d24ecf33SMark           ierr = KSPSetType(eksp, KSPCG);CHKERRQ(ierr);
1108d24ecf33SMark         }
1109d24ecf33SMark       }
1110c0decd05SBarry Smith       ierr = KSPSetErrorIfNotConverged(eksp,pc->erroriffailure);CHKERRQ(ierr);
11112e68589bSMark F. Adams       ierr = KSPSetNormType(eksp, KSP_NORM_NONE);CHKERRQ(ierr);
1112c2436cacSMark F. Adams 
1113c2436cacSMark F. Adams       ierr = KSPSetInitialGuessNonzero(eksp, PETSC_FALSE);CHKERRQ(ierr);
111423ee1639SBarry Smith       ierr = KSPSetOperators(eksp, Amat, Amat);CHKERRQ(ierr);
11152e68589bSMark F. Adams 
1116da509fc8SJed Brown       ierr = KSPGetPC(eksp, &epc);CHKERRQ(ierr);
1117da509fc8SJed Brown       ierr = PCSetType(epc, PCJACOBI);CHKERRQ(ierr);  /* smoother in smoothed agg. */
1118c2436cacSMark F. Adams 
1119*074aec55SMark Adams       ierr = KSPSetTolerances(eksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, 10);CHKERRQ(ierr); // 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
1120*074aec55SMark Adams 
11210ed2132dSStefano Zampini       ierr = KSPSetFromOptions(eksp);CHKERRQ(ierr);
11220ed2132dSStefano Zampini       ierr = KSPSetComputeSingularValues(eksp,PETSC_TRUE);CHKERRQ(ierr);
11232e68589bSMark F. Adams       ierr = KSPSolve(eksp, bb, xx);CHKERRQ(ierr);
1124c0decd05SBarry Smith       ierr = KSPCheckSolve(eksp,pc,xx);CHKERRQ(ierr);
11252e68589bSMark F. Adams 
11262e68589bSMark F. Adams       ierr = KSPComputeExtremeSingularValues(eksp, &emax, &emin);CHKERRQ(ierr);
11277d3de750SJacob Faibussowitsch       ierr = PetscInfo(pc,"Smooth P0: max eigen=%e min=%e PC=%s\n",(double)emax,(double)emin,PCJACOBI);CHKERRQ(ierr);
11282e68589bSMark F. Adams       ierr = VecDestroy(&xx);CHKERRQ(ierr);
11292e68589bSMark F. Adams       ierr = VecDestroy(&bb);CHKERRQ(ierr);
11302e68589bSMark F. Adams       ierr = KSPDestroy(&eksp);CHKERRQ(ierr);
11312e68589bSMark F. Adams     }
113218c3aa7eSMark     if (pc_gamg->use_sa_esteig) {
113318c3aa7eSMark       mg->min_eigen_DinvA[pc_gamg->current_level] = emin;
113418c3aa7eSMark       mg->max_eigen_DinvA[pc_gamg->current_level] = emax;
11357d3de750SJacob Faibussowitsch       ierr = PetscInfo(pc,"Smooth P0: level %D, cache spectra %g %g\n",pc_gamg->current_level,(double)emin,(double)emax);CHKERRQ(ierr);
113618c3aa7eSMark     } else {
113718c3aa7eSMark       mg->min_eigen_DinvA[pc_gamg->current_level] = 0;
113818c3aa7eSMark       mg->max_eigen_DinvA[pc_gamg->current_level] = 0;
113918c3aa7eSMark     }
114018c3aa7eSMark   } else {
114118c3aa7eSMark     mg->min_eigen_DinvA[pc_gamg->current_level] = 0;
114218c3aa7eSMark     mg->max_eigen_DinvA[pc_gamg->current_level] = 0;
114318c3aa7eSMark   }
11442e68589bSMark F. Adams 
11452a808120SBarry Smith   /* smooth P0 */
11462a808120SBarry Smith   for (jj = 0; jj < pc_gamg_agg->nsmooths; jj++) {
11472a808120SBarry Smith     Mat tMat;
11482a808120SBarry Smith     Vec diag;
11492a808120SBarry Smith 
11502a808120SBarry Smith     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr);
11512a808120SBarry Smith 
11522e68589bSMark F. Adams     /* smooth P1 := (I - omega/lam D^{-1}A)P0 */
11534555aa8cSStefano Zampini     ierr = PetscLogEventBegin(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2],0,0,0,0);CHKERRQ(ierr);
11542e68589bSMark F. Adams     ierr = MatMatMult(Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat);CHKERRQ(ierr);
11554555aa8cSStefano Zampini     ierr = PetscLogEventEnd(petsc_gamg_setup_matmat_events[pc_gamg->current_level][2],0,0,0,0);CHKERRQ(ierr);
1156b4da3a1bSStefano Zampini     ierr = MatProductClear(tMat);CHKERRQ(ierr);
11570a545947SLisandro Dalcin     ierr = MatCreateVecs(Amat, &diag, NULL);CHKERRQ(ierr);
11582e68589bSMark F. Adams     ierr = MatGetDiagonal(Amat, diag);CHKERRQ(ierr); /* effectively PCJACOBI */
11592e68589bSMark F. Adams     ierr = VecReciprocal(diag);CHKERRQ(ierr);
11600a545947SLisandro Dalcin     ierr = MatDiagonalScale(tMat, diag, NULL);CHKERRQ(ierr);
11612e68589bSMark F. Adams     ierr = VecDestroy(&diag);CHKERRQ(ierr);
1162b4da3a1bSStefano Zampini 
1163d5d25401SBarry Smith     /* TODO: Set a PCFailedReason and exit the building of the AMG preconditioner */
11642c71b3e2SJacob Faibussowitsch     PetscCheckFalse(emax == 0.0,PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"Computed maximum singular value as zero");
1165d5d25401SBarry Smith     /* TODO: Document the 1.4 and don't hardwire it in this routine */
1166b4ec6429SMark F. Adams     alpha = -1.4/emax;
1167b4da3a1bSStefano Zampini 
11682e68589bSMark F. Adams     ierr = MatAYPX(tMat, alpha, Prol, SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
11692e68589bSMark F. Adams     ierr = MatDestroy(&Prol);CHKERRQ(ierr);
11702e68589bSMark F. Adams     Prol = tMat;
11710cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr);
11722e68589bSMark F. Adams   }
1173fd1112cbSBarry Smith   ierr = PetscLogEventEnd(PC_GAMGOptProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
1174c8b0795cSMark F. Adams   *a_P = Prol;
1175c8b0795cSMark F. Adams   PetscFunctionReturn(0);
11762e68589bSMark F. Adams }
11772e68589bSMark F. Adams 
1178c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
1179c8b0795cSMark F. Adams /*
1180c8b0795cSMark F. Adams    PCCreateGAMG_AGG
11812e68589bSMark F. Adams 
1182c8b0795cSMark F. Adams   Input Parameter:
1183c8b0795cSMark F. Adams    . pc -
1184c8b0795cSMark F. Adams */
1185c8b0795cSMark F. Adams PetscErrorCode  PCCreateGAMG_AGG(PC pc)
1186c8b0795cSMark F. Adams {
1187c8b0795cSMark F. Adams   PetscErrorCode ierr;
1188c8b0795cSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
1189c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1190c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg;
11912e68589bSMark F. Adams 
1192c8b0795cSMark F. Adams   PetscFunctionBegin;
1193c8b0795cSMark F. Adams   /* create sub context for SA */
1194b00a9115SJed Brown   ierr            = PetscNewLog(pc,&pc_gamg_agg);CHKERRQ(ierr);
1195c8b0795cSMark F. Adams   pc_gamg->subctx = pc_gamg_agg;
1196c8b0795cSMark F. Adams 
11971ab5ffc9SJed Brown   pc_gamg->ops->setfromoptions = PCSetFromOptions_GAMG_AGG;
11989b8ffb57SJed Brown   pc_gamg->ops->destroy        = PCDestroy_GAMG_AGG;
1199c8b0795cSMark F. Adams   /* reset does not do anything; setup not virtual */
1200c8b0795cSMark F. Adams 
1201c8b0795cSMark F. Adams   /* set internal function pointers */
1202fd1112cbSBarry Smith   pc_gamg->ops->graph             = PCGAMGGraph_AGG;
12031ab5ffc9SJed Brown   pc_gamg->ops->coarsen           = PCGAMGCoarsen_AGG;
12041ab5ffc9SJed Brown   pc_gamg->ops->prolongator       = PCGAMGProlongator_AGG;
1205fd1112cbSBarry Smith   pc_gamg->ops->optprolongator    = PCGAMGOptProlongator_AGG;
12061ab5ffc9SJed Brown   pc_gamg->ops->createdefaultdata = PCSetData_AGG;
12075adeb434SBarry Smith   pc_gamg->ops->view              = PCView_GAMG_AGG;
1208c8b0795cSMark F. Adams 
1209cab9ed1eSBarry Smith   pc_gamg_agg->square_graph = 1;
1210cab9ed1eSBarry Smith   pc_gamg_agg->sym_graph    = PETSC_FALSE;
1211cab9ed1eSBarry Smith   pc_gamg_agg->nsmooths     = 1;
1212cab9ed1eSBarry Smith 
1213e0877f53SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNSmooths_C",PCGAMGSetNSmooths_AGG);CHKERRQ(ierr);
1214e0877f53SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetSymGraph_C",PCGAMGSetSymGraph_AGG);CHKERRQ(ierr);
1215e0877f53SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetSquareGraph_C",PCGAMGSetSquareGraph_AGG);CHKERRQ(ierr);
1216bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_AGG);CHKERRQ(ierr);
12172e68589bSMark F. Adams   PetscFunctionReturn(0);
12182e68589bSMark F. Adams }
1219