xref: /petsc/src/ksp/pc/impls/gamg/agg.c (revision 484f0a72670a35713b4330055ba2ddc33a7e6cf6)
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*/
6af0996ceSBarry Smith #include <petsc/private/kspimpl.h>
766f2ef4bSMark Adams #include <petsc/private/dmimpl.h>
82e68589bSMark F. Adams 
92e68589bSMark F. Adams #include <petscblaslapack.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:
282e68589bSMark F. Adams .  -pc_gamg_agg_nsmooths
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 /*@
62c8b0795cSMark F. Adams    PCGAMGSetSymGraph -
63c8b0795cSMark F. Adams 
64c8b0795cSMark F. Adams    Not Collective on PC
65c8b0795cSMark F. Adams 
66c8b0795cSMark F. Adams    Input Parameters:
67c8b0795cSMark F. Adams .  pc - the preconditioner context
68c8b0795cSMark F. Adams 
69c8b0795cSMark F. Adams    Options Database Key:
70c8b0795cSMark F. Adams .  -pc_gamg_sym_graph
71c8b0795cSMark F. Adams 
72c8b0795cSMark F. Adams    Level: intermediate
73c8b0795cSMark F. Adams 
74c8b0795cSMark F. Adams    Concepts: Aggregation AMG preconditioner
75c8b0795cSMark F. Adams 
76c8b0795cSMark F. Adams .seealso: ()
77c8b0795cSMark F. Adams @*/
78c8b0795cSMark F. Adams PetscErrorCode PCGAMGSetSymGraph(PC pc, PetscBool n)
79c8b0795cSMark F. Adams {
80c8b0795cSMark F. Adams   PetscErrorCode ierr;
81c8b0795cSMark F. Adams 
82c8b0795cSMark F. Adams   PetscFunctionBegin;
83c8b0795cSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
84c8b0795cSMark F. Adams   ierr = PetscTryMethod(pc,"PCGAMGSetSymGraph_C",(PC,PetscBool),(pc,n));CHKERRQ(ierr);
85c8b0795cSMark F. Adams   PetscFunctionReturn(0);
86c8b0795cSMark F. Adams }
87c8b0795cSMark F. Adams 
88c8b0795cSMark F. Adams #undef __FUNCT__
89e0877f53SBarry Smith #define __FUNCT__ "PCGAMGSetSymGraph_AGG"
90e0877f53SBarry Smith static PetscErrorCode PCGAMGSetSymGraph_AGG(PC pc, PetscBool n)
91c8b0795cSMark F. Adams {
92c8b0795cSMark F. Adams   PC_MG       *mg          = (PC_MG*)pc->data;
93c8b0795cSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG*)mg->innerctx;
94c8b0795cSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
95c8b0795cSMark F. Adams 
96c8b0795cSMark F. Adams   PetscFunctionBegin;
97c8b0795cSMark F. Adams   pc_gamg_agg->sym_graph = n;
982e68589bSMark F. Adams   PetscFunctionReturn(0);
992e68589bSMark F. Adams }
1002e68589bSMark F. Adams 
101ef4ad70eSMark F. Adams #undef __FUNCT__
102ef4ad70eSMark F. Adams #define __FUNCT__ "PCGAMGSetSquareGraph"
103ef4ad70eSMark F. Adams /*@
104ef4ad70eSMark F. Adams    PCGAMGSetSquareGraph -
105ef4ad70eSMark F. Adams 
106ef4ad70eSMark F. Adams    Not Collective on PC
107ef4ad70eSMark F. Adams 
108ef4ad70eSMark F. Adams    Input Parameters:
109ef4ad70eSMark F. Adams .  pc - the preconditioner context
110ef4ad70eSMark F. Adams 
111ef4ad70eSMark F. Adams    Options Database Key:
112ef4ad70eSMark F. Adams .  -pc_gamg_square_graph
113ef4ad70eSMark F. Adams 
114ef4ad70eSMark F. Adams    Level: intermediate
115ef4ad70eSMark F. Adams 
116ef4ad70eSMark F. Adams    Concepts: Aggregation AMG preconditioner
117ef4ad70eSMark F. Adams 
118ef4ad70eSMark F. Adams .seealso: ()
119ef4ad70eSMark F. Adams @*/
1209ab59c8bSMark Adams PetscErrorCode PCGAMGSetSquareGraph(PC pc, PetscInt n)
121ef4ad70eSMark F. Adams {
122ef4ad70eSMark F. Adams   PetscErrorCode ierr;
123ef4ad70eSMark F. Adams 
124ef4ad70eSMark F. Adams   PetscFunctionBegin;
125ef4ad70eSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1269ab59c8bSMark Adams   ierr = PetscTryMethod(pc,"PCGAMGSetSquareGraph_C",(PC,PetscInt),(pc,n));CHKERRQ(ierr);
127ef4ad70eSMark F. Adams   PetscFunctionReturn(0);
128ef4ad70eSMark F. Adams }
129ef4ad70eSMark F. Adams 
130ef4ad70eSMark F. Adams #undef __FUNCT__
131e0877f53SBarry Smith #define __FUNCT__ "PCGAMGSetSquareGraph_AGG"
132e0877f53SBarry Smith static PetscErrorCode PCGAMGSetSquareGraph_AGG(PC pc, PetscInt n)
133ef4ad70eSMark F. Adams {
134ef4ad70eSMark F. Adams   PC_MG       *mg          = (PC_MG*)pc->data;
135ef4ad70eSMark F. Adams   PC_GAMG     *pc_gamg     = (PC_GAMG*)mg->innerctx;
136ef4ad70eSMark F. Adams   PC_GAMG_AGG *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
137ef4ad70eSMark F. Adams 
138ef4ad70eSMark F. Adams   PetscFunctionBegin;
139ef4ad70eSMark F. Adams   pc_gamg_agg->square_graph = n;
140ef4ad70eSMark F. Adams   PetscFunctionReturn(0);
141ef4ad70eSMark F. Adams }
142ef4ad70eSMark F. Adams 
1432e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
1442e68589bSMark F. Adams /*
1452e68589bSMark F. Adams    PCSetFromOptions_GAMG_AGG
1462e68589bSMark F. Adams 
1472e68589bSMark F. Adams   Input Parameter:
1482e68589bSMark F. Adams    . pc -
1492e68589bSMark F. Adams */
1502e68589bSMark F. Adams #undef __FUNCT__
1512e68589bSMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG_AGG"
152e0877f53SBarry Smith static PetscErrorCode PCSetFromOptions_GAMG_AGG(PetscOptionItems *PetscOptionsObject,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;
157c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
1582e68589bSMark F. Adams 
1592e68589bSMark F. Adams   PetscFunctionBegin;
160e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"GAMG-AGG options");CHKERRQ(ierr);
1612e68589bSMark F. Adams   {
1622e68589bSMark F. Adams     /* -pc_gamg_agg_nsmooths */
163d3042614SMark Adams     pc_gamg_agg->nsmooths = 1;
1642fa5cd67SKarl Rupp 
1658afaa268SBarry 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);
166c8b0795cSMark F. Adams 
167c8b0795cSMark F. Adams     /* -pc_gamg_sym_graph */
168c8b0795cSMark F. Adams     pc_gamg_agg->sym_graph = PETSC_FALSE;
1692fa5cd67SKarl Rupp 
1708afaa268SBarry 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);
171ef4ad70eSMark F. Adams 
172ef4ad70eSMark F. Adams     /* -pc_gamg_square_graph */
1739ab59c8bSMark Adams     pc_gamg_agg->square_graph = 1;
1742fa5cd67SKarl Rupp 
1759ab59c8bSMark 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);
1762e68589bSMark F. Adams   }
1772e68589bSMark F. Adams   ierr = PetscOptionsTail();CHKERRQ(ierr);
1782e68589bSMark F. Adams   PetscFunctionReturn(0);
1792e68589bSMark F. Adams }
1802e68589bSMark F. Adams 
1812e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
1822e68589bSMark F. Adams /*
1832e68589bSMark F. Adams    PCDestroy_AGG
1842e68589bSMark F. Adams 
1852e68589bSMark F. Adams   Input Parameter:
1862e68589bSMark F. Adams    . pc -
1872e68589bSMark F. Adams */
1882e68589bSMark F. Adams #undef __FUNCT__
1899b8ffb57SJed Brown #define __FUNCT__ "PCDestroy_GAMG_AGG"
190e0877f53SBarry Smith static PetscErrorCode PCDestroy_GAMG_AGG(PC pc)
1912e68589bSMark F. Adams {
1922e68589bSMark F. Adams   PetscErrorCode ierr;
1932e68589bSMark F. Adams   PC_MG          *mg          = (PC_MG*)pc->data;
1942e68589bSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
1952e68589bSMark F. Adams 
1962e68589bSMark F. Adams   PetscFunctionBegin;
1979b8ffb57SJed Brown   ierr = PetscFree(pc_gamg->subctx);CHKERRQ(ierr);
1983ae0bb68SMark Adams   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",NULL);CHKERRQ(ierr);
1992e68589bSMark F. Adams   PetscFunctionReturn(0);
2002e68589bSMark F. Adams }
2012e68589bSMark F. Adams 
2022e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
2032e68589bSMark F. Adams /*
2042e68589bSMark F. Adams    PCSetCoordinates_AGG
205302f38e8SMark F. Adams      - collective
2062e68589bSMark F. Adams 
2072e68589bSMark F. Adams    Input Parameter:
2082e68589bSMark F. Adams    . pc - the preconditioner context
209a2f3521dSMark F. Adams    . ndm - dimesion of data (used for dof/vertex for Stokes)
210302f38e8SMark F. Adams    . a_nloc - number of vertices local
211302f38e8SMark F. Adams    . coords - [a_nloc][ndm] - interleaved coordinate data: {x_0, y_0, z_0, x_1, y_1, ...}
2122e68589bSMark F. Adams */
2131e6b0712SBarry Smith 
2142e68589bSMark F. Adams #undef __FUNCT__
2152e68589bSMark F. Adams #define __FUNCT__ "PCSetCoordinates_AGG"
2161e6b0712SBarry Smith static PetscErrorCode PCSetCoordinates_AGG(PC pc, PetscInt ndm, PetscInt a_nloc, PetscReal *coords)
2172e68589bSMark F. Adams {
2182e68589bSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
2192e68589bSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
2202e68589bSMark F. Adams   PetscErrorCode ierr;
22169344418SMark F. Adams   PetscInt       arrsz,kk,ii,jj,nloc,ndatarows,ndf;
222a2f3521dSMark F. Adams   Mat            mat = pc->pmat;
2232e68589bSMark F. Adams 
2242e68589bSMark F. Adams   PetscFunctionBegin;
225a2f3521dSMark F. Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
226a2f3521dSMark F. Adams   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
227302f38e8SMark F. Adams   nloc = a_nloc;
2282e68589bSMark F. Adams 
2292e68589bSMark F. Adams   /* SA: null space vectors */
23069344418SMark F. Adams   ierr = MatGetBlockSize(mat, &ndf);CHKERRQ(ierr); /* this does not work for Stokes */
23169344418SMark F. Adams   if (coords && ndf==1) pc_gamg->data_cell_cols = 1; /* scalar w/ coords and SA (not needed) */
232a2f3521dSMark F. Adams   else if (coords) {
233c666adbfSMark F. Adams     if (ndm > ndf) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"degrees of motion %d > block size %d",ndm,ndf);
23469344418SMark F. Adams     pc_gamg->data_cell_cols = (ndm==2 ? 3 : 6); /* displacement elasticity */
23569344418SMark F. Adams     if (ndm != ndf) {
236c666adbfSMark 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);
237a2f3521dSMark F. Adams     }
238806fa848SBarry Smith   } else pc_gamg->data_cell_cols = ndf; /* no data, force SA with constant null space vectors */
23971959b99SBarry Smith   pc_gamg->data_cell_rows = ndatarows = ndf;
24071959b99SBarry 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);
241c8b0795cSMark F. Adams   arrsz = nloc*pc_gamg->data_cell_rows*pc_gamg->data_cell_cols;
2422e68589bSMark F. Adams 
2432e68589bSMark F. Adams   /* create data - syntactic sugar that should be refactored at some point */
2447f66b68fSMark Adams   if (!pc_gamg->data || (pc_gamg->data_sz != arrsz)) {
2452e68589bSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
246854ce69bSBarry Smith     ierr = PetscMalloc1(arrsz+1, &pc_gamg->data);CHKERRQ(ierr);
2472e68589bSMark F. Adams   }
2482e68589bSMark F. Adams   /* copy data in - column oriented */
2492e68589bSMark F. Adams   for (kk=0; kk<nloc; kk++) {
25069344418SMark F. Adams     const PetscInt M     = nloc*pc_gamg->data_cell_rows; /* stride into data */
25169344418SMark F. Adams     PetscReal      *data = &pc_gamg->data[kk*ndatarows]; /* start of cell */
252c8b0795cSMark F. Adams     if (pc_gamg->data_cell_cols==1) *data = 1.0;
2532e68589bSMark F. Adams     else {
25469344418SMark F. Adams       /* translational modes */
2552fa5cd67SKarl Rupp       for (ii=0;ii<ndatarows;ii++) {
2562fa5cd67SKarl Rupp         for (jj=0;jj<ndatarows;jj++) {
25769344418SMark F. Adams           if (ii==jj)data[ii*M + jj] = 1.0;
2582e68589bSMark F. Adams           else data[ii*M + jj] = 0.0;
2592fa5cd67SKarl Rupp         }
2602fa5cd67SKarl Rupp       }
26169344418SMark F. Adams 
26269344418SMark F. Adams       /* rotational modes */
2632e68589bSMark F. Adams       if (coords) {
26469344418SMark F. Adams         if (ndm == 2) {
2652e68589bSMark F. Adams           data   += 2*M;
2662e68589bSMark F. Adams           data[0] = -coords[2*kk+1];
2672e68589bSMark F. Adams           data[1] =  coords[2*kk];
268806fa848SBarry Smith         } else {
2692e68589bSMark F. Adams           data   += 3*M;
2702e68589bSMark F. Adams           data[0] = 0.0;             data[M+0] =  coords[3*kk+2]; data[2*M+0] = -coords[3*kk+1];
2712e68589bSMark F. Adams           data[1] = -coords[3*kk+2]; data[M+1] = 0.0;             data[2*M+1] =  coords[3*kk];
2722e68589bSMark F. Adams           data[2] =  coords[3*kk+1]; data[M+2] = -coords[3*kk];   data[2*M+2] = 0.0;
2732e68589bSMark F. Adams         }
2742e68589bSMark F. Adams       }
2752e68589bSMark F. Adams     }
2762e68589bSMark F. Adams   }
2772e68589bSMark F. Adams 
2782e68589bSMark F. Adams   pc_gamg->data_sz = arrsz;
2792e68589bSMark F. Adams   PetscFunctionReturn(0);
2802e68589bSMark F. Adams }
2812e68589bSMark F. Adams 
282b43b03e9SMark F. Adams typedef PetscInt NState;
283b43b03e9SMark F. Adams static const NState NOT_DONE=-2;
284b43b03e9SMark F. Adams static const NState DELETED =-1;
285b43b03e9SMark F. Adams static const NState REMOVED =-3;
286b43b03e9SMark F. Adams #define IS_SELECTED(s) (s!=DELETED && s!=NOT_DONE && s!=REMOVED)
287b43b03e9SMark F. Adams 
288c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
289c8b0795cSMark F. Adams /*
290b43b03e9SMark F. Adams    smoothAggs - greedy grab of with G1 (unsquared graph) -- AIJ specific
291b43b03e9SMark F. Adams      - AGG-MG specific: clears singletons out of 'selected_2'
292c8b0795cSMark F. Adams 
293c8b0795cSMark F. Adams    Input Parameter:
2942adfe9d3SBarry Smith    . Gmat_2 - glabal matrix of graph (data not defined)   base (squared) graph
2952adfe9d3SBarry Smith    . Gmat_1 - base graph to grab with                 base graph
296c8b0795cSMark F. Adams    Input/Output Parameter:
2970cbbd2e1SMark F. Adams    . aggs_2 - linked list of aggs with gids)
298c8b0795cSMark F. Adams */
299c8b0795cSMark F. Adams #undef __FUNCT__
300c8b0795cSMark F. Adams #define __FUNCT__ "smoothAggs"
3012adfe9d3SBarry Smith static PetscErrorCode smoothAggs(Mat Gmat_2, Mat Gmat_1,PetscCoarsenData *aggs_2)
302c8b0795cSMark F. Adams {
303c8b0795cSMark F. Adams   PetscErrorCode ierr;
304c8b0795cSMark F. Adams   PetscBool      isMPI;
305c8b0795cSMark F. Adams   Mat_SeqAIJ     *matA_1, *matB_1=0, *matA_2, *matB_2=0;
3063b4367a7SBarry Smith   MPI_Comm       comm;
3070cbbd2e1SMark F. Adams   PetscInt       lid,*ii,*idx,ix,Iend,my0,kk,n,j;
308c8b0795cSMark F. Adams   Mat_MPIAIJ     *mpimat_2 = 0, *mpimat_1=0;
309c8b0795cSMark F. Adams   const PetscInt nloc      = Gmat_2->rmap->n;
3100cbbd2e1SMark F. Adams   PetscScalar    *cpcol_1_state,*cpcol_2_state,*cpcol_2_par_orig,*lid_parent_gid;
3110cbbd2e1SMark F. Adams   PetscInt       *lid_cprowID_1;
312c8b0795cSMark F. Adams   NState         *lid_state;
3130cbbd2e1SMark F. Adams   Vec            ghost_par_orig2;
314c8b0795cSMark F. Adams 
315c8b0795cSMark F. Adams   PetscFunctionBegin;
3163b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Gmat_2,&comm);CHKERRQ(ierr);
317c8b0795cSMark F. Adams   ierr = MatGetOwnershipRange(Gmat_1,&my0,&Iend);CHKERRQ(ierr);
318c8b0795cSMark F. Adams 
3190cbbd2e1SMark F. Adams   if (PETSC_FALSE) {
320c8b0795cSMark F. Adams     PetscViewer viewer; char fname[32]; static int llev=0;
321c8b0795cSMark F. Adams     sprintf(fname,"Gmat2_%d.m",llev++);
3223b4367a7SBarry Smith     PetscViewerASCIIOpen(comm,fname,&viewer);
3236a9046bcSBarry Smith     ierr = PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
324c8b0795cSMark F. Adams     ierr = MatView(Gmat_2, viewer);CHKERRQ(ierr);
325f159cad9SBarry Smith     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
326c8b0795cSMark F. Adams     ierr = PetscViewerDestroy(&viewer);
327c8b0795cSMark F. Adams   }
328c8b0795cSMark F. Adams 
329c8b0795cSMark F. Adams   /* get submatrices */
330251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)Gmat_1, MATMPIAIJ, &isMPI);CHKERRQ(ierr);
331c8b0795cSMark F. Adams   if (isMPI) {
332c8b0795cSMark F. Adams     /* grab matrix objects */
333c8b0795cSMark F. Adams     mpimat_2 = (Mat_MPIAIJ*)Gmat_2->data;
334c8b0795cSMark F. Adams     mpimat_1 = (Mat_MPIAIJ*)Gmat_1->data;
335c8b0795cSMark F. Adams     matA_1   = (Mat_SeqAIJ*)mpimat_1->A->data;
336c8b0795cSMark F. Adams     matB_1   = (Mat_SeqAIJ*)mpimat_1->B->data;
337c8b0795cSMark F. Adams     matA_2   = (Mat_SeqAIJ*)mpimat_2->A->data;
338c8b0795cSMark F. Adams     matB_2   = (Mat_SeqAIJ*)mpimat_2->B->data;
339c8b0795cSMark F. Adams 
340c8b0795cSMark F. Adams     /* force compressed row storage for B matrix in AuxMat */
34111e456e1SBarry Smith     ierr = MatCheckCompressedRow(mpimat_1->B,matB_1->nonzerorowcnt,&matB_1->compressedrow,matB_1->i,Gmat_1->rmap->n,-1.0);CHKERRQ(ierr);
342c8b0795cSMark F. Adams 
343785e854fSJed Brown     ierr = PetscMalloc1(nloc, &lid_cprowID_1);CHKERRQ(ierr);
3440cbbd2e1SMark F. Adams     for (lid = 0; lid < nloc; lid++) lid_cprowID_1[lid] = -1;
345c8b0795cSMark F. Adams     for (ix=0; ix<matB_1->compressedrow.nrows; ix++) {
346c8b0795cSMark F. Adams       PetscInt lid = matB_1->compressedrow.rindex[ix];
347c8b0795cSMark F. Adams       lid_cprowID_1[lid] = ix;
348c8b0795cSMark F. Adams     }
349806fa848SBarry Smith   } else {
350c8b0795cSMark F. Adams     matA_1        = (Mat_SeqAIJ*)Gmat_1->data;
351c8b0795cSMark F. Adams     matA_2        = (Mat_SeqAIJ*)Gmat_2->data;
3520298fd71SBarry Smith     lid_cprowID_1 = NULL;
353c8b0795cSMark F. Adams   }
35478a438d6SMark Adams   if (nloc>0) {
35571959b99SBarry Smith     if (!(matA_1 && !matA_1->compressedrow.use)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"!(matA_1 && !matA_1->compressedrow.use)");
3567f66b68fSMark Adams     if (!(!matB_1 || matB_1->compressedrow.use)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"!(matB_1==0 || matB_1->compressedrow.use)");
35771959b99SBarry Smith     if (!(matA_2 && !matA_2->compressedrow.use)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"!(matA_2 && !matA_2->compressedrow.use)");
3587f66b68fSMark Adams     if (!(!matB_2 || matB_2->compressedrow.use)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"!(matB_2==0 || matB_2->compressedrow.use)");
35978a438d6SMark Adams   }
360c8b0795cSMark F. Adams   /* get state of locals and selected gid for deleted */
361e632b94dSBarry Smith   ierr = PetscMalloc2(nloc, &lid_state,nloc, &lid_parent_gid);CHKERRQ(ierr);
362c8b0795cSMark F. Adams   for (lid = 0; lid < nloc; lid++) {
3630cbbd2e1SMark F. Adams     lid_parent_gid[lid] = -1.0;
364c8b0795cSMark F. Adams     lid_state[lid]      = DELETED;
365c8b0795cSMark F. Adams   }
3660cbbd2e1SMark F. Adams 
3670cbbd2e1SMark F. Adams   /* set lid_state */
3680cbbd2e1SMark F. Adams   for (lid = 0; lid < nloc; lid++) {
36941b27cdeSMark F. Adams     PetscCDPos pos;
370e78576d6SMark F. Adams     ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr);
371e78576d6SMark F. Adams     if (pos) {
372e78576d6SMark F. Adams       PetscInt gid1;
3732fa5cd67SKarl Rupp 
374*484f0a72SBarry Smith       ierr = PetscCDIntNdGetID(pos, &gid1);CHKERRQ(ierr);
37571959b99SBarry Smith       if (gid1 != lid+my0) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gid1 %D != lid %D + my0 %D",gid1,lid,my0);
3760cbbd2e1SMark F. Adams       lid_state[lid] = gid1;
377b43b03e9SMark F. Adams     }
378b43b03e9SMark F. Adams   }
3790cbbd2e1SMark F. Adams 
3800cbbd2e1SMark F. Adams   /* map local to selected local, DELETED means a ghost owns it */
381c8b0795cSMark F. Adams   for (lid=kk=0; lid<nloc; lid++) {
382c8b0795cSMark F. Adams     NState state = lid_state[lid];
383c8b0795cSMark F. Adams     if (IS_SELECTED(state)) {
38441b27cdeSMark F. Adams       PetscCDPos pos;
385e78576d6SMark F. Adams       ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr);
386e78576d6SMark F. Adams       while (pos) {
387e78576d6SMark F. Adams         PetscInt gid1;
388*484f0a72SBarry Smith         ierr = PetscCDIntNdGetID(pos, &gid1);CHKERRQ(ierr);
389e78576d6SMark F. Adams         ierr = PetscCDGetNextPos(aggs_2,lid,&pos);CHKERRQ(ierr);
390e78576d6SMark F. Adams 
3912fa5cd67SKarl Rupp         if (gid1 >= my0 && gid1 < Iend) lid_parent_gid[gid1-my0] = (PetscScalar)(lid + my0);
392c8b0795cSMark F. Adams       }
3930cbbd2e1SMark F. Adams     }
3940cbbd2e1SMark F. Adams   }
3950cbbd2e1SMark F. Adams   /* get 'cpcol_1/2_state' & cpcol_2_par_orig - uses mpimat_1/2->lvec for temp space */
396c8b0795cSMark F. Adams   if (isMPI) {
397c8b0795cSMark F. Adams     Vec tempVec;
398c8b0795cSMark F. Adams     /* get 'cpcol_1_state' */
3992a7a6963SBarry Smith     ierr = MatCreateVecs(Gmat_1, &tempVec, 0);CHKERRQ(ierr);
400c8b0795cSMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
401c8b0795cSMark F. Adams       PetscScalar v = (PetscScalar)lid_state[kk];
402c8b0795cSMark F. Adams       ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr);
403c8b0795cSMark F. Adams     }
404c8b0795cSMark F. Adams     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
405c8b0795cSMark F. Adams     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
406806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
407806fa848SBarry Smith     ierr = VecScatterEnd(mpimat_1->Mvctx,tempVec, mpimat_1->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
408c8b0795cSMark F. Adams     ierr = VecGetArray(mpimat_1->lvec, &cpcol_1_state);CHKERRQ(ierr);
409c8b0795cSMark F. Adams     /* get 'cpcol_2_state' */
410806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
411806fa848SBarry Smith     ierr =   VecScatterEnd(mpimat_2->Mvctx,tempVec, mpimat_2->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
412c8b0795cSMark F. Adams     ierr = VecGetArray(mpimat_2->lvec, &cpcol_2_state);CHKERRQ(ierr);
4130cbbd2e1SMark F. Adams     /* get 'cpcol_2_par_orig' */
4140cbbd2e1SMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
4150cbbd2e1SMark F. Adams       PetscScalar v = (PetscScalar)lid_parent_gid[kk];
4160cbbd2e1SMark F. Adams       ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr);
4170cbbd2e1SMark F. Adams     }
4180cbbd2e1SMark F. Adams     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
4190cbbd2e1SMark F. Adams     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
4200cbbd2e1SMark F. Adams     ierr = VecDuplicate(mpimat_2->lvec, &ghost_par_orig2);CHKERRQ(ierr);
421806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
422806fa848SBarry Smith     ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghost_par_orig2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4230cbbd2e1SMark F. Adams     ierr = VecGetArray(ghost_par_orig2, &cpcol_2_par_orig);CHKERRQ(ierr);
4240cbbd2e1SMark F. Adams 
425c8b0795cSMark F. Adams     ierr = VecDestroy(&tempVec);CHKERRQ(ierr);
426c8b0795cSMark F. Adams   } /* ismpi */
427c8b0795cSMark F. Adams 
428c8b0795cSMark F. Adams   /* doit */
429c8b0795cSMark F. Adams   for (lid=0; lid<nloc; lid++) {
430c8b0795cSMark F. Adams     NState state = lid_state[lid];
4310cbbd2e1SMark F. Adams     if (IS_SELECTED(state)) {
4320cbbd2e1SMark F. Adams       /* steal locals */
433c8b0795cSMark F. Adams       ii  = matA_1->i; n = ii[lid+1] - ii[lid];
434c8b0795cSMark F. Adams       idx = matA_1->j + ii[lid];
435c8b0795cSMark F. Adams       for (j=0; j<n; j++) {
4360cbbd2e1SMark F. Adams         PetscInt lidj   = idx[j], sgid;
437c8b0795cSMark F. Adams         NState   statej = lid_state[lidj];
4380cbbd2e1SMark F. Adams         if (statej==DELETED && (sgid=(PetscInt)PetscRealPart(lid_parent_gid[lidj])) != lid+my0) { /* steal local */
4390cbbd2e1SMark F. Adams           lid_parent_gid[lidj] = (PetscScalar)(lid+my0); /* send this if sgid is not local */
4400cbbd2e1SMark F. Adams           if (sgid >= my0 && sgid < Iend) {       /* I'm stealing this local from a local sgid */
4410cbbd2e1SMark F. Adams             PetscInt   hav=0,slid=sgid-my0,gidj=lidj+my0;
4420298fd71SBarry Smith             PetscCDPos pos,last=NULL;
443c8b0795cSMark F. Adams             /* looking for local from local so id_llist_2 works */
444e78576d6SMark F. Adams             ierr = PetscCDGetHeadPos(aggs_2,slid,&pos);CHKERRQ(ierr);
445e78576d6SMark F. Adams             while (pos) {
446e78576d6SMark F. Adams               PetscInt gid;
447*484f0a72SBarry Smith               ierr = PetscCDIntNdGetID(pos, &gid);CHKERRQ(ierr);
4480cbbd2e1SMark F. Adams               if (gid == gidj) {
44971959b99SBarry Smith                 if (!last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"last cannot be null");
45041b27cdeSMark F. Adams                 ierr = PetscCDRemoveNextNode(aggs_2, slid, last);CHKERRQ(ierr);
45141b27cdeSMark F. Adams                 ierr = PetscCDAppendNode(aggs_2, lid, pos);CHKERRQ(ierr);
4520cbbd2e1SMark F. Adams                 hav  = 1;
453c8b0795cSMark F. Adams                 break;
454806fa848SBarry Smith               } else last = pos;
455e78576d6SMark F. Adams 
456e78576d6SMark F. Adams               ierr = PetscCDGetNextPos(aggs_2,slid,&pos);CHKERRQ(ierr);
457c8b0795cSMark F. Adams             }
458c8b0795cSMark F. Adams             if (hav!=1) {
45971959b99SBarry Smith               if (!hav) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix");
460c666adbfSMark F. Adams               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %d times???",hav);
461c8b0795cSMark F. Adams             }
462806fa848SBarry Smith           } else {            /* I'm stealing this local, owned by a ghost */
4631692a7d8SBarry Smith             if (sgid != -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Have un-symmetric graph (apparently). Use '-pc_gamg_sym_graph true' to symetrize the graph or '-pc_gamg_threshold -1.0' if the matrix is structurally symmetric.");
46441b27cdeSMark F. Adams             ierr = PetscCDAppendID(aggs_2, lid, lidj+my0);CHKERRQ(ierr);
465c8b0795cSMark F. Adams           }
466c8b0795cSMark F. Adams         }
4670cbbd2e1SMark F. Adams       } /* local neighbors */
468806fa848SBarry Smith     } else if (state == DELETED && lid_cprowID_1) {
4690cbbd2e1SMark F. Adams       PetscInt sgidold = (PetscInt)PetscRealPart(lid_parent_gid[lid]);
470c8b0795cSMark F. Adams       /* see if I have a selected ghost neighbor that will steal me */
471c8b0795cSMark F. Adams       if ((ix=lid_cprowID_1[lid]) != -1) {
472c8b0795cSMark F. Adams         ii  = matB_1->compressedrow.i; n = ii[ix+1] - ii[ix];
473c8b0795cSMark F. Adams         idx = matB_1->j + ii[ix];
474c8b0795cSMark F. Adams         for (j=0; j<n; j++) {
475c8b0795cSMark F. Adams           PetscInt cpid   = idx[j];
476c8b0795cSMark F. Adams           NState   statej = (NState)PetscRealPart(cpcol_1_state[cpid]);
477c8b0795cSMark F. Adams           if (IS_SELECTED(statej) && sgidold != (PetscInt)statej) { /* ghost will steal this, remove from my list */
4780cbbd2e1SMark F. Adams             lid_parent_gid[lid] = (PetscScalar)statej; /* send who selected */
4790cbbd2e1SMark F. Adams             if (sgidold>=my0 && sgidold<Iend) { /* this was mine */
4800cbbd2e1SMark F. Adams               PetscInt   hav=0,oldslidj=sgidold-my0;
4810298fd71SBarry Smith               PetscCDPos pos,last=NULL;
4820cbbd2e1SMark F. Adams               /* remove from 'oldslidj' list */
483e78576d6SMark F. Adams               ierr = PetscCDGetHeadPos(aggs_2,oldslidj,&pos);CHKERRQ(ierr);
484e78576d6SMark F. Adams               while (pos) {
485e78576d6SMark F. Adams                 PetscInt gid;
486*484f0a72SBarry Smith                 ierr = PetscCDIntNdGetID(pos, &gid);CHKERRQ(ierr);
4870cbbd2e1SMark F. Adams                 if (lid+my0 == gid) {
4880cbbd2e1SMark F. Adams                   /* id_llist_2[lastid] = id_llist_2[flid];   /\* remove lid from oldslidj list *\/ */
48971959b99SBarry Smith                   if (!last) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"last cannot be null");
49041b27cdeSMark F. Adams                   ierr = PetscCDRemoveNextNode(aggs_2, oldslidj, last);CHKERRQ(ierr);
4910cbbd2e1SMark F. Adams                   /* ghost (PetscScalar)statej will add this later */
4920cbbd2e1SMark F. Adams                   hav = 1;
493c8b0795cSMark F. Adams                   break;
494806fa848SBarry Smith                 } else last = pos;
495e78576d6SMark F. Adams 
496e78576d6SMark F. Adams                 ierr = PetscCDGetNextPos(aggs_2,oldslidj,&pos);CHKERRQ(ierr);
497c8b0795cSMark F. Adams               }
498c8b0795cSMark F. Adams               if (hav!=1) {
4997f66b68fSMark Adams                 if (!hav) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"failed to find adj in 'selected' lists - structurally unsymmetric matrix");
500c666adbfSMark F. Adams                 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"found node %d times???",hav);
501c8b0795cSMark F. Adams               }
502806fa848SBarry Smith             } else {
5030cbbd2e1SMark F. Adams               /* ghosts remove this later */
5040cbbd2e1SMark F. Adams             }
505c8b0795cSMark F. Adams           }
506c8b0795cSMark F. Adams         }
507c8b0795cSMark F. Adams       }
508c8b0795cSMark F. Adams     } /* selected/deleted */
509c8b0795cSMark F. Adams   } /* node loop */
510c8b0795cSMark F. Adams 
511c8b0795cSMark F. Adams   if (isMPI) {
5120cbbd2e1SMark F. Adams     PetscScalar     *cpcol_2_parent,*cpcol_2_gid;
5130cbbd2e1SMark F. Adams     Vec             tempVec,ghostgids2,ghostparents2;
5140cbbd2e1SMark F. Adams     PetscInt        cpid,nghost_2;
5151943db53SBarry Smith     PCGAMGHashTable gid_cpid;
516c8b0795cSMark F. Adams 
5170cbbd2e1SMark F. Adams     ierr = VecGetSize(mpimat_2->lvec, &nghost_2);CHKERRQ(ierr);
5182a7a6963SBarry Smith     ierr = MatCreateVecs(Gmat_2, &tempVec, 0);CHKERRQ(ierr);
5190cbbd2e1SMark F. Adams 
5200cbbd2e1SMark F. Adams     /* get 'cpcol_2_parent' */
521c8b0795cSMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
5220cbbd2e1SMark F. Adams       ierr = VecSetValues(tempVec, 1, &j, &lid_parent_gid[kk], INSERT_VALUES);CHKERRQ(ierr);
523c8b0795cSMark F. Adams     }
524c8b0795cSMark F. Adams     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
525c8b0795cSMark F. Adams     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
5260cbbd2e1SMark F. Adams     ierr = VecDuplicate(mpimat_2->lvec, &ghostparents2);CHKERRQ(ierr);
527806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
528806fa848SBarry Smith     ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostparents2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5290cbbd2e1SMark F. Adams     ierr = VecGetArray(ghostparents2, &cpcol_2_parent);CHKERRQ(ierr);
5300cbbd2e1SMark F. Adams 
5310cbbd2e1SMark F. Adams     /* get 'cpcol_2_gid' */
5320cbbd2e1SMark F. Adams     for (kk=0,j=my0; kk<nloc; kk++,j++) {
5330cbbd2e1SMark F. Adams       PetscScalar v = (PetscScalar)j;
5340cbbd2e1SMark F. Adams       ierr = VecSetValues(tempVec, 1, &j, &v, INSERT_VALUES);CHKERRQ(ierr);
5350cbbd2e1SMark F. Adams     }
5360cbbd2e1SMark F. Adams     ierr = VecAssemblyBegin(tempVec);CHKERRQ(ierr);
5370cbbd2e1SMark F. Adams     ierr = VecAssemblyEnd(tempVec);CHKERRQ(ierr);
5380cbbd2e1SMark F. Adams     ierr = VecDuplicate(mpimat_2->lvec, &ghostgids2);CHKERRQ(ierr);
539806fa848SBarry Smith     ierr = VecScatterBegin(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
540806fa848SBarry Smith     ierr = VecScatterEnd(mpimat_2->Mvctx,tempVec, ghostgids2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5410cbbd2e1SMark F. Adams     ierr = VecGetArray(ghostgids2, &cpcol_2_gid);CHKERRQ(ierr);
5420cbbd2e1SMark F. Adams 
543c8b0795cSMark F. Adams     ierr = VecDestroy(&tempVec);CHKERRQ(ierr);
544c8b0795cSMark F. Adams 
5450cbbd2e1SMark F. Adams     /* look for deleted ghosts and add to table */
5461943db53SBarry Smith     ierr = PCGAMGHashTableCreate(2*nghost_2+1, &gid_cpid);CHKERRQ(ierr);
5470cbbd2e1SMark F. Adams     for (cpid = 0; cpid < nghost_2; cpid++) {
5480cbbd2e1SMark F. Adams       NState state = (NState)PetscRealPart(cpcol_2_state[cpid]);
5490cbbd2e1SMark F. Adams       if (state==DELETED) {
5500cbbd2e1SMark F. Adams         PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]);
5510cbbd2e1SMark F. Adams         PetscInt sgid_old = (PetscInt)PetscRealPart(cpcol_2_par_orig[cpid]);
5520cbbd2e1SMark F. Adams         if (sgid_old == -1 && sgid_new != -1) {
5530cbbd2e1SMark F. Adams           PetscInt gid = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]);
5541943db53SBarry Smith           ierr = PCGAMGHashTableAdd(&gid_cpid, gid, cpid);CHKERRQ(ierr);
5550cbbd2e1SMark F. Adams         }
5560cbbd2e1SMark F. Adams       }
5570cbbd2e1SMark F. Adams     }
558c8b0795cSMark F. Adams 
5590cbbd2e1SMark F. Adams     /* look for deleted ghosts and see if they moved - remove it */
560c8b0795cSMark F. Adams     for (lid=0; lid<nloc; lid++) {
561c8b0795cSMark F. Adams       NState state = lid_state[lid];
562c8b0795cSMark F. Adams       if (IS_SELECTED(state)) {
5630298fd71SBarry Smith         PetscCDPos pos,last=NULL;
564c8b0795cSMark F. Adams         /* look for deleted ghosts and see if they moved */
565e78576d6SMark F. Adams         ierr = PetscCDGetHeadPos(aggs_2,lid,&pos);CHKERRQ(ierr);
566e78576d6SMark F. Adams         while (pos) {
567e78576d6SMark F. Adams           PetscInt gid;
568*484f0a72SBarry Smith           ierr = PetscCDIntNdGetID(pos, &gid);CHKERRQ(ierr);
569e78576d6SMark F. Adams 
5700cbbd2e1SMark F. Adams           if (gid < my0 || gid >= Iend) {
5711943db53SBarry Smith             ierr = PCGAMGHashTableFind(&gid_cpid, gid, &cpid);CHKERRQ(ierr);
5720cbbd2e1SMark F. Adams             if (cpid != -1) {
5730cbbd2e1SMark F. Adams               /* a moved ghost - */
5740cbbd2e1SMark F. Adams               /* id_llist_2[lastid] = id_llist_2[flid];    /\* remove 'flid' from list *\/ */
57541b27cdeSMark F. Adams               ierr = PetscCDRemoveNextNode(aggs_2, lid, last);CHKERRQ(ierr);
576806fa848SBarry Smith             } else last = pos;
577806fa848SBarry Smith           } else last = pos;
578e78576d6SMark F. Adams 
579e78576d6SMark F. Adams           ierr = PetscCDGetNextPos(aggs_2,lid,&pos);CHKERRQ(ierr);
580c8b0795cSMark F. Adams         } /* loop over list of deleted */
581c8b0795cSMark F. Adams       } /* selected */
582c8b0795cSMark F. Adams     }
5831943db53SBarry Smith     ierr = PCGAMGHashTableDestroy(&gid_cpid);CHKERRQ(ierr);
584c8b0795cSMark F. Adams 
5850cbbd2e1SMark F. Adams     /* look at ghosts, see if they changed - and it */
5860cbbd2e1SMark F. Adams     for (cpid = 0; cpid < nghost_2; cpid++) {
5870cbbd2e1SMark F. Adams       PetscInt sgid_new = (PetscInt)PetscRealPart(cpcol_2_parent[cpid]);
5880cbbd2e1SMark F. Adams       if (sgid_new >= my0 && sgid_new < Iend) { /* this is mine */
5890cbbd2e1SMark F. Adams         PetscInt   gid     = (PetscInt)PetscRealPart(cpcol_2_gid[cpid]);
5900cbbd2e1SMark F. Adams         PetscInt   slid_new=sgid_new-my0,hav=0;
59141b27cdeSMark F. Adams         PetscCDPos pos;
5920cbbd2e1SMark F. Adams         /* search for this gid to see if I have it */
593e78576d6SMark F. Adams         ierr = PetscCDGetHeadPos(aggs_2,slid_new,&pos);CHKERRQ(ierr);
594e78576d6SMark F. Adams         while (pos) {
595e78576d6SMark F. Adams           PetscInt gidj;
596*484f0a72SBarry Smith           ierr = PetscCDIntNdGetID(pos, &gidj);CHKERRQ(ierr);
597e78576d6SMark F. Adams           ierr = PetscCDGetNextPos(aggs_2,slid_new,&pos);CHKERRQ(ierr);
598e78576d6SMark F. Adams 
5990cbbd2e1SMark F. Adams           if (gidj == gid) { hav = 1; break; }
600c8b0795cSMark F. Adams         }
601c8b0795cSMark F. Adams         if (hav != 1) {
602ffc955d6SMark F. Adams           /* insert 'flidj' into head of llist */
60341b27cdeSMark F. Adams           ierr = PetscCDAppendID(aggs_2, slid_new, gid);CHKERRQ(ierr);
604c8b0795cSMark F. Adams         }
605c8b0795cSMark F. Adams       }
606c8b0795cSMark F. Adams     }
607c8b0795cSMark F. Adams 
6080cbbd2e1SMark F. Adams     ierr = VecRestoreArray(mpimat_1->lvec, &cpcol_1_state);CHKERRQ(ierr);
6090cbbd2e1SMark F. Adams     ierr = VecRestoreArray(mpimat_2->lvec, &cpcol_2_state);CHKERRQ(ierr);
6100cbbd2e1SMark F. Adams     ierr = VecRestoreArray(ghostparents2, &cpcol_2_parent);CHKERRQ(ierr);
6110cbbd2e1SMark F. Adams     ierr = VecRestoreArray(ghostgids2, &cpcol_2_gid);CHKERRQ(ierr);
612c8b0795cSMark F. Adams     ierr = PetscFree(lid_cprowID_1);CHKERRQ(ierr);
6130cbbd2e1SMark F. Adams     ierr = VecDestroy(&ghostgids2);CHKERRQ(ierr);
6140cbbd2e1SMark F. Adams     ierr = VecDestroy(&ghostparents2);CHKERRQ(ierr);
6150cbbd2e1SMark F. Adams     ierr = VecDestroy(&ghost_par_orig2);CHKERRQ(ierr);
616c8b0795cSMark F. Adams   }
617c8b0795cSMark F. Adams 
618e632b94dSBarry Smith   ierr = PetscFree2(lid_state,lid_parent_gid);CHKERRQ(ierr);
619c8b0795cSMark F. Adams   PetscFunctionReturn(0);
620c8b0795cSMark F. Adams }
6212e68589bSMark F. Adams 
6222e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
6232e68589bSMark F. Adams /*
624a2f3521dSMark F. Adams    PCSetData_AGG - called if data is not set with PCSetCoordinates.
625a2f3521dSMark F. Adams       Looks in Mat for near null space.
626a2f3521dSMark F. Adams       Does not work for Stokes
6272e68589bSMark F. Adams 
6282e68589bSMark F. Adams   Input Parameter:
6292e68589bSMark F. Adams    . pc -
630a2f3521dSMark F. Adams    . a_A - matrix to get (near) null space out of.
6312e68589bSMark F. Adams */
6322e68589bSMark F. Adams #undef __FUNCT__
6332e68589bSMark F. Adams #define __FUNCT__ "PCSetData_AGG"
634e0877f53SBarry Smith static PetscErrorCode PCSetData_AGG(PC pc, Mat a_A)
6352e68589bSMark F. Adams {
6362e68589bSMark F. Adams   PetscErrorCode ierr;
637b8cd405aSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
638b8cd405aSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
639b8cd405aSMark F. Adams   MatNullSpace   mnull;
6402e68589bSMark F. Adams   PetscFunctionBegin;
64166f2ef4bSMark Adams 
642b8cd405aSMark F. Adams   ierr = MatGetNearNullSpace(a_A, &mnull);CHKERRQ(ierr);
643b8cd405aSMark F. Adams   if (!mnull) {
64466f2ef4bSMark Adams     DM dm;
64566f2ef4bSMark Adams     ierr = PCGetDM(pc, &dm);CHKERRQ(ierr);
64666f2ef4bSMark Adams     if (!dm) {
64766f2ef4bSMark Adams       ierr = MatGetDM(a_A, &dm);CHKERRQ(ierr);
64866f2ef4bSMark Adams     }
64966f2ef4bSMark Adams     if (dm) {
65066f2ef4bSMark Adams       PetscObject deformation;
651b0d30dd6SMatthew G. Knepley       PetscInt    Nf;
652b0d30dd6SMatthew G. Knepley 
653b0d30dd6SMatthew G. Knepley       ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr);
654b0d30dd6SMatthew G. Knepley       if (Nf) {
65566f2ef4bSMark Adams         ierr = DMGetField(dm, 0, &deformation);CHKERRQ(ierr);
65666f2ef4bSMark Adams         ierr = PetscObjectQuery((PetscObject)deformation,"nearnullspace",(PetscObject*)&mnull);CHKERRQ(ierr);
65766f2ef4bSMark Adams         if (!mnull) {
65866f2ef4bSMark Adams           ierr = PetscObjectQuery((PetscObject)deformation,"nullspace",(PetscObject*)&mnull);CHKERRQ(ierr);
65966f2ef4bSMark Adams         }
66066f2ef4bSMark Adams       }
66166f2ef4bSMark Adams     }
662b0d30dd6SMatthew G. Knepley   }
66366f2ef4bSMark Adams 
66466f2ef4bSMark Adams   if (!mnull) {
665a2f3521dSMark F. Adams     PetscInt bs,NN,MM;
6669d1b15c3SMark F. Adams     ierr = MatGetBlockSize(a_A, &bs);CHKERRQ(ierr);
66771959b99SBarry Smith     ierr = MatGetLocalSize(a_A, &MM, &NN);CHKERRQ(ierr);
66871959b99SBarry Smith     if (MM % bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MM %D must be divisible by bs %D",MM,bs);
6690298fd71SBarry Smith     ierr = PCSetCoordinates_AGG(pc, bs, MM/bs, NULL);CHKERRQ(ierr);
670806fa848SBarry Smith   } else {
671b8cd405aSMark F. Adams     PetscReal *nullvec;
672b8cd405aSMark F. Adams     PetscBool has_const;
673b8cd405aSMark F. Adams     PetscInt  i,j,mlocal,nvec,bs;
674b8cd405aSMark F. Adams     const Vec *vecs; const PetscScalar *v;
6750298fd71SBarry Smith     ierr = MatGetLocalSize(a_A,&mlocal,NULL);CHKERRQ(ierr);
676b8cd405aSMark F. Adams     ierr = MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs);CHKERRQ(ierr);
677a0dea87dSMark Adams     pc_gamg->data_sz = (nvec+!!has_const)*mlocal;
678785e854fSJed Brown     ierr = PetscMalloc1((nvec+!!has_const)*mlocal,&nullvec);CHKERRQ(ierr);
679b8cd405aSMark F. Adams     if (has_const) for (i=0; i<mlocal; i++) nullvec[i] = 1.0;
680b8cd405aSMark F. Adams     for (i=0; i<nvec; i++) {
681b8cd405aSMark F. Adams       ierr = VecGetArrayRead(vecs[i],&v);CHKERRQ(ierr);
682b8cd405aSMark F. Adams       for (j=0; j<mlocal; j++) nullvec[(i+!!has_const)*mlocal + j] = PetscRealPart(v[j]);
683b8cd405aSMark F. Adams       ierr = VecRestoreArrayRead(vecs[i],&v);CHKERRQ(ierr);
684b8cd405aSMark F. Adams     }
685b8cd405aSMark F. Adams     pc_gamg->data           = nullvec;
686b8cd405aSMark F. Adams     pc_gamg->data_cell_cols = (nvec+!!has_const);
6872fa5cd67SKarl Rupp 
68806e133e7SMark Adams     ierr = MatGetBlockSize(a_A, &bs);CHKERRQ(ierr);
6892fa5cd67SKarl Rupp 
690b8cd405aSMark F. Adams     pc_gamg->data_cell_rows = bs;
691b8cd405aSMark F. Adams   }
6922e68589bSMark F. Adams   PetscFunctionReturn(0);
6932e68589bSMark F. Adams }
6942e68589bSMark F. Adams 
6952e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
6962e68589bSMark F. Adams /*
6972e68589bSMark F. Adams  formProl0
6982e68589bSMark F. Adams 
6992e68589bSMark F. Adams    Input Parameter:
7002adfe9d3SBarry Smith    . agg_llists - list of arrays with aggregates -- list from selected vertices of aggregate unselected vertices
7012adfe9d3SBarry Smith    . bs - row block size
7020cbbd2e1SMark F. Adams    . nSAvec - column bs of new P
7030cbbd2e1SMark F. Adams    . my0crs - global index of start of locals
7042adfe9d3SBarry Smith    . data_stride - bs*(nloc nodes + ghost nodes) [data_stride][nSAvec]
7052e68589bSMark F. Adams    . data_in[data_stride*nSAvec] - local data on fine grid
7062e68589bSMark F. Adams    . flid_fgid[data_stride/bs] - make local to global IDs, includes ghosts in 'locals_llist'
7072e68589bSMark F. Adams   Output Parameter:
7082e68589bSMark F. Adams    . a_data_out - in with fine grid data (w/ghosts), out with coarse grid data
7092e68589bSMark F. Adams    . a_Prol - prolongation operator
7102e68589bSMark F. Adams */
7112e68589bSMark F. Adams #undef __FUNCT__
7122e68589bSMark F. Adams #define __FUNCT__ "formProl0"
7132adfe9d3SBarry 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)
7142e68589bSMark F. Adams {
7152e68589bSMark F. Adams   PetscErrorCode  ierr;
716ac187d40SBarry Smith   PetscInt        Istart,my0,Iend,nloc,clid,flid = 0,aggID,kk,jj,ii,mm,ndone,nSelected,minsz,nghosts,out_data_stride;
7173b4367a7SBarry Smith   MPI_Comm        comm;
71873911c69SBarry Smith   PetscMPIInt     rank;
7192e68589bSMark F. Adams   PetscReal       *out_data;
72041b27cdeSMark F. Adams   PetscCDPos      pos;
7211943db53SBarry Smith   PCGAMGHashTable fgid_flid;
7220cbbd2e1SMark F. Adams 
723797e13b7SMark F. Adams /* #define OUT_AGGS */
724519f805aSKarl Rupp #if defined(OUT_AGGS)
7250298fd71SBarry Smith   static PetscInt llev = 0; char fname[32]; FILE *file = NULL; PetscInt pM;
7269057884aSMark F. Adams #endif
7272e68589bSMark F. Adams 
7282e68589bSMark F. Adams   PetscFunctionBegin;
7293b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)a_Prol,&comm);CHKERRQ(ierr);
7303b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
7312e68589bSMark F. Adams   ierr = MatGetOwnershipRange(a_Prol, &Istart, &Iend);CHKERRQ(ierr);
73271959b99SBarry Smith   nloc = (Iend-Istart)/bs; my0 = Istart/bs;
73371959b99SBarry 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);
7340cbbd2e1SMark F. Adams   Iend   /= bs;
7350cbbd2e1SMark F. Adams   nghosts = data_stride/bs - nloc;
7362e68589bSMark F. Adams 
7371943db53SBarry Smith   ierr = PCGAMGHashTableCreate(2*nghosts+1, &fgid_flid);CHKERRQ(ierr);
7380cbbd2e1SMark F. Adams   for (kk=0; kk<nghosts; kk++) {
7391943db53SBarry Smith     ierr = PCGAMGHashTableAdd(&fgid_flid, flid_fgid[nloc+kk], nloc+kk);CHKERRQ(ierr);
7402e68589bSMark F. Adams   }
7412e68589bSMark F. Adams 
742519f805aSKarl Rupp #if defined(OUT_AGGS)
743c5df96a5SBarry Smith   sprintf(fname,"aggs_%d_%d.m",llev++,rank);
7442fa5cd67SKarl Rupp   if (llev==1) file = fopen(fname,"w");
7450cbbd2e1SMark F. Adams   MatGetSize(a_Prol, &pM, &jj);
7460cbbd2e1SMark F. Adams #endif
7470cbbd2e1SMark F. Adams 
7480cbbd2e1SMark F. Adams   /* count selected -- same as number of cols of P */
7490cbbd2e1SMark F. Adams   for (nSelected=mm=0; mm<nloc; mm++) {
750e78576d6SMark F. Adams     PetscBool ise;
751e78576d6SMark F. Adams     ierr = PetscCDEmptyAt(agg_llists, mm, &ise);CHKERRQ(ierr);
752e78576d6SMark F. Adams     if (!ise) nSelected++;
7530cbbd2e1SMark F. Adams   }
7540cbbd2e1SMark F. Adams   ierr = MatGetOwnershipRangeColumn(a_Prol, &ii, &jj);CHKERRQ(ierr);
75571959b99SBarry Smith   if ((ii/nSAvec) != my0crs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ii %D /nSAvec %D  != my0crs %D",ii,nSAvec,my0crs);
75671959b99SBarry 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);
7570cbbd2e1SMark F. Adams 
7582e68589bSMark F. Adams   /* aloc space for coarse point data (output) */
7590cbbd2e1SMark F. Adams   out_data_stride = nSelected*nSAvec;
7602fa5cd67SKarl Rupp 
761785e854fSJed Brown   ierr = PetscMalloc1(out_data_stride*nSAvec, &out_data);CHKERRQ(ierr);
76233812677SSatish Balay   for (ii=0;ii<out_data_stride*nSAvec;ii++) out_data[ii]=PETSC_MAX_REAL;
7630cbbd2e1SMark F. Adams   *a_data_out = out_data; /* output - stride nSelected*nSAvec */
7642e68589bSMark F. Adams 
7652e68589bSMark F. Adams   /* find points and set prolongation */
766c8b0795cSMark F. Adams   minsz = 100;
7672e68589bSMark F. Adams   ndone = 0;
7680cbbd2e1SMark F. Adams   for (mm = clid = 0; mm < nloc; mm++) {
769e78576d6SMark F. Adams     ierr = PetscCDSizeAt(agg_llists, mm, &jj);CHKERRQ(ierr);
770e78576d6SMark F. Adams     if (jj > 0) {
7710cbbd2e1SMark F. Adams       const PetscInt lid = mm, cgid = my0crs + clid;
7720cbbd2e1SMark F. Adams       PetscInt       cids[100]; /* max bs */
7730cbbd2e1SMark F. Adams       PetscBLASInt   asz  =jj,M=asz*bs,N=nSAvec,INFO;
7742e68589bSMark F. Adams       PetscBLASInt   Mdata=M+((N-M>0) ? N-M : 0),LDA=Mdata,LWORK=N*bs;
7752e68589bSMark F. Adams       PetscScalar    *qqc,*qqr,*TAU,*WORK;
7762e68589bSMark F. Adams       PetscInt       *fids;
77765d7b583SSatish Balay       PetscReal      *data;
7780cbbd2e1SMark F. Adams       /* count agg */
7790cbbd2e1SMark F. Adams       if (asz<minsz) minsz = asz;
7800cbbd2e1SMark F. Adams 
7810cbbd2e1SMark F. Adams       /* get block */
782e632b94dSBarry Smith       ierr = PetscMalloc5(Mdata*N, &qqc,M*N, &qqr,N, &TAU,LWORK, &WORK,M, &fids);CHKERRQ(ierr);
7832e68589bSMark F. Adams 
7842e68589bSMark F. Adams       aggID = 0;
785e78576d6SMark F. Adams       ierr  = PetscCDGetHeadPos(agg_llists,lid,&pos);CHKERRQ(ierr);
786e78576d6SMark F. Adams       while (pos) {
787e78576d6SMark F. Adams         PetscInt gid1;
788*484f0a72SBarry Smith         ierr = PetscCDIntNdGetID(pos, &gid1);CHKERRQ(ierr);
789e78576d6SMark F. Adams         ierr = PetscCDGetNextPos(agg_llists,lid,&pos);CHKERRQ(ierr);
790e78576d6SMark F. Adams 
7910cbbd2e1SMark F. Adams         if (gid1 >= my0 && gid1 < Iend) flid = gid1 - my0;
7920cbbd2e1SMark F. Adams         else {
7931943db53SBarry Smith           ierr = PCGAMGHashTableFind(&fgid_flid, gid1, &flid);CHKERRQ(ierr);
79471959b99SBarry Smith           if (flid < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot find gid1 in table");
7950cbbd2e1SMark F. Adams         }
7962e68589bSMark F. Adams         /* copy in B_i matrix - column oriented */
79765d7b583SSatish Balay         data = &data_in[flid*bs];
7983d3eaba7SBarry Smith         for (ii = 0; ii < bs; ii++) {
7992e68589bSMark F. Adams           for (jj = 0; jj < N; jj++) {
8000cbbd2e1SMark F. Adams             PetscReal d = data[jj*data_stride + ii];
8010cbbd2e1SMark F. Adams             qqc[jj*Mdata + aggID*bs + ii] = d;
8022e68589bSMark F. Adams           }
8032e68589bSMark F. Adams         }
804519f805aSKarl Rupp #if defined(OUT_AGGS)
805b2a4f308SMark F. Adams         if (llev==1) {
8069057884aSMark F. Adams           char     str[] = "plot(%e,%e,'r*'), hold on,\n", col[] = "rgbkmc", sim[] = "*os+h>d<vx^";
8070cbbd2e1SMark F. Adams           PetscInt MM,pi,pj;
808c5df96a5SBarry Smith           str[12] = col[(clid+17*rank)%6]; str[13] = sim[(clid+17*rank)%11];
809f7620de1SMatthew G Knepley           MM      = (PetscInt)(PetscSqrtReal((PetscReal)pM));
8100cbbd2e1SMark F. Adams           pj      = gid1/MM; pi = gid1%MM;
811b2a4f308SMark F. Adams           fprintf(file,str,(double)pi,(double)pj);
812b2a4f308SMark F. Adams           /* fprintf(file,str,data[2*data_stride+1],-data[2*data_stride]); */
8139057884aSMark F. Adams         }
8149057884aSMark F. Adams #endif
8152e68589bSMark F. Adams         /* set fine IDs */
8162e68589bSMark F. Adams         for (kk=0; kk<bs; kk++) fids[aggID*bs + kk] = flid_fgid[flid]*bs + kk;
8172e68589bSMark F. Adams 
8182e68589bSMark F. Adams         aggID++;
8190cbbd2e1SMark F. Adams       }
8202e68589bSMark F. Adams 
8212e68589bSMark F. Adams       /* pad with zeros */
8222e68589bSMark F. Adams       for (ii = asz*bs; ii < Mdata; ii++) {
8232e68589bSMark F. Adams         for (jj = 0; jj < N; jj++, kk++) {
8242e68589bSMark F. Adams           qqc[jj*Mdata + ii] = .0;
8252e68589bSMark F. Adams         }
8262e68589bSMark F. Adams       }
8272e68589bSMark F. Adams 
8282e68589bSMark F. Adams       ndone += aggID;
8292e68589bSMark F. Adams       /* QR */
83084df3f90SSatish Balay       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
8318b83055fSJed Brown       PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&Mdata, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
83284df3f90SSatish Balay       ierr = PetscFPTrapPop();CHKERRQ(ierr);
833d23427c4SJed Brown       if (INFO != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xGEQRF error");
8342e68589bSMark F. Adams       /* get R - column oriented - output B_{i+1} */
8352e68589bSMark F. Adams       {
8362e68589bSMark F. Adams         PetscReal *data = &out_data[clid*nSAvec];
8372e68589bSMark F. Adams         for (jj = 0; jj < nSAvec; jj++) {
8382e68589bSMark F. Adams           for (ii = 0; ii < nSAvec; ii++) {
83933812677SSatish 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);
8400cbbd2e1SMark F. Adams            if (ii <= jj) data[jj*out_data_stride + ii] = PetscRealPart(qqc[jj*Mdata + ii]);
8410cbbd2e1SMark F. Adams            else data[jj*out_data_stride + ii] = 0.;
8422e68589bSMark F. Adams           }
8432e68589bSMark F. Adams         }
8442e68589bSMark F. Adams       }
8452e68589bSMark F. Adams 
8462e68589bSMark F. Adams       /* get Q - row oriented */
8478b83055fSJed Brown       PetscStackCallBLAS("LAPACKungqr",LAPACKungqr_(&Mdata, &N, &N, qqc, &LDA, TAU, WORK, &LWORK, &INFO));
848c666adbfSMark F. Adams       if (INFO != 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xORGQR error arg %d",-INFO);
8492e68589bSMark F. Adams 
8502e68589bSMark F. Adams       for (ii = 0; ii < M; ii++) {
8512e68589bSMark F. Adams         for (jj = 0; jj < N; jj++) {
8522e68589bSMark F. Adams           qqr[N*ii + jj] = qqc[jj*Mdata + ii];
8532e68589bSMark F. Adams         }
8542e68589bSMark F. Adams       }
8552e68589bSMark F. Adams 
8562e68589bSMark F. Adams       /* add diagonal block of P0 */
857c8b0795cSMark F. Adams       for (kk=0; kk<N; kk++) {
858c8b0795cSMark F. Adams         cids[kk] = N*cgid + kk; /* global col IDs in P0 */
859c8b0795cSMark F. Adams       }
8602e68589bSMark F. Adams       ierr = MatSetValues(a_Prol,M,fids,N,cids,qqr,INSERT_VALUES);CHKERRQ(ierr);
8612e68589bSMark F. Adams 
862e632b94dSBarry Smith       ierr = PetscFree5(qqc,qqr,TAU,WORK,fids);CHKERRQ(ierr);
863b43b03e9SMark F. Adams       clid++;
8640cbbd2e1SMark F. Adams     } /* coarse agg */
8650cbbd2e1SMark F. Adams   } /* for all fine nodes */
8660cbbd2e1SMark F. Adams   ierr = MatAssemblyBegin(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8670cbbd2e1SMark F. Adams   ierr = MatAssemblyEnd(a_Prol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
8682e68589bSMark F. Adams 
869b2566f29SBarry Smith /* ierr = MPIU_Allreduce(&ndone, &ii, 1, MPIU_INT, MPIU_SUM, comm); */
8702e68589bSMark F. Adams /* MatGetSize(a_Prol, &kk, &jj); */
871b2566f29SBarry Smith /* ierr = MPIU_Allreduce(&minsz, &jj, 1, MPIU_INT, MPIU_MIN, comm); */
8723b4367a7SBarry 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); */
8732e68589bSMark F. Adams 
874519f805aSKarl Rupp #if defined(OUT_AGGS)
875b2a4f308SMark F. Adams   if (llev==1) fclose(file);
8769057884aSMark F. Adams #endif
8771943db53SBarry Smith   ierr = PCGAMGHashTableDestroy(&fgid_flid);CHKERRQ(ierr);
8782e68589bSMark F. Adams   PetscFunctionReturn(0);
8792e68589bSMark F. Adams }
8802e68589bSMark F. Adams 
8815adeb434SBarry Smith #undef __FUNCT__
8825adeb434SBarry Smith #define __FUNCT__ "PCView_GAMG_AGG"
8835adeb434SBarry Smith static PetscErrorCode PCView_GAMG_AGG(PC pc,PetscViewer viewer)
8845adeb434SBarry Smith {
8855adeb434SBarry Smith   PetscErrorCode ierr;
8865adeb434SBarry Smith   PC_MG          *mg      = (PC_MG*)pc->data;
8875adeb434SBarry Smith   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
8885adeb434SBarry Smith   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
8895adeb434SBarry Smith 
8905adeb434SBarry Smith   PetscFunctionBegin;
8915adeb434SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"      AGG specific options\n");CHKERRQ(ierr);
8925adeb434SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"        Symmetric graph %s\n",pc_gamg_agg->sym_graph ? "true" : "false");CHKERRQ(ierr);
8935adeb434SBarry Smith   PetscFunctionReturn(0);
8945adeb434SBarry Smith }
8955adeb434SBarry Smith 
8962e68589bSMark F. Adams /* -------------------------------------------------------------------------- */
8972e68589bSMark F. Adams /*
898fd1112cbSBarry Smith    PCGAMGGraph_AGG
8992e68589bSMark F. Adams 
9002e68589bSMark F. Adams   Input Parameter:
9012e68589bSMark F. Adams    . pc - this
9022e68589bSMark F. Adams    . Amat - matrix on this fine level
9032e68589bSMark F. Adams   Output Parameter:
904c8b0795cSMark F. Adams    . a_Gmat -
9052e68589bSMark F. Adams */
9062e68589bSMark F. Adams #undef __FUNCT__
907fd1112cbSBarry Smith #define __FUNCT__ "PCGAMGGraph_AGG"
908e0877f53SBarry Smith static PetscErrorCode PCGAMGGraph_AGG(PC pc,Mat Amat,Mat *a_Gmat)
909c8b0795cSMark F. Adams {
910c8b0795cSMark F. Adams   PetscErrorCode            ierr;
911c8b0795cSMark F. Adams   PC_MG                     *mg          = (PC_MG*)pc->data;
912c8b0795cSMark F. Adams   PC_GAMG                   *pc_gamg     = (PC_GAMG*)mg->innerctx;
913c8b0795cSMark F. Adams   const PetscReal           vfilter      = pc_gamg->threshold;
914c8b0795cSMark F. Adams   PC_GAMG_AGG               *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
915e0940f08SMark F. Adams   Mat                       Gmat;
9163b4367a7SBarry Smith   MPI_Comm                  comm;
91767744fedSMark F. Adams   PetscBool /* set,flg , */ symm;
918c8b0795cSMark F. Adams 
919c8b0795cSMark F. Adams   PetscFunctionBegin;
9203b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr);
921fd1112cbSBarry Smith   ierr = PetscLogEventBegin(PC_GAMGGraph_AGG,0,0,0,0);CHKERRQ(ierr);
922c8b0795cSMark F. Adams 
92367744fedSMark F. Adams   /* ierr = MatIsSymmetricKnown(Amat, &set, &flg);CHKERRQ(ierr); || !(set && flg) -- this causes lot of symm calls */
924c666adbfSMark F. Adams   symm = (PetscBool)(pc_gamg_agg->sym_graph); /* && !pc_gamg_agg->square_graph; */
9250cbbd2e1SMark F. Adams 
9262d7fac45SMark F. Adams   ierr = PCGAMGCreateGraph(Amat, &Gmat);CHKERRQ(ierr);
927bf4339c2SBarry Smith   ierr = PCGAMGFilterGraph(&Gmat, vfilter, symm);CHKERRQ(ierr);
928c8b0795cSMark F. Adams 
929e0940f08SMark F. Adams   *a_Gmat = Gmat;
930c8b0795cSMark F. Adams 
931fd1112cbSBarry Smith   ierr = PetscLogEventEnd(PC_GAMGGraph_AGG,0,0,0,0);CHKERRQ(ierr);
932c8b0795cSMark F. Adams   PetscFunctionReturn(0);
933c8b0795cSMark F. Adams }
934c8b0795cSMark F. Adams 
935c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
936c8b0795cSMark F. Adams /*
937b43b03e9SMark F. Adams    PCGAMGCoarsen_AGG
938c8b0795cSMark F. Adams 
939c8b0795cSMark F. Adams   Input Parameter:
940e0940f08SMark F. Adams    . a_pc - this
941e0940f08SMark F. Adams   Input/Output Parameter:
9420cbbd2e1SMark F. Adams    . a_Gmat1 - graph on this fine level - coarsening can change this (squares it)
943c8b0795cSMark F. Adams   Output Parameter:
9440cbbd2e1SMark F. Adams    . agg_lists - list of aggregates
945c8b0795cSMark F. Adams */
946c8b0795cSMark F. Adams #undef __FUNCT__
947b43b03e9SMark F. Adams #define __FUNCT__ "PCGAMGCoarsen_AGG"
948e0877f53SBarry Smith static PetscErrorCode PCGAMGCoarsen_AGG(PC a_pc,Mat *a_Gmat1,PetscCoarsenData **agg_lists)
949c8b0795cSMark F. Adams {
950c8b0795cSMark F. Adams   PetscErrorCode ierr;
951e0940f08SMark F. Adams   PC_MG          *mg          = (PC_MG*)a_pc->data;
952c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
953c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
9540cbbd2e1SMark F. Adams   Mat            mat,Gmat2, Gmat1 = *a_Gmat1;  /* squared graph */
9550cbbd2e1SMark F. Adams   IS             perm;
9565cfd4bd9SMark Adams   PetscInt       Istart,Iend,Ii,nloc,bs,n,m;
957c8b0795cSMark F. Adams   PetscInt       *permute;
958c8b0795cSMark F. Adams   PetscBool      *bIndexSet;
959b43b03e9SMark F. Adams   MatCoarsen     crs;
9603b4367a7SBarry Smith   MPI_Comm       comm;
96173911c69SBarry Smith   PetscMPIInt    rank;
962e2c00dcbSBarry Smith   PetscReal      rr;
963e2c00dcbSBarry Smith   PetscInt       iSwapIndex;
964c8b0795cSMark F. Adams 
965c8b0795cSMark F. Adams   PetscFunctionBegin;
9660cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr);
9673b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Gmat1,&comm);CHKERRQ(ierr);
9683b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
969e0940f08SMark F. Adams   ierr = MatGetLocalSize(Gmat1, &n, &m);CHKERRQ(ierr);
97071959b99SBarry Smith   ierr = MatGetBlockSize(Gmat1, &bs);CHKERRQ(ierr);
97171959b99SBarry Smith   if (bs != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"bs %D must be 1",bs);
972c8b0795cSMark F. Adams   nloc = n/bs;
973c8b0795cSMark F. Adams 
9749ab59c8bSMark Adams   if (pc_gamg->current_level < pc_gamg_agg->square_graph) {
975e2c00dcbSBarry 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);
976806fa848SBarry Smith     ierr = MatTransposeMatMult(Gmat1, Gmat1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Gmat2);CHKERRQ(ierr);
977806fa848SBarry Smith   } else Gmat2 = Gmat1;
978c8b0795cSMark F. Adams 
9795cfd4bd9SMark Adams   /* get MIS aggs - randomize */
980785e854fSJed Brown   ierr = PetscMalloc1(nloc, &permute);CHKERRQ(ierr);
981e2c00dcbSBarry Smith   ierr = PetscCalloc1(nloc, &bIndexSet);CHKERRQ(ierr);
982c8b0795cSMark F. Adams   for (Ii = 0; Ii < nloc; Ii++) {
983c8b0795cSMark F. Adams     permute[Ii]   = Ii;
984c8b0795cSMark F. Adams   }
9855cfd4bd9SMark Adams   ierr = MatGetOwnershipRange(Gmat1, &Istart, &Iend);CHKERRQ(ierr);
986c8b0795cSMark F. Adams   for (Ii = 0; Ii < nloc; Ii++) {
98714a9496bSBarry Smith     ierr = PetscRandomGetValueReal(pc_gamg->random,&rr);CHKERRQ(ierr);
988e2c00dcbSBarry Smith     iSwapIndex = (PetscInt) (rr*nloc);
989c8b0795cSMark F. Adams     if (!bIndexSet[iSwapIndex] && iSwapIndex != Ii) {
990c8b0795cSMark F. Adams       PetscInt iTemp = permute[iSwapIndex];
991c8b0795cSMark F. Adams       permute[iSwapIndex]   = permute[Ii];
992c8b0795cSMark F. Adams       permute[Ii]           = iTemp;
993c8b0795cSMark F. Adams       bIndexSet[iSwapIndex] = PETSC_TRUE;
994c8b0795cSMark F. Adams     }
995c8b0795cSMark F. Adams   }
996c8b0795cSMark F. Adams   ierr = PetscFree(bIndexSet);CHKERRQ(ierr);
997c8b0795cSMark F. Adams 
998806fa848SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF, nloc, permute, PETSC_USE_POINTER, &perm);CHKERRQ(ierr);
9990cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
10000cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr);
1001b43b03e9SMark F. Adams #endif
10023b4367a7SBarry Smith   ierr = MatCoarsenCreate(comm, &crs);CHKERRQ(ierr);
10039057884aSMark F. Adams   ierr = MatCoarsenSetFromOptions(crs);CHKERRQ(ierr);
1004b43b03e9SMark F. Adams   ierr = MatCoarsenSetGreedyOrdering(crs, perm);CHKERRQ(ierr);
1005b43b03e9SMark F. Adams   ierr = MatCoarsenSetAdjacency(crs, Gmat2);CHKERRQ(ierr);
1006b43b03e9SMark F. Adams   ierr = MatCoarsenSetStrictAggs(crs, PETSC_TRUE);CHKERRQ(ierr);
1007b43b03e9SMark F. Adams   ierr = MatCoarsenApply(crs);CHKERRQ(ierr);
10080cbbd2e1SMark F. Adams   ierr = MatCoarsenGetData(crs, agg_lists);CHKERRQ(ierr); /* output */
1009b43b03e9SMark F. Adams   ierr = MatCoarsenDestroy(&crs);CHKERRQ(ierr);
1010b43b03e9SMark F. Adams 
1011c8b0795cSMark F. Adams   ierr = ISDestroy(&perm);CHKERRQ(ierr);
1012c8b0795cSMark F. Adams   ierr = PetscFree(permute);CHKERRQ(ierr);
10130cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
10140cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET4],0,0,0,0);CHKERRQ(ierr);
1015b43b03e9SMark F. Adams #endif
10169f3f12c8SMark F. Adams 
1017c8b0795cSMark F. Adams   /* smooth aggs */
1018e0940f08SMark F. Adams   if (Gmat2 != Gmat1) {
10190cbbd2e1SMark F. Adams     const PetscCoarsenData *llist = *agg_lists;
10200cbbd2e1SMark F. Adams     ierr     = smoothAggs(Gmat2, Gmat1, *agg_lists);CHKERRQ(ierr);
1021c8b0795cSMark F. Adams     ierr     = MatDestroy(&Gmat1);CHKERRQ(ierr);
1022e0940f08SMark F. Adams     *a_Gmat1 = Gmat2; /* output */
102341b27cdeSMark F. Adams     ierr     = PetscCDGetMat(llist, &mat);CHKERRQ(ierr);
10243b4367a7SBarry Smith     if (mat) SETERRQ(comm,PETSC_ERR_ARG_WRONG, "Auxilary matrix with squared graph????");
1025806fa848SBarry Smith   } else {
10260cbbd2e1SMark F. Adams     const PetscCoarsenData *llist = *agg_lists;
10279ab59c8bSMark Adams     /* see if we have a matrix that takes precedence (returned from MatCoarsenApply) */
102841b27cdeSMark F. Adams     ierr = PetscCDGetMat(llist, &mat);CHKERRQ(ierr);
10290cbbd2e1SMark F. Adams     if (mat) {
10300cbbd2e1SMark F. Adams       ierr     = MatDestroy(&Gmat1);CHKERRQ(ierr);
10310cbbd2e1SMark F. Adams       *a_Gmat1 = mat; /* output */
10320cbbd2e1SMark F. Adams     }
10330cbbd2e1SMark F. Adams   }
10340cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGCoarsen_AGG,0,0,0,0);CHKERRQ(ierr);
1035c8b0795cSMark F. Adams   PetscFunctionReturn(0);
1036c8b0795cSMark F. Adams }
1037c8b0795cSMark F. Adams 
1038c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
1039c8b0795cSMark F. Adams /*
10400cbbd2e1SMark F. Adams  PCGAMGProlongator_AGG
1041c8b0795cSMark F. Adams 
1042c8b0795cSMark F. Adams  Input Parameter:
1043c8b0795cSMark F. Adams  . pc - this
1044c8b0795cSMark F. Adams  . Amat - matrix on this fine level
1045c8b0795cSMark F. Adams  . Graph - used to get ghost data for nodes in
10460cbbd2e1SMark F. Adams  . agg_lists - list of aggregates
1047c8b0795cSMark F. Adams  Output Parameter:
1048c8b0795cSMark F. Adams  . a_P_out - prolongation operator to the next level
1049c8b0795cSMark F. Adams  */
1050c8b0795cSMark F. Adams #undef __FUNCT__
10510cbbd2e1SMark F. Adams #define __FUNCT__ "PCGAMGProlongator_AGG"
1052e0877f53SBarry Smith static PetscErrorCode PCGAMGProlongator_AGG(PC pc,Mat Amat,Mat Gmat,PetscCoarsenData *agg_lists,Mat *a_P_out)
10532e68589bSMark F. Adams {
10542e68589bSMark F. Adams   PC_MG          *mg       = (PC_MG*)pc->data;
10552e68589bSMark F. Adams   PC_GAMG        *pc_gamg  = (PC_GAMG*)mg->innerctx;
10564a2f8832SBarry Smith   const PetscInt col_bs = pc_gamg->data_cell_cols;
10572e68589bSMark F. Adams   PetscErrorCode ierr;
1058c8b0795cSMark F. Adams   PetscInt       Istart,Iend,nloc,ii,jj,kk,my0,nLocalSelected,bs;
1059c8b0795cSMark F. Adams   Mat            Prol;
1060c5df96a5SBarry Smith   PetscMPIInt    rank, size;
10613b4367a7SBarry Smith   MPI_Comm       comm;
1062c8b0795cSMark F. Adams   PetscReal      *data_w_ghost;
1063c8b0795cSMark F. Adams   PetscInt       myCrs0, nbnodes=0, *flid_fgid;
1064d9558ea9SBarry Smith   MatType        mtype;
10652e68589bSMark F. Adams 
10662e68589bSMark F. Adams   PetscFunctionBegin;
10673b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr);
10684a2f8832SBarry Smith   if (col_bs < 1) SETERRQ(comm,PETSC_ERR_PLIB,"Column bs cannot be less than 1");
10690cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
10703b4367a7SBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
10713b4367a7SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
10722e68589bSMark F. Adams   ierr = MatGetOwnershipRange(Amat, &Istart, &Iend);CHKERRQ(ierr);
1073c8b0795cSMark F. Adams   ierr = MatGetBlockSize(Amat, &bs);CHKERRQ(ierr);
107471959b99SBarry Smith   nloc = (Iend-Istart)/bs; my0 = Istart/bs;
107571959b99SBarry Smith   if ((Iend-Istart) % bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"(Iend %D - Istart %D) not divisible by bs %D",Iend,Istart,bs);
10762e68589bSMark F. Adams 
10772e68589bSMark F. Adams   /* get 'nLocalSelected' */
10780cbbd2e1SMark F. Adams   for (ii=0, nLocalSelected = 0; ii < nloc; ii++) {
1079e78576d6SMark F. Adams     PetscBool ise;
10800cbbd2e1SMark F. Adams     /* filter out singletons 0 or 1? */
1081e78576d6SMark F. Adams     ierr = PetscCDEmptyAt(agg_lists, ii, &ise);CHKERRQ(ierr);
1082e78576d6SMark F. Adams     if (!ise) nLocalSelected++;
10832e68589bSMark F. Adams   }
10842e68589bSMark F. Adams 
10852e68589bSMark F. Adams   /* create prolongator, create P matrix */
1086d9558ea9SBarry Smith   ierr = MatGetType(Amat,&mtype);CHKERRQ(ierr);
10873b4367a7SBarry Smith   ierr = MatCreate(comm, &Prol);CHKERRQ(ierr);
1088806fa848SBarry Smith   ierr = MatSetSizes(Prol,nloc*bs,nLocalSelected*col_bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1089a2f3521dSMark F. Adams   ierr = MatSetBlockSizes(Prol, bs, col_bs);CHKERRQ(ierr);
1090d9558ea9SBarry Smith   ierr = MatSetType(Prol, mtype);CHKERRQ(ierr);
10914a2f8832SBarry Smith   ierr = MatSeqAIJSetPreallocation(Prol, col_bs, NULL);CHKERRQ(ierr);
10924a2f8832SBarry Smith   ierr = MatMPIAIJSetPreallocation(Prol,col_bs, NULL,col_bs, NULL);CHKERRQ(ierr);
10932e68589bSMark F. Adams 
10942e68589bSMark F. Adams   /* can get all points "removed" */
1095c8b0795cSMark F. Adams   ierr =  MatGetSize(Prol, &kk, &ii);CHKERRQ(ierr);
10967f66b68fSMark Adams   if (!ii) {
1097569f4572SMark Adams     ierr = PetscInfo(pc,"No selected points on coarse grid\n");CHKERRQ(ierr);
10982e68589bSMark F. Adams     ierr = MatDestroy(&Prol);CHKERRQ(ierr);
10990298fd71SBarry Smith     *a_P_out = NULL;  /* out */
1100e87675ddSBarry Smith     ierr = PetscLogEventEnd(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
11012e68589bSMark F. Adams     PetscFunctionReturn(0);
11022e68589bSMark F. Adams   }
1103bf4339c2SBarry Smith   ierr = PetscInfo1(pc,"New grid %D nodes\n",ii/col_bs);CHKERRQ(ierr);
1104c8b0795cSMark F. Adams   ierr = MatGetOwnershipRangeColumn(Prol, &myCrs0, &kk);CHKERRQ(ierr);
11050cbbd2e1SMark F. Adams 
110671959b99SBarry 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);
1107c8b0795cSMark F. Adams   myCrs0 = myCrs0/col_bs;
110871959b99SBarry 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);
11092e68589bSMark F. Adams 
11102e68589bSMark F. Adams   /* create global vector of data in 'data_w_ghost' */
11110cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
11120cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr);
11132e68589bSMark F. Adams #endif
1114c5df96a5SBarry Smith   if (size > 1) { /*  */
11152e68589bSMark F. Adams     PetscReal *tmp_gdata,*tmp_ldata,*tp2;
1116785e854fSJed Brown     ierr = PetscMalloc1(nloc, &tmp_ldata);CHKERRQ(ierr);
11174a2f8832SBarry Smith     for (jj = 0; jj < col_bs; jj++) {
1118c8b0795cSMark F. Adams       for (kk = 0; kk < bs; kk++) {
1119a2f3521dSMark F. Adams         PetscInt        ii,stride;
1120c8b0795cSMark F. Adams         const PetscReal *tp = pc_gamg->data + jj*bs*nloc + kk;
11212fa5cd67SKarl Rupp         for (ii = 0; ii < nloc; ii++, tp += bs) tmp_ldata[ii] = *tp;
11222fa5cd67SKarl Rupp 
1123806fa848SBarry Smith         ierr = PCGAMGGetDataWithGhosts(Gmat, 1, tmp_ldata, &stride, &tmp_gdata);CHKERRQ(ierr);
1124a2f3521dSMark F. Adams 
11257f66b68fSMark Adams         if (!jj && !kk) { /* now I know how many todal nodes - allocate */
11264a2f8832SBarry Smith           ierr    = PetscMalloc1(stride*bs*col_bs, &data_w_ghost);CHKERRQ(ierr);
1127a2f3521dSMark F. Adams           nbnodes = bs*stride;
11282e68589bSMark F. Adams         }
1129a2f3521dSMark F. Adams         tp2 = data_w_ghost + jj*bs*stride + kk;
11302fa5cd67SKarl Rupp         for (ii = 0; ii < stride; ii++, tp2 += bs) *tp2 = tmp_gdata[ii];
11312e68589bSMark F. Adams         ierr = PetscFree(tmp_gdata);CHKERRQ(ierr);
11322e68589bSMark F. Adams       }
11332e68589bSMark F. Adams     }
11342e68589bSMark F. Adams     ierr = PetscFree(tmp_ldata);CHKERRQ(ierr);
1135806fa848SBarry Smith   } else {
1136c8b0795cSMark F. Adams     nbnodes      = bs*nloc;
1137c8b0795cSMark F. Adams     data_w_ghost = (PetscReal*)pc_gamg->data;
11382e68589bSMark F. Adams   }
11392e68589bSMark F. Adams 
11402e68589bSMark F. Adams   /* get P0 */
1141c5df96a5SBarry Smith   if (size > 1) {
11422e68589bSMark F. Adams     PetscReal *fid_glid_loc,*fiddata;
1143a2f3521dSMark F. Adams     PetscInt  stride;
11442e68589bSMark F. Adams 
1145785e854fSJed Brown     ierr = PetscMalloc1(nloc, &fid_glid_loc);CHKERRQ(ierr);
11462e68589bSMark F. Adams     for (kk=0; kk<nloc; kk++) fid_glid_loc[kk] = (PetscReal)(my0+kk);
1147806fa848SBarry Smith     ierr = PCGAMGGetDataWithGhosts(Gmat, 1, fid_glid_loc, &stride, &fiddata);CHKERRQ(ierr);
1148785e854fSJed Brown     ierr = PetscMalloc1(stride, &flid_fgid);CHKERRQ(ierr);
1149a2f3521dSMark F. Adams     for (kk=0; kk<stride; kk++) flid_fgid[kk] = (PetscInt)fiddata[kk];
11502e68589bSMark F. Adams     ierr = PetscFree(fiddata);CHKERRQ(ierr);
1151a2f3521dSMark F. Adams 
115271959b99SBarry Smith     if (stride != nbnodes/bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"stride %D != nbnodes %D/bs %D",stride,nbnodes,bs);
11532e68589bSMark F. Adams     ierr = PetscFree(fid_glid_loc);CHKERRQ(ierr);
1154806fa848SBarry Smith   } else {
1155785e854fSJed Brown     ierr = PetscMalloc1(nloc, &flid_fgid);CHKERRQ(ierr);
11562e68589bSMark F. Adams     for (kk=0; kk<nloc; kk++) flid_fgid[kk] = my0 + kk;
11572e68589bSMark F. Adams   }
11580cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
11590cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET7],0,0,0,0);CHKERRQ(ierr);
11600cbbd2e1SMark F. Adams   ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr);
11612e68589bSMark F. Adams #endif
1162c8b0795cSMark F. Adams   {
11630298fd71SBarry Smith     PetscReal *data_out = NULL;
11644a2f8832SBarry Smith     ierr = formProl0(agg_lists, bs, col_bs, myCrs0, nbnodes,data_w_ghost, flid_fgid, &data_out, Prol);CHKERRQ(ierr);
1165c8b0795cSMark F. Adams     ierr = PetscFree(pc_gamg->data);CHKERRQ(ierr);
11662fa5cd67SKarl Rupp 
1167c8b0795cSMark F. Adams     pc_gamg->data           = data_out;
11684a2f8832SBarry Smith     pc_gamg->data_cell_rows = col_bs;
11694a2f8832SBarry Smith     pc_gamg->data_sz        = col_bs*col_bs*nLocalSelected;
1170c8b0795cSMark F. Adams   }
11710cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
11720cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET8],0,0,0,0);CHKERRQ(ierr);
1173c8b0795cSMark F. Adams #endif
117461ba4676SBarry Smith   if (size > 1) {ierr = PetscFree(data_w_ghost);CHKERRQ(ierr);}
11752e68589bSMark F. Adams   ierr = PetscFree(flid_fgid);CHKERRQ(ierr);
11762e68589bSMark F. Adams 
1177c8b0795cSMark F. Adams   *a_P_out = Prol;  /* out */
1178ed4e82eaSPeter Brune 
11790cbbd2e1SMark F. Adams   ierr = PetscLogEventEnd(PC_GAMGProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
1180c8b0795cSMark F. Adams   PetscFunctionReturn(0);
1181c8b0795cSMark F. Adams }
1182c8b0795cSMark F. Adams 
1183c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
1184c8b0795cSMark F. Adams /*
1185fd1112cbSBarry Smith    PCGAMGOptProlongator_AGG
1186c8b0795cSMark F. Adams 
1187c8b0795cSMark F. Adams   Input Parameter:
1188c8b0795cSMark F. Adams    . pc - this
1189c8b0795cSMark F. Adams    . Amat - matrix on this fine level
1190c8b0795cSMark F. Adams  In/Output Parameter:
11912a808120SBarry Smith    . a_P - prolongation operator to the next level
1192c8b0795cSMark F. Adams */
1193c8b0795cSMark F. Adams #undef __FUNCT__
1194fd1112cbSBarry Smith #define __FUNCT__ "PCGAMGOptProlongator_AGG"
1195e0877f53SBarry Smith static PetscErrorCode PCGAMGOptProlongator_AGG(PC pc,Mat Amat,Mat *a_P)
1196c8b0795cSMark F. Adams {
1197c8b0795cSMark F. Adams   PetscErrorCode ierr;
1198c8b0795cSMark F. Adams   PC_MG          *mg          = (PC_MG*)pc->data;
1199c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg     = (PC_GAMG*)mg->innerctx;
1200c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg = (PC_GAMG_AGG*)pc_gamg->subctx;
1201c8b0795cSMark F. Adams   PetscInt       jj;
1202c8b0795cSMark F. Adams   Mat            Prol  = *a_P;
12033b4367a7SBarry Smith   MPI_Comm       comm;
12042a808120SBarry Smith   KSP            eksp;
12052a808120SBarry Smith   Vec            bb, xx;
12062a808120SBarry Smith   PC             epc;
12072a808120SBarry Smith   PetscReal      alpha, emax, emin;
1208c8b0795cSMark F. Adams 
1209c8b0795cSMark F. Adams   PetscFunctionBegin;
12103b4367a7SBarry Smith   ierr = PetscObjectGetComm((PetscObject)Amat,&comm);CHKERRQ(ierr);
1211fd1112cbSBarry Smith   ierr = PetscLogEventBegin(PC_GAMGOptProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
1212c8b0795cSMark F. Adams 
12132a808120SBarry Smith   /* compute maximum value of operator to be used in smoother */
12142a808120SBarry Smith   if (0 < pc_gamg_agg->nsmooths) {
1215e7e70905SMark Adams     ierr = MatCreateVecs(Amat, &bb, 0);CHKERRQ(ierr);
1216e7e70905SMark Adams     ierr = MatCreateVecs(Amat, &xx, 0);CHKERRQ(ierr);
121714a9496bSBarry Smith     ierr = VecSetRandom(bb,pc_gamg->random);CHKERRQ(ierr);
1218e696ed0bSMark F. Adams 
12193b4367a7SBarry Smith     ierr = KSPCreate(comm,&eksp);CHKERRQ(ierr);
1220422a814eSBarry Smith     ierr = KSPSetErrorIfNotConverged(eksp,pc->erroriffailure);CHKERRQ(ierr);
1221806fa848SBarry Smith     ierr = KSPSetTolerances(eksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,10);CHKERRQ(ierr);
12222e68589bSMark F. Adams     ierr = KSPSetNormType(eksp, KSP_NORM_NONE);CHKERRQ(ierr);
1223da509fc8SJed Brown     ierr = KSPSetOptionsPrefix(eksp,((PetscObject)pc)->prefix);CHKERRQ(ierr);
1224c2436cacSMark F. Adams     ierr = KSPAppendOptionsPrefix(eksp, "gamg_est_");CHKERRQ(ierr);
1225c2436cacSMark F. Adams     ierr = KSPSetFromOptions(eksp);CHKERRQ(ierr);
1226c2436cacSMark F. Adams 
1227c2436cacSMark F. Adams     ierr = KSPSetInitialGuessNonzero(eksp, PETSC_FALSE);CHKERRQ(ierr);
122823ee1639SBarry Smith     ierr = KSPSetOperators(eksp, Amat, Amat);CHKERRQ(ierr);
12292e68589bSMark F. Adams     ierr = KSPSetComputeSingularValues(eksp,PETSC_TRUE);CHKERRQ(ierr);
12302e68589bSMark F. Adams 
1231da509fc8SJed Brown     ierr = KSPGetPC(eksp, &epc);CHKERRQ(ierr);
1232da509fc8SJed Brown     ierr = PCSetType(epc, PCJACOBI);CHKERRQ(ierr);  /* smoother in smoothed agg. */
1233c2436cacSMark F. Adams 
12342e68589bSMark F. Adams     /* solve - keep stuff out of logging */
12352e68589bSMark F. Adams     ierr = PetscLogEventDeactivate(KSP_Solve);CHKERRQ(ierr);
12362e68589bSMark F. Adams     ierr = PetscLogEventDeactivate(PC_Apply);CHKERRQ(ierr);
12372e68589bSMark F. Adams     ierr = KSPSolve(eksp, bb, xx);CHKERRQ(ierr);
12382e68589bSMark F. Adams     ierr = PetscLogEventActivate(KSP_Solve);CHKERRQ(ierr);
12392e68589bSMark F. Adams     ierr = PetscLogEventActivate(PC_Apply);CHKERRQ(ierr);
12402e68589bSMark F. Adams 
12412e68589bSMark F. Adams     ierr = KSPComputeExtremeSingularValues(eksp, &emax, &emin);CHKERRQ(ierr);
1242569f4572SMark Adams     ierr = PetscInfo3(pc,"Smooth P0: max eigen=%e min=%e PC=%s\n",emax,emin,PCJACOBI);CHKERRQ(ierr);
12432e68589bSMark F. Adams     ierr = VecDestroy(&xx);CHKERRQ(ierr);
12442e68589bSMark F. Adams     ierr = VecDestroy(&bb);CHKERRQ(ierr);
12452e68589bSMark F. Adams     ierr = KSPDestroy(&eksp);CHKERRQ(ierr);
12462e68589bSMark F. Adams   }
12472e68589bSMark F. Adams 
12482a808120SBarry Smith   /* smooth P0 */
12492a808120SBarry Smith   for (jj = 0; jj < pc_gamg_agg->nsmooths; jj++) {
12502a808120SBarry Smith     Mat       tMat;
12512a808120SBarry Smith     Vec       diag;
12522a808120SBarry Smith 
12532a808120SBarry Smith #if defined PETSC_GAMG_USE_LOG
12542a808120SBarry Smith     ierr = PetscLogEventBegin(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr);
12552a808120SBarry Smith #endif
12562a808120SBarry Smith 
12572e68589bSMark F. Adams     /* smooth P1 := (I - omega/lam D^{-1}A)P0 */
12582e68589bSMark F. Adams     ierr  = MatMatMult(Amat, Prol, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &tMat);CHKERRQ(ierr);
12592a7a6963SBarry Smith     ierr  = MatCreateVecs(Amat, &diag, 0);CHKERRQ(ierr);
12602e68589bSMark F. Adams     ierr  = MatGetDiagonal(Amat, diag);CHKERRQ(ierr); /* effectively PCJACOBI */
12612e68589bSMark F. Adams     ierr  = VecReciprocal(diag);CHKERRQ(ierr);
12622e68589bSMark F. Adams     ierr  = MatDiagonalScale(tMat, diag, 0);CHKERRQ(ierr);
12632e68589bSMark F. Adams     ierr  = VecDestroy(&diag);CHKERRQ(ierr);
1264b4ec6429SMark F. Adams     alpha = -1.4/emax;
12652e68589bSMark F. Adams     ierr  = MatAYPX(tMat, alpha, Prol, SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
12662e68589bSMark F. Adams     ierr  = MatDestroy(&Prol);CHKERRQ(ierr);
12672e68589bSMark F. Adams     Prol  = tMat;
12680cbbd2e1SMark F. Adams #if defined PETSC_GAMG_USE_LOG
12690cbbd2e1SMark F. Adams     ierr = PetscLogEventEnd(petsc_gamg_setup_events[SET9],0,0,0,0);CHKERRQ(ierr);
12702e68589bSMark F. Adams #endif
12712e68589bSMark F. Adams   }
1272fd1112cbSBarry Smith   ierr = PetscLogEventEnd(PC_GAMGOptProlongator_AGG,0,0,0,0);CHKERRQ(ierr);
1273c8b0795cSMark F. Adams   *a_P = Prol;
1274c8b0795cSMark F. Adams   PetscFunctionReturn(0);
12752e68589bSMark F. Adams }
12762e68589bSMark F. Adams 
1277c8b0795cSMark F. Adams /* -------------------------------------------------------------------------- */
1278c8b0795cSMark F. Adams /*
1279c8b0795cSMark F. Adams    PCCreateGAMG_AGG
12802e68589bSMark F. Adams 
1281c8b0795cSMark F. Adams   Input Parameter:
1282c8b0795cSMark F. Adams    . pc -
1283c8b0795cSMark F. Adams */
1284c8b0795cSMark F. Adams #undef __FUNCT__
1285c8b0795cSMark F. Adams #define __FUNCT__ "PCCreateGAMG_AGG"
1286c8b0795cSMark F. Adams PetscErrorCode  PCCreateGAMG_AGG(PC pc)
1287c8b0795cSMark F. Adams {
1288c8b0795cSMark F. Adams   PetscErrorCode ierr;
1289c8b0795cSMark F. Adams   PC_MG          *mg      = (PC_MG*)pc->data;
1290c8b0795cSMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
1291c8b0795cSMark F. Adams   PC_GAMG_AGG    *pc_gamg_agg;
12922e68589bSMark F. Adams 
1293c8b0795cSMark F. Adams   PetscFunctionBegin;
1294c8b0795cSMark F. Adams   /* create sub context for SA */
1295b00a9115SJed Brown   ierr            = PetscNewLog(pc,&pc_gamg_agg);CHKERRQ(ierr);
1296c8b0795cSMark F. Adams   pc_gamg->subctx = pc_gamg_agg;
1297c8b0795cSMark F. Adams 
12981ab5ffc9SJed Brown   pc_gamg->ops->setfromoptions = PCSetFromOptions_GAMG_AGG;
12999b8ffb57SJed Brown   pc_gamg->ops->destroy        = PCDestroy_GAMG_AGG;
1300c8b0795cSMark F. Adams   /* reset does not do anything; setup not virtual */
1301c8b0795cSMark F. Adams 
1302c8b0795cSMark F. Adams   /* set internal function pointers */
1303fd1112cbSBarry Smith   pc_gamg->ops->graph             = PCGAMGGraph_AGG;
13041ab5ffc9SJed Brown   pc_gamg->ops->coarsen           = PCGAMGCoarsen_AGG;
13051ab5ffc9SJed Brown   pc_gamg->ops->prolongator       = PCGAMGProlongator_AGG;
1306fd1112cbSBarry Smith   pc_gamg->ops->optprolongator    = PCGAMGOptProlongator_AGG;
13071ab5ffc9SJed Brown   pc_gamg->ops->createdefaultdata = PCSetData_AGG;
13085adeb434SBarry Smith   pc_gamg->ops->view              = PCView_GAMG_AGG;
1309c8b0795cSMark F. Adams 
1310e0877f53SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetNSmooths_C",PCGAMGSetNSmooths_AGG);CHKERRQ(ierr);
1311e0877f53SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetSymGraph_C",PCGAMGSetSymGraph_AGG);CHKERRQ(ierr);
1312e0877f53SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGAMGSetSquareGraph_C",PCGAMGSetSquareGraph_AGG);CHKERRQ(ierr);
1313bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_AGG);CHKERRQ(ierr);
13142e68589bSMark F. Adams   PetscFunctionReturn(0);
13152e68589bSMark F. Adams }
1316