xref: /petsc/src/ksp/pc/impls/gamg/agg.c (revision b817416efb44ab96f162295e4cd7ec72f89365b8)
12e68589bSMark F. Adams /*
2*b817416eSBarry 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*/
6539c167fSBarry Smith /*  Next line needed to deactivate KSP_Solve logging */
7af0996ceSBarry Smith #include <petsc/private/kspimpl.h>
82e68589bSMark F. Adams #include <petscblaslapack.h>
9539c167fSBarry Smith #include <petscdm.h>
102e68589bSMark F. Adams 
112e68589bSMark F. Adams typedef struct {
12c8b0795cSMark F. Adams   PetscInt  nsmooths;
13c8b0795cSMark F. Adams   PetscBool sym_graph;
149ab59c8bSMark Adams   PetscInt  square_graph;
152e68589bSMark F. Adams } PC_GAMG_AGG;
162e68589bSMark F. Adams 
172e68589bSMark F. Adams /*@
182e68589bSMark F. Adams    PCGAMGSetNSmooths - Set number of smoothing steps (1 is typical)
192e68589bSMark F. Adams 
202e68589bSMark F. Adams    Not Collective on PC
212e68589bSMark F. Adams 
222e68589bSMark F. Adams    Input Parameters:
232e68589bSMark F. Adams .  pc - the preconditioner context
242e68589bSMark F. Adams 
252e68589bSMark F. Adams    Options Database Key:
26cab9ed1eSBarry Smith .  -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation
272e68589bSMark F. Adams 
282e68589bSMark F. Adams    Level: intermediate
292e68589bSMark F. Adams 
302e68589bSMark F. Adams    Concepts: Aggregation AMG preconditioner
312e68589bSMark F. Adams 
322e68589bSMark F. Adams .seealso: ()
332e68589bSMark F. Adams @*/
342e68589bSMark F. Adams PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n)
352e68589bSMark F. Adams {
362e68589bSMark F. Adams   PetscErrorCode ierr;
372e68589bSMark F. Adams 
382e68589bSMark F. Adams   PetscFunctionBegin;
392e68589bSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
402e68589bSMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetNSmooths_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
412e68589bSMark F. Adams   PetscFunctionReturn(0);
422e68589bSMark F. Adams }
432e68589bSMark F. Adams 
44e0877f53SBarry Smith static PetscErrorCode PCGAMGSetNSmooths_AGG(PC pc, PetscInt n)
452e68589bSMark F. Adams {
462e68589bSMark F. Adams   PC_MG       *mg          = (PC_MG*)pc->data;
472e68589bSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG*)mg->innerctx;
48c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
492e68589bSMark F. Adams 
502e68589bSMark F. Adams   PetscFunctionBegin;
51c8b0795cSMark F. Adams   pc_gamg_agg->nsmooths = n;
52c8b0795cSMark F. Adams   PetscFunctionReturn(0);
53c8b0795cSMark F. Adams }
54c8b0795cSMark F. Adams 
55c8b0795cSMark F. Adams /*@
56cab9ed1eSBarry Smith    PCGAMGSetSymGraph - Symmetrize the graph before computing the aggregation. Some algorithms require the graph be symmetric
57c8b0795cSMark F. Adams 
58c8b0795cSMark F. Adams    Not Collective on PC
59c8b0795cSMark F. Adams 
60c8b0795cSMark F. Adams    Input Parameters:
61cab9ed1eSBarry Smith +  pc - the preconditioner context
62cab9ed1eSBarry Smith .  n - PETSC_TRUE or PETSC_FALSE
63c8b0795cSMark F. Adams 
64c8b0795cSMark F. Adams    Options Database Key:
65cab9ed1eSBarry Smith .  -pc_gamg_sym_graph <true,default=false> - symmetrize the graph before computing the aggregation
66c8b0795cSMark F. Adams 
67c8b0795cSMark F. Adams    Level: intermediate
68c8b0795cSMark F. Adams 
69c8b0795cSMark F. Adams    Concepts: Aggregation AMG preconditioner
70c8b0795cSMark F. Adams 
71cab9ed1eSBarry Smith .seealso: PCGAMGSetSquareGraph()
72c8b0795cSMark F. Adams @*/
73c8b0795cSMark F. Adams PetscErrorCode PCGAMGSetSymGraph(PC pc, PetscBool n)
74c8b0795cSMark F. Adams {
75c8b0795cSMark F. Adams   PetscErrorCode ierr;
76c8b0795cSMark F. Adams 
77c8b0795cSMark F. Adams   PetscFunctionBegin;
78c8b0795cSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
79c8b0795cSMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetSymGraph_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
80c8b0795cSMark F. Adams   PetscFunctionReturn(0);
81c8b0795cSMark F. Adams }
82c8b0795cSMark F. Adams 
83e0877f53SBarry Smith static PetscErrorCode PCGAMGSetSymGraph_AGG(PC pc, PetscBool n)
84c8b0795cSMark F. Adams {
85c8b0795cSMark F. Adams   PC_MG       *mg          = (PC_MG*)pc->data;
86c8b0795cSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG*)mg->innerctx;
87c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
88c8b0795cSMark F. Adams 
89c8b0795cSMark F. Adams   PetscFunctionBegin;
90c8b0795cSMark F. Adams   pc_gamg_agg->sym_graph = n;
912e68589bSMark F. Adams   PetscFunctionReturn(0);
922e68589bSMark F. Adams }
932e68589bSMark F. Adams 
94ef4ad70eSMark F. Adams /*@
95cab9ed1eSBarry Smith    PCGAMGSetSquareGraph -  Square the graph, ie. compute A'*A before aggregating it
96ef4ad70eSMark F. Adams 
97ef4ad70eSMark F. Adams    Not Collective on PC
98ef4ad70eSMark F. Adams 
99ef4ad70eSMark F. Adams    Input Parameters:
100cab9ed1eSBarry Smith +  pc - the preconditioner context
101cab9ed1eSBarry Smith -  n - PETSC_TRUE or PETSC_FALSE
102ef4ad70eSMark F. Adams 
103ef4ad70eSMark F. Adams    Options Database Key:
104cab9ed1eSBarry Smith .  -pc_gamg_square_graph <n,default = 1> - number of levels to square the graph on before aggregating it
105ef4ad70eSMark F. Adams 
106ef4ad70eSMark F. Adams    Level: intermediate
107ef4ad70eSMark F. Adams 
108ef4ad70eSMark F. Adams    Concepts: Aggregation AMG preconditioner
109ef4ad70eSMark F. Adams 
110cab9ed1eSBarry Smith .seealso: PCGAMGSetSymGraph()
111ef4ad70eSMark F. Adams @*/
1129ab59c8bSMark Adams PetscErrorCode PCGAMGSetSquareGraph(PC pc, PetscInt n)
113ef4ad70eSMark F. Adams {
114ef4ad70eSMark F. Adams   PetscErrorCode ierr;
115ef4ad70eSMark F. Adams 
116ef4ad70eSMark F. Adams   PetscFunctionBegin;
117ef4ad70eSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1189ab59c8bSMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetSquareGraph_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
119ef4ad70eSMark F. Adams   PetscFunctionReturn(0);
120ef4ad70eSMark F. Adams }
121ef4ad70eSMark F. Adams 
122e0877f53SBarry Smith static PetscErrorCode PCGAMGSetSquareGraph_AGG(PC pc, PetscInt n)
123ef4ad70eSMark F. Adams {
124ef4ad70eSMark F. Adams   PC_MG       *mg          = (PC_MG*)pc->data;
125ef4ad70eSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG*)mg->innerctx;
126ef4ad70eSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
127ef4ad70eSMark F. Adams 
128ef4ad70eSMark F. Adams   PetscFunctionBegin;
129ef4ad70eSMark F. Adams   pc_gamg_agg->square_graph = n;
130ef4ad70eSMark F. Adams   PetscFunctionReturn(0);
131ef4ad70eSMark F. Adams }
132ef4ad70eSMark F. Adams 
133e0877f53SBarry Smith static PetscErrorCode PCSetFromOptions_GAMG_AGG(PetscOptionItems *PetscOptionsObject,PC pc)
1342e68589bSMark F. Adams {
1352e68589bSMark F. Adams   PetscErrorCode ierr;
1362e68589bSMark F. Adams   PC_MG          *mg          = (PC_MG*)pc->data;
1372e68589bSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
138c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
1392e68589bSMark F. Adams 
1402e68589bSMark F. Adams   PetscFunctionBegin;
141e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"GAMG-AGG options");CHKERRQ(ierr);
1422e68589bSMark F. Adams   {
1438afaa268SBarry 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);
1448afaa268SBarry 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);
1459ab59c8bSMark 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);
1462e68589bSMark F. Adams   }
1472e68589bSMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
1482e68589bSMark F. Adams   PetscFunctionReturn(0);
1492e68589bSMark F. Adams }
1502e68589bSMark F. Adams 
1512e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
152e0877f53SBarry Smith static PetscErrorCode PCDestroy_GAMG_AGG(PC pc)
1532e68589bSMark F. Adams {
1542e68589bSMark F. Adams   PetscErrorCode ierr;
1552e68589bSMark F. Adams   PC_MG          *mg          = (PC_MG*)pc->data;
1562e68589bSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
1572e68589bSMark F. Adams 
1582e68589bSMark F. Adams   PetscFunctionBegin;
1599b8ffb57SJed Brown   ierr = PetscFree(pc_gamg->subctx);CHKERRQ(ierr);
1603ae0bb68SMark Adams   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",NULL);CHKERRQ(ierr);
1612e68589bSMark F. Adams   PetscFunctionReturn(0);
1622e68589bSMark F. Adams }
1632e68589bSMark F. Adams 
1642e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
1652e68589bSMark F. Adams /*
1662e68589bSMark F. Adams    PCSetCoordinates_AGG
167302f38e8SMark F. Adams      - collective
1682e68589bSMark F. Adams 
1692e68589bSMark F. Adams    Input Parameter:
1702e68589bSMark F. Adams    . pc - the preconditioner context
171a2f3521dSMark F. Adams    . ndm - dimesion of data (used for dof/vertex for Stokes)
172302f38e8SMark F. Adams    . a_nloc - number of vertices local
173302f38e8SMark F. Adams    . coords - [a_nloc][ndm] - interleaved coordinate data: {x_0, y_0, z_0, x_1, y_1, ...}
1742e68589bSMark F. Adams */
1751e6b0712SBarry Smith 
1761e6b0712SBarry Smith static PetscErrorCode PCSetCoordinates_AGG(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords)
1772e68589bSMark F. Adams {
1782e68589bSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
1792e68589bSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1802e68589bSMark F. Adams   PetscErrorCode ierr;
18169344418SMark F. Adams   PetscInt       arrsz,kk,ii,jj,nloc,ndatarows,ndf;
182a2f3521dSMark F. Adams   Mat            mat = pc->pmat;
1832e68589bSMark F. Adams 
1842e68589bSMark F. Adams   PetscFunctionBegin;
185a2f3521dSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
186a2f3521dSMark F. Adams   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
187302f38e8SMark F. Adams   nloc = a_nloc;
1882e68589bSMark F. Adams 
1892e68589bSMark F. Adams   /* SA: null space vectors */
19069344418SMark F. Adams   ierr = MatGetBlockSize(mat, &ndf);CHKERRQ(ierr); /* this does not work for Stokes */
19169344418SMark F. Adams   if (coords && ndf==1) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */
192a2f3521dSMark F. Adams   else if (coords) {
193*b817416eSBarry Smith     if (ndm > ndf) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"degrees of motion %D > block size %D",ndm,ndf);
19469344418SMark F. Adams     pc_gamg->data_cell_cols = (ndm==2 ? 3 : 6); /* displacement elasticity */
19569344418SMark F. Adams     if (ndm != ndf) {
196*b817416eSBarry Smith       if (pc_gamg->data_cell_cols != ndf) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Don't know how to create null space for ndm=%D, ndf=%D.  Use MatSetNearNullSpace.",ndm,ndf);
197a2f3521dSMark F. Adams     }
198806fa848SBarry Smith   } else pc_gamg->data_cell_cols = ndf; /* no data, force SA with constant null space vectors */
19971959b99SBarry Smith   pc_gamg->data_cell_rows = ndatarows = ndf;
20071959b99SBarry Smith   if (pc_gamg->data_cell_cols <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"pc_gamg->data_cell_cols %D <= 0",pc_gamg->data_cell_cols);
201c8b0795cSMark F. Adams   arrsz = nloc*pc_gamg->data_cell_rows*pc_gamg->data_cell_cols;
2022e68589bSMark F. Adams 
2032e68589bSMark F. Adams   /* create data - syntactic sugar that should be refactored at some point */
2047f66b68fSMark Adams   if (!pc_gamg->data || (pc_gamg->data_sz != arrsz)) {
2052e68589bSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
206854ce69bSBarry Smith     ierr = PetscMalloc1(arrsz+1, &pc_gamg->data);CHKERRQ(ierr);
2072e68589bSMark F. Adams   }
2082e68589bSMark F. Adams   /* copy data in - column oriented */
2092e68589bSMark F. Adams   for (kk=0; kk<nloc; kk++) {
21069344418SMark F. Adams     const PetscInt M     = nloc*pc_gamg->data_cell_rows; /* stride into data */
21169344418SMark F. Adams     PetscReal      *data = &pc_gamg->data[kk*ndatarows]; /* start of cell */
212c8b0795cSMark F. Adams     if (pc_gamg->data_cell_cols==1) *data = 1.0;
2132e68589bSMark F. Adams     else {
21469344418SMark F. Adams       /* translational modes */
2152fa5cd67SKarl Rupp       for (ii=0;ii<ndatarows;ii++) {
2162fa5cd67SKarl Rupp         for (jj=0;jj<ndatarows;jj++) {
21769344418SMark F. Adams           if (ii==jj)data[ii*M + jj] = 1.0;
2182e68589bSMark F. Adams           else data[ii*M + jj] = 0.0;
2192fa5cd67SKarl Rupp         }
2202fa5cd67SKarl Rupp       }
22169344418SMark F. Adams 
22269344418SMark F. Adams       /* rotational modes */
2232e68589bSMark F. Adams       if (coords) {
22469344418SMark F. Adams         if (ndm == 2) {
2252e68589bSMark F. Adams           data   += 2*M;
2262e68589bSMark F. Adams           data[0] = -coords[2*kk+1];
2272e68589bSMark F. Adams           data[1] =  coords[2*kk];
228806fa848SBarry Smith         } else {
2292e68589bSMark F. Adams           data   += 3*M;
2302e68589bSMark F. Adams           data[0] = 0.0;             data[M+0] =  coords[3*kk+2]; data[2*M+0] = -coords[3*kk+1];
2312e68589bSMark F. Adams           data[1] = -coords[3*kk+2]; data[M+1] = 0.0;             data[2*M+1] =  coords[3*kk];
2322e68589bSMark F. Adams           data[2] =  coords[3*kk+1]; data[M+2] = -coords[3*kk];   data[2*M+2] = 0.0;
2332e68589bSMark F. Adams         }
2342e68589bSMark F. Adams       }
2352e68589bSMark F. Adams     }
2362e68589bSMark F. Adams   }
2372e68589bSMark F. Adams 
2382e68589bSMark F. Adams   pc_gamg->data_sz = arrsz;
2392e68589bSMark F. Adams   PetscFunctionReturn(0);
2402e68589bSMark F. Adams }
2412e68589bSMark F. Adams 
242b43b03e9SMark F. Adams typedef PetscInt NState;
243b43b03e9SMark F. Adams static const NState NOT_DONE=-2;
244b43b03e9SMark F. Adams static const NState DELETED =-1;
245b43b03e9SMark F. Adams static const NState REMOVED =-3;
246b43b03e9SMark F. Adams #define IS_SELECTED(s) (s!=DELETED && s!=NOT_DONE && s!=REMOVED)
247b43b03e9SMark F. Adams 
248c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
249c8b0795cSMark F. Adams /*
250b43b03e9SMark F. Adams    smoothAggs - greedy grab of with G1 (unsquared graph) -- AIJ specific
251b43b03e9SMark F. Adams      - AGG-MG specific: clears singletons out of 'selected_2'
252c8b0795cSMark F. Adams 
253c8b0795cSMark F. Adams    Input Parameter:
2542adfe9d3SBarry Smith    . Gmat_2 - glabal matrix of graph (data not defined)   base (squared) graph
2552adfe9d3SBarry Smith    . Gmat_1 - base graph to grab with                 base graph
256c8b0795cSMark F. Adams    Input/Output Parameter:
2570cbbd2e1SMark F. Adams    . aggs_2 - linked list of aggs with gids)
258c8b0795cSMark F. Adams */
259e3df7ea3SBarry Smith static PetscErrorCode smoothAggs(PC pc,Mat Gmat_2, Mat Gmat_1,PetscCoarsenData *aggs_2)
260c8b0795cSMark F. Adams {
261c8b0795cSMark F. Adams   PetscErrorCode ierr;
262c8b0795cSMark F. Adams   PetscBool      isMPI;
2633fa9c902SMark Adams   Mat_SeqAIJ     *matA_1, *matB_1=0;
2643b4367a7SBarry Smith   MPI_Comm       comm;
2650cbbd2e1SMark F. Adams   PetscInt       lid,*ii,*idx,ix,Iend,my0,kk,n,j;
266c8b0795cSMark F. Adams   Mat_MPIAIJ     *mpimat_2 = 0, *mpimat_1=0;
267c8b0795cSMark F. Adams   const PetscInt nloc      = Gmat_2->rmap->n;
2680cbbd2e1SMark F. Adams   PetscScalar    *cpcol_1_state,*cpcol_2_state,*cpcol_2_par_orig,*lid_parent_gid;
2690cbbd2e1SMark F. Adams   PetscInt       *lid_cprowID_1;
270c8b0795cSMark F. Adams   NState         *lid_state;
2710cbbd2e1SMark F. Adams   Vec            ghost_par_orig2;
272c8b0795cSMark F. Adams 
273c8b0795cSMark F. Adams   PetscFunctionBegin;
2743b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Gmat_2,&comm);CHKERRQ(ierr);
275c8b0795cSMark F. Adams   ierr = MatGetOwnershipRange(Gmat_1,&my0,&Iend);CHKERRQ(ierr);
276c8b0795cSMark F. Adams 
277c8b0795cSMark F. Adams   /* get submatrices */
278251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)Gmat_1, MATMPIAIJ, &isMPI);CHKERRQ(ierr);
279c8b0795cSMark F. Adams   if (isMPI) {
280c8b0795cSMark F. Adams     /* grab matrix objects */
281c8b0795cSMark F. Adams     mpimat_2 = (Mat_MPIAIJ*)Gmat_2->data;
282c8b0795cSMark F. Adams     mpimat_1 = (Mat_MPIAIJ*)Gmat_1->data;
283c8b0795cSMark F. Adams     matA_1   = (Mat_SeqAIJ*)mpimat_1->A->data;
284c8b0795cSMark F. Adams     matB_1   = (Mat_SeqAIJ*)mpimat_1->B->data;
285c8b0795cSMark F. Adams 
286c8b0795cSMark F. Adams     /* force compressed row storage for B matrix in AuxMat */
28711e456e1SBarry Smith     ierr = MatCheckCompressedRow(mpimat_1->B,matB_1->nonzerorowcnt,&matB_1->compressedrow,matB_1->i,Gmat_1->rmap->n,-1.0);CHKERRQ(ierr);
288c8b0795cSMark F. Adams 
289785e854fSJed Brown     ierr = PetscMalloc1(nloc, &lid_cprowID_1);CHKERRQ(ierr);
2900cbbd2e1SMark F. Adams     for (lid = 0; lid < nloc; lid++) lid_cprowID_1[lid] = -1;
291c8b0795cSMark F. Adams     for (ix=0; ix<matB_1->compressedrow.nrows; ix++) {
292c8b0795cSMark F. Adams       PetscInt lid = matB_1->compressedrow.rindex[ix];
293c8b0795cSMark F. Adams       lid_cprowID_1[lid] = ix;
294c8b0795cSMark F. Adams     }
295806fa848SBarry Smith   } else {
29615687449SMark Adams     PetscBool        isAIJ;
2974099cc6bSBarry Smith     ierr = PetscObjectBaseTypeCompare((PetscObject)Gmat_1,MATSEQAIJ,&isAIJ);CHKERRQ(ierr);
29815687449SMark Adams     if (!isAIJ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Require AIJ matrix.");
299c8b0795cSMark F. Adams     matA_1        = (Mat_SeqAIJ*)Gmat_1->data;
3000298fd71SBarry Smith     lid_cprowID_1 = NULL;
301c8b0795cSMark F. Adams   }
30278a438d6SMark Adams   if (nloc>0) {
303359038b3SMark Adams     if (matB_1 && !matB_1->compressedrow.use) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"matB_1 && !matB_1->compressedrow.use: PETSc bug???");
30478a438d6SMark Adams   }
305c8b0795cSMark F. Adams   /* get state of locals and selected gid for deleted */
306e632b94dSBarry Smith   ierr = PetscMalloc2(nloc, &lid_state,nloc, &lid_parent_gid);CHKERRQ(ierr);
307c8b0795cSMark F. Adams   for (lid = 0; lid < nloc; lid++) {
3080cbbd2e1SMark F. Adams     lid_parent_gid[lid] = -1.0;
309c8b0795cSMark F. Adams     lid_state[lid]      = DELETED;
310c8b0795cSMark F. Adams   }
3110cbbd2e1SMark F. Adams 
3120cbbd2e1SMark F. Adams   /* set lid_state */
3130cbbd2e1SMark F. Adams   for (lid = 0; lid < nloc; lid++) {
314539c167fSBarry Smith     PetscCDIntNd *pos;
315e78576d6SMark F. Adams     ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr);
316e78576d6SMark F. Adams     if (pos) {
317e78576d6SMark F. Adams       PetscInt gid1;
3182fa5cd67SKarl Rupp 
319484f0a72SBarry Smith       ierr = PetscCDIntNdGetID(pos, &gid1);CHKERRQ(ierr);
32071959b99SBarry Smith       if (gid1 != lid+my0) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gid1 %D != lid %D + my0 %D",gid1,lid,my0);
3210cbbd2e1SMark F. Adams       lid_state[lid] = gid1;
322b43b03e9SMark F. Adams     }
323b43b03e9SMark F. Adams   }
3240cbbd2e1SMark F. Adams 
3250cbbd2e1SMark F. Adams   /* map local to selected local, DELETED means a ghost owns it */
326c8b0795cSMark F. Adams   for (lid=kk=0; lid<nloc; lid++) {
327c8b0795cSMark F. Adams     NState state = lid_state[lid];
328c8b0795cSMark F. Adams     if (IS_SELECTED(state)) {
329539c167fSBarry Smith       PetscCDIntNd *pos;
330e78576d6SMark F. Adams       ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr);
331e78576d6SMark F. Adams       while (pos) {
332e78576d6SMark F. Adams         PetscInt gid1;
333484f0a72SBarry Smith         ierr = PetscCDIntNdGetID(pos, &gid1);CHKERRQ(ierr);
334e78576d6SMark F. Adams         ierr = PetscCDGetNextPos(aggs_2,lid,&pos);CHKERRQ(ierr);
335e78576d6SMark F. Adams 
3362fa5cd67SKarl Rupp         if (gid1 >= my0 && gid1 < Iend) lid_parent_gid[gid1-my0] = (PetscScalar)(lid + my0);
337c8b0795cSMark F. Adams       }
3380cbbd2e1SMark F. Adams     }
3390cbbd2e1SMark F. Adams   }
3400cbbd2e1SMark F. Adams   /* get 'cpcol_1/2_state' & cpcol_2_par_orig - uses mpimat_1/2->lvec for temp space */
341c8b0795cSMark F. Adams   if (isMPI) {
342c8b0795cSMark F. Adams     Vec tempVec;
343c8b0795cSMark F. Adams     /* get 'cpcol_1_state' */
3442a7a6963SBarry Smith     ierr = MatCreateVecs(Gmat_1, &tempVec, 0);CHKERRQ(ierr);
345c8b0795cSMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
346c8b0795cSMark F. Adams       PetscScalar v = (PetscScalar)lid_state[kk];
347c8b0795cSMark F. Adams       ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr);
348c8b0795cSMark F. Adams     }
349c8b0795cSMark F. Adams     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
350c8b0795cSMark F. Adams     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
351806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
352806fa848SBarry Smith     ierr = VecScatterEnd(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
353c8b0795cSMark F. Adams     ierr = VecGetArray(mpimat_1->lvec, &cpcol_1_state);CHKERRQ(ierr);
354c8b0795cSMark F. Adams     /* get 'cpcol_2_state' */
355806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
356806fa848SBarry Smith     ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
357c8b0795cSMark F. Adams     ierr = VecGetArray(mpimat_2->lvec, &cpcol_2_state);CHKERRQ(ierr);
3580cbbd2e1SMark F. Adams     /* get 'cpcol_2_par_orig' */
3590cbbd2e1SMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
3600cbbd2e1SMark F. Adams       PetscScalar v = (PetscScalar)lid_parent_gid[kk];
3610cbbd2e1SMark F. Adams       ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr);
3620cbbd2e1SMark F. Adams     }
3630cbbd2e1SMark F. Adams     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
3640cbbd2e1SMark F. Adams     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
3650cbbd2e1SMark F. Adams     ierr = VecDuplicate(mpimat_2->lvec, &ghost_par_orig2);CHKERRQ(ierr);
366806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
367806fa848SBarry Smith     ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
3680cbbd2e1SMark F. Adams     ierr = VecGetArray(ghost_par_orig2, &cpcol_2_par_orig);CHKERRQ(ierr);
3690cbbd2e1SMark F. Adams 
370c8b0795cSMark F. Adams     ierr = VecDestroy(&tempVec);CHKERRQ(ierr);
371c8b0795cSMark F. Adams   } /* ismpi */
372c8b0795cSMark F. Adams 
373c8b0795cSMark F. Adams   /* doit */
374c8b0795cSMark F. Adams   for (lid=0; lid<nloc; lid++) {
375c8b0795cSMark F. Adams     NState state = lid_state[lid];
3760cbbd2e1SMark F. Adams     if (IS_SELECTED(state)) {
3770cbbd2e1SMark F. Adams       /* steal locals */
378c8b0795cSMark F. Adams       ii  = matA_1->i; n = ii[lid+1] - ii[lid];
379c8b0795cSMark F. Adams       idx = matA_1->j + ii[lid];
380c8b0795cSMark F. Adams       for (j=0; j<n; j++) {
3810cbbd2e1SMark F. Adams         PetscInt lidj   = idx[j], sgid;
382c8b0795cSMark F. Adams         NState   statej = lid_state[lidj];
3830cbbd2e1SMark F. Adams         if (statej==DELETED && (sgid=(PetscInt)PetscRealPart(lid_parent_gid[lidj])) != lid+my0) { /* steal local */
3840cbbd2e1SMark F. Adams           lid_parent_gid[lidj] = (PetscScalar)(lid+my0); /* send this if sgid is not local */
3850cbbd2e1SMark F. Adams           if (sgid >= my0 && sgid < Iend) {       /* I'm stealing this local from a local sgid */
3860cbbd2e1SMark F. Adams             PetscInt     hav=0,slid=sgid-my0,gidj=lidj+my0;
387539c167fSBarry Smith             PetscCDIntNd *pos,*last=NULL;
388c8b0795cSMark F. Adams             /* looking for local from local so id_llist_2 works */
389e78576d6SMark F. Adams             ierr = PetscCDGetHeadPos(aggs_2,slid,&pos);CHKERRQ(ierr);
390e78576d6SMark F. Adams             while (pos) {
391e78576d6SMark F. Adams               PetscInt gid;
392484f0a72SBarry Smith               ierr = PetscCDIntNdGetID(pos, &gid);CHKERRQ(ierr);
3930cbbd2e1SMark F. Adams               if (gid == gidj) {
39471959b99SBarry Smith                 if (!last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"last cannot be null");
39541b27cdeSMark F. Adams                 ierr = PetscCDRemoveNextNode(aggs_2, slid, last);CHKERRQ(ierr);
39641b27cdeSMark F. Adams                 ierr = PetscCDAppendNode(aggs_2, lid, pos);CHKERRQ(ierr);
3970cbbd2e1SMark F. Adams                 hav  = 1;
398c8b0795cSMark F. Adams                 break;
399806fa848SBarry Smith               } else last = pos;
400e78576d6SMark F. Adams 
401e78576d6SMark F. Adams               ierr = PetscCDGetNextPos(aggs_2,slid,&pos);CHKERRQ(ierr);
402c8b0795cSMark F. Adams             }
403c8b0795cSMark F. Adams             if (hav!=1) {
40471959b99SBarry Smith               if (!hav) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix");
405*b817416eSBarry Smith               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %D times???",hav);
406c8b0795cSMark F. Adams             }
407806fa848SBarry Smith           } else {            /* I'm stealing this local, owned by a ghost */
408abf98081SSatish Balay             if (sgid != -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Have un-symmetric graph (apparently). Use '-%spc_gamg_sym_graph true' to symetrize the graph or '-%spc_gamg_threshold -1' if the matrix is structurally symmetric.",((PetscObject)pc)->prefix,((PetscObject)pc)->prefix);
40941b27cdeSMark F. Adams             ierr = PetscCDAppendID(aggs_2, lid, lidj+my0);CHKERRQ(ierr);
410c8b0795cSMark F. Adams           }
411c8b0795cSMark F. Adams         }
4120cbbd2e1SMark F. Adams       } /* local neighbors */
413806fa848SBarry Smith     } else if (state == DELETED && lid_cprowID_1) {
4140cbbd2e1SMark F. Adams       PetscInt sgidold = (PetscInt)PetscRealPart(lid_parent_gid[lid]);
415c8b0795cSMark F. Adams       /* see if I have a selected ghost neighbor that will steal me */
416c8b0795cSMark F. Adams       if ((ix=lid_cprowID_1[lid]) != -1) {
417c8b0795cSMark F. Adams         ii  = matB_1->compressedrow.i; n = ii[ix+1] - ii[ix];
418c8b0795cSMark F. Adams         idx = matB_1->j + ii[ix];
419c8b0795cSMark F. Adams         for (j=0; j<n; j++) {
420c8b0795cSMark F. Adams           PetscInt cpid   = idx[j];
421c8b0795cSMark F. Adams           NState   statej = (NState)PetscRealPart(cpcol_1_state[cpid]);
422c8b0795cSMark F. Adams           if (IS_SELECTED(statej) && sgidold != (PetscInt)statej) { /* ghost will steal this, remove from my list */
4230cbbd2e1SMark F. Adams             lid_parent_gid[lid] = (PetscScalar)statej; /* send who selected */
4240cbbd2e1SMark F. Adams             if (sgidold>=my0 && sgidold<Iend) { /* this was mine */
4250cbbd2e1SMark F. Adams               PetscInt     hav=0,oldslidj=sgidold-my0;
426539c167fSBarry Smith               PetscCDIntNd *pos,*last=NULL;
4270cbbd2e1SMark F. Adams               /* remove from 'oldslidj' list */
428e78576d6SMark F. Adams               ierr = PetscCDGetHeadPos(aggs_2,oldslidj,&pos);CHKERRQ(ierr);
429e78576d6SMark F. Adams               while (pos) {
430e78576d6SMark F. Adams                 PetscInt gid;
431484f0a72SBarry Smith                 ierr = PetscCDIntNdGetID(pos, &gid);CHKERRQ(ierr);
4320cbbd2e1SMark F. Adams                 if (lid+my0 == gid) {
4330cbbd2e1SMark F. Adams                   /* id_llist_2[lastid] = id_llist_2[flid];   /\* remove lid from oldslidj list *\/ */
43471959b99SBarry Smith                   if (!last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"last cannot be null");
43541b27cdeSMark F. Adams                   ierr = PetscCDRemoveNextNode(aggs_2, oldslidj, last);CHKERRQ(ierr);
4360cbbd2e1SMark F. Adams                   /* ghost (PetscScalar)statej will add this later */
4370cbbd2e1SMark F. Adams                   hav = 1;
438c8b0795cSMark F. Adams                   break;
439806fa848SBarry Smith                 } else last = pos;
440e78576d6SMark F. Adams 
441e78576d6SMark F. Adams                 ierr = PetscCDGetNextPos(aggs_2,oldslidj,&pos);CHKERRQ(ierr);
442c8b0795cSMark F. Adams               }
443c8b0795cSMark F. Adams               if (hav!=1) {
4447f66b68fSMark Adams                 if (!hav) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix");
445*b817416eSBarry Smith                 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %D times???",hav);
446c8b0795cSMark F. Adams               }
447806fa848SBarry Smith             } else {
4480cbbd2e1SMark F. Adams               /* ghosts remove this later */
4490cbbd2e1SMark F. Adams             }
450c8b0795cSMark F. Adams           }
451c8b0795cSMark F. Adams         }
452c8b0795cSMark F. Adams       }
453c8b0795cSMark F. Adams     } /* selected/deleted */
454c8b0795cSMark F. Adams   } /* node loop */
455c8b0795cSMark F. Adams 
456c8b0795cSMark F. Adams   if (isMPI) {
4570cbbd2e1SMark F. Adams     PetscScalar     *cpcol_2_parent,*cpcol_2_gid;
4580cbbd2e1SMark F. Adams     Vec             tempVec,ghostgids2,ghostparents2;
4590cbbd2e1SMark F. Adams     PetscInt        cpid,nghost_2;
4601943db53SBarry Smith     PCGAMGHashTable gid_cpid;
461c8b0795cSMark F. Adams 
4620cbbd2e1SMark F. Adams     ierr = VecGetSize(mpimat_2->lvec, &nghost_2);CHKERRQ(ierr);
4632a7a6963SBarry Smith     ierr = MatCreateVecs(Gmat_2, &tempVec, 0);CHKERRQ(ierr);
4640cbbd2e1SMark F. Adams 
4650cbbd2e1SMark F. Adams     /* get 'cpcol_2_parent' */
466c8b0795cSMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
4670cbbd2e1SMark F. Adams       ierr = VecSetValues(tempVec, 1, &j, &lid_parent_gid[kk], INSERT_VALUES);CHKERRQ(ierr);
468c8b0795cSMark F. Adams     }
469c8b0795cSMark F. Adams     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
470c8b0795cSMark F. Adams     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
4710cbbd2e1SMark F. Adams     ierr = VecDuplicate(mpimat_2->lvec, &ghostparents2);CHKERRQ(ierr);
472806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
473806fa848SBarry Smith     ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4740cbbd2e1SMark F. Adams     ierr = VecGetArray(ghostparents2, &cpcol_2_parent);CHKERRQ(ierr);
4750cbbd2e1SMark F. Adams 
4760cbbd2e1SMark F. Adams     /* get 'cpcol_2_gid' */
4770cbbd2e1SMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
4780cbbd2e1SMark F. Adams       PetscScalar v = (PetscScalar)j;
4790cbbd2e1SMark F. Adams       ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr);
4800cbbd2e1SMark F. Adams     }
4810cbbd2e1SMark F. Adams     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
4820cbbd2e1SMark F. Adams     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
4830cbbd2e1SMark F. Adams     ierr = VecDuplicate(mpimat_2->lvec, &ghostgids2);CHKERRQ(ierr);
484806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
485806fa848SBarry Smith     ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4860cbbd2e1SMark F. Adams     ierr = VecGetArray(ghostgids2, &cpcol_2_gid);CHKERRQ(ierr);
487c8b0795cSMark F. Adams     ierr = VecDestroy(&tempVec);CHKERRQ(ierr);
488c8b0795cSMark F. Adams 
4890cbbd2e1SMark F. Adams     /* look for deleted ghosts and add to table */
4901943db53SBarry Smith     ierr = PCGAMGHashTableCreate(2*nghost_2+1, &gid_cpid);CHKERRQ(ierr);
4910cbbd2e1SMark F. Adams     for (cpid = 0; cpid < nghost_2; cpid++) {
4920cbbd2e1SMark F. Adams       NState state = (NState)PetscRealPart(cpcol_2_state[cpid]);
4930cbbd2e1SMark F. Adams       if (state==DELETED) {
4940cbbd2e1SMark F. Adams         PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]);
4950cbbd2e1SMark F. Adams         PetscInt sgid_old = (PetscInt)PetscRealPart(cpcol_2_par_orig[cpid]);
4960cbbd2e1SMark F. Adams         if (sgid_old == -1 && sgid_new != -1) {
4970cbbd2e1SMark F. Adams           PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]);
4981943db53SBarry Smith           ierr = PCGAMGHashTableAdd(&gid_cpid, gid, cpid);CHKERRQ(ierr);
4990cbbd2e1SMark F. Adams         }
5000cbbd2e1SMark F. Adams       }
5010cbbd2e1SMark F. Adams     }
502c8b0795cSMark F. Adams 
5030cbbd2e1SMark F. Adams     /* look for deleted ghosts and see if they moved - remove it */
504c8b0795cSMark F. Adams     for (lid=0; lid<nloc; lid++) {
505c8b0795cSMark F. Adams       NState state = lid_state[lid];
506c8b0795cSMark F. Adams       if (IS_SELECTED(state)) {
507539c167fSBarry Smith         PetscCDIntNd *pos,*last=NULL;
508c8b0795cSMark F. Adams         /* look for deleted ghosts and see if they moved */
509e78576d6SMark F. Adams         ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr);
510e78576d6SMark F. Adams         while (pos) {
511e78576d6SMark F. Adams           PetscInt gid;
512484f0a72SBarry Smith           ierr = PetscCDIntNdGetID(pos, &gid);CHKERRQ(ierr);
513e78576d6SMark F. Adams 
5140cbbd2e1SMark F. Adams           if (gid < my0 || gid >= Iend) {
5151943db53SBarry Smith             ierr = PCGAMGHashTableFind(&gid_cpid, gid, &cpid);CHKERRQ(ierr);
5160cbbd2e1SMark F. Adams             if (cpid != -1) {
5170cbbd2e1SMark F. Adams               /* a moved ghost - */
5180cbbd2e1SMark F. Adams               /* id_llist_2[lastid] = id_llist_2[flid];    /\* remove 'flid' from list *\/ */
51941b27cdeSMark F. Adams               ierr = PetscCDRemoveNextNode(aggs_2, lid, last);CHKERRQ(ierr);
520806fa848SBarry Smith             } else last = pos;
521806fa848SBarry Smith           } else last = pos;
522e78576d6SMark F. Adams 
523e78576d6SMark F. Adams           ierr = PetscCDGetNextPos(aggs_2,lid,&pos);CHKERRQ(ierr);
524c8b0795cSMark F. Adams         } /* loop over list of deleted */
525c8b0795cSMark F. Adams       } /* selected */
526c8b0795cSMark F. Adams     }
5271943db53SBarry Smith     ierr = PCGAMGHashTableDestroy(&gid_cpid);CHKERRQ(ierr);
528c8b0795cSMark F. Adams 
5290cbbd2e1SMark F. Adams     /* look at ghosts, see if they changed - and it */
5300cbbd2e1SMark F. Adams     for (cpid = 0; cpid < nghost_2; cpid++) {
5310cbbd2e1SMark F. Adams       PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]);
5320cbbd2e1SMark F. Adams       if (sgid_new >= my0 && sgid_new < Iend) { /* this is mine */
5330cbbd2e1SMark F. Adams         PetscInt     gid     = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]);
5340cbbd2e1SMark F. Adams         PetscInt     slid_new=sgid_new-my0,hav=0;
535539c167fSBarry Smith         PetscCDIntNd *pos;
536539c167fSBarry Smith 
5370cbbd2e1SMark F. Adams         /* search for this gid to see if I have it */
538e78576d6SMark F. Adams         ierr = PetscCDGetHeadPos(aggs_2,slid_new,&pos);CHKERRQ(ierr);
539e78576d6SMark F. Adams         while (pos) {
540e78576d6SMark F. Adams           PetscInt gidj;
541484f0a72SBarry Smith           ierr = PetscCDIntNdGetID(pos, &gidj);CHKERRQ(ierr);
542e78576d6SMark F. Adams           ierr = PetscCDGetNextPos(aggs_2,slid_new,&pos);CHKERRQ(ierr);
543e78576d6SMark F. Adams 
5440cbbd2e1SMark F. Adams           if (gidj == gid) { hav = 1; break; }
545c8b0795cSMark F. Adams         }
546c8b0795cSMark F. Adams         if (hav != 1) {
547ffc955d6SMark F. Adams           /* insert 'flidj' into head of llist */
54841b27cdeSMark F. Adams           ierr = PetscCDAppendID(aggs_2, slid_new, gid);CHKERRQ(ierr);
549c8b0795cSMark F. Adams         }
550c8b0795cSMark F. Adams       }
551c8b0795cSMark F. Adams     }
552c8b0795cSMark F. Adams 
5530cbbd2e1SMark F. Adams     ierr = VecRestoreArray(mpimat_1->lvec, &cpcol_1_state);CHKERRQ(ierr);
5540cbbd2e1SMark F. Adams     ierr = VecRestoreArray(mpimat_2->lvec, &cpcol_2_state);CHKERRQ(ierr);
5550cbbd2e1SMark F. Adams     ierr = VecRestoreArray(ghostparents2, &cpcol_2_parent);CHKERRQ(ierr);
5560cbbd2e1SMark F. Adams     ierr = VecRestoreArray(ghostgids2, &cpcol_2_gid);CHKERRQ(ierr);
557c8b0795cSMark F. Adams     ierr = PetscFree(lid_cprowID_1);CHKERRQ(ierr);
5580cbbd2e1SMark F. Adams     ierr = VecDestroy(&ghostgids2);CHKERRQ(ierr);
5590cbbd2e1SMark F. Adams     ierr = VecDestroy(&ghostparents2);CHKERRQ(ierr);
5600cbbd2e1SMark F. Adams     ierr = VecDestroy(&ghost_par_orig2);CHKERRQ(ierr);
561c8b0795cSMark F. Adams   }
562c8b0795cSMark F. Adams 
563e632b94dSBarry Smith   ierr = PetscFree2(lid_state,lid_parent_gid);CHKERRQ(ierr);
564c8b0795cSMark F. Adams   PetscFunctionReturn(0);
565c8b0795cSMark F. Adams }
5662e68589bSMark F. Adams 
5672e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
5682e68589bSMark F. Adams /*
569a2f3521dSMark F. Adams    PCSetData_AGG - called if data is not set with PCSetCoordinates.
570a2f3521dSMark F. Adams       Looks in Mat for near null space.
571a2f3521dSMark F. Adams       Does not work for Stokes
5722e68589bSMark F. Adams 
5732e68589bSMark F. Adams   Input Parameter:
5742e68589bSMark F. Adams    . pc -
575a2f3521dSMark F. Adams    . a_A - matrix to get (near) null space out of.
5762e68589bSMark F. Adams */
577e0877f53SBarry Smith static PetscErrorCode PCSetData_AGG(PC pc, Mat a_A)
5782e68589bSMark F. Adams {
5792e68589bSMark F. Adams   PetscErrorCode ierr;
580b8cd405aSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
581b8cd405aSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
582b8cd405aSMark F. Adams   MatNullSpace   mnull;
58366f2ef4bSMark Adams 
584*b817416eSBarry Smith   PetscFunctionBegin;
585b8cd405aSMark F. Adams   ierr = MatGetNearNullSpace(a_A, &mnull);CHKERRQ(ierr);
586b8cd405aSMark F. Adams   if (!mnull) {
58766f2ef4bSMark Adams     DM dm;
58866f2ef4bSMark Adams     ierr = PCGetDM(pc, &dm);CHKERRQ(ierr);
58966f2ef4bSMark Adams     if (!dm) {
59066f2ef4bSMark Adams       ierr = MatGetDM(a_A, &dm);CHKERRQ(ierr);
59166f2ef4bSMark Adams     }
59266f2ef4bSMark Adams     if (dm) {
59366f2ef4bSMark Adams       PetscObject deformation;
594b0d30dd6SMatthew G. Knepley       PetscInt    Nf;
595b0d30dd6SMatthew G. Knepley 
596b0d30dd6SMatthew G. Knepley       ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr);
597b0d30dd6SMatthew G. Knepley       if (Nf) {
59866f2ef4bSMark Adams         ierr = DMGetField(dm, 0, &deformation);CHKERRQ(ierr);
59966f2ef4bSMark Adams         ierr = PetscObjectQuery((PetscObject)deformation,"nearnullspace",(PetscObject*)&mnull);CHKERRQ(ierr);
60066f2ef4bSMark Adams         if (!mnull) {
60166f2ef4bSMark Adams           ierr = PetscObjectQuery((PetscObject)deformation,"nullspace",(PetscObject*)&mnull);CHKERRQ(ierr);
60266f2ef4bSMark Adams         }
60366f2ef4bSMark Adams       }
60466f2ef4bSMark Adams     }
605b0d30dd6SMatthew G. Knepley   }
60666f2ef4bSMark Adams 
60766f2ef4bSMark Adams   if (!mnull) {
608a2f3521dSMark F. Adams     PetscInt bs,NN,MM;
6099d1b15c3SMark F. Adams     ierr = MatGetBlockSize(a_A, &bs);CHKERRQ(ierr);
61071959b99SBarry Smith     ierr = MatGetLocalSize(a_A, &MM, &NN);CHKERRQ(ierr);
61171959b99SBarry Smith     if (MM % bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MM %D must be divisible by bs %D",MM,bs);
6120298fd71SBarry Smith     ierr = PCSetCoordinates_AGG(pc, bs, MM/bs, NULL);CHKERRQ(ierr);
613806fa848SBarry Smith   } else {
614b8cd405aSMark F. Adams     PetscReal *nullvec;
615b8cd405aSMark F. Adams     PetscBool has_const;
616b8cd405aSMark F. Adams     PetscInt  i,j,mlocal,nvec,bs;
617b8cd405aSMark F. Adams     const Vec *vecs; const PetscScalar *v;
618*b817416eSBarry Smith 
6190298fd71SBarry Smith     ierr = MatGetLocalSize(a_A,&mlocal,NULL);CHKERRQ(ierr);
620b8cd405aSMark F. Adams     ierr = MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs);CHKERRQ(ierr);
621a0dea87dSMark Adams     pc_gamg->data_sz = (nvec+!!has_const)*mlocal;
622785e854fSJed Brown     ierr = PetscMalloc1((nvec+!!has_const)*mlocal,&nullvec);CHKERRQ(ierr);
623b8cd405aSMark F. Adams     if (has_const) for (i=0; i<mlocal; i++) nullvec[i] = 1.0;
624b8cd405aSMark F. Adams     for (i=0; i<nvec; i++) {
625b8cd405aSMark F. Adams       ierr = VecGetArrayRead(vecs[i],&v);CHKERRQ(ierr);
626b8cd405aSMark F. Adams       for (j=0; j<mlocal; j++) nullvec[(i+!!has_const)*mlocal + j] = PetscRealPart(v[j]);
627b8cd405aSMark F. Adams       ierr = VecRestoreArrayRead(vecs[i],&v);CHKERRQ(ierr);
628b8cd405aSMark F. Adams     }
629b8cd405aSMark F. Adams     pc_gamg->data           = nullvec;
630b8cd405aSMark F. Adams     pc_gamg->data_cell_cols = (nvec+!!has_const);
63106e133e7SMark Adams     ierr = MatGetBlockSize(a_A, &bs);CHKERRQ(ierr);
632b8cd405aSMark F. Adams     pc_gamg->data_cell_rows = bs;
633b8cd405aSMark F. Adams   }
6342e68589bSMark F. Adams   PetscFunctionReturn(0);
6352e68589bSMark F. Adams }
6362e68589bSMark F. Adams 
6372e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
6382e68589bSMark F. Adams /*
6392e68589bSMark F. Adams  formProl0
6402e68589bSMark F. Adams 
6412e68589bSMark F. Adams    Input Parameter:
6422adfe9d3SBarry Smith    . agg_llists - list of arrays with aggregates -- list from selected vertices of aggregate unselected vertices
6432adfe9d3SBarry Smith    . bs - row block size
6440cbbd2e1SMark F. Adams    . nSAvec - column bs of new P
6450cbbd2e1SMark F. Adams    . my0crs - global index of start of locals
6462adfe9d3SBarry Smith    . data_stride - bs*(nloc nodes + ghost nodes) [data_stride][nSAvec]
6472e68589bSMark F. Adams    . data_in[data_stride*nSAvec] - local data on fine grid
6482e68589bSMark F. Adams    . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist'
6492e68589bSMark F. Adams   Output Parameter:
6502e68589bSMark F. Adams    . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data
6512e68589bSMark F. Adams    . a_Prol - prolongation operator
6522e68589bSMark F. Adams */
6532adfe9d3SBarry 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)
6542e68589bSMark F. Adams {
6552e68589bSMark F. Adams   PetscErrorCode  ierr;
656ac187d40SBarry Smith   PetscInt        Istart,my0,Iend,nloc,clid,flid = 0,aggID,kk,jj,ii,mm,ndone,nSelected,minsz,nghosts,out_data_stride;
6573b4367a7SBarry Smith   MPI_Comm        comm;
65873911c69SBarry Smith   PetscMPIInt     rank;
6592e68589bSMark F. Adams   PetscReal       *out_data;
660539c167fSBarry Smith   PetscCDIntNd    *pos;
6611943db53SBarry Smith   PCGAMGHashTable fgid_flid;
6620cbbd2e1SMark F. Adams 
6632e68589bSMark F. Adams   PetscFunctionBegin;
6643b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)a_Prol,&comm);CHKERRQ(ierr);
6653b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
6662e68589bSMark F. Adams   ierr = MatGetOwnershipRange(a_Prol, &Istart, &Iend);CHKERRQ(ierr);
66771959b99SBarry Smith   nloc = (Iend-Istart)/bs; my0 = Istart/bs;
668*b817416eSBarry Smith   if ((Iend-Istart) % bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Iend %D - Istart %D must be divisible by bs %D",Iend,Istart,bs);
6690cbbd2e1SMark F. Adams   Iend   /= bs;
6700cbbd2e1SMark F. Adams   nghosts = data_stride/bs - nloc;
6712e68589bSMark F. Adams 
6721943db53SBarry Smith   ierr = PCGAMGHashTableCreate(2*nghosts+1, &fgid_flid);CHKERRQ(ierr);
6730cbbd2e1SMark F. Adams   for (kk=0; kk<nghosts; kk++) {
6741943db53SBarry Smith     ierr = PCGAMGHashTableAdd(&fgid_flid, flid_fgid[nloc+kk], nloc+kk);CHKERRQ(ierr);
6752e68589bSMark F. Adams   }
6762e68589bSMark F. Adams 
6770cbbd2e1SMark F. Adams   /* count selected -- same as number of cols of P */
6780cbbd2e1SMark F. Adams   for (nSelected=mm=0; mm<nloc; mm++) {
679e78576d6SMark F. Adams     PetscBool ise;
680e78576d6SMark F. Adams     ierr = PetscCDEmptyAt(agg_llists, mm, &ise);CHKERRQ(ierr);
681e78576d6SMark F. Adams     if (!ise) nSelected++;
6820cbbd2e1SMark F. Adams   }
6830cbbd2e1SMark F. Adams   ierr = MatGetOwnershipRangeColumn(a_Prol, &ii, &jj);CHKERRQ(ierr);
68471959b99SBarry Smith   if ((ii/nSAvec) != my0crs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D /nSAvec %D  != my0crs %D",ii,nSAvec,my0crs);
68571959b99SBarry Smith   if (nSelected != (jj-ii)/nSAvec) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_PLIB,"nSelected %D != (jj %D - ii %D)/nSAvec %D",nSelected,jj,ii,nSAvec);
6860cbbd2e1SMark F. Adams 
6872e68589bSMark F. Adams   /* aloc space for coarse point data (output) */
6880cbbd2e1SMark F. Adams   out_data_stride = nSelected*nSAvec;
6892fa5cd67SKarl Rupp 
690785e854fSJed Brown   ierr = PetscMalloc1(out_data_stride*nSAvec, &out_data);CHKERRQ(ierr);
69133812677SSatish Balay   for (ii=0;ii<out_data_stride*nSAvec;ii++) out_data[ii]=PETSC_MAX_REAL;
6920cbbd2e1SMark F. Adams   *a_data_out = out_data; /* output - stride nSelected*nSAvec */
6932e68589bSMark F. Adams 
6942e68589bSMark F. Adams   /* find points and set prolongation */
695c8b0795cSMark F. Adams   minsz = 100;
6962e68589bSMark F. Adams   ndone = 0;
6970cbbd2e1SMark F. Adams   for (mm = clid = 0; mm < nloc; mm++) {
698e78576d6SMark F. Adams     ierr = PetscCDSizeAt(agg_llists, mm, &jj);CHKERRQ(ierr);
699e78576d6SMark F. Adams     if (jj > 0) {
7000cbbd2e1SMark F. Adams       const PetscInt lid = mm, cgid = my0crs + clid;
7010cbbd2e1SMark F. Adams       PetscInt       cids[100]; /* max bs */
7020cbbd2e1SMark F. Adams       PetscBLASInt   asz  =jj,M=asz*bs,N=nSAvec,INFO;
7032e68589bSMark F. Adams       PetscBLASInt   Mdata=M+((N-M>0) ? N-M : 0),LDA=Mdata,LWORK=N*bs;
7042e68589bSMark F. Adams       PetscScalar    *qqc,*qqr,*TAU,*WORK;
7052e68589bSMark F. Adams       PetscInt       *fids;
70665d7b583SSatish Balay       PetscReal      *data;
707*b817416eSBarry Smith 
7080cbbd2e1SMark F. Adams       /* count agg */
7090cbbd2e1SMark F. Adams       if (asz<minsz) minsz = asz;
7100cbbd2e1SMark F. Adams 
7110cbbd2e1SMark F. Adams       /* get block */
712e632b94dSBarry Smith       ierr = PetscMalloc5(Mdata*N, &qqc,M*N, &qqr,N, &TAU,LWORK, &WORK,M, &fids);CHKERRQ(ierr);
7132e68589bSMark F. Adams 
7142e68589bSMark F. Adams       aggID = 0;
715e78576d6SMark F. Adams       ierr  = PetscCDGetHeadPos(agg_llists,lid,&pos);CHKERRQ(ierr);
716e78576d6SMark F. Adams       while (pos) {
717e78576d6SMark F. Adams         PetscInt gid1;
718484f0a72SBarry Smith         ierr = PetscCDIntNdGetID(pos, &gid1);CHKERRQ(ierr);
719e78576d6SMark F. Adams         ierr = PetscCDGetNextPos(agg_llists,lid,&pos);CHKERRQ(ierr);
720e78576d6SMark F. Adams 
7210cbbd2e1SMark F. Adams         if (gid1 >= my0 && gid1 < Iend) flid = gid1 - my0;
7220cbbd2e1SMark F. Adams         else {
7231943db53SBarry Smith           ierr = PCGAMGHashTableFind(&fgid_flid, gid1, &flid);CHKERRQ(ierr);
72471959b99SBarry Smith           if (flid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot find gid1 in table");
7250cbbd2e1SMark F. Adams         }
7262e68589bSMark F. Adams         /* copy in B_i matrix - column oriented */
72765d7b583SSatish Balay         data = &data_in[flid*bs];
7283d3eaba7SBarry Smith         for (ii = 0; ii < bs; ii++) {
7292e68589bSMark F. Adams           for (jj = 0; jj < N; jj++) {
7300cbbd2e1SMark F. Adams             PetscReal d = data[jj*data_stride + ii];
7310cbbd2e1SMark F. Adams             qqc[jj*Mdata + aggID*bs + ii] = d;
7322e68589bSMark F. Adams           }
7332e68589bSMark F. Adams         }
7342e68589bSMark F. Adams         /* set fine IDs */
7352e68589bSMark F. Adams         for (kk=0; kk<bs; kk++) fids[aggID*bs + kk] = flid_fgid[flid]*bs + kk;
7362e68589bSMark F. Adams         aggID++;
7370cbbd2e1SMark F. Adams       }
7382e68589bSMark F. Adams 
7392e68589bSMark F. Adams       /* pad with zeros */
7402e68589bSMark F. Adams       for (ii = asz*bs; ii < Mdata; ii++) {
7412e68589bSMark F. Adams         for (jj = 0; jj < N; jj++, kk++) {
7422e68589bSMark F. Adams           qqc[jj*Mdata + ii] = .0;
7432e68589bSMark F. Adams         }
7442e68589bSMark F. Adams       }
7452e68589bSMark F. Adams 
7462e68589bSMark F. Adams       ndone += aggID;
7472e68589bSMark F. Adams       /* QR */
74884df3f90SSatish Balay       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
7498b83055fSJed Brown       PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
75084df3f90SSatish Balay       ierr = PetscFPTrapPop();CHKERRQ(ierr);
751d23427c4SJed Brown       if (INFO != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xGEQRF error");
7522e68589bSMark F. Adams       /* get R - column oriented - output B_{i+1} */
7532e68589bSMark F. Adams       {
7542e68589bSMark F. Adams         PetscReal *data = &out_data[clid*nSAvec];
7552e68589bSMark F. Adams         for (jj = 0; jj < nSAvec; jj++) {
7562e68589bSMark F. Adams           for (ii = 0; ii < nSAvec; ii++) {
75733812677SSatish Balay             if (data[jj*out_data_stride + ii] != PETSC_MAX_REAL) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"data[jj*out_data_stride + ii] != %e",PETSC_MAX_REAL);
7580cbbd2e1SMark F. Adams            if (ii <= jj) data[jj*out_data_stride + ii] = PetscRealPart(qqc[jj*Mdata + ii]);
7590cbbd2e1SMark F. Adams            else data[jj*out_data_stride + ii] = 0.;
7602e68589bSMark F. Adams           }
7612e68589bSMark F. Adams         }
7622e68589bSMark F. Adams       }
7632e68589bSMark F. Adams 
7642e68589bSMark F. Adams       /* get Q - row oriented */
765c964aadfSJose E. Roman       PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
766c666adbfSMark F. Adams       if (INFO != 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xORGQR error arg %d",-INFO);
7672e68589bSMark F. Adams 
7682e68589bSMark F. Adams       for (ii = 0; ii < M; ii++) {
7692e68589bSMark F. Adams         for (jj = 0; jj < N; jj++) {
7702e68589bSMark F. Adams           qqr[N*ii + jj] = qqc[jj*Mdata + ii];
7712e68589bSMark F. Adams         }
7722e68589bSMark F. Adams       }
7732e68589bSMark F. Adams 
7742e68589bSMark F. Adams       /* add diagonal block of P0 */
775c8b0795cSMark F. Adams       for (kk=0; kk<N; kk++) {
776c8b0795cSMark F. Adams         cids[kk] = N*cgid + kk; /* global col IDs in P0 */
777c8b0795cSMark F. Adams       }
7782e68589bSMark F. Adams       ierr = MatSetValues(a_Prol,M,fids,N,cids,qqr,INSERT_VALUES);CHKERRQ(ierr);
779e632b94dSBarry Smith       ierr = PetscFree5(qqc,qqr,TAU,WORK,fids);CHKERRQ(ierr);
780b43b03e9SMark F. Adams       clid++;
7810cbbd2e1SMark F. Adams     } /* coarse agg */
7820cbbd2e1SMark F. Adams   } /* for all fine nodes */
7830cbbd2e1SMark F. Adams   ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7840cbbd2e1SMark F. Adams   ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7851943db53SBarry Smith   ierr = PCGAMGHashTableDestroy(&fgid_flid);CHKERRQ(ierr);
7862e68589bSMark F. Adams   PetscFunctionReturn(0);
7872e68589bSMark F. Adams }
7882e68589bSMark F. Adams 
7895adeb434SBarry Smith static PetscErrorCode PCView_GAMG_AGG(PC pc,PetscViewer viewer)
7905adeb434SBarry Smith {
7915adeb434SBarry Smith   PetscErrorCode ierr;
7925adeb434SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
7935adeb434SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
7945adeb434SBarry Smith   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
7955adeb434SBarry Smith 
7965adeb434SBarry Smith   PetscFunctionBegin;
7975adeb434SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"      AGG specific options\n");CHKERRQ(ierr);
7985adeb434SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"        Symmetric graph %s\n",pc_gamg_agg->sym_graph ? "true" : "false");CHKERRQ(ierr);
799cab9ed1eSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"        Number of levels to square graph %D\n",pc_gamg_agg->square_graph);CHKERRQ(ierr);
800cab9ed1eSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"        Number smoothing steps %D\n",pc_gamg_agg->nsmooths);CHKERRQ(ierr);
8015adeb434SBarry Smith   PetscFunctionReturn(0);
8025adeb434SBarry Smith }
8035adeb434SBarry Smith 
8042e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
8052e68589bSMark F. Adams /*
806fd1112cbSBarry Smith    PCGAMGGraph_AGG
8072e68589bSMark F. Adams 
8082e68589bSMark F. Adams   Input Parameter:
8092e68589bSMark F. Adams    . pc - this
8102e68589bSMark F. Adams    . Amat - matrix on this fine level
8112e68589bSMark F. Adams   Output Parameter:
812c8b0795cSMark F. Adams    . a_Gmat -
8132e68589bSMark F. Adams */
814e0877f53SBarry Smith static PetscErrorCode PCGAMGGraph_AGG(PC pc,Mat Amat,Mat *a_Gmat)
815c8b0795cSMark F. Adams {
816c8b0795cSMark F. Adams   PetscErrorCode            ierr;
817c8b0795cSMark F. Adams   PC_MG                     *mg          = (PC_MG*)pc->data;
818c8b0795cSMark F. Adams   PC_GAMG                   *pc_gamg     = (PC_GAMG*)mg->innerctx;
819c1eae691SMark Adams   const PetscReal           vfilter      = pc_gamg->threshold[pc_gamg->current_level];
820c8b0795cSMark F. Adams   PC_GAMG_AGG               *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
821e0940f08SMark F. Adams   Mat                       Gmat;
8223b4367a7SBarry Smith   MPI_Comm                  comm;
82367744fedSMark F. Adams   PetscBool /* set,flg , */ symm;
824c8b0795cSMark F. Adams 
825c8b0795cSMark F. Adams   PetscFunctionBegin;
8263b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr);
827fd1112cbSBarry Smith   ierr = PetscLogEventBegin(PC_GAMGGraph_AGG,0,0,0,0);CHKERRQ(ierr);
828c8b0795cSMark F. Adams 
82967744fedSMark F. Adams   /* ierr = MatIsSymmetricKnown(Amat, &set, &flg);CHKERRQ(ierr); || !(set && flg) -- this causes lot of symm calls */
830c666adbfSMark F. Adams   symm = (PetscBool)(pc_gamg_agg->sym_graph); /* && !pc_gamg_agg->square_graph; */
8310cbbd2e1SMark F. Adams 
8322d7fac45SMark F. Adams   ierr = PCGAMGCreateGraph(Amat, &Gmat);CHKERRQ(ierr);
833bf4339c2SBarry Smith   ierr = PCGAMGFilterGraph(&Gmat, vfilter, symm);CHKERRQ(ierr);
834e0940f08SMark F. Adams   *a_Gmat = Gmat;
835fd1112cbSBarry Smith   ierr = PetscLogEventEnd(PC_GAMGGraph_AGG,0,0,0,0);CHKERRQ(ierr);
836c8b0795cSMark F. Adams   PetscFunctionReturn(0);
837c8b0795cSMark F. Adams }
838c8b0795cSMark F. Adams 
839c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
840c8b0795cSMark F. Adams /*
841b43b03e9SMark F. Adams    PCGAMGCoarsen_AGG
842c8b0795cSMark F. Adams 
843c8b0795cSMark F. Adams   Input Parameter:
844e0940f08SMark F. Adams    . a_pc - this
845e0940f08SMark F. Adams   Input/Output Parameter:
8460cbbd2e1SMark F. Adams    . a_Gmat1 - graph on this fine level - coarsening can change this (squares it)
847c8b0795cSMark F. Adams   Output Parameter:
8480cbbd2e1SMark F. Adams    . agg_lists - list of aggregates
849c8b0795cSMark F. Adams */
850e0877f53SBarry Smith static PetscErrorCode PCGAMGCoarsen_AGG(PC a_pc,Mat *a_Gmat1,PetscCoarsenData **agg_lists)
851c8b0795cSMark F. Adams {
852c8b0795cSMark F. Adams   PetscErrorCode ierr;
853e0940f08SMark F. Adams   PC_MG          *mg          = (PC_MG*)a_pc->data;
854c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
855c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
8560cbbd2e1SMark F. Adams   Mat            mat,Gmat2, Gmat1 = *a_Gmat1;  /* squared graph */
8570cbbd2e1SMark F. Adams   IS             perm;
8585cfd4bd9SMark Adams   PetscInt       Istart,Iend,Ii,nloc,bs,n,m;
859c8b0795cSMark F. Adams   PetscInt       *permute;
860c8b0795cSMark F. Adams   PetscBool      *bIndexSet;
861b43b03e9SMark F. Adams   MatCoarsen     crs;
8623b4367a7SBarry Smith   MPI_Comm       comm;
86373911c69SBarry Smith   PetscMPIInt    rank;
864aea53286SMark Adams   PetscReal      hashfact;
865e2c00dcbSBarry Smith   PetscInt       iSwapIndex;
8663b50add6SMark Adams   PetscRandom    random;
867c8b0795cSMark F. Adams 
868c8b0795cSMark F. Adams   PetscFunctionBegin;
8690cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr);
8703b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Gmat1,&comm);CHKERRQ(ierr);
8713b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
872e0940f08SMark F. Adams   ierr = MatGetLocalSize(Gmat1, &n, &m);CHKERRQ(ierr);
87371959b99SBarry Smith   ierr = MatGetBlockSize(Gmat1, &bs);CHKERRQ(ierr);
87471959b99SBarry Smith   if (bs != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"bs %D must be 1",bs);
875c8b0795cSMark F. Adams   nloc = n/bs;
876c8b0795cSMark F. Adams 
8779ab59c8bSMark Adams   if (pc_gamg->current_level < pc_gamg_agg->square_graph) {
878*b817416eSBarry Smith     ierr = PetscInfo2(a_pc,"Square Graph on level %D of %D to square\n",pc_gamg->current_level+1,pc_gamg_agg->square_graph);CHKERRQ(ierr);
879806fa848SBarry Smith     ierr = MatTransposeMatMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2);CHKERRQ(ierr);
880806fa848SBarry Smith   } else Gmat2 = Gmat1;
881c8b0795cSMark F. Adams 
8825cfd4bd9SMark Adams   /* get MIS aggs - randomize */
883785e854fSJed Brown   ierr = PetscMalloc1(nloc, &permute);CHKERRQ(ierr);
884e2c00dcbSBarry Smith   ierr = PetscCalloc1(nloc, &bIndexSet);CHKERRQ(ierr);
885c8b0795cSMark F. Adams   for (Ii = 0; Ii < nloc; Ii++) {
886c8b0795cSMark F. Adams     permute[Ii]   = Ii;
887c8b0795cSMark F. Adams   }
8883b50add6SMark Adams   ierr = PetscRandomCreate(PETSC_COMM_SELF,&random);CHKERRQ(ierr);
8895cfd4bd9SMark Adams   ierr = MatGetOwnershipRange(Gmat1, &Istart, &Iend);CHKERRQ(ierr);
890c8b0795cSMark F. Adams   for (Ii = 0; Ii < nloc; Ii++) {
8913b50add6SMark Adams     ierr = PetscRandomGetValueReal(random,&hashfact);CHKERRQ(ierr);
892aea53286SMark Adams     iSwapIndex = (PetscInt) (hashfact*nloc)%nloc;
893c8b0795cSMark F. Adams     if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) {
894c8b0795cSMark F. Adams       PetscInt iTemp = permute[iSwapIndex];
895c8b0795cSMark F. Adams       permute[iSwapIndex]   = permute[Ii];
896c8b0795cSMark F. Adams       permute[Ii]           = iTemp;
897c8b0795cSMark F. Adams       bIndexSet[iSwapIndex] = PETSC_TRUE;
898c8b0795cSMark F. Adams     }
899c8b0795cSMark F. Adams   }
900c8b0795cSMark F. Adams   ierr = PetscFree(bIndexSet);CHKERRQ(ierr);
9013b50add6SMark Adams   ierr = PetscRandomDestroy(&random);CHKERRQ(ierr);
902806fa848SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm);CHKERRQ(ierr);
9030cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
9040cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr);
905b43b03e9SMark F. Adams #endif
9063b4367a7SBarry Smith   ierr = MatCoarsenCreate(comm, &crs);CHKERRQ(ierr);
9079057884aSMark F. Adams   ierr = MatCoarsenSetFromOptions(crs);CHKERRQ(ierr);
908b43b03e9SMark F. Adams   ierr = MatCoarsenSetGreedyOrdering(crs, perm);CHKERRQ(ierr);
909b43b03e9SMark F. Adams   ierr = MatCoarsenSetAdjacency(crs, Gmat2);CHKERRQ(ierr);
910b43b03e9SMark F. Adams   ierr = MatCoarsenSetStrictAggs(crs, PETSC_TRUE);CHKERRQ(ierr);
911b43b03e9SMark F. Adams   ierr = MatCoarsenApply(crs);CHKERRQ(ierr);
9120cbbd2e1SMark F. Adams   ierr = MatCoarsenGetData(crs, agg_lists);CHKERRQ(ierr); /* output */
913b43b03e9SMark F. Adams   ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr);
914b43b03e9SMark F. Adams 
915c8b0795cSMark F. Adams   ierr = ISDestroy(&perm);CHKERRQ(ierr);
916c8b0795cSMark F. Adams   ierr = PetscFree(permute);CHKERRQ(ierr);
9170cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
9180cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr);
919b43b03e9SMark F. Adams #endif
9209f3f12c8SMark F. Adams 
921c8b0795cSMark F. Adams   /* smooth aggs */
922e0940f08SMark F. Adams   if (Gmat2 != Gmat1) {
9230cbbd2e1SMark F. Adams     const PetscCoarsenData *llist = *agg_lists;
924e3df7ea3SBarry Smith     ierr     = smoothAggs(a_pc,Gmat2, Gmat1, *agg_lists);CHKERRQ(ierr);
925c8b0795cSMark F. Adams     ierr     = MatDestroy(&Gmat1);CHKERRQ(ierr);
926e0940f08SMark F. Adams     *a_Gmat1 = Gmat2; /* output */
92741b27cdeSMark F. Adams     ierr     = PetscCDGetMat(llist, &mat);CHKERRQ(ierr);
9283b4367a7SBarry Smith     if (mat) SETERRQ(comm,PETSC_ERR_ARG_WRONG, "Auxilary matrix with squared graph????");
929806fa848SBarry Smith   } else {
9300cbbd2e1SMark F. Adams     const PetscCoarsenData *llist = *agg_lists;
9319ab59c8bSMark Adams     /* see if we have a matrix that takes precedence (returned from MatCoarsenApply) */
93241b27cdeSMark F. Adams     ierr = PetscCDGetMat(llist, &mat);CHKERRQ(ierr);
9330cbbd2e1SMark F. Adams     if (mat) {
9340cbbd2e1SMark F. Adams       ierr     = MatDestroy(&Gmat1);CHKERRQ(ierr);
9350cbbd2e1SMark F. Adams       *a_Gmat1 = mat; /* output */
9360cbbd2e1SMark F. Adams     }
9370cbbd2e1SMark F. Adams   }
9380cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr);
939c8b0795cSMark F. Adams   PetscFunctionReturn(0);
940c8b0795cSMark F. Adams }
941c8b0795cSMark F. Adams 
942c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
943c8b0795cSMark F. Adams /*
9440cbbd2e1SMark F. Adams  PCGAMGProlongator_AGG
945c8b0795cSMark F. Adams 
946c8b0795cSMark F. Adams  Input Parameter:
947c8b0795cSMark F. Adams  . pc - this
948c8b0795cSMark F. Adams  . Amat - matrix on this fine level
949c8b0795cSMark F. Adams  . Graph - used to get ghost data for nodes in
9500cbbd2e1SMark F. Adams  . agg_lists - list of aggregates
951c8b0795cSMark F. Adams  Output Parameter:
952c8b0795cSMark F. Adams  . a_P_out - prolongation operator to the next level
953c8b0795cSMark F. Adams  */
954e0877f53SBarry Smith static PetscErrorCode PCGAMGProlongator_AGG(PC pc,Mat Amat,Mat Gmat,PetscCoarsenData *agg_lists,Mat *a_P_out)
9552e68589bSMark F. Adams {
9562e68589bSMark F. Adams   PC_MG          *mg       = (PC_MG*)pc->data;
9572e68589bSMark F. Adams   PC_GAMG        *pc_gamg  = (PC_GAMG*)mg->innerctx;
9584a2f8832SBarry Smith   const PetscInt col_bs = pc_gamg->data_cell_cols;
9592e68589bSMark F. Adams   PetscErrorCode ierr;
960c8b0795cSMark F. Adams   PetscInt       Istart,Iend,nloc,ii,jj,kk,my0,nLocalSelected,bs;
961c8b0795cSMark F. Adams   Mat            Prol;
962c5df96a5SBarry Smith   PetscMPIInt    rank, size;
9633b4367a7SBarry Smith   MPI_Comm       comm;
964c8b0795cSMark F. Adams   PetscReal      *data_w_ghost;
965c8b0795cSMark F. Adams   PetscInt       myCrs0, nbnodes=0, *flid_fgid;
966d9558ea9SBarry Smith   MatType        mtype;
9672e68589bSMark F. Adams 
9682e68589bSMark F. Adams   PetscFunctionBegin;
9693b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr);
9704a2f8832SBarry Smith   if (col_bs < 1) SETERRQ(comm,PETSC_ERR_PLIB,"Column bs cannot be less than 1");
9710cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
9723b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
9733b4367a7SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
9742e68589bSMark F. Adams   ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr);
975c8b0795cSMark F. Adams   ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr);
97671959b99SBarry Smith   nloc = (Iend-Istart)/bs; my0 = Istart/bs;
97771959b99SBarry Smith   if ((Iend-Istart) % bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"(Iend %D - Istart %D) not divisible by bs %D",Iend,Istart,bs);
9782e68589bSMark F. Adams 
9792e68589bSMark F. Adams   /* get 'nLocalSelected' */
9800cbbd2e1SMark F. Adams   for (ii=0, nLocalSelected = 0; ii < nloc; ii++) {
981e78576d6SMark F. Adams     PetscBool ise;
9820cbbd2e1SMark F. Adams     /* filter out singletons 0 or 1? */
983e78576d6SMark F. Adams     ierr = PetscCDEmptyAt(agg_lists, ii, &ise);CHKERRQ(ierr);
984e78576d6SMark F. Adams     if (!ise) nLocalSelected++;
9852e68589bSMark F. Adams   }
9862e68589bSMark F. Adams 
9872e68589bSMark F. Adams   /* create prolongator, create P matrix */
988d9558ea9SBarry Smith   ierr = MatGetType(Amat,&mtype);CHKERRQ(ierr);
9893b4367a7SBarry Smith   ierr = MatCreate(comm, &Prol);CHKERRQ(ierr);
990806fa848SBarry Smith   ierr = MatSetSizes(Prol,nloc*bs,nLocalSelected*col_bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
991a2f3521dSMark F. Adams   ierr = MatSetBlockSizes(Prol, bs, col_bs);CHKERRQ(ierr);
992d9558ea9SBarry Smith   ierr = MatSetType(Prol, mtype);CHKERRQ(ierr);
9934a2f8832SBarry Smith   ierr = MatSeqAIJSetPreallocation(Prol, col_bs, NULL);CHKERRQ(ierr);
9944a2f8832SBarry Smith   ierr = MatMPIAIJSetPreallocation(Prol,col_bs, NULL,col_bs, NULL);CHKERRQ(ierr);
9952e68589bSMark F. Adams 
9962e68589bSMark F. Adams   /* can get all points "removed" */
997c8b0795cSMark F. Adams   ierr =  MatGetSize(Prol, &kk, &ii);CHKERRQ(ierr);
9987f66b68fSMark Adams   if (!ii) {
999569f4572SMark Adams     ierr = PetscInfo(pc,"No selected points on coarse grid\n");CHKERRQ(ierr);
10002e68589bSMark F. Adams     ierr = MatDestroy(&Prol);CHKERRQ(ierr);
10010298fd71SBarry Smith     *a_P_out = NULL;  /* out */
1002e87675ddSBarry Smith     ierr = PetscLogEventEnd(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
10032e68589bSMark F. Adams     PetscFunctionReturn(0);
10042e68589bSMark F. Adams   }
1005bf4339c2SBarry Smith   ierr = PetscInfo1(pc,"New grid %D nodes\n",ii/col_bs);CHKERRQ(ierr);
1006c8b0795cSMark F. Adams   ierr = MatGetOwnershipRangeColumn(Prol, &myCrs0, &kk);CHKERRQ(ierr);
10070cbbd2e1SMark F. Adams 
100871959b99SBarry Smith   if ((kk-myCrs0) % col_bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"(kk %D -myCrs0 %D) not divisible by col_bs %D",kk,myCrs0,col_bs);
1009c8b0795cSMark F. Adams   myCrs0 = myCrs0/col_bs;
101071959b99SBarry Smith   if ((kk/col_bs-myCrs0) != nLocalSelected) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_PLIB,"(kk %D/col_bs %D - myCrs0 %D) != nLocalSelected %D)",kk,col_bs,myCrs0,nLocalSelected);
10112e68589bSMark F. Adams 
10122e68589bSMark F. Adams   /* create global vector of data in 'data_w_ghost' */
10130cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
10140cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr);
10152e68589bSMark F. Adams #endif
1016c5df96a5SBarry Smith   if (size > 1) { /*  */
10172e68589bSMark F. Adams     PetscReal *tmp_gdata,*tmp_ldata,*tp2;
1018785e854fSJed Brown     ierr = PetscMalloc1(nloc, &tmp_ldata);CHKERRQ(ierr);
10194a2f8832SBarry Smith     for (jj = 0; jj < col_bs; jj++) {
1020c8b0795cSMark F. Adams       for (kk = 0; kk < bs; kk++) {
1021a2f3521dSMark F. Adams         PetscInt        ii,stride;
1022c8b0795cSMark F. Adams         const PetscReal *tp = pc_gamg->data + jj*bs*nloc + kk;
10232fa5cd67SKarl Rupp         for (ii = 0; ii < nloc; ii++, tp += bs) tmp_ldata[ii] = *tp;
10242fa5cd67SKarl Rupp 
1025806fa848SBarry Smith         ierr = PCGAMGGetDataWithGhosts(Gmat, 1, tmp_ldata, &stride, &tmp_gdata);CHKERRQ(ierr);
1026a2f3521dSMark F. Adams 
10277f66b68fSMark Adams         if (!jj && !kk) { /* now I know how many todal nodes - allocate */
10284a2f8832SBarry Smith           ierr    = PetscMalloc1(stride*bs*col_bs, &data_w_ghost);CHKERRQ(ierr);
1029a2f3521dSMark F. Adams           nbnodes = bs*stride;
10302e68589bSMark F. Adams         }
1031a2f3521dSMark F. Adams         tp2 = data_w_ghost + jj*bs*stride + kk;
10322fa5cd67SKarl Rupp         for (ii = 0; ii < stride; ii++, tp2 += bs) *tp2 = tmp_gdata[ii];
10332e68589bSMark F. Adams         ierr = PetscFree(tmp_gdata);CHKERRQ(ierr);
10342e68589bSMark F. Adams       }
10352e68589bSMark F. Adams     }
10362e68589bSMark F. Adams     ierr = PetscFree(tmp_ldata);CHKERRQ(ierr);
1037806fa848SBarry Smith   } else {
1038c8b0795cSMark F. Adams     nbnodes      = bs*nloc;
1039c8b0795cSMark F. Adams     data_w_ghost = (PetscReal*)pc_gamg->data;
10402e68589bSMark F. Adams   }
10412e68589bSMark F. Adams 
10422e68589bSMark F. Adams   /* get P0 */
1043c5df96a5SBarry Smith   if (size > 1) {
10442e68589bSMark F. Adams     PetscReal *fid_glid_loc,*fiddata;
1045a2f3521dSMark F. Adams     PetscInt  stride;
10462e68589bSMark F. Adams 
1047785e854fSJed Brown     ierr = PetscMalloc1(nloc, &fid_glid_loc);CHKERRQ(ierr);
10482e68589bSMark F. Adams     for (kk=0; kk<nloc; kk++) fid_glid_loc[kk] = (PetscReal)(my0+kk);
1049806fa848SBarry Smith     ierr = PCGAMGGetDataWithGhosts(Gmat, 1, fid_glid_loc, &stride, &fiddata);CHKERRQ(ierr);
1050785e854fSJed Brown     ierr = PetscMalloc1(stride, &flid_fgid);CHKERRQ(ierr);
1051a2f3521dSMark F. Adams     for (kk=0; kk<stride; kk++) flid_fgid[kk] = (PetscInt)fiddata[kk];
10522e68589bSMark F. Adams     ierr = PetscFree(fiddata);CHKERRQ(ierr);
1053a2f3521dSMark F. Adams 
105471959b99SBarry Smith     if (stride != nbnodes/bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"stride %D != nbnodes %D/bs %D",stride,nbnodes,bs);
10552e68589bSMark F. Adams     ierr = PetscFree(fid_glid_loc);CHKERRQ(ierr);
1056806fa848SBarry Smith   } else {
1057785e854fSJed Brown     ierr = PetscMalloc1(nloc, &flid_fgid);CHKERRQ(ierr);
10582e68589bSMark F. Adams     for (kk=0; kk<nloc; kk++) flid_fgid[kk] = my0 + kk;
10592e68589bSMark F. Adams   }
10600cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
10610cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr);
10620cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr);
10632e68589bSMark F. Adams #endif
1064c8b0795cSMark F. Adams   {
10650298fd71SBarry Smith     PetscReal *data_out = NULL;
10664a2f8832SBarry Smith     ierr = formProl0(agg_lists, bs, col_bs, myCrs0, nbnodes,data_w_ghost, flid_fgid, &data_out, Prol);CHKERRQ(ierr);
1067c8b0795cSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
10682fa5cd67SKarl Rupp 
1069c8b0795cSMark F. Adams     pc_gamg->data           = data_out;
10704a2f8832SBarry Smith     pc_gamg->data_cell_rows = col_bs;
10714a2f8832SBarry Smith     pc_gamg->data_sz        = col_bs*col_bs*nLocalSelected;
1072c8b0795cSMark F. Adams   }
10730cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
10740cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr);
1075c8b0795cSMark F. Adams #endif
107661ba4676SBarry Smith   if (size > 1) {ierr = PetscFree(data_w_ghost);CHKERRQ(ierr);}
10772e68589bSMark F. Adams   ierr = PetscFree(flid_fgid);CHKERRQ(ierr);
10782e68589bSMark F. Adams 
1079c8b0795cSMark F. Adams   *a_P_out = Prol;  /* out */
1080ed4e82eaSPeter Brune 
10810cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
1082c8b0795cSMark F. Adams   PetscFunctionReturn(0);
1083c8b0795cSMark F. Adams }
1084c8b0795cSMark F. Adams 
1085c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
1086c8b0795cSMark F. Adams /*
1087fd1112cbSBarry Smith    PCGAMGOptProlongator_AGG
1088c8b0795cSMark F. Adams 
1089c8b0795cSMark F. Adams   Input Parameter:
1090c8b0795cSMark F. Adams    . pc - this
1091c8b0795cSMark F. Adams    . Amat - matrix on this fine level
1092c8b0795cSMark F. Adams  In/Output Parameter:
10932a808120SBarry Smith    . a_P - prolongation operator to the next level
1094c8b0795cSMark F. Adams */
1095e0877f53SBarry Smith static PetscErrorCode PCGAMGOptProlongator_AGG(PC pc,Mat Amat,Mat *a_P)
1096c8b0795cSMark F. Adams {
1097c8b0795cSMark F. Adams   PetscErrorCode ierr;
1098c8b0795cSMark F. Adams   PC_MG          *mg          = (PC_MG*)pc->data;
1099c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
1100c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
1101c8b0795cSMark F. Adams   PetscInt       jj;
1102c8b0795cSMark F. Adams   Mat            Prol  = *a_P;
11033b4367a7SBarry Smith   MPI_Comm       comm;
11042a808120SBarry Smith   KSP            eksp;
11052a808120SBarry Smith   Vec            bb, xx;
11062a808120SBarry Smith   PC             epc;
11072a808120SBarry Smith   PetscReal      alpha, emax, emin;
11083b50add6SMark Adams   PetscRandom    random;
1109c8b0795cSMark F. Adams 
1110c8b0795cSMark F. Adams   PetscFunctionBegin;
11113b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr);
1112fd1112cbSBarry Smith   ierr = PetscLogEventBegin(PC_GAMGOptProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
1113c8b0795cSMark F. Adams 
11142a808120SBarry Smith   /* compute maximum value of operator to be used in smoother */
11152a808120SBarry Smith   if (0 < pc_gamg_agg->nsmooths) {
1116e7e70905SMark Adams     ierr = MatCreateVecs(Amat, &bb, 0);CHKERRQ(ierr);
1117e7e70905SMark Adams     ierr = MatCreateVecs(Amat, &xx, 0);CHKERRQ(ierr);
11183b50add6SMark Adams     ierr = PetscRandomCreate(PETSC_COMM_SELF,&random);CHKERRQ(ierr);
11193b50add6SMark Adams     ierr = VecSetRandom(bb,random);CHKERRQ(ierr);
11203b50add6SMark Adams     ierr = PetscRandomDestroy(&random);CHKERRQ(ierr);
1121e696ed0bSMark F. Adams 
11223b4367a7SBarry Smith     ierr = KSPCreate(comm,&eksp);CHKERRQ(ierr);
1123806fa848SBarry Smith     ierr = KSPSetTolerances(eksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,10);CHKERRQ(ierr);
11242e68589bSMark F. Adams     ierr = KSPSetNormType(eksp, KSP_NORM_NONE);CHKERRQ(ierr);
11253ca41b3bSMark Adams     ierr = KSPSetErrorIfNotConverged(eksp,PETSC_FALSE);CHKERRQ(ierr);
1126c2436cacSMark F. Adams 
1127c2436cacSMark F. Adams     ierr = KSPSetInitialGuessNonzero(eksp, PETSC_FALSE);CHKERRQ(ierr);
112823ee1639SBarry Smith     ierr = KSPSetOperators(eksp, Amat, Amat);CHKERRQ(ierr);
11292e68589bSMark F. Adams     ierr = KSPSetComputeSingularValues(eksp,PETSC_TRUE);CHKERRQ(ierr);
11302e68589bSMark F. Adams 
1131da509fc8SJed Brown     ierr = KSPGetPC(eksp, &epc);CHKERRQ(ierr);
1132da509fc8SJed Brown     ierr = PCSetType(epc, PCJACOBI);CHKERRQ(ierr);  /* smoother in smoothed agg. */
1133c2436cacSMark F. Adams 
1134da7ca8b8SBarry Smith     ierr = KSPSetOptionsPrefix(eksp,((PetscObject)pc)->prefix);CHKERRQ(ierr);
11353ca41b3bSMark Adams     ierr = KSPAppendOptionsPrefix(eksp, "gamg_est_");CHKERRQ(ierr);
11363ca41b3bSMark Adams     ierr = KSPSetFromOptions(eksp);CHKERRQ(ierr);
11373ca41b3bSMark Adams 
11382e68589bSMark F. Adams     /* solve - keep stuff out of logging */
11392e68589bSMark F. Adams     ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr);
11402e68589bSMark F. Adams     ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr);
11412e68589bSMark F. Adams     ierr = KSPSolve(eksp, bb, xx);CHKERRQ(ierr);
11422e68589bSMark F. Adams     ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr);
11432e68589bSMark F. Adams     ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr);
11442e68589bSMark F. Adams 
11452e68589bSMark F. Adams     ierr = KSPComputeExtremeSingularValues(eksp, &emax, &emin);CHKERRQ(ierr);
1146569f4572SMark Adams     ierr = PetscInfo3(pc,"Smooth P0: max eigen=%e min=%e PC=%s\n",emax,emin,PCJACOBI);CHKERRQ(ierr);
11472e68589bSMark F. Adams     ierr = VecDestroy(&xx);CHKERRQ(ierr);
11482e68589bSMark F. Adams     ierr = VecDestroy(&bb);CHKERRQ(ierr);
11492e68589bSMark F. Adams     ierr = KSPDestroy(&eksp);CHKERRQ(ierr);
11502e68589bSMark F. Adams   }
11512e68589bSMark F. Adams 
11522a808120SBarry Smith   /* smooth P0 */
11532a808120SBarry Smith   for (jj = 0; jj < pc_gamg_agg->nsmooths; jj++) {
11542a808120SBarry Smith     Mat       tMat;
11552a808120SBarry Smith     Vec       diag;
11562a808120SBarry Smith 
11572a808120SBarry Smith #if defined PETSC_GAMG_USE_LOG
11582a808120SBarry Smith     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr);
11592a808120SBarry Smith #endif
11602a808120SBarry Smith 
11612e68589bSMark F. Adams     /* smooth P1 := (I - omega/lam D^{-1}A)P0 */
11622e68589bSMark F. Adams     ierr  = MatMatMult(Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat);CHKERRQ(ierr);
11632a7a6963SBarry Smith     ierr  = MatCreateVecs(Amat, &diag, 0);CHKERRQ(ierr);
11642e68589bSMark F. Adams     ierr  = MatGetDiagonal(Amat, diag);CHKERRQ(ierr); /* effectively PCJACOBI */
11652e68589bSMark F. Adams     ierr  = VecReciprocal(diag);CHKERRQ(ierr);
11662e68589bSMark F. Adams     ierr  = MatDiagonalScale(tMat, diag, 0);CHKERRQ(ierr);
11672e68589bSMark F. Adams     ierr  = VecDestroy(&diag);CHKERRQ(ierr);
1168b4ec6429SMark F. Adams     alpha = -1.4/emax;
11692e68589bSMark F. Adams     ierr  = MatAYPX(tMat, alpha, Prol, SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
11702e68589bSMark F. Adams     ierr  = MatDestroy(&Prol);CHKERRQ(ierr);
11712e68589bSMark F. Adams     Prol  = tMat;
11720cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
11730cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr);
11742e68589bSMark F. Adams #endif
11752e68589bSMark F. Adams   }
1176fd1112cbSBarry Smith   ierr = PetscLogEventEnd(PC_GAMGOptProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
1177c8b0795cSMark F. Adams   *a_P = Prol;
1178c8b0795cSMark F. Adams   PetscFunctionReturn(0);
11792e68589bSMark F. Adams }
11802e68589bSMark F. Adams 
1181c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
1182c8b0795cSMark F. Adams /*
1183c8b0795cSMark F. Adams    PCCreateGAMG_AGG
11842e68589bSMark F. Adams 
1185c8b0795cSMark F. Adams   Input Parameter:
1186c8b0795cSMark F. Adams    . pc -
1187c8b0795cSMark F. Adams */
1188c8b0795cSMark F. Adams PetscErrorCode  PCCreateGAMG_AGG(PC pc)
1189c8b0795cSMark F. Adams {
1190c8b0795cSMark F. Adams   PetscErrorCode ierr;
1191c8b0795cSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
1192c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1193c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg;
11942e68589bSMark F. Adams 
1195c8b0795cSMark F. Adams   PetscFunctionBegin;
1196c8b0795cSMark F. Adams   /* create sub context for SA */
1197b00a9115SJed Brown   ierr            = PetscNewLog(pc,&pc_gamg_agg);CHKERRQ(ierr);
1198c8b0795cSMark F. Adams   pc_gamg->subctx = pc_gamg_agg;
1199c8b0795cSMark F. Adams 
12001ab5ffc9SJed Brown   pc_gamg->ops->setfromoptions = PCSetFromOptions_GAMG_AGG;
12019b8ffb57SJed Brown   pc_gamg->ops->destroy        = PCDestroy_GAMG_AGG;
1202c8b0795cSMark F. Adams   /* reset does not do anything; setup not virtual */
1203c8b0795cSMark F. Adams 
1204c8b0795cSMark F. Adams   /* set internal function pointers */
1205fd1112cbSBarry Smith   pc_gamg->ops->graph             = PCGAMGGraph_AGG;
12061ab5ffc9SJed Brown   pc_gamg->ops->coarsen           = PCGAMGCoarsen_AGG;
12071ab5ffc9SJed Brown   pc_gamg->ops->prolongator       = PCGAMGProlongator_AGG;
1208fd1112cbSBarry Smith   pc_gamg->ops->optprolongator    = PCGAMGOptProlongator_AGG;
12091ab5ffc9SJed Brown   pc_gamg->ops->createdefaultdata = PCSetData_AGG;
12105adeb434SBarry Smith   pc_gamg->ops->view              = PCView_GAMG_AGG;
1211c8b0795cSMark F. Adams 
1212cab9ed1eSBarry Smith   pc_gamg_agg->square_graph = 1;
1213cab9ed1eSBarry Smith   pc_gamg_agg->sym_graph    = PETSC_FALSE;
1214cab9ed1eSBarry Smith   pc_gamg_agg->nsmooths     = 1;
1215cab9ed1eSBarry Smith 
1216e0877f53SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNSmooths_C",PCGAMGSetNSmooths_AGG);CHKERRQ(ierr);
1217e0877f53SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetSymGraph_C",PCGAMGSetSymGraph_AGG);CHKERRQ(ierr);
1218e0877f53SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetSquareGraph_C",PCGAMGSetSquareGraph_AGG);CHKERRQ(ierr);
1219bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_AGG);CHKERRQ(ierr);
12202e68589bSMark F. Adams   PetscFunctionReturn(0);
12212e68589bSMark F. Adams }
1222