xref: /petsc/src/ksp/pc/impls/gamg/agg.c (revision e3df7ea3ec46b014a8ed04aa00925be89b358ba7)
12e68589bSMark F. Adams /*
22e68589bSMark F. Adams  GAMG geometric-algebric multiogrid 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 #undef __FUNCT__
182e68589bSMark F. Adams #define __FUNCT__ "PCGAMGSetNSmooths"
192e68589bSMark F. Adams /*@
202e68589bSMark F. Adams    PCGAMGSetNSmooths - Set number of smoothing steps (1 is typical)
212e68589bSMark F. Adams 
222e68589bSMark F. Adams    Not Collective on PC
232e68589bSMark F. Adams 
242e68589bSMark F. Adams    Input Parameters:
252e68589bSMark F. Adams .  pc - the preconditioner context
262e68589bSMark F. Adams 
272e68589bSMark F. Adams    Options Database Key:
28cab9ed1eSBarry Smith .  -pc_gamg_agg_nsmooths <nsmooth, default=1> - number of smoothing steps to use with smooth aggregation
292e68589bSMark F. Adams 
302e68589bSMark F. Adams    Level: intermediate
312e68589bSMark F. Adams 
322e68589bSMark F. Adams    Concepts: Aggregation AMG preconditioner
332e68589bSMark F. Adams 
342e68589bSMark F. Adams .seealso: ()
352e68589bSMark F. Adams @*/
362e68589bSMark F. Adams PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n)
372e68589bSMark F. Adams {
382e68589bSMark F. Adams   PetscErrorCode ierr;
392e68589bSMark F. Adams 
402e68589bSMark F. Adams   PetscFunctionBegin;
412e68589bSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
422e68589bSMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetNSmooths_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
432e68589bSMark F. Adams   PetscFunctionReturn(0);
442e68589bSMark F. Adams }
452e68589bSMark F. Adams 
462e68589bSMark F. Adams #undef __FUNCT__
47e0877f53SBarry Smith #define __FUNCT__ "PCGAMGSetNSmooths_AGG"
48e0877f53SBarry Smith static PetscErrorCode PCGAMGSetNSmooths_AGG(PC pc, PetscInt n)
492e68589bSMark F. Adams {
502e68589bSMark F. Adams   PC_MG       *mg          = (PC_MG*)pc->data;
512e68589bSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG*)mg->innerctx;
52c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
532e68589bSMark F. Adams 
542e68589bSMark F. Adams   PetscFunctionBegin;
55c8b0795cSMark F. Adams   pc_gamg_agg->nsmooths = n;
56c8b0795cSMark F. Adams   PetscFunctionReturn(0);
57c8b0795cSMark F. Adams }
58c8b0795cSMark F. Adams 
59c8b0795cSMark F. Adams #undef __FUNCT__
60c8b0795cSMark F. Adams #define __FUNCT__ "PCGAMGSetSymGraph"
61c8b0795cSMark F. Adams /*@
62cab9ed1eSBarry Smith    PCGAMGSetSymGraph - Symmetrize the graph before computing the aggregation. Some algorithms require the graph be symmetric
63c8b0795cSMark F. Adams 
64c8b0795cSMark F. Adams    Not Collective on PC
65c8b0795cSMark F. Adams 
66c8b0795cSMark F. Adams    Input Parameters:
67cab9ed1eSBarry Smith +  pc - the preconditioner context
68cab9ed1eSBarry Smith .  n - PETSC_TRUE or PETSC_FALSE
69c8b0795cSMark F. Adams 
70c8b0795cSMark F. Adams    Options Database Key:
71cab9ed1eSBarry Smith .  -pc_gamg_sym_graph <true,default=false> - symmetrize the graph before computing the aggregation
72c8b0795cSMark F. Adams 
73c8b0795cSMark F. Adams    Level: intermediate
74c8b0795cSMark F. Adams 
75c8b0795cSMark F. Adams    Concepts: Aggregation AMG preconditioner
76c8b0795cSMark F. Adams 
77cab9ed1eSBarry Smith .seealso: PCGAMGSetSquareGraph()
78c8b0795cSMark F. Adams @*/
79c8b0795cSMark F. Adams PetscErrorCode PCGAMGSetSymGraph(PC pc, PetscBool n)
80c8b0795cSMark F. Adams {
81c8b0795cSMark F. Adams   PetscErrorCode ierr;
82c8b0795cSMark F. Adams 
83c8b0795cSMark F. Adams   PetscFunctionBegin;
84c8b0795cSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
85c8b0795cSMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetSymGraph_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
86c8b0795cSMark F. Adams   PetscFunctionReturn(0);
87c8b0795cSMark F. Adams }
88c8b0795cSMark F. Adams 
89c8b0795cSMark F. Adams #undef __FUNCT__
90e0877f53SBarry Smith #define __FUNCT__ "PCGAMGSetSymGraph_AGG"
91e0877f53SBarry Smith static PetscErrorCode PCGAMGSetSymGraph_AGG(PC pc, PetscBool n)
92c8b0795cSMark F. Adams {
93c8b0795cSMark F. Adams   PC_MG       *mg          = (PC_MG*)pc->data;
94c8b0795cSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG*)mg->innerctx;
95c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
96c8b0795cSMark F. Adams 
97c8b0795cSMark F. Adams   PetscFunctionBegin;
98c8b0795cSMark F. Adams   pc_gamg_agg->sym_graph = n;
992e68589bSMark F. Adams   PetscFunctionReturn(0);
1002e68589bSMark F. Adams }
1012e68589bSMark F. Adams 
102ef4ad70eSMark F. Adams #undef __FUNCT__
103ef4ad70eSMark F. Adams #define __FUNCT__ "PCGAMGSetSquareGraph"
104ef4ad70eSMark F. Adams /*@
105cab9ed1eSBarry Smith    PCGAMGSetSquareGraph -  Square the graph, ie. compute A'*A before aggregating it
106ef4ad70eSMark F. Adams 
107ef4ad70eSMark F. Adams    Not Collective on PC
108ef4ad70eSMark F. Adams 
109ef4ad70eSMark F. Adams    Input Parameters:
110cab9ed1eSBarry Smith +  pc - the preconditioner context
111cab9ed1eSBarry Smith -  n - PETSC_TRUE or PETSC_FALSE
112ef4ad70eSMark F. Adams 
113ef4ad70eSMark F. Adams    Options Database Key:
114cab9ed1eSBarry Smith .  -pc_gamg_square_graph <n,default = 1> - number of levels to square the graph on before aggregating it
115ef4ad70eSMark F. Adams 
116ef4ad70eSMark F. Adams    Level: intermediate
117ef4ad70eSMark F. Adams 
118ef4ad70eSMark F. Adams    Concepts: Aggregation AMG preconditioner
119ef4ad70eSMark F. Adams 
120cab9ed1eSBarry Smith .seealso: PCGAMGSetSymGraph()
121ef4ad70eSMark F. Adams @*/
1229ab59c8bSMark Adams PetscErrorCode PCGAMGSetSquareGraph(PC pc, PetscInt n)
123ef4ad70eSMark F. Adams {
124ef4ad70eSMark F. Adams   PetscErrorCode ierr;
125ef4ad70eSMark F. Adams 
126ef4ad70eSMark F. Adams   PetscFunctionBegin;
127ef4ad70eSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1289ab59c8bSMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetSquareGraph_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
129ef4ad70eSMark F. Adams   PetscFunctionReturn(0);
130ef4ad70eSMark F. Adams }
131ef4ad70eSMark F. Adams 
132ef4ad70eSMark F. Adams #undef __FUNCT__
133e0877f53SBarry Smith #define __FUNCT__ "PCGAMGSetSquareGraph_AGG"
134e0877f53SBarry Smith static PetscErrorCode PCGAMGSetSquareGraph_AGG(PC pc, PetscInt n)
135ef4ad70eSMark F. Adams {
136ef4ad70eSMark F. Adams   PC_MG       *mg          = (PC_MG*)pc->data;
137ef4ad70eSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG*)mg->innerctx;
138ef4ad70eSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
139ef4ad70eSMark F. Adams 
140ef4ad70eSMark F. Adams   PetscFunctionBegin;
141ef4ad70eSMark F. Adams   pc_gamg_agg->square_graph = n;
142ef4ad70eSMark F. Adams   PetscFunctionReturn(0);
143ef4ad70eSMark F. Adams }
144ef4ad70eSMark F. Adams 
1452e68589bSMark F. Adams #undef __FUNCT__
1462e68589bSMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG_AGG"
147e0877f53SBarry Smith static PetscErrorCode PCSetFromOptions_GAMG_AGG(PetscOptionItems *PetscOptionsObject,PC pc)
1482e68589bSMark F. Adams {
1492e68589bSMark F. Adams   PetscErrorCode ierr;
1502e68589bSMark F. Adams   PC_MG          *mg          = (PC_MG*)pc->data;
1512e68589bSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
152c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
1532e68589bSMark F. Adams 
1542e68589bSMark F. Adams   PetscFunctionBegin;
155e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"GAMG-AGG options");CHKERRQ(ierr);
1562e68589bSMark F. Adams   {
1578afaa268SBarry 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);
1588afaa268SBarry 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);
1599ab59c8bSMark 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);
1602e68589bSMark F. Adams   }
1612e68589bSMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
1622e68589bSMark F. Adams   PetscFunctionReturn(0);
1632e68589bSMark F. Adams }
1642e68589bSMark F. Adams 
1652e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
1662e68589bSMark F. Adams #undef __FUNCT__
1679b8ffb57SJed Brown #define __FUNCT__ "PCDestroy_GAMG_AGG"
168e0877f53SBarry Smith static PetscErrorCode PCDestroy_GAMG_AGG(PC pc)
1692e68589bSMark F. Adams {
1702e68589bSMark F. Adams   PetscErrorCode ierr;
1712e68589bSMark F. Adams   PC_MG          *mg          = (PC_MG*)pc->data;
1722e68589bSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
1732e68589bSMark F. Adams 
1742e68589bSMark F. Adams   PetscFunctionBegin;
1759b8ffb57SJed Brown   ierr = PetscFree(pc_gamg->subctx);CHKERRQ(ierr);
1763ae0bb68SMark Adams   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",NULL);CHKERRQ(ierr);
1772e68589bSMark F. Adams   PetscFunctionReturn(0);
1782e68589bSMark F. Adams }
1792e68589bSMark F. Adams 
1802e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
1812e68589bSMark F. Adams /*
1822e68589bSMark F. Adams    PCSetCoordinates_AGG
183302f38e8SMark F. Adams      - collective
1842e68589bSMark F. Adams 
1852e68589bSMark F. Adams    Input Parameter:
1862e68589bSMark F. Adams    . pc - the preconditioner context
187a2f3521dSMark F. Adams    . ndm - dimesion of data (used for dof/vertex for Stokes)
188302f38e8SMark F. Adams    . a_nloc - number of vertices local
189302f38e8SMark F. Adams    . coords - [a_nloc][ndm] - interleaved coordinate data: {x_0, y_0, z_0, x_1, y_1, ...}
1902e68589bSMark F. Adams */
1911e6b0712SBarry Smith 
1922e68589bSMark F. Adams #undef __FUNCT__
1932e68589bSMark F. Adams #define __FUNCT__ "PCSetCoordinates_AGG"
1941e6b0712SBarry Smith static PetscErrorCode PCSetCoordinates_AGG(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords)
1952e68589bSMark F. Adams {
1962e68589bSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
1972e68589bSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1982e68589bSMark F. Adams   PetscErrorCode ierr;
19969344418SMark F. Adams   PetscInt       arrsz,kk,ii,jj,nloc,ndatarows,ndf;
200a2f3521dSMark F. Adams   Mat            mat = pc->pmat;
2012e68589bSMark F. Adams 
2022e68589bSMark F. Adams   PetscFunctionBegin;
203a2f3521dSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
204a2f3521dSMark F. Adams   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
205302f38e8SMark F. Adams   nloc = a_nloc;
2062e68589bSMark F. Adams 
2072e68589bSMark F. Adams   /* SA: null space vectors */
20869344418SMark F. Adams   ierr = MatGetBlockSize(mat, &ndf);CHKERRQ(ierr); /* this does not work for Stokes */
20969344418SMark F. Adams   if (coords && ndf==1) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */
210a2f3521dSMark F. Adams   else if (coords) {
211c666adbfSMark F. Adams     if (ndm > ndf) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"degrees of motion %d > block size %d",ndm,ndf);
21269344418SMark F. Adams     pc_gamg->data_cell_cols = (ndm==2 ? 3 : 6); /* displacement elasticity */
21369344418SMark F. Adams     if (ndm != ndf) {
214c666adbfSMark F. Adams       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);
215a2f3521dSMark F. Adams     }
216806fa848SBarry Smith   } else pc_gamg->data_cell_cols = ndf; /* no data, force SA with constant null space vectors */
21771959b99SBarry Smith   pc_gamg->data_cell_rows = ndatarows = ndf;
21871959b99SBarry 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);
219c8b0795cSMark F. Adams   arrsz = nloc*pc_gamg->data_cell_rows*pc_gamg->data_cell_cols;
2202e68589bSMark F. Adams 
2212e68589bSMark F. Adams   /* create data - syntactic sugar that should be refactored at some point */
2227f66b68fSMark Adams   if (!pc_gamg->data || (pc_gamg->data_sz != arrsz)) {
2232e68589bSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
224854ce69bSBarry Smith     ierr = PetscMalloc1(arrsz+1, &pc_gamg->data);CHKERRQ(ierr);
2252e68589bSMark F. Adams   }
2262e68589bSMark F. Adams   /* copy data in - column oriented */
2272e68589bSMark F. Adams   for (kk=0; kk<nloc; kk++) {
22869344418SMark F. Adams     const PetscInt M     = nloc*pc_gamg->data_cell_rows; /* stride into data */
22969344418SMark F. Adams     PetscReal      *data = &pc_gamg->data[kk*ndatarows]; /* start of cell */
230c8b0795cSMark F. Adams     if (pc_gamg->data_cell_cols==1) *data = 1.0;
2312e68589bSMark F. Adams     else {
23269344418SMark F. Adams       /* translational modes */
2332fa5cd67SKarl Rupp       for (ii=0;ii<ndatarows;ii++) {
2342fa5cd67SKarl Rupp         for (jj=0;jj<ndatarows;jj++) {
23569344418SMark F. Adams           if (ii==jj)data[ii*M + jj] = 1.0;
2362e68589bSMark F. Adams           else data[ii*M + jj] = 0.0;
2372fa5cd67SKarl Rupp         }
2382fa5cd67SKarl Rupp       }
23969344418SMark F. Adams 
24069344418SMark F. Adams       /* rotational modes */
2412e68589bSMark F. Adams       if (coords) {
24269344418SMark F. Adams         if (ndm == 2) {
2432e68589bSMark F. Adams           data   += 2*M;
2442e68589bSMark F. Adams           data[0] = -coords[2*kk+1];
2452e68589bSMark F. Adams           data[1] =  coords[2*kk];
246806fa848SBarry Smith         } else {
2472e68589bSMark F. Adams           data   += 3*M;
2482e68589bSMark F. Adams           data[0] = 0.0;             data[M+0] =  coords[3*kk+2]; data[2*M+0] = -coords[3*kk+1];
2492e68589bSMark F. Adams           data[1] = -coords[3*kk+2]; data[M+1] = 0.0;             data[2*M+1] =  coords[3*kk];
2502e68589bSMark F. Adams           data[2] =  coords[3*kk+1]; data[M+2] = -coords[3*kk];   data[2*M+2] = 0.0;
2512e68589bSMark F. Adams         }
2522e68589bSMark F. Adams       }
2532e68589bSMark F. Adams     }
2542e68589bSMark F. Adams   }
2552e68589bSMark F. Adams 
2562e68589bSMark F. Adams   pc_gamg->data_sz = arrsz;
2572e68589bSMark F. Adams   PetscFunctionReturn(0);
2582e68589bSMark F. Adams }
2592e68589bSMark F. Adams 
260b43b03e9SMark F. Adams typedef PetscInt NState;
261b43b03e9SMark F. Adams static const NState NOT_DONE=-2;
262b43b03e9SMark F. Adams static const NState DELETED =-1;
263b43b03e9SMark F. Adams static const NState REMOVED =-3;
264b43b03e9SMark F. Adams #define IS_SELECTED(s) (s!=DELETED && s!=NOT_DONE && s!=REMOVED)
265b43b03e9SMark F. Adams 
266c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
267c8b0795cSMark F. Adams /*
268b43b03e9SMark F. Adams    smoothAggs - greedy grab of with G1 (unsquared graph) -- AIJ specific
269b43b03e9SMark F. Adams      - AGG-MG specific: clears singletons out of 'selected_2'
270c8b0795cSMark F. Adams 
271c8b0795cSMark F. Adams    Input Parameter:
2722adfe9d3SBarry Smith    . Gmat_2 - glabal matrix of graph (data not defined)   base (squared) graph
2732adfe9d3SBarry Smith    . Gmat_1 - base graph to grab with                 base graph
274c8b0795cSMark F. Adams    Input/Output Parameter:
2750cbbd2e1SMark F. Adams    . aggs_2 - linked list of aggs with gids)
276c8b0795cSMark F. Adams */
277c8b0795cSMark F. Adams #undef __FUNCT__
278c8b0795cSMark F. Adams #define __FUNCT__ "smoothAggs"
279*e3df7ea3SBarry Smith static PetscErrorCode smoothAggs(PC pc,Mat Gmat_2, Mat Gmat_1,PetscCoarsenData *aggs_2)
280c8b0795cSMark F. Adams {
281c8b0795cSMark F. Adams   PetscErrorCode ierr;
282c8b0795cSMark F. Adams   PetscBool      isMPI;
2833fa9c902SMark Adams   Mat_SeqAIJ     *matA_1, *matB_1=0;
2843b4367a7SBarry Smith   MPI_Comm       comm;
2850cbbd2e1SMark F. Adams   PetscInt       lid,*ii,*idx,ix,Iend,my0,kk,n,j;
286c8b0795cSMark F. Adams   Mat_MPIAIJ     *mpimat_2 = 0, *mpimat_1=0;
287c8b0795cSMark F. Adams   const PetscInt nloc      = Gmat_2->rmap->n;
2880cbbd2e1SMark F. Adams   PetscScalar    *cpcol_1_state,*cpcol_2_state,*cpcol_2_par_orig,*lid_parent_gid;
2890cbbd2e1SMark F. Adams   PetscInt       *lid_cprowID_1;
290c8b0795cSMark F. Adams   NState         *lid_state;
2910cbbd2e1SMark F. Adams   Vec            ghost_par_orig2;
292c8b0795cSMark F. Adams 
293c8b0795cSMark F. Adams   PetscFunctionBegin;
2943b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Gmat_2,&comm);CHKERRQ(ierr);
295c8b0795cSMark F. Adams   ierr = MatGetOwnershipRange(Gmat_1,&my0,&Iend);CHKERRQ(ierr);
296c8b0795cSMark F. Adams 
2970cbbd2e1SMark F. Adams   if (PETSC_FALSE) {
298c8b0795cSMark F. Adams     PetscViewer viewer; char fname[32]; static int llev=0;
299c8b0795cSMark F. Adams     sprintf(fname,"Gmat2_%d.m",llev++);
3003b4367a7SBarry Smith     PetscViewerASCIIOpen(comm,fname,&viewer);
3016a9046bcSBarry Smith     ierr = PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
302c8b0795cSMark F. Adams     ierr = MatView(Gmat_2, viewer);CHKERRQ(ierr);
303f159cad9SBarry Smith     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
304c8b0795cSMark F. Adams     ierr = PetscViewerDestroy(&viewer);
305c8b0795cSMark F. Adams   }
306c8b0795cSMark F. Adams 
307c8b0795cSMark F. Adams   /* get submatrices */
308251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)Gmat_1, MATMPIAIJ, &isMPI);CHKERRQ(ierr);
309c8b0795cSMark F. Adams   if (isMPI) {
310c8b0795cSMark F. Adams     /* grab matrix objects */
311c8b0795cSMark F. Adams     mpimat_2 = (Mat_MPIAIJ*)Gmat_2->data;
312c8b0795cSMark F. Adams     mpimat_1 = (Mat_MPIAIJ*)Gmat_1->data;
313c8b0795cSMark F. Adams     matA_1   = (Mat_SeqAIJ*)mpimat_1->A->data;
314c8b0795cSMark F. Adams     matB_1   = (Mat_SeqAIJ*)mpimat_1->B->data;
315c8b0795cSMark F. Adams 
316c8b0795cSMark F. Adams     /* force compressed row storage for B matrix in AuxMat */
31711e456e1SBarry Smith     ierr = MatCheckCompressedRow(mpimat_1->B,matB_1->nonzerorowcnt,&matB_1->compressedrow,matB_1->i,Gmat_1->rmap->n,-1.0);CHKERRQ(ierr);
318c8b0795cSMark F. Adams 
319785e854fSJed Brown     ierr = PetscMalloc1(nloc, &lid_cprowID_1);CHKERRQ(ierr);
3200cbbd2e1SMark F. Adams     for (lid = 0; lid < nloc; lid++) lid_cprowID_1[lid] = -1;
321c8b0795cSMark F. Adams     for (ix=0; ix<matB_1->compressedrow.nrows; ix++) {
322c8b0795cSMark F. Adams       PetscInt lid = matB_1->compressedrow.rindex[ix];
323c8b0795cSMark F. Adams       lid_cprowID_1[lid] = ix;
324c8b0795cSMark F. Adams     }
325806fa848SBarry Smith   } else {
32615687449SMark Adams     PetscBool        isAIJ;
32715687449SMark Adams     ierr = PetscObjectTypeCompare((PetscObject)Gmat_1,MATSEQAIJ,&isAIJ);CHKERRQ(ierr);
32815687449SMark Adams     if (!isAIJ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Require AIJ matrix.");
329c8b0795cSMark F. Adams     matA_1        = (Mat_SeqAIJ*)Gmat_1->data;
3300298fd71SBarry Smith     lid_cprowID_1 = NULL;
331c8b0795cSMark F. Adams   }
33278a438d6SMark Adams   if (nloc>0) {
333359038b3SMark Adams     if (matB_1 && !matB_1->compressedrow.use) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"matB_1 && !matB_1->compressedrow.use: PETSc bug???");
33478a438d6SMark Adams   }
335c8b0795cSMark F. Adams   /* get state of locals and selected gid for deleted */
336e632b94dSBarry Smith   ierr = PetscMalloc2(nloc, &lid_state,nloc, &lid_parent_gid);CHKERRQ(ierr);
337c8b0795cSMark F. Adams   for (lid = 0; lid < nloc; lid++) {
3380cbbd2e1SMark F. Adams     lid_parent_gid[lid] = -1.0;
339c8b0795cSMark F. Adams     lid_state[lid]      = DELETED;
340c8b0795cSMark F. Adams   }
3410cbbd2e1SMark F. Adams 
3420cbbd2e1SMark F. Adams   /* set lid_state */
3430cbbd2e1SMark F. Adams   for (lid = 0; lid < nloc; lid++) {
344539c167fSBarry Smith     PetscCDIntNd *pos;
345e78576d6SMark F. Adams     ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr);
346e78576d6SMark F. Adams     if (pos) {
347e78576d6SMark F. Adams       PetscInt gid1;
3482fa5cd67SKarl Rupp 
349484f0a72SBarry Smith       ierr = PetscCDIntNdGetID(pos, &gid1);CHKERRQ(ierr);
35071959b99SBarry Smith       if (gid1 != lid+my0) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gid1 %D != lid %D + my0 %D",gid1,lid,my0);
3510cbbd2e1SMark F. Adams       lid_state[lid] = gid1;
352b43b03e9SMark F. Adams     }
353b43b03e9SMark F. Adams   }
3540cbbd2e1SMark F. Adams 
3550cbbd2e1SMark F. Adams   /* map local to selected local, DELETED means a ghost owns it */
356c8b0795cSMark F. Adams   for (lid=kk=0; lid<nloc; lid++) {
357c8b0795cSMark F. Adams     NState state = lid_state[lid];
358c8b0795cSMark F. Adams     if (IS_SELECTED(state)) {
359539c167fSBarry Smith       PetscCDIntNd *pos;
360e78576d6SMark F. Adams       ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr);
361e78576d6SMark F. Adams       while (pos) {
362e78576d6SMark F. Adams         PetscInt gid1;
363484f0a72SBarry Smith         ierr = PetscCDIntNdGetID(pos, &gid1);CHKERRQ(ierr);
364e78576d6SMark F. Adams         ierr = PetscCDGetNextPos(aggs_2,lid,&pos);CHKERRQ(ierr);
365e78576d6SMark F. Adams 
3662fa5cd67SKarl Rupp         if (gid1 >= my0 && gid1 < Iend) lid_parent_gid[gid1-my0] = (PetscScalar)(lid + my0);
367c8b0795cSMark F. Adams       }
3680cbbd2e1SMark F. Adams     }
3690cbbd2e1SMark F. Adams   }
3700cbbd2e1SMark F. Adams   /* get 'cpcol_1/2_state' & cpcol_2_par_orig - uses mpimat_1/2->lvec for temp space */
371c8b0795cSMark F. Adams   if (isMPI) {
372c8b0795cSMark F. Adams     Vec tempVec;
373c8b0795cSMark F. Adams     /* get 'cpcol_1_state' */
3742a7a6963SBarry Smith     ierr = MatCreateVecs(Gmat_1, &tempVec, 0);CHKERRQ(ierr);
375c8b0795cSMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
376c8b0795cSMark F. Adams       PetscScalar v = (PetscScalar)lid_state[kk];
377c8b0795cSMark F. Adams       ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr);
378c8b0795cSMark F. Adams     }
379c8b0795cSMark F. Adams     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
380c8b0795cSMark F. Adams     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
381806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
382806fa848SBarry Smith     ierr = VecScatterEnd(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
383c8b0795cSMark F. Adams     ierr = VecGetArray(mpimat_1->lvec, &cpcol_1_state);CHKERRQ(ierr);
384c8b0795cSMark F. Adams     /* get 'cpcol_2_state' */
385806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
386806fa848SBarry Smith     ierr =   VecScatterEnd(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
387c8b0795cSMark F. Adams     ierr = VecGetArray(mpimat_2->lvec, &cpcol_2_state);CHKERRQ(ierr);
3880cbbd2e1SMark F. Adams     /* get 'cpcol_2_par_orig' */
3890cbbd2e1SMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
3900cbbd2e1SMark F. Adams       PetscScalar v = (PetscScalar)lid_parent_gid[kk];
3910cbbd2e1SMark F. Adams       ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr);
3920cbbd2e1SMark F. Adams     }
3930cbbd2e1SMark F. Adams     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
3940cbbd2e1SMark F. Adams     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
3950cbbd2e1SMark F. Adams     ierr = VecDuplicate(mpimat_2->lvec, &ghost_par_orig2);CHKERRQ(ierr);
396806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
397806fa848SBarry Smith     ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
3980cbbd2e1SMark F. Adams     ierr = VecGetArray(ghost_par_orig2, &cpcol_2_par_orig);CHKERRQ(ierr);
3990cbbd2e1SMark F. Adams 
400c8b0795cSMark F. Adams     ierr = VecDestroy(&tempVec);CHKERRQ(ierr);
401c8b0795cSMark F. Adams   } /* ismpi */
402c8b0795cSMark F. Adams 
403c8b0795cSMark F. Adams   /* doit */
404c8b0795cSMark F. Adams   for (lid=0; lid<nloc; lid++) {
405c8b0795cSMark F. Adams     NState state = lid_state[lid];
4060cbbd2e1SMark F. Adams     if (IS_SELECTED(state)) {
4070cbbd2e1SMark F. Adams       /* steal locals */
408c8b0795cSMark F. Adams       ii  = matA_1->i; n = ii[lid+1] - ii[lid];
409c8b0795cSMark F. Adams       idx = matA_1->j + ii[lid];
410c8b0795cSMark F. Adams       for (j=0; j<n; j++) {
4110cbbd2e1SMark F. Adams         PetscInt lidj   = idx[j], sgid;
412c8b0795cSMark F. Adams         NState   statej = lid_state[lidj];
4130cbbd2e1SMark F. Adams         if (statej==DELETED && (sgid=(PetscInt)PetscRealPart(lid_parent_gid[lidj])) != lid+my0) { /* steal local */
4140cbbd2e1SMark F. Adams           lid_parent_gid[lidj] = (PetscScalar)(lid+my0); /* send this if sgid is not local */
4150cbbd2e1SMark F. Adams           if (sgid >= my0 && sgid < Iend) {       /* I'm stealing this local from a local sgid */
4160cbbd2e1SMark F. Adams             PetscInt   hav=0,slid=sgid-my0,gidj=lidj+my0;
417539c167fSBarry Smith             PetscCDIntNd *pos,*last=NULL;
418c8b0795cSMark F. Adams             /* looking for local from local so id_llist_2 works */
419e78576d6SMark F. Adams             ierr = PetscCDGetHeadPos(aggs_2,slid,&pos);CHKERRQ(ierr);
420e78576d6SMark F. Adams             while (pos) {
421e78576d6SMark F. Adams               PetscInt gid;
422484f0a72SBarry Smith               ierr = PetscCDIntNdGetID(pos, &gid);CHKERRQ(ierr);
4230cbbd2e1SMark F. Adams               if (gid == gidj) {
42471959b99SBarry Smith                 if (!last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"last cannot be null");
42541b27cdeSMark F. Adams                 ierr = PetscCDRemoveNextNode(aggs_2, slid, last);CHKERRQ(ierr);
42641b27cdeSMark F. Adams                 ierr = PetscCDAppendNode(aggs_2, lid, pos);CHKERRQ(ierr);
4270cbbd2e1SMark F. Adams                 hav  = 1;
428c8b0795cSMark F. Adams                 break;
429806fa848SBarry Smith               } else last = pos;
430e78576d6SMark F. Adams 
431e78576d6SMark F. Adams               ierr = PetscCDGetNextPos(aggs_2,slid,&pos);CHKERRQ(ierr);
432c8b0795cSMark F. Adams             }
433c8b0795cSMark F. Adams             if (hav!=1) {
43471959b99SBarry Smith               if (!hav) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix");
435c666adbfSMark F. Adams               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %d times???",hav);
436c8b0795cSMark F. Adams             }
437806fa848SBarry Smith           } else {            /* I'm stealing this local, owned by a ghost */
438*e3df7ea3SBarry Smith             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.0' if the matrix is structurally symmetric.",((PetscObject)pc)->prefix,((PetscObject)pc)->prefix);
43941b27cdeSMark F. Adams             ierr = PetscCDAppendID(aggs_2, lid, lidj+my0);CHKERRQ(ierr);
440c8b0795cSMark F. Adams           }
441c8b0795cSMark F. Adams         }
4420cbbd2e1SMark F. Adams       } /* local neighbors */
443806fa848SBarry Smith     } else if (state == DELETED && lid_cprowID_1) {
4440cbbd2e1SMark F. Adams       PetscInt sgidold = (PetscInt)PetscRealPart(lid_parent_gid[lid]);
445c8b0795cSMark F. Adams       /* see if I have a selected ghost neighbor that will steal me */
446c8b0795cSMark F. Adams       if ((ix=lid_cprowID_1[lid]) != -1) {
447c8b0795cSMark F. Adams         ii  = matB_1->compressedrow.i; n = ii[ix+1] - ii[ix];
448c8b0795cSMark F. Adams         idx = matB_1->j + ii[ix];
449c8b0795cSMark F. Adams         for (j=0; j<n; j++) {
450c8b0795cSMark F. Adams           PetscInt cpid   = idx[j];
451c8b0795cSMark F. Adams           NState   statej = (NState)PetscRealPart(cpcol_1_state[cpid]);
452c8b0795cSMark F. Adams           if (IS_SELECTED(statej) && sgidold != (PetscInt)statej) { /* ghost will steal this, remove from my list */
4530cbbd2e1SMark F. Adams             lid_parent_gid[lid] = (PetscScalar)statej; /* send who selected */
4540cbbd2e1SMark F. Adams             if (sgidold>=my0 && sgidold<Iend) { /* this was mine */
4550cbbd2e1SMark F. Adams               PetscInt   hav=0,oldslidj=sgidold-my0;
456539c167fSBarry Smith               PetscCDIntNd *pos,*last=NULL;
4570cbbd2e1SMark F. Adams               /* remove from 'oldslidj' list */
458e78576d6SMark F. Adams               ierr = PetscCDGetHeadPos(aggs_2,oldslidj,&pos);CHKERRQ(ierr);
459e78576d6SMark F. Adams               while (pos) {
460e78576d6SMark F. Adams                 PetscInt gid;
461484f0a72SBarry Smith                 ierr = PetscCDIntNdGetID(pos, &gid);CHKERRQ(ierr);
4620cbbd2e1SMark F. Adams                 if (lid+my0 == gid) {
4630cbbd2e1SMark F. Adams                   /* id_llist_2[lastid] = id_llist_2[flid];   /\* remove lid from oldslidj list *\/ */
46471959b99SBarry Smith                   if (!last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"last cannot be null");
46541b27cdeSMark F. Adams                   ierr = PetscCDRemoveNextNode(aggs_2, oldslidj, last);CHKERRQ(ierr);
4660cbbd2e1SMark F. Adams                   /* ghost (PetscScalar)statej will add this later */
4670cbbd2e1SMark F. Adams                   hav = 1;
468c8b0795cSMark F. Adams                   break;
469806fa848SBarry Smith                 } else last = pos;
470e78576d6SMark F. Adams 
471e78576d6SMark F. Adams                 ierr = PetscCDGetNextPos(aggs_2,oldslidj,&pos);CHKERRQ(ierr);
472c8b0795cSMark F. Adams               }
473c8b0795cSMark F. Adams               if (hav!=1) {
4747f66b68fSMark Adams                 if (!hav) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix");
475c666adbfSMark F. Adams                 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %d times???",hav);
476c8b0795cSMark F. Adams               }
477806fa848SBarry Smith             } else {
4780cbbd2e1SMark F. Adams               /* ghosts remove this later */
4790cbbd2e1SMark F. Adams             }
480c8b0795cSMark F. Adams           }
481c8b0795cSMark F. Adams         }
482c8b0795cSMark F. Adams       }
483c8b0795cSMark F. Adams     } /* selected/deleted */
484c8b0795cSMark F. Adams   } /* node loop */
485c8b0795cSMark F. Adams 
486c8b0795cSMark F. Adams   if (isMPI) {
4870cbbd2e1SMark F. Adams     PetscScalar     *cpcol_2_parent,*cpcol_2_gid;
4880cbbd2e1SMark F. Adams     Vec             tempVec,ghostgids2,ghostparents2;
4890cbbd2e1SMark F. Adams     PetscInt        cpid,nghost_2;
4901943db53SBarry Smith     PCGAMGHashTable gid_cpid;
491c8b0795cSMark F. Adams 
4920cbbd2e1SMark F. Adams     ierr = VecGetSize(mpimat_2->lvec, &nghost_2);CHKERRQ(ierr);
4932a7a6963SBarry Smith     ierr = MatCreateVecs(Gmat_2, &tempVec, 0);CHKERRQ(ierr);
4940cbbd2e1SMark F. Adams 
4950cbbd2e1SMark F. Adams     /* get 'cpcol_2_parent' */
496c8b0795cSMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
4970cbbd2e1SMark F. Adams       ierr = VecSetValues(tempVec, 1, &j, &lid_parent_gid[kk], INSERT_VALUES);CHKERRQ(ierr);
498c8b0795cSMark F. Adams     }
499c8b0795cSMark F. Adams     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
500c8b0795cSMark F. Adams     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
5010cbbd2e1SMark F. Adams     ierr = VecDuplicate(mpimat_2->lvec, &ghostparents2);CHKERRQ(ierr);
502806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
503806fa848SBarry Smith     ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5040cbbd2e1SMark F. Adams     ierr = VecGetArray(ghostparents2, &cpcol_2_parent);CHKERRQ(ierr);
5050cbbd2e1SMark F. Adams 
5060cbbd2e1SMark F. Adams     /* get 'cpcol_2_gid' */
5070cbbd2e1SMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
5080cbbd2e1SMark F. Adams       PetscScalar v = (PetscScalar)j;
5090cbbd2e1SMark F. Adams       ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr);
5100cbbd2e1SMark F. Adams     }
5110cbbd2e1SMark F. Adams     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
5120cbbd2e1SMark F. Adams     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
5130cbbd2e1SMark F. Adams     ierr = VecDuplicate(mpimat_2->lvec, &ghostgids2);CHKERRQ(ierr);
514806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
515806fa848SBarry Smith     ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5160cbbd2e1SMark F. Adams     ierr = VecGetArray(ghostgids2, &cpcol_2_gid);CHKERRQ(ierr);
5170cbbd2e1SMark F. Adams 
518c8b0795cSMark F. Adams     ierr = VecDestroy(&tempVec);CHKERRQ(ierr);
519c8b0795cSMark F. Adams 
5200cbbd2e1SMark F. Adams     /* look for deleted ghosts and add to table */
5211943db53SBarry Smith     ierr = PCGAMGHashTableCreate(2*nghost_2+1, &gid_cpid);CHKERRQ(ierr);
5220cbbd2e1SMark F. Adams     for (cpid = 0; cpid < nghost_2; cpid++) {
5230cbbd2e1SMark F. Adams       NState state = (NState)PetscRealPart(cpcol_2_state[cpid]);
5240cbbd2e1SMark F. Adams       if (state==DELETED) {
5250cbbd2e1SMark F. Adams         PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]);
5260cbbd2e1SMark F. Adams         PetscInt sgid_old = (PetscInt)PetscRealPart(cpcol_2_par_orig[cpid]);
5270cbbd2e1SMark F. Adams         if (sgid_old == -1 && sgid_new != -1) {
5280cbbd2e1SMark F. Adams           PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]);
5291943db53SBarry Smith           ierr = PCGAMGHashTableAdd(&gid_cpid, gid, cpid);CHKERRQ(ierr);
5300cbbd2e1SMark F. Adams         }
5310cbbd2e1SMark F. Adams       }
5320cbbd2e1SMark F. Adams     }
533c8b0795cSMark F. Adams 
5340cbbd2e1SMark F. Adams     /* look for deleted ghosts and see if they moved - remove it */
535c8b0795cSMark F. Adams     for (lid=0; lid<nloc; lid++) {
536c8b0795cSMark F. Adams       NState state = lid_state[lid];
537c8b0795cSMark F. Adams       if (IS_SELECTED(state)) {
538539c167fSBarry Smith         PetscCDIntNd *pos,*last=NULL;
539c8b0795cSMark F. Adams         /* look for deleted ghosts and see if they moved */
540e78576d6SMark F. Adams         ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr);
541e78576d6SMark F. Adams         while (pos) {
542e78576d6SMark F. Adams           PetscInt gid;
543484f0a72SBarry Smith           ierr = PetscCDIntNdGetID(pos, &gid);CHKERRQ(ierr);
544e78576d6SMark F. Adams 
5450cbbd2e1SMark F. Adams           if (gid < my0 || gid >= Iend) {
5461943db53SBarry Smith             ierr = PCGAMGHashTableFind(&gid_cpid, gid, &cpid);CHKERRQ(ierr);
5470cbbd2e1SMark F. Adams             if (cpid != -1) {
5480cbbd2e1SMark F. Adams               /* a moved ghost - */
5490cbbd2e1SMark F. Adams               /* id_llist_2[lastid] = id_llist_2[flid];    /\* remove 'flid' from list *\/ */
55041b27cdeSMark F. Adams               ierr = PetscCDRemoveNextNode(aggs_2, lid, last);CHKERRQ(ierr);
551806fa848SBarry Smith             } else last = pos;
552806fa848SBarry Smith           } else last = pos;
553e78576d6SMark F. Adams 
554e78576d6SMark F. Adams           ierr = PetscCDGetNextPos(aggs_2,lid,&pos);CHKERRQ(ierr);
555c8b0795cSMark F. Adams         } /* loop over list of deleted */
556c8b0795cSMark F. Adams       } /* selected */
557c8b0795cSMark F. Adams     }
5581943db53SBarry Smith     ierr = PCGAMGHashTableDestroy(&gid_cpid);CHKERRQ(ierr);
559c8b0795cSMark F. Adams 
5600cbbd2e1SMark F. Adams     /* look at ghosts, see if they changed - and it */
5610cbbd2e1SMark F. Adams     for (cpid = 0; cpid < nghost_2; cpid++) {
5620cbbd2e1SMark F. Adams       PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]);
5630cbbd2e1SMark F. Adams       if (sgid_new >= my0 && sgid_new < Iend) { /* this is mine */
5640cbbd2e1SMark F. Adams         PetscInt     gid     = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]);
5650cbbd2e1SMark F. Adams         PetscInt     slid_new=sgid_new-my0,hav=0;
566539c167fSBarry Smith         PetscCDIntNd *pos;
567539c167fSBarry Smith 
5680cbbd2e1SMark F. Adams         /* search for this gid to see if I have it */
569e78576d6SMark F. Adams         ierr = PetscCDGetHeadPos(aggs_2,slid_new,&pos);CHKERRQ(ierr);
570e78576d6SMark F. Adams         while (pos) {
571e78576d6SMark F. Adams           PetscInt gidj;
572484f0a72SBarry Smith           ierr = PetscCDIntNdGetID(pos, &gidj);CHKERRQ(ierr);
573e78576d6SMark F. Adams           ierr = PetscCDGetNextPos(aggs_2,slid_new,&pos);CHKERRQ(ierr);
574e78576d6SMark F. Adams 
5750cbbd2e1SMark F. Adams           if (gidj == gid) { hav = 1; break; }
576c8b0795cSMark F. Adams         }
577c8b0795cSMark F. Adams         if (hav != 1) {
578ffc955d6SMark F. Adams           /* insert 'flidj' into head of llist */
57941b27cdeSMark F. Adams           ierr = PetscCDAppendID(aggs_2, slid_new, gid);CHKERRQ(ierr);
580c8b0795cSMark F. Adams         }
581c8b0795cSMark F. Adams       }
582c8b0795cSMark F. Adams     }
583c8b0795cSMark F. Adams 
5840cbbd2e1SMark F. Adams     ierr = VecRestoreArray(mpimat_1->lvec, &cpcol_1_state);CHKERRQ(ierr);
5850cbbd2e1SMark F. Adams     ierr = VecRestoreArray(mpimat_2->lvec, &cpcol_2_state);CHKERRQ(ierr);
5860cbbd2e1SMark F. Adams     ierr = VecRestoreArray(ghostparents2, &cpcol_2_parent);CHKERRQ(ierr);
5870cbbd2e1SMark F. Adams     ierr = VecRestoreArray(ghostgids2, &cpcol_2_gid);CHKERRQ(ierr);
588c8b0795cSMark F. Adams     ierr = PetscFree(lid_cprowID_1);CHKERRQ(ierr);
5890cbbd2e1SMark F. Adams     ierr = VecDestroy(&ghostgids2);CHKERRQ(ierr);
5900cbbd2e1SMark F. Adams     ierr = VecDestroy(&ghostparents2);CHKERRQ(ierr);
5910cbbd2e1SMark F. Adams     ierr = VecDestroy(&ghost_par_orig2);CHKERRQ(ierr);
592c8b0795cSMark F. Adams   }
593c8b0795cSMark F. Adams 
594e632b94dSBarry Smith   ierr = PetscFree2(lid_state,lid_parent_gid);CHKERRQ(ierr);
595c8b0795cSMark F. Adams   PetscFunctionReturn(0);
596c8b0795cSMark F. Adams }
5972e68589bSMark F. Adams 
5982e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
5992e68589bSMark F. Adams /*
600a2f3521dSMark F. Adams    PCSetData_AGG - called if data is not set with PCSetCoordinates.
601a2f3521dSMark F. Adams       Looks in Mat for near null space.
602a2f3521dSMark F. Adams       Does not work for Stokes
6032e68589bSMark F. Adams 
6042e68589bSMark F. Adams   Input Parameter:
6052e68589bSMark F. Adams    . pc -
606a2f3521dSMark F. Adams    . a_A - matrix to get (near) null space out of.
6072e68589bSMark F. Adams */
6082e68589bSMark F. Adams #undef __FUNCT__
6092e68589bSMark F. Adams #define __FUNCT__ "PCSetData_AGG"
610e0877f53SBarry Smith static PetscErrorCode PCSetData_AGG(PC pc, Mat a_A)
6112e68589bSMark F. Adams {
6122e68589bSMark F. Adams   PetscErrorCode ierr;
613b8cd405aSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
614b8cd405aSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
615b8cd405aSMark F. Adams   MatNullSpace   mnull;
6162e68589bSMark F. Adams   PetscFunctionBegin;
61766f2ef4bSMark Adams 
618b8cd405aSMark F. Adams   ierr = MatGetNearNullSpace(a_A, &mnull);CHKERRQ(ierr);
619b8cd405aSMark F. Adams   if (!mnull) {
62066f2ef4bSMark Adams     DM dm;
62166f2ef4bSMark Adams     ierr = PCGetDM(pc, &dm);CHKERRQ(ierr);
62266f2ef4bSMark Adams     if (!dm) {
62366f2ef4bSMark Adams       ierr = MatGetDM(a_A, &dm);CHKERRQ(ierr);
62466f2ef4bSMark Adams     }
62566f2ef4bSMark Adams     if (dm) {
62666f2ef4bSMark Adams       PetscObject deformation;
627b0d30dd6SMatthew G. Knepley       PetscInt    Nf;
628b0d30dd6SMatthew G. Knepley 
629b0d30dd6SMatthew G. Knepley       ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr);
630b0d30dd6SMatthew G. Knepley       if (Nf) {
63166f2ef4bSMark Adams         ierr = DMGetField(dm, 0, &deformation);CHKERRQ(ierr);
63266f2ef4bSMark Adams         ierr = PetscObjectQuery((PetscObject)deformation,"nearnullspace",(PetscObject*)&mnull);CHKERRQ(ierr);
63366f2ef4bSMark Adams         if (!mnull) {
63466f2ef4bSMark Adams           ierr = PetscObjectQuery((PetscObject)deformation,"nullspace",(PetscObject*)&mnull);CHKERRQ(ierr);
63566f2ef4bSMark Adams         }
63666f2ef4bSMark Adams       }
63766f2ef4bSMark Adams     }
638b0d30dd6SMatthew G. Knepley   }
63966f2ef4bSMark Adams 
64066f2ef4bSMark Adams   if (!mnull) {
641a2f3521dSMark F. Adams     PetscInt bs,NN,MM;
6429d1b15c3SMark F. Adams     ierr = MatGetBlockSize(a_A, &bs);CHKERRQ(ierr);
64371959b99SBarry Smith     ierr = MatGetLocalSize(a_A, &MM, &NN);CHKERRQ(ierr);
64471959b99SBarry Smith     if (MM % bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MM %D must be divisible by bs %D",MM,bs);
6450298fd71SBarry Smith     ierr = PCSetCoordinates_AGG(pc, bs, MM/bs, NULL);CHKERRQ(ierr);
646806fa848SBarry Smith   } else {
647b8cd405aSMark F. Adams     PetscReal *nullvec;
648b8cd405aSMark F. Adams     PetscBool has_const;
649b8cd405aSMark F. Adams     PetscInt  i,j,mlocal,nvec,bs;
650b8cd405aSMark F. Adams     const Vec *vecs; const PetscScalar *v;
6510298fd71SBarry Smith     ierr = MatGetLocalSize(a_A,&mlocal,NULL);CHKERRQ(ierr);
652b8cd405aSMark F. Adams     ierr = MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs);CHKERRQ(ierr);
653a0dea87dSMark Adams     pc_gamg->data_sz = (nvec+!!has_const)*mlocal;
654785e854fSJed Brown     ierr = PetscMalloc1((nvec+!!has_const)*mlocal,&nullvec);CHKERRQ(ierr);
655b8cd405aSMark F. Adams     if (has_const) for (i=0; i<mlocal; i++) nullvec[i] = 1.0;
656b8cd405aSMark F. Adams     for (i=0; i<nvec; i++) {
657b8cd405aSMark F. Adams       ierr = VecGetArrayRead(vecs[i],&v);CHKERRQ(ierr);
658b8cd405aSMark F. Adams       for (j=0; j<mlocal; j++) nullvec[(i+!!has_const)*mlocal + j] = PetscRealPart(v[j]);
659b8cd405aSMark F. Adams       ierr = VecRestoreArrayRead(vecs[i],&v);CHKERRQ(ierr);
660b8cd405aSMark F. Adams     }
661b8cd405aSMark F. Adams     pc_gamg->data           = nullvec;
662b8cd405aSMark F. Adams     pc_gamg->data_cell_cols = (nvec+!!has_const);
6632fa5cd67SKarl Rupp 
66406e133e7SMark Adams     ierr = MatGetBlockSize(a_A, &bs);CHKERRQ(ierr);
6652fa5cd67SKarl Rupp 
666b8cd405aSMark F. Adams     pc_gamg->data_cell_rows = bs;
667b8cd405aSMark F. Adams   }
6682e68589bSMark F. Adams   PetscFunctionReturn(0);
6692e68589bSMark F. Adams }
6702e68589bSMark F. Adams 
6712e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
6722e68589bSMark F. Adams /*
6732e68589bSMark F. Adams  formProl0
6742e68589bSMark F. Adams 
6752e68589bSMark F. Adams    Input Parameter:
6762adfe9d3SBarry Smith    . agg_llists - list of arrays with aggregates -- list from selected vertices of aggregate unselected vertices
6772adfe9d3SBarry Smith    . bs - row block size
6780cbbd2e1SMark F. Adams    . nSAvec - column bs of new P
6790cbbd2e1SMark F. Adams    . my0crs - global index of start of locals
6802adfe9d3SBarry Smith    . data_stride - bs*(nloc nodes + ghost nodes) [data_stride][nSAvec]
6812e68589bSMark F. Adams    . data_in[data_stride*nSAvec] - local data on fine grid
6822e68589bSMark F. Adams    . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist'
6832e68589bSMark F. Adams   Output Parameter:
6842e68589bSMark F. Adams    . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data
6852e68589bSMark F. Adams    . a_Prol - prolongation operator
6862e68589bSMark F. Adams */
6872e68589bSMark F. Adams #undef __FUNCT__
6882e68589bSMark F. Adams #define __FUNCT__ "formProl0"
6892adfe9d3SBarry 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)
6902e68589bSMark F. Adams {
6912e68589bSMark F. Adams   PetscErrorCode  ierr;
692ac187d40SBarry Smith   PetscInt        Istart,my0,Iend,nloc,clid,flid = 0,aggID,kk,jj,ii,mm,ndone,nSelected,minsz,nghosts,out_data_stride;
6933b4367a7SBarry Smith   MPI_Comm        comm;
69473911c69SBarry Smith   PetscMPIInt     rank;
6952e68589bSMark F. Adams   PetscReal       *out_data;
696539c167fSBarry Smith   PetscCDIntNd    *pos;
6971943db53SBarry Smith   PCGAMGHashTable fgid_flid;
6980cbbd2e1SMark F. Adams 
699797e13b7SMark F. Adams /* #define OUT_AGGS */
700519f805aSKarl Rupp #if defined(OUT_AGGS)
7010298fd71SBarry Smith   static PetscInt llev = 0; char fname[32]; FILE *file = NULL; PetscInt pM;
7029057884aSMark F. Adams #endif
7032e68589bSMark F. Adams 
7042e68589bSMark F. Adams   PetscFunctionBegin;
7053b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)a_Prol,&comm);CHKERRQ(ierr);
7063b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
7072e68589bSMark F. Adams   ierr = MatGetOwnershipRange(a_Prol, &Istart, &Iend);CHKERRQ(ierr);
70871959b99SBarry Smith   nloc = (Iend-Istart)/bs; my0 = Istart/bs;
70971959b99SBarry 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);
7100cbbd2e1SMark F. Adams   Iend   /= bs;
7110cbbd2e1SMark F. Adams   nghosts = data_stride/bs - nloc;
7122e68589bSMark F. Adams 
7131943db53SBarry Smith   ierr = PCGAMGHashTableCreate(2*nghosts+1, &fgid_flid);CHKERRQ(ierr);
7140cbbd2e1SMark F. Adams   for (kk=0; kk<nghosts; kk++) {
7151943db53SBarry Smith     ierr = PCGAMGHashTableAdd(&fgid_flid, flid_fgid[nloc+kk], nloc+kk);CHKERRQ(ierr);
7162e68589bSMark F. Adams   }
7172e68589bSMark F. Adams 
718519f805aSKarl Rupp #if defined(OUT_AGGS)
719c5df96a5SBarry Smith   sprintf(fname,"aggs_%d_%d.m",llev++,rank);
7202fa5cd67SKarl Rupp   if (llev==1) file = fopen(fname,"w");
7210cbbd2e1SMark F. Adams   MatGetSize(a_Prol, &pM, &jj);
7220cbbd2e1SMark F. Adams #endif
7230cbbd2e1SMark F. Adams 
7240cbbd2e1SMark F. Adams   /* count selected -- same as number of cols of P */
7250cbbd2e1SMark F. Adams   for (nSelected=mm=0; mm<nloc; mm++) {
726e78576d6SMark F. Adams     PetscBool ise;
727e78576d6SMark F. Adams     ierr = PetscCDEmptyAt(agg_llists, mm, &ise);CHKERRQ(ierr);
728e78576d6SMark F. Adams     if (!ise) nSelected++;
7290cbbd2e1SMark F. Adams   }
7300cbbd2e1SMark F. Adams   ierr = MatGetOwnershipRangeColumn(a_Prol, &ii, &jj);CHKERRQ(ierr);
73171959b99SBarry Smith   if ((ii/nSAvec) != my0crs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D /nSAvec %D  != my0crs %D",ii,nSAvec,my0crs);
73271959b99SBarry 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);
7330cbbd2e1SMark F. Adams 
7342e68589bSMark F. Adams   /* aloc space for coarse point data (output) */
7350cbbd2e1SMark F. Adams   out_data_stride = nSelected*nSAvec;
7362fa5cd67SKarl Rupp 
737785e854fSJed Brown   ierr = PetscMalloc1(out_data_stride*nSAvec, &out_data);CHKERRQ(ierr);
73833812677SSatish Balay   for (ii=0;ii<out_data_stride*nSAvec;ii++) out_data[ii]=PETSC_MAX_REAL;
7390cbbd2e1SMark F. Adams   *a_data_out = out_data; /* output - stride nSelected*nSAvec */
7402e68589bSMark F. Adams 
7412e68589bSMark F. Adams   /* find points and set prolongation */
742c8b0795cSMark F. Adams   minsz = 100;
7432e68589bSMark F. Adams   ndone = 0;
7440cbbd2e1SMark F. Adams   for (mm = clid = 0; mm < nloc; mm++) {
745e78576d6SMark F. Adams     ierr = PetscCDSizeAt(agg_llists, mm, &jj);CHKERRQ(ierr);
746e78576d6SMark F. Adams     if (jj > 0) {
7470cbbd2e1SMark F. Adams       const PetscInt lid = mm, cgid = my0crs + clid;
7480cbbd2e1SMark F. Adams       PetscInt       cids[100]; /* max bs */
7490cbbd2e1SMark F. Adams       PetscBLASInt   asz  =jj,M=asz*bs,N=nSAvec,INFO;
7502e68589bSMark F. Adams       PetscBLASInt   Mdata=M+((N-M>0) ? N-M : 0),LDA=Mdata,LWORK=N*bs;
7512e68589bSMark F. Adams       PetscScalar    *qqc,*qqr,*TAU,*WORK;
7522e68589bSMark F. Adams       PetscInt       *fids;
75365d7b583SSatish Balay       PetscReal      *data;
7540cbbd2e1SMark F. Adams       /* count agg */
7550cbbd2e1SMark F. Adams       if (asz<minsz) minsz = asz;
7560cbbd2e1SMark F. Adams 
7570cbbd2e1SMark F. Adams       /* get block */
758e632b94dSBarry Smith       ierr = PetscMalloc5(Mdata*N, &qqc,M*N, &qqr,N, &TAU,LWORK, &WORK,M, &fids);CHKERRQ(ierr);
7592e68589bSMark F. Adams 
7602e68589bSMark F. Adams       aggID = 0;
761e78576d6SMark F. Adams       ierr  = PetscCDGetHeadPos(agg_llists,lid,&pos);CHKERRQ(ierr);
762e78576d6SMark F. Adams       while (pos) {
763e78576d6SMark F. Adams         PetscInt gid1;
764484f0a72SBarry Smith         ierr = PetscCDIntNdGetID(pos, &gid1);CHKERRQ(ierr);
765e78576d6SMark F. Adams         ierr = PetscCDGetNextPos(agg_llists,lid,&pos);CHKERRQ(ierr);
766e78576d6SMark F. Adams 
7670cbbd2e1SMark F. Adams         if (gid1 >= my0 && gid1 < Iend) flid = gid1 - my0;
7680cbbd2e1SMark F. Adams         else {
7691943db53SBarry Smith           ierr = PCGAMGHashTableFind(&fgid_flid, gid1, &flid);CHKERRQ(ierr);
77071959b99SBarry Smith           if (flid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot find gid1 in table");
7710cbbd2e1SMark F. Adams         }
7722e68589bSMark F. Adams         /* copy in B_i matrix - column oriented */
77365d7b583SSatish Balay         data = &data_in[flid*bs];
7743d3eaba7SBarry Smith         for (ii = 0; ii < bs; ii++) {
7752e68589bSMark F. Adams           for (jj = 0; jj < N; jj++) {
7760cbbd2e1SMark F. Adams             PetscReal d = data[jj*data_stride + ii];
7770cbbd2e1SMark F. Adams             qqc[jj*Mdata + aggID*bs + ii] = d;
7782e68589bSMark F. Adams           }
7792e68589bSMark F. Adams         }
780519f805aSKarl Rupp #if defined(OUT_AGGS)
781b2a4f308SMark F. Adams         if (llev==1) {
7829057884aSMark F. Adams           char     str[] = "plot(%e,%e,'r*'), hold on,\n", col[] = "rgbkmc", sim[] = "*os+h>d<vx^";
7830cbbd2e1SMark F. Adams           PetscInt MM,pi,pj;
784c5df96a5SBarry Smith           str[12] = col[(clid+17*rank)%6]; str[13] = sim[(clid+17*rank)%11];
785f7620de1SMatthew G Knepley           MM      = (PetscInt)(PetscSqrtReal((PetscReal)pM));
7860cbbd2e1SMark F. Adams           pj      = gid1/MM; pi = gid1%MM;
787b2a4f308SMark F. Adams           fprintf(file,str,(double)pi,(double)pj);
788b2a4f308SMark F. Adams           /* fprintf(file,str,data[2*data_stride+1],-data[2*data_stride]); */
7899057884aSMark F. Adams         }
7909057884aSMark F. Adams #endif
7912e68589bSMark F. Adams         /* set fine IDs */
7922e68589bSMark F. Adams         for (kk=0; kk<bs; kk++) fids[aggID*bs + kk] = flid_fgid[flid]*bs + kk;
7932e68589bSMark F. Adams 
7942e68589bSMark F. Adams         aggID++;
7950cbbd2e1SMark F. Adams       }
7962e68589bSMark F. Adams 
7972e68589bSMark F. Adams       /* pad with zeros */
7982e68589bSMark F. Adams       for (ii = asz*bs; ii < Mdata; ii++) {
7992e68589bSMark F. Adams         for (jj = 0; jj < N; jj++, kk++) {
8002e68589bSMark F. Adams           qqc[jj*Mdata + ii] = .0;
8012e68589bSMark F. Adams         }
8022e68589bSMark F. Adams       }
8032e68589bSMark F. Adams 
8042e68589bSMark F. Adams       ndone += aggID;
8052e68589bSMark F. Adams       /* QR */
80684df3f90SSatish Balay       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
8078b83055fSJed Brown       PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
80884df3f90SSatish Balay       ierr = PetscFPTrapPop();CHKERRQ(ierr);
809d23427c4SJed Brown       if (INFO != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xGEQRF error");
8102e68589bSMark F. Adams       /* get R - column oriented - output B_{i+1} */
8112e68589bSMark F. Adams       {
8122e68589bSMark F. Adams         PetscReal *data = &out_data[clid*nSAvec];
8132e68589bSMark F. Adams         for (jj = 0; jj < nSAvec; jj++) {
8142e68589bSMark F. Adams           for (ii = 0; ii < nSAvec; ii++) {
81533812677SSatish 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);
8160cbbd2e1SMark F. Adams            if (ii <= jj) data[jj*out_data_stride + ii] = PetscRealPart(qqc[jj*Mdata + ii]);
8170cbbd2e1SMark F. Adams            else data[jj*out_data_stride + ii] = 0.;
8182e68589bSMark F. Adams           }
8192e68589bSMark F. Adams         }
8202e68589bSMark F. Adams       }
8212e68589bSMark F. Adams 
8222e68589bSMark F. Adams       /* get Q - row oriented */
8238b83055fSJed Brown       PetscStackCallBLAS("LAPACKungqr",LAPACKungqr_(&Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
824c666adbfSMark F. Adams       if (INFO != 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xORGQR error arg %d",-INFO);
8252e68589bSMark F. Adams 
8262e68589bSMark F. Adams       for (ii = 0; ii < M; ii++) {
8272e68589bSMark F. Adams         for (jj = 0; jj < N; jj++) {
8282e68589bSMark F. Adams           qqr[N*ii + jj] = qqc[jj*Mdata + ii];
8292e68589bSMark F. Adams         }
8302e68589bSMark F. Adams       }
8312e68589bSMark F. Adams 
8322e68589bSMark F. Adams       /* add diagonal block of P0 */
833c8b0795cSMark F. Adams       for (kk=0; kk<N; kk++) {
834c8b0795cSMark F. Adams         cids[kk] = N*cgid + kk; /* global col IDs in P0 */
835c8b0795cSMark F. Adams       }
8362e68589bSMark F. Adams       ierr = MatSetValues(a_Prol,M,fids,N,cids,qqr,INSERT_VALUES);CHKERRQ(ierr);
8372e68589bSMark F. Adams 
838e632b94dSBarry Smith       ierr = PetscFree5(qqc,qqr,TAU,WORK,fids);CHKERRQ(ierr);
839b43b03e9SMark F. Adams       clid++;
8400cbbd2e1SMark F. Adams     } /* coarse agg */
8410cbbd2e1SMark F. Adams   } /* for all fine nodes */
8420cbbd2e1SMark F. Adams   ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8430cbbd2e1SMark F. Adams   ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8442e68589bSMark F. Adams 
845b2566f29SBarry Smith /* ierr = MPIU_Allreduce(&ndone, &ii, 1, MPIU_INT, MPIU_SUM, comm); */
8462e68589bSMark F. Adams /* MatGetSize(a_Prol, &kk, &jj); */
847b2566f29SBarry Smith /* ierr = MPIU_Allreduce(&minsz, &jj, 1, MPIU_INT, MPIU_MIN, comm); */
8483b4367a7SBarry Smith /* PetscPrintf(comm," **** [%d]%s %d total done, %d nodes (%d local done), min agg. size = %d\n",rank,__FUNCT__,ii,kk/bs,ndone,jj); */
8492e68589bSMark F. Adams 
850519f805aSKarl Rupp #if defined(OUT_AGGS)
851b2a4f308SMark F. Adams   if (llev==1) fclose(file);
8529057884aSMark F. Adams #endif
8531943db53SBarry Smith   ierr = PCGAMGHashTableDestroy(&fgid_flid);CHKERRQ(ierr);
8542e68589bSMark F. Adams   PetscFunctionReturn(0);
8552e68589bSMark F. Adams }
8562e68589bSMark F. Adams 
8575adeb434SBarry Smith #undef __FUNCT__
8585adeb434SBarry Smith #define __FUNCT__ "PCView_GAMG_AGG"
8595adeb434SBarry Smith static PetscErrorCode PCView_GAMG_AGG(PC pc,PetscViewer viewer)
8605adeb434SBarry Smith {
8615adeb434SBarry Smith   PetscErrorCode ierr;
8625adeb434SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
8635adeb434SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
8645adeb434SBarry Smith   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
8655adeb434SBarry Smith 
8665adeb434SBarry Smith   PetscFunctionBegin;
8675adeb434SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"      AGG specific options\n");CHKERRQ(ierr);
8685adeb434SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"        Symmetric graph %s\n",pc_gamg_agg->sym_graph ? "true" : "false");CHKERRQ(ierr);
869cab9ed1eSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"        Number of levels to square graph %D\n",pc_gamg_agg->square_graph);CHKERRQ(ierr);
870cab9ed1eSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"        Number smoothing steps %D\n",pc_gamg_agg->nsmooths);CHKERRQ(ierr);
8715adeb434SBarry Smith   PetscFunctionReturn(0);
8725adeb434SBarry Smith }
8735adeb434SBarry Smith 
8742e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
8752e68589bSMark F. Adams /*
876fd1112cbSBarry Smith    PCGAMGGraph_AGG
8772e68589bSMark F. Adams 
8782e68589bSMark F. Adams   Input Parameter:
8792e68589bSMark F. Adams    . pc - this
8802e68589bSMark F. Adams    . Amat - matrix on this fine level
8812e68589bSMark F. Adams   Output Parameter:
882c8b0795cSMark F. Adams    . a_Gmat -
8832e68589bSMark F. Adams */
8842e68589bSMark F. Adams #undef __FUNCT__
885fd1112cbSBarry Smith #define __FUNCT__ "PCGAMGGraph_AGG"
886e0877f53SBarry Smith static PetscErrorCode PCGAMGGraph_AGG(PC pc,Mat Amat,Mat *a_Gmat)
887c8b0795cSMark F. Adams {
888c8b0795cSMark F. Adams   PetscErrorCode            ierr;
889c8b0795cSMark F. Adams   PC_MG                     *mg          = (PC_MG*)pc->data;
890c8b0795cSMark F. Adams   PC_GAMG                   *pc_gamg     = (PC_GAMG*)mg->innerctx;
891c8b0795cSMark F. Adams   const PetscReal           vfilter      = pc_gamg->threshold;
892c8b0795cSMark F. Adams   PC_GAMG_AGG               *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
893e0940f08SMark F. Adams   Mat                       Gmat;
8943b4367a7SBarry Smith   MPI_Comm                  comm;
89567744fedSMark F. Adams   PetscBool /* set,flg , */ symm;
896c8b0795cSMark F. Adams 
897c8b0795cSMark F. Adams   PetscFunctionBegin;
8983b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr);
899fd1112cbSBarry Smith   ierr = PetscLogEventBegin(PC_GAMGGraph_AGG,0,0,0,0);CHKERRQ(ierr);
900c8b0795cSMark F. Adams 
90167744fedSMark F. Adams   /* ierr = MatIsSymmetricKnown(Amat, &set, &flg);CHKERRQ(ierr); || !(set && flg) -- this causes lot of symm calls */
902c666adbfSMark F. Adams   symm = (PetscBool)(pc_gamg_agg->sym_graph); /* && !pc_gamg_agg->square_graph; */
9030cbbd2e1SMark F. Adams 
9042d7fac45SMark F. Adams   ierr = PCGAMGCreateGraph(Amat, &Gmat);CHKERRQ(ierr);
905bf4339c2SBarry Smith   ierr = PCGAMGFilterGraph(&Gmat, vfilter, symm);CHKERRQ(ierr);
906c8b0795cSMark F. Adams 
907e0940f08SMark F. Adams   *a_Gmat = Gmat;
908c8b0795cSMark F. Adams 
909fd1112cbSBarry Smith   ierr = PetscLogEventEnd(PC_GAMGGraph_AGG,0,0,0,0);CHKERRQ(ierr);
910c8b0795cSMark F. Adams   PetscFunctionReturn(0);
911c8b0795cSMark F. Adams }
912c8b0795cSMark F. Adams 
913c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
914c8b0795cSMark F. Adams /*
915b43b03e9SMark F. Adams    PCGAMGCoarsen_AGG
916c8b0795cSMark F. Adams 
917c8b0795cSMark F. Adams   Input Parameter:
918e0940f08SMark F. Adams    . a_pc - this
919e0940f08SMark F. Adams   Input/Output Parameter:
9200cbbd2e1SMark F. Adams    . a_Gmat1 - graph on this fine level - coarsening can change this (squares it)
921c8b0795cSMark F. Adams   Output Parameter:
9220cbbd2e1SMark F. Adams    . agg_lists - list of aggregates
923c8b0795cSMark F. Adams */
924c8b0795cSMark F. Adams #undef __FUNCT__
925b43b03e9SMark F. Adams #define __FUNCT__ "PCGAMGCoarsen_AGG"
926e0877f53SBarry Smith static PetscErrorCode PCGAMGCoarsen_AGG(PC a_pc,Mat *a_Gmat1,PetscCoarsenData **agg_lists)
927c8b0795cSMark F. Adams {
928c8b0795cSMark F. Adams   PetscErrorCode ierr;
929e0940f08SMark F. Adams   PC_MG          *mg          = (PC_MG*)a_pc->data;
930c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
931c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
9320cbbd2e1SMark F. Adams   Mat            mat,Gmat2, Gmat1 = *a_Gmat1;  /* squared graph */
9330cbbd2e1SMark F. Adams   IS             perm;
9345cfd4bd9SMark Adams   PetscInt       Istart,Iend,Ii,nloc,bs,n,m;
935c8b0795cSMark F. Adams   PetscInt       *permute;
936c8b0795cSMark F. Adams   PetscBool      *bIndexSet;
937b43b03e9SMark F. Adams   MatCoarsen     crs;
9383b4367a7SBarry Smith   MPI_Comm       comm;
93973911c69SBarry Smith   PetscMPIInt    rank;
940aea53286SMark Adams   PetscReal      hashfact;
941e2c00dcbSBarry Smith   PetscInt       iSwapIndex;
9423b50add6SMark Adams   PetscRandom    random;
943c8b0795cSMark F. Adams 
944c8b0795cSMark F. Adams   PetscFunctionBegin;
9450cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr);
9463b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Gmat1,&comm);CHKERRQ(ierr);
9473b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
948e0940f08SMark F. Adams   ierr = MatGetLocalSize(Gmat1, &n, &m);CHKERRQ(ierr);
94971959b99SBarry Smith   ierr = MatGetBlockSize(Gmat1, &bs);CHKERRQ(ierr);
95071959b99SBarry Smith   if (bs != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"bs %D must be 1",bs);
951c8b0795cSMark F. Adams   nloc = n/bs;
952c8b0795cSMark F. Adams 
9539ab59c8bSMark Adams   if (pc_gamg->current_level < pc_gamg_agg->square_graph) {
954e2c00dcbSBarry 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);
955806fa848SBarry Smith     ierr = MatTransposeMatMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2);CHKERRQ(ierr);
956806fa848SBarry Smith   } else Gmat2 = Gmat1;
957c8b0795cSMark F. Adams 
9585cfd4bd9SMark Adams   /* get MIS aggs - randomize */
959785e854fSJed Brown   ierr = PetscMalloc1(nloc, &permute);CHKERRQ(ierr);
960e2c00dcbSBarry Smith   ierr = PetscCalloc1(nloc, &bIndexSet);CHKERRQ(ierr);
961c8b0795cSMark F. Adams   for (Ii = 0; Ii < nloc; Ii++) {
962c8b0795cSMark F. Adams     permute[Ii]   = Ii;
963c8b0795cSMark F. Adams   }
9643b50add6SMark Adams   ierr = PetscRandomCreate(PETSC_COMM_SELF,&random);CHKERRQ(ierr);
9655cfd4bd9SMark Adams   ierr = MatGetOwnershipRange(Gmat1, &Istart, &Iend);CHKERRQ(ierr);
966c8b0795cSMark F. Adams   for (Ii = 0; Ii < nloc; Ii++) {
9673b50add6SMark Adams     ierr = PetscRandomGetValueReal(random,&hashfact);CHKERRQ(ierr);
968aea53286SMark Adams     iSwapIndex = (PetscInt) (hashfact*nloc)%nloc;
969c8b0795cSMark F. Adams     if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) {
970c8b0795cSMark F. Adams       PetscInt iTemp = permute[iSwapIndex];
971c8b0795cSMark F. Adams       permute[iSwapIndex]   = permute[Ii];
972c8b0795cSMark F. Adams       permute[Ii]           = iTemp;
973c8b0795cSMark F. Adams       bIndexSet[iSwapIndex] = PETSC_TRUE;
974c8b0795cSMark F. Adams     }
975c8b0795cSMark F. Adams   }
976c8b0795cSMark F. Adams   ierr = PetscFree(bIndexSet);CHKERRQ(ierr);
9773b50add6SMark Adams   ierr = PetscRandomDestroy(&random);CHKERRQ(ierr);
978806fa848SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm);CHKERRQ(ierr);
9790cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
9800cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr);
981b43b03e9SMark F. Adams #endif
9823b4367a7SBarry Smith   ierr = MatCoarsenCreate(comm, &crs);CHKERRQ(ierr);
9839057884aSMark F. Adams   ierr = MatCoarsenSetFromOptions(crs);CHKERRQ(ierr);
984b43b03e9SMark F. Adams   ierr = MatCoarsenSetGreedyOrdering(crs, perm);CHKERRQ(ierr);
985b43b03e9SMark F. Adams   ierr = MatCoarsenSetAdjacency(crs, Gmat2);CHKERRQ(ierr);
986b43b03e9SMark F. Adams   ierr = MatCoarsenSetStrictAggs(crs, PETSC_TRUE);CHKERRQ(ierr);
987b43b03e9SMark F. Adams   ierr = MatCoarsenApply(crs);CHKERRQ(ierr);
9880cbbd2e1SMark F. Adams   ierr = MatCoarsenGetData(crs, agg_lists);CHKERRQ(ierr); /* output */
989b43b03e9SMark F. Adams   ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr);
990b43b03e9SMark F. Adams 
991c8b0795cSMark F. Adams   ierr = ISDestroy(&perm);CHKERRQ(ierr);
992c8b0795cSMark F. Adams   ierr = PetscFree(permute);CHKERRQ(ierr);
9930cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
9940cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr);
995b43b03e9SMark F. Adams #endif
9969f3f12c8SMark F. Adams 
997c8b0795cSMark F. Adams   /* smooth aggs */
998e0940f08SMark F. Adams   if (Gmat2 != Gmat1) {
9990cbbd2e1SMark F. Adams     const PetscCoarsenData *llist = *agg_lists;
1000*e3df7ea3SBarry Smith     ierr     = smoothAggs(a_pc,Gmat2, Gmat1, *agg_lists);CHKERRQ(ierr);
1001c8b0795cSMark F. Adams     ierr     = MatDestroy(&Gmat1);CHKERRQ(ierr);
1002e0940f08SMark F. Adams     *a_Gmat1 = Gmat2; /* output */
100341b27cdeSMark F. Adams     ierr     = PetscCDGetMat(llist, &mat);CHKERRQ(ierr);
10043b4367a7SBarry Smith     if (mat) SETERRQ(comm,PETSC_ERR_ARG_WRONG, "Auxilary matrix with squared graph????");
1005806fa848SBarry Smith   } else {
10060cbbd2e1SMark F. Adams     const PetscCoarsenData *llist = *agg_lists;
10079ab59c8bSMark Adams     /* see if we have a matrix that takes precedence (returned from MatCoarsenApply) */
100841b27cdeSMark F. Adams     ierr = PetscCDGetMat(llist, &mat);CHKERRQ(ierr);
10090cbbd2e1SMark F. Adams     if (mat) {
10100cbbd2e1SMark F. Adams       ierr     = MatDestroy(&Gmat1);CHKERRQ(ierr);
10110cbbd2e1SMark F. Adams       *a_Gmat1 = mat; /* output */
10120cbbd2e1SMark F. Adams     }
10130cbbd2e1SMark F. Adams   }
10140cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr);
1015c8b0795cSMark F. Adams   PetscFunctionReturn(0);
1016c8b0795cSMark F. Adams }
1017c8b0795cSMark F. Adams 
1018c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
1019c8b0795cSMark F. Adams /*
10200cbbd2e1SMark F. Adams  PCGAMGProlongator_AGG
1021c8b0795cSMark F. Adams 
1022c8b0795cSMark F. Adams  Input Parameter:
1023c8b0795cSMark F. Adams  . pc - this
1024c8b0795cSMark F. Adams  . Amat - matrix on this fine level
1025c8b0795cSMark F. Adams  . Graph - used to get ghost data for nodes in
10260cbbd2e1SMark F. Adams  . agg_lists - list of aggregates
1027c8b0795cSMark F. Adams  Output Parameter:
1028c8b0795cSMark F. Adams  . a_P_out - prolongation operator to the next level
1029c8b0795cSMark F. Adams  */
1030c8b0795cSMark F. Adams #undef __FUNCT__
10310cbbd2e1SMark F. Adams #define __FUNCT__ "PCGAMGProlongator_AGG"
1032e0877f53SBarry Smith static PetscErrorCode PCGAMGProlongator_AGG(PC pc,Mat Amat,Mat Gmat,PetscCoarsenData *agg_lists,Mat *a_P_out)
10332e68589bSMark F. Adams {
10342e68589bSMark F. Adams   PC_MG          *mg       = (PC_MG*)pc->data;
10352e68589bSMark F. Adams   PC_GAMG        *pc_gamg  = (PC_GAMG*)mg->innerctx;
10364a2f8832SBarry Smith   const PetscInt col_bs = pc_gamg->data_cell_cols;
10372e68589bSMark F. Adams   PetscErrorCode ierr;
1038c8b0795cSMark F. Adams   PetscInt       Istart,Iend,nloc,ii,jj,kk,my0,nLocalSelected,bs;
1039c8b0795cSMark F. Adams   Mat            Prol;
1040c5df96a5SBarry Smith   PetscMPIInt    rank, size;
10413b4367a7SBarry Smith   MPI_Comm       comm;
1042c8b0795cSMark F. Adams   PetscReal      *data_w_ghost;
1043c8b0795cSMark F. Adams   PetscInt       myCrs0, nbnodes=0, *flid_fgid;
1044d9558ea9SBarry Smith   MatType        mtype;
10452e68589bSMark F. Adams 
10462e68589bSMark F. Adams   PetscFunctionBegin;
10473b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr);
10484a2f8832SBarry Smith   if (col_bs < 1) SETERRQ(comm,PETSC_ERR_PLIB,"Column bs cannot be less than 1");
10490cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
10503b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
10513b4367a7SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
10522e68589bSMark F. Adams   ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr);
1053c8b0795cSMark F. Adams   ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr);
105471959b99SBarry Smith   nloc = (Iend-Istart)/bs; my0 = Istart/bs;
105571959b99SBarry Smith   if ((Iend-Istart) % bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"(Iend %D - Istart %D) not divisible by bs %D",Iend,Istart,bs);
10562e68589bSMark F. Adams 
10572e68589bSMark F. Adams   /* get 'nLocalSelected' */
10580cbbd2e1SMark F. Adams   for (ii=0, nLocalSelected = 0; ii < nloc; ii++) {
1059e78576d6SMark F. Adams     PetscBool ise;
10600cbbd2e1SMark F. Adams     /* filter out singletons 0 or 1? */
1061e78576d6SMark F. Adams     ierr = PetscCDEmptyAt(agg_lists, ii, &ise);CHKERRQ(ierr);
1062e78576d6SMark F. Adams     if (!ise) nLocalSelected++;
10632e68589bSMark F. Adams   }
10642e68589bSMark F. Adams 
10652e68589bSMark F. Adams   /* create prolongator, create P matrix */
1066d9558ea9SBarry Smith   ierr = MatGetType(Amat,&mtype);CHKERRQ(ierr);
10673b4367a7SBarry Smith   ierr = MatCreate(comm, &Prol);CHKERRQ(ierr);
1068806fa848SBarry Smith   ierr = MatSetSizes(Prol,nloc*bs,nLocalSelected*col_bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1069a2f3521dSMark F. Adams   ierr = MatSetBlockSizes(Prol, bs, col_bs);CHKERRQ(ierr);
1070d9558ea9SBarry Smith   ierr = MatSetType(Prol, mtype);CHKERRQ(ierr);
10714a2f8832SBarry Smith   ierr = MatSeqAIJSetPreallocation(Prol, col_bs, NULL);CHKERRQ(ierr);
10724a2f8832SBarry Smith   ierr = MatMPIAIJSetPreallocation(Prol,col_bs, NULL,col_bs, NULL);CHKERRQ(ierr);
10732e68589bSMark F. Adams 
10742e68589bSMark F. Adams   /* can get all points "removed" */
1075c8b0795cSMark F. Adams   ierr =  MatGetSize(Prol, &kk, &ii);CHKERRQ(ierr);
10767f66b68fSMark Adams   if (!ii) {
1077569f4572SMark Adams     ierr = PetscInfo(pc,"No selected points on coarse grid\n");CHKERRQ(ierr);
10782e68589bSMark F. Adams     ierr = MatDestroy(&Prol);CHKERRQ(ierr);
10790298fd71SBarry Smith     *a_P_out = NULL;  /* out */
1080e87675ddSBarry Smith     ierr = PetscLogEventEnd(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
10812e68589bSMark F. Adams     PetscFunctionReturn(0);
10822e68589bSMark F. Adams   }
1083bf4339c2SBarry Smith   ierr = PetscInfo1(pc,"New grid %D nodes\n",ii/col_bs);CHKERRQ(ierr);
1084c8b0795cSMark F. Adams   ierr = MatGetOwnershipRangeColumn(Prol, &myCrs0, &kk);CHKERRQ(ierr);
10850cbbd2e1SMark F. Adams 
108671959b99SBarry 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);
1087c8b0795cSMark F. Adams   myCrs0 = myCrs0/col_bs;
108871959b99SBarry 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);
10892e68589bSMark F. Adams 
10902e68589bSMark F. Adams   /* create global vector of data in 'data_w_ghost' */
10910cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
10920cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr);
10932e68589bSMark F. Adams #endif
1094c5df96a5SBarry Smith   if (size > 1) { /*  */
10952e68589bSMark F. Adams     PetscReal *tmp_gdata,*tmp_ldata,*tp2;
1096785e854fSJed Brown     ierr = PetscMalloc1(nloc, &tmp_ldata);CHKERRQ(ierr);
10974a2f8832SBarry Smith     for (jj = 0; jj < col_bs; jj++) {
1098c8b0795cSMark F. Adams       for (kk = 0; kk < bs; kk++) {
1099a2f3521dSMark F. Adams         PetscInt        ii,stride;
1100c8b0795cSMark F. Adams         const PetscReal *tp = pc_gamg->data + jj*bs*nloc + kk;
11012fa5cd67SKarl Rupp         for (ii = 0; ii < nloc; ii++, tp += bs) tmp_ldata[ii] = *tp;
11022fa5cd67SKarl Rupp 
1103806fa848SBarry Smith         ierr = PCGAMGGetDataWithGhosts(Gmat, 1, tmp_ldata, &stride, &tmp_gdata);CHKERRQ(ierr);
1104a2f3521dSMark F. Adams 
11057f66b68fSMark Adams         if (!jj && !kk) { /* now I know how many todal nodes - allocate */
11064a2f8832SBarry Smith           ierr    = PetscMalloc1(stride*bs*col_bs, &data_w_ghost);CHKERRQ(ierr);
1107a2f3521dSMark F. Adams           nbnodes = bs*stride;
11082e68589bSMark F. Adams         }
1109a2f3521dSMark F. Adams         tp2 = data_w_ghost + jj*bs*stride + kk;
11102fa5cd67SKarl Rupp         for (ii = 0; ii < stride; ii++, tp2 += bs) *tp2 = tmp_gdata[ii];
11112e68589bSMark F. Adams         ierr = PetscFree(tmp_gdata);CHKERRQ(ierr);
11122e68589bSMark F. Adams       }
11132e68589bSMark F. Adams     }
11142e68589bSMark F. Adams     ierr = PetscFree(tmp_ldata);CHKERRQ(ierr);
1115806fa848SBarry Smith   } else {
1116c8b0795cSMark F. Adams     nbnodes      = bs*nloc;
1117c8b0795cSMark F. Adams     data_w_ghost = (PetscReal*)pc_gamg->data;
11182e68589bSMark F. Adams   }
11192e68589bSMark F. Adams 
11202e68589bSMark F. Adams   /* get P0 */
1121c5df96a5SBarry Smith   if (size > 1) {
11222e68589bSMark F. Adams     PetscReal *fid_glid_loc,*fiddata;
1123a2f3521dSMark F. Adams     PetscInt  stride;
11242e68589bSMark F. Adams 
1125785e854fSJed Brown     ierr = PetscMalloc1(nloc, &fid_glid_loc);CHKERRQ(ierr);
11262e68589bSMark F. Adams     for (kk=0; kk<nloc; kk++) fid_glid_loc[kk] = (PetscReal)(my0+kk);
1127806fa848SBarry Smith     ierr = PCGAMGGetDataWithGhosts(Gmat, 1, fid_glid_loc, &stride, &fiddata);CHKERRQ(ierr);
1128785e854fSJed Brown     ierr = PetscMalloc1(stride, &flid_fgid);CHKERRQ(ierr);
1129a2f3521dSMark F. Adams     for (kk=0; kk<stride; kk++) flid_fgid[kk] = (PetscInt)fiddata[kk];
11302e68589bSMark F. Adams     ierr = PetscFree(fiddata);CHKERRQ(ierr);
1131a2f3521dSMark F. Adams 
113271959b99SBarry Smith     if (stride != nbnodes/bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"stride %D != nbnodes %D/bs %D",stride,nbnodes,bs);
11332e68589bSMark F. Adams     ierr = PetscFree(fid_glid_loc);CHKERRQ(ierr);
1134806fa848SBarry Smith   } else {
1135785e854fSJed Brown     ierr = PetscMalloc1(nloc, &flid_fgid);CHKERRQ(ierr);
11362e68589bSMark F. Adams     for (kk=0; kk<nloc; kk++) flid_fgid[kk] = my0 + kk;
11372e68589bSMark F. Adams   }
11380cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
11390cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr);
11400cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr);
11412e68589bSMark F. Adams #endif
1142c8b0795cSMark F. Adams   {
11430298fd71SBarry Smith     PetscReal *data_out = NULL;
11444a2f8832SBarry Smith     ierr = formProl0(agg_lists, bs, col_bs, myCrs0, nbnodes,data_w_ghost, flid_fgid, &data_out, Prol);CHKERRQ(ierr);
1145c8b0795cSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
11462fa5cd67SKarl Rupp 
1147c8b0795cSMark F. Adams     pc_gamg->data           = data_out;
11484a2f8832SBarry Smith     pc_gamg->data_cell_rows = col_bs;
11494a2f8832SBarry Smith     pc_gamg->data_sz        = col_bs*col_bs*nLocalSelected;
1150c8b0795cSMark F. Adams   }
11510cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
11520cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr);
1153c8b0795cSMark F. Adams #endif
115461ba4676SBarry Smith   if (size > 1) {ierr = PetscFree(data_w_ghost);CHKERRQ(ierr);}
11552e68589bSMark F. Adams   ierr = PetscFree(flid_fgid);CHKERRQ(ierr);
11562e68589bSMark F. Adams 
1157c8b0795cSMark F. Adams   *a_P_out = Prol;  /* out */
1158ed4e82eaSPeter Brune 
11590cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
1160c8b0795cSMark F. Adams   PetscFunctionReturn(0);
1161c8b0795cSMark F. Adams }
1162c8b0795cSMark F. Adams 
1163c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
1164c8b0795cSMark F. Adams /*
1165fd1112cbSBarry Smith    PCGAMGOptProlongator_AGG
1166c8b0795cSMark F. Adams 
1167c8b0795cSMark F. Adams   Input Parameter:
1168c8b0795cSMark F. Adams    . pc - this
1169c8b0795cSMark F. Adams    . Amat - matrix on this fine level
1170c8b0795cSMark F. Adams  In/Output Parameter:
11712a808120SBarry Smith    . a_P - prolongation operator to the next level
1172c8b0795cSMark F. Adams */
1173c8b0795cSMark F. Adams #undef __FUNCT__
1174fd1112cbSBarry Smith #define __FUNCT__ "PCGAMGOptProlongator_AGG"
1175e0877f53SBarry Smith static PetscErrorCode PCGAMGOptProlongator_AGG(PC pc,Mat Amat,Mat *a_P)
1176c8b0795cSMark F. Adams {
1177c8b0795cSMark F. Adams   PetscErrorCode ierr;
1178c8b0795cSMark F. Adams   PC_MG          *mg          = (PC_MG*)pc->data;
1179c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
1180c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
1181c8b0795cSMark F. Adams   PetscInt       jj;
1182c8b0795cSMark F. Adams   Mat            Prol  = *a_P;
11833b4367a7SBarry Smith   MPI_Comm       comm;
11842a808120SBarry Smith   KSP            eksp;
11852a808120SBarry Smith   Vec            bb, xx;
11862a808120SBarry Smith   PC             epc;
11872a808120SBarry Smith   PetscReal      alpha, emax, emin;
11883b50add6SMark Adams   PetscRandom    random;
1189c8b0795cSMark F. Adams 
1190c8b0795cSMark F. Adams   PetscFunctionBegin;
11913b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr);
1192fd1112cbSBarry Smith   ierr = PetscLogEventBegin(PC_GAMGOptProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
1193c8b0795cSMark F. Adams 
11942a808120SBarry Smith   /* compute maximum value of operator to be used in smoother */
11952a808120SBarry Smith   if (0 < pc_gamg_agg->nsmooths) {
1196e7e70905SMark Adams     ierr = MatCreateVecs(Amat, &bb, 0);CHKERRQ(ierr);
1197e7e70905SMark Adams     ierr = MatCreateVecs(Amat, &xx, 0);CHKERRQ(ierr);
11983b50add6SMark Adams     ierr = PetscRandomCreate(PETSC_COMM_SELF,&random);CHKERRQ(ierr);
11993b50add6SMark Adams     ierr = VecSetRandom(bb,random);CHKERRQ(ierr);
12003b50add6SMark Adams     ierr = PetscRandomDestroy(&random);CHKERRQ(ierr);
1201e696ed0bSMark F. Adams 
12023b4367a7SBarry Smith     ierr = KSPCreate(comm,&eksp);CHKERRQ(ierr);
1203422a814eSBarry Smith     ierr = KSPSetErrorIfNotConverged(eksp,pc->erroriffailure);CHKERRQ(ierr);
1204806fa848SBarry Smith     ierr = KSPSetTolerances(eksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,10);CHKERRQ(ierr);
12052e68589bSMark F. Adams     ierr = KSPSetNormType(eksp, KSP_NORM_NONE);CHKERRQ(ierr);
1206da509fc8SJed Brown     ierr = KSPSetOptionsPrefix(eksp,((PetscObject)pc)->prefix);CHKERRQ(ierr);
1207c2436cacSMark F. Adams     ierr = KSPAppendOptionsPrefix(eksp, "gamg_est_");CHKERRQ(ierr);
1208c2436cacSMark F. Adams     ierr = KSPSetFromOptions(eksp);CHKERRQ(ierr);
1209c2436cacSMark F. Adams 
1210c2436cacSMark F. Adams     ierr = KSPSetInitialGuessNonzero(eksp, PETSC_FALSE);CHKERRQ(ierr);
121123ee1639SBarry Smith     ierr = KSPSetOperators(eksp, Amat, Amat);CHKERRQ(ierr);
12122e68589bSMark F. Adams     ierr = KSPSetComputeSingularValues(eksp,PETSC_TRUE);CHKERRQ(ierr);
12132e68589bSMark F. Adams 
1214da509fc8SJed Brown     ierr = KSPGetPC(eksp, &epc);CHKERRQ(ierr);
1215da509fc8SJed Brown     ierr = PCSetType(epc, PCJACOBI);CHKERRQ(ierr);  /* smoother in smoothed agg. */
1216c2436cacSMark F. Adams 
12172e68589bSMark F. Adams     /* solve - keep stuff out of logging */
12182e68589bSMark F. Adams     ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr);
12192e68589bSMark F. Adams     ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr);
12202e68589bSMark F. Adams     ierr = KSPSolve(eksp, bb, xx);CHKERRQ(ierr);
12212e68589bSMark F. Adams     ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr);
12222e68589bSMark F. Adams     ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr);
12232e68589bSMark F. Adams 
12242e68589bSMark F. Adams     ierr = KSPComputeExtremeSingularValues(eksp, &emax, &emin);CHKERRQ(ierr);
1225569f4572SMark Adams     ierr = PetscInfo3(pc,"Smooth P0: max eigen=%e min=%e PC=%s\n",emax,emin,PCJACOBI);CHKERRQ(ierr);
12262e68589bSMark F. Adams     ierr = VecDestroy(&xx);CHKERRQ(ierr);
12272e68589bSMark F. Adams     ierr = VecDestroy(&bb);CHKERRQ(ierr);
12282e68589bSMark F. Adams     ierr = KSPDestroy(&eksp);CHKERRQ(ierr);
12292e68589bSMark F. Adams   }
12302e68589bSMark F. Adams 
12312a808120SBarry Smith   /* smooth P0 */
12322a808120SBarry Smith   for (jj = 0; jj < pc_gamg_agg->nsmooths; jj++) {
12332a808120SBarry Smith     Mat       tMat;
12342a808120SBarry Smith     Vec       diag;
12352a808120SBarry Smith 
12362a808120SBarry Smith #if defined PETSC_GAMG_USE_LOG
12372a808120SBarry Smith     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr);
12382a808120SBarry Smith #endif
12392a808120SBarry Smith 
12402e68589bSMark F. Adams     /* smooth P1 := (I - omega/lam D^{-1}A)P0 */
12412e68589bSMark F. Adams     ierr  = MatMatMult(Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat);CHKERRQ(ierr);
12422a7a6963SBarry Smith     ierr  = MatCreateVecs(Amat, &diag, 0);CHKERRQ(ierr);
12432e68589bSMark F. Adams     ierr  = MatGetDiagonal(Amat, diag);CHKERRQ(ierr); /* effectively PCJACOBI */
12442e68589bSMark F. Adams     ierr  = VecReciprocal(diag);CHKERRQ(ierr);
12452e68589bSMark F. Adams     ierr  = MatDiagonalScale(tMat, diag, 0);CHKERRQ(ierr);
12462e68589bSMark F. Adams     ierr  = VecDestroy(&diag);CHKERRQ(ierr);
1247b4ec6429SMark F. Adams     alpha = -1.4/emax;
12482e68589bSMark F. Adams     ierr  = MatAYPX(tMat, alpha, Prol, SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
12492e68589bSMark F. Adams     ierr  = MatDestroy(&Prol);CHKERRQ(ierr);
12502e68589bSMark F. Adams     Prol  = tMat;
12510cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
12520cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr);
12532e68589bSMark F. Adams #endif
12542e68589bSMark F. Adams   }
1255fd1112cbSBarry Smith   ierr = PetscLogEventEnd(PC_GAMGOptProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
1256c8b0795cSMark F. Adams   *a_P = Prol;
1257c8b0795cSMark F. Adams   PetscFunctionReturn(0);
12582e68589bSMark F. Adams }
12592e68589bSMark F. Adams 
1260c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
1261c8b0795cSMark F. Adams /*
1262c8b0795cSMark F. Adams    PCCreateGAMG_AGG
12632e68589bSMark F. Adams 
1264c8b0795cSMark F. Adams   Input Parameter:
1265c8b0795cSMark F. Adams    . pc -
1266c8b0795cSMark F. Adams */
1267c8b0795cSMark F. Adams #undef __FUNCT__
1268c8b0795cSMark F. Adams #define __FUNCT__ "PCCreateGAMG_AGG"
1269c8b0795cSMark F. Adams PetscErrorCode  PCCreateGAMG_AGG(PC pc)
1270c8b0795cSMark F. Adams {
1271c8b0795cSMark F. Adams   PetscErrorCode ierr;
1272c8b0795cSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
1273c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1274c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg;
12752e68589bSMark F. Adams 
1276c8b0795cSMark F. Adams   PetscFunctionBegin;
1277c8b0795cSMark F. Adams   /* create sub context for SA */
1278b00a9115SJed Brown   ierr            = PetscNewLog(pc,&pc_gamg_agg);CHKERRQ(ierr);
1279c8b0795cSMark F. Adams   pc_gamg->subctx = pc_gamg_agg;
1280c8b0795cSMark F. Adams 
12811ab5ffc9SJed Brown   pc_gamg->ops->setfromoptions = PCSetFromOptions_GAMG_AGG;
12829b8ffb57SJed Brown   pc_gamg->ops->destroy        = PCDestroy_GAMG_AGG;
1283c8b0795cSMark F. Adams   /* reset does not do anything; setup not virtual */
1284c8b0795cSMark F. Adams 
1285c8b0795cSMark F. Adams   /* set internal function pointers */
1286fd1112cbSBarry Smith   pc_gamg->ops->graph             = PCGAMGGraph_AGG;
12871ab5ffc9SJed Brown   pc_gamg->ops->coarsen           = PCGAMGCoarsen_AGG;
12881ab5ffc9SJed Brown   pc_gamg->ops->prolongator       = PCGAMGProlongator_AGG;
1289fd1112cbSBarry Smith   pc_gamg->ops->optprolongator    = PCGAMGOptProlongator_AGG;
12901ab5ffc9SJed Brown   pc_gamg->ops->createdefaultdata = PCSetData_AGG;
12915adeb434SBarry Smith   pc_gamg->ops->view              = PCView_GAMG_AGG;
1292c8b0795cSMark F. Adams 
1293cab9ed1eSBarry Smith   pc_gamg_agg->square_graph = 1;
1294cab9ed1eSBarry Smith   pc_gamg_agg->sym_graph    = PETSC_FALSE;
1295cab9ed1eSBarry Smith   pc_gamg_agg->nsmooths     = 1;
1296cab9ed1eSBarry Smith 
1297cab9ed1eSBarry Smith 
1298e0877f53SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNSmooths_C",PCGAMGSetNSmooths_AGG);CHKERRQ(ierr);
1299e0877f53SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetSymGraph_C",PCGAMGSetSymGraph_AGG);CHKERRQ(ierr);
1300e0877f53SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetSquareGraph_C",PCGAMGSetSquareGraph_AGG);CHKERRQ(ierr);
1301bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_AGG);CHKERRQ(ierr);
13022e68589bSMark F. Adams   PetscFunctionReturn(0);
13032e68589bSMark F. Adams }
1304